MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
plor-transfer.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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 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// The grid transfer operators are represented using either
22// InterpolationGridTransfer or L2ProjectionGridTransfer (depending on the
23// options requested by the user). The two transfer operators are then:
24//
25// 1. R: HO -> LOR, defined by GridTransfer::ForwardOperator
26// 2. P: LOR -> HO, defined by GridTransfer::BackwardOperator
27//
28// While defined generally, these operators have some nice properties for
29// particular finite element spaces. For example they satisfy PR=I, plus mass
30// conservation in both directions for L2 fields.
31//
32// Compile with: make plor-transfer
33//
34// Sample runs: plor-transfer
35// plor-transfer -h1
36// plor-transfer -t
37// plor-transfer -m ../../data/star-q2.mesh -lref 5 -p 4
38// plor-transfer -m ../../data/star-mixed.mesh -lref 3 -p 2
39// plor-transfer -lref 4 -o 4 -lo 0 -p 1
40// plor-transfer -lref 5 -o 4 -lo 0 -p 1
41// plor-transfer -lref 5 -o 4 -lo 3 -p 2
42// plor-transfer -lref 5 -o 4 -lo 0 -p 3
43
44#include "mfem.hpp"
45#include <fstream>
46#include <iostream>
47
48using namespace std;
49using namespace mfem;
50
51int problem = 1; // problem type
52
53int Wx = 0, Wy = 0; // window position
54int Ww = 350, Wh = 350; // window size
55int offx = Ww+5, offy = Wh+25; // window offsets
56
57string space;
58string direction;
59
60// Exact functions to project
61real_t RHO_exact(const Vector &x);
62
63// Helper functions
64void visualize(VisItDataCollection &, string, int, int);
66 string);
67
68int main(int argc, char *argv[])
69{
70 // Initialize MPI and HYPRE.
71 Mpi::Init(argc, argv);
73
74 // Parse command-line options.
75 const char *mesh_file = "../../data/star.mesh";
76 int order = 3;
77 int lref = order+1;
78 int lorder = 0;
79 bool vis = true;
80 bool useH1 = false;
81 bool use_pointwise_transfer = false;
82
83 OptionsParser args(argc, argv);
84 args.AddOption(&mesh_file, "-m", "--mesh",
85 "Mesh file to use.");
86 args.AddOption(&problem, "-p", "--problem",
87 "Problem type (see the RHO_exact function).");
88 args.AddOption(&order, "-o", "--order",
89 "Finite element order (polynomial degree) or -1 for"
90 " isoparametric space.");
91 args.AddOption(&lref, "-lref", "--lor-ref-level", "LOR refinement level.");
92 args.AddOption(&lorder, "-lo", "--lor-order",
93 "LOR space order (polynomial degree, zero by default).");
94 args.AddOption(&vis, "-vis", "--visualization", "-no-vis",
95 "--no-visualization",
96 "Enable or disable GLVis visualization.");
97 args.AddOption(&useH1, "-h1", "--use-h1", "-l2", "--use-l2",
98 "Use H1 spaces instead of L2.");
99 args.AddOption(&use_pointwise_transfer, "-t", "--use-pointwise-transfer",
100 "-no-t", "--dont-use-pointwise-transfer",
101 "Use pointwise transfer operators instead of L2 projection.");
102 args.ParseCheck();
103
104 // Read the mesh from the given mesh file.
105 Mesh serial_mesh(mesh_file, 1, 1);
106 ParMesh mesh(MPI_COMM_WORLD, serial_mesh);
107 serial_mesh.Clear();
108 int dim = mesh.Dimension();
109
110 // Create the low-order refined mesh
111 int basis_lor = BasisType::GaussLobatto; // BasisType::ClosedUniform;
112 ParMesh mesh_lor = ParMesh::MakeRefined(mesh, lref, basis_lor);
113
114 // Create spaces
115 FiniteElementCollection *fec, *fec_lor;
116 if (useH1)
117 {
118 space = "H1";
119 if (lorder == 0)
120 {
121 lorder = 1;
122 if (Mpi::Root())
123 {
124 cerr << "Switching the H1 LOR space order from 0 to 1\n";
125 }
126 }
127 fec = new H1_FECollection(order, dim);
128 fec_lor = new H1_FECollection(lorder, dim);
129 }
130 else
131 {
132 space = "L2";
133 fec = new L2_FECollection(order, dim);
134 fec_lor = new L2_FECollection(lorder, dim);
135 }
136
137 ParFiniteElementSpace fespace(&mesh, fec);
138 ParFiniteElementSpace fespace_lor(&mesh_lor, fec_lor);
139
140 ParGridFunction rho(&fespace);
141 ParGridFunction rho_lor(&fespace_lor);
142
143 // Data collections for vis/analysis
144 VisItDataCollection HO_dc(MPI_COMM_WORLD, "HO", &mesh);
145 HO_dc.RegisterField("density", &rho);
146 VisItDataCollection LOR_dc(MPI_COMM_WORLD, "LOR", &mesh_lor);
147 LOR_dc.RegisterField("density", &rho_lor);
148
149 ParBilinearForm M_ho(&fespace);
151 M_ho.Assemble();
152 M_ho.Finalize();
153 HypreParMatrix* M_ho_tdof = M_ho.ParallelAssemble();
154
155 ParBilinearForm M_lor(&fespace_lor);
157 M_lor.Assemble();
158 M_lor.Finalize();
159 HypreParMatrix* M_lor_tdof = M_lor.ParallelAssemble();
160
161 // HO projections
162 direction = "HO -> LOR @ HO";
164 rho.ProjectCoefficient(RHO);
165 // Make sure AMR constraints are satisfied
166 rho.SetTrueVector();
167 rho.SetFromTrueVector();
168
169 real_t ho_mass = compute_mass(&fespace, -1.0, HO_dc, "HO ");
170 if (vis) { visualize(HO_dc, "HO", Wx, Wy); Wx += offx; }
171
172 GridTransfer *gt;
173 if (use_pointwise_transfer)
174 {
175 gt = new InterpolationGridTransfer(fespace, fespace_lor);
176 }
177 else
178 {
179 gt = new L2ProjectionGridTransfer(fespace, fespace_lor);
180 }
181 const Operator &R = gt->ForwardOperator();
182
183 // HO->LOR restriction
184 direction = "HO -> LOR @ LOR";
185 R.Mult(rho, rho_lor);
186 compute_mass(&fespace_lor, ho_mass, LOR_dc, "R(HO) ");
187 if (vis) { visualize(LOR_dc, "R(HO)", Wx, Wy); Wx += offx; }
188 auto global_max = [](const Vector& v)
189 {
190 real_t max = v.Normlinf();
191 MPI_Allreduce(MPI_IN_PLACE, &max, 1, MPITypeMap<real_t>::mpi_type,
192 MPI_MAX, MPI_COMM_WORLD);
193 return max;
194 };
195
197 {
198 const Operator &P = gt->BackwardOperator();
199 // LOR->HO prolongation
200 direction = "HO -> LOR @ HO";
201 ParGridFunction rho_prev = rho;
202 P.Mult(rho_lor, rho);
203 compute_mass(&fespace, ho_mass, HO_dc, "P(R(HO)) ");
204 if (vis) { visualize(HO_dc, "P(R(HO))", Wx, Wy); Wx = 0; Wy += offy; }
205
206 rho_prev -= rho;
207 Vector rho_prev_true(fespace.GetTrueVSize());
208 rho_prev.GetTrueDofs(rho_prev_true);
209 real_t l_inf = global_max(rho_prev_true);
210 if (Mpi::Root())
211 {
212 cout.precision(12);
213 cout << "|HO - P(R(HO))|_∞ = " << l_inf << endl;
214 }
215 }
216
217 // HO* to LOR* dual fields
218 ParLinearForm M_rho(&fespace), M_rho_lor(&fespace_lor);
219 auto global_sum = [](const Vector& v)
220 {
221 real_t sum = v.Sum();
222 MPI_Allreduce(MPI_IN_PLACE, &sum, 1, MPITypeMap<real_t>::mpi_type,
223 MPI_SUM, MPI_COMM_WORLD);
224 return sum;
225 };
226 if (!use_pointwise_transfer && gt->SupportsBackwardsOperator())
227 {
228 Vector M_rho_true(fespace.GetTrueVSize());
229 M_ho_tdof->Mult(rho.GetTrueVector(), M_rho_true);
230 fespace.GetRestrictionOperator()->MultTranspose(M_rho_true, M_rho);
231 const Operator &P = gt->BackwardOperator();
232 P.MultTranspose(M_rho, M_rho_lor);
233 real_t ho_dual_mass = global_sum(M_rho);
234 real_t lor_dual_mass = global_sum(M_rho_lor);
235 if (Mpi::Root())
236 {
237 cout << "HO -> LOR dual field: " << abs(ho_dual_mass - lor_dual_mass) << "\n\n";
238 }
239 }
240
241 // LOR projections
242 direction = "LOR -> HO @ LOR";
243 rho_lor.ProjectCoefficient(RHO);
244 ParGridFunction rho_lor_prev = rho_lor;
245 real_t lor_mass = compute_mass(&fespace_lor, -1.0, LOR_dc, "LOR ");
246 if (vis) { visualize(LOR_dc, "LOR", Wx, Wy); Wx += offx; }
247
249 {
250 const Operator &P = gt->BackwardOperator();
251 // Prolongate to HO space
252 direction = "LOR -> HO @ HO";
253 P.Mult(rho_lor, rho);
254 compute_mass(&fespace, lor_mass, HO_dc, "P(LOR) ");
255 if (vis) { visualize(HO_dc, "P(LOR)", Wx, Wy); Wx += offx; }
256
257 // Restrict back to LOR space. This won't give the original function because
258 // the rho_lor doesn't necessarily live in the range of R.
259 direction = "LOR -> HO @ LOR";
260 R.Mult(rho, rho_lor);
261 compute_mass(&fespace_lor, lor_mass, LOR_dc, "R(P(LOR))");
262 if (vis) { visualize(LOR_dc, "R(P(LOR))", Wx, Wy); }
263
264 rho_lor_prev -= rho_lor;
265 Vector rho_lor_prev_true(fespace_lor.GetTrueVSize());
266 rho_lor_prev.GetTrueDofs(rho_lor_prev_true);
267 real_t l_inf = global_max(rho_lor_prev_true);
268 if (Mpi::Root())
269 {
270 cout.precision(12);
271 cout << "|LOR - R(P(LOR))|_∞ = " << l_inf << endl;
272 }
273 }
274
275 // LOR* to HO* dual fields
276 if (!use_pointwise_transfer)
277 {
278 Vector M_rho_lor_true(fespace_lor.GetTrueVSize());
279 M_lor_tdof->Mult(rho_lor.GetTrueVector(), M_rho_lor_true);
280 fespace_lor.GetRestrictionOperator()->MultTranspose(M_rho_lor_true,
281 M_rho_lor);
282 R.MultTranspose(M_rho_lor, M_rho);
283 real_t ho_dual_mass = global_sum(M_rho);
284 real_t lor_dual_mass = global_sum(M_rho_lor);
285
286 cout << lor_dual_mass << '\n';
287 cout << ho_dual_mass << '\n';
288
289 if (Mpi::Root())
290 {
291 cout << "LOR -> HO dual field: " << abs(ho_dual_mass - lor_dual_mass) << '\n';
292 }
293 }
294
295 delete fec;
296 delete fec_lor;
297 delete M_ho_tdof;
298 delete M_lor_tdof;
299 delete gt;
300
301 return 0;
302}
303
304
306{
307 switch (problem)
308 {
309 case 1: // smooth field
310 return x(1)+0.25*cos(2*M_PI*x.Norml2());
311 case 2: // cubic function
312 return x(1)*x(1)*x(1) + 2*x(0)*x(1) + x(0);
313 case 3: // sharp gradient
314 return M_PI/2-atan(5*(2*x.Norml2()-1));
315 case 4: // basis function
316 return (x.Norml2() < 0.1) ? 1 : 0;
317 default:
318 return 1.0;
319 }
320}
321
322
323void visualize(VisItDataCollection &dc, string prefix, int x, int y)
324{
325 int w = Ww, h = Wh;
326
327 char vishost[] = "localhost";
328 int visport = 19916;
329
330 socketstream sol_sockL2(vishost, visport);
331 sol_sockL2 << "parallel " << Mpi::WorldSize() << " " << Mpi::WorldRank() <<
332 "\n";
333 sol_sockL2.precision(8);
334 sol_sockL2 << "solution\n" << *dc.GetMesh() << *dc.GetField("density")
335 << "window_geometry " << x << " " << y << " " << w << " " << h
336 << "plot_caption '" << space << " " << prefix << " Density'"
337 << "window_title '" << direction << "'" << flush;
338}
339
340
342 VisItDataCollection &dc, string prefix)
343{
344 ConstantCoefficient one(1.0);
345 ParLinearForm lf(L2);
347 lf.Assemble();
348
349 real_t newmass = lf(*dc.GetParField("density"));
350 if (Mpi::Root())
351 {
352 cout.precision(18);
353 cout << space << " " << prefix << " mass = " << newmass;
354 if (massL2 >= 0)
355 {
356 cout.precision(4);
357 cout << " (" << fabs(newmass-massL2)*100/massL2 << "%)";
358 }
359 cout << endl;
360 }
361 return newmass;
362}
@ GaussLobatto
Closed type.
Definition fe_base.hpp:32
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY....
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
A coefficient that is constant across space and time.
GridFunction * GetField(const std::string &field_name)
Get a pointer to a grid function in the collection.
ParGridFunction * GetParField(const std::string &field_name)
Get a pointer to a parallel grid function in the collection.
Mesh * GetMesh()
Get a pointer to the mesh in the collection.
Class for domain integration .
Definition lininteg.hpp:109
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
A general function coefficient.
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition gridfunc.hpp:144
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition gridfunc.hpp:150
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition gridfunc.hpp:130
Base class for transfer algorithms that construct transfer Operators between two finite element (FE) ...
Definition transfer.hpp:30
virtual bool SupportsBackwardsOperator() const
Definition transfer.hpp:102
virtual const Operator & ForwardOperator()=0
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
virtual const Operator & BackwardOperator()=0
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, real_t alpha=1.0, real_t beta=0.0) const
Computes y = alpha * A * x + beta * y.
Definition hypre.cpp:1815
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.hpp:74
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
Definition transfer.hpp:124
Transfer data in L2 and H1 finite element spaces between a coarse mesh and an embedded refined mesh u...
Definition transfer.hpp:171
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Mesh data type.
Definition mesh.hpp:56
void Clear()
Clear the contents of the Mesh.
Definition mesh.hpp:730
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
static bool Root()
Return true if the rank in MPI_COMM_WORLD is zero.
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
Abstract operator.
Definition operator.hpp:25
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition operator.hpp:93
void ParseCheck(std::ostream &out=mfem::out)
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 'var' to receive the value. Enable/disable tags are used to set the bool...
Definition optparser.hpp:82
Class for parallel bilinear form.
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void Assemble(int skip_zeros=1)
Assemble the local matrix.
Abstract parallel finite element space.
Definition pfespace.hpp:29
const Operator * GetRestrictionOperator() const override
int GetTrueVSize() const override
Return the number of local vector true dofs.
Definition pfespace.hpp:289
Class for parallel grid function.
Definition pgridfunc.hpp:33
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
void ProjectCoefficient(Coefficient &coeff) override
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Class for parallel linear form.
void Assemble()
Assembles the ParLinearForm i.e. sums over all domain/bdr integrators.
Class for parallel meshes.
Definition pmesh.hpp:34
static ParMesh MakeRefined(ParMesh &orig_mesh, int ref_factor, int ref_type)
Create a uniformly refined (by any factor) version of orig_mesh.
Definition pmesh.cpp:1375
Vector data type.
Definition vector.hpp:80
real_t Norml2() const
Returns the l2 norm of the vector.
Definition vector.cpp:832
Data collection with VisIt I/O routines.
virtual void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
int dim
Definition ex24.cpp:53
int main()
const int visport
float real_t
Definition config.hpp:43
const char vishost[]
int Ww
int Wy
void visualize(VisItDataCollection &, string, int, int)
int Wx
real_t compute_mass(ParFiniteElementSpace *, real_t, VisItDataCollection &, string)
int problem
int offx
int Wh
string direction
real_t RHO_exact(const Vector &x)
string space
int offy
Helper struct to convert a C++ type to an MPI type.