MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex26.cpp
Go to the documentation of this file.
1// MFEM Example 26
2//
3// Compile with: make ex26
4//
5// Sample runs: ex26 -m ../data/star.mesh
6// ex26 -m ../data/fichera.mesh
7// ex26 -m ../data/beam-hex.mesh
8//
9// Device sample runs:
10// ex26 -d cuda
11// ex26 -d raja-cuda
12// ex26 -d occa-cuda
13// ex26 -d raja-omp
14// ex26 -d occa-omp
15// ex26 -d ceed-cpu
16// ex26 -d ceed-cuda
17// ex26 -m ../data/beam-hex.mesh -d cuda
18//
19// Description: This example code demonstrates the use of MFEM to define a
20// simple finite element discretization of the Laplace problem
21// -Delta u = 1 with homogeneous Dirichlet boundary conditions
22// as in Example 1.
23//
24// It highlights on the creation of a hierarchy of discretization
25// spaces with partial assembly and the construction of an
26// efficient multigrid preconditioner for the iterative solver.
27//
28// We recommend viewing Example 1 before viewing this example.
29
30#include "mfem.hpp"
31#include <fstream>
32#include <iostream>
33
34using namespace std;
35using namespace mfem;
36
37// Class for constructing a multigrid preconditioner for the diffusion operator.
38// This example multigrid preconditioner class demonstrates the creation of the
39// diffusion bilinear forms and operators using partial assembly for all spaces
40// in the FiniteElementSpaceHierarchy. The preconditioner uses a CG solver on
41// the coarsest level and second order Chebyshev accelerated smoothers on the
42// other levels.
43class DiffusionMultigrid : public GeometricMultigrid
44{
45private:
47
48public:
49 // Constructs a diffusion multigrid for the given FiniteElementSpaceHierarchy
50 // and the array of essential boundaries
51 DiffusionMultigrid(FiniteElementSpaceHierarchy& fespaces, Array<int>& ess_bdr)
52 : GeometricMultigrid(fespaces, ess_bdr), coeff(1.0)
53 {
54 ConstructCoarseOperatorAndSolver(fespaces.GetFESpaceAtLevel(0));
55 for (int level = 1; level < fespaces.GetNumLevels(); ++level)
56 {
57 ConstructOperatorAndSmoother(fespaces.GetFESpaceAtLevel(level), level);
58 }
59 }
60
61private:
62 void ConstructBilinearForm(FiniteElementSpace& fespace)
63 {
64 BilinearForm* form = new BilinearForm(&fespace);
65 form->SetAssemblyLevel(AssemblyLevel::PARTIAL);
67 form->Assemble();
68 bfs.Append(form);
69 }
70
71 void ConstructCoarseOperatorAndSolver(FiniteElementSpace& coarse_fespace)
72 {
73 ConstructBilinearForm(coarse_fespace);
74
75 OperatorPtr opr;
77 bfs[0]->FormSystemMatrix(*essentialTrueDofs[0], opr);
78 opr.SetOperatorOwner(false);
79
80 CGSolver* pcg = new CGSolver();
81 pcg->SetPrintLevel(-1);
82 pcg->SetMaxIter(200);
83 pcg->SetRelTol(sqrt(1e-4));
84 pcg->SetAbsTol(0.0);
85 pcg->SetOperator(*opr.Ptr());
86
87 AddLevel(opr.Ptr(), pcg, true, true);
88 }
89
90 void ConstructOperatorAndSmoother(FiniteElementSpace& fespace, int level)
91 {
92 const Array<int> &ess_tdof_list = *essentialTrueDofs[level];
93 ConstructBilinearForm(fespace);
94
95 OperatorPtr opr;
97 bfs[level]->FormSystemMatrix(ess_tdof_list, opr);
98 opr.SetOperatorOwner(false);
99
100 Vector diag(fespace.GetTrueVSize());
101 bfs[level]->AssembleDiagonal(diag);
102
103 Solver* smoother = new OperatorChebyshevSmoother(*opr, diag, ess_tdof_list, 2);
104 AddLevel(opr.Ptr(), smoother, true, true);
105 }
106};
107
108
109int main(int argc, char *argv[])
110{
111 // 1. Parse command-line options.
112 const char *mesh_file = "../data/star.mesh";
113 int geometric_refinements = 0;
114 int order_refinements = 2;
115 const char *device_config = "cpu";
116 bool visualization = true;
117
118 OptionsParser args(argc, argv);
119 args.AddOption(&mesh_file, "-m", "--mesh",
120 "Mesh file to use.");
121 args.AddOption(&geometric_refinements, "-gr", "--geometric-refinements",
122 "Number of geometric refinements done prior to order refinements.");
123 args.AddOption(&order_refinements, "-or", "--order-refinements",
124 "Number of order refinements. Finest level in the hierarchy has order 2^{or}.");
125 args.AddOption(&device_config, "-d", "--device",
126 "Device configuration string, see Device::Configure().");
127 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
128 "--no-visualization",
129 "Enable or disable GLVis visualization.");
130 args.Parse();
131 if (!args.Good())
132 {
133 args.PrintUsage(cout);
134 return 1;
135 }
136 args.PrintOptions(cout);
137
138 // 2. Enable hardware devices such as GPUs, and programming models such as
139 // CUDA, OCCA, RAJA and OpenMP based on command line options.
140 Device device(device_config);
141 device.Print();
142
143 // 3. Read the mesh from the given mesh file. We can handle triangular,
144 // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
145 // the same code.
146 Mesh *mesh = new Mesh(mesh_file, 1, 1);
147 int dim = mesh->Dimension();
148
149 // 4. Refine the mesh to increase the resolution. In this example we do
150 // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
151 // largest number that gives a final mesh with no more than 50,000
152 // elements.
153 {
154 int ref_levels =
155 (int)floor(log(5000./mesh->GetNE())/log(2.)/dim);
156 for (int l = 0; l < ref_levels; l++)
157 {
158 mesh->UniformRefinement();
159 }
160 }
161
162 // 5. Define a finite element space hierarchy on the mesh. Here we use
163 // continuous Lagrange finite elements. We start with order 1 on the
164 // coarse level and geometrically refine the spaces by the specified
165 // amount. Afterwards, we increase the order of the finite elements
166 // by a factor of 2 for each additional level.
168 FiniteElementSpace *coarse_fespace = new FiniteElementSpace(mesh, fec);
169 FiniteElementSpaceHierarchy fespaces(mesh, coarse_fespace, true, true);
170
172 collections.Append(fec);
173 for (int level = 0; level < geometric_refinements; ++level)
174 {
175 fespaces.AddUniformlyRefinedLevel();
176 }
177 for (int level = 0; level < order_refinements; ++level)
178 {
179 collections.Append(new H1_FECollection((int)std::pow(2, level+1), dim));
180 fespaces.AddOrderRefinedLevel(collections.Last());
181 }
182
183 cout << "Number of finite element unknowns: "
184 << fespaces.GetFinestFESpace().GetTrueVSize() << endl;
185
186 // 6. Set up the linear form b(.) which corresponds to the right-hand side of
187 // the FEM linear system, which in this case is (1,phi_i) where phi_i are
188 // the basis functions in the finite element fespace.
189 LinearForm *b = new LinearForm(&fespaces.GetFinestFESpace());
190 ConstantCoefficient one(1.0);
191 b->AddDomainIntegrator(new DomainLFIntegrator(one));
192 b->Assemble();
193
194 // 7. Define the solution vector x as a finite element grid function
195 // corresponding to fespace. Initialize x with initial guess of zero,
196 // which satisfies the boundary conditions.
197 GridFunction x(&fespaces.GetFinestFESpace());
198 x = 0.0;
199
200 // 8. Create the multigrid operator using the previously created
201 // FiniteElementSpaceHierarchy and additional boundary information. This
202 // operator is then used to create the MultigridSolver as a preconditioner
203 // in the iterative solver.
204 Array<int> ess_bdr(mesh->bdr_attributes.Max());
205 ess_bdr = 1;
206
207 DiffusionMultigrid M(fespaces, ess_bdr);
208 M.SetCycleType(Multigrid::CycleType::VCYCLE, 1, 1);
209
210 OperatorPtr A;
211 Vector B, X;
212 M.FormFineLinearSystem(x, *b, A, X, B);
213 cout << "Size of linear system: " << A->Height() << endl;
214
215 // 9. Solve the linear system A X = B.
216 PCG(*A, M, B, X, 1, 2000, 1e-12, 0.0);
217
218 // 10. Recover the solution as a finite element grid function.
219 M.RecoverFineFEMSolution(X, *b, x);
220
221 // 11. Save the refined mesh and the solution. This output can be viewed
222 // later using GLVis: "glvis -m refined.mesh -g sol.gf".
223 ofstream mesh_ofs("refined.mesh");
224 mesh_ofs.precision(8);
225 fespaces.GetFinestFESpace().GetMesh()->Print(mesh_ofs);
226 ofstream sol_ofs("sol.gf");
227 sol_ofs.precision(8);
228 x.Save(sol_ofs);
229
230 // 12. Send the solution by socket to a GLVis server.
231 if (visualization)
232 {
233 char vishost[] = "localhost";
234 int visport = 19916;
235 socketstream sol_sock(vishost, visport);
236 sol_sock.precision(8);
237 sol_sock << "solution\n" << *fespaces.GetFinestFESpace().GetMesh() << x <<
238 flush;
239 }
240
241 // 13. Free the used memory.
242 delete b;
243 for (int level = 0; level < collections.Size(); ++level)
244 {
245 delete collections[level];
246 }
247
248 return 0;
249}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
int Size() const
Return the logical size of the array.
Definition array.hpp:144
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:769
T & Last()
Return the last element in the array.
Definition array.hpp:802
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
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
A coefficient that is constant across space and time.
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
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
int GetNumLevels() const
Returns the number of levels in the hierarchy.
virtual const FiniteElementSpace & GetFinestFESpace() const
Returns the finite element space at the finest level.
virtual void AddUniformlyRefinedLevel(int dim=1, int ordering=Ordering::byVDIM)
Adds one level to the hierarchy by uniformly refining the mesh on the previous level.
virtual const FiniteElementSpace & GetFESpaceAtLevel(int level) const
Returns the finite element space at the given level.
virtual void AddOrderRefinedLevel(FiniteElementCollection *fec, int dim=1, int ordering=Ordering::byVDIM)
Adds one level to the hierarchy by using a different finite element order defined through FiniteEleme...
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
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:559
Geometric multigrid associated with a hierarchy of finite element spaces.
void RecoverFineFEMSolution(const Vector &X, const Vector &b, Vector &x)
Recover the solution of a linear system formed with FormFineLinearSystem()
Array< Array< int > * > essentialTrueDofs
Array< BilinearForm * > bfs
const FiniteElementSpaceHierarchy & fespaces
void FormFineLinearSystem(Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)
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.
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
void SetRelTol(real_t rtol)
Definition solvers.hpp:209
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
Vector with associated FE space and LinearFormIntegrators.
Mesh data type.
Definition mesh.hpp:56
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:282
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Definition mesh.hpp:2288
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
void AddLevel(Operator *op, Solver *smoother, bool ownOperator, bool ownSmoother)
Adds a level to the multigrid operator hierarchy.
Definition multigrid.cpp:87
void SetCycleType(CycleType cycleType_, int preSmoothingSteps_, int postSmoothingSteps_)
Set cycle type and number of pre- and post-smoothing steps used by Mult.
Definition multigrid.cpp:98
Chebyshev accelerated smoothing with given vector, no matrix necessary.
Definition solvers.hpp:394
Pointer to an Operator of a specified type.
Definition handle.hpp:34
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
Definition handle.hpp:120
void SetType(Operator::Type tid)
Invoke Clear() and set a new type id.
Definition handle.hpp:132
Operator * Ptr() const
Access the underlying Operator pointer.
Definition handle.hpp:87
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
@ ANY_TYPE
ID for the base class Operator, i.e. any type.
Definition operator.hpp:285
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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
bool Good() const
Return true if the command line options were parsed successfully.
Base class for solvers.
Definition operator.hpp:683
Vector data type.
Definition vector.hpp:80
int dim
Definition ex24.cpp:53
int main()
real_t b
Definition lissajous.cpp:42
const int visport
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
Definition solvers.cpp:913
const char vishost[]