MFEM  v4.3.0
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-2021, 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  IntegrationRule *ir_simplex; // IntegrationRule to split quads/hex -> simplex
55  struct gslib::findpts_data_2 *fdata2D; // gslib's internal data
56  struct gslib::findpts_data_3 *fdata3D; // gslib's internal data
57  struct gslib::crystal *cr; // gslib's internal data
58  struct gslib::comm *gsl_comm; // gslib's internal data
62  bool setupflag; // flag to indicate whether gslib data has been setup
63  double default_interp_value; // used for points that are not found in the mesh
64  AvgType avgtype; // average type used for L2 functions
65 
66  /// Get GridFunction from MFEM format to GSLIB format
67  virtual void GetNodeValues(const GridFunction &gf_in, Vector &node_vals);
68  /// Get nodal coordinates from mesh to the format expected by GSLIB for quads
69  /// and hexes
70  virtual void GetQuadHexNodalCoordinates();
71  /// Convert simplices to quad/hexes and then get nodal coordinates for each
72  /// split element into format expected by GSLIB
73  virtual void GetSimplexNodalCoordinates();
74 
75  /// Use GSLIB for communication and interpolation
76  virtual void InterpolateH1(const GridFunction &field_in, Vector &field_out);
77  /// Uses GSLIB Crystal Router for communication followed by MFEM's
78  /// interpolation functions
79  virtual void InterpolateGeneral(const GridFunction &field_in,
80  Vector &field_out);
81  /// Map {r,s,t} coordinates from [-1,1] to [0,1] for MFEM. For simplices mesh
82  /// find the original element number (that was split into micro quads/hexes
83  /// by GetSimplexNodalCoordinates())
84  virtual void MapRefPosAndElemIndices();
85 
86 public:
88 
89 #ifdef MFEM_USE_MPI
90  FindPointsGSLIB(MPI_Comm comm_);
91 #endif
92 
93  virtual ~FindPointsGSLIB();
94 
95  /** Initializes the internal mesh in gslib, by sending the positions of the
96  Gauss-Lobatto nodes of the input Mesh object @a m.
97  Note: not tested with periodic (L2).
98  Note: the input mesh @a m must have Nodes set.
99 
100  @param[in] m Input mesh.
101  @param[in] bb_t (Optional) Relative size of bounding box around
102  each element.
103  @param[in] newt_tol (Optional) Newton tolerance for the gslib
104  search methods.
105  @param[in] npt_max (Optional) Number of points for simultaneous
106  iteration. This alters performance and
107  memory footprint.*/
108  void Setup(Mesh &m, const double bb_t = 0.1,
109  const double newt_tol = 1.0e-12,
110  const int npt_max = 256);
111  /** Searches positions given in physical space by @a point_pos. These positions
112  must by ordered by nodes: (XXX...,YYY...,ZZZ).
113  This function populates the following member variables:
114  #gsl_code Return codes for each point: inside element (0),
115  element boundary (1), not found (2).
116  #gsl_proc MPI proc ids where the points were found.
117  #gsl_elem Element ids where the points were found.
118  Defaults to 0 for points that were not found.
119  #gsl_mfem_elem Element ids corresponding to MFEM-mesh where the points
120  were found. #gsl_mfem_elem != #gsl_elem for simplices
121  Defaults to 0 for points that were not found.
122  #gsl_ref Reference coordinates of the found point.
123  Ordered by vdim (XYZ,XYZ,XYZ...). Defaults to -1 for
124  points that were not found. Note: the gslib reference
125  frame is [-1,1].
126  #gsl_mfem_ref Reference coordinates #gsl_ref mapped to [0,1].
127  Defaults to 0 for points that were not found.
128  #gsl_dist Distance between the sought and the found point
129  in physical space. */
130  void FindPoints(const Vector &point_pos);
131  /// Setup FindPoints and search positions
132  void FindPoints(Mesh &m, const Vector &point_pos,
133  const double bb_t = 0.1,
134  const double newt_tol = 1.0e-12, const int npt_max = 256);
135 
136  /** Interpolation of field values at prescribed reference space positions.
137  @param[in] field_in Function values that will be interpolated on the
138  reference positions. Note: it is assumed that
139  @a field_in is in H1 and in the same space as the
140  mesh that was given to Setup().
141  @param[out] field_out Interpolated values. For points that are not found
142  the value is set to #default_interp_value. */
143  virtual void Interpolate(const GridFunction &field_in, Vector &field_out);
144  /** Search positions and interpolate */
145  void Interpolate(const Vector &point_pos, const GridFunction &field_in,
146  Vector &field_out);
147  /** Setup FindPoints, search positions and interpolate */
148  void Interpolate(Mesh &m, const Vector &point_pos,
149  const GridFunction &field_in, Vector &field_out);
150 
151  /// Average type to be used for L2 functions in-case a point is located at
152  /// an element boundary where the function might be multi-valued.
153  virtual void SetL2AvgType(AvgType avgtype_) { avgtype = avgtype_; }
154 
155  /// Set the default interpolation value for points that are not found in the
156  /// mesh.
157  virtual void SetDefaultInterpolationValue(double interp_value_)
158  {
159  default_interp_value = interp_value_;
160  }
161 
162  /** Cleans up memory allocated internally by gslib.
163  Note that in parallel, this must be called before MPI_Finalize(), as it
164  calls MPI_Comm_free() for internal gslib communicators. */
165  virtual void FreeData();
166 
167  /// Return code for each point searched by FindPoints: inside element (0), on
168  /// element boundary (1), or not found (2).
169  virtual const Array<unsigned int> &GetCode() const { return gsl_code; }
170  /// Return element number for each point found by FindPoints.
171  virtual const Array<unsigned int> &GetElem() const { return gsl_mfem_elem; }
172  /// Return MPI rank on which each point was found by FindPoints.
173  virtual const Array<unsigned int> &GetProc() const { return gsl_proc; }
174  /// Return reference coordinates for each point found by FindPoints.
175  virtual const Vector &GetReferencePosition() const { return gsl_mfem_ref; }
176  /// Return distance between the sought and the found point in physical space,
177  /// for each point found by FindPoints.
178  virtual const Vector &GetDist() const { return gsl_dist; }
179 
180  /// Return element number for each point found by FindPoints corresponding to
181  /// GSLIB mesh. gsl_mfem_elem != gsl_elem for mesh with simplices.
182  virtual const Array<unsigned int> &GetGSLIBElem() const { return gsl_elem; }
183  /// Return reference coordinates in [-1,1] (internal range in GSLIB) for each
184  /// point found by FindPoints.
185  virtual const Vector &GetGSLIBReferencePosition() const { return gsl_ref; }
186 };
187 
188 /** \brief OversetFindPointsGSLIB enables use of findpts for arbitrary number of
189  overlapping grids.
190  The parameters in this class are the same as FindPointsGSLIB with the
191  difference of additional inputs required to account for more than 1 mesh. */
193 {
194 protected:
195  bool overset;
196  unsigned int u_meshid;
197  Vector distfint; // Used to store nodal vals of grid func. passed to findpts
198 
199 public:
201  overset(true) { }
202 
203 #ifdef MFEM_USE_MPI
204  OversetFindPointsGSLIB(MPI_Comm comm_) : FindPointsGSLIB(comm_),
205  overset(true) { }
206 #endif
207 
208  /** Initializes the internal mesh in gslib, by sending the positions of the
209  Gauss-Lobatto nodes of the input Mesh object @a m.
210  Note: not tested with periodic meshes (L2).
211  Note: the input mesh @a m must have Nodes set.
212 
213  @param[in] m Input mesh.
214  @param[in] meshid A unique # for each overlapping mesh. This id is
215  used to make sure that points being searched are not
216  looked for in the mesh that they belong to.
217  @param[in] gfmax (Optional) GridFunction in H1 that is used as a
218  discriminator when one point is located in multiple
219  meshes. The mesh that maximizes gfmax is chosen.
220  For example, using the distance field based on the
221  overlapping boundaries is helpful for convergence
222  during Schwarz iterations.
223  @param[in] bb_t (Optional) Relative size of bounding box around
224  each element.
225  @param[in] newt_tol (Optional) Newton tolerance for the gslib
226  search methods.
227  @param[in] npt_max (Optional) Number of points for simultaneous
228  iteration. This alters performance and
229  memory footprint.*/
230  void Setup(Mesh &m, const int meshid, GridFunction *gfmax = NULL,
231  const double bb_t = 0.1, const double newt_tol = 1.0e-12,
232  const int npt_max = 256);
233 
234  /** Searches positions given in physical space by @a point_pos. All output
235  Arrays and Vectors are expected to have the correct size.
236 
237  @param[in] point_pos Positions to be found. Must by ordered by nodes
238  (XXX...,YYY...,ZZZ).
239  @param[in] point_id Index of the mesh that the point belongs to
240  (corresponding to @a meshid in Setup). */
241  void FindPoints(const Vector &point_pos, Array<unsigned int> &point_id);
242 
243  /** Search positions and interpolate */
244  void Interpolate(const Vector &point_pos, Array<unsigned int> &point_id,
245  const GridFunction &field_in, Vector &field_out);
247 };
248 
249 } // namespace mfem
250 
251 #endif // MFEM_USE_GSLIB
252 
253 #endif // MFEM_GSLIB
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:582
Array< unsigned int > gsl_elem
Definition: gslib.hpp:60
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:178
virtual ~FindPointsGSLIB()
Definition: gslib.cpp:54
Array< unsigned int > gsl_mfem_elem
Definition: gslib.hpp:60
virtual const Vector & GetReferencePosition() const
Return reference coordinates for each point found by FindPoints.
Definition: gslib.hpp:175
OversetFindPointsGSLIB(MPI_Comm comm_)
Definition: gslib.hpp:204
Array< unsigned int > gsl_proc
Definition: gslib.hpp:60
virtual void GetNodeValues(const GridFunction &gf_in, Vector &node_vals)
Get GridFunction from MFEM format to GSLIB format.
Definition: gslib.cpp:233
virtual void SetL2AvgType(AvgType avgtype_)
Definition: gslib.hpp:153
virtual void MapRefPosAndElemIndices()
Definition: gslib.cpp:492
struct gslib::findpts_data_3 * fdata3D
Definition: gslib.hpp:56
OversetFindPointsGSLIB enables use of findpts for arbitrary number of overlapping grids...
Definition: gslib.hpp:192
virtual void FreeData()
Definition: gslib.cpp:212
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:75
double default_interp_value
Definition: gslib.hpp:63
struct gslib::crystal * cr
Definition: gslib.hpp:57
void FindPoints(const Vector &point_pos, Array< unsigned int > &point_id)
Definition: gslib.cpp:919
Array< unsigned int > gsl_code
Definition: gslib.hpp:60
virtual void GetSimplexNodalCoordinates()
Definition: gslib.cpp:328
void FindPoints(const Vector &point_pos)
Definition: gslib.cpp:129
virtual void InterpolateGeneral(const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:699
virtual void GetQuadHexNodalCoordinates()
Definition: gslib.cpp:290
virtual const Array< unsigned int > & GetProc() const
Return MPI rank on which each point was found by FindPoints.
Definition: gslib.hpp:173
struct gslib::findpts_data_2 * fdata2D
Definition: gslib.hpp:55
virtual void SetDefaultInterpolationValue(double interp_value_)
Definition: gslib.hpp:157
virtual const Array< unsigned int > & GetElem() const
Return element number for each point found by FindPoints.
Definition: gslib.hpp:171
virtual const Array< unsigned int > & GetCode() const
Definition: gslib.hpp:169
struct gslib::comm * gsl_comm
Definition: gslib.hpp:58
void Interpolate(const Vector &point_pos, Array< unsigned int > &point_id, const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:986
IntegrationRule * ir_simplex
Definition: gslib.hpp:54
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:182
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:847
virtual const Vector & GetGSLIBReferencePosition() const
Definition: gslib.hpp:185
virtual void InterpolateH1(const GridFunction &field_in, Vector &field_out)
Use GSLIB for communication and interpolation.
Definition: gslib.cpp:657