Chaste Commit::baa90ac2819b962188b7562f2326be23c47859a7
VertexBasedCellPopulation.cpp
1/*
2
3Copyright (c) 2005-2024, University of Oxford.
4All rights reserved.
5
6University of Oxford means the Chancellor, Masters and Scholars of the
7University of Oxford, having an administrative office at Wellington
8Square, Oxford OX1 2JD, UK.
9
10This file is part of Chaste.
11
12Redistribution and use in source and binary forms, with or without
13modification, 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
23THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33
34*/
35
36#include "VertexBasedCellPopulation.hpp"
37#include "Warnings.hpp"
38#include "ShortAxisVertexBasedDivisionRule.hpp"
39#include "StepSizeException.hpp"
40#include "WildTypeCellMutationState.hpp"
41#include "Cylindrical2dVertexMesh.hpp"
42#include "SmartPointers.hpp"
43#include "T2SwapCellKiller.hpp"
44#include "ApoptoticCellProperty.hpp"
45#include "CellSrnModel.hpp"
46#include "CellPopulationElementWriter.hpp"
47#include "VertexT1SwapLocationsWriter.hpp"
48#include "VertexT2SwapLocationsWriter.hpp"
49#include "VertexT3SwapLocationsWriter.hpp"
50#include "VertexIntersectionSwapLocationsWriter.hpp"
51#include "AbstractCellBasedSimulation.hpp"
52
53template<unsigned DIM>
55 std::vector<CellPtr>& rCells,
56 bool deleteMesh,
57 bool validate,
58 const std::vector<unsigned> locationIndices)
59 : AbstractOffLatticeCellPopulation<DIM>(rMesh, rCells, locationIndices),
60 mDeleteMesh(deleteMesh),
61 mOutputCellRearrangementLocations(true),
62 mRestrictVertexMovement(true)
63{
66
67 // If no location indices are specified, associate with elements from the mesh (assumed to be sequentially ordered).
68 std::list<CellPtr>::iterator it = this->mCells.begin();
69 for (unsigned i=0; it != this->mCells.end(); ++it, ++i)
70 {
71 unsigned index = locationIndices.empty() ? i : locationIndices[i]; // assume that the ordering matches
73 }
74
75 // Check each element has only one cell attached
76 if (validate)
77 {
78 Validate();
79 }
80
81 // If cells contain an SRN model, then we need to track mesh operations
82 // and update SRNs accordingly. Here we assume the first cell is a representative of other cells
83 if ((*this->mCells.begin())->HasSrnModel())
84 {
85 mPopulationSrn.SetVertexCellPopulation(this);
87 }
88}
89
90template<unsigned DIM>
94 mDeleteMesh(true),
95 mOutputCellRearrangementLocations(true),
96 mRestrictVertexMovement(true),
97 mPopulationSrn(rPopSrn)
98{
100 mPopulationSrn.SetVertexCellPopulation(this);
101}
102
103template<unsigned DIM>
105{
106 if (mDeleteMesh)
107 {
108 delete &this->mrMesh;
109 }
110}
111
112template<unsigned DIM>
114{
115 // Take the average of the cells containing this vertex
116 double average_damping_constant = 0.0;
117
118 std::set<unsigned> containing_elements = GetNode(nodeIndex)->rGetContainingElementIndices();
119
120 unsigned num_containing_elements = containing_elements.size();
121 if (num_containing_elements == 0)
122 {
123 EXCEPTION("At time " << SimulationTime::Instance()->GetTime() << ", Node " << nodeIndex << " is not contained in any elements, so GetDampingConstant() returns zero");
124 }
125
126 double temp = 1.0/((double) num_containing_elements);
127 for (std::set<unsigned>::iterator iter = containing_elements.begin();
128 iter != containing_elements.end();
129 ++iter)
130 {
131 CellPtr p_cell = this->GetCellUsingLocationIndex(*iter);
132 bool cell_is_wild_type = p_cell->GetMutationState()->IsType<WildTypeCellMutationState>();
133
134 if (cell_is_wild_type)
135 {
136 average_damping_constant += this->GetDampingConstantNormal()*temp;
137 }
138 else
139 {
140 average_damping_constant += this->GetDampingConstantMutant()*temp;
141 }
142 }
143
144 return average_damping_constant;
145}
146
147template<unsigned DIM>
149{
150 return *mpMutableVertexMesh;
151}
152
153template<unsigned DIM>
155{
156 return *mpMutableVertexMesh;
157}
158
159template<unsigned DIM>
161{
162 return mpMutableVertexMesh->GetElement(elementIndex);
163}
164
165template<unsigned DIM>
167{
168 return this->mrMesh.GetNumNodes();
169}
170
171template<unsigned DIM>
173{
174 return mpMutableVertexMesh->GetCentroidOfElement(this->mCellLocationMap[pCell.get()]);
175}
176
177template<unsigned DIM>
179{
180 return this->mrMesh.GetNode(index);
181}
182
183template<unsigned DIM>
185{
186 unsigned elem_index = this->GetLocationIndexUsingCell(pCell);
187 return this->rGetMesh().GetNeighbouringElementIndices(elem_index);
188}
189
190template<unsigned DIM>
191std::set<std::pair<unsigned, unsigned>>
193{
194 std::set<std::pair<unsigned, unsigned>> neighbours;
195 auto cellLocationIndex = this->GetLocationIndexUsingCell(pCell);
196 auto p_element = this->GetElement(cellLocationIndex);
197 auto global_edge_index = p_element->GetEdgeGlobalIndex(edgeLocalIndex);
198 auto neighbour_element_indices = p_element->GetNeighbouringElementAtEdgeIndex(edgeLocalIndex);
199
200 // Normally there is only one neighbouring element
201 for (auto neighbour_element_index : neighbour_element_indices)
202 {
203 auto p_neighbour_element = this->GetElement(neighbour_element_index);
204
205 // Iterate over neighbouring element indices
206 for (unsigned elem_index = 0; elem_index < p_neighbour_element->GetNumEdges(); elem_index++)
207 {
208 // If the neighbours edge matches EdgeLocalIndex
209 if (p_neighbour_element->GetEdge(elem_index)->GetIndex() == global_edge_index)
210 {
211 neighbours.insert(std::pair<unsigned, unsigned>(neighbour_element_index, elem_index));
212 }
213 }
214 }
215
216 return neighbours;
217}
218
219template<unsigned DIM>
221{
222 return mpMutableVertexMesh->AddNode(pNewNode);
223}
224
225template<unsigned DIM>
226void VertexBasedCellPopulation<DIM>::SetNode(unsigned nodeIndex, ChastePoint<DIM>& rNewLocation)
227{
228 mpMutableVertexMesh->SetNode(nodeIndex, rNewLocation);
229}
230
231template<unsigned DIM>
233{
234 return mpMutableVertexMesh->GetElement(this->GetLocationIndexUsingCell(pCell));
235}
236
237template<unsigned DIM>
239{
240 return mpMutableVertexMesh->GetNumElements();
241}
242
243template<unsigned DIM>
244CellPtr VertexBasedCellPopulation<DIM>::AddCell(CellPtr pNewCell, CellPtr pParentCell)
245{
246 // Get the element associated with this cell
247 VertexElement<DIM, DIM>* p_element = GetElementCorrespondingToCell(pParentCell);
248
249 // Get the orientation of division
250 c_vector<double, DIM> division_vector = mpVertexBasedDivisionRule->CalculateCellDivisionVector(pParentCell, *this);
251
252 // Divide the element
253 unsigned new_element_index = mpMutableVertexMesh->DivideElementAlongGivenAxis(p_element, division_vector, true);
254
255 // Associate the new cell with the element
256 this->mCells.push_back(pNewCell);
257
258 // Update location cell map
259 CellPtr p_created_cell = this->mCells.back();
260 this->SetCellUsingLocationIndex(new_element_index,p_created_cell);
261 this->mCellLocationMap[p_created_cell.get()] = new_element_index;
262
263 return p_created_cell;
264}
265
266template<unsigned DIM>
268{
269 unsigned num_removed = 0;
270
271 for (std::list<CellPtr>::iterator it = this->mCells.begin();
272 it != this->mCells.end();
273 )
274 {
275 if ((*it)->IsDead())
276 {
277 // Count the cell as dead
278 num_removed++;
279
280 // Remove the element from the mesh if it is not deleted yet
282 if (!(this->GetElement(this->GetLocationIndexUsingCell((*it)))->IsDeleted()))
283 {
284 // This warning relies on the fact that there is only one other possibility for
285 // vertex elements to be marked as deleted: a T2 swap
286 WARN_ONCE_ONLY("A Cell is removed without performing a T2 swap. This could leave a void in the mesh.");
287 mpMutableVertexMesh->DeleteElementPriorToReMesh(this->GetLocationIndexUsingCell((*it)));
288 }
289
290 // Delete the cell
291 it = this->mCells.erase(it);
292 }
293 else
294 {
295 ++it;
296 }
297 }
298 return num_removed;
299}
300
301template<unsigned DIM>
302void VertexBasedCellPopulation<DIM>::CheckForStepSizeException(unsigned nodeIndex, c_vector<double,DIM>& rDisplacement, double dt)
303{
304 double length = norm_2(rDisplacement);
305
306 if(mRestrictVertexMovement)
307 {
308 if (length > 0.5*mpMutableVertexMesh->GetCellRearrangementThreshold())
309 {
310 rDisplacement *= 0.5*mpMutableVertexMesh->GetCellRearrangementThreshold()/length;
311
312 std::ostringstream message;
313 message << "Vertices are moving more than half the CellRearrangementThreshold. This could cause elements to become inverted ";
314 message << "so the motion has been restricted. Use a smaller timestep to avoid these warnings.";
315
316 double suggested_step = 0.95*dt*((0.5*mpMutableVertexMesh->GetCellRearrangementThreshold())/length);
317
318 // The first time we see this behaviour, throw a StepSizeException, but not more than once
319 if(mThrowStepSizeException)
320 {
321 mThrowStepSizeException = false;
322 throw StepSizeException(suggested_step, message.str(), false);
323 }
324 }
325 }
326}
327
328template<unsigned DIM>
330{
331 return GetElementCorrespondingToCell(pCell)->IsDeleted();
332}
333
334template<unsigned DIM>
335void VertexBasedCellPopulation<DIM>::Update(bool hasHadBirthsOrDeaths)
336{
337 VertexElementMap element_map(mpMutableVertexMesh->GetNumAllElements());
338 mpMutableVertexMesh->ReMesh(element_map);
339
340 if (!element_map.IsIdentityMap())
341 {
342 // Fix up the mappings between CellPtrs and VertexElements
344 std::map<Cell*, unsigned> old_map = this->mCellLocationMap;
345
346 this->mCellLocationMap.clear();
347 this->mLocationCellMap.clear();
348
349 for (std::list<CellPtr>::iterator cell_iter = this->mCells.begin();
350 cell_iter != this->mCells.end();
351 ++cell_iter)
352 {
353 // The cell vector should only ever contain living cells
354 unsigned old_elem_index = old_map[(*cell_iter).get()];
355 assert(!element_map.IsDeleted(old_elem_index));
356
357 unsigned new_elem_index = element_map.GetNewIndex(old_elem_index);
358 this->SetCellUsingLocationIndex(new_elem_index, *cell_iter);
359 }
360
361 // Check that each VertexElement has only one CellPtr associated with it in the updated cell population
362 Validate();
363 }
364
365 if (this->GetNumAllCells()>0)
366 {
367 //First cell is representative of other cells
368 bool EdgeModelOrNot = (*this->mCells.begin())->GetSrnModel()->HasEdgeModel();
369
370 if (EdgeModelOrNot)
371 {
372 // Note that SRN initialisation in daughter cell is handled through Cell::Divide() method
373 mPopulationSrn.UpdateSrnAfterBirthOrDeath(element_map);
374 }
375 }
376
377 element_map.ResetToIdentity();
378}
379
380template<unsigned DIM>
382{
383 // Check each element has only one cell attached
384 std::vector<unsigned> validated_element = std::vector<unsigned>(this->GetNumElements(), 0);
385 for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = this->Begin();
386 cell_iter != this->End();
387 ++cell_iter)
388 {
389 unsigned elem_index = this->GetLocationIndexUsingCell(*cell_iter);
390 validated_element[elem_index]++;
391 }
392
393 for (unsigned i=0; i<validated_element.size(); i++)
394 {
395 if (validated_element[i] == 0)
396 {
397 EXCEPTION("At time " << SimulationTime::Instance()->GetTime() <<", Element " << i << " does not appear to have a cell associated with it");
398 }
399
400 if (validated_element[i] > 1)
401 {
402 // This should never be reached as you can only set one cell per element index
403 EXCEPTION("At time " << SimulationTime::Instance()->GetTime() <<", Element " << i << " appears to have " << validated_element[i] << " cells associated with it");
404 }
405 }
406}
407
408template<unsigned DIM>
410{
411 pPopulationWriter->Visit(this);
412}
413
414template<unsigned DIM>
416{
417 pPopulationCountWriter->Visit(this);
418}
419
420template<unsigned DIM>
422{
423 pPopulationEventWriter->Visit(this);
424}
425
426template<unsigned DIM>
427void VertexBasedCellPopulation<DIM>::AcceptCellWriter(boost::shared_ptr<AbstractCellWriter<DIM, DIM> > pCellWriter, CellPtr pCell)
428{
429 pCellWriter->VisitCell(pCell, this);
430}
431
432template<unsigned DIM>
434{
435 // Get the vertex element index corresponding to this cell
436 unsigned elem_index = this->GetLocationIndexUsingCell(pCell);
437
438 // Get the element rosette rank from the vertex mesh
439 unsigned rosette_rank = mpMutableVertexMesh->GetRosetteRankOfElement(elem_index);
440
441 return rosette_rank;
442}
443
444template<unsigned DIM>
446{
447 // Get the vertex element index corresponding to this cell
448 unsigned elem_index = this->GetLocationIndexUsingCell(pCell);
449
450 // Get the cell's volume from the vertex mesh
451 double cell_volume = mpMutableVertexMesh->GetVolumeOfElement(elem_index);
452
453 return cell_volume;
454}
455
456template<unsigned DIM>
458{
459 unsigned num_cells = this->GetNumAllCells();
460 if (num_cells>0)
461 {
462 auto cells = this->rGetCells();
463 boost::shared_ptr<CellEdgeData> p_cell_edge_data = (*cells.begin())->GetCellEdgeData();
464 //If edge SRNs are specified, then write vtk results into a mesh where quantities
465 //associated with each edge are taken into account. We assume that the first cell is
466 //representative of all cells.
467 //Edge VTKs are also written if cells contain CellEdgeData
468
469 //If cells contain edge data
470 if (p_cell_edge_data->GetNumItems() != 0&&mWriteEdgeVtkResults)
471 {
472 this->WriteCellEdgeVtkResultsToFile(rDirectory);
473 return;
474 }
475
476 //If cells don't contain edge data, output CellData only
477 if (p_cell_edge_data->GetNumItems() == 0&&mWriteCellVtkResults)
478 {
479 this->WriteCellVtkResultsToFile(rDirectory);
480 }
481 }
482}
483
484template<unsigned DIM>
486{
487#ifdef CHASTE_VTK
488
489 // Create mesh writer for VTK output
490 VertexMeshWriter<DIM, DIM> mesh_writer(rDirectory, "results", false);
491
492 // We avoid writing out CellData if the population is empty (i.e. no cells).
493 unsigned num_cells = this->GetNumAllCells();
494
495 if (num_cells > 0)
496 {
497 // Iterate over any cell writers that are present
498 for (auto cell_writer_iter = this->mCellWriters.begin();
499 cell_writer_iter != this->mCellWriters.end();
500 ++cell_writer_iter)
501 {
502 // Create vector to store VTK cell data
503 std::vector<double> vtk_cell_data(num_cells);
504
505 // Iterate over vertex elements ///\todo #2512 - replace with loop over cells
506 for (auto elem_iter = mpMutableVertexMesh->GetElementIteratorBegin();
507 elem_iter != mpMutableVertexMesh->GetElementIteratorEnd();
508 ++elem_iter)
509 {
510 // Get index of this element in the vertex mesh
511 unsigned elem_index = elem_iter->GetIndex();
512
513 // Get the cell corresponding to this element
514 CellPtr p_cell = this->GetCellUsingLocationIndex(elem_index);
515 assert(p_cell);
516
517 // Populate the vector of VTK cell data
518 vtk_cell_data[elem_index] = (*cell_writer_iter)->GetCellDataForVtkOutput(p_cell, this);
519 }
520
521 mesh_writer.AddCellData((*cell_writer_iter)->GetVtkCellDataName(), vtk_cell_data);
522 }
523
524 // When outputting any CellData, we assume that the first cell is representative of all cells
525 unsigned num_cell_data_items = this->Begin()->GetCellData()->GetNumItems();
526 std::vector<std::string> cell_data_names = this->Begin()->GetCellData()->GetKeys();
527
528 std::vector<std::vector<double> > cell_data;
529 for (unsigned var=0; var<num_cell_data_items; var++)
530 {
531 std::vector<double> cell_data_var(num_cells);
532 cell_data.push_back(cell_data_var);
533 }
534
535 // Loop over vertex elements ///\todo #2512 - replace with loop over cells
536 for (auto elem_iter = mpMutableVertexMesh->GetElementIteratorBegin();
537 elem_iter != mpMutableVertexMesh->GetElementIteratorEnd();
538 ++elem_iter)
539 {
540 // Get index of this element in the vertex mesh
541 unsigned elem_index = elem_iter->GetIndex();
542
543 // Get the cell corresponding to this element
544 CellPtr p_cell = this->GetCellUsingLocationIndex(elem_index);
545 assert(p_cell);
546
547 for (unsigned var=0; var<num_cell_data_items; var++)
548 {
549 cell_data[var][elem_index] = p_cell->GetCellData()->GetItem(cell_data_names[var]);
550 }
551 }
552 for (unsigned var=0; var<num_cell_data_items; var++)
553 {
554 mesh_writer.AddCellData(cell_data_names[var], cell_data[var]);
555 }
556 }
557
558 unsigned num_timesteps = SimulationTime::Instance()->GetTimeStepsElapsed();
559 std::stringstream time;
560 time << num_timesteps;
561
562 mesh_writer.WriteVtkUsingMesh(*mpMutableVertexMesh, time.str());
563
564 *(this->mpVtkMetaFile) << " <DataSet timestep=\"";
565 *(this->mpVtkMetaFile) << num_timesteps;
566 *(this->mpVtkMetaFile) << "\" group=\"\" part=\"0\" file=\"results_";
567 *(this->mpVtkMetaFile) << num_timesteps;
568 *(this->mpVtkMetaFile) << ".vtu\"/>\n";
569
570#endif //CHASTE_VTK
571}
572
573template<unsigned DIM>
575{
576#ifdef CHASTE_VTK
577 //Writes cell only data
578 {
579 // Create mesh writer for VTK output
580 VertexMeshWriter<DIM, DIM> mesh_writer(rDirectory, "cell_results", false);
581
582 // Iterate over any cell writers that are present
583 unsigned num_cells = this->GetNumAllCells();
584 for (typename std::vector<boost::shared_ptr<AbstractCellWriter<DIM, DIM> > >::iterator cell_writer_iter = this->mCellWriters.begin();
585 cell_writer_iter != this->mCellWriters.end();
586 ++cell_writer_iter)
587 {
588 // Create vector to store VTK cell data
589 std::vector<double> vtk_cell_data(num_cells);
590
591 // Iterate over vertex elements ///\todo #2512 - replace with loop over cells
592 for (typename VertexMesh<DIM,DIM>::VertexElementIterator elem_iter = mpMutableVertexMesh->GetElementIteratorBegin();
593 elem_iter != mpMutableVertexMesh->GetElementIteratorEnd();
594 ++elem_iter)
595 {
596 // Get index of this element in the vertex mesh
597 unsigned elem_index = elem_iter->GetIndex();
598
599 // Get the cell corresponding to this element
600 CellPtr p_cell = this->GetCellUsingLocationIndex(elem_index);
601 assert(p_cell);
602
603 // Populate the vector of VTK cell data
604 vtk_cell_data[elem_index] = (*cell_writer_iter)->GetCellDataForVtkOutput(p_cell, this);
605 }
606
607 mesh_writer.AddCellData((*cell_writer_iter)->GetVtkCellDataName(), vtk_cell_data);
608 }
609
610 // When outputting any CellData, we assume that the first cell is representative of all cells
611 unsigned num_cell_data_items = this->Begin()->GetCellData()->GetNumItems();
612 std::vector<std::string> cell_data_names = this->Begin()->GetCellData()->GetKeys();
613
614 std::vector<std::vector<double> > cell_data;
615 for (unsigned var=0; var<num_cell_data_items; var++)
616 {
617 std::vector<double> cell_data_var(num_cells);
618 cell_data.push_back(cell_data_var);
619 }
620
621 // Loop over vertex elements ///\todo #2512 - replace with loop over cells
622 for (typename VertexMesh<DIM,DIM>::VertexElementIterator elem_iter = mpMutableVertexMesh->GetElementIteratorBegin();
623 elem_iter != mpMutableVertexMesh->GetElementIteratorEnd();
624 ++elem_iter)
625 {
626 // Get index of this element in the vertex mesh
627 unsigned elem_index = elem_iter->GetIndex();
628
629 // Get the cell corresponding to this element
630 CellPtr p_cell = this->GetCellUsingLocationIndex(elem_index);
631 assert(p_cell);
632
633 for (unsigned var=0; var<num_cell_data_items; var++)
634 {
635 cell_data[var][elem_index] = p_cell->GetCellData()->GetItem(cell_data_names[var]);
636 }
637 }
638 for (unsigned var=0; var<num_cell_data_items; var++)
639 {
640 mesh_writer.AddCellData(cell_data_names[var], cell_data[var]);
641 }
642
643 unsigned num_timesteps = SimulationTime::Instance()->GetTimeStepsElapsed();
644 std::stringstream time;
645 time << num_timesteps;
646
647 mesh_writer.WriteVtkUsingMesh(*mpMutableVertexMesh, time.str());
648
649 *(this->mpVtkMetaFile) << " <DataSet timestep=\"";
650 *(this->mpVtkMetaFile) << num_timesteps;
651 *(this->mpVtkMetaFile) << "\" group=\"\" part=\"0\" file=\"cell_results_";
652 *(this->mpVtkMetaFile) << num_timesteps;
653 *(this->mpVtkMetaFile) << ".vtu\"/>\n";
654 }
655 // Create mesh writer for VTK output
656 TrapezoidEdgeVertexMeshWriter<DIM, DIM> mesh_writer(rDirectory, "results", false);
657 unsigned num_edges = 0;
658
659 // Here elements are synonymous with cells
660 const unsigned num_cells = this->GetNumElements();
661 //Similarly as in TrapEdgeVerteMeshWriter, but instead of nodes
662 //we fill edge arrays
663 //The first value MUST be zero
664 std::vector<unsigned> cell_offset_dist(num_cells);
665 //The order of stored data is illustrated below:
666 // [_____|_][____|_]
667 // ^^^^ ^
668 // edge cell interior
669 for (unsigned i=1; i<num_cells; ++i)
670 {
671 cell_offset_dist[i] = cell_offset_dist[i-1]+this->GetElement(i-1)->GetNumEdges()+1;
672 num_edges += this->GetElement(i)->GetNumEdges();
673 }
674 // Total number of edges
675 num_edges += this->GetElement(0)->GetNumEdges();
676
677 // Iterate over any cell writers that are present
678 // The data that is written below is associated with the entire cell
679 // E.g. the cell age. Thus, edges also share the same characteristic
680 for (auto cell_writer : this->mCellWriters)
681 {
682 // Create vector to store VTK cell data
683 std::vector<double> vtk_cell_data(num_edges+num_cells);
684
685 // Iterate over vertex elements ///\todo #2512 - replace with loop over cells
686 for (auto elem_iter = mpMutableVertexMesh->GetElementIteratorBegin();
687 elem_iter != mpMutableVertexMesh->GetElementIteratorEnd();
688 ++elem_iter)
689 {
690 // Get index of this element in the vertex mesh
691 unsigned elem_index = elem_iter->GetIndex();
692
693 // Get the cell corresponding to this element
694 CellPtr p_cell = this->GetCellUsingLocationIndex(elem_index);
695 assert(p_cell);
696
697 // Write the same property in all the edges
698 unsigned elem_num_edges = elem_iter->GetNumEdges();
699
700 // Edge data
701 for (unsigned e = 0; e < elem_num_edges; ++e)
702 {
703 // Populate the vector of VTK cell data
704 vtk_cell_data[cell_offset_dist[elem_index]+e] = cell_writer->GetCellDataForVtkOutput(p_cell, this);
705 }
706 // Internal cell data
707 vtk_cell_data[cell_offset_dist[elem_index]+elem_num_edges] = cell_writer->GetCellDataForVtkOutput(p_cell, this);
708 }
709 mesh_writer.AddCellData(cell_writer->GetVtkCellDataName(), vtk_cell_data);
710 }
711 // When outputting CellData and CellEdgeData, we assume that the first cell
712 // and its edges are representative of the population
713
714 // Get cell/edge data names
715 const unsigned num_cell_data_items = this->Begin()->GetCellData()->GetNumItems();
716 std::vector<std::string> cell_data_names = this->Begin()->GetCellData()->GetKeys();
717
718 const unsigned num_edge_data_items = this->Begin()->GetCellEdgeData()->GetNumItems();
719 std::vector<std::string> edge_data_names = this->Begin()->GetCellEdgeData()->GetKeys();
720
721 // Total number of data items. Each data item (edge+interior) has its own value
722 const unsigned num_data_items = num_edge_data_items + num_cell_data_items;
723 std::vector<std::string> data_names(num_data_items);
724 for (unsigned var=0; var<num_edge_data_items; ++var)
725 {
726 data_names[var] = edge_data_names[var];
727 }
728 for (unsigned var=num_edge_data_items; var<num_data_items; ++var)
729 {
730 data_names[var] = cell_data_names[var-num_edge_data_items];
731 }
732 std::vector<std::vector<double> > data_values(num_data_items,
733 std::vector<double>(num_edges + num_cells));
734
735 // Writing CellEdgeData values to the edges of the cells
736 // Loop over vertex elements ///\todo #2512 - replace with loop over cells
737 for (auto elem_iter = mpMutableVertexMesh->GetElementIteratorBegin();
738 elem_iter != mpMutableVertexMesh->GetElementIteratorEnd();
739 ++elem_iter)
740 {
741 const unsigned elem_index = elem_iter->GetIndex();
742 CellPtr p_cell = this->GetCellUsingLocationIndex(elem_index);
743 assert(p_cell);
744
745 // Edge data for interior data is set to zero.
746 unsigned elem_num_edges = elem_iter->GetNumEdges();
747 for (unsigned var = 0; var < num_edge_data_items; ++var)
748 {
749 // Write the same property in all the edges
750 for (unsigned e = 0; e < elem_num_edges; ++e)
751 {
752 data_values[var][cell_offset_dist[elem_index]+e] = p_cell->GetCellEdgeData()->GetItem(data_names[var])[e];
753 }
754 // Cell interior is set to zero
755 data_values[var][cell_offset_dist[elem_index]+elem_num_edges] = 0.0;
756 }
757 }
758
759 // Writing CellData values to the interior of the cells
760 // Loop over vertex elements ///\todo #2512 - replace with loop over cells
761 for (auto elem_iter = mpMutableVertexMesh->GetElementIteratorBegin();
762 elem_iter != mpMutableVertexMesh->GetElementIteratorEnd();
763 ++elem_iter)
764 {
765 // Get index of this element in the vertex mesh
766 unsigned elem_index = elem_iter->GetIndex();
767 unsigned elem_num_edges = elem_iter->GetNumEdges();
768
769 // Get the cell corresponding to this element
770 CellPtr p_cell = this->GetCellUsingLocationIndex(elem_index);
771 assert(p_cell);
772
773 for (unsigned var = num_edge_data_items; var<num_data_items; ++var)
774 {
775 // Data in the edges are set to 0
776 for (unsigned e = 0; e < elem_num_edges; ++e)
777 {
778 data_values[var][cell_offset_dist[elem_index]+e] = 0.0;
779 }
780 // Filling cell interior data
781 data_values[var][cell_offset_dist[elem_index]+elem_num_edges] = p_cell->GetCellData()->GetItem(data_names[var]);
782 }
783 }
784
785 for (unsigned var=0; var<num_data_items; var++)
786 {
787 mesh_writer.AddCellData(data_names[var], data_values[var]);
788 }
789
790 unsigned num_timesteps = SimulationTime::Instance()->GetTimeStepsElapsed();
791 std::stringstream time;
792 time << num_timesteps;
793
794 mesh_writer.WriteVtkUsingMesh(*mpMutableVertexMesh, time.str());
795
796 *(this->mpVtkMetaFile) << " <DataSet timestep=\"";
797 *(this->mpVtkMetaFile) << num_timesteps;
798 *(this->mpVtkMetaFile) << "\" group=\"\" part=\"0\" file=\"results_";
799 *(this->mpVtkMetaFile) << num_timesteps;
800 *(this->mpVtkMetaFile) << ".vtu\"/>\n";
801#endif //CHASTE_VTK
802}
803
804template<unsigned DIM>
806{
807 if (this->mOutputResultsForChasteVisualizer)
808 {
809 if (!this-> template HasWriter<CellPopulationElementWriter>())
810 {
811 this-> template AddPopulationWriter<CellPopulationElementWriter>();
812 }
813 }
814
815 if (mOutputCellRearrangementLocations)
816 {
817 if (!this-> template HasWriter<VertexT1SwapLocationsWriter>())
818 {
819 this-> template AddPopulationWriter<VertexT1SwapLocationsWriter>();
820 }
821 if (!this-> template HasWriter<VertexT2SwapLocationsWriter>())
822 {
823 this-> template AddPopulationWriter<VertexT2SwapLocationsWriter>();
824 }
825 if (!this-> template HasWriter<VertexT3SwapLocationsWriter>())
826 {
827 this-> template AddPopulationWriter<VertexT3SwapLocationsWriter>();
828 }
829 if (!this-> template HasWriter<VertexIntersectionSwapLocationsWriter>())
830 {
831 this-> template AddPopulationWriter<VertexIntersectionSwapLocationsWriter>();
832 }
833 }
834
836}
837
838template<unsigned DIM>
840{
841 return mOutputCellRearrangementLocations;
842}
843
844template<unsigned DIM>
846{
847 mOutputCellRearrangementLocations = outputCellRearrangementLocations;
848}
849
850template<unsigned DIM>
852{
853 *rParamsFile << "\t\t<CellRearrangementThreshold>" << mpMutableVertexMesh->GetCellRearrangementThreshold() << "</CellRearrangementThreshold>\n";
854 *rParamsFile << "\t\t<T2Threshold>" << mpMutableVertexMesh->GetT2Threshold() << "</T2Threshold>\n";
855 *rParamsFile << "\t\t<CellRearrangementRatio>" << mpMutableVertexMesh->GetCellRearrangementRatio() << "</CellRearrangementRatio>\n";
856 *rParamsFile << "\t\t<OutputCellRearrangementLocations>" << mOutputCellRearrangementLocations << "</OutputCellRearrangementLocations>\n";
857
858 // Add the division rule parameters
859 *rParamsFile << "\t\t<VertexBasedDivisionRule>\n";
860 mpVertexBasedDivisionRule->OutputCellVertexBasedDivisionRuleInfo(rParamsFile);
861 *rParamsFile << "\t\t</VertexBasedDivisionRule>\n";
862
863 // Call method on direct parent class
865}
866
867template<unsigned DIM>
868double VertexBasedCellPopulation<DIM>::GetWidth(const unsigned& rDimension)
869{
870 // Call GetWidth() on the mesh
871 double width = this->mrMesh.GetWidth(rDimension);
872
873 return width;
874}
875
876template<unsigned DIM>
878{
879 return mpMutableVertexMesh->GetNeighbouringNodeIndices(index);
880}
881
882template<unsigned DIM>
883boost::shared_ptr<AbstractVertexBasedDivisionRule<DIM> > VertexBasedCellPopulation<DIM>::GetVertexBasedDivisionRule()
884{
885 return mpVertexBasedDivisionRule;
886}
887
888template<unsigned DIM>
890{
891 mpVertexBasedDivisionRule = pVertexBasedDivisionRule;
892}
893
894template<unsigned DIM>
896{
897 // This method only works in 2D sequential
898 if (DIM != 2)
899 {
900 EXCEPTION("This function is only valid in 2D"); // LCOV_EXCL_LINE
901 }
902 assert(PetscTools::IsSequential());
903
904 unsigned num_vertex_nodes = mpMutableVertexMesh->GetNumNodes();
905 unsigned num_vertex_elements = mpMutableVertexMesh->GetNumElements();
906
907 std::string mesh_file_name = "mesh";
908
909 // Get a unique temporary foldername
910 std::stringstream pid;
911 pid << getpid();
912 OutputFileHandler output_file_handler("2D_temporary_tetrahedral_mesh_" + pid.str());
913 std::string output_dir = output_file_handler.GetOutputDirectoryFullPath();
914
915 // Compute the number of nodes in the TetrahedralMesh
916 unsigned num_tetrahedral_nodes = num_vertex_nodes + num_vertex_elements;
917
918 // Write node file
919 out_stream p_node_file = output_file_handler.OpenOutputFile(mesh_file_name+".node");
920 (*p_node_file) << std::scientific;
921 (*p_node_file) << std::setprecision(20);
922 (*p_node_file) << num_tetrahedral_nodes << "\t2\t0\t1" << std::endl;
923
924 // Begin by writing each node in the VertexMesh
925 for (unsigned node_index=0; node_index<num_vertex_nodes; node_index++)
926 {
927 Node<DIM>* p_node = mpMutableVertexMesh->GetNode(node_index);
928
930 unsigned index = p_node->GetIndex();
931 const c_vector<double, DIM>& r_location = p_node->rGetLocation();
932 unsigned is_boundary_node = p_node->IsBoundaryNode() ? 1 : 0;
933
934 (*p_node_file) << index << "\t" << r_location[0] << "\t" << r_location[1] << "\t" << is_boundary_node << std::endl;
935 }
936
937 // Now write an additional node at each VertexElement's centroid
938 unsigned num_tetrahedral_elements = 0;
939 for (unsigned vertex_elem_index=0; vertex_elem_index<num_vertex_elements; vertex_elem_index++)
940 {
941 unsigned index = num_vertex_nodes + vertex_elem_index;
942
943 c_vector<double, DIM> location = mpMutableVertexMesh->GetCentroidOfElement(vertex_elem_index);
944
945 // Any node located at a VertexElement's centroid will not be a boundary node
946 unsigned is_boundary_node = 0;
947 (*p_node_file) << index << "\t" << location[0] << "\t" << location[1] << "\t" << is_boundary_node << std::endl;
948
949 // Also keep track of how many tetrahedral elements there will be
950 num_tetrahedral_elements += mpMutableVertexMesh->GetElement(vertex_elem_index)->GetNumNodes();
951 }
952 p_node_file->close();
953
954 // Write element file
955 out_stream p_elem_file = output_file_handler.OpenOutputFile(mesh_file_name+".ele");
956 (*p_elem_file) << std::scientific;
957 (*p_elem_file) << num_tetrahedral_elements << "\t3\t0" << std::endl;
958
959 std::set<std::pair<unsigned, unsigned> > tetrahedral_edges;
960
961 unsigned tetrahedral_elem_index = 0;
962 for (unsigned vertex_elem_index=0; vertex_elem_index<num_vertex_elements; vertex_elem_index++)
963 {
964 VertexElement<DIM, DIM>* p_vertex_element = mpMutableVertexMesh->GetElement(vertex_elem_index);
965
966 // Iterate over nodes owned by this VertexElement
967 unsigned num_nodes_in_vertex_element = p_vertex_element->GetNumNodes();
968 for (unsigned local_index=0; local_index<num_nodes_in_vertex_element; local_index++)
969 {
970 unsigned node_0_index = p_vertex_element->GetNodeGlobalIndex(local_index);
971 unsigned node_1_index = p_vertex_element->GetNodeGlobalIndex((local_index+1)%num_nodes_in_vertex_element);
972 unsigned node_2_index = num_vertex_nodes + vertex_elem_index;
973
974 (*p_elem_file) << tetrahedral_elem_index++ << "\t" << node_0_index << "\t" << node_1_index << "\t" << node_2_index << std::endl;
975
976 // Add edges to the set if they are not already present
977 std::pair<unsigned, unsigned> edge_0 = this->CreateOrderedPair(node_0_index, node_1_index);
978 std::pair<unsigned, unsigned> edge_1 = this->CreateOrderedPair(node_1_index, node_2_index);
979 std::pair<unsigned, unsigned> edge_2 = this->CreateOrderedPair(node_2_index, node_0_index);
980
981 tetrahedral_edges.insert(edge_0);
982 tetrahedral_edges.insert(edge_1);
983 tetrahedral_edges.insert(edge_2);
984 }
985 }
986 p_elem_file->close();
987
988 // Write edge file
989 out_stream p_edge_file = output_file_handler.OpenOutputFile(mesh_file_name+".edge");
990 (*p_edge_file) << std::scientific;
991 (*p_edge_file) << tetrahedral_edges.size() << "\t1" << std::endl;
992
993 unsigned edge_index = 0;
994 for (std::set<std::pair<unsigned, unsigned> >::iterator edge_iter = tetrahedral_edges.begin();
995 edge_iter != tetrahedral_edges.end();
996 ++edge_iter)
997 {
998 std::pair<unsigned, unsigned> this_edge = *edge_iter;
999
1000 // To be a boundary edge both nodes need to be boundary nodes.
1001 bool is_boundary_edge = false;
1002 if (this_edge.first < mpMutableVertexMesh->GetNumNodes() &&
1003 this_edge.second < mpMutableVertexMesh->GetNumNodes())
1004 {
1005 is_boundary_edge = (mpMutableVertexMesh->GetNode(this_edge.first)->IsBoundaryNode() &&
1006 mpMutableVertexMesh->GetNode(this_edge.second)->IsBoundaryNode() );
1007 }
1008 unsigned is_boundary_edge_unsigned = is_boundary_edge ? 1 : 0;
1009
1010 (*p_edge_file) << edge_index++ << "\t" << this_edge.first << "\t" << this_edge.second << "\t" << is_boundary_edge_unsigned << std::endl;
1011 }
1012 p_edge_file->close();
1013
1014 // Having written the mesh to file, now construct it using TrianglesMeshReader
1016
1017 // Nested scope so reader is destroyed before we remove the temporary files
1018 {
1019 TrianglesMeshReader<DIM, DIM> mesh_reader(output_dir + mesh_file_name);
1020 p_mesh->ConstructFromMeshReader(mesh_reader);
1021 }
1022
1023 // Delete the temporary files
1024 output_file_handler.FindFile("").Remove();
1025
1026 // The original files have been deleted, it is better if the mesh object forgets about them
1028
1029 return p_mesh;
1030}
1031
1032template<unsigned DIM>
1033std::vector< c_vector< double, DIM > > VertexBasedCellPopulation<DIM>::GetLocationsOfT2Swaps()
1034{
1035 return mLocationsOfT2Swaps;
1036}
1037
1038template<unsigned DIM>
1040{
1041 return mCellIdsOfT2Swaps;
1042}
1043
1044template<unsigned DIM>
1045void VertexBasedCellPopulation<DIM>::AddLocationOfT2Swap(c_vector< double, DIM> locationOfT2Swap)
1046{
1047 mLocationsOfT2Swaps.push_back(locationOfT2Swap);
1048}
1049
1050template<unsigned DIM>
1052{
1053 mCellIdsOfT2Swaps.push_back(idOfT2Swap);
1054}
1055
1056template<unsigned DIM>
1058{
1059 mCellIdsOfT2Swaps.clear();
1060 mLocationsOfT2Swaps.clear();
1061}
1062
1063template<unsigned DIM>
1065{
1066 bool non_apoptotic_cell_present = true;
1067
1068 if (pdeNodeIndex < this->GetNumNodes())
1069 {
1070 std::set<unsigned> containing_element_indices = this->GetNode(pdeNodeIndex)->rGetContainingElementIndices();
1071
1072 for (std::set<unsigned>::iterator iter = containing_element_indices.begin();
1073 iter != containing_element_indices.end();
1074 iter++)
1075 {
1076 if (this->GetCellUsingLocationIndex(*iter)->template HasCellProperty<ApoptoticCellProperty>() )
1077 {
1078 non_apoptotic_cell_present = false;
1079 break;
1080 }
1081 }
1082 }
1083 else
1084 {
1085 /*
1086 * This node of the tetrahedral finite element mesh is in the centre of the element of the
1087 * vertex-based cell population, so we can use an offset to compute which cell to interrogate.
1088 */
1089 non_apoptotic_cell_present = !(this->GetCellUsingLocationIndex(pdeNodeIndex - this->GetNumNodes())->template HasCellProperty<ApoptoticCellProperty>());
1090 }
1091
1092 return non_apoptotic_cell_present;
1093}
1094
1095template<unsigned DIM>
1097 unsigned pdeNodeIndex,
1098 std::string& rVariableName,
1099 bool dirichletBoundaryConditionApplies,
1100 double dirichletBoundaryValue)
1101{
1102 unsigned num_nodes = this->GetNumNodes();
1103 double value = 0.0;
1104
1105 // Cells correspond to nodes in the centre of the vertex element; nodes on vertices have averaged values from containing cells
1106
1107 if (pdeNodeIndex >= num_nodes)
1108 {
1109 // Offset to relate elements in vertex mesh to nodes in tetrahedral mesh
1110 assert(pdeNodeIndex-num_nodes < num_nodes);
1111
1112 CellPtr p_cell = this->GetCellUsingLocationIndex(pdeNodeIndex - num_nodes);
1113 value = p_cell->GetCellData()->GetItem(rVariableName);
1114 }
1115 else
1116 {
1118 if (dirichletBoundaryConditionApplies)
1119 {
1120 // We need to impose the Dirichlet boundaries again here as not represented in cell data
1121 value = dirichletBoundaryValue;
1122 }
1123 else
1124 {
1125 assert(pdeNodeIndex < num_nodes);
1126 Node<DIM>* p_node = this->GetNode(pdeNodeIndex);
1127
1128 // Average over data from containing elements (cells)
1129 std::set<unsigned> containing_elements = p_node->rGetContainingElementIndices();
1130 for (std::set<unsigned>::iterator index_iter = containing_elements.begin();
1131 index_iter != containing_elements.end();
1132 ++index_iter)
1133 {
1134 assert(*index_iter < num_nodes);
1135 CellPtr p_cell = this->GetCellUsingLocationIndex(*index_iter);
1136 value += p_cell->GetCellData()->GetItem(rVariableName);
1137 }
1138 value /= containing_elements.size();
1139 }
1140 }
1141
1142 return value;
1143}
1144
1145template<unsigned DIM>
1147{
1148 return 0.002;
1149}
1150
1151template<unsigned DIM>
1153{
1154 if (bool(dynamic_cast<Cylindrical2dVertexMesh*>(&(this->mrMesh))))
1155 {
1156 *pVizSetupFile << "MeshWidth\t" << this->GetWidth(0) << "\n";
1157 }
1158}
1159
1160template<unsigned DIM>
1162{
1163 MAKE_PTR_ARGS(T2SwapCellKiller<DIM>, p_t2_swap_cell_killer, (this));
1164 pSimulation->AddCellKiller(p_t2_swap_cell_killer);
1165}
1166
1167
1168template<unsigned DIM>
1170{
1171 return mRestrictVertexMovement;
1172}
1173
1174template<unsigned DIM>
1176{
1177 mRestrictVertexMovement = restrictMovement;
1178}
1179
1180template<unsigned DIM>
1185
1186template<unsigned DIM>
1191
1192template<unsigned DIM>
1194{
1195 mWriteCellVtkResults = new_val;
1196}
1197
1198template<unsigned DIM>
1200{
1201 mWriteEdgeVtkResults = new_val;
1202}
1203
1204// Explicit instantiation
1205template class VertexBasedCellPopulation<1>;
1206template class VertexBasedCellPopulation<2>;
1207template class VertexBasedCellPopulation<3>;
1208
1209// Serialization for Boost >= 1.36
#define EXCEPTION(message)
#define EXPORT_TEMPLATE_CLASS_SAME_DIMS(CLASS)
#define MAKE_PTR_ARGS(TYPE, NAME, ARGS)
void AddCellKiller(boost::shared_ptr< AbstractCellKiller< SPACE_DIM > > pCellKiller)
virtual void AddCellUsingLocationIndex(unsigned index, CellPtr pCell)
virtual void OpenWritersFiles(OutputFileHandler &rOutputFileHandler)
AbstractMesh< ELEMENT_DIM, SPACE_DIM > & mrMesh
unsigned GetNumNodes() const
unsigned GetNodeGlobalIndex(unsigned localIndex) const
void SetMeshHasChangedSinceLoading()
virtual void OutputCellPopulationParameters(out_stream &rParamsFile)
void Remove() const
void SetMeshOperationTracking(const bool track)
Definition Node.hpp:59
std::set< unsigned > & rGetContainingElementIndices()
Definition Node.cpp:300
const c_vector< double, SPACE_DIM > & rGetLocation() const
Definition Node.cpp:139
bool IsBoundaryNode() const
Definition Node.cpp:164
unsigned GetIndex() const
Definition Node.cpp:158
std::string GetOutputDirectoryFullPath() const
FileFinder FindFile(std::string leafName) const
out_stream OpenOutputFile(const std::string &rFileName, std::ios_base::openmode mode=std::ios::out|std::ios::trunc) const
static bool IsSequential()
static SimulationTime * Instance()
unsigned GetTimeStepsElapsed() const
void ConstructFromMeshReader(AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &rMeshReader)
void WriteVtkUsingMesh(VertexMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, const std::string &stamp="")
void AddCellData(std::string dataName, std::vector< double > dataPayload)
double GetDampingConstant(unsigned nodeIndex)
virtual void AcceptCellWriter(boost::shared_ptr< AbstractCellWriter< DIM, DIM > > pCellWriter, CellPtr pCell)
virtual void AcceptPopulationWriter(boost::shared_ptr< AbstractCellPopulationWriter< DIM, DIM > > pPopulationWriter)
void OutputCellPopulationParameters(out_stream &rParamsFile)
std::set< unsigned > GetNeighbouringNodeIndices(unsigned index)
c_vector< double, DIM > GetLocationOfCellCentre(CellPtr pCell)
virtual void WriteCellVtkResultsToFile(const std::string &rDirectory)
void SetWriteCellVtkResults(const bool new_val)
boost::shared_ptr< AbstractVertexBasedDivisionRule< DIM > > GetVertexBasedDivisionRule()
virtual void WriteVtkResultsToFile(const std::string &rDirectory)
virtual void OpenWritersFiles(OutputFileHandler &rOutputFileHandler)
virtual void AcceptPopulationEventWriter(boost::shared_ptr< AbstractCellPopulationEventWriter< DIM, DIM > > pPopulationEventWriter)
VertexElement< DIM, DIM > * GetElementCorrespondingToCell(CellPtr pCell)
virtual void SimulationSetupHook(AbstractCellBasedSimulation< DIM, DIM > *pSimulation)
bool IsCellAssociatedWithADeletedLocation(CellPtr pCell)
unsigned AddNode(Node< DIM > *pNewNode)
CellPtr AddCell(CellPtr pNewCell, CellPtr pParentCell=CellPtr())
VertexBasedCellPopulation(MutableVertexMesh< DIM, DIM > &rMesh, std::vector< CellPtr > &rCells, bool deleteMesh=false, bool validate=true, const std::vector< unsigned > locationIndices=std::vector< unsigned >())
unsigned GetRosetteRankOfCell(CellPtr pCell)
VertexElement< DIM, DIM > * GetElement(unsigned elementIndex)
virtual bool IsPdeNodeAssociatedWithNonApoptoticCell(unsigned pdeNodeIndex)
void AddLocationOfT2Swap(c_vector< double, DIM > locationOfT2Swap)
void SetWriteEdgeVtkResults(const bool new_val)
virtual TetrahedralMesh< DIM, DIM > * GetTetrahedralMeshForPdeModifier()
VertexBasedPopulationSrn< DIM > & rGetVertexBasedPopulationSrn()
void SetOutputCellRearrangementLocations(bool outputCellRearrangementLocations)
MutableVertexMesh< DIM, DIM > & rGetMesh()
void SetVertexBasedDivisionRule(boost::shared_ptr< AbstractVertexBasedDivisionRule< DIM > > pVertexBasedDivisionRule)
MutableVertexMesh< DIM, DIM > * mpMutableVertexMesh
VertexBasedPopulationSrn< DIM > mPopulationSrn
virtual void WriteCellEdgeVtkResultsToFile(const std::string &rDirectory)
std::set< std::pair< unsigned, unsigned > > GetNeighbouringEdgeIndices(CellPtr pCell, unsigned edgeLocalIndex)
virtual void CheckForStepSizeException(unsigned nodeIndex, c_vector< double, DIM > &rDisplacement, double dt)
std::vector< c_vector< double, DIM > > GetLocationsOfT2Swaps()
void AddCellIdOfT2Swap(unsigned idOfT2Swap)
virtual double GetCellDataItemAtPdeNode(unsigned pdeNodeIndex, std::string &rVariableName, bool dirichletBoundaryConditionApplies=false, double dirichletBoundaryValue=0.0)
virtual void WriteDataToVisualizerSetupFile(out_stream &pVizSetupFile)
std::set< unsigned > GetNeighbouringLocationIndices(CellPtr pCell)
void Update(bool hasHadBirthsOrDeaths=true)
void SetRestrictVertexMovementBoolean(bool restrictVertexMovement)
boost::shared_ptr< AbstractVertexBasedDivisionRule< DIM > > mpVertexBasedDivisionRule
std::vector< unsigned > GetCellIdsOfT2Swaps()
void SetNode(unsigned index, ChastePoint< DIM > &rNewLocation)
virtual void AcceptPopulationCountWriter(boost::shared_ptr< AbstractCellPopulationCountWriter< DIM, DIM > > pPopulationCountWriter)
double GetWidth(const unsigned &rDimension)
Node< DIM > * GetNode(unsigned index)
unsigned GetNewIndex(unsigned oldIndex) const
bool IsDeleted(unsigned index)
void WriteVtkUsingMesh(VertexMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, std::string stamp="")
void AddCellData(std::string dataName, std::vector< double > dataPayload)