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