MOOSE - Multiscale Object Oriented Simulation Environment
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
Ksolve.h
Go to the documentation of this file.
1 /**********************************************************************
2 ** This program is part of 'MOOSE', the
3 ** Messaging Object Oriented Simulation Environment.
4 ** Copyright (C) 2003-2014 Upinder S. Bhalla. and NCBS
5 ** It is made available under the terms of the
6 ** GNU Lesser General Public License version 2.1
7 ** See the file COPYING.LIB for the full notice.
8 **********************************************************************/
9 
10 #ifndef _KSOLVE_H
11 #define _KSOLVE_H
12 
13 class Stoich;
14 
16 {
17 public:
18  Ksolve();
19  ~Ksolve();
20 
22  // Field assignment stuff
25  string getMethod() const;
26  void setMethod( string method );
27 
29  double getEpsAbs() const;
30  void setEpsAbs( double val );
31 
33  double getEpsRel() const;
34  void setEpsRel( double val );
35 
37  Id getStoich() const;
38  void setStoich( Id stoich );
39 
41  Id getDsolve() const;
42  void setDsolve( Id dsolve );
43 
44 
45  unsigned int getNumLocalVoxels() const;
46  unsigned int getNumAllVoxels() const;
51  void setNumAllVoxels( unsigned int num );
52 
54  vector< double > getNvec( unsigned int voxel) const;
55  void setNvec( unsigned int voxel, vector< double > vec );
56 
57 #if PARALLELIZE_KSOLVE_WITH_CPP11_ASYNC
58  // Set number of threads to use (for deterministic case only).
59  unsigned int getNumThreads( ) const;
60  void setNumThreads( unsigned int x );
61 
62  // Parallel advance().
63  void parallel_advance(int begin, int end, size_t nWorkers, ProcPtr p);
64 #endif
65 
70  double getEstimatedDt() const;
71 
73  // Dest Finfos
75  void process( const Eref& e, ProcPtr p );
76  void reinit( const Eref& e, ProcPtr p );
77  void initProc( const Eref& e, ProcPtr p );
78  void initReinit( const Eref& e, ProcPtr p );
84  void updateVoxelVol( vector< double > vols );
86  // Utility for SrcFinfo
88 
90  // Solver interface functions
92  unsigned int getPoolIndex( const Eref& e ) const;
93  unsigned int getVoxelIndex( const Eref& e ) const;
94 
96  // ZombiePoolInterface inherited functions
98 
99  void setN( const Eref& e, double v );
100  double getN( const Eref& e ) const;
101  void setNinit( const Eref& e, double v );
102  double getNinit( const Eref& e ) const;
103  void setDiffConst( const Eref& e, double v );
104  double getDiffConst( const Eref& e ) const;
105 
111  void setNumPools( unsigned int num );
112  unsigned int getNumPools() const;
113 
114  VoxelPoolsBase* pools( unsigned int i );
115  double volume( unsigned int i ) const;
116 
117  void getBlock( vector< double >& values ) const;
118  void setBlock( const vector< double >& values );
119 
120  void matchJunctionVols( vector< double >& vols, Id otherCompt )
121  const;
122 
127  void updateRateTerms( unsigned int index );
128 
130  // for debugging
131  void print() const;
132 
134  static const Cinfo* initCinfo();
135 
136 private:
137 
138  string method_;
139  double epsAbs_;
140  double epsRel_;
141 
142 #if PARALLELIZE_KSOLVE_WITH_CPP11_ASYNC
143 
146  unsigned int numThreads_;
147 #endif
148 
156  vector< VoxelPools > pools_;
157 
159  unsigned int startVoxel_;
160 
163 
168 
171 
172 };
173 
174 #endif // _KSOLVE_H
void updateRateTerms(unsigned int index)
Definition: Ksolve.cpp:641
vector< VoxelPools > pools_
Definition: Ksolve.h:156
void setDsolve(Id dsolve)
Assigns the diffusion solver. Used by the reac solvers.
Definition: Ksolve.cpp:406
static const Cinfo * initCinfo()
Definition: Ksolve.cpp:40
void getBlock(vector< double > &values) const
Definition: Ksolve.cpp:753
void setNinit(const Eref &e, double v)
Set initial # of molecules in given pool and voxel. Bdry cond.
Definition: Ksolve.cpp:698
VoxelPoolsBase * pools(unsigned int i)
Return a pointer to the specified VoxelPool.
Definition: Ksolve.cpp:739
void setNumPools(unsigned int num)
Definition: Ksolve.cpp:723
string getMethod() const
Assigns integration method.
Definition: Ksolve.cpp:245
void reinit(const Eref &e, ProcPtr p)
Definition: Ksolve.cpp:600
Id dsolve_
Definition: Ksolve.h:167
~Ksolve()
Definition: Ksolve.cpp:236
Stoich * stoichPtr_
Utility ptr used to help Pool Id lookups by the Ksolve.
Definition: Ksolve.h:162
void updateVoxelVol(vector< double > vols)
Definition: Ksolve.cpp:800
double volume(unsigned int i) const
Return volume of voxel i.
Definition: Ksolve.cpp:746
void print() const
Definition: Ksolve.cpp:823
unsigned int getPoolIndex(const Eref &e) const
Return pool index, using Stoich ptr to do lookup.
Definition: Ksolve.cpp:666
Definition: Stoich.h:49
string method_
Definition: Ksolve.h:138
unsigned int getNumPools() const
gets number of pools (species) handled by system.
Definition: Ksolve.cpp:732
ZombiePoolInterface * dsolvePtr_
Pointer to diffusion solver.
Definition: Ksolve.h:170
double getDiffConst(const Eref &e) const
Diffusion constant: Only one per pool, voxel number is ignored.
Definition: Ksolve.cpp:718
void setEpsRel(double val)
Definition: Ksolve.cpp:293
void setStoich(Id stoich)
Definition: Ksolve.cpp:353
vector< double > getNvec(unsigned int voxel) const
Returns the vector of pool Num at the specified voxel.
Definition: Ksolve.cpp:447
void matchJunctionVols(vector< double > &vols, Id otherCompt) const
unsigned int startVoxel_
First voxel indexed on the current node.
Definition: Ksolve.h:159
double getEpsAbs() const
Assigns Absolute tolerance for integration.
Definition: Ksolve.cpp:274
Definition: Eref.h:26
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
Definition: Ksolve.h:15
void setN(const Eref &e, double v)
Set # of molecules in given pool and voxel. Varies with time.
Definition: Ksolve.cpp:683
void setNvec(unsigned int voxel, vector< double > vec)
Definition: Ksolve.cpp:457
Definition: Id.h:17
void initReinit(const Eref &e, ProcPtr p)
Definition: Ksolve.cpp:630
double getNinit(const Eref &e) const
get initial # of molecules in given pool and voxel. Bdry cond.
Definition: Ksolve.cpp:705
double epsRel_
Definition: Ksolve.h:140
void setDiffConst(const Eref &e, double v)
Diffusion constant: Only one per pool, voxel number is ignored.
Definition: Ksolve.cpp:713
Id getDsolve() const
Inherited from ZombiePoolInterface.
Definition: Ksolve.cpp:401
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 epsAbs_
Definition: Ksolve.h:139
double getEstimatedDt() const
Definition: Ksolve.cpp:474
void setEpsAbs(double val)
Definition: Ksolve.cpp:279
void setBlock(const vector< double > &values)
Definition: Ksolve.cpp:776
Definition: Cinfo.h:18
unsigned int getVoxelIndex(const Eref &e) const
Definition: Ksolve.cpp:671
Ksolve()
Definition: Ksolve.cpp:216
double getN(const Eref &e) const
Get # of molecules in given pool and voxel. Varies with time.
Definition: Ksolve.cpp:690