Tomographer  v4.0
Tomographer C++ Framework Documentation
tspacellhwalker.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_TSPACELLHWALKER_H
29 #define TOMOGRAPHER_DENSEDM_TSPACELLHWALKER_H
30 
31 #include <cstddef>
32 #include <cmath>
33 
34 #include <random>
35 
36 #include <boost/math/constants/constants.hpp>
37 
43 #include <tomographer/mhrw.h>
44 
53 namespace Tomographer {
54 namespace DenseDM {
55 namespace TSpace {
56 
71 template<typename DenseLLH_, typename Rng, typename LoggerType>
72 TOMOGRAPHER_EXPORT class LLHMHWalker
73  : public Tools::NeedOwnOperatorNew<typename DenseLLH_::DMTypes::MatrixType>::ProviderType
74 {
75 public:
77  typedef DenseLLH_ DenseLLH;
79  typedef typename DenseLLH::DMTypes DMTypes;
81  typedef typename DenseLLH::LLHValueType LLHValueType;
83  typedef typename DMTypes::MatrixType MatrixType;
87  typedef typename DMTypes::RealScalar RealScalar;
90 
94  typedef RealScalar StepRealType;
95 
97  typedef MatrixType PointType;
99  typedef LLHValueType FnValueType;
101  enum {
106  };
107 
108 private:
109 
110  const DenseLLH & _llh;
111  Rng & _rng;
112  std::normal_distribution<RealScalar> _normal_distr_rnd;
113  const ParamX<DMTypes> _px;
114 
115  LoggerType & _log;
116 
117  MatrixType _startpt;
118 
119 public:
120 
131  LLHMHWalker(const MatrixType & startpt, const DenseLLH & llh, Rng & rng,
132  LoggerType & log_)
133  : _llh(llh),
134  _rng(rng),
135  _normal_distr_rnd(0.0, 1.0),
136  _px(llh.dmt),
137  _log(log_),
138  _startpt(startpt)
139  {
140  }
141 
146  LLHMHWalker(LLHMHWalker && moveother)
147  : _llh(std::move(moveother._llh)),
148  _rng(moveother._rng),
149  _normal_distr_rnd(moveother._normal_distr_rnd),
150  _px(std::move(moveother._px)),
151  _log(moveother._log),
152  _startpt(moveother._startpt)
153  {
154  }
155 
156 
166  inline void init()
167  {
168  _log.debug("TSpace::LLHMHWalker", "Starting random walk");
169  }
170 
172  inline const MatrixType & startPoint()
173  {
174  // It's fine to hard-code "1e-3" because for any type, valid T-matrices have norm == 1
175  if (_startpt.norm() > 1e-3) {
176  // nonzero matrix given: that's the starting point.
177  return _startpt;
178  }
179 
180  // zero matrix given: means to choose random starting point
181  MatrixType T(_llh.dmt.initMatrixType());
182  T = Tools::denseRandom<MatrixType>(
183  _rng, _normal_distr_rnd, _llh.dmt.dim(), _llh.dmt.dim()
184  );
185  _startpt = T/T.norm(); // normalize to be on surface of the sphere
186 
187  _log.debug("TSpace::LLHMHWalker", [&](std::ostream & str) {
188  str << "Chosen random start point T = \n" << _startpt;
189  });
190 
191  // return start point
192  return _startpt;
193  }
194 
196  inline void thermalizingDone()
197  {
198  }
199 
201  inline void done()
202  {
203  }
204 
212  TOMOGRAPHER_ENABLED_IF(DenseLLH::LLHCalcType == LLHCalcTypeX)
213  inline LLHValueType fnLogVal(const MatrixType & T) const
214  {
215  MatrixType rho(T*T.adjoint());
216  VectorParamType x = _px.HermToX(rho);
217  LLHValueType llhval = _llh.logLikelihoodX(x);
218  // _log.longdebug("fnlogval(%s) = %g\n", streamcstr(x.transpose()), llhval);
219  return llhval;
220  }
221 
229  TOMOGRAPHER_ENABLED_IF(DenseLLH::LLHCalcType == LLHCalcTypeRho)
230  inline LLHValueType fnLogVal(const MatrixType & T) const
231  {
232  LLHValueType llhval = _llh.logLikelihoodRho(T*T.adjoint());
233  // _log.longdebug("fnlogval( %s ) = %g\n", streamcstr(T*T.adjoint()), llhval);
234  return llhval;
235  }
236 
237 
239  inline MatrixType jumpFn(const MatrixType& cur_T, RealScalar step_size)
240  {
241  MatrixType DeltaT(Tools::denseRandom<MatrixType>(
242  _rng, _normal_distr_rnd, _llh.dmt.dim(), _llh.dmt.dim()
243  ));
244 
245  MatrixType new_T(cur_T + step_size * DeltaT);
246 
247  // renormalize to "project" onto the large T-space sphere
248  new_T /= new_T.norm(); //Matrix<>.norm() is Frobenius norm.
249 
250  // _log.longdebug("jump_fn(): step_size=%g, cur_T =\n%s\nDeltaT = \n%s\nnew_T = \n%s",
251  // step_size, streamstr(cur_T).c_str(), streamstr(DeltaT).c_str(),
252  // streamstr(new_T).c_str());
253 
254  // hope for Return Value Optimization by compiler
255  return new_T;
256  }
257 
258 };
259 
260 
261 
262 } // namespace TSpace
263 } // namespace DenseDM
264 } // namespace Tomographer
265 
266 
267 #endif
void init()
Provided for MHRandomWalk. Initializes some fields and prepares for a random walk.
Base namespace for the Tomographer project.
Definition: densellh.h:45
const MatrixType & startPoint()
Return the starting point given in the constructor, or a random start point.
DMTypes::VectorParamType VectorParamType
Type of an X-parameterization of a density operator (see X Parameterization)
C++ types for describing dense density matrices in various parameterizations.
DenseLLH::DMTypes DMTypes
The data types of our problem.
Provide appropriate operator new() definitions for a structure which has a member of the given stored...
STL namespace.
LLHMHWalker(LLHMHWalker &&moveother)
A MHWalker Interface needs to have a move constructor.
LLHMHWalker(const MatrixType &startpt, const DenseLLH &llh, Rng &rng, LoggerType &log_)
Constructor which just initializes the given fields.
STL class.
DMTypes::RealScalar RealScalar
The real scalar corresponding to our data types. Usually a double.
DenseLLH::LLHValueType LLHValueType
The loglikelihood function value type (see DenseLLH Interface e.g. IndepMeasLLH)
RealScalar StepRealType
The real scalar type for the step size — this is our scalar type (see MHWalker Interface) ...
MatrixType PointType
Provided for MHRandomWalk. A point in our random walk = a density matrix.
Tools for parameterizing hermitian matrices with the X Parameterization.
Some declarations for the DenseLLH type interface.
VectorParamType HermToX(MatrixTypeConstRef Herm) const
Get the X-parameterization corresponding to a given hermitian matrix.
Definition: param_herm_x.h:81
A random walk in the density matrix space of a Hilbert state space of a quantum system.
Provides the logarithm MH function value at each point (see Role of UseFnSyntaxType) ...
Definition: mhrw.h:68
void done()
Callback for after random walk is finished. No-op.
DenseLLH_ DenseLLH
The DenseLLH interface object type.
Managing the need for specific overrides to operator new() for some types (especially Eigen types) ...
void thermalizingDone()
Callback for after thermalizing is done. No-op.
LLHValueType FnValueType
Provided for MHRandomWalk. The function value type is the loglikelihood value type.
The DenseLLH-compatible object exposes a method logLikelihoodRho(), taking as argument a (const ref t...
Definition: densellh.h:61
MatrixType jumpFn(const MatrixType &cur_T, RealScalar step_size)
Decides of a new point to jump to for the random walk.
We will calculate the log-likelihood function, which is the logarithm of the Metropolis-Hastings func...
LLHValueType fnLogVal(const MatrixType &T) const
Calculate the logarithm of the Metropolis-Hastings function value.
Routines for performing a Metropolis-Hastings random walk.
DMTypes::MatrixType MatrixType
The matrix type for a density operator on our quantum system.
STL class.
DMTypes::ComplexScalar ComplexScalar
The complex real scalar corresponding to our data types. Usually a std::complex<double>.
The DenseLLH-compatible object exposes a method logLikelihoodX(), taking as argument a (const ref to ...
Definition: densellh.h:66
Utilities for logging messages.
RealScalar_ RealScalar
Real scalar type, given in template parameter. Usually double is fine.
Definition: dmtypes.h:121