MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
pfindpts.cpp
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 // ----------------------------------------------------------------------
13 // Parallel Find Points Miniapp: Evaluate grid function in physical space
14 // ----------------------------------------------------------------------
15 //
16 // This miniapp demonstrates the interpolation of a high-order grid function on
17 // a set of points in physical-space. The miniapp is based on GSLIB-FindPoints,
18 // which provides two key functionalities. First, for a given set of points in
19 // the physical-space, it determines the computational coordinates (element
20 // number, reference-space coordinates inside the element, and processor number
21 // [in parallel]) for each point. Second, based on computational coordinates, it
22 // interpolates a grid function in the given points. Inside GSLIB, computation
23 // of the coordinates requires use of a Hash Table to identify the candidate
24 // processor and element for each point, followed by the Newton's method to
25 // determine the reference-space coordinates inside the candidate element.
26 //
27 // Compile with: make pfindpts
28 //
29 // Sample runs:
30 // mpirun -np 2 pfindpts -m ../../data/rt-2d-p4-tri.mesh -o 4
31 // mpirun -np 2 pfindpts -m ../../data/inline-tri.mesh -o 3
32 // mpirun -np 2 pfindpts -m ../../data/inline-quad.mesh -o 3
33 // mpirun -np 2 pfindpts -m ../../data/inline-quad.mesh -o 3 -hr
34 // mpirun -np 2 pfindpts -m ../../data/inline-tet.mesh -o 3
35 // mpirun -np 2 pfindpts -m ../../data/inline-hex.mesh -o 3
36 // mpirun -np 2 pfindpts -m ../../data/inline-wedge.mesh -o 3
37 // mpirun -np 2 pfindpts -m ../../data/amr-quad.mesh -o 2
38 // mpirun -np 2 pfindpts -m ../../data/rt-2d-q3.mesh -o 3 -mo 4 -ft 2
39 // mpirun -np 2 pfindpts -m ../../data/inline-quad.mesh -ft 1 -no-vis -sr0
40 
41 #include "mfem.hpp"
42 
43 using namespace mfem;
44 using namespace std;
45 
46 // Scalar function to project
47 double field_func(const Vector &x)
48 {
49  const int dim = x.Size();
50  double res = 0.0;
51  for (int d = 0; d < dim; d++) { res += x(d) * x(d); }
52  return res;
53 }
54 
55 void F_exact(const Vector &p, Vector &F)
56 {
57  F(0) = field_func(p);
58  for (int i = 1; i < F.Size(); i++) { F(i) = (i+1)*F(0); }
59 }
60 
61 int main (int argc, char *argv[])
62 {
63  // Initialize MPI and HYPRE.
64  Mpi::Init(argc, argv);
65  int num_procs = Mpi::WorldSize();
66  int myid = Mpi::WorldRank();
67  Hypre::Init();
68 
69  // Set the method's default parameters.
70  const char *mesh_file = "../../data/rt-2d-q3.mesh";
71  int order = 3;
72  int mesh_poly_deg = 3;
73  int rs_levels = 0;
74  int rp_levels = 0;
75  bool visualization = true;
76  int fieldtype = 0;
77  int ncomp = 1;
78  bool search_on_rank_0 = false;
79  bool hrefinement = false;
80 
81  // Parse command-line options.
82  OptionsParser args(argc, argv);
83  args.AddOption(&mesh_file, "-m", "--mesh",
84  "Mesh file to use.");
85  args.AddOption(&order, "-o", "--order",
86  "Finite element order (polynomial degree).");
87  args.AddOption(&mesh_poly_deg, "-mo", "--mesh-order",
88  "Polynomial degree of mesh finite element space.");
89  args.AddOption(&rs_levels, "-rs", "--refine-serial",
90  "Number of times to refine the mesh uniformly in serial.");
91  args.AddOption(&rp_levels, "-rp", "--refine-parallel",
92  "Number of times to refine the mesh uniformly in parallel.");
93  args.AddOption(&fieldtype, "-ft", "--field-type",
94  "Field type: 0 - H1, 1 - L2, 2 - H(div), 3 - H(curl).");
95  args.AddOption(&ncomp, "-nc", "--ncomp",
96  "Number of components for H1 or L2 GridFunctions");
97  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
98  "--no-visualization",
99  "Enable or disable GLVis visualization.");
100  args.AddOption(&search_on_rank_0, "-sr0", "--search-on-r0", "-no-sr0",
101  "--no-search-on-r0",
102  "Enable search only on rank 0 (disable to search points on all tasks).");
103  args.AddOption(&hrefinement, "-hr", "--h-refinement", "-no-hr",
104  "--no-h-refinement",
105  "Do random h refinements to mesh.");
106 
107  args.Parse();
108  if (!args.Good())
109  {
110  args.PrintUsage(cout);
111  return 1;
112  }
113  if (myid == 0) { args.PrintOptions(cout); }
114 
115  // Initialize and refine the starting mesh.
116  Mesh *mesh = new Mesh(mesh_file, 1, 1, false);
117  for (int lev = 0; lev < rs_levels; lev++) { mesh->UniformRefinement(); }
118  const int dim = mesh->Dimension();
119  if (myid == 0)
120  {
121  cout << "Mesh curvature of the original mesh: ";
122  if (mesh->GetNodes()) { cout << mesh->GetNodes()->OwnFEC()->Name(); }
123  else { cout << "(NONE)"; }
124  cout << endl;
125  }
126 
127  // Mesh bounding box (for the full serial mesh).
128  Vector pos_min, pos_max;
129  MFEM_VERIFY(mesh_poly_deg > 0, "The order of the mesh must be positive.");
130  mesh->GetBoundingBox(pos_min, pos_max, mesh_poly_deg);
131  if (myid == 0)
132  {
133  cout << "--- Generating equidistant point for:\n"
134  << "x in [" << pos_min(0) << ", " << pos_max(0) << "]\n"
135  << "y in [" << pos_min(1) << ", " << pos_max(1) << "]\n";
136  if (dim == 3)
137  {
138  cout << "z in [" << pos_min(2) << ", " << pos_max(2) << "]\n";
139  }
140  }
141 
142  // Distribute the mesh.
143  if (hrefinement) { mesh->EnsureNCMesh(); }
144  ParMesh pmesh(MPI_COMM_WORLD, *mesh);
145  delete mesh;
146  for (int lev = 0; lev < rp_levels; lev++) { pmesh.UniformRefinement(); }
147 
148  // Random h-refinements to mesh
149  if (hrefinement) { pmesh.RandomRefinement(0.5); }
150 
151  // Curve the mesh based on the chosen polynomial degree.
152  H1_FECollection fecm(mesh_poly_deg, dim);
153  ParFiniteElementSpace pfespace(&pmesh, &fecm, dim);
154  pmesh.SetNodalFESpace(&pfespace);
155  if (myid == 0)
156  {
157  cout << "Mesh curvature of the curved mesh: " << fecm.Name() << endl;
158  }
159 
160  MFEM_VERIFY(ncomp > 0, "Invalid number of components.");
161  int vec_dim = ncomp;
162  FiniteElementCollection *fec = NULL;
163  if (fieldtype == 0)
164  {
165  fec = new H1_FECollection(order, dim);
166  if (myid == 0) { cout << "H1-GridFunction\n"; }
167  }
168  else if (fieldtype == 1)
169  {
170  fec = new L2_FECollection(order, dim);
171  if (myid == 0) { cout << "L2-GridFunction\n"; }
172  }
173  else if (fieldtype == 2)
174  {
175  fec = new RT_FECollection(order, dim);
176  ncomp = 1;
177  vec_dim = dim;
178  if (myid == 0) { cout << "H(div)-GridFunction\n"; }
179  }
180  else if (fieldtype == 3)
181  {
182  fec = new ND_FECollection(order, dim);
183  ncomp = 1;
184  vec_dim = dim;
185  if (myid == 0) { cout << "H(curl)-GridFunction\n"; }
186  }
187  else
188  {
189  if (myid == 0) { MFEM_ABORT("Invalid FECollection type."); }
190  }
191  ParFiniteElementSpace sc_fes(&pmesh, fec, ncomp);
192  ParGridFunction field_vals(&sc_fes);
193 
194  // Project the GridFunction using VectorFunctionCoefficient.
196  field_vals.ProjectCoefficient(F);
197 
198  // Display the mesh and the field through glvis.
199  if (visualization)
200  {
201  char vishost[] = "localhost";
202  int visport = 19916;
203  socketstream sout;
204  sout.open(vishost, visport);
205  if (!sout)
206  {
207  if (myid == 0)
208  {
209  cout << "Unable to connect to GLVis server at "
210  << vishost << ':' << visport << endl;
211  }
212  }
213  else
214  {
215  sout << "parallel " << num_procs << " " << myid << "\n";
216  sout.precision(8);
217  sout << "solution\n" << pmesh << field_vals;
218  if (dim == 2) { sout << "keys RmjA*****\n"; }
219  if (dim == 3) { sout << "keys mA\n"; }
220  sout << flush;
221  }
222  }
223 
224  // Generate equidistant points in physical coordinates over the whole mesh.
225  // Note that some points might be outside, if the mesh is not a box. Note
226  // also that all tasks search the same points (not mandatory).
227  const int pts_cnt_1D = 10;
228  int pts_cnt = pow(pts_cnt_1D, dim);
229  Vector vxyz(pts_cnt * dim);
230  if (dim == 2)
231  {
233  const IntegrationRule &ir = el.GetNodes();
234  for (int i = 0; i < ir.GetNPoints(); i++)
235  {
236  const IntegrationPoint &ip = ir.IntPoint(i);
237  vxyz(i) = pos_min(0) + ip.x * (pos_max(0)-pos_min(0));
238  vxyz(pts_cnt + i) = pos_min(1) + ip.y * (pos_max(1)-pos_min(1));
239  }
240  }
241  else
242  {
244  const IntegrationRule &ir = el.GetNodes();
245  for (int i = 0; i < ir.GetNPoints(); i++)
246  {
247  const IntegrationPoint &ip = ir.IntPoint(i);
248  vxyz(i) = pos_min(0) + ip.x * (pos_max(0)-pos_min(0));
249  vxyz(pts_cnt + i) = pos_min(1) + ip.y * (pos_max(1)-pos_min(1));
250  vxyz(2*pts_cnt + i) = pos_min(2) + ip.z * (pos_max(2)-pos_min(2));
251  }
252  }
253 
254  if ( (myid != 0) && (search_on_rank_0) )
255  {
256  pts_cnt = 0;
257  vxyz.Destroy();
258  }
259 
260  // Find and Interpolate FE function values on the desired points.
261  Vector interp_vals(pts_cnt*vec_dim);
262  FindPointsGSLIB finder(MPI_COMM_WORLD);
263  finder.Setup(pmesh);
264  finder.Interpolate(vxyz, field_vals, interp_vals);
265  Array<unsigned int> code_out = finder.GetCode();
266  Array<unsigned int> task_id_out = finder.GetProc();
267  Vector dist_p_out = finder.GetDist();
268 
269  // Print the results for task 0 since either 1) all tasks have the
270  // same set of points or 2) only task 0 has any points.
271  if (myid == 0 )
272  {
273  int face_pts = 0, not_found = 0, found_loc = 0, found_away = 0;
274  double max_err = 0.0, max_dist = 0.0;
275  Vector pos(dim);
276  int npt = 0;
277  for (int j = 0; j < vec_dim; j++)
278  {
279  for (int i = 0; i < pts_cnt; i++)
280  {
281  if (j == 0)
282  {
283  (task_id_out[i] == (unsigned)myid) ? found_loc++ : found_away++;
284  }
285 
286  if (code_out[i] < 2)
287  {
288  for (int d = 0; d < dim; d++) { pos(d) = vxyz(d * pts_cnt + i); }
289  Vector exact_val(vec_dim);
290  F_exact(pos, exact_val);
291  max_err = std::max(max_err, fabs(exact_val(j) - interp_vals(npt)));
292  max_dist = std::max(max_dist, dist_p_out(i));
293  if (code_out[i] == 1 && j == 0) { face_pts++; }
294  }
295  else { if (j == 0) { not_found++; } }
296  npt++;
297  }
298  }
299 
300  cout << setprecision(16)
301  << "Searched unique points: " << pts_cnt
302  << "\nFound on local mesh: " << found_loc
303  << "\nFound on other tasks: " << found_away
304  << "\nMax interp error: " << max_err
305  << "\nMax dist (of found): " << max_dist
306  << "\nPoints not found: " << not_found
307  << "\nPoints on faces: " << face_pts << endl;
308  }
309 
310  // Free the internal gslib data.
311  finder.FreeData();
312 
313  delete fec;
314 
315  return 0;
316 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
const char vishost[]
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:584
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(), hypre will be finalized automatically at program exit.
Definition: hypre.hpp:67
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
virtual const Vector & GetDist() const
Definition: gslib.hpp:178
Arbitrary order L2 elements in 3D on a cube.
Definition: fe_l2.hpp:62
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
Definition: mesh.cpp:129
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:525
static int WorldSize()
Return the size of MPI_COMM_WORLD.
double field_func(const Vector &x)
Definition: findpts.cpp:103
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
virtual void FreeData()
Definition: gslib.cpp:214
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition: mesh.cpp:9394
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9737
void F_exact(const Vector &, Vector &)
Definition: ex4.cpp:268
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
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition: fe_base.hpp:390
void RandomRefinement(double prob, bool aniso=false, int nonconforming=-1, int nc_limit=0)
Refine each element with given probability. Uses GeneralRefinement.
Definition: mesh.cpp:9418
static void Init()
Singleton creation with Mpi::Init();.
const int visport
virtual const char * Name() const
Definition: fe_coll.hpp:237
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:5273
int Dimension() const
Definition: mesh.hpp:999
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:337
FDualNumber< tbase > pow(const FDualNumber< tbase > &a, const FDualNumber< tbase > &b)
pow([dual number],[dual number])
Definition: fdual.hpp:543
A general vector function coefficient.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
virtual const Array< unsigned int > & GetProc() const
Return MPI rank on which each point was found by FindPoints.
Definition: gslib.hpp:173
Class for integration point with weight.
Definition: intrules.hpp:25
int dim
Definition: ex24.cpp:53
virtual const Array< unsigned int > & GetCode() const
Definition: gslib.hpp:169
Arbitrary order L2 elements in 2D on a square.
Definition: fe_l2.hpp:39
void Destroy()
Destroy a vector.
Definition: vector.hpp:598
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:324
int open(const char hostname[], int port)
Open the socket stream on &#39;port&#39; at &#39;hostname&#39;.
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:411
FindPointsGSLIB can robustly evaluate a GridFunction on an arbitrary collection of points...
Definition: gslib.hpp:47
Vector data type.
Definition: vector.hpp:60
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7841
Nodes: x_i = i/(n-1), i=0,...,n-1.
Definition: fe_base.hpp:34
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Class for parallel meshes.
Definition: pmesh.hpp:32
int main()
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:284
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150