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