MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex5.cpp
Go to the documentation of this file.
1 // MFEM Example 5
2 //
3 // Compile with: make ex5
4 //
5 // Sample runs: ex5 -m ../data/square-disc.mesh
6 // ex5 -m ../data/star.mesh
7 // ex5 -m ../data/star.mesh -pa
8 // ex5 -m ../data/beam-tet.mesh
9 // ex5 -m ../data/beam-hex.mesh
10 // ex5 -m ../data/beam-hex.mesh -pa
11 // ex5 -m ../data/escher.mesh
12 // ex5 -m ../data/fichera.mesh
13 //
14 // Device sample runs:
15 // ex5 -m ../data/star.mesh -pa -d cuda
16 // ex5 -m ../data/star.mesh -pa -d raja-cuda
17 // ex5 -m ../data/star.mesh -pa -d raja-omp
18 // ex5 -m ../data/beam-hex.mesh -pa -d cuda
19 //
20 // Description: This example code solves a simple 2D/3D mixed Darcy problem
21 // corresponding to the saddle point system
22 //
23 // k*u + grad p = f
24 // - div u = g
25 //
26 // with natural boundary condition -p = <given pressure>.
27 // Here, we use a given exact solution (u,p) and compute the
28 // corresponding r.h.s. (f,g). We discretize with Raviart-Thomas
29 // finite elements (velocity u) and piecewise discontinuous
30 // polynomials (pressure p).
31 //
32 // The example demonstrates the use of the BlockOperator class, as
33 // well as the collective saving of several grid functions in
34 // VisIt (visit.llnl.gov) and ParaView (paraview.org) formats.
35 //
36 // We recommend viewing examples 1-4 before viewing this example.
37 
38 #include "mfem.hpp"
39 #include <fstream>
40 #include <iostream>
41 #include <algorithm>
42 
43 using namespace std;
44 using namespace mfem;
45 
46 // Define the analytical solution and forcing terms / boundary conditions
47 void uFun_ex(const Vector & x, Vector & u);
48 double pFun_ex(const Vector & x);
49 void fFun(const Vector & x, Vector & f);
50 double gFun(const Vector & x);
51 double f_natural(const Vector & x);
52 
53 int main(int argc, char *argv[])
54 {
55  StopWatch chrono;
56 
57  // 1. Parse command-line options.
58  const char *mesh_file = "../data/star.mesh";
59  int order = 1;
60  bool pa = false;
61  const char *device_config = "cpu";
62  bool visualization = 1;
63 
64  OptionsParser args(argc, argv);
65  args.AddOption(&mesh_file, "-m", "--mesh",
66  "Mesh file to use.");
67  args.AddOption(&order, "-o", "--order",
68  "Finite element order (polynomial degree).");
69  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
70  "--no-partial-assembly", "Enable Partial Assembly.");
71  args.AddOption(&device_config, "-d", "--device",
72  "Device configuration string, see Device::Configure().");
73  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
74  "--no-visualization",
75  "Enable or disable GLVis visualization.");
76  args.Parse();
77  if (!args.Good())
78  {
79  args.PrintUsage(cout);
80  return 1;
81  }
82  args.PrintOptions(cout);
83 
84  // 2. Enable hardware devices such as GPUs, and programming models such as
85  // CUDA, OCCA, RAJA and OpenMP based on command line options.
86  Device device(device_config);
87  device.Print();
88 
89  // 3. Read the mesh from the given mesh file. We can handle triangular,
90  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
91  // the same code.
92  Mesh *mesh = new Mesh(mesh_file, 1, 1);
93  int dim = mesh->Dimension();
94 
95  // 4. Refine the mesh to increase the resolution. In this example we do
96  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
97  // largest number that gives a final mesh with no more than 10,000
98  // elements.
99  {
100  int ref_levels =
101  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
102  for (int l = 0; l < ref_levels; l++)
103  {
104  mesh->UniformRefinement();
105  }
106  }
107 
108  // 5. Define a finite element space on the mesh. Here we use the
109  // Raviart-Thomas finite elements of the specified order.
110  FiniteElementCollection *hdiv_coll(new RT_FECollection(order, dim));
111  FiniteElementCollection *l2_coll(new L2_FECollection(order, dim));
112 
113  FiniteElementSpace *R_space = new FiniteElementSpace(mesh, hdiv_coll);
114  FiniteElementSpace *W_space = new FiniteElementSpace(mesh, l2_coll);
115 
116  // 6. Define the BlockStructure of the problem, i.e. define the array of
117  // offsets for each variable. The last component of the Array is the sum
118  // of the dimensions of each block.
119  Array<int> block_offsets(3); // number of variables + 1
120  block_offsets[0] = 0;
121  block_offsets[1] = R_space->GetVSize();
122  block_offsets[2] = W_space->GetVSize();
123  block_offsets.PartialSum();
124 
125  std::cout << "***********************************************************\n";
126  std::cout << "dim(R) = " << block_offsets[1] - block_offsets[0] << "\n";
127  std::cout << "dim(W) = " << block_offsets[2] - block_offsets[1] << "\n";
128  std::cout << "dim(R+W) = " << block_offsets.Last() << "\n";
129  std::cout << "***********************************************************\n";
130 
131  // 7. Define the coefficients, analytical solution, and rhs of the PDE.
132  ConstantCoefficient k(1.0);
133 
134  VectorFunctionCoefficient fcoeff(dim, fFun);
135  FunctionCoefficient fnatcoeff(f_natural);
136  FunctionCoefficient gcoeff(gFun);
137 
138  VectorFunctionCoefficient ucoeff(dim, uFun_ex);
140 
141  // 8. Allocate memory (x, rhs) for the analytical solution and the right hand
142  // side. Define the GridFunction u,p for the finite element solution and
143  // linear forms fform and gform for the right hand side. The data
144  // allocated by x and rhs are passed as a reference to the grid functions
145  // (u,p) and the linear forms (fform, gform).
146  MemoryType mt = device.GetMemoryType();
147  BlockVector x(block_offsets, mt), rhs(block_offsets, mt);
148 
149  LinearForm *fform(new LinearForm);
150  fform->Update(R_space, rhs.GetBlock(0), 0);
153  fform->Assemble();
154  fform->SyncAliasMemory(rhs);
155 
156  LinearForm *gform(new LinearForm);
157  gform->Update(W_space, rhs.GetBlock(1), 0);
158  gform->AddDomainIntegrator(new DomainLFIntegrator(gcoeff));
159  gform->Assemble();
160  gform->SyncAliasMemory(rhs);
161 
162  // 9. Assemble the finite element matrices for the Darcy operator
163  //
164  // D = [ M B^T ]
165  // [ B 0 ]
166  // where:
167  //
168  // M = \int_\Omega k u_h \cdot v_h d\Omega u_h, v_h \in R_h
169  // B = -\int_\Omega \div u_h q_h d\Omega u_h \in R_h, q_h \in W_h
170  BilinearForm *mVarf(new BilinearForm(R_space));
171  MixedBilinearForm *bVarf(new MixedBilinearForm(R_space, W_space));
172 
173  if (pa) { mVarf->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
175  mVarf->Assemble();
176  if (!pa) { mVarf->Finalize(); }
177 
178  if (pa) { bVarf->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
180  bVarf->Assemble();
181  if (!pa) { bVarf->Finalize(); }
182 
183  BlockOperator darcyOp(block_offsets);
184 
185  TransposeOperator *Bt = NULL;
186 
187  if (pa)
188  {
189  Bt = new TransposeOperator(bVarf);
190 
191  darcyOp.SetBlock(0,0, mVarf);
192  darcyOp.SetBlock(0,1, Bt, -1.0);
193  darcyOp.SetBlock(1,0, bVarf, -1.0);
194  }
195  else
196  {
197  SparseMatrix &M(mVarf->SpMat());
198  SparseMatrix &B(bVarf->SpMat());
199  B *= -1.;
200  if (Device::IsEnabled()) { B.BuildTranspose(); }
201  Bt = new TransposeOperator(&B);
202 
203  darcyOp.SetBlock(0,0, &M);
204  darcyOp.SetBlock(0,1, Bt);
205  darcyOp.SetBlock(1,0, &B);
206  }
207 
208  // 10. Construct the operators for preconditioner
209  //
210  // P = [ diag(M) 0 ]
211  // [ 0 B diag(M)^-1 B^T ]
212  //
213  // Here we use Symmetric Gauss-Seidel to approximate the inverse of the
214  // pressure Schur Complement
215  SparseMatrix *MinvBt = NULL;
216  Vector Md(mVarf->Height());
217 
218  BlockDiagonalPreconditioner darcyPrec(block_offsets);
219  Solver *invM, *invS;
220  SparseMatrix *S = NULL;
221 
222  if (pa)
223  {
224  mVarf->AssembleDiagonal(Md);
225  auto Md_host = Md.HostRead();
226  Vector invMd(mVarf->Height());
227  for (int i=0; i<mVarf->Height(); ++i)
228  {
229  invMd(i) = 1.0 / Md_host[i];
230  }
231 
232  Vector BMBt_diag(bVarf->Height());
233  bVarf->AssembleDiagonal_ADAt(invMd, BMBt_diag);
234 
235  Array<int> ess_tdof_list; // empty
236 
237  invM = new OperatorJacobiSmoother(Md, ess_tdof_list);
238  invS = new OperatorJacobiSmoother(BMBt_diag, ess_tdof_list);
239  }
240  else
241  {
242  SparseMatrix &M(mVarf->SpMat());
243  M.GetDiag(Md);
244  Md.HostReadWrite();
245 
246  SparseMatrix &B(bVarf->SpMat());
247  MinvBt = Transpose(B);
248 
249  for (int i = 0; i < Md.Size(); i++)
250  {
251  MinvBt->ScaleRow(i, 1./Md(i));
252  }
253 
254  S = Mult(B, *MinvBt);
255 
256  invM = new DSmoother(M);
257 
258 #ifndef MFEM_USE_SUITESPARSE
259  invS = new GSSmoother(*S);
260 #else
261  invS = new UMFPackSolver(*S);
262 #endif
263  }
264 
265  invM->iterative_mode = false;
266  invS->iterative_mode = false;
267 
268  darcyPrec.SetDiagonalBlock(0, invM);
269  darcyPrec.SetDiagonalBlock(1, invS);
270 
271  // 11. Solve the linear system with MINRES.
272  // Check the norm of the unpreconditioned residual.
273  int maxIter(1000);
274  double rtol(1.e-6);
275  double atol(1.e-10);
276 
277  chrono.Clear();
278  chrono.Start();
279  MINRESSolver solver;
280  solver.SetAbsTol(atol);
281  solver.SetRelTol(rtol);
282  solver.SetMaxIter(maxIter);
283  solver.SetOperator(darcyOp);
284  solver.SetPreconditioner(darcyPrec);
285  solver.SetPrintLevel(1);
286  x = 0.0;
287  solver.Mult(rhs, x);
288  if (device.IsEnabled()) { x.HostRead(); }
289  chrono.Stop();
290 
291  if (solver.GetConverged())
292  {
293  std::cout << "MINRES converged in " << solver.GetNumIterations()
294  << " iterations with a residual norm of "
295  << solver.GetFinalNorm() << ".\n";
296  }
297  else
298  {
299  std::cout << "MINRES did not converge in " << solver.GetNumIterations()
300  << " iterations. Residual norm is " << solver.GetFinalNorm()
301  << ".\n";
302  }
303  std::cout << "MINRES solver took " << chrono.RealTime() << "s.\n";
304 
305  // 12. Create the grid functions u and p. Compute the L2 error norms.
306  GridFunction u, p;
307  u.MakeRef(R_space, x.GetBlock(0), 0);
308  p.MakeRef(W_space, x.GetBlock(1), 0);
309 
310  int order_quad = max(2, 2*order+1);
311  const IntegrationRule *irs[Geometry::NumGeom];
312  for (int i=0; i < Geometry::NumGeom; ++i)
313  {
314  irs[i] = &(IntRules.Get(i, order_quad));
315  }
316 
317  double err_u = u.ComputeL2Error(ucoeff, irs);
318  double norm_u = ComputeLpNorm(2., ucoeff, *mesh, irs);
319  double err_p = p.ComputeL2Error(pcoeff, irs);
320  double norm_p = ComputeLpNorm(2., pcoeff, *mesh, irs);
321 
322  std::cout << "|| u_h - u_ex || / || u_ex || = " << err_u / norm_u << "\n";
323  std::cout << "|| p_h - p_ex || / || p_ex || = " << err_p / norm_p << "\n";
324 
325  // 13. Save the mesh and the solution. This output can be viewed later using
326  // GLVis: "glvis -m ex5.mesh -g sol_u.gf" or "glvis -m ex5.mesh -g
327  // sol_p.gf".
328  {
329  ofstream mesh_ofs("ex5.mesh");
330  mesh_ofs.precision(8);
331  mesh->Print(mesh_ofs);
332 
333  ofstream u_ofs("sol_u.gf");
334  u_ofs.precision(8);
335  u.Save(u_ofs);
336 
337  ofstream p_ofs("sol_p.gf");
338  p_ofs.precision(8);
339  p.Save(p_ofs);
340  }
341 
342  // 14. Save data in the VisIt format
343  VisItDataCollection visit_dc("Example5", mesh);
344  visit_dc.RegisterField("velocity", &u);
345  visit_dc.RegisterField("pressure", &p);
346  visit_dc.Save();
347 
348  // 15. Save data in the ParaView format
349  ParaViewDataCollection paraview_dc("Example5", mesh);
350  paraview_dc.SetPrefixPath("ParaView");
351  paraview_dc.SetLevelsOfDetail(order);
352  paraview_dc.SetCycle(0);
353  paraview_dc.SetDataFormat(VTKFormat::BINARY);
354  paraview_dc.SetHighOrderOutput(true);
355  paraview_dc.SetTime(0.0); // set the time
356  paraview_dc.RegisterField("velocity",&u);
357  paraview_dc.RegisterField("pressure",&p);
358  paraview_dc.Save();
359 
360  // 16. Send the solution by socket to a GLVis server.
361  if (visualization)
362  {
363  char vishost[] = "localhost";
364  int visport = 19916;
365  socketstream u_sock(vishost, visport);
366  u_sock.precision(8);
367  u_sock << "solution\n" << *mesh << u << "window_title 'Velocity'" << endl;
368  socketstream p_sock(vishost, visport);
369  p_sock.precision(8);
370  p_sock << "solution\n" << *mesh << p << "window_title 'Pressure'" << endl;
371  }
372 
373  // 17. Free the used memory.
374  delete fform;
375  delete gform;
376  delete invM;
377  delete invS;
378  delete S;
379  delete Bt;
380  delete MinvBt;
381  delete mVarf;
382  delete bVarf;
383  delete W_space;
384  delete R_space;
385  delete l2_coll;
386  delete hdiv_coll;
387  delete mesh;
388 
389  return 0;
390 }
391 
392 
393 void uFun_ex(const Vector & x, Vector & u)
394 {
395  double xi(x(0));
396  double yi(x(1));
397  double zi(0.0);
398  if (x.Size() == 3)
399  {
400  zi = x(2);
401  }
402 
403  u(0) = - exp(xi)*sin(yi)*cos(zi);
404  u(1) = - exp(xi)*cos(yi)*cos(zi);
405 
406  if (x.Size() == 3)
407  {
408  u(2) = exp(xi)*sin(yi)*sin(zi);
409  }
410 }
411 
412 // Change if needed
413 double pFun_ex(const Vector & x)
414 {
415  double xi(x(0));
416  double yi(x(1));
417  double zi(0.0);
418 
419  if (x.Size() == 3)
420  {
421  zi = x(2);
422  }
423 
424  return exp(xi)*sin(yi)*cos(zi);
425 }
426 
427 void fFun(const Vector & x, Vector & f)
428 {
429  f = 0.0;
430 }
431 
432 double gFun(const Vector & x)
433 {
434  if (x.Size() == 3)
435  {
436  return -pFun_ex(x);
437  }
438  else
439  {
440  return 0;
441  }
442 }
443 
444 double f_natural(const Vector & x)
445 {
446  return (-pFun_ex(x));
447 }
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:97
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:540
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1387
void BuildTranspose() const
Build and store internally the transpose of this matrix which will be used in the methods AddMultTran...
Definition: sparsemat.cpp:822
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
int GetNumIterations() const
Definition: solvers.hpp:103
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
Data type for scaled Jacobi-type smoother of sparse matrix.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:920
void fFun(const Vector &x, Vector &f)
Definition: ex5.cpp:427
A class to handle Vectors in a block fashion.
Definition: blockvector.hpp:30
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
void SetDataFormat(VTKFormat fmt)
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:456
Helper class for ParaView visualization data.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:102
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
virtual void Save()
Save the collection and a VisIt root file.
virtual void AssembleDiagonal(Vector &diag) const
Assemble the diagonal of the bilinear form into diag. Note that diag is a tdof Vector.
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:279
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
double RealTime()
Definition: tic_toc.cpp:426
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:652
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
MINRES method.
Definition: solvers.hpp:426
void SyncAliasMemory(const Vector &v) const
Update the alias memory location of the vector to match v.
Definition: vector.hpp:235
static bool IsEnabled()
Return true if any backend other than Backend::CPU is enabled.
Definition: device.hpp:246
double GetFinalNorm() const
Definition: solvers.hpp:105
Data type for Gauss-Seidel smoother of sparse matrix.
void Stop()
Stop the stopwatch.
Definition: tic_toc.cpp:416
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds a domain integrator. Assumes ownership of bfi.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:839
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3484
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1453
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
Compute the Lp norm of a function f. .
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:439
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
Data type sparse matrix.
Definition: sparsemat.hpp:41
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
constexpr char vishost[]
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:130
double f_natural(const Vector &x)
Definition: ex5.cpp:444
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
constexpr int visport
Timing object.
Definition: tic_toc.hpp:34
static MemoryType GetMemoryType()
(DEPRECATED) Equivalent to GetDeviceMemoryType().
Definition: device.hpp:277
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
Definition: solvers.hpp:100
void Assemble(int skip_zeros=1)
double pFun_ex(const Vector &x)
Definition: ex5.cpp:413
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:201
void SetHighOrderOutput(bool high_order_output_)
void AssembleDiagonal_ADAt(const Vector &D, Vector &diag) const
Assemble the diagonal of ADA^T into diag, where A is this mixed bilinear form and D is a diagonal...
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1439
int GetConverged() const
Definition: solvers.hpp:104
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:70
int Dimension() const
Definition: mesh.hpp:911
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
void SetTime(double t)
Set physical time (for time-dependent simulations)
void Start()
Clear the elapsed time and start the stopwatch.
Definition: tic_toc.cpp:411
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:341
A general vector function coefficient.
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACY.
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
void SetAbsTol(double atol)
Definition: solvers.hpp:99
void SetRelTol(double rtol)
Definition: solvers.hpp:98
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
MemoryType
Memory types supported by MFEM.
Definition: mem_manager.hpp:31
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 &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
Definition: array.cpp:103
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
Definition: table.cpp:414
The transpose of a given operator. Switches the roles of the methods Mult() and MultTranspose().
Definition: operator.hpp:711
void Update()
Update the object according to the associated FE space fes.
Definition: linearform.hpp:187
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix: .
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
void GetDiag(Vector &d) const
Returns the Diagonal of A.
Definition: sparsemat.cpp:549
void uFun_ex(const Vector &x, Vector &u)
Definition: ex5.cpp:393
double gFun(const Vector &x)
Definition: ex5.cpp:432
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
int dim
Definition: ex24.cpp:53
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
T & Last()
Return the last element in the array.
Definition: array.hpp:779
void SetLevelsOfDetail(int levels_of_detail_)
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:266
A general function coefficient.
Vector data type.
Definition: vector.hpp:60
Vector with associated FE space and LinearFormIntegrators.
Definition: linearform.hpp:23
virtual void Save() override
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix.
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
Base class for solvers.
Definition: operator.hpp:648
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:121
A class to handle Block systems in a matrix-free implementation.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:377
int main()
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
void Clear()
Clear the elapsed time on the stopwatch and restart it if it&#39;s running.
Definition: tic_toc.cpp:406
void SetBlock(int iRow, int iCol, Operator *op, double c=1.0)
Add a block op in the block-entry (iblock, jblock).
Vector & GetBlock(int i)
Get the i-th vector in the block.
Definition: blockvector.hpp:90
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:285
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150