Chaste  Build::
DiscreteContinuumMesh.cpp
1 /*
2 
3 Copyright (c) 2005-2016, University of Oxford.
4  All rights reserved.
5 
6  University of Oxford means the Chancellor, Masters and Scholars of the
7  University of Oxford, having an administrative office at Wellington
8  Square, Oxford OX1 2JD, UK.
9 
10  This file is part of Chaste.
11 
12  Redistribution and use in source and binary forms, with or without
13  modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24  AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25  IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26  ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27  LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28  CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29  GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30  HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31  LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32  OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34  */
35 
36 #include <boost/lexical_cast.hpp>
37 #define _BACKWARD_BACKWARD_WARNING_H 1 //Cut out the vtk deprecated warning
38 #include <vtkIdList.h>
39 #include <vtkCellArray.h>
40 #include <vtkLine.h>
41 #include <vtkPoints.h>
42 #include <vtkTriangle.h>
43 #include <vtkPoints.h>
44 #include <vtkTetra.h>
45 #include <vtkGenericCell.h>
46 #include "Exception.hpp"
47 #include "Warnings.hpp"
48 #include "Facet.hpp"
49 #include "Polygon.hpp"
50 #include "DiscreteContinuumMesh.hpp"
51 #include "Element.hpp"
52 #include "UblasVectorInclude.hpp"
53 #include "AbstractTetrahedralMesh.hpp"
54 #include "VtkMeshWriter.hpp"
55 #include "BaseUnits.hpp"
56 
57 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
59  mAttributes(),
60  mReferenceLength(BaseUnits::Instance()->GetReferenceLengthScale()),
61  mpVtkMesh(),
62  mpVtkCellLocator(vtkSmartPointer<vtkCellLocator>::New()),
63  mVtkRepresentationUpToDate(false),
64  mPointElementMap(),
65  mSegmentElementMap(),
66  mCellElementMap(),
67  mpNetwork(),
68  mpCellPopulation(),
69  mCellPopulationReferenceLength(BaseUnits::Instance()->GetReferenceLengthScale())
70 {
71 
72 }
73 
74 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
75 boost::shared_ptr<DiscreteContinuumMesh<ELEMENT_DIM, SPACE_DIM> > DiscreteContinuumMesh<ELEMENT_DIM, SPACE_DIM>::Create()
76 {
77  MAKE_PTR(DiscreteContinuumMesh<ELEMENT_DIM>, pSelf);
78  return pSelf;
79 }
80 
81 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
83 {
84 
85 }
86 
87 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
89 {
90  return mReferenceLength;
91 }
92 
93 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
95 {
97  {
98  mpVtkMesh = vtkSmartPointer<vtkUnstructuredGrid>::New();
99  vtkSmartPointer<vtkPoints> p_vtk_points = vtkSmartPointer<vtkPoints>::New();
100  std::vector<c_vector<double, SPACE_DIM> > node_locations = GetNodeLocations();
101  p_vtk_points->SetNumberOfPoints(node_locations.size());
102  for(unsigned idx=0; idx<node_locations.size(); idx++)
103  {
104  if(SPACE_DIM==3)
105  {
106  p_vtk_points->InsertPoint(idx, node_locations[idx][0], node_locations[idx][1], node_locations[idx][2]);
107  }
108  else
109  {
110  p_vtk_points->InsertPoint(idx, node_locations[idx][0], node_locations[idx][1], 0.0);
111  }
112  }
113  mpVtkMesh->SetPoints(p_vtk_points);
114 
115  // Add vtk tets or triangles
116  std::vector<std::vector<unsigned> > element_connectivity = GetConnectivity();
117  unsigned num_elements = element_connectivity.size();
118  mpVtkMesh->Allocate(num_elements, num_elements);
119 
120  for(unsigned idx=0; idx<num_elements; idx++)
121  {
122  if(ELEMENT_DIM==3)
123  {
124  vtkSmartPointer<vtkTetra> p_vtk_element = vtkSmartPointer<vtkTetra>::New();
125  unsigned num_nodes = element_connectivity[idx].size();
126  for(unsigned jdx=0; jdx<num_nodes; jdx++)
127  {
128  p_vtk_element->GetPointIds()->SetId(jdx, element_connectivity[idx][jdx]);
129  }
130  mpVtkMesh->InsertNextCell(p_vtk_element->GetCellType(), p_vtk_element->GetPointIds());
131  }
132  else
133  {
134  vtkSmartPointer<vtkTriangle> p_vtk_element = vtkSmartPointer<vtkTriangle>::New();
135  unsigned num_nodes = element_connectivity[idx].size();
136  for(unsigned jdx=0; jdx<num_nodes; jdx++)
137  {
138  p_vtk_element->GetPointIds()->SetId(jdx, element_connectivity[idx][jdx]);
139  }
140  mpVtkMesh->InsertNextCell(p_vtk_element->GetCellType(), p_vtk_element->GetPointIds());
141  }
142  }
143  mpVtkCellLocator = vtkSmartPointer<vtkCellLocator>::New();
144  mpVtkCellLocator->SetDataSet(mpVtkMesh);
145  mpVtkCellLocator->BuildLocator();
146 
148  }
149  return mpVtkMesh;
150 }
151 
152 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
154 {
155  return mAttributes;
156 }
157 
158 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
160 {
161  std::vector<c_vector<double, SPACE_DIM> > centroids(this->GetNumElements());
162  for(unsigned idx=0; idx<this->GetNumElements(); idx++)
163  {
164  centroids[idx] = this->GetElement(idx)->CalculateCentroid();
165  }
166  return centroids;
167 }
168 
169 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
171 {
173  {
175  }
176 
177  mPointElementMap.clear();
178  mPointElementMap = std::vector<std::vector<unsigned> >(this->GetNumElements());
179 
180  std::vector<int> cell_indices(points.size());
181  for(unsigned idx=0; idx<points.size(); idx++)
182  {
183  double x_coords[3];
184  x_coords[0] = points[idx].GetLocation(mReferenceLength)[0];
185  x_coords[1] = points[idx].GetLocation(mReferenceLength)[1];
186  if(SPACE_DIM == 3)
187  {
188  x_coords[2] = points[idx].GetLocation(mReferenceLength)[2];
189  }
190  else
191  {
192  x_coords[2] = 0.0;
193  }
194 
195  int cell_id=-1;
196  cell_id = mpVtkCellLocator->FindCell(x_coords);
197 
198  if(cell_id>=0)
199  {
200  mPointElementMap[cell_id].push_back(idx);
201  }
202  }
203  return mPointElementMap;
204 }
205 
206 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
207 const std::vector<std::vector<CellPtr> >& DiscreteContinuumMesh<ELEMENT_DIM, SPACE_DIM>::GetElementCellMap(bool update)
208 {
209  if (!update)
210  {
211  return mCellElementMap;
212  }
213 
215  {
217  }
218 
219  if (!mpCellPopulation)
220  {
221  EXCEPTION("A cell population has not been set. Can not create a cell element map.");
222  }
223 
224  // Loop over all cells and associate cells with the points
225  mCellElementMap.clear();
226  mCellElementMap = std::vector<std::vector<CellPtr> >(this->GetNumElements());
227 
228  double cell_mesh_length_scaling = mCellPopulationReferenceLength/mReferenceLength;
229  for (typename AbstractCellPopulation<SPACE_DIM>::Iterator cell_iter = mpCellPopulation->Begin();
230  cell_iter != mpCellPopulation->End(); ++cell_iter)
231  {
232  c_vector<double, SPACE_DIM> location = mpCellPopulation->GetLocationOfCellCentre(*cell_iter);
233  double x_coords[3];
234  x_coords[0] = location[0]*cell_mesh_length_scaling;
235  x_coords[1] = location[1]*cell_mesh_length_scaling;
236  if(SPACE_DIM == 3)
237  {
238  x_coords[2] = location[2]*cell_mesh_length_scaling;
239  }
240  else
241  {
242  x_coords[2] = 0.0;
243  }
244  int cell_id=-1;
245  cell_id = mpVtkCellLocator->FindCell(x_coords);
246 
247  if(cell_id>=0)
248  {
249  mCellElementMap[cell_id].push_back(*cell_iter);
250  }
251  }
252  return mCellElementMap;
253 }
254 
255 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
256 std::vector<std::vector<boost::shared_ptr<VesselSegment<SPACE_DIM> > > > DiscreteContinuumMesh<ELEMENT_DIM, SPACE_DIM>::GetElementSegmentMap(bool update,
257  bool useVesselSurface)
258 {
259  if (!update)
260  {
261  return mSegmentElementMap;
262  }
263 
264  if (!mpNetwork)
265  {
266  EXCEPTION("A vessel network has not been set. Can not create a vessel element map.");
267  }
268 
269  mSegmentElementMap.clear();
270  unsigned num_elements = this->GetNumElements();
271  mSegmentElementMap = std::vector<std::vector<boost::shared_ptr<VesselSegment<SPACE_DIM> > > >(num_elements);
272 
273  std::vector<boost::shared_ptr<VesselSegment<SPACE_DIM> > > segments = mpNetwork->GetVesselSegments();
274 
275  for (unsigned jdx = 0; jdx < segments.size(); jdx++)
276  {
277  DimensionalChastePoint<SPACE_DIM> loc1 = segments[jdx]->GetNode(0)->rGetLocation();
278  DimensionalChastePoint<SPACE_DIM> loc2 = segments[jdx]->GetNode(1)->rGetLocation();
279  double x_coords1[3];
280  x_coords1[0] = loc1.GetLocation(mReferenceLength)[0];
281  x_coords1[1] = loc1.GetLocation(mReferenceLength)[1];
282  if(SPACE_DIM == 3)
283  {
284  x_coords1[2] = loc1.GetLocation(mReferenceLength)[2];
285  }
286  else
287  {
288  x_coords1[2] = 0.0;
289  }
290  double x_coords2[3];
291  x_coords2[0] = loc2.GetLocation(mReferenceLength)[0];
292  x_coords2[1] = loc2.GetLocation(mReferenceLength)[1];
293  if(SPACE_DIM == 3)
294  {
295  x_coords2[2] = loc2.GetLocation(mReferenceLength)[2];
296  }
297  else
298  {
299  x_coords2[2] = 0.0;
300  }
301 
302  vtkSmartPointer<vtkIdList> p_id_list = vtkSmartPointer<vtkIdList>::New();
303  mpVtkCellLocator->FindCellsAlongLine(x_coords1, x_coords2, 1.e-8, p_id_list);
304  unsigned num_intersections = p_id_list->GetNumberOfIds();
305  for(unsigned idx=0; idx<num_intersections; idx++)
306  {
307  mSegmentElementMap[p_id_list->GetId(idx)].push_back(segments[jdx]);
308  }
309  }
310 
311  return mSegmentElementMap;
312 }
313 
314 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
316 {
317  std::vector<std::vector<unsigned> > connectivity;
318  unsigned num_elements = AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetNumElements();
319  for (unsigned idx = 0; idx < num_elements; idx++)
320  {
321  std::vector<unsigned> node_indexes;
322  unsigned num_nodes = AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetElement(idx)->GetNumNodes();
323  for (unsigned jdx = 0; jdx < num_nodes; jdx++)
324  {
325  node_indexes.push_back(AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetElement(idx)->GetNodeGlobalIndex(jdx));
326  }
327  connectivity.push_back(node_indexes);
328  }
329  return connectivity;
330 }
331 
332 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
333 std::vector<c_vector<double, SPACE_DIM> > DiscreteContinuumMesh<ELEMENT_DIM, SPACE_DIM>::GetNodeLocations()
334 {
335  unsigned num_nodes = AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetNumNodes();
336  std::vector<c_vector<double, SPACE_DIM> > locations(num_nodes);
337  for (unsigned idx = 0; idx < num_nodes; idx++)
338  {
339  locations[idx] = AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetNode(idx)->rGetLocation();
340  }
341  return locations;
342 }
343 
344 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
345 std::vector<DimensionalChastePoint<SPACE_DIM> > DiscreteContinuumMesh<ELEMENT_DIM, SPACE_DIM>::GetNodeLocationsAsPoints()
346 {
347  unsigned num_nodes = AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetNumNodes();
348  std::vector<DimensionalChastePoint<SPACE_DIM> > locations(num_nodes);
349  for (unsigned idx = 0; idx < num_nodes; idx++)
350  {
351  locations[idx] = DimensionalChastePoint<SPACE_DIM>(AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetNode(idx)->rGetLocation(), mReferenceLength);
352  }
353  return locations;
354 }
355 
356 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
357 void DiscreteContinuumMesh<ELEMENT_DIM, SPACE_DIM>::ImportDiscreteContinuumMeshFromTetgen(tetgen::tetgenio& mesherOutput, unsigned numberOfElements,
358  int *elementList, unsigned numberOfFaces, int *faceList,
359  int *edgeMarkerList)
360 {
361  unsigned nodes_per_element = mesherOutput.numberofcorners;
362 
363  assert(nodes_per_element == ELEMENT_DIM + 1 || nodes_per_element == (ELEMENT_DIM + 1) * (ELEMENT_DIM + 2) / 2);
364 
365  for (unsigned i = 0; i < this->mBoundaryElements.size(); i++)
366  {
367  delete this->mBoundaryElements[i];
368  }
369  for (unsigned i = 0; i < this->mElements.size(); i++)
370  {
371  delete this->mElements[i];
372  }
373  for (unsigned i = 0; i < this->mNodes.size(); i++)
374  {
375  delete this->mNodes[i];
376  }
377 
378  this->mNodes.clear();
379  this->mElements.clear();
380  this->mBoundaryElements.clear();
381  this->mBoundaryNodes.clear();
382 
383  // Construct the nodes
384  for (unsigned node_index = 0; node_index < (unsigned) mesherOutput.numberofpoints; node_index++)
385  {
386  this->mNodes.push_back(new Node<SPACE_DIM>(node_index, &mesherOutput.pointlist[node_index * SPACE_DIM], false));
387  }
388 
389  // Construct the elements
390  this->mElements.reserve(numberOfElements);
391 
392  unsigned real_element_index = 0;
393  for (unsigned element_index = 0; element_index < numberOfElements; element_index++)
394  {
395  std::vector<Node<SPACE_DIM>*> nodes;
396  for (unsigned j = 0; j < ELEMENT_DIM + 1; j++)
397  {
398  unsigned global_node_index = elementList[element_index * (nodes_per_element) + j];
399  assert(global_node_index < this->mNodes.size());
400  nodes.push_back(this->mNodes[global_node_index]);
401 
402  }
403 
404  /*
405  * For some reason, tetgen in library mode makes its initial Delaunay mesh
406  * with very thin slivers. Hence we expect to ignore some of the elements!
407  */
408  Element<ELEMENT_DIM, SPACE_DIM>* p_element;
409  try
410  {
411  p_element = new Element<ELEMENT_DIM, SPACE_DIM>(real_element_index, nodes);
412 
413  // Shouldn't throw after this point
414  this->mElements.push_back(p_element);
415 
416  // Add the internals to quadratics
417  for (unsigned j = ELEMENT_DIM + 1; j < nodes_per_element; j++)
418  {
419  unsigned global_node_index = elementList[element_index * nodes_per_element + j];
420  assert(global_node_index < this->mNodes.size());
421  this->mElements[real_element_index]->AddNode(this->mNodes[global_node_index]);
422  this->mNodes[global_node_index]->AddElement(real_element_index);
423  this->mNodes[global_node_index]->MarkAsInternal();
424  }
425  real_element_index++;
426  }
427  catch (Exception &)
428  {
429  if (SPACE_DIM == 2)
430  {
431  WARNING("Triangle has produced a zero area (collinear) element");
432  }
433  else
434  {
435  WARNING("Tetgen has produced a zero volume (coplanar) element");
436  }
437  }
438  }
439 
440  // Construct the BoundaryElements (and mark boundary nodes)
441  unsigned next_boundary_element_index = 0;
442  for (unsigned boundary_element_index = 0; boundary_element_index < numberOfFaces; boundary_element_index++)
443  {
444  /*
445  * Tetgen produces only boundary faces (set edgeMarkerList to NULL).
446  * Triangle marks which edges are on the boundary.
447  */
448  if (edgeMarkerList == NULL || edgeMarkerList[boundary_element_index] == 1)
449  {
450  std::vector<Node<SPACE_DIM>*> nodes;
451  for (unsigned j = 0; j < ELEMENT_DIM; j++)
452  {
453  unsigned global_node_index = faceList[boundary_element_index * ELEMENT_DIM + j];
454  assert(global_node_index < this->mNodes.size());
455  nodes.push_back(this->mNodes[global_node_index]);
456  if (!nodes[j]->IsBoundaryNode())
457  {
458  nodes[j]->SetAsBoundaryNode();
459  this->mBoundaryNodes.push_back(nodes[j]);
460  }
461  }
462 
463  /*
464  * For some reason, tetgen in library mode makes its initial Delaunay mesh
465  * with very thin slivers. Hence we expect to ignore some of the elements!
466  */
467  BoundaryElement<ELEMENT_DIM - 1, SPACE_DIM>* p_b_element;
468  try
469  {
470  p_b_element = new BoundaryElement<ELEMENT_DIM - 1, SPACE_DIM>(next_boundary_element_index, nodes);
471  this->mBoundaryElements.push_back(p_b_element);
472  next_boundary_element_index++;
473  }
474  catch (Exception &)
475  {
476  // Tetgen is feeding us lies //Watch this space for coverage
477  assert(SPACE_DIM == 3);
478  }
479  }
480  }
481 
482  this->RefreshJacobianCachedData();
483 }
484 
485 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
487 {
488  mAttributes = attributes;
489 }
490 
491 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
492 void DiscreteContinuumMesh<ELEMENT_DIM, SPACE_DIM>::SetCellPopulation(AbstractCellPopulation<SPACE_DIM>& rCellPopulation,
493  units::quantity<unit::length> cellLengthScale)
494 {
495  mpCellPopulation = &rCellPopulation;
496  mCellPopulationReferenceLength = cellLengthScale;
497 }
498 
499 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
501 {
502  mpNetwork = pNetwork;
503 }
504 
505 // Explicit instantiation
506 template class DiscreteContinuumMesh<2> ;
507 template class DiscreteContinuumMesh<3> ;
vtkSmartPointer< vtkUnstructuredGrid > mpVtkMesh
std::vector< c_vector< double, SPACE_DIM > > GetElementCentroids()
std::vector< c_vector< double, SPACE_DIM > > GetNodeLocations()
std::vector< unsigned > GetElementRegionMarkers()
units::quantity< unit::length > mCellPopulationReferenceLength
std::vector< std::vector< unsigned > > GetConnectivity()
std::vector< std::vector< CellPtr > > mCellElementMap
boost::shared_ptr< VesselNetwork< SPACE_DIM > > mpNetwork
const std::vector< std::vector< CellPtr > > & GetElementCellMap(bool update=true)
std::vector< std::vector< boost::shared_ptr< VesselSegment< SPACE_DIM > > > > mSegmentElementMap
std::vector< std::vector< boost::shared_ptr< VesselSegment< SPACE_DIM > > > > GetElementSegmentMap(bool update=true, bool useVesselSurface=false)
void ImportDiscreteContinuumMeshFromTetgen(tetgen::tetgenio &mesherOutput, unsigned numberOfElements, int *elementList, unsigned numberOfFaces, int *faceList, int *edgeMarkerList)
vtkSmartPointer< vtkUnstructuredGrid > GetAsVtkUnstructuredGrid()
units::quantity< unit::length > mReferenceLength
c_vector< double, DIM > GetLocation(units::quantity< unit::length > scale)
std::vector< std::vector< unsigned > > GetPointElementMap(std::vector< DimensionalChastePoint< SPACE_DIM > > points)
void SetCellPopulation(AbstractCellPopulation< SPACE_DIM > &rCellPopulation, units::quantity< unit::length > cellLengthScale)
vtkSmartPointer< vtkCellLocator > mpVtkCellLocator
std::vector< DimensionalChastePoint< SPACE_DIM > > GetNodeLocationsAsPoints()
void SetAttributes(std::vector< unsigned > attributes)
units::quantity< unit::length > GetReferenceLengthScale()
static boost::shared_ptr< DiscreteContinuumMesh< ELEMENT_DIM, SPACE_DIM > > Create()
void SetVesselNetwork(boost::shared_ptr< VesselNetwork< SPACE_DIM > > pNetwork)
std::vector< std::vector< unsigned > > mPointElementMap
std::vector< unsigned > mAttributes
AbstractCellPopulation< SPACE_DIM > * mpCellPopulation