AbstractCardiacTissue.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2011
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00026 
00027 */
00028 
00029 #include "AbstractCardiacTissue.hpp"
00030 
00031 #include "DistributedVector.hpp"
00032 #include "AxisymmetricConductivityTensors.hpp"
00033 #include "OrthotropicConductivityTensors.hpp"
00034 #include "Exception.hpp"
00035 #include "ChastePoint.hpp"
00036 #include "AbstractChasteRegion.hpp"
00037 #include "HeartEventHandler.hpp"
00038 #include "PetscTools.hpp"
00039 #include "PetscVecTools.hpp"
00040 
00041 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00042 AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::AbstractCardiacTissue(
00043             AbstractCardiacCellFactory<ELEMENT_DIM,SPACE_DIM>* pCellFactory,
00044             bool exchangeHalos)
00045     : mpMesh(pCellFactory->GetMesh()),
00046       mpDistributedVectorFactory(mpMesh->GetDistributedVectorFactory()),
00047       mpConductivityModifier(NULL),
00048       mHasPurkinje(false),
00049       mDoCacheReplication(true),
00050       mMeshUnarchived(false),
00051       mExchangeHalos(exchangeHalos)
00052 {
00053     //This constructor is called from the Initialise() method of the CardiacProblem class
00054     assert(pCellFactory != NULL);
00055     assert(pCellFactory->GetMesh() != NULL);
00056 
00057     if (PetscTools::IsSequential())
00058     {
00059         //Remove the request for a halo exchange
00060         mExchangeHalos = false;
00061     }
00062 
00063     unsigned num_local_nodes = mpDistributedVectorFactory->GetLocalOwnership();
00064     unsigned ownership_range_low = mpDistributedVectorFactory->GetLow();
00065     mCellsDistributed.resize(num_local_nodes);
00066 
00067     // Figure out if we're dealing with Purkinje
00068     AbstractPurkinjeCellFactory<ELEMENT_DIM,SPACE_DIM>* p_purkinje_cell_factory
00069        = dynamic_cast<AbstractPurkinjeCellFactory<ELEMENT_DIM,SPACE_DIM>*>(pCellFactory);
00070     if (p_purkinje_cell_factory)
00071     {
00072         mHasPurkinje = true;
00073         mPurkinjeCellsDistributed.resize(num_local_nodes);
00074     }
00075 
00077     // Set up cells
00079     try
00080     {
00081         for (unsigned local_index = 0; local_index < num_local_nodes; local_index++)
00082         {
00083             unsigned global_index = ownership_range_low + local_index;
00084             mCellsDistributed[local_index] = pCellFactory->CreateCardiacCellForNode(global_index);
00085             mCellsDistributed[local_index]->SetUsedInTissueSimulation();
00086 
00087             if (mHasPurkinje)
00088             {
00089                 mPurkinjeCellsDistributed[local_index] = p_purkinje_cell_factory->CreatePurkinjeCellForNode(global_index);
00090                 mPurkinjeCellsDistributed[local_index]->SetUsedInTissueSimulation();
00091             }
00092         }
00093 
00094         pCellFactory->FinaliseCellCreation(&mCellsDistributed,
00095                                            mpDistributedVectorFactory->GetLow(),
00096                                            mpDistributedVectorFactory->GetHigh());
00097         if (mHasPurkinje)
00098         {
00099             p_purkinje_cell_factory->FinalisePurkinjeCellCreation(&mPurkinjeCellsDistributed,
00100                                                                   mpDistributedVectorFactory->GetLow(),
00101                                                                   mpDistributedVectorFactory->GetHigh());
00102         }
00103     }
00104     catch (const Exception& e)
00105     {
00106         // This catch statement is quite tricky to cover, but it is actually done in TestCardiacSimulation::TestMono1dSodiumBlockBySettingNamedParameter
00107 
00108         // Errors thrown creating cells will often be process-specific
00109         PetscTools::ReplicateException(true);
00110 
00111         // Delete cells
00112         // Should really do this for other processes too, but this is all we need
00113         // to get memory testing to pass, and leaking when we're about to die isn't
00114         // that bad!
00115         for (std::vector<AbstractCardiacCell*>::iterator cell_iterator = mCellsDistributed.begin();
00116              cell_iterator != mCellsDistributed.end();
00117              ++cell_iterator)
00118         {
00119             delete (*cell_iterator);
00120         }
00121 
00122         throw e;
00123     }
00124     PetscTools::ReplicateException(false);
00125 
00126     // Halo nodes (if required)
00127     SetUpHaloCells(pCellFactory);
00128 
00129     HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
00130     mIionicCacheReplicated.Resize( pCellFactory->GetNumberOfCells() );
00131     mIntracellularStimulusCacheReplicated.Resize( pCellFactory->GetNumberOfCells() );
00132     if (mHasPurkinje)
00133     {
00134         mPurkinjeIionicCacheReplicated.Resize( pCellFactory->GetNumberOfCells() );
00135     }
00136     HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
00137 
00138     if(HeartConfig::Instance()->IsMeshProvided() && HeartConfig::Instance()->GetLoadMesh())
00139     {
00140         mFibreFilePathNoExtension = "./" + HeartConfig::Instance()->GetMeshName();
00141     }
00142     else
00143     {
00144         // As of 10671 fibre orientation can only be defined when loading a mesh from disc.
00145         mFibreFilePathNoExtension = "";
00146     }
00147     CreateIntracellularConductivityTensor();
00148 }
00149 
00150 // Constructor used for archiving
00151 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00152 AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::AbstractCardiacTissue(AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* pMesh)
00153     : mpMesh(pMesh),
00154       mpDistributedVectorFactory(mpMesh->GetDistributedVectorFactory()),
00155       mHasPurkinje(false),
00156       mDoCacheReplication(true),
00157       mMeshUnarchived(true),
00158       mExchangeHalos(false)
00159 {
00160     mIionicCacheReplicated.Resize(mpDistributedVectorFactory->GetProblemSize());
00161     mIntracellularStimulusCacheReplicated.Resize(mpDistributedVectorFactory->GetProblemSize());
00162 
00163     mFibreFilePathNoExtension = ArchiveLocationInfo::GetArchiveDirectory() + ArchiveLocationInfo::GetMeshFilename();
00164     CreateIntracellularConductivityTensor();
00165 }
00166 
00167 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00168 AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::~AbstractCardiacTissue()
00169 {
00170     // Delete cells
00171     for (std::vector<AbstractCardiacCell*>::iterator iter = mCellsDistributed.begin();
00172          iter != mCellsDistributed.end();
00173          ++iter)
00174     {
00175         delete (*iter);
00176     }
00177 
00178     // Delete cells for halo nodes
00179     for (std::vector<AbstractCardiacCell*>::iterator iter = mHaloCellsDistributed.begin();
00180          iter != mHaloCellsDistributed.end();
00181          ++iter)
00182     {
00183         delete (*iter);
00184     }
00185 
00186     delete mpIntracellularConductivityTensors;
00187 
00188     // Delete Purkinje cells
00189     for (std::vector<AbstractCardiacCell*>::iterator iter = mPurkinjeCellsDistributed.begin();
00190          iter != mPurkinjeCellsDistributed.end();
00191          ++iter)
00192     {
00193         delete (*iter);
00194     }
00195 
00196     // If the mesh was unarchived we need to free it explicitly.
00197     if (mMeshUnarchived)
00198     {
00199         delete mpMesh;
00200     }
00201 }
00202 
00203 
00204 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00205 bool AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::HasPurkinje()
00206 {
00207     return mHasPurkinje;
00208 }
00209 
00210 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00211 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::CreateIntracellularConductivityTensor()
00212 {
00213     HeartEventHandler::BeginEvent(HeartEventHandler::READ_MESH);
00214     mpConfig = HeartConfig::Instance();
00215 
00216     if (mpConfig->IsMeshProvided() && mpConfig->GetLoadMesh())
00217     {
00218         assert(mFibreFilePathNoExtension != "");
00219 
00220         switch (mpConfig->GetConductivityMedia())
00221         {
00222             case cp::media_type::Orthotropic:
00223             {
00224                 mpIntracellularConductivityTensors = new OrthotropicConductivityTensors<ELEMENT_DIM,SPACE_DIM>;
00225                 FileFinder ortho_file(mFibreFilePathNoExtension + ".ortho", RelativeTo::AbsoluteOrCwd);
00226                 assert(ortho_file.Exists());
00227                 mpIntracellularConductivityTensors->SetFibreOrientationFile(ortho_file);
00228                 break;
00229             }
00230 
00231             case cp::media_type::Axisymmetric:
00232             {
00233                 mpIntracellularConductivityTensors = new AxisymmetricConductivityTensors<ELEMENT_DIM,SPACE_DIM>;
00234                 FileFinder axi_file(mFibreFilePathNoExtension + ".axi", RelativeTo::AbsoluteOrCwd);
00235                 assert(axi_file.Exists());
00236                 mpIntracellularConductivityTensors->SetFibreOrientationFile(axi_file);
00237                 break;
00238             }
00239 
00240             case cp::media_type::NoFibreOrientation:
00242                 mpIntracellularConductivityTensors = new OrthotropicConductivityTensors<ELEMENT_DIM,SPACE_DIM>;
00243                 break;
00244 
00245             default :
00246                 NEVER_REACHED;
00247         }
00248     }
00249     else // Slab defined in config file or SetMesh() called; no fibre orientation assumed
00250     {
00252         mpIntracellularConductivityTensors = new OrthotropicConductivityTensors<ELEMENT_DIM,SPACE_DIM>;
00253     }
00254 
00255     c_vector<double, SPACE_DIM> intra_conductivities;
00256     mpConfig->GetIntracellularConductivities(intra_conductivities);
00257 
00258     // this definition must be here (and not inside the if statement) because SetNonConstantConductivities() will keep
00259     // a pointer to it and we don't want it to go out of scope before Init() is called
00260     unsigned num_local_elements = mpMesh->GetNumLocalElements();
00261     std::vector<c_vector<double, SPACE_DIM> > hetero_intra_conductivities;
00262 
00263     if (mpConfig->GetConductivityHeterogeneitiesProvided())
00264     {
00265         try
00266         {
00267             assert(hetero_intra_conductivities.size()==0);
00268             hetero_intra_conductivities.resize(num_local_elements, intra_conductivities);
00269         }
00270         catch(std::bad_alloc &r_bad_alloc)
00271         {
00272 #define COVERAGE_IGNORE
00273             std::cout << "Failed to allocate std::vector of size " << num_local_elements << std::endl;
00274             PetscTools::ReplicateException(true);
00275             throw r_bad_alloc;
00276 #undef COVERAGE_IGNORE
00277         }
00278         PetscTools::ReplicateException(false);
00279 
00280         std::vector<boost::shared_ptr<AbstractChasteRegion<SPACE_DIM> > > conductivities_heterogeneity_areas;
00281         std::vector< c_vector<double,3> > intra_h_conductivities;
00282         std::vector< c_vector<double,3> > extra_h_conductivities;
00283         HeartConfig::Instance()->GetConductivityHeterogeneities(conductivities_heterogeneity_areas,
00284                                                                 intra_h_conductivities,
00285                                                                 extra_h_conductivities);
00286 
00287         unsigned local_element_index = 0;
00288 
00289         for (typename AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>::ElementIterator it = mpMesh->GetElementIteratorBegin();
00290              it != mpMesh->GetElementIteratorEnd();
00291              ++it)
00292         {
00293 //            unsigned element_index = it->GetIndex();
00294             // if element centroid is contained in the region
00295             ChastePoint<SPACE_DIM> element_centroid(it->CalculateCentroid());
00296             for (unsigned region_index=0; region_index< conductivities_heterogeneity_areas.size(); region_index++)
00297             {
00298                 if ( conductivities_heterogeneity_areas[region_index]->DoesContain(element_centroid) )
00299                 {
00300                     //We don't use ublas vector assignment here, because we might be getting a subvector of a 3-vector
00301                     for (unsigned i=0; i<SPACE_DIM; i++)
00302                     {
00303                         hetero_intra_conductivities[local_element_index][i] = intra_h_conductivities[region_index][i];
00304                     }
00305                 }
00306             }
00307             local_element_index++;
00308         }
00309 
00310         mpIntracellularConductivityTensors->SetNonConstantConductivities(&hetero_intra_conductivities);
00311     }
00312     else
00313     {
00314         mpIntracellularConductivityTensors->SetConstantConductivities(intra_conductivities);
00315     }
00316 
00317     mpIntracellularConductivityTensors->Init(this->mpMesh);
00318     HeartEventHandler::EndEvent(HeartEventHandler::READ_MESH);
00319 }
00320 
00321 
00322 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00323 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::SetCacheReplication(bool doCacheReplication)
00324 {
00325     mDoCacheReplication = doCacheReplication;
00326 }
00327 
00328 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00329 bool AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::GetDoCacheReplication()
00330 {
00331     return mDoCacheReplication;
00332 }
00333 
00334 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00335 const c_matrix<double, SPACE_DIM, SPACE_DIM>& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetIntracellularConductivityTensor(unsigned elementIndex)
00336 {
00337     assert( mpIntracellularConductivityTensors);
00338     if(mpConductivityModifier==NULL)
00339     {
00340         return (*mpIntracellularConductivityTensors)[elementIndex];
00341     }
00342     else
00343     {
00344         return mpConductivityModifier->rGetModifiedConductivityTensor(elementIndex, (*mpIntracellularConductivityTensors)[elementIndex]);
00345     }
00346 }
00347 
00348 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00349 const c_matrix<double, SPACE_DIM, SPACE_DIM>& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetExtracellularConductivityTensor(unsigned elementIndex)
00350 {
00351      EXCEPTION("Monodomain tissues do not have extracellular conductivity tensors.");
00352 }
00353 
00354 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00355 AbstractCardiacCell* AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::GetCardiacCell( unsigned globalIndex )
00356 {
00357     assert(mpDistributedVectorFactory->GetLow() <= globalIndex &&
00358            globalIndex < mpDistributedVectorFactory->GetHigh());
00359     return mCellsDistributed[globalIndex - mpDistributedVectorFactory->GetLow()];
00360 }
00361 
00362 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00363 AbstractCardiacCell* AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::GetPurkinjeCell( unsigned globalIndex )
00364 {
00365     assert(mpDistributedVectorFactory->GetLow() <= globalIndex &&
00366            globalIndex < mpDistributedVectorFactory->GetHigh());
00367     EXCEPT_IF_NOT(mHasPurkinje);
00368     return mPurkinjeCellsDistributed[globalIndex - mpDistributedVectorFactory->GetLow()];
00369 }
00370 
00371 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00372 AbstractCardiacCell* AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::GetCardiacCellOrHaloCell( unsigned globalIndex )
00373 {
00374     std::map<unsigned, unsigned>::const_iterator node_position;
00375     // First search the halo
00376     if ((node_position=mHaloGlobalToLocalIndexMap.find(globalIndex)) != mHaloGlobalToLocalIndexMap.end())
00377     {
00378         //Found a halo node
00379         return mHaloCellsDistributed[node_position->second];
00380     }
00381     // Then search the owned node
00382     if ( mpDistributedVectorFactory->IsGlobalIndexLocal(globalIndex)  )
00383     {
00384         //Found an owned node
00385         return mCellsDistributed[globalIndex - mpDistributedVectorFactory->GetLow()];
00386     }
00387     // Not here
00388     EXCEPTION("Requested node/halo " << globalIndex << " does not belong to processor " << PetscTools::GetMyRank());
00389 }
00390 
00391 
00392 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00393 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::CalculateHaloNodesFromNodeExchange()
00394 {
00395     std::set<unsigned> halos_as_set;
00396     for (unsigned proc=0; proc<PetscTools::GetNumProcs(); proc++)
00397     {
00398         halos_as_set.insert(mNodesToReceivePerProcess[proc].begin(), mNodesToReceivePerProcess[proc].end());
00399     }
00400     mHaloNodes = std::vector<unsigned>(halos_as_set.begin(), halos_as_set.end());
00401     //PRINT_VECTOR(mHaloNodes);
00402 }
00403 
00404 
00405 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00406 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::SetUpHaloCells(AbstractCardiacCellFactory<ELEMENT_DIM,SPACE_DIM>* pCellFactory)
00407 {
00408     if (mExchangeHalos)
00409     {
00410         mpMesh->CalculateNodeExchange(mNodesToSendPerProcess, mNodesToReceivePerProcess);
00411         //Note that the following call will not work for a TetrahedralMesh which has no concept of halo nodes.
00412         //mpMesh->GetHaloNodeIndices( mHaloNodes );
00413         CalculateHaloNodesFromNodeExchange();
00414         unsigned num_halo_nodes = mHaloNodes.size();
00415         mHaloCellsDistributed.resize( num_halo_nodes );
00416 
00417         try
00418         {
00419             for (unsigned local_index = 0; local_index < num_halo_nodes; local_index++)
00420             {
00421                 unsigned global_index = mHaloNodes[local_index];
00422                 mHaloCellsDistributed[local_index] = pCellFactory->CreateCardiacCellForNode(global_index);
00423                 mHaloCellsDistributed[local_index]->SetUsedInTissueSimulation();
00424                 mHaloGlobalToLocalIndexMap[global_index] = local_index;
00425             }
00426 
00427             // No need to call FinaliseCellCreation() as halo node cardiac cells will
00428             // never be stimulated (their values are communicated from the process that
00429             // owns them.
00430         }
00431         catch (const Exception& e)
00432         {
00433             // Errors thrown creating cells will often be process-specific
00434             PetscTools::ReplicateException(true);
00435 
00436             // Delete cells
00437             // Should really do this for other processes too, but this is all we need
00438             // to get memory testing to pass, and leaking when we're about to die isn't
00439             // that bad!
00440             for (std::vector<AbstractCardiacCell*>::iterator cell_iterator = mHaloCellsDistributed.begin();
00441                  cell_iterator != mHaloCellsDistributed.end();
00442                  ++cell_iterator)
00443             {
00444                 delete (*cell_iterator);
00445             }
00446 
00447             throw e;
00448         }
00449         PetscTools::ReplicateException(false);
00450     }
00451 }
00452 
00453 
00454 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00455 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::SolveCellSystems(Vec existingSolution, double time, double nextTime, bool updateVoltage)
00456 {
00457     HeartEventHandler::BeginEvent(HeartEventHandler::SOLVE_ODES);
00458 
00459     DistributedVector dist_solution = mpDistributedVectorFactory->CreateDistributedVector(existingSolution);
00460     DistributedVector::Stripe voltage(dist_solution, 0);
00461 
00462     try
00463     {
00464         for (DistributedVector::Iterator index = dist_solution.Begin();
00465              index != dist_solution.End();
00466              ++index)
00467         {
00468             // overwrite the voltage with the input value
00469             mCellsDistributed[index.Local]->SetVoltage( voltage[index] );
00470 
00471             if(!updateVoltage)
00472             {
00473                 // solve
00474                 // Note: Voltage is not be updated. The voltage is updated in the PDE solve.
00475                 mCellsDistributed[index.Local]->ComputeExceptVoltage(time, nextTime);
00476             }
00477             else
00478             {
00479                 // solve, including updating the voltage (for the operator-splitting implementation of the monodomain solver)
00480                 mCellsDistributed[index.Local]->SolveAndUpdateState(time, nextTime);
00481                 voltage[index] = mCellsDistributed[index.Local]->GetVoltage();
00482             }
00483 
00484             // update the Iionic and stimulus caches
00485             UpdateCaches(index.Global, index.Local, nextTime);
00486         }
00487 
00488         if(updateVoltage)
00489         {
00490             dist_solution.Restore();
00491         }
00492     }
00493     catch (Exception &e)
00494     {
00495         PetscTools::ReplicateException(true);
00496         throw e;
00497     }
00498 
00499     PetscTools::ReplicateException(false);
00500     HeartEventHandler::EndEvent(HeartEventHandler::SOLVE_ODES);
00501 
00502     // Communicate new state variable values to halo nodes
00503     if (mExchangeHalos)
00504     {
00505         for ( unsigned rank_offset = 1; rank_offset < PetscTools::GetNumProcs(); rank_offset++ )
00506         {
00507             unsigned send_to      = (PetscTools::GetMyRank() + rank_offset) % (PetscTools::GetNumProcs());
00508             unsigned receive_from = (PetscTools::GetMyRank() + PetscTools::GetNumProcs()- rank_offset ) % (PetscTools::GetNumProcs());
00509 
00510             unsigned number_of_cells_to_send    = mNodesToSendPerProcess[send_to].size();
00511             unsigned number_of_cells_to_receive = mNodesToReceivePerProcess[receive_from].size();
00512 
00513             // Pack
00514             if ( number_of_cells_to_send > 0 )
00515             {
00516                 unsigned send_size = 0;
00517                 for (unsigned i=0; i<number_of_cells_to_send; i++)
00518                 {
00519                     unsigned global_cell_index = mNodesToSendPerProcess[send_to][i];
00520                     send_size += mCellsDistributed[global_cell_index - mpDistributedVectorFactory->GetLow()]->GetNumberOfStateVariables();
00521                 }
00522 
00523                 double send_data[send_size];
00524 
00525                 unsigned send_index = 0;
00526                 for (unsigned cell = 0; cell < number_of_cells_to_send; cell++)
00527                 {
00528                     unsigned global_cell_index = mNodesToSendPerProcess[send_to][cell];
00529                     AbstractCardiacCell* p_cell = mCellsDistributed[global_cell_index - mpDistributedVectorFactory->GetLow()];
00530                     std::vector<double>& cell_data = p_cell->rGetStateVariables();
00531                     const unsigned num_state_vars = p_cell->GetNumberOfStateVariables();
00532                     for (unsigned state_variable = 0; state_variable < num_state_vars; state_variable++)
00533                     {
00534                         send_data[send_index++] = cell_data[state_variable];
00535                     }
00536                 }
00537 
00538                 // Send
00539                 int ret;
00540                 ret = MPI_Send( send_data,
00541                                 send_size,
00542                                 MPI_DOUBLE,
00543                                 send_to,
00544                                 0,
00545                                 PETSC_COMM_WORLD );
00546                 assert ( ret == MPI_SUCCESS );
00547             }
00548 
00549             if ( number_of_cells_to_receive > 0 )
00550             {
00551                 // Receive
00552                 unsigned receive_size = 0;
00553                 for (unsigned i=0; i<number_of_cells_to_receive; i++)
00554                 {
00555                     unsigned halo_cell_index = mHaloGlobalToLocalIndexMap[mNodesToReceivePerProcess[receive_from][i]];
00556                     receive_size += mHaloCellsDistributed[halo_cell_index]->GetNumberOfStateVariables();
00557                 }
00558 
00559                 double receive_data[receive_size];
00560                 MPI_Status status;
00561 
00562                 int ret;
00563                 ret = MPI_Recv( receive_data,
00564                                 receive_size,
00565                                 MPI_DOUBLE,
00566                                 receive_from,
00567                                 0,
00568                                 PETSC_COMM_WORLD,
00569                                 &status );
00570                 assert ( ret == MPI_SUCCESS);
00571 
00572                 // Unpack
00573                 unsigned receive_index = 0;
00574                 for ( unsigned cell = 0; cell < number_of_cells_to_receive; cell++ )
00575                 {
00576                     AbstractCardiacCell* p_cell = mHaloCellsDistributed[mHaloGlobalToLocalIndexMap[mNodesToReceivePerProcess[receive_from][cell]]];
00577                     std::vector<double> cell_data;
00578                     cell_data.resize(p_cell->GetNumberOfStateVariables());
00579                     const unsigned number_of_state_variables = p_cell->GetNumberOfStateVariables();
00580                     for (unsigned state_variable = 0; state_variable < number_of_state_variables; state_variable++)
00581                     {
00582                         cell_data[state_variable] = receive_data[receive_index++];
00583                     }
00584                     p_cell->SetStateVariables(cell_data);
00585                 }
00586             }
00587         }
00588     }
00589 
00590     HeartEventHandler::BeginEvent(HeartEventHandler::COMMUNICATION);
00591     if ( mDoCacheReplication )
00592     {
00593         ReplicateCaches();
00594     }
00595     HeartEventHandler::EndEvent(HeartEventHandler::COMMUNICATION);
00596 }
00597 
00598 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00599 ReplicatableVector& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetIionicCacheReplicated()
00600 {
00601     return mIionicCacheReplicated;
00602 }
00603 
00604 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00605 ReplicatableVector& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetIntracellularStimulusCacheReplicated()
00606 {
00607     return mIntracellularStimulusCacheReplicated;
00608 }
00609 
00610 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00611 ReplicatableVector& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetPurkinjeIionicCacheReplicated()
00612 {
00613     EXCEPT_IF_NOT(mHasPurkinje);
00614     return mPurkinjeIionicCacheReplicated;
00615 }
00616 
00617 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00618 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::UpdateCaches(unsigned globalIndex, unsigned localIndex, double nextTime)
00619 {
00620     mIionicCacheReplicated[globalIndex] = mCellsDistributed[localIndex]->GetIIonic();
00621     mIntracellularStimulusCacheReplicated[globalIndex] = mCellsDistributed[localIndex]->GetIntracellularStimulus(nextTime);
00622     if (mHasPurkinje)
00623     {
00625         NEVER_REACHED;
00626 //        mPurkinjeIionicCacheReplicated[globalIndex] = mPurkinjeCellsDistributed[localIndex]->GetIIonic();
00627     }
00628 }
00629 
00630 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00631 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::ReplicateCaches()
00632 {
00633     mIionicCacheReplicated.Replicate(mpDistributedVectorFactory->GetLow(), mpDistributedVectorFactory->GetHigh());
00634     mIntracellularStimulusCacheReplicated.Replicate(mpDistributedVectorFactory->GetLow(), mpDistributedVectorFactory->GetHigh());
00635     if (mHasPurkinje)
00636     {
00638         NEVER_REACHED;
00639 //        mPurkinjeIionicCacheReplicated.Replicate(mpDistributedVectorFactory->GetLow(), mpDistributedVectorFactory->GetHigh());
00640     }
00641 }
00642 
00643 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00644 const std::vector<AbstractCardiacCell*>& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetCellsDistributed() const
00645 {
00646     return mCellsDistributed;
00647 }
00648 
00649 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00650 const std::vector<AbstractCardiacCell*>& AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::rGetPurkinjeCellsDistributed() const
00651 {
00652     EXCEPT_IF_NOT(mHasPurkinje);
00653     return mPurkinjeCellsDistributed;
00654 }
00655 
00656 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00657 const AbstractTetrahedralMesh<ELEMENT_DIM,SPACE_DIM>* AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::pGetMesh() const
00658 {
00659     return mpMesh;
00660 }
00661 
00662 template <unsigned ELEMENT_DIM,unsigned SPACE_DIM>
00663 void AbstractCardiacTissue<ELEMENT_DIM,SPACE_DIM>::SetConductivityModifier(AbstractConductivityModifier<ELEMENT_DIM,SPACE_DIM>* pModifier)
00664 {
00665     assert(pModifier!=NULL);
00666     assert(mpConductivityModifier==NULL); // shouldn't be called twice for example, or with two different modifiers (remove this assert
00667                                           // if for whatever reason want to be able to overwrite modifiers
00668     mpConductivityModifier = pModifier;
00669 }
00670 
00671 
00673 // Explicit instantiation
00675 
00676 template class AbstractCardiacTissue<1,1>;
00677 template class AbstractCardiacTissue<1,2>;
00678 template class AbstractCardiacTissue<1,3>;
00679 template class AbstractCardiacTissue<2,2>;
00680 template class AbstractCardiacTissue<3,3>;
Generated on Thu Dec 22 13:00:07 2011 for Chaste by  doxygen 1.6.3