CArl
Code Arlequin / C++ implementation
carl::Intersection_Search Class Reference

Class containing the structure needed to find all the intersections between two meshes, inside the coupling region mesh. More...

#include <intersection_search.h>

Public Member Functions

 Intersection_Search (libMesh::Mesh &mesh_A, libMesh::Mesh &mesh_B, libMesh::Mesh &mesh_Coupling, libMesh::Mesh &mesh_I, const std::string &output_base=std::string("test"), IntersectionMeshingMethod MeshingMethod=IntersectionMeshingMethod::CGAL, double Min_Inter_Volume=1E-15, bool bDoPerf_log=true, bool bDebugOutput=false)
 Constructor. More...
 
libMesh::Mesh & mesh_A ()
 Reference to the mesh A. More...
 
libMesh::Mesh & mesh_B ()
 Reference to the mesh B. More...
 
libMesh::Mesh & mesh_Coupling ()
 Reference to the coupling mesh. More...
 
void PreparePreallocationAndLoad (SearchMethod search_type=BRUTE)
 Do a preliminary search, to optimize the intersection search and construction. More...
 
void PreallocateAndPartitionCoupling ()
 Preallocate Intersection_Search::m_Intersection_Pairs_multimap and repartition the coupling mesh. More...
 
void BuildIntersections (SearchMethod search_type=BRUTE)
 Find and build all the intersections. More...
 
void SetScalingFiles (const std::string &timing_data_file_base)
 Set up timing output file. More...
 
void SkipIntersectionConstruction (bool bSkipIntersectionConstruction)
 Set the "Skip intersection construction" flag. More...
 
void SkipIntersectionPartitioning (bool bSkipIntersectionPartitioning)
 Set the "Skip the coupling mesh repartition" flag. More...
 
void CalculateGlobalVolume ()
 Calculate the total intersection volume over all the processors. More...
 

Protected Member Functions

void BuildCoupledPatches (const libMesh::Elem *Query_elem, int patch_counter=0)
 Build both patches associated a given query element from the coupling region mesh. More...
 
void FindPatchIntersections_Brute (const libMesh::Elem *Query_elem)
 Find all the intersections between the patches associated to the Query_elem, using a brute force method. All elements from a patch are tested against all the elements from the other patch. More...
 
void FindPatchIntersections_Front (const libMesh::Elem *Query_elem)
 Find all the intersections between the patches associated to the Query_elem, using an advancing front method. More...
 
void FindFirstPair (Patch_construction *Patch_guide, Patch_construction *Patch_probed, std::pair< unsigned int, unsigned int > &First_intersection)
 Find the first intersecting pair from a patch. More...
 
void BruteForce_FindFirstPair (Patch_construction *Patch_guide, Patch_construction *Patch_probed, std::pair< unsigned int, unsigned int > &First_intersection)
 Find the first intersecting pair from a patch, doing a full scan of the patches. More...
 
void FindAndBuildIntersections_Brute ()
 Find and build all the intersections, using the brute force method. More...
 
void FindIntersections_Brute ()
 Find all the intersections, using the brute force method. More...
 
void FindAndBuildIntersections_Front ()
 Find and build all the intersections, using the advancing front method. More...
 
void FindIntersections_Front ()
 Find all the intersections, using the advancing front method. More...
 
void CalculateIntersectionVolume (const libMesh::Elem *Query_elem)
 Calculate the volume of the intersections associated to Query_elem without updating the intersection mesh. Currently unused More...
 
void UpdateCouplingIntersection (const libMesh::Elem *Query_elem)
 Build the intersections associated to the Query_elem and update the intersection mesh. More...
 

Protected Attributes

libMesh::Mesh & m_Mesh_A
 
libMesh::Mesh & m_Mesh_B
 
libMesh::Mesh & m_Mesh_Coupling
 Mesh representing the coupling region. More...
 
const libMesh::Parallel::Communicator & m_comm
 Global communicator. More...
 
const unsigned int m_nodes
 Number of MPI processors. More...
 
const unsigned int m_rank
 This processor's rank (or ID in the communicator) More...
 
const libMesh::Parallel::Communicator & m_local_comm
 Local communicator, used for mesh patches. More...
 
Patch_construction m_Patch_Constructor_A
 Patch_construction object for mesh A. More...
 
Patch_construction m_Patch_Constructor_B
 Patch_construction object for mesh B. More...
 
IntersectionMeshingMethod m_MeshingMethod
 Flag defining the intersection meshing algorithm. More...
 
Mesh_Intersection m_Mesh_Intersection
 Object containing the intersection mesh. More...
 
std::unordered_multimap< unsigned int, unsigned int > m_Intersection_Pairs_multimap
 Multimap containing the intersection pairs. More...
 
Intersection_Tools m_Intersection_test
 Intersection operations for the main intersection tests. More...
 
Intersection_Tools m_Intersection_test_neighbors
 Intersection operations for the neighbor intersection tests. (why do we need both?) More...
 
std::vector< unsigned int > m_Nb_Of_Intersections_Elem_C
 Vector containing the number of intersections found inside each of the coupling mesh elements. More...
 
libMesh::ErrorVector m_coupling_weights
 libMesh::ErrorVector used to repartition the coupling mesh. More...
 
bool m_bSaveInterData
 Boolean flag determining if we should save the intersection data or not. Default: true. More...
 
bool m_bDidPreliminarySearch
 Boolean flag determining if we did a preliminary intersection search, used for the coupling mesh repartition and memory preallocations. More...
 
bool m_bHavePreallocData
 Boolean flag determining if we have the data for a proper preallocation of m_Intersection_Pairs_multimap. More...
 
bool m_bIntersectionsBuilt
 Boolean flag indicating if the intersections were built. (A bit useless right now ...) More...
 
double m_Min_Inter_Volume
 Volume cutoff for the intersection polyhedrons. Default: 1E-15. More...
 
bool m_bSkipIntersectionConstruction
 Boolean indicating if we should skip the intersection construction. Default: false. More...
 
bool m_bSkipIntersectionPartitioning
 Boolean indicating if we should skip the intersection repartitioning. Default: false. More...
 
std::string m_Output_filename_base
 Output filenames base. More...
 
bool MASTER_bPerfLog_intersection_search
 Do performance log? Default: true. More...
 
libMesh::PerfLog m_perf_log
 libMesh::PerfLog object. More...
 
bool m_bPrintDebug
 Print debug data. Default: false. More...
 
bool m_bPrintTimingData
 Print timing data. Default: false. More...
 
bool m_bPrintIntersectionsPerPartData
 Print intersections per partition. Default: false. More...
 
std::string m_timing_data_file_base
 Output filenames base for the timing data. More...
 

Detailed Description

Class containing the structure needed to find all the intersections between two meshes, inside the coupling region mesh.

Definition at line 28 of file intersection_search.h.

Constructor & Destructor Documentation

carl::Intersection_Search::Intersection_Search ( libMesh::Mesh &  mesh_A,
libMesh::Mesh &  mesh_B,
libMesh::Mesh &  mesh_Coupling,
libMesh::Mesh &  mesh_I,
const std::string &  output_base = std::string("test"),
IntersectionMeshingMethod  MeshingMethod = IntersectionMeshingMethod::CGAL,
double  Min_Inter_Volume = 1E-15,
bool  bDoPerf_log = true,
bool  bDebugOutput = false 
)
inline

Constructor.

Definition at line 179 of file intersection_search.h.

