MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
block-solvers.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11//
12// ----------------------------------------------------------
13// Block Solvers Miniapp: Compare Saddle Point System Solvers
14// ----------------------------------------------------------
15//
16// This miniapp compares various linear solvers for the saddle point system
17// obtained from mixed finite element discretization of the simple mixed Darcy
18// problem in ex5p
19//
20// k*u + grad p = f
21// - div u = g
22//
23// with natural boundary condition -p = <given pressure>. We use a given exact
24// solution (u,p) and compute the corresponding r.h.s. (f,g). We discretize
25// with Raviart-Thomas finite elements (velocity u) and piecewise discontinuous
26// polynomials (pressure p).
27//
28// The solvers being compared include:
29// 1. MINRES preconditioned by a block diagonal preconditioner
30// 2. The divergence free solver (couple and decoupled modes)
31// 3. The Bramble-Pasciak solver (using BPCG or regular PCG)
32//
33// We recommend viewing example 5 before viewing this miniapp.
34//
35// Sample runs:
36//
37// mpirun -np 8 block-solvers -r 2 -o 0
38// mpirun -np 8 block-solvers -m anisotropic.mesh -c anisotropic.coeff -eb anisotropic.bdr
39//
40//
41// NOTE: The coefficient file (provided through -c) defines a piecewise constant
42// scalar coefficient k. The number of entries in this file should equal
43// to the number of "element attributes" in the mesh file. The value of
44// the coefficient in elements with the i-th attribute is given by the
45// i-th entry of the coefficient file.
46//
47//
48// NOTE: The essential boundary attribute file (provided through -eb) defines
49// which attributes to impose essential boundary condition (on u). The
50// number of entries in this file should equal to the number of "boundary
51// attributes" in the mesh file. If the i-th entry of the file is nonzero
52// (respectively 0), essential (respectively natural) boundary condition
53// will be imposed on boundary with the i-th attribute.
54
55#include "mfem.hpp"
56#include "bramble_pasciak.hpp"
57#include "div_free_solver.hpp"
58#include <fstream>
59#include <iostream>
60#include <memory>
61
62using namespace std;
63using namespace mfem;
64using namespace blocksolvers;
65
66// Exact solution, u and p, and r.h.s., f and g.
67void u_exact(const Vector & x, Vector & u);
68real_t p_exact(const Vector & x);
69void f_exact(const Vector & x, Vector & f);
70real_t g_exact(const Vector & x);
71real_t natural_bc(const Vector & x);
72
73/** Wrapper for assembling the discrete Darcy problem (ex5p)
74 [ M B^T ] [u] = [f]
75 [ B 0 ] [p] = [g]
76 where:
77 M = \int_\Omega (k u_h) \cdot v_h dx,
78 B = -\int_\Omega (div_h u_h) q_h dx,
79 f = \int_\Omega f_exact v_h dx + \int_D natural_bc v_h dS,
80 g = \int_\Omega g_exact q_h dx,
81 u_h, v_h \in R_h (Raviart-Thomas finite element space),
82 q_h \in W_h (piecewise discontinuous polynomials),
83 D: subset of the boundary where natural boundary condition is imposed. */
84class DarcyProblem
85{
86 OperatorPtr M_;
87 OperatorPtr B_;
88 Vector rhs_;
89 Vector ess_data_;
92 ParMesh mesh_;
93 ParBilinearForm *mVarf_;
96 FunctionCoefficient pcoeff_;
97 DFSSpaces dfs_spaces_;
98 PWConstCoefficient mass_coeff;
100public:
101 DarcyProblem(Mesh &mesh, int num_refines, int order, const char *coef_file,
102 Array<int> &ess_bdr, DFSParameters param);
103
104 HypreParMatrix& GetM() { return *M_.As<HypreParMatrix>(); }
105 HypreParMatrix& GetB() { return *B_.As<HypreParMatrix>(); }
106 const Vector& GetRHS() { return rhs_; }
107 const Vector& GetEssentialBC() { return ess_data_; }
108 const DFSData& GetDFSData() const { return dfs_spaces_.GetDFSData(); }
109 void ShowError(const Vector &sol, bool verbose);
110 void VisualizeSolution(const Vector &sol, std::string tag);
111 ParBilinearForm* GetMform() const { return mVarf_; }
112 ParMixedBilinearForm* GetBform() const { return bVarf_; }
113};
114
115DarcyProblem::DarcyProblem(Mesh &mesh, int num_refs, int order,
116 const char *coef_file, Array<int> &ess_bdr,
117 DFSParameters dfs_param)
118 : mesh_(MPI_COMM_WORLD, mesh), ucoeff_(mesh.Dimension(), u_exact),
119 pcoeff_(p_exact), dfs_spaces_(order, num_refs, &mesh_, ess_bdr, dfs_param),
120 mass_coeff()
121{
122 for (int l = 0; l < num_refs; l++)
123 {
124 mesh_.UniformRefinement();
125 dfs_spaces_.CollectDFSData();
126 }
127
128 Vector coef_vector(mesh.GetNE());
129 coef_vector = 1.0;
130 if (std::strcmp(coef_file, ""))
131 {
132 ifstream coef_str(coef_file);
133 coef_vector.Load(coef_str, mesh.GetNE());
134 }
135
136 mass_coeff.UpdateConstants(coef_vector);
140
141 u_.SetSpace(dfs_spaces_.GetHdivFES());
142 p_.SetSpace(dfs_spaces_.GetL2FES());
143 p_ = 0.0;
144 u_ = 0.0;
145 u_.ProjectBdrCoefficientNormal(ucoeff_, ess_bdr);
146
147 ParLinearForm fform(dfs_spaces_.GetHdivFES());
148 fform.AddDomainIntegrator(new VectorFEDomainLFIntegrator(fcoeff));
149 fform.AddBoundaryIntegrator(new VectorFEBoundaryFluxLFIntegrator(natcoeff));
150 fform.Assemble();
151
152 ParLinearForm gform(dfs_spaces_.GetL2FES());
153 gform.AddDomainIntegrator(new DomainLFIntegrator(gcoeff));
154 gform.Assemble();
155
156 mVarf_ = new ParBilinearForm(dfs_spaces_.GetHdivFES());
157 bVarf_ = new ParMixedBilinearForm(dfs_spaces_.GetHdivFES(),
158 dfs_spaces_.GetL2FES());
159
160 mVarf_->AddDomainIntegrator(new VectorFEMassIntegrator(mass_coeff));
161 mVarf_->ComputeElementMatrices();
162 mVarf_->Assemble();
163 mVarf_->EliminateEssentialBC(ess_bdr, u_, fform);
164
165 mVarf_->Finalize();
166 M_.Reset(mVarf_->ParallelAssemble());
167
169 bVarf_->Assemble();
170 bVarf_->SpMat() *= -1.0;
171 bVarf_->EliminateTrialDofs(ess_bdr, u_, gform);
172 bVarf_->Finalize();
173 B_.Reset(bVarf_->ParallelAssemble());
174
175 rhs_.SetSize(M_->NumRows() + B_->NumRows());
176 Vector rhs_block0(rhs_.GetData(), M_->NumRows());
177 Vector rhs_block1(rhs_.GetData()+M_->NumRows(), B_->NumRows());
178 fform.ParallelAssemble(rhs_block0);
179 gform.ParallelAssemble(rhs_block1);
180
181 ess_data_.SetSize(M_->NumRows() + B_->NumRows());
182 ess_data_ = 0.0;
183 Vector ess_data_block0(ess_data_.GetData(), M_->NumRows());
184 u_.ParallelProject(ess_data_block0);
185
186 int order_quad = max(2, 2*order+1);
187 for (int i=0; i < Geometry::NumGeom; ++i)
188 {
189 irs_[i] = &(IntRules.Get(i, order_quad));
190 }
191}
192
193void DarcyProblem::ShowError(const Vector& sol, bool verbose)
194{
195 u_.Distribute(Vector(sol.GetData(), M_->NumRows()));
196 p_.Distribute(Vector(sol.GetData()+M_->NumRows(), B_->NumRows()));
197
198 real_t err_u = u_.ComputeL2Error(ucoeff_, irs_);
199 real_t norm_u = ComputeGlobalLpNorm(2, ucoeff_, mesh_, irs_);
200 real_t err_p = p_.ComputeL2Error(pcoeff_, irs_);
201 real_t norm_p = ComputeGlobalLpNorm(2, pcoeff_, mesh_, irs_);
202
203 if (!verbose) { return; }
204 mfem::out << "|| u_h - u_ex || / || u_ex || = " << err_u / norm_u << "\n";
205 mfem::out << "|| p_h - p_ex || / || p_ex || = " << err_p / norm_p << "\n";
206}
207
208void DarcyProblem::VisualizeSolution(const Vector& sol, string tag)
209{
210 int num_procs, myid;
211 MPI_Comm_size(mesh_.GetComm(), &num_procs);
212 MPI_Comm_rank(mesh_.GetComm(), &myid);
213
214 u_.Distribute(Vector(sol.GetData(), M_->NumRows()));
215 p_.Distribute(Vector(sol.GetData()+M_->NumRows(), B_->NumRows()));
216
217 const char vishost[] = "localhost";
218 const int visport = 19916;
220 u_sock << "parallel " << num_procs << " " << myid << "\n";
221 u_sock.precision(8);
222 u_sock << "solution\n" << mesh_ << u_ << "window_title 'Velocity ("
223 << tag << " solver)'" << endl;
224 MPI_Barrier(mesh_.GetComm());
226 p_sock << "parallel " << num_procs << " " << myid << "\n";
227 p_sock.precision(8);
228 p_sock << "solution\n" << mesh_ << p_ << "window_title 'Pressure ("
229 << tag << " solver)'" << endl;
230}
231
232bool IsAllNeumannBoundary(const Array<int>& ess_bdr_attr)
233{
234 for (int attr : ess_bdr_attr) { if (attr == 0) { return false; } }
235 return true;
236}
237
238int main(int argc, char *argv[])
239{
240#ifdef HYPRE_USING_GPU
241 mfem::out << "\nAs of mfem-4.3 and hypre-2.22.0 (July 2021) this miniapp\n"
242 << "is NOT supported with the GPU version of hypre.\n\n";
243 return MFEM_SKIP_RETURN_VALUE;
244#endif
245
246 // Initialize MPI and HYPRE.
247 Mpi::Init(argc, argv);
248 Hypre::Init();
249
250 StopWatch chrono;
251
252 // Parse command-line options.
253 const char *mesh_file = "../../data/beam-hex.mesh";
254 const char *coef_file = "";
255 const char *ess_bdr_attr_file = "";
256 int order = 0;
257 int ser_ref_levels = 1;
258 int par_ref_levels = 1;
259 bool show_error = false;
260 bool visualization = false;
261
262 DFSParameters param;
263 BPSParameters bps_param;
264
265 OptionsParser args(argc, argv);
266 args.AddOption(&mesh_file, "-m", "--mesh",
267 "Mesh file to use.");
268 args.AddOption(&order, "-o", "--order",
269 "Finite element order (polynomial degree).");
270 args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
271 "Number of serial refinement steps.");
272 args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
273 "Number of parallel refinement steps.");
274 args.AddOption(&coef_file, "-c", "--coef",
275 "Coefficient file to use.");
276 args.AddOption(&ess_bdr_attr_file, "-eb", "--ess-bdr",
277 "Essential boundary attribute file to use.");
278 args.AddOption(&show_error, "-se", "--show-error", "-no-se",
279 "--no-show-error",
280 "Show or not show approximation error.");
281 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
282 "--no-visualization",
283 "Enable or disable GLVis visualization.");
284 args.Parse();
285 if (Mpi::Root()) { args.ParseCheck(); }
286
287 if (Mpi::Root() && par_ref_levels == 0)
288 {
289 mfem::out << "WARNING: DivFree solver is equivalent to BDPMinresSolver "
290 << "when par_ref_levels == 0.\n";
291 }
292
293 // Initialize the mesh, boundary attributes, and solver parameters
294 Mesh *mesh = new Mesh(mesh_file, 1, 1);
295
296 int dim = mesh->Dimension();
297
298 if (Mpi::Root())
299 {
300 mfem::out << "Number of serial refinements: " << ser_ref_levels << "\n"
301 << "Number of parallel refinements: " << par_ref_levels << "\n";
302 }
303
304 for (int i = 0; i < ser_ref_levels; ++i)
305 {
306 mesh->UniformRefinement();
307 }
308
309 if (Mpi::Root() && Mpi::WorldSize() > mesh->GetNE())
310 {
311 mfem::out << "\nWARNING: Number of processors is greater than the number of "
312 << "elements in the mesh.\n"
313 << "Number of processors: " << Mpi::WorldSize() << "\n"
314 << "Number of elements: " << mesh->GetNE() << "\n\n";
315 }
316
317 Array<int> ess_bdr(mesh->bdr_attributes.Max());
318 ess_bdr = 0;
319 if (std::strcmp(ess_bdr_attr_file, ""))
320 {
321 ifstream ess_bdr_attr_str(ess_bdr_attr_file);
322 ess_bdr.Load(mesh->bdr_attributes.Max(), ess_bdr_attr_str);
323 }
324 if (IsAllNeumannBoundary(ess_bdr))
325 {
326 if (Mpi::Root())
327 {
328 mfem::out << "\nSolution is not unique when Neumann boundary condition is "
329 << "imposed on the entire boundary. \nPlease provide a different "
330 << "boundary condition.\n";
331 }
332 delete mesh;
333 return 0;
334 }
335
336 string line = "**********************************************************\n";
337
338 chrono.Restart();
339
340 // Generate components of the saddle point problem
341 DarcyProblem darcy(*mesh, par_ref_levels, order, coef_file, ess_bdr, param);
342 HypreParMatrix& M = darcy.GetM();
343 HypreParMatrix& B = darcy.GetB();
344 const DFSData& DFS_data = darcy.GetDFSData();
345 delete mesh;
346
347 if (Mpi::Root())
348 {
349 mfem::out << line << "System assembled in " << chrono.RealTime() << "s.\n";
350 mfem::out << "Dimension of the physical space: " << dim << "\n";
351 mfem::out << "Size of the discrete Darcy system: " << M.M() + B.M() << "\n";
352 if (par_ref_levels > 0)
353 {
354 mfem::out << "Dimension of the divergence free subspace: "
355 << DFS_data.C.back().Ptr()->NumCols() << "\n\n";
356 }
357 }
358
359 // Setup various solvers for the discrete problem
360 std::map<const DarcySolver*, real_t> setup_time;
361 chrono.Restart();
362 BDPMinresSolver bdp(M, B, param);
363 setup_time[&bdp] = chrono.RealTime();
364
365 chrono.Restart();
366 DivFreeSolver dfs_dm(M, B, DFS_data);
367 setup_time[&dfs_dm] = chrono.RealTime();
368
369 chrono.Restart();
370 const_cast<bool&>(DFS_data.param.coupled_solve) = true;
371 DivFreeSolver dfs_cm(M, B, DFS_data);
372 setup_time[&dfs_cm] = chrono.RealTime();
373
374 chrono.Restart();
375 BramblePasciakSolver bp_bpcg(darcy.GetMform(), darcy.GetBform(), bps_param);
376 setup_time[&bp_bpcg] = chrono.RealTime();
377
378 chrono.Restart();
379 bps_param.use_bpcg = false;
380 BramblePasciakSolver bp_pcg(darcy.GetMform(), darcy.GetBform(), bps_param);
381 setup_time[&bp_pcg] = chrono.RealTime();
382
383 std::map<const DarcySolver*, std::string> solver_to_name;
384 solver_to_name[&bdp] = "Block-diagonal-preconditioned MINRES";
385 solver_to_name[&dfs_dm] = "Divergence free (decoupled mode)";
386 solver_to_name[&dfs_cm] = "Divergence free (coupled mode)";
387 solver_to_name[&bp_bpcg] = "Bramble Pasciak CG (using BPCG)";
388 solver_to_name[&bp_pcg] = "Bramble Pasciak CG (using regular PCG)";
389
390 // Solve the problem using all solvers
391 for (const auto& solver_pair : solver_to_name)
392 {
393 auto& solver = solver_pair.first;
394 auto& name = solver_pair.second;
395
396 Vector sol = darcy.GetEssentialBC();
397
398 chrono.Restart();
399 solver->Mult(darcy.GetRHS(), sol);
400 chrono.Stop();
401
402 if (Mpi::Root())
403 {
404 mfem::out << line << name << " solver:\n Setup time: "
405 << setup_time[solver] << "s.\n Solve time: "
406 << chrono.RealTime() << "s.\n Total time: "
407 << setup_time[solver] + chrono.RealTime() << "s.\n"
408 << " Iteration count: " << solver->GetNumIterations() <<"\n\n";
409 }
410 if (show_error && std::strcmp(coef_file, "") == 0)
411 {
412 darcy.ShowError(sol, Mpi::Root());
413 }
414 else if (show_error && Mpi::Root())
415 {
416 mfem::out << "Exact solution is unknown for coefficient '" << coef_file
417 << "'.\nApproximation error is computed in this case!\n\n";
418 }
419
420 if (visualization) { darcy.VisualizeSolution(sol, name); }
421
422 }
423
424 return 0;
425}
426
427void u_exact(const Vector & x, Vector & u)
428{
429 real_t xi(x(0));
430 real_t yi(x(1));
431 real_t zi(x.Size() == 3 ? x(2) : 0.0);
432
433 u(0) = - exp(xi)*sin(yi)*cos(zi);
434 u(1) = - exp(xi)*cos(yi)*cos(zi);
435 if (x.Size() == 3)
436 {
437 u(2) = exp(xi)*sin(yi)*sin(zi);
438 }
439}
440
442{
443 real_t xi(x(0));
444 real_t yi(x(1));
445 real_t zi(x.Size() == 3 ? x(2) : 0.0);
446 return exp(xi)*sin(yi)*cos(zi);
447}
448
449void f_exact(const Vector & x, Vector & f)
450{
451 f = 0.0;
452}
453
455{
456 if (x.Size() == 3) { return -p_exact(x); }
457 return 0;
458}
459
461{
462 return (-p_exact(x));
463}
real_t p_exact(const Vector &x)
real_t g_exact(const Vector &x)
real_t natural_bc(const Vector &x)
void u_exact(const Vector &x, Vector &u)
void f_exact(const Vector &x, Vector &f)
bool IsAllNeumannBoundary(const Array< int > &ess_bdr_attr)
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
void Load(std::istream &in, int fmt=0)
Read an Array from the stream in using format fmt. The format fmt can be:
Definition array.cpp:53
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY....
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void ComputeElementMatrices()
Compute and store internally all element matrices.
void EliminateEssentialBC(const Array< int > &bdr_attr_is_ess, const Vector &sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Eliminate essential boundary DOFs from the system.
Class for domain integration .
Definition lininteg.hpp:109
A general function coefficient.
static const int NumGeom
Definition geom.hpp:42
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, const Array< int > &bdr_attr)
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
HYPRE_BigInt M() const
Returns the global number of rows.
Definition hypre.hpp:627
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.hpp:74
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
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
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 UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY.
void Assemble(int skip_zeros=1)
void EliminateTrialDofs(const Array< int > &bdr_attr_is_ess, const Vector &sol, Vector &rhs)
Eliminate essential boundary DOFs from the columns of the system.
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix: .
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds a domain integrator. Assumes ownership of bfi.
static bool Root()
Return true if the rank in MPI_COMM_WORLD is zero.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
Pointer to an Operator of a specified type.
Definition handle.hpp:34
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get().
Definition handle.hpp:104
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition handle.hpp:145
int NumRows() const
Get the number of rows (size of output) of the Operator. Synonym with Height().
Definition operator.hpp:69
void ParseCheck(std::ostream &out=mfem::out)
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
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
A piecewise constant coefficient with the constants keyed off the element attribute numbers.
void UpdateConstants(Vector &c)
Update the constants with vector c.
Class for parallel bilinear form.
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void Assemble(int skip_zeros=1)
Assemble the local matrix.
Class for parallel grid function.
Definition pgridfunc.hpp:33
real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
void SetSpace(FiniteElementSpace *f) override
Associate a new FiniteElementSpace with the ParGridFunction.
Definition pgridfunc.cpp:96
void Distribute(const Vector *tv)
Class for parallel linear form.
Class for parallel meshes.
Definition pmesh.hpp:34
MPI_Comm GetComm() const
Definition pmesh.hpp:402
Class for parallel bilinear form using different test and trial FE spaces.
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P_test^t A P_trial.
Timing object.
Definition tic_toc.hpp:36
double RealTime()
Return the number of real seconds elapsed since the stopwatch was started.
Definition tic_toc.cpp:432
void Restart()
Clears and restarts the stopwatch. Equivalent to Clear() followed by Start().
Definition tic_toc.cpp:416
void Stop()
Stop the stopwatch.
Definition tic_toc.cpp:422
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition lininteg.hpp:347
A general vector function coefficient.
Vector data type.
Definition vector.hpp:80
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
Definition vector.hpp:227
Wrapper for the block-diagonal-preconditioned MINRES employed in ex5p.cpp.
Bramble-Pasciak Solver for Darcy equation.
const DFSData & GetDFSData() const
ParFiniteElementSpace * GetHdivFES() const
ParFiniteElementSpace * GetL2FES() const
int dim
Definition ex24.cpp:53
int main()
const int visport
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
real_t ComputeGlobalLpNorm(real_t p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition globals.hpp:66
float real_t
Definition config.hpp:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:486
const char vishost[]
Parameters for the BramblePasciakSolver method.
Data for the divergence free solver.
std::vector< OperatorPtr > C
Parameters for the divergence free solver.