![]() |
MFEM v4.7.0
Finite element discretization library
|
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points. There are three key functions in FindPointsGSLIB: More...
#include <gslib.hpp>
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. | |
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. | |
virtual const Array< unsigned int > & | GetProc () const |
Return MPI rank on which each point was found by FindPoints. | |
virtual const Vector & | GetReferencePosition () const |
Return reference coordinates for each point found by FindPoints. | |
virtual const Vector & | GetDist () const |
virtual const Array< unsigned int > & | GetGSLIBElem () const |
virtual const Vector & | GetGSLIBReferencePosition () const |
Protected Member Functions | |
virtual void | InterpolateH1 (const GridFunction &field_in, Vector &field_out) |
Use GSLIB for communication and interpolation. | |
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. | |
virtual void | MapRefPosAndElemIndices () |
Protected Attributes | |
Mesh * | mesh |
Array< Mesh * > | mesh_split |
Array< IntegrationRule * > | ir_split |
Array< FiniteElementSpace * > | fes_rst_map |
Array< GridFunction * > | gf_rst_map |
FiniteElementCollection * | fec_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 |
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points. There are three key functions in FindPointsGSLIB:
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.FindPointsGSLIB provides interface to use these functions individually or using a single call.
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.
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 ) |
|
virtual |
|
inlinevirtual |
|
inlinevirtual |
|
inlinevirtual |
|
inlinevirtual |
|
inlinevirtual |
|
protectedvirtual |
Get GridFunction value at the points expected by GSLIB.
|
inlinevirtual |
|
inlinevirtual |
|
virtual |
Interpolation of field values at prescribed reference space positions.
[in] | field_in | Function 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_out | Interpolated values. For points that are not found the value is set to default_interp_value. |
Reimplemented in mfem::OversetFindPointsGSLIB.
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.
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.
|
protectedvirtual |
|
protectedvirtual |
|
protectedvirtual |
|
inlinevirtual |
|
inlinevirtual |
|
inlinevirtual |
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.
[in] | m | Input 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. |
|
protectedvirtual |
|
protectedvirtual |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |