SamplePointContainer Class Referenceabstract

Base class for all sample point containers. More...

#include <sample_point_container.h>

+ Inheritance diagram for SamplePointContainer:

Public Member Functions

 SamplePointContainer (Mesh *mesh_pt, const Vector< std::pair< double, double >> &min_and_max_coordinates, const bool &use_eulerian_coordinates_during_setup, const bool &ignore_halo_elements_during_locate_zeta_search, const unsigned &nsample_points_generated_per_element)
 Constructor. More...
 
 SamplePointContainer ()
 
 SamplePointContainer (const SamplePointContainer &data)=delete
 Broken copy constructor. More...
 
void operator= (const SamplePointContainer &)=delete
 Broken assignment operator. More...
 
virtual ~SamplePointContainer ()
 Virtual destructor. More...
 
virtual void locate_zeta (const Vector< double > &zeta, GeomObject *&sub_geom_object_pt, Vector< double > &s)=0
 
virtual unsignedtotal_number_of_sample_points_visited_during_locate_zeta_from_top_level ()
 
virtual unsigned total_number_of_sample_points_computed_recursively () const =0
 
virtual unsigned ndim_zeta () const =0
 Dimension of the zeta ( = dim of local coordinate of elements) More...
 
Mesh * mesh_pt () const
 Pointer to mesh from whose FiniteElements sample points are created. More...
 
const std::pair< double, double > & min_and_max_coordinates (const unsigned &i) const
 
const Vector< std::pair< double, double > > & min_and_max_coordinates () const
 
bool use_eulerian_coordinates_during_setup () const
 
unsignednsample_points_generated_per_element ()
 "Measure of" number of sample points generated in each element More...
 
doublemax_search_radius ()
 

Static Public Attributes

static std::ofstream Visited_sample_points_file
 File to record sequence of visited sample points in. More...
 
static bool Always_fail_elemental_locate_zeta = false
 Boolean flag to make to make locate zeta fail. More...
 
static bool Use_equally_spaced_interior_sample_points = true
 
static bool Enable_timing_of_setup = false
 Time setup? More...
 
static double Percentage_offset = 5.0
 Offset of sample point container boundaries beyond max/min coords. More...
 

Protected Member Functions

void setup_min_and_max_coordinates ()
 Setup the min and max coordinates for the mesh, in each dimension. More...
 

Protected Attributes

Mesh * Mesh_pt
 Pointer to mesh from whose FiniteElements sample points are created. More...
 
Vector< std::pair< double, double > > Min_and_max_coordinates
 
bool Use_eulerian_coordinates_during_setup
 
unsigned Nsample_points_generated_per_element
 "Measure of" number of sample points generated in each element More...
 
unsigned Total_number_of_sample_points_visited_during_locate_zeta_from_top_level
 
double Max_search_radius
 

Detailed Description

Base class for all sample point containers.

///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////////

Constructor & Destructor Documentation

◆ SamplePointContainer() [1/3]

SamplePointContainer::SamplePointContainer ( Mesh *  mesh_pt,
const Vector< std::pair< double, double >> &  min_and_max_coordinates,
const bool use_eulerian_coordinates_during_setup,
const bool ignore_halo_elements_during_locate_zeta_search,
const unsigned nsample_points_generated_per_element 
)
inline

Constructor.

215  : Mesh_pt(mesh_pt),
219 #ifdef OOMPH_HAS_MPI
220  Ignore_halo_elements_during_locate_zeta_search(
221  ignore_halo_elements_during_locate_zeta_search),
222 #endif
226  {
227  // Don't limit max. search radius
228  Max_search_radius = DBL_MAX;
229  }
Vector< std::pair< double, double > > Min_and_max_coordinates
Definition: sample_point_container.h:363
unsigned & nsample_points_generated_per_element()
"Measure of" number of sample points generated in each element
Definition: sample_point_container.h:314
const Vector< std::pair< double, double > > & min_and_max_coordinates() const
Definition: sample_point_container.h:290
bool use_eulerian_coordinates_during_setup() const
Definition: sample_point_container.h:308
unsigned Nsample_points_generated_per_element
"Measure of" number of sample points generated in each element
Definition: sample_point_container.h:377
unsigned Total_number_of_sample_points_visited_during_locate_zeta_from_top_level
Definition: sample_point_container.h:382
bool Use_eulerian_coordinates_during_setup
Definition: sample_point_container.h:367
Mesh * Mesh_pt
Pointer to mesh from whose FiniteElements sample points are created.
Definition: sample_point_container.h:357
double Max_search_radius
Definition: sample_point_container.h:389
Mesh * mesh_pt() const
Pointer to mesh from whose FiniteElements sample points are created.
Definition: sample_point_container.h:275