187  :
188  m_Mesh_A { mesh_A },
189  m_Mesh_B { mesh_B },
191  m_comm { m_Mesh_Coupling.comm() },
192  m_nodes { m_comm.size() },
193  m_rank { m_comm.rank() },
194  m_local_comm { mesh_I.comm() },
195  m_Patch_Constructor_A { Patch_construction(m_Mesh_A,m_local_comm)},
196  m_Patch_Constructor_B { Patch_construction(m_Mesh_B,m_local_comm)},
197  m_MeshingMethod { MeshingMethod },
198  m_Mesh_Intersection { Mesh_Intersection(mesh_I,m_Mesh_A,m_Mesh_B,m_MeshingMethod)},
199  m_bSaveInterData { true },
200  m_bDidPreliminarySearch { false },
201  m_bHavePreallocData { false },
202  m_bIntersectionsBuilt { false },
203  m_Min_Inter_Volume { Min_Inter_Volume },
206  m_Output_filename_base { output_base + "_r_" + std::to_string(m_rank) + "_n_" + std::to_string(m_nodes)},
208  m_perf_log { libMesh::PerfLog("Intersection search", MASTER_bPerfLog_intersection_search) },
209  m_bPrintDebug { bDebugOutput },
210  m_bPrintTimingData { false },
212  {
213  // Reserve space for the intersection multimap
214  m_Nb_Of_Intersections_Elem_C.resize(mesh_Coupling.n_elem(),0);
215  };
bool m_bDidPreliminarySearch
Boolean flag determining if we did a preliminary intersection search, used for the coupling mesh repa...
bool m_bIntersectionsBuilt
Boolean flag indicating if the intersections were built. (A bit useless right now ...
double m_Min_Inter_Volume
Volume cutoff for the intersection polyhedrons. Default: 1E-15.
bool m_bSkipIntersectionConstruction
Boolean indicating if we should skip the intersection construction. Default: false.
const libMesh::Parallel::Communicator & m_local_comm
Local communicator, used for mesh patches.
std::string m_Output_filename_base
Output filenames base.
bool m_bPrintIntersectionsPerPartData
Print intersections per partition. Default: false.
libMesh::Mesh & m_Mesh_Coupling
Mesh representing the coupling region.
const libMesh::Parallel::Communicator & m_comm
Global communicator.
const unsigned int m_nodes
Number of MPI processors.
IntersectionMeshingMethod m_MeshingMethod
Flag defining the intersection meshing algorithm.
bool m_bSkipIntersectionPartitioning
Boolean indicating if we should skip the intersection repartitioning. Default: false.
Patch_construction m_Patch_Constructor_B
Patch_construction object for mesh B.
libMesh::PerfLog m_perf_log
libMesh::PerfLog object.
Patch_construction m_Patch_Constructor_A
Patch_construction object for mesh A.
std::vector< unsigned int > m_Nb_Of_Intersections_Elem_C
Vector containing the number of intersections found inside each of the coupling mesh elements...
bool m_bHavePreallocData
Boolean flag determining if we have the data for a proper preallocation of m_Intersection_Pairs_multi...
Mesh_Intersection m_Mesh_Intersection
Object containing the intersection mesh.
bool m_bPrintDebug
Print debug data. Default: false.
libMesh::Mesh & mesh_A()
Reference to the mesh A.
libMesh::Mesh & mesh_Coupling()
Reference to the coupling mesh.
bool MASTER_bPerfLog_intersection_search
Do performance log? Default: true.
bool m_bSaveInterData
Boolean flag determining if we should save the intersection data or not. Default: true...
const unsigned int m_rank
This processor's rank (or ID in the communicator)
libMesh::Mesh & mesh_B()
Reference to the mesh B.
bool m_bPrintTimingData
Print timing data. Default: false.

Member Function Documentation

void carl::Intersection_Search::BruteForce_FindFirstPair ( Patch_construction Patch_guide,
Patch_construction Patch_probed,
std::pair< unsigned int, unsigned int > &  First_intersection 
)
protected

Find the first intersecting pair from a patch, doing a full scan of the patches.

This is, essentially, a brute force algorithm set to stop after the first positive test).

Definition at line 197 of file intersection_search.cpp.

200  {
201  bool bFoundFirstInter = false;
202  std::unordered_set<unsigned int> & Patch_Set_Guide = Patch_guide->elem_indexes();
203  std::unordered_set<unsigned int> & Patch_Set_Probed = Patch_probed->elem_indexes();
204 
205  libMesh::ReplicatedMesh& Mesh_patch_guide = Patch_guide->patch_mesh();
206  libMesh::ReplicatedMesh& Mesh_patch_probed = Patch_probed->patch_mesh();
207 
208  std::unordered_set<unsigned int>::iterator it_patch_Guide;
209  std::unordered_set<unsigned int>::iterator it_patch_Probed;
210 
211  for( it_patch_Guide = Patch_Set_Guide.begin();
212  it_patch_Guide != Patch_Set_Guide.end();
213  ++it_patch_Guide)
214  {
215  const libMesh::Elem * elem_Guide = Mesh_patch_guide.elem(*it_patch_Guide);
216 
217  for( it_patch_Probed = Patch_Set_Probed.begin();
218  it_patch_Probed != Patch_Set_Probed.end();
219  ++it_patch_Probed)
220  {
221  const libMesh::Elem * elem_Probed = Mesh_patch_probed.elem(*it_patch_Probed);
222 
223  bFoundFirstInter = m_Intersection_test.libMesh_exact_do_intersect_inside_coupling(elem_Guide,elem_Probed);
224 
225  if(bFoundFirstInter)
226  {
227  First_intersection.first = *it_patch_Guide;
228  First_intersection.second = *it_patch_Probed;
229  return;
230  }
231  }
232  }
233 
234  homemade_assert_msg(bFoundFirstInter,"Could't find a first intersecting pair. Are you sure that the patches do intersect?\n");
235  }
Intersection_Tools m_Intersection_test
Intersection operations for the main intersection tests.
bool libMesh_exact_do_intersect_inside_coupling(const libMesh::Elem *elem_A, const libMesh::Elem *elem_B, bool bCreateNewNefForA=true)
Determinate if two elements intersect inside the coupling region.
#define homemade_assert_msg(asserted, msg)
Definition: common_header.h:63
void carl::Intersection_Search::BuildCoupledPatches ( const libMesh::Elem *  Query_elem,
int  patch_counter = 0 
)
protected

Build both patches associated a given query element from the coupling region mesh.

Calls the method Patch_construction::BuildPatch for both patches, and export the patch meshes if Intersection_Search::m_bPrintDebug == true.

Definition at line 30 of file intersection_search.cpp.

31  {
32  // Unbreakable Patches!
34 
35  // Trusty Patches!
37 
38  if(m_bPrintDebug)
39  {
40  std::string filename = "meshes/3D/tests/output/patch_mesh_A_" + std::to_string(patch_counter) + "_" + std::to_string(m_rank);
42  filename = "meshes/3D/tests/output/patch_mesh_B_" + std::to_string(patch_counter) + "_" + std::to_string(m_rank);
44  }
45  }
void export_patch_mesh(std::string &filename_base)
Export the patch mesh to a file.
Patch_construction m_Patch_Constructor_B
Patch_construction object for mesh B.
Patch_construction m_Patch_Constructor_A
Patch_construction object for mesh A.
void BuildPatch(const libMesh::Elem *Query_elem)
Build the patch mesh covering a given "Query_elem".
bool m_bPrintDebug
Print debug data. Default: false.
const unsigned int m_rank
This processor's rank (or ID in the communicator)
void carl::Intersection_Search::BuildIntersections ( SearchMethod  search_type = BRUTE)

Find and build all the intersections.

Definition at line 750 of file intersection_search.cpp.

