MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
gslib.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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
20namespace gslib
21{
22struct comm;
23struct crystal;
24struct hash_data_3;
25struct hash_data_2;
26struct gs_data;
27}
28
29namespace mfem
30{
31
32/** \brief FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary
33 * collection of points.
34 *
35 * There are three key functions in FindPointsGSLIB:
36 *
37 * 1. Setup - constructs the internal data structures of gslib. See \ref Setup.
38 *
39 * 2. FindPoints - for any given arbitrary set of points in physical space,
40 * gslib finds the element number, MPI rank, and the reference space
41 * coordinates inside the element that each point is located in. gslib also
42 * returns a code that indicates whether the point was found inside an
43 * element, on element border, or not found in the domain.
44 * For points returned as found on `element border`, the point is either
45 * on an element edge/face or near the domain boundary, and gslib also
46 * returns a distance to the border. Points near (but outside) the domain
47 * boundary must then be marked as not found using the distance returned
48 * by gslib. See \ref FindPoints.
49 *
50 * 3. Interpolate - Interpolates any grid function at the points found using 2.
51 * For functions in L2 finite element space, use \ref SetL2AvgType to
52 * specify how to interpolate values at points located at element boundaries
53 * where the function might be multi-valued. See \ref Interpolate.
54 *
55 * FindPointsGSLIB also provides interface to use these functions through a
56 * single call.
57 *
58 * For custom interpolation (e.g., evaluating strain rate tensor), we provide
59 * functions that use gslib to send element index and corresponding
60 * reference-space coordinates for each point to the mpi rank that the element
61 * is located on. Then, custom interpolation can be defined locally by the user
62 * before sending the values back to mpi ranks where the query originated from.
63 * See \ref DistributePointInfoToOwningMPIRanks and
64 * \ref DistributeInterpolatedValues.
65 */
67{
68public:
69 enum AvgType {NONE, ARITHMETIC, HARMONIC}; // Average type for L2 functions
70
71protected:
73 Array<Mesh *> mesh_split; // Meshes used to split simplices.
74 // IntegrationRules for simplex->Quad/Hex and to project to p_max in-case of
75 // p-refinement.
77 Array<FiniteElementSpace *> fes_rst_map; //FESpaces to map Quad/Hex->Simplex
78 Array<GridFunction *> gf_rst_map; // GridFunctions to map Quad/Hex->Simplex
80 void *fdataD;
81 struct gslib::crystal *cr; // gslib's internal data
82 struct gslib::comm *gsl_comm; // gslib's internal data
83 int dim, points_cnt; // mesh dimension and number of points
86 Array<unsigned int> recv_proc, recv_index; // data for custom interpolation
87 bool setupflag; // flag to indicate if gslib data has been setup
88 double default_interp_value; // used for points that are not found in the mesh
89 AvgType avgtype; // average type used for L2 functions
92 int NE_split_total; // total number of elements after mesh splitting
93 int mesh_points_cnt; // number of mesh nodes
94 // Tolerance to ignore points found beyond the mesh boundary.
95 // i.e. if ||x*-x(r)||_2^2 > bdr_tol, we mark point as not found.
96 double bdr_tol;
97 // Use CPU functions for Mesh/GridFunction on device for gslib1.0.7
98 bool gpu_to_cpu_fallback = false;
99
100 // Device specific data used for FindPoints
101 struct
102 {
103 bool setup_device = false;
104 bool find_device = false;
106 double newt_tol; // Tolerance specified during setup for Newton solve
107 struct gslib::crystal *cr;
108 struct gslib::hash_data_3 *hash3;
109 struct gslib::hash_data_2 *hash2;
114
115 /// Use GSLIB for communication and interpolation
116 virtual void InterpolateH1(const GridFunction &field_in, Vector &field_out);
117 /// Uses GSLIB Crystal Router for communication followed by MFEM's
118 /// interpolation functions
119 virtual void InterpolateGeneral(const GridFunction &field_in,
120 Vector &field_out);
121
122 /// Since GSLIB is designed to work with quads/hexes, we split every
123 /// triangle/tet/prism/pyramid element into quads/hexes.
124 virtual void SetupSplitMeshes();
125
126 /// Setup integration points that will be used to interpolate the nodal
127 /// location at points expected by GSLIB.
129 IntegrationRule *irule,
130 int order);
131
132 /// Get GridFunction value at the points expected by GSLIB.
133 virtual void GetNodalValues(const GridFunction *gf_in, Vector &node_vals);
134
135 /// Map {r,s,t} coordinates from [-1,1] to [0,1] for MFEM. For simplices,
136 /// find the original element number (that was split into micro quads/hexes)
137 /// during the setup phase.
138 virtual void MapRefPosAndElemIndices();
139
140 // Device functions
141 // FindPoints locally on device for 3D.
142 void FindPointsLocal3(const Vector &point_pos, int point_pos_ordering,
143 Array<unsigned int> &gsl_code_dev_l,
144 Array<unsigned int> &gsl_elem_dev_l, Vector &gsl_ref_l,
145 Vector &gsl_dist_l, int npt);
146
147 // FindPoints locally on device for 2D.
148 void FindPointsLocal2(const Vector &point_pos, int point_pos_ordering,
149 Array<unsigned int> &gsl_code_dev_l,
150 Array<unsigned int> &gsl_elem_dev_l, Vector &gsl_ref_l,
151 Vector &gsl_dist_l, int npt);
152
153 // Interpolate on device for 3D.
154 void InterpolateLocal3(const Vector &field_in,
155 Array<int> &gsl_elem_dev_l,
156 Vector &gsl_ref_l,
157 Vector &field_out,
158 int npt, int ncomp,
159 int nel, int dof1dsol);
160 // Interpolate on device for 2D.
161 void InterpolateLocal2(const Vector &field_in,
162 Array<int> &gsl_elem_dev_l,
163 Vector &gsl_ref_l,
164 Vector &field_out,
165 int npt, int ncomp,
166 int nel, int dof1dsol);
167
168 // Prepare data for device functions.
169 void SetupDevice();
170
171 /** Searches positions given in physical space by @a point_pos.
172 These positions can be ordered byNodes: (XXX...,YYY...,ZZZ) or
173 byVDim: (XYZ,XYZ,....XYZ) specified by @a point_pos_ordering. */
174 void FindPointsOnDevice(const Vector &point_pos,
175 int point_pos_ordering = Ordering::byNODES);
176
177 /** Interpolation of field values at prescribed reference space positions.
178 @param[in] field_in_evec E-vector of grid function to be interpolated.
179 Assumed ordering is NDOFSxVDIMxNEL
180 @param[in] nel Number of elements in the mesh.
181 @param[in] ncomp Number of components in the field.
182 @param[in] dof1dsol Number of degrees of freedom in each reference
183 space direction.
184 @param[in] ordering Ordering of the out field values: byNodes/byVDIM
185
186 @param[out] field_out Interpolated values. For points that are not found
187 the value is set to #default_interp_value. */
188 void InterpolateOnDevice(const Vector &field_in_evec, Vector &field_out,
189 const int nel, const int ncomp,
190 const int dof1dsol, const int ordering);
191public:
193
194#ifdef MFEM_USE_MPI
195 FindPointsGSLIB(MPI_Comm comm_);
196#endif
197
198 virtual ~FindPointsGSLIB();
199
200 /** Initializes the internal mesh in gslib, by sending the positions of the
201 Gauss-Lobatto nodes of the input Mesh object \p m.
202 Note: not tested with periodic (L2).
203 Note: the input mesh \p m must have Nodes set.
204
205 @param[in] m Input mesh.
206 @param[in] bb_t (Optional) Relative size of bounding box around
207 each element.
208 @param[in] newt_tol (Optional) Newton tolerance for the gslib
209 search methods.
210 @param[in] npt_max (Optional) Number of points for simultaneous
211 iteration. This alters performance and
212 memory footprint.*/
213 void Setup(Mesh &m, const double bb_t = 0.1,
214 const double newt_tol = 1.0e-12,
215 const int npt_max = 256);
216 /** Searches positions given in physical space by \p point_pos.
217 These positions can be ordered byNodes: (XXX...,YYY...,ZZZ) or
218 byVDim: (XYZ,XYZ,....XYZ) specified by \p point_pos_ordering.
219 This function populates the following member variables:
220 #gsl_code Return codes for each point: inside element (0),
221 element boundary (1), not found (2).
222 #gsl_proc MPI proc ids where the points were found.
223 #gsl_elem Element ids where the points were found.
224 Defaults to 0 for points that were not found.
225 #gsl_mfem_elem Element ids corresponding to MFEM-mesh where the points
226 were found. #gsl_mfem_elem != #gsl_elem for simplices
227 Defaults to 0 for points that were not found.
228 #gsl_ref Reference coordinates of the found point.
229 Ordered by vdim (XYZ,XYZ,XYZ...). Defaults to -1 for
230 points that were not found. Note: the gslib reference
231 frame is [-1,1].
232 #gsl_mfem_ref Reference coordinates #gsl_ref mapped to [0,1].
233 Defaults to 0 for points that were not found.
234 #gsl_dist Distance between the sought and the found point
235 in physical space. */
236 void FindPoints(const Vector &point_pos,
237 int point_pos_ordering = Ordering::byNODES);
238 /// Setup FindPoints and search positions
239 void FindPoints(Mesh &m, const Vector &point_pos,
240 int point_pos_ordering = Ordering::byNODES,
241 const double bb_t = 0.1, const double newt_tol = 1.0e-12,
242 const int npt_max = 256);
243
244 /** Interpolation of field values at prescribed reference space positions.
245 @param[in] field_in Function values that will be interpolated on the
246 reference positions. Note: it is assumed that
247 \p field_in is in H1 and in the same space as the
248 mesh that was given to Setup().
249 @param[out] field_out Interpolated values. For points that are not found
250 the value is set to #default_interp_value. */
251 virtual void Interpolate(const GridFunction &field_in, Vector &field_out);
252 /** Search positions and interpolate. The ordering (byNODES or byVDIM) of
253 the output values in \p field_out corresponds to the ordering used
254 in the input GridFunction \p field_in. */
255 void Interpolate(const Vector &point_pos, const GridFunction &field_in,
256 Vector &field_out,
257 int point_pos_ordering = Ordering::byNODES);
258 /** Setup FindPoints, search positions and interpolate. The ordering (byNODES
259 or byVDIM) of the output values in \p field_out corresponds to the
260 ordering used in the input GridFunction \p field_in. */
261 void Interpolate(Mesh &m, const Vector &point_pos,
262 const GridFunction &field_in, Vector &field_out,
263 int point_pos_ordering = Ordering::byNODES);
264
265 /// Average type to be used for L2 functions in-case a point is located at
266 /// an element boundary where the function might be multi-valued.
267 virtual void SetL2AvgType(AvgType avgtype_) { avgtype = avgtype_; }
268
269 /// Set the default interpolation value for points that are not found in the
270 /// mesh.
271 virtual void SetDefaultInterpolationValue(double interp_value_)
272 {
273 default_interp_value = interp_value_;
274 }
275
276 /// Set the tolerance for detecting points outside the 'curvilinear' boundary
277 /// that gslib may return as found on the boundary. Points found on boundary
278 /// with distance greater than @ bdr_tol are marked as not found.
280 {
281 bdr_tol = bdr_tol_;
282 }
283
284 /// Enable/Disable use of CPU functions for GPU data if the gslib version
285 /// is older.
286 virtual void SetGPUtoCPUFallback(bool mode) { gpu_to_cpu_fallback = mode; }
287
288 /** Cleans up memory allocated internally by gslib.
289 Note that in parallel, this must be called before MPI_Finalize(), as it
290 calls MPI_Comm_free() for internal gslib communicators. */
291 virtual void FreeData();
292
293 /// Return code for each point searched by FindPoints: inside element (0), on
294 /// element boundary (1), or not found (2).
295 virtual const Array<unsigned int> &GetCode() const { return gsl_code; }
296 /// Return element number for each point found by FindPoints.
297 virtual const Array<unsigned int> &GetElem() const { return gsl_mfem_elem; }
298 /// Return MPI rank on which each point was found by FindPoints.
299 virtual const Array<unsigned int> &GetProc() const { return gsl_proc; }
300 /// Return reference coordinates for each point found by FindPoints.
301 virtual const Vector &GetReferencePosition() const { return gsl_mfem_ref; }
302 /// Return distance between the sought and the found point in physical space,
303 /// for each point found by FindPoints.
304 virtual const Vector &GetDist() const { return gsl_dist; }
305
306 /// Return element number for each point found by FindPoints corresponding to
307 /// GSLIB mesh. gsl_mfem_elem != gsl_elem for mesh with simplices.
308 virtual const Array<unsigned int> &GetGSLIBElem() const { return gsl_elem; }
309 /// Return reference coordinates in [-1,1] (internal range in GSLIB) for each
310 /// point found by FindPoints.
311 virtual const Vector &GetGSLIBReferencePosition() const { return gsl_ref; }
312
313 /** @name Methods to support a custom interpolation procedure.
314 \brief The physical-space point that the user seeks to interpolate at
315 could be located inside an element on another mpi rank.
316 To enable a custom interpolation procedure (e.g., strain tensor computation)
317 we need a mechanism to first send element indices and reference-space
318 coordinates to the mpi-ranks where each point is found. Then the custom
319 interpolation can be done locally by the user before sending the
320 interpolated values back to the mpi-ranks that the query originated from.
321 Example usage looks something like this:
322
323 FindPoints() -> DistributePointInfoToOwningMPIRanks() -> Computation by
324 user -> DistributeInterpolatedValues().
325 */
326 ///@{
327 /// Distribute element indices in #gsl_mfem_elem, the reference coordinates
328 /// #gsl_mfem_ref, and the code #gsl_code to the corresponding mpi-rank
329 /// #gsl_proc for each point. The received information is provided locally
330 /// in \p recv_elem, \p recv_ref (ordered by vdim), and \p recv_code.
331 /// Note: The user can send empty Array/Vectors to the method as they are
332 /// appropriately sized and filled internally.
334 Array<unsigned int> &recv_elem, Vector &recv_ref,
335 Array<unsigned int> &recv_code);
336 /// Return interpolated values back to the mpi-ranks #recv_proc that had
337 /// sent the element indices and corresponding reference-space coordinates.
338 /// Specify \p vdim and \p ordering (by nodes or by vdim) based on how the
339 /// \p int_vals are structured. The received values are filled in
340 /// \p field_out consistent with the original ordering of the points that
341 /// were used in \ref FindPoints.
342 virtual void DistributeInterpolatedValues(const Vector &int_vals,
343 const int vdim,
344 const int ordering,
345 Vector &field_out) const;
346 ///@}
347
348 /// Return the axis-aligned bounding boxes (AABB) computed during \ref Setup.
349 /// The size of the returned vector is (nel x nverts x dim), where nel is the
350 /// number of elements (after splitting for simplcies), nverts is number of
351 /// vertices (4 in 2D, 8 in 3D), and dim is the spatial dimension.
353
354 /// Return the oriented bounding boxes (OBB) computed during \ref Setup.
355 /// Each OBB is represented using the inverse transformation (A^{-1}) and
356 /// its center (x_c), such that a point x is inside the OBB if:
357 /// -1 <= A^{-1}(x-x_c) <= 1.
358 /// The inverse transformation is returned in \p obbA, a DenseTensor of
359 /// size (dim x dim x nel), and the OBB centers are returned in \p obbC,
360 /// a vector of size (nel x dim). The vertices of the OBBs are returned in
361 /// \p obbV, a vector of size (nel x nverts x dim) .
362 void GetOrientedBoundingBoxes(DenseTensor &obbA, Vector &obbC, Vector &obbV);
363};
364
365/** \brief OversetFindPointsGSLIB enables use of findpts for arbitrary number of
366 overlapping grids.
367
368 The parameters in this class are the same as FindPointsGSLIB with the
369 difference of additional inputs required to account for more than 1 mesh. */
371{
372protected:
374 unsigned int u_meshid;
375 Vector distfint; // Used to store nodal vals of grid func. passed to findpts
376
377public:
380
381#ifdef MFEM_USE_MPI
382 OversetFindPointsGSLIB(MPI_Comm comm_) : FindPointsGSLIB(comm_),
383 overset(true) { }
384#endif
385
386 /** Initializes the internal mesh in gslib, by sending the positions of the
387 Gauss-Lobatto nodes of the input Mesh object \p m.
388 Note: not tested with periodic meshes (L2).
389 Note: the input mesh \p m must have Nodes set.
390
391 @param[in] m Input mesh.
392 @param[in] meshid A unique # for each overlapping mesh. This id is
393 used to make sure that points being searched are not
394 looked for in the mesh that they belong to.
395 @param[in] gfmax (Optional) GridFunction in H1 that is used as a
396 discriminator when one point is located in multiple
397 meshes. The mesh that maximizes gfmax is chosen.
398 For example, using the distance field based on the
399 overlapping boundaries is helpful for convergence
400 during Schwarz iterations.
401 @param[in] bb_t (Optional) Relative size of bounding box around
402 each element.
403 @param[in] newt_tol (Optional) Newton tolerance for the gslib
404 search methods.
405 @param[in] npt_max (Optional) Number of points for simultaneous
406 iteration. This alters performance and
407 memory footprint.*/
408 void Setup(Mesh &m, const int meshid, GridFunction *gfmax = NULL,
409 const double bb_t = 0.1, const double newt_tol = 1.0e-12,
410 const int npt_max = 256);
411
412 /** Searches positions given in physical space by \p point_pos. All output
413 Arrays and Vectors are expected to have the correct size.
414
415 @param[in] point_pos Positions to be found.
416 @param[in] point_id Index of the mesh that the point belongs
417 to (corresponding to \p meshid in Setup).
418 @param[in] point_pos_ordering Ordering of the points:
419 byNodes: (XXX...,YYY...,ZZZ) or
420 byVDim: (XYZ,XYZ,....XYZ) */
421 void FindPoints(const Vector &point_pos,
422 Array<unsigned int> &point_id,
423 int point_pos_ordering = Ordering::byNODES);
424
425 /** Search positions and interpolate */
426 void Interpolate(const Vector &point_pos, Array<unsigned int> &point_id,
427 const GridFunction &field_in, Vector &field_out,
428 int point_pos_ordering = Ordering::byNODES);
430};
431
432/** \brief Class for gather-scatter (gs) operations on Vectors based on
433 corresponding global identifiers.
434
435 This functionality is useful for gs-ops on DOF values across processor
436 boundary, where the global identifier would be the corresponding true DOF
437 index. Operations currently supported are min, max, sum, and multiplication.
438 Note: identifier 0 does not participate in the gather-scatter operation and
439 a given identifier can be included multiple times on a given rank.
440 For example, consider a vector, v:
441 - v = [0.3, 0.4, 0.25, 0.7] on rank1,
442 - v = [0.6, 0.1] on rank 2,
443 - v = [-0.2, 0.3, 0.7, 0.] on rank 3.
444
445 Consider a corresponding Array<int>, a:
446 - a = [1, 2, 3, 1] on rank 1,
447 - a = [3, 2] on rank 2,
448 - a = [1, 2, 0, 3] on rank 3.
449
450 A gather-scatter "minimum" operation, done as follows:
451 GSOPGSLIB gs = GSOPGSLIB(MPI_COMM_WORLD, a);
452 gs.GS(v, GSOp::MIN);
453 would return into v:
454 - v = [-0.2, 0.1, 0., -0.2] on rank 1,
455 - v = [0., 0.1] on rank 2,
456 - v = [-0.2, 0.1, 0.7, 0.] on rank 3,
457 where the values have been compared across all processors based on the
458 integer identifier. */
460{
461protected:
462 struct gslib::crystal *cr; // gslib's internal data
463 struct gslib::comm *gsl_comm; // gslib's internal data
464 struct gslib::gs_data *gsl_data = NULL;
466
467public:
469
470#ifdef MFEM_USE_MPI
471 GSOPGSLIB(MPI_Comm comm_, Array<long long> &ids);
472#endif
473
474 virtual ~GSOPGSLIB();
475
476 /// Supported operation types. See class description.
477 enum GSOp {ADD, MUL, MIN, MAX};
478
479 /// Update the identifiers used for the gather-scatter operator.
480 /// Same \p ids get grouped together and id == 0 does not participate.
481 /// See class description.
482 void UpdateIdentifiers(const Array<long long> &ids);
483
484 /// Gather-Scatter operation on senddata. Must match length of unique
485 /// identifiers used in the constructor. See class description.
486 void GS(Vector &senddata, GSOp op);
487};
488
489} // namespace mfem
490
491#endif // MFEM_USE_GSLIB
492
493#endif // MFEM_GSLIB
Rank 3 tensor (array of matrices)
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points.
Definition gslib.hpp:67
virtual void DistributePointInfoToOwningMPIRanks(Array< unsigned int > &recv_elem, Vector &recv_ref, Array< unsigned int > &recv_code)
Definition gslib.cpp:2084
struct mfem::FindPointsGSLIB::@24 DEV
virtual ~FindPointsGSLIB()
Definition gslib.cpp:120
void FindPointsOnDevice(const Vector &point_pos, int point_pos_ordering=Ordering::byNODES)
Definition gslib.cpp:489
Array< unsigned int > gsl_code
Definition gslib.hpp:84
Array< Mesh * > mesh_split
Definition gslib.hpp:73
virtual const Vector & GetDist() const
Definition gslib.hpp:304
void InterpolateLocal2(const Vector &field_in, Array< int > &gsl_elem_dev_l, Vector &gsl_ref_l, Vector &field_out, int npt, int ncomp, int nel, int dof1dsol)
void FindPoints(const Vector &point_pos, int point_pos_ordering=Ordering::byNODES)
Definition gslib.cpp:242
void GetAxisAlignedBoundingBoxes(Vector &aabb)
Definition gslib.cpp:2198
void FindPointsLocal3(const Vector &point_pos, int point_pos_ordering, Array< unsigned int > &gsl_code_dev_l, Array< unsigned int > &gsl_elem_dev_l, Vector &gsl_ref_l, Vector &gsl_dist_l, int npt)
virtual const Vector & GetReferencePosition() const
Return reference coordinates for each point found by FindPoints.
Definition gslib.hpp:301
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:163
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
Definition gslib.cpp:1721
virtual void InterpolateH1(const GridFunction &field_in, Vector &field_out)
Use GSLIB for communication and interpolation.
Definition gslib.cpp:1847
Array< FiniteElementSpace * > fes_rst_map
Definition gslib.hpp:77
virtual void SetGPUtoCPUFallback(bool mode)
Definition gslib.hpp:286
virtual void GetNodalValues(const GridFunction *gf_in, Vector &node_vals)
Get GridFunction value at the points expected by GSLIB.
Definition gslib.cpp:1434
struct gslib::hash_data_3 * hash3
Definition gslib.hpp:108
virtual const Array< unsigned int > & GetCode() const
Definition gslib.hpp:295
Array< int > split_element_index
Definition gslib.hpp:91
FiniteElementCollection * fec_map_lin
Definition gslib.hpp:79
virtual const Array< unsigned int > & GetElem() const
Return element number for each point found by FindPoints.
Definition gslib.hpp:297
virtual void DistributeInterpolatedValues(const Vector &int_vals, const int vdim, const int ordering, Vector &field_out) const
Definition gslib.cpp:2143
struct gslib::comm * gsl_comm
Definition gslib.hpp:82
void InterpolateLocal3(const Vector &field_in, Array< int > &gsl_elem_dev_l, Vector &gsl_ref_l, Vector &field_out, int npt, int ncomp, int nel, int dof1dsol)
Array< unsigned int > gsl_elem
Definition gslib.hpp:84
Array< unsigned int > recv_proc
Definition gslib.hpp:86
virtual const Array< unsigned int > & GetGSLIBElem() const
Definition gslib.hpp:308
Array< unsigned int > gsl_proc
Definition gslib.hpp:84
virtual void SetupIntegrationRuleForSplitMesh(Mesh *mesh, IntegrationRule *irule, int order)
Definition gslib.cpp:1386
virtual const Vector & GetGSLIBReferencePosition() const
Definition gslib.hpp:311
Array< unsigned int > gsl_mfem_elem
Definition gslib.hpp:84
Array< IntegrationRule * > ir_split
Definition gslib.hpp:76
Array< unsigned int > loc_hash_offset
Definition gslib.hpp:111
double default_interp_value
Definition gslib.hpp:88
virtual void SetupSplitMeshes()
Definition gslib.cpp:1159
virtual void MapRefPosAndElemIndices()
Definition gslib.cpp:1544
virtual void FreeData()
Definition gslib.cpp:1128
virtual void SetDefaultInterpolationValue(double interp_value_)
Definition gslib.hpp:271
virtual void InterpolateGeneral(const GridFunction &field_in, Vector &field_out)
Definition gslib.cpp:1923
Array< GridFunction * > gf_rst_map
Definition gslib.hpp:78
virtual void SetL2AvgType(AvgType avgtype_)
Definition gslib.hpp:267
virtual const Array< unsigned int > & GetProc() const
Return MPI rank on which each point was found by FindPoints.
Definition gslib.hpp:299
Array< int > split_element_map
Definition gslib.hpp:90
void GetOrientedBoundingBoxes(DenseTensor &obbA, Vector &obbC, Vector &obbV)
Definition gslib.cpp:2268
void InterpolateOnDevice(const Vector &field_in_evec, Vector &field_out, const int nel, const int ncomp, const int dof1dsol, const int ordering)
Definition gslib.cpp:856
struct gslib::hash_data_2 * hash2
Definition gslib.hpp:109
Array< unsigned int > recv_index
Definition gslib.hpp:86
void FindPointsLocal2(const Vector &point_pos, int point_pos_ordering, Array< unsigned int > &gsl_code_dev_l, Array< unsigned int > &gsl_elem_dev_l, Vector &gsl_ref_l, Vector &gsl_dist_l, int npt)
virtual void SetDistanceToleranceForPointsFoundOnBoundary(double bdr_tol_)
Definition gslib.hpp:279
struct gslib::crystal * cr
Definition gslib.hpp:81
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
Class for gather-scatter (gs) operations on Vectors based on corresponding global identifiers.
Definition gslib.hpp:460
struct gslib::comm * gsl_comm
Definition gslib.hpp:463
void GS(Vector &senddata, GSOp op)
Definition gslib.cpp:2621
GSOPGSLIB(Array< long long > &ids)
Definition gslib.cpp:2567
virtual ~GSOPGSLIB()
Definition gslib.cpp:2596
GSOp
Supported operation types. See class description.
Definition gslib.hpp:477
struct gslib::crystal * cr
Definition gslib.hpp:462
struct gslib::gs_data * gsl_data
Definition gslib.hpp:464
void UpdateIdentifiers(const Array< long long > &ids)
Definition gslib.cpp:2605
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Mesh data type.
Definition mesh.hpp:64
OversetFindPointsGSLIB enables use of findpts for arbitrary number of overlapping grids.
Definition gslib.hpp:371
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:2382
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:2557
void FindPoints(const Vector &point_pos, Array< unsigned int > &point_id, int point_pos_ordering=Ordering::byNODES)
Definition gslib.cpp:2477
OversetFindPointsGSLIB(MPI_Comm comm_)
Definition gslib.hpp:382
Vector data type.
Definition vector.hpp:82