Chaste  Release::2017.1
TetrahedralMesh.cpp
1 /*
2 
3 Copyright (c) 2005-2017, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34 */
35 
36 #include "TetrahedralMesh.hpp"
37 
38 #include <iostream>
39 #include <cassert>
40 #include <sstream>
41 #include <map>
42 #include <limits>
43 
44 #include "BoundaryElement.hpp"
45 #include "Element.hpp"
46 #include "Exception.hpp"
47 #include "Node.hpp"
48 #include "OutputFileHandler.hpp"
49 #include "PetscTools.hpp"
50 #include "RandomNumberGenerator.hpp"
51 #include "Warnings.hpp"
52 
53 // Jonathan Shewchuk's triangle and Hang Si's tetgen
54 #define REAL double
55 #define VOID void
56 #include "triangle.h"
57 #include "tetgen.h"
58 #undef REAL
59 #undef VOID
60 
61 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
63 {
64  Clear();
65 }
66 
67 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
70 {
71  assert(rMeshReader.HasNodePermutation() == false);
72  this->mMeshFileBaseName = rMeshReader.GetMeshFileBaseName();
73 
74  // Record number of corner nodes
75  unsigned num_nodes = rMeshReader.GetNumNodes();
76 
77  /*
78  * Reserve memory for nodes, so we don't have problems with
79  * pointers stored in elements becoming invalid.
80  */
81  this->mNodes.reserve(num_nodes);
82 
83  rMeshReader.Reset();
84 
85  //typename std::map<std::pair<unsigned,unsigned>,unsigned>::const_iterator iterator;
86  //std::map<std::pair<unsigned,unsigned>,unsigned> internal_nodes_map;
87 
88  // Add nodes
89  std::vector<double> coords;
90  for (unsigned node_index=0; node_index < num_nodes; node_index++)
91  {
92  coords = rMeshReader.GetNextNode();
93  Node<SPACE_DIM>* p_node = new Node<SPACE_DIM>(node_index, coords, false);
94 
95  for (unsigned i=0; i<rMeshReader.GetNodeAttributes().size(); i++)
96  {
97  double attribute = rMeshReader.GetNodeAttributes()[i];
98  p_node->AddNodeAttribute(attribute);
99 
100  }
101  this->mNodes.push_back(p_node);
102  }
103 
104  //unsigned new_node_index = mNumCornerNodes;
105 
106  rMeshReader.Reset();
107  // Add elements
108  //new_node_index = mNumCornerNodes;
109  this->mElements.reserve(rMeshReader.GetNumElements());
110 
111  for (unsigned element_index=0; element_index < (unsigned) rMeshReader.GetNumElements(); element_index++)
112  {
113  ElementData element_data = rMeshReader.GetNextElementData();
114  std::vector<Node<SPACE_DIM>*> nodes;
115 
116  /*
117  * NOTE: currently just reading element vertices from mesh reader - even if it
118  * does contain information about internal nodes (ie for quadratics) this is
119  * ignored here and used elsewhere: ie don't do this:
120  * unsigned nodes_size = node_indices.size();
121  */
122  for (unsigned j=0; j<ELEMENT_DIM+1; j++) // num vertices=ELEMENT_DIM+1, may not be equal to nodes_size.
123  {
124  assert(element_data.NodeIndices[j] < this->mNodes.size());
125  nodes.push_back(this->mNodes[element_data.NodeIndices[j]]);
126  }
127 
128  Element<ELEMENT_DIM,SPACE_DIM>* p_element = new Element<ELEMENT_DIM,SPACE_DIM>(element_index, nodes);
129 
130  this->mElements.push_back(p_element);
131 
132  if (rMeshReader.GetNumElementAttributes() > 0)
133  {
134  assert(rMeshReader.GetNumElementAttributes() == 1);
135  double attribute_value = element_data.AttributeValue;
136  p_element->SetAttribute(attribute_value);
137  }
138  }
139 
140  // Add boundary elements and nodes
141  for (unsigned face_index=0; face_index<(unsigned)rMeshReader.GetNumFaces(); face_index++)
142  {
143  ElementData face_data = rMeshReader.GetNextFaceData();
144  std::vector<unsigned> node_indices = face_data.NodeIndices;
145 
146  /*
147  * NOTE: unlike the above where we just read element *vertices* from mesh reader, here we are
148  * going to read a quadratic mesh with internal elements.
149  * (There are only a few meshes with internals in the face file that we might as well use them.)
150  *
151  */
152  std::vector<Node<SPACE_DIM>*> nodes;
153  for (unsigned node_index=0; node_index<node_indices.size(); node_index++)
154  {
155  assert(node_indices[node_index] < this->mNodes.size());
156  // Add Node pointer to list for creating an element
157  nodes.push_back(this->mNodes[node_indices[node_index]]);
158  }
159 
160  // This is a boundary face, so ensure all its nodes are marked as boundary nodes
161  for (unsigned j=0; j<nodes.size(); j++)
162  {
163  if (!nodes[j]->IsBoundaryNode())
164  {
165  nodes[j]->SetAsBoundaryNode();
166  this->mBoundaryNodes.push_back(nodes[j]);
167  }
168 
169  // Register the index that this bounday element will have with the node
170  nodes[j]->AddBoundaryElement(face_index);
171  }
172 
173  // The added elements will be deleted in our destructor
174  BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>* p_boundary_element = new BoundaryElement<ELEMENT_DIM-1,SPACE_DIM>(face_index, nodes);
175  this->mBoundaryElements.push_back(p_boundary_element);
176 
177  if (rMeshReader.GetNumFaceAttributes() > 0)
178  {
179  assert(rMeshReader.GetNumFaceAttributes() == 1);
180  double attribute_value = face_data.AttributeValue;
181  p_boundary_element->SetAttribute(attribute_value);
182  }
183  }
184 
185  RefreshJacobianCachedData();
186 
187  rMeshReader.Reset();
188 }
189 
190 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
191 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ReadNodesPerProcessorFile(const std::string& rNodesPerProcessorFile)
192 {
193  std::vector<unsigned> nodes_per_processor_vec;
194 
195  std::ifstream file_stream(rNodesPerProcessorFile.c_str());
196  if (file_stream.is_open())
197  {
198  while (file_stream)
199  {
200  unsigned nodes_per_processor;
201  file_stream >> nodes_per_processor;
202 
203  if (file_stream)
204  {
205  nodes_per_processor_vec.push_back(nodes_per_processor);
206  }
207  }
208  }
209  else
210  {
211  EXCEPTION("Unable to read nodes per processor file " + rNodesPerProcessorFile);
212  }
213 
214  unsigned sum = 0;
215  for (unsigned i=0; i<nodes_per_processor_vec.size(); i++)
216  {
217  sum += nodes_per_processor_vec[i];
218  }
219 
220  if (sum != this->GetNumNodes())
221  {
222  EXCEPTION("Sum of nodes per processor, " << sum
223  << ", not equal to number of nodes in mesh, " << this->GetNumNodes());
224  }
225 
226  unsigned num_owned=nodes_per_processor_vec[PetscTools::GetMyRank()];
227 
228  if (nodes_per_processor_vec.size() != PetscTools::GetNumProcs())
229  {
230  EXCEPTION("Number of processes doesn't match the size of the nodes-per-processor file");
231  }
232  delete this->mpDistributedVectorFactory;
233  this->mpDistributedVectorFactory = new DistributedVectorFactory(this->GetNumNodes(), num_owned);
234 }
235 
236 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
238 {
239  /*
240  * Each face of each element is a set of node indices.
241  * We form a set of these in order to get their parity:
242  * all faces which appear once are inserted into the set;
243  * all faces which appear twice are inserted and then removed from the set;
244  * we're assuming that faces never appear more than twice.
245  */
246  std::set< std::set<unsigned> > odd_parity_faces;
247 
248  for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = this->GetElementIteratorBegin();
249  iter != this->GetElementIteratorEnd();
250  ++iter)
251  {
252  for (unsigned face_index=0; face_index<=ELEMENT_DIM; face_index++)
253  {
254  std::set<unsigned> face_info;
255  for (unsigned node_index=0; node_index<=ELEMENT_DIM; node_index++)
256  {
257  // Leave one index out each time
258  if (node_index != face_index)
259  {
260  face_info.insert(iter->GetNodeGlobalIndex(node_index));
261  }
262  }
263  // Face is now formed - attempt to find it
264  std::set< std::set<unsigned> >::iterator find_face=odd_parity_faces.find(face_info);
265  if (find_face != odd_parity_faces.end())
266  {
267  // Face was in set, so it now has even parity.
268  // Remove it via the iterator
269  odd_parity_faces.erase(find_face);
270  }
271  else
272  {
273  // Face is not in set so it now has odd parity. Insert it
274  odd_parity_faces.insert(face_info);
275  }
276 
277  }
278  }
279 
280  /*
281  * At this point the odd parity faces should be the same as the
282  * boundary elements. We could check this explicitly or we could
283  * just count them.
284  */
285  return(odd_parity_faces.size() == this->GetNumBoundaryElements());
286 }
287 
288 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
290 {
291  double mesh_volume = 0.0;
292 
293  for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = this->GetElementIteratorBegin();
294  iter != this->GetElementIteratorEnd();
295  ++iter)
296  {
297  mesh_volume += iter->GetVolume(mElementJacobianDeterminants[iter->GetIndex()]);
298  }
299 
300  return mesh_volume;
301 }
302 
303 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
305 {
306  // ELEMENT_DIM-1 is the dimension of the boundary element
307  assert(ELEMENT_DIM >= 1);
308  const unsigned bound_element_dim = ELEMENT_DIM-1;
309  assert(bound_element_dim < 3);
310  if (bound_element_dim == 0)
311  {
312  return 0.0;
313  }
314 
315  double mesh_surface = 0.0;
316  typename TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::BoundaryElementIterator it = this->GetBoundaryElementIteratorBegin();
317 
318  while (it != this->GetBoundaryElementIteratorEnd())
319  {
320  mesh_surface += mBoundaryElementJacobianDeterminants[(*it)->GetIndex()];
321  it++;
322  }
323 
324  if (bound_element_dim == 2)
325  {
326  mesh_surface /= 2.0;
327  }
328 
329  return mesh_surface;
330 }
331 
332 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
334 {
335  // Make a permutation vector of the identity
337  std::vector<unsigned> perm;
338  p_rng->Shuffle(this->mNodes.size(), perm);
339 
340  // Call the non-random version
341  PermuteNodes(perm);
342 }
343 
344 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
345 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::PermuteNodes(const std::vector<unsigned>& perm)
346 {
347  // Let's not do this if there are any deleted nodes
348  assert( this->GetNumAllNodes() == this->GetNumNodes());
349 
350  assert(perm.size() == this->mNodes.size());
351 
352  // Copy the node pointers
353  std::vector< Node<SPACE_DIM>* > copy_m_nodes;
354  copy_m_nodes.assign(this->mNodes.begin(), this->mNodes.end());
355 
356  for (unsigned original_index=0; original_index<this->mNodes.size(); original_index++)
357  {
358  assert(perm[original_index] < this->mNodes.size());
359  //perm[original_index] holds the new assigned index of that node
360  this->mNodes[ perm[original_index] ] = copy_m_nodes[original_index];
361  }
362 
363  // Update indices
364  for (unsigned index=0; index<this->mNodes.size(); index++)
365  {
366  this->mNodes[index]->SetIndex(index);
367  }
368 
369  // Copy the permutation vector into the mesh
370  this->mNodePermutation = perm;
371 }
372 
373 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
374 unsigned TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetContainingElementIndexWithInitialGuess(const ChastePoint<SPACE_DIM>& rTestPoint, unsigned startingElementGuess, bool strict)
375 {
376  assert(startingElementGuess<this->GetNumElements());
377 
378  /*
379  * Let m=startingElementGuess, N=num_elem-1.
380  * We search from in this order: m, m+1, m+2, .. , N, 0, 1, .., m-1.
381  */
382  unsigned i = startingElementGuess;
383  bool reached_end = false;
384 
385  while (!reached_end)
386  {
387  if (this->mElements[i]->IncludesPoint(rTestPoint, strict))
388  {
389  assert(!this->mElements[i]->IsDeleted());
390  return i;
391  }
392 
393  // Increment
394  i++;
395  if (i==this->GetNumElements())
396  {
397  i=0;
398  }
399 
400  // Back to the beginning yet?
401  if (i==startingElementGuess)
402  {
403  reached_end = true;
404  }
405  }
406 
407  // If it's in none of the elements, then throw
408  std::stringstream ss;
409  ss << "Point [";
410  for (unsigned j=0; (int)j<(int)SPACE_DIM-1; j++)
411  {
412  ss << rTestPoint[j] << ",";
413  }
414  ss << rTestPoint[SPACE_DIM-1] << "] is not in mesh - all elements tested";
415  EXCEPTION(ss.str());
416 }
417 
418 
419 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
421 {
422  EXCEPT_IF_NOT(ELEMENT_DIM == SPACE_DIM); // LCOV_EXCL_LINE // CalculateInterpolationWeights hits an assertion otherwise
423  double max_min_weight = -std::numeric_limits<double>::infinity();
424  unsigned closest_index = 0;
425  for (unsigned i=0; i<this->mElements.size(); i++)
426  {
427  c_vector<double, ELEMENT_DIM+1> weight = this->mElements[i]->CalculateInterpolationWeights(rTestPoint);
428  double neg_weight_sum = 0.0;
429  for (unsigned j=0; j<=ELEMENT_DIM; j++)
430  {
431  if (weight[j] < 0.0)
432  {
433  neg_weight_sum += weight[j];
434  }
435  }
436  if (neg_weight_sum > max_min_weight)
437  {
438  max_min_weight = neg_weight_sum;
439  closest_index = i;
440  }
441  }
442  assert(!this->mElements[closest_index]->IsDeleted());
443  return closest_index;
444 }
445 
446 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
448 {
449  std::vector<unsigned> element_indices;
450  for (unsigned i=0; i<this->mElements.size(); i++)
451  {
452  if (this->mElements[i]->IncludesPoint(rTestPoint))
453  {
454  assert(!this->mElements[i]->IsDeleted());
455  element_indices.push_back(i);
456  }
457  }
458  return element_indices;
459 }
460 
461 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
463 {
464  // Three loops, just like the destructor. note we don't delete boundary nodes.
465  for (unsigned i=0; i<this->mBoundaryElements.size(); i++)
466  {
467  delete this->mBoundaryElements[i];
468  }
469  for (unsigned i=0; i<this->mElements.size(); i++)
470  {
471  delete this->mElements[i];
472  }
473  for (unsigned i=0; i<this->mNodes.size(); i++)
474  {
475  delete this->mNodes[i];
476  }
477 
478  this->mNodes.clear();
479  this->mElements.clear();
480  this->mBoundaryElements.clear();
481  this->mBoundaryNodes.clear();
482 }
483 
484 
485 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
487 {
488  assert(SPACE_DIM == 2); // LCOV_EXCL_LINE
489  assert(SPACE_DIM == ELEMENT_DIM); // LCOV_EXCL_LINE
490 
491  double x_difference = this->mNodes[indexB]->rGetLocation()[0] - this->mNodes[indexA]->rGetLocation()[0];
492  double y_difference = this->mNodes[indexB]->rGetLocation()[1] - this->mNodes[indexA]->rGetLocation()[1];
493 
494  if (x_difference == 0)
495  {
496  if (y_difference > 0)
497  {
498  return M_PI/2.0;
499  }
500  else if (y_difference < 0)
501  {
502  return -M_PI/2.0;
503  }
504  else
505  {
506  EXCEPTION("Tried to compute polar angle of (0,0)");
507  }
508  }
509 
510  double angle = atan2(y_difference,x_difference);
511  return angle;
512 }
513 
515 // Edge iterator class //
517 
518 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
520 {
521  assert((*this) != mrMesh.EdgesEnd());
523  return p_element->GetNode(mNodeALocalIndex);
524 }
525 
526 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
528 {
529  assert((*this) != mrMesh.EdgesEnd());
531  return p_element->GetNode(mNodeBLocalIndex);
532 }
533 
534 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
536 {
537  return (mElemIndex != rOther.mElemIndex ||
540 }
541 
542 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
544 {
545  bool already_seen_this_edge;
546 
547  unsigned num_elements = mrMesh.GetNumAllElements();
548  std::pair<unsigned, unsigned> current_node_pair;
549  do
550  {
551  /*
552  * Advance to the next edge in the mesh.
553  * Node indices are incremented modulo #nodes_per_elem.
554  */
555  mNodeBLocalIndex = (mNodeBLocalIndex + 1) % (ELEMENT_DIM+1);
557  {
558  mNodeALocalIndex = (mNodeALocalIndex + 1) % (ELEMENT_DIM+1);
559  mNodeBLocalIndex = (mNodeALocalIndex + 1) % (ELEMENT_DIM+1);
560  }
561 
562  if (mNodeALocalIndex == 0 && mNodeBLocalIndex == 1) // advance to next element...
563  {
564  // ...skipping deleted ones
565  do
566  {
567  mElemIndex++;
568  }
569  while (mElemIndex!=num_elements && mrMesh.GetElement(mElemIndex)->IsDeleted());
570  }
571 
572  if (mElemIndex != num_elements)
573  {
575  unsigned node_a_global_index = p_element->GetNodeGlobalIndex(mNodeALocalIndex);
576  unsigned node_b_global_index = p_element->GetNodeGlobalIndex(mNodeBLocalIndex);
577  if (node_b_global_index < node_a_global_index)
578  {
579  // Swap them over
580  unsigned temp = node_a_global_index;
581  node_a_global_index = node_b_global_index;
582  node_b_global_index = temp;
583  }
584 
585  // Check we haven't seen it before
586  current_node_pair = std::pair<unsigned, unsigned>(node_a_global_index, node_b_global_index);
587  already_seen_this_edge = (mEdgesVisited.count(current_node_pair) != 0);
588  }
589  else
590  {
591  already_seen_this_edge = false;
592  }
593  }
594 
595  while (already_seen_this_edge);
596  mEdgesVisited.insert(current_node_pair);
597 
598  return (*this);
599 }
600 
601 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
603  : mrMesh(rMesh),
604  mElemIndex(elemIndex),
605  mNodeALocalIndex(0),
607 {
608  if (elemIndex == mrMesh.GetNumAllElements())
609  {
610  return;
611  }
612 
613  mEdgesVisited.clear();
614 
615  // Add the current node pair to the store
616  unsigned node_a_global_index = mrMesh.GetElement(mElemIndex)->GetNodeGlobalIndex(mNodeALocalIndex);
617  unsigned node_b_global_index = mrMesh.GetElement(mElemIndex)->GetNodeGlobalIndex(mNodeBLocalIndex);
618  if (node_b_global_index < node_a_global_index)
619  {
620  // Swap them over
621  unsigned temp = node_a_global_index;
622  node_a_global_index = node_b_global_index;
623  node_b_global_index = temp;
624  }
625 
626  // Check we haven't seen it before
627  std::pair<unsigned, unsigned> current_node_pair = std::pair<unsigned, unsigned>(node_a_global_index, node_b_global_index);
628  mEdgesVisited.insert(current_node_pair);
629 }
630 
631 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
633 {
634  unsigned first_element_index = 0;
635  while (first_element_index!=this->GetNumAllElements() && this->GetElement(first_element_index)->IsDeleted())
636  {
637  first_element_index++;
638  }
639  return EdgeIterator(*this, first_element_index);
640 }
641 
642 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
644 {
645  return EdgeIterator(*this, this->GetNumAllElements());
646 }
647 
648 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
650 {
652 }
653 
654 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
656 {
657  assert(index < this->mNodes.size() );
658  return index;
659 }
660 
661 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
663 {
664  assert(index < this->mElements.size() );
665  return index;
666 }
667 
668 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
670 {
671  assert(index < this->mBoundaryElements.size() );
672  return index;
673 }
674 
675 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
677 {
678  unsigned num_elements = this->GetNumAllElements();
679  unsigned num_boundary_elements = this->GetNumAllBoundaryElements();
680 
681  // Make sure we have enough space
682  this->mElementJacobians.resize(num_elements);
683  this->mElementInverseJacobians.resize(num_elements);
684 
685  if (ELEMENT_DIM < SPACE_DIM)
686  {
687  this->mElementWeightedDirections.resize(num_elements);
688  }
689 
690  this->mBoundaryElementWeightedDirections.resize(num_boundary_elements);
691 
692  this->mElementJacobianDeterminants.resize(num_elements);
693  this->mBoundaryElementJacobianDeterminants.resize(num_boundary_elements);
694 
695  // Update caches
697  iter != this->GetElementIteratorEnd();
698  ++iter)
699  {
700  unsigned index = iter->GetIndex();
701  iter->CalculateInverseJacobian(this->mElementJacobians[index], this->mElementJacobianDeterminants[index], this->mElementInverseJacobians[index]);
702  }
703 
704  if (ELEMENT_DIM < SPACE_DIM)
705  {
707  iter != this->GetElementIteratorEnd();
708  ++iter)
709  {
710  unsigned index = iter->GetIndex();
711  iter->CalculateWeightedDirection(this->mElementWeightedDirections[index], this->mElementJacobianDeterminants[index]);
712  }
713  }
714 
716  itb != this->GetBoundaryElementIteratorEnd();
717  itb++)
718  {
719  unsigned index = (*itb)->GetIndex();
720  (*itb)->CalculateWeightedDirection(this->mBoundaryElementWeightedDirections[index], this->mBoundaryElementJacobianDeterminants[index]);
721  }
722 }
723 
724 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
725 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetJacobianForElement(unsigned elementIndex, c_matrix<double, SPACE_DIM, SPACE_DIM>& rJacobian, double& rJacobianDeterminant) const
726 {
727  assert(ELEMENT_DIM <= SPACE_DIM);
728  assert(elementIndex < this->mElementJacobians.size());
729  rJacobian = this->mElementJacobians[elementIndex];
730  rJacobianDeterminant = this->mElementJacobianDeterminants[elementIndex];
731 }
732 
733 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
734 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetInverseJacobianForElement(unsigned elementIndex, c_matrix<double, SPACE_DIM, ELEMENT_DIM>& rJacobian, double& rJacobianDeterminant, c_matrix<double, ELEMENT_DIM, SPACE_DIM>& rInverseJacobian) const
735 {
736  assert(ELEMENT_DIM <= SPACE_DIM); // LCOV_EXCL_LINE
737  assert(elementIndex < this->mElementInverseJacobians.size());
738  rInverseJacobian = this->mElementInverseJacobians[elementIndex];
739  rJacobian = this->mElementJacobians[elementIndex];
740  rJacobianDeterminant = this->mElementJacobianDeterminants[elementIndex];
741 }
742 
743 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
744 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetWeightedDirectionForElement(unsigned elementIndex, c_vector<double, SPACE_DIM>& rWeightedDirection, double& rJacobianDeterminant) const
745 {
746  assert(ELEMENT_DIM < SPACE_DIM); // LCOV_EXCL_LINE
747  assert(elementIndex < this->mElementWeightedDirections.size());
748  rWeightedDirection = this->mElementWeightedDirections[elementIndex];
749  rJacobianDeterminant = this->mElementJacobianDeterminants[elementIndex];
750 }
751 
752 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
753 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetWeightedDirectionForBoundaryElement(unsigned elementIndex, c_vector<double, SPACE_DIM>& rWeightedDirection, double& rJacobianDeterminant) const
754 {
755  assert(elementIndex < this->mBoundaryElementWeightedDirections.size());
756  rWeightedDirection = this->mBoundaryElementWeightedDirections[elementIndex];
757  rJacobianDeterminant = this->mBoundaryElementJacobianDeterminants[elementIndex];
758 }
759 
760 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
762 {
763  mesherIo.numberofpoints = 0;
764  mesherIo.pointlist = nullptr;
765  mesherIo.numberofpointattributes = 0;
766  mesherIo.pointattributelist = (double *) nullptr;
767  mesherIo.pointmarkerlist = (int *) nullptr;
768  mesherIo.numberofsegments = 0;
769  mesherIo.numberofholes = 0;
770  mesherIo.numberofregions = 0;
771  mesherIo.trianglelist = (int *) nullptr;
772  mesherIo.triangleattributelist = (double *) nullptr;
773  mesherIo.numberoftriangleattributes = 0;
774  mesherIo.edgelist = (int *) nullptr;
775  mesherIo.edgemarkerlist = (int *) nullptr;
776 }
777 
778 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
780 {
781  if (mesherIo.numberofpoints != 0)
782  {
783  mesherIo.numberofpoints=0;
784  free(mesherIo.pointlist);
785  }
786 
787  // These (and the above) should actually be safe since we explicity set to NULL above
788  free(mesherIo.pointattributelist);
789  free(mesherIo.pointmarkerlist);
790  free(mesherIo.trianglelist);
791  free(mesherIo.triangleattributelist);
792  free(mesherIo.edgelist);
793  free(mesherIo.edgemarkerlist);
794 }
795 
796 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
797 template <class MESHER_IO>
798 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ExportToMesher(NodeMap& map, MESHER_IO& mesherInput, int *elementList)
799 {
800  if (SPACE_DIM == 2)
801  {
802  mesherInput.pointlist = (double *) malloc(this->GetNumNodes() * SPACE_DIM * sizeof(double));
803  }
804  else
805  {
806  mesherInput.pointlist = new double[this->GetNumNodes() * SPACE_DIM];
807  }
808 
809  mesherInput.numberofpoints = this->GetNumNodes();
810  unsigned new_index = 0;
811  for (unsigned i=0; i<this->GetNumAllNodes(); i++)
812  {
813  if (this->mNodes[i]->IsDeleted())
814  {
815  map.SetDeleted(i);
816  }
817  else
818  {
819  map.SetNewIndex(i, new_index);
820  for (unsigned j=0; j<SPACE_DIM; j++)
821  {
822  mesherInput.pointlist[SPACE_DIM*new_index + j] = this->mNodes[i]->rGetLocation()[j];
823  }
824  new_index++;
825  }
826  }
827  if (elementList != nullptr)
828  {
829  unsigned element_index = 0;
830 
831  // Assume there is enough space for this
832  mesherInput.numberofcorners=ELEMENT_DIM+1;
834  elem_iter != this->GetElementIteratorEnd();
835  ++elem_iter)
836  {
837 
838  for (unsigned j=0; j<=ELEMENT_DIM; j++)
839  {
840  elementList[element_index*(ELEMENT_DIM+1) + j] = (*elem_iter).GetNodeGlobalIndex(j);
841  }
842  element_index++;
843  }
844  }
845 }
846 
847 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
848 template <class MESHER_IO>
849 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ImportFromMesher(MESHER_IO& mesherOutput, unsigned numberOfElements, int *elementList, unsigned numberOfFaces, int *faceList, int *edgeMarkerList)
850 {
851  unsigned nodes_per_element = mesherOutput.numberofcorners;
852 
853  assert( nodes_per_element == ELEMENT_DIM+1 || nodes_per_element == (ELEMENT_DIM+1)*(ELEMENT_DIM+2)/2 );
854 
855  Clear();
856 
857  // Construct the nodes
858  for (unsigned node_index=0; node_index<(unsigned)mesherOutput.numberofpoints; node_index++)
859  {
860  this->mNodes.push_back(new Node<SPACE_DIM>(node_index, &mesherOutput.pointlist[node_index * SPACE_DIM], false));
861  }
862 
863  // Construct the elements
864  this->mElements.reserve(numberOfElements);
865 
866  unsigned real_element_index=0;
867  for (unsigned element_index=0; element_index<numberOfElements; element_index++)
868  {
869  std::vector<Node<SPACE_DIM>*> nodes;
870  for (unsigned j=0; j<ELEMENT_DIM+1; j++)
871  {
872  unsigned global_node_index = elementList[element_index*(nodes_per_element) + j];
873  assert(global_node_index < this->mNodes.size());
874  nodes.push_back(this->mNodes[global_node_index]);
875 
876  }
877 
878  /*
879  * For some reason, tetgen in library mode makes its initial Delaunay mesh
880  * with very thin slivers. Hence we expect to ignore some of the elements!
881  */
883  try
884  {
885  p_element = new Element<ELEMENT_DIM, SPACE_DIM>(real_element_index, nodes);
886 
887  // Shouldn't throw after this point
888  this->mElements.push_back(p_element);
889 
890  // Add the internals to quadratics
891  for (unsigned j=ELEMENT_DIM+1; j<nodes_per_element; j++)
892  {
893  unsigned global_node_index = elementList[element_index*nodes_per_element + j];
894  assert(global_node_index < this->mNodes.size());
895  this->mElements[real_element_index]->AddNode( this->mNodes[global_node_index] );
896  this->mNodes[global_node_index]->AddElement(real_element_index);
897  this->mNodes[global_node_index]->MarkAsInternal();
898  }
899  real_element_index++;
900  }
901  catch (Exception &)
902  {
903  if (SPACE_DIM == 2)
904  {
905  WARNING("Triangle has produced a zero area (collinear) element");
906  }
907  else
908  {
909  WARNING("Tetgen has produced a zero volume (coplanar) element");
910  }
911  }
912  }
913 
914  // Construct the BoundaryElements (and mark boundary nodes)
915  unsigned next_boundary_element_index = 0;
916  for (unsigned boundary_element_index=0; boundary_element_index<numberOfFaces; boundary_element_index++)
917  {
918  /*
919  * Tetgen produces only boundary faces (set edgeMarkerList to NULL).
920  * Triangle marks which edges are on the boundary.
921  */
922  if (edgeMarkerList == nullptr || edgeMarkerList[boundary_element_index] == 1)
923  {
924  std::vector<Node<SPACE_DIM>*> nodes;
925  for (unsigned j=0; j<ELEMENT_DIM; j++)
926  {
927  unsigned global_node_index = faceList[boundary_element_index*ELEMENT_DIM + j];
928  assert(global_node_index < this->mNodes.size());
929  nodes.push_back(this->mNodes[global_node_index]);
930  if (!nodes[j]->IsBoundaryNode())
931  {
932  nodes[j]->SetAsBoundaryNode();
933  this->mBoundaryNodes.push_back(nodes[j]);
934  }
935  }
936 
937  /*
938  * For some reason, tetgen in library mode makes its initial Delaunay mesh
939  * with very thin slivers. Hence we expect to ignore some of the elements!
940  */
941  BoundaryElement<ELEMENT_DIM-1, SPACE_DIM>* p_b_element;
942  try
943  {
944  p_b_element = new BoundaryElement<ELEMENT_DIM-1, SPACE_DIM>(next_boundary_element_index, nodes);
945  this->mBoundaryElements.push_back(p_b_element);
946  next_boundary_element_index++;
947  }
948  // LCOV_EXCL_START
949  catch (Exception &)
950  {
951  // Tetgen is feeding us lies
952  /*
953  * Note: this code is covered in profiling (Test3dOffLatticeRepresentativeSimulation).
954  * It's hard to replicate Tetgen's behaviour with a unit test.
955  */
956  assert(SPACE_DIM == 3);
957  }
958  // LCOV_EXCL_STOP
959  }
960  }
961 
963 }
964 
965 // Explicit instantiation
966 template class TetrahedralMesh<1,1>;
967 template class TetrahedralMesh<1,2>;
968 template class TetrahedralMesh<1,3>;
969 template class TetrahedralMesh<2,2>;
970 template class TetrahedralMesh<2,3>;
971 template class TetrahedralMesh<3,3>;
972 
977 template void TetrahedralMesh<2,2>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
978 template void TetrahedralMesh<2,2>::ImportFromMesher<triangulateio>(triangulateio&, unsigned, int *, unsigned, int *, int *);
979 
980 template void TetrahedralMesh<3,3>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
981 template void TetrahedralMesh<3,3>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned, int *, unsigned, int *, int *);
982 
983 //The following don't ever need to be instantiated, but are needed to keep some compilers happy
984 template void TetrahedralMesh<1,2>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
985 template void TetrahedralMesh<1,2>::ImportFromMesher<triangulateio>(triangulateio&, unsigned, int *, unsigned, int *, int *);
986 
987 template void TetrahedralMesh<1,3>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
988 template void TetrahedralMesh<1,3>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned, int *, unsigned, int *, int *);
989 template void TetrahedralMesh<2,3>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
990 template void TetrahedralMesh<2,3>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned, int *, unsigned, int *, int *);
991 
992 //Intel compilation with IPO thinks that it's missing some bizarre instantiations
993 template void TetrahedralMesh<3u, 3u>::ImportFromMesher<triangulateio>(triangulateio&, unsigned int, int*, unsigned int, int*, int*);
994 template void TetrahedralMesh<1u, 1u>::ImportFromMesher<triangulateio>(triangulateio&, unsigned int, int*, unsigned int, int*, int*);
995 template void TetrahedralMesh<1u, 1u>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned int, int*, unsigned int, int*, int*);
996 template void TetrahedralMesh<2u, 2u>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned int, int*, unsigned int, int*, int*);
997 template void TetrahedralMesh<1u, 1u>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
998 
999 // Compiler on ARC cluser HAL requires the following
1000 template void TetrahedralMesh<3u, 3u>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1001 template void TetrahedralMesh<1u, 1u>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1002 template void TetrahedralMesh<2u, 2u>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
1003 
1004 // Intel v11 compilation thinks that it's missing even more bizarre instantiations
1005 //template void TetrahedralMesh<2,2>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
1006 //template void TetrahedralMesh<3,3>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1007 //template void TetrahedralMesh<1,3>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1008 //template void TetrahedralMesh<1,1>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1009 //template void TetrahedralMesh<1,2>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
1010 //template void TetrahedralMesh<2,3>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1011 //template void TetrahedralMesh<1,3>::ImportFromMesher<triangulateio>(triangulateio&, unsigned, int *, unsigned, int *, int *);
1012 //template void TetrahedralMesh<2,3>::ImportFromMesher<triangulateio>(triangulateio&, unsigned, int *, unsigned, int *, int *);
1013 //template void TetrahedralMesh<1,2>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned, int *, unsigned, int *, int *);
1014 #ifdef _MSC_VER
1015 //MSVC requires these
1016 template void TetrahedralMesh<2,2>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
1017 template void TetrahedralMesh<3,3>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1018 template void TetrahedralMesh<1,3>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1019 template void TetrahedralMesh<1,1>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1020 template void TetrahedralMesh<1,2>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
1021 template void TetrahedralMesh<2,3>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1022 template void TetrahedralMesh<1,3>::ImportFromMesher<triangulateio>(triangulateio&, unsigned, int *, unsigned, int *, int *);
1023 template void TetrahedralMesh<2,3>::ImportFromMesher<triangulateio>(triangulateio&, unsigned, int *, unsigned, int *, int *);
1024 template void TetrahedralMesh<1,2>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned, int *, unsigned, int *, int *);
1025 #endif
1026 
1031 // Serialization for Boost >= 1.36
1032 #define CHASTE_SERIALIZATION_CPP
std::vector< c_matrix< double, SPACE_DIM, ELEMENT_DIM > > mElementJacobians
ElementIterator GetElementIteratorBegin(bool skipDeletedElements=true)
double GetAngleBetweenNodes(unsigned indexA, unsigned indexB)
std::vector< Element< ELEMENT_DIM, SPACE_DIM > * > mElements
void SetAttribute(double attribute)
std::vector< Node< SPACE_DIM > * > mBoundaryNodes
void ImportFromMesher(MESHER_IO &mesherOutput, unsigned numberOfElements, int *elementList, unsigned numberOfFaces, int *faceList, int *edgeMarkerList)
virtual ElementData GetNextElementData()=0
Definition: Node.hpp:58
std::vector< unsigned > GetContainingElementIndices(const ChastePoint< SPACE_DIM > &rTestPoint)
BoundaryElementIterator GetBoundaryElementIteratorBegin() const
unsigned GetNodeGlobalIndex(unsigned localIndex) const
virtual void GetWeightedDirectionForBoundaryElement(unsigned elementIndex, c_vector< double, SPACE_DIM > &rWeightedDirection, double &rJacobianDeterminant) const
void AddNodeAttribute(double attribute)
Definition: Node.cpp:170
unsigned GetContainingElementIndexWithInitialGuess(const ChastePoint< SPACE_DIM > &rTestPoint, unsigned startingElementGuess, bool strict=false)
virtual void GetInverseJacobianForElement(unsigned elementIndex, c_matrix< double, SPACE_DIM, ELEMENT_DIM > &rJacobian, double &rJacobianDeterminant, c_matrix< double, ELEMENT_DIM, SPACE_DIM > &rInverseJacobian) const
#define EXCEPTION(message)
Definition: Exception.hpp:143
virtual void RefreshJacobianCachedData()
std::vector< c_matrix< double, ELEMENT_DIM, SPACE_DIM > > mElementInverseJacobians
Element< ELEMENT_DIM, SPACE_DIM > * GetElement(unsigned index) const
void ConstructFromMeshReader(AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &rMeshReader)
virtual unsigned GetNumNodes() const
virtual void GetJacobianForElement(unsigned elementIndex, c_matrix< double, SPACE_DIM, SPACE_DIM > &rJacobian, double &rJacobianDeterminant) const
virtual unsigned GetNumAllNodes() const
Node< SPACE_DIM > * GetNodeA()
void ExportToMesher(NodeMap &map, MESHER_IO &mesherInput, int *elementList=nullptr)
void Shuffle(std::vector< boost::shared_ptr< T > > &rValues)
virtual ElementData GetNextFaceData()=0
std::vector< unsigned > NodeIndices
void SetNewIndex(unsigned oldIndex, unsigned newIndex)
Definition: NodeMap.cpp:66
EdgeIterator(TetrahedralMesh &rMesh, unsigned elemIndex)
virtual bool HasNodePermutation()
Node< SPACE_DIM > * GetNodeB()
void ReadNodesPerProcessorFile(const std::string &rNodesPerProcessorFile)
unsigned SolveBoundaryElementMapping(unsigned index) const
std::vector< Node< SPACE_DIM > * > mNodes
virtual unsigned GetNumElements() const =0
void InitialiseTriangulateIo(triangulateio &mesherIo)
std::vector< BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > * > mBoundaryElements
unsigned SolveElementMapping(unsigned index) const
virtual void Reset()=0
void SetDeleted(unsigned index)
Definition: NodeMap.cpp:76
std::vector< double > mElementJacobianDeterminants
unsigned GetNearestElementIndex(const ChastePoint< SPACE_DIM > &rTestPoint)
std::vector< c_vector< double, SPACE_DIM > > mBoundaryElementWeightedDirections
#define EXCEPT_IF_NOT(test)
Definition: Exception.hpp:158
static RandomNumberGenerator * Instance()
#define EXPORT_TEMPLATE_CLASS_ALL_DIMS(CLASS)
virtual unsigned GetNumFaceAttributes() const
bool operator!=(const typename TetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::EdgeIterator &rOther)
virtual unsigned GetNumElementAttributes() const
virtual void Clear()
EdgeIterator EdgesBegin()
unsigned SolveNodeMapping(unsigned index) const
std::vector< double > mBoundaryElementJacobianDeterminants
unsigned GetNumAllBoundaryElements() const
virtual void GetWeightedDirectionForElement(unsigned elementIndex, c_vector< double, SPACE_DIM > &rWeightedDirection, double &rJacobianDeterminant) const
virtual std::string GetMeshFileBaseName()
static unsigned GetNumProcs()
Definition: PetscTools.cpp:108
virtual std::vector< double > GetNextNode()=0
std::vector< BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > * >::const_iterator BoundaryElementIterator
std::vector< c_vector< double, SPACE_DIM > > mElementWeightedDirections
static unsigned GetMyRank()
Definition: PetscTools.cpp:114
EdgeIterator EdgesEnd()
std::set< std::pair< unsigned, unsigned > > mEdgesVisited
gcov doesn&#39;t like this file...
virtual unsigned GetNumNodes() const =0
void FreeTriangulateIo(triangulateio &mesherIo)
virtual std::vector< double > GetNodeAttributes()
virtual unsigned GetNumFaces() const =0
BoundaryElementIterator GetBoundaryElementIteratorEnd() const