Tomographer
v5.4
Tomographer C++ Framework Documentation
|
A random walk in the density matrix space of a Hilbert state space of a quantum system (uniform random unitary jumps; slower) More...
#include <tomographer/densedm/tspacellhwalker.h>
Public Types | |
enum | |
see MHWalker Interface More... | |
typedef DenseLLHType_ | DenseLLHType |
The DenseLLH interface object type. | |
typedef RngType_ | RngType |
The random number generator type. | |
typedef LoggerType_ | LoggerType |
The logger type. | |
typedef DenseLLHType::DMTypes | DMTypes |
The data types of our problem. | |
typedef DenseLLHType::LLHValueType | LLHValueType |
The loglikelihood function value type (see DenseLLH Interface e.g. IndepMeasLLH) | |
typedef DMTypes::MatrixType | MatrixType |
The matrix type for a density operator on our quantum system. | |
typedef DMTypes::VectorParamType | VectorParamType |
Type of an X-parameterization of a density operator (see X Parameterization) | |
typedef DMTypes::RealScalar | RealScalar |
The real scalar corresponding to our data types. Usually a double . | |
typedef DMTypes::ComplexScalar | ComplexScalar |
The complex real scalar corresponding to our data types. Usually a std::complex<double> . | |
typedef MHWalkerParamsStepSize< RealScalar > | WalkerParams |
The real scalar type for the step size — this scalar type is the only parameter we need (see MHWalker Interface) | |
typedef MatrixType | PointType |
Provided for MHRandomWalk. A point in our random walk = a density matrix. | |
typedef LLHValueType | FnValueType |
Provided for MHRandomWalk. The function value type is the loglikelihood value type. | |
Public Member Functions | |
LLHMHWalker (const MatrixType &startpt, const DenseLLHType &llh, RngType &rng, LoggerType &log_) | |
Constructor which just initializes the given fields. More... | |
void | init () |
Provided for MHRandomWalk. Initializes some fields and prepares for a random walk. More... | |
const MatrixType & | startPoint () |
Return the starting point given in the constructor, or a random start point. | |
void | thermalizingDone () |
Callback for after thermalizing is done. No-op. | |
void | done () |
Callback for after random walk is finished. No-op. | |
LLHValueType | fnLogVal (const MatrixType &T) const |
Calculate the logarithm of the Metropolis-Hastings function value. More... | |
MatrixType | jumpFn (const MatrixType &cur_T, WalkerParams params) |
Decides of a new point to jump to for the random walk. | |
A random walk in the density matrix space of a Hilbert state space of a quantum system (uniform random unitary jumps; slower)
The random walk explores the density operators on a Hilbert space of a given dimension, according to the distribution given by the likelihood function on the Hilbert-schmidt uniform prior.
DenseLLHType | A type satisfying the DenseLLH Interface |
RngType | A std::random random number generator (such as std::mt19937) |
LoggerType | A logger type (see Logging and Loggers) |
Definition at line 120 of file tspacellhwalker.h.
anonymous enum |
Enumerator | |
---|---|
UseFnSyntaxType | We will calculate the log-likelihood function, which is the logarithm of the Metropolis-Hastings function we should be calculating. |
Definition at line 154 of file tspacellhwalker.h.
|
inline |
Constructor which just initializes the given fields.
The llh object is responsible for calculating the log-likelihood function of the tomography experiment (see class documentation).
The startpt argument specifies the starting point for the random walk (matrix in T Parameterization). If you provide a zero startpt here, then a random starting point will be chosen using the rng random number generator to generate a random point on the sphere.
Definition at line 184 of file tspacellhwalker.h.
|
inline |
Calculate the logarithm of the Metropolis-Hastings function value.
Definition at line 248 of file tspacellhwalker.h.
|
inline |
Provided for MHRandomWalk. Initializes some fields and prepares for a random walk.
Definition at line 204 of file tspacellhwalker.h.