MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
block-solvers.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2021, 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. The divergence free solver (couple and decoupled modes)
30 // 2. MINRES preconditioned by a block diagonal preconditioner
31 //
32 // We recommend viewing example 5 before viewing this miniapp.
33 //
34 // Sample runs:
35 //
36 // mpirun -np 8 block-solvers-compare -r 2 -o 0
37 // mpirun -np 8 block-solvers-compare -m anisotropic.mesh -c anisotropic.coeff -be anisotropic.bdr
38 //
39 //
40 // NOTE: The coefficient file (provided through -c) defines a piecewise constant
41 // scalar coefficient k. The number of entries in this file should equal
42 // to the number of "element attributes" in the mesh file. The value of
43 // the coefficient in elements with the i-th attribute is given by the
44 // i-th entry of the coefficient file.
45 //
46 //
47 // NOTE: The essential boundary attribute file (provided through -eb) defines
48 // which attributes to impose essential boundary condition (on u). The
49 // number of entries in this file should equal to the number of "boundary
50 // attributes" in the mesh file. If the i-th entry of the file is nonzero
51 // (respectively 0), essential (respectively natural) boundary condition
52 // will be imposed on boundary with the i-th attribute.
53 
54 #include "mfem.hpp"
55 #include "div_free_solver.hpp"
56 #include <fstream>
57 #include <iostream>
58 #include <memory>
59 
60 using namespace std;
61 using namespace mfem;
62 using namespace blocksolvers;
63 
64 // Exact solution, u and p, and r.h.s., f and g.
65 void u_exact(const Vector & x, Vector & u);
66 double p_exact(const Vector & x);
67 void f_exact(const Vector & x, Vector & f);
68 double g_exact(const Vector & x);
69 double natural_bc(const Vector & x);
70 
71 /** Wrapper for assembling the discrete Darcy problem (ex5p)
72  [ M B^T ] [u] = [f]
73  [ B 0 ] [p] = [g]
74  where:
75  M = \int_\Omega (k u_h) \cdot v_h dx,
76  B = -\int_\Omega (div_h u_h) q_h dx,
77  f = \int_\Omega f_exact v_h dx + \int_D natural_bc v_h dS,
78  g = \int_\Omega g_exact q_h dx,
79  u_h, v_h \in R_h (Raviart-Thomas finite element space),
80  q_h \in W_h (piecewise discontinuous polynomials),
81  D: subset of the boundary where natural boundary condition is imposed. */
82 class DarcyProblem
83 {
84  OperatorPtr M_;
85  OperatorPtr B_;
86  Vector rhs_;
87  Vector ess_data_;
88  ParGridFunction u_;
89  ParGridFunction p_;
90  ParMesh mesh_;
92  FunctionCoefficient pcoeff_;
93  DFSSpaces dfs_spaces_;
94  const IntegrationRule *irs_[Geometry::NumGeom];
95 public:
96  DarcyProblem(Mesh &mesh, int num_refines, int order, const char *coef_file,
97  Array<int> &ess_bdr, DFSParameters param);
98 
99  HypreParMatrix& GetM() { return *M_.As<HypreParMatrix>(); }
100  HypreParMatrix& GetB() { return *B_.As<HypreParMatrix>(); }
101  const Vector& GetRHS() { return rhs_; }
102  const Vector& GetEssentialBC() { return ess_data_; }
103  const DFSData& GetDFSData() const { return dfs_spaces_.GetDFSData(); }
104  void ShowError(const Vector &sol, bool verbose);
105  void VisualizeSolution(const Vector &sol, string tag);
106 };
107 
108 DarcyProblem::DarcyProblem(Mesh &mesh, int num_refs, int order,
109  const char *coef_file, Array<int> &ess_bdr,
110  DFSParameters dfs_param)
111  : mesh_(MPI_COMM_WORLD, mesh), ucoeff_(mesh.Dimension(), u_exact),
112  pcoeff_(p_exact), dfs_spaces_(order, num_refs, &mesh_, ess_bdr, dfs_param)
113 {
114  for (int l = 0; l < num_refs; l++)
115  {
116  mesh_.UniformRefinement();
117  dfs_spaces_.CollectDFSData();
118  }
119 
120  Vector coef_vector(mesh.GetNE());
121  coef_vector = 1.0;
122  if (std::strcmp(coef_file, ""))
123  {
124  ifstream coef_str(coef_file);
125  coef_vector.Load(coef_str, mesh.GetNE());
126  }
127  PWConstCoefficient mass_coeff(coef_vector);
128  VectorFunctionCoefficient fcoeff(mesh_.Dimension(), f_exact);
131 
132  u_.SetSpace(dfs_spaces_.GetHdivFES());
133  p_.SetSpace(dfs_spaces_.GetL2FES());
134  p_ = 0.0;
135  u_ = 0.0;
136  u_.ProjectBdrCoefficientNormal(ucoeff_, ess_bdr);
137 
138  ParLinearForm fform(dfs_spaces_.GetHdivFES());
140  fform.AddBoundaryIntegrator(new VectorFEBoundaryFluxLFIntegrator(natcoeff));
141  fform.Assemble();
142 
143  ParLinearForm gform(dfs_spaces_.GetL2FES());
144  gform.AddDomainIntegrator(new DomainLFIntegrator(gcoeff));
145  gform.Assemble();
146 
147  ParBilinearForm mVarf(dfs_spaces_.GetHdivFES());
148  ParMixedBilinearForm bVarf(dfs_spaces_.GetHdivFES(), dfs_spaces_.GetL2FES());
149 
150  mVarf.AddDomainIntegrator(new VectorFEMassIntegrator(mass_coeff));
151  mVarf.Assemble();
152  mVarf.EliminateEssentialBC(ess_bdr, u_, fform);
153  mVarf.Finalize();
154  M_.Reset(mVarf.ParallelAssemble());
155 
156  bVarf.AddDomainIntegrator(new VectorFEDivergenceIntegrator);
157  bVarf.Assemble();
158  bVarf.SpMat() *= -1.0;
159  bVarf.EliminateTrialDofs(ess_bdr, u_, gform);
160  bVarf.Finalize();
161  B_.Reset(bVarf.ParallelAssemble());
162 
163  rhs_.SetSize(M_->NumRows() + B_->NumRows());
164  Vector rhs_block0(rhs_.GetData(), M_->NumRows());
165  Vector rhs_block1(rhs_.GetData()+M_->NumRows(), B_->NumRows());
166  fform.ParallelAssemble(rhs_block0);
167  gform.ParallelAssemble(rhs_block1);
168 
169  ess_data_.SetSize(M_->NumRows() + B_->NumRows());
170  ess_data_ = 0.0;
171  Vector ess_data_block0(ess_data_.GetData(), M_->NumRows());
172  u_.ParallelProject(ess_data_block0);
173 
174  int order_quad = max(2, 2*order+1);
175  for (int i=0; i < Geometry::NumGeom; ++i)
176  {
177  irs_[i] = &(IntRules.Get(i, order_quad));
178  }
179 }
180 
181 void DarcyProblem::ShowError(const Vector& sol, bool verbose)
182 {
183  u_.Distribute(Vector(sol.GetData(), M_->NumRows()));
184  p_.Distribute(Vector(sol.GetData()+M_->NumRows(), B_->NumRows()));
185 
186  double err_u = u_.ComputeL2Error(ucoeff_, irs_);
187  double norm_u = ComputeGlobalLpNorm(2, ucoeff_, mesh_, irs_);
188  double err_p = p_.ComputeL2Error(pcoeff_, irs_);
189  double norm_p = ComputeGlobalLpNorm(2, pcoeff_, mesh_, irs_);
190 
191  if (!verbose) { return; }
192  cout << "|| u_h - u_ex || / || u_ex || = " << err_u / norm_u << "\n";
193  cout << "|| p_h - p_ex || / || p_ex || = " << err_p / norm_p << "\n";
194 }
195 
196 void DarcyProblem::VisualizeSolution(const Vector& sol, string tag)
197 {
198  int num_procs, myid;
199  MPI_Comm_size(mesh_.GetComm(), &num_procs);
200  MPI_Comm_rank(mesh_.GetComm(), &myid);
201 
202  u_.Distribute(Vector(sol.GetData(), M_->NumRows()));
203  p_.Distribute(Vector(sol.GetData()+M_->NumRows(), B_->NumRows()));
204 
205  const char vishost[] = "localhost";
206  const int visport = 19916;
207  socketstream u_sock(vishost, visport);
208  u_sock << "parallel " << num_procs << " " << myid << "\n";
209  u_sock.precision(8);
210  u_sock << "solution\n" << mesh_ << u_ << "window_title 'Velocity ("
211  << tag << " solver)'" << endl;
212  MPI_Barrier(mesh_.GetComm());
213  socketstream p_sock(vishost, visport);
214  p_sock << "parallel " << num_procs << " " << myid << "\n";
215  p_sock.precision(8);
216  p_sock << "solution\n" << mesh_ << p_ << "window_title 'Pressure ("
217  << tag << " solver)'" << endl;
218 }
219 
220 bool IsAllNeumannBoundary(const Array<int>& ess_bdr_attr)
221 {
222  for (int attr : ess_bdr_attr) { if (attr == 0) { return false; } }
223  return true;
224 }
225 
226 int main(int argc, char *argv[])
227 {
228  // Initialize MPI.
229  MPI_Session mpi(argc, argv);
230 
231  StopWatch chrono;
232  auto ResetTimer = [&chrono]() { chrono.Clear(); chrono.Start(); };
233 
234  // Parse command-line options.
235  const char *mesh_file = "../../data/beam-hex.mesh";
236  const char *coef_file = "";
237  const char *ess_bdr_attr_file = "";
238  int order = 0;
239  int par_ref_levels = 2;
240  bool show_error = false;
241  bool visualization = false;
242  DFSParameters param;
243  OptionsParser args(argc, argv);
244  args.AddOption(&mesh_file, "-m", "--mesh",
245  "Mesh file to use.");
246  args.AddOption(&order, "-o", "--order",
247  "Finite element order (polynomial degree).");
248  args.AddOption(&par_ref_levels, "-r", "--ref",
249  "Number of parallel refinement steps.");
250  args.AddOption(&coef_file, "-c", "--coef",
251  "Coefficient file to use.");
252  args.AddOption(&ess_bdr_attr_file, "-eb", "--ess-bdr",
253  "Essential boundary attribute file to use.");
254  args.AddOption(&show_error, "-se", "--show-error", "-no-se",
255  "--no-show-error",
256  "Show or not show approximation error.");
257  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
258  "--no-visualization",
259  "Enable or disable GLVis visualization.");
260  args.Parse();
261  if (!args.Good())
262  {
263  if (mpi.Root()) { args.PrintUsage(cout); }
264  return 1;
265  }
266  if (mpi.Root()) { args.PrintOptions(cout); }
267 
268  if (mpi.Root() && par_ref_levels == 0)
269  {
270  std::cout << "WARNING: DivFree solver is equivalent to BDPMinresSolver "
271  << "when par_ref_levels == 0.\n";
272  }
273 
274  // Initialize the mesh, boundary attributes, and solver parameters
275  Mesh *mesh = new Mesh(mesh_file, 1, 1);
276  int dim = mesh->Dimension();
277  int ser_ref_lvls = (int)ceil(log(mpi.WorldSize()/mesh->GetNE())/log(2.)/dim);
278  for (int i = 0; i < ser_ref_lvls; ++i)
279  {
280  mesh->UniformRefinement();
281  }
282 
283  Array<int> ess_bdr(mesh->bdr_attributes.Max());
284  ess_bdr = 0;
285  if (std::strcmp(ess_bdr_attr_file, ""))
286  {
287  ifstream ess_bdr_attr_str(ess_bdr_attr_file);
288  ess_bdr.Load(mesh->bdr_attributes.Max(), ess_bdr_attr_str);
289  }
290  if (IsAllNeumannBoundary(ess_bdr))
291  {
292  if (mpi.Root())
293  {
294  cout << "\nSolution is not unique when Neumann boundary condition is "
295  << "imposed on the entire boundary. \nPlease provide a different "
296  << "boundary condition.\n";
297  }
298  delete mesh;
299  return 0;
300  }
301 
302  string line = "**********************************************************\n";
303 
304  ResetTimer();
305 
306  // Generate components of the saddle point problem
307  DarcyProblem darcy(*mesh, par_ref_levels, order, coef_file, ess_bdr, param);
308  HypreParMatrix& M = darcy.GetM();
309  HypreParMatrix& B = darcy.GetB();
310  const DFSData& DFS_data = darcy.GetDFSData();
311  delete mesh;
312 
313  if (mpi.Root())
314  {
315  cout << line << "System assembled in " << chrono.RealTime() << "s.\n";
316  cout << "Dimension of the physical space: " << dim << "\n";
317  cout << "Size of the discrete Darcy system: " << M.M() + B.M() << "\n";
318  if (par_ref_levels > 0)
319  {
320  cout << "Dimension of the divergence free subspace: "
321  << DFS_data.C.back().Ptr()->NumCols() << "\n\n";
322  }
323  }
324 
325  // Setup various solvers for the discrete problem
326  std::map<const DarcySolver*, double> setup_time;
327  ResetTimer();
328  BDPMinresSolver bdp(M, B, param);
329  setup_time[&bdp] = chrono.RealTime();
330 
331  ResetTimer();
332  DivFreeSolver dfs_dm(M, B, DFS_data);
333  setup_time[&dfs_dm] = chrono.RealTime();
334 
335  ResetTimer();
336  const_cast<bool&>(DFS_data.param.coupled_solve) = true;
337  DivFreeSolver dfs_cm(M, B, DFS_data);
338  setup_time[&dfs_cm] = chrono.RealTime();
339 
340  std::map<const DarcySolver*, std::string> solver_to_name;
341  solver_to_name[&bdp] = "Block-diagonal-preconditioned MINRES";
342  solver_to_name[&dfs_dm] = "Divergence free (decoupled mode)";
343  solver_to_name[&dfs_cm] = "Divergence free (coupled mode)";
344 
345  // Solve the problem using all solvers
346  for (const auto& solver_pair : solver_to_name)
347  {
348  auto& solver = solver_pair.first;
349  auto& name = solver_pair.second;
350 
351  Vector sol = darcy.GetEssentialBC();
352 
353  ResetTimer();
354  solver->Mult(darcy.GetRHS(), sol);
355  chrono.Stop();
356 
357  if (mpi.Root())
358  {
359  cout << line << name << " solver:\n Setup time: "
360  << setup_time[solver] << "s.\n Solve time: "
361  << chrono.RealTime() << "s.\n Total time: "
362  << setup_time[solver] + chrono.RealTime() << "s.\n"
363  << " Iteration count: " << solver->GetNumIterations() <<"\n\n";
364  }
365  if (show_error && std::strcmp(coef_file, "") == 0)
366  {
367  darcy.ShowError(sol, mpi.Root());
368  }
369  else if (show_error && mpi.Root())
370  {
371  cout << "Exact solution is unknown for coefficient '" << coef_file
372  << "'.\nApproximation error is computed in this case!\n\n";
373  }
374 
375  if (visualization) { darcy.VisualizeSolution(sol, name); }
376 
377  }
378 
379  return 0;
380 }
381 
382 void u_exact(const Vector & x, Vector & u)
383 {
384  double xi(x(0));
385  double yi(x(1));
386  double zi(x.Size() == 3 ? x(2) : 0.0);
387 
388  u(0) = - exp(xi)*sin(yi)*cos(zi);
389  u(1) = - exp(xi)*cos(yi)*cos(zi);
390  if (x.Size() == 3)
391  {
392  u(2) = exp(xi)*sin(yi)*sin(zi);
393  }
394 }
395 
396 double p_exact(const Vector & x)
397 {
398  double xi(x(0));
399  double yi(x(1));
400  double zi(x.Size() == 3 ? x(2) : 0.0);
401  return exp(xi)*sin(yi)*cos(zi);
402 }
403 
404 void f_exact(const Vector & x, Vector & f)
405 {
406  f = 0.0;
407 }
408 
409 double g_exact(const Vector & x)
410 {
411  if (x.Size() == 3) { return -p_exact(x); }
412  return 0;
413 }
414 
415 double natural_bc(const Vector & x)
416 {
417  return (-p_exact(x));
418 }
int WorldSize() const
Return MPI_COMM_WORLD&#39;s size.
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:97
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
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:920
static const int NumGeom
Definition: geom.hpp:41
bool IsAllNeumannBoundary(const Array< int > &ess_bdr_attr)
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
double RealTime()
Definition: tic_toc.cpp:426
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:199
HYPRE_BigInt M() const
Returns the global number of rows.
Definition: hypre.hpp:521
void Stop()
Stop the stopwatch.
Definition: tic_toc.cpp:416
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds a domain integrator. Assumes ownership of bfi.
Data for the divergence free solver.
Class for parallel linear form.
Definition: plinearform.hpp:26
A simple convenience class that calls MPI_Init() at construction and MPI_Finalize() at destruction...
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
Wrapper for the block-diagonal-preconditioned MINRES defined in ex5p.cpp.
constexpr char vishost[]
double natural_bc(const Vector &x)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
constexpr int visport
Timing object.
Definition: tic_toc.hpp:34
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
bool Root() const
Return true if WorldRank() == 0.
int Dimension() const
Definition: mesh.hpp:911
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
void Start()
Clear the elapsed time and start the stopwatch.
Definition: tic_toc.cpp:411
double p_exact(const Vector &x)
Definition: ex24.cpp:380
A general vector function coefficient.
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:204
void u_exact(const Vector &x, Vector &u)
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 f_exact(const Vector &, Vector &)
Definition: ex3.cpp:257
Parameters for the divergence free solver.
double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
A piecewise constant coefficient with the constants keyed off the element attribute numbers...
Definition: coefficient.hpp:94
std::vector< OperatorPtr > C
Class for parallel bilinear form using different test and trial FE spaces.
int dim
Definition: ex24.cpp:53
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
Class for parallel bilinear form.
double g_exact(const Vector &x)
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:266
A general function coefficient.
Vector data type.
Definition: vector.hpp:60
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:277
Class for parallel meshes.
Definition: pmesh.hpp:32
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:377
int main()
void Clear()
Clear the elapsed time on the stopwatch and restart it if it&#39;s running.
Definition: tic_toc.cpp:406
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150