36#include "FineCoarseMeshPair.hpp"
40 : mrFineMesh(rFineMesh),
41 mrCoarseMesh(rCoarseMesh),
42 mpFineMeshBoxCollection(nullptr),
43 mpCoarseMeshBoxCollection(nullptr)
63 DeleteFineBoxCollection();
64 DeleteCoarseBoxCollection();
70 if (mpFineMeshBoxCollection !=
nullptr)
72 delete mpFineMeshBoxCollection;
73 mpFineMeshBoxCollection =
nullptr;
80 if (mpCoarseMeshBoxCollection !=
nullptr)
82 delete mpCoarseMeshBoxCollection;
83 mpCoarseMeshBoxCollection =
nullptr;
94 SetUpBoxes(mrFineMesh, boxWidth, mpFineMeshBoxCollection);
100 SetUpBoxes(mrCoarseMesh, boxWidth, mpCoarseMeshBoxCollection);
103template<
unsigned DIM>
110 delete rpBoxCollection;
111 rpBoxCollection =
nullptr;
118 c_vector<double,2*DIM> extended_min_and_max;
119 for (
unsigned i=0; i<DIM; i++)
121 double width = bounding_box.
GetWidth(i);
124 extended_min_and_max(2*i) = bounding_box.
rGetLowerCorner()[i] - 0.05*width;
127 extended_min_and_max(2*i+1) = bounding_box.
rGetUpperCorner()[i] + 0.05*width;
139 boxWidth = (extended_min_and_max(1) - extended_min_and_max(0))/19.000000001;
144 if (boxWidth < min_max_edge_length[1])
146 boxWidth = 1.1*min_max_edge_length[1];
158 std::set<unsigned> box_indices_each_node_this_elem;
159 for (
unsigned j=0; j<DIM+1; j++)
163 if (rpBoxCollection->
IsOwned(p_node))
166 box_indices_each_node_this_elem.insert(box_index);
169 for (std::set<unsigned>::iterator iter = box_indices_each_node_this_elem.begin();
170 iter != box_indices_each_node_this_elem.end();
174 rpBoxCollection->
rGetBox( *iter ).AddElement(p_element);
187template<
unsigned DIM>
191 if (mpFineMeshBoxCollection ==
nullptr)
193 EXCEPTION(
"Call SetUpBoxesOnFineMesh() before ComputeFineElementsAndWeightsForCoarseQuadPoints()");
200 mFineMeshElementsAndWeights.resize(quad_point_posns.
Size());
202 for (
unsigned i=0; i<mFineMeshElementsAndWeights.size(); i++)
204 mFineMeshElementsAndWeights[i].ElementNum = 0u;
205 mFineMeshElementsAndWeights[i].Weights = zero_vector<double>(DIM+1);
212 std::cout <<
"\nComputing fine elements and weights for coarse quad points\n";
217 ResetStatisticsVariables();
218 for (
unsigned i=0; i<quad_point_posns.
Size(); i++)
223 std::cout <<
"\t" << i <<
" of " << quad_point_posns.
Size() << std::flush;
228 unsigned box_for_this_point = mpFineMeshBoxCollection->CalculateContainingBox( quad_point_posns.
rGet(i) );
229 if (mpFineMeshBoxCollection->IsBoxOwned(box_for_this_point))
234 ComputeFineElementAndWeightForGivenPoint(point, safeMode, box_for_this_point, i);
238 assert(mFineMeshElementsAndWeights[i].ElementNum == 0u);
239 assert(norm_2(mFineMeshElementsAndWeights[i].Weights) == 0.0 );
242 ShareFineElementData();
243 if (mStatisticsCounters[1] > 0)
245 WARNING(mStatisticsCounters[1] <<
" of " << quad_point_posns.
Size() <<
" coarse-mesh quadrature points were outside the fine mesh");
249template<
unsigned DIM>
252 if (mpFineMeshBoxCollection==
nullptr)
254 EXCEPTION(
"Call SetUpBoxesOnFineMesh() before ComputeFineElementsAndWeightsForCoarseNodes()");
258 mFineMeshElementsAndWeights.resize(mrCoarseMesh.GetNumNodes());
260 for (
unsigned i=0; i<mFineMeshElementsAndWeights.size(); i++)
262 mFineMeshElementsAndWeights[i].ElementNum = 0u;
263 mFineMeshElementsAndWeights[i].Weights = zero_vector<double>(DIM+1);
269 std::cout <<
"\nComputing fine elements and weights for coarse nodes\n";
274 ResetStatisticsVariables();
275 for (
unsigned i=0; i<mrCoarseMesh.GetNumNodes(); i++)
280 std::cout <<
"\t" << i <<
" of " << mrCoarseMesh.GetNumNodes() << std::flush;
284 Node<DIM>* p_node = mrCoarseMesh.GetNode(i);
287 unsigned box_for_this_point = mpFineMeshBoxCollection->CalculateContainingBox( p_node->
rGetModifiableLocation() );
288 if (mpFineMeshBoxCollection->IsBoxOwned(box_for_this_point))
293 ComputeFineElementAndWeightForGivenPoint(point, safeMode, box_for_this_point, i);
296 ShareFineElementData();
305template<
unsigned DIM>
308 unsigned boxForThisPoint,
311 std::set<unsigned> test_element_indices;
319 CollectElementsInContainingBox(mpFineMeshBoxCollection, boxForThisPoint, test_element_indices);
322 c_vector<double,DIM+1> weight;
328 elem_index = mrFineMesh.GetContainingElementIndex(rPoint,
330 test_element_indices,
332 weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
334 mStatisticsCounters[0]++;
339 test_element_indices.clear();
341 CollectElementsInLocalBoxes(mpFineMeshBoxCollection, boxForThisPoint, test_element_indices);
345 elem_index = mrFineMesh.GetContainingElementIndex(rPoint,
false, test_element_indices,
true);
346 weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
347 mStatisticsCounters[0]++;
356 elem_index = mrFineMesh.GetContainingElementIndex(rPoint,
false);
357 weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
358 mStatisticsCounters[0]++;
364 elem_index = mrFineMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
365 weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
367 mNotInMesh.push_back(index);
368 mNotInMeshNearestElementWeights.push_back(weight);
369 mStatisticsCounters[1]++;
374 assert(test_element_indices.size() > 0);
381 elem_index = mrFineMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
382 weight = mrFineMesh.GetElement(elem_index)->CalculateInterpolationWeights(rPoint);
384 mNotInMesh.push_back(index);
385 mNotInMeshNearestElementWeights.push_back(weight);
386 mStatisticsCounters[1]++;
391 mFineMeshElementsAndWeights[index].ElementNum = elem_index;
392 mFineMeshElementsAndWeights[index].Weights = weight;
403template<
unsigned DIM>
406 if (mpCoarseMeshBoxCollection==
nullptr)
408 EXCEPTION(
"Call SetUpBoxesOnCoarseMesh() before ComputeCoarseElementsForFineNodes()");
414 std::cout <<
"\nComputing coarse elements for fine nodes\n";
417 mCoarseElementsForFineNodes.clear();
418 mCoarseElementsForFineNodes.resize(mrFineMesh.GetNumNodes(), 0.0);
420 ResetStatisticsVariables();
421 for (
unsigned i=0; i<mCoarseElementsForFineNodes.size(); i++)
426 std::cout <<
"\t" << i <<
" of " << mCoarseElementsForFineNodes.size() << std::flush;
433 unsigned box_for_this_point = mpCoarseMeshBoxCollection->CalculateContainingBox(mrFineMesh.GetNode(i)->rGetModifiableLocation());
434 if (mpCoarseMeshBoxCollection->IsBoxOwned(box_for_this_point))
436 mCoarseElementsForFineNodes[i] = ComputeCoarseElementForGivenPoint(point, safeMode, box_for_this_point);
439 ShareCoarseElementData();
442template<
unsigned DIM>
445 if (mpCoarseMeshBoxCollection==
nullptr)
447 EXCEPTION(
"Call SetUpBoxesOnCoarseMesh() before ComputeCoarseElementsForFineElementCentroids()");
453 std::cout <<
"\nComputing coarse elements for fine element centroids\n";
457 mCoarseElementsForFineElementCentroids.clear();
458 mCoarseElementsForFineElementCentroids.resize(mrFineMesh.GetNumElements(), 0.0);
460 ResetStatisticsVariables();
461 for (
unsigned i=0; i<mrFineMesh.GetNumElements(); i++)
466 std::cout <<
"\t" << i <<
" of " << mrFineMesh.GetNumElements() << std::flush;
470 c_vector<double,DIM> point_cvec = mrFineMesh.GetElement(i)->CalculateCentroid();
474 unsigned box_for_this_point = mpCoarseMeshBoxCollection->CalculateContainingBox( point_cvec );
476 if (mpCoarseMeshBoxCollection->IsBoxOwned(box_for_this_point))
478 mCoarseElementsForFineElementCentroids[i] = ComputeCoarseElementForGivenPoint(point, safeMode, box_for_this_point);
481 ShareCoarseElementData();
484template<
unsigned DIM>
487 unsigned boxForThisPoint)
495 std::set<unsigned> test_element_indices;
496 CollectElementsInContainingBox(mpCoarseMeshBoxCollection, boxForThisPoint, test_element_indices);
502 elem_index = mrCoarseMesh.GetContainingElementIndex(rPoint,
504 test_element_indices,
507 mStatisticsCounters[0]++;
512 test_element_indices.clear();
513 CollectElementsInLocalBoxes(mpCoarseMeshBoxCollection, boxForThisPoint, test_element_indices);
517 elem_index = mrCoarseMesh.GetContainingElementIndex(rPoint,
false, test_element_indices,
true);
518 mStatisticsCounters[0]++;
527 elem_index = mrCoarseMesh.GetContainingElementIndex(rPoint,
false);
529 mStatisticsCounters[0]++;
534 elem_index = mrCoarseMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
535 mStatisticsCounters[1]++;
540 assert(test_element_indices.size() > 0);
547 elem_index = mrCoarseMesh.GetNearestElementIndexFromTestElements(rPoint,test_element_indices);
548 mStatisticsCounters[1]++;
560template<
unsigned DIM>
563 std::set<unsigned>& rElementIndices)
565 for (
typename std::set<
Element<DIM,DIM>*>::iterator elem_iter = rpBoxCollection->
rGetBox(boxIndex).rGetElementsContained().begin();
566 elem_iter != rpBoxCollection->
rGetBox(boxIndex).rGetElementsContained().end();
569 rElementIndices.insert((*elem_iter)->GetIndex());
573template<
unsigned DIM>
576 std::set<unsigned>& rElementIndices)
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();
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();
587 rElementIndices.insert((*elem_iter)->GetIndex());
596template<
unsigned DIM>
600 mNotInMeshNearestElementWeights.clear();
601 mStatisticsCounters.resize(2, 0u);
604template<
unsigned DIM>
614 std::vector<unsigned> local_counters = mStatisticsCounters;
615 MPI_Allreduce(&local_counters[0], &mStatisticsCounters[0], 2u, MPI_UNSIGNED, MPI_SUM, PETSC_COMM_WORLD);
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++)
625 all_element_indices[index] = mFineMeshElementsAndWeights[index].ElementNum;
626 for (
unsigned j=0; j<DIM+1; j++)
628 all_weights[index*(DIM+1)+j] = mFineMeshElementsAndWeights[index].Weights[j];
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);
639 for (
unsigned index=0; index<mFineMeshElementsAndWeights.size(); index++)
641 mFineMeshElementsAndWeights[index].ElementNum = all_element_indices[index];
642 for (
unsigned j=0; j<DIM+1; j++)
644 mFineMeshElementsAndWeights[index].Weights[j] = all_weights[index*(DIM+1)+j] ;
649template<
unsigned DIM>
659 std::vector<unsigned> local_counters = mStatisticsCounters;
660 MPI_Allreduce(&local_counters[0], &mStatisticsCounters[0], 2u, MPI_UNSIGNED, MPI_SUM, PETSC_COMM_WORLD);
663 if (mCoarseElementsForFineNodes.empty() ==
false)
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);
668 if (mCoarseElementsForFineElementCentroids.empty() ==
false)
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);
675template<
unsigned DIM>
678 std::cout <<
"\nFineCoarseMeshPair statistics for the last-called method:\n";
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";
688 if (mNotInMesh.size() > 0)
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++)
693 std::cout <<
"\t\t" << mNotInMesh[i] <<
", " << mNotInMeshNearestElementWeights[i] <<
"\n";
#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)
bool IsOwned(Node< DIM > *pNode)
unsigned CalculateContainingBox(Node< DIM > *pNode)
bool IsBoxOwned(unsigned globalIndex)
void SetupAllLocalBoxes()
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 DeleteCoarseBoxCollection()
void ComputeCoarseElementsForFineNodes(bool safeMode)
void ComputeCoarseElementsForFineElementCentroids(bool safeMode)
void CollectElementsInLocalBoxes(DistributedBoxCollection< DIM > *&rpBoxCollection, unsigned boxIndex, std::set< unsigned > &rElementIndices)
void ResetStatisticsVariables()
void DeleteFineBoxCollection()
void ComputeFineElementAndWeightForGivenPoint(ChastePoint< DIM > &rPoint, bool safeMode, unsigned boxForThisPoint, unsigned index)
void ShareFineElementData()
void ShareCoarseElementData()
c_vector< double, SPACE_DIM > & rGetModifiableLocation()
const c_vector< double, SPACE_DIM > & rGetLocation() const
c_vector< double, DIM > & rGet(unsigned elementIndex, unsigned quadIndex)