Chaste Release::3.1
CellBasedPdeHandler.cpp
00001 /*
00002 
00003 Copyright (c) 2005-2012, University of Oxford.
00004 All rights reserved.
00005 
00006 University of Oxford means the Chancellor, Masters and Scholars of the
00007 University of Oxford, having an administrative office at Wellington
00008 Square, Oxford OX1 2JD, UK.
00009 
00010 This file is part of Chaste.
00011 
00012 Redistribution and use in source and binary forms, with or without
00013 modification, are permitted provided that the following conditions are met:
00014  * Redistributions of source code must retain the above copyright notice,
00015    this list of conditions and the following disclaimer.
00016  * Redistributions in binary form must reproduce the above copyright notice,
00017    this list of conditions and the following disclaimer in the documentation
00018    and/or other materials provided with the distribution.
00019  * Neither the name of the University of Oxford nor the names of its
00020    contributors may be used to endorse or promote products derived from this
00021    software without specific prior written permission.
00022 
00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00033 
00034 */
00035 
00036 #include "CellBasedPdeHandler.hpp"
00037 #include "MeshBasedCellPopulationWithGhostNodes.hpp"
00038 #include "NodeBasedCellPopulation.hpp"
00039 #include "PottsBasedCellPopulation.hpp"
00040 #include "BoundaryConditionsContainer.hpp"
00041 #include "SimpleLinearEllipticSolver.hpp"
00042 #include "CellBasedPdeSolver.hpp"
00043 #include "Exception.hpp"
00044 
00045 template<unsigned DIM>
00046 CellBasedPdeHandler<DIM>::CellBasedPdeHandler(AbstractCellPopulation<DIM>* pCellPopulation,
00047                                               bool deleteMemberPointersInDestructor)
00048     : mpCellPopulation(pCellPopulation),
00049       mWriteAverageRadialPdeSolution(false),
00050       mWriteDailyAverageRadialPdeSolution(false),
00051       mAverageRadialSolutionVariableName(""),
00052       mSetBcsOnCoarseBoundary(true),
00053       mNumRadialIntervals(UNSIGNED_UNSET),
00054       mpCoarsePdeMesh(NULL),
00055       mDeleteMemberPointersInDestructor(deleteMemberPointersInDestructor)
00056 {
00057     // We must be using a NodeBasedCellPopulation or MeshBasedCellPopulation, with at least one cell
00059     assert(mpCellPopulation->GetNumRealCells() != 0);
00060 }
00061 
00062 template<unsigned DIM>
00063 CellBasedPdeHandler<DIM>::~CellBasedPdeHandler()
00064 {
00065     /*
00066      * Avoid memory leaks. Note that we do not take responsibility for
00067      * deleting mpCellPopulation, as this object is usually owned by a
00068      * subclass of AbstractCellBasedSimulation, which deletes the cell
00069      * population upon destruction if restored from an archive.
00070      */
00071     if (mDeleteMemberPointersInDestructor)
00072     {
00073         for (unsigned i=0; i<mPdeAndBcCollection.size(); i++)
00074         {
00075             delete mPdeAndBcCollection[i];
00076         }
00077     }
00078     if (mpCoarsePdeMesh)
00079     {
00080         delete mpCoarsePdeMesh;
00081     }
00082 }
00083 
00084 template<unsigned DIM>
00085 const AbstractCellPopulation<DIM>* CellBasedPdeHandler<DIM>::GetCellPopulation() const
00086 {
00087     return mpCellPopulation;
00088 }
00089 
00090 template<unsigned DIM>
00091 TetrahedralMesh<DIM,DIM>* CellBasedPdeHandler<DIM>::GetCoarsePdeMesh()
00092 {
00093     return mpCoarsePdeMesh;
00094 }
00095 
00096 template<unsigned DIM>
00097 void CellBasedPdeHandler<DIM>::AddPdeAndBc(PdeAndBoundaryConditions<DIM>* pPdeAndBc)
00098 {
00099     if (!mPdeAndBcCollection.empty() && (pPdeAndBc->rGetDependentVariableName()==""))
00100     {
00101         EXCEPTION("When adding more than one PDE to CellBasedPdeHandler set the dependent variable name using SetDependentVariableName(name).");
00102     }
00103     for (unsigned i=0; i< mPdeAndBcCollection.size(); i++)
00104     {
00105         if (pPdeAndBc->rGetDependentVariableName() == mPdeAndBcCollection[i]->rGetDependentVariableName())
00106         {
00107             EXCEPTION("The name "+ pPdeAndBc->rGetDependentVariableName()+ " has already been used in the PDE collection");
00108         }
00109     }
00110     mPdeAndBcCollection.push_back(pPdeAndBc);
00111 }
00112 
00113 template<unsigned DIM>
00114 Vec CellBasedPdeHandler<DIM>::GetPdeSolution(const std::string& rName)
00115 {
00116     if (rName != "")
00117     {
00118         for (unsigned i=0; i<mPdeAndBcCollection.size(); i++)
00119         {
00120             if (mPdeAndBcCollection[i]->rGetDependentVariableName() == rName)
00121             {
00122                 return mPdeAndBcCollection[i]->GetSolution();
00123             }
00124         }
00125 
00126         EXCEPTION("The PDE collection does not contain a PDE named " + rName);
00127     }
00128     else
00129     {
00130         assert(mPdeAndBcCollection.size()==1);
00131         return mPdeAndBcCollection[0]->GetSolution();
00132     }
00133 }
00134 
00135 template<unsigned DIM>
00136 void CellBasedPdeHandler<DIM>::InitialiseCellPdeElementMap()
00137 {
00138     if (mpCoarsePdeMesh == NULL)
00139     {
00140         EXCEPTION("InitialiseCellPdeElementMap() should only be called if mpCoarsePdeMesh is set up.");
00141     }
00142 
00143     mCellPdeElementMap.clear();
00144 
00145     // Find the element of mpCoarsePdeMesh that contains each cell and populate mCellPdeElementMap
00146     for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
00147          cell_iter != mpCellPopulation->End();
00148          ++cell_iter)
00149     {
00150         const ChastePoint<DIM>& r_position_of_cell = mpCellPopulation->GetLocationOfCellCentre(*cell_iter);
00151         unsigned elem_index = mpCoarsePdeMesh->GetContainingElementIndex(r_position_of_cell);
00152         mCellPdeElementMap[*cell_iter] = elem_index;
00153     }
00154 }
00155 
00156 template<unsigned DIM>
00157 void CellBasedPdeHandler<DIM>::UpdateCellPdeElementMap()
00158 {
00159     // Find the element of mpCoarsePdeMesh that contains each cell and populate mCellPdeElementMap
00160     for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
00161          cell_iter != mpCellPopulation->End();
00162          ++cell_iter)
00163     {
00164         const ChastePoint<DIM>& r_position_of_cell = mpCellPopulation->GetLocationOfCellCentre(*cell_iter);
00165         unsigned elem_index = mpCoarsePdeMesh->GetContainingElementIndexWithInitialGuess(r_position_of_cell, mCellPdeElementMap[*cell_iter]);
00166         mCellPdeElementMap[*cell_iter] = elem_index;
00167     }
00168 }
00169 
00170 template<unsigned DIM>
00171 void CellBasedPdeHandler<DIM>::OpenResultsFiles(std::string outputDirectory)
00172 {
00173     // If using a NodeBasedCellPopulation or a PottsBasedCellPopulation, mpCoarsePdeMesh must be set up
00174     if (PdeSolveNeedsCoarseMesh() && mpCoarsePdeMesh==NULL)
00175     {
00176         EXCEPTION("Trying to solve a PDE on a cell population that doesn't have a mesh. Try calling UseCoarsePdeMesh().");
00177     }
00178 
00179     if (mpCoarsePdeMesh != NULL)
00180     {
00181         InitialiseCellPdeElementMap();
00182 
00183         // Write mesh to file
00184         TrianglesMeshWriter<DIM,DIM> mesh_writer(outputDirectory+"/coarse_mesh_output", "coarse_mesh",false);
00185         mesh_writer.WriteFilesUsingMesh(*mpCoarsePdeMesh);
00186     }
00187 
00188     if (PetscTools::AmMaster())
00189     {
00190         OutputFileHandler output_file_handler(outputDirectory+"/", false);
00191 
00192         if (mpCoarsePdeMesh != NULL)
00193         {
00194             mpVizPdeSolutionResultsFile = output_file_handler.OpenOutputFile("results.vizcoarsepdesolution");
00195         }
00196         else
00197         {
00198             mpVizPdeSolutionResultsFile = output_file_handler.OpenOutputFile("results.vizpdesolution");
00199         }
00200 
00201         if (mWriteAverageRadialPdeSolution)
00202         {
00203             mpAverageRadialPdeSolutionResultsFile = output_file_handler.OpenOutputFile("radial_dist.dat");
00204         }
00205     }
00206 
00207     double current_time = SimulationTime::Instance()->GetTime();
00208     WritePdeSolution(current_time);
00209 }
00210 
00211 template<unsigned DIM>
00212 void CellBasedPdeHandler<DIM>::CloseResultsFiles()
00213 {
00214     // Close results files
00215     if (PetscTools::AmMaster())
00216     {
00217         mpVizPdeSolutionResultsFile->close();
00218         if (mWriteAverageRadialPdeSolution)
00219         {
00220             WriteAverageRadialPdeSolution(SimulationTime::Instance()->GetTime());
00221             mpAverageRadialPdeSolutionResultsFile->close();
00222         }
00223     }
00224 }
00225 
00226 template<unsigned DIM>
00227 void CellBasedPdeHandler<DIM>::UseCoarsePdeMesh(double stepSize, ChasteCuboid<DIM> meshCuboid, bool centreOnCellPopulation)
00228 {
00229     // If solving PDEs on a coarse mesh, each PDE must have an averaged source term
00230     if (mPdeAndBcCollection.empty())
00231     {
00232         EXCEPTION("mPdeAndBcCollection should be populated prior to calling UseCoarsePdeMesh().");
00233     }
00234     for (unsigned pde_index=0; pde_index<mPdeAndBcCollection.size(); pde_index++)
00235     {
00236         if (mPdeAndBcCollection[pde_index]->HasAveragedSourcePde() == false)
00237         {
00238             EXCEPTION("UseCoarsePdeMesh() should only be called if averaged-source PDEs are specified.");
00239         }
00240     }
00241 
00242     // Create a regular coarse tetrahedral mesh
00243     mpCoarsePdeMesh = new TetrahedralMesh<DIM,DIM>;
00244     switch (DIM)
00245     {
00246         case 1:
00247             mpCoarsePdeMesh->ConstructRegularSlabMesh(stepSize, meshCuboid.GetWidth(0));
00248             break;
00249         case 2:
00250             mpCoarsePdeMesh->ConstructRegularSlabMesh(stepSize, meshCuboid.GetWidth(0), meshCuboid.GetWidth(1));
00251             break;
00252         case 3:
00253             mpCoarsePdeMesh->ConstructRegularSlabMesh(stepSize, meshCuboid.GetWidth(0), meshCuboid.GetWidth(1), meshCuboid.GetWidth(2));
00254             break;
00255         default:
00256             NEVER_REACHED;
00257     }
00258 
00259     if (centreOnCellPopulation)
00260     {
00261         // Find the centre of the coarse PDE mesh
00262         c_vector<double,DIM> centre_of_coarse_mesh = zero_vector<double>(DIM);
00263         for (unsigned i=0; i<mpCoarsePdeMesh->GetNumNodes(); i++)
00264         {
00265             centre_of_coarse_mesh += mpCoarsePdeMesh->GetNode(i)->rGetLocation();
00266         }
00267         centre_of_coarse_mesh /= mpCoarsePdeMesh->GetNumNodes();
00268 
00269         // Translate the centre of coarse PDE mesh to the centre of the cell population
00270         c_vector<double,DIM> centre_of_cell_population = mpCellPopulation->GetCentroidOfCellPopulation();
00271         mpCoarsePdeMesh->Translate(centre_of_cell_population - centre_of_coarse_mesh);
00272     }
00273     else
00274     {
00275         // Get centroid of meshCuboid
00276         ChastePoint<DIM> upper = meshCuboid.rGetUpperCorner();
00277         ChastePoint<DIM> lower = meshCuboid.rGetLowerCorner();
00278         c_vector<double,DIM> centre_of_cuboid = 0.5*(upper.rGetLocation() + lower.rGetLocation());
00279 
00280         // Find the centre of the coarse PDE mesh
00281         c_vector<double,DIM> centre_of_coarse_mesh = zero_vector<double>(DIM);
00282         for (unsigned i=0; i<mpCoarsePdeMesh->GetNumNodes(); i++)
00283         {
00284             centre_of_coarse_mesh += mpCoarsePdeMesh->GetNode(i)->rGetLocation();
00285         }
00286         centre_of_coarse_mesh /= mpCoarsePdeMesh->GetNumNodes();
00287 
00288         mpCoarsePdeMesh->Translate(centre_of_cuboid - centre_of_coarse_mesh);
00289     }
00290 }
00291 
00292 template<unsigned DIM>
00293 void CellBasedPdeHandler<DIM>::SolvePdeAndWriteResultsToFile(unsigned samplingTimestepMultiple)
00294 {
00295     // Record whether we are solving PDEs on a coarse mesh
00296     bool using_coarse_pde_mesh = (mpCoarsePdeMesh != NULL);
00297 
00298     // If solving PDEs on a coarse mesh, each PDE should have an averaged source term; otherwise none should
00299     assert(!mPdeAndBcCollection.empty());
00300     for (unsigned pde_index=0; pde_index<mPdeAndBcCollection.size(); pde_index++)
00301     {
00302         assert(mPdeAndBcCollection[pde_index]);
00303         assert(mPdeAndBcCollection[pde_index]->HasAveragedSourcePde() == using_coarse_pde_mesh);
00304     }
00305 
00306     // Make sure the cell population is in a nice state
00307     mpCellPopulation->Update();
00308 
00309     // Store a pointer to the (population-level or coarse) mesh
00310     TetrahedralMesh<DIM,DIM>* p_mesh;
00311     if (using_coarse_pde_mesh)
00312     {
00313         p_mesh = mpCoarsePdeMesh;
00314     }
00315     else
00316     {
00317         // If not using a coarse PDE mesh, we must be using a MeshBasedCellPopulation
00318         p_mesh = &(static_cast<MeshBasedCellPopulation<DIM>*>(mpCellPopulation)->rGetMesh());
00319     }
00320 
00321     // Loop over elements of mPdeAndBcCollection
00322     for (unsigned pde_index=0; pde_index<mPdeAndBcCollection.size(); pde_index++)
00323     {
00324         // Get pointer to this PdeAndBoundaryConditions object
00325         PdeAndBoundaryConditions<DIM>* p_pde_and_bc = mPdeAndBcCollection[pde_index];
00326 
00327         // Set up boundary conditions
00328         AbstractBoundaryCondition<DIM>* p_bc = p_pde_and_bc->GetBoundaryCondition();
00329         BoundaryConditionsContainer<DIM,DIM,1> bcc(false);
00330 
00331         if (p_pde_and_bc->IsNeumannBoundaryCondition()) // this BC is of Neumann type
00332         {
00333             if (using_coarse_pde_mesh)
00334             {
00336                 EXCEPTION("Neumann BCs not yet implemented when using a coarse PDE mesh");
00337             }
00338             else
00339             {
00340                 for (typename TetrahedralMesh<DIM,DIM>::BoundaryElementIterator elem_iter = p_mesh->GetBoundaryElementIteratorBegin();
00341                      elem_iter != p_mesh->GetBoundaryElementIteratorEnd();
00342                      ++elem_iter)
00343                 {
00344                     bcc.AddNeumannBoundaryCondition(*elem_iter, p_bc);
00345                 }
00346             }
00347         }
00348         else // assume that if the BC is of Neumann type, then it is Dirichlet
00349         {
00350             if (using_coarse_pde_mesh && !mSetBcsOnCoarseBoundary)
00351             {
00352                 // Get the set of coarse element indices that contain cells
00353                 std::set<unsigned> coarse_element_indices_in_map;
00354                 for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
00355                      cell_iter != mpCellPopulation->End();
00356                      ++cell_iter)
00357                 {
00358                     coarse_element_indices_in_map.insert(mCellPdeElementMap[*cell_iter]);
00359                 }
00360 
00361                 // Find the node indices associated with elements whose indices are NOT in the set coarse_element_indices_in_map
00362                 std::set<unsigned> coarse_mesh_boundary_node_indices;
00363                 for (unsigned i=0; i<p_mesh->GetNumElements(); i++)
00364                 {
00365                     if (coarse_element_indices_in_map.find(i) == coarse_element_indices_in_map.end())
00366                     {
00367                         Element<DIM,DIM>* p_element = p_mesh->GetElement(i);
00368                         for (unsigned j=0; j<DIM+1; j++)
00369                         {
00370                             unsigned node_index = p_element->GetNodeGlobalIndex(j);
00371                             coarse_mesh_boundary_node_indices.insert(node_index);
00372                         }
00373                     }
00374                 }
00375 
00376                 // Apply boundary condition to the nodes in the set coarse_mesh_boundary_node_indices
00377                 for (std::set<unsigned>::iterator iter = coarse_mesh_boundary_node_indices.begin();
00378                      iter != coarse_mesh_boundary_node_indices.end();
00379                      ++iter)
00380                 {
00381                     bcc.AddDirichletBoundaryCondition(p_mesh->GetNode(*iter), p_bc, 0, false);
00382                 }
00383             }
00384             else // apply BC at boundary nodes of (population-level or coarse) mesh
00385             {
00386                 for (typename TetrahedralMesh<DIM,DIM>::BoundaryNodeIterator node_iter = p_mesh->GetBoundaryNodeIteratorBegin();
00387                      node_iter != p_mesh->GetBoundaryNodeIteratorEnd();
00388                      ++node_iter)
00389                 {
00390                     bcc.AddDirichletBoundaryCondition(*node_iter, p_bc);
00391                 }
00392             }
00393         }
00394 
00395         // If the solution at the previous timestep exists...
00396         PetscInt previous_solution_size = 0;
00397         if (p_pde_and_bc->GetSolution())
00398         {
00399             VecGetSize(p_pde_and_bc->GetSolution(), &previous_solution_size);
00400         }
00401 
00402         // ...then record whether it is the correct size...
00403         bool is_previous_solution_size_correct = (previous_solution_size == (int)p_mesh->GetNumNodes());
00404 
00405         // ...and if it is, store it as an initial guess for the PDE solver
00406         Vec initial_guess;
00407         if (is_previous_solution_size_correct)
00408         {
00409             // This Vec is copied by the solver's Solve() method, so must be deleted here too
00410             VecDuplicate(p_pde_and_bc->GetSolution(), &initial_guess);
00411             VecCopy(p_pde_and_bc->GetSolution(), initial_guess);
00412             p_pde_and_bc->DestroySolution();
00413         }
00414         else
00415         {
00417             if (!using_coarse_pde_mesh && p_pde_and_bc->GetSolution())
00418             {
00419                 assert(previous_solution_size != 0);
00420                 p_pde_and_bc->DestroySolution();
00421             }
00422         }
00423 
00424         // Create a PDE solver and solve the PDE on the (population-level or coarse) mesh
00425         if (using_coarse_pde_mesh)
00426         {
00427             // When using a coarse PDE mesh, we must set up the source terms before solving the PDE.
00428             // pass in mCellPdeElementMap to speed up finding cells.
00429             this->UpdateCellPdeElementMap();
00430             p_pde_and_bc->SetUpSourceTermsForAveragedSourcePde(p_mesh, &mCellPdeElementMap);
00431 
00432             SimpleLinearEllipticSolver<DIM,DIM> solver(p_mesh, p_pde_and_bc->GetPde(), &bcc);
00433 
00434             // If we have an initial guess, use this when solving the system...
00435             if (is_previous_solution_size_correct)
00436             {
00437                 p_pde_and_bc->SetSolution(solver.Solve(initial_guess));
00438                 PetscTools::Destroy(initial_guess);
00439             }
00440             else // ...otherwise do not supply one
00441             {
00442                 p_pde_and_bc->SetSolution(solver.Solve());
00443             }
00444         }
00445         else
00446         {
00447             CellBasedPdeSolver<DIM> solver(p_mesh, p_pde_and_bc->GetPde(), &bcc);
00448 
00449             // If we have an initial guess, use this...
00450             if (is_previous_solution_size_correct)
00451             {
00452                 p_pde_and_bc->SetSolution(solver.Solve(initial_guess));
00453                 PetscTools::Destroy(initial_guess);
00454             }
00455             else // ...otherwise do not supply one
00456             {
00457                 p_pde_and_bc->SetSolution(solver.Solve());
00458             }
00459         }
00460 
00461         // Store the PDE solution in an accessible form
00462         ReplicatableVector solution_repl(p_pde_and_bc->GetSolution());
00463 
00464         // Having solved the PDE, now update CellData
00465         for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
00466              cell_iter != mpCellPopulation->End();
00467              ++cell_iter)
00468         {
00469             unsigned node_index = mpCellPopulation->GetLocationIndexUsingCell(*cell_iter);
00470             double solution_at_node = 0.0;
00471 
00472             if (using_coarse_pde_mesh)
00473             {
00474                 // When using a coarse PDE mesh, the cells are not nodes of the mesh, so we must interpolate
00475 
00476                 // Find the element in the coarse mesh that contains this cell. CellElementMap has been updated so use this.
00477                 unsigned elem_index = mCellPdeElementMap[*cell_iter];
00478                 Element<DIM,DIM>* p_element = mpCoarsePdeMesh->GetElement(elem_index);
00479 
00480                 const ChastePoint<DIM>& node_location = mpCellPopulation->GetLocationOfCellCentre(*cell_iter);
00481 
00482                 c_vector<double,DIM+1> weights = p_element->CalculateInterpolationWeights(node_location);
00483                 for (unsigned i=0; i<DIM+1; i++)
00484                 {
00485                     double nodal_value = solution_repl[p_element->GetNodeGlobalIndex(i)];
00486                     solution_at_node += nodal_value * weights(i);
00487                 }
00488             }
00489             else
00490             {
00491                 solution_at_node = solution_repl[node_index];
00492             }
00493             cell_iter->GetCellData()->SetItem(mPdeAndBcCollection[pde_index]->rGetDependentVariableName(), solution_at_node);
00494         }
00495     }
00496 
00497     // Write results to file if required
00498     SimulationTime* p_time = SimulationTime::Instance();
00499     double time_next_step = p_time->GetTime() + p_time->GetTimeStep();
00500     if ((p_time->GetTimeStepsElapsed()+1)%samplingTimestepMultiple == 0)
00501     {
00502         WritePdeSolution(time_next_step);
00503     }
00504 
00505 #define COVERAGE_IGNORE
00506 
00507     if (!using_coarse_pde_mesh)
00508     {
00509         if (mWriteDailyAverageRadialPdeSolution)
00510         {
00512             p_time = SimulationTime::Instance();
00513             time_next_step = p_time->GetTime() + p_time->GetTimeStep();
00514             unsigned num_timesteps_per_day = (unsigned) (DBL_EPSILON + 24/SimulationTime::Instance()->GetTimeStep());
00515             if ((p_time->GetTimeStepsElapsed()+1) % num_timesteps_per_day == 0)
00516             {
00517                 WriteAverageRadialPdeSolution(time_next_step);
00518             }
00519         }
00520     }
00521 #undef COVERAGE_IGNORE
00522 }
00523 
00524 template<unsigned DIM>
00525 unsigned CellBasedPdeHandler<DIM>::FindCoarseElementContainingCell(CellPtr pCell)
00526 {
00527     // Get containing element at last timestep from mCellPdeElementMap
00528     unsigned old_element_index = mCellPdeElementMap[pCell];
00529 
00530     // Create a std::set of guesses for the current containing element
00531     std::set<unsigned> test_elements;
00532     test_elements.insert(old_element_index);
00533 
00534     Element<DIM,DIM>* p_element = mpCoarsePdeMesh->GetElement(old_element_index);
00535     for (unsigned local_index=0; local_index<DIM+1; local_index++)
00536     {
00537         std::set<unsigned> element_indices = p_element->GetNode(local_index)->rGetContainingElementIndices();
00538         for (std::set<unsigned>::iterator iter = element_indices.begin();
00539              iter != element_indices.end();
00540              ++iter)
00541         {
00542             test_elements.insert(*iter);
00543         }
00544     }
00545 
00546     // Find new element, using the previous one as a guess
00547     const ChastePoint<DIM>& r_cell_position = mpCellPopulation->GetLocationOfCellCentre(pCell);
00548     unsigned new_element_index = mpCoarsePdeMesh->GetContainingElementIndex(r_cell_position, false, test_elements);
00549 
00550     // Update mCellPdeElementMap
00551     mCellPdeElementMap[pCell] = new_element_index;
00552 
00553     return new_element_index;
00554 }
00555 
00556 template<unsigned DIM>
00557 void CellBasedPdeHandler<DIM>::WritePdeSolution(double time)
00558 {
00559     if (PetscTools::AmMaster())
00560     {
00561         (*mpVizPdeSolutionResultsFile) << time << "\t";
00562 
00563         for (unsigned pde_index=0; pde_index<mPdeAndBcCollection.size(); pde_index++)
00564         {
00565             if (mpCoarsePdeMesh != NULL)
00566             {
00567                 PdeAndBoundaryConditions<DIM>* p_pde_and_bc = mPdeAndBcCollection[pde_index];
00568 
00569                 for (unsigned i=0; i<mpCoarsePdeMesh->GetNumNodes(); i++)
00570                 {
00571                     (*mpVizPdeSolutionResultsFile) << i << " ";
00572                     c_vector<double,DIM> location = mpCoarsePdeMesh->GetNode(i)->rGetLocation();
00573                     for (unsigned k=0; k<DIM; k++)
00574                     {
00575                         (*mpVizPdeSolutionResultsFile) << location[k] << " ";
00576                     }
00577 
00578                     if (p_pde_and_bc->GetSolution())
00579                     {
00580                         ReplicatableVector solution_repl(p_pde_and_bc->GetSolution());
00581                         (*mpVizPdeSolutionResultsFile) << solution_repl[i] << " ";
00582                     }
00583                     else
00584                     {
00586 
00587                         // Find the nearest cell to this coarse mesh node
00588                         unsigned nearest_node_index = 0;
00589                         double nearest_node_distance = DBL_MAX;
00590                         for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
00591                              cell_iter != mpCellPopulation->End();
00592                              ++cell_iter)
00593                         {
00594                             c_vector<double, DIM> node_location = mpCellPopulation->GetLocationOfCellCentre(*cell_iter);
00595                             if (norm_2(node_location - location) < nearest_node_distance)
00596                             {
00597                                 nearest_node_index = mpCellPopulation->GetLocationIndexUsingCell(*cell_iter);
00598                             }
00599                         }
00600 
00601                         CellPtr p_cell = mpCellPopulation->GetCellUsingLocationIndex(nearest_node_index);
00602                         double solution = p_cell->GetCellData()->GetItem(mPdeAndBcCollection[pde_index]->rGetDependentVariableName());
00603                         (*mpVizPdeSolutionResultsFile) << solution << " ";
00604                     }
00605                 }
00606             }
00607             else
00608             {
00609                 for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
00610                      cell_iter != mpCellPopulation->End();
00611                      ++cell_iter)
00612                 {
00613                     unsigned node_index = mpCellPopulation->GetLocationIndexUsingCell(*cell_iter);
00614                     (*mpVizPdeSolutionResultsFile) << node_index << " ";
00615                     const c_vector<double,DIM>& position = mpCellPopulation->GetLocationOfCellCentre(*cell_iter);
00616                     for (unsigned i=0; i<DIM; i++)
00617                     {
00618                         (*mpVizPdeSolutionResultsFile) << position[i] << " ";
00619                     }
00620                     double solution = cell_iter->GetCellData()->GetItem(mPdeAndBcCollection[pde_index]->rGetDependentVariableName());
00621                     (*mpVizPdeSolutionResultsFile) << solution << " ";
00622                 }
00623             }
00624         }
00625         (*mpVizPdeSolutionResultsFile) << "\n";
00626     }
00627 }
00628 
00629 template<unsigned DIM>
00630 void CellBasedPdeHandler<DIM>::SetWriteAverageRadialPdeSolution(const std::string& rName, unsigned numRadialIntervals, bool writeDailyResults)
00631 {
00632     mAverageRadialSolutionVariableName = rName;
00633     mWriteAverageRadialPdeSolution = true;
00634     mNumRadialIntervals = numRadialIntervals;
00635     mWriteDailyAverageRadialPdeSolution = writeDailyResults;
00636 }
00637 
00638 template<unsigned DIM>
00639 void CellBasedPdeHandler<DIM>::SetImposeBcsOnCoarseBoundary(bool setBcsOnCoarseBoundary)
00640 {
00641     mSetBcsOnCoarseBoundary = setBcsOnCoarseBoundary;
00642 }
00643 
00644 template<unsigned DIM>
00645 void CellBasedPdeHandler<DIM>::WriteAverageRadialPdeSolution(double time)
00646 {
00647     (*mpAverageRadialPdeSolutionResultsFile) << time << " ";
00648 
00649     // Calculate the centre of the cell population
00650     c_vector<double,DIM> centre = mpCellPopulation->GetCentroidOfCellPopulation();
00651 
00652     // Calculate the distance between each node and the centre of the cell population, as well as the maximum of these
00653     std::map<double, CellPtr> radius_cell_map;
00654     double max_distance_from_centre = 0.0;
00655     for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
00656          cell_iter != mpCellPopulation->End();
00657          ++cell_iter)
00658     {
00659         double distance = norm_2(mpCellPopulation->GetLocationOfCellCentre(*cell_iter) - centre);
00660         radius_cell_map[distance] = *cell_iter;
00661 
00662         if (distance > max_distance_from_centre)
00663         {
00664             max_distance_from_centre = distance;
00665         }
00666     }
00667 
00668     // Create vector of radius intervals
00669     std::vector<double> radius_intervals;
00670     for (unsigned i=0; i<mNumRadialIntervals; i++)
00671     {
00672         double upper_radius = max_distance_from_centre*((double) i+1)/((double) mNumRadialIntervals);
00673         radius_intervals.push_back(upper_radius);
00674     }
00675 
00676     // Calculate PDE solution in each radial interval
00677     double lower_radius = 0.0;
00678     for (unsigned i=0; i<mNumRadialIntervals; i++)
00679     {
00680         unsigned counter = 0;
00681         double average_solution = 0.0;
00682 
00683         for (std::map<double, CellPtr>::iterator iter = radius_cell_map.begin(); iter != radius_cell_map.end(); ++iter)
00684         {
00685             if (iter->first > lower_radius && iter->first <= radius_intervals[i])
00686             {
00687                 average_solution += (iter->second)->GetCellData()->GetItem(mAverageRadialSolutionVariableName);
00688                 counter++;
00689             }
00690         }
00691         if (counter != 0)
00692         {
00693             average_solution /= (double) counter;
00694         }
00695 
00696         // Write results to file
00697         (*mpAverageRadialPdeSolutionResultsFile) << radius_intervals[i] << " " << average_solution << " ";
00698         lower_radius = radius_intervals[i];
00699     }
00700     (*mpAverageRadialPdeSolutionResultsFile) << "\n";
00701 }
00702 
00703 template<unsigned DIM>
00704 bool CellBasedPdeHandler<DIM>::GetWriteAverageRadialPdeSolution()
00705 {
00706     return mWriteAverageRadialPdeSolution;
00707 }
00708 
00709 template<unsigned DIM>
00710 bool CellBasedPdeHandler<DIM>::GetWriteDailyAverageRadialPdeSolution()
00711 {
00712     return mWriteDailyAverageRadialPdeSolution;
00713 }
00714 
00715 template<unsigned DIM>
00716 bool CellBasedPdeHandler<DIM>::GetImposeBcsOnCoarseBoundary()
00717 {
00718     return mSetBcsOnCoarseBoundary;
00719 }
00720 
00721 template<unsigned DIM>
00722 unsigned CellBasedPdeHandler<DIM>::GetNumRadialIntervals()
00723 {
00724     return mNumRadialIntervals;
00725 }
00726 
00727 template<unsigned DIM>
00728 void CellBasedPdeHandler<DIM>::OutputParameters(out_stream& rParamsFile)
00729 {
00730     std::string type = GetIdentifier();
00731 
00732     *rParamsFile << "\t\t<" << type << ">\n";
00733     *rParamsFile << "\t\t<WriteAverageRadialPdeSolution>" << mWriteAverageRadialPdeSolution << "</WriteAverageRadialPdeSolution>\n";
00734     *rParamsFile << "\t\t<WriteDailyAverageRadialPdeSolution>" << mWriteDailyAverageRadialPdeSolution << "</WriteDailyAverageRadialPdeSolution>\n";
00735     *rParamsFile << "\t\t<SetBcsOnCoarseBoundary>" << mSetBcsOnCoarseBoundary << "</SetBcsOnCoarseBoundary>\n";
00736     *rParamsFile << "\t\t<NumRadialIntervals>" << mNumRadialIntervals << "</NumRadialIntervals>\n";
00737     *rParamsFile << "\t\t</" << type << ">\n";
00738 }
00739 
00740 template<unsigned DIM>
00741 bool CellBasedPdeHandler<DIM>::PdeSolveNeedsCoarseMesh()
00742 {
00743     return ((dynamic_cast<NodeBasedCellPopulation<DIM>*>(mpCellPopulation) != NULL) || (dynamic_cast<PottsBasedCellPopulation<DIM>*>(mpCellPopulation) != NULL));
00744 }
00745 
00746 // Serialization for Boost >= 1.36
00747 #include "SerializationExportWrapperForCpp.hpp"
00748 EXPORT_TEMPLATE_CLASS_SAME_DIMS(CellBasedPdeHandler)
00749 
00750 
00751 // Explicit instantiation
00753 
00754 template class CellBasedPdeHandler<1>;
00755 template class CellBasedPdeHandler<2>;
00756 template class CellBasedPdeHandler<3>;