Tomographer  v1.0a
Tomographer C++ Framework Documentation
param_herm_x.h
Go to the documentation of this file.
1 /* This file is part of the Tomographer project, which is distributed under the
2  * terms of the MIT license.
3  *
4  * The MIT License (MIT)
5  *
6  * Copyright (c) 2015 ETH Zurich, Institute for Theoretical Physics, Philippe Faist
7  *
8  * Permission is hereby granted, free of charge, to any person obtaining a copy
9  * of this software and associated documentation files (the "Software"), to deal
10  * in the Software without restriction, including without limitation the rights
11  * to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
12  * copies of the Software, and to permit persons to whom the Software is
13  * furnished to do so, subject to the following conditions:
14  *
15  * The above copyright notice and this permission notice shall be included in
16  * all copies or substantial portions of the Software.
17  *
18  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
19  * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
20  * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
21  * AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
22  * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
23  * OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
24  * SOFTWARE.
25  */
26 
27 #ifndef HERM_PARAM_X_H
28 #define HERM_PARAM_X_H
29 
36 #include <cmath>
37 #include <complex>
38 
39 #include <boost/math/constants/constants.hpp>
40 
41 #include <Eigen/Core>
42 
43 
44 namespace Tomographer {
45 
46 
52 template<bool OnlyLowerTri=false, typename Derived1=Eigen::MatrixXd, typename Derived2=Eigen::MatrixXd>
54 {
55  { using namespace Eigen; EIGEN_STATIC_ASSERT_LVALUE(Derived1); }
56 
57  const int dim = Herm.rows();
58  const int dimtri = dim*(dim-1)/2;
59  eigen_assert(dim == Herm.cols()); // assert Herm is (dim x dim)
60  eigen_assert(x.rows() == dim*dim && x.cols() == 1); // assert x is (dim*dim x 1)
61 
62  typedef typename Derived1::Scalar Scalar;
63  typedef typename Eigen::NumTraits<Scalar>::Real RealScalar;
64 
65  Herm.diagonal().real() = x.block(0,0,dim,1);
66  Herm.diagonal().imag().setZero();
67 
68  int n, m;
69  for (n = 1; n < dim; ++n) {
70  for (m = 0; m < n; ++m) {
71  const int k = dim + n*(n-1)/2 + m;
72  const int l = dimtri + k;
73  Herm(n,m) = boost::math::constants::half_root_two<RealScalar>() * Scalar(x(k), x(l));
74  if (!OnlyLowerTri) {
75  // complex conj. on opposite triangular part
76  Herm(m,n) = boost::math::constants::half_root_two<RealScalar>() * Scalar(x(k), -x(l));
77  }
78  }
79  }
80 };
81 
82 
92 template<typename Derived1, typename Derived2>
94 {
95  { using namespace Eigen; EIGEN_STATIC_ASSERT_LVALUE(Derived1); }
96 
97  const int dim = Herm.rows();
98  const int dimtri = dim*(dim-1)/2;
99  eigen_assert(dim == Herm.cols()); // assert Herm is (dim x dim)
100  eigen_assert(x.rows() == dim*dim && x.cols() == 1); // assert x is (dim*dim x 1)
101 
102  typedef typename Derived1::Scalar Scalar;
103  typedef typename Eigen::NumTraits<Scalar>::Real RealScalar;
104 
105  x.block(0,0,dim,1) = Herm.real().diagonal();
106 
107  int n, m;
108  for (n = 1; n < dim; ++n) {
109  for (m = 0; m < n; ++m) {
110  const int k = dim + n*(n-1)/2 + m;
111  const int l = dimtri + k;
112  x(k) = Herm(n,m).real() * boost::math::constants::root_two<RealScalar>();
113  x(l) = Herm(n,m).imag() * boost::math::constants::root_two<RealScalar>();
114  }
115  }
116 }
117 
118 
119 } // namespace Tomographer
120 
121 
122 #endif
Base namespace for the Tomographer project.
Definition: dmmhrw.h:51
const ImagReturnType imag() const
RealReturnType real() const
Block< Derived > block(Index startRow, Index startCol, Index blockRows, Index blockCols)
DiagonalReturnType diagonal()
void param_x_to_herm(Eigen::MatrixBase< Derived1 > &Herm, const Eigen::DenseBase< Derived2 > &x)
Get the Hermitian matrix parameterized by the "X-parameter" vector x.
Definition: param_herm_x.h:53
void param_herm_to_x(Eigen::DenseBase< Derived1 > &x, const Eigen::MatrixBase< Derived2 > &Herm)
Get the X-parameterization corresponding to a given hermitian matrix.
Definition: param_herm_x.h:93