Chaste  Release::3.4
CellBasedPdeHandler.cpp
1 /*
2 
3 Copyright (c) 2005-2016, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34 */
35 
36 #include "CellBasedPdeHandler.hpp"
37 #include "MeshBasedCellPopulationWithGhostNodes.hpp"
38 #include "NodeBasedCellPopulation.hpp"
39 #include "PottsBasedCellPopulation.hpp"
40 #include "VertexBasedCellPopulation.hpp"
41 #include "CaBasedCellPopulation.hpp"
42 #include "SimpleLinearEllipticSolver.hpp"
43 #include "CellBasedPdeSolver.hpp"
44 #include "Exception.hpp"
45 #include "VtkMeshWriter.hpp"
46 
47 template<unsigned DIM>
49  bool deleteMemberPointersInDestructor)
50  : mpCellPopulation(pCellPopulation),
51  mWriteAverageRadialPdeSolution(false),
52  mWriteDailyAverageRadialPdeSolution(false),
53  mAverageRadialSolutionVariableName(""),
54  mSetBcsOnCoarseBoundary(true),
55  mNumRadialIntervals(UNSIGNED_UNSET),
56  mpCoarsePdeMesh(NULL),
57  mDeleteMemberPointersInDestructor(deleteMemberPointersInDestructor)
58 {
59  // We must be using a CellPopulation with at least one cell
61  assert(mpCellPopulation->GetNumRealCells() != 0);
62 }
63 
64 template<unsigned DIM>
66 {
67  /*
68  * Avoid memory leaks. Note that we do not take responsibility for
69  * deleting mpCellPopulation, as this object is usually owned by a
70  * subclass of AbstractCellBasedSimulation, which deletes the cell
71  * population upon destruction if restored from an archive.
72  */
73  if (mDeleteMemberPointersInDestructor)
74  {
75  for (unsigned i=0; i<mPdeAndBcCollection.size(); i++)
76  {
77  delete mPdeAndBcCollection[i];
78  }
79  }
80  if (mpCoarsePdeMesh)
81  {
82  delete mpCoarsePdeMesh;
83  }
84 }
85 
86 template<unsigned DIM>
88 {
89  return mpCellPopulation;
90 }
91 
92 template<unsigned DIM>
94 {
95  return mpCoarsePdeMesh;
96 }
97 
98 template<unsigned DIM>
100 {
101  if (!mPdeAndBcCollection.empty() && (pPdeAndBc->rGetDependentVariableName()==""))
102  {
103  EXCEPTION("When adding more than one PDE to CellBasedPdeHandler set the dependent variable name using SetDependentVariableName(name).");
104  }
105  for (unsigned i=0; i< mPdeAndBcCollection.size(); i++)
106  {
107  if (pPdeAndBc->rGetDependentVariableName() == mPdeAndBcCollection[i]->rGetDependentVariableName())
108  {
109  EXCEPTION("The name "+ pPdeAndBc->rGetDependentVariableName()+ " has already been used in the PDE collection");
110  }
111  }
112  mPdeAndBcCollection.push_back(pPdeAndBc);
113 }
114 
115 template<unsigned DIM>
117 {
118  if (rName != "")
119  {
120  for (unsigned i=0; i<mPdeAndBcCollection.size(); i++)
121  {
122  if (mPdeAndBcCollection[i]->rGetDependentVariableName() == rName)
123  {
124  return mPdeAndBcCollection[i]->GetSolution();
125  }
126  }
127 
128  EXCEPTION("The PDE collection does not contain a PDE named " + rName);
129  }
130  else
131  {
132  assert(mPdeAndBcCollection.size()==1);
133  return mPdeAndBcCollection[0]->GetSolution();
134  }
135 }
136 
137 template<unsigned DIM>
139 {
140  if (mpCoarsePdeMesh == NULL)
141  {
142  EXCEPTION("InitialiseCellPdeElementMap() should only be called if mpCoarsePdeMesh is set up.");
143  }
144 
145  mCellPdeElementMap.clear();
146 
147  // Find the element of mpCoarsePdeMesh that contains each cell and populate mCellPdeElementMap
148  for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
149  cell_iter != mpCellPopulation->End();
150  ++cell_iter)
151  {
152  const ChastePoint<DIM>& r_position_of_cell = mpCellPopulation->GetLocationOfCellCentre(*cell_iter);
153  unsigned elem_index = mpCoarsePdeMesh->GetContainingElementIndex(r_position_of_cell);
154  mCellPdeElementMap[*cell_iter] = elem_index;
155  }
156 }
157 
158 template<unsigned DIM>
160 {
161  // Find the element of mpCoarsePdeMesh that contains each cell and populate mCellPdeElementMap
162  for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
163  cell_iter != mpCellPopulation->End();
164  ++cell_iter)
165  {
166  const ChastePoint<DIM>& r_position_of_cell = mpCellPopulation->GetLocationOfCellCentre(*cell_iter);
167  unsigned elem_index = mpCoarsePdeMesh->GetContainingElementIndexWithInitialGuess(r_position_of_cell, mCellPdeElementMap[*cell_iter]);
168  mCellPdeElementMap[*cell_iter] = elem_index;
169  }
170 }
171 
172 template<unsigned DIM>
173 void CellBasedPdeHandler<DIM>::OpenResultsFiles(std::string outputDirectory)
174 {
175  // If appropriate, make a coarse mesh which exactly overlays the lattice sites of a PottsMesh (used for all OnLattice simulations)
176  if ((dynamic_cast<CaBasedCellPopulation<DIM>*>(mpCellPopulation) != NULL) && mpCoarsePdeMesh==NULL)
177  {
178  assert(DIM ==2);
179  ChasteCuboid<DIM> cuboid = mpCellPopulation->rGetMesh().CalculateBoundingBox();
180 
181  // Currently only works with square meshes
182  assert(cuboid.GetWidth(0) == cuboid.GetWidth(1));
183 
184  UseCoarsePdeMesh(1, cuboid, false);
185  }
186 
187  // If using a NodeBasedCellPopulation a VertexBasedCellPopulation, a CaBasedCellPopulation or a PottsBasedCellPopulation, mpCoarsePdeMesh must be set up
188  if (PdeSolveNeedsCoarseMesh() && mpCoarsePdeMesh==NULL)
189  {
190  EXCEPTION("Trying to solve a PDE on a cell population that doesn't have a mesh. Try calling UseCoarsePdeMesh().");
191  }
192 
193  if (mpCoarsePdeMesh != NULL)
194  {
195  InitialiseCellPdeElementMap();
196 
197  // Write mesh to file
198  TrianglesMeshWriter<DIM,DIM> mesh_writer(outputDirectory+"/coarse_mesh_output", "coarse_mesh",false);
199  mesh_writer.WriteFilesUsingMesh(*mpCoarsePdeMesh);
200  }
201 
202  if (PetscTools::AmMaster())
203  {
204  OutputFileHandler output_file_handler(outputDirectory+"/", false);
205 
206  if (mpCoarsePdeMesh != NULL)
207  {
208  mpVizPdeSolutionResultsFile = output_file_handler.OpenOutputFile("results.vizcoarsepdesolution");
209  }
210  else
211  {
212  mpVizPdeSolutionResultsFile = output_file_handler.OpenOutputFile("results.vizpdesolution");
213  }
214 
215  if (mWriteAverageRadialPdeSolution)
216  {
217  mpAverageRadialPdeSolutionResultsFile = output_file_handler.OpenOutputFile("radial_dist.dat");
218  }
219  }
220 
221  mDirPath = outputDirectory; // caching the path to the output directory for VTK
222  //double current_time = SimulationTime::Instance()->GetTime();
223  //WritePdeSolution(current_time);
224 }
225 
226 template<unsigned DIM>
228 {
229  // Close results files
230  if (PetscTools::AmMaster())
231  {
232  mpVizPdeSolutionResultsFile->close();
233  if (mWriteAverageRadialPdeSolution)
234  {
235  WriteAverageRadialPdeSolution(SimulationTime::Instance()->GetTime());
236  mpAverageRadialPdeSolutionResultsFile->close();
237  }
238  }
239 }
240 
241 template<unsigned DIM>
242 void CellBasedPdeHandler<DIM>::UseCoarsePdeMesh(double stepSize, ChasteCuboid<DIM> meshCuboid, bool centreOnCellPopulation)
243 {
244  if (mPdeAndBcCollection.empty())
245  {
246  EXCEPTION("mPdeAndBcCollection should be populated prior to calling UseCoarsePdeMesh().");
247  }
248  // If solving PDEs on a coarse mesh, each PDE must have an averaged source term
249  for (unsigned pde_index=0; pde_index<mPdeAndBcCollection.size(); pde_index++)
250  {
251  if (mPdeAndBcCollection[pde_index]->HasAveragedSourcePde() == false && !dynamic_cast<CaBasedCellPopulation<DIM>*>(mpCellPopulation))
252  {
253  EXCEPTION("UseCoarsePdeMesh() should only be called if averaged-source PDEs are specified.");
254  }
255  }
256 
257  // Create a regular coarse tetrahedral mesh
258  mpCoarsePdeMesh = new TetrahedralMesh<DIM,DIM>;
259  switch (DIM)
260  {
261  case 1:
262  mpCoarsePdeMesh->ConstructRegularSlabMesh(stepSize, meshCuboid.GetWidth(0));
263  break;
264  case 2:
265  mpCoarsePdeMesh->ConstructRegularSlabMesh(stepSize, meshCuboid.GetWidth(0), meshCuboid.GetWidth(1));
266  break;
267  case 3:
268  mpCoarsePdeMesh->ConstructRegularSlabMesh(stepSize, meshCuboid.GetWidth(0), meshCuboid.GetWidth(1), meshCuboid.GetWidth(2));
269  break;
270  default:
272  }
273 
274  if (centreOnCellPopulation)
275  {
276  // Find the centre of the coarse PDE mesh
277  c_vector<double,DIM> centre_of_coarse_mesh = zero_vector<double>(DIM);
278  for (unsigned i=0; i<mpCoarsePdeMesh->GetNumNodes(); i++)
279  {
280  centre_of_coarse_mesh += mpCoarsePdeMesh->GetNode(i)->rGetLocation();
281  }
282  centre_of_coarse_mesh /= mpCoarsePdeMesh->GetNumNodes();
283 
284  // Translate the centre of coarse PDE mesh to the centre of the cell population
285  c_vector<double,DIM> centre_of_cell_population = mpCellPopulation->GetCentroidOfCellPopulation();
286  mpCoarsePdeMesh->Translate(centre_of_cell_population - centre_of_coarse_mesh);
287  }
288  else
289  {
290  // Get centroid of meshCuboid
291  ChastePoint<DIM> upper = meshCuboid.rGetUpperCorner();
292  ChastePoint<DIM> lower = meshCuboid.rGetLowerCorner();
293  c_vector<double,DIM> centre_of_cuboid = 0.5*(upper.rGetLocation() + lower.rGetLocation());
294 
295  // Find the centre of the coarse PDE mesh
296  c_vector<double,DIM> centre_of_coarse_mesh = zero_vector<double>(DIM);
297  for (unsigned i=0; i<mpCoarsePdeMesh->GetNumNodes(); i++)
298  {
299  centre_of_coarse_mesh += mpCoarsePdeMesh->GetNode(i)->rGetLocation();
300  }
301  centre_of_coarse_mesh /= mpCoarsePdeMesh->GetNumNodes();
302 
303  mpCoarsePdeMesh->Translate(centre_of_cuboid - centre_of_coarse_mesh);
304  }
305 }
306 
307 template<unsigned DIM>
308 void CellBasedPdeHandler<DIM>::SolvePdeAndWriteResultsToFile(unsigned samplingTimestepMultiple)
309 {
310  // Record whether we are solving PDEs on a coarse mesh
311  bool using_coarse_pde_mesh = (mpCoarsePdeMesh != NULL);
312 
313  // If solving PDEs on a coarse mesh, each PDE should have an averaged source term; otherwise none should
314  assert(!mPdeAndBcCollection.empty());
315  for (unsigned pde_index=0; pde_index<mPdeAndBcCollection.size(); pde_index++)
316  {
317  assert(mPdeAndBcCollection[pde_index]);
318  assert(mPdeAndBcCollection[pde_index]->HasAveragedSourcePde() == using_coarse_pde_mesh || dynamic_cast<CaBasedCellPopulation<DIM>*>(mpCellPopulation));
319  }
320 
321  // Make sure the cell population is in a nice state
322  mpCellPopulation->Update();
323 
324  // Store a pointer to the (population-level or coarse) mesh
325  TetrahedralMesh<DIM,DIM>* p_mesh;
326  if (using_coarse_pde_mesh)
327  {
328  p_mesh = mpCoarsePdeMesh;
329  }
330  else
331  {
332  // If not using a coarse PDE mesh, we must be using a MeshBasedCellPopulation
333  p_mesh = &(static_cast<MeshBasedCellPopulation<DIM>*>(mpCellPopulation)->rGetMesh());
334  }
335 
336  // Loop over elements of mPdeAndBcCollection
337  for (unsigned pde_index=0; pde_index<mPdeAndBcCollection.size(); pde_index++)
338  {
339  // Get pointer to this PdeAndBoundaryConditions object
340  PdeAndBoundaryConditions<DIM>* p_pde_and_bc = mPdeAndBcCollection[pde_index];
341 
342  // Set up boundary conditions
343  std::auto_ptr<BoundaryConditionsContainer<DIM,DIM,1> > p_bcc = ConstructBoundaryConditionsContainer(p_pde_and_bc, p_mesh);
344 
345  // If the solution at the previous timestep exists...
346  PetscInt previous_solution_size = 0;
347  if (p_pde_and_bc->GetSolution())
348  {
349  VecGetSize(p_pde_and_bc->GetSolution(), &previous_solution_size);
350  }
351 
352  // ...then record whether it is the correct size...
353  bool is_previous_solution_size_correct = (previous_solution_size == (int)p_mesh->GetNumNodes());
354 
355  // ...and if it is, store it as an initial guess for the PDE solver
356  Vec initial_guess;
357  if (is_previous_solution_size_correct)
358  {
359  // This Vec is copied by the solver's Solve() method, so must be deleted here too
360  VecDuplicate(p_pde_and_bc->GetSolution(), &initial_guess);
361  VecCopy(p_pde_and_bc->GetSolution(), initial_guess);
362  p_pde_and_bc->DestroySolution();
363  }
364  else
365  {
367  if (!using_coarse_pde_mesh && p_pde_and_bc->GetSolution())
368  {
369  assert(previous_solution_size != 0);
370  p_pde_and_bc->DestroySolution();
371  }
372  }
373 
374  if (using_coarse_pde_mesh)
375  {
376  // We update mCellPdeElementMap before setting up source terms to speed up
377  // finding cells in the case of AveragedSourcePdes. This is also used with
378  // non-AveragedSourcePdes when using a CaBasedCellPopulation.
379  this->UpdateCellPdeElementMap();
380  }
381 
382  // Create a PDE solver and solve the PDE on the (population-level or coarse) mesh
383  if (p_pde_and_bc->HasAveragedSourcePde())
384  {
385  // When using a coarse PDE mesh, we must set up the source terms before solving the PDE.
386  // Pass in already updated CellPdeElementMap to speed up finding cells.
387  p_pde_and_bc->SetUpSourceTermsForAveragedSourcePde(p_mesh, &mCellPdeElementMap);
388 
389  SimpleLinearEllipticSolver<DIM,DIM> solver(p_mesh, p_pde_and_bc->GetPde(), p_bcc.get());
390 
391  // If we have an initial guess, use this when solving the system...
392  if (is_previous_solution_size_correct)
393  {
394  p_pde_and_bc->SetSolution(solver.Solve(initial_guess));
395  PetscTools::Destroy(initial_guess);
396  }
397  else // ...otherwise do not supply one
398  {
399  p_pde_and_bc->SetSolution(solver.Solve());
400  }
401  }
402  else
403  {
404  CellBasedPdeSolver<DIM> solver(p_mesh, p_pde_and_bc->GetPde(), p_bcc.get());
405 
406  // If we have an initial guess, use this...
407  if (is_previous_solution_size_correct)
408  {
409  p_pde_and_bc->SetSolution(solver.Solve(initial_guess));
410  PetscTools::Destroy(initial_guess);
411  }
412  else // ...otherwise do not supply one
413  {
414  p_pde_and_bc->SetSolution(solver.Solve());
415  }
416  }
417 
418  // Store the PDE solution in an accessible form
419  ReplicatableVector solution_repl(p_pde_and_bc->GetSolution());
420 
421  // Having solved the PDE, now update CellData
422  for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
423  cell_iter != mpCellPopulation->End();
424  ++cell_iter)
425  {
426  unsigned node_index = mpCellPopulation->GetLocationIndexUsingCell(*cell_iter);
427  double solution_at_node = 0.0;
428 
429  if (using_coarse_pde_mesh)
430  {
431  // When using a coarse PDE mesh, the cells are not nodes of the mesh, so we must interpolate
432 
433  // Find the element in the coarse mesh that contains this cell. CellElementMap has been updated so use this.
434  unsigned elem_index = mCellPdeElementMap[*cell_iter];
435  Element<DIM,DIM>* p_element = mpCoarsePdeMesh->GetElement(elem_index);
436 
437  const ChastePoint<DIM>& node_location = mpCellPopulation->GetLocationOfCellCentre(*cell_iter);
438 
439  c_vector<double,DIM+1> weights = p_element->CalculateInterpolationWeights(node_location);
440  for (unsigned i=0; i<DIM+1; i++)
441  {
442  double nodal_value = solution_repl[p_element->GetNodeGlobalIndex(i)];
443  solution_at_node += nodal_value * weights(i);
444  }
445  }
446  else
447  {
448  solution_at_node = solution_repl[node_index];
449  }
450  cell_iter->GetCellData()->SetItem(mPdeAndBcCollection[pde_index]->rGetDependentVariableName(), solution_at_node);
451  }
452  }
453 
454  // Write results to file if required
456  if ((p_time->GetTimeStepsElapsed())%samplingTimestepMultiple == 0)
457  {
458  WritePdeSolution(p_time->GetTime());
459  }
460 #define COVERAGE_IGNORE
461  if (!using_coarse_pde_mesh)
463  {
464  if (mWriteDailyAverageRadialPdeSolution)
465  {
467  p_time = SimulationTime::Instance();
468  unsigned num_timesteps_per_day = (unsigned) (DBL_EPSILON + 24/SimulationTime::Instance()->GetTimeStep());
469  if ((p_time->GetTimeStepsElapsed()) % num_timesteps_per_day == 0)
470  {
471  WriteAverageRadialPdeSolution(p_time->GetTime());
472  }
473  }
474  }
475 #undef COVERAGE_IGNORE
476 }
477 
478 template<unsigned DIM>
479 std::auto_ptr<BoundaryConditionsContainer<DIM,DIM,1> > CellBasedPdeHandler<DIM>::ConstructBoundaryConditionsContainer(
482 {
483  std::auto_ptr<BoundaryConditionsContainer<DIM,DIM,1> > p_bcc(new BoundaryConditionsContainer<DIM,DIM,1>(false));
484 
486 
487  if (pPdeAndBc->IsNeumannBoundaryCondition()) // this BC is of Neumann type
488  {
489  // Note p_mesh is the coarse mesh or the natural mesh as appropriate
491  elem_iter != pMesh->GetBoundaryElementIteratorEnd();
492  ++elem_iter)
493  {
494  p_bcc->AddNeumannBoundaryCondition(*elem_iter, p_bc);
495  }
496  }
497  else // assume that if the BC is not of Neumann type, then it is Dirichlet
498  {
499  bool using_coarse_pde_mesh = (mpCoarsePdeMesh != NULL);
500 
501  if (using_coarse_pde_mesh && !mSetBcsOnCoarseBoundary)
502  {
503  // Get the set of coarse element indices that contain cells
504  std::set<unsigned> coarse_element_indices_in_map;
505  for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
506  cell_iter != mpCellPopulation->End();
507  ++cell_iter)
508  {
509  coarse_element_indices_in_map.insert(mCellPdeElementMap[*cell_iter]);
510  }
511 
512  // Find the node indices associated with elements whose indices are NOT in the set coarse_element_indices_in_map
513  std::set<unsigned> coarse_mesh_boundary_node_indices;
514  for (unsigned i=0; i<pMesh->GetNumElements(); i++)
515  {
516  if (coarse_element_indices_in_map.find(i) == coarse_element_indices_in_map.end())
517  {
518  Element<DIM,DIM>* p_element = pMesh->GetElement(i);
519  for (unsigned j=0; j<DIM+1; j++)
520  {
521  unsigned node_index = p_element->GetNodeGlobalIndex(j);
522  coarse_mesh_boundary_node_indices.insert(node_index);
523  }
524  }
525  }
526 
527  // Apply boundary condition to the nodes in the set coarse_mesh_boundary_node_indices
528  for (std::set<unsigned>::iterator iter = coarse_mesh_boundary_node_indices.begin();
529  iter != coarse_mesh_boundary_node_indices.end();
530  ++iter)
531  {
532  p_bcc->AddDirichletBoundaryCondition(pMesh->GetNode(*iter), p_bc, 0, false);
533  }
534  }
535  else // apply BC at boundary nodes of (population-level or coarse) mesh
536  {
538  node_iter != pMesh->GetBoundaryNodeIteratorEnd();
539  ++node_iter)
540  {
541  p_bcc->AddDirichletBoundaryCondition(*node_iter, p_bc);
542  }
543  }
544  }
545 
546  return p_bcc;
547 }
548 
549 template<unsigned DIM>
551 {
552  // Get containing element at last timestep from mCellPdeElementMap
553  unsigned old_element_index = mCellPdeElementMap[pCell];
554 
555  // Create a std::set of guesses for the current containing element
556  std::set<unsigned> test_elements;
557  test_elements.insert(old_element_index);
558 
559  Element<DIM,DIM>* p_element = mpCoarsePdeMesh->GetElement(old_element_index);
560  for (unsigned local_index=0; local_index<DIM+1; local_index++)
561  {
562  std::set<unsigned> element_indices = p_element->GetNode(local_index)->rGetContainingElementIndices();
563  for (std::set<unsigned>::iterator iter = element_indices.begin();
564  iter != element_indices.end();
565  ++iter)
566  {
567  test_elements.insert(*iter);
568  }
569  }
570 
571  // Find new element, using the previous one as a guess
572  const ChastePoint<DIM>& r_cell_position = mpCellPopulation->GetLocationOfCellCentre(pCell);
573  unsigned new_element_index = mpCoarsePdeMesh->GetContainingElementIndex(r_cell_position, false, test_elements);
574 
575  // Update mCellPdeElementMap
576  mCellPdeElementMap[pCell] = new_element_index;
577 
578  return new_element_index;
579 }
580 
581 template<unsigned DIM>
583 {
584  if (PetscTools::AmMaster())
585  {
586  (*mpVizPdeSolutionResultsFile) << time << "\t";
587 
588 #ifdef CHASTE_VTK
589  // Note that this mesh writer is only constructed and used if mpCoarsePdeMesh exists
590  VtkMeshWriter<DIM,DIM>* p_vtk_mesh_writer = NULL;
591  if (DIM>1 && mpCoarsePdeMesh != NULL )
592  {
593  std::ostringstream time_string;
594  time_string << SimulationTime::Instance()->GetTimeStepsElapsed()+1;
595  std::string results_file = "pde_results_"+time_string.str();
596  // Note that this mesh writer is always constructed, but is only used if mpCoarsePdeMesh exists
597  p_vtk_mesh_writer = new VtkMeshWriter<DIM,DIM>(mDirPath, results_file, false);
598  }
599 #endif //CHASTE_VTK
600  for (unsigned pde_index=0; pde_index<mPdeAndBcCollection.size(); pde_index++)
601  {
602  if (mpCoarsePdeMesh != NULL)
603  {
604  PdeAndBoundaryConditions<DIM>* p_pde_and_bc = mPdeAndBcCollection[pde_index];
605  assert(p_pde_and_bc->rGetDependentVariableName() != "");
606 
607 #ifdef CHASTE_VTK
608  if (p_pde_and_bc->GetSolution())
609  {
610  if (DIM>1)
611  {
612  ReplicatableVector solution_repl(p_pde_and_bc->GetSolution());
613  std::vector<double> pde_solution;
614  for (unsigned i=0; i<mpCoarsePdeMesh->GetNumNodes(); i++)
615  {
616  pde_solution.push_back(solution_repl[i]);
617  }
618 
619  p_vtk_mesh_writer->AddPointData(p_pde_and_bc->rGetDependentVariableName(),pde_solution);
620  }
621  }
622 
623 #endif //CHASTE_VTK
624 
625  for (unsigned i=0; i<mpCoarsePdeMesh->GetNumNodes(); i++)
626  {
627  (*mpVizPdeSolutionResultsFile) << i << " ";
628  c_vector<double,DIM> location = mpCoarsePdeMesh->GetNode(i)->rGetLocation();
629  for (unsigned k=0; k<DIM; k++)
630  {
631  (*mpVizPdeSolutionResultsFile) << location[k] << " ";
632  }
633 
634  if (p_pde_and_bc->GetSolution())
635  {
636  ReplicatableVector solution_repl(p_pde_and_bc->GetSolution());
637  (*mpVizPdeSolutionResultsFile) << solution_repl[i] << " ";
638  }
639  else
640  {
641  // should only come into this method AFTER solving the PDE
643  }
644  }
645  }
646  else // Not coarse mesh
647  {
648  for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
649  cell_iter != mpCellPopulation->End();
650  ++cell_iter)
651  {
652  unsigned node_index = mpCellPopulation->GetLocationIndexUsingCell(*cell_iter);
653  (*mpVizPdeSolutionResultsFile) << node_index << " ";
654  const c_vector<double,DIM>& position = mpCellPopulation->GetLocationOfCellCentre(*cell_iter);
655  for (unsigned i=0; i<DIM; i++)
656  {
657  (*mpVizPdeSolutionResultsFile) << position[i] << " ";
658  }
659  double solution = cell_iter->GetCellData()->GetItem(mPdeAndBcCollection[pde_index]->rGetDependentVariableName());
660  (*mpVizPdeSolutionResultsFile) << solution << " ";
661  }
662  }
663  }
664  (*mpVizPdeSolutionResultsFile) << "\n";
665 #ifdef CHASTE_VTK
666  if (p_vtk_mesh_writer != NULL)
667  {
668  p_vtk_mesh_writer->WriteFilesUsingMesh(*mpCoarsePdeMesh);
669  delete p_vtk_mesh_writer;
670  }
671 #endif //CHASTE_VTK
672  }
673 }
674 
675 template<unsigned DIM>
676 void CellBasedPdeHandler<DIM>::SetWriteAverageRadialPdeSolution(const std::string& rName, unsigned numRadialIntervals, bool writeDailyResults)
677 {
678  mAverageRadialSolutionVariableName = rName;
679  mWriteAverageRadialPdeSolution = true;
680  mNumRadialIntervals = numRadialIntervals;
681  mWriteDailyAverageRadialPdeSolution = writeDailyResults;
682 }
683 
684 template<unsigned DIM>
686 {
687  mSetBcsOnCoarseBoundary = setBcsOnCoarseBoundary;
688 }
689 
690 template<unsigned DIM>
692 {
693  (*mpAverageRadialPdeSolutionResultsFile) << time << " ";
694 
695  // Calculate the centre of the cell population
696  c_vector<double,DIM> centre = mpCellPopulation->GetCentroidOfCellPopulation();
697 
698  // Calculate the distance between each node and the centre of the cell population, as well as the maximum of these
699  std::map<double, CellPtr> radius_cell_map;
700  double max_distance_from_centre = 0.0;
701  for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = mpCellPopulation->Begin();
702  cell_iter != mpCellPopulation->End();
703  ++cell_iter)
704  {
705  double distance = norm_2(mpCellPopulation->GetLocationOfCellCentre(*cell_iter) - centre);
706  radius_cell_map[distance] = *cell_iter;
707 
708  if (distance > max_distance_from_centre)
709  {
710  max_distance_from_centre = distance;
711  }
712  }
713 
714  // Create vector of radius intervals
715  std::vector<double> radius_intervals;
716  for (unsigned i=0; i<mNumRadialIntervals; i++)
717  {
718  double upper_radius = max_distance_from_centre*((double) i+1)/((double) mNumRadialIntervals);
719  radius_intervals.push_back(upper_radius);
720  }
721 
722  // Calculate PDE solution in each radial interval
723  double lower_radius = 0.0;
724  for (unsigned i=0; i<mNumRadialIntervals; i++)
725  {
726  unsigned counter = 0;
727  double average_solution = 0.0;
728 
729  for (std::map<double, CellPtr>::iterator iter = radius_cell_map.begin(); iter != radius_cell_map.end(); ++iter)
730  {
731  if (iter->first > lower_radius && iter->first <= radius_intervals[i])
732  {
733  average_solution += (iter->second)->GetCellData()->GetItem(mAverageRadialSolutionVariableName);
734  counter++;
735  }
736  }
737  if (counter != 0)
738  {
739  average_solution /= (double) counter;
740  }
741 
742  // Write results to file
743  (*mpAverageRadialPdeSolutionResultsFile) << radius_intervals[i] << " " << average_solution << " ";
744  lower_radius = radius_intervals[i];
745  }
746  (*mpAverageRadialPdeSolutionResultsFile) << "\n";
747 }
748 
749 template<unsigned DIM>
750 double CellBasedPdeHandler<DIM>::GetPdeSolutionAtPoint(const c_vector<double,DIM>& rPoint, const std::string& rVariable)
751 {
752  double solution_at_point = 0.0;
753 
754  unsigned pde_index = UINT_MAX;
755 
756  // Loop over elements of mPdeAndBcCollection to find correct PDE
757  for (unsigned i=0; i<mPdeAndBcCollection.size(); i++)
758  {
759  if (mPdeAndBcCollection[i]->rGetDependentVariableName() == rVariable)
760  {
761  pde_index = i;
762  break;
763  }
764  }
765  if (pde_index == UINT_MAX)
766  {
767  EXCEPTION("Tried to get the solution of a variable name: " + rVariable + ". There is no PDE with that variable.");
768  }
769  PdeAndBoundaryConditions<DIM>* p_pde_and_bc = mPdeAndBcCollection[pde_index];
770 
771  Element<DIM,DIM>* p_containing_element;
772 
773  if (mpCoarsePdeMesh != NULL)
774  {
775  // Find PDE element containing point
776  unsigned elem_index = mpCoarsePdeMesh->GetContainingElementIndex(ChastePoint<DIM>(rPoint));
777  p_containing_element = mpCoarsePdeMesh->GetElement(elem_index);
778  }
779  else // Tetrahedral mesh
780  {
781  // If not using a coarse PDE mesh, we must be using a MeshBasedCellPopulation
782  TetrahedralMesh<DIM,DIM>* p_tetrahedral_mesh = &(static_cast<MeshBasedCellPopulation<DIM>*>(mpCellPopulation)->rGetMesh());
783 
784  unsigned elem_index = p_tetrahedral_mesh->GetContainingElementIndex(ChastePoint<DIM>(rPoint));
785  p_containing_element = p_tetrahedral_mesh->GetElement(elem_index);
786  }
787 
788  // Interpolate solution
789  if (p_pde_and_bc->GetSolution())
790  {
791  ReplicatableVector solution_repl(p_pde_and_bc->GetSolution());
792  c_vector<double,DIM+1> weights = p_containing_element->CalculateInterpolationWeights(rPoint);
793  for (unsigned i=0; i<DIM+1; i++)
794  {
795  double nodal_value = solution_repl[p_containing_element->GetNodeGlobalIndex(i)];
796  solution_at_point += nodal_value * weights(i);
797  }
798  }
799  else
800  {
801  // should only come into this method AFTER solving the PDE
803  }
804 
805  return solution_at_point;
806 }
807 
808 template<unsigned DIM>
810 {
811  return mWriteAverageRadialPdeSolution;
812 }
813 
814 template<unsigned DIM>
816 {
817  return mWriteDailyAverageRadialPdeSolution;
818 }
819 
820 template<unsigned DIM>
822 {
823  return mSetBcsOnCoarseBoundary;
824 }
825 
826 template<unsigned DIM>
828 {
829  return mNumRadialIntervals;
830 }
831 
832 template<unsigned DIM>
833 void CellBasedPdeHandler<DIM>::OutputParameters(out_stream& rParamsFile)
834 {
835  std::string type = GetIdentifier();
836 
837  *rParamsFile << "\t\t<" << type << ">\n";
838  *rParamsFile << "\t\t<WriteAverageRadialPdeSolution>" << mWriteAverageRadialPdeSolution << "</WriteAverageRadialPdeSolution>\n";
839  *rParamsFile << "\t\t<WriteDailyAverageRadialPdeSolution>" << mWriteDailyAverageRadialPdeSolution << "</WriteDailyAverageRadialPdeSolution>\n";
840  *rParamsFile << "\t\t<SetBcsOnCoarseBoundary>" << mSetBcsOnCoarseBoundary << "</SetBcsOnCoarseBoundary>\n";
841  *rParamsFile << "\t\t<NumRadialIntervals>" << mNumRadialIntervals << "</NumRadialIntervals>\n";
842  *rParamsFile << "\t\t</" << type << ">\n";
843 }
844 
845 template<unsigned DIM>
847 {
848  return ((dynamic_cast<NodeBasedCellPopulation<DIM>*>(mpCellPopulation) != NULL)
849  || (dynamic_cast<PottsBasedCellPopulation<DIM>*>(mpCellPopulation) != NULL)
850  || (dynamic_cast<CaBasedCellPopulation<DIM>*>(mpCellPopulation) != NULL)
851  || (dynamic_cast<VertexBasedCellPopulation<DIM>*>(mpCellPopulation) != NULL));
852 }
853 
854 // Serialization for Boost >= 1.36
857 
858 template class CellBasedPdeHandler<1>;
860 template class CellBasedPdeHandler<2>;
861 template class CellBasedPdeHandler<3>;
double GetPdeSolutionAtPoint(const c_vector< double, DIM > &rPoint, const std::string &rVariable)
c_vector< double, DIM > & rGetLocation()
Definition: ChastePoint.cpp:76
BoundaryElementIterator GetBoundaryElementIteratorBegin() const
unsigned GetNodeGlobalIndex(unsigned localIndex) const
void AddPdeAndBc(PdeAndBoundaryConditions< DIM > *pPdeAndBc)
void ConstructRegularSlabMesh(double spaceStep, double width, double height=0, double depth=0)
CellBasedPdeHandler(AbstractCellPopulation< DIM > *pCellPopulation, bool deleteMemberPointersInDestructor=false)
Node< SPACE_DIM > * GetNode(unsigned index) const
virtual unsigned GetNumElements() const
#define EXCEPTION(message)
Definition: Exception.hpp:143
static SimulationTime * Instance()
Element< ELEMENT_DIM, SPACE_DIM > * GetElement(unsigned index) const
static bool AmMaster()
Definition: PetscTools.cpp:120
virtual unsigned GetNumNodes() const
unsigned GetTimeStepsElapsed() const
unsigned FindCoarseElementContainingCell(CellPtr pCell)
BoundaryNodeIterator GetBoundaryNodeIteratorBegin() const
double GetTimeStep() const
#define NEVER_REACHED
Definition: Exception.hpp:206
virtual void UseCoarsePdeMesh(double stepSize, ChasteCuboid< DIM > meshCuboid, bool centreOnCellPopulation=false)
void SetImposeBcsOnCoarseBoundary(bool setBcsOnCoarseBoundary)
Node< SPACE_DIM > * GetNode(unsigned localIndex) const
void WriteFilesUsingMesh(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, bool keepOriginalElementIndexing=true)
out_stream OpenOutputFile(const std::string &rFileName, std::ios_base::openmode mode=std::ios::out|std::ios::trunc) const
const ChastePoint< SPACE_DIM > & rGetUpperCorner() const
#define EXPORT_TEMPLATE_CLASS_SAME_DIMS(CLASS)
virtual void SolvePdeAndWriteResultsToFile(unsigned samplingTimestepMultiple)
void WriteAverageRadialPdeSolution(double time)
void OpenResultsFiles(std::string outputDirectory)
const unsigned UNSIGNED_UNSET
Definition: Exception.hpp:52
unsigned GetContainingElementIndex(const ChastePoint< SPACE_DIM > &rTestPoint, bool strict=false, std::set< unsigned > testElements=std::set< unsigned >(), bool onlyTryWithTestElements=false)
TetrahedralMesh< DIM, DIM > * GetCoarsePdeMesh()
AbstractCellPopulation< DIM > * mpCellPopulation
virtual void OutputParameters(out_stream &rParamsFile)
double GetWidth(unsigned rDimension) const
static void Destroy(Vec &rVec)
Definition: PetscTools.hpp:351
virtual std::auto_ptr< BoundaryConditionsContainer< DIM, DIM, 1 > > ConstructBoundaryConditionsContainer(PdeAndBoundaryConditions< DIM > *pPdeAndBc, TetrahedralMesh< DIM, DIM > *pMesh)
double GetTime() const
void SetWriteAverageRadialPdeSolution(const std::string &rName, unsigned numRadialIntervals=10, bool writeDailyResults=false)
AbstractBoundaryCondition< DIM > * GetBoundaryCondition() const
virtual void WritePdeSolution(double time)
BoundaryNodeIterator GetBoundaryNodeIteratorEnd() const
AbstractLinearEllipticPde< DIM, DIM > * GetPde()
virtual void WriteFilesUsingMesh(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, bool keepOriginalElementIndexing=true)
void SetUpSourceTermsForAveragedSourcePde(TetrahedralMesh< DIM, DIM > *pMesh, std::map< CellPtr, unsigned > *pCellPdeElementMap=NULL)
c_vector< double, SPACE_DIM+1 > CalculateInterpolationWeights(const ChastePoint< SPACE_DIM > &rTestPoint)
Definition: Element.cpp:228
void AddPointData(std::string name, std::vector< double > data)
const AbstractCellPopulation< DIM > * GetCellPopulation() const
BoundaryElementIterator GetBoundaryElementIteratorEnd() const
virtual Vec GetPdeSolution(const std::string &rName="")
const ChastePoint< SPACE_DIM > & rGetLowerCorner() const