Chaste  Release::3.4
TetrahedralMesh.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 "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 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
237 {
238  /*
239  * Each face of each element is a set of node indices.
240  * We form a set of these in order to get their parity:
241  * all faces which appear once are inserted into the set;
242  * all faces which appear twice are inserted and then removed from the set;
243  * we're assuming that faces never appear more than twice.
244  */
245  std::set< std::set<unsigned> > odd_parity_faces;
246 
247  for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = this->GetElementIteratorBegin();
248  iter != this->GetElementIteratorEnd();
249  ++iter)
250  {
251  for (unsigned face_index=0; face_index<=ELEMENT_DIM; face_index++)
252  {
253  std::set<unsigned> face_info;
254  for (unsigned node_index=0; node_index<=ELEMENT_DIM; node_index++)
255  {
256  // Leave one index out each time
257  if (node_index != face_index)
258  {
259  face_info.insert(iter->GetNodeGlobalIndex(node_index));
260  }
261  }
262  // Face is now formed - attempt to find it
263  std::set< std::set<unsigned> >::iterator find_face=odd_parity_faces.find(face_info);
264  if (find_face != odd_parity_faces.end())
265  {
266  // Face was in set, so it now has even parity.
267  // Remove it via the iterator
268  odd_parity_faces.erase(find_face);
269  }
270  else
271  {
272  // Face is not in set so it now has odd parity. Insert it
273  odd_parity_faces.insert(face_info);
274  }
275 
276  }
277  }
278 
279  /*
280  * At this point the odd parity faces should be the same as the
281  * boundary elements. We could check this explicitly or we could
282  * just count them.
283  */
284  return(odd_parity_faces.size() == this->GetNumBoundaryElements());
285 }
286 
287 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
289 {
290  double mesh_volume = 0.0;
291 
292  for (typename AbstractTetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ElementIterator iter = this->GetElementIteratorBegin();
293  iter != this->GetElementIteratorEnd();
294  ++iter)
295  {
296  mesh_volume += iter->GetVolume(mElementJacobianDeterminants[iter->GetIndex()]);
297  }
298 
299  return mesh_volume;
300 }
301 
302 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
304 {
305  // ELEMENT_DIM-1 is the dimension of the boundary element
306  assert(ELEMENT_DIM >= 1);
307  const unsigned bound_element_dim = ELEMENT_DIM-1;
308  assert(bound_element_dim < 3);
309  if ( bound_element_dim == 0)
310  {
311  return 0.0;
312  }
313 
314  double mesh_surface = 0.0;
315  typename TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::BoundaryElementIterator it = this->GetBoundaryElementIteratorBegin();
316 
317  while (it != this->GetBoundaryElementIteratorEnd())
318  {
319  mesh_surface += mBoundaryElementJacobianDeterminants[(*it)->GetIndex()];
320  it++;
321  }
322 
323  if ( bound_element_dim == 2)
324  {
325  mesh_surface /= 2.0;
326  }
327 
328  return mesh_surface;
329 }
330 
331 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
333 {
334  // Make a permutation vector of the identity
336  std::vector<unsigned> perm;
337  p_rng->Shuffle(this->mNodes.size(), perm);
338 
339  // Call the non-random version
340  PermuteNodes(perm);
341 }
342 
343 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
344 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::PermuteNodes(const std::vector<unsigned>& perm)
345 {
346  // Let's not do this if there are any deleted nodes
347  assert( this->GetNumAllNodes() == this->GetNumNodes());
348 
349  assert(perm.size() == this->mNodes.size());
350 
351  // Copy the node pointers
352  std::vector< Node<SPACE_DIM>* > copy_m_nodes;
353  copy_m_nodes.assign(this->mNodes.begin(), this->mNodes.end());
354 
355  for (unsigned original_index=0; original_index<this->mNodes.size(); original_index++)
356  {
357  assert(perm[original_index] < this->mNodes.size());
358  //perm[original_index] holds the new assigned index of that node
359  this->mNodes[ perm[original_index] ] = copy_m_nodes[original_index];
360  }
361 
362  // Update indices
363  for (unsigned index=0; index<this->mNodes.size(); index++)
364  {
365  this->mNodes[index]->SetIndex(index);
366  }
367 
368  // Copy the permutation vector into the mesh
369  this->mNodePermutation = perm;
370 }
371 
372 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
373 unsigned TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetContainingElementIndexWithInitialGuess(const ChastePoint<SPACE_DIM>& rTestPoint, unsigned startingElementGuess, bool strict)
374 {
375  assert(startingElementGuess<this->GetNumElements());
376 
377  /*
378  * Let m=startingElementGuess, N=num_elem-1.
379  * We search from in this order: m, m+1, m+2, .. , N, 0, 1, .., m-1.
380  */
381  unsigned i = startingElementGuess;
382  bool reached_end = false;
383 
384  while (!reached_end)
385  {
386  if (this->mElements[i]->IncludesPoint(rTestPoint, strict))
387  {
388  assert(!this->mElements[i]->IsDeleted());
389  return i;
390  }
391 
392  // Increment
393  i++;
394  if (i==this->GetNumElements())
395  {
396  i=0;
397  }
398 
399  // Back to the beginning yet?
400  if (i==startingElementGuess)
401  {
402  reached_end = true;
403  }
404  }
405 
406  // If it's in none of the elements, then throw
407  std::stringstream ss;
408  ss << "Point [";
409  for (unsigned j=0; (int)j<(int)SPACE_DIM-1; j++)
410  {
411  ss << rTestPoint[j] << ",";
412  }
413  ss << rTestPoint[SPACE_DIM-1] << "] is not in mesh - all elements tested";
414  EXCEPTION(ss.str());
415 }
416 
417 
418 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
420 {
421  double max_min_weight = -std::numeric_limits<double>::infinity();
422  unsigned closest_index = 0;
423  for (unsigned i=0; i<this->mElements.size(); i++)
424  {
425  c_vector<double, ELEMENT_DIM+1> weight=this->mElements[i]->CalculateInterpolationWeights(rTestPoint);
426  double neg_weight_sum=0.0;
427  for (unsigned j=0; j<=ELEMENT_DIM; j++)
428  {
429  if (weight[j] < 0.0)
430  {
431  neg_weight_sum += weight[j];
432  }
433  }
434  if (neg_weight_sum > max_min_weight)
435  {
436  max_min_weight = neg_weight_sum;
437  closest_index = i;
438  }
439  }
440  assert(!this->mElements[closest_index]->IsDeleted());
441  return closest_index;
442 }
443 
444 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
446 {
447  std::vector<unsigned> element_indices;
448  for (unsigned i=0; i<this->mElements.size(); i++)
449  {
450  if (this->mElements[i]->IncludesPoint(rTestPoint))
451  {
452  assert(!this->mElements[i]->IsDeleted());
453  element_indices.push_back(i);
454  }
455  }
456  return element_indices;
457 }
458 
459 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
461 {
462  // Three loops, just like the destructor. note we don't delete boundary nodes.
463  for (unsigned i=0; i<this->mBoundaryElements.size(); i++)
464  {
465  delete this->mBoundaryElements[i];
466  }
467  for (unsigned i=0; i<this->mElements.size(); i++)
468  {
469  delete this->mElements[i];
470  }
471  for (unsigned i=0; i<this->mNodes.size(); i++)
472  {
473  delete this->mNodes[i];
474  }
475 
476  this->mNodes.clear();
477  this->mElements.clear();
478  this->mBoundaryElements.clear();
479  this->mBoundaryNodes.clear();
480 }
481 
482 
483 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
485 {
486  assert(SPACE_DIM == 2);
487  assert(SPACE_DIM == ELEMENT_DIM);
488 
489  double x_difference = this->mNodes[indexB]->rGetLocation()[0] - this->mNodes[indexA]->rGetLocation()[0];
490  double y_difference = this->mNodes[indexB]->rGetLocation()[1] - this->mNodes[indexA]->rGetLocation()[1];
491 
492  if (x_difference == 0)
493  {
494  if (y_difference > 0)
495  {
496  return M_PI/2.0;
497  }
498  else if (y_difference < 0)
499  {
500  return -M_PI/2.0;
501  }
502  else
503  {
504  EXCEPTION("Tried to compute polar angle of (0,0)");
505  }
506  }
507 
508  double angle = atan2(y_difference,x_difference);
509  return angle;
510 }
511 
513 // Edge iterator class //
515 
516 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
518 {
519  assert((*this) != mrMesh.EdgesEnd());
521  return p_element->GetNode(mNodeALocalIndex);
522 }
523 
524 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
526 {
527  assert((*this) != mrMesh.EdgesEnd());
528  Element<ELEMENT_DIM,SPACE_DIM>* p_element = mrMesh.GetElement(mElemIndex);
529  return p_element->GetNode(mNodeBLocalIndex);
530 }
531 
532 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
534 {
535  return (mElemIndex != rOther.mElemIndex ||
536  mNodeALocalIndex != rOther.mNodeALocalIndex ||
537  mNodeBLocalIndex != rOther.mNodeBLocalIndex);
538 }
539 
540 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
542 {
543  bool already_seen_this_edge;
544 
545  unsigned num_elements = mrMesh.GetNumAllElements();
546  std::pair<unsigned, unsigned> current_node_pair;
547  do
548  {
549  /*
550  * Advance to the next edge in the mesh.
551  * Node indices are incremented modulo #nodes_per_elem.
552  */
553  mNodeBLocalIndex = (mNodeBLocalIndex + 1) % (ELEMENT_DIM+1);
554  if (mNodeBLocalIndex == mNodeALocalIndex)
555  {
556  mNodeALocalIndex = (mNodeALocalIndex + 1) % (ELEMENT_DIM+1);
557  mNodeBLocalIndex = (mNodeALocalIndex + 1) % (ELEMENT_DIM+1);
558  }
559 
560  if (mNodeALocalIndex == 0 && mNodeBLocalIndex == 1) // advance to next element...
561  {
562  // ...skipping deleted ones
563  do
564  {
565  mElemIndex++;
566  }
567  while (mElemIndex!=num_elements && mrMesh.GetElement(mElemIndex)->IsDeleted());
568  }
569 
570  if (mElemIndex != num_elements)
571  {
572  Element<ELEMENT_DIM, SPACE_DIM>* p_element = mrMesh.GetElement(mElemIndex);
573  unsigned node_a_global_index = p_element->GetNodeGlobalIndex(mNodeALocalIndex);
574  unsigned node_b_global_index = p_element->GetNodeGlobalIndex(mNodeBLocalIndex);
575  if (node_b_global_index < node_a_global_index)
576  {
577  // Swap them over
578  unsigned temp = node_a_global_index;
579  node_a_global_index = node_b_global_index;
580  node_b_global_index = temp;
581  }
582 
583  // Check we haven't seen it before
584  current_node_pair = std::pair<unsigned, unsigned>(node_a_global_index, node_b_global_index);
585  already_seen_this_edge = (mEdgesVisited.count(current_node_pair) != 0);
586  }
587  else
588  {
589  already_seen_this_edge = false;
590  }
591  }
592 
593  while (already_seen_this_edge);
594  mEdgesVisited.insert(current_node_pair);
595 
596  return (*this);
597 }
598 
599 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
601  : mrMesh(rMesh),
602  mElemIndex(elemIndex),
603  mNodeALocalIndex(0),
604  mNodeBLocalIndex(1)
605 {
606  if (elemIndex == mrMesh.GetNumAllElements())
607  {
608  return;
609  }
610 
611  mEdgesVisited.clear();
612 
613  // Add the current node pair to the store
614  unsigned node_a_global_index = mrMesh.GetElement(mElemIndex)->GetNodeGlobalIndex(mNodeALocalIndex);
615  unsigned node_b_global_index = mrMesh.GetElement(mElemIndex)->GetNodeGlobalIndex(mNodeBLocalIndex);
616  if (node_b_global_index < node_a_global_index)
617  {
618  // Swap them over
619  unsigned temp = node_a_global_index;
620  node_a_global_index = node_b_global_index;
621  node_b_global_index = temp;
622  }
623 
624  // Check we haven't seen it before
625  std::pair<unsigned, unsigned> current_node_pair = std::pair<unsigned, unsigned>(node_a_global_index, node_b_global_index);
626  mEdgesVisited.insert(current_node_pair);
627 }
628 
629 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
631 {
632  unsigned first_element_index = 0;
633  while (first_element_index!=this->GetNumAllElements() && this->GetElement(first_element_index)->IsDeleted())
634  {
635  first_element_index++;
636  }
637  return EdgeIterator(*this, first_element_index);
638 }
639 
640 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
642 {
643  return EdgeIterator(*this, this->GetNumAllElements());
644 }
645 
646 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
648 {
650 }
651 
652 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
654 {
655  assert(index < this->mNodes.size() );
656  return index;
657 }
658 
659 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
661 {
662  assert(index < this->mElements.size() );
663  return index;
664 }
665 
666 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
668 {
669  assert(index < this->mBoundaryElements.size() );
670  return index;
671 }
672 
673 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
675 {
676  unsigned num_elements = this->GetNumAllElements();
677  unsigned num_boundary_elements = this->GetNumAllBoundaryElements();
678 
679  // Make sure we have enough space
680  this->mElementJacobians.resize(num_elements);
681  this->mElementInverseJacobians.resize(num_elements);
682 
683  if (ELEMENT_DIM < SPACE_DIM)
684  {
685  this->mElementWeightedDirections.resize(num_elements);
686  }
687 
688  this->mBoundaryElementWeightedDirections.resize(num_boundary_elements);
689 
690  this->mElementJacobianDeterminants.resize(num_elements);
691  this->mBoundaryElementJacobianDeterminants.resize(num_boundary_elements);
692 
693  // Update caches
695  iter != this->GetElementIteratorEnd();
696  ++iter)
697  {
698  unsigned index = iter->GetIndex();
699  iter->CalculateInverseJacobian(this->mElementJacobians[index], this->mElementJacobianDeterminants[index], this->mElementInverseJacobians[index]);
700  }
701 
702  if (ELEMENT_DIM < SPACE_DIM)
703  {
705  iter != this->GetElementIteratorEnd();
706  ++iter)
707  {
708  unsigned index = iter->GetIndex();
709  iter->CalculateWeightedDirection(this->mElementWeightedDirections[index], this->mElementJacobianDeterminants[index]);
710  }
711  }
712 
714  itb != this->GetBoundaryElementIteratorEnd();
715  itb++)
716  {
717  unsigned index = (*itb)->GetIndex();
718  (*itb)->CalculateWeightedDirection(this->mBoundaryElementWeightedDirections[index], this->mBoundaryElementJacobianDeterminants[index]);
719  }
720 }
721 
722 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
723 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetJacobianForElement(unsigned elementIndex, c_matrix<double, SPACE_DIM, SPACE_DIM>& rJacobian, double& rJacobianDeterminant) const
724 {
725  assert(ELEMENT_DIM <= SPACE_DIM);
726  assert(elementIndex < this->mElementJacobians.size());
727  rJacobian = this->mElementJacobians[elementIndex];
728  rJacobianDeterminant = this->mElementJacobianDeterminants[elementIndex];
729 }
730 
731 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
732 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
733 {
734  assert(ELEMENT_DIM <= SPACE_DIM);
735  assert(elementIndex < this->mElementInverseJacobians.size());
736  rInverseJacobian = this->mElementInverseJacobians[elementIndex];
737  rJacobian = this->mElementJacobians[elementIndex];
738  rJacobianDeterminant = this->mElementJacobianDeterminants[elementIndex];
739 }
740 
741 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
742 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetWeightedDirectionForElement(unsigned elementIndex, c_vector<double, SPACE_DIM>& rWeightedDirection, double& rJacobianDeterminant) const
743 {
744  assert(ELEMENT_DIM < SPACE_DIM);
745  assert(elementIndex < this->mElementWeightedDirections.size());
746  rWeightedDirection = this->mElementWeightedDirections[elementIndex];
747  rJacobianDeterminant = this->mElementJacobianDeterminants[elementIndex];
748 }
749 
750 template<unsigned ELEMENT_DIM, unsigned SPACE_DIM>
751 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::GetWeightedDirectionForBoundaryElement(unsigned elementIndex, c_vector<double, SPACE_DIM>& rWeightedDirection, double& rJacobianDeterminant) const
752 {
753  assert(elementIndex < this->mBoundaryElementWeightedDirections.size());
754  rWeightedDirection = this->mBoundaryElementWeightedDirections[elementIndex];
755  rJacobianDeterminant = this->mBoundaryElementJacobianDeterminants[elementIndex];
756 }
757 
758 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
760 {
761  mesherIo.numberofpoints = 0;
762  mesherIo.pointlist = NULL;
763  mesherIo.numberofpointattributes = 0;
764  mesherIo.pointattributelist = (double *) NULL;
765  mesherIo.pointmarkerlist = (int *) NULL;
766  mesherIo.numberofsegments = 0;
767  mesherIo.numberofholes = 0;
768  mesherIo.numberofregions = 0;
769  mesherIo.trianglelist = (int *) NULL;
770  mesherIo.triangleattributelist = (double *) NULL;
771  mesherIo.numberoftriangleattributes = 0;
772  mesherIo.edgelist = (int *) NULL;
773  mesherIo.edgemarkerlist = (int *) NULL;
774 }
775 
776 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
778 {
779  if (mesherIo.numberofpoints != 0)
780  {
781  mesherIo.numberofpoints=0;
782  free(mesherIo.pointlist);
783  }
784 
785  // These (and the above) should actually be safe since we explicity set to NULL above
786  free(mesherIo.pointattributelist);
787  free(mesherIo.pointmarkerlist);
788  free(mesherIo.trianglelist);
789  free(mesherIo.triangleattributelist);
790  free(mesherIo.edgelist);
791  free(mesherIo.edgemarkerlist);
792 }
793 
794 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
795 template <class MESHER_IO>
796 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ExportToMesher(NodeMap& map, MESHER_IO& mesherInput, int *elementList)
797 {
798  if (SPACE_DIM == 2)
799  {
800  mesherInput.pointlist = (double *) malloc(this->GetNumNodes() * SPACE_DIM * sizeof(double));
801  }
802  else
803  {
804  mesherInput.pointlist = new double[this->GetNumNodes() * SPACE_DIM];
805  }
806 
807  mesherInput.numberofpoints = this->GetNumNodes();
808  unsigned new_index = 0;
809  for (unsigned i=0; i<this->GetNumAllNodes(); i++)
810  {
811  if (this->mNodes[i]->IsDeleted())
812  {
813  map.SetDeleted(i);
814  }
815  else
816  {
817  map.SetNewIndex(i, new_index);
818  for (unsigned j=0; j<SPACE_DIM; j++)
819  {
820  mesherInput.pointlist[SPACE_DIM*new_index + j] = this->mNodes[i]->rGetLocation()[j];
821  }
822  new_index++;
823  }
824  }
825  if (elementList != NULL)
826  {
827  unsigned element_index = 0;
828 
829  // Assume there is enough space for this
830  mesherInput.numberofcorners=ELEMENT_DIM+1;
832  elem_iter != this->GetElementIteratorEnd();
833  ++elem_iter)
834  {
835 
836  for (unsigned j=0; j<=ELEMENT_DIM; j++)
837  {
838  elementList[element_index*(ELEMENT_DIM+1) + j] = (*elem_iter).GetNodeGlobalIndex(j);
839  }
840  element_index++;
841  }
842  }
843 }
844 
845 template <unsigned ELEMENT_DIM, unsigned SPACE_DIM>
846 template <class MESHER_IO>
847 void TetrahedralMesh<ELEMENT_DIM, SPACE_DIM>::ImportFromMesher(MESHER_IO& mesherOutput, unsigned numberOfElements, int *elementList, unsigned numberOfFaces, int *faceList, int *edgeMarkerList)
848 {
849  unsigned nodes_per_element = mesherOutput.numberofcorners;
850 
851  assert( nodes_per_element == ELEMENT_DIM+1 || nodes_per_element == (ELEMENT_DIM+1)*(ELEMENT_DIM+2)/2 );
852 
853  Clear();
854 
855  // Construct the nodes
856  for (unsigned node_index=0; node_index<(unsigned)mesherOutput.numberofpoints; node_index++)
857  {
858  this->mNodes.push_back(new Node<SPACE_DIM>(node_index, &mesherOutput.pointlist[node_index * SPACE_DIM], false));
859  }
860 
861  // Construct the elements
862  this->mElements.reserve(numberOfElements);
863 
864  unsigned real_element_index=0;
865  for (unsigned element_index=0; element_index<numberOfElements; element_index++)
866  {
867  std::vector<Node<SPACE_DIM>*> nodes;
868  for (unsigned j=0; j<ELEMENT_DIM+1; j++)
869  {
870  unsigned global_node_index = elementList[element_index*(nodes_per_element) + j];
871  assert(global_node_index < this->mNodes.size());
872  nodes.push_back(this->mNodes[global_node_index]);
873 
874  }
875 
876  /*
877  * For some reason, tetgen in library mode makes its initial Delaunay mesh
878  * with very thin slivers. Hence we expect to ignore some of the elements!
879  */
881  try
882  {
883  p_element = new Element<ELEMENT_DIM, SPACE_DIM>(real_element_index, nodes);
884 
885  // Shouldn't throw after this point
886  this->mElements.push_back(p_element);
887 
888  // Add the internals to quadratics
889  for (unsigned j=ELEMENT_DIM+1; j<nodes_per_element; j++)
890  {
891  unsigned global_node_index = elementList[element_index*nodes_per_element + j];
892  assert(global_node_index < this->mNodes.size());
893  this->mElements[real_element_index]->AddNode( this->mNodes[global_node_index] );
894  this->mNodes[global_node_index]->AddElement(real_element_index);
895  this->mNodes[global_node_index]->MarkAsInternal();
896  }
897  real_element_index++;
898  }
899  catch (Exception &)
900  {
901  if (SPACE_DIM == 2)
902  {
903  WARNING("Triangle has produced a zero area (collinear) element");
904  }
905  else
906  {
907  WARNING("Tetgen has produced a zero volume (coplanar) element");
908  }
909  }
910  }
911 
912  // Construct the BoundaryElements (and mark boundary nodes)
913  unsigned next_boundary_element_index = 0;
914  for (unsigned boundary_element_index=0; boundary_element_index<numberOfFaces; boundary_element_index++)
915  {
916  /*
917  * Tetgen produces only boundary faces (set edgeMarkerList to NULL).
918  * Triangle marks which edges are on the boundary.
919  */
920  if (edgeMarkerList == NULL || edgeMarkerList[boundary_element_index] == 1)
921  {
922  std::vector<Node<SPACE_DIM>*> nodes;
923  for (unsigned j=0; j<ELEMENT_DIM; j++)
924  {
925  unsigned global_node_index = faceList[boundary_element_index*ELEMENT_DIM + j];
926  assert(global_node_index < this->mNodes.size());
927  nodes.push_back(this->mNodes[global_node_index]);
928  if (!nodes[j]->IsBoundaryNode())
929  {
930  nodes[j]->SetAsBoundaryNode();
931  this->mBoundaryNodes.push_back(nodes[j]);
932  }
933  }
934 
935  /*
936  * For some reason, tetgen in library mode makes its initial Delaunay mesh
937  * with very thin slivers. Hence we expect to ignore some of the elements!
938  */
939  BoundaryElement<ELEMENT_DIM-1, SPACE_DIM>* p_b_element;
940  try
941  {
942  p_b_element = new BoundaryElement<ELEMENT_DIM-1, SPACE_DIM>(next_boundary_element_index, nodes);
943  this->mBoundaryElements.push_back(p_b_element);
944  next_boundary_element_index++;
945  }
946  catch (Exception &)
947  {
948  // Tetgen is feeding us lies //Watch this space for coverage
949  assert(SPACE_DIM == 3);
950  }
951  }
952  }
953 
955 }
956 
958 // Explicit instantiation
960 
961 template class TetrahedralMesh<1,1>;
962 template class TetrahedralMesh<1,2>;
963 template class TetrahedralMesh<1,3>;
964 template class TetrahedralMesh<2,2>;
965 template class TetrahedralMesh<2,3>;
966 template class TetrahedralMesh<3,3>;
967 
972 template void TetrahedralMesh<2,2>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
973 template void TetrahedralMesh<2,2>::ImportFromMesher<triangulateio>(triangulateio&, unsigned, int *, unsigned, int *, int *);
974 
975 template void TetrahedralMesh<3,3>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
976 template void TetrahedralMesh<3,3>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned, int *, unsigned, int *, int *);
977 
978 //The following don't ever need to be instantiated, but are needed to keep some compilers happy
979 template void TetrahedralMesh<1,2>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
980 template void TetrahedralMesh<1,2>::ImportFromMesher<triangulateio>(triangulateio&, unsigned, int *, unsigned, int *, int *);
981 
982 template void TetrahedralMesh<1,3>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
983 template void TetrahedralMesh<1,3>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned, int *, unsigned, int *, int *);
984 template void TetrahedralMesh<2,3>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
985 template void TetrahedralMesh<2,3>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned, int *, unsigned, int *, int *);
986 
987 //Intel compilation with IPO thinks that it's missing some bizarre instantiations
988 template void TetrahedralMesh<3u, 3u>::ImportFromMesher<triangulateio>(triangulateio&, unsigned int, int*, unsigned int, int*, int*);
989 template void TetrahedralMesh<1u, 1u>::ImportFromMesher<triangulateio>(triangulateio&, unsigned int, int*, unsigned int, int*, int*);
990 template void TetrahedralMesh<1u, 1u>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned int, int*, unsigned int, int*, int*);
991 template void TetrahedralMesh<2u, 2u>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned int, int*, unsigned int, int*, int*);
992 template void TetrahedralMesh<1u, 1u>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
993 
994 // Compiler on ARC cluser HAL requires the following
995 template void TetrahedralMesh<3u, 3u>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
996 template void TetrahedralMesh<1u, 1u>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
997 template void TetrahedralMesh<2u, 2u>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
998 
999 // Intel v11 compilation thinks that it's missing even more bizarre instantiations
1000 //template void TetrahedralMesh<2,2>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
1001 //template void TetrahedralMesh<3,3>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1002 //template void TetrahedralMesh<1,3>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1003 //template void TetrahedralMesh<1,1>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1004 //template void TetrahedralMesh<1,2>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
1005 //template void TetrahedralMesh<2,3>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1006 //template void TetrahedralMesh<1,3>::ImportFromMesher<triangulateio>(triangulateio&, unsigned, int *, unsigned, int *, int *);
1007 //template void TetrahedralMesh<2,3>::ImportFromMesher<triangulateio>(triangulateio&, unsigned, int *, unsigned, int *, int *);
1008 //template void TetrahedralMesh<1,2>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned, int *, unsigned, int *, int *);
1009 #ifdef _MSC_VER
1010 //MSVC requires these
1011 template void TetrahedralMesh<2,2>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
1012 template void TetrahedralMesh<3,3>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1013 template void TetrahedralMesh<1,3>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1014 template void TetrahedralMesh<1,1>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1015 template void TetrahedralMesh<1,2>::ExportToMesher<tetgen::tetgenio>(NodeMap&, tetgen::tetgenio&, int*);
1016 template void TetrahedralMesh<2,3>::ExportToMesher<triangulateio>(NodeMap&, triangulateio&, int*);
1017 template void TetrahedralMesh<1,3>::ImportFromMesher<triangulateio>(triangulateio&, unsigned, int *, unsigned, int *, int *);
1018 template void TetrahedralMesh<2,3>::ImportFromMesher<triangulateio>(triangulateio&, unsigned, int *, unsigned, int *, int *);
1019 template void TetrahedralMesh<1,2>::ImportFromMesher<tetgen::tetgenio>(tetgen::tetgenio&, unsigned, int *, unsigned, int *, int *);
1020 #endif
1021 
1026 // Serialization for Boost >= 1.36
1027 #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:171
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 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)
Node< SPACE_DIM > * GetNode(unsigned localIndex) const
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:75
std::vector< double > mElementJacobianDeterminants
unsigned GetNearestElementIndex(const ChastePoint< SPACE_DIM > &rTestPoint)
std::vector< c_vector< double, SPACE_DIM > > mBoundaryElementWeightedDirections
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
void ExportToMesher(NodeMap &map, MESHER_IO &mesherInput, int *elementList=NULL)
EdgeIterator EdgesEnd()
std::set< std::pair< unsigned, unsigned > > mEdgesVisited
virtual unsigned GetNumNodes() const =0
void FreeTriangulateIo(triangulateio &mesherIo)
virtual std::vector< double > GetNodeAttributes()
virtual unsigned GetNumFaces() const =0
BoundaryElementIterator GetBoundaryElementIteratorEnd() const