MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
gslib.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, 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  *
42  * 3. Interpolate - Interpolates any grid function at the points found using 2.
43  *
44  * FindPointsGSLIB provides interface to use these functions individually or
45  * using a single call.
46  */
48 {
49 public:
50  enum AvgType {NONE, ARITHMETIC, HARMONIC}; // Average type for L2 functions
51 
52 protected:
54  Array<Mesh *> mesh_split; // Meshes used to split simplices.
55  Array<IntegrationRule *> ir_split; // IntegrationRules for simplex->Quad/Hex
57  fes_rst_map; // FESpaces to map info Quad/Hex->Simplex
58  Array<GridFunction *> gf_rst_map; // GridFunctions to map info Quad/Hex->Simplex
60  struct gslib::findpts_data_2 *fdata2D; // gslib's internal data
61  struct gslib::findpts_data_3 *fdata3D; // gslib's internal data
62  struct gslib::crystal *cr; // gslib's internal data
63  struct gslib::comm *gsl_comm; // gslib's internal data
67  bool setupflag; // flag to indicate whether gslib data has been setup
68  double default_interp_value; // used for points that are not found in the mesh
69  AvgType avgtype; // average type used for L2 functions
73 
74  /// Use GSLIB for communication and interpolation
75  virtual void InterpolateH1(const GridFunction &field_in, Vector &field_out);
76  /// Uses GSLIB Crystal Router for communication followed by MFEM's
77  /// interpolation functions
78  virtual void InterpolateGeneral(const GridFunction &field_in,
79  Vector &field_out);
80 
81  /// Since GSLIB is designed to work with quads/hexes, we split every
82  /// triangle/tet/prism/pyramid element into quads/hexes.
83  virtual void SetupSplitMeshes();
84 
85  /// Setup integration points that will be used to interpolate the nodal
86  /// location at points expected by GSLIB.
88  IntegrationRule *irule,
89  int order);
90 
91  /// Get GridFunction value at the points expected by GSLIB.
92  virtual void GetNodalValues(const GridFunction *gf_in, Vector &node_vals);
93 
94  /// Map {r,s,t} coordinates from [-1,1] to [0,1] for MFEM. For simplices,
95  /// find the original element number (that was split into micro quads/hexes)
96  /// during the setup phase.
97  virtual void MapRefPosAndElemIndices();
98 
99 public:
100  FindPointsGSLIB();
101 
102 #ifdef MFEM_USE_MPI
103  FindPointsGSLIB(MPI_Comm comm_);
104 #endif
105 
106  virtual ~FindPointsGSLIB();
107 
108  /** Initializes the internal mesh in gslib, by sending the positions of the
109  Gauss-Lobatto nodes of the input Mesh object @a m.
110  Note: not tested with periodic (L2).
111  Note: the input mesh @a m must have Nodes set.
112 
113  @param[in] m Input mesh.
114  @param[in] bb_t (Optional) Relative size of bounding box around
115  each element.
116  @param[in] newt_tol (Optional) Newton tolerance for the gslib
117  search methods.
118  @param[in] npt_max (Optional) Number of points for simultaneous
119  iteration. This alters performance and
120  memory footprint.*/
121  void Setup(Mesh &m, const double bb_t = 0.1,
122  const double newt_tol = 1.0e-12,
123  const int npt_max = 256);
124  /** Searches positions given in physical space by @a point_pos.
125  These positions can be ordered byNodes: (XXX...,YYY...,ZZZ) or
126  byVDim: (XYZ,XYZ,....XYZ) specified by @a point_pos_ordering.
127  This function populates the following member variables:
128  #gsl_code Return codes for each point: inside element (0),
129  element boundary (1), not found (2).
130  #gsl_proc MPI proc ids where the points were found.
131  #gsl_elem Element ids where the points were found.
132  Defaults to 0 for points that were not found.
133  #gsl_mfem_elem Element ids corresponding to MFEM-mesh where the points
134  were found. #gsl_mfem_elem != #gsl_elem for simplices
135  Defaults to 0 for points that were not found.
136  #gsl_ref Reference coordinates of the found point.
137  Ordered by vdim (XYZ,XYZ,XYZ...). Defaults to -1 for
138  points that were not found. Note: the gslib reference
139  frame is [-1,1].
140  #gsl_mfem_ref Reference coordinates #gsl_ref mapped to [0,1].
141  Defaults to 0 for points that were not found.
142  #gsl_dist Distance between the sought and the found point
143  in physical space. */
144  void FindPoints(const Vector &point_pos,
145  int point_pos_ordering = Ordering::byNODES);
146  /// Setup FindPoints and search positions
147  void FindPoints(Mesh &m, const Vector &point_pos,
148  int point_pos_ordering = Ordering::byNODES,
149  const double bb_t = 0.1,
150  const double newt_tol = 1.0e-12, const int npt_max = 256);
151 
152  /** Interpolation of field values at prescribed reference space positions.
153  @param[in] field_in Function values that will be interpolated on the
154  reference positions. Note: it is assumed that
155  @a field_in is in H1 and in the same space as the
156  mesh that was given to Setup().
157  @param[out] field_out Interpolated values. For points that are not found
158  the value is set to #default_interp_value. */
159  virtual void Interpolate(const GridFunction &field_in, Vector &field_out);
160  /** Search positions and interpolate. The ordering (byNODES or byVDIM) of
161  the output values in @a field_out corresponds to the ordering used
162  in the input GridFunction @a field_in. */
163  void Interpolate(const Vector &point_pos, const GridFunction &field_in,
164  Vector &field_out,
165  int point_pos_ordering = Ordering::byNODES);
166  /** Setup FindPoints, search positions and interpolate. The ordering (byNODES
167  or byVDIM) of the output values in @a field_out corresponds to the
168  ordering used in the input GridFunction @a field_in. */
169  void Interpolate(Mesh &m, const Vector &point_pos,
170  const GridFunction &field_in, Vector &field_out,
171  int point_pos_ordering = Ordering::byNODES);
172 
173  /// Average type to be used for L2 functions in-case a point is located at
174  /// an element boundary where the function might be multi-valued.
175  virtual void SetL2AvgType(AvgType avgtype_) { avgtype = avgtype_; }
176 
177  /// Set the default interpolation value for points that are not found in the
178  /// mesh.
179  virtual void SetDefaultInterpolationValue(double interp_value_)
180  {
181  default_interp_value = interp_value_;
182  }
183 
184  /** Cleans up memory allocated internally by gslib.
185  Note that in parallel, this must be called before MPI_Finalize(), as it
186  calls MPI_Comm_free() for internal gslib communicators. */
187  virtual void FreeData();
188 
189  /// Return code for each point searched by FindPoints: inside element (0), on
190  /// element boundary (1), or not found (2).
191  virtual const Array<unsigned int> &GetCode() const { return gsl_code; }
192  /// Return element number for each point found by FindPoints.
193  virtual const Array<unsigned int> &GetElem() const { return gsl_mfem_elem; }
194  /// Return MPI rank on which each point was found by FindPoints.
195  virtual const Array<unsigned int> &GetProc() const { return gsl_proc; }
196  /// Return reference coordinates for each point found by FindPoints.
197  virtual const Vector &GetReferencePosition() const { return gsl_mfem_ref; }
198  /// Return distance between the sought and the found point in physical space,
199  /// for each point found by FindPoints.
200  virtual const Vector &GetDist() const { return gsl_dist; }
201 
202  /// Return element number for each point found by FindPoints corresponding to
203  /// GSLIB mesh. gsl_mfem_elem != gsl_elem for mesh with simplices.
204  virtual const Array<unsigned int> &GetGSLIBElem() const { return gsl_elem; }
205  /// Return reference coordinates in [-1,1] (internal range in GSLIB) for each
206  /// point found by FindPoints.
207  virtual const Vector &GetGSLIBReferencePosition() const { return gsl_ref; }
208 };
209 
210 /** \brief OversetFindPointsGSLIB enables use of findpts for arbitrary number of
211  overlapping grids.
212  The parameters in this class are the same as FindPointsGSLIB with the
213  difference of additional inputs required to account for more than 1 mesh. */
215 {
216 protected:
217  bool overset;
218  unsigned int u_meshid;
219  Vector distfint; // Used to store nodal vals of grid func. passed to findpts
220 
221 public:
223  overset(true) { }
224 
225 #ifdef MFEM_USE_MPI
226  OversetFindPointsGSLIB(MPI_Comm comm_) : FindPointsGSLIB(comm_),
227  overset(true) { }
228 #endif
229 
230  /** Initializes the internal mesh in gslib, by sending the positions of the
231  Gauss-Lobatto nodes of the input Mesh object @a m.
232  Note: not tested with periodic meshes (L2).
233  Note: the input mesh @a m must have Nodes set.
234 
235  @param[in] m Input mesh.
236  @param[in] meshid A unique # for each overlapping mesh. This id is
237  used to make sure that points being searched are not
238  looked for in the mesh that they belong to.
239  @param[in] gfmax (Optional) GridFunction in H1 that is used as a
240  discriminator when one point is located in multiple
241  meshes. The mesh that maximizes gfmax is chosen.
242  For example, using the distance field based on the
243  overlapping boundaries is helpful for convergence
244  during Schwarz iterations.
245  @param[in] bb_t (Optional) Relative size of bounding box around
246  each element.
247  @param[in] newt_tol (Optional) Newton tolerance for the gslib
248  search methods.
249  @param[in] npt_max (Optional) Number of points for simultaneous
250  iteration. This alters performance and
251  memory footprint.*/
252  void Setup(Mesh &m, const int meshid, GridFunction *gfmax = NULL,
253  const double bb_t = 0.1, const double newt_tol = 1.0e-12,
254  const int npt_max = 256);
255 
256  /** Searches positions given in physical space by @a point_pos. All output
257  Arrays and Vectors are expected to have the correct size.
258 
259  @param[in] point_pos Positions to be found.
260  @param[in] point_id Index of the mesh that the point belongs
261  to (corresponding to @a meshid in Setup).
262  @param[in] point_pos_ordering Ordering of the points:
263  byNodes: (XXX...,YYY...,ZZZ) or
264  byVDim: (XYZ,XYZ,....XYZ) */
265  void FindPoints(const Vector &point_pos,
266  Array<unsigned int> &point_id,
267  int point_pos_ordering = Ordering::byNODES);
268 
269  /** Search positions and interpolate */
270  void Interpolate(const Vector &point_pos, Array<unsigned int> &point_id,
271  const GridFunction &field_in, Vector &field_out,
272  int point_pos_ordering = Ordering::byNODES);
274 };
275 
276 } // namespace mfem
277 
278 #endif // MFEM_USE_GSLIB
279 
280 #endif // MFEM_GSLIB
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:801
Array< unsigned int > gsl_elem
Definition: gslib.hpp:65
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
virtual const Vector & GetDist() const
Definition: gslib.hpp:200
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:58
Array< unsigned int > gsl_mfem_elem
Definition: gslib.hpp:65
virtual void SetupSplitMeshes()
Definition: gslib.cpp:294
virtual const Vector & GetReferencePosition() const
Return reference coordinates for each point found by FindPoints.
Definition: gslib.hpp:197
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:1263
OversetFindPointsGSLIB(MPI_Comm comm_)
Definition: gslib.hpp:226
Array< unsigned int > gsl_proc
Definition: gslib.hpp:65
Array< Mesh * > mesh_split
Definition: gslib.hpp:54
virtual void SetL2AvgType(AvgType avgtype_)
Definition: gslib.hpp:175
virtual void MapRefPosAndElemIndices()
Definition: gslib.cpp:654
struct gslib::findpts_data_3 * fdata3D
Definition: gslib.hpp:61
virtual void GetNodalValues(const GridFunction *gf_in, Vector &node_vals)
Get GridFunction value at the points expected by GSLIB.
Definition: gslib.cpp:560
Array< IntegrationRule * > ir_split
Definition: gslib.hpp:55
OversetFindPointsGSLIB enables use of findpts for arbitrary number of overlapping grids...
Definition: gslib.hpp:214
virtual void FreeData()
Definition: gslib.cpp:265
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:68
struct gslib::crystal * cr
Definition: gslib.hpp:62
Array< int > split_element_index
Definition: gslib.hpp:71
Array< FiniteElementSpace * > fes_rst_map
Definition: gslib.hpp:57
Array< unsigned int > gsl_code
Definition: gslib.hpp:65
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:946
virtual void SetupIntegrationRuleForSplitMesh(Mesh *mesh, IntegrationRule *irule, int order)
Definition: gslib.cpp:512
virtual const Array< unsigned int > & GetProc() const
Return MPI rank on which each point was found by FindPoints.
Definition: gslib.hpp:195
struct gslib::findpts_data_2 * fdata2D
Definition: gslib.hpp:60
virtual void SetDefaultInterpolationValue(double interp_value_)
Definition: gslib.hpp:179
virtual const Array< unsigned int > & GetElem() const
Return element number for each point found by FindPoints.
Definition: gslib.hpp:193
Array< int > split_element_map
Definition: gslib.hpp:70
virtual const Array< unsigned int > & GetCode() const
Definition: gslib.hpp:191
struct gslib::comm * gsl_comm
Definition: gslib.hpp:63
FiniteElementCollection * fec_map_lin
Definition: gslib.hpp:59
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points...
Definition: gslib.hpp:47
Vector data type.
Definition: vector.hpp:60
virtual const Array< unsigned int > & GetGSLIBElem() const
Definition: gslib.hpp:204
void FindPoints(const Vector &point_pos, Array< unsigned int > &point_id, int point_pos_ordering=Ordering::byNODES)
Definition: gslib.cpp:1187
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:1106
virtual const Vector & GetGSLIBReferencePosition() const
Definition: gslib.hpp:207
virtual void InterpolateH1(const GridFunction &field_in, Vector &field_out)
Use GSLIB for communication and interpolation.
Definition: gslib.cpp:883