MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
lor_solvers.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// Low-Order Refined Solvers Miniapp
14// ---------------------------------
15//
16// This miniapp illustrates the use of low-order refined preconditioners for
17// finite element problems defined using H1, H(curl), H(div), or L2 finite
18// element spaces. The following problems are solved, depending on the chosen
19// finite element space:
20//
21// H1 and L2: definite Helmholtz problem, u - Delta u = f
22// (in L2 discretized using the symmetric interior penalty DG method)
23//
24// H(curl): definite Maxwell problem, u + curl curl u = f
25//
26// H(div): grad-div problem, u - grad(div u) = f
27//
28// In each case, the high-order finite element problem is preconditioned using a
29// low-order finite element discretization defined on a Gauss-Lobatto refined
30// mesh. The low-order problem is solved using a direct solver if MFEM is
31// compiled with SuiteSparse enabled, or with one iteration of symmetric
32// Gauss-Seidel smoothing otherwise.
33//
34// For vector finite element spaces, the special "Integrated" basis type is used
35// to obtain spectral equivalence between the high-order and low-order refined
36// discretizations. This basis is defined in reference [1] and spectral
37// equivalence is shown in [2]:
38//
39// [1]. M. Gerritsma. Edge functions for spectral element methods. Spectral and
40// High Order Methods for Partial Differential Equations. (2010)
41// [2]. C. Dohrmann. Spectral equivalence properties of higher-order tensor
42// product finite elements and applications to preconditioning. (2021)
43//
44// The action of the high-order operator is computed using MFEM's partial
45// assembly/matrix-free algorithms (except in the case of L2, which remains
46// future work).
47//
48// Compile with: make lor_solvers
49//
50// Sample runs:
51//
52// lor_solvers -fe h
53// lor_solvers -fe n
54// lor_solvers -fe r
55// lor_solvers -fe l
56// lor_solvers -m ../../data/amr-quad.mesh -fe h
57// lor_solvers -m ../../data/amr-quad.mesh -fe n
58// lor_solvers -m ../../data/amr-quad.mesh -fe r
59// lor_solvers -m ../../data/amr-quad.mesh -fe l
60// lor_solvers -m ../../data/star-surf.mesh -fe h
61// lor_solvers -m ../../data/star-surf.mesh -fe n
62// lor_solvers -m ../../data/star-surf.mesh -fe r
63//
64// Device sample runs:
65// lor_solvers -fe h -d cuda
66// lor_solvers -fe n -d cuda
67// lor_solvers -fe r -d cuda
68// lor_solvers -fe l -d cuda
69
70#include "mfem.hpp"
71#include <fstream>
72#include <iostream>
73#include <memory>
74
75#include "lor_mms.hpp"
76
77using namespace std;
78using namespace mfem;
79
80int main(int argc, char *argv[])
81{
82 const char *mesh_file = "../../data/star.mesh";
83 int ref_levels = 1;
84 int order = 3;
85 const char *fe = "h";
86 const char *device_config = "cpu";
87 bool visualization = true;
88
89 OptionsParser args(argc, argv);
90 args.AddOption(&mesh_file, "-m", "--mesh", "Mesh file to use.");
91 args.AddOption(&ref_levels, "-r", "--refine",
92 "Number of times to refine the mesh uniformly.");
93 args.AddOption(&order, "-o", "--order", "Polynomial degree.");
94 args.AddOption(&fe, "-fe", "--fe-type",
95 "FE type. h for H1, n for Hcurl, r for Hdiv, l for L2");
96 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
97 "--no-visualization",
98 "Enable or disable GLVis visualization.");
99 args.AddOption(&device_config, "-d", "--device",
100 "Device configuration string, see Device::Configure().");
101 args.ParseCheck();
102
103 Device device(device_config);
104 device.Print();
105
106 bool H1 = false, ND = false, RT = false, L2 = false;
107 if (string(fe) == "h") { H1 = true; }
108 else if (string(fe) == "n") { ND = true; }
109 else if (string(fe) == "r") { RT = true; }
110 else if (string(fe) == "l") { L2 = true; }
111 else { MFEM_ABORT("Bad FE type. Must be 'h', 'n', 'r', or 'l'."); }
112
113 real_t kappa = (order+1)*(order+1); // Penalty used for DG discretizations
114
115 Mesh mesh(mesh_file, 1, 1);
116 const int dim = mesh.Dimension();
117 const int sdim = mesh.SpaceDimension();
118 MFEM_VERIFY(dim == 2 || dim == 3, "Mesh dimension must be 2 or 3.");
119 MFEM_VERIFY(!L2 || dim == sdim, "DG surface meshes not supported.");
120 for (int l = 0; l < ref_levels; l++) { mesh.UniformRefinement(); }
121
122 FunctionCoefficient f_coeff(f(1.0)), u_coeff(u);
123 VectorFunctionCoefficient f_vec_coeff(sdim, f_vec(RT)),
124 u_vec_coeff(sdim, u_vec);
125
127 unique_ptr<FiniteElementCollection> fec;
128 if (H1) { fec.reset(new H1_FECollection(order, dim, b1)); }
129 else if (ND) { fec.reset(new ND_FECollection(order, dim, b1, b2)); }
130 else if (RT) { fec.reset(new RT_FECollection(order-1, dim, b1, b2)); }
131 else { fec.reset(new L2_FECollection(order, dim, b1)); }
132
133 FiniteElementSpace fes(&mesh, fec.get());
134 cout << "Number of DOFs: " << fes.GetTrueVSize() << endl;
135
136 Array<int> ess_dofs;
137 // In DG, boundary conditions are enforced weakly, so no essential DOFs.
138 if (!L2) { fes.GetBoundaryTrueDofs(ess_dofs); }
139
140 BilinearForm a(&fes);
141 if (H1 || L2)
142 {
143 a.AddDomainIntegrator(new MassIntegrator);
144 a.AddDomainIntegrator(new DiffusionIntegrator);
145 }
146 else
147 {
148 a.AddDomainIntegrator(new VectorFEMassIntegrator);
149 }
150
151 if (ND) { a.AddDomainIntegrator(new CurlCurlIntegrator); }
152 else if (RT) { a.AddDomainIntegrator(new DivDivIntegrator); }
153 else if (L2)
154 {
155 a.AddInteriorFaceIntegrator(new DGDiffusionIntegrator(-1.0, kappa));
156 a.AddBdrFaceIntegrator(new DGDiffusionIntegrator(-1.0, kappa));
157 }
158 // Partial assembly not currently supported for DG or for surface meshes with
159 // vector finite elements (ND or RT).
160 if (!L2 && (H1 || sdim == dim)) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
161 a.Assemble();
162
163 LinearForm b(&fes);
164 if (H1 || L2) { b.AddDomainIntegrator(new DomainLFIntegrator(f_coeff)); }
165 else { b.AddDomainIntegrator(new VectorFEDomainLFIntegrator(f_vec_coeff)); }
166 if (L2)
167 {
168 // DG boundary conditions are enforced weakly with this integrator.
169 b.AddBdrFaceIntegrator(new DGDirichletLFIntegrator(u_coeff, -1.0, kappa));
170 }
171 b.Assemble();
172
173 GridFunction x(&fes);
174 if (H1 || L2) { x.ProjectCoefficient(u_coeff);}
175 else { x.ProjectCoefficient(u_vec_coeff); }
176
177 Vector X, B;
179 a.FormLinearSystem(ess_dofs, x, b, A, X, B);
180
181#ifdef MFEM_USE_SUITESPARSE
182 LORSolver<UMFPackSolver> solv_lor(a, ess_dofs);
183#else
184 LORSolver<GSSmoother> solv_lor(a, ess_dofs);
185#endif
186
187 CGSolver cg;
188 cg.SetAbsTol(0.0);
189 cg.SetRelTol(1e-12);
190 cg.SetMaxIter(500);
191 cg.SetPrintLevel(1);
192 cg.SetOperator(*A);
193 cg.SetPreconditioner(solv_lor);
194 cg.Mult(B, X);
195
196 a.RecoverFEMSolution(X, b, x);
197
198 if (sdim == dim)
199 {
200 real_t er =
201 (H1 || L2) ? x.ComputeL2Error(u_coeff) : x.ComputeL2Error(u_vec_coeff);
202 cout << "L2 error: " << er << endl;
203 }
204
205 if (visualization)
206 {
207 // Save the solution and mesh to disk. The output can be viewed using
208 // GLVis as follows: "glvis -m mesh.mesh -g sol.gf"
209 x.Save("sol.gf");
210 mesh.Save("mesh.mesh");
211
212 // Also save the solution for visualization using ParaView
213 ParaViewDataCollection dc("LOR", &mesh);
214 dc.SetPrefixPath("ParaView");
215 dc.SetHighOrderOutput(true);
216 dc.SetLevelsOfDetail(order);
217 dc.RegisterField("u", &x);
218 dc.SetCycle(0);
219 dc.SetTime(0.0);
220 dc.Save();
221 }
222
223 return 0;
224}
@ GaussLobatto
Closed type.
Definition fe_base.hpp:32
@ IntegratedGLL
Integrated GLL indicator functions.
Definition fe_base.hpp:39
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
Conjugate gradient method.
Definition solvers.hpp:513
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition solvers.hpp:526
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the Conjugate Gradient method.
Definition solvers.cpp:718
Integrator for for Nedelec elements.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetTime(real_t t)
Set physical time (for time-dependent simulations)
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition device.hpp:123
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition device.cpp:286
for Raviart-Thomas elements
Class for domain integration .
Definition lininteg.hpp:109
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:716
void GetBoundaryTrueDofs(Array< int > &boundary_dofs, int component=-1)
Get a list of all boundary true dofs, boundary_dofs. For spaces with 'vdim' > 1, the 'component' para...
Definition fespace.cpp:632
A general function coefficient.
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
virtual real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
void SetRelTol(real_t rtol)
Definition solvers.hpp:209
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition solvers.cpp:173
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition solvers.cpp:71
void SetMaxIter(int max_it)
Definition solvers.hpp:211
void SetAbsTol(real_t atol)
Definition solvers.hpp:210
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
Represents a solver of type SolverType created using the low-order refined version of the given Bilin...
Definition lor.hpp:204
Vector with associated FE space and LinearFormIntegrators.
Mesh data type.
Definition mesh.hpp:56
virtual void Save(const std::string &fname, int precision=16) const
Definition mesh.cpp:11481
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1163
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition fe_coll.hpp:465
Pointer to an Operator of a specified type.
Definition handle.hpp:34
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
Helper class for ParaView visualization data.
void SetHighOrderOutput(bool high_order_output_)
void SetLevelsOfDetail(int levels_of_detail_)
virtual void Save() override
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition fe_coll.hpp:386
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition lininteg.hpp:347
A general vector function coefficient.
Vector data type.
Definition vector.hpp:80
real_t kappa
Definition ex24.cpp:54
int dim
Definition ex24.cpp:53
int main()
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
void u_vec(const Vector &xvec, Vector &u)
Definition lor_mms.hpp:50
std::function< void(const Vector &, Vector &)> f_vec(bool grad_div_problem)
Definition lor_mms.hpp:68
float real_t
Definition config.hpp:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30