|
| 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) |
|
void | FindPoints (Mesh &m, const Vector &point_pos, 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) |
|
void | Interpolate (Mesh &m, const Vector &point_pos, const GridFunction &field_in, Vector &field_out) |
|
virtual void | SetL2AvgType (AvgType avgtype_) |
|
virtual void | SetDefaultInterpolationValue (double interp_value_) |
|
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 Vector & | GetReferencePosition () const |
| Return reference coordinates for each point found by FindPoints. More...
|
|
virtual const Vector & | GetDist () const |
|
virtual const Array< unsigned
int > & | GetGSLIBElem () const |
|
virtual const Vector & | GetGSLIBReferencePosition () const |
|
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points. There are three key functions in FindPointsGSLIB:
- Setup - constructs the internal data structures of gslib.
- 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.
- 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 47 of file gslib.hpp.
void mfem::FindPointsGSLIB::FindPoints |
( |
const Vector & |
point_pos | ) |
|
Searches positions given in physical space by point_pos. These positions must by ordered by nodes: (XXX...,YYY...,ZZZ). 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 131 of file gslib.cpp.