751  {
753  {
754  // Must do an preeeety expensive preallocation. Ouch ...
755  m_Intersection_Pairs_multimap.reserve(m_Mesh_A.n_elem()*m_Mesh_B.n_elem());
756  }
757 
758  m_bSaveInterData = true;
759  m_bIntersectionsBuilt = false;
760 
761  switch (search_type)
762  {
764  break;
765 
767  break;
768 
771  break;
772  }
773 
775  {
777  m_bIntersectionsBuilt = true;
778  }
779  }
bool m_bIntersectionsBuilt
Boolean flag indicating if the intersections were built. (A bit useless right now ...
bool m_bSkipIntersectionConstruction
Boolean indicating if we should skip the intersection construction. Default: false.
void FindAndBuildIntersections_Front()
Find and build all the intersections, using the advancing front method.
std::string m_Output_filename_base
Output filenames base.
void export_intersection_data(const std::string &filename_base, const std::string &mesh_format=std::string(".e"))
Export both the intersection mesh its data structures.
std::unordered_multimap< unsigned int, unsigned int > m_Intersection_Pairs_multimap
Multimap containing the intersection pairs.
bool m_bHavePreallocData
Boolean flag determining if we have the data for a proper preallocation of m_Intersection_Pairs_multi...
Mesh_Intersection m_Mesh_Intersection
Object containing the intersection mesh.
bool m_bSaveInterData
Boolean flag determining if we should save the intersection data or not. Default: true...
void FindAndBuildIntersections_Brute()
Find and build all the intersections, using the brute force method.
void carl::Intersection_Search::CalculateGlobalVolume ( )

Calculate the total intersection volume over all the processors.

Definition at line 784 of file intersection_search.cpp.

785  {
786  // Use a barrier to guarantee that all procs are in the same position
787  m_comm.barrier();
788 
789  // Calculate the volume of the intersection on each processor
790  double global_volume = m_Mesh_Intersection.get_total_volume();
791 
792  // Add it!
793  m_comm.sum(global_volume);
794 
795  // Calculate the volume of the coupling region on each processor
796  libMesh::Mesh::const_element_iterator it_coupl = m_Mesh_Coupling.local_elements_begin();
797  libMesh::Mesh::const_element_iterator it_coupl_end = m_Mesh_Coupling.local_elements_end();
798  double real_volume = 0;
799 
800  for( ; it_coupl != it_coupl_end; ++it_coupl)
801  {
802  const libMesh::Elem * Query_elem = * it_coupl;
803  real_volume += Query_elem->volume();
804  }
805 
806  // Add it!
807  m_comm.sum(real_volume);
808 
809  std::cout << " TOTAL volume, proc. " << m_rank << std::endl;
810  std::cout << " -> Intersection volume / real : " << global_volume << " / " << real_volume << std::endl << std::endl;
811 
812  }
double get_total_volume()
Calculate the intersection mesh's total volume.
libMesh::Mesh & m_Mesh_Coupling
Mesh representing the coupling region.
const libMesh::Parallel::Communicator & m_comm
Global communicator.
Mesh_Intersection m_Mesh_Intersection
Object containing the intersection mesh.
const unsigned int m_rank
This processor's rank (or ID in the communicator)
void carl::Intersection_Search::CalculateIntersectionVolume ( const libMesh::Elem *  Query_elem)
protected

Calculate the volume of the intersections associated to Query_elem without updating the intersection mesh. Currently unused

Definition at line 817 of file intersection_search.cpp.

818  {
819  std::unordered_set<unsigned int> & Patch_Set_A = m_Patch_Constructor_A.elem_indexes();
820 
821  std::unordered_set<unsigned int>::iterator it_patch_A;
822 
823  bool bDoIntersect = true;
824  bool bCreateNewNefForA = true;
825 
827 
828  std::set<libMesh::Point> points_out;
829  double total_volume = 0;
830 
831  // Debug vars
832  int nbOfTests = 0;
833  int nbOfPositiveTests = 0;
834  for( it_patch_A = Patch_Set_A.begin();
835  it_patch_A != Patch_Set_A.end();
836  ++it_patch_A)
837  {
838  auto iterator_pair = m_Intersection_Pairs_multimap.equal_range(*it_patch_A);
839 
840  const libMesh::Elem * elem_A = m_Mesh_A.elem(*it_patch_A);
841  bCreateNewNefForA = true;
842 
843  for(auto it_patch_B = iterator_pair.first ;
844  it_patch_B != iterator_pair.second ;
845  ++it_patch_B )
846  {
847  const libMesh::Elem * elem_B = m_Mesh_B.elem(it_patch_B->second);
848  points_out.clear();
849 
850  m_perf_log.push("Build intersection polyhedron","Build intersections");
851  bDoIntersect = m_Intersection_test.libMesh_exact_intersection_inside_coupling(elem_A,elem_B,points_out,bCreateNewNefForA,true,false);
852  m_perf_log.pop("Build intersection polyhedron","Build intersections");
853  ++nbOfTests;
854 
855  if(bDoIntersect)
856  {
857  bCreateNewNefForA = false;
858  m_perf_log.push("Calculate volume","Build intersections");
859  total_volume += m_Mesh_Intersection.get_intersection_volume(points_out);
860  m_perf_log.pop("Calculate volume","Build intersections");
861  ++nbOfPositiveTests;
862  }
863  }
864  }
865 
866  if(m_bPrintDebug)
867  {
868  std::cout << " DEBUG: calculate the volume" << std::endl;
869  std::cout << " -> Intersection volume / real : " << total_volume << " / " << Query_elem->volume() << std::endl << std::endl;
870  std::cout << " -> Positives / tests : " << nbOfPositiveTests << " / " << nbOfTests
871  << " (" << 100.*nbOfPositiveTests/nbOfTests << "%)" << std::endl << std::endl;
872  }
873  }
void libmesh_set_coupling_nef_polyhedron(const libMesh::Elem *elem_C)
Set the Nef polyhedron for the coupling mesh element.
std::unordered_set< unsigned int > & elem_indexes()
Returns the set of elements inside the patch.
libMesh::PerfLog m_perf_log
libMesh::PerfLog object.
Patch_construction m_Patch_Constructor_A
Patch_construction object for mesh A.
std::unordered_multimap< unsigned int, unsigned int > m_Intersection_Pairs_multimap
Multimap containing the intersection pairs.
Mesh_Intersection m_Mesh_Intersection
Object containing the intersection mesh.
Intersection_Tools m_Intersection_test
Intersection operations for the main intersection tests.
bool m_bPrintDebug
Print debug data. Default: false.
double get_intersection_volume(std::set< libMesh::Point > &input_points)
Calculate the volume of an polyhedron defined by the point set.
bool libMesh_exact_intersection_inside_coupling(const libMesh::Elem *elem_A, const libMesh::Elem *elem_B, std::set< libMesh::Point > &points_out, bool bCreateNewNefForA=true, bool bConvertPoints=true, bool bTestNeeded=false)
Build the intersection between two elements intersect inside the coupling region. ...
void carl::Intersection_Search::FindAndBuildIntersections_Brute ( )
protected

Find and build all the intersections, using the brute force method.

Skips build step if Intersection_Search::m_bSkipIntersectionConstruction == true.

Definition at line 438 of file intersection_search.cpp.

439  {
440  // Prepare iterators
441  libMesh::Mesh::const_element_iterator it_coupl = m_Mesh_Coupling.local_elements_begin();
442  libMesh::Mesh::const_element_iterator it_coupl_end = m_Mesh_Coupling.local_elements_end();
443 
444  int patch_counter = 0;
445  double real_volume = 0;
446 
448 
449  for( ; it_coupl != it_coupl_end; ++it_coupl)
450  {
451  const libMesh::Elem * Query_elem = * it_coupl;
452 
453  m_perf_log.push("Build patches","Brute force");
454  BuildCoupledPatches(Query_elem,patch_counter);
455  m_perf_log.pop("Build patches","Brute force");
456 
457  ++patch_counter;
458 
459  if(m_bPrintDebug)
460  {
461  real_volume += Query_elem->volume();
462  }
463 
464  m_perf_log.push("Find intersections","Brute force");
465  FindPatchIntersections_Brute(Query_elem);
466  m_perf_log.pop("Find intersections","Brute force");
468  {
469  m_perf_log.push("Build intersections","Brute force");
470  UpdateCouplingIntersection(Query_elem);
471  m_perf_log.pop("Build intersections","Brute force");
472  }
473  };
474 
476  {
477  m_perf_log.push("Prepare mesh","Brute force");
479  m_perf_log.pop("Prepare mesh","Brute force");
480  }
481 
483  {
484  std::vector<double> timing_find_intersections(m_nodes,0);
485  std::vector<double> timing_build_intersections(m_nodes,0);
486 
487  libMesh::PerfData performance_data = m_perf_log.get_perf_data("Find intersections","Brute force");
488  timing_find_intersections[m_rank] = performance_data.tot_time_incl_sub;
489  performance_data = m_perf_log.get_perf_data("Build intersections","Brute force");
490  timing_build_intersections[m_rank] = performance_data.tot_time_incl_sub;
491 
492  m_comm.sum(timing_find_intersections);
493  m_comm.sum(timing_build_intersections);
494 
495  if(m_rank == 0)
496  {
497  print_stats_to_file(timing_find_intersections,m_timing_data_file_base + "_search_brute.dat");
499  {
500  print_stats_to_file(timing_build_intersections,m_timing_data_file_base + "_build_brute.dat");
501  }
502  }
503  }
504 
506  {
507  m_perf_log.push("Calculate volume","Brute force");
508  double total_volume = m_Mesh_Intersection.get_total_volume();
509  m_perf_log.push("Calculate volume","Brute force");
510 
511  std::cout << " DEBUG volume on proc. " << m_rank << "(BRUTE)" << std::endl;
512  std::cout << " -> Mesh elems, nodes : " << m_Mesh_Intersection.mesh().n_elem() << " , " << m_Mesh_Intersection.mesh().n_nodes() << std::endl;
513  std::cout << " -> Intersection volume / real : " << total_volume << " / " << real_volume << std::endl << std::endl;
514  }
515  }
std::string m_timing_data_file_base
Output filenames base for the timing data.
double get_total_volume()
Calculate the intersection mesh's total volume.
bool m_bSkipIntersectionConstruction
Boolean indicating if we should skip the intersection construction. Default: false.
const libMesh::ReplicatedMesh & mesh()
Returns the address of the intersection mesh.
void FindPatchIntersections_Brute(const libMesh::Elem *Query_elem)
Find all the intersections between the patches associated to the Query_elem, using a brute force meth...
libMesh::Mesh & m_Mesh_Coupling
Mesh representing the coupling region.
const libMesh::Parallel::Communicator & m_comm
Global communicator.
const unsigned int m_nodes
Number of MPI processors.
void BuildCoupledPatches(const libMesh::Elem *Query_elem, int patch_counter=0)
Build both patches associated a given query element from the coupling region mesh.
void print_stats_to_file(std::vector< double > &vec_data, const std::string filename)
void UpdateCouplingIntersection(const libMesh::Elem *Query_elem)
Build the intersections associated to the Query_elem and update the intersection mesh.
libMesh::PerfLog m_perf_log
libMesh::PerfLog object.
Mesh_Intersection m_Mesh_Intersection
Object containing the intersection mesh.
bool m_bPrintDebug
Print debug data. Default: false.
void initialize()
Initialize the data structures.
void prepare_for_use()
Prepare the intersection mesh for use.
const unsigned int m_rank
This processor's rank (or ID in the communicator)
bool m_bPrintTimingData
Print timing data. Default: false.
void carl::Intersection_Search::FindAndBuildIntersections_Front ( )
protected

Find and build all the intersections, using the advancing front method.

Skips build step if Intersection_Search::m_bSkipIntersectionConstruction == true.

Definition at line 545 of file intersection_search.cpp.

546  {
547  // Prepare iterators
548  libMesh::Mesh::const_element_iterator it_coupl = m_Mesh_Coupling.local_elements_begin();
549  libMesh::Mesh::const_element_iterator it_coupl_end = m_Mesh_Coupling.local_elements_end();
550 
551  int patch_counter = 0;
552  double real_volume = 0;
553 
555 
556  for( ; it_coupl != it_coupl_end; ++it_coupl)
557  {
558  const libMesh::Elem * Query_elem = * it_coupl;
559 
560  m_perf_log.push("Build patches","Advancing front");
561  BuildCoupledPatches(Query_elem,patch_counter);
562  m_perf_log.pop("Build patches","Advancing front");
563 
564  ++patch_counter;
565 
566  if(m_bPrintDebug)
567  {
568  real_volume += Query_elem->volume();
569  }
570  m_perf_log.push("Find intersections","Advancing front");
571  FindPatchIntersections_Front(Query_elem);
572  m_perf_log.pop("Find intersections","Advancing front");
574  {
575  m_perf_log.push("Build intersections","Advancing front");
576  UpdateCouplingIntersection(Query_elem);
577  m_perf_log.pop("Build intersections","Advancing front");
578  }
579  };
580 
582  {
583  m_perf_log.push("Prepare mesh","Advancing front");
585  m_perf_log.pop("Prepare mesh","Advancing front");
586  }
587 
589  {
590  std::vector<double> timing_find_intersections(m_nodes,0);
591  std::vector<double> timing_build_intersections(m_nodes,0);
592 
593  libMesh::PerfData performance_data = m_perf_log.get_perf_data("Find intersections","Advancing front");
594  timing_find_intersections[m_rank] = performance_data.tot_time_incl_sub;
595  performance_data = m_perf_log.get_perf_data("Build intersections","Advancing front");
596  timing_build_intersections[m_rank] = performance_data.tot_time_incl_sub;
597 
598  m_comm.sum(timing_find_intersections);
599  m_comm.sum(timing_build_intersections);
600 
601  if(m_rank == 0)
602  {
603  print_stats_to_file(timing_find_intersections,m_timing_data_file_base + "_search_advancing.dat");
605  {
606  print_stats_to_file(timing_build_intersections,m_timing_data_file_base + "_build_advancing.dat");
607  }
608  }
609  }
610 
612  {
613  m_perf_log.push("Calculate volume","Advancing front");
614  double total_volume = m_Mesh_Intersection.get_total_volume();
615  m_perf_log.push("Calculate volume","Advancing front");
616 
617  std::cout << " DEBUG volume on proc. " << m_rank << "(FRONT)" << std::endl;
618  std::cout << " -> Mesh elems, nodes : " << m_Mesh_Intersection.mesh().n_elem() << " , " << m_Mesh_Intersection.mesh().n_nodes() << std::endl;
619  std::cout << " -> Intersection volume / real : " << total_volume << " / " << real_volume << std::endl << std::endl;
620  }
621  }
std::string m_timing_data_file_base
Output filenames base for the timing data.
double get_total_volume()
Calculate the intersection mesh's total volume.
bool m_bSkipIntersectionConstruction
Boolean indicating if we should skip the intersection construction. Default: false.
const libMesh::ReplicatedMesh & mesh()
Returns the address of the intersection mesh.
libMesh::Mesh & m_Mesh_Coupling
Mesh representing the coupling region.
const libMesh::Parallel::Communicator & m_comm
Global communicator.
const unsigned int m_nodes
Number of MPI processors.
void BuildCoupledPatches(const libMesh::Elem *Query_elem, int patch_counter=0)
Build both patches associated a given query element from the coupling region mesh.
void print_stats_to_file(std::vector< double > &vec_data, const std::string filename)
void UpdateCouplingIntersection(const libMesh::Elem *Query_elem)
Build the intersections associated to the Query_elem and update the intersection mesh.
libMesh::PerfLog m_perf_log
libMesh::PerfLog object.
Mesh_Intersection m_Mesh_Intersection
Object containing the intersection mesh.
bool m_bPrintDebug
Print debug data. Default: false.
void initialize()
Initialize the data structures.
void FindPatchIntersections_Front(const libMesh::Elem *Query_elem)
Find all the intersections between the patches associated to the Query_elem, using an advancing front...
void prepare_for_use()
Prepare the intersection mesh for use.
const unsigned int m_rank
This processor's rank (or ID in the communicator)
bool m_bPrintTimingData
Print timing data. Default: false.
void carl::Intersection_Search::FindFirstPair ( Patch_construction Patch_guide,
Patch_construction Patch_probed,
std::pair< unsigned int, unsigned int > &  First_intersection 
)
protected

Find the first intersecting pair from a patch.

Do so by :

  1. locating the elements from the guide patch that contain the vertices of an element from the probed patch.
  2. for each one of these guide elements, find the one that intersects the probed element inside the coupling.
  3. if this fails, use a brute force search method (call Intersection_Search::BruteForce_FindFirstPair).

Definition at line 126 of file intersection_search.cpp.

129  {
130  // Set up some references for a simpler code
131  libMesh::ReplicatedMesh& Mesh_patch_guide = Patch_guide->patch_mesh();
132  libMesh::ReplicatedMesh& Mesh_patch_probed = Patch_probed->patch_mesh();
133 
134  // Set up locator for the first intersecting pair
135  std::unique_ptr<libMesh::PointLocatorBase> Patch_guide_Locator = Mesh_patch_guide.sub_point_locator();
136 
137  // Instruction needed to avoid the code from crashing if a query is
138  // outside the mesh
139  Patch_guide_Locator->enable_out_of_mesh_mode();
140 
141  // Find the first pair
142  libMesh::Elem * Patch_probed_first_elem = NULL;
143  unsigned int Patch_probed_first_elem_id = 0;
144 
145  // Set of intersecting terms:
146  std::set<unsigned int> Intersecting_guide_elems;
147  bool bFoundIntersection = false;
148  for(libMesh::ReplicatedMesh::element_iterator element_probed_it = Mesh_patch_probed.elements_begin();
149  element_probed_it != Mesh_patch_probed.elements_end();
150  ++element_probed_it)
151  {
152  Patch_probed_first_elem = * element_probed_it;
153  Patch_probed_first_elem_id = Patch_probed->convert_patch_to_parent_elem_id(Patch_probed_first_elem->id());
154  bFoundIntersection = m_Intersection_test.FindAllIntersection(Patch_probed_first_elem,Patch_guide_Locator,Intersecting_guide_elems);
155  if(bFoundIntersection)
156  {
157  break;
158  }
159  }
160 
161  // Search for one intersecting term inside the guide patch
162  libMesh::Elem * Patch_guide_first_elem = NULL;
163  unsigned int Patch_guide_first_elem_id;
164  bool bFoundFirstInter = false;
165 
166  for(std::set<unsigned int>::iterator it_inter = Intersecting_guide_elems.begin();
167  it_inter != Intersecting_guide_elems.end();
168  ++it_inter)
169  {
170  Patch_guide_first_elem = Mesh_patch_guide.elem(*it_inter);
171  Patch_guide_first_elem_id = Patch_guide->convert_patch_to_parent_elem_id(*it_inter);
172 
173  // Must test if the intersection is inside the coupling region
174  bFoundFirstInter = m_Intersection_test_neighbors.libMesh_exact_do_intersect_inside_coupling(Patch_probed_first_elem,Patch_guide_first_elem);
175 
176  if(bFoundFirstInter)
177  {
178  First_intersection.first = Patch_guide_first_elem_id;
179  First_intersection.second = Patch_probed_first_elem_id;
180  break;
181  }
182  }
183 
184  if(!bFoundFirstInter)
185  {
186  // Couldn't find the first intersection using a point search ...
187  // Will have to do the things the hard way ...
188  BruteForce_FindFirstPair(Patch_guide,Patch_probed,First_intersection);
189  }
190  };
void BruteForce_FindFirstPair(Patch_construction *Patch_guide, Patch_construction *Patch_probed, std::pair< unsigned int, unsigned int > &First_intersection)
Find the first intersecting pair from a patch, doing a full scan of the patches.
bool FindAllIntersection(const libMesh::Elem *Query_elem, std::unique_ptr< libMesh::PointLocatorBase > &point_locator, std::set< unsigned int > &Intersecting_elems)
Find all elements from the mesh intersecting Query_elem, using the mesh's libMesh::PointLocatorBase.
Intersection_Tools m_Intersection_test
Intersection operations for the main intersection tests.
bool libMesh_exact_do_intersect_inside_coupling(const libMesh::Elem *elem_A, const libMesh::Elem *elem_B, bool bCreateNewNefForA=true)
Determinate if two elements intersect inside the coupling region.
Intersection_Tools m_Intersection_test_neighbors
Intersection operations for the neighbor intersection tests. (why do we need both?)
void carl::Intersection_Search::FindIntersections_Brute ( )
protected

Find all the intersections, using the brute force method.

Skips saving intersecting element pairs if Intersection_Search::m_bSaveInterData == true.

Definition at line 517 of file intersection_search.cpp.

518  {
519  // Prepare iterators
520  libMesh::Mesh::const_element_iterator it_coupl = m_Mesh_Coupling.local_elements_begin();
521  libMesh::Mesh::const_element_iterator it_coupl_end = m_Mesh_Coupling.local_elements_end();
522 
523  int patch_counter = 0;
524 
525  for( ; it_coupl != it_coupl_end; ++it_coupl)
526  {
527  const libMesh::Elem * Query_elem = * it_coupl;
528 
529  m_perf_log.push("Build patches","Brute force (preamble run)");
530  BuildCoupledPatches(Query_elem,patch_counter);
531  m_perf_log.pop("Build patches","Brute force (preamble run)");
532 
533  ++patch_counter;
534 
535  m_perf_log.push("Find intersections","Brute force (preamble run)");
536  FindPatchIntersections_Brute(Query_elem);
537  m_perf_log.pop("Find intersections","Brute force (preamble run)");
538  };
539  }
void FindPatchIntersections_Brute(const libMesh::Elem *Query_elem)
Find all the intersections between the patches associated to the Query_elem, using a brute force meth...
libMesh::Mesh & m_Mesh_Coupling
Mesh representing the coupling region.
void BuildCoupledPatches(const libMesh::Elem *Query_elem, int patch_counter=0)
Build both patches associated a given query element from the coupling region mesh.
libMesh::PerfLog m_perf_log
libMesh::PerfLog object.
void carl::Intersection_Search::FindIntersections_Front ( )
protected

Find all the intersections, using the advancing front method.

Skips saving intersecting element pairs if Intersection_Search::m_bSaveInterData == true.

Definition at line 623 of file intersection_search.cpp.

624  {
625  // Prepare iterators
626  libMesh::Mesh::const_element_iterator it_coupl = m_Mesh_Coupling.local_elements_begin();
627  libMesh::Mesh::const_element_iterator it_coupl_end = m_Mesh_Coupling.local_elements_end();
628 
629  int patch_counter = 0;
630 
631  for( ; it_coupl != it_coupl_end; ++it_coupl)
632  {
633  const libMesh::Elem * Query_elem = * it_coupl;
634 
635  m_perf_log.push("Build patches (preamble)","Advancing front (preamble run)");
636  BuildCoupledPatches(Query_elem,patch_counter);
637  m_perf_log.pop("Build patches","Advancing front (preamble run)");
638 
639  ++patch_counter;
640 
641  m_perf_log.push("Find intersections","Advancing front (preamble run)");
642  FindPatchIntersections_Front(Query_elem);
643  m_perf_log.pop("Find intersections","Advancing front (preamble run)");
644  };
645  }
libMesh::Mesh & m_Mesh_Coupling
Mesh representing the coupling region.
void BuildCoupledPatches(const libMesh::Elem *Query_elem, int patch_counter=0)
Build both patches associated a given query element from the coupling region mesh.
libMesh::PerfLog m_perf_log
libMesh::PerfLog object.
void FindPatchIntersections_Front(const libMesh::Elem *Query_elem)
Find all the intersections between the patches associated to the Query_elem, using an advancing front...
void carl::Intersection_Search::FindPatchIntersections_Brute ( const libMesh::Elem *  Query_elem)
protected

Find all the intersections between the patches associated to the Query_elem, using a brute force method. All elements from a patch are tested against all the elements from the other patch.

Definition at line 52 of file intersection_search.cpp.

53  {
54  m_perf_log.push("Preamble","Brute force algorithm");
55 
56  // Code for the brute force intersection tests
58 
59  // Intersection_Tools
61 
62  // Boolean: do they intersect?
63  bool bDoIntersect = false;
64 
65  std::unordered_set<unsigned int> & Patch_Set_A = m_Patch_Constructor_A.elem_indexes();
66  std::unordered_set<unsigned int> & Patch_Set_B = m_Patch_Constructor_B.elem_indexes();
67 
68  std::unordered_set<unsigned int>::iterator it_patch_A;
69  std::unordered_set<unsigned int>::iterator it_patch_B;
70 
71  // Debug vars
72  int nbOfTests = 0;
73  int nbOfPositiveTests = 0;
74 
75  m_perf_log.pop("Preamble","Brute force algorithm");
76 
77  for( it_patch_A = Patch_Set_A.begin();
78  it_patch_A != Patch_Set_A.end();
79  ++it_patch_A)
80  {
81  const libMesh::Elem * elem_A = m_Mesh_A.elem(*it_patch_A);
82 
83  for( it_patch_B = Patch_Set_B.begin();
84  it_patch_B != Patch_Set_B.end();
85  ++it_patch_B)
86  {
87  m_perf_log.push("Find intersection","Brute force algorithm");
88  ++nbOfTests;
89  const libMesh::Elem * elem_B = m_Mesh_B.elem(*it_patch_B);
90 
91  bDoIntersect = m_Intersection_test.libMesh_exact_do_intersect(elem_A,elem_B);
92  m_perf_log.pop("Find intersection","Brute force algorithm");
93 
94  if(bDoIntersect)
95  {
97  {
98  m_Intersection_Pairs_multimap.insert(std::pair<unsigned int, unsigned int>(*it_patch_A,*it_patch_B));
99  }
100  ++nbOfPositiveTests;
101  }
102  }
103  }
104 
105  // Partitioning vars
106  unsigned int query_idx = Query_elem->id();
107  m_Nb_Of_Intersections_Elem_C[query_idx] = nbOfPositiveTests;
108 
109  if(m_bPrintDebug)
110  {
111  std::cout << " DEBUG: brute force search results" << std::endl;
112  std::cout << " -> Positives / tests : " << nbOfPositiveTests << " / " << nbOfTests
113  << " (" << 100.*nbOfPositiveTests/nbOfTests << "%)" << std::endl << std::endl;
114  }
115  }
void libmesh_set_coupling_nef_polyhedron(const libMesh::Elem *elem_C)
Set the Nef polyhedron for the coupling mesh element.
Patch_construction m_Patch_Constructor_B
Patch_construction object for mesh B.
std::unordered_set< unsigned int > & elem_indexes()
Returns the set of elements inside the patch.
libMesh::PerfLog m_perf_log
libMesh::PerfLog object.
Patch_construction m_Patch_Constructor_A
Patch_construction object for mesh A.
std::unordered_multimap< unsigned int, unsigned int > m_Intersection_Pairs_multimap
Multimap containing the intersection pairs.
std::vector< unsigned int > m_Nb_Of_Intersections_Elem_C
Vector containing the number of intersections found inside each of the coupling mesh elements...
Intersection_Tools m_Intersection_test
Intersection operations for the main intersection tests.
bool m_bPrintDebug
Print debug data. Default: false.
bool libMesh_exact_do_intersect(const libMesh::Elem *elem_A, const libMesh::Elem *elem_B)
Determinate if two elements intersect.
bool m_bSaveInterData
Boolean flag determining if we should save the intersection data or not. Default: true...
void carl::Intersection_Search::FindPatchIntersections_Front ( const libMesh::Elem *  Query_elem)
protected

Find all the intersections between the patches associated to the Query_elem, using an advancing front method.

Definition at line 241 of file intersection_search.cpp.

242  {
243  m_perf_log.push("Preamble","Advancing front algorithm");
244  // Code for the front intersection tests
246 
247  // Intersection_Tools
250 
251  // --- Set up variables needed by Gander's algorithm
252 
253  // First, determinate which patch will be the guide, and which will be
254  // the probed
255  Patch_construction * Patch_guide = &m_Patch_Constructor_A;
256  Patch_construction * Patch_probed = &m_Patch_Constructor_B;
257 
258  // Exchange them if the guide is smaller
259  bool bGuidedByB = Patch_guide->size() > Patch_probed->size();
260  if(bGuidedByB)
261  {
262  if(m_bPrintDebug)
263  {
264  std::cout << " DEBUG: Probe: A | Guide: B" << std::endl;
265  }
266  Patch_guide = &m_Patch_Constructor_B;
267  Patch_probed = &m_Patch_Constructor_A;
268  }
269  else
270  {
271  if(m_bPrintDebug)
272  {
273  std::cout << " DEBUG: Probe: B | Guide: A" << std::endl;
274  }
275  }
276 
277  // Ids of the working elements
278  unsigned int Guide_working_elem_id;
279  unsigned int Probed_working_elem_id;
280 
281  // Find the first pair, (guide, probed)
282  std::pair<unsigned int,unsigned int> First_intersection;
283  BruteForce_FindFirstPair(Patch_guide,Patch_probed,First_intersection);
284 
285  // Index and iterators guide neighbor elements to be tested
286  unsigned int guide_neighbor_index = 0;
287  std::unordered_set<unsigned int>::iterator it_neigh, it_neigh_end;
288 
289  // Boolean saying if two elements intersect (exactly)
290  bool bDoIntersect = false;
291 
292  // Boolean saying if a neighbor element intersects (exactly)
293  bool bDoNeighborIntersect = false;
294 
295  // Boolean used to short-circuit the neighbor search test if all
296  // the concerned elements already have neighbours.
297  bool bDoTestGuideNeighbours = true;
298 
299  // --- Preamble - initializations
300  // Insert the first elements in the queues
301  Patch_guide->FrontSearch_initialize();
302  Patch_probed->FrontSearch_initialize();
303 
304  Patch_guide->intersection_queue_push_back(First_intersection.first);
305  Patch_probed->intersection_queue_push_back(First_intersection.second);
306  Patch_guide->set_elem_as_inside_queue(First_intersection.first);
307 
308  // Debug vars
309  int nbOfGuideElemsTested = 0;
310  int nbOfTests = 0;
311  int nbOfPositiveTests = 0;
312 
313  bool bCreateNewNefForGuideNeigh = true;
314 
315  m_perf_log.pop("Preamble","Advancing front algorithm");
316 
317  while(!Patch_guide->intersection_queue_empty())
318  {
319  m_perf_log.push("Set up new guide element","Advancing front algorithm");
320 
321  // Pop out the first element from the guide intersection queue
322  Guide_working_elem_id = Patch_guide->intersection_queue_extract_front_elem();
323  Patch_guide->set_elem_as_treated(Guide_working_elem_id);
324 
325  // Get its pointer
326  const libMesh::Elem * elem_guide = Patch_guide->current_elem_pointer();
327 
328  ++nbOfGuideElemsTested;
329 
330  // Set up the neighbor test short-circuits
331  bDoTestGuideNeighbours = Patch_guide->set_neighbors_to_search_next_pairs();
332 
333  // Clear and initialize the probed patch lists
334  Patch_probed->FrontSearch_reset();
335 
336  // Pop out the first element from the probed intersection queue, and
337  // insert it inside the test queue
338  Patch_probed->FrontSearch_prepare_for_probed_test();
339 
340  m_perf_log.pop("Set up new guide element","Advancing front algorithm");
341 
342  while(!Patch_probed->test_queue_empty())
343  {
344  m_perf_log.push("Find intersection","Advancing front algorithm");
345  // Pop out the first elem to be tested
346  Probed_working_elem_id = Patch_probed->test_queue_extract_front_elem();
347  Patch_probed->set_elem_as_treated(Probed_working_elem_id);
348 
349  const libMesh::Elem * elem_probed = Patch_probed->current_elem_pointer();
350 
351  // Test the intersection
352  bDoIntersect = m_Intersection_test.libMesh_exact_do_intersect(elem_guide,elem_probed);
353  ++nbOfTests;
354  m_perf_log.pop("Find intersection","Advancing front algorithm");
355 
356  // If they do intersect, we have to ...
357  if(bDoIntersect)
358  {
359  // 1) add elements to intersection multimap, watching out
360  // for the element order
361  if(m_bSaveInterData)
362  {
363  if(bGuidedByB)
364  {
365  m_Intersection_Pairs_multimap.insert(std::pair<unsigned int, unsigned int>(Probed_working_elem_id,Guide_working_elem_id));
366  }
367  else
368  {
369  m_Intersection_Pairs_multimap.insert(std::pair<unsigned int, unsigned int>(Guide_working_elem_id,Probed_working_elem_id));
370  }
371  }
372  ++nbOfPositiveTests;
373 
374  m_perf_log.push("Update test queue","Advancing front algorithm");
375  // 2) add elem_probed's neighbors, if not treated yet
376  Patch_probed->add_neighbors_to_test_list();
377  m_perf_log.pop("Update test queue","Advancing front algorithm");
378 
379  m_perf_log.push("Update intersection queue","Advancing front algorithm");
380  // 3) determinate if any of the guide neighbors are good
381  // candidates with the probed element.
382  if(bDoTestGuideNeighbours)
383  {
384  it_neigh = Patch_guide->neighbors_to_search_next_pair().begin();
385  it_neigh_end = Patch_guide->neighbors_to_search_next_pair().end();
386 
387  bCreateNewNefForGuideNeigh = true;
388 
389  for( ; it_neigh != it_neigh_end; ++it_neigh )
390  {
391  // This element doesn't have an intersecting pair
392  // yet. Test it, but only save it if it is inside
393  // the coupling region
394  guide_neighbor_index = *it_neigh;
395  const libMesh::Elem * elem_neigh = Patch_guide->elem(guide_neighbor_index);
396  bDoNeighborIntersect = m_Intersection_test_neighbors.libMesh_exact_do_intersect_inside_coupling(elem_probed,elem_neigh,bCreateNewNefForGuideNeigh);
397 
398  if(bDoNeighborIntersect)
399  {
400  bCreateNewNefForGuideNeigh = false;
401  }
402 
403  if(bDoNeighborIntersect)
404  {
405  // Now it has an intersecting pair!
406  Patch_guide->intersection_queue_push_back(guide_neighbor_index);
407  Patch_probed->intersection_queue_push_back(Probed_working_elem_id);
408 
409  // Set the guide elements as "already inside the queue"
410  Patch_guide->set_elem_as_inside_queue(guide_neighbor_index);
411  }
412  }
413 
414  bDoTestGuideNeighbours = Patch_guide->set_neighbors_to_search_next_pairs();
415  }
416  m_perf_log.pop("Update intersection queue","Advancing front algorithm");
417  }
418  }
419  }
420 
421  // Partitioning vars
422  unsigned int query_idx = Query_elem->id();
423  m_Nb_Of_Intersections_Elem_C[query_idx] = nbOfPositiveTests;
424 
425  if(m_bPrintDebug)
426  {
427  std::cout << " DEBUG: advancing front search results" << std::endl;
428  std::cout << " -> Guide elements tested / all : " << nbOfGuideElemsTested << " / " << Patch_guide->size() << std::endl;
429  std::cout << " -> Positives / tests : " << nbOfPositiveTests << " / " << nbOfTests
430  << " (" << 100.*nbOfPositiveTests/nbOfTests << "%)" << std::endl << std::endl;
431  }
432  }
void libmesh_set_coupling_nef_polyhedron(const libMesh::Elem *elem_C)
Set the Nef polyhedron for the coupling mesh element.
void BruteForce_FindFirstPair(Patch_construction *Patch_guide, Patch_construction *Patch_probed, std::pair< unsigned int, unsigned int > &First_intersection)
Find the first intersecting pair from a patch, doing a full scan of the patches.
Patch_construction m_Patch_Constructor_B
Patch_construction object for mesh B.
libMesh::PerfLog m_perf_log
libMesh::PerfLog object.
Patch_construction m_Patch_Constructor_A
Patch_construction object for mesh A.
std::unordered_multimap< unsigned int, unsigned int > m_Intersection_Pairs_multimap
Multimap containing the intersection pairs.
std::vector< unsigned int > m_Nb_Of_Intersections_Elem_C
Vector containing the number of intersections found inside each of the coupling mesh elements...
Intersection_Tools m_Intersection_test
Intersection operations for the main intersection tests.
bool m_bPrintDebug
Print debug data. Default: false.
bool libMesh_exact_do_intersect_inside_coupling(const libMesh::Elem *elem_A, const libMesh::Elem *elem_B, bool bCreateNewNefForA=true)
Determinate if two elements intersect inside the coupling region.
bool libMesh_exact_do_intersect(const libMesh::Elem *elem_A, const libMesh::Elem *elem_B)
Determinate if two elements intersect.
unsigned int size()
Returns the number of elements inside the patch.
bool m_bSaveInterData
Boolean flag determining if we should save the intersection data or not. Default: true...
Intersection_Tools m_Intersection_test_neighbors
Intersection operations for the neighbor intersection tests. (why do we need both?)
libMesh::Mesh & carl::Intersection_Search::mesh_A ( )

Reference to the mesh A.

Definition at line 12 of file intersection_search.cpp.

13  {
14  return m_Mesh_A;
15  }
libMesh::Mesh & carl::Intersection_Search::mesh_B ( )

Reference to the mesh B.

Definition at line 17 of file intersection_search.cpp.

18  {
19  return m_Mesh_B;
20  }
libMesh::Mesh & carl::Intersection_Search::mesh_Coupling ( )

Reference to the coupling mesh.

Definition at line 22 of file intersection_search.cpp.

23  {
24  return m_Mesh_Coupling;
25  }
libMesh::Mesh & m_Mesh_Coupling
Mesh representing the coupling region.
void carl::Intersection_Search::PreallocateAndPartitionCoupling ( )

Preallocate Intersection_Search::m_Intersection_Pairs_multimap and repartition the coupling mesh.

Definition at line 673 of file intersection_search.cpp.

674  {
676  {
677  // Redo the partitioning
679 
680  for(unsigned int iii = 0; iii < m_Nb_Of_Intersections_Elem_C.size(); ++iii)
681  {
683  }
684 
686  {
687  libMesh::Partitioner * dummy_partitioner = m_Mesh_Coupling.partitioner().get();
688  dummy_partitioner->attach_weights(&m_coupling_weights);
689  m_Mesh_Coupling.partition(m_nodes);
690  }
691 
692  // Allocate the intersection maps ...
693  std::vector<unsigned int> nb_of_inters_per_rank(m_nodes,0);
694  unsigned int dummy_nb_of_inters = 0;
695  libMesh::Mesh::const_element_iterator it_local = m_Mesh_Coupling.local_elements_begin();
696  libMesh::Mesh::const_element_iterator it_local_end = m_Mesh_Coupling.local_elements_end();
697  for( ; it_local != it_local_end; ++ it_local)
698  {
699  const libMesh::Elem * dummy_elem = * it_local;
700  nb_of_inters_per_rank[m_rank] += m_coupling_weights[dummy_elem->id()];
701  }
702  m_comm.sum(nb_of_inters_per_rank);
703  m_Intersection_Pairs_multimap.reserve(2*dummy_nb_of_inters);
704 
705  // ... and the intersection grid
706  m_Mesh_Intersection.preallocate_grid(192*dummy_nb_of_inters);
707  m_bHavePreallocData = true;
708 
710  {
711  if(m_rank == 0)
712  {
713  std::string intersection_statistics_filename;
715  {
716  intersection_statistics_filename = m_timing_data_file_base + "_inters_per_partition__new_partitioning.dat";
717  }
718  else
719  {
720  intersection_statistics_filename = m_timing_data_file_base + "_inters_per_partition__original.dat";
721  }
722 
723  libMesh::StatisticsVector<int> intersection_statistics(m_nodes,0);
724  for(unsigned int iii = 0; iii < m_nodes; ++iii)
725  {
726  intersection_statistics[iii] = nb_of_inters_per_rank[iii];
727  }
728 
729  std::ofstream inter_statistics_output(intersection_statistics_filename,std::ofstream::app);
730 
731  inter_statistics_output << intersection_statistics.minimum() << " "
732  << intersection_statistics.maximum() << " "
733  << intersection_statistics.mean() << " "
734  << intersection_statistics.median() << " "
735  << intersection_statistics.stddev() << std::endl;
736 
737  inter_statistics_output.close();
738  }
739 
740  }
741  }
742  }
bool m_bDidPreliminarySearch
Boolean flag determining if we did a preliminary intersection search, used for the coupling mesh repa...
std::string m_timing_data_file_base
Output filenames base for the timing data.
bool m_bPrintIntersectionsPerPartData
Print intersections per partition. Default: false.
libMesh::Mesh & m_Mesh_Coupling
Mesh representing the coupling region.
const libMesh::Parallel::Communicator & m_comm
Global communicator.
const unsigned int m_nodes
Number of MPI processors.
bool m_bSkipIntersectionPartitioning
Boolean indicating if we should skip the intersection repartitioning. Default: false.
std::unordered_multimap< unsigned int, unsigned int > m_Intersection_Pairs_multimap
Multimap containing the intersection pairs.
std::vector< unsigned int > m_Nb_Of_Intersections_Elem_C
Vector containing the number of intersections found inside each of the coupling mesh elements...
bool m_bHavePreallocData
Boolean flag determining if we have the data for a proper preallocation of m_Intersection_Pairs_multi...
Mesh_Intersection m_Mesh_Intersection
Object containing the intersection mesh.
libMesh::ErrorVector m_coupling_weights
libMesh::ErrorVector used to repartition the coupling mesh.
void preallocate_grid(int map_preallocation)
Preallocate the discrete points grid.
const unsigned int m_rank
This processor's rank (or ID in the communicator)
void carl::Intersection_Search::PreparePreallocationAndLoad ( SearchMethod  search_type = BRUTE)

Do a preliminary search, to optimize the intersection search and construction.

This method finds all the intersecting element pairs, only saving their number for each element of the coupling mesh. This information is used to repartition the coupling mesh and to preallocate the Intersection_Search::m_Intersection_Pairs_multimap.

Definition at line 652 of file intersection_search.cpp.

653  {
654  m_bSaveInterData = false;
655  switch (search_type)
656  {
658  break;
659 
661  break;
662 
665  break;
666  }
667 
669 
671  }
bool m_bDidPreliminarySearch
Boolean flag determining if we did a preliminary intersection search, used for the coupling mesh repa...
void FindIntersections_Front()
Find all the intersections, using the advancing front method.
void FindIntersections_Brute()
Find all the intersections, using the brute force method.
const libMesh::Parallel::Communicator & m_comm
Global communicator.
std::vector< unsigned int > m_Nb_Of_Intersections_Elem_C
Vector containing the number of intersections found inside each of the coupling mesh elements...
bool m_bSaveInterData
Boolean flag determining if we should save the intersection data or not. Default: true...
void carl::Intersection_Search::SetScalingFiles ( const std::string &  timing_data_file_base)

Set up timing output file.

Definition at line 938 of file intersection_search.cpp.

939  {
940  m_bPrintTimingData = true;
942  m_timing_data_file_base = timing_data_file_base;
943  }
std::string m_timing_data_file_base
Output filenames base for the timing data.
bool m_bPrintIntersectionsPerPartData
Print intersections per partition. Default: false.
bool m_bPrintTimingData
Print timing data. Default: false.
void carl::Intersection_Search::SkipIntersectionConstruction ( bool  bSkipIntersectionConstruction)
inline

Set the "Skip intersection construction" flag.

Definition at line 245 of file intersection_search.h.

246  {
247  m_bSkipIntersectionConstruction = bSkipIntersectionConstruction;
248  }
bool m_bSkipIntersectionConstruction
Boolean indicating if we should skip the intersection construction. Default: false.
void carl::Intersection_Search::SkipIntersectionPartitioning ( bool  bSkipIntersectionPartitioning)
inline

Set the "Skip the coupling mesh repartition" flag.

Definition at line 251 of file intersection_search.h.

252  {
253  m_bSkipIntersectionPartitioning = bSkipIntersectionPartitioning;
254  }
bool m_bSkipIntersectionPartitioning
Boolean indicating if we should skip the intersection repartitioning. Default: false.
void carl::Intersection_Search::UpdateCouplingIntersection ( const libMesh::Elem *  Query_elem)
protected

Build the intersections associated to the Query_elem and update the intersection mesh.

Definition at line 878 of file intersection_search.cpp.

879  {
880  // Preamble
881  std::unordered_set<unsigned int> & Patch_Set_A = m_Patch_Constructor_A.elem_indexes();
882  std::unordered_set<unsigned int>::iterator it_patch_A;
883 
884  bool bDoIntersect = true;
885  bool bCreateNewNefForA = true;
887 
888  std::set<libMesh::Point> points_out;
889 
890  // Debug vars
891  int nbOfTests = 0;
892  int nbOfPositiveTests = 0;
893 
894  // For each element inside patch A ...
895  for( it_patch_A = Patch_Set_A.begin();
896  it_patch_A != Patch_Set_A.end();
897  ++it_patch_A)
898  {
899  auto iterator_pair = m_Intersection_Pairs_multimap.equal_range(*it_patch_A);
900 
901  const libMesh::Elem * elem_A = m_Mesh_A.elem(*it_patch_A);
902  bCreateNewNefForA = true;
903 
904  // ... get the intersections with patch B ...
905  for(auto it_patch_B = iterator_pair.first ;
906  it_patch_B != iterator_pair.second ;
907  ++it_patch_B )
908  {
909  const libMesh::Elem * elem_B = m_Mesh_B.elem(it_patch_B->second);
910  points_out.clear();
911 
912  // ... test if they intersect inside the coupling region ...
913  m_perf_log.push("Build intersection polyhedron","Build intersections");
914  bDoIntersect = m_Intersection_test.libMesh_exact_intersection_inside_coupling(elem_A,elem_B,points_out,bCreateNewNefForA,true,false);
915  m_perf_log.pop("Build intersection polyhedron","Build intersections");
916  ++nbOfTests;
917 
918  if(bDoIntersect)
919  {
920  // ... and, if they do, update the intersection mesh!
921  bCreateNewNefForA = false;
922  m_perf_log.push("Update intersection","Build intersections");
923  m_Mesh_Intersection.increase_intersection_mesh(points_out,*it_patch_A,it_patch_B->second,Query_elem->id());
924  m_perf_log.pop("Update intersection","Build intersections");
925  ++nbOfPositiveTests;
926  }
927  }
928  }
929 
930  if(m_bPrintDebug)
931  {
932  std::cout << " DEBUG: update intersection mesh" << std::endl;
933  std::cout << " -> Positives / tests : " << nbOfPositiveTests << " / " << nbOfTests
934  << " (" << 100.*nbOfPositiveTests/nbOfTests << "%)" << std::endl << std::endl;
935  }
936  }
void increase_intersection_mesh(const std::set< libMesh::Point > &input_points, unsigned int elem_idx_A, unsigned int elem_idx_B, unsigned int elem_idx_C)
Increase the intersection mesh with the current polyhedron mesh and update its data structures...
void libmesh_set_coupling_nef_polyhedron(const libMesh::Elem *elem_C)
Set the Nef polyhedron for the coupling mesh element.
std::unordered_set< unsigned int > & elem_indexes()
Returns the set of elements inside the patch.
libMesh::PerfLog m_perf_log
libMesh::PerfLog object.
Patch_construction m_Patch_Constructor_A
Patch_construction object for mesh A.
std::unordered_multimap< unsigned int, unsigned int > m_Intersection_Pairs_multimap
Multimap containing the intersection pairs.
Mesh_Intersection m_Mesh_Intersection
Object containing the intersection mesh.
Intersection_Tools m_Intersection_test
Intersection operations for the main intersection tests.
bool m_bPrintDebug
Print debug data. Default: false.
bool libMesh_exact_intersection_inside_coupling(const libMesh::Elem *elem_A, const libMesh::Elem *elem_B, std::set< libMesh::Point > &points_out, bool bCreateNewNefForA=true, bool bConvertPoints=true, bool bTestNeeded=false)
Build the intersection between two elements intersect inside the coupling region. ...

Member Data Documentation

bool carl::Intersection_Search::m_bDidPreliminarySearch
protected

Boolean flag determining if we did a preliminary intersection search, used for the coupling mesh repartition and memory preallocations.

Definition at line 79 of file intersection_search.h.

bool carl::Intersection_Search::m_bHavePreallocData
protected

Boolean flag determining if we have the data for a proper preallocation of m_Intersection_Pairs_multimap.

Definition at line 82 of file intersection_search.h.

bool carl::Intersection_Search::m_bIntersectionsBuilt
protected

Boolean flag indicating if the intersections were built. (A bit useless right now ...)

Definition at line 85 of file intersection_search.h.

bool carl::Intersection_Search::m_bPrintDebug
protected

Print debug data. Default: false.

Definition at line 100 of file intersection_search.h.

bool carl::Intersection_Search::m_bPrintIntersectionsPerPartData
protected

Print intersections per partition. Default: false.

Definition at line 102 of file intersection_search.h.

bool carl::Intersection_Search::m_bPrintTimingData
protected

Print timing data. Default: false.

Definition at line 101 of file intersection_search.h.

bool carl::Intersection_Search::m_bSaveInterData
protected

Boolean flag determining if we should save the intersection data or not. Default: true.

Definition at line 76 of file intersection_search.h.

bool carl::Intersection_Search::m_bSkipIntersectionConstruction
protected

Boolean indicating if we should skip the intersection construction. Default: false.

Definition at line 91 of file intersection_search.h.

bool carl::Intersection_Search::m_bSkipIntersectionPartitioning
protected

Boolean indicating if we should skip the intersection repartitioning. Default: false.

Definition at line 94 of file intersection_search.h.

const libMesh::Parallel::Communicator& carl::Intersection_Search::m_comm
protected

Global communicator.

Definition at line 38 of file intersection_search.h.

libMesh::ErrorVector carl::Intersection_Search::m_coupling_weights
protected

libMesh::ErrorVector used to repartition the coupling mesh.

Definition at line 73 of file intersection_search.h.

std::unordered_multimap<unsigned int,unsigned int> carl::Intersection_Search::m_Intersection_Pairs_multimap
protected

Multimap containing the intersection pairs.

The element indexes of the mesh A are used as the keys of the multimap

Definition at line 61 of file intersection_search.h.

Intersection_Tools carl::Intersection_Search::m_Intersection_test
protected

Intersection operations for the main intersection tests.

Definition at line 64 of file intersection_search.h.

Intersection_Tools carl::Intersection_Search::m_Intersection_test_neighbors
protected

Intersection operations for the neighbor intersection tests. (why do we need both?)

Definition at line 67 of file intersection_search.h.

const libMesh::Parallel::Communicator& carl::Intersection_Search::m_local_comm
protected

Local communicator, used for mesh patches.

Definition at line 41 of file intersection_search.h.

libMesh::Mesh& carl::Intersection_Search::m_Mesh_A
protected

Definition at line 33 of file intersection_search.h.

libMesh::Mesh& carl::Intersection_Search::m_Mesh_B
protected

Definition at line 34 of file intersection_search.h.

libMesh::Mesh& carl::Intersection_Search::m_Mesh_Coupling
protected

Mesh representing the coupling region.

Definition at line 35 of file intersection_search.h.

Mesh_Intersection carl::Intersection_Search::m_Mesh_Intersection
protected

Object containing the intersection mesh.

Definition at line 55 of file intersection_search.h.

IntersectionMeshingMethod carl::Intersection_Search::m_MeshingMethod
protected

Flag defining the intersection meshing algorithm.

Can be either carl::LIBMESH_TETGEN (use LibMesh's tetgen module, problematic with Intel compilers) or carl::CGAL (use CGAL's Triangulation_3). Default: carl::CGAL.

Definition at line 52 of file intersection_search.h.

double carl::Intersection_Search::m_Min_Inter_Volume
protected

Volume cutoff for the intersection polyhedrons. Default: 1E-15.

Definition at line 88 of file intersection_search.h.

std::vector<unsigned int> carl::Intersection_Search::m_Nb_Of_Intersections_Elem_C
protected

Vector containing the number of intersections found inside each of the coupling mesh elements.

Definition at line 70 of file intersection_search.h.

const unsigned int carl::Intersection_Search::m_nodes
protected

Number of MPI processors.

Definition at line 39 of file intersection_search.h.

std::string carl::Intersection_Search::m_Output_filename_base
protected

Output filenames base.

Definition at line 96 of file intersection_search.h.

Patch_construction carl::Intersection_Search::m_Patch_Constructor_A
protected

Patch_construction object for mesh A.

Definition at line 44 of file intersection_search.h.

Patch_construction carl::Intersection_Search::m_Patch_Constructor_B
protected

Patch_construction object for mesh B.

Definition at line 45 of file intersection_search.h.

libMesh::PerfLog carl::Intersection_Search::m_perf_log
protected

libMesh::PerfLog object.

Definition at line 99 of file intersection_search.h.

const unsigned int carl::Intersection_Search::m_rank
protected

This processor's rank (or ID in the communicator)

Definition at line 40 of file intersection_search.h.

std::string carl::Intersection_Search::m_timing_data_file_base
protected

Output filenames base for the timing data.

Definition at line 103 of file intersection_search.h.

bool carl::Intersection_Search::MASTER_bPerfLog_intersection_search
protected

Do performance log? Default: true.

Definition at line 98 of file intersection_search.h.


The documentation for this class was generated from the following files: