Chaste Commit::1fd4e48e3990e67db148bc1bc4cf6991a0049d0c
FineCoarseMeshPair.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 "FineCoarseMeshPair.hpp"
37
38template<unsigned DIM>
40 : mrFineMesh(rFineMesh),
41 mrCoarseMesh(rCoarseMesh),
42 mpFineMeshBoxCollection(nullptr),
43 mpCoarseMeshBoxCollection(nullptr)
44{
46}
47
48template<unsigned DIM>
53
54template<unsigned DIM>
56{
57 return mrCoarseMesh;
58}
59
60template<unsigned DIM>
62{
63 DeleteFineBoxCollection();
64 DeleteCoarseBoxCollection();
65}
66
67template<unsigned DIM>
69{
70 if (mpFineMeshBoxCollection != nullptr)
71 {
72 delete mpFineMeshBoxCollection;
73 mpFineMeshBoxCollection = nullptr;
74 }
75}
76
77template<unsigned DIM>
79{
80 if (mpCoarseMeshBoxCollection != nullptr)
81 {
82 delete mpCoarseMeshBoxCollection;
83 mpCoarseMeshBoxCollection = nullptr;
84 }
85}
86
88// Setting up boxes methods
90
91template<unsigned DIM>
93{
94 SetUpBoxes(mrFineMesh, boxWidth, mpFineMeshBoxCollection);
95}
96
97template<unsigned DIM>
99{
100 SetUpBoxes(mrCoarseMesh, boxWidth, mpCoarseMeshBoxCollection);
101}
102
103template<unsigned DIM>
105 double boxWidth,
106 DistributedBoxCollection<DIM>*& rpBoxCollection)
107{
108 if (rpBoxCollection)
109 {
110 delete rpBoxCollection;
111 rpBoxCollection = nullptr;
112 }
113
114 // Compute min and max values for the fine mesh nodes
115 ChasteCuboid<DIM> bounding_box = rMesh.CalculateBoundingBox();
116
117 // Set up the boxes using a domain that is slightly larger than the fine mesh
118 c_vector<double,2*DIM> extended_min_and_max;
119 for (unsigned i=0; i<DIM; i++)
120 {
121 double width = bounding_box.GetWidth(i);
122
123 // Subtract from the minima
124 extended_min_and_max(2*i) = bounding_box.rGetLowerCorner()[i] - 0.05*width;
125
126 // Add to the maxima
127 extended_min_and_max(2*i+1) = bounding_box.rGetUpperCorner()[i] + 0.05*width;
128 }
129
130 if (boxWidth < 0)
131 {
132 /*
133 * Use default value = max(max_edge_length, w20), where w20 is the width
134 * corresponding to 20 boxes in the x-direction.
135 *
136 * BoxCollection creates an extra box so divide by 19 not 20. Add a little
137 * bit on to ensure minor numerical fluctuations don't change the answer.
138 */
139 boxWidth = (extended_min_and_max(1) - extended_min_and_max(0))/19.000000001;
140
141 // Determine the maximum edge length
142 c_vector<double, 2> min_max_edge_length = rMesh.CalculateMinMaxEdgeLengths();
143
144 if (boxWidth < min_max_edge_length[1])
145 {
146 boxWidth = 1.1*min_max_edge_length[1];
147 }
148 }
149
150 rpBoxCollection = new DistributedBoxCollection<DIM>(boxWidth, extended_min_and_max);
151 rpBoxCollection->SetupAllLocalBoxes();
152
153 // For each element, if ANY of its nodes are physically in a box, put that element in that box
154 for (unsigned i=0; i<rMesh.GetNumElements(); i++)
155 {
156 Element<DIM,DIM>* p_element = rMesh.GetElement(i);
157
158 std::set<unsigned> box_indices_each_node_this_elem;
159 for (unsigned j=0; j<DIM+1; j++) // num vertices per element
160 {
161 Node<DIM>* p_node = p_element->GetNode(j);
162 // Only take note of box inclusions which are in our domain
163 if (rpBoxCollection->IsOwned(p_node))
164 {
165 unsigned box_index = rpBoxCollection->CalculateContainingBox(p_node);
166 box_indices_each_node_this_elem.insert(box_index);
167 }
168 }
169 for (std::set<unsigned>::iterator iter = box_indices_each_node_this_elem.begin();
170 iter != box_indices_each_node_this_elem.end();
171 ++iter)
172 {
173 assert(rpBoxCollection->IsBoxOwned( *iter ));
174 rpBoxCollection->rGetBox( *iter ).AddElement(p_element);
175 }
176 }
177}
178
180// ComputeFineElementsAndWeightsForCoarseQuadPoints()
181// and
182// ComputeFineElementsAndWeightsForCoarseNodes()
183// and
184// common method
186
187template<unsigned DIM>
189 bool safeMode)
190{
191 if (mpFineMeshBoxCollection == nullptr)
192 {
193 EXCEPTION("Call SetUpBoxesOnFineMesh() before ComputeFineElementsAndWeightsForCoarseQuadPoints()");
194 }
195
196 // Get the quad point (physical) positions
197 QuadraturePointsGroup<DIM> quad_point_posns(mrCoarseMesh, rQuadRule);
198
199 // Resize the elements and weights vector.
200 mFineMeshElementsAndWeights.resize(quad_point_posns.Size());
201 // Make sure that, in parallel, silent processes have their structs initialised to zero values
202 for (unsigned i=0; i<mFineMeshElementsAndWeights.size(); i++)
203 {
204 mFineMeshElementsAndWeights[i].ElementNum = 0u;
205 mFineMeshElementsAndWeights[i].Weights = zero_vector<double>(DIM+1);
206 }
207
208
209 // LCOV_EXCL_START
210 if (CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
211 {
212 std::cout << "\nComputing fine elements and weights for coarse quad points\n";
213 }
214 // LCOV_EXCL_STOP
215
216
217 ResetStatisticsVariables();
218 for (unsigned i=0; i<quad_point_posns.Size(); i++)
219 {
220 // LCOV_EXCL_START
221 if (CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
222 {
223 std::cout << "\t" << i << " of " << quad_point_posns.Size() << std::flush;
224 }
225 // LCOV_EXCL_STOP
226
227 // Get the box this point is in
228 unsigned box_for_this_point = mpFineMeshBoxCollection->CalculateContainingBox( quad_point_posns.rGet(i) );
229 if (mpFineMeshBoxCollection->IsBoxOwned(box_for_this_point))
230 {
231 // A chaste point version of the c-vector is needed for the GetContainingElement call.
232 ChastePoint<DIM> point(quad_point_posns.rGet(i));
233
234 ComputeFineElementAndWeightForGivenPoint(point, safeMode, box_for_this_point, i);
235 }
236 else
237 {
238 assert(mFineMeshElementsAndWeights[i].ElementNum == 0u);
239 assert(norm_2(mFineMeshElementsAndWeights[i].Weights) == 0.0 );
240 }
241 }
242 ShareFineElementData();
243 if (mStatisticsCounters[1] > 0)
244 {
245 WARNING(mStatisticsCounters[1] << " of " << quad_point_posns.Size() << " coarse-mesh quadrature points were outside the fine mesh");
246 }
247}
248
249template<unsigned DIM>
251{
252 if (mpFineMeshBoxCollection==nullptr)
253 {
254 EXCEPTION("Call SetUpBoxesOnFineMesh() before ComputeFineElementsAndWeightsForCoarseNodes()");
255 }
256
257 // Resize the elements and weights vector.
258 mFineMeshElementsAndWeights.resize(mrCoarseMesh.GetNumNodes());
259 // Make sure that, in parallel, silent processes have their structs initialised to zero values
260 for (unsigned i=0; i<mFineMeshElementsAndWeights.size(); i++)
261 {
262 mFineMeshElementsAndWeights[i].ElementNum = 0u;
263 mFineMeshElementsAndWeights[i].Weights = zero_vector<double>(DIM+1);
264 }
265
266 // LCOV_EXCL_START
267 if (CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
268 {
269 std::cout << "\nComputing fine elements and weights for coarse nodes\n";
270 }
271 // LCOV_EXCL_STOP
272
273
274 ResetStatisticsVariables();
275 for (unsigned i=0; i<mrCoarseMesh.GetNumNodes(); i++)
276 {
277 // LCOV_EXCL_START
278 if (CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
279 {
280 std::cout << "\t" << i << " of " << mrCoarseMesh.GetNumNodes() << std::flush;
281 }
282 // LCOV_EXCL_STOP
283
284 Node<DIM>* p_node = mrCoarseMesh.GetNode(i);
285
286 // Get the box this point is in
287 unsigned box_for_this_point = mpFineMeshBoxCollection->CalculateContainingBox( p_node->rGetModifiableLocation() );
288 if (mpFineMeshBoxCollection->IsBoxOwned(box_for_this_point))
289 {
290 // A chaste point version of the c-vector is needed for the GetContainingElement call
291 ChastePoint<DIM> point(p_node->rGetLocation());
292
293 ComputeFineElementAndWeightForGivenPoint(point, safeMode, box_for_this_point, i);
294 }
295 }
296 ShareFineElementData();
297}
298
305template<unsigned DIM>
307 bool safeMode,
308 unsigned boxForThisPoint,
309 unsigned index)
310{
311 std::set<unsigned> test_element_indices;
312
313 /*
314 * The elements to try (initially) are those contained in the box the point is in.
315 *
316 * Note: it is possible the point to be in an element not 'in' this box, as it is
317 * possible for all element nodes to be in different boxes.
318 */
319 CollectElementsInContainingBox(mpFineMeshBoxCollection, boxForThisPoint, test_element_indices);
320
321 unsigned elem_index;
322 c_vector<double,DIM+1> weight;
323
324 try
325 {
326 //std::cout << "\n" << "# test elements initially " << test_element_indices.size() << "\n";
327 // Try these elements only, initially
328 elem_index = mrFineMesh.GetContainingElementIndex(rPoint,
329 false,
330 test_element_indices,
331 true /* quit if not in test_elements */);
332 weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
333
334 mStatisticsCounters[0]++;
335 }
336 catch(Exception&) //not_in_box
337 {
338 // Now try all elements, trying the elements contained in the boxes local to this element first
339 test_element_indices.clear();
340
341 CollectElementsInLocalBoxes(mpFineMeshBoxCollection, boxForThisPoint, test_element_indices);
342
343 try
344 {
345 elem_index = mrFineMesh.GetContainingElementIndex(rPoint, false, test_element_indices, true);
346 weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
347 mStatisticsCounters[0]++;
348 }
349 catch(Exception&) //not_in_local_boxes
350 {
351 if (safeMode)
352 {
353 // Try the remaining elements
354 try
355 {
356 elem_index = mrFineMesh.GetContainingElementIndex(rPoint, false);
357 weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
358 mStatisticsCounters[0]++;
359
360 }
361 catch (Exception&) // not_in_mesh
362 {
363 // The point is not in ANY element, so store the nearest element and corresponding weights
364 elem_index = mrFineMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
365 weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
366
367 mNotInMesh.push_back(index);
368 mNotInMeshNearestElementWeights.push_back(weight);
369 mStatisticsCounters[1]++;
370 }
371 }
372 else
373 {
374 assert(test_element_indices.size() > 0); // boxes probably too small if this fails
375
376 /*
377 * Immediately assume it isn't in the rest of the mesh - this should be the
378 * case assuming the box width was chosen suitably. Store the nearest element
379 * and corresponding weights.
380 */
381 elem_index = mrFineMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
382 weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
383
384 mNotInMesh.push_back(index);
385 mNotInMeshNearestElementWeights.push_back(weight);
386 mStatisticsCounters[1]++;
387 }
388 }
389 }
390
391 mFineMeshElementsAndWeights[index].ElementNum = elem_index;
392 mFineMeshElementsAndWeights[index].Weights = weight;
393}
394
396// ComputeCoarseElementsForFineNodes
397// and
398// ComputeCoarseElementsForFineElementCentroids
399// and
400// common method
402
403template<unsigned DIM>
405{
406 if (mpCoarseMeshBoxCollection==nullptr)
407 {
408 EXCEPTION("Call SetUpBoxesOnCoarseMesh() before ComputeCoarseElementsForFineNodes()");
409 }
410
411 // LCOV_EXCL_START
412 if (CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
413 {
414 std::cout << "\nComputing coarse elements for fine nodes\n";
415 }
416 // LCOV_EXCL_STOP
417 mCoarseElementsForFineNodes.clear();
418 mCoarseElementsForFineNodes.resize(mrFineMesh.GetNumNodes(), 0.0);
419
420 ResetStatisticsVariables();
421 for (unsigned i=0; i<mCoarseElementsForFineNodes.size(); i++)
422 {
423 // LCOV_EXCL_START
424 if (CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
425 {
426 std::cout << "\t" << i << " of " << mCoarseElementsForFineNodes.size() << std::flush;
427 }
428 // LCOV_EXCL_STOP
429
430 ChastePoint<DIM> point = mrFineMesh.GetNode(i)->GetPoint();
431
432 // Get the box this point is in
433 unsigned box_for_this_point = mpCoarseMeshBoxCollection->CalculateContainingBox(mrFineMesh.GetNode(i)->rGetModifiableLocation());
434 if (mpCoarseMeshBoxCollection->IsBoxOwned(box_for_this_point))
435 {
436 mCoarseElementsForFineNodes[i] = ComputeCoarseElementForGivenPoint(point, safeMode, box_for_this_point);
437 }
438 }
439 ShareCoarseElementData();
440}
441
442template<unsigned DIM>
444{
445 if (mpCoarseMeshBoxCollection==nullptr)
446 {
447 EXCEPTION("Call SetUpBoxesOnCoarseMesh() before ComputeCoarseElementsForFineElementCentroids()");
448 }
449
450 // LCOV_EXCL_START
451 if (CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
452 {
453 std::cout << "\nComputing coarse elements for fine element centroids\n";
454 }
455 // LCOV_EXCL_STOP
456
457 mCoarseElementsForFineElementCentroids.clear();
458 mCoarseElementsForFineElementCentroids.resize(mrFineMesh.GetNumElements(), 0.0);
459
460 ResetStatisticsVariables();
461 for (unsigned i=0; i<mrFineMesh.GetNumElements(); i++)
462 {
463 // LCOV_EXCL_START
464 if (CommandLineArguments::Instance()->OptionExists("-mesh_pair_verbose"))
465 {
466 std::cout << "\t" << i << " of " << mrFineMesh.GetNumElements() << std::flush;
467 }
468 // LCOV_EXCL_STOP
469
470 c_vector<double,DIM> point_cvec = mrFineMesh.GetElement(i)->CalculateCentroid();
471 ChastePoint<DIM> point(point_cvec);
472
473 // Get the box this point is in
474 unsigned box_for_this_point = mpCoarseMeshBoxCollection->CalculateContainingBox( point_cvec );
475
476 if (mpCoarseMeshBoxCollection->IsBoxOwned(box_for_this_point))
477 {
478 mCoarseElementsForFineElementCentroids[i] = ComputeCoarseElementForGivenPoint(point, safeMode, box_for_this_point);
479 }
480 }
481 ShareCoarseElementData();
482}
483
484template<unsigned DIM>
486 bool safeMode,
487 unsigned boxForThisPoint)
488{
495 std::set<unsigned> test_element_indices;
496 CollectElementsInContainingBox(mpCoarseMeshBoxCollection, boxForThisPoint, test_element_indices);
497
498 unsigned elem_index;
499
500 try
501 {
502 elem_index = mrCoarseMesh.GetContainingElementIndex(rPoint,
503 false,
504 test_element_indices,
505 true /* quit if not in test_elements */);
506
507 mStatisticsCounters[0]++;
508 }
509 catch(Exception&) // not_in_box
510 {
511 // Now try all elements, trying the elements contained in the boxes local to this element first
512 test_element_indices.clear();
513 CollectElementsInLocalBoxes(mpCoarseMeshBoxCollection, boxForThisPoint, test_element_indices);
514
515 try
516 {
517 elem_index = mrCoarseMesh.GetContainingElementIndex(rPoint, false, test_element_indices, true);
518 mStatisticsCounters[0]++;
519 }
520 catch(Exception&) // not_in_local_boxes
521 {
522 if (safeMode)
523 {
524 // Try the remaining elements
525 try
526 {
527 elem_index = mrCoarseMesh.GetContainingElementIndex(rPoint, false);
528
529 mStatisticsCounters[0]++;
530 }
531 catch (Exception&) // not_in_mesh
532 {
533 // The point is not in ANY element, so store the nearest element and corresponding weights
534 elem_index = mrCoarseMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
535 mStatisticsCounters[1]++;
536 }
537 }
538 else
539 {
540 assert(test_element_indices.size() > 0); // boxes probably too small if this fails
541
542 /*
543 * Immediately assume it isn't in the rest of the mesh - this should be the
544 * case assuming the box width was chosen suitably. Store the nearest element
545 * and corresponding weights.
546 */
547 elem_index = mrCoarseMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
548 mStatisticsCounters[1]++;
549 }
550 }
551 }
552
553 return elem_index;
554}
555
557// Helper methods for code
559
560template<unsigned DIM>
562 unsigned boxIndex,
563 std::set<unsigned>& rElementIndices)
564{
565 for (typename std::set<Element<DIM,DIM>*>::iterator elem_iter = rpBoxCollection->rGetBox(boxIndex).rGetElementsContained().begin();
566 elem_iter != rpBoxCollection->rGetBox(boxIndex).rGetElementsContained().end();
567 ++elem_iter)
568 {
569 rElementIndices.insert((*elem_iter)->GetIndex());
570 }
571}
572
573template<unsigned DIM>
575 unsigned boxIndex,
576 std::set<unsigned>& rElementIndices)
577{
578 std::set<unsigned> local_boxes = rpBoxCollection->rGetLocalBoxes(boxIndex);
579 for (std::set<unsigned>::iterator local_box_iter = local_boxes.begin();
580 local_box_iter != local_boxes.end();
581 ++local_box_iter)
582 {
583 for (typename std::set<Element<DIM,DIM>*>::iterator elem_iter = rpBoxCollection->rGetBox(*local_box_iter).rGetElementsContained().begin();
584 elem_iter != rpBoxCollection->rGetBox(*local_box_iter).rGetElementsContained().end();
585 ++elem_iter)
586 {
587 rElementIndices.insert((*elem_iter)->GetIndex());
588 }
589 }
590}
591
593// Statistics related methods
595
596template<unsigned DIM>
598{
599 mNotInMesh.clear();
600 mNotInMeshNearestElementWeights.clear();
601 mStatisticsCounters.resize(2, 0u);
602}
603
604template<unsigned DIM>
606{
608 {
609 return;
610 }
611
612 // This sums the results so it isn't idempotent: you get a different result if you call this method twice
613 // Should not matter: the methods which call this helper method have reset everything which is about to be shared
614 std::vector<unsigned> local_counters = mStatisticsCounters;
615 MPI_Allreduce(&local_counters[0], &mStatisticsCounters[0], 2u, MPI_UNSIGNED, MPI_SUM, PETSC_COMM_WORLD);
616
617
618 // Get all the element number and weights into a contiguous format
619 unsigned elements_size = mFineMeshElementsAndWeights.size();
620 std::vector<unsigned> all_element_indices(elements_size);
621 unsigned weights_size = elements_size*(DIM+1);
622 std::vector<double> all_weights(weights_size);
623 for (unsigned index=0; index<mFineMeshElementsAndWeights.size(); index++)
624 {
625 all_element_indices[index] = mFineMeshElementsAndWeights[index].ElementNum;
626 for (unsigned j=0; j<DIM+1; j++)
627 {
628 all_weights[index*(DIM+1)+j] = mFineMeshElementsAndWeights[index].Weights[j];
629 }
630 }
631
632 // Copy and share
633 std::vector<unsigned> local_all_element_indices = all_element_indices;
634 std::vector<double> local_all_weights = all_weights;
635 MPI_Allreduce(&local_all_element_indices[0], &all_element_indices[0], elements_size, MPI_UNSIGNED, MPI_SUM, PETSC_COMM_WORLD);
636 MPI_Allreduce( &local_all_weights[0], &all_weights[0], weights_size, MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
637
638 // Put back into the regular data structure
639 for (unsigned index=0; index<mFineMeshElementsAndWeights.size(); index++)
640 {
641 mFineMeshElementsAndWeights[index].ElementNum = all_element_indices[index];
642 for (unsigned j=0; j<DIM+1; j++)
643 {
644 mFineMeshElementsAndWeights[index].Weights[j] = all_weights[index*(DIM+1)+j] ;
645 }
646 }
647}
648
649template<unsigned DIM>
651{
653 {
654 return;
655 }
656
657 // This sums the results so it isn't idempotent: you get a different result if you call this method twice
658 // Should not matter: the methods which call this helper method have reset #mStatisticsCounters
659 std::vector<unsigned> local_counters = mStatisticsCounters;
660 MPI_Allreduce(&local_counters[0], &mStatisticsCounters[0], 2u, MPI_UNSIGNED, MPI_SUM, PETSC_COMM_WORLD);
661
662 // The rest uses "max" so it is idempotent. You can safely re-share results between processes without them changing.
663 if (mCoarseElementsForFineNodes.empty() == false)
664 {
665 std::vector<unsigned> temp_coarse_elements = mCoarseElementsForFineNodes;
666 MPI_Allreduce( &temp_coarse_elements[0], &mCoarseElementsForFineNodes[0], mCoarseElementsForFineNodes.size(), MPI_UNSIGNED, MPI_MAX, PETSC_COMM_WORLD);
667 }
668 if (mCoarseElementsForFineElementCentroids.empty() == false)
669 {
670 std::vector<unsigned> temp_coarse_elements = mCoarseElementsForFineElementCentroids;
671 MPI_Allreduce( &temp_coarse_elements[0], &mCoarseElementsForFineElementCentroids[0], mCoarseElementsForFineElementCentroids.size(), MPI_UNSIGNED, MPI_MAX, PETSC_COMM_WORLD);
672 }
673}
674
675template<unsigned DIM>
677{
678 std::cout << "\nFineCoarseMeshPair statistics for the last-called method:\n";
679
680// std::cout << "\tNum points for which containing element was found, using box containing that point = " << mStatisticsCounters[0] << "\n";
681// std::cout << "\tNum points for which containing element was in local box = " << mStatisticsCounters[1] << "\n";
682// std::cout << "\tNum points for which containing element was in an element in a non-local box = " << mStatisticsCounters[2] << "\n";
683// std::cout << "\tNum points for which no containing element was found = " << mStatisticsCounters[3] << "\n";
684
685 std::cout << "\tNum points for which containing element was found: " << mStatisticsCounters[0] << "\n";
686 std::cout << "\tNum points for which no containing element was found = " << mStatisticsCounters[1] << "\n";
687
688 if (mNotInMesh.size() > 0)
689 {
690 std::cout << "\tIndices and weights for points (nodes/quad points) for which no containing element was found:\n";
691 for (unsigned i=0; i<mNotInMesh.size(); i++)
692 {
693 std::cout << "\t\t" << mNotInMesh[i] << ", " << mNotInMeshNearestElementWeights[i] << "\n";
694 }
695 }
696}
697
699
700template class FineCoarseMeshPair<1>;
701template class FineCoarseMeshPair<2>;
702template class FineCoarseMeshPair<3>;
#define EXCEPTION(message)
Node< SPACE_DIM > * GetNode(unsigned localIndex) const
ChasteCuboid< SPACE_DIM > CalculateBoundingBox(const std::vector< Node< SPACE_DIM > * > &rNodes) const
Element< ELEMENT_DIM, SPACE_DIM > * GetElement(unsigned index) const
virtual c_vector< double, 2 > CalculateMinMaxEdgeLengths()
virtual unsigned GetNumElements() const
const ChastePoint< SPACE_DIM > & rGetUpperCorner() const
const ChastePoint< SPACE_DIM > & rGetLowerCorner() const
double GetWidth(unsigned rDimension) const
static CommandLineArguments * Instance()
std::set< unsigned > & rGetLocalBoxes(unsigned boxIndex)
unsigned CalculateContainingBox(Node< DIM > *pNode)
bool IsBoxOwned(unsigned globalIndex)
Box< DIM > & rGetBox(unsigned boxIndex)
void SetUpBoxesOnFineMesh(double boxWidth=-1)
void CollectElementsInContainingBox(DistributedBoxCollection< DIM > *&rpBoxCollection, unsigned boxIndex, std::set< unsigned > &rElementIndices)
void ComputeFineElementsAndWeightsForCoarseQuadPoints(GaussianQuadratureRule< DIM > &rQuadRule, bool safeMode)
void SetUpBoxes(AbstractTetrahedralMesh< DIM, DIM > &rMesh, double boxWidth, DistributedBoxCollection< DIM > *&rpBoxCollection)
const AbstractTetrahedralMesh< DIM, DIM > & GetCoarseMesh() const
FineCoarseMeshPair(AbstractTetrahedralMesh< DIM, DIM > &rFineMesh, AbstractTetrahedralMesh< DIM, DIM > &rCoarseMesh)
void SetUpBoxesOnCoarseMesh(double boxWidth=-1)
const AbstractTetrahedralMesh< DIM, DIM > & GetFineMesh() const
void ComputeFineElementsAndWeightsForCoarseNodes(bool safeMode)
unsigned ComputeCoarseElementForGivenPoint(ChastePoint< DIM > &rPoint, bool safeMode, unsigned boxForThisPoint)
void ComputeCoarseElementsForFineNodes(bool safeMode)
void ComputeCoarseElementsForFineElementCentroids(bool safeMode)
void CollectElementsInLocalBoxes(DistributedBoxCollection< DIM > *&rpBoxCollection, unsigned boxIndex, std::set< unsigned > &rElementIndices)
void ComputeFineElementAndWeightForGivenPoint(ChastePoint< DIM > &rPoint, bool safeMode, unsigned boxForThisPoint, unsigned index)
Definition Node.hpp:59
c_vector< double, SPACE_DIM > & rGetModifiableLocation()
Definition Node.cpp:151
const c_vector< double, SPACE_DIM > & rGetLocation() const
Definition Node.cpp:139
static bool IsSequential()
c_vector< double, DIM > & rGet(unsigned elementIndex, unsigned quadIndex)