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>
|
| 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) |
|
void | Interpolate (const Vector &point_pos, Array< unsigned int > &point_id, const GridFunction &field_in, Vector &field_out) |
|
| 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 |
|
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.
Definition at line 192 of file gslib.hpp.
mfem::OversetFindPointsGSLIB::OversetFindPointsGSLIB |
( |
| ) |
|
|
inline |
mfem::OversetFindPointsGSLIB::OversetFindPointsGSLIB |
( |
MPI_Comm |
comm_ | ) |
|
|
inline |
void mfem::OversetFindPointsGSLIB::FindPoints |
( |
const Vector & |
point_pos, |
|
|
Array< unsigned int > & |
point_id |
|
) |
| |
Searches positions given in physical space by point_pos. All output Arrays and Vectors are expected to have the correct size.
- Parameters
-
[in] | point_pos | Positions to be found. Must by ordered by nodes (XXX...,YYY...,ZZZ). |
[in] | point_id | Index of the mesh that the point belongs to (corresponding to meshid in Setup). |
Definition at line 919 of file gslib.cpp.
void mfem::OversetFindPointsGSLIB::Interpolate |
( |
const Vector & |
point_pos, |
|
|
Array< unsigned int > & |
point_id, |
|
|
const GridFunction & |
field_in, |
|
|
Vector & |
field_out |
|
) |
| |
Search positions and interpolate
Definition at line 986 of file gslib.cpp.
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.
- Parameters
-
[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. |
Definition at line 847 of file gslib.cpp.
Vector mfem::OversetFindPointsGSLIB::distfint |
|
protected |
bool mfem::OversetFindPointsGSLIB::overset |
|
protected |
unsigned int mfem::OversetFindPointsGSLIB::u_meshid |
|
protected |
The documentation for this class was generated from the following files: