Toast++  1.0.2 (r.539)
Forward and inverse modelling in optical tomography
Public Member Functions | Protected Member Functions | List of all members
TFwdSolverMW< T > Class Template Reference

Templated forward solver class for multi-wavelength problems. More...

#include <fwdsolver_mw.h>

Inheritance diagram for TFwdSolverMW< T >:
TFwdSolver< T >

Public Member Functions

 TFwdSolverMW (const QMMesh *mesh, LSOLVER linsolver, double tol=1e-10)
 Constructor. Creates a forward solver instance. More...
 
 TFwdSolverMW (const QMMesh *mesh, char *solver, double tol=1e-10)
 Constructor. Creates a forward solver instance. More...
 
 TFwdSolverMW (const QMMesh *mesh, ParamParser &pp)
 Constructor. Creates a forward solver instance. More...
 
 ~TFwdSolverMW ()
 Destructor. Destroys the forward solver instance.
 
TVector< T > ProjectAll_wavel (const TCompRowMatrix< T > &qvec, const TCompRowMatrix< T > &mvec, const MWsolution &sol, double omega, DataScale scl=DATA_DEFAULT)
 Return boundary data for all sources and all detectors at all wavelengths. More...
 
RVector ProjectAll_wavel_real (const TCompRowMatrix< T > &qvec, const TCompRowMatrix< T > &mvec, const MWsolution &sol, double omega, DataScale scl=DATA_DEFAULT)
 Return boundary data for the complex case in a real vector. More...
 
- Public Member Functions inherited from TFwdSolver< T >
 TFwdSolver (const QMMesh *mesh, LSOLVER linsolver, double tol=1e-10)
 Constructor. Creates a forward solver instance. More...
 
 TFwdSolver (const QMMesh *mesh, const char *solver, double tol=1e-10)
 Constructor. Creates a forward solver instance. More...
 
 TFwdSolver (const QMMesh *mesh, ParamParser &pp)
 Constructor. Creates a forward solver instance. More...
 
 ~TFwdSolver ()
 Destructor. Destroys the forward solver instance.
 
void SetDataScaling (DataScale scl)
 Set the default scaling for data returned by the projection methods. More...
 
DataScale GetDataScaling () const
 Returns the current default setting for data scaling. More...
 
void SetPrecon (PreconType type)
 Set the preconditioning method. More...
 
void SetPhaseUnwrap (bool unwrap)
 
void ReadParams (ParamParser &pp)
 Read solver parameters from a parameter file. More...
 
void WriteParams (ParamParser &pp)
 Write the current solver settings to a parameter file. More...
 
void SetLinSolver (const char *solver, double tol=1e-10)
 Set the solver type and tolerance. More...
 
LSOLVER LinSolver () const
 Returns the current solver type. More...
 
double GetLinSolverTol () const
 Returns the solver tolerance. More...
 
void SetLinSolverMaxit (int maxit)
 Set max iteration count for iterative solver. More...
 
const QMMeshMeshPtr () const
 Returns a pointer to the FEM mesh associated with the forward solver instance. More...
 
void Allocate ()
 Evaluates fill structure of system matrices and allocates dynamic memory for them. More...
 
void AssembleSystemMatrix (const Solution &sol, double omega=0, bool elbasis=false)
 Construct the FEM system matrix from a set of parameter distributions for a given modulation frequency. More...
 
void AssembleSystemMatrixComponent (const RVector &prm, int type)
 
void AssembleMassMatrix (const Mesh *mesh=0)
 Construct the FEM mass matrix for time-dependent problems. More...
 
void Reset (const Solution &sol, double omega=0, bool elbasis=false)
 Reset the forward solver by re-building the system matrices from a new set of parameters. More...
 
void CalcField (const TVector< T > &qvec, TVector< T > &phi, IterativeSolverResult *res=0) const
 Calculate photon density field for a given source vector. More...
 
void CalcFields (const TCompRowMatrix< T > &qvec, TVector< T > *phi, IterativeSolverResult *res=0) const
 Calculate photon density fields for all sources. More...
 
void CalcFields (int q0, int q1, const TCompRowMatrix< T > &qvec, TVector< T > *phi, IterativeSolverResult *res=0) const
 Calculate photon density fields for a range of sources. More...
 
TVector< T > ProjectSingle (int q, const TCompRowMatrix< T > &mvec, const TVector< T > &phi, DataScale scl=DATA_DEFAULT) const
 Return boundary data for a single source, given the corresponding photon density field. More...
 
TVector< T > ProjectAll (const TCompRowMatrix< T > &mvec, const TVector< T > *dphi, DataScale scl=DATA_DEFAULT)
 Return boundary data for all sources and all detectors, given photon density fields for all sources. More...
 
TVector< T > ProjectAll (const TCompRowMatrix< T > &qvec, const TCompRowMatrix< T > &mvec, const Solution &sol, double omega, DataScale scl=DATA_DEFAULT)
 Return boundary data for all sources and all detectors, given parameter distributions. More...
 
RVector ProjectAll_real (const TCompRowMatrix< T > &mvec, const TVector< T > *dphi, DataScale scl=DATA_DEFAULT)
 Return boundary data for the complex case in a real vector. More...
 
FVector ProjectAll_singlereal (const TCompRowMatrix< T > &mvec, const TVector< T > *dphi, DataScale scl=DATA_DEFAULT)
 
