MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex39.cpp
Go to the documentation of this file.
1// MFEM Example 39
2//
3// Compile with: make ex39
4//
5// Sample runs: ex39
6// ex39 -ess "Southern Boundary"
7// ex39 -src Base
8//
9// Description: This example code demonstrates the use of named attribute
10// sets in MFEM to specify material regions, boundary regions,
11// or source regions by name rather than attribute numbers. It
12// also demonstrates how new named attribute sets may be created
13// from arbitrary groupings of attribute numbers and used as a
14// convenient shorthand to refer to those groupings in other
15// portions of the application or through the command line.
16//
17// The particular problem being solved here is nearly the same
18// as that in example 1 i.e. a simple finite element
19// discretization of the Laplace problem -Delta u = 1 with
20// homogeneous Dirichlet boundary conditions and, in this case,
21// an inhomogeneous diffusion coefficient. The diffusion
22// coefficient is given a small default value throughout the
23// domain which is increased by two separate amounts in two named
24// regions.
25//
26// This example makes use of a specific input mesh, "compass.msh",
27// containing named domain and boundary regions generated by Gmsh
28// and stored in their "msh" format (version 2.2). This file
29// defines eight boundary regions corresponding to eight compass
30// headings; "ENE", "NNE", "NNW", "WSW", "SSW", "SSE", and "ESE".
31// It also defines nine domain regions; "Base", "N Even", "N Odd",
32// "W Even", "W Odd", "S Even", "S Odd", "E Even", and "E Odd".
33// These regions split the four compass pointers into two halves
34// each and also label the remaining elements as "Base". Starting
35// with these named regions we test the construction of named
36// sets as well as reading and writing these named groupings from
37// and to mesh files.
38//
39// The example highlights the use of named attribute sets for
40// both subdomains and boundaries in different contexts as well
41// as basic methods to create named sets from existing attributes.
42
43#include "mfem.hpp"
44#include <fstream>
45#include <iostream>
46
47using namespace std;
48using namespace mfem;
49
50int main(int argc, char *argv[])
51{
52 // 1. Parse command-line options.
53 const char *mesh_file = "../data/compass.msh";
54 int order = 1;
55 string source_name = "Rose Even";
56 string ess_name = "Boundary";
57 bool visualization = true;
58
59 OptionsParser args(argc, argv);
60 args.AddOption(&mesh_file, "-m", "--mesh",
61 "Mesh file to use.");
62 args.AddOption(&order, "-o", "--order",
63 "Finite element order (polynomial degree) or -1 for"
64 " isoparametric space.");
65 args.AddOption(&source_name,"-src","--source-attr-name",
66 "Name of attribute set containing source.");
67 args.AddOption(&ess_name,"-ess","--ess-attr-name",
68 "Name of attribute set containing essential BC.");
69 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
70 "--no-visualization",
71 "Enable or disable GLVis visualization.");
72 args.ParseCheck();
73
74 // 2. Read the mesh from the given mesh file. We can handle triangular,
75 // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
76 // the same code.
77 Mesh mesh(mesh_file, 1, 1);
78 int dim = mesh.Dimension();
79
80 // 3. Refine the mesh to increase the resolution. In this example we do
81 // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
82 // largest number that gives a final mesh with no more than 50,000
83 // elements.
84 {
85 int ref_levels =
86 (int)floor(log(50000./mesh.GetNE())/log(2.)/dim);
87 for (int l = 0; l < ref_levels; l++)
88 {
89 mesh.UniformRefinement();
90 }
91 }
92
93 // 4a. Display attribute set names contained in the initial mesh
94 AttributeSets &attr_sets = mesh.attribute_sets;
95 AttributeSets &bdr_attr_sets = mesh.bdr_attribute_sets;
96 {
97 std::set<string> names = attr_sets.GetAttributeSetNames();
98 cout << "Element Attribute Set Names: ";
99 for (auto const &set_name : names)
100 {
101 cout << " \"" << set_name << "\"";
102 }
103 cout << endl;
104
105 std::set<string> bdr_names = bdr_attr_sets.GetAttributeSetNames();
106 cout << "Boundary Attribute Set Names: ";
107 for (auto const &bdr_set_name : bdr_names)
108 {
109 cout << " \"" << bdr_set_name << "\"";
110 }
111 cout << endl;
112 }
113
114 // 4b. Define new regions based on existing attribute sets
115 {
116 Array<int> & Na = attr_sets.GetAttributeSet("N Even");
117 Array<int> & Nb = attr_sets.GetAttributeSet("N Odd");
118 Array<int> & Sa = attr_sets.GetAttributeSet("S Even");
119 Array<int> & Sb = attr_sets.GetAttributeSet("S Odd");
120 Array<int> & Ea = attr_sets.GetAttributeSet("E Even");
121 Array<int> & Eb = attr_sets.GetAttributeSet("E Odd");
122 Array<int> & Wa = attr_sets.GetAttributeSet("W Even");
123 Array<int> & Wb = attr_sets.GetAttributeSet("W Odd");
124
125 // Create a new set spanning the North point
126 attr_sets.SetAttributeSet("North", Na);
127 attr_sets.AddToAttributeSet("North", Nb);
128
129 // Create a new set spanning the South point
130 attr_sets.SetAttributeSet("South", Sa);
131 attr_sets.AddToAttributeSet("South", Sb);
132
133 // Create a new set spanning the East point
134 attr_sets.SetAttributeSet("East", Ea);
135 attr_sets.AddToAttributeSet("East", Eb);
136
137 // Create a new set spanning the West point
138 attr_sets.SetAttributeSet("West", Wa);
139 attr_sets.AddToAttributeSet("West", Wb);
140
141 // Create a new set consisting of the "a" sides of the compass rose
142 attr_sets.SetAttributeSet("Rose Even", Na);
143 attr_sets.AddToAttributeSet("Rose Even", Sa);
144 attr_sets.AddToAttributeSet("Rose Even", Ea);
145 attr_sets.AddToAttributeSet("Rose Even", Wa);
146
147 // Create a new set consisting of the "b" sides of the compass rose
148 attr_sets.SetAttributeSet("Rose Odd", Nb);
149 attr_sets.AddToAttributeSet("Rose Odd", Sb);
150 attr_sets.AddToAttributeSet("Rose Odd", Eb);
151 attr_sets.AddToAttributeSet("Rose Odd", Wb);
152
153 // Create a new set consisting of the full compass rose
154 Array<int> & Ra = attr_sets.GetAttributeSet("Rose Even");
155 Array<int> & Rb = attr_sets.GetAttributeSet("Rose Odd");
156 attr_sets.SetAttributeSet("Rose", Ra);
157 attr_sets.AddToAttributeSet("Rose", Rb);
158 }
159 // 4c. Define new boundary regions based on existing boundary attribute sets
160 {
161 Array<int> & NNE = bdr_attr_sets.GetAttributeSet("NNE");
162 Array<int> & NNW = bdr_attr_sets.GetAttributeSet("NNW");
163 Array<int> & ENE = bdr_attr_sets.GetAttributeSet("ENE");
164 Array<int> & ESE = bdr_attr_sets.GetAttributeSet("ESE");
165 Array<int> & SSE = bdr_attr_sets.GetAttributeSet("SSE");
166 Array<int> & SSW = bdr_attr_sets.GetAttributeSet("SSW");
167 Array<int> & WNW = bdr_attr_sets.GetAttributeSet("WNW");
168 Array<int> & WSW = bdr_attr_sets.GetAttributeSet("WSW");
169
170 bdr_attr_sets.SetAttributeSet("Northern Boundary", NNE);
171 bdr_attr_sets.AddToAttributeSet("Northern Boundary", NNW);
172
173 bdr_attr_sets.SetAttributeSet("Southern Boundary", SSE);
174 bdr_attr_sets.AddToAttributeSet("Southern Boundary", SSW);
175
176 bdr_attr_sets.SetAttributeSet("Eastern Boundary", ENE);
177 bdr_attr_sets.AddToAttributeSet("Eastern Boundary", ESE);
178
179 bdr_attr_sets.SetAttributeSet("Western Boundary", WNW);
180 bdr_attr_sets.AddToAttributeSet("Western Boundary", WSW);
181
182 bdr_attr_sets.SetAttributeSet("Boundary",
183 bdr_attr_sets.GetAttributeSet
184 ("Northern Boundary"));
185 bdr_attr_sets.AddToAttributeSet("Boundary",
186 bdr_attr_sets.GetAttributeSet
187 ("Southern Boundary"));
188 bdr_attr_sets.AddToAttributeSet("Boundary",
189 bdr_attr_sets.GetAttributeSet
190 ("Eastern Boundary"));
191 bdr_attr_sets.AddToAttributeSet("Boundary",
192 bdr_attr_sets.GetAttributeSet
193 ("Western Boundary"));
194 }
195
196 // 5. Define a finite element space on the mesh. Here we use continuous
197 // Lagrange finite elements of the specified order.
198 H1_FECollection fec(order, mesh.Dimension());
199 FiniteElementSpace fespace(&mesh, &fec);
200 cout << "Number of finite element unknowns: "
201 << fespace.GetTrueVSize() << endl;
202
203 // 6. Determine the list of true (i.e. conforming) essential boundary dofs.
204 // In this example, the boundary conditions are defined by marking all
205 // the boundary regions corresponding to the boundary attributes
206 // contained in the set named "ess_name" as essential (Dirichlet) and
207 // converting them to a list of true dofs.
208 Array<int> ess_tdof_list;
209 if (bdr_attr_sets.AttributeSetExists(ess_name))
210 {
211 Array<int> ess_bdr_marker = bdr_attr_sets.GetAttributeSetMarker(ess_name);
212 fespace.GetEssentialTrueDofs(ess_bdr_marker, ess_tdof_list);
213 }
214
215 // 7. Set up the linear form b(.) which corresponds to the right-hand side of
216 // the FEM linear system, which in this case is (1_s,phi_i) where phi_i
217 // are the basis functions in fespace and 1_s is an indicator function
218 // equal to 1 on the region defined by the named set "source_name" and
219 // zero elsewhere.
220 Array<int> source_marker = attr_sets.GetAttributeSetMarker(source_name);
221
222 LinearForm b(&fespace);
223 ConstantCoefficient one(1.0);
224 b.AddDomainIntegrator(new DomainLFIntegrator(one), source_marker);
225 b.Assemble();
226
227 // 8. Define the solution vector x as a finite element grid function
228 // corresponding to fespace. Initialize x with initial guess of zero,
229 // which satisfies the boundary conditions.
230 GridFunction x(&fespace);
231 x = 0.0;
232
233 // 9. Set up the bilinear form a(.,.) on the finite element space
234 // corresponding to the Laplacian operator -Delta, by adding the
235 // Diffusion domain integrator.
236 BilinearForm a(&fespace);
237
238 ConstantCoefficient defaultCoef(1.0e-6);
239 ConstantCoefficient baseCoef(1.0);
240 ConstantCoefficient roseCoef(2.0);
241
242 Array<int> base_marker = attr_sets.GetAttributeSetMarker("Base");
243 Array<int> rose_marker = attr_sets.GetAttributeSetMarker("Rose Even");
244
245 // Impose a very small diffusion coefficient across the entire mesh
246 a.AddDomainIntegrator(new DiffusionIntegrator(defaultCoef));
247
248 // Impose an additional, stronger diffusion coefficient in select regions
249 a.AddDomainIntegrator(new DiffusionIntegrator(baseCoef), base_marker);
250 a.AddDomainIntegrator(new DiffusionIntegrator(roseCoef), rose_marker);
251
252 // 10. Assemble the bilinear form and the corresponding linear system,
253 // applying any necessary transformations.
254 a.Assemble();
255
256 SparseMatrix A;
257 Vector B, X;
258 a.FormLinearSystem(ess_tdof_list, x, b, A, X, B);
259
260 cout << "Size of linear system: " << A.Height() << endl;
261
262 // 11. Solve the system using PCG with symmetric Gauss-Seidel preconditioner.
263 GSSmoother M(A);
264 PCG(A, M, B, X, 1, 800, 1e-12, 0.0);
265
266 // 12. Recover the solution as a finite element grid function.
267 a.RecoverFEMSolution(X, b, x);
268
269 // 13. Save the refined mesh and the solution. This output can be viewed
270 // later using GLVis: "glvis -m refined.mesh -g sol.gf".
271 mesh.Save("refined.mesh");
272 x.Save("sol.gf");
273
274 // 14. Send the solution by socket to a GLVis server.
275 if (visualization)
276 {
277 char vishost[] = "localhost";
278 int visport = 19916;
279 socketstream sol_sock(vishost, visport);
280 sol_sock.precision(8);
281 sol_sock << "solution\n" << mesh << x << "keys Rjmm" << flush;
282 }
283
284 return 0;
285}
void SetAttributeSet(const std::string &set_name, const Array< int > &attr)
Create a new attribute set.
Array< int > GetAttributeSetMarker(const std::string &set_name)
Return a marker array corresponding to a named attribute set.
bool AttributeSetExists(const std::string &name) const
Return true is the named attribute set is present.
Array< int > & GetAttributeSet(const std::string &set_name)
Access a named attribute set.
void AddToAttributeSet(const std::string &set_name, int attr)
Add a single entry to an existing attribute set.
std::set< std::string > GetAttributeSetNames() const
Return all attribute set names as an STL set.
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.
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
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
AttributeSets bdr_attribute_sets
Named sets of boundary element attributes.
Definition mesh.hpp:288
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
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
AttributeSets attribute_sets
Named sets of element attributes.
Definition mesh.hpp:285
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
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
Data type sparse matrix.
Definition sparsemat.hpp:51
Vector data type.
Definition vector.hpp:80
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
const char vishost[]