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