MFEM
v4.5.2
Finite element discretization library
|
OversetFindPointsGSLIB enables use of findpts for arbitrary number of overlapping grids. The parameters in this class are the same as FindPointsGSLIB with the difference of additional inputs required to account for more than 1 mesh. More...
#include <gslib.hpp>
Public Member Functions | |
OversetFindPointsGSLIB () | |
OversetFindPointsGSLIB (MPI_Comm comm_) | |
void | Setup (Mesh &m, const int meshid, GridFunction *gfmax=NULL, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256) |
void | FindPoints (const Vector &point_pos, Array< unsigned int > &point_id, int point_pos_ordering=Ordering::byNODES) |
void | Interpolate (const Vector &point_pos, Array< unsigned int > &point_id, const GridFunction &field_in, Vector &field_out, int point_pos_ordering=Ordering::byNODES) |
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) |
Public Member Functions inherited from mfem::FindPointsGSLIB | |
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 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 |
Protected Attributes | |
bool | overset |
unsigned int | u_meshid |
Vector | distfint |
Protected Attributes inherited from mfem::FindPointsGSLIB | |
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 |
Additional Inherited Members | |
Public Types inherited from mfem::FindPointsGSLIB | |
enum | AvgType { NONE, ARITHMETIC, HARMONIC } |
Protected Member Functions inherited from mfem::FindPointsGSLIB | |
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 () |
OversetFindPointsGSLIB enables use of findpts for arbitrary number of overlapping grids. The parameters in this class are the same as FindPointsGSLIB with the difference of additional inputs required to account for more than 1 mesh.
|
inline |
|
inline |
void mfem::OversetFindPointsGSLIB::FindPoints | ( | const Vector & | point_pos, |
Array< unsigned int > & | point_id, | ||
int | point_pos_ordering = Ordering::byNODES |
||
) |
Searches positions given in physical space by point_pos. All output Arrays and Vectors are expected to have the correct size.
[in] | point_pos | Positions to be found. |
[in] | point_id | Index of the mesh that the point belongs to (corresponding to meshid in Setup). |
[in] | point_pos_ordering | Ordering of the points: byNodes: (XXX...,YYY...,ZZZ) or byVDim: (XYZ,XYZ,....XYZ) |
void mfem::OversetFindPointsGSLIB::Interpolate | ( | const Vector & | point_pos, |
Array< unsigned int > & | point_id, | ||
const GridFunction & | field_in, | ||
Vector & | field_out, | ||
int | point_pos_ordering = Ordering::byNODES |
||
) |
void mfem::FindPointsGSLIB::Interpolate |
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 |
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. |
void mfem::FindPointsGSLIB::Interpolate |
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.
void mfem::OversetFindPointsGSLIB::Setup | ( | Mesh & | m, |
const int | meshid, | ||
GridFunction * | gfmax = NULL , |
||
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 meshes (L2). Note: the input mesh m must have Nodes set.
[in] | m | Input mesh. |
[in] | meshid | A unique # for each overlapping mesh. This id is used to make sure that points being searched are not looked for in the mesh that they belong to. |
[in] | gfmax | (Optional) GridFunction in H1 that is used as a discriminator when one point is located in multiple meshes. The mesh that maximizes gfmax is chosen. For example, using the distance field based on the overlapping boundaries is helpful for convergence during Schwarz iterations. |
[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. |
|
protected |