MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
lor-transfer.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, 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 // LOR Transfer Miniapp: Map functions between HO and LOR spaces
14 // --------------------------------------------------------------
15 //
16 // This miniapp visualizes the maps between a high-order (HO) finite element
17 // space, typically using high-order functions on a high-order mesh, and a
18 // low-order refined (LOR) finite element space, typically defined by 0th or 1st
19 // order functions on a low-order refinement of the HO mesh.
20 //
21 // Two main operators are illustrated:
22 //
23 // 1. R: HO -> LOR, defined by FiniteElementSpace::GetTransferOperator
24 // 2. P: LOR -> HO, defined by FiniteElementSpace::GetReverseTransferOperator
25 //
26 // While defined generally, these operators have some nice properties for
27 // particular finite element spaces. For example they satisfy PR=I, plus mass
28 // conservation in both directions for L2 fields.
29 //
30 // Compile with: make lor-transfer
31 //
32 // Sample runs: lor-transfer
33 // lor-transfer -h1
34 // lor-transfer -t
35 // lor-transfer -m ../../data/star-q2.mesh -lref 5 -p 4
36 // lor-transfer -lref 4 -o 4 -lo 0 -p 1
37 // lor-transfer -lref 5 -o 4 -lo 0 -p 1
38 // lor-transfer -lref 5 -o 4 -lo 3 -p 2
39 // lor-transfer -lref 5 -o 4 -lo 0 -p 3
40 
41 #include "mfem.hpp"
42 #include <fstream>
43 #include <iostream>
44 
45 using namespace std;
46 using namespace mfem;
47 
48 int problem = 1; // problem type
49 
50 int Wx = 0, Wy = 0; // window position
51 int Ww = 350, Wh = 350; // window size
52 int offx = Ww+5, offy = Wh+25; // window offsets
53 
54 string space;
55 string direction;
56 
57 // Exact functions to project
58 double RHO_exact(const Vector &x);
59 
60 // Helper functions
61 void visualize(VisItDataCollection &, string, int, int);
63  string);
64 
65 int main(int argc, char *argv[])
66 {
67  // Parse command-line options.
68  const char *mesh_file = "../../data/star.mesh";
69  int order = 4;
70  int lref = order;
71  int lorder = 0;
72  bool vis = true;
73  bool useH1 = false;
74  bool use_transfer = false;
75 
76  OptionsParser args(argc, argv);
77  args.AddOption(&mesh_file, "-m", "--mesh",
78  "Mesh file to use.");
79  args.AddOption(&problem, "-p", "--problem",
80  "Problem type (see the RHO_exact function).");
81  args.AddOption(&order, "-o", "--order",
82  "Finite element order (polynomial degree) or -1 for"
83  " isoparametric space.");
84  args.AddOption(&lref, "-lref", "--lor-ref-level", "LOR refinement level.");
85  args.AddOption(&lorder, "-lo", "--lor-order",
86  "LOR space order (polynomial degree, zero by default).");
87  args.AddOption(&vis, "-vis", "--visualization", "-no-vis",
88  "--no-visualization",
89  "Enable or disable GLVis visualization.");
90  args.AddOption(&useH1, "-h1", "--use-h1", "-l2", "--use-l2",
91  "Use H1 spaces instead of L2.");
92  args.AddOption(&use_transfer, "-t", "--use-pointwise-transfer", "-no-t",
93  "--dont-use-pointwise-transfer",
94  "Use pointwise transfer operators instead of L2 projection.");
95  args.Parse();
96  if (!args.Good())
97  {
98  args.PrintUsage(cout);
99  return 1;
100  }
101  args.PrintOptions(cout);
102 
103  // Read the mesh from the given mesh file.
104  Mesh mesh(mesh_file, 1, 1);
105  int dim = mesh.Dimension();
106 
107  // Create the low-order refined mesh
108  int basis_lor = BasisType::GaussLobatto; // BasisType::ClosedUniform;
109  Mesh mesh_lor(&mesh, lref, basis_lor);
110 
111  // Create spaces
112  FiniteElementCollection *fec, *fec_lor;
113  if (useH1)
114  {
115  space = "H1";
116  if (lorder == 0)
117  {
118  lorder = 1;
119  cerr << "Switching the H1 LOR space order from 0 to 1\n";
120  }
121  fec = new H1_FECollection(order-1, dim);
122  fec_lor = new H1_FECollection(lorder, dim);
123  }
124  else
125  {
126  space = "L2";
127  fec = new L2_FECollection(order-1, dim);
128  fec_lor = new L2_FECollection(lorder, dim);
129  }
130 
131  FiniteElementSpace fespace(&mesh, fec);
132  FiniteElementSpace fespace_lor(&mesh_lor, fec_lor);
133 
134  GridFunction rho(&fespace);
135  GridFunction rho_lor(&fespace_lor);
136 
137  // Data collections for vis/analysis
138  VisItDataCollection HO_dc("HO", &mesh);
139  HO_dc.RegisterField("density", &rho);
140  VisItDataCollection LOR_dc("LOR", &mesh_lor);
141  LOR_dc.RegisterField("density", &rho_lor);
142 
143 
144  // HO projections
145  direction = "HO -> LOR @ HO";
147  rho.ProjectCoefficient(RHO);
148  double ho_mass = compute_mass(&fespace, -1.0, HO_dc, "HO ");
149  if (vis) { visualize(HO_dc, "HO", Wx, Wy); Wx += offx; }
150 
151  GridTransfer *gt;
152  if (use_transfer)
153  {
154  gt = new InterpolationGridTransfer(fespace, fespace_lor);
155  }
156  else
157  {
158  gt = new L2ProjectionGridTransfer(fespace, fespace_lor);
159  }
160  const Operator &R = gt->ForwardOperator();
161  const Operator &P = gt->BackwardOperator();
162 
163  // HO->LOR restriction
164  direction = "HO -> LOR @ LOR";
165  R.Mult(rho, rho_lor);
166  compute_mass(&fespace_lor, ho_mass, LOR_dc, "R(HO) ");
167  if (vis) { visualize(LOR_dc, "R(HO)", Wx, Wy); Wx += offx; }
168 
169  // LOR->HO prolongation
170  direction = "HO -> LOR @ HO";
171  GridFunction rho_prev = rho;
172  P.Mult(rho_lor, rho);
173  compute_mass(&fespace, ho_mass, HO_dc, "P(R(HO)) ");
174  if (vis) { visualize(HO_dc, "P(R(HO))", Wx, Wy); Wx = 0; Wy += offy; }
175 
176  rho_prev -= rho;
177  cout.precision(12);
178  cout << "|HO - P(R(HO))|_∞ = " << rho_prev.Normlinf() << endl << endl;
179 
180  // LOR projections
181  direction = "LOR -> HO @ LOR";
182  rho_lor.ProjectCoefficient(RHO);
183  GridFunction rho_lor_prev = rho_lor;
184  double lor_mass = compute_mass(&fespace_lor, -1.0, LOR_dc, "LOR ");
185  if (vis) { visualize(LOR_dc, "LOR", Wx, Wy); Wx += offx; }
186 
187  // Prolongate to HO space
188  direction = "LOR -> HO @ HO";
189  P.Mult(rho_lor, rho);
190  compute_mass(&fespace, lor_mass, HO_dc, "P(LOR) ");
191  if (vis) { visualize(HO_dc, "P(LOR)", Wx, Wy); Wx += offx; }
192 
193  // Restrict back to LOR space. This won't give the original function because
194  // the rho_lor doesn't necessarily live in the range of R.
195  direction = "LOR -> HO @ LOR";
196  R.Mult(rho, rho_lor);
197  compute_mass(&fespace_lor, lor_mass, LOR_dc, "R(P(LOR))");
198  if (vis) { visualize(LOR_dc, "R(P(LOR))", Wx, Wy); }
199 
200  rho_lor_prev -= rho_lor;
201  cout.precision(12);
202  cout << "|LOR - R(P(LOR))|_∞ = " << rho_lor_prev.Normlinf() << endl;
203 
204  delete fec;
205  delete fec_lor;
206 
207  return 0;
208 }
209 
210 
211 double RHO_exact(const Vector &x)
212 {
213  switch (problem)
214  {
215  case 1: // smooth field
216  return x(1)+0.25*cos(2*M_PI*x.Norml2());
217  case 2: // cubic function
218  return x(1)*x(1)*x(1) + 2*x(0)*x(1) + x(0);
219  case 3: // sharp gradient
220  return M_PI/2-atan(5*(2*x.Norml2()-1));
221  case 4: // basis function
222  return (x.Norml2() < 0.1) ? 1 : 0;
223  default:
224  return 1.0;
225  }
226 }
227 
228 
229 void visualize(VisItDataCollection &dc, string prefix, int x, int y)
230 {
231  int w = Ww, h = Wh;
232 
233  char vishost[] = "localhost";
234  int visport = 19916;
235 
236  socketstream sol_sockL2(vishost, visport);
237  sol_sockL2.precision(8);
238  sol_sockL2 << "solution\n" << *dc.GetMesh() << *dc.GetField("density")
239  << "window_geometry " << x << " " << y << " " << w << " " << h
240  << "plot_caption '" << space << " " << prefix << " Density'"
241  << "window_title '" << direction << "'" << flush;
242 }
243 
244 
245 double compute_mass(FiniteElementSpace *L2, double massL2,
246  VisItDataCollection &dc, string prefix)
247 {
248  ConstantCoefficient one(1.0);
249  BilinearForm ML2(L2);
250  ML2.AddDomainIntegrator(new MassIntegrator(one));
251  ML2.Assemble();
252 
253  GridFunction rhoone(L2);
254  rhoone = 1.0;
255 
256  double newmass = ML2.InnerProduct(*dc.GetField("density"),rhoone);
257  cout.precision(18);
258  cout << space << " " << prefix << " mass = " << newmass;
259  if (massL2 >= 0)
260  {
261  cout.precision(4);
262  cout << " (" << fabs(newmass-massL2)*100/massL2 << "%)";
263  }
264  cout << endl;
265  return newmass;
266 }
void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
Definition: ex10.cpp:379
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
GridFunction * GetField(const std::string &field_name)
Get a pointer to a grid function in the collection.
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:744
int offx
int Wx
virtual const Operator & BackwardOperator()=0
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
double Normlinf() const
Returns the l_infinity norm of the vector.
Definition: vector.cpp:762
double RHO_exact(const Vector &x)
int main(int argc, char *argv[])
Definition: ex1.cpp:66
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 const Operator & ForwardOperator()=0
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
constexpr int visport
Data collection with VisIt I/O routines.
int Dimension() const
Definition: mesh.hpp:788
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
Base class for transfer algorithms that construct transfer Operators between two finite element (FE) ...
Definition: fespace.hpp:778
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
Definition: fespace.hpp:870
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
int problem
Definition: ex15.cpp:54
string space
double InnerProduct(const Vector &x, const Vector &y) const
Compute .
string direction
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
int Wh
int dim
Definition: ex24.cpp:53
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:304
Transfer data between a coarse mesh and an embedded refined mesh using L2 projection.
Definition: fespace.hpp:922
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:2252
Mesh * GetMesh()
Get a pointer to the mesh in the collection.
A general function coefficient.
Vector data type.
Definition: vector.hpp:51
int offy
int Wy
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:159
int Ww
Abstract operator.
Definition: operator.hpp:24
double compute_mass(FiniteElementSpace *, double, VisItDataCollection &, string)
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:221
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145