MOOSE - Multiscale Object Oriented Simulation Environment
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
Ksolve Class Reference

#include <Ksolve.h>

+ Inheritance diagram for Ksolve:
+ Collaboration diagram for Ksolve:

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
 
VoxelPoolsBasepools (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 ()
 
- Public Member Functions inherited from ZombiePoolInterface
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 CinfoinitCinfo ()
 

Private Attributes

Id dsolve_
 
ZombiePoolInterfacedsolvePtr_
 Pointer to diffusion solver. More...
 
double epsAbs_
 
double epsRel_
 
string method_
 
vector< VoxelPoolspools_
 
unsigned int startVoxel_
 First voxel indexed on the current node. More...
 
StoichstoichPtr_
 Utility ptr used to help Pool Id lookups by the Ksolve. More...
 

Additional Inherited Members

- Protected Attributes inherited from ZombiePoolInterface
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_
 

Detailed Description

Definition at line 15 of file Ksolve.h.

Constructor & Destructor Documentation

Ksolve::Ksolve ( )

Definition at line 216 of file Ksolve.cpp.

Referenced by initCinfo().

217  :
218 #if USE_GSL
219  method_( "rk5" ),
220 #elif USE_BOOST_ODE
221  method_( "rk5a" ),
222 #endif
223  epsAbs_( 1e-7 ),
224  epsRel_( 1e-7 ),
225 #if PARALLELIZE_KSOLVE_WITH_CPP11_ASYNC
226  numThreads_( 3 ),
227 #endif
228  pools_( 1 ),
229  startVoxel_( 0 ),
230  dsolve_(),
231  dsolvePtr_( 0 )
232 {
233  ;
234 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
Id dsolve_
Definition: Ksolve.h:167
string method_
Definition: Ksolve.h:138
ZombiePoolInterface * dsolvePtr_
Pointer to diffusion solver.
Definition: Ksolve.h:170
unsigned int startVoxel_
First voxel indexed on the current node.
Definition: Ksolve.h:159
double epsRel_
Definition: Ksolve.h:140
double epsAbs_
Definition: Ksolve.h:139

+ Here is the caller graph for this function:

Ksolve::~Ksolve ( )

Definition at line 236 of file Ksolve.cpp.

237 {
238  ;
239 }

Member Function Documentation

void Ksolve::getBlock ( vector< double > &  values) const
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().

754 {
755  unsigned int startVoxel = values[0];
756  unsigned int numVoxels = values[1];
757  unsigned int startPool = values[2];
758  unsigned int numPools = values[3];
759 
760  assert( startVoxel >= startVoxel_ );
761  assert( numVoxels <= pools_.size() );
762  assert( pools_.size() > 0 );
763  assert( numPools + startPool <= pools_[0].size() );
764  values.resize( 4 + numVoxels * numPools );
765 
766  for ( unsigned int i = 0; i < numVoxels; ++i )
767  {
768  const double* v = pools_[ startVoxel + i ].S();
769  for ( unsigned int j = 0; j < numPools; ++j )
770  {
771  values[ 4 + j * numVoxels + i] = v[ j + startPool ];
772  }
773  }
774 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
unsigned int startVoxel_
First voxel indexed on the current node.
Definition: Ksolve.h:159

+ Here is the caller graph for this function:

double Ksolve::getDiffConst ( const Eref e) const
virtual

Diffusion constant: Only one per pool, voxel number is ignored.

Implements ZombiePoolInterface.

Definition at line 718 of file Ksolve.cpp.

719 {
720  return 0;
721 }
Id Ksolve::getDsolve ( ) const

Inherited from ZombiePoolInterface.

Assigns Dsolve object to Ksolve.

Definition at line 401 of file Ksolve.cpp.

References dsolve_.

402 {
403  return dsolve_;
404 }
Id dsolve_
Definition: Ksolve.h:167
double Ksolve::getEpsAbs ( ) const

Assigns Absolute tolerance for integration.

Definition at line 274 of file Ksolve.cpp.

References epsAbs_.

Referenced by initCinfo().

275 {
276  return epsAbs_;
277 }
double epsAbs_
Definition: Ksolve.h:139

+ Here is the caller graph for this function:

double Ksolve::getEpsRel ( ) const

Assigns Relative tolerance for integration.

Definition at line 288 of file Ksolve.cpp.

References epsRel_.

Referenced by initCinfo().

289 {
290  return epsRel_;
291 }
double epsRel_
Definition: Ksolve.h:140

+ Here is the caller graph for this function:

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().

475 {
476  static const double EPSILON = 1e-15;
477  vector< double > s( stoichPtr_->getNumAllPools(), 1.0 );
478  vector< double > v( stoichPtr_->getNumRates(), 0.0 );
479  double maxVel = 0.0;
480  if ( pools_.size() > 0.0 )
481  {
482  pools_[0].updateReacVelocities( &s[0], v );
483  for ( vector< double >::iterator
484  i = v.begin(); i != v.end(); ++i )
485  if ( maxVel < *i )
486  maxVel = *i;
487  }
488  if ( maxVel < EPSILON )
489  return 0.1; // Based on typical sig pathway reac rates.
490  // Heuristic: the largest velocity times dt should be 10% of mol conc.
491  return 0.1 / maxVel;
492 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
Stoich * stoichPtr_
Utility ptr used to help Pool Id lookups by the Ksolve.
Definition: Ksolve.h:162
unsigned int getNumRates() const
Definition: Stoich.cpp:568
unsigned int getNumAllPools() const
Definition: Stoich.cpp:520
#define EPSILON
Definition: MatrixOps.h:28

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

string Ksolve::getMethod ( ) const

Assigns integration method.

Definition at line 245 of file Ksolve.cpp.

References method_.

Referenced by initCinfo().

246 {
247  return method_;
248 }
string method_
Definition: Ksolve.h:138

+ Here is the caller graph for this function:

double Ksolve::getN ( const Eref e) const
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_.

691 {
692  unsigned int vox = getVoxelIndex( e );
693  if ( vox != OFFNODE )
694  return pools_[vox].getN( getPoolIndex( e ) );
695  return 0.0;
696 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
unsigned int getPoolIndex(const Eref &e) const
Return pool index, using Stoich ptr to do lookup.
Definition: Ksolve.cpp:666
const unsigned int OFFNODE
Definition: Ksolve.cpp:38
unsigned int getVoxelIndex(const Eref &e) const
Definition: Ksolve.cpp:671

+ Here is the call graph for this function:

double Ksolve::getNinit ( const Eref e) const
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_.

706 {
707  unsigned int vox = getVoxelIndex( e );
708  if ( vox != OFFNODE )
709  return pools_[vox].getNinit( getPoolIndex( e ) );
710  return 0.0;
711 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
unsigned int getPoolIndex(const Eref &e) const
Return pool index, using Stoich ptr to do lookup.
Definition: Ksolve.cpp:666
const unsigned int OFFNODE
Definition: Ksolve.cpp:38
unsigned int getVoxelIndex(const Eref &e) const
Definition: Ksolve.cpp:671

+ Here is the call graph for this function:

unsigned int Ksolve::getNumAllVoxels ( ) const

Definition at line 432 of file Ksolve.cpp.

References pools_.

Referenced by initCinfo().

433 {
434  return pools_.size(); // Need to redo.
435 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156

+ Here is the caller graph for this function:

unsigned int Ksolve::getNumLocalVoxels ( ) const
virtual

Inherited from ZombiePoolInterface.

Implements ZombiePoolInterface.

Definition at line 427 of file Ksolve.cpp.

References pools_.

Referenced by initCinfo(), and process().

428 {
429  return pools_.size();
430 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156

+ Here is the caller graph for this function:

unsigned int Ksolve::getNumPools ( ) const
virtual

gets number of pools (species) handled by system.

Implements ZombiePoolInterface.

Definition at line 732 of file Ksolve.cpp.

References pools_.

Referenced by initCinfo().

733 {
734  if ( pools_.size() > 0 )
735  return pools_[0].size();
736  return 0;
737 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156

+ Here is the caller graph for this function:

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().

448 {
449  static vector< double > dummy;
450  if ( voxel < pools_.size() )
451  {
452  return const_cast< VoxelPools* >( &( pools_[ voxel ] ) )->Svec();
453  }
454  return dummy;
455 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
static DestFinfo dummy("dummy","This Finfo is a dummy. If you are reading this you have used an invalid index", 0)
vector< double > & Svec()
Returns a handle to the mol # vector.

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

unsigned int Ksolve::getPoolIndex ( const Eref er) const
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().

667 {
668  return stoichPtr_->convertIdToPoolIndex( e.id() );
669 }
Stoich * stoichPtr_
Utility ptr used to help Pool Id lookups by the Ksolve.
Definition: Ksolve.h:162
unsigned int convertIdToPoolIndex(Id id) const
Definition: Stoich.cpp:1464

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

Id Ksolve::getStoich ( ) const

Assigns Stoich object to Ksolve.

Definition at line 317 of file Ksolve.cpp.

References ZombiePoolInterface::stoich_.

Referenced by initCinfo().

318 {
319  return stoich_;
320 }

+ Here is the caller graph for this function:

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().

672 {
673  unsigned int ret = e.dataIndex();
674  if ( ret < startVoxel_ || ret >= startVoxel_ + pools_.size() )
675  return OFFNODE;
676  return ret - startVoxel_;
677 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
unsigned int dataIndex() const
Definition: Eref.h:50
unsigned int startVoxel_
First voxel indexed on the current node.
Definition: Ksolve.h:159
const unsigned int OFFNODE
Definition: Ksolve.cpp:38

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

const Cinfo * Ksolve::initCinfo ( )
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().

41 {
43  // Field definitions
45 
46  static ValueFinfo< Ksolve, string > method (
47  "method",
48  "Integration method, using GSL. So far only explict. Options are:"
49  "rk5: The default Runge-Kutta-Fehlberg 5th order adaptive dt method"
50  "gsl: alias for the above"
51  "rk4: The Runge-Kutta 4th order fixed dt method"
52  "rk2: The Runge-Kutta 2,3 embedded fixed dt method"
53  "rkck: The Runge-Kutta Cash-Karp (4,5) method"
54  "rk8: The Runge-Kutta Prince-Dormand (8,9) method" ,
57  );
58 
59  static ValueFinfo< Ksolve, double > epsAbs (
60  "epsAbs",
61  "Absolute permissible integration error range.",
64  );
65 
66  static ValueFinfo< Ksolve, double > epsRel (
67  "epsRel",
68  "Relative permissible integration error range.",
71  );
72 
73 
74  static ValueFinfo< Ksolve, Id > compartment(
75  "compartment",
76  "Compartment in which the Ksolve reaction system lives.",
79  );
80 
81  static ReadOnlyValueFinfo< Ksolve, unsigned int > numLocalVoxels(
82  "numLocalVoxels",
83  "Number of voxels in the core reac-diff system, on the "
84  "current solver. ",
86  );
87  static LookupValueFinfo<
88  Ksolve, unsigned int, vector< double > > nVec(
89  "nVec",
90  "vector of pool counts. Index specifies which voxel.",
93  );
94  static ValueFinfo< Ksolve, unsigned int > numAllVoxels(
95  "numAllVoxels",
96  "Number of voxels in the entire reac-diff system, "
97  "including proxy voxels to represent abutting compartments.",
100  );
101 
102 #if PARALLELIZE_KSOLVE_WITH_CPP11_ASYNC
103  static ValueFinfo< Ksolve, unsigned int > numThreads (
104  "numThreads",
105  "Number of threads to use (applicable in deterministic case)",
106  &Ksolve::setNumThreads,
107  &Ksolve::getNumThreads
108  );
109 #endif
110 
111  static ValueFinfo< Ksolve, unsigned int > numPools(
112  "numPools",
113  "Number of molecular pools in the entire reac-diff system, "
114  "including variable, function and buffered.",
117  );
118 
119  static ReadOnlyValueFinfo< Ksolve, double > estimatedDt(
120  "estimatedDt",
121  "Estimated timestep for reac system based on Euler error",
123  );
124  static ReadOnlyValueFinfo< Ksolve, Id > stoich(
125  "stoich",
126  "Id for stoichiometry object tied to this Ksolve",
128  );
129 
130 
132  // DestFinfo definitions
134 
135  static DestFinfo process( "process",
136  "Handles process call from Clock",
138  static DestFinfo reinit( "reinit",
139  "Handles reinit call from Clock",
141 
142  static DestFinfo initProc( "initProc",
143  "Handles initProc call from Clock",
145  static DestFinfo initReinit( "initReinit",
146  "Handles initReinit call from Clock",
148 
149  static DestFinfo voxelVol( "voxelVol",
150  "Handles updates to all voxels. Comes from parent "
151  "ChemCompt object.",
152  new OpFunc1< Ksolve, vector< double > >(
154  );
156  // Shared definitions
158  static Finfo* procShared[] =
159  {
160  &process, &reinit
161  };
162  static SharedFinfo proc( "proc",
163  "Shared message for process and reinit. These are used for "
164  "all regular Ksolve calculations including interfacing with "
165  "the diffusion calculations by a Dsolve.",
166  procShared, sizeof( procShared ) / sizeof( const Finfo* )
167  );
168  static Finfo* initShared[] =
169  {
171  };
172  static SharedFinfo init( "init",
173  "Shared message for initProc and initReinit. This is used"
174  " when the system has cross-compartment reactions. ",
175  initShared, sizeof( initShared ) / sizeof( const Finfo* )
176  );
177 
178  static Finfo* ksolveFinfos[] =
179  {
180  &method, // Value
181  &epsAbs, // Value
182  &epsRel , // Value
183 #if PARALLELIZE_KSOLVE_WITH_CPP11_ASYNC
184  &numThreads, // Value
185 #endif
186  &compartment, // Value
187  &numLocalVoxels, // ReadOnlyValue
188  &nVec, // LookupValue
189  &numAllVoxels, // ReadOnlyValue
190  &numPools, // Value
191  &estimatedDt, // ReadOnlyValue
192  &stoich, // ReadOnlyValue
193  &voxelVol, // DestFinfo
194  &proc, // SharedFinfo
195  &init, // SharedFinfo
196  };
197 
198  static Dinfo< Ksolve > dinfo;
199  static Cinfo ksolveCinfo(
200  "Ksolve",
202  ksolveFinfos,
203  sizeof(ksolveFinfos)/sizeof(Finfo *),
204  &dinfo
205  );
206 
207  return &ksolveCinfo;
208 }
Id init(int argc, char **argv, bool &doUnitTests, bool &doRegressionTests, unsigned int &benchmark)
Definition: main.cpp:150
void setNumPools(unsigned int num)
Definition: Ksolve.cpp:723
string getMethod() const
Assigns integration method.
Definition: Ksolve.cpp:245
Definition: Dinfo.h:60
void reinit(const Eref &e, ProcPtr p)
Definition: Ksolve.cpp:600
void updateVoxelVol(vector< double > vols)
Definition: Ksolve.cpp:800
static const Cinfo * ksolveCinfo
Definition: Ksolve.cpp:210
virtual void setCompartment(Id compartment)
Assigns compartment.
unsigned int getNumPools() const
gets number of pools (species) handled by system.
Definition: Ksolve.cpp:732
void setEpsRel(double val)
Definition: Ksolve.cpp:293
vector< double > getNvec(unsigned int voxel) const
Returns the vector of pool Num at the specified voxel.
Definition: Ksolve.cpp:447
double getEpsAbs() const
Assigns Absolute tolerance for integration.
Definition: Ksolve.cpp:274
Definition: OpFunc.h:27
void setNumAllVoxels(unsigned int num)
Definition: Ksolve.cpp:438
void process(const Eref &e, ProcPtr p)
Definition: Ksolve.cpp:497
void setMethod(string method)
Definition: Ksolve.cpp:250
unsigned int getNumLocalVoxels() const
Inherited from ZombiePoolInterface.
Definition: Ksolve.cpp:427
Id getStoich() const
Assigns Stoich object to Ksolve.
Definition: Ksolve.cpp:317
void setNvec(unsigned int voxel, vector< double > vec)
Definition: Ksolve.cpp:457
void initReinit(const Eref &e, ProcPtr p)
Definition: Ksolve.cpp:630
static const Cinfo * initCinfo()
Definition: Neutral.cpp:16
void initProc(const Eref &e, ProcPtr p)
Definition: Ksolve.cpp:626
double getEpsRel() const
Assigns Relative tolerance for integration.
Definition: Ksolve.cpp:288
unsigned int getNumAllVoxels() const
Definition: Ksolve.cpp:432
double getEstimatedDt() const
Definition: Ksolve.cpp:474
void setEpsAbs(double val)
Definition: Ksolve.cpp:279
Definition: Cinfo.h:18
Ksolve()
Definition: Ksolve.cpp:216
Definition: Finfo.h:12

+ Here is the call graph for this function:

void Ksolve::initProc ( const Eref e,
ProcPtr  p 
)

Definition at line 626 of file Ksolve.cpp.

Referenced by initCinfo().

627 {
628 }

+ Here is the caller graph for this function:

void Ksolve::initReinit ( const Eref e,
ProcPtr  p 
)

Definition at line 630 of file Ksolve.cpp.

References ProcInfo::dt, pools_, and reinit().

Referenced by initCinfo().

631 {
632  for ( unsigned int i = 0 ; i < pools_.size(); ++i )
633  pools_[i].reinit( p->dt );
634 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
void reinit(const Eref &e, ProcPtr p)
Definition: Ksolve.cpp:600
double dt
Definition: ProcInfo.h:18

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void Ksolve::matchJunctionVols ( vector< double > &  vols,
Id  otherCompt 
) const
VoxelPoolsBase * Ksolve::pools ( unsigned int  i)
virtual

Return a pointer to the specified VoxelPool.

Implements ZombiePoolInterface.

Definition at line 739 of file Ksolve.cpp.

References pools_.

740 {
741  if ( pools_.size() > i )
742  return &pools_[i];
743  return 0;
744 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
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_.

824 {
825  cout << "path = " << stoichPtr_->getKsolve().path() <<
826  ", numPools = " << pools_.size() << "\n";
827  for ( unsigned int i = 0; i < pools_.size(); ++i )
828  {
829  cout << "pools[" << i << "] contents = ";
830  pools_[i].print();
831  }
832  cout << "method = " << method_ << ", stoich=" << stoich_.path() <<endl;
833  cout << "dsolve = " << dsolve_.path() << endl;
834  cout << "compartment = " << compartment_.path() << endl;
835 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
std::string path(const std::string &separator="/") const
Definition: Id.cpp:76
Id dsolve_
Definition: Ksolve.h:167
Id getKsolve() const
Definition: Stoich.cpp:446
Stoich * stoichPtr_
Utility ptr used to help Pool Id lookups by the Ksolve.
Definition: Ksolve.h:162
string method_
Definition: Ksolve.h:138
Id compartment_
Id of Chem compartment used to figure out volumes of voxels.

+ Here is the call graph for this function:

void Ksolve::process ( const Eref e,
ProcPtr  p 
)

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().

498 {
499  if ( isBuilt_ == false )
500  return;
501 
502  // First, handle incoming diffusion values, update S with those.
503  if ( dsolvePtr_ )
504  {
505  vector< double > dvalues( 4 );
506  dvalues[0] = 0;
507  dvalues[1] = getNumLocalVoxels();
508  dvalues[2] = 0;
509  dvalues[3] = stoichPtr_->getNumVarPools();
510 
511  dsolvePtr_->getBlock( dvalues );
512  // Second, set the prev_ value in DiffPoolVec
513  dsolvePtr_->setPrev();
514  setBlock( dvalues );
515  }
516 
517  size_t nvPools = pools_.size( );
518 
519 #ifdef PARALLELIZE_KSOLVE_WITH_CPP11_ASYNC
520  // Third, do the numerical integration for all reactions.
521  size_t grainSize = min( nvPools, 1 + (nvPools / numThreads_ ) );
522  size_t nWorkers = nvPools / grainSize;
523 
524  if( 1 == nWorkers || 1 == nvPools )
525  {
526  if( numThreads_ > 1 )
527  {
528 #ifndef NDEBUG
529  cout << "Debug: Reset to 1 threads " << endl;
530 #endif
531  numThreads_ = 1;
532  }
533 
534  for ( size_t i = 0; i < nvPools; i++ )
535  pools_[i].advance( p );
536  }
537  else
538  {
539  /*-----------------------------------------------------------------------------
540  * Somewhat complicated computation to compute the number of threads. 1
541  * thread per (at least) voxel pool is ideal situation.
542  *-----------------------------------------------------------------------------*/
543  //cout << "Grain size " << grainSize << " Workers : " << nWorkers << endl;
544  for (size_t i = 0; i < nWorkers; i++)
545  parallel_advance( i * grainSize, (i+1) * grainSize, nWorkers, p );
546  }
547 #else
548  for ( size_t i = 0; i < nvPools; i++ )
549  pools_[i].advance( p );
550 #endif
551 
552 
553  // Assemble and send the integrated values off for the Dsolve.
554  if ( dsolvePtr_ )
555  {
556  vector< double > kvalues( 4 );
557  kvalues[0] = 0;
558  kvalues[1] = getNumLocalVoxels();
559  kvalues[2] = 0;
560  kvalues[3] = stoichPtr_->getNumVarPools();
561  getBlock( kvalues );
562  dsolvePtr_->setBlock( kvalues );
563 
564  // Now use the values in the Dsolve to update junction fluxes
565  // for diffusion, channels, and xreacs
567  }
568 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
void getBlock(vector< double > &values) const
Definition: Ksolve.cpp:753
Stoich * stoichPtr_
Utility ptr used to help Pool Id lookups by the Ksolve.
Definition: Ksolve.h:162
virtual void updateJunctions(double dt)
Used for telling Dsolver to handle all ops across Junctions.
void advance(vector< double > &y, const vector< Triplet< double > > &ops, const vector< double > &diagVal)
ZombiePoolInterface * dsolvePtr_
Pointer to diffusion solver.
Definition: Ksolve.h:170
double dt
Definition: ProcInfo.h:18
bool isBuilt_
Flag: True when solver setup has been completed.
virtual void setPrev()
Used to tell Dsolver to assign 'prev' values.
unsigned int getNumVarPools() const
Returns number of local pools that are updated by solver.
Definition: Stoich.cpp:510
unsigned int getNumLocalVoxels() const
Inherited from ZombiePoolInterface.
Definition: Ksolve.cpp:427
virtual void setBlock(const vector< double > &values)=0
void setBlock(const vector< double > &values)
Definition: Ksolve.cpp:776
virtual void getBlock(vector< double > &values) const =0

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void Ksolve::reinit ( const Eref e,
ProcPtr  p 
)

Definition at line 600 of file Ksolve.cpp.

References ProcInfo::dt, ZombiePoolInterface::isBuilt_, pools_, and stoichPtr_.

Referenced by initCinfo(), and initReinit().

601 {
602  if ( !stoichPtr_ )
603  return;
604 
605  if ( isBuilt_ )
606  {
607  for ( unsigned int i = 0 ; i < pools_.size(); ++i )
608  pools_[i].reinit( p->dt );
609  }
610  else
611  {
612  cout << "Warning:Ksolve::reinit: Reaction system not initialized\n";
613  return;
614  }
615 
616 #if PARALLELIZE_KSOLVE_WITH_CPP11_ASYNC
617  if( 1 < getNumThreads( ) )
618  cout << "Debug: User wants Ksolve with " << numThreads_ << " threads" << endl;
619 #endif
620 
621 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
void reinit(const Eref &e, ProcPtr p)
Definition: Ksolve.cpp:600
Stoich * stoichPtr_
Utility ptr used to help Pool Id lookups by the Ksolve.
Definition: Ksolve.h:162
double dt
Definition: ProcInfo.h:18
bool isBuilt_
Flag: True when solver setup has been completed.

+ Here is the caller graph for this function:

void Ksolve::setBlock ( const vector< double > &  values)
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().

777 {
778  unsigned int startVoxel = values[0];
779  unsigned int numVoxels = values[1];
780  unsigned int startPool = values[2];
781  unsigned int numPools = values[3];
782 
783  assert( startVoxel >= startVoxel_ );
784  assert( numVoxels <= pools_.size() );
785  assert( pools_.size() > 0 );
786  assert( numPools + startPool <= pools_[0].size() );
787  assert( values.size() == 4 + numVoxels * numPools );
788 
789  for ( unsigned int i = 0; i < numVoxels; ++i )
790  {
791  double* v = pools_[ startVoxel + i ].varS();
792  for ( unsigned int j = 0; j < numPools; ++j )
793  {
794  v[ j + startPool ] = values[ 4 + j * numVoxels + i ];
795  }
796  }
797 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
unsigned int startVoxel_
First voxel indexed on the current node.
Definition: Ksolve.h:159

+ Here is the caller graph for this function:

void Ksolve::setDiffConst ( const Eref e,
double  val 
)
virtual

Diffusion constant: Only one per pool, voxel number is ignored.

Implements ZombiePoolInterface.

Definition at line 713 of file Ksolve.cpp.

714 {
715  ; // Do nothing.
716 }
void Ksolve::setDsolve ( Id  dsolve)
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().

407 {
408  if ( dsolve == Id () )
409  {
410  dsolvePtr_ = 0;
411  dsolve_ = Id();
412  }
413  else if ( dsolve.element()->cinfo()->isA( "Dsolve" ) )
414  {
415  dsolve_ = dsolve;
416  dsolvePtr_ = reinterpret_cast< ZombiePoolInterface* >(
417  dsolve.eref().data() );
418  }
419  else
420  {
421  cout << "Warning: Ksolve::setDsolve: Object '" << dsolve.path() <<
422  "' should be class Dsolve, is: " <<
423  dsolve.element()->cinfo()->name() << endl;
424  }
425 }
char * data() const
Definition: Eref.cpp:41
Element * element() const
Synonym for Id::operator()()
Definition: Id.cpp:113
std::string path(const std::string &separator="/") const
Definition: Id.cpp:76
Id dsolve_
Definition: Ksolve.h:167
Eref eref() const
Definition: Id.cpp:125
const std::string & name() const
Definition: Cinfo.cpp:260
ZombiePoolInterface * dsolvePtr_
Pointer to diffusion solver.
Definition: Ksolve.h:170
bool isA(const string &ancestor) const
Definition: Cinfo.cpp:280
const Cinfo * cinfo() const
Definition: Element.cpp:66
Definition: Id.h:17

+ Here is the call graph for this function:

void Ksolve::setEpsAbs ( double  val)

Definition at line 279 of file Ksolve.cpp.

References epsAbs_.

Referenced by initCinfo().

280 {
281  if ( epsAbs < 0 )
282  epsAbs_ = 1.0e-4;
283  else
284  epsAbs_ = epsAbs;
285 }
double epsAbs_
Definition: Ksolve.h:139

+ Here is the caller graph for this function:

void Ksolve::setEpsRel ( double  val)

Definition at line 293 of file Ksolve.cpp.

References epsRel_.

Referenced by initCinfo().

294 {
295  if ( epsRel < 0 )
296  {
297  epsRel_ = 1.0e-6;
298  }
299  else
300  {
301  epsRel_ = epsRel;
302  }
303 }
double epsRel_
Definition: Ksolve.h:140

+ Here is the caller graph for this function:

void Ksolve::setMethod ( string  method)

Definition at line 250 of file Ksolve.cpp.

References method_.

Referenced by initCinfo().

251 {
252 #if USE_GSL
253  if ( method == "rk5" || method == "gsl" )
254  {
255  method_ = "rk5";
256  }
257  else if ( method == "rk4" || method == "rk2" ||
258  method == "rk8" || method == "rkck" )
259  {
260  method_ = method;
261  }
262  else
263  {
264  cout << "Warning: Ksolve::setMethod: '" << method <<
265  "' not known, using rk5\n";
266  method_ = "rk5";
267  }
268 #elif USE_BOOST_ODE
269  // TODO: Check for boost related methods.
270  method_ = method;
271 #endif
272 }
string method_
Definition: Ksolve.h:138

+ Here is the caller graph for this function:

void Ksolve::setN ( const Eref e,
double  val 
)
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_.

684 {
685  unsigned int vox = getVoxelIndex( e );
686  if ( vox != OFFNODE )
687  pools_[vox].setN( getPoolIndex( e ), v );
688 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
unsigned int getPoolIndex(const Eref &e) const
Return pool index, using Stoich ptr to do lookup.
Definition: Ksolve.cpp:666
const unsigned int OFFNODE
Definition: Ksolve.cpp:38
unsigned int getVoxelIndex(const Eref &e) const
Definition: Ksolve.cpp:671

+ Here is the call graph for this function:

void Ksolve::setNinit ( const Eref e,
double  val 
)
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_.

699 {
700  unsigned int vox = getVoxelIndex( e );
701  if ( vox != OFFNODE )
702  pools_[vox].setNinit( getPoolIndex( e ), v );
703 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
unsigned int getPoolIndex(const Eref &e) const
Return pool index, using Stoich ptr to do lookup.
Definition: Ksolve.cpp:666
const unsigned int OFFNODE
Definition: Ksolve.cpp:38
unsigned int getVoxelIndex(const Eref &e) const
Definition: Ksolve.cpp:671

+ Here is the call graph for this function:

void Ksolve::setNumAllVoxels ( unsigned int  num)
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().

439 {
440  if ( numVoxels == 0 )
441  {
442  return;
443  }
444  pools_.resize( numVoxels );
445 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156

+ Here is the caller graph for this function:

void Ksolve::setNumPools ( unsigned int  num)
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().

724 {
725  unsigned int numVoxels = pools_.size();
726  for ( unsigned int i = 0 ; i < numVoxels; ++i )
727  {
728  pools_[i].resizeArrays( numPoolSpecies );
729  }
730 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156

+ Here is the caller graph for this function:

void Ksolve::setNvec ( unsigned int  voxel,
vector< double >  vec 
)

Definition at line 457 of file Ksolve.cpp.

References pools_.

Referenced by initCinfo().

458 {
459  if ( voxel < pools_.size() )
460  {
461  if ( nVec.size() != pools_[voxel].size() )
462  {
463  cout << "Warning: Ksolve::setNvec: size mismatch ( " <<
464  nVec.size() << ", " << pools_[voxel].size() << ")\n";
465  return;
466  }
467  double* s = pools_[voxel].varS();
468  for ( unsigned int i = 0; i < nVec.size(); ++i )
469  s[i] = nVec[i];
470  }
471 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156

+ Here is the caller graph for this function:

void Ksolve::setStoich ( Id  stoich)
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_.

354 {
355  assert( stoich.element()->cinfo()->isA( "Stoich" ) );
356  stoich_ = stoich;
357  stoichPtr_ = reinterpret_cast< Stoich* >( stoich.eref().data() );
358  if ( !isBuilt_ )
359  {
360  OdeSystem ode;
361  ode.epsAbs = epsAbs_;
362  ode.epsRel = epsRel_;
363  // ode.initStepSize = getEstimatedDt();
364  ode.initStepSize = 0.01; // This will be overridden at reinit.
365  ode.method = method_;
366 #ifdef USE_GSL
367  ode.gslSys.dimension = stoichPtr_->getNumAllPools();
368  if ( ode.gslSys.dimension == 0 ) {
369  stoichPtr_ = 0;
370  return; // No pools, so don't bother.
371  }
372  innerSetMethod( ode, method_ );
373  ode.gslSys.function = &VoxelPools::gslFunc;
374  ode.gslSys.jacobian = 0;
375  innerSetMethod( ode, method_ );
376  unsigned int numVoxels = pools_.size();
377  for ( unsigned int i = 0 ; i < numVoxels; ++i )
378  {
379  ode.gslSys.params = &pools_[i];
380  pools_[i].setStoich( stoichPtr_, &ode );
381  // pools_[i].setIntDt( ode.initStepSize ); // We're setting it up anyway
382  }
383 #elif USE_BOOST_ODE
384  ode.dimension = stoichPtr_->getNumAllPools();
385  ode.boostSys.epsAbs = epsAbs_;
386  ode.boostSys.epsRel = epsRel_;
387  ode.boostSys.method = method_;
388  if ( ode.dimension == 0 )
389  return; // No pools, so don't bother.
390  unsigned int numVoxels = pools_.size();
391  for ( unsigned int i = 0 ; i < numVoxels; ++i )
392  {
393  ode.boostSys.params = &pools_[i];
394  pools_[i].setStoich( stoichPtr_, &ode );
395  }
396 #endif
397  isBuilt_ = true;
398  }
399 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156
char * data() const
Definition: Eref.cpp:41
Element * element() const
Synonym for Id::operator()()
Definition: Id.cpp:113
std::string method
Definition: OdeSystem.h:28
double epsAbs
Definition: OdeSystem.h:36
Stoich * stoichPtr_
Utility ptr used to help Pool Id lookups by the Ksolve.
Definition: Ksolve.h:162
double epsRel
Definition: OdeSystem.h:37
Eref eref() const
Definition: Id.cpp:125
Definition: Stoich.h:49
string method_
Definition: Ksolve.h:138
double initStepSize
Definition: OdeSystem.h:35
bool isBuilt_
Flag: True when solver setup has been completed.
bool isA(const string &ancestor) const
Definition: Cinfo.cpp:280
const Cinfo * cinfo() const
Definition: Element.cpp:66
unsigned int getNumAllPools() const
Definition: Stoich.cpp:520
double epsRel_
Definition: Ksolve.h:140
double epsAbs_
Definition: Ksolve.h:139

+ Here is the call graph for this function:

void Ksolve::updateRateTerms ( unsigned int  index)
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().

642 {
643  if ( index == ~0U )
644  {
645  // unsigned int numCrossRates = stoichPtr_->getNumRates() - stoichPtr_->getNumCoreRates();
646  for ( unsigned int i = 0 ; i < pools_.size(); ++i )
647  {
648  // pools_[i].resetXreacScale( numCrossRates );
649  pools_[i].updateAllRateTerms( stoichPtr_->getRateTerms(),
651  }
652  }
653  else if ( index < stoichPtr_->getNumRates() )
654  {
655  for ( unsigned int i = 0 ; i < pools_.size(); ++i )
657  stoichPtr_->getNumCoreRates(), index );
658  }
659 }
void updateRateTerms(unsigned int index)
Definition: Ksolve.cpp:641
vector< VoxelPools > pools_
Definition: Ksolve.h:156
Stoich * stoichPtr_
Utility ptr used to help Pool Id lookups by the Ksolve.
Definition: Ksolve.h:162
unsigned int getNumCoreRates() const
Number of rate terms for reactions purely on this compartment.
Definition: Stoich.cpp:574
const vector< RateTerm * > & getRateTerms() const
Returns a reference to the entire rates_ vector.
Definition: Stoich.cpp:589

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

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().

801 {
802  // For now we assume identical numbers of voxels. Also assume
803  // identical voxel junctions. But it should not be too hard to
804  // update those too.
805  if ( vols.size() == pools_.size() )
806  {
807  for ( unsigned int i = 0; i < vols.size(); ++i )
808  {
809  pools_[i].setVolumeAndDependencies( vols[i] );
810  }
811  updateRateTerms( ~0U );
812  }
813 }
void updateRateTerms(unsigned int index)
Definition: Ksolve.cpp:641
vector< VoxelPools > pools_
Definition: Ksolve.h:156

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

double Ksolve::volume ( unsigned int  i) const
virtual

Return volume of voxel i.

Implements ZombiePoolInterface.

Definition at line 746 of file Ksolve.cpp.

References pools_.

747 {
748  if ( pools_.size() > i )
749  return pools_[i].getVolume();
750  return 0.0;
751 }
vector< VoxelPools > pools_
Definition: Ksolve.h:156

Member Data Documentation

Id Ksolve::dsolve_
private

Id of diffusion solver, needed for coordinating numerics.

Definition at line 167 of file Ksolve.h.

Referenced by getDsolve(), print(), and setDsolve().

ZombiePoolInterface* Ksolve::dsolvePtr_
private

Pointer to diffusion solver.

Definition at line 170 of file Ksolve.h.

Referenced by process(), and setDsolve().

double Ksolve::epsAbs_
private

Definition at line 139 of file Ksolve.h.

Referenced by getEpsAbs(), setEpsAbs(), and setStoich().

double Ksolve::epsRel_
private

Definition at line 140 of file Ksolve.h.

Referenced by getEpsRel(), setEpsRel(), and setStoich().

string Ksolve::method_
private

Definition at line 138 of file Ksolve.h.

Referenced by getMethod(), print(), setMethod(), and setStoich().

vector< VoxelPools > Ksolve::pools_
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().

unsigned int Ksolve::startVoxel_
private

First voxel indexed on the current node.

Definition at line 159 of file Ksolve.h.

Referenced by getBlock(), getVoxelIndex(), and setBlock().

Stoich* Ksolve::stoichPtr_
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().


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