Tomographer  v1.0a
Tomographer C++ Framework Documentation
tomoproblem.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 TOMOPROBLEM_H
28 #define TOMOPROBLEM_H
29 
37 namespace Tomographer {
38 
39 
56 template<typename MatrQ_, typename LLHValueType_ = double, bool UseCLoopInstead = false>
58 {
59  typedef MatrQ_ MatrQ;
60  typedef LLHValueType_ LLHValueType;
61 
63  const MatrQ matq;
64 
66  const int dim;
68  const int dim2;
70  const int Ndof;
71 
73  LLHValueType NMeasAmplifyFactor;
74 
79 
83 
89 
96  IndepMeasTomoProblem(MatrQ matq_)
97  : matq(matq_),
98  dim(matq_.dim()),
99  dim2(dim*dim),
100  Ndof(dim2 - 1),
101  NMeasAmplifyFactor(1),
102  rho_MLE(matq.initMatrixType()),
103  T_MLE(matq.initMatrixType()),
104  x_MLE(matq.initVectorParamType())
105  {
106  // NOTE: Exn & N are left uninitialized, because we don't yet know how many POVM
107  // effects there will be. (That will probably be read from a MAT file or such...)
108  assert(MatrQ::FixedDim == Eigen::Dynamic || MatrQ::FixedDim == (int)matq.dim());
109  }
110 
118  inline LLHValueType calc_llh(const typename MatrQ::VectorParamType & x) const;
119 };
120 
121 
122 
123 namespace tomo_internal
124 {
132  template<typename MatrQ, typename LLHValueType, bool UseCLoopInstead>
133  struct LLH_Calculator
134  {
144  static inline LLHValueType calc_llh(const IndepMeasTomoProblem<MatrQ,LLHValueType,UseCLoopInstead> *data,
145  const typename MatrQ::VectorParamType & x)
146  {
147  return -2 * data->NMeasAmplifyFactor * (
148  data->Nx.template cast<LLHValueType>() * (data->Exn * x).array().log()
149  ).sum();
150  }
151  };
152 
153  //
154  // Specialization of helper for implementing the C for loop instead.
155  //
156  // no doxygen doc needed.
157  //
158  template<typename MatrQ, typename LLHValueType>
159  struct LLH_Calculator<MatrQ, LLHValueType, true>
160  {
161  static inline LLHValueType calc_llh(const IndepMeasTomoProblem<MatrQ,LLHValueType,true> *data,
162  const typename MatrQ::VectorParamType & x)
163  {
164  std::size_t k;
165  typename MatrQ::RealScalar val = 0;
166  for (k = 0; k < (std::size_t)data->Exn.rows(); ++k) {
167  const typename MatrQ::RealScalar thisval = (data->Exn.row(k) * x);
168  //std::cout << "k = "<<k<<"; thisval = "<<thisval << "\n";
169  val += data->Nx[k] * std::log(thisval);
170  }
171  //std::cout << "val = " << val << "\n\n";
172  return -2 * data->NMeasAmplifyFactor * val;
173  }
174  };
175 }
176 
177 
178 template<typename MatrQ, typename LLHValueType, bool UseCLoopInstead>
180  const typename MatrQ::VectorParamType & x
181  ) const
182 {
183  return tomo_internal::LLH_Calculator<MatrQ, LLHValueType, UseCLoopInstead>::calc_llh(this, x);
184 }
185 
186 
187 
188 
189 
190 } // namespace Tomographer
191 
192 
193 
194 #endif
Eigen::Matrix< std::complex< RealScalar >, FixedDim, FixedDim > MatrixType
Complex dim x dim Matrix.
Definition: matrq.h:132
Base namespace for the Tomographer project.
Definition: dmmhrw.h:51
T log(T...args)
IndepMeasTomoProblem(MatrQ matq_)
Constructs an IndepMeasTomoProblem instance.
Definition: tomoproblem.h:96
const MatrQ matq
The data types for this problem.
Definition: tomoproblem.h:63
Eigen::Array< IntFreqType, Eigen::Dynamic, 1, 0, FixedMaxParamList, 1 > FreqListType
dynamic Array of integers [maximum FixedMaxParamList entries or Dynamic]
Definition: matrq.h:149
const int dim2
The square of the dimension of the Hilbert space, dim2=dim*dim.
Definition: tomoproblem.h:68
MatrQ::MatrixType T_MLE
Maximum likelihood estimate as T-parameterized density matrix.
Definition: tomoproblem.h:86
const int Ndof
The number of degrees of freedom, .
Definition: tomoproblem.h:70
MatrQ::FreqListType Nx
The frequency list, i.e. number of times each POVM effect was observed.
Definition: tomoproblem.h:78
tomo_internal::matrq_traits< MatrQ< FixedDim_, FixedMaxParamList_, RealScalar_, IntFreqType_ > >::RealScalar RealScalar
Real scalar type (usually double)
Definition: matrq.h:125
MatrQ::VectorParamType x_MLE
X-Parameterized version of rho_MLE (param_herm_to_x())
Definition: tomoproblem.h:88
Eigen::Matrix< RealScalar, Eigen::Dynamic, FixedDim2, Eigen::RowMajor, FixedMaxParamList, FixedDim2 > VectorParamListType
dynamic Matrix with rows = dim*dim Vectors (row-major) [maximum FixedMaxParamList rows...
Definition: matrq.h:144
MatrQ::MatrixType rho_MLE
Maximum likelihood estimate as density matrix.
Definition: tomoproblem.h:81
LLHValueType calc_llh(const typename MatrQ::VectorParamType &x) const
Calculates the log-likelihood function.
Definition: tomoproblem.h:179
const int dim
The dimension of the Hilbert space.
Definition: tomoproblem.h:66
Eigen::Matrix< RealScalar, FixedDim2, 1 > VectorParamType
Real dim*dim Vector.
Definition: matrq.h:135
Stores data for a tomography problem with independent measurements.
Definition: tomoproblem.h:57
LLHValueType NMeasAmplifyFactor
A factor to artificially amplify the number of measurements by.
Definition: tomoproblem.h:73
MatrQ::VectorParamListType Exn
The POVM Entries, parameterized with X-param (param_herm_to_x())
Definition: tomoproblem.h:76