Chaste  Release::2017.1
MeshBasedCellPopulation.cpp
1 /*
2 
3 Copyright (c) 2005-2017, 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 "MeshBasedCellPopulation.hpp"
37 #include "VtkMeshWriter.hpp"
38 #include "CellBasedEventHandler.hpp"
39 #include "Cylindrical2dMesh.hpp"
40 #include "Cylindrical2dVertexMesh.hpp"
41 #include "NodesOnlyMesh.hpp"
42 #include "CellId.hpp"
43 #include "CellVolumesWriter.hpp"
44 #include "CellPopulationElementWriter.hpp"
45 #include "VoronoiDataWriter.hpp"
46 #include "NodeVelocityWriter.hpp"
47 #include "CellPopulationAreaWriter.hpp"
48 
49 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
51  std::vector<CellPtr>& rCells,
52  const std::vector<unsigned> locationIndices,
53  bool deleteMesh,
54  bool validate)
55  : AbstractCentreBasedCellPopulation<ELEMENT_DIM,SPACE_DIM>(rMesh, rCells, locationIndices),
56  mpVoronoiTessellation(nullptr),
57  mDeleteMesh(deleteMesh),
58  mUseAreaBasedDampingConstant(false),
59  mAreaBasedDampingConstantParameter(0.1),
60  mWriteVtkAsPoints(false),
61  mOutputMeshInVtk(false),
62  mHasVariableRestLength(false)
63 {
64  mpMutableMesh = static_cast<MutableMesh<ELEMENT_DIM,SPACE_DIM>* >(&(this->mrMesh));
65 
66  assert(this->mCells.size() <= this->mrMesh.GetNumNodes());
67 
68  if (validate)
69  {
70  Validate();
71  }
72 
73  // Initialise the applied force at each node to zero
74  for (typename AbstractMesh<ELEMENT_DIM, SPACE_DIM>::NodeIterator node_iter = this->rGetMesh().GetNodeIteratorBegin();
75  node_iter != this->rGetMesh().GetNodeIteratorEnd();
76  ++node_iter)
77  {
78  node_iter->ClearAppliedForce();
79  }
80 }
81 
82 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
84  : AbstractCentreBasedCellPopulation<ELEMENT_DIM,SPACE_DIM>(rMesh)
85 {
86  mpMutableMesh = static_cast<MutableMesh<ELEMENT_DIM,SPACE_DIM>* >(&(this->mrMesh));
87  mpVoronoiTessellation = nullptr;
88  mDeleteMesh = true;
89 }
90 
91 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
93 {
94  delete mpVoronoiTessellation;
95 
96  if (mDeleteMesh)
97  {
98  delete &this->mrMesh;
99  }
100 }
101 
102 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
104 {
106 }
107 
108 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
110 {
111  assert(SPACE_DIM == 2); // LCOV_EXCL_LINE
112  mUseAreaBasedDampingConstant = useAreaBasedDampingConstant;
113 }
114 
115 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
117 {
118  return mpMutableMesh->AddNode(pNewNode);
119 }
120 
121 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
123 {
124  static_cast<MutableMesh<ELEMENT_DIM,SPACE_DIM>&>((this->mrMesh)).SetNode(nodeIndex, rNewLocation, false);
125 }
126 
127 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
129 {
131 
132  /*
133  * The next code block computes the area-dependent damping constant as given by equation
134  * (5) in the following reference: van Leeuwen et al. 2009. An integrative computational model
135  * for intestinal tissue renewal. Cell Prolif. 42(5):617-636. doi:10.1111/j.1365-2184.2009.00627.x
136  */
138  {
147  assert(SPACE_DIM == 2); // LCOV_EXCL_LINE
148 
149  double rest_length = 1.0;
151 
157  double d1 = 2.0*(1.0 - d0)/(sqrt(3.0)*rest_length*rest_length);
158 
159  double area_cell = GetVolumeOfVoronoiElement(nodeIndex);
160 
166  assert(area_cell < 1000);
167 
168  damping_multiplier = d0 + area_cell*d1;
169  }
170 
171  return damping_multiplier;
172 }
173 
174 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
176 {
177  std::vector<bool> validated_node = std::vector<bool>(this->GetNumNodes(), false);
178 
179  for (typename AbstractCellPopulation<ELEMENT_DIM,SPACE_DIM>::Iterator cell_iter=this->Begin(); cell_iter!=this->End(); ++cell_iter)
180  {
181  unsigned node_index = this->GetLocationIndexUsingCell(*cell_iter);
182  validated_node[node_index] = true;
183  }
184 
185  for (unsigned i=0; i<validated_node.size(); i++)
186  {
187  if (!validated_node[i])
188  {
189  EXCEPTION("At time " << SimulationTime::Instance()->GetTime() << ", Node " << i << " does not appear to have a cell associated with it");
190  }
191  }
192 }
193 
194 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
196 {
197  return *mpMutableMesh;
198 }
199 
200 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
202 {
203  return *mpMutableMesh;
204 }
205 
206 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
208 {
209  return mpMutableMesh;
210 }
211 
212 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
214 {
215  unsigned num_removed = 0;
216  for (std::list<CellPtr>::iterator it = this->mCells.begin();
217  it != this->mCells.end();
218  )
219  {
220  if ((*it)->IsDead())
221  {
222  // Check if this cell is in a marked spring
223  std::vector<const std::pair<CellPtr,CellPtr>*> pairs_to_remove; // Pairs that must be purged
224  for (std::set<std::pair<CellPtr,CellPtr> >::iterator it1 = this->mMarkedSprings.begin();
225  it1 != this->mMarkedSprings.end();
226  ++it1)
227  {
228  const std::pair<CellPtr,CellPtr>& r_pair = *it1;
229 
230  for (unsigned i=0; i<2; i++)
231  {
232  CellPtr p_cell = (i==0 ? r_pair.first : r_pair.second);
233 
234  if (p_cell == *it)
235  {
236  // Remember to purge this spring
237  pairs_to_remove.push_back(&r_pair);
238  break;
239  }
240  }
241  }
242 
243  // Purge any marked springs that contained this cell
244  for (std::vector<const std::pair<CellPtr,CellPtr>* >::iterator pair_it = pairs_to_remove.begin();
245  pair_it != pairs_to_remove.end();
246  ++pair_it)
247  {
248  this->mMarkedSprings.erase(**pair_it);
249  }
250 
251  // Remove the node from the mesh
252  num_removed++;
254 
255  // Update mappings between cells and location indices
256  unsigned location_index_of_removed_node = this->GetLocationIndexUsingCell((*it));
257  this->RemoveCellUsingLocationIndex(location_index_of_removed_node, (*it));
258 
259  // Update vector of cells
260  it = this->mCells.erase(it);
261  }
262  else
263  {
264  ++it;
265  }
266  }
267 
268  return num_removed;
269 }
270 
271 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
273 {
275  bool output_node_velocities = (this-> template HasWriter<NodeVelocityWriter>());
276 
287  std::map<unsigned, double> old_node_radius_map;
288  old_node_radius_map.clear();
289  if (this->mrMesh.GetNodeIteratorBegin()->HasNodeAttributes())
290  {
291  if (this->mrMesh.GetNodeIteratorBegin()->GetRadius() > 0.0)
292  {
293  for (typename AbstractMesh<ELEMENT_DIM, SPACE_DIM>::NodeIterator node_iter = this->mrMesh.GetNodeIteratorBegin();
294  node_iter != this->mrMesh.GetNodeIteratorEnd();
295  ++node_iter)
296  {
297  unsigned node_index = node_iter->GetIndex();
298  old_node_radius_map[node_index] = node_iter->GetRadius();
299  }
300  }
301  }
302 
303  std::map<unsigned, c_vector<double, SPACE_DIM> > old_node_applied_force_map;
304  old_node_applied_force_map.clear();
305  if (output_node_velocities)
306  {
307  /*
308  * If outputting node velocities, we must keep a record of the applied force at each
309  * node, since this will be cleared during the remeshing process. We then restore
310  * these attributes to the nodes after calling ReMesh().
311  */
312  for (typename AbstractMesh<ELEMENT_DIM, SPACE_DIM>::NodeIterator node_iter = this->mrMesh.GetNodeIteratorBegin();
313  node_iter != this->mrMesh.GetNodeIteratorEnd();
314  ++node_iter)
315  {
316  unsigned node_index = node_iter->GetIndex();
317  old_node_applied_force_map[node_index] = node_iter->rGetAppliedForce();
318  }
319  }
320 
321  NodeMap node_map(this->mrMesh.GetNumAllNodes());
322 
323  // We must use a static_cast to call ReMesh() as this method is not defined in parent mesh classes
324  static_cast<MutableMesh<ELEMENT_DIM,SPACE_DIM>&>((this->mrMesh)).ReMesh(node_map);
325 
326  if (!node_map.IsIdentityMap())
327  {
328  UpdateGhostNodesAfterReMesh(node_map);
329 
330  // Update the mappings between cells and location indices
331  std::map<Cell*, unsigned> old_cell_location_map = this->mCellLocationMap;
332 
333  // Remove any dead pointers from the maps (needed to avoid archiving errors)
334  this->mLocationCellMap.clear();
335  this->mCellLocationMap.clear();
336 
337  for (std::list<CellPtr>::iterator it = this->mCells.begin(); it != this->mCells.end(); ++it)
338  {
339  unsigned old_node_index = old_cell_location_map[(*it).get()];
340 
341  // This shouldn't ever happen, as the cell vector only contains living cells
342  assert(!node_map.IsDeleted(old_node_index));
343 
344  unsigned new_node_index = node_map.GetNewIndex(old_node_index);
345  this->SetCellUsingLocationIndex(new_node_index,*it);
346 
347  if (old_node_radius_map[old_node_index] > 0.0)
348  {
349  this->GetNode(new_node_index)->SetRadius(old_node_radius_map[old_node_index]);
350  }
351  if (output_node_velocities)
352  {
353  this->GetNode(new_node_index)->AddAppliedForceContribution(old_node_applied_force_map[old_node_index]);
354  }
355  }
356 
357  this->Validate();
358  }
359  else
360  {
361  if (old_node_radius_map[this->mCellLocationMap[(*(this->mCells.begin())).get()]] > 0.0)
362  {
363  for (std::list<CellPtr>::iterator it = this->mCells.begin(); it != this->mCells.end(); ++it)
364  {
365  unsigned node_index = this->mCellLocationMap[(*it).get()];
366  this->GetNode(node_index)->SetRadius(old_node_radius_map[node_index]);
367  }
368  }
369  if (output_node_velocities)
370  {
371  for (std::list<CellPtr>::iterator it = this->mCells.begin(); it != this->mCells.end(); ++it)
372  {
373  unsigned node_index = this->mCellLocationMap[(*it).get()];
374  this->GetNode(node_index)->AddAppliedForceContribution(old_node_applied_force_map[node_index]);
375  }
376  }
377  }
378 
379  // Purge any marked springs that are no longer springs
380  std::vector<const std::pair<CellPtr,CellPtr>*> springs_to_remove;
381  for (std::set<std::pair<CellPtr,CellPtr> >::iterator spring_it = this->mMarkedSprings.begin();
382  spring_it != this->mMarkedSprings.end();
383  ++spring_it)
384  {
385  CellPtr p_cell_1 = spring_it->first;
386  CellPtr p_cell_2 = spring_it->second;
387  Node<SPACE_DIM>* p_node_1 = this->GetNodeCorrespondingToCell(p_cell_1);
388  Node<SPACE_DIM>* p_node_2 = this->GetNodeCorrespondingToCell(p_cell_2);
389 
390  bool joined = false;
391 
392  // For each element containing node1, if it also contains node2 then the cells are joined
393  std::set<unsigned> node2_elements = p_node_2->rGetContainingElementIndices();
394  for (typename Node<SPACE_DIM>::ContainingElementIterator elem_iter = p_node_1->ContainingElementsBegin();
395  elem_iter != p_node_1->ContainingElementsEnd();
396  ++elem_iter)
397  {
398  if (node2_elements.find(*elem_iter) != node2_elements.end())
399  {
400  joined = true;
401  break;
402  }
403  }
404 
405  // If no longer joined, remove this spring from the set
406  if (!joined)
407  {
408  springs_to_remove.push_back(&(*spring_it));
409  }
410  }
411 
412  // Remove any springs necessary
413  for (std::vector<const std::pair<CellPtr,CellPtr>* >::iterator spring_it = springs_to_remove.begin();
414  spring_it != springs_to_remove.end();
415  ++spring_it)
416  {
417  this->mMarkedSprings.erase(**spring_it);
418  }
419 
420  // Tessellate if needed
422 
424 }
425 
426 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
428 {
429  if ((SPACE_DIM==2 || SPACE_DIM==3)&&(ELEMENT_DIM==SPACE_DIM))
430  {
431  CellBasedEventHandler::BeginEvent(CellBasedEventHandler::TESSELLATION);
433  this-> template HasWriter<VoronoiDataWriter>() ||
434  this-> template HasWriter<CellPopulationAreaWriter>() ||
435  this-> template HasWriter<CellVolumesWriter>())
436  {
438  }
439  CellBasedEventHandler::EndEvent(CellBasedEventHandler::TESSELLATION);
440  }
441 }
442 
443 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
445 {
446  // Only implemented for 2D elements
447  assert(ELEMENT_DIM == 2); // LCOV_EXCL_LINE
448 
449  std::vector<c_vector<unsigned, 5> > new_nodes;
450  new_nodes = rGetMesh().SplitLongEdges(springDivisionThreshold);
451 
452  // Add new cells onto new nodes
453  for (unsigned index=0; index<new_nodes.size(); index++)
454  {
455  // Copy the cell attached to one of the neighbouring nodes onto the new node
456  unsigned new_node_index = new_nodes[index][0];
457  unsigned node_a_index = new_nodes[index][1];
458  unsigned node_b_index = new_nodes[index][2];
459 
460  CellPtr p_neighbour_cell = this->GetCellUsingLocationIndex(node_a_index);
461 
462  // Create copy of cell property collection to modify for daughter cell
463  CellPropertyCollection daughter_property_collection = p_neighbour_cell->rGetCellPropertyCollection();
464 
465  // Remove the CellId from the daughter cell a new one will be assigned in the constructor
466  daughter_property_collection.RemoveProperty<CellId>();
467 
468  CellPtr p_new_cell(new Cell(p_neighbour_cell->GetMutationState(),
469  p_neighbour_cell->GetCellCycleModel()->CreateCellCycleModel(),
470  p_neighbour_cell->GetSrnModel()->CreateSrnModel(),
471  false,
472  daughter_property_collection));
473 
474  // Add new cell to cell population
475  this->mCells.push_back(p_new_cell);
476  this->AddCellUsingLocationIndex(new_node_index,p_new_cell);
477 
478  // Update rest lengths
479 
480  // Remove old node pair // note node_a_index < node_b_index
481  std::pair<unsigned,unsigned> node_pair = this->CreateOrderedPair(node_a_index, node_b_index);
482  double old_rest_length = mSpringRestLengths[node_pair];
483 
484  std::map<std::pair<unsigned,unsigned>, double>::iterator iter = mSpringRestLengths.find(node_pair);
485  mSpringRestLengths.erase(iter);
486 
487  // Add new pairs
488  node_pair = this->CreateOrderedPair(node_a_index, new_node_index);
489  mSpringRestLengths[node_pair] = 0.5*old_rest_length;
490 
491  node_pair = this->CreateOrderedPair(node_b_index, new_node_index);
492  mSpringRestLengths[node_pair] = 0.5*old_rest_length;
493 
494  // If necessary add other new spring rest lengths
495  for (unsigned pair_index=3; pair_index<5; pair_index++)
496  {
497  unsigned other_node_index = new_nodes[index][pair_index];
498 
499  if (other_node_index != UNSIGNED_UNSET)
500  {
501  node_pair = this->CreateOrderedPair(other_node_index, new_node_index);
502  double new_rest_length = rGetMesh().GetDistanceBetweenNodes(new_node_index, other_node_index);
503  mSpringRestLengths[node_pair] = new_rest_length;
504  }
505  }
506  }
507 }
508 
509 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
511 {
512  return this->mrMesh.GetNode(index);
513 }
514 
515 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
517 {
518  return this->mrMesh.GetNumAllNodes();
519 }
520 
521 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
523 {
524 }
525 
526 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
527 CellPtr MeshBasedCellPopulation<ELEMENT_DIM,SPACE_DIM>::AddCell(CellPtr pNewCell, CellPtr pParentCell)
528 {
529  assert(pNewCell);
530  assert(pParentCell);
531 
532  // Add new cell to population
533  CellPtr p_created_cell = AbstractCentreBasedCellPopulation<ELEMENT_DIM,SPACE_DIM>::AddCell(pNewCell, pParentCell);
534  assert(p_created_cell == pNewCell);
535 
536  // Mark spring between parent cell and new cell
537  std::pair<CellPtr,CellPtr> cell_pair = this->CreateCellPair(pParentCell, p_created_cell);
538  this->MarkSpring(cell_pair);
539 
540  // Return pointer to new cell
541  return p_created_cell;
542 }
543 
545 // Output methods //
547 
548 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
550 {
552  {
553  if (!this-> template HasWriter<CellPopulationElementWriter>())
554  {
555  this-> template AddPopulationWriter<CellPopulationElementWriter>();
556  }
557  }
558 
560 }
561 
562 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
564 {
565  if (SimulationTime::Instance()->GetTimeStepsElapsed() == 0 && this->mpVoronoiTessellation == nullptr)
566  {
567  TessellateIfNeeded(); // Update isn't run on time-step zero
568  }
569 
571 }
572 
573 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
575 {
576  pPopulationWriter->Visit(this);
577 }
578 
579 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
581 {
582  pPopulationCountWriter->Visit(this);
583 }
584 
585 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
587 {
588  pCellWriter->VisitCell(pCell, this);
589 }
590 
591 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
593 {
594 #ifdef CHASTE_VTK
595  // Store the present time as a string
596  unsigned num_timesteps = SimulationTime::Instance()->GetTimeStepsElapsed();
597  std::stringstream time;
598  time << num_timesteps;
599 
600  // Store the number of cells for which to output data to VTK
601  unsigned num_cells_from_mesh = GetNumNodes();
602  if (!mWriteVtkAsPoints && (mpVoronoiTessellation != nullptr))
603  {
604  num_cells_from_mesh = mpVoronoiTessellation->GetNumElements();
605  }
606 
607  // When outputting any CellData, we assume that the first cell is representative of all cells
608  unsigned num_cell_data_items = this->Begin()->GetCellData()->GetNumItems();
609  std::vector<std::string> cell_data_names = this->Begin()->GetCellData()->GetKeys();
610 
611  std::vector<std::vector<double> > cell_data;
612  for (unsigned var=0; var<num_cell_data_items; var++)
613  {
614  std::vector<double> cell_data_var(num_cells_from_mesh);
615  cell_data.push_back(cell_data_var);
616  }
617 
618  if (mOutputMeshInVtk)
619  {
620  // Create mesh writer for VTK output
621  VtkMeshWriter<ELEMENT_DIM, SPACE_DIM> mesh_writer(rDirectory, "mesh_"+time.str(), false);
622  mesh_writer.WriteFilesUsingMesh(rGetMesh());
623  }
624 
625  if (mWriteVtkAsPoints)
626  {
627  // Create mesh writer for VTK output
628  VtkMeshWriter<SPACE_DIM, SPACE_DIM> cells_writer(rDirectory, "results_"+time.str(), false);
629 
630  // Iterate over any cell writers that are present
631  unsigned num_cells = this->GetNumAllCells();
632  for (typename std::vector<boost::shared_ptr<AbstractCellWriter<ELEMENT_DIM, SPACE_DIM> > >::iterator cell_writer_iter = this->mCellWriters.begin();
633  cell_writer_iter != this->mCellWriters.end();
634  ++cell_writer_iter)
635  {
636  // Create vector to store VTK cell data
637  std::vector<double> vtk_cell_data(num_cells);
638 
639  // Loop over cells
640  for (typename AbstractCellPopulation<ELEMENT_DIM,SPACE_DIM>::Iterator cell_iter = this->Begin();
641  cell_iter != this->End();
642  ++cell_iter)
643  {
644  // Get the node index corresponding to this cell
645  unsigned node_index = this->GetLocationIndexUsingCell(*cell_iter);
646 
647  // Populate the vector of VTK cell data
648  vtk_cell_data[node_index] = (*cell_writer_iter)->GetCellDataForVtkOutput(*cell_iter, this);
649  }
650 
651  cells_writer.AddPointData((*cell_writer_iter)->GetVtkCellDataName(), vtk_cell_data);
652  }
653 
654  // Loop over cells
655  for (typename AbstractCellPopulation<ELEMENT_DIM,SPACE_DIM>::Iterator cell_iter = this->Begin();
656  cell_iter != this->End();
657  ++cell_iter)
658  {
659  // Get the node index corresponding to this cell
660  unsigned node_index = this->GetLocationIndexUsingCell(*cell_iter);
661 
662  for (unsigned var=0; var<num_cell_data_items; var++)
663  {
664  cell_data[var][node_index] = cell_iter->GetCellData()->GetItem(cell_data_names[var]);
665  }
666  }
667  for (unsigned var=0; var<num_cell_data_items; var++)
668  {
669  cells_writer.AddPointData(cell_data_names[var], cell_data[var]);
670  }
671 
672  // Make a copy of the nodes in a disposable mesh for writing
673  {
674  std::vector<Node<SPACE_DIM>* > nodes;
675  for (unsigned index=0; index<this->mrMesh.GetNumNodes(); index++)
676  {
677  Node<SPACE_DIM>* p_node = this->mrMesh.GetNode(index);
678  nodes.push_back(p_node);
679  }
680 
682  mesh.ConstructNodesWithoutMesh(nodes, 1.5); // Arbitrary cut off as connectivity not used.
683  cells_writer.WriteFilesUsingMesh(mesh);
684  }
685 
686  *(this->mpVtkMetaFile) << " <DataSet timestep=\"";
687  *(this->mpVtkMetaFile) << num_timesteps;
688  *(this->mpVtkMetaFile) << "\" group=\"\" part=\"0\" file=\"results_";
689  *(this->mpVtkMetaFile) << num_timesteps;
690  *(this->mpVtkMetaFile) << ".vtu\"/>\n";
691  }
692  else if (mpVoronoiTessellation != nullptr)
693  {
694  // Create mesh writer for VTK output
695  VertexMeshWriter<ELEMENT_DIM, SPACE_DIM> mesh_writer(rDirectory, "results", false);
696  std::vector<double> cell_volumes(num_cells_from_mesh);
697 
698  // Iterate over any cell writers that are present
699  unsigned num_cells = this->GetNumAllCells();
700  for (typename std::vector<boost::shared_ptr<AbstractCellWriter<ELEMENT_DIM, SPACE_DIM> > >::iterator cell_writer_iter = this->mCellWriters.begin();
701  cell_writer_iter != this->mCellWriters.end();
702  ++cell_writer_iter)
703  {
704  // Create vector to store VTK cell data
705  std::vector<double> vtk_cell_data(num_cells);
706 
707  // Loop over elements of mpVoronoiTessellation
708  for (typename VertexMesh<ELEMENT_DIM, SPACE_DIM>::VertexElementIterator elem_iter = mpVoronoiTessellation->GetElementIteratorBegin();
709  elem_iter != mpVoronoiTessellation->GetElementIteratorEnd();
710  ++elem_iter)
711  {
712  // Get index of this element in mpVoronoiTessellation
713  unsigned elem_index = elem_iter->GetIndex();
714 
715  // Get the cell corresponding to this element, via the index of the corresponding node in mrMesh
716  unsigned node_index = mpVoronoiTessellation->GetDelaunayNodeIndexCorrespondingToVoronoiElementIndex(elem_index);
717  CellPtr p_cell = this->GetCellUsingLocationIndex(node_index);
718 
719  // Populate the vector of VTK cell data
720  vtk_cell_data[elem_index] = (*cell_writer_iter)->GetCellDataForVtkOutput(p_cell, this);
721  }
722 
723  mesh_writer.AddCellData((*cell_writer_iter)->GetVtkCellDataName(), vtk_cell_data);
724  }
725 
726  // Loop over elements of mpVoronoiTessellation
727  for (typename VertexMesh<ELEMENT_DIM, SPACE_DIM>::VertexElementIterator elem_iter = mpVoronoiTessellation->GetElementIteratorBegin();
728  elem_iter != mpVoronoiTessellation->GetElementIteratorEnd();
729  ++elem_iter)
730  {
731  // Get index of this element in mpVoronoiTessellation
732  unsigned elem_index = elem_iter->GetIndex();
733 
734  // Get the cell corresponding to this element, via the index of the corresponding node in mrMesh
735  unsigned node_index = mpVoronoiTessellation->GetDelaunayNodeIndexCorrespondingToVoronoiElementIndex(elem_index);
736  CellPtr p_cell = this->GetCellUsingLocationIndex(node_index);
737 
738  for (unsigned var=0; var<num_cell_data_items; var++)
739  {
740  cell_data[var][elem_index] = p_cell->GetCellData()->GetItem(cell_data_names[var]);
741  }
742  }
743 
744  for (unsigned var=0; var<cell_data.size(); var++)
745  {
746  mesh_writer.AddCellData(cell_data_names[var], cell_data[var]);
747  }
748 
749  mesh_writer.WriteVtkUsingMesh(*mpVoronoiTessellation, time.str());
750  *(this->mpVtkMetaFile) << " <DataSet timestep=\"";
751  *(this->mpVtkMetaFile) << num_timesteps;
752  *(this->mpVtkMetaFile) << "\" group=\"\" part=\"0\" file=\"results_";
753  *(this->mpVtkMetaFile) << num_timesteps;
754  *(this->mpVtkMetaFile) << ".vtu\"/>\n";
755  }
756 #endif //CHASTE_VTK
757 }
758 
759 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
761 {
762  double cell_volume = 0;
763 
764  if (ELEMENT_DIM == SPACE_DIM)
765  {
766  // Ensure that the Voronoi tessellation exists
767  if (mpVoronoiTessellation == nullptr)
768  {
770  }
771 
772  // Get the node index corresponding to this cell
773  unsigned node_index = this->GetLocationIndexUsingCell(pCell);
774 
775  // Try to get the element index of the Voronoi tessellation corresponding to this node index
776  try
777  {
778  unsigned element_index = mpVoronoiTessellation->GetVoronoiElementIndexCorrespondingToDelaunayNodeIndex(node_index);
779 
780  // Get the cell's volume from the Voronoi tessellation
781  cell_volume = mpVoronoiTessellation->GetVolumeOfElement(element_index);
782  }
783  catch (Exception&)
784  {
785  // If it doesn't exist this must be a boundary cell, so return infinite volume
786  cell_volume = DBL_MAX;
787  }
788  }
789  else if (SPACE_DIM==3 && ELEMENT_DIM==2)
790  {
791  unsigned node_index = this->GetLocationIndexUsingCell(pCell);
792 
793  Node<SPACE_DIM>* p_node = rGetMesh().GetNode(node_index);
794 
795  assert(!(p_node->rGetContainingElementIndices().empty()));
796 
797  for (typename Node<SPACE_DIM>::ContainingElementIterator elem_iter = p_node->ContainingElementsBegin();
798  elem_iter != p_node->ContainingElementsEnd();
799  ++elem_iter)
800  {
801  Element<ELEMENT_DIM,SPACE_DIM>* p_element = rGetMesh().GetElement(*elem_iter);
802 
803  c_matrix<double, SPACE_DIM, ELEMENT_DIM> jacob;
804  double det;
805 
806  p_element->CalculateJacobian(jacob, det);
807 
808  cell_volume += fabs(p_element->GetVolume(det));
809  }
810 
811  // This calculation adds a third of each element to the total area
812  cell_volume /= 3.0;
813  }
814  else
815  {
816  // Not implemented for other dimensions
818  }
819 
820  return cell_volume;
821 }
822 
823 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
825 {
826  mWriteVtkAsPoints = writeVtkAsPoints;
827 }
828 
829 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
831 {
832  return mWriteVtkAsPoints;
833 }
834 
835 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
837 {
838  mOutputMeshInVtk = outputMeshInVtk;
839 }
840 
841 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
843 {
844  return mOutputMeshInVtk;
845 }
846 
847 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
849 {
850  if (bool(dynamic_cast<Cylindrical2dMesh*>(&(this->mrMesh))))
851  {
852  *pVizSetupFile << "MeshWidth\t" << this->GetWidth(0) << "\n";
853  }
854 }
855 
857 // Spring iterator class //
859 
860 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
862 {
863  return mEdgeIter.GetNodeA();
864 }
865 
866 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
868 {
869  return mEdgeIter.GetNodeB();
870 }
871 
872 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
874 {
875  assert((*this) != mrCellPopulation.SpringsEnd());
876  return mrCellPopulation.GetCellUsingLocationIndex(mEdgeIter.GetNodeA()->GetIndex());
877 }
878 
879 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
881 {
882  assert((*this) != mrCellPopulation.SpringsEnd());
883  return mrCellPopulation.GetCellUsingLocationIndex(mEdgeIter.GetNodeB()->GetIndex());
884 }
885 
886 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
888 {
889  return (mEdgeIter != rOther.mEdgeIter);
890 }
891 
892 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
894 {
895  bool edge_is_ghost = false;
896 
897  do
898  {
899  ++mEdgeIter;
900  if (*this != mrCellPopulation.SpringsEnd())
901  {
902  bool a_is_ghost = mrCellPopulation.IsGhostNode(mEdgeIter.GetNodeA()->GetIndex());
903  bool b_is_ghost = mrCellPopulation.IsGhostNode(mEdgeIter.GetNodeB()->GetIndex());
904 
905  edge_is_ghost = (a_is_ghost || b_is_ghost);
906  }
907  }
908  while (*this!=mrCellPopulation.SpringsEnd() && edge_is_ghost);
909 
910  return (*this);
911 }
912 
913 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
917  : mrCellPopulation(rCellPopulation),
918  mEdgeIter(edgeIter)
919 {
920  if (mEdgeIter!=static_cast<MutableMesh<ELEMENT_DIM,SPACE_DIM>*>(&(this->mrCellPopulation.mrMesh))->EdgesEnd())
921  {
922  bool a_is_ghost = mrCellPopulation.IsGhostNode(mEdgeIter.GetNodeA()->GetIndex());
923  bool b_is_ghost = mrCellPopulation.IsGhostNode(mEdgeIter.GetNodeB()->GetIndex());
924 
925  if (a_is_ghost || b_is_ghost)
926  {
927  ++(*this);
928  }
929  }
930 }
931 
932 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
934 {
935  return SpringIterator(*this, static_cast<MutableMesh<ELEMENT_DIM,SPACE_DIM>&>((this->mrMesh)).EdgesBegin());
936 }
937 
938 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
940 {
941  return SpringIterator(*this, static_cast<MutableMesh<ELEMENT_DIM,SPACE_DIM>&>((this->mrMesh)).EdgesEnd());
942 }
943 
947 template<>
949 {
950  delete mpVoronoiTessellation;
951 
952  // Check if the mesh associated with this cell population is periodic
953  bool is_mesh_periodic = false;
954  if (bool(dynamic_cast<Cylindrical2dMesh*>(&mrMesh)))
955  {
956  is_mesh_periodic = true;
957  mpVoronoiTessellation = new Cylindrical2dVertexMesh(static_cast<Cylindrical2dMesh &>(this->mrMesh));
958  }
959  else
960  {
961  mpVoronoiTessellation = new VertexMesh<2, 2>(static_cast<MutableMesh<2, 2> &>((this->mrMesh)), is_mesh_periodic);
962  }
963 }
964 
968 // LCOV_EXCL_START
969 template<>
971 {
972  // We don't allow tessellation yet.
974 }
975 // LCOV_EXCL_STOP
976 
982 template<>
984 {
985  delete mpVoronoiTessellation;
986  mpVoronoiTessellation = new VertexMesh<3, 3>(static_cast<MutableMesh<3, 3> &>((this->mrMesh)));
987 }
988 
993 // LCOV_EXCL_START
994 template<>
996 {
997  // No 1D Voronoi tessellation
999 }
1000 // LCOV_EXCL_STOP
1001 
1002 
1007 // LCOV_EXCL_START
1008 template<>
1010 {
1011  // No 1D Voronoi tessellation
1012  NEVER_REACHED;
1013 }
1014 // LCOV_EXCL_STOP
1015 
1020 // LCOV_EXCL_START
1021 template<>
1023 {
1024  // No 1D Voronoi tessellation
1025  NEVER_REACHED;
1026 }
1027 // LCOV_EXCL_STOP
1028 
1029 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
1031 {
1032  assert(mpVoronoiTessellation!=nullptr);
1033  return mpVoronoiTessellation;
1034 }
1035 
1036 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
1038 {
1039  unsigned element_index = mpVoronoiTessellation->GetVoronoiElementIndexCorrespondingToDelaunayNodeIndex(index);
1040  double volume = mpVoronoiTessellation->GetVolumeOfElement(element_index);
1041  return volume;
1042 }
1043 
1044 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
1046 {
1047  unsigned element_index = mpVoronoiTessellation->GetVoronoiElementIndexCorrespondingToDelaunayNodeIndex(index);
1048  double surface_area = mpVoronoiTessellation->GetSurfaceAreaOfElement(element_index);
1049  return surface_area;
1050 }
1051 
1052 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
1054 {
1055  unsigned element_index1 = mpVoronoiTessellation->GetVoronoiElementIndexCorrespondingToDelaunayNodeIndex(index1);
1056  unsigned element_index2 = mpVoronoiTessellation->GetVoronoiElementIndexCorrespondingToDelaunayNodeIndex(index2);
1057  try
1058  {
1059  double edge_length = mpVoronoiTessellation->GetEdgeLength(element_index1, element_index2);
1060  return edge_length;
1061  }
1062  catch (Exception&)
1063  {
1064  // The edge was between two (potentially infinite) cells on the boundary of the mesh
1065  EXCEPTION("Spring iterator tried to calculate interaction between degenerate cells on the boundary of the mesh. Have you set ghost layers correctly?");
1066  }
1067 }
1068 
1069 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
1071 {
1072  bool res = true;
1073  for (std::list<CellPtr>::iterator it=this->mCells.begin();
1074  it!=this->mCells.end();
1075  ++it)
1076  {
1077  CellPtr p_cell = *it;
1078  assert(p_cell);
1079  AbstractCellCycleModel* p_model = p_cell->GetCellCycleModel();
1080  assert(p_model);
1081 
1082  // Check cell exists in cell population
1083  unsigned node_index = this->GetLocationIndexUsingCell(p_cell);
1084  std::cout << "Cell at node " << node_index << " addr " << p_cell << std::endl << std::flush;
1085  CellPtr p_cell_in_cell_population = this->GetCellUsingLocationIndex(node_index);
1086 // LCOV_EXCL_START //Debugging code. Shouldn't fail under normal conditions
1087  if (p_cell_in_cell_population != p_cell)
1088  {
1089  std::cout << " Mismatch with cell population" << std::endl << std::flush;
1090  res = false;
1091  }
1092 
1093  // Check model links back to cell
1094  if (p_model->GetCell() != p_cell)
1095  {
1096  std::cout << " Mismatch with cycle model" << std::endl << std::flush;
1097  res = false;
1098  }
1099  }
1100  UNUSED_OPT(res);
1101  assert(res);
1102 // LCOV_EXCL_STOP
1103 
1104  res = true;
1105  for (std::set<std::pair<CellPtr,CellPtr> >::iterator it1 = this->mMarkedSprings.begin();
1106  it1 != this->mMarkedSprings.end();
1107  ++it1)
1108  {
1109  const std::pair<CellPtr,CellPtr>& r_pair = *it1;
1110 
1111  for (unsigned i=0; i<2; i++)
1112  {
1113  CellPtr p_cell = (i==0 ? r_pair.first : r_pair.second);
1114 
1115  assert(p_cell);
1116  AbstractCellCycleModel* p_model = p_cell->GetCellCycleModel();
1117  assert(p_model);
1118  unsigned node_index = this->GetLocationIndexUsingCell(p_cell);
1119  std::cout << "Cell at node " << node_index << " addr " << p_cell << std::endl << std::flush;
1120 
1121 // LCOV_EXCL_START //Debugging code. Shouldn't fail under normal conditions
1122  // Check cell is alive
1123  if (p_cell->IsDead())
1124  {
1125  std::cout << " Cell is dead" << std::endl << std::flush;
1126  res = false;
1127  }
1128 
1129  // Check cell exists in cell population
1130  CellPtr p_cell_in_cell_population = this->GetCellUsingLocationIndex(node_index);
1131  if (p_cell_in_cell_population != p_cell)
1132  {
1133  std::cout << " Mismatch with cell population" << std::endl << std::flush;
1134  res = false;
1135  }
1136 
1137  // Check model links back to cell
1138  if (p_model->GetCell() != p_cell)
1139  {
1140  std::cout << " Mismatch with cycle model" << std::endl << std::flush;
1141  res = false;
1142  }
1143  }
1144 // LCOV_EXCL_STOP
1145  }
1146  assert(res);
1147 }
1148 
1149 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
1151 {
1153 }
1154 
1155 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
1157 {
1158  assert(areaBasedDampingConstantParameter >= 0.0);
1159  mAreaBasedDampingConstantParameter = areaBasedDampingConstantParameter;
1160 }
1161 
1162 // LCOV_EXCL_START
1163 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
1165 {
1166  //mNodePairs.Clear();
1167  NEVER_REACHED;
1168  return mNodePairs;
1169 }
1170 // LCOV_EXCL_STOP
1171 
1172 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
1174 {
1175  *rParamsFile << "\t\t<UseAreaBasedDampingConstant>" << mUseAreaBasedDampingConstant << "</UseAreaBasedDampingConstant>\n";
1176  *rParamsFile << "\t\t<AreaBasedDampingConstantParameter>" << mAreaBasedDampingConstantParameter << "</AreaBasedDampingConstantParameter>\n";
1177  *rParamsFile << "\t\t<WriteVtkAsPoints>" << mWriteVtkAsPoints << "</WriteVtkAsPoints>\n";
1178  *rParamsFile << "\t\t<OutputMeshInVtk>" << mOutputMeshInVtk << "</OutputMeshInVtk>\n";
1179  *rParamsFile << "\t\t<HasVariableRestLength>" << mHasVariableRestLength << "</HasVariableRestLength>\n";
1180 
1181  // Call method on direct parent class
1183 }
1184 
1185 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
1187 {
1188  // Call GetWidth() on the mesh
1189  double width = this->mrMesh.GetWidth(rDimension);
1190  return width;
1191 }
1192 
1193 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
1195 {
1196  // Get pointer to this node
1197  Node<SPACE_DIM>* p_node = this->mrMesh.GetNode(index);
1198 
1199  // Loop over containing elements
1200  std::set<unsigned> neighbouring_node_indices;
1201  for (typename Node<SPACE_DIM>::ContainingElementIterator elem_iter = p_node->ContainingElementsBegin();
1202  elem_iter != p_node->ContainingElementsEnd();
1203  ++elem_iter)
1204  {
1205  // Get pointer to this containing element
1206  Element<ELEMENT_DIM,SPACE_DIM>* p_element = static_cast<MutableMesh<ELEMENT_DIM,SPACE_DIM>&>((this->mrMesh)).GetElement(*elem_iter);
1207 
1208  // Loop over nodes contained in this element
1209  for (unsigned i=0; i<p_element->GetNumNodes(); i++)
1210  {
1211  // Get index of this node and add its index to the set if not the original node
1212  unsigned node_index = p_element->GetNodeGlobalIndex(i);
1213  if (node_index != index)
1214  {
1215  neighbouring_node_indices.insert(node_index);
1216  }
1217  }
1218  }
1219  return neighbouring_node_indices;
1220 }
1221 
1222 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
1224 {
1225  mSpringRestLengths.clear();
1226 
1227  // Iterate over all springs and add calculate separation of adjacent node pairs
1228  for (SpringIterator spring_iterator = SpringsBegin();
1229  spring_iterator != SpringsEnd();
1230  ++spring_iterator)
1231  {
1232  // Note that nodeA_global_index is always less than nodeB_global_index
1233  Node<SPACE_DIM>* p_nodeA = spring_iterator.GetNodeA();
1234  Node<SPACE_DIM>* p_nodeB = spring_iterator.GetNodeB();
1235 
1236  unsigned nodeA_global_index = p_nodeA->GetIndex();
1237  unsigned nodeB_global_index = p_nodeB->GetIndex();
1238 
1239  // Calculate the distance between nodes
1240  double separation = rGetMesh().GetDistanceBetweenNodes(nodeA_global_index, nodeB_global_index);
1241 
1242  // Order node indices
1243  std::pair<unsigned,unsigned> node_pair = this->CreateOrderedPair(nodeA_global_index, nodeB_global_index);
1244 
1245  mSpringRestLengths[node_pair] = separation;
1246  }
1247  mHasVariableRestLength = true;
1248 }
1249 
1250 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
1252 {
1254  {
1255  std::pair<unsigned,unsigned> node_pair = this->CreateOrderedPair(indexA, indexB);
1256  std::map<std::pair<unsigned,unsigned>, double>::const_iterator iter = mSpringRestLengths.find(node_pair);
1257 
1258  if (iter != mSpringRestLengths.end())
1259  {
1260  // Return the stored rest length
1261  return iter->second;
1262  }
1263  else
1264  {
1265  EXCEPTION("Tried to get a rest length of an edge that doesn't exist. You can only use variable rest lengths if SetUpdateCellPopulationRule is set on the simulation.");
1266  }
1267  }
1268  else
1269  {
1270  return 1.0;
1271  }
1272 }
1273 
1274 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
1275 void MeshBasedCellPopulation<ELEMENT_DIM,SPACE_DIM>::SetRestLength(unsigned indexA, unsigned indexB, double restLength)
1276 {
1278  {
1279  std::pair<unsigned,unsigned> node_pair = this->CreateOrderedPair(indexA, indexB);
1280  std::map<std::pair<unsigned,unsigned>, double>::iterator iter = mSpringRestLengths.find(node_pair);
1281 
1282  if (iter != mSpringRestLengths.end())
1283  {
1284  // modify the stored rest length
1285  iter->second = restLength;
1286  }
1287  else
1288  {
1289  EXCEPTION("Tried to set the rest length of an edge not in the mesh.");
1290  }
1291  }
1292  else
1293  {
1294  EXCEPTION("Tried to set a rest length in a simulation with fixed rest length. You can only use variable rest lengths if SetUpdateCellPopulationRule is set on the simulation.");
1295  }
1296 }
1297 
1298 // Explicit instantiation
1299 template class MeshBasedCellPopulation<1,1>;
1300 template class MeshBasedCellPopulation<1,2>;
1301 template class MeshBasedCellPopulation<2,2>;
1302 template class MeshBasedCellPopulation<1,3>;
1303 template class MeshBasedCellPopulation<2,3>;
1304 template class MeshBasedCellPopulation<3,3>;
1305 
1306 // Serialization for Boost >= 1.36
double GetSurfaceAreaOfVoronoiElement(unsigned index)
double GetVolumeOfVoronoiElement(unsigned index)
void AddCellData(std::string dataName, std::vector< double > dataPayload)
double GetWidth(const unsigned &rDimension)
double GetVoronoiEdgeLength(unsigned index1, unsigned index2)
virtual CellPtr AddCell(CellPtr pNewCell, CellPtr pParentCell)
Definition: Cell.hpp:89
unsigned AddNode(Node< SPACE_DIM > *pNewNode)
virtual CellPtr GetCellUsingLocationIndex(unsigned index)
Definition: Node.hpp:58
unsigned GetLocationIndexUsingCell(CellPtr pCell)
double GetVolume(double determinant) const
double GetDampingConstant(unsigned nodeIndex)
std::vector< boost::shared_ptr< AbstractCellWriter< ELEMENT_DIM, SPACE_DIM > > > mCellWriters
unsigned GetNodeGlobalIndex(unsigned localIndex) const
void DeleteNodePriorToReMesh(unsigned index)
virtual void RemoveCellUsingLocationIndex(unsigned index, CellPtr pCell)
virtual void AcceptCellWriter(boost::shared_ptr< AbstractCellWriter< ELEMENT_DIM, SPACE_DIM > > pCellWriter, CellPtr pCell)
std::vector< std::pair< Node< SPACE_DIM > *, Node< SPACE_DIM > * > > & rGetNodePairs()
virtual void WriteDataToVisualizerSetupFile(out_stream &pVizSetupFile)
#define EXCEPTION(message)
Definition: Exception.hpp:143
virtual void AcceptPopulationCountWriter(boost::shared_ptr< AbstractCellPopulationCountWriter< ELEMENT_DIM, SPACE_DIM > > pPopulationCountWriter)
static SimulationTime * Instance()
void SetNode(unsigned nodeIndex, ChastePoint< SPACE_DIM > &rNewLocation)
Element< ELEMENT_DIM, SPACE_DIM > * GetElement(unsigned index) const
std::set< unsigned > & rGetContainingElementIndices()
Definition: Node.cpp:300
virtual void WriteVtkResultsToFile(const std::string &rDirectory)
unsigned GetTimeStepsElapsed() const
#define NEVER_REACHED
Definition: Exception.hpp:206
virtual void WriteResultsToFiles(const std::string &rDirectory)
virtual void SetNode(unsigned index, ChastePoint< SPACE_DIM > point, bool concreteMove=true)
std::pair< CellPtr, CellPtr > CreateCellPair(CellPtr pCell1, CellPtr pCell2)
MutableMesh< ELEMENT_DIM, SPACE_DIM > & rGetMesh()
bool operator!=(const typename MeshBasedCellPopulation< ELEMENT_DIM, SPACE_DIM >::SpringIterator &rOther)
virtual void OutputCellPopulationParameters(out_stream &rParamsFile)
void WriteFilesUsingMesh(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, bool keepOriginalElementIndexing=true)
double GetRestLength(unsigned indexA, unsigned indexB)
void SetWriteVtkAsPoints(bool writeVtkAsPoints)
std::pair< unsigned, unsigned > CreateOrderedPair(unsigned index1, unsigned index2)
MutableMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator mEdgeIter
ContainingElementIterator ContainingElementsBegin() const
Definition: Node.hpp:485
MutableMesh< ELEMENT_DIM, SPACE_DIM > * mpMutableMesh
std::map< unsigned, std::set< CellPtr > > mLocationCellMap
void SetCellUsingLocationIndex(unsigned index, CellPtr pCell)
std::map< Cell *, unsigned > mCellLocationMap
void SetOutputMeshInVtk(bool outputMeshInVtk)
void DivideLongSprings(double springDivisionThreshold)
virtual void AcceptPopulationWriter(boost::shared_ptr< AbstractCellPopulationWriter< ELEMENT_DIM, SPACE_DIM > > pPopulationWriter)
Node< SPACE_DIM > * GetNodeCorrespondingToCell(CellPtr pCell)
VertexMesh< ELEMENT_DIM, SPACE_DIM > * mpVoronoiTessellation
unsigned GetNumNodes() const
SpringIterator(MeshBasedCellPopulation< ELEMENT_DIM, SPACE_DIM > &rCellPopulation, typename MutableMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator edgeIter)
virtual void Update(bool hasHadBirthsOrDeaths=true)
const unsigned UNSIGNED_UNSET
Definition: Exception.hpp:52
void WriteVtkUsingMesh(VertexMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, std::string stamp="")
virtual void ReMesh(NodeMap &map)
CellPtr AddCell(CellPtr pNewCell, CellPtr pParentCell=CellPtr())
#define EXPORT_TEMPLATE_CLASS_ALL_DIMS(CLASS)
void MarkSpring(std::pair< CellPtr, CellPtr > &rCellPair)
virtual void UpdateGhostNodesAfterReMesh(NodeMap &rMap)
void SetRestLength(unsigned indexA, unsigned indexB, double restLength)
void SetAreaBasedDampingConstant(bool useAreaBasedDampingConstant)
virtual void OpenWritersFiles(OutputFileHandler &rOutputFileHandler)
virtual void WriteResultsToFiles(const std::string &rDirectory)
void ConstructNodesWithoutMesh(const std::vector< Node< SPACE_DIM > * > &rNodes, double maxInteractionDistance)
std::set< unsigned > GetNeighbouringNodeIndices(unsigned index)
void SetMeshHasChangedSinceLoading()
void CalculateJacobian(c_matrix< double, SPACE_DIM, ELEMENT_DIM > &rJacobian, double &rJacobianDeterminant)
Node< SPACE_DIM > * GetNode(unsigned index)
#define UNUSED_OPT(var)
Definition: Exception.hpp:216
VertexMesh< ELEMENT_DIM, SPACE_DIM > * GetVoronoiTessellation()
MeshBasedCellPopulation(MutableMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, std::vector< CellPtr > &rCells, const std::vector< unsigned > locationIndices=std::vector< unsigned >(), bool deleteMesh=false, bool validate=true)
unsigned GetIndex() const
Definition: Node.cpp:158
double GetVolumeOfCell(CellPtr pCell)
void OutputCellPopulationParameters(out_stream &rParamsFile)
std::vector< std::pair< Node< SPACE_DIM > *, Node< SPACE_DIM > * > > mNodePairs
EdgeIterator EdgesEnd()
void SetAreaBasedDampingConstantParameter(double areaBasedDampingConstantParameter)
virtual void AddCellUsingLocationIndex(unsigned index, CellPtr pCell)
AbstractMesh< ELEMENT_DIM, SPACE_DIM > & mrMesh
virtual double GetDampingConstant(unsigned nodeIndex)
std::map< std::pair< unsigned, unsigned >, double > mSpringRestLengths
ContainingElementIterator ContainingElementsEnd() const
Definition: Node.hpp:493
virtual TetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * GetTetrahedralMeshForPdeModifier()
std::set< std::pair< CellPtr, CellPtr > > mMarkedSprings
virtual void OpenWritersFiles(OutputFileHandler &rOutputFileHandler)
MeshBasedCellPopulation< ELEMENT_DIM, SPACE_DIM > & mrCellPopulation