#include <gslib.hpp>
|
| FindPointsGSLIB () |
|
| FindPointsGSLIB (MPI_Comm _comm) |
|
| ~FindPointsGSLIB () |
|
void | Setup (Mesh &m, double bb_t, double newt_tol, int npt_max) |
|
void | FindPoints (Vector &point_pos, Array< unsigned int > &codes, Array< unsigned int > &proc_ids, Array< unsigned int > &elem_ids, Vector &ref_pos, Vector &dist) |
|
void | Interpolate (Array< unsigned int > &codes, Array< unsigned int > &proc_ids, Array< unsigned int > &elem_ids, Vector &ref_pos, const GridFunction &field_in, Vector &field_out) |
|
void | FreeData () |
|
Definition at line 27 of file gslib.hpp.
mfem::FindPointsGSLIB::FindPointsGSLIB |
( |
| ) |
|
mfem::FindPointsGSLIB::FindPointsGSLIB |
( |
MPI_Comm |
_comm | ) |
|
mfem::FindPointsGSLIB::~FindPointsGSLIB |
( |
| ) |
|
void mfem::FindPointsGSLIB::FindPoints |
( |
Vector & |
point_pos, |
|
|
Array< unsigned int > & |
codes, |
|
|
Array< unsigned int > & |
proc_ids, |
|
|
Array< unsigned int > & |
elem_ids, |
|
|
Vector & |
ref_pos, |
|
|
Vector & |
dist |
|
) |
| |
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). |
[out] | codes | Return codes for each point: inside element (0), element boundary (1), not found (2). |
[out] | proc_ids | MPI proc ids where the points were found. |
[out] | elem_ids | Element ids where the points were found. |
[out] | ref_pos | Reference coordinates of the found point. Ordered by vdim (XYZ,XYZ,XYZ...). Note: the gslib reference frame is [-1,1]. |
[out] | dist | Distance between the seeked and the found point in physical space. |
Definition at line 115 of file gslib.cpp.
void mfem::FindPointsGSLIB::FreeData |
( |
| ) |
|
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 185 of file gslib.cpp.
void mfem::FindPointsGSLIB::GetNodeValues |
( |
const GridFunction & |
gf_in, |
|
|
Vector & |
node_vals |
|
) |
| |
|
protected |
void mfem::FindPointsGSLIB::Interpolate |
( |
Array< unsigned int > & |
codes, |
|
|
Array< unsigned int > & |
proc_ids, |
|
|
Array< unsigned int > & |
elem_ids, |
|
|
Vector & |
ref_pos, |
|
|
const GridFunction & |
field_in, |
|
|
Vector & |
field_out |
|
) |
| |
Interpolation of field values at prescribed reference space positions.
- Parameters
-
[in] | codes | Return codes for each point: inside element (0), element boundary (1), not found (2). |
[in] | proc_ids | MPI proc ids where the points were found. |
[in] | elem_ids | Element ids where the points were found. |
[in] | ref_pos | Reference coordinates of the found point. Ordered by vdim (XYZ,XYZ,XYZ...). Note: the gslib reference frame is [-1,1]. |
[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. |
Definition at line 155 of file gslib.cpp.
void mfem::FindPointsGSLIB::Setup |
( |
Mesh & |
m, |
|
|
double |
bb_t, |
|
|
double |
newt_tol, |
|
|
int |
npt_max |
|
) |
| |
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 (DG meshes). Note: the input mesh m must have Nodes set.
- Parameters
-
[in] | m | Input mesh. |
[in] | bb_t | Relative size of bounding box around each element. |
[in] | newt_tol | Newton tolerance for the gslib search methods. |
[in] | npt_max | Number of points for simultaneous iteration. This alters performance and memory footprint. |
Definition at line 58 of file gslib.cpp.
int mfem::FindPointsGSLIB::dim |
|
protected |
struct findpts_data_2* mfem::FindPointsGSLIB::fdata2D |
|
protected |
struct findpts_data_3* mfem::FindPointsGSLIB::fdata3D |
|
protected |
struct comm* mfem::FindPointsGSLIB::gsl_comm |
|
protected |
Vector mfem::FindPointsGSLIB::gsl_mesh |
|
protected |
Mesh* mfem::FindPointsGSLIB::mesh |
|
protected |
The documentation for this class was generated from the following files: