Tomographer  v5.3
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 #include <boost/serialization/serialization.hpp>
44 
45 #include <Eigen/Core>
46 
47 #include <tomographer/tools/cxxutil.h> // static_or_dynamic, tomographer_assert()
48 
49 
50 
51 namespace Tomographer {
52 namespace DenseDM {
53 
57 template<typename DMTypes_>
58 class TOMOGRAPHER_EXPORT ParamX
59 {
60 public:
61  typedef DMTypes_ DMTypes;
62  typedef typename DMTypes::MatrixType MatrixType;
66  typedef typename DMTypes::RealScalar RealScalar;
67  typedef typename DMTypes::ComplexScalar ComplexScalar;
68  typedef Eigen::Index IndexType;
69 
73  ParamX(DMTypes dmt) : dim(dmt.dim()) { }
74 
82  inline VectorParamType HermToX(MatrixTypeConstRef Herm) const
83  {
84  // don't construct a DMTypes object, because it could have a special constructor...
85  //DMTypes dmt(dim);
86  const Eigen::Index dim2 = dim*dim;
87 
88  // hope RVO kicks in
89  VectorParamType x = VectorParamType::Zero(dim2,1);
90  const Eigen::Index dimtri = (dim2 - dim)/2;
91 
92  tomographer_assert(dim == Herm.cols()); // assert Herm is (dim x dim)
93 
94  x.block(0,0,dim,1) = Herm.real().diagonal();
95 
96  IndexType k = dim;
97  IndexType n, m;
98  for (n = 1; n < dim; ++n) {
99  for (m = 0; m < n; ++m) {
100  x(k) = Herm(n,m).real() * boost::math::constants::root_two<RealScalar>();
101  x(dimtri + k) = Herm(n,m).imag() * boost::math::constants::root_two<RealScalar>();
102  ++k;
103  }
104  }
105 
106  return x;
107  }
108 
114  template<bool OnlyLowerTri = false>
115  inline MatrixType XToHerm(VectorParamTypeConstRef x) const
116  {
117  const Eigen::Index dim2 = dim*dim;
118 
119  // should be optimized by compiler via RVO
120  MatrixType Herm = MatrixType::Zero(dim,dim);
121 
122  const IndexType dimtri = (dim2-dim)/2;
123  tomographer_assert(x.rows() == dim2 && x.cols() == 1); // assert x is (dim*dim x 1)
124 
125  Herm.diagonal().real() = x.block(0,0,dim,1);
126  Herm.diagonal().imag().setZero();
127 
128  IndexType k = dim;
129  IndexType n, m;
130  for (n = 1; n < dim; ++n) {
131  for (m = 0; m < n; ++m) {
132  Herm(n,m) = boost::math::constants::half_root_two<RealScalar>() * ComplexScalar(x(k), x(dimtri + k));
133  if (!OnlyLowerTri) {
134  // complex conj. on opposite triangular part
135  Herm(m,n) = boost::math::constants::half_root_two<RealScalar>() * ComplexScalar(x(k), -x(dimtri + k));
136  }
137  ++k;
138  }
139  }
140  return Herm;
141  }
142 
143 
146  ParamX() : dim(-1) { }
147 
148 protected:
149  Eigen::Index dim;
150 
151 private:
152  friend boost::serialization::access;
153  template<typename Archive>
154  void serialize(Archive & a, unsigned int /* version */)
155  {
156  a & dim;
157  }
158 };
159 
160 
161 
162 
163 } // namespace DenseDM
164 } // namespace Tomographer
165 
166 
167 #endif
Base namespace for the Tomographer project.
Definition: densellh.h:45
const ImagReturnType imag() const
Convert hermitian matrices to vectors via their X Parameterization.
Definition: param_herm_x.h:58
STL class.
RealReturnType real() const
MatrixType XToHerm(VectorParamTypeConstRef x) const
Get the Hermitian matrix parameterized by the "X-parameter" vector x.
Definition: param_herm_x.h:115
Block< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > > block(Index startRow, Index startCol, Index blockRows, Index blockCols)
ParamX()
Construct an invalid ParamX object – only for use with Boost.Serialization.
Definition: param_herm_x.h:146
VectorParamType HermToX(MatrixTypeConstRef Herm) const
Get the X-parameterization corresponding to a given hermitian matrix.
Definition: param_herm_x.h:82
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:73
#define tomographer_assert(...)
Assertion test macro.
Definition: cxxdefs.h:84
RealScalar_ RealScalar
Real scalar type, given in template parameter. Usually double is fine.
Definition: dmtypes.h:121