MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::FindPointsGSLIB Class Reference

FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points. 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.
 
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 SetGPUtoCPUFallback (bool mode)
 
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 VectorGetReferencePosition () const
 Return reference coordinates for each point found by FindPoints.
 
virtual const VectorGetDist () const
 
virtual const Array< unsigned int > & GetGSLIBElem () const
 
virtual const VectorGetGSLIBReferencePosition () const
 
void GetAxisAlignedBoundingBoxes (Vector &aabb)
 
void GetOrientedBoundingBoxes (DenseTensor &obbA, Vector &obbC, Vector &obbV)
 
Methods to support a custom interpolation procedure.

The physical-space point that the user seeks to interpolate at could be located inside an element on another mpi rank. To enable a custom interpolation procedure (e.g., strain tensor computation) we need a mechanism to first send element indices and reference-space coordinates to the mpi-ranks where each point is found. Then the custom interpolation can be done locally by the user before sending the interpolated values back to the mpi-ranks that the query originated from. Example usage looks something like this:

FindPoints() -> DistributePointInfoToOwningMPIRanks() -> Computation by user -> DistributeInterpolatedValues().

virtual void DistributePointInfoToOwningMPIRanks (Array< unsigned int > &recv_elem, Vector &recv_ref, Array< unsigned int > &recv_code)
 
virtual void DistributeInterpolatedValues (const Vector &int_vals, const int vdim, const int ordering, Vector &field_out) 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 ()
 
void FindPointsLocal3 (const Vector &point_pos, int point_pos_ordering, Array< unsigned int > &gsl_code_dev_l, Array< unsigned int > &gsl_elem_dev_l, Vector &gsl_ref_l, Vector &gsl_dist_l, int npt)
 
void FindPointsLocal2 (const Vector &point_pos, int point_pos_ordering, Array< unsigned int > &gsl_code_dev_l, Array< unsigned int > &gsl_elem_dev_l, Vector &gsl_ref_l, Vector &gsl_dist_l, int npt)
 
void InterpolateLocal3 (const Vector &field_in, Array< int > &gsl_elem_dev_l, Vector &gsl_ref_l, Vector &field_out, int npt, int ncomp, int nel, int dof1dsol)
 
void InterpolateLocal2 (const Vector &field_in, Array< int > &gsl_elem_dev_l, Vector &gsl_ref_l, Vector &field_out, int npt, int ncomp, int nel, int dof1dsol)
 
void SetupDevice ()
 
void FindPointsOnDevice (const Vector &point_pos, int point_pos_ordering=Ordering::byNODES)
 
void InterpolateOnDevice (const Vector &field_in_evec, Vector &field_out, const int nel, const int ncomp, const int dof1dsol, const int ordering)
 

Protected Attributes

Meshmesh
 
Array< Mesh * > mesh_split
 
Array< IntegrationRule * > ir_split
 
Array< FiniteElementSpace * > fes_rst_map
 
Array< GridFunction * > gf_rst_map
 
FiniteElementCollectionfec_map_lin
 
void * fdataD
 
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
 
Array< unsigned int > recv_proc
 
Array< unsigned int > recv_index
 
bool setupflag
 
double default_interp_value
 
AvgType avgtype
 
Array< int > split_element_map
 
Array< int > split_element_index
 
int NE_split_total
 
int mesh_points_cnt
 
double bdr_tol
 
bool gpu_to_cpu_fallback = false
 
struct { 
 
   bool   setup_device = false 
 
   bool   find_device = false 
 
   int   local_hash_size 
 
   int   dof1d 
 
   int   dof1d_sol 
 
   int   h_o_size 
 
   int   h_nx 
 
   double   newt_tol 
 
   struct gslib::crystal *   cr 
 
   struct gslib::hash_data_3 *   hash3 
 
   struct gslib::hash_data_2 *   hash2 
 
   Vector   bb 
 
   Vector   wtend 
 
   Vector   gll1d 
 
   Vector   lagcoeff 
 
   Vector   gll1d_sol 
 
   Vector   lagcoeff_sol 
 
   Array< unsigned int >   loc_hash_offset 
 
   Vector   loc_hash_min 
 
