MOOSE - Multiscale Object Oriented Simulation Environment
|
#include <Ksolve.h>
Public Member Functions | |
void | getBlock (vector< double > &values) const |
double | getDiffConst (const Eref &e) const |
Diffusion constant: Only one per pool, voxel number is ignored. More... | |
Id | getDsolve () const |
Inherited from ZombiePoolInterface. More... | |
double | getEpsAbs () const |
Assigns Absolute tolerance for integration. More... | |
double | getEpsRel () const |
Assigns Relative tolerance for integration. More... | |
double | getEstimatedDt () const |
string | getMethod () const |
Assigns integration method. More... | |
double | getN (const Eref &e) const |
Get # of molecules in given pool and voxel. Varies with time. More... | |
double | getNinit (const Eref &e) const |
get initial # of molecules in given pool and voxel. Bdry cond. More... | |
unsigned int | getNumAllVoxels () const |
unsigned int | getNumLocalVoxels () const |
Inherited from ZombiePoolInterface. More... | |
unsigned int | getNumPools () const |
gets number of pools (species) handled by system. More... | |
vector< double > | getNvec (unsigned int voxel) const |
Returns the vector of pool Num at the specified voxel. More... | |
unsigned int | getPoolIndex (const Eref &e) const |
Return pool index, using Stoich ptr to do lookup. More... | |
Id | getStoich () const |
Assigns Stoich object to Ksolve. More... | |
unsigned int | getVoxelIndex (const Eref &e) const |
void | initProc (const Eref &e, ProcPtr p) |
void | initReinit (const Eref &e, ProcPtr p) |
Ksolve () | |
void | matchJunctionVols (vector< double > &vols, Id otherCompt) const |
VoxelPoolsBase * | pools (unsigned int i) |
Return a pointer to the specified VoxelPool. More... | |
void | print () const |
void | process (const Eref &e, ProcPtr p) |
void | reinit (const Eref &e, ProcPtr p) |
void | setBlock (const vector< double > &values) |
void | setDiffConst (const Eref &e, double v) |
Diffusion constant: Only one per pool, voxel number is ignored. More... | |
void | setDsolve (Id dsolve) |
Assigns the diffusion solver. Used by the reac solvers. More... | |
void | setEpsAbs (double val) |
void | setEpsRel (double val) |
void | setMethod (string method) |
void | setN (const Eref &e, double v) |
Set # of molecules in given pool and voxel. Varies with time. More... | |
void | setNinit (const Eref &e, double v) |
Set initial # of molecules in given pool and voxel. Bdry cond. More... | |
void | setNumAllVoxels (unsigned int num) |
void | setNumPools (unsigned int num) |
void | setNvec (unsigned int voxel, vector< double > vec) |
void | setStoich (Id stoich) |
void | updateRateTerms (unsigned int index) |
void | updateVoxelVol (vector< double > vols) |
double | volume (unsigned int i) const |
Return volume of voxel i. More... | |
~Ksolve () | |
![]() | |
Id | getCompartment () const |
virtual void | setCompartment (Id compartment) |
Assigns compartment. More... | |
virtual void | setMotorConst (const Eref &e, double val) |
virtual void | setPrev () |
Used to tell Dsolver to assign 'prev' values. More... | |
virtual void | updateJunctions (double dt) |
Used for telling Dsolver to handle all ops across Junctions. More... | |
ZombiePoolInterface () | |
Static Public Member Functions | |
static const Cinfo * | initCinfo () |
Private Attributes | |
Id | dsolve_ |
ZombiePoolInterface * | dsolvePtr_ |
Pointer to diffusion solver. More... | |
double | epsAbs_ |
double | epsRel_ |
string | method_ |
vector< VoxelPools > | pools_ |
unsigned int | startVoxel_ |
First voxel indexed on the current node. More... | |
Stoich * | stoichPtr_ |
Utility ptr used to help Pool Id lookups by the Ksolve. More... | |
Additional Inherited Members | |
![]() | |
Id | compartment_ |
Id of Chem compartment used to figure out volumes of voxels. More... | |
bool | isBuilt_ |
Flag: True when solver setup has been completed. More... | |
Id | stoich_ |
Ksolve::Ksolve | ( | ) |
Definition at line 216 of file Ksolve.cpp.
Referenced by initCinfo().
Ksolve::~Ksolve | ( | ) |
Definition at line 236 of file Ksolve.cpp.
|
virtual |
Gets block of data. The first 4 entries are passed in on the 'values' vector: the start voxel, numVoxels, start pool#, numPools. These are followed by numVoxels * numPools of data values which are filled in by the function. We assert that the entire requested block is present in this ZombiePoolInterface. The block is organized as an array of arrays of voxels; values[pool#][voxel#]
Note that numVoxels and numPools are the number in the current block, not the upper limit of the block. So values.size() == 4 + numPools * numVoxels.
Implements ZombiePoolInterface.
Definition at line 753 of file Ksolve.cpp.
References pools_, and startVoxel_.
Referenced by process().
|
virtual |
Diffusion constant: Only one per pool, voxel number is ignored.
Implements ZombiePoolInterface.
Definition at line 718 of file Ksolve.cpp.
Id Ksolve::getDsolve | ( | ) | const |
Inherited from ZombiePoolInterface.
Assigns Dsolve object to Ksolve.
Definition at line 401 of file Ksolve.cpp.
References dsolve_.
double Ksolve::getEpsAbs | ( | ) | const |
Assigns Absolute tolerance for integration.
Definition at line 274 of file Ksolve.cpp.
References epsAbs_.
Referenced by initCinfo().
double Ksolve::getEpsRel | ( | ) | const |
Assigns Relative tolerance for integration.
Definition at line 288 of file Ksolve.cpp.
References epsRel_.
Referenced by initCinfo().
double Ksolve::getEstimatedDt | ( | ) | const |
This does a quick and dirty estimate of the timestep suitable for this sytem
Definition at line 474 of file Ksolve.cpp.
References EPSILON, Stoich::getNumAllPools(), Stoich::getNumRates(), pools_, and stoichPtr_.
Referenced by initCinfo().
string Ksolve::getMethod | ( | ) | const |
Assigns integration method.
Definition at line 245 of file Ksolve.cpp.
References method_.
Referenced by initCinfo().
|
virtual |
Get # of molecules in given pool and voxel. Varies with time.
Implements ZombiePoolInterface.
Definition at line 690 of file Ksolve.cpp.
References getPoolIndex(), getVoxelIndex(), OFFNODE, and pools_.
|
virtual |
get initial # of molecules in given pool and voxel. Bdry cond.
Implements ZombiePoolInterface.
Definition at line 705 of file Ksolve.cpp.
References getPoolIndex(), getVoxelIndex(), OFFNODE, and pools_.
unsigned int Ksolve::getNumAllVoxels | ( | ) | const |
Definition at line 432 of file Ksolve.cpp.
References pools_.
Referenced by initCinfo().
|
virtual |
Inherited from ZombiePoolInterface.
Implements ZombiePoolInterface.
Definition at line 427 of file Ksolve.cpp.
References pools_.
Referenced by initCinfo(), and process().
|
virtual |
gets number of pools (species) handled by system.
Implements ZombiePoolInterface.
Definition at line 732 of file Ksolve.cpp.
References pools_.
Referenced by initCinfo().
vector< double > Ksolve::getNvec | ( | unsigned int | voxel | ) | const |
Returns the vector of pool Num at the specified voxel.
Definition at line 447 of file Ksolve.cpp.
References dummy, pools_, and VoxelPoolsBase::Svec().
Referenced by initCinfo().
|
virtual |
Return pool index, using Stoich ptr to do lookup.
Implements ZombiePoolInterface.
Definition at line 666 of file Ksolve.cpp.
References Stoich::convertIdToPoolIndex(), Eref::id(), and stoichPtr_.
Referenced by getN(), getNinit(), setN(), and setNinit().
Id Ksolve::getStoich | ( | ) | const |
Assigns Stoich object to Ksolve.
Definition at line 317 of file Ksolve.cpp.
References ZombiePoolInterface::stoich_.
Referenced by initCinfo().
unsigned int Ksolve::getVoxelIndex | ( | const Eref & | e | ) | const |
Definition at line 671 of file Ksolve.cpp.
References Eref::dataIndex(), OFFNODE, pools_, and startVoxel_.
Referenced by getN(), getNinit(), setN(), and setNinit().
|
static |
Definition at line 40 of file Ksolve.cpp.
References ZombiePoolInterface::getCompartment(), getEpsAbs(), getEpsRel(), getEstimatedDt(), getMethod(), getNumAllVoxels(), getNumLocalVoxels(), getNumPools(), getNvec(), getStoich(), init(), Neutral::initCinfo(), initProc(), initReinit(), Ksolve(), ksolveCinfo, process(), reinit(), ZombiePoolInterface::setCompartment(), setEpsAbs(), setEpsRel(), setMethod(), setNumAllVoxels(), setNumPools(), setNvec(), and updateVoxelVol().
Definition at line 626 of file Ksolve.cpp.
Referenced by initCinfo().
Definition at line 630 of file Ksolve.cpp.
References ProcInfo::dt, pools_, and reinit().
Referenced by initCinfo().
void Ksolve::matchJunctionVols | ( | vector< double > & | vols, |
Id | otherCompt | ||
) | const |
|
virtual |
Return a pointer to the specified VoxelPool.
Implements ZombiePoolInterface.
Definition at line 739 of file Ksolve.cpp.
References pools_.
void Ksolve::print | ( | ) | const |
Definition at line 823 of file Ksolve.cpp.
References ZombiePoolInterface::compartment_, dsolve_, Stoich::getKsolve(), method_, Id::path(), pools_, ZombiePoolInterface::stoich_, and stoichPtr_.
Definition at line 497 of file Ksolve.cpp.
References advance(), dsolvePtr_, ProcInfo::dt, ZombiePoolInterface::getBlock(), getBlock(), getNumLocalVoxels(), Stoich::getNumVarPools(), ZombiePoolInterface::isBuilt_, pools_, ZombiePoolInterface::setBlock(), setBlock(), ZombiePoolInterface::setPrev(), stoichPtr_, and ZombiePoolInterface::updateJunctions().
Referenced by initCinfo().
Definition at line 600 of file Ksolve.cpp.
References ProcInfo::dt, ZombiePoolInterface::isBuilt_, pools_, and stoichPtr_.
Referenced by initCinfo(), and initReinit().
|
virtual |
Sets block of data. The first 4 entries on the 'values' vector are the start voxel, numVoxels, start pool#, numPools. These are followed by numVoxels * numPools of data values.
Implements ZombiePoolInterface.
Definition at line 776 of file Ksolve.cpp.
References pools_, and startVoxel_.
Referenced by process().
|
virtual |
Diffusion constant: Only one per pool, voxel number is ignored.
Implements ZombiePoolInterface.
Definition at line 713 of file Ksolve.cpp.
|
virtual |
Assigns the diffusion solver. Used by the reac solvers.
Implements ZombiePoolInterface.
Definition at line 406 of file Ksolve.cpp.
References Element::cinfo(), Eref::data(), dsolve_, dsolvePtr_, Id::element(), Id::eref(), Cinfo::isA(), Cinfo::name(), and Id::path().
void Ksolve::setEpsAbs | ( | double | val | ) |
Definition at line 279 of file Ksolve.cpp.
References epsAbs_.
Referenced by initCinfo().
void Ksolve::setEpsRel | ( | double | val | ) |
Definition at line 293 of file Ksolve.cpp.
References epsRel_.
Referenced by initCinfo().
void Ksolve::setMethod | ( | string | method | ) |
Definition at line 250 of file Ksolve.cpp.
References method_.
Referenced by initCinfo().
|
virtual |
Set # of molecules in given pool and voxel. Varies with time.
Implements ZombiePoolInterface.
Definition at line 683 of file Ksolve.cpp.
References getPoolIndex(), getVoxelIndex(), OFFNODE, and pools_.
|
virtual |
Set initial # of molecules in given pool and voxel. Bdry cond.
Implements ZombiePoolInterface.
Definition at line 698 of file Ksolve.cpp.
References getPoolIndex(), getVoxelIndex(), OFFNODE, and pools_.
|
virtual |
Assigns the number of voxels used in the entire reac-diff system. Note that fewer than this may be used on any given node.
Implements ZombiePoolInterface.
Definition at line 438 of file Ksolve.cpp.
References pools_.
Referenced by initCinfo().
|
virtual |
Assigns number of different pools (chemical species) present in each voxel. Inherited.
Implements ZombiePoolInterface.
Definition at line 723 of file Ksolve.cpp.
References pools_.
Referenced by initCinfo().
void Ksolve::setNvec | ( | unsigned int | voxel, |
vector< double > | vec | ||
) |
Definition at line 457 of file Ksolve.cpp.
References pools_.
Referenced by initCinfo().
|
virtual |
Informs the ZPI about the stoich, used during subsequent computations. Called to wrap up the model building. The Stoich does this call after it has set up its own path.
Implements ZombiePoolInterface.
Definition at line 353 of file Ksolve.cpp.
References Element::cinfo(), Eref::data(), Id::element(), OdeSystem::epsAbs, epsAbs_, OdeSystem::epsRel, epsRel_, Id::eref(), Stoich::getNumAllPools(), OdeSystem::initStepSize, Cinfo::isA(), ZombiePoolInterface::isBuilt_, OdeSystem::method, method_, pools_, ZombiePoolInterface::stoich_, and stoichPtr_.
|
virtual |
Rescale specified voxel rate term following rate constant change or volume change. If index == ~0U then does all terms.
updateRateTerms obtains the latest parameters for the rates_ vector, and has each of the pools update its parameters including rescaling for volumes.
Implements ZombiePoolInterface.
Definition at line 641 of file Ksolve.cpp.
References Stoich::getNumCoreRates(), Stoich::getRateTerms(), pools_, and stoichPtr_.
Referenced by updateVoxelVol().
void Ksolve::updateVoxelVol | ( | vector< double > | vols | ) |
Handles request to change volumes of voxels in this Ksolve, and all cascading effects of this. At this point it won't handle change in size of voxel array.
Definition at line 800 of file Ksolve.cpp.
References pools_, and updateRateTerms().
Referenced by initCinfo().
|
virtual |
Return volume of voxel i.
Implements ZombiePoolInterface.
Definition at line 746 of file Ksolve.cpp.
References pools_.
|
private |
Id of diffusion solver, needed for coordinating numerics.
Definition at line 167 of file Ksolve.h.
Referenced by getDsolve(), print(), and setDsolve().
|
private |
Pointer to diffusion solver.
Definition at line 170 of file Ksolve.h.
Referenced by process(), and setDsolve().
|
private |
Definition at line 139 of file Ksolve.h.
Referenced by getEpsAbs(), setEpsAbs(), and setStoich().
|
private |
Definition at line 140 of file Ksolve.h.
Referenced by getEpsRel(), setEpsRel(), and setStoich().
|
private |
Definition at line 138 of file Ksolve.h.
Referenced by getMethod(), print(), setMethod(), and setStoich().
|
private |
Each VoxelPools entry handles all the pools in a single voxel. Each entry knows how to update itself in order to complete the kinetic calculations for that voxel. The ksolver does multinode management by indexing only the subset of entries present on this node.
Definition at line 156 of file Ksolve.h.
Referenced by getBlock(), getEstimatedDt(), getN(), getNinit(), getNumAllVoxels(), getNumLocalVoxels(), getNumPools(), getNvec(), getVoxelIndex(), initReinit(), pools(), print(), process(), reinit(), setBlock(), setN(), setNinit(), setNumAllVoxels(), setNumPools(), setNvec(), setStoich(), updateRateTerms(), updateVoxelVol(), and volume().
|
private |
First voxel indexed on the current node.
Definition at line 159 of file Ksolve.h.
Referenced by getBlock(), getVoxelIndex(), and setBlock().
|
private |
Utility ptr used to help Pool Id lookups by the Ksolve.
Definition at line 162 of file Ksolve.h.
Referenced by getEstimatedDt(), getPoolIndex(), print(), process(), reinit(), setStoich(), and updateRateTerms().