References Max_search_radius.

◆ SamplePointContainer() [2/3]

SamplePointContainer::SamplePointContainer ( )
inline

Broken default constructor; needed for broken copy constructors. Don't call. It will die.

234  {
235  // Throw the error
236  throw OomphLibError("Broken default constructor. Don't call this!",
239  }
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ SamplePointContainer() [3/3]

SamplePointContainer::SamplePointContainer ( const SamplePointContainer data)
delete

Broken copy constructor.

◆ ~SamplePointContainer()

virtual SamplePointContainer::~SamplePointContainer ( )
inlinevirtual

Virtual destructor.

248 {}

Member Function Documentation

◆ locate_zeta()

virtual void SamplePointContainer::locate_zeta ( const Vector< double > &  zeta,
GeomObject *&  sub_geom_object_pt,
Vector< double > &  s 
)
pure virtual

Find sub-GeomObject (finite element) and the local coordinate s within it that contains point with global coordinate zeta. sub_geom_object_pt=0 if point can't be found.

Implemented in NonRefineableBinArray, and RefineableBinArray.

Referenced by oomph::MeshAsGeomObject::locate_zeta().

◆ max_search_radius()

double& SamplePointContainer::max_search_radius ( )
inline

Set maximum search radius for locate zeta. This is initialised do DBL_MAX so we brutally search through the entire bin structure, no matter how big it is until we've found the required point (or failed to do so. This can be VERY costly with fine meshes. Here the user takes full responsibility and states that we have no chance in hell to find the required point in a bin whose closest vertex is further than the specified max search radius.

328  {
329  return Max_search_radius;
330  }

References Max_search_radius.

Referenced by RefineableBinArray::locate_zeta(), and oomph::LineVisualiser::setup().

◆ mesh_pt()

Mesh* SamplePointContainer::mesh_pt ( ) const
inline

Pointer to mesh from whose FiniteElements sample points are created.

276  {
277  return Mesh_pt;
278  }

References Mesh_pt.

◆ min_and_max_coordinates() [1/2]

const Vector<std::pair<double, double> >& SamplePointContainer::min_and_max_coordinates ( ) const
inline

Vector of pair of doubles for min and maximum coordinates. min (first) and max. (second) coordinates

291  {
293  }

References Min_and_max_coordinates.

Referenced by RefineableBinArray::output_neighbouring_bins().

◆ min_and_max_coordinates() [2/2]

const std::pair<double, double>& SamplePointContainer::min_and_max_coordinates ( const unsigned i) const
inline

Pair of doubles for min and maximum coordinates in i-th direction: min (first) and max. (second) coordinates

284  {
285  return Min_and_max_coordinates[i];
286  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9

References i, and Min_and_max_coordinates.

Referenced by oomph::RefineableGmshTetMesh< ELEMENT >::adapt().

◆ ndim_zeta()

virtual unsigned SamplePointContainer::ndim_zeta ( ) const
pure virtual

Dimension of the zeta ( = dim of local coordinate of elements)

Implemented in BinArray.

◆ nsample_points_generated_per_element()

unsigned& SamplePointContainer::nsample_points_generated_per_element ( )
inline

"Measure of" number of sample points generated in each element

315  {
317  }

References Nsample_points_generated_per_element.

◆ operator=()

void SamplePointContainer::operator= ( const SamplePointContainer )
delete

Broken assignment operator.

◆ setup_min_and_max_coordinates()

void SamplePointContainer::setup_min_and_max_coordinates ( )
protected

Setup the min and max coordinates for the mesh, in each dimension.

Helper function to compute the min and max coordinates for the mesh, in each dimension

685  {
686  // Get the lagrangian dimension
687  int n_lagrangian = ndim_zeta();
688 
689  // Storage locally (i.e. in parallel on each processor)
690  // for the minimum and maximum coordinates
691  double zeta_min_local[n_lagrangian];
692  double zeta_max_local[n_lagrangian];
693  for (int i = 0; i < n_lagrangian; i++)
694  {
695  zeta_min_local[i] = DBL_MAX;
696  zeta_max_local[i] = -DBL_MAX;
697  }
698 
699  // Loop over the elements of the mesh
700  unsigned n_el = Mesh_pt->nelement();
701  for (unsigned e = 0; e < n_el; e++)
702  {
703  FiniteElement* el_pt = Mesh_pt->finite_element_pt(e);
704 
705  // Get the number of vertices (nplot=2 does the trick)
706  unsigned n_plot = 2;
707  unsigned n_plot_points = el_pt->nplot_points(n_plot);
708 
709  // Loop over the number of plot points
710  for (unsigned iplot = 0; iplot < n_plot_points; iplot++)
711  {
712  Vector<double> s_local(n_lagrangian);
713  Vector<double> zeta_global(n_lagrangian);
714 
715  // Get the local s -- need to sample over the entire range
716  // of the elements!
717  bool use_equally_spaced_interior_sample_points = false;
718  el_pt->get_s_plot(
719  iplot, n_plot, s_local, use_equally_spaced_interior_sample_points);
720 
721  // Now interpolate to global coordinates
723  {
724  el_pt->interpolated_x(s_local, zeta_global);
725  }
726  else
727  {
728  el_pt->interpolated_zeta(s_local, zeta_global);
729  }
730 
731  // Check the max and min in each direction
732  for (int i = 0; i < n_lagrangian; i++)
733  {
734  // Is the coordinate less than the minimum?
735  if (zeta_global[i] < zeta_min_local[i])
736  {
737  zeta_min_local[i] = zeta_global[i];
738  }
739  // Is the coordinate bigger than the maximum?
740  if (zeta_global[i] > zeta_max_local[i])
741  {
742  zeta_max_local[i] = zeta_global[i];
743  }
744  }
745  }
746  }
747 
748  // Global extrema - in parallel, need to get max/min across all processors
749  double zeta_min[n_lagrangian];
750  double zeta_max[n_lagrangian];
751  for (int i = 0; i < n_lagrangian; i++)
752  {
753  zeta_min[i] = 0.0;
754  zeta_max[i] = 0.0;
755  }
756 
757 #ifdef OOMPH_HAS_MPI
758  // If the mesh has been distributed and we want consistent bins
759  // across all processors
760  if (Mesh_pt->is_mesh_distributed())
761  {
762  // .. we need a non-null communicator!
763  if (Mesh_pt->communicator_pt() != 0)
764  {
765  int n_proc = Mesh_pt->communicator_pt()->nproc();
766  if (n_proc > 1)
767  {
768  // Get the minima and maxima over all processors
769  MPI_Allreduce(zeta_min_local,
770  zeta_min,
771  n_lagrangian,
772  MPI_DOUBLE,
773  MPI_MIN,
774  Mesh_pt->communicator_pt()->mpi_comm());
775  MPI_Allreduce(zeta_max_local,
776  zeta_max,
777  n_lagrangian,
778  MPI_DOUBLE,
779  MPI_MAX,
780  Mesh_pt->communicator_pt()->mpi_comm());
781  }
782  }
783  else // Null communicator - throw an error
784  {
785  std::ostringstream error_message_stream;
786  error_message_stream << "Communicator not set for a Mesh\n"
787  << "that was created from a distributed Mesh";
788  throw OomphLibError(error_message_stream.str(),
791  }
792  }
793  else // If the mesh hasn't been distributed then the
794  // max and min are the same on all processors
795  {
796  for (int i = 0; i < n_lagrangian; i++)
797  {
798  zeta_min[i] = zeta_min_local[i];
799  zeta_max[i] = zeta_max_local[i];
800  }
801  }
802 #else // If we're not using MPI then the mesh can't be distributed
803  for (int i = 0; i < n_lagrangian; i++)
804  {
805  zeta_min[i] = zeta_min_local[i];
806  zeta_max[i] = zeta_max_local[i];
807  }
808 #endif
809 
810  // Decrease/increase min and max to allow for any overshoot in
811  // meshes that may move around
812  // There's no point in doing this for DIM_LAGRANGIAN==1
813  for (int i = 0; i < n_lagrangian; i++)
814  {
815  double length = zeta_max[i] - zeta_min[i];
816  zeta_min[i] -= ((Percentage_offset / 100.0) * length);
817  zeta_max[i] += ((Percentage_offset / 100.0) * length);
818  }
819 
820  // Set the entries as the Min/Max_coords vector
821  Min_and_max_coordinates.resize(n_lagrangian);
822  for (int i = 0; i < n_lagrangian; i++)
823  {
824  Min_and_max_coordinates[i].first = zeta_min[i];
825  Min_and_max_coordinates[i].second = zeta_max[i];
826  }
827  }
Array< double, 1, 3 > e(1./3., 0.5, 2.)
static double Percentage_offset
Offset of sample point container boundaries beyond max/min coords.
Definition: sample_point_container.h:349
virtual unsigned ndim_zeta() const =0
Dimension of the zeta ( = dim of local coordinate of elements)

References e(), i, OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by NonRefineableBinArray::NonRefineableBinArray(), and RefineableBinArray::RefineableBinArray().

◆ total_number_of_sample_points_computed_recursively()

virtual unsigned SamplePointContainer::total_number_of_sample_points_computed_recursively ( ) const
pure virtual

Total number of sample points in sample point container, possibly computed recursively.

Implemented in NonRefineableBinArray, and RefineableBinArray.

◆ total_number_of_sample_points_visited_during_locate_zeta_from_top_level()

virtual unsigned& SamplePointContainer::total_number_of_sample_points_visited_during_locate_zeta_from_top_level ( )
inlinevirtual

Counter to keep track of how many sample points we've visited during top level call to locate_zeta. Virtual so it can be overloaded for different versions.

Reimplemented in RefineableBinArray.

References Total_number_of_sample_points_visited_during_locate_zeta_from_top_level.

Referenced by NonRefineableBinArray::locate_zeta().

◆ use_eulerian_coordinates_during_setup()

bool SamplePointContainer::use_eulerian_coordinates_during_setup ( ) const
inline

Use Eulerian coordinates (i.e. interpolated_x) rather than zeta itself (i.e. interpolated_zeta) to identify point.

309  {
311  }

References Use_eulerian_coordinates_during_setup.

Referenced by NonRefineableBinArray::locate_zeta().

Member Data Documentation

◆ Always_fail_elemental_locate_zeta

bool SamplePointContainer::Always_fail_elemental_locate_zeta = false
static

Boolean flag to make to make locate zeta fail.

Boolean flag to make to make locate zeta fail. Used for debugging/ illustration of search procedures.

Referenced by check_locate_zeta(), RefineableBin::locate_zeta(), RefineableBinArray::locate_zeta(), and NonRefineableBinArray::locate_zeta().

◆ Enable_timing_of_setup

bool SamplePointContainer::Enable_timing_of_setup = false
static

◆ Max_search_radius

double SamplePointContainer::Max_search_radius
protected

Max radius beyond which we stop searching the bin. Initialised to DBL_MAX so keep going until the point is found or until we've searched every single bin. Overwriting this means we won't search in bins whose closest vertex is at a distance greater than Max_search_radius from the point to be located.

Referenced by NonRefineableBinArray::locate_zeta(), max_search_radius(), and SamplePointContainer().

◆ Mesh_pt

Mesh* SamplePointContainer::Mesh_pt
protected

◆ Min_and_max_coordinates

◆ Nsample_points_generated_per_element

unsigned SamplePointContainer::Nsample_points_generated_per_element
protected

"Measure of" number of sample points generated in each element

Referenced by RefineableBinArray::fill_bin_array(), NonRefineableBinArray::fill_bin_array(), and nsample_points_generated_per_element().

◆ Percentage_offset

double SamplePointContainer::Percentage_offset = 5.0
static

Offset of sample point container boundaries beyond max/min coords.

Referenced by NonRefineableBinArray::fill_bin_array().

◆ Total_number_of_sample_points_visited_during_locate_zeta_from_top_level

unsigned SamplePointContainer::Total_number_of_sample_points_visited_during_locate_zeta_from_top_level
protected

◆ Use_equally_spaced_interior_sample_points

bool SamplePointContainer::Use_equally_spaced_interior_sample_points = true
static

Use equally spaced sample points? (otherwise vertices are sampled repeatedly

Referenced by RefineableBinArray::fill_bin_array(), NonRefineableBinArray::fill_bin_array(), and RefineableBin::locate_zeta().

◆ Use_eulerian_coordinates_during_setup

bool SamplePointContainer::Use_eulerian_coordinates_during_setup
protected

Use Eulerian coordinates (i.e. interpolated_x) rather than zeta itself (i.e. interpolated_zeta) to identify point.

Referenced by RefineableBinArray::fill_bin_array(), NonRefineableBinArray::fill_bin_array(), NonRefineableBinArray::output_bins(), and use_eulerian_coordinates_during_setup().

◆ Visited_sample_points_file

std::ofstream SamplePointContainer::Visited_sample_points_file
static

File to record sequence of visited sample points in.

File to record sequence of visited sample points in. Used for debugging/ illustration of search procedures.

/////////////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////////////// SamplePointContainer base class /////////////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////////////

Referenced by check_locate_zeta(), RefineableBin::locate_zeta(), and NonRefineableBinArray::locate_zeta().


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