MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
field-interp.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 // Field Interp Miniapp: Transfer a grid function between meshes
14 // -------------------------------------------------------------
15 //
16 // This miniapp provides the capability to transfer a grid function (H1, L2,
17 // H(div), and H(curl)) from one mesh onto another using GSLIB-FindPoints. Using
18 // FindPoints, we identify the nodal positions of the target mesh with respect
19 // to the source mesh and then interpolate the source grid function. The
20 // interpolated values are then projected onto the desired finite element space
21 // on the target mesh. Finally, the transferred solution is visualized using
22 // GLVis. Note that the source grid function can be a user-defined vector
23 // function or a grid function file that is compatible with the source mesh.
24 //
25 // Compile with: make field-interp
26 //
27 // Sample runs:
28 // field-interp
29 // field-interp -o 1
30 // field-interp -fts 3 -ft 0
31 // field-interp -m1 triple-pt-1.mesh -s1 triple-pt-1.gf -m2 triple-pt-2.mesh -ft 1
32 // field-interp -m2 ../meshing/amr-quad-q2.mesh -ft 0 -r 1
33 
34 #include "mfem.hpp"
35 #include <fstream>
36 
37 using namespace mfem;
38 using namespace std;
39 
40 // Scalar function to project
41 double scalar_func(const Vector &x)
42 {
43  const int dim = x.Size();
44  double res = 0.0;
45  for (int d = 0; d < dim; d++) { res += x(d) * x(d); }
46  return res;
47 }
48 
49 void vector_func(const Vector &p, Vector &F)
50 {
51  F(0) = scalar_func(p);
52  for (int i = 1; i < F.Size(); i++) { F(i) = (i+1)*pow(-1, i)*F(0); }
53 }
54 
55 int main (int argc, char *argv[])
56 {
57  // Set the method's default parameters.
58  const char *src_mesh_file = "../meshing/square01.mesh";
59  const char *tar_mesh_file = "../../data/inline-tri.mesh";
60  const char *src_sltn_file = "must_be_provided_by_the_user.gf";
61  int src_fieldtype = 0;
62  int src_ncomp = 1;
63  int src_gf_ordering = 0;
64  int ref_levels = 0;
65  int fieldtype = -1;
66  int order = 3;
67  bool visualization = true;
68 
69  // Parse command-line options.
70  OptionsParser args(argc, argv);
71  args.AddOption(&src_mesh_file, "-m1", "--mesh1",
72  "Mesh file for the starting solution.");
73  args.AddOption(&tar_mesh_file, "-m2", "--mesh2",
74  "Mesh file for interpolation.");
75  args.AddOption(&src_sltn_file, "-s1", "--solution1",
76  "(optional) GridFunction file compatible with src_mesh_file."
77  "Set src_fieldtype to -1 if this option is used.");
78  args.AddOption(&src_fieldtype, "-fts", "--field-type-src",
79  "Source GridFunction type:"
80  "0 - H1 (default), 1 - L2, 2 - H(div), 3 - H(curl).");
81  args.AddOption(&src_ncomp, "-nc", "--ncomp",
82  "Number of components for H1 or L2 GridFunctions.");
83  args.AddOption(&src_gf_ordering, "-gfo", "--gfo",
84  "Node ordering: 0 (byNodes) or 1 (byVDim)");
85  args.AddOption(&ref_levels, "-r", "--refine",
86  "Number of refinements of the interpolation mesh.");
87  args.AddOption(&fieldtype, "-ft", "--field-type",
88  "Target GridFunction type: -1 - source GridFunction type (default),"
89  "0 - H1, 1 - L2, 2 - H(div), 3 - H(curl).");
90  args.AddOption(&order, "-o", "--order",
91  "Order of the interpolated solution.");
92  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
93  "--no-visualization",
94  "Enable or disable GLVis visualization.");
95  args.Parse();
96  if (!args.Good())
97  {
98  args.PrintUsage(cout);
99  return 1;
100  }
101  args.PrintOptions(cout);
102 
103  // Input meshes.
104  Mesh mesh_1(src_mesh_file, 1, 1, false);
105  Mesh mesh_2(tar_mesh_file, 1, 1, false);
106  const int dim = mesh_1.Dimension();
107  MFEM_ASSERT(dim == mesh_2.Dimension(), "Source and target meshes "
108  "must be in the same dimension.");
109  MFEM_VERIFY(dim > 1, "GSLIB requires a 2D or a 3D mesh" );
110 
111  for (int lev = 0; lev < ref_levels; lev++)
112  {
113  mesh_2.UniformRefinement();
114  }
115 
116  if (mesh_1.GetNodes() == NULL) { mesh_1.SetCurvature(1); }
117  if (mesh_2.GetNodes() == NULL) { mesh_2.SetCurvature(1); }
118  const int mesh_poly_deg =
119  mesh_2.GetNodes()->FESpace()->GetElementOrder(0);
120  cout << "Source mesh curvature: "
121  << mesh_1.GetNodes()->OwnFEC()->Name() << endl
122  << "Target mesh curvature: "
123  << mesh_2.GetNodes()->OwnFEC()->Name() << endl;
124 
125  int src_vdim = src_ncomp;
126  FiniteElementCollection *src_fec = NULL;
127  FiniteElementSpace *src_fes = NULL;
128  GridFunction *func_source = NULL;
129  if (src_fieldtype < 0) // use src_sltn_file
130  {
131  ifstream mat_stream_1(src_sltn_file);
132  func_source = new GridFunction(&mesh_1, mat_stream_1);
133  src_vdim = func_source->FESpace()->GetVDim();
134  }
135  else if (src_fieldtype == 0)
136  {
137  src_fec = new H1_FECollection(order, dim);
138  }
139  else if (src_fieldtype == 1)
140  {
141  src_fec = new L2_FECollection(order, dim);
142  }
143  else if (src_fieldtype == 2)
144  {
145  src_fec = new RT_FECollection(order, dim);
146  src_ncomp = 1;
147  src_vdim = dim;
148  }
149  else if (src_fieldtype == 3)
150  {
151  src_fec = new ND_FECollection(order, dim);
152  src_ncomp = 1;
153  src_vdim = dim;
154  }
155  else
156  {
157  MFEM_ABORT("Invalid FECollection type.");
158  }
159 
160  if (src_fieldtype > -1)
161  {
162  MFEM_VERIFY(src_gf_ordering >= 0 && src_gf_ordering <= 1,
163  " Source grid function ordering must be 0 (byNodes)"
164  " or 1 (byVDIM.");
165  src_fes = new FiniteElementSpace(&mesh_1, src_fec, src_ncomp, src_gf_ordering);
166  func_source = new GridFunction(src_fes);
167  // Project the grid function using VectorFunctionCoefficient.
169  func_source->ProjectCoefficient(F);
170  }
171 
172  // Display the starting mesh and the field.
173  if (visualization)
174  {
175  char vishost[] = "localhost";
176  int visport = 19916;
177  socketstream sout1;
178  sout1.open(vishost, visport);
179  if (!sout1)
180  {
181  cout << "Unable to connect to GLVis server at "
182  << vishost << ':' << visport << endl;
183  }
184  else
185  {
186  sout1.precision(8);
187  sout1 << "solution\n" << mesh_1 << *func_source
188  << "window_title 'Source mesh and solution'"
189  << "window_geometry 0 0 600 600";
190  if (dim == 2) { sout1 << "keys RmjAc"; }
191  if (dim == 3) { sout1 << "keys mA\n"; }
192  sout1 << flush;
193  }
194  }
195 
196  const Geometry::Type gt = mesh_2.GetNodalFESpace()->GetFE(0)->GetGeomType();
197  MFEM_VERIFY(gt != Geometry::PRISM, "Wedge elements are not currently "
198  "supported.");
199  MFEM_VERIFY(mesh_2.GetNumGeometries(mesh_2.Dimension()) == 1, "Mixed meshes"
200  "are not currently supported.");
201 
202  // Ensure the source grid function can be transferred using GSLIB-FindPoints.
203  const FiniteElementCollection *fec_in = func_source->FESpace()->FEColl();
204  std::cout << "Source FE collection: " << fec_in->Name() << std::endl;
205 
206  if (src_fieldtype < 0)
207  {
208  const H1_FECollection *fec_h1 = dynamic_cast<const H1_FECollection *>(fec_in);
209  const L2_FECollection *fec_l2 = dynamic_cast<const L2_FECollection *>(fec_in);
210  const RT_FECollection *fec_rt = dynamic_cast<const RT_FECollection *>(fec_in);
211  const ND_FECollection *fec_nd = dynamic_cast<const ND_FECollection *>(fec_in);
212  if (fec_h1) { src_fieldtype = 0; }
213  else if (fec_l2) { src_fieldtype = 1; }
214  else if (fec_rt) { src_fieldtype = 2; }
215  else if (fec_nd) { src_fieldtype = 3; }
216  else { MFEM_ABORT("GridFunction type not supported yet."); }
217  }
218  if (fieldtype < 0) { fieldtype = src_fieldtype; }
219 
220  // Setup the FiniteElementSpace and GridFunction on the target mesh.
221  FiniteElementCollection *tar_fec = NULL;
222  FiniteElementSpace *tar_fes = NULL;
223 
224  int tar_vdim = src_vdim;
225  if (fieldtype == 0)
226  {
227  tar_fec = new H1_FECollection(order, dim);
228  tar_vdim = (src_fieldtype > 1) ? dim : src_vdim;
229  }
230  else if (fieldtype == 1)
231  {
232  tar_fec = new L2_FECollection(order, dim);
233  tar_vdim = (src_fieldtype > 1) ? dim : src_vdim;
234  }
235  else if (fieldtype == 2)
236  {
237  tar_fec = new RT_FECollection(order, dim);
238  tar_vdim = 1;
239  MFEM_VERIFY(src_fieldtype > 1, "Cannot interpolate a scalar "
240  "grid function to a vector");
241 
242  }
243  else if (fieldtype == 3)
244  {
245  tar_fec = new ND_FECollection(order, dim);
246  tar_vdim = 1;
247  MFEM_VERIFY(src_fieldtype > 1, "Cannot interpolate a scalar "
248  "grid function to a vector");
249  }
250  else
251  {
252  MFEM_ABORT("GridFunction type not supported.");
253  }
254  std::cout << "Target FE collection: " << tar_fec->Name() << std::endl;
255  tar_fes = new FiniteElementSpace(&mesh_2, tar_fec, tar_vdim,
256  src_fes->GetOrdering());
257  GridFunction func_target(tar_fes);
258 
259  const int NE = mesh_2.GetNE(),
260  nsp = tar_fes->GetFE(0)->GetNodes().GetNPoints(),
261  tar_ncomp = func_target.VectorDim();
262 
263  // Generate list of points where the grid function will be evaluated.
264  Vector vxyz;
265  int point_ordering;
266  if (fieldtype == 0 && order == mesh_poly_deg)
267  {
268  vxyz = *mesh_2.GetNodes();
269  point_ordering = mesh_2.GetNodes()->FESpace()->GetOrdering();
270  }
271  else
272  {
273  vxyz.SetSize(nsp*NE*dim);
274  for (int i = 0; i < NE; i++)
275  {
276  const FiniteElement *fe = tar_fes->GetFE(i);
277  const IntegrationRule ir = fe->GetNodes();
279 
280  DenseMatrix pos;
281  et->Transform(ir, pos);
282  Vector rowx(vxyz.GetData() + i*nsp, nsp),
283  rowy(vxyz.GetData() + i*nsp + NE*nsp, nsp),
284  rowz;
285  if (dim == 3)
286  {
287  rowz.SetDataAndSize(vxyz.GetData() + i*nsp + 2*NE*nsp, nsp);
288  }
289  pos.GetRow(0, rowx);
290  pos.GetRow(1, rowy);
291  if (dim == 3) { pos.GetRow(2, rowz); }
292  }
293  point_ordering = Ordering::byNODES;
294  }
295  const int nodes_cnt = vxyz.Size() / dim;
296 
297  // Evaluate source grid function.
298  Vector interp_vals(nodes_cnt*tar_ncomp);
299  FindPointsGSLIB finder;
300  finder.Setup(mesh_1);
301  finder.Interpolate(vxyz, *func_source, interp_vals, point_ordering);
302 
303  // Project the interpolated values to the target FiniteElementSpace.
304  if (fieldtype <= 1) // H1 or L2
305  {
306  if ((fieldtype == 0 && order == mesh_poly_deg) || fieldtype == 1)
307  {
308  func_target = interp_vals;
309  }
310  else // H1 - but mesh order != GridFunction order
311  {
312  Array<int> vdofs;
313  Vector vals;
314  Vector elem_dof_vals(nsp*tar_ncomp);
315 
316  for (int i = 0; i < mesh_2.GetNE(); i++)
317  {
318  tar_fes->GetElementVDofs(i, vdofs);
319  vals.SetSize(vdofs.Size());
320  for (int j = 0; j < nsp; j++)
321  {
322  for (int d = 0; d < tar_ncomp; d++)
323  {
324  // Arrange values byNodes
325  int idx = src_fes->GetOrdering() == Ordering::byNODES ?
326  d*nsp*NE + i*nsp + j : i*nsp*dim + d + j*dim;
327  elem_dof_vals(j + d*nsp) = interp_vals(idx);
328  }
329  }
330  func_target.SetSubVector(vdofs, elem_dof_vals);
331  }
332  }
333  }
334  else // H(div) or H(curl)
335  {
336  Array<int> vdofs;
337  Vector vals;
338  Vector elem_dof_vals(nsp*tar_ncomp);
339 
340  for (int i = 0; i < mesh_2.GetNE(); i++)
341  {
342  tar_fes->GetElementVDofs(i, vdofs);
343  vals.SetSize(vdofs.Size());
344  for (int j = 0; j < nsp; j++)
345  {
346  for (int d = 0; d < tar_ncomp; d++)
347  {
348  // Arrange values byVDim
349  int idx = src_fes->GetOrdering() == Ordering::byNODES ?
350  d*nsp*NE + i*nsp + j : i*nsp*dim + d + j*dim;
351  elem_dof_vals(j*tar_ncomp+d) = interp_vals(idx);
352  }
353  }
354  tar_fes->GetFE(i)->ProjectFromNodes(elem_dof_vals,
355  *tar_fes->GetElementTransformation(i),
356  vals);
357  func_target.SetSubVector(vdofs, vals);
358  }
359  }
360 
361  // Visualize the transferred solution.
362  if (visualization)
363  {
364  char vishost[] = "localhost";
365  int visport = 19916;
366  socketstream sout1;
367  sout1.open(vishost, visport);
368  if (!sout1)
369  {
370  cout << "Unable to connect to GLVis server at "
371  << vishost << ':' << visport << endl;
372  }
373  else
374  {
375  sout1.precision(8);
376  sout1 << "solution\n" << mesh_2 << func_target
377  << "window_title 'Target mesh and solution'"
378  << "window_geometry 600 0 600 600";
379  if (dim == 2) { sout1 << "keys RmjAc"; }
380  if (dim == 3) { sout1 << "keys mA\n"; }
381  sout1 << flush;
382  }
383  }
384 
385  // Output the target mesh with the interpolated solution.
386  ostringstream rho_name;
387  rho_name << "interpolated.gf";
388  ofstream rho_ofs(rho_name.str().c_str());
389  rho_ofs.precision(8);
390  func_target.Save(rho_ofs);
391  rho_ofs.close();
392 
393  // Free the internal gslib data.
394  finder.FreeData();
395 
396  // Delete remaining memory.
397  delete func_source;
398  delete src_fes;
399  delete src_fec;
400  delete tar_fes;
401  delete tar_fec;
402 
403  return 0;
404 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
Abstract class for all finite elements.
Definition: fe_base.hpp:235
const char vishost[]
virtual void Interpolate(const GridFunction &field_in, Vector &field_out)
Definition: gslib.cpp:801
Ordering::Type GetOrdering() const
Return the ordering method.
Definition: fespace.hpp:599
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
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
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition: mesh.cpp:5908
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:200
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:923
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:209
int close()
Close the socketstream.
virtual void ProjectFromNodes(Vector &vc, ElementTransformation &Trans, Vector &dofs) const
Given a vector of values at the finite element nodes and a transformation, compute its projection (ap...
Definition: fe_base.cpp:137
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
double scalar_func(const Vector &x)
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:281
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:265
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9816
void GetRow(int r, Vector &row) const
Definition: densemat.cpp:1276
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:5343
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
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition: fe_base.hpp:390
const int visport
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:581
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:668
int Dimension() const
Definition: mesh.hpp:1006
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:341
A general vector function coefficient.
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:647
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
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
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
Definition: vector.hpp:157
virtual const char * Name() const
Definition: fe_coll.hpp:65
void vector_func(const Vector &p, Vector &F)
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:5338
int dim
Definition: ex24.cpp:53
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2781
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:2399
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:415
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:601
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7908
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:220
int main()
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:288
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150