MURE
Loading...
Searching...
No Matches
Public Member Functions | List of all members
CRAMSolver Class Reference

Base class of Bateman-type solvers like Runge-Kutta or CRAM or any other. More...

#include <CRAMSolver.hxx>

Inheritance diagram for CRAMSolver:
BatemanSolver

Public Member Functions

 CRAMSolver (string type="IPF", int order=48, bool aut=true)
 Normal constructor.
 
Public methods
void Solve (double *N0, double **Matrix, int N, double t1, double t2) override
 Abstract method to Solve Bateman equation : must be overriden.
 
void SetTheMatrix (double **f)
 
void printCoeffs (int order=0)
 
void SetNSubsteps (int n)
 
void SetOrder (int n)
 
- Public Member Functions inherited from BatemanSolver
 BatemanSolver ()
 Normal constructor.
 
 BatemanSolver (const BatemanSolver &BS)
 Copy constructor.
 
virtual ~BatemanSolver ()
 Destructor.
 
int GetNumberOfEquationSize ()
 
void SetNumberOfEquations (int N)
 
void SetPrecision (double eps=1e-5)
 
void SetForbidNegativeValue ()
 Forbid negative value during integration.
 
void SetXe135EquilibriumIndex (int Te135, int I135, int Xe135)
 Set index (position) in EvolvingSystem::fEvolvingAtoms vector for Te-135, I-135 and Xe-135.
 

Private methods

t_SpMat_SolverVector CRAMSpMatSolver
 
t_CRAMCoeff CRAMCoefficients
 
string CRAMType
 
string SparseSolver
 
int CRAMOrder
 
double ** fTheMatrix
 The evolution Matrix.
 
double CRAMDeltaT
 
double CRAMDeltaTsubstep
 
int nSubsteps
 
Timer HPTimer
 
const double cutoff = 1e-20
 
const double significantDigits = 10
 
const int maxorder = 48
 
bool automatic = false
 
bool fBuildEqns
 
static const vector< int > fCramOrder
 
static const vector< t_CRAMCoeffCRAMCoeff
 
void BuildEqns (double *N0, double **Matrix)
 Setup Eigen3 solver, compute inverse matrices.
 
t_CRAMCoeff GetCRAMCoefficients (int order=0)
 
t_SpMat_SolverVector GetSparseLUSolvers ()
 
int GetNumberSubsteps (double &N0_i, double &A_ii, double dT, int order, double sigDigits)
 
double GetSignificantDigits (double &N0_i, double &A_ii, double dT, int order, int niter)
 
int GetNeededCRAMOrder (double &N0_i, double &A_ii, double dT, double sigDigits)
 
void SetOptimalStrategy (double *N0, double **fMatrix)
 
int GetNearestCRAMOrder (int order)
 
double GetLG10 (double &n0, double &A_ii, double dT, int nT)
 
void CleanUp ()
 

Additional Inherited Members

- Protected Attributes inherited from BatemanSolver
int fNVar
 The size of the composition vector and /or number of ZAIs involved.
 
double fPrecision
 Precision of the solution.
 
double fTmin
 The lowest value of the decay constant considered.
 
bool fIsNegativeValueAllowed
 
int fIdxTe135
 index of Te-135 in the fEvolvingAtoms array
 
int fIdxI135
 index of I-135 in the fEvolvingAtoms array
 
int fIdxXe135
 index of Xe-135 in the fEvolvingAtoms array
 

Detailed Description

Base class of Bateman-type solvers like Runge-Kutta or CRAM or any other.

This class only holds common variables and defines template functions for inheriting classes

Author
kernie
Version
1.0

Constructor & Destructor Documentation

◆ CRAMSolver()

CRAMSolver::CRAMSolver ( string  type = "IPF",
int  order = 48,
bool  aut = true 
)

Normal constructor.

Member Function Documentation

◆ BuildEqns()

void CRAMSolver::BuildEqns ( double *  N0,
double **  Matrix 
)
private

