00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036 #include "AbstractCardiacTissue.hpp"
00037
00038 #include <boost/scoped_array.hpp>
00039
00040 #include "DistributedVector.hpp"
00041 #include "AxisymmetricConductivityTensors.hpp"
00042 #include "OrthotropicConductivityTensors.hpp"
00043 #include "Exception.hpp"
00044 #include "ChastePoint.hpp"
00045 #include "AbstractChasteRegion.hpp"
00046 #include "HeartEventHandler.hpp"
00047 #include "PetscTools.hpp"
00048 #include "PetscVecTools.hpp"
00049 #include "AbstractCvodeCell.hpp"
00050 #include "Warnings.hpp"
00051
00052 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00053 AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::AbstractCardiacTissue(
00054 AbstractCardiacCellFactory<ELEMENT_DIM,SPACE_DIM>* pCellFactory,
00055 bool exchangeHalos)
00056 : mpMesh(pCellFactory->GetMesh()),
00057 mpDistributedVectorFactory(mpMesh->GetDistributedVectorFactory()),
00058 mpConductivityModifier(NULL),
00059 mHasPurkinje(false),
00060 mDoCacheReplication(true),
00061 mMeshUnarchived(false),
00062 mExchangeHalos(exchangeHalos)
00063 {
00064
00065 assert(pCellFactory != NULL);
00066 assert(pCellFactory->GetMesh() != NULL);
00067
00068 if (PetscTools::IsSequential())
00069 {
00070
00071 mExchangeHalos = false;
00072 }
00073
00074 unsigned num_local_nodes = mpDistributedVectorFactory->GetLocalOwnership();
00075 bool process_has_no_nodes = (num_local_nodes == 0u);
00076 if (PetscTools::ReplicateBool(process_has_no_nodes))
00077 {
00078
00079
00080
00081
00082
00083 #define COVERAGE_IGNORE
00084
00085 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");
00086 #undef COVERAGE_IGNORE
00087 }
00088 unsigned ownership_range_low = mpDistributedVectorFactory->GetLow();
00089 mCellsDistributed.resize(num_local_nodes);
00090
00091
00092 AbstractPurkinjeCellFactory<ELEMENT_DIM,SPACE_DIM>* p_purkinje_cell_factory
00093 = dynamic_cast<AbstractPurkinjeCellFactory<ELEMENT_DIM,SPACE_DIM>*>(pCellFactory);
00094 if (p_purkinje_cell_factory)
00095 {
00096 mHasPurkinje = true;
00097 mPurkinjeCellsDistributed.resize(num_local_nodes);
00098 }
00099
00101
00103 try
00104 {
00105 for (unsigned local_index = 0; local_index < num_local_nodes; local_index++)
00106 {
00107 unsigned global_index = ownership_range_low + local_index;
00108 Node<SPACE_DIM>* p_node = mpMesh->GetNode(global_index);
00109 mCellsDistributed[local_index] = pCellFactory->CreateCardiacCellForNode(p_node);
00110 mCellsDistributed[local_index]->SetUsedInTissueSimulation();
00111
00112 if (mHasPurkinje)
00113 {
00114 mPurkinjeCellsDistributed[local_index]
00115 = p_purkinje_cell_factory->CreatePurkinjeCellForNode(p_node, mCellsDistributed[local_index]);
00116 mPurkinjeCellsDistributed[local_index]->SetUsedInTissueSimulation();
00117 }
00118 }
00119
00120 pCellFactory->FinaliseCellCreation(&mCellsDistributed,
00121 mpDistributedVectorFactory->GetLow(),
00122 mpDistributedVectorFactory->GetHigh());
00123 if (mHasPurkinje)
00124 {
00125 p_purkinje_cell_factory->FinalisePurkinjeCellCreation(&mPurkinjeCellsDistributed,
00126 mpDistributedVectorFactory->GetLow(),
00127 mpDistributedVectorFactory->GetHigh());
00128 }
00129 }
00130 catch (const Exception& e)
00131 {
00132
00133
00134
00135
00136 PetscTools::ReplicateException(true);
00137
00138
00139
00140
00141
00142 for (std::vector<AbstractCardiacCellInterface*>::iterator cell_iterator = mCellsDistributed.begin();
00143 cell_iterator != mCellsDistributed.end();
00144 ++cell_iterator)
00145 {
00146 delete (*cell_iterator);
00147 }
00148
00149 throw e;
00150 }
00151 PetscTools::ReplicateException(false);
00152
00153
00154 SetUpHaloCells(pCellFactory);
00155
00156 HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
00157 mIionicCacheReplicated.Resize( pCellFactory->GetNumberOfCells() );
00158 mIntracellularStimulusCacheReplicated.Resize( pCellFactory->GetNumberOfCells() );
00159
00160 if (mHasPurkinje)
00161 {
00162 mPurkinjeIionicCacheReplicated.Resize( pCellFactory->GetNumberOfCells() );
00163 mPurkinjeIntracellularStimulusCacheReplicated.Resize( pCellFactory->GetNumberOfCells() );
00164 }
00165 HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
00166
00167 if(HeartConfig::Instance()->IsMeshProvided() && HeartConfig::Instance()->GetLoadMesh())
00168 {
00169 mFibreFilePathNoExtension = HeartConfig::Instance()->GetMeshName();
00170 }
00171 else
00172 {
00173
00174 mFibreFilePathNoExtension = "";
00175 }
00176 CreateIntracellularConductivityTensor();
00177 }
00178
00179
00180 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00181 AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::AbstractCardiacTissue(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh)
00182 : mpMesh(pMesh),
00183 mpDistributedVectorFactory(mpMesh->GetDistributedVectorFactory()),
00184 mHasPurkinje(false),
00185 mDoCacheReplication(true),
00186 mMeshUnarchived(true),
00187 mExchangeHalos(false)
00188 {
00189 mIionicCacheReplicated.Resize(mpDistributedVectorFactory->GetProblemSize());
00190 mIntracellularStimulusCacheReplicated.Resize(mpDistributedVectorFactory->GetProblemSize());
00191
00192 mFibreFilePathNoExtension = ArchiveLocationInfo::GetArchiveDirectory() + ArchiveLocationInfo::GetMeshFilename();
00193 CreateIntracellularConductivityTensor();
00194 }
00195
00196 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00197 AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::~AbstractCardiacTissue()
00198 {
00199
00200 for (std::vector<AbstractCardiacCellInterface*>::iterator iter = mCellsDistributed.begin();
00201 iter != mCellsDistributed.end();
00202 ++iter)
00203 {
00204 delete (*iter);
00205 }
00206
00207
00208 for (std::vector<AbstractCardiacCellInterface*>::iterator iter = mHaloCellsDistributed.begin();
00209 iter != mHaloCellsDistributed.end();
00210 ++iter)
00211 {
00212 delete (*iter);
00213 }
00214
00215 delete mpIntracellularConductivityTensors;
00216
00217
00218 for (std::vector<AbstractCardiacCellInterface*>::iterator iter = mPurkinjeCellsDistributed.begin();
00219 iter != mPurkinjeCellsDistributed.end();
00220 ++iter)
00221 {
00222 delete (*iter);
00223 }
00224
00225
00226 if (mMeshUnarchived)
00227 {
00228 delete mpMesh;
00229 }
00230 }
00231
00232
00233 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00234 bool AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::HasPurkinje()
00235 {
00236 return mHasPurkinje;
00237 }
00238
00239 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00240 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::CreateIntracellularConductivityTensor()
00241 {
00242 HeartEventHandler::BeginEvent(HeartEventHandler::READ_MESH);
00243 mpConfig = HeartConfig::Instance();
00244
00245 if (mpConfig->IsMeshProvided() && mpConfig->GetLoadMesh())
00246 {
00247 assert(mFibreFilePathNoExtension != "");
00248
00249 switch (mpConfig->GetConductivityMedia())
00250 {
00251 case cp::media_type::Orthotropic:
00252 {
00253 mpIntracellularConductivityTensors = new OrthotropicConductivityTensors<ELEMENT_DIM,SPACE_DIM>;
00254 FileFinder ortho_file(mFibreFilePathNoExtension + ".ortho", RelativeTo::AbsoluteOrCwd);
00255 assert(ortho_file.Exists());
00256 mpIntracellularConductivityTensors->SetFibreOrientationFile(ortho_file);
00257 break;
00258 }
00259
00260 case cp::media_type::Axisymmetric:
00261 {
00262 mpIntracellularConductivityTensors = new AxisymmetricConductivityTensors<ELEMENT_DIM,SPACE_DIM>;
00263 FileFinder axi_file(mFibreFilePathNoExtension + ".axi", RelativeTo::AbsoluteOrCwd);
00264 assert(axi_file.Exists());
00265 mpIntracellularConductivityTensors->SetFibreOrientationFile(axi_file);
00266 break;
00267 }
00268
00269 case cp::media_type::NoFibreOrientation:
00271 mpIntracellularConductivityTensors = new OrthotropicConductivityTensors<ELEMENT_DIM,SPACE_DIM>;
00272 break;
00273
00274 default :
00275 NEVER_REACHED;
00276 }
00277 }
00278 else
00279 {
00281 mpIntracellularConductivityTensors = new OrthotropicConductivityTensors<ELEMENT_DIM,SPACE_DIM>;
00282 }
00283
00284 c_vector<double, SPACE_DIM> intra_conductivities;
00285 mpConfig->GetIntracellularConductivities(intra_conductivities);
00286
00287
00288
00289 unsigned num_local_elements = mpMesh->GetNumLocalElements();
00290 std::vector<c_vector<double, SPACE_DIM> > hetero_intra_conductivities;
00291
00292 if (mpConfig->GetConductivityHeterogeneitiesProvided())
00293 {
00294 try
00295 {
00296 assert(hetero_intra_conductivities.size()==0);
00297 hetero_intra_conductivities.resize(num_local_elements, intra_conductivities);
00298 }
00299 catch(std::bad_alloc &r_bad_alloc)
00300 {
00301 #define COVERAGE_IGNORE
00302 std::cout << "Failed to allocate std::vector of size " << num_local_elements << std::endl;
00303 PetscTools::ReplicateException(true);
00304 throw r_bad_alloc;
00305 #undef COVERAGE_IGNORE
00306 }
00307 PetscTools::ReplicateException(false);
00308
00309 std::vector<boost::shared_ptr<AbstractChasteRegion<SPACE_DIM> > > conductivities_heterogeneity_areas;
00310 std::vector< c_vector<double,3> > intra_h_conductivities;
00311 std::vector< c_vector<double,3> > extra_h_conductivities;
00312 HeartConfig::Instance()->GetConductivityHeterogeneities(conductivities_heterogeneity_areas,
00313 intra_h_conductivities,
00314 extra_h_conductivities);
00315
00316 unsigned local_element_index = 0;
00317
00318 for (typename AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>::ElementIterator it = mpMesh->GetElementIteratorBegin();
00319 it != mpMesh->GetElementIteratorEnd();
00320 ++it)
00321 {
00322
00323
00324 ChastePoint<SPACE_DIM> element_centroid(it->CalculateCentroid());
00325 for (unsigned region_index=0; region_index< conductivities_heterogeneity_areas.size(); region_index++)
00326 {
00327 if ( conductivities_heterogeneity_areas[region_index]->DoesContain(element_centroid) )
00328 {
00329
00330 for (unsigned i=0; i<SPACE_DIM; i++)
00331 {
00332 hetero_intra_conductivities[local_element_index][i] = intra_h_conductivities[region_index][i];
00333 }
00334 }
00335 }
00336 local_element_index++;
00337 }
00338
00339 mpIntracellularConductivityTensors->SetNonConstantConductivities(&hetero_intra_conductivities);
00340 }
00341 else
00342 {
00343 mpIntracellularConductivityTensors->SetConstantConductivities(intra_conductivities);
00344 }
00345
00346 mpIntracellularConductivityTensors->Init(this->mpMesh);
00347 HeartEventHandler::EndEvent(HeartEventHandler::READ_MESH);
00348 }
00349
00350
00351 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00352 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::SetCacheReplication(bool doCacheReplication)
00353 {
00354 mDoCacheReplication = doCacheReplication;
00355 }
00356
00357 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00358 bool AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::GetDoCacheReplication()
00359 {
00360 return mDoCacheReplication;
00361 }
00362
00363 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00364 const c_matrix<double, SPACE_DIM, SPACE_DIM>& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetIntracellularConductivityTensor(unsigned elementIndex)
00365 {
00366 assert( mpIntracellularConductivityTensors);
00367 if (mpConductivityModifier==NULL)
00368 {
00369 return (*mpIntracellularConductivityTensors)[elementIndex];
00370 }
00371 else
00372 {
00373 return mpConductivityModifier->rGetModifiedConductivityTensor(elementIndex, (*mpIntracellularConductivityTensors)[elementIndex], 0u);
00374 }
00375 }
00376
00377 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00378 const c_matrix<double, SPACE_DIM, SPACE_DIM>& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetExtracellularConductivityTensor(unsigned elementIndex)
00379 {
00380 EXCEPTION("Monodomain tissues do not have extracellular conductivity tensors.");
00381 }
00382
00383 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00384 AbstractCardiacCellInterface* AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::GetCardiacCell( unsigned globalIndex )
00385 {
00386 assert(mpDistributedVectorFactory->GetLow() <= globalIndex &&
00387 globalIndex < mpDistributedVectorFactory->GetHigh());
00388 return mCellsDistributed[globalIndex - mpDistributedVectorFactory->GetLow()];
00389 }
00390
00391 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00392 AbstractCardiacCellInterface* AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::GetPurkinjeCell( unsigned globalIndex )
00393 {
00394 assert(mpDistributedVectorFactory->GetLow() <= globalIndex &&
00395 globalIndex < mpDistributedVectorFactory->GetHigh());
00396 EXCEPT_IF_NOT(mHasPurkinje);
00397 return mPurkinjeCellsDistributed[globalIndex - mpDistributedVectorFactory->GetLow()];
00398 }
00399
00400 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00401 AbstractCardiacCellInterface* AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::GetCardiacCellOrHaloCell( unsigned globalIndex )
00402 {
00403 std::map<unsigned, unsigned>::const_iterator node_position;
00404
00405 if ((node_position=mHaloGlobalToLocalIndexMap.find(globalIndex)) != mHaloGlobalToLocalIndexMap.end())
00406 {
00407
00408 return mHaloCellsDistributed[node_position->second];
00409 }
00410
00411 if ( mpDistributedVectorFactory->IsGlobalIndexLocal(globalIndex) )
00412 {
00413
00414 return mCellsDistributed[globalIndex - mpDistributedVectorFactory->GetLow()];
00415 }
00416
00417 EXCEPTION("Requested node/halo " << globalIndex << " does not belong to processor " << PetscTools::GetMyRank());
00418 }
00419
00420
00421 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00422 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::CalculateHaloNodesFromNodeExchange()
00423 {
00424 std::set<unsigned> halos_as_set;
00425 for (unsigned proc=0; proc<PetscTools::GetNumProcs(); proc++)
00426 {
00427 halos_as_set.insert(mNodesToReceivePerProcess[proc].begin(), mNodesToReceivePerProcess[proc].end());
00428 }
00429 mHaloNodes = std::vector<unsigned>(halos_as_set.begin(), halos_as_set.end());
00430
00431 }
00432
00433
00434 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00435 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::SetUpHaloCells(AbstractCardiacCellFactory<ELEMENT_DIM,SPACE_DIM>* pCellFactory)
00436 {
00437 if (mExchangeHalos)
00438 {
00439 mpMesh->CalculateNodeExchange(mNodesToSendPerProcess, mNodesToReceivePerProcess);
00440
00441
00442 CalculateHaloNodesFromNodeExchange();
00443 unsigned num_halo_nodes = mHaloNodes.size();
00444 mHaloCellsDistributed.resize( num_halo_nodes );
00445 for (unsigned local_index = 0; local_index < num_halo_nodes; local_index++)
00446 {
00447 unsigned global_index = mHaloNodes[local_index];
00448
00449 Node<SPACE_DIM>* p_node = mpMesh->GetNodeOrHaloNode(global_index);
00450 mHaloCellsDistributed[local_index] = pCellFactory->CreateCardiacCellForNode(p_node);
00451 mHaloCellsDistributed[local_index]->SetUsedInTissueSimulation();
00452 mHaloGlobalToLocalIndexMap[global_index] = local_index;
00453 }
00454
00455
00456
00457 }
00458 }
00459
00460
00461 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00462 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::SolveCellSystems(Vec existingSolution, double time, double nextTime, bool updateVoltage)
00463 {
00464 if(mHasPurkinje)
00465 {
00466
00467 assert(!updateVoltage);
00468
00469
00470
00471 assert(PetscVecTools::GetSize(existingSolution)==2*mpMesh->GetNumNodes());
00472 }
00473
00474 HeartEventHandler::BeginEvent(HeartEventHandler::SOLVE_ODES);
00475
00476 DistributedVector dist_solution = mpDistributedVectorFactory->CreateDistributedVector(existingSolution);
00477
00479
00481 DistributedVector::Stripe voltage(dist_solution, 0);
00482 try
00483 {
00484 double voltage_before_update;
00485 for (DistributedVector::Iterator index = dist_solution.Begin();
00486 index != dist_solution.End();
00487 ++index)
00488 {
00489 voltage_before_update = voltage[index];
00490 mCellsDistributed[index.Local]->SetVoltage( voltage_before_update );
00491
00492
00494 try
00495 {
00496 if (!updateVoltage)
00497 {
00498
00499
00500 #ifndef CHASTE_CVODE
00501 mCellsDistributed[index.Local]->ComputeExceptVoltage(time, nextTime);
00502 #else
00503
00504
00505 try
00506 {
00507 mCellsDistributed[index.Local]->ComputeExceptVoltage(time, nextTime);
00508 }
00509 catch (Exception &e)
00510 {
00511
00512
00513 if(dynamic_cast<AbstractCvodeCell*>(mCellsDistributed[index.Local]))
00514 {
00515
00516 static_cast<AbstractCvodeCell*>(mCellsDistributed[index.Local])->ResetSolver();
00517 mCellsDistributed[index.Local]->ComputeExceptVoltage(time, nextTime);
00518 WARNING("Global node " << index.Global << " had an ODE solving problem in t = [" << time <<
00519 ", " << nextTime << "] ms. This was fixed by a reset of CVODE, but may suggest PDE time"
00520 " step should be reduced, or CVODE tolerances relaxed.");
00521 }
00522 else
00523 {
00524 throw e;
00525 }
00526 }
00527 #endif // CHASTE_CVODE
00528 }
00529 else
00530 {
00531
00532 mCellsDistributed[index.Local]->SolveAndUpdateState(time, nextTime);
00533 voltage[index] = mCellsDistributed[index.Local]->GetVoltage();
00534 }
00535 }
00536 catch (Exception &e)
00537 {
00538 std::cout << std::setprecision(16);
00539 std::cout << "Global node " << index.Global << " had problems with ODE solve between "
00540 "t = " << time << " and " << nextTime << "ms.\n";
00541
00542 std::cout << "Voltage at this node before solve was " << voltage_before_update << "mV\n"
00543 "(this SHOULD NOT necessarily be the same as the one in the state variables,\n"
00544 "which can be ignored and stay at the initial condition - the voltage is dictated by PDE instead of state variable.)\n";
00545
00546 std::cout << "Stimulus current (NB converted to micro-Amps per cm^3) applied here is equal to:\n\t"
00547 << mCellsDistributed[index.Local]->GetIntracellularStimulus(time) << " at t = " << time << "ms,\n\t"
00548 << mCellsDistributed[index.Local]->GetIntracellularStimulus(nextTime) << " at t = " << nextTime << "ms.\n";
00549
00550 std::cout << "Cell model: " << dynamic_cast<AbstractUntemplatedParameterisedSystem*>(mCellsDistributed[index.Local])->GetSystemName() << "\n";
00551
00552 std::cout << "All state variables are now:\n";
00553 std::vector<double> state_vars = mCellsDistributed[index.Local]->GetStdVecStateVariables();
00554 std::vector<std::string> state_var_names = mCellsDistributed[index.Local]->rGetStateVariableNames();
00555 for (unsigned i=0; i<state_vars.size(); i++)
00556 {
00557 std::cout << "\t" << state_var_names[i] << "\t:\t" << state_vars[i] << "\n";
00558 }
00559 std::cout << std::flush;
00560
00561 throw e;
00562 }
00563
00564 UpdateCaches(index.Global, index.Local, nextTime);
00565 }
00566
00567 if (updateVoltage)
00568 {
00569 dist_solution.Restore();
00570 }
00571 }
00572 catch (Exception &e)
00573 {
00574 PetscTools::ReplicateException(true);
00575 throw e;
00576 }
00577
00579
00581 if (mHasPurkinje)
00582 {
00583 DistributedVector::Stripe purkinje_voltage(dist_solution, 1);
00584 try
00585 {
00586 for (DistributedVector::Iterator index = dist_solution.Begin();
00587 index != dist_solution.End();
00588 ++index)
00589 {
00590
00591 mPurkinjeCellsDistributed[index.Local]->SetVoltage( purkinje_voltage[index] );
00592
00593
00594
00595 mPurkinjeCellsDistributed[index.Local]->ComputeExceptVoltage(time, nextTime);
00596
00597
00598 UpdatePurkinjeCaches(index.Global, index.Local, nextTime);
00599 }
00600 }
00601 catch (Exception&)
00602 {
00606
00607 NEVER_REACHED;
00608
00609
00610 }
00611 }
00612
00613
00614
00615 PetscTools::ReplicateException(false);
00616 HeartEventHandler::EndEvent(HeartEventHandler::SOLVE_ODES);
00617
00618
00619 if (mExchangeHalos)
00620 {
00621 assert(!mHasPurkinje);
00622
00623 for ( unsigned rank_offset = 1; rank_offset < PetscTools::GetNumProcs(); rank_offset++ )
00624 {
00625 unsigned send_to = (PetscTools::GetMyRank() + rank_offset) % (PetscTools::GetNumProcs());
00626 unsigned receive_from = (PetscTools::GetMyRank() + PetscTools::GetNumProcs()- rank_offset ) % (PetscTools::GetNumProcs());
00627
00628 unsigned number_of_cells_to_send = mNodesToSendPerProcess[send_to].size();
00629 unsigned number_of_cells_to_receive = mNodesToReceivePerProcess[receive_from].size();
00630
00631
00632 unsigned send_size = 0;
00633 for (unsigned i=0; i<number_of_cells_to_send; i++)
00634 {
00635 unsigned global_cell_index = mNodesToSendPerProcess[send_to][i];
00636 send_size += mCellsDistributed[global_cell_index - mpDistributedVectorFactory->GetLow()]->GetNumberOfStateVariables();
00637 }
00638
00639 boost::scoped_array<double> send_data(new double[send_size]);
00640
00641 unsigned send_index = 0;
00642 for (unsigned cell = 0; cell < number_of_cells_to_send; cell++)
00643 {
00644 unsigned global_cell_index = mNodesToSendPerProcess[send_to][cell];
00645 AbstractCardiacCellInterface* p_cell = mCellsDistributed[global_cell_index - mpDistributedVectorFactory->GetLow()];
00646 std::vector<double> cell_data = p_cell->GetStdVecStateVariables();
00647 const unsigned num_state_vars = p_cell->GetNumberOfStateVariables();
00648 for (unsigned state_variable = 0; state_variable < num_state_vars; state_variable++)
00649 {
00650 send_data[send_index++] = cell_data[state_variable];
00651 }
00652 }
00653
00654 unsigned receive_size = 0;
00655 for (unsigned i=0; i<number_of_cells_to_receive; i++)
00656 {
00657 unsigned halo_cell_index = mHaloGlobalToLocalIndexMap[mNodesToReceivePerProcess[receive_from][i]];
00658 receive_size += mHaloCellsDistributed[halo_cell_index]->GetNumberOfStateVariables();
00659 }
00660
00661 boost::scoped_array<double> receive_data(new double[receive_size]);
00662
00663
00664 int ret;
00665 MPI_Status status;
00666 ret = MPI_Sendrecv(send_data.get(), send_size,
00667 MPI_DOUBLE,
00668 send_to, 0,
00669 receive_data.get(), receive_size,
00670 MPI_DOUBLE,
00671 receive_from, 0,
00672 PETSC_COMM_WORLD, &status);
00673 UNUSED_OPT(ret);
00674 assert ( ret == MPI_SUCCESS);
00675
00676
00677 unsigned receive_index = 0;
00678 for ( unsigned cell = 0; cell < number_of_cells_to_receive; cell++ )
00679 {
00680 AbstractCardiacCellInterface* p_cell = mHaloCellsDistributed[mHaloGlobalToLocalIndexMap[mNodesToReceivePerProcess[receive_from][cell]]];
00681 const unsigned number_of_state_variables = p_cell->GetNumberOfStateVariables();
00682
00683 std::vector<double> cell_data(number_of_state_variables);
00684 for (unsigned state_variable = 0; state_variable < number_of_state_variables; state_variable++)
00685 {
00686 cell_data[state_variable] = receive_data[receive_index++];
00687 }
00688 p_cell->SetStateVariables(cell_data);
00689 }
00690 }
00691 }
00692
00693 HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
00694 if ( mDoCacheReplication )
00695 {
00696 ReplicateCaches();
00697 }
00698 HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
00699 }
00700
00701 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00702 ReplicatableVector& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetIionicCacheReplicated()
00703 {
00704 return mIionicCacheReplicated;
00705 }
00706
00707 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00708 ReplicatableVector& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetIntracellularStimulusCacheReplicated()
00709 {
00710 return mIntracellularStimulusCacheReplicated;
00711 }
00712
00713 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00714 ReplicatableVector& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetPurkinjeIionicCacheReplicated()
00715 {
00716 EXCEPT_IF_NOT(mHasPurkinje);
00717 return mPurkinjeIionicCacheReplicated;
00718 }
00719
00720 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00721 ReplicatableVector& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetPurkinjeIntracellularStimulusCacheReplicated()
00722 {
00723 EXCEPT_IF_NOT(mHasPurkinje);
00724 return mPurkinjeIntracellularStimulusCacheReplicated;
00725 }
00726
00727 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00728 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::UpdateCaches(unsigned globalIndex, unsigned localIndex, double nextTime)
00729 {
00730 mIionicCacheReplicated[globalIndex] = mCellsDistributed[localIndex]->GetIIonic();
00731 mIntracellularStimulusCacheReplicated[globalIndex] = mCellsDistributed[localIndex]->GetIntracellularStimulus(nextTime);
00732 }
00733
00734 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00735 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::UpdatePurkinjeCaches(unsigned globalIndex, unsigned localIndex, double nextTime)
00736 {
00737 assert(mHasPurkinje);
00738 mPurkinjeIionicCacheReplicated[globalIndex] = mPurkinjeCellsDistributed[localIndex]->GetIIonic();
00739 mPurkinjeIntracellularStimulusCacheReplicated[globalIndex] = mPurkinjeCellsDistributed[localIndex]->GetIntracellularStimulus(nextTime);
00740 }
00741
00742 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00743 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::ReplicateCaches()
00744 {
00745
00746
00747 assert(!mHasPurkinje);
00748
00749 mIionicCacheReplicated.Replicate(mpDistributedVectorFactory->GetLow(), mpDistributedVectorFactory->GetHigh());
00750 mIntracellularStimulusCacheReplicated.Replicate(mpDistributedVectorFactory->GetLow(), mpDistributedVectorFactory->GetHigh());
00751
00752
00753
00754
00755
00756
00757 }
00758
00759 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00760 const std::vector<AbstractCardiacCellInterface*>& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetCellsDistributed() const
00761 {
00762 return mCellsDistributed;
00763 }
00764
00765 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00766 const std::vector<AbstractCardiacCellInterface*>& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetPurkinjeCellsDistributed() const
00767 {
00768 EXCEPT_IF_NOT(mHasPurkinje);
00769 return mPurkinjeCellsDistributed;
00770 }
00771
00772 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00773 const AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::pGetMesh() const
00774 {
00775 return mpMesh;
00776 }
00777
00778 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00779 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::SetConductivityModifier(AbstractConductivityModifier<ELEMENT_DIM,SPACE_DIM>* pModifier)
00780 {
00781 assert(pModifier!=NULL);
00782 assert(mpConductivityModifier==NULL);
00783
00784 mpConductivityModifier = pModifier;
00785 }
00786
00787
00789
00791
00792 template class AbstractCardiacTissue<1,1>;
00793 template class AbstractCardiacTissue<1,2>;
00794 template class AbstractCardiacTissue<1,3>;
00795 template class AbstractCardiacTissue<2,2>;
00796 template class AbstractCardiacTissue<3,3>;