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