Chaste  Release::3.4
AbstractCardiacTissue.cpp
1 /*
2 
3 Copyright (c) 2005-2016, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34 */
35 
36 #include "AbstractCardiacTissue.hpp"
37 
38 #include <boost/scoped_array.hpp>
39 
40 #include "DistributedVector.hpp"
41 #include "AxisymmetricConductivityTensors.hpp"
42 #include "OrthotropicConductivityTensors.hpp"
43 #include "Exception.hpp"
44 #include "ChastePoint.hpp"
45 #include "AbstractChasteRegion.hpp"
46 #include "HeartEventHandler.hpp"
47 #include "PetscTools.hpp"
48 #include "PetscVecTools.hpp"
49 #include "AbstractCvodeCell.hpp"
50 #include "Warnings.hpp"
51 
52 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
55  bool exchangeHalos)
56  : mpMesh(pCellFactory->GetMesh()),
57  mpDistributedVectorFactory(mpMesh->GetDistributedVectorFactory()),
58  mpConductivityModifier(NULL),
59  mHasPurkinje(false),
60  mDoCacheReplication(true),
61  mMeshUnarchived(false),
62  mExchangeHalos(exchangeHalos)
63 {
64  //This constructor is called from the Initialise() method of the CardiacProblem class
65  assert(pCellFactory != NULL);
66  assert(pCellFactory->GetMesh() != NULL);
67 
69  {
70  //Remove the request for a halo exchange
71  mExchangeHalos = false;
72  }
73 
74  unsigned num_local_nodes = mpDistributedVectorFactory->GetLocalOwnership();
75  bool process_has_no_nodes = (num_local_nodes == 0u);
76  if (PetscTools::ReplicateBool(process_has_no_nodes))
77  {
78  /* If there are no nodes on a process then there is a potential for an error in preconditioning.
79  * This is dangerous because the process without nodes may segfault and not propagate the error
80  * to the other processes. Therefore, to avoid deadlock, we share this potential for error between
81  * processes and throw an exception.
82  */
83 #define COVERAGE_IGNORE
84  // This problem normally occurs on 3 or more processes, so we can't cover it - coverage only runs with 1 and 2 processes.
85  EXCEPTION("No cells were assigned some process in AbstractCardiacTissue constructor. Advice: Make total number of processors no greater than number of nodes in the mesh");
86 #undef COVERAGE_IGNORE
87  }
88  unsigned ownership_range_low = mpDistributedVectorFactory->GetLow();
89  mCellsDistributed.resize(num_local_nodes);
90 
91  // Figure out if we're dealing with Purkinje
93  = dynamic_cast<AbstractPurkinjeCellFactory<ELEMENT_DIM,SPACE_DIM>*>(pCellFactory);
94  if (p_purkinje_cell_factory)
95  {
96  mHasPurkinje = true;
97  mPurkinjeCellsDistributed.resize(num_local_nodes);
98  }
99 
101  // Set up cells
103  try
104  {
105  for (unsigned local_index = 0; local_index < num_local_nodes; local_index++)
106  {
107  unsigned global_index = ownership_range_low + local_index;
108  Node<SPACE_DIM>* p_node = mpMesh->GetNode(global_index);
109  mCellsDistributed[local_index] = pCellFactory->CreateCardiacCellForNode(p_node);
110  mCellsDistributed[local_index]->SetUsedInTissueSimulation();
111 
112  if (mHasPurkinje)
113  {
114  mPurkinjeCellsDistributed[local_index]
115  = p_purkinje_cell_factory->CreatePurkinjeCellForNode(p_node, mCellsDistributed[local_index]);
116  mPurkinjeCellsDistributed[local_index]->SetUsedInTissueSimulation();
117  }
118  }
119 
123  if (mHasPurkinje)
124  {
128  }
129  }
130  catch (const Exception& e)
131  {
132  // This catch statement is quite tricky to cover, but it is actually done
133  // in TestCardiacSimulation::TestMono1dSodiumBlockBySettingNamedParameter()
134 
135  // Errors thrown creating cells will often be process-specific
137 
138  // Delete cells
139  // Should really do this for other processes too, but this is all we need
140  // to get memory testing to pass, and leaking when we're about to die isn't
141  // that bad!
142  for (std::vector<AbstractCardiacCellInterface*>::iterator cell_iterator = mCellsDistributed.begin();
143  cell_iterator != mCellsDistributed.end();
144  ++cell_iterator)
145  {
146  delete (*cell_iterator);
147  }
148 
149  throw e;
150  }
152 
153  // Halo nodes (if required)
154  SetUpHaloCells(pCellFactory);
155 
156  HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
159 
160  if (mHasPurkinje)
161  {
164  }
165  HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
166 
167  if(HeartConfig::Instance()->IsMeshProvided() && HeartConfig::Instance()->GetLoadMesh())
168  {
170  }
171  else
172  {
173  // As of r10671 fibre orientation can only be defined when loading a mesh from disc.
175  }
177 }
178 
179 // Constructor used for archiving
180 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
182  : mpMesh(pMesh),
183  mpDistributedVectorFactory(mpMesh->GetDistributedVectorFactory()),
184  mHasPurkinje(false),
185  mDoCacheReplication(true),
186  mMeshUnarchived(true),
187  mExchangeHalos(false)
188 {
191 
194 }
195 
196 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
198 {
199  // Delete cells
200  for (std::vector<AbstractCardiacCellInterface*>::iterator iter = mCellsDistributed.begin();
201  iter != mCellsDistributed.end();
202  ++iter)
203  {
204  delete (*iter);
205  }
206 
207  // Delete cells for halo nodes
208  for (std::vector<AbstractCardiacCellInterface*>::iterator iter = mHaloCellsDistributed.begin();
209  iter != mHaloCellsDistributed.end();
210  ++iter)
211  {
212  delete (*iter);
213  }
214 
215  delete mpIntracellularConductivityTensors;
216 
217  // Delete Purkinje cells
218  for (std::vector<AbstractCardiacCellInterface*>::iterator iter = mPurkinjeCellsDistributed.begin();
219  iter != mPurkinjeCellsDistributed.end();
220  ++iter)
221  {
222  delete (*iter);
223  }
224 
225  // If the mesh was unarchived we need to free it explicitly.
226  if (mMeshUnarchived)
227  {
228  delete mpMesh;
229  }
230 }
231 
232 
233 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
235 {
236  return mHasPurkinje;
237 }
238 
239 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
241 {
242  HeartEventHandler::BeginEvent(HeartEventHandler::READ_MESH);
243  mpConfig = HeartConfig::Instance();
244 
245  if (mpConfig->IsMeshProvided() && mpConfig->GetLoadMesh())
246  {
247  assert(mFibreFilePathNoExtension != "");
248 
249  switch (mpConfig->GetConductivityMedia())
250  {
251  case cp::media_type::Orthotropic:
252  {
253  mpIntracellularConductivityTensors = new OrthotropicConductivityTensors<ELEMENT_DIM,SPACE_DIM>;
254  FileFinder ortho_file(mFibreFilePathNoExtension + ".ortho", RelativeTo::AbsoluteOrCwd);
255  assert(ortho_file.Exists());
256  mpIntracellularConductivityTensors->SetFibreOrientationFile(ortho_file);
257  break;
258  }
259 
260  case cp::media_type::Axisymmetric:
261  {
262  mpIntracellularConductivityTensors = new AxisymmetricConductivityTensors<ELEMENT_DIM,SPACE_DIM>;
263  FileFinder axi_file(mFibreFilePathNoExtension + ".axi", RelativeTo::AbsoluteOrCwd);
264  assert(axi_file.Exists());
265  mpIntracellularConductivityTensors->SetFibreOrientationFile(axi_file);
266  break;
267  }
268 
269  case cp::media_type::NoFibreOrientation:
271  mpIntracellularConductivityTensors = new OrthotropicConductivityTensors<ELEMENT_DIM,SPACE_DIM>;
272  break;
273 
274  default :
276  }
277  }
278  else // Slab defined in config file or SetMesh() called; no fibre orientation assumed
279  {
281  mpIntracellularConductivityTensors = new OrthotropicConductivityTensors<ELEMENT_DIM,SPACE_DIM>;
282  }
283 
284  c_vector<double, SPACE_DIM> intra_conductivities;
285  mpConfig->GetIntracellularConductivities(intra_conductivities);
286 
287  // this definition must be here (and not inside the if statement) because SetNonConstantConductivities() will keep
288  // a pointer to it and we don't want it to go out of scope before Init() is called
289  unsigned num_local_elements = mpMesh->GetNumLocalElements();
290  std::vector<c_vector<double, SPACE_DIM> > hetero_intra_conductivities;
291 
292  if (mpConfig->GetConductivityHeterogeneitiesProvided())
293  {
294  try
295  {
296  assert(hetero_intra_conductivities.size()==0);
297  hetero_intra_conductivities.resize(num_local_elements, intra_conductivities);
298  }
299  catch(std::bad_alloc &r_bad_alloc)
300  {
301 #define COVERAGE_IGNORE
302  std::cout << "Failed to allocate std::vector of size " << num_local_elements << std::endl;
304  throw r_bad_alloc;
305 #undef COVERAGE_IGNORE
306  }
308 
309  std::vector<boost::shared_ptr<AbstractChasteRegion<SPACE_DIM> > > conductivities_heterogeneity_areas;
310  std::vector< c_vector<double,3> > intra_h_conductivities;
311  std::vector< c_vector<double,3> > extra_h_conductivities;
312  HeartConfig::Instance()->GetConductivityHeterogeneities(conductivities_heterogeneity_areas,
313  intra_h_conductivities,
314  extra_h_conductivities);
315 
316  unsigned local_element_index = 0;
317 
318  for (typename AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>::ElementIterator it = mpMesh->GetElementIteratorBegin();
319  it != mpMesh->GetElementIteratorEnd();
320  ++it)
321  {
322 // unsigned element_index = it->GetIndex();
323  // if element centroid is contained in the region
324  ChastePoint<SPACE_DIM> element_centroid(it->CalculateCentroid());
325  for (unsigned region_index=0; region_index< conductivities_heterogeneity_areas.size(); region_index++)
326  {
327  if ( conductivities_heterogeneity_areas[region_index]->DoesContain(element_centroid) )
328  {
329  //We don't use ublas vector assignment here, because we might be getting a subvector of a 3-vector
330  for (unsigned i=0; i<SPACE_DIM; i++)
331  {
332  hetero_intra_conductivities[local_element_index][i] = intra_h_conductivities[region_index][i];
333  }
334  }
335  }
336  local_element_index++;
337  }
338 
339  mpIntracellularConductivityTensors->SetNonConstantConductivities(&hetero_intra_conductivities);
340  }
341  else
342  {
343  mpIntracellularConductivityTensors->SetConstantConductivities(intra_conductivities);
344  }
345 
346  mpIntracellularConductivityTensors->Init(this->mpMesh);
347  HeartEventHandler::EndEvent(HeartEventHandler::READ_MESH);
348 }
349 
350 
351 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
353 {
354  mDoCacheReplication = doCacheReplication;
355 }
356 
357 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
359 {
360  return mDoCacheReplication;
361 }
362 
363 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
364 const c_matrix<double, SPACE_DIM, SPACE_DIM>& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetIntracellularConductivityTensor(unsigned elementIndex)
365 {
366  assert( mpIntracellularConductivityTensors);
367  if (mpConductivityModifier==NULL)
368  {
369  return (*mpIntracellularConductivityTensors)[elementIndex];
370  }
371  else
372  {
373  return mpConductivityModifier->rGetModifiedConductivityTensor(elementIndex, (*mpIntracellularConductivityTensors)[elementIndex], 0u);
374  }
375 }
376 
377 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
378 const c_matrix<double, SPACE_DIM, SPACE_DIM>& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetExtracellularConductivityTensor(unsigned elementIndex)
379 {
380  EXCEPTION("Monodomain tissues do not have extracellular conductivity tensors.");
381 }
382 
383 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
385 {
386  assert(mpDistributedVectorFactory->GetLow() <= globalIndex &&
387  globalIndex < mpDistributedVectorFactory->GetHigh());
388  return mCellsDistributed[globalIndex - mpDistributedVectorFactory->GetLow()];
389 }
390 
391 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
393 {
394  assert(mpDistributedVectorFactory->GetLow() <= globalIndex &&
395  globalIndex < mpDistributedVectorFactory->GetHigh());
396  EXCEPT_IF_NOT(mHasPurkinje);
397  return mPurkinjeCellsDistributed[globalIndex - mpDistributedVectorFactory->GetLow()];
398 }
399 
400 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
402 {
403  std::map<unsigned, unsigned>::const_iterator node_position;
404  // First search the halo
405  if ((node_position=mHaloGlobalToLocalIndexMap.find(globalIndex)) != mHaloGlobalToLocalIndexMap.end())
406  {
407  //Found a halo node
408  return mHaloCellsDistributed[node_position->second];
409  }
410  // Then search the owned node
411  if ( mpDistributedVectorFactory->IsGlobalIndexLocal(globalIndex) )
412  {
413  //Found an owned node
414  return mCellsDistributed[globalIndex - mpDistributedVectorFactory->GetLow()];
415  }
416  // Not here
417  EXCEPTION("Requested node/halo " << globalIndex << " does not belong to processor " << PetscTools::GetMyRank());
418 }
419 
420 
421 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
423 {
424  std::set<unsigned> halos_as_set;
425  for (unsigned proc=0; proc<PetscTools::GetNumProcs(); proc++)
426  {
427  halos_as_set.insert(mNodesToReceivePerProcess[proc].begin(), mNodesToReceivePerProcess[proc].end());
428  }
429  mHaloNodes = std::vector<unsigned>(halos_as_set.begin(), halos_as_set.end());
430  //PRINT_VECTOR(mHaloNodes);
431 }
432 
433 
434 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
436 {
437  if (mExchangeHalos)
438  {
439  mpMesh->CalculateNodeExchange(mNodesToSendPerProcess, mNodesToReceivePerProcess);
440  // Note that the following call will not work for a TetrahedralMesh which has
441  // no concept of halo nodes.
442  //mpMesh->GetHaloNodeIndices( mHaloNodes );
443  CalculateHaloNodesFromNodeExchange();
444  unsigned num_halo_nodes = mHaloNodes.size();
445  mHaloCellsDistributed.resize( num_halo_nodes );
446  for (unsigned local_index = 0; local_index < num_halo_nodes; local_index++)
447  {
448  unsigned global_index = mHaloNodes[local_index];
449  // These are all halo nodes, so we use the "GetNodeOrHaloNode" variety of GetNode
450  Node<SPACE_DIM>* p_node = mpMesh->GetNodeOrHaloNode(global_index);
451  mHaloCellsDistributed[local_index] = pCellFactory->CreateCardiacCellForNode(p_node);
452  mHaloCellsDistributed[local_index]->SetUsedInTissueSimulation();
453  mHaloGlobalToLocalIndexMap[global_index] = local_index;
454  }
455  // No need to call FinaliseCellCreation() as halo node cardiac cells will
456  // never be stimulated (their values are communicated from the process that
457  // owns them).
458  }
459 }
460 
461 
462 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
463 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::SolveCellSystems(Vec existingSolution, double time, double nextTime, bool updateVoltage)
464 {
465  if(mHasPurkinje)
466  {
467  // can't do Purkinje and operator splitting
468  assert(!updateVoltage);
469  // The code below assumes Purkinje is are monodomain, so the vector has two stripes.
470  // The assert will fail the first time bidomain purkinje is coded - need to decide what
471  // ordering the three stripes (V, V_purk, phi_e) are in
472  assert(PetscVecTools::GetSize(existingSolution)==2*mpMesh->GetNumNodes());
473  }
474 
475  HeartEventHandler::BeginEvent(HeartEventHandler::SOLVE_ODES);
476 
477  DistributedVector dist_solution = mpDistributedVectorFactory->CreateDistributedVector(existingSolution);
478 
480  // Solve cell models (except purkinje cell models)
482  DistributedVector::Stripe voltage(dist_solution, 0);
483  try
484  {
485  double voltage_before_update;
486  for (DistributedVector::Iterator index = dist_solution.Begin();
487  index != dist_solution.End();
488  ++index)
489  {
490  voltage_before_update = voltage[index];
491  mCellsDistributed[index.Local]->SetVoltage( voltage_before_update );
492 
493  // Added a try-catch here to provide more output to screen when an error occurs.
495  try
496  {
497  if (!updateVoltage)
498  {
499  // solve ODE system at this node.
500  // Note: Voltage is not being updated. The voltage is updated in the PDE solve.
501 #ifndef CHASTE_CVODE
502  mCellsDistributed[index.Local]->ComputeExceptVoltage(time, nextTime);
503 #else
504  // If CVODE is enabled, and this is a CVODE cell
505  // there's a chance we can recover this by doing a reset so put the above call in a try...catch.
506  try
507  {
508  mCellsDistributed[index.Local]->ComputeExceptVoltage(time, nextTime);
509  }
510  catch (Exception &e)
511  {
512  // Try an 'emergency' reset if this is a CVODE cell.
513  // See #2594 for why we think this may be necessary.
514  if(dynamic_cast<AbstractCvodeCell*>(mCellsDistributed[index.Local]))
515  {
516  // Reset the CVODE cell, this leads to a call to CVodeReInit.
517  static_cast<AbstractCvodeCell*>(mCellsDistributed[index.Local])->ResetSolver();
518  mCellsDistributed[index.Local]->ComputeExceptVoltage(time, nextTime);
519  WARNING("Global node " << index.Global << " had an ODE solving problem in t = [" << time <<
520  ", " << nextTime << "] ms. This was fixed by a reset of CVODE, but may suggest PDE time"
521  " step should be reduced, or CVODE tolerances relaxed.");
522  }
523  else
524  {
525  throw e;
526  }
527  }
528 #endif // CHASTE_CVODE
529  }
530  else
531  {
532  // solve, including updating the voltage (for the operator-splitting implementation of the monodomain solver)
533  mCellsDistributed[index.Local]->SolveAndUpdateState(time, nextTime);
534  voltage[index] = mCellsDistributed[index.Local]->GetVoltage();
535  }
536  }
537  catch (Exception &e)
538  {
539  std::cout << std::setprecision(16);
540  std::cout << "Global node " << index.Global << " had problems with ODE solve between "
541  "t = " << time << " and " << nextTime << "ms.\n";
542 
543  std::cout << "Voltage at this node before solve was " << voltage_before_update << "mV\n"
544  "(this SHOULD NOT necessarily be the same as the one in the state variables,\n"
545  "which can be ignored and stay at the initial condition - the voltage is dictated by PDE instead of state variable.)\n";
546 
547  std::cout << "Stimulus current (NB converted to micro-Amps per cm^3) applied here is equal to:\n\t"
548  << mCellsDistributed[index.Local]->GetIntracellularStimulus(time) << " at t = " << time << "ms,\n\t"
549  << mCellsDistributed[index.Local]->GetIntracellularStimulus(nextTime) << " at t = " << nextTime << "ms.\n";
550 
551  std::cout << "Cell model: " << dynamic_cast<AbstractUntemplatedParameterisedSystem*>(mCellsDistributed[index.Local])->GetSystemName() << "\n";
552 
553  std::cout << "All state variables are now:\n";
554  std::vector<double> state_vars = mCellsDistributed[index.Local]->GetStdVecStateVariables();
555  std::vector<std::string> state_var_names = mCellsDistributed[index.Local]->rGetStateVariableNames();
556  for (unsigned i=0; i<state_vars.size(); i++)
557  {
558  std::cout << "\t" << state_var_names[i] << "\t:\t" << state_vars[i] << "\n";
559  }
560  std::cout << std::flush;
561 
562  throw e;
563  }
564  // update the Iionic and stimulus caches
565  UpdateCaches(index.Global, index.Local, nextTime);
566  }
567 
568  if (updateVoltage)
569  {
570  dist_solution.Restore();
571  }
572  }
573  catch (Exception &e)
574  {
576  throw e;
577  }
578 
580  // Solve purkinje cell models
582  if (mHasPurkinje)
583  {
584  DistributedVector::Stripe purkinje_voltage(dist_solution, 1);
585  try
586  {
587  for (DistributedVector::Iterator index = dist_solution.Begin();
588  index != dist_solution.End();
589  ++index)
590  {
591  // overwrite the voltage with the input value
592  mPurkinjeCellsDistributed[index.Local]->SetVoltage( purkinje_voltage[index] );
593 
594  // solve
595  // Note: Voltage is not being updated. The voltage is updated in the PDE solve.
596  mPurkinjeCellsDistributed[index.Local]->ComputeExceptVoltage(time, nextTime);
597 
598  // update the Iionic and stimulus caches
599  UpdatePurkinjeCaches(index.Global, index.Local, nextTime);
600  }
601  }
602  catch (Exception&)
603  {
607 
609  //PetscTools::ReplicateException(true);
610  //throw e;
611  }
612  }
613 
614 
615 
617  HeartEventHandler::EndEvent(HeartEventHandler::SOLVE_ODES);
618 
619  // Communicate new state variable values to halo nodes
620  if (mExchangeHalos)
621  {
622  assert(!mHasPurkinje);
623 
624  for ( unsigned rank_offset = 1; rank_offset < PetscTools::GetNumProcs(); rank_offset++ )
625  {
626  unsigned send_to = (PetscTools::GetMyRank() + rank_offset) % (PetscTools::GetNumProcs());
627  unsigned receive_from = (PetscTools::GetMyRank() + PetscTools::GetNumProcs()- rank_offset ) % (PetscTools::GetNumProcs());
628 
629  unsigned number_of_cells_to_send = mNodesToSendPerProcess[send_to].size();
630  unsigned number_of_cells_to_receive = mNodesToReceivePerProcess[receive_from].size();
631 
632  // Pack send buffer
633  unsigned send_size = 0;
634  for (unsigned i=0; i<number_of_cells_to_send; i++)
635  {
636  unsigned global_cell_index = mNodesToSendPerProcess[send_to][i];
637  send_size += mCellsDistributed[global_cell_index - mpDistributedVectorFactory->GetLow()]->GetNumberOfStateVariables();
638  }
639 
640  boost::scoped_array<double> send_data(new double[send_size]);
641 
642  unsigned send_index = 0;
643  for (unsigned cell = 0; cell < number_of_cells_to_send; cell++)
644  {
645  unsigned global_cell_index = mNodesToSendPerProcess[send_to][cell];
646  AbstractCardiacCellInterface* p_cell = mCellsDistributed[global_cell_index - mpDistributedVectorFactory->GetLow()];
647  std::vector<double> cell_data = p_cell->GetStdVecStateVariables();
648  const unsigned num_state_vars = p_cell->GetNumberOfStateVariables();
649  for (unsigned state_variable = 0; state_variable < num_state_vars; state_variable++)
650  {
651  send_data[send_index++] = cell_data[state_variable];
652  }
653  }
654  // Receive buffer
655  unsigned receive_size = 0;
656  for (unsigned i=0; i<number_of_cells_to_receive; i++)
657  {
658  unsigned halo_cell_index = mHaloGlobalToLocalIndexMap[mNodesToReceivePerProcess[receive_from][i]];
659  receive_size += mHaloCellsDistributed[halo_cell_index]->GetNumberOfStateVariables();
660  }
661 
662  boost::scoped_array<double> receive_data(new double[receive_size]);
663 
664  // Send and receive
665  int ret;
666  MPI_Status status;
667  ret = MPI_Sendrecv(send_data.get(), send_size,
668  MPI_DOUBLE,
669  send_to, 0,
670  receive_data.get(), receive_size,
671  MPI_DOUBLE,
672  receive_from, 0,
673  PETSC_COMM_WORLD, &status);
674  UNUSED_OPT(ret);
675  assert ( ret == MPI_SUCCESS);
676 
677  // Unpack
678  unsigned receive_index = 0;
679  for ( unsigned cell = 0; cell < number_of_cells_to_receive; cell++ )
680  {
681  AbstractCardiacCellInterface* p_cell = mHaloCellsDistributed[mHaloGlobalToLocalIndexMap[mNodesToReceivePerProcess[receive_from][cell]]];
682  const unsigned number_of_state_variables = p_cell->GetNumberOfStateVariables();
683 
684  std::vector<double> cell_data(number_of_state_variables);
685  for (unsigned state_variable = 0; state_variable < number_of_state_variables; state_variable++)
686  {
687  cell_data[state_variable] = receive_data[receive_index++];
688  }
689  p_cell->SetStateVariables(cell_data);
690  }
691  }
692  }
693 
694  HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
695  if ( mDoCacheReplication )
696  {
697  ReplicateCaches();
698  }
699  HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
700 }
701 
702 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
704 {
705  return mIionicCacheReplicated;
706 }
707 
708 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
710 {
711  return mIntracellularStimulusCacheReplicated;
712 }
713 
714 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
716 {
717  EXCEPT_IF_NOT(mHasPurkinje);
718  return mPurkinjeIionicCacheReplicated;
719 }
720 
721 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
723 {
724  EXCEPT_IF_NOT(mHasPurkinje);
725  return mPurkinjeIntracellularStimulusCacheReplicated;
726 }
727 
728 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
729 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::UpdateCaches(unsigned globalIndex, unsigned localIndex, double nextTime)
730 {
731  mIionicCacheReplicated[globalIndex] = mCellsDistributed[localIndex]->GetIIonic();
732  mIntracellularStimulusCacheReplicated[globalIndex] = mCellsDistributed[localIndex]->GetIntracellularStimulus(nextTime);
733 }
734 
735 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
736 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::UpdatePurkinjeCaches(unsigned globalIndex, unsigned localIndex, double nextTime)
737 {
738  assert(mHasPurkinje);
739  mPurkinjeIionicCacheReplicated[globalIndex] = mPurkinjeCellsDistributed[localIndex]->GetIIonic();
740  mPurkinjeIntracellularStimulusCacheReplicated[globalIndex] = mPurkinjeCellsDistributed[localIndex]->GetIntracellularStimulus(nextTime);
741 }
742 
743 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
745 {
746  // ReplicateCaches only needed for SVI (and non-matrix based assembly which is no longer in code)
747  // which is not implemented with Purkinje. See commented code below if introducing this.
748  assert(!mHasPurkinje);
749 
750  mIionicCacheReplicated.Replicate(mpDistributedVectorFactory->GetLow(), mpDistributedVectorFactory->GetHigh());
751  mIntracellularStimulusCacheReplicated.Replicate(mpDistributedVectorFactory->GetLow(), mpDistributedVectorFactory->GetHigh());
752 
753  //if (mHasPurkinje)
754  //{
755  // mPurkinjeIionicCacheReplicated.Replicate(mpDistributedVectorFactory->GetLow(), mpDistributedVectorFactory->GetHigh());
756  // mPurkinjeIntracellularStimulusCacheReplicated.Replicate(mpDistributedVectorFactory->GetLow(), mpDistributedVectorFactory->GetHigh());
757  //}
758 }
759 
760 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
761 const std::vector<AbstractCardiacCellInterface*>& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetCellsDistributed() const
762 {
763  return mCellsDistributed;
764 }
765 
766 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
767 const std::vector<AbstractCardiacCellInterface*>& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetPurkinjeCellsDistributed() const
768 {
769  EXCEPT_IF_NOT(mHasPurkinje);
770  return mPurkinjeCellsDistributed;
771 }
772 
773 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
775 {
776  return mpMesh;
777 }
778 
779 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
781 {
782  assert(pModifier!=NULL);
783  assert(mpConductivityModifier==NULL); // shouldn't be called twice for example, or with two different modifiers (remove this assert
784  // if for whatever reason want to be able to overwrite modifiers)
785  mpConductivityModifier = pModifier;
786 }
787 
788 
790 // Explicit instantiation
792 
793 template class AbstractCardiacTissue<1,1>;
794 template class AbstractCardiacTissue<1,2>;
795 template class AbstractCardiacTissue<1,3>;
796 template class AbstractCardiacTissue<2,2>;
797 template class AbstractCardiacTissue<3,3>;
virtual void FinaliseCellCreation(std::vector< AbstractCardiacCellInterface * > *pCellsDistributed, unsigned lo, unsigned hi)
ReplicatableVector mIntracellularStimulusCacheReplicated
virtual AbstractCardiacCellInterface * CreateCardiacCellForNode(Node< SPACE_DIM > *pNode)
static unsigned GetSize(Vec vector)
void SetConductivityModifier(AbstractConductivityModifier< ELEMENT_DIM, SPACE_DIM > *pModifier)
DistributedVectorFactory * mpDistributedVectorFactory
AbstractCardiacCellInterface * CreatePurkinjeCellForNode(Node< SPACE_DIM > *pNode, AbstractCardiacCellInterface *pCardiacCell)
Definition: Node.hpp:58
ReplicatableVector & rGetPurkinjeIionicCacheReplicated()
static bool ReplicateBool(bool flag)
Definition: PetscTools.cpp:186
AbstractCardiacTissue(AbstractCardiacCellFactory< ELEMENT_DIM, SPACE_DIM > *pCellFactory, bool exchangeHalos=false)
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * GetMesh()
virtual void SolveCellSystems(Vec existingSolution, double time, double nextTime, bool updateVoltage=false)
static std::string GetMeshFilename()
#define EXCEPTION(message)
Definition: Exception.hpp:143
void UpdateCaches(unsigned globalIndex, unsigned localIndex, double nextTime)
ReplicatableVector mPurkinjeIionicCacheReplicated
std::vector< AbstractCardiacCellInterface * > mPurkinjeCellsDistributed
virtual std::vector< double > GetStdVecStateVariables()=0
const std::vector< std::string > & rGetStateVariableNames() const
#define NEVER_REACHED
Definition: Exception.hpp:206
ReplicatableVector mPurkinjeIntracellularStimulusCacheReplicated
std::string GetMeshName() const
void GetConductivityHeterogeneities(std::vector< boost::shared_ptr< AbstractChasteRegion< DIM > > > &conductivitiesHeterogeneityAreas, std::vector< c_vector< double, 3 > > &intraConductivities, std::vector< c_vector< double, 3 > > &extraConductivities) const
virtual void SetStateVariables(const std::vector< double > &rVariables)=0
virtual unsigned GetNumberOfStateVariables() const =0
std::vector< AbstractCardiacCellInterface * > mCellsDistributed
static bool IsSequential()
Definition: PetscTools.cpp:91
AbstractCardiacCellInterface * GetCardiacCellOrHaloCell(unsigned globalIndex)
ReplicatableVector & rGetIntracellularStimulusCacheReplicated()
AbstractCardiacCellInterface * GetPurkinjeCell(unsigned globalIndex)
void SetUpHaloCells(AbstractCardiacCellFactory< ELEMENT_DIM, SPACE_DIM > *pCellFactory)
virtual void FinalisePurkinjeCellCreation(std::vector< AbstractCardiacCellInterface * > *pPurkinjeCellsDistributed, unsigned lo, unsigned hi)
static void ReplicateException(bool flag)
Definition: PetscTools.cpp:198
ReplicatableVector mIionicCacheReplicated
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
void ComputeExceptVoltage(double tStart, double tEnd)
#define EXCEPT_IF_NOT(test)
Definition: Exception.hpp:158
const std::vector< AbstractCardiacCellInterface * > & rGetCellsDistributed() const
virtual const c_matrix< double, SPACE_DIM, SPACE_DIM > & rGetExtracellularConductivityTensor(unsigned elementIndex)
bool Exists() const
Definition: FileFinder.cpp:180
#define UNUSED_OPT(var)
Definition: Exception.hpp:216
void Resize(unsigned size)
const std::vector< AbstractCardiacCellInterface * > & rGetPurkinjeCellsDistributed() const
void SetCacheReplication(bool doCacheReplication)
ReplicatableVector & rGetPurkinjeIntracellularStimulusCacheReplicated()
static unsigned GetNumProcs()
Definition: PetscTools.cpp:108
const c_matrix< double, SPACE_DIM, SPACE_DIM > & rGetIntracellularConductivityTensor(unsigned elementIndex)
static unsigned GetMyRank()
Definition: PetscTools.cpp:114
static std::string GetArchiveDirectory()
AbstractCardiacCellInterface * GetCardiacCell(unsigned globalIndex)
void UpdatePurkinjeCaches(unsigned globalIndex, unsigned localIndex, double nextTime)
static HeartConfig * Instance()
const AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * pGetMesh() const
ReplicatableVector & rGetIionicCacheReplicated()