Setup Eigen3 solver, compute inverse matrices.

◆ CleanUp()

void CRAMSolver::CleanUp ( )
private

◆ GetCRAMCoefficients()

t_CRAMCoeff CRAMSolver::GetCRAMCoefficients ( int  order = 0)
private

◆ GetLG10()

double CRAMSolver::GetLG10 ( double &  n0,
double &  A_ii,
double  dT,
int  nT 
)
private

◆ GetNearestCRAMOrder()

int CRAMSolver::GetNearestCRAMOrder ( int  order)
private

◆ GetNeededCRAMOrder()

int CRAMSolver::GetNeededCRAMOrder ( double &  N0_i,
double &  A_ii,
double  dT,
double  sigDigits 
)
private

◆ GetNumberSubsteps()

int CRAMSolver::GetNumberSubsteps ( double &  N0_i,
double &  A_ii,
double  dT,
int  order,
double  sigDigits 
)
private

◆ GetSignificantDigits()

double CRAMSolver::GetSignificantDigits ( double &  N0_i,
double &  A_ii,
double  dT,
int  order,
int  niter 
)
private

◆ GetSparseLUSolvers()

t_SpMat_SolverVector CRAMSolver::GetSparseLUSolvers ( )
private

◆ printCoeffs()

void CRAMSolver::printCoeffs ( int  order = 0)

◆ SetNSubsteps()

void CRAMSolver::SetNSubsteps ( int  n)
inline

◆ SetOptimalStrategy()

void CRAMSolver::SetOptimalStrategy ( double *  N0,
double **  fMatrix 
)
private

◆ SetOrder()

void CRAMSolver::SetOrder ( int  n)
inline

◆ SetTheMatrix()

void CRAMSolver::SetTheMatrix ( double **  f)
inline

◆ Solve()

void CRAMSolver::Solve ( double *  N0,
double **  Matrix,
int  N,
double  t1,
double  t2 
)
overridevirtual

Abstract method to Solve Bateman equation : must be overriden.

For a given vector

\[\vec{N_0}\]

of initial composition and a given matrix NxN, integrate Bateman equation from t1 to t2.

Parameters
N0: initial composition vector
Matrix: matrix of Bateman eqs
N: initial composition vector
t1: initial time
t2: final time

Implements BatemanSolver.

Member Data Documentation

◆ automatic

bool CRAMSolver::automatic = false
private

◆ CRAMCoeff

const vector< t_CRAMCoeff > CRAMSolver::CRAMCoeff
staticprivate

◆ CRAMCoefficients

t_CRAMCoeff CRAMSolver::CRAMCoefficients
private

◆ CRAMDeltaT

double CRAMSolver::CRAMDeltaT
private

◆ CRAMDeltaTsubstep

double CRAMSolver::CRAMDeltaTsubstep
private

◆ CRAMOrder

int CRAMSolver::CRAMOrder
private

◆ CRAMSpMatSolver

t_SpMat_SolverVector CRAMSolver::CRAMSpMatSolver
private

◆ CRAMType

string CRAMSolver::CRAMType
private

◆ cutoff

const double CRAMSolver::cutoff = 1e-20
private

◆ fBuildEqns

bool CRAMSolver::fBuildEqns
private

◆ fCramOrder

const vector< int > CRAMSolver::fCramOrder
staticprivate

◆ fTheMatrix

double** CRAMSolver::fTheMatrix
private

The evolution Matrix.

◆ HPTimer

Timer CRAMSolver::HPTimer
private

◆ maxorder

const int CRAMSolver::maxorder = 48
private

◆ nSubsteps

int CRAMSolver::nSubsteps
private

◆ significantDigits

const double CRAMSolver::significantDigits = 10
private

◆ SparseSolver

string CRAMSolver::SparseSolver
private

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

MURE Project, documentation generated by Doxygen 1.9.7 - Fri Jan 19 2024