MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
findpts.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 // 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 findpts
28 //
29 // Sample runs:
30 // findpts -m ../../data/rt-2d-p4-tri.mesh -o 4
31 // findpts -m ../../data/inline-tri.mesh -o 3
32 // findpts -m ../../data/inline-quad.mesh -o 3
33 // findpts -m ../../data/inline-quad.mesh -o 3 -hr -pr
34 // findpts -m ../../data/inline-tet.mesh -o 3
35 // findpts -m ../../data/inline-hex.mesh -o 3
36 // findpts -m ../../data/inline-wedge.mesh -o 3
37 // findpts -m ../../data/amr-quad.mesh -o 2
38 // findpts -m ../../data/rt-2d-q3.mesh -o 3 -mo 4 -ft 2
39 
40 #include "mfem.hpp"
41 
42 using namespace mfem;
43 using namespace std;
44 
45 // Experimental - required for visualizing functions on p-refined spaces.
46 GridFunction* ProlongToMaxOrder(const GridFunction *x, const int fieldtype)
47 {
48  const FiniteElementSpace *fespace = x->FESpace();
49  Mesh *mesh = fespace->GetMesh();
50  const FiniteElementCollection *fec = fespace->FEColl();
51 
52  // find the max order in the space
53  int max_order = 1;
54  for (int i = 0; i < mesh->GetNE(); i++)
55  {
56  max_order = std::max(fespace->GetElementOrder(i), max_order);
57  }
58 
59  // create a visualization space of max order for all elements
60  FiniteElementCollection *fecInt = NULL;
61  if (fieldtype == 0)
62  {
63  fecInt = new H1_FECollection(max_order, mesh->Dimension());
64  }
65  else if (fieldtype == 1)
66  {
67  fecInt = new L2_FECollection(max_order, mesh->Dimension());
68  }
69  FiniteElementSpace *spaceInt = new FiniteElementSpace(mesh, fecInt);
70 
72  DenseMatrix I;
73 
74  GridFunction *xInt = new GridFunction(spaceInt);
75 
76  // interpolate solution vector in the larger space
77  for (int i = 0; i < mesh->GetNE(); i++)
78  {
79  Geometry::Type geom = mesh->GetElementGeometry(i);
81 
82  Array<int> dofs;
83  fespace->GetElementDofs(i, dofs);
84  Vector elemvect, vectInt;
85  x->GetSubVector(dofs, elemvect);
86 
87  const auto *fe = fec->GetFE(geom, fespace->GetElementOrder(i));
88  const auto *feInt = fecInt->GetFE(geom, max_order);
89 
90  feInt->GetTransferMatrix(*fe, T, I);
91  spaceInt->GetElementDofs(i, dofs);
92  vectInt.SetSize(dofs.Size());
93 
94  I.Mult(elemvect, vectInt);
95  xInt->SetSubVector(dofs, vectInt);
96  }
97 
98  xInt->MakeOwner(fecInt);
99  return xInt;
100 }
101 
102 // Scalar function to project
103 double field_func(const Vector &x)
104 {
105  const int dim = x.Size();
106  double res = 0.0;
107  for (int d = 0; d < dim; d++) { res += x(d) * x(d); }
108  return res;
109 }
110 
111 void F_exact(const Vector &p, Vector &F)
112 {
113  F(0) = field_func(p);
114  for (int i = 1; i < F.Size(); i++) { F(i) = (i+1)*F(0); }
115 }
116 
117 int main (int argc, char *argv[])
118 {
119  // Set the method's default parameters.
120  const char *mesh_file = "../../data/rt-2d-q3.mesh";
121  int order = 3;
122  int mesh_poly_deg = 3;
123  int rs_levels = 0;
124  bool visualization = true;
125  int fieldtype = 0;
126  int ncomp = 1;
127  bool hrefinement = false;
128  bool prefinement = false;
129 
130  // Parse command-line options.
131  OptionsParser args(argc, argv);
132  args.AddOption(&mesh_file, "-m", "--mesh",
133  "Mesh file to use.");
134  args.AddOption(&order, "-o", "--order",
135  "Finite element order (polynomial degree).");
136  args.AddOption(&mesh_poly_deg, "-mo", "--mesh-order",
137  "Polynomial degree of mesh finite element space.");
138  args.AddOption(&rs_levels, "-rs", "--refine-serial",
139  "Number of times to refine the mesh uniformly in serial.");
140  args.AddOption(&fieldtype, "-ft", "--field-type",
141  "Field type: 0 - H1, 1 - L2, 2 - H(div), 3 - H(curl).");
142  args.AddOption(&ncomp, "-nc", "--ncomp",
143  "Number of components for H1 or L2 GridFunctions");
144  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
145  "--no-visualization",
146  "Enable or disable GLVis visualization.");
147  args.AddOption(&hrefinement, "-hr", "--h-refinement", "-no-hr",
148  "--no-h-refinement",
149  "Do random h refinements to mesh.");
150  args.AddOption(&prefinement, "-pr", "--p-refinement", "-no-pr",
151  "--no-p-refinement",
152  "Do random p refinements to solution field.");
153 
154  args.Parse();
155  if (!args.Good())
156  {
157  args.PrintUsage(cout);
158  return 1;
159  }
160  args.PrintOptions(cout);
161 
162  // Initialize and refine the starting mesh.
163  Mesh mesh(mesh_file, 1, 1, false);
164  for (int lev = 0; lev < rs_levels; lev++) { mesh.UniformRefinement(); }
165  const int dim = mesh.Dimension();
166  cout << "Mesh curvature of the original mesh: ";
167  if (mesh.GetNodes()) { cout << mesh.GetNodes()->OwnFEC()->Name(); }
168  else { cout << "(NONE)"; }
169  cout << endl;
170 
171  // Mesh bounding box.
172  Vector pos_min, pos_max;
173  MFEM_VERIFY(mesh_poly_deg > 0, "The order of the mesh must be positive.");
174  mesh.GetBoundingBox(pos_min, pos_max, mesh_poly_deg);
175  if (hrefinement || prefinement) { mesh.EnsureNCMesh(); }
176  cout << "--- Generating equidistant point for:\n"
177  << "x in [" << pos_min(0) << ", " << pos_max(0) << "]\n"
178  << "y in [" << pos_min(1) << ", " << pos_max(1) << "]\n";
179  if (dim == 3)
180  {
181  cout << "z in [" << pos_min(2) << ", " << pos_max(2) << "]\n";
182  }
183 
184  // Random h-refinements to mesh
185  if (hrefinement) { mesh.RandomRefinement(0.5); }
186 
187  // Curve the mesh based on the chosen polynomial degree.
188  H1_FECollection fecm(mesh_poly_deg, dim);
189  FiniteElementSpace fespace(&mesh, &fecm, dim);
190  mesh.SetNodalFESpace(&fespace);
191  cout << "Mesh curvature of the curved mesh: " << fecm.Name() << endl;
192 
193  MFEM_VERIFY(ncomp > 0, "Invalid number of components.");
194  int vec_dim = ncomp;
195  FiniteElementCollection *fec = NULL;
196  if (fieldtype == 0)
197  {
198  fec = new H1_FECollection(order, dim);
199  cout << "H1-GridFunction\n";
200  }
201  else if (fieldtype == 1)
202  {
203  fec = new L2_FECollection(order, dim);
204  cout << "L2-GridFunction\n";
205  }
206  else if (fieldtype == 2)
207  {
208  fec = new RT_FECollection(order, dim);
209  ncomp = 1;
210  vec_dim = dim;
211  cout << "H(div)-GridFunction\n";
212  }
213  else if (fieldtype == 3)
214  {
215  fec = new ND_FECollection(order, dim);
216  ncomp = 1;
217  vec_dim = dim;
218  cout << "H(curl)-GridFunction\n";
219  }
220  else
221  {
222  MFEM_ABORT("Invalid field type.");
223  }
224  FiniteElementSpace sc_fes(&mesh, fec, ncomp);
225  GridFunction field_vals(&sc_fes);
226 
227  // Random p-refinements to the solution field
228  if (prefinement)
229  {
230  for (int e = 0; e < mesh.GetNE(); e++)
231  {
232  if (rand() % 2 == 0)
233  {
234  int element_order = sc_fes.GetElementOrder(e);
235  sc_fes.SetElementOrder(e, element_order + 1);
236  }
237  }
238  sc_fes.Update(false);
239  field_vals.Update();
240  }
241 
242  // Project the GridFunction using VectorFunctionCoefficient.
244  field_vals.ProjectCoefficient(F);
245 
246  GridFunction *field_vals_pref = prefinement ?
247  ProlongToMaxOrder(&field_vals, fieldtype) :
248  &field_vals;
249 
250  // Display the mesh and the field through glvis.
251  if (visualization)
252  {
253  char vishost[] = "localhost";
254  int visport = 19916;
255  socketstream sout;
256  sout.open(vishost, visport);
257  if (!sout)
258  {
259  cout << "Unable to connect to GLVis server at "
260  << vishost << ':' << visport << endl;
261  }
262  else
263  {
264  sout.precision(8);
265  sout << "solution\n" << mesh << *field_vals_pref;
266  if (dim == 2) { sout << "keys RmjA*****\n"; }
267  if (dim == 3) { sout << "keys mA\n"; }
268  sout << flush;
269  }
270  }
271 
272  // Generate equidistant points in physical coordinates over the whole mesh.
273  // Note that some points might be outside, if the mesh is not a box. Note
274  // also that all tasks search the same points (not mandatory).
275  const int pts_cnt_1D = 25;
276  int pts_cnt = pow(pts_cnt_1D, dim);
277  Vector vxyz(pts_cnt * dim);
278  if (dim == 2)
279  {
281  const IntegrationRule &ir = el.GetNodes();
282  for (int i = 0; i < ir.GetNPoints(); i++)
283  {
284  const IntegrationPoint &ip = ir.IntPoint(i);
285  vxyz(i) = pos_min(0) + ip.x * (pos_max(0)-pos_min(0));
286  vxyz(pts_cnt + i) = pos_min(1) + ip.y * (pos_max(1)-pos_min(1));
287  }
288  }
289  else
290  {
292  const IntegrationRule &ir = el.GetNodes();
293  for (int i = 0; i < ir.GetNPoints(); i++)
294  {
295  const IntegrationPoint &ip = ir.IntPoint(i);
296  vxyz(i) = pos_min(0) + ip.x * (pos_max(0)-pos_min(0));
297  vxyz(pts_cnt + i) = pos_min(1) + ip.y * (pos_max(1)-pos_min(1));
298  vxyz(2*pts_cnt + i) = pos_min(2) + ip.z * (pos_max(2)-pos_min(2));
299  }
300  }
301 
302  // Find and Interpolate FE function values on the desired points.
303  Vector interp_vals(pts_cnt*vec_dim);
304  FindPointsGSLIB finder;
305  finder.Setup(mesh);
307  finder.Interpolate(vxyz, field_vals, interp_vals);
308  Array<unsigned int> code_out = finder.GetCode();
309  Vector dist_p_out = finder.GetDist();
310 
311  int face_pts = 0, not_found = 0, found = 0;
312  double max_err = 0.0, max_dist = 0.0;
313  Vector pos(dim);
314  int npt = 0;
315  for (int j = 0; j < vec_dim; j++)
316  {
317  for (int i = 0; i < pts_cnt; i++)
318  {
319  if (code_out[i] < 2)
320  {
321  if (j == 0) { found++; }
322  for (int d = 0; d < dim; d++) { pos(d) = vxyz(d * pts_cnt + i); }
323  Vector exact_val(vec_dim);
324  F_exact(pos, exact_val);
325  max_err = std::max(max_err, fabs(exact_val(j) - interp_vals[npt]));
326  max_dist = std::max(max_dist, dist_p_out(i));
327  if (code_out[i] == 1 && j == 0) { face_pts++; }
328  }
329  else { if (j == 0) { not_found++; } }
330  npt++;
331  }
332  }
333 
334  cout << setprecision(16)
335  << "Searched points: " << pts_cnt
336  << "\nFound points: " << found
337  << "\nMax interp error: " << max_err
338  << "\nMax dist (of found): " << max_dist
339  << "\nPoints not found: " << not_found
340  << "\nPoints on faces: " << face_pts << endl;
341 
342  // Free the internal gslib data.
343  finder.FreeData();
344 
345  if (prefinement) { delete field_vals_pref; }
346  delete fec;
347 
348  return 0;
349 }
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
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:560
const FiniteElement * GetFE(Geometry::Type geom, int p) const
Variable order version of FiniteElementForGeometry().
Definition: fe_coll.hpp:161
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 void Update(bool want_transform=true)
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation ...
Definition: fespace.cpp:3336
virtual const Vector & GetDist() const
Definition: gslib.hpp:178
Arbitrary order L2 elements in 3D on a cube.
Definition: fe_l2.hpp:62
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:521
void SetElementOrder(int i, int p)
Sets the order of the i&#39;th finite element.
Definition: fespace.cpp:151
void GetBoundingBox(Vector &min, Vector &max, int ref=2)
Returns the minimum and maximum corners of the mesh bounding box.
Definition: mesh.cpp:129
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition: vector.cpp:534
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
int GetElementOrder(int i) const
Returns the order of the i&#39;th finite element.
Definition: fespace.cpp:178
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:928
void MakeOwner(FiniteElementCollection *fec_)
Make the GridFunction the owner of fec and fes.
Definition: gridfunc.hpp:124
double field_func(const Vector &x)
Definition: findpts.cpp:103
virtual void SetL2AvgType(AvgType avgtype_)
Definition: gslib.hpp:153
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
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:433
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9737
GridFunction * ProlongToMaxOrder(const GridFunction *x, const int fieldtype)
Definition: findpts.cpp:46
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
const int visport
virtual const char * Name() const
Definition: fe_coll.hpp:237
void SetNodalFESpace(FiniteElementSpace *nfes)
Definition: mesh.cpp:5273
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:652
int Dimension() const
Definition: mesh.hpp:999
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
virtual DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element &#39;elem&#39;.
Definition: fespace.cpp:2680
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.
void SetIdentityTransformation(Geometry::Type GeomType)
Set the FiniteElement Geometry for the reference elements being used.
Definition: eltrans.cpp:383
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
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
Geometry::Type GetElementGeometry(int i) const
Definition: mesh.hpp:1048
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: gridfunc.cpp:164
Class for integration point with weight.
Definition: intrules.hpp:25
A standard isoparametric element transformation.
Definition: eltrans.hpp:361
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 PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:324
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2395
int open(const char hostname[], int port)
Open the socket stream on &#39;port&#39; at &#39;hostname&#39;.
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
const FiniteElementCollection * FEColl() const
Definition: fespace.hpp:577
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:160
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
virtual void GetTransferMatrix(const FiniteElement &fe, ElementTransformation &Trans, DenseMatrix &I) const
Return interpolation matrix, I, which maps dofs from a coarse element, fe, to the fine dofs on this f...
Definition: fe_base.cpp:118
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