RVector ProjectAll_real (const TCompRowMatrix< T > &qvec, const TCompRowMatrix< T > &mvec, const Solution &sol, double omega, DataScale scl=DATA_DEFAULT)
 Return boundary data for all sources and all detectors, given parameters defined via a Solution instance. More...
 
FVector ProjectAll_singlereal (const TCompRowMatrix< T > &qvec, const TCompRowMatrix< T > &mvec, const Solution &sol, double omega, DataScale scl=DATA_DEFAULT)
 

Protected Member Functions

void Setup ()
 Forward solver initialisation routines.
 
void Cleanup ()
 Clean-up routines.
 
- Protected Member Functions inherited from TFwdSolver< T >
void Setup ()
 
void SetupType ()
 
void DeleteType ()
 
RVector UnfoldComplex (const TVector< T > &vec) const
 Unfold real and imaginary parts of a complex vector. More...
 
FVector UnfoldSComplex (const TVector< T > &vec) const
 

Additional Inherited Members

- Public Attributes inherited from TFwdSolver< T >
const QMMeshmeshptr
 pointer to the associated FEM mesh
 
LSOLVER solvertp
 linear solver type
 
IterativeMethod method
 iterative solver method, if applicable
 
TCompRowMatrix< T > * F
 FEM system matrix.
 
TCompRowMatrix< T > * FL
 lower triangle of system matrix decomposition
 
TVector< T > * Fd
 diagonal of Cholesky factorisation
 
TPreconditioner< T > * precon
 preconditioner instance
 
TCompRowMatrix< T > * B
 mass matrix; only used for time-domain problems
 
void * SuperLU
 SuperLU solver engine.
 
double iterative_tol
 iterative solver tolerance
 
int iterative_maxit
 iterative solver max iterations (0 for auto)
 
PreconType precontp
 preconditioner
 
- Protected Attributes inherited from TFwdSolver< T >
DataScale dscale
 default data scaling: DATA_LIN or DATA_LOG
 
TVector< T > * pphi
 work buffer for field calculation
 
bool unwrap_phase
 use phase unwrapping?
 

Detailed Description

template<class T>
class TFwdSolverMW< T >

Templated forward solver class for multi-wavelength problems.

TFwdSolverMW is an extension to TFwdSolver which allows to calculate the boundary projections for multiple wavelengths with a single call to ProjectAll_wavel.

Constructor & Destructor Documentation

template<class T >
TFwdSolverMW< T >::TFwdSolverMW ( const QMMesh mesh,
LSOLVER  linsolver,
double  tol = 1e-10 
)

Constructor. Creates a forward solver instance.

Parameters
meshpointer to associated QMMesh object
linsolverlinear solver type (see Linear solver methods)
tollinear solver tolerance (only used for iterative methods)
template<class T >
TFwdSolverMW< T >::TFwdSolverMW ( const QMMesh mesh,
char *  solver,
double  tol = 1e-10 
)

Constructor. Creates a forward solver instance.

Parameters
meshpointer to associated QMMesh object
solverlinear solver type, provided as a string (see notes)
tollinear solver tolerance (only used for iterative methods)
Note
For valid strings to define the linear solver, see SetLinSolver.
template<class T >
TFwdSolverMW< T >::TFwdSolverMW ( const QMMesh mesh,
ParamParser pp 
)

Constructor. Creates a forward solver instance.

Parameters
meshpointer to associated QMMesh object
ppparser to read solver options from.
Note
For recognised items in the configuration file used by the parser, see the ReadParams method.

Member Function Documentation

template<class T >
TVector<T> TFwdSolverMW< T >::ProjectAll_wavel ( const TCompRowMatrix< T > &  qvec,
const TCompRowMatrix< T > &  mvec,
const MWsolution sol,
double  omega,
DataScale  scl = DATA_DEFAULT 
)

Return boundary data for all sources and all detectors at all wavelengths.

Parameters
qvecarray of source column vectors
mvecarray of measurement column vectors
solmulti-wavelength solution instance
omegamodulation frequency (cycles/ps)
scloutput data scaling: linear or logarithmic
Returns
data vector (size nw*nQM, with nw=number of wavelengths)
Note
This method sequentially resets the system matrix for all wavelengths according to the optical parameters stored in sol, then calculates the fields and projections.
The returned vector is arranged in blocks of data for each wavelength. Each block is structured in the same way as the vector returned by TFwdSolver::ProjectAll.
template<class T >
RVector TFwdSolverMW< T >::ProjectAll_wavel_real ( const TCompRowMatrix< T > &  qvec,
const TCompRowMatrix< T > &  mvec,
const MWsolution sol,
double  omega,
DataScale  scl = DATA_DEFAULT 
)

Return boundary data for the complex case in a real vector.

Parameters
qvecarray of source column vectors
mvecarray of measurement column vectors
solmulti-wavelength solution instance
omegamodulation frequency (cycles/ps)
scloutput data scaling: linear or logarithmic
Returns
data vector (size nw*nQM*2, with nw=number of wavelengths)
Note
This method sequentially resets the system matrix for all wavelengths according to the optical parameters stored in sol, then calculates the fields and projections.
The returned vector consists of a real and imaginary block. Both blocks consist of sub-blocks for each wavelength, and each sub-block contains further blocks for each source, consisting of the measurements for that source.

The documentation for this class was generated from the following file: