MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex31.cpp
Go to the documentation of this file.
1// MFEM Example 31
2//
3// Compile with: make ex31
4//
5// Sample runs: ex31 -m ../data/inline-segment.mesh -o 2
6// ex31 -m ../data/hexagon.mesh -o 2
7// ex31 -m ../data/star.mesh -o 2
8// ex31 -m ../data/fichera.mesh -o 3 -r 1
9// ex31 -m ../data/square-disc-nurbs.mesh -o 3
10// ex31 -m ../data/amr-quad.mesh -o 2 -r 1
11// ex31 -m ../data/amr-hex.mesh -r 1
12//
13// Description: This example code solves a simple electromagnetic diffusion
14// problem corresponding to the second order definite Maxwell
15// equation curl curl E + sigma E = f with boundary condition
16// E x n = <given tangential field>. In this example sigma is an
17// anisotropic 3x3 tensor. Here, we use a given exact solution E
18// and compute the corresponding r.h.s. f. We discretize with
19// Nedelec finite elements in 1D, 2D, or 3D.
20//
21// The example demonstrates the use of restricted H(curl) finite
22// element spaces with the curl-curl and the (vector finite
23// element) mass bilinear form, as well as the computation of
24// discretization error when the exact solution is known. These
25// restricted spaces allow the solution of 1D or 2D
26// electromagnetic problems which involve 3D field vectors. Such
27// problems arise in plasma physics and crystallography.
28//
29// We recommend viewing example 3 before viewing this example.
30
31#include "mfem.hpp"
32#include <fstream>
33#include <iostream>
34
35using namespace std;
36using namespace mfem;
37
38// Exact solution, E, and r.h.s., f. See below for implementation.
39void E_exact(const Vector &, Vector &);
40void CurlE_exact(const Vector &, Vector &);
41void f_exact(const Vector &, Vector &);
43int dim;
44
45int main(int argc, char *argv[])
46{
47 // 1. Parse command-line options.
48 const char *mesh_file = "../data/inline-quad.mesh";
49 int ref_levels = 2;
50 int order = 1;
51 bool visualization = 1;
52
53 OptionsParser args(argc, argv);
54 args.AddOption(&mesh_file, "-m", "--mesh",
55 "Mesh file to use.");
56 args.AddOption(&ref_levels, "-r", "--refine",
57 "Number of times to refine the mesh uniformly.");
58 args.AddOption(&order, "-o", "--order",
59 "Finite element order (polynomial degree).");
60 args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
61 " solution.");
62 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
63 "--no-visualization",
64 "Enable or disable GLVis visualization.");
65 args.ParseCheck();
66
67 kappa = freq * M_PI;
68
69 // 2. Read the mesh from the given mesh file. We can handle triangular,
70 // quadrilateral, or mixed meshes with the same code.
71 Mesh mesh(mesh_file, 1, 1);
72 dim = mesh.Dimension();
73
74 // 3. Refine the mesh to increase the resolution. In this example we do
75 // 'ref_levels' of uniform refinement (2 by default, or specified on
76 // the command line with -r).
77 for (int lev = 0; lev < ref_levels; lev++)
78 {
79 mesh.UniformRefinement();
80 }
81
82 // 4. Define a finite element space on the mesh. Here we use the Nedelec
83 // finite elements of the specified order restricted to 1D, 2D, or 3D
84 // depending on the dimension of the given mesh file.
85 FiniteElementCollection *fec = NULL;
86 if (dim == 1)
87 {
88 fec = new ND_R1D_FECollection(order, dim);
89 }
90 else if (dim == 2)
91 {
92 fec = new ND_R2D_FECollection(order, dim);
93 }
94 else
95 {
96 fec = new ND_FECollection(order, dim);
97 }
98 FiniteElementSpace fespace(&mesh, fec);
99 int size = fespace.GetTrueVSize();
100 cout << "Number of H(Curl) unknowns: " << size << endl;
101
102 // 5. Determine the list of true essential boundary dofs. In this example,
103 // the boundary conditions are defined by marking all the boundary
104 // attributes from the mesh as essential (Dirichlet) and converting them
105 // to a list of true dofs.
106 Array<int> ess_tdof_list;
107 if (mesh.bdr_attributes.Size())
108 {
109 Array<int> ess_bdr(mesh.bdr_attributes.Max());
110 ess_bdr = 1;
111 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
112 }
113
114 // 6. Set up the linear form b(.) which corresponds to the right-hand side
115 // of the FEM linear system, which in this case is (f,phi_i) where f is
116 // given by the function f_exact and phi_i are the basis functions in
117 // the finite element fespace.
119 LinearForm b(&fespace);
120 b.AddDomainIntegrator(new VectorFEDomainLFIntegrator(f));
121 b.Assemble();
122
123 // 7. Define the solution vector x as a finite element grid function
124 // corresponding to fespace. Initialize x by projecting the exact
125 // solution. Note that only values from the boundary edges will be used
126 // when eliminating the non-homogeneous boundary condition to modify the
127 // r.h.s. vector b.
128 GridFunction sol(&fespace);
131 sol.ProjectCoefficient(E);
132
133 // 8. Set up the bilinear form corresponding to the EM diffusion operator
134 // curl muinv curl + sigma I, by adding the curl-curl and the mass domain
135 // integrators.
136 DenseMatrix sigmaMat(3);
137 sigmaMat(0,0) = 2.0; sigmaMat(1,1) = 2.0; sigmaMat(2,2) = 2.0;
138 sigmaMat(0,2) = 0.0; sigmaMat(2,0) = 0.0;
139 sigmaMat(0,1) = M_SQRT1_2; sigmaMat(1,0) = M_SQRT1_2; // 1/sqrt(2) in cmath
140 sigmaMat(1,2) = M_SQRT1_2; sigmaMat(2,1) = M_SQRT1_2;
141
142 ConstantCoefficient muinv(1.0);
144 BilinearForm a(&fespace);
145 a.AddDomainIntegrator(new CurlCurlIntegrator(muinv));
146 a.AddDomainIntegrator(new VectorFEMassIntegrator(sigma));
147
148 // 9. Assemble the bilinear form and the corresponding linear system,
149 // applying any necessary transformations such as: eliminating boundary
150 // conditions, applying conforming constraints for non-conforming AMR,
151 // etc.
152 a.Assemble();
153
154 OperatorPtr A;
155 Vector B, X;
156
157 a.FormLinearSystem(ess_tdof_list, sol, b, A, X, B);
158
159 // 10. Solve the system A X = B.
160
161#ifndef MFEM_USE_SUITESPARSE
162 // 11. Define a simple symmetric Gauss-Seidel preconditioner and use it to
163 // solve the system Ax=b with PCG.
164 GSSmoother M((SparseMatrix&)(*A));
165 PCG(*A, M, B, X, 1, 500, 1e-12, 0.0);
166#else
167 // 11. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the
168 // system.
169 UMFPackSolver umf_solver;
170 umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
171 umf_solver.SetOperator(*A);
172 umf_solver.Mult(B, X);
173#endif
174
175 // 12. Recover the solution as a finite element grid function.
176 a.RecoverFEMSolution(X, b, sol);
177
178 // 13. Compute and print the H(Curl) norm of the error.
179 {
180 real_t error = sol.ComputeHCurlError(&E, &CurlE);
181 cout << "\n|| E_h - E ||_{H(Curl)} = " << error << '\n' << endl;
182 }
183
184
185 // 14. Save the refined mesh and the solution. This output can be viewed
186 // later using GLVis: "glvis -m refined.mesh -g sol.gf".
187 {
188 ofstream mesh_ofs("refined.mesh");
189 mesh_ofs.precision(8);
190 mesh.Print(mesh_ofs);
191 ofstream sol_ofs("sol.gf");
192 sol_ofs.precision(8);
193 sol.Save(sol_ofs);
194 }
195
196 // 15. Send the solution by socket to a GLVis server.
197 if (visualization)
198 {
199 char vishost[] = "localhost";
200 int visport = 19916;
201
202 VectorGridFunctionCoefficient solCoef(&sol);
203 CurlGridFunctionCoefficient dsolCoef(&sol);
204
205 if (dim ==1)
206 {
210 socketstream dy_sock(vishost, visport);
211 socketstream dz_sock(vishost, visport);
212 x_sock.precision(8);
213 y_sock.precision(8);
214 z_sock.precision(8);
215 dy_sock.precision(8);
216 dz_sock.precision(8);
217
218 Vector xVec(3); xVec = 0.0; xVec(0) = 1;
219 Vector yVec(3); yVec = 0.0; yVec(1) = 1;
220 Vector zVec(3); zVec = 0.0; zVec(2) = 1;
221 VectorConstantCoefficient xVecCoef(xVec);
222 VectorConstantCoefficient yVecCoef(yVec);
223 VectorConstantCoefficient zVecCoef(zVec);
224
225 H1_FECollection fec_h1(order, dim);
226 L2_FECollection fec_l2(order-1, dim);
227
228 FiniteElementSpace fes_h1(&mesh, &fec_h1);
229 FiniteElementSpace fes_l2(&mesh, &fec_l2);
230
231 GridFunction xComp(&fes_l2);
232 GridFunction yComp(&fes_h1);
233 GridFunction zComp(&fes_h1);
234
235 GridFunction dyComp(&fes_l2);
236 GridFunction dzComp(&fes_l2);
237
238 InnerProductCoefficient xCoef(xVecCoef, solCoef);
239 InnerProductCoefficient yCoef(yVecCoef, solCoef);
240 InnerProductCoefficient zCoef(zVecCoef, solCoef);
241
242 xComp.ProjectCoefficient(xCoef);
243 yComp.ProjectCoefficient(yCoef);
244 zComp.ProjectCoefficient(zCoef);
245
246 x_sock << "solution\n" << mesh << xComp << flush
247 << "window_title 'X component'" << endl;
248 y_sock << "solution\n" << mesh << yComp << flush
249 << "window_geometry 403 0 400 350 "
250 << "window_title 'Y component'" << endl;
251 z_sock << "solution\n" << mesh << zComp << flush
252 << "window_geometry 806 0 400 350 "
253 << "window_title 'Z component'" << endl;
254
255 InnerProductCoefficient dyCoef(yVecCoef, dsolCoef);
256 InnerProductCoefficient dzCoef(zVecCoef, dsolCoef);
257
258 dyComp.ProjectCoefficient(dyCoef);
259 dzComp.ProjectCoefficient(dzCoef);
260
261 dy_sock << "solution\n" << mesh << dyComp << flush
262 << "window_geometry 403 375 400 350 "
263 << "window_title 'Y component of Curl'" << endl;
264 dz_sock << "solution\n" << mesh << dzComp << flush
265 << "window_geometry 806 375 400 350 "
266 << "window_title 'Z component of Curl'" << endl;
267 }
268 else if (dim == 2)
269 {
270 socketstream xy_sock(vishost, visport);
272 socketstream dxy_sock(vishost, visport);
273 socketstream dz_sock(vishost, visport);
274
275 DenseMatrix xyMat(2,3); xyMat = 0.0;
276 xyMat(0,0) = 1.0; xyMat(1,1) = 1.0;
277 MatrixConstantCoefficient xyMatCoef(xyMat);
278 Vector zVec(3); zVec = 0.0; zVec(2) = 1;
279 VectorConstantCoefficient zVecCoef(zVec);
280
281 MatrixVectorProductCoefficient xyCoef(xyMatCoef, solCoef);
282 InnerProductCoefficient zCoef(zVecCoef, solCoef);
283
284 H1_FECollection fec_h1(order, dim);
285 ND_FECollection fec_nd(order, dim);
286 RT_FECollection fec_rt(order-1, dim);
287 L2_FECollection fec_l2(order-1, dim);
288
289 FiniteElementSpace fes_h1(&mesh, &fec_h1);
290 FiniteElementSpace fes_nd(&mesh, &fec_nd);
291 FiniteElementSpace fes_rt(&mesh, &fec_rt);
292 FiniteElementSpace fes_l2(&mesh, &fec_l2);
293
294 GridFunction xyComp(&fes_nd);
295 GridFunction zComp(&fes_h1);
296
297 GridFunction dxyComp(&fes_rt);
298 GridFunction dzComp(&fes_l2);
299
300 xyComp.ProjectCoefficient(xyCoef);
301 zComp.ProjectCoefficient(zCoef);
302
303 xy_sock.precision(8);
304 xy_sock << "solution\n" << mesh << xyComp
305 << "window_title 'XY components'\n" << flush;
306 z_sock << "solution\n" << mesh << zComp << flush
307 << "window_geometry 403 0 400 350 "
308 << "window_title 'Z component'" << endl;
309
310 MatrixVectorProductCoefficient dxyCoef(xyMatCoef, dsolCoef);
311 InnerProductCoefficient dzCoef(zVecCoef, dsolCoef);
312
313 dxyComp.ProjectCoefficient(dxyCoef);
314 dzComp.ProjectCoefficient(dzCoef);
315
316 dxy_sock << "solution\n" << mesh << dxyComp << flush
317 << "window_geometry 0 375 400 350 "
318 << "window_title 'XY components of Curl'" << endl;
319 dz_sock << "solution\n" << mesh << dzComp << flush
320 << "window_geometry 403 375 400 350 "
321 << "window_title 'Z component of Curl'" << endl;
322 }
323 else
324 {
325 socketstream sol_sock(vishost, visport);
326 socketstream dsol_sock(vishost, visport);
327
328 RT_FECollection fec_rt(order-1, dim);
329
330 FiniteElementSpace fes_rt(&mesh, &fec_rt);
331
332 GridFunction dsol(&fes_rt);
333
334 dsol.ProjectCoefficient(dsolCoef);
335
336 sol_sock.precision(8);
337 sol_sock << "solution\n" << mesh << sol
338 << "window_title 'Solution'" << flush << endl;
339 dsol_sock << "solution\n" << mesh << dsol << flush
340 << "window_geometry 0 375 400 350 "
341 << "window_title 'Curl of solution'" << endl;
342 }
343 }
344
345 // 16. Free the used memory.
346 delete fec;
347
348 return 0;
349}
350
351
352void E_exact(const Vector &x, Vector &E)
353{
354 if (dim == 1)
355 {
356 E(0) = 1.1 * sin(kappa * x(0) + 0.0 * M_PI);
357 E(1) = 1.2 * sin(kappa * x(0) + 0.4 * M_PI);
358 E(2) = 1.3 * sin(kappa * x(0) + 0.9 * M_PI);
359 }
360 else if (dim == 2)
361 {
362 E(0) = 1.1 * sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
363 E(1) = 1.2 * sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
364 E(2) = 1.3 * sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
365 }
366 else
367 {
368 E(0) = 1.1 * sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
369 E(1) = 1.2 * sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
370 E(2) = 1.3 * sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
371 E *= cos(kappa * x(2));
372 }
373}
374
375void CurlE_exact(const Vector &x, Vector &dE)
376{
377 if (dim == 1)
378 {
379 real_t c4 = cos(kappa * x(0) + 0.4 * M_PI);
380 real_t c9 = cos(kappa * x(0) + 0.9 * M_PI);
381
382 dE(0) = 0.0;
383 dE(1) = -1.3 * c9;
384 dE(2) = 1.2 * c4;
385 dE *= kappa;
386 }
387 else if (dim == 2)
388 {
389 real_t c0 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
390 real_t c4 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
391 real_t c9 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
392
393 dE(0) = 1.3 * c9;
394 dE(1) = -1.3 * c9;
395 dE(2) = 1.2 * c4 - 1.1 * c0;
396 dE *= kappa * M_SQRT1_2;
397 }
398 else
399 {
400 real_t s0 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
401 real_t c0 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
402 real_t s4 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
403 real_t c4 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
404 real_t c9 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
405 real_t sk = sin(kappa * x(2));
406 real_t ck = cos(kappa * x(2));
407
408 dE(0) = 1.2 * s4 * sk + 1.3 * M_SQRT1_2 * c9 * ck;
409 dE(1) = -1.1 * s0 * sk - 1.3 * M_SQRT1_2 * c9 * ck;
410 dE(2) = -M_SQRT1_2 * (1.1 * c0 - 1.2 * c4) * ck;
411 dE *= kappa;
412 }
413}
414
415void f_exact(const Vector &x, Vector &f)
416{
417 if (dim == 1)
418 {
419 real_t s0 = sin(kappa * x(0) + 0.0 * M_PI);
420 real_t s4 = sin(kappa * x(0) + 0.4 * M_PI);
421 real_t s9 = sin(kappa * x(0) + 0.9 * M_PI);
422
423 f(0) = 2.2 * s0 + 1.2 * M_SQRT1_2 * s4;
424 f(1) = 1.2 * (2.0 + kappa * kappa) * s4 +
425 M_SQRT1_2 * (1.1 * s0 + 1.3 * s9);
426 f(2) = 1.3 * (2.0 + kappa * kappa) * s9 + 1.2 * M_SQRT1_2 * s4;
427 }
428 else if (dim == 2)
429 {
430 real_t s0 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
431 real_t s4 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
432 real_t s9 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
433
434 f(0) = 0.55 * (4.0 + kappa * kappa) * s0 +
435 0.6 * (M_SQRT2 - kappa * kappa) * s4;
436 f(1) = 0.55 * (M_SQRT2 - kappa * kappa) * s0 +
437 0.6 * (4.0 + kappa * kappa) * s4 +
438 0.65 * M_SQRT2 * s9;
439 f(2) = 0.6 * M_SQRT2 * s4 + 1.3 * (2.0 + kappa * kappa) * s9;
440 }
441 else
442 {
443 real_t s0 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
444 real_t c0 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
445 real_t s4 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
446 real_t c4 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
447 real_t s9 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
448 real_t c9 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
449 real_t sk = sin(kappa * x(2));
450 real_t ck = cos(kappa * x(2));
451
452 f(0) = 0.55 * (4.0 + 3.0 * kappa * kappa) * s0 * ck +
453 0.6 * (M_SQRT2 - kappa * kappa) * s4 * ck -
454 0.65 * M_SQRT2 * kappa * kappa * c9 * sk;
455
456 f(1) = 0.55 * (M_SQRT2 - kappa * kappa) * s0 * ck +
457 0.6 * (4.0 + 3.0 * kappa * kappa) * s4 * ck +
458 0.65 * M_SQRT2 * s9 * ck -
459 0.65 * M_SQRT2 * kappa * kappa * c9 * sk;
460
461 f(2) = 0.6 * M_SQRT2 * s4 * ck -
462 M_SQRT2 * kappa * kappa * (0.55 * c0 + 0.6 * c4) * sk
463 + 1.3 * (2.0 + kappa * kappa) * s9 * ck;
464 }
465}
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.
Definition: coefficient.hpp:85
Integrator for for Nedelec elements.
Vector coefficient defined as the Curl of a vector GridFunction.
Data type dense matrix using column-major storage.
Definition: densemat.hpp:24
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:27
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 real_t ComputeHCurlError(VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns the error measured in H(curl)-norm for ND elements.
Definition: gridfunc.cpp:3127
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3628
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2346
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:260
Scalar coefficient defined as the inner product of two vector coefficients.
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition: fe_coll.hpp:330
Vector with associated FE space and LinearFormIntegrators.
Definition: linearform.hpp:25
A matrix coefficient that is constant in space and time.
Vector coefficient defined as a product of a matrix coefficient and a vector coefficient.
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 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
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:465
Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 1D.
Definition: fe_coll.hpp:520
Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 2D.
Definition: fe_coll.hpp:584
Pointer to an Operator of a specified type.
Definition: handle.hpp:34
void ParseCheck(std::ostream &out=mfem::out)
Definition: optparser.cpp:255
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
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:386
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 coefficient that is constant in space and time.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:347
A general vector function coefficient.
Vector coefficient defined by a vector GridFunction.
Vector data type.
Definition: vector.hpp:80
real_t sigma(const Vector &x)
Definition: maxwell.cpp:102
real_t kappa
Definition: ex31.cpp:42
void CurlE_exact(const Vector &, Vector &)
Definition: ex31.cpp:375
int dim
Definition: ex31.cpp:43
void E_exact(const Vector &, Vector &)
Definition: ex31.cpp:352
void f_exact(const Vector &, Vector &)
Definition: ex31.cpp:415
real_t freq
Definition: ex31.cpp:42
int visport
char vishost[]
int main()
real_t b
Definition: lissajous.cpp:42
real_t a
Definition: lissajous.cpp:41
real_t f(const Vector &p)
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
float real_t
Definition: config.hpp:43