MFEM  v4.5.2
Finite element discretization library
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
mfem::FindPointsGSLIB Class Reference

FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points. There are three key functions in FindPointsGSLIB: More...

#include <gslib.hpp>

Inheritance diagram for mfem::FindPointsGSLIB:
[legend]
Collaboration diagram for mfem::FindPointsGSLIB:
[legend]

Public Types

enum  AvgType { NONE, ARITHMETIC, HARMONIC }
 

Public Member Functions

 FindPointsGSLIB ()
 
 FindPointsGSLIB (MPI_Comm comm_)
 
virtual ~FindPointsGSLIB ()
 
void Setup (Mesh &m, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
 
void FindPoints (const Vector &point_pos, int point_pos_ordering=Ordering::byNODES)
 
void FindPoints (Mesh &m, const Vector &point_pos, int point_pos_ordering=Ordering::byNODES, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
 Setup FindPoints and search positions. More...
 
virtual void Interpolate (const GridFunction &field_in, Vector &field_out)
 
void Interpolate (const Vector &point_pos, const GridFunction &field_in, Vector &field_out, int point_pos_ordering=Ordering::byNODES)
 
void Interpolate (Mesh &m, const Vector &point_pos, const GridFunction &field_in, Vector &field_out, int point_pos_ordering=Ordering::byNODES)
 
virtual void SetL2AvgType (AvgType avgtype_)
 
virtual void SetDefaultInterpolationValue (double interp_value_)
 
virtual void SetDistanceToleranceForPointsFoundOnBoundary (double bdr_tol_)
 
virtual void FreeData ()
 
virtual const Array< unsigned int > & GetCode () const
 
virtual const Array< unsigned int > & GetElem () const
 Return element number for each point found by FindPoints. More...
 
virtual const Array< unsigned int > & GetProc () const
 Return MPI rank on which each point was found by FindPoints. More...
 
virtual const VectorGetReferencePosition () const
 Return reference coordinates for each point found by FindPoints. More...
 
virtual const VectorGetDist () const
 
virtual const Array< unsigned int > & GetGSLIBElem () const
 
virtual const VectorGetGSLIBReferencePosition () const
 

Protected Member Functions

virtual void InterpolateH1 (const GridFunction &field_in, Vector &field_out)
 Use GSLIB for communication and interpolation. More...
 
virtual void InterpolateGeneral (const GridFunction &field_in, Vector &field_out)
 
virtual void SetupSplitMeshes ()
 
virtual void SetupIntegrationRuleForSplitMesh (Mesh *mesh, IntegrationRule *irule, int order)
 
virtual void GetNodalValues (const GridFunction *gf_in, Vector &node_vals)
 Get GridFunction value at the points expected by GSLIB. More...
 
virtual void MapRefPosAndElemIndices ()
 

Protected Attributes

Meshmesh
 
Array< Mesh * > mesh_split
 
Array< IntegrationRule * > ir_split
 
Array< FiniteElementSpace * > fes_rst_map
 
Array< GridFunction * > gf_rst_map
 
FiniteElementCollectionfec_map_lin
 
struct gslib::findpts_data_2 * fdata2D
 
struct gslib::findpts_data_3 * fdata3D
 
struct gslib::crystal * cr
 
struct gslib::comm * gsl_comm
 
int dim
 
int points_cnt
 
Array< unsigned int > gsl_code
 
Array< unsigned int > gsl_proc
 
Array< unsigned int > gsl_elem
 
Array< unsigned int > gsl_mfem_elem
 
Vector gsl_mesh
 
Vector gsl_ref
 
Vector gsl_dist
 
Vector gsl_mfem_ref
 
bool setupflag
 
double default_interp_value
 
AvgType avgtype
 
Array< int > split_element_map
 
Array< int > split_element_index
 
int NE_split_total
 
double bdr_tol
 

Detailed Description

FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points. There are three key functions in FindPointsGSLIB:

  1. Setup - constructs the internal data structures of gslib.
  2. FindPoints - for any given arbitrary set of points in physical space, gslib finds the element number, MPI rank, and the reference space coordinates inside the element that each point is located in. gslib also returns a code that indicates whether the point was found inside an element, on element border, or not found in the domain. For points returned as found on element border, the point is either on an element edge/face or near the domain boundary, and gslib also returns a distance to the border. Points near (but outside) the domain boundary must then be marked as not found using the distance returned by gslib.
  3. Interpolate - Interpolates any grid function at the points found using 2.

FindPointsGSLIB provides interface to use these functions individually or using a single call.

Definition at line 52 of file gslib.hpp.

Member Enumeration Documentation

◆ AvgType

Enumerator
NONE 
ARITHMETIC 
HARMONIC 

Definition at line 55 of file gslib.hpp.

Constructor & Destructor Documentation

◆ FindPointsGSLIB() [1/2]

mfem::FindPointsGSLIB::FindPointsGSLIB ( )

Definition at line 35 of file gslib.cpp.

◆ FindPointsGSLIB() [2/2]

mfem::FindPointsGSLIB::FindPointsGSLIB ( MPI_Comm  comm_)

Definition at line 82 of file gslib.cpp.

◆ ~FindPointsGSLIB()

mfem::FindPointsGSLIB::~FindPointsGSLIB ( )
virtual

Definition at line 67 of file gslib.cpp.

Member Function Documentation

◆ FindPoints() [1/2]

void mfem::FindPointsGSLIB::FindPoints ( const Vector point_pos,
int  point_pos_ordering = Ordering::byNODES 
)

Searches positions given in physical space by point_pos. These positions can be ordered byNodes: (XXX...,YYY...,ZZZ) or byVDim: (XYZ,XYZ,....XYZ) specified by point_pos_ordering. This function populates the following member variables: gsl_code Return codes for each point: inside element (0), element boundary (1), not found (2). gsl_proc MPI proc ids where the points were found. gsl_elem Element ids where the points were found. Defaults to 0 for points that were not found. gsl_mfem_elem Element ids corresponding to MFEM-mesh where the points were found. gsl_mfem_elem != gsl_elem for simplices Defaults to 0 for points that were not found. gsl_ref Reference coordinates of the found point. Ordered by vdim (XYZ,XYZ,XYZ...). Defaults to -1 for points that were not found. Note: the gslib reference frame is [-1,1]. gsl_mfem_ref Reference coordinates gsl_ref mapped to [0,1]. Defaults to 0 for points that were not found. gsl_dist Distance between the sought and the found point in physical space.

Definition at line 171 of file gslib.cpp.

◆ FindPoints() [2/2]

void mfem::FindPointsGSLIB::FindPoints ( Mesh m,
const Vector point_pos,
int  point_pos_ordering = Ordering::byNODES,
const double  bb_t = 0.1,
const double  newt_tol = 1.0e-12,
const int  npt_max = 256 
)

Setup FindPoints and search positions.

Definition at line 240 of file gslib.cpp.

◆ FreeData()

void mfem::FindPointsGSLIB::FreeData ( )
virtual

Cleans up memory allocated internally by gslib. Note that in parallel, this must be called before MPI_Finalize(), as it calls MPI_Comm_free() for internal gslib communicators.

Definition at line 267 of file gslib.cpp.

◆ GetCode()

virtual const Array<unsigned int>& mfem::FindPointsGSLIB::GetCode ( ) const
inlinevirtual

Return code for each point searched by FindPoints: inside element (0), on element boundary (1), or not found (2).

Definition at line 206 of file gslib.hpp.

◆ GetDist()

virtual const Vector& mfem::FindPointsGSLIB::GetDist ( ) const
inlinevirtual

Return distance between the sought and the found point in physical space, for each point found by FindPoints.

Definition at line 215 of file gslib.hpp.

◆ GetElem()

virtual const Array<unsigned int>& mfem::FindPointsGSLIB::GetElem ( ) const
inlinevirtual

Return element number for each point found by FindPoints.

Definition at line 208 of file gslib.hpp.

◆ GetGSLIBElem()

virtual const Array<unsigned int>& mfem::FindPointsGSLIB::GetGSLIBElem ( ) const
inlinevirtual

Return element number for each point found by FindPoints corresponding to GSLIB mesh. gsl_mfem_elem != gsl_elem for mesh with simplices.

Definition at line 219 of file gslib.hpp.

◆ GetGSLIBReferencePosition()

virtual const Vector& mfem::FindPointsGSLIB::GetGSLIBReferencePosition ( ) const
inlinevirtual

Return reference coordinates in [-1,1] (internal range in GSLIB) for each point found by FindPoints.

Definition at line 222 of file gslib.hpp.

◆ GetNodalValues()

void mfem::FindPointsGSLIB::GetNodalValues ( const GridFunction gf_in,
Vector node_vals 
)
protectedvirtual

Get GridFunction value at the points expected by GSLIB.

Definition at line 562 of file gslib.cpp.

◆ GetProc()

virtual const Array<unsigned int>& mfem::FindPointsGSLIB::GetProc ( ) const
inlinevirtual

Return MPI rank on which each point was found by FindPoints.

Definition at line 210 of file gslib.hpp.

◆ GetReferencePosition()

virtual const Vector& mfem::FindPointsGSLIB::GetReferencePosition ( ) const
inlinevirtual

Return reference coordinates for each point found by FindPoints.

Definition at line 212 of file gslib.hpp.

◆ Interpolate() [1/3]

void mfem::FindPointsGSLIB::Interpolate ( const GridFunction field_in,
Vector field_out 
)
virtual

Interpolation of field values at prescribed reference space positions.

Parameters
[in]field_inFunction values that will be interpolated on the reference positions. Note: it is assumed that field_in is in H1 and in the same space as the mesh that was given to Setup().
[out]field_outInterpolated values. For points that are not found the value is set to default_interp_value.

Definition at line 803 of file gslib.cpp.

◆ Interpolate() [2/3]

void mfem::FindPointsGSLIB::Interpolate ( const Vector point_pos,
const GridFunction field_in,
Vector field_out,
int  point_pos_ordering = Ordering::byNODES 
)

Search positions and interpolate. The ordering (byNODES or byVDIM) of the output values in field_out corresponds to the ordering used in the input GridFunction field_in.

Definition at line 251 of file gslib.cpp.

◆ Interpolate() [3/3]

void mfem::FindPointsGSLIB::Interpolate ( Mesh m,
const Vector point_pos,
const GridFunction field_in,
Vector field_out,
int  point_pos_ordering = Ordering::byNODES 
)

Setup FindPoints, search positions and interpolate. The ordering (byNODES or byVDIM) of the output values in field_out corresponds to the ordering used in the input GridFunction field_in.

Definition at line 259 of file gslib.cpp.

◆ InterpolateGeneral()

void mfem::FindPointsGSLIB::InterpolateGeneral ( const GridFunction field_in,
Vector field_out 
)
protectedvirtual

Uses GSLIB Crystal Router for communication followed by MFEM's interpolation functions

Definition at line 948 of file gslib.cpp.

◆ InterpolateH1()

void mfem::FindPointsGSLIB::InterpolateH1 ( const GridFunction field_in,
Vector field_out 
)
protectedvirtual

Use GSLIB for communication and interpolation.

Definition at line 885 of file gslib.cpp.

◆ MapRefPosAndElemIndices()

void mfem::FindPointsGSLIB::MapRefPosAndElemIndices ( )
protectedvirtual

Map {r,s,t} coordinates from [-1,1] to [0,1] for MFEM. For simplices, find the original element number (that was split into micro quads/hexes) during the setup phase.

Definition at line 656 of file gslib.cpp.

◆ SetDefaultInterpolationValue()

virtual void mfem::FindPointsGSLIB::SetDefaultInterpolationValue ( double  interp_value_)
inlinevirtual

Set the default interpolation value for points that are not found in the mesh.

Definition at line 186 of file gslib.hpp.

◆ SetDistanceToleranceForPointsFoundOnBoundary()

virtual void mfem::FindPointsGSLIB::SetDistanceToleranceForPointsFoundOnBoundary ( double  bdr_tol_)
inlinevirtual

Set the tolerance for detecting points outside the 'curvilinear' boundary that gslib may return as found on the boundary. Points found on boundary with distance greater than @ bdr_tol are marked as not found.

Definition at line 194 of file gslib.hpp.

◆ SetL2AvgType()

virtual void mfem::FindPointsGSLIB::SetL2AvgType ( AvgType  avgtype_)
inlinevirtual

Average type to be used for L2 functions in-case a point is located at an element boundary where the function might be multi-valued.

Definition at line 182 of file gslib.hpp.

◆ Setup()

void mfem::FindPointsGSLIB::Setup ( Mesh m,
const double  bb_t = 0.1,
const double  newt_tol = 1.0e-12,
const int  npt_max = 256 
)

Initializes the internal mesh in gslib, by sending the positions of the Gauss-Lobatto nodes of the input Mesh object m. Note: not tested with periodic (L2). Note: the input mesh m must have Nodes set.

Parameters
[in]mInput mesh.
[in]bb_t(Optional) Relative size of bounding box around each element.
[in]newt_tol(Optional) Newton tolerance for the gslib search methods.
[in]npt_max(Optional) Number of points for simultaneous iteration. This alters performance and memory footprint.

Definition at line 107 of file gslib.cpp.

◆ SetupIntegrationRuleForSplitMesh()

void mfem::FindPointsGSLIB::SetupIntegrationRuleForSplitMesh ( Mesh mesh,
IntegrationRule irule,
int  order 
)
protectedvirtual

Setup integration points that will be used to interpolate the nodal location at points expected by GSLIB.

Definition at line 514 of file gslib.cpp.

◆ SetupSplitMeshes()

void mfem::FindPointsGSLIB::SetupSplitMeshes ( )
protectedvirtual

Since GSLIB is designed to work with quads/hexes, we split every triangle/tet/prism/pyramid element into quads/hexes.

Definition at line 296 of file gslib.cpp.

Member Data Documentation

◆ avgtype

AvgType mfem::FindPointsGSLIB::avgtype
protected

Definition at line 74 of file gslib.hpp.

◆ bdr_tol

double mfem::FindPointsGSLIB::bdr_tol
protected

Definition at line 79 of file gslib.hpp.

◆ cr

struct gslib::crystal* mfem::FindPointsGSLIB::cr
protected

Definition at line 67 of file gslib.hpp.

◆ default_interp_value

double mfem::FindPointsGSLIB::default_interp_value
protected

Definition at line 73 of file gslib.hpp.

◆ dim

int mfem::FindPointsGSLIB::dim
protected

Definition at line 69 of file gslib.hpp.

◆ fdata2D

struct gslib::findpts_data_2* mfem::FindPointsGSLIB::fdata2D
protected

Definition at line 65 of file gslib.hpp.

◆ fdata3D

struct gslib::findpts_data_3* mfem::FindPointsGSLIB::fdata3D
protected

Definition at line 66 of file gslib.hpp.

◆ fec_map_lin

FiniteElementCollection* mfem::FindPointsGSLIB::fec_map_lin
protected

Definition at line 64 of file gslib.hpp.

◆ fes_rst_map

Array<FiniteElementSpace *> mfem::FindPointsGSLIB::fes_rst_map
protected

Definition at line 62 of file gslib.hpp.

◆ gf_rst_map

Array<GridFunction *> mfem::FindPointsGSLIB::gf_rst_map
protected

Definition at line 63 of file gslib.hpp.

◆ gsl_code

Array<unsigned int> mfem::FindPointsGSLIB::gsl_code
protected

Definition at line 70 of file gslib.hpp.

◆ gsl_comm

struct gslib::comm* mfem::FindPointsGSLIB::gsl_comm
protected

Definition at line 68 of file gslib.hpp.

◆ gsl_dist

Vector mfem::FindPointsGSLIB::gsl_dist
protected

Definition at line 71 of file gslib.hpp.

◆ gsl_elem

Array<unsigned int> mfem::FindPointsGSLIB::gsl_elem
protected

Definition at line 70 of file gslib.hpp.

◆ gsl_mesh

Vector mfem::FindPointsGSLIB::gsl_mesh
protected

Definition at line 71 of file gslib.hpp.

◆ gsl_mfem_elem

Array<unsigned int> mfem::FindPointsGSLIB::gsl_mfem_elem
protected

Definition at line 70 of file gslib.hpp.

◆ gsl_mfem_ref

Vector mfem::FindPointsGSLIB::gsl_mfem_ref
protected

Definition at line 71 of file gslib.hpp.

◆ gsl_proc

Array<unsigned int> mfem::FindPointsGSLIB::gsl_proc
protected

Definition at line 70 of file gslib.hpp.

◆ gsl_ref

Vector mfem::FindPointsGSLIB::gsl_ref
protected

Definition at line 71 of file gslib.hpp.

◆ ir_split

Array<IntegrationRule *> mfem::FindPointsGSLIB::ir_split
protected

Definition at line 60 of file gslib.hpp.

◆ mesh

Mesh* mfem::FindPointsGSLIB::mesh
protected

Definition at line 58 of file gslib.hpp.

◆ mesh_split

Array<Mesh *> mfem::FindPointsGSLIB::mesh_split
protected

Definition at line 59 of file gslib.hpp.

◆ NE_split_total

int mfem::FindPointsGSLIB::NE_split_total
protected

Definition at line 77 of file gslib.hpp.

◆ points_cnt

int mfem::FindPointsGSLIB::points_cnt
protected

Definition at line 69 of file gslib.hpp.

◆ setupflag

bool mfem::FindPointsGSLIB::setupflag
protected

Definition at line 72 of file gslib.hpp.

◆ split_element_index

Array<int> mfem::FindPointsGSLIB::split_element_index
protected

Definition at line 76 of file gslib.hpp.

◆ split_element_map

Array<int> mfem::FindPointsGSLIB::split_element_map
protected

Definition at line 75 of file gslib.hpp.


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