Chaste Commit::1fd4e48e3990e67db148bc1bc4cf6991a0049d0c
PottsBasedCellPopulation.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 "PottsBasedCellPopulation.hpp"
37#include "RandomNumberGenerator.hpp"
38#include "AbstractPottsUpdateRule.hpp"
39#include "NodesOnlyMesh.hpp"
40#include "CellPopulationElementWriter.hpp"
41#include "CellIdWriter.hpp"
42
43// Needed to convert mesh in order to write nodes to VTK (visualize as glyphs)
44#include "VtkMeshWriter.hpp"
45
46template<unsigned DIM>
48{
49 // Check each element has only one cell associated with it
50 std::vector<unsigned> validated_element = std::vector<unsigned>(this->GetNumElements(), 0);
51
52 for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = this->Begin();
53 cell_iter != this->End();
54 ++cell_iter)
55 {
56 unsigned elem_index = this->GetLocationIndexUsingCell(*cell_iter);
57 validated_element[elem_index]++;
58 }
59
60 for (unsigned i=0; i<validated_element.size(); i++)
61 {
62 if (validated_element[i] == 0)
63 {
64 EXCEPTION("At time " << SimulationTime::Instance()->GetTime() << ", Element " << i << " does not appear to have a cell associated with it");
65 }
66
67 if (validated_element[i] > 1)
68 {
69 EXCEPTION("At time " << SimulationTime::Instance()->GetTime() << ", Element " << i << " appears to have " << validated_element[i] << " cells associated with it");
70 }
71 }
72}
73
74template<unsigned DIM>
76 std::vector<CellPtr>& rCells,
77 bool deleteMesh,
78 bool validate,
79 const std::vector<unsigned> locationIndices)
80 : AbstractOnLatticeCellPopulation<DIM>(rMesh, rCells, locationIndices, deleteMesh),
81 mpElementTessellation(nullptr),
82 mpMutableMesh(nullptr),
83 mTemperature(0.1),
84 mNumSweepsPerTimestep(1)
85{
86 mpPottsMesh = static_cast<PottsMesh<DIM>* >(&(this->mrMesh));
87 // Check each element has only one cell associated with it
88 if (validate)
89 {
90 Validate();
91 }
92}
93
94template<unsigned DIM>
97 mpElementTessellation(nullptr),
98 mpMutableMesh(nullptr),
99 mTemperature(0.1),
100 mNumSweepsPerTimestep(1)
101{
102 mpPottsMesh = static_cast<PottsMesh<DIM>* >(&(this->mrMesh));
103}
104
105template<unsigned DIM>
107{
108 // This pointer is always null because PottsBasedCellPopulation::CreateElementTessellation
109 // is not implemented. See #1666 in the trac ticket archive for more information.
110 assert(mpElementTessellation == nullptr);
111
112 delete mpMutableMesh;
113
114 if (this->mDeleteMesh)
115 {
116 delete &this->mrMesh;
117 }
118}
119
120template<unsigned DIM>
122{
123 return *mpPottsMesh;
124}
125
126template<unsigned DIM>
128{
129 return *mpPottsMesh;
130}
131
132template<unsigned DIM>
134{
135 std::vector<Node<DIM>*> temp_nodes;
136
137 // Create nodes at the centre of the cells
138 for (typename AbstractCellPopulation<DIM>::Iterator cell_iter = this->Begin();
139 cell_iter != this->End();
140 ++cell_iter)
141 {
142 unsigned index = this->GetLocationIndexUsingCell(*cell_iter);
143 c_vector<double, DIM> location = this->GetLocationOfCellCentre(*cell_iter);
144 temp_nodes.push_back(new Node<DIM>(index, location));
145 }
146
147 return new MutableMesh<DIM, DIM>(temp_nodes);
148}
149
150template<unsigned DIM>
152{
153 return mpPottsMesh->GetElement(elementIndex);
154}
155
156template<unsigned DIM>
158{
159 return mpPottsMesh->GetNumElements();
160}
161
162template<unsigned DIM>
164{
165 return this->mrMesh.GetNode(index);
166}
167
168template<unsigned DIM>
170{
171 return this->mrMesh.GetNumNodes();
172}
173
174template<unsigned DIM>
176{
177 unsigned elem_index = this->GetLocationIndexUsingCell(pCell);
178 return mpPottsMesh->GetNeighbouringElementIndices(elem_index);
179}
180
181template<unsigned DIM>
183{
184 return mpPottsMesh->GetCentroidOfElement(this->GetLocationIndexUsingCell(pCell));
185}
186
187template<unsigned DIM>
189{
190 return mpPottsMesh->GetElement(this->GetLocationIndexUsingCell(pCell));
191}
192
193template<unsigned DIM>
194CellPtr PottsBasedCellPopulation<DIM>::AddCell(CellPtr pNewCell, CellPtr pParentCell)
195{
196 // Get the element associated with this cell
197 PottsElement<DIM>* p_element = GetElementCorrespondingToCell(pParentCell);
198
199 // Divide the element
200 unsigned new_element_index = mpPottsMesh->DivideElement(p_element, true); // new element will be below the existing element
201
202 // Associate the new cell with the element
203 this->mCells.push_back(pNewCell);
204
205 // Update location cell map
206 CellPtr p_created_cell = this->mCells.back();
207 this->SetCellUsingLocationIndex(new_element_index,p_created_cell);
208 return p_created_cell;
209}
210
211template<unsigned DIM>
213{
214 unsigned num_removed = 0;
215
216 for (std::list<CellPtr>::iterator cell_iter = this->mCells.begin();
217 cell_iter != this->mCells.end();
218 )
219 {
220 if ((*cell_iter)->IsDead())
221 {
222 // Get the location index corresponding to this cell
223 unsigned location_index = this->GetLocationIndexUsingCell(*cell_iter);
224
225 // Use this to remove the cell from the population
226 mpPottsMesh->DeleteElement(location_index);
227
228 // Erase cell and update counter
229 cell_iter = this->mCells.erase(cell_iter);
230 num_removed++;
231 }
232 else
233 {
234 ++cell_iter;
235 }
236 }
237 return num_removed;
238}
239
240template<unsigned DIM>
242{
243 /*
244 * This method implements a Monte Carlo method to update the cell population.
245 * We sample randomly from all nodes in the mesh. Once we have selected a target
246 * node we randomly select a neighbour. The Hamiltonian is evaluated in the
247 * current configuration (H_0) and with the target node added to the same
248 * element as the neighbour (H_1). Based on the vale of deltaH = H_1 - H_0,
249 * the switch is either made or not.
250 *
251 * For each time step (i.e. each time this method is called) we sample
252 * mrMesh.GetNumNodes() nodes. This is known as a Monte Carlo Step (MCS).
253 */
254
256 unsigned num_nodes = this->mrMesh.GetNumNodes();
257
258 // Randomly permute mUpdateRuleCollection if specified
259 if (this->mIterateRandomlyOverUpdateRuleCollection)
260 {
261 // Randomly permute mUpdateRuleCollection
262 p_gen->Shuffle(this->mUpdateRuleCollection);
263 }
264
265 for (unsigned i=0; i<num_nodes*mNumSweepsPerTimestep; i++)
266 {
267 unsigned node_index;
268
269 if (this->mUpdateNodesInRandomOrder)
270 {
271 node_index = p_gen->randMod(num_nodes);
272 }
273 else
274 {
275 // Loop over nodes in index order.
276 node_index = i%num_nodes;
277 }
278
279 Node<DIM>* p_node = this->mrMesh.GetNode(node_index);
280
281 // Each node in the mesh must be in at most one element
282 assert(p_node->GetNumContainingElements() <= 1);
283
284 // Find a random available neighbouring node to overwrite current site
285 std::set<unsigned> neighbouring_node_indices = mpPottsMesh->GetMooreNeighbouringNodeIndices(node_index);
286 unsigned neighbour_location_index;
287
288 if (!neighbouring_node_indices.empty())
289 {
290 unsigned num_neighbours = neighbouring_node_indices.size();
291 unsigned chosen_neighbour = p_gen->randMod(num_neighbours);
292
293 std::set<unsigned>::iterator neighbour_iter = neighbouring_node_indices.begin();
294 for (unsigned j=0; j<chosen_neighbour; j++)
295 {
296 neighbour_iter++;
297 }
298
299 neighbour_location_index = *neighbour_iter;
300
301 std::set<unsigned> containing_elements = p_node->rGetContainingElementIndices();
302 std::set<unsigned> neighbour_containing_elements = GetNode(neighbour_location_index)->rGetContainingElementIndices();
303 // Only calculate Hamiltonian and update elements if the nodes are from different elements, or one is from the medium
304 if ((!containing_elements.empty() && neighbour_containing_elements.empty())
305 || (containing_elements.empty() && !neighbour_containing_elements.empty())
306 || (!containing_elements.empty() && !neighbour_containing_elements.empty() && *containing_elements.begin() != *neighbour_containing_elements.begin()))
307 {
308 double delta_H = 0.0; // This is H_1-H_0.
309
310 // Now add contributions to the Hamiltonian from each AbstractPottsUpdateRule
311 for (typename std::vector<boost::shared_ptr<AbstractUpdateRule<DIM> > >::iterator iter = this->mUpdateRuleCollection.begin();
312 iter != this->mUpdateRuleCollection.end();
313 ++iter)
314 {
315 // This static cast is fine, since we assert the update rule must be a Potts update rule in AddUpdateRule()
316 double dH = (boost::static_pointer_cast<AbstractPottsUpdateRule<DIM> >(*iter))->EvaluateHamiltonianContribution(neighbour_location_index, p_node->GetIndex(), *this);
317 delta_H += dH;
318 }
319
320 // Generate a uniform random number to do the random motion
321 double random_number = p_gen->ranf();
322
323 double p = exp(-delta_H/mTemperature);
324 if (delta_H <= 0 || random_number < p)
325 {
326 // Do swap
327
328 // Remove the current node from any elements containing it (there should be at most one such element)
329 for (std::set<unsigned>::iterator iter = containing_elements.begin();
330 iter != containing_elements.end();
331 ++iter)
332 {
333 GetElement(*iter)->DeleteNode(GetElement(*iter)->GetNodeLocalIndex(node_index));
334
336 }
337
338 // Next add the current node to any elements containing the neighbouring node (there should be at most one such element)
339 for (std::set<unsigned>::iterator iter = neighbour_containing_elements.begin();
340 iter != neighbour_containing_elements.end();
341 ++iter)
342 {
343 GetElement(*iter)->AddNode(this->mrMesh.GetNode(node_index));
344 }
345 }
346 }
347 }
348 }
349}
350
351template<unsigned DIM>
353{
354 return GetElementCorrespondingToCell(pCell)->IsDeleted();
355}
356
357template<unsigned DIM>
358void PottsBasedCellPopulation<DIM>::Update(bool hasHadBirthsOrDeaths)
359{
360}
361
362template<unsigned DIM>
364{
365 if (this->mOutputResultsForChasteVisualizer)
366 {
367 if (!this-> template HasWriter<CellPopulationElementWriter>())
368 {
369 this-> template AddPopulationWriter<CellPopulationElementWriter>();
370 }
371 }
372 // Add a CellID writer so that a VTK file will contain IDs for visualisation. (It will also dump a "loggedcell.dat" file as a side-effect.)
373 this-> template AddCellWriter<CellIdWriter>();
374
376}
377
378template<unsigned DIM>
379void PottsBasedCellPopulation<DIM>::WriteResultsToFiles(const std::string& rDirectory)
380{
381 CreateElementTessellation(); // To be used to output to the visualizer
382
384}
385
386template<unsigned DIM>
388{
389 pPopulationWriter->Visit(this);
390}
391
392template<unsigned DIM>
394{
395 pPopulationCountWriter->Visit(this);
396}
397
398template<unsigned DIM>
400{
401 pPopulationEventWriter->Visit(this);
402}
403
404template<unsigned DIM>
405void PottsBasedCellPopulation<DIM>::AcceptCellWriter(boost::shared_ptr<AbstractCellWriter<DIM, DIM> > pCellWriter, CellPtr pCell)
406{
407 pCellWriter->VisitCell(pCell, this);
408}
409
410template<unsigned DIM>
412{
413 // Get element index corresponding to this cell
414 unsigned elem_index = this->GetLocationIndexUsingCell(pCell);
415
416 // Get volume of this element in the Potts mesh
417 double cell_volume = mpPottsMesh->GetVolumeOfElement(elem_index);
418
419 return cell_volume;
420}
421
422template<unsigned DIM>
423double PottsBasedCellPopulation<DIM>::GetWidth(const unsigned& rDimension)
424{
425 // Call GetWidth() on the mesh
426 double width = this->mrMesh.GetWidth(rDimension);
427
428 return width;
429}
430
431template<unsigned DIM>
433{
434 assert(bool(dynamic_cast<AbstractPottsUpdateRule<DIM>*>(pUpdateRule.get())));
435 this->mUpdateRuleCollection.push_back(pUpdateRule);
436}
437
438template<unsigned DIM>
440{
442// delete mpElementTessellation;
443//
444// ///\todo this code would need to be extended if the domain were required to be periodic
445//
446// std::vector<Node<2>*> nodes;
447// for (unsigned node_index=0; node_index<mrMesh.GetNumNodes(); node_index++)
448// {
449// Node<2>* p_temp_node = mrMesh.GetNode(node_index);
450// nodes.push_back(p_temp_node);
451// }
452// MutableMesh<2,2> mesh(nodes);
453// mpElementTessellation = new VertexMesh<2,2>(mesh);
454}
455
456template<unsigned DIM>
458{
459// assert(mpElementTessellation != NULL);
460 return mpElementTessellation;
461}
462
463template<unsigned DIM>
465{
466 delete mpMutableMesh;
467
468 // Get the nodes of the PottsMesh
469 std::vector<Node<DIM>*> nodes;
470 for (unsigned node_index=0; node_index<this->mrMesh.GetNumNodes(); node_index++)
471 {
472 const c_vector<double, DIM>& r_location = this->mrMesh.GetNode(node_index)->rGetLocation();
473 nodes.push_back(new Node<DIM>(node_index, r_location));
474 }
475
476 mpMutableMesh = new MutableMesh<DIM,DIM>(nodes);
477}
478
479template<unsigned DIM>
481{
482 assert(mpMutableMesh);
483 return mpMutableMesh;
484}
485
486template<unsigned DIM>
488{
489 *rParamsFile << "\t\t<Temperature>" << mTemperature << "</Temperature>\n";
490 *rParamsFile << "\t\t<NumSweepsPerTimestep>" << mNumSweepsPerTimestep << "</NumSweepsPerTimestep>\n";
491
492 // Call method on direct parent class
494}
495
496template<unsigned DIM>
498{
499 mTemperature = temperature;
500}
501
502template<unsigned DIM>
504{
505 return mTemperature;
506}
507
508template<unsigned DIM>
510{
511 mNumSweepsPerTimestep = numSweepsPerTimestep;
512}
513
514template<unsigned DIM>
516{
517 return mNumSweepsPerTimestep;
518}
519
520template<unsigned DIM>
522{
523#ifdef CHASTE_VTK
524 unsigned num_timesteps = SimulationTime::Instance()->GetTimeStepsElapsed();
525 std::stringstream time;
526 time << num_timesteps;
527
528 // Create mesh writer for VTK output
529 VtkMeshWriter<DIM, DIM> mesh_writer(rDirectory, "results_"+time.str(), false);
530
531 // Iterate over any cell writers that are present
532 unsigned num_nodes = GetNumNodes();
533 for (typename std::vector<boost::shared_ptr<AbstractCellWriter<DIM, DIM> > >::iterator cell_writer_iter = this->mCellWriters.begin();
534 cell_writer_iter != this->mCellWriters.end();
535 ++cell_writer_iter)
536 {
537 // Create vector to store VTK cell data
538 std::vector<double> vtk_cell_data(num_nodes);
539
540 // Iterate over nodes in the mesh
541 for (typename AbstractMesh<DIM,DIM>::NodeIterator iter = mpPottsMesh->GetNodeIteratorBegin();
542 iter != mpPottsMesh->GetNodeIteratorEnd();
543 ++iter)
544 {
545 // Get the index of this node in the mesh and those elements (i.e. cells) that contain this node
546 unsigned node_index = iter->GetIndex();
547 std::set<unsigned> element_indices = iter->rGetContainingElementIndices();
548
549 // If there are no elements associated with this node, then we set the value of any VTK cell data to be -1 at this node...
550 if (element_indices.empty())
551 {
552 // Populate the vector of VTK cell data
553 vtk_cell_data[node_index] = -1.0;
554 }
555 else
556 {
557 // ... otherwise there should be exactly one element (i.e. cell) containing this node
558 assert(element_indices.size() == 1);
559 unsigned elem_index = *(element_indices.begin());
560 CellPtr p_cell = this->GetCellUsingLocationIndex(elem_index);
561
562 // Populate the vector of VTK cell data
563 vtk_cell_data[node_index] = (*cell_writer_iter)->GetCellDataForVtkOutput(p_cell, this);
564 }
565 }
566
567 mesh_writer.AddPointData((*cell_writer_iter)->GetVtkCellDataName(), vtk_cell_data);
568 }
569
570 // When outputting any CellData, we assume that the first cell is representative of all cells
571 unsigned num_cell_data_items = this->Begin()->GetCellData()->GetNumItems();
572 std::vector<std::string> cell_data_names = this->Begin()->GetCellData()->GetKeys();
573
574 std::vector<std::vector<double> > cell_data;
575 for (unsigned var=0; var<num_cell_data_items; var++)
576 {
577 std::vector<double> cell_data_var(num_nodes);
578 cell_data.push_back(cell_data_var);
579 }
580
581 for (typename AbstractMesh<DIM,DIM>::NodeIterator iter = mpPottsMesh->GetNodeIteratorBegin();
582 iter != mpPottsMesh->GetNodeIteratorEnd();
583 ++iter)
584 {
585 // Get the index of this node in the mesh and those elements (i.e. cells) that contain this node
586 unsigned node_index = iter->GetIndex();
587 std::set<unsigned> element_indices = iter->rGetContainingElementIndices();
588
589 // If there are no elements associated with this node, then we set the value of any VTK cell data to be -1 at this node...
590 if (element_indices.empty())
591 {
592 for (unsigned var=0; var<num_cell_data_items; var++)
593 {
594 cell_data[var][node_index] = -1.0;
595 }
596 }
597 else
598 {
599 // ... otherwise there should be exactly one element (i.e. cell) containing this node
600 assert(element_indices.size() == 1);
601 unsigned elem_index = *(element_indices.begin());
602 CellPtr p_cell = this->GetCellUsingLocationIndex(elem_index);
603
604 for (unsigned var=0; var<num_cell_data_items; var++)
605 {
606 cell_data[var][node_index] = p_cell->GetCellData()->GetItem(cell_data_names[var]);
607 }
608 }
609 }
610 for (unsigned var=0; var<cell_data.size(); var++)
611 {
612 mesh_writer.AddPointData(cell_data_names[var], cell_data[var]);
613 }
614
615 /*
616 * At present, the VTK writer can only write things which inherit from AbstractTetrahedralMeshWriter.
617 * For now, we do an explicit conversion to NodesOnlyMesh. This can be written to VTK, then visualized
618 * as glyphs in Paraview.
619 */
620 NodesOnlyMesh<DIM> temp_mesh;
621 temp_mesh.ConstructNodesWithoutMesh(*mpPottsMesh, 1.5); // This cut-off is arbitrary, as node connectivity is not used here
622 mesh_writer.WriteFilesUsingMesh(temp_mesh);
623
624 *(this->mpVtkMetaFile) << " <DataSet timestep=\"";
625 *(this->mpVtkMetaFile) << num_timesteps;
626 *(this->mpVtkMetaFile) << "\" group=\"\" part=\"0\" file=\"results_";
627 *(this->mpVtkMetaFile) << num_timesteps;
628 *(this->mpVtkMetaFile) << ".vtu\"/>\n";
629
630 // Extra part to output the outlines of cells
631 if (DIM ==2)
632 {
633 std::vector<Node<2>*> outline_nodes;
634 std::vector<VertexElement<2,2>*> outline_elements;
635
636 unsigned outline_node_index = 0;
637 unsigned outline_element_index = 0;
638 for (typename AbstractMesh<DIM,DIM>::NodeIterator iter = mpPottsMesh->GetNodeIteratorBegin();
639 iter != mpPottsMesh->GetNodeIteratorEnd();
640 ++iter)
641 {
642 // Get the index of this node in the mesh and those elements (i.e. cells) that contain this node
643 unsigned node_index = iter->GetIndex();
644 std::set<unsigned> element_indices = iter->rGetContainingElementIndices();
645
646 std::set<unsigned> target_neighbouring_node_indices = this->rGetMesh().GetVonNeumannNeighbouringNodeIndices(node_index);
647
648 for (std::set<unsigned>::iterator neighbour_iter = target_neighbouring_node_indices.begin();
649 neighbour_iter != target_neighbouring_node_indices.end();
650 ++neighbour_iter)
651 {
652 std::set<unsigned> neighbouring_element_indices = this->rGetMesh().GetNode(*neighbour_iter)->rGetContainingElementIndices();
653
654 // If different cells add a line
655 if (element_indices != neighbouring_element_indices)
656 {
657 std::vector<Node<2>*> element_nodes;
658
659 const c_vector<double, 2>& r_node_location = this->mrMesh.GetNode(node_index)->rGetLocation();
660 const c_vector<double, 2>& r_neighbour_node_location = this->mrMesh.GetNode(*neighbour_iter)->rGetLocation();
661
662 c_vector<double, 2> unit_tangent = r_neighbour_node_location - r_node_location;
663
664 if (norm_2(unit_tangent) > 1.0) // It's a periodic neighbour
665 {
666 c_vector<double, 2> mid_point = 0.5*(r_node_location + r_neighbour_node_location);
667 if (unit_tangent(0) == 0)
668 {
669 if (r_node_location(1) < r_neighbour_node_location(1))
670 {
671 mid_point(1) = r_node_location(1) - 0.5;
672 }
673 else
674 {
675 mid_point(1) = r_node_location(1) + 0.5;
676 }
677 }
678 else
679 {
680 assert(unit_tangent(1) == 0);
681
682 if (r_node_location(0) < r_neighbour_node_location(0))
683 {
684 mid_point(0) = r_node_location(0) - 0.5;
685 }
686 else
687 {
688 mid_point(0) = r_node_location(0) + 0.5;
689 }
690 }
691 assert(norm_2(unit_tangent) > 0);
692 unit_tangent = unit_tangent/norm_2(unit_tangent);
693
694 c_vector<double, DIM> unit_normal;
695 unit_normal(0) = -unit_tangent(1);
696 unit_normal(1) = unit_tangent(0);
697
698 std::vector<Node<2>*> element_nodes;
699
700 // Need at least three points to visualise in Paraview
701 Node<2>* p_node_1 = new Node<2>(outline_node_index, mid_point - 0.5*unit_normal);
702 outline_nodes.push_back(p_node_1);
703 element_nodes.push_back(outline_nodes[outline_node_index]);
704 outline_node_index++;
705
706 Node<2>* p_node_2 = new Node<2>(outline_node_index, mid_point);
707 outline_nodes.push_back(p_node_2);
708 element_nodes.push_back(outline_nodes[outline_node_index]);
709 outline_node_index++;
710
711 Node<2>* p_node_3 = new Node<2>(outline_node_index, mid_point + 0.5*unit_normal);
712 outline_nodes.push_back(p_node_3);
713 element_nodes.push_back(outline_nodes[outline_node_index]);
714 outline_node_index++;
715
716 VertexElement<2,2>* p_element = new VertexElement<2,2>(outline_element_index, element_nodes);
717 outline_elements.push_back(p_element);
718 outline_element_index++;
719
720 }
721 else // Standard Neighbour
722 {
723 c_vector<double, 2> mid_point = 0.5*(r_node_location + r_neighbour_node_location);
724
725 assert(norm_2(unit_tangent) > 0);
726 unit_tangent = unit_tangent/norm_2(unit_tangent);
727
728 c_vector<double, 2> unit_normal;
729 unit_normal(0) = -unit_tangent(1);
730 unit_normal(1) = unit_tangent(0);
731
732 std::vector<Node<2>*> element_nodes;
733
734 // Need at least three points to visualise in Paraview
735 Node<2>* p_node_1 = new Node<2>(outline_node_index, mid_point - 0.5*unit_normal);
736 outline_nodes.push_back(p_node_1);
737 element_nodes.push_back(outline_nodes[outline_node_index]);
738 outline_node_index++;
739
740 Node<2>* p_node_2 = new Node<2>(outline_node_index, mid_point);
741 outline_nodes.push_back(p_node_2);
742 element_nodes.push_back(outline_nodes[outline_node_index]);
743 outline_node_index++;
744
745 Node<2>* p_node_3 = new Node<2>(outline_node_index, mid_point + 0.5*unit_normal);
746 outline_nodes.push_back(p_node_3);
747 element_nodes.push_back(outline_nodes[outline_node_index]);
748 outline_node_index++;
749
750 VertexElement<2,2>* p_element = new VertexElement<2,2>(outline_element_index, element_nodes);
751 outline_elements.push_back(p_element);
752 outline_element_index++;
753 }
754 }
755 }
756 }
757 VertexMesh<2,2> cell_outline_mesh(outline_nodes,outline_elements);
758
759 VertexMeshWriter<2, 2> outline_mesh_writer(rDirectory, "outlines", false);
760 outline_mesh_writer.WriteVtkUsingMesh(cell_outline_mesh, time.str());
761 outline_mesh_writer.WriteFilesUsingMesh(cell_outline_mesh);
762 }
763#endif //CHASTE_VTK
764}
765
766template<unsigned DIM>
768 unsigned pdeNodeIndex,
769 std::string& rVariableName,
770 bool dirichletBoundaryConditionApplies,
771 double dirichletBoundaryValue)
772{
773 CellPtr p_cell = this->GetCellUsingLocationIndex(pdeNodeIndex);
774 double value = p_cell->GetCellData()->GetItem(rVariableName);
775 return value;
776}
777
778// Explicit instantiation
779template class PottsBasedCellPopulation<1>;
780template class PottsBasedCellPopulation<2>;
781template class PottsBasedCellPopulation<3>;
782
783// Serialization for Boost >= 1.36
#define EXCEPTION(message)
#define EXPORT_TEMPLATE_CLASS_SAME_DIMS(CLASS)
virtual void OpenWritersFiles(OutputFileHandler &rOutputFileHandler)
AbstractMesh< ELEMENT_DIM, SPACE_DIM > & mrMesh
virtual void WriteResultsToFiles(const std::string &rDirectory)
Node< SPACE_DIM > * GetNode(unsigned index) const
virtual void OutputCellPopulationParameters(out_stream &rParamsFile)
Definition Node.hpp:59
std::set< unsigned > & rGetContainingElementIndices()
Definition Node.cpp:300
unsigned GetNumContainingElements() const
Definition Node.cpp:312
unsigned GetIndex() const
Definition Node.cpp:158
void ConstructNodesWithoutMesh(const std::vector< Node< SPACE_DIM > * > &rNodes, double maxInteractionDistance)
virtual void OpenWritersFiles(OutputFileHandler &rOutputFileHandler)
void Update(bool hasHadBirthsOrDeaths=true)
virtual void AcceptPopulationEventWriter(boost::shared_ptr< AbstractCellPopulationEventWriter< DIM, DIM > > pPopulationEventWriter)
VertexMesh< DIM, DIM > * GetElementTessellation()
virtual void WriteResultsToFiles(const std::string &rDirectory)
virtual double GetCellDataItemAtPdeNode(unsigned pdeNodeIndex, std::string &rVariableName, bool dirichletBoundaryConditionApplies=false, double dirichletBoundaryValue=0.0)
virtual void AcceptCellWriter(boost::shared_ptr< AbstractCellWriter< DIM, DIM > > pCellWriter, CellPtr pCell)
PottsElement< DIM > * GetElementCorrespondingToCell(CellPtr pCell)
virtual void AcceptPopulationWriter(boost::shared_ptr< AbstractCellPopulationWriter< DIM, DIM > > pPopulationWriter)
c_vector< double, DIM > GetLocationOfCellCentre(CellPtr pCell)
PottsBasedCellPopulation(PottsMesh< DIM > &rMesh, std::vector< CellPtr > &rCells, bool deleteMesh=false, bool validate=true, const std::vector< unsigned > locationIndices=std::vector< unsigned >())
std::set< unsigned > GetNeighbouringLocationIndices(CellPtr pCell)
PottsElement< DIM > * GetElement(unsigned elementIndex)
bool IsCellAssociatedWithADeletedLocation(CellPtr pCell)
void SetNumSweepsPerTimestep(unsigned numSweepsPerTimestep)
CellPtr AddCell(CellPtr pNewCell, CellPtr pParentCell=CellPtr())
Node< DIM > * GetNode(unsigned index)
virtual void AcceptPopulationCountWriter(boost::shared_ptr< AbstractCellPopulationCountWriter< DIM, DIM > > pPopulationCountWriter)
virtual void WriteVtkResultsToFile(const std::string &rDirectory)
virtual TetrahedralMesh< DIM, DIM > * GetTetrahedralMeshForPdeModifier()
double GetWidth(const unsigned &rDimension)
MutableMesh< DIM, DIM > * GetMutableMesh()
virtual void AddUpdateRule(boost::shared_ptr< AbstractUpdateRule< DIM > > pUpdateRule)
void OutputCellPopulationParameters(out_stream &rParamsFile)
void SetTemperature(double temperature)
static RandomNumberGenerator * Instance()
void Shuffle(std::vector< boost::shared_ptr< T > > &rValues)
unsigned randMod(unsigned base)
static SimulationTime * Instance()
unsigned GetTimeStepsElapsed() const
void WriteVtkUsingMesh(VertexMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, std::string stamp="")
void WriteFilesUsingMesh(VertexMesh< ELEMENT_DIM, SPACE_DIM > &rMesh)
void AddPointData(std::string name, std::vector< double > data)
void WriteFilesUsingMesh(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, bool keepOriginalElementIndexing=true)