MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
ex1.cpp
Go to the documentation of this file.
1// MFEM Example 1
2//
3// Compile with: make ex1
4//
5// Sample runs: ex1 -m ../data/square-disc.mesh
6// ex1 -m ../data/star.mesh
7// ex1 -m ../data/star-mixed.mesh
8// ex1 -m ../data/escher.mesh
9// ex1 -m ../data/fichera.mesh
10// ex1 -m ../data/fichera-mixed.mesh
11// ex1 -m ../data/toroid-wedge.mesh
12// ex1 -m ../data/octahedron.mesh -o 1
13// ex1 -m ../data/periodic-annulus-sector.msh
14// ex1 -m ../data/periodic-torus-sector.msh
15// ex1 -m ../data/square-disc-p2.vtk -o 2
16// ex1 -m ../data/square-disc-p3.mesh -o 3
17// ex1 -m ../data/square-disc-nurbs.mesh -o -1
18// ex1 -m ../data/star-mixed-p2.mesh -o 2
19// ex1 -m ../data/disc-nurbs.mesh -o -1
20// ex1 -m ../data/pipe-nurbs.mesh -o -1
21// ex1 -m ../data/fichera-mixed-p2.mesh -o 2
22// ex1 -m ../data/star-surf.mesh
23// ex1 -m ../data/square-disc-surf.mesh
24// ex1 -m ../data/inline-segment.mesh
25// ex1 -m ../data/amr-quad.mesh
26// ex1 -m ../data/amr-hex.mesh
27// ex1 -m ../data/fichera-amr.mesh
28// ex1 -m ../data/mobius-strip.mesh
29// ex1 -m ../data/mobius-strip.mesh -o -1 -sc
30// ex1 -m ../data/nc3-nurbs.mesh -o -1
31//
32// Device sample runs:
33// ex1 -pa -d cuda
34// ex1 -fa -d cuda
35// ex1 -pa -d raja-cuda
36// * ex1 -pa -d raja-hip
37// ex1 -pa -d occa-cuda
38// ex1 -pa -d raja-omp
39// ex1 -pa -d occa-omp
40// ex1 -pa -d ceed-cpu
41// ex1 -pa -d ceed-cpu -o 4 -a
42// ex1 -pa -d ceed-cpu -m ../data/square-mixed.mesh
43// ex1 -pa -d ceed-cpu -m ../data/fichera-mixed.mesh
44// * ex1 -pa -d ceed-cuda
45// * ex1 -pa -d ceed-hip
46// ex1 -pa -d ceed-cuda:/gpu/cuda/shared
47// ex1 -pa -d ceed-cuda:/gpu/cuda/shared -m ../data/square-mixed.mesh
48// ex1 -pa -d ceed-cuda:/gpu/cuda/shared -m ../data/fichera-mixed.mesh
49// ex1 -m ../data/beam-hex.mesh -pa -d cuda
50// ex1 -m ../data/beam-tet.mesh -pa -d ceed-cpu
51// ex1 -m ../data/beam-tet.mesh -pa -d ceed-cuda:/gpu/cuda/ref
52//
53// Description: This example code demonstrates the use of MFEM to define a
54// simple finite element discretization of the Poisson problem
55// -Delta u = 1 with homogeneous Dirichlet boundary conditions.
56// Specifically, we discretize using a FE space of the specified
57// order, or if order < 1 using an isoparametric/isogeometric
58// space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
59// NURBS mesh, etc.)
60//
61// The example highlights the use of mesh refinement, finite
62// element grid functions, as well as linear and bilinear forms
63// corresponding to the left-hand side and right-hand side of the
64// discrete linear system. We also cover the explicit elimination
65// of essential boundary conditions, static condensation, and the
66// optional connection to the GLVis tool for visualization.
67
68#include "mfem.hpp"
69#include <fstream>
70#include <iostream>
71
72using namespace std;
73using namespace mfem;
74
75int main(int argc, char *argv[])
76{
77 // 1. Parse command-line options.
78 const char *mesh_file = "../data/star.mesh";
79 int order = 1;
80 bool static_cond = false;
81 bool pa = false;
82 bool fa = false;
83 const char *device_config = "cpu";
84 bool visualization = true;
85 bool algebraic_ceed = false;
86
87 OptionsParser args(argc, argv);
88 args.AddOption(&mesh_file, "-m", "--mesh",
89 "Mesh file to use.");
90 args.AddOption(&order, "-o", "--order",
91 "Finite element order (polynomial degree) or -1 for"
92 " isoparametric space.");
93 args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
94 "--no-static-condensation", "Enable static condensation.");
95 args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
96 "--no-partial-assembly", "Enable Partial Assembly.");
97 args.AddOption(&fa, "-fa", "--full-assembly", "-no-fa",
98 "--no-full-assembly", "Enable Full Assembly.");
99 args.AddOption(&device_config, "-d", "--device",
100 "Device configuration string, see Device::Configure().");
101#ifdef MFEM_USE_CEED
102 args.AddOption(&algebraic_ceed, "-a", "--algebraic", "-no-a", "--no-algebraic",
103 "Use algebraic Ceed solver");
104#endif
105 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
106 "--no-visualization",
107 "Enable or disable GLVis visualization.");
108 args.Parse();
109 if (!args.Good())
110 {
111 args.PrintUsage(cout);
112 return 1;
113 }
114 args.PrintOptions(cout);
115
116 // 2. Enable hardware devices such as GPUs, and programming models such as
117 // CUDA, OCCA, RAJA and OpenMP based on command line options.
118 Device device(device_config);
119 device.Print();
120
121 // 3. Read the mesh from the given mesh file. We can handle triangular,
122 // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
123 // the same code.
124 Mesh mesh(mesh_file, 1, 1);
125 int dim = mesh.Dimension();
126
127 // 4. Refine the mesh to increase the resolution. In this example we do
128 // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
129 // largest number that gives a final mesh with no more than 50,000
130 // elements.
131 {
132 int ref_levels =
133 (int)floor(log(50000./mesh.GetNE())/log(2.)/dim);
134 for (int l = 0; l < ref_levels; l++)
135 {
136 mesh.UniformRefinement();
137 }
138 }
139
140 // 5. Define a finite element space on the mesh. Here we use continuous
141 // Lagrange finite elements of the specified order. If order < 1, we
142 // instead use an isoparametric/isogeometric space.
144 bool delete_fec;
145 if (order > 0)
146 {
147 fec = new H1_FECollection(order, dim);
148 delete_fec = true;
149 }
150 else if (mesh.GetNodes())
151 {
152 fec = mesh.GetNodes()->OwnFEC();
153 delete_fec = false;
154 cout << "Using isoparametric FEs: " << fec->Name() << endl;
155 }
156 else
157 {
158 fec = new H1_FECollection(order = 1, dim);
159 delete_fec = true;
160 }
161 FiniteElementSpace fespace(&mesh, fec);
162 cout << "Number of finite element unknowns: "
163 << fespace.GetTrueVSize() << endl;
164
165 // 6. Determine the list of true (i.e. conforming) essential boundary dofs.
166 // In this example, the boundary conditions are defined by marking all
167 // the external boundary attributes from the mesh as essential (Dirichlet)
168 // and converting them to a list of true dofs.
169 Array<int> ess_tdof_list;
170 if (mesh.bdr_attributes.Size())
171 {
172 Array<int> ess_bdr(mesh.bdr_attributes.Max());
173 ess_bdr = 0;
174 // Apply boundary conditions on all external boundaries:
175 mesh.MarkExternalBoundaries(ess_bdr);
176 // Boundary conditions can also be applied based on named attributes:
177 // mesh.MarkNamedBoundaries(set_name, ess_bdr)
178
179 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
180 }
181
182 // 7. Set up the linear form b(.) which corresponds to the right-hand side of
183 // the FEM linear system, which in this case is (1,phi_i) where phi_i are
184 // the basis functions in the finite element fespace.
185 LinearForm b(&fespace);
186 ConstantCoefficient one(1.0);
187 b.AddDomainIntegrator(new DomainLFIntegrator(one));
188 b.Assemble();
189
190 // 8. Define the solution vector x as a finite element grid function
191 // corresponding to fespace. Initialize x with initial guess of zero,
192 // which satisfies the boundary conditions.
193 GridFunction x(&fespace);
194 x = 0.0;
195
196 // 9. Set up the bilinear form a(.,.) on the finite element space
197 // corresponding to the Laplacian operator -Delta, by adding the Diffusion
198 // domain integrator.
199 BilinearForm a(&fespace);
200 if (pa) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
201 if (fa)
202 {
203 a.SetAssemblyLevel(AssemblyLevel::FULL);
204 // Sort the matrix column indices when running on GPU or with OpenMP (i.e.
205 // when Device::IsEnabled() returns true). This makes the results
206 // bit-for-bit deterministic at the cost of somewhat longer run time.
207 a.EnableSparseMatrixSorting(Device::IsEnabled());
208 }
209 a.AddDomainIntegrator(new DiffusionIntegrator(one));
210
211 // 10. Assemble the bilinear form and the corresponding linear system,
212 // applying any necessary transformations such as: eliminating boundary
213 // conditions, applying conforming constraints for non-conforming AMR,
214 // static condensation, etc.
215 if (static_cond) { a.EnableStaticCondensation(); }
216 a.Assemble();
217
218 OperatorPtr A;
219 Vector B, X;
220 a.FormLinearSystem(ess_tdof_list, x, b, A, X, B);
221
222 cout << "Size of linear system: " << A->Height() << endl;
223
224 // 11. Solve the linear system A X = B.
225 if (!pa)
226 {
227#ifndef MFEM_USE_SUITESPARSE
228 // Use a simple symmetric Gauss-Seidel preconditioner with PCG.
229 GSSmoother M((SparseMatrix&)(*A));
230 PCG(*A, M, B, X, 1, 200, 1e-12, 0.0);
231#else
232 // If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
233 UMFPackSolver umf_solver;
234 umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
235 umf_solver.SetOperator(*A);
236 umf_solver.Mult(B, X);
237#endif
238 }
239 else
240 {
241 if (UsesTensorBasis(fespace))
242 {
243 if (algebraic_ceed)
244 {
245 ceed::AlgebraicSolver M(a, ess_tdof_list);
246 PCG(*A, M, B, X, 1, 400, 1e-12, 0.0);
247 }
248 else
249 {
250 OperatorJacobiSmoother M(a, ess_tdof_list);
251 PCG(*A, M, B, X, 1, 400, 1e-12, 0.0);
252 }
253 }
254 else
255 {
256 CG(*A, B, X, 1, 400, 1e-12, 0.0);
257 }
258 }
259
260 // 12. Recover the solution as a finite element grid function.
261 a.RecoverFEMSolution(X, b, x);
262
263 // 13. Save the refined mesh and the solution. This output can be viewed later
264 // using GLVis: "glvis -m refined.mesh -g sol.gf".
265 ofstream mesh_ofs("refined.mesh");
266 mesh_ofs.precision(8);
267 mesh.Print(mesh_ofs);
268 ofstream sol_ofs("sol.gf");
269 sol_ofs.precision(8);
270 x.Save(sol_ofs);
271
272 // 14. Send the solution by socket to a GLVis server.
273 if (visualization)
274 {
275 char vishost[] = "localhost";
276 int visport = 19916;
277 socketstream sol_sock(vishost, visport);
278 sol_sock.precision(8);
279 sol_sock << "solution\n" << mesh << x << flush;
280 }
281
282 // 15. Free the used memory.
283 if (delete_fec)
284 {
285 delete fec;
286 }
287
288 return 0;
289}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:69
int Size() const
Return the logical size of the array.
Definition array.hpp:166
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
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:124
void Print(std::ostream &os=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition device.cpp:317
static bool IsEnabled()
Return true if any backend other than Backend::CPU is enabled.
Definition device.hpp:247
Class for domain integration .
Definition lininteg.hpp:108
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
virtual const char * Name() const
Definition fe_coll.hpp:79
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:208
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:829
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition fespace.cpp:628
Data type for Gauss-Seidel smoother of sparse matrix.
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:32
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:279
Vector with associated FE space and LinearFormIntegrators.
Mesh data type.
Definition mesh.hpp:65
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:304
virtual void MarkExternalBoundaries(Array< int > &bdr_marker, bool excl=true) const
Mark boundary attributes of external boundaries.
Definition mesh.cpp:1741
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Print the mesh to the given stream using the default MFEM mesh format.
Definition mesh.hpp:2567
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1377
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1306
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:9627
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:11637
Pointer to an Operator of a specified type.
Definition handle.hpp:34
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition solvers.hpp:422
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
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.
Data type sparse matrix.
Definition sparsemat.hpp:51
Direct sparse solver using UMFPACK.
Definition solvers.hpp:1210
real_t Control[UMFPACK_CONTROL]
Definition solvers.hpp:1220
void SetOperator(const Operator &op) override
Factorize the given Operator op which must be a SparseMatrix.
Definition solvers.cpp:3368
void Mult(const Vector &b, Vector &x) const override
Direct solution of the linear system using UMFPACK.
Definition solvers.cpp:3463
Vector data type.
Definition vector.hpp:82
Wrapper for AlgebraicMultigrid object.
int dim
Definition ex24.cpp:53
int main()
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
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:1067
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
Definition solvers.cpp:1052
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
Definition fespace.hpp:1548
const char vishost[]