   Vector   loc_hash_fac 
 
DEV 
 

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. See Setup.
  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. See FindPoints.
  3. Interpolate - Interpolates any grid function at the points found using 2. For functions in L2 finite element space, use SetL2AvgType to specify how to interpolate values at points located at element boundaries where the function might be multi-valued. See Interpolate.

FindPointsGSLIB also provides interface to use these functions through a single call.

For custom interpolation (e.g., evaluating strain rate tensor), we provide functions that use gslib to send element index and corresponding reference-space coordinates for each point to the mpi rank that the element is located on. Then, custom interpolation can be defined locally by the user before sending the values back to mpi ranks where the query originated from. See DistributePointInfoToOwningMPIRanks and DistributeInterpolatedValues.

Definition at line 66 of file gslib.hpp.

Member Enumeration Documentation

◆ AvgType

Enumerator
NONE 
ARITHMETIC 
HARMONIC 

Definition at line 69 of file gslib.hpp.

Constructor & Destructor Documentation

◆ FindPointsGSLIB() [1/2]

mfem::FindPointsGSLIB::FindPointsGSLIB ( )

Definition at line 87 of file gslib.cpp.

◆ FindPointsGSLIB() [2/2]

mfem::FindPointsGSLIB::FindPointsGSLIB ( MPI_Comm comm_)

Definition at line 137 of file gslib.cpp.

◆ ~FindPointsGSLIB()

mfem::FindPointsGSLIB::~FindPointsGSLIB ( )
virtual

Definition at line 120 of file gslib.cpp.

Member Function Documentation

◆ DistributeInterpolatedValues()

void mfem::FindPointsGSLIB::DistributeInterpolatedValues ( const Vector & int_vals,
const int vdim,
const int ordering,
Vector & field_out ) const
virtual

Return interpolated values back to the mpi-ranks recv_proc that had sent the element indices and corresponding reference-space coordinates. Specify vdim and ordering (by nodes or by vdim) based on how the int_vals are structured. The received values are filled in field_out consistent with the original ordering of the points that were used in FindPoints.

Definition at line 2143 of file gslib.cpp.

◆ DistributePointInfoToOwningMPIRanks()

void mfem::FindPointsGSLIB::DistributePointInfoToOwningMPIRanks ( Array< unsigned int > & recv_elem,
Vector & recv_ref,
Array< unsigned int > & recv_code )
virtual

Distribute element indices in gsl_mfem_elem, the reference coordinates gsl_mfem_ref, and the code gsl_code to the corresponding mpi-rank gsl_proc for each point. The received information is provided locally in recv_elem, recv_ref (ordered by vdim), and recv_code. Note: The user can send empty Array/Vectors to the method as they are appropriately sized and filled internally.

Definition at line 2084 of file gslib.cpp.

◆ 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 242 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 1101 of file gslib.cpp.

◆ FindPointsLocal2()

void mfem::FindPointsGSLIB::FindPointsLocal2 ( const Vector & point_pos,
int point_pos_ordering,
Array< unsigned int > & gsl_code_dev_l,
Array< unsigned int > & gsl_elem_dev_l,
Vector & gsl_ref_l,
Vector & gsl_dist_l,
int npt )
protected

Definition at line 1151 of file findpts_local_2.cpp.

◆ FindPointsLocal3()

void mfem::FindPointsGSLIB::FindPointsLocal3 ( const Vector & point_pos,
int point_pos_ordering,
Array< unsigned int > & gsl_code_dev_l,
Array< unsigned int > & gsl_elem_dev_l,
Vector & gsl_ref_l,
Vector & gsl_dist_l,
int npt )
protected

Definition at line 1785 of file findpts_local_3.cpp.

◆ FindPointsOnDevice()

void mfem::FindPointsGSLIB::FindPointsOnDevice ( const Vector & point_pos,
int point_pos_ordering = Ordering::byNODES )
protected

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.

Definition at line 489 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 1128 of file gslib.cpp.

◆ GetAxisAlignedBoundingBoxes()

void mfem::FindPointsGSLIB::GetAxisAlignedBoundingBoxes ( Vector & aabb)

Return the axis-aligned bounding boxes (AABB) computed during Setup. The size of the returned vector is (nel x nverts x dim), where nel is the number of elements (after splitting for simplcies), nverts is number of vertices (4 in 2D, 8 in 3D), and dim is the spatial dimension.

Definition at line 2198 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 295 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 304 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 297 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 308 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 311 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 1434 of file gslib.cpp.

◆ GetOrientedBoundingBoxes()

void mfem::FindPointsGSLIB::GetOrientedBoundingBoxes ( DenseTensor & obbA,
Vector & obbC,
Vector & obbV )

Return the oriented bounding boxes (OBB) computed during Setup. Each OBB is represented using the inverse transformation (A^{-1}) and its center (x_c), such that a point x is inside the OBB if: -1 <= A^{-1}(x-x_c) <= 1. The inverse transformation is returned in obbA, a DenseTensor of size (dim x dim x nel), and the OBB centers are returned in obbC, a vector of size (nel x dim). The vertices of the OBBs are returned in obbV, a vector of size (nel x nverts x dim) .

Definition at line 2268 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 299 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 301 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.

Reimplemented in mfem::OversetFindPointsGSLIB.

Definition at line 1721 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 1112 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 1120 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 1923 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 1847 of file gslib.cpp.

◆ InterpolateLocal2()

void mfem::FindPointsGSLIB::InterpolateLocal2 ( const Vector & field_in,
Array< int > & gsl_elem_dev_l,
Vector & gsl_ref_l,
Vector & field_out,
int npt,
int ncomp,
int nel,
int dof1dsol )
protected

Definition at line 118 of file interpolate_local_2.cpp.

◆ InterpolateLocal3()

void mfem::FindPointsGSLIB::InterpolateLocal3 ( const Vector & field_in,
Array< int > & gsl_elem_dev_l,
Vector & gsl_ref_l,
Vector & field_out,
int npt,
int ncomp,
int nel,
int dof1dsol )
protected

Definition at line 123 of file interpolate_local_3.cpp.

◆ InterpolateOnDevice()

void mfem::FindPointsGSLIB::InterpolateOnDevice ( const Vector & field_in_evec,
Vector & field_out,
const int nel,
const int ncomp,
const int dof1dsol,
const int ordering )
protected

Interpolation of field values at prescribed reference space positions.

Parameters
[in]field_in_evecE-vector of grid function to be interpolated. Assumed ordering is NDOFSxVDIMxNEL
[in]nelNumber of elements in the mesh.
[in]ncompNumber of components in the field.
[in]dof1dsolNumber of degrees of freedom in each reference space direction.
[in]orderingOrdering of the out field values: byNodes/byVDIM
[out]field_outInterpolated values. For points that are not found the value is set to default_interp_value.

Definition at line 856 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 1544 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 271 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 279 of file gslib.hpp.

◆ SetGPUtoCPUFallback()

virtual void mfem::FindPointsGSLIB::SetGPUtoCPUFallback ( bool mode)
inlinevirtual

Enable/Disable use of CPU functions for GPU data if the gslib version is older.

Definition at line 286 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 267 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 163 of file gslib.cpp.

◆ SetupDevice()

void mfem::FindPointsGSLIB::SetupDevice ( )
protected

Definition at line 370 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 1386 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 1159 of file gslib.cpp.

Member Data Documentation

◆ avgtype

AvgType mfem::FindPointsGSLIB::avgtype
protected

Definition at line 89 of file gslib.hpp.

◆ bb

Vector mfem::FindPointsGSLIB::bb
mutable

Definition at line 110 of file gslib.hpp.

◆ bdr_tol

double mfem::FindPointsGSLIB::bdr_tol
protected

Definition at line 96 of file gslib.hpp.

◆ cr

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

Definition at line 81 of file gslib.hpp.

◆ default_interp_value

double mfem::FindPointsGSLIB::default_interp_value
protected

Definition at line 88 of file gslib.hpp.

◆ [struct]

struct { ... } mfem::FindPointsGSLIB::DEV

◆ dim

int mfem::FindPointsGSLIB::dim
protected

Definition at line 83 of file gslib.hpp.

◆ dof1d

int mfem::FindPointsGSLIB::dof1d

Definition at line 105 of file gslib.hpp.

◆ dof1d_sol

int mfem::FindPointsGSLIB::dof1d_sol

Definition at line 105 of file gslib.hpp.

◆ fdataD

void* mfem::FindPointsGSLIB::fdataD
protected

Definition at line 80 of file gslib.hpp.

◆ fec_map_lin

FiniteElementCollection* mfem::FindPointsGSLIB::fec_map_lin
protected

Definition at line 79 of file gslib.hpp.

◆ fes_rst_map

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

Definition at line 77 of file gslib.hpp.

◆ find_device

bool mfem::FindPointsGSLIB::find_device = false

Definition at line 104 of file gslib.hpp.

◆ gf_rst_map

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

Definition at line 78 of file gslib.hpp.

◆ gll1d

Vector mfem::FindPointsGSLIB::gll1d

Definition at line 110 of file gslib.hpp.

◆ gll1d_sol

Vector mfem::FindPointsGSLIB::gll1d_sol

Definition at line 110 of file gslib.hpp.

◆ gpu_to_cpu_fallback

bool mfem::FindPointsGSLIB::gpu_to_cpu_fallback = false
protected

Definition at line 98 of file gslib.hpp.

◆ gsl_code

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

Definition at line 84 of file gslib.hpp.

◆ gsl_comm

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

Definition at line 82 of file gslib.hpp.

◆ gsl_dist

Vector mfem::FindPointsGSLIB::gsl_dist
protected

Definition at line 85 of file gslib.hpp.

◆ gsl_elem

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

Definition at line 84 of file gslib.hpp.

◆ gsl_mesh

Vector mfem::FindPointsGSLIB::gsl_mesh
protected

Definition at line 85 of file gslib.hpp.

◆ gsl_mfem_elem

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

Definition at line 84 of file gslib.hpp.

◆ gsl_mfem_ref

Vector mfem::FindPointsGSLIB::gsl_mfem_ref
protected

Definition at line 85 of file gslib.hpp.

◆ gsl_proc

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

Definition at line 84 of file gslib.hpp.

◆ gsl_ref

Vector mfem::FindPointsGSLIB::gsl_ref
protected

Definition at line 85 of file gslib.hpp.

◆ h_nx

int mfem::FindPointsGSLIB::h_nx

Definition at line 105 of file gslib.hpp.

◆ h_o_size

int mfem::FindPointsGSLIB::h_o_size

Definition at line 105 of file gslib.hpp.

◆ hash2

struct gslib::hash_data_2* mfem::FindPointsGSLIB::hash2

Definition at line 109 of file gslib.hpp.

◆ hash3

struct gslib::hash_data_3* mfem::FindPointsGSLIB::hash3

Definition at line 108 of file gslib.hpp.

◆ ir_split

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

Definition at line 76 of file gslib.hpp.

◆ lagcoeff

Vector mfem::FindPointsGSLIB::lagcoeff

Definition at line 110 of file gslib.hpp.

◆ lagcoeff_sol

Vector mfem::FindPointsGSLIB::lagcoeff_sol

Definition at line 110 of file gslib.hpp.

◆ loc_hash_fac

Vector mfem::FindPointsGSLIB::loc_hash_fac

Definition at line 112 of file gslib.hpp.

◆ loc_hash_min

Vector mfem::FindPointsGSLIB::loc_hash_min
mutable

Definition at line 112 of file gslib.hpp.

◆ loc_hash_offset

Array<unsigned int> mfem::FindPointsGSLIB::loc_hash_offset
mutable

Definition at line 111 of file gslib.hpp.

◆ local_hash_size

int mfem::FindPointsGSLIB::local_hash_size

Definition at line 105 of file gslib.hpp.

◆ mesh

Mesh* mfem::FindPointsGSLIB::mesh
protected

Definition at line 72 of file gslib.hpp.

◆ mesh_points_cnt

int mfem::FindPointsGSLIB::mesh_points_cnt
protected

Definition at line 93 of file gslib.hpp.

◆ mesh_split

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

Definition at line 73 of file gslib.hpp.

◆ NE_split_total

int mfem::FindPointsGSLIB::NE_split_total
protected

Definition at line 92 of file gslib.hpp.

◆ newt_tol

double mfem::FindPointsGSLIB::newt_tol

Definition at line 106 of file gslib.hpp.

◆ points_cnt

int mfem::FindPointsGSLIB::points_cnt
protected

Definition at line 83 of file gslib.hpp.

◆ recv_index

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

Definition at line 86 of file gslib.hpp.

◆ recv_proc

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

Definition at line 86 of file gslib.hpp.

◆ setup_device

bool mfem::FindPointsGSLIB::setup_device = false

Definition at line 103 of file gslib.hpp.

◆ setupflag

bool mfem::FindPointsGSLIB::setupflag
protected

Definition at line 87 of file gslib.hpp.

◆ split_element_index

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

Definition at line 91 of file gslib.hpp.

◆ split_element_map

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

Definition at line 90 of file gslib.hpp.

◆ wtend

Vector mfem::FindPointsGSLIB::wtend

Definition at line 110 of file gslib.hpp.


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