MFEM  v4.5.2
Finite element discretization library
gslib.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #ifndef MFEM_GSLIB
13 #define MFEM_GSLIB
14 
15 #include "../config/config.hpp"
16 #include "gridfunc.hpp"
17 
18 #ifdef MFEM_USE_GSLIB
19 
20 namespace gslib
21 {
22 struct comm;
23 struct findpts_data_2;
24 struct findpts_data_3;
25 struct crystal;
26 }
27 
28 namespace mfem
29 {
30 
31 /** \brief FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary
32  * collection of points. There are three key functions in FindPointsGSLIB:
33  *
34  * 1. Setup - constructs the internal data structures of gslib.
35  *
36  * 2. FindPoints - for any given arbitrary set of points in physical space,
37  * gslib finds the element number, MPI rank, and the reference space
38  * coordinates inside the element that each point is located in. gslib also
39  * returns a code that indicates whether the point was found inside an
40  * element, on element border, or not found in the domain.
41  * For points returned as found on `element border`, the point is either
42  * on an element edge/face or near the domain boundary, and gslib also
43  * returns a distance to the border. Points near (but outside) the domain
44  * boundary must then be marked as not found using the distance returned
45  * by gslib.
46  *
47  * 3. Interpolate - Interpolates any grid function at the points found using 2.
48  *
49  * FindPointsGSLIB provides interface to use these functions individually or
50  * using a single call.
51  */
53 {
54 public:
55  enum AvgType {NONE, ARITHMETIC, HARMONIC}; // Average type for L2 functions
56 
57 protected:
59  Array<Mesh *> mesh_split; // Meshes used to split simplices.
60  Array<IntegrationRule *> ir_split; // IntegrationRules for simplex->Quad/Hex
62  fes_rst_map; // FESpaces to map info Quad/Hex->Simplex
63  Array<GridFunction *> gf_rst_map; // GridFunctions to map info Quad/Hex->Simplex
65  struct gslib::findpts_data_2 *fdata2D; // gslib's internal data
66  struct gslib::findpts_data_3 *fdata3D; // gslib's internal data
67  struct gslib::crystal *cr; // gslib's internal data
68  struct gslib::comm *gsl_comm; // gslib's internal data
72  bool setupflag; // flag to indicate whether gslib data has been setup
73  double default_interp_value; // used for points that are not found in the mesh
74  AvgType avgtype; // average type used for L2 functions
78  // Tolerance to ignore points just outside elements at the boundary.
79  double bdr_tol;
80 
81  /// Use GSLIB for communication and interpolation
82  virtual void InterpolateH1(const GridFunction &field_in, Vector &field_out);
83  /// Uses GSLIB Crystal Router for communication followed by MFEM's
84  /// interpolation functions
85  virtual void InterpolateGeneral(const GridFunction &field_in,
86  Vector &field_out);
87 
88  /// Since GSLIB is designed to work with quads/hexes, we split every
89  /// triangle/tet/prism/pyramid element into quads/hexes.
90  virtual void SetupSplitMeshes();
91 
92  /// Setup integration points that will be used to interpolate the nodal
93  /// location at points expected by GSLIB.
95  IntegrationRule *irule,
96  int order);
97 
98  /// Get GridFunction value at the points expected by GSLIB.
99  virtual void GetNodalValues(const GridFunction *gf_in, Vector &node_vals);
100 
101  /// Map {r,s,t} coordinates from [-1,1] to [0,1] for MFEM. For simplices,
102  /// find the original element number (that was split into micro quads/hexes)
103  /// during the setup phase.
104  virtual void MapRefPosAndElemIndices();
105 
106 public:
107  FindPointsGSLIB();
108 
109 #ifdef MFEM_USE_MPI
110  FindPointsGSLIB(MPI_Comm comm_);
111 #endif
112 
113  virtual ~FindPointsGSLIB();
114 
115  /** Initializes the internal mesh in gslib, by sending the positions of the
116  Gauss-Lobatto nodes of the input Mesh object @a m.
117  Note: not tested with periodic (L2).
118  Note: the input mesh @a m must have Nodes set.
119 
120  @param[in] m Input mesh.
121  @param[in] bb_t (Optional) Relative size of bounding box around
122  each element.
123  @param[in] newt_tol (Optional) Newton tolerance for the gslib
124  search methods.
125  @param[in] npt_max (Optional) Number of points for simultaneous
126  iteration. This alters performance and
127  memory footprint.*/
128  void Setup(Mesh &m, const double bb_t = 0.1,
129  const double newt_tol = 1.0e-12,
130  const int npt_max = 256);
131  /** Searches positions given in physical space by @a point_pos.
132  These positions can be ordered byNodes: (XXX...,YYY...,ZZZ) or
133  byVDim: (XYZ,XYZ,....XYZ) specified by @a point_pos_ordering.
134  This function populates the following member variables:
135  #gsl_code Return codes for each point: inside element (0),
136  element boundary (1), not found (2).
137  #gsl_proc MPI proc ids where the points were found.
138  #gsl_elem Element ids where the points were found.
139  Defaults to 0 for points that were not found.
140  #gsl_mfem_elem Element ids corresponding to MFEM-mesh where the points
141  were found. #gsl_mfem_elem != #gsl_elem for simplices
142  Defaults to 0 for points that were not found.
143  #gsl_ref Reference coordinates of the found point.
144  Ordered by vdim (XYZ,XYZ,XYZ...). Defaults to -1 for
145  points that were not found. Note: the gslib reference
146  frame is [-1,1].
147  #gsl_mfem_ref Reference coordinates #gsl_ref mapped to [0,1].
148  Defaults to 0 for points that were not found.
149  #gsl_dist Distance between the sought and the found point
150  in physical space. */
151  void FindPoints(const Vector &point_pos,
152  int point_pos_ordering = Ordering::byNODES);
153  /// Setup FindPoints and search positions
154  void FindPoints(Mesh &m, const Vector &point_pos,
155  int point_pos_ordering = Ordering::byNODES,
156  const double bb_t = 0.1,
157  const double newt_tol = 1.0e-12, const int npt_max = 256);
158 
159  /** Interpolation of field values at prescribed reference space positions.
160  @param[in] field_in Function values that will be interpolated on the
161  reference positions. Note: it is assumed that
162  @a field_in is in H1 and in the same space as the
163  mesh that was given to Setup().
164  @param[out] field_out Interpolated values. For points that are not found
165  the value is set to #default_interp_value. */
166  virtual void Interpolate(const GridFunction &field_in, Vector &field_out);
167  /** Search positions and interpolate. The ordering (byNODES or byVDIM) of
168  the output values in @a field_out corresponds to the ordering used
169  in the input GridFunction @a field_in. */
170  void Interpolate(const Vector &point_pos, const GridFunction &field_in,
171  Vector &field_out,
172  int point_pos_ordering = Ordering::byNODES);
173  /** Setup FindPoints, search positions and interpolate. The ordering (byNODES
174  or byVDIM) of the output values in @a field_out corresponds to the
175  ordering used in the input GridFunction @a field_in. */
176  void Interpolate(Mesh &m, const Vector &point_pos,
177  const GridFunction &field_in, Vector &field_out,
178  int point_pos_ordering = Ordering::byNODES);
179 
180  /// Average type to be used for L2 functions in-case a point is located at
181  /// an element boundary where the function might be multi-valued.
182  virtual void SetL2AvgType(AvgType avgtype_) { avgtype = avgtype_; }
183 
184  /// Set the default interpolation value for points that are not found in the
185  /// mesh.
186  virtual void SetDefaultInterpolationValue(double interp_value_)
187  {
188  default_interp_value = interp_value_;
189  }
190 
191  /// Set the tolerance for detecting points outside the 'curvilinear' boundary
192  /// that gslib may return as found on the boundary. Points found on boundary
193  /// with distance greater than @ bdr_tol are marked as not found.
194  virtual void SetDistanceToleranceForPointsFoundOnBoundary(double bdr_tol_)
195  {
196  bdr_tol = bdr_tol_;
197  }
198 
199  /** Cleans up memory allocated internally by gslib.
200  Note that in parallel, this must be called before MPI_Finalize(), as it
201  calls MPI_Comm_free() for internal gslib communicators. */
202  virtual void FreeData();
203 
204  /// Return code for each point searched by FindPoints: inside element (0), on
205  /// element boundary (1), or not found (2).
206  virtual const Array<unsigned int> &GetCode() const { return gsl_code; }
207  /// Return element number for each point found by FindPoints.
208  virtual const Array<unsigned int> &GetElem() const { return gsl_mfem_elem; }
209  /// Return MPI rank on which each point was found by FindPoints.
210  virtual const Array<unsigned int> &GetProc() const { return gsl_proc; }
211  /// Return reference coordinates for each point found by FindPoints.
212  virtual const Vector &GetReferencePosition() const { return gsl_mfem_ref; }
213  /// Return distance between the sought and the found point in physical space,
214  /// for each point found by FindPoints.
215  virtual const Vector &GetDist() const { return gsl_dist; }
216 
217  /// Return element number for each point found by FindPoints corresponding to
218  /// GSLIB mesh. gsl_mfem_elem != gsl_elem for mesh with simplices.
219  virtual const Array<unsigned int> &GetGSLIBElem() const { return gsl_elem; }
220  /// Return reference coordinates in [-1,1] (internal range in GSLIB) for each
221  /// point found by FindPoints.
222  virtual const Vector &GetGSLIBReferencePosition() const { return gsl_ref; }
223 };
224 
225 /** \brief OversetFindPointsGSLIB enables use of findpts for arbitrary number of
226  overlapping grids.
227  The parameters in this class are the same as FindPointsGSLIB with the
228  difference of additional inputs required to account for more than 1 mesh. */
230 {
231 protected:
232  bool overset;
233  unsigned int u_meshid;
234  Vector distfint; // Used to store nodal vals of grid func. passed to findpts
235 
236 public:
238  overset(true) { }
239 
240 #ifdef MFEM_USE_MPI
241  OversetFindPointsGSLIB(MPI_Comm comm_) : FindPointsGSLIB(comm_),
242  overset(true) { }
243 #endif
244 
245  /** Initializes the internal mesh in gslib, by sending the positions of the
246  Gauss-Lobatto nodes of the input Mesh object @a m.
247  Note: not tested with periodic meshes (L2).
248  Note: the input mesh @a m must have Nodes set.
249 
250  @param[in] m Input mesh.
251  @param[in] meshid A unique # for each overlapping mesh. This id is
252  used to make sure that points being searched are not
253  looked for in the mesh that they belong to.
254  @param[in] gfmax (Optional) GridFunction in H1 that is used as a
255  discriminator when one point is located in multiple
256  meshes. The mesh that maximizes gfmax is chosen.
257  For example, using the distance field based on the
258  overlapping boundaries is helpful for convergence
259  during Schwarz iterations.
260  @param[in] bb_t (Optional) Relative size of bounding box around
261  each element.
262  @param[in] newt_tol (Optional) Newton tolerance for the gslib
263  search methods.
264  @param[in] npt_max (Optional) Number of points for simultaneous
265  iteration. This alters performance and
266  memory footprint.*/
267  void Setup(Mesh &m, const int meshid, GridFunction *gfmax = NULL,
268  const double bb_t = 0.1, const double newt_tol = 1.0e-12,
269  const int npt_max = 256);
270 
271  /** Searches positions given in physical space by @a point_pos. All output
272  Arrays and Vectors are expected to have the correct size.
273 
274  @param[in] point_pos Positions to be found.
275  @param[in] point_id Index of the mesh that the point belongs
276  to (corresponding to @a meshid in Setup).
277  @param[in] point_pos_ordering Ordering of the points:
278  byNodes: (XXX...,YYY...,ZZZ) or
279  byVDim: (XYZ,XYZ,....XYZ) */
280  void FindPoints(const Vector &point_pos,
281  Array<unsigned int> &point_id,
282  int point_pos_ordering = Ordering::byNODES);
283 
284  /** Search positions and interpolate */
285  void Interpolate(const Vector &point_pos, Array<unsigned int> &point_id,
286  const GridFunction &field_in, Vector &field_out,
287  int point_pos_ordering = Ordering::byNODES);
289 };
290 
291 } // namespace mfem
292 
293 #endif // MFEM_USE_GSLIB
294 
295 #endif // MFEM_GSLIB
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:803
virtual const Array< unsigned int > & GetProc() const
Return MPI rank on which each point was found by FindPoints.
Definition: gslib.hpp:210
Array< unsigned int > gsl_elem
Definition: gslib.hpp:70
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
virtual const Array< unsigned int > & GetGSLIBElem() const
Definition: gslib.hpp:219
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
virtual ~FindPointsGSLIB()
Definition: gslib.cpp:67
void FindPoints(const Vector &point_pos, int point_pos_ordering=Ordering::byNODES)
Definition: gslib.cpp:171
Array< GridFunction * > gf_rst_map
Definition: gslib.hpp:63
Array< unsigned int > gsl_mfem_elem
Definition: gslib.hpp:70
virtual void SetupSplitMeshes()
Definition: gslib.cpp:296
void Interpolate(const Vector &point_pos, Array< unsigned int > &point_id, const GridFunction &field_in, Vector &field_out, int point_pos_ordering=Ordering::byNODES)
Definition: gslib.cpp:1267
OversetFindPointsGSLIB(MPI_Comm comm_)
Definition: gslib.hpp:241
Array< unsigned int > gsl_proc
Definition: gslib.hpp:70
Array< Mesh * > mesh_split
Definition: gslib.hpp:59
virtual void SetL2AvgType(AvgType avgtype_)
Definition: gslib.hpp:182
virtual const Array< unsigned int > & GetElem() const
Return element number for each point found by FindPoints.
Definition: gslib.hpp:208
virtual void MapRefPosAndElemIndices()
Definition: gslib.cpp:656
struct gslib::findpts_data_3 * fdata3D
Definition: gslib.hpp:66
virtual void GetNodalValues(const GridFunction *gf_in, Vector &node_vals)
Get GridFunction value at the points expected by GSLIB.
Definition: gslib.cpp:562
Array< IntegrationRule * > ir_split
Definition: gslib.hpp:60
virtual const Vector & GetGSLIBReferencePosition() const
Definition: gslib.hpp:222
OversetFindPointsGSLIB enables use of findpts for arbitrary number of overlapping grids...
Definition: gslib.hpp:229
virtual void FreeData()
Definition: gslib.cpp:267
void Setup(Mesh &m, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
Definition: gslib.cpp:107
double default_interp_value
Definition: gslib.hpp:73
struct gslib::crystal * cr
Definition: gslib.hpp:67
Array< int > split_element_index
Definition: gslib.hpp:76
Array< FiniteElementSpace * > fes_rst_map
Definition: gslib.hpp:62
Array< unsigned int > gsl_code
Definition: gslib.hpp:70
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
virtual void InterpolateGeneral(const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:948
virtual void SetupIntegrationRuleForSplitMesh(Mesh *mesh, IntegrationRule *irule, int order)
Definition: gslib.cpp:514
struct gslib::findpts_data_2 * fdata2D
Definition: gslib.hpp:65
virtual void SetDefaultInterpolationValue(double interp_value_)
Definition: gslib.hpp:186
virtual void SetDistanceToleranceForPointsFoundOnBoundary(double bdr_tol_)
Definition: gslib.hpp:194
Array< int > split_element_map
Definition: gslib.hpp:75
virtual const Vector & GetDist() const
Definition: gslib.hpp:215
struct gslib::comm * gsl_comm
Definition: gslib.hpp:68
FiniteElementCollection * fec_map_lin
Definition: gslib.hpp:64
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points...
Definition: gslib.hpp:52
Vector data type.
Definition: vector.hpp:60
Definition: gslib.cpp:23
virtual const Array< unsigned int > & GetCode() const
Definition: gslib.hpp:206
void FindPoints(const Vector &point_pos, Array< unsigned int > &point_id, int point_pos_ordering=Ordering::byNODES)
Definition: gslib.cpp:1189
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)
Definition: gslib.cpp:1108
virtual void InterpolateH1(const GridFunction &field_in, Vector &field_out)
Use GSLIB for communication and interpolation.
Definition: gslib.cpp:885
virtual const Vector & GetReferencePosition() const
Return reference coordinates for each point found by FindPoints.
Definition: gslib.hpp:212