MOOSE - Multiscale Object Oriented Simulation Environment
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
VoxelPoolsBase.cpp
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-2010 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 #include "header.h"
10 
11 #include "VoxelPoolsBase.h"
12 #include "../mesh/VoxelJunction.h"
13 #include "XferInfo.h"
14 #include "ZombiePoolInterface.h"
15 #include "RateTerm.h"
16 #include "FuncTerm.h"
17 #include "SparseMatrix.h"
18 #include "KinSparseMatrix.h"
19 #include "Stoich.h"
20 
22 // Class definitions
24 
26  stoichPtr_( 0 ),
27  S_(1),
28  Sinit_(1),
29  volume_(1.0)
30 {
31  ;
32 }
33 
35 {}
36 
38 // Array ops
41 void VoxelPoolsBase::resizeArrays( unsigned int totNumPools )
42 {
43  S_.resize( totNumPools, 0.0 );
44  Sinit_.resize( totNumPools, 0.0);
45 }
46 
48 {
49  S_ = Sinit_;
50 }
51 
53 // Access functions
55 const double* VoxelPoolsBase::S() const
56 {
57  return &S_[0];
58 }
59 
60 vector< double >& VoxelPoolsBase::Svec()
61 {
62  return S_;
63 }
64 
66 {
67  return &S_[0];
68 }
69 
70 const double* VoxelPoolsBase::Sinit() const
71 {
72  return &Sinit_[0];
73 }
74 
76 {
77  return &Sinit_[0];
78 }
79 
80 unsigned int VoxelPoolsBase::size() const
81 {
82  return Sinit_.size();
83 }
84 
85 void VoxelPoolsBase::setVolume( double vol )
86 {
87  volume_ = vol;
88 }
89 
91 {
92  return volume_;
93 }
94 
96 {
97  double ratio = vol / volume_;
98  volume_ = vol;
99  for ( vector< double >::iterator
100  i = Sinit_.begin(); i != Sinit_.end(); ++i )
101  *i *= ratio;
102 
103  for ( vector< double >::iterator i = S_.begin(); i != S_.end(); ++i )
104  *i *= ratio;
105 
106  // I would like to update the xReacScaleSubstreates and Products here,
107  // but I don't know the order of their reactions. So leave it to
108  // a subsequent call via Ksolve or Stoich.
109 }
110 
112  double ratio, const Stoich* stoichPtr )
113 {
114  volume_ *= ratio; // Scale vol
115  for ( vector< double >::iterator
116  i = Sinit_.begin(); i != Sinit_.end(); ++i )
117  *i *= ratio; // Scale Bufs
118  // Here we also need to set the Ns for the buffered pools.
119  unsigned int start = stoichPtr_->getNumVarPools();
120  unsigned int end = start + stoichPtr_->getNumBufPools();
121  assert( end == Sinit_.size() );
122  for ( unsigned int i = start; i < end; ++i )
123  {
124  // Must not reassign pools that are controlled by functions.
125  if ( !stoichPtr->isFuncTarget(i) )
126  S_[i] = Sinit_[i];
127  }
128 
129  // Scale rates. Start by clearing out old rates if any
130  for ( unsigned int i = 0; i < rates_.size(); ++i )
131  delete( rates_[i] );
132 
133  unsigned int numCoreRates = stoichPtr->getNumCoreRates();
134  const vector< RateTerm* >& rates = stoichPtr->getRateTerms();
135  rates_.resize( rates.size() );
136  for ( unsigned int i = 0; i < numCoreRates; ++i )
137  rates_[i] = rates[i]->copyWithVolScaling( getVolume(), 1, 1 );
138  for ( unsigned int i = numCoreRates; i < rates.size(); ++i )
139  {
140  rates_[i] = rates[i]->copyWithVolScaling( getVolume(),
141  getXreacScaleSubstrates(i - numCoreRates),
142  getXreacScaleProducts(i - numCoreRates ) );
143  }
144 }
145 
147 // Zombie Pool Access functions
149 
150 void VoxelPoolsBase::setN( unsigned int i, double v )
151 {
152  S_[i] = v;
153  if ( S_[i] < 0.0 )
154  S_[i] = 0.0;
155 }
156 
157 double VoxelPoolsBase::getN( unsigned int i ) const
158 {
159  return S_[i];
160 }
161 
162 void VoxelPoolsBase::setNinit( unsigned int i, double v )
163 {
164  Sinit_[i] = v;
165  if ( Sinit_[i] < 0.0 )
166  Sinit_[i] = 0.0;
167 }
168 
169 double VoxelPoolsBase::getNinit( unsigned int i ) const
170 {
171  return Sinit_[i];
172 }
173 
174 void VoxelPoolsBase::setDiffConst( unsigned int i, double v )
175 {
176  ; // Do nothing.
177 }
178 
179 double VoxelPoolsBase::getDiffConst( unsigned int i ) const
180 {
181  return 0;
182 }
183 
185 // Handle cross compartment reactions
188  const vector< unsigned int >& poolIndex,
189  const vector< double >& values,
190  const vector< double >& lastValues,
191  unsigned int voxelIndex )
192 {
193  unsigned int offset = voxelIndex * poolIndex.size();
194  vector< double >::const_iterator i = values.begin() + offset;
195  vector< double >::const_iterator j = lastValues.begin() + offset;
196  for ( vector< unsigned int >::const_iterator
197  k = poolIndex.begin(); k != poolIndex.end(); ++k )
198  {
199  S_[*k] += *i++ - *j++;
200  }
201 }
202 
204  const vector< unsigned int >& poolIndex,
205  const vector< double >& values,
206  unsigned int numProxyPools,
207  unsigned int voxelIndex )
208 {
209  unsigned int offset = voxelIndex * poolIndex.size();
210  vector< double >::const_iterator i = values.begin() + offset;
211  unsigned int proxyEndIndex = stoichPtr_->getNumVarPools() +
213  for ( vector< unsigned int >::const_iterator
214  k = poolIndex.begin(); k != poolIndex.end(); ++k )
215  {
216  // if ( *k >= S_.size() - numProxyPools )
217  if ( *k >= stoichPtr_->getNumVarPools() && *k < proxyEndIndex )
218  {
219  // cout << S_[*k] << ", " << Sinit_[*k] << ", " << *i << endl;
220  Sinit_[*k] = *i;
221  S_[*k] = *i;
222  }
223  i++;
224  }
225 }
226 
228  unsigned int voxelIndex,
229  vector< double >& values,
230  const vector< unsigned int >& poolIndex)
231 {
232  unsigned int offset = voxelIndex * poolIndex.size();
233  vector< double >::iterator i = values.begin() + offset;
234  for ( vector< unsigned int >::const_iterator
235  k = poolIndex.begin(); k != poolIndex.end(); ++k )
236  {
237  *i++ = S_[*k];
238  }
239 }
240 
242  unsigned int comptIndex, Id otherComptId, unsigned int voxel )
243 {
244  if ( comptIndex >= proxyPoolVoxels_.size() )
245  {
246  proxyPoolVoxels_.resize( comptIndex + 1 );
247  }
248 
249  proxyPoolVoxels_[comptIndex].push_back( voxel );
250  proxyComptMap_[otherComptId] = comptIndex;
251 }
252 
254  unsigned int comptIndex, unsigned int transferIndex )
255 {
256  if ( comptIndex >= proxyTransferIndex_.size() )
257  proxyTransferIndex_.resize( comptIndex + 1 );
258  proxyTransferIndex_[comptIndex].push_back( transferIndex );
259 }
260 
261 bool VoxelPoolsBase::hasXfer( unsigned int comptIndex ) const
262 {
263  if ( comptIndex >= proxyPoolVoxels_.size() )
264  return false;
265  return (proxyPoolVoxels_[ comptIndex ].size() > 0);
266 }
267 
269 {
270  if ( i1 == Id () )
271  return false;
272  if ( proxyComptMap_.find( i1 ) == proxyComptMap_.end() )
273  return false;
274  if ( i2 == Id() ) // This is intentionally blank, only one jn.
275  return true;
276  // If there is an i2 but it isn't on the map, then not connected.
277  if ( proxyComptMap_.find( i2 ) == proxyComptMap_.end() )
278  return false;
279  return true;
280 }
281 
283 // Cross reaction stuff.
285 
286 void VoxelPoolsBase::resetXreacScale( unsigned int size )
287 {
288  xReacScaleSubstrates_.assign( size, 1.0 );
289  xReacScaleProducts_.assign( size, 1.0 );
290 }
291 
292 void VoxelPoolsBase::forwardReacVolumeFactor( unsigned int i, double volume )
293 {
294  assert( i < xReacScaleSubstrates_.size() );
295  xReacScaleSubstrates_[i] *= volume / getVolume();
296  // cout << "forwardReacVolumeFactor[" << i << "] = " << xReacScaleSubstrates_[i] <<endl;
297 }
298 
299 void VoxelPoolsBase::backwardReacVolumeFactor( unsigned int i, double volume )
300 {
301  assert( i < xReacScaleProducts_.size() );
302  xReacScaleProducts_[i] *= volume / getVolume();
303  // cout << "backwardReacVolumeFactor[" << i << "] = "<< xReacScaleProducts_[i] <<endl;
304 }
305 
306 double VoxelPoolsBase::getXreacScaleSubstrates( unsigned int i ) const
307 {
308  return xReacScaleSubstrates_[i];
309 }
310 
311 double VoxelPoolsBase::getXreacScaleProducts( unsigned int i ) const
312 {
313  return xReacScaleProducts_[i];
314 }
315 
321  const vector< Id >& offSolverReacs,
322  const vector< pair< Id, Id > >& offSolverReacCompts )
323 {
324  assert (offSolverReacs.size() == offSolverReacCompts.size() );
325  // unsigned int numCoreRates = stoichPtr_->getNumCoreRates();
326  for ( unsigned int i = 0; i < offSolverReacCompts.size(); ++i )
327  {
328  const pair< Id, Id >& p = offSolverReacCompts[i];
329  if ( !isVoxelJunctionPresent( p.first, p.second) )
330  {
331  Id reacId = offSolverReacs[i];
332  const Cinfo* reacCinfo = reacId.element()->cinfo();
333  unsigned int k = stoichPtr_->convertIdToReacIndex( offSolverReacs[i] );
334  // Start by replacing the immediate cross reaction term.
335  if ( rates_[k] )
336  delete rates_[k];
337  rates_[k] = new ExternReac;
338  if ( stoichPtr_->getOneWay() )
339  {
340  k++; // Delete the next entry too, it is the reverse reacn.
341  assert( k < rates_.size() );
342  if ( reacCinfo->isA( "ReacBase" ) )
343  {
344  if ( rates_[k] )
345  delete rates_[k];
346  rates_[k] = new ExternReac;
347  }
348  if ( reacCinfo->isA( "CplxEnzBase" ) ) // Delete next two.
349  {
350  if ( rates_[k] )
351  delete rates_[k];
352  rates_[k] = new ExternReac;
353  k++;
354  assert( k < rates_.size() );
355  if ( rates_[k] )
356  delete rates_[k];
357  rates_[k] = new ExternReac;
358  }
359  }
360  else
361  {
362  if ( reacCinfo->isA( "CplxEnzBase" ) ) // Delete next one.
363  {
364  k++;
365  assert( k < rates_.size() );
366  if ( rates_[k] )
367  delete rates_[k];
368  rates_[k] = new ExternReac;
369  }
370  }
371  }
372  }
373 }
374 
377 {
378  cout << "S_.size=" << S_.size() << ", volume = " << volume_ << endl;
379  cout << "proxyPoolsVoxels.size()=" << proxyPoolVoxels_.size() <<
380  ", proxyTransferIndex.size()=" << proxyTransferIndex_.size() <<
381  endl;
382  assert( proxyPoolVoxels_.size() == proxyTransferIndex_.size() );
383  for ( unsigned int i = 0; i < proxyPoolVoxels_.size(); ++i )
384  {
385  cout << "ppv[" << i << "]=";
386  const vector< unsigned int >& ppv = proxyPoolVoxels_[i];
387  for ( unsigned int j = 0; j < ppv.size(); ++j )
388  {
389  cout << " " << ppv[j];
390  }
391  cout << endl;
392  }
393  for ( unsigned int i = 0; i < proxyTransferIndex_.size(); ++i )
394  {
395  cout << "pti[" << i << "]=";
396  const vector< unsigned int >& pti = proxyTransferIndex_[i];
397  for ( unsigned int j = 0; j < pti.size(); ++j )
398  {
399  cout << " " << pti[j];
400  }
401  cout << endl;
402  }
403  cout <<
404  "xReacScaleSubstrates.size()=" << xReacScaleSubstrates_.size() <<
405  ", xReacScaleProducts.size()=" << xReacScaleProducts_.size() <<
406  endl;
407  assert( xReacScaleSubstrates_.size() == xReacScaleProducts_.size() );
408  for ( unsigned int i = 0; i < xReacScaleSubstrates_.size(); ++i )
409  {
410  cout << i << " " << xReacScaleSubstrates_[i] << " " <<
411  xReacScaleProducts_[i] << endl;
412  }
413  cout << "############## RATES ######################\n";
414  for ( unsigned int i = 0; i < rates_.size(); ++i )
415  {
416  cout << i << " : " << rates_[i]->getR1() << ", " <<
417  rates_[i]->getR2() << endl;
418  }
419 }
double getXreacScaleProducts(unsigned int i) const
const Stoich * stoichPtr_
double getDiffConst(unsigned int) const
void print() const
Debugging utility.
const double * S() const
void addProxyTransferIndex(unsigned int comptIndex, unsigned int transferIndex)
Element * element() const
Synonym for Id::operator()()
Definition: Id.cpp:113
vector< double > xReacScaleProducts_
double getNinit(unsigned int) const
vector< double > S_
double getVolume() const
Return the volume of the voxel.
double * varSinit()
vector< double > Sinit_
vector< vector< unsigned int > > proxyTransferIndex_
void xferIn(const vector< unsigned int > &poolIndex, const vector< double > &values, const vector< double > &lastValues, unsigned int voxelIndex)
virtual ~VoxelPoolsBase()
bool isVoxelJunctionPresent(Id i1, Id i2) const
Definition: Stoich.h:49
bool hasXfer(unsigned int comptIndex) const
True when this voxel has data to be transferred.
double getN(unsigned int) const
vector< double > & Svec()
Returns a handle to the mol # vector.
void setN(unsigned int i, double v)
void setDiffConst(unsigned int, double v)
void filterCrossRateTerms(const vector< Id > &offSolverReacs, const vector< pair< Id, Id > > &offSolverReacCompts)
double getXreacScaleSubstrates(unsigned int i) const
void resizeArrays(unsigned int totNumPools)
Allocate # of pools.
void xferOut(unsigned int voxelIndex, vector< double > &values, const vector< unsigned int > &poolIndex)
Assembles data values for sending out for x-compt reacs.
void xferInOnlyProxies(const vector< unsigned int > &poolIndex, const vector< double > &values, unsigned int numProxyPools, unsigned int voxelIndex)
unsigned int size() const
const double * Sinit() const
bool isA(const string &ancestor) const
Definition: Cinfo.cpp:280
void scaleVolsBufsRates(double ratio, const Stoich *stoichPtr)
const Cinfo * cinfo() const
Definition: Element.cpp:66
void resetXreacScale(unsigned int size)
vector< double > xReacScaleSubstrates_
unsigned int getNumBufPools() const
Returns number of local buffered pools.
Definition: Stoich.cpp:515
unsigned int getNumVarPools() const
Returns number of local pools that are updated by solver.
Definition: Stoich.cpp:510
void setVolume(double vol)
Just assigns the volume without any cascading to other values.
bool isFuncTarget(unsigned int poolIndex) const
Returns true if the specified pool is controlled by a func.
Definition: Stoich.cpp:606
void setNinit(unsigned int, double v)
unsigned int convertIdToReacIndex(Id id) const
Definition: Stoich.cpp:1474
static const Cinfo * reacCinfo
Definition: Reac.cpp:41
Definition: Id.h:17
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
virtual void setVolumeAndDependencies(double vol)
vector< vector< unsigned int > > proxyPoolVoxels_
bool getOneWay() const
Definition: Stoich.cpp:305
void forwardReacVolumeFactor(unsigned int i, double volume)
unsigned int getNumProxyPools() const
Definition: Stoich.cpp:525
map< Id, unsigned int > proxyComptMap_
vector< RateTerm * > rates_
Definition: Cinfo.h:18
void addProxyVoxy(unsigned int comptIndex, Id comptId, unsigned int voxel)
Adds the index of a voxel to which proxy data should go.
void backwardReacVolumeFactor(unsigned int i, double volume)