Tomographer  v3.0
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) 2016 ETH Zurich, Institute for Theoretical Physics, Philippe Faist
7  * Copyright (c) 2017 Caltech, Institute for Quantum Information and Matter, Philippe Faist
8  *
9  * Permission is hereby granted, free of charge, to any person obtaining a copy
10  * of this software and associated documentation files (the "Software"), to deal
11  * in the Software without restriction, including without limitation the rights
12  * to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
13  * copies of the Software, and to permit persons to whom the Software is
14  * furnished to do so, subject to the following conditions:
15  *
16  * The above copyright notice and this permission notice shall be included in
17  * all copies or substantial portions of the Software.
18  *
19  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
20  * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
21  * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
22  * AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
23  * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
24  * OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
25  * SOFTWARE.
26  */
27 
28 #ifndef TOMOGRAPHER_DENSEDM_PARAM_HERM_X_H
29 #define TOMOGRAPHER_DENSEDM_PARAM_HERM_X_H
30 
31 
32 
39 #include <cmath>
40 #include <complex>
41 
42 #include <boost/math/constants/constants.hpp>
43 
44 #include <Eigen/Core>
45 
46 #include <tomographer/tools/cxxutil.h> // static_or_dynamic, tomographer_assert()
47 
48 
49 
50 namespace Tomographer {
51 namespace DenseDM {
52 
56 template<typename DMTypes_>
57 class ParamX {
58 public:
59  typedef DMTypes_ DMTypes;
60  typedef typename DMTypes::MatrixType MatrixType;
64  typedef typename DMTypes::RealScalar RealScalar;
65  typedef typename DMTypes::ComplexScalar ComplexScalar;
66  typedef typename MatrixType::Index IndexType;
67 
71  ParamX(DMTypes dmt) : _dmt(dmt) { }
72 
80  inline VectorParamType HermToX(MatrixTypeConstRef Herm) const
81  {
82  // hope RVO kicks in
83  VectorParamType x(_dmt.initVectorParamType());
84  const IndexType dimtri = (_dmt.dim2() - _dmt.dim())/2;
85 
86  tomographer_assert((IndexType)_dmt.dim() == Herm.cols()); // assert Herm is (dim x dim)
87 
88  x.block(0,0,_dmt.dim(),1) = Herm.real().diagonal();
89 
90  IndexType k = _dmt.dim();
91  IndexType n, m;
92  for (n = 1; n < (IndexType)_dmt.dim(); ++n) {
93  for (m = 0; m < n; ++m) {
94  x(k) = Herm(n,m).real() * boost::math::constants::root_two<RealScalar>();
95  x(dimtri + k) = Herm(n,m).imag() * boost::math::constants::root_two<RealScalar>();
96  ++k;
97  }
98  }
99 
100  return x;
101  }
102 
108  template<bool OnlyLowerTri = false>
109  inline MatrixType XToHerm(VectorParamTypeConstRef x) const
110  {
111  // should be optimized by compiler via RVO
112  MatrixType Herm(_dmt.initMatrixType());
113 
114  const IndexType dimtri = (_dmt.dim2()-_dmt.dim())/2;
115  tomographer_assert(x.rows() == (IndexType)_dmt.dim2() && x.cols() == 1); // assert x is (dim*dim x 1)
116 
117  Herm.diagonal().real() = x.block(0,0,_dmt.dim(),1);
118  Herm.diagonal().imag().setZero();
119 
120  IndexType k = _dmt.dim();
121  IndexType n, m;
122  for (n = 1; n < (IndexType)_dmt.dim(); ++n) {
123  for (m = 0; m < n; ++m) {
124  Herm(n,m) = boost::math::constants::half_root_two<RealScalar>() * ComplexScalar(x(k), x(dimtri + k));
125  if (!OnlyLowerTri) {
126  // complex conj. on opposite triangular part
127  Herm(m,n) = boost::math::constants::half_root_two<RealScalar>() * ComplexScalar(x(k), -x(dimtri + k));
128  }
129  ++k;
130  }
131  }
132  return Herm;
133  }
134 
135 protected:
136  const DMTypes _dmt;
137 };
138 
139 
140 
141 
142 } // namespace DenseDM
143 } // namespace Tomographer
144 
145 
146 #endif
Base namespace for the Tomographer project.
Definition: densellh.h:45
internal::traits< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > >::Index Index
Convert hermitian matrices to vectors via their X Parameterization.
Definition: param_herm_x.h:57
STL class.
MatrixType XToHerm(VectorParamTypeConstRef x) const
Get the Hermitian matrix parameterized by the "X-parameter" vector x.
Definition: param_herm_x.h:109
VectorParamType HermToX(MatrixTypeConstRef Herm) const
Get the X-parameterization corresponding to a given hermitian matrix.
Definition: param_herm_x.h:80
Some C++ utilities, with a tad of C++11 tricks.
ParamX(DMTypes dmt)
Constructor. Just give it the DMTypes instance.
Definition: param_herm_x.h:71
RealScalar_ RealScalar
Real scalar type, given in template parameter. Usually double is fine.
Definition: dmtypes.h:120