MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
bramble_pasciak.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#include "bramble_pasciak.hpp"
14
15using namespace std;
16
17namespace mfem
18{
19namespace blocksolvers
20{
21/// Bramble-Pasciak Solver
23 ParBilinearForm *mVarf,
25 const BPSParameters &param)
26 : DarcySolver(mVarf->ParFESpace()->GetTrueVSize(),
27 bVarf->TestFESpace()->GetTrueVSize())
28{
29 M_.reset(mVarf->ParallelAssemble());
30 B_.reset(bVarf->ParallelAssemble());
31 Q_.reset(ConstructMassPreconditioner(*mVarf, param.q_scaling));
32
33 Vector diagM;
34 M_->GetDiag(diagM);
35 auto BT = B_->Transpose();
36 auto invDBt = new HypreParMatrix(*BT);
37 invDBt->InvScaleRows(diagM);
38 auto S = ParMult(B_.get(), invDBt);
39 M0_.Reset(new HypreDiagScale(*M_));
40 M1_.Reset(new HypreBoomerAMG(*S));
41 M1_.As<HypreBoomerAMG>()->SetPrintLevel(0);
42
43 Init(*M_, *B_, *Q_, *M0_.As<Solver>(), *M1_.As<Solver>(), param);
44}
45
48 Solver &M0, Solver &M1,
49 const BPSParameters &param)
50 : DarcySolver(M.NumRows(), B.NumRows())
51{
52 Init(M, B, Q, M0, M1, param);
53}
54
55void BramblePasciakSolver::Init(
57 Solver &M0, Solver &M1,
58 const BPSParameters &param)
59{
60 auto Bt = new TransposeOperator(&B);
61 auto invQ = new HypreDiagScale(Q);
62
63 use_bpcg = param.use_bpcg;
64
65 if (use_bpcg)
66 {
67 oop_ = new BlockOperator(offsets_);
68 oop_->owns_blocks = false;
69 oop_->SetBlock(0, 0, &M);
70 oop_->SetBlock(0, 1, Bt);
71 oop_->SetBlock(1, 0, &B);
72
73 // cpc_ unused in bpcg
74 auto temp_cpc = new BlockDiagonalPreconditioner(offsets_);
75 temp_cpc->owns_blocks = true;
76 temp_cpc->SetDiagonalBlock(0, invQ);
77 temp_cpc->SetDiagonalBlock(1, &M1);
78 // tri(1,0) = B M0 = B invQ
79 auto id_m = new IdentityOperator(M.NumRows());
80 auto id_b = new IdentityOperator(B.NumRows());
81 auto BinvQ = new ProductOperator(&B, invQ, false, false);
82 // tri
83 auto temp_tri = new BlockOperator(offsets_);
84 temp_tri->owns_blocks = true;
85 temp_tri->SetBlock(0, 0, id_m);
86 temp_tri->SetBlock(1, 1, id_b, -1.0);
87 temp_tri->SetBlock(1, 0, BinvQ);
88
89 ppc_ = new ProductOperator(temp_cpc, temp_tri, true, true);
90
91 ipc_ = new BlockOperator(offsets_);
92 ipc_->owns_blocks = false;
93 ipc_->SetDiagonalBlock(0, invQ);
94
95 // bpcg
96 solver_.reset(new BPCGSolver(M.GetComm(), *ipc_, *ppc_));
97 solver_->SetOperator(*oop_);
98 }
99 else
100 {
101 // oop_ unused in cg
102 auto temp_oop = new BlockOperator(offsets_);
103 temp_oop->owns_blocks = false;
104 temp_oop->SetBlock(0, 0, &M);
105 temp_oop->SetBlock(0, 1, Bt);
106 temp_oop->SetBlock(1, 0, &B);
107
108 // ipc_ unused in cg
109 auto temp_ipc = new BlockOperator(offsets_);
110 temp_ipc->owns_blocks = false;
111 temp_ipc->SetDiagonalBlock(0, invQ);
112
113 // temp_AN = temp_oop * temp_ipc
114 auto temp_AN = new ProductOperator(temp_oop, temp_ipc, true, true);
115
116 // Required for updating the RHS
117 auto id = new IdentityOperator(M.NumRows()+B.NumRows());
118 map_ = new SumOperator(temp_AN, 1.0, id, -1.0, true, true);
119
120 mop_ = new ProductOperator(map_, temp_oop, false, true);
121
123 cpc_->owns_blocks = true;
124 cpc_->SetDiagonalBlock(0, &M0);
125 cpc_->SetDiagonalBlock(1, &M1);
126
127 // (P)CG
128 solver_.reset(new CGSolver(M.GetComm()));
129 solver_->SetOperator(*mop_);
130 solver_->SetPreconditioner(*cpc_);
131 }
132 SetOptions(*solver_, param);
133}
134
136 ParBilinearForm &mVarf, real_t q_scaling)
137{
138 MFEM_ASSERT((q_scaling > 0.0) && (q_scaling < 1.0),
139 "Invalid Q-scaling factor: q_scaling = " << q_scaling );
140 ParBilinearForm qVarf(mVarf.ParFESpace());
141 qVarf.AllocateMatrix();
142#ifndef MFEM_USE_LAPACK
143 if (Mpi::Root())
144 {
145 mfem::out << "Warning: Using inverse power method to compute the minimum "
146 << "eigenvalue of the small eigenvalue problem.\n";
147 mfem::out << " Consider compiling MFEM with LAPACK support.\n";
148 }
149#endif
150 for (int i = 0; i < mVarf.ParFESpace()->GetNE(); ++i)
151 {
152 DenseMatrix M_i, Q_i;
153 Vector diag_i;
154 real_t scaling = 0.0, eval_i = 0.0;
155 mVarf.ComputeElementMatrix(i, M_i);
156 M_i.GetDiag(diag_i);
157 // M_i <- D^{-1/2} M_i D^{-1/2}, where D = diag(M_i)
158 M_i.InvSymmetricScaling(diag_i);
159 // M_i x = ev diag(M_i) x
160#ifdef MFEM_USE_LAPACK
161 DenseMatrix evec;
162 Vector eval;
163 M_i.Eigenvalues(eval, evec);
164 eval_i = eval.Min();
165#else
166 // Inverse power method
167 Vector x(M_i.Height()), Mx(M_i.Height()), diff(M_i.Height());
168 real_t eval_prev = 0.0;
169 int iter = 0;
170 x.Randomize(696383552+779345*i);
171#if defined(MFEM_USE_DOUBLE)
172 const real_t rel_tol = 1e-12;
173#elif defined(MFEM_USE_SINGLE)
174 const real_t rel_tol = 1e-6;
175#else
176#error "Only single and double precision are supported!"
177 const real_t rel_tol = 1e-12;
178#endif
179 DenseMatrixInverse M_i_inv(M_i);
180 do
181 {
182 eval_prev = eval_i;
183 M_i_inv.Mult(x, Mx);
184 eval_i = Mx.Norml2();
185 x.Set(1.0/eval_i, Mx);
186 ++iter;
187 }
188 while ((iter < 1000) && (fabs(eval_i - eval_prev)/fabs(eval_i) > rel_tol));
189 MFEM_VERIFY(fabs(eval_i - eval_prev)/fabs(eval_i) <= rel_tol,
190 "Inverse power method did not converge."
191 << "\n\t iter = " << iter
192 << "\n\t eval_i = " << eval_i
193 << "\n\t eval_prev = " << eval_prev
194 << "\n\t fabs(eval_i - eval_prev)/fabs(eval_i) = "
195 << fabs(eval_i - eval_prev)/fabs(eval_i));
196 eval_i = 1.0/eval_i;
197#endif
198 scaling = q_scaling*eval_i;
199 diag_i.Set(scaling, diag_i);
200 Q_i.Diag(diag_i.GetData(), diag_i.Size());
201 qVarf.AssembleElementMatrix(i, Q_i, 1);
202 }
203 qVarf.Finalize();
204 return qVarf.ParallelAssemble();
205}
206
207void BramblePasciakSolver::Mult(const Vector & x, Vector & y) const
208{
209 if (!use_bpcg)
210 {
211 Vector transformed_rhs(x.Size());
212 map_->Mult(x, transformed_rhs);
213 solver_->Mult(transformed_rhs, y);
214 }
215 else
216 {
217 solver_->Mult(x, y);
218 }
219 for (int dof : ess_zero_dofs_) { y[dof] = 0.0; }
220}
221
222/// Bramble-Pasciak CG
224{
226
227 r.SetSize(width, mt); r.UseDevice(true);
228 p.SetSize(width, mt); p.UseDevice(true);
229 g.SetSize(width, mt); g.UseDevice(true);
230 t.SetSize(width, mt); t.UseDevice(true);
231 r_bar.SetSize(width, mt); r_bar.UseDevice(true);
232 r_red.SetSize(width, mt); r_red.UseDevice(true);
233 g_red.SetSize(width, mt); g_red.UseDevice(true);
234}
235
236void BPCGSolver::Mult(const Vector &b, Vector &x) const
237{
238 int i;
239 real_t delta, delta0, del0;
240 real_t alpha, beta, gamma;
241
242 // Initialization
243 x.UseDevice(true);
244 if (iterative_mode)
245 {
246 oper->Mult(x, r);
247 subtract(b, r, r); // r = b - A x
248 }
249 else
250 {
251 r = b;
252 x = 0.0;
253 }
254
255 pprec->Mult(r,r_bar); // r_bar = P r
256 p = r_bar;
257 oper->Mult(p, g); // g = A p
258 oper->Mult(r_bar, t); // t = A r_bar
259 iprec->Mult(r, r_red); // r_red = N r
260
261 delta = delta0 = Dot(t, r_red) - Dot(r_bar, r); // Dot(Pr, r)
262 if (delta0 >= 0.0) { initial_norm = sqrt(delta0); }
263 MFEM_ASSERT(IsFinite(delta), "norm = " << delta);
265 {
266 mfem::out << " Iteration : " << setw(3) << 0 << " (P r, r) = "
267 << delta << (print_options.first_and_last ? " ...\n" : "\n");
268 }
269 Monitor(0, delta, r, x);
270
271 if (delta < 0.0)
272 {
274 {
275 mfem::out << "BPCG: The preconditioner is not positive definite. (Pr, r) = "
276 << delta << '\n';
277 }
278 converged = false;
279 final_iter = 0;
282 return;
283 }
284 del0 = std::max(delta*rel_tol*rel_tol, abs_tol*abs_tol);
285 if (delta <= del0)
286 {
287 converged = true;
288 final_iter = 0;
289 final_norm = sqrt(delta);
290 return;
291 }
292
293 iprec->Mult(g, g_red);
294 gamma = Dot(g, g_red) - Dot(g,p); // Dot(Ap, p)
295 MFEM_ASSERT(IsFinite(gamma), "den (gamma) = " << gamma);
296 if (gamma <= 0.0)
297 {
298 if (Dot(r_bar, r_bar) > 0.0 && print_options.warnings)
299 {
300 mfem::out << "BPCG: The operator is not positive definite. (Ar, r) = "
301 << gamma << '\n';
302 }
303 if (gamma == 0.0)
304 {
305 converged = false;
306 final_iter = 0;
307 final_norm = sqrt(delta);
308 return;
309 }
310 }
311
312 // Start iteration
313 converged = false;
315 for (i = 1; true; )
316 {
317 alpha = delta0/gamma;
318 add(x, alpha, p, x); // x = x + alpha p
319 add(r, -alpha, g, r); // r = r - alpha g
320
321 pprec->Mult(r, r_bar); // r_bar = P r
322 iprec->Mult(r, r_red); // r_red = N r
323 oper->Mult(r_bar, t); // t = A r_bar
324 delta = Dot(t, r_red) - Dot(r_bar,r);
325
326 // Check
327 MFEM_ASSERT(IsFinite(delta), "norm = " << delta);
328 if (delta < 0.0)
329 {
331 {
332 mfem::out << "BPCG: The preconditioner is not positive definite. (Pr, r) = "
333 << delta << '\n';
334 }
335 converged = false;
336 final_iter = i;
337 break;
338 }
340 {
341 mfem::out << " Iteration : " << setw(3) << i << " (Pr, r) = "
342 << delta << std::endl;
343 }
344 Monitor(i, delta, r, x);
345 if (delta <= del0)
346 {
347 converged = true;
348 final_iter = i;
349 break;
350 }
351 if (++i > max_iter)
352 {
353 break;
354 }
355 // End check
356
357 beta = delta/delta0;
358 add(r_bar, beta, p, p); // p = r_bar + beta p
359 add(t, beta, g, g); // g = t + beta g
360
361 delta0 = delta;
362 iprec->Mult(g, g_red);
363 gamma = Dot(g, g_red) - Dot(g,p); // Dot(Ap, p)
364 MFEM_ASSERT(IsFinite(gamma), "den (gamma) = " << gamma);
365 if (gamma <= 0.0)
366 {
367 if (Dot(r_bar, r_bar) > 0.0 && print_options.warnings)
368 {
369 mfem::out << "BPCG: The operator is not positive definite. (Ar, r) = "
370 << gamma << '\n';
371 }
372 if (gamma == 0.0)
373 {
374 final_iter = i;
375 break;
376 }
377 }
378 }
379
381 {
382 mfem::out << " Iteration : " << setw(3) << final_iter << " (Pr, r) = "
383 << delta << '\n';
384 }
386 {
387 mfem::out << "BPCG: Number of iterations: " << final_iter << '\n';
388 }
391 {
392 const auto arf = pow (gamma/delta0, 0.5/final_iter);
393 mfem::out << "Average reduction factor = " << arf << '\n';
394 }
396 {
397 mfem::out << "BPCG: No convergence!" << '\n';
398 }
399
400 final_norm = sqrt(delta);
401 Monitor(final_iter, final_norm, r, x, true);
402}
403} // namespace blocksolvers
404} // namespace mfem
void AllocateMatrix()
Pre-allocate the internal SparseMatrix before assembly. If the internal flag precompute_sparsity is s...
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY....
void AssembleElementMatrix(int i, const DenseMatrix &elmat, int skip_zeros=1)
Assemble the given element matrix.
void ComputeElementMatrix(int i, DenseMatrix &elmat)
Compute the element matrix of the given element.
A class to handle Block diagonal preconditioners in a matrix-free implementation.
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
A class to handle Block systems in a matrix-free implementation.
void SetDiagonalBlock(int iblock, Operator *op, real_t c=1.0)
Add block op in the block-entry (iblock, iblock).
void SetBlock(int iRow, int iCol, Operator *op, real_t c=1.0)
Add a block op in the block-entry (iblock, jblock).
Conjugate gradient method.
Definition solvers.hpp:513
void Mult(const real_t *x, real_t *y) const
Matrix vector multiplication with the inverse of dense matrix.
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void GetDiag(Vector &d) const
Returns the diagonal of the matrix.
void Eigenvalues(Vector &ev)
Compute eigenvalues of A x = ev x where A = *this.
Definition densemat.hpp:272
void Diag(real_t c, int n)
Creates n x n diagonal matrix with diagonal elements c.
void InvSymmetricScaling(const Vector &s)
InvSymmetricScaling this = diag(sqrt(1./s)) * this * diag(sqrt(1./s))
Definition densemat.cpp:484
The BoomerAMG solver in hypre.
Definition hypre.hpp:1691
Jacobi preconditioner in hypre.
Definition hypre.hpp:1481
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
MPI_Comm GetComm() const
MPI communicator.
Definition hypre.hpp:578
Identity Operator I: x -> x.
Definition operator.hpp:706
A class to initialize the size of a Tensor.
Definition dtensor.hpp:55
real_t abs_tol
Absolute tolerance.
Definition solvers.hpp:158
void Monitor(int it, real_t norm, const Vector &r, const Vector &x, bool final=false) const
Monitor both the residual r and the solution x.
Definition solvers.cpp:190
real_t rel_tol
Relative tolerance.
Definition solvers.hpp:155
PrintLevel print_options
Output behavior for the iterative solver.
Definition solvers.hpp:138
virtual real_t Dot(const Vector &x, const Vector &y) const
Return the standard (l2, i.e., Euclidean) inner product of x and y.
Definition solvers.cpp:55
const Operator * oper
Definition solvers.hpp:121
int max_iter
Limit for the number of iterations the solver is allowed to do.
Definition solvers.hpp:152
static bool Root()
Return true if the rank in MPI_COMM_WORLD is zero.
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
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
Definition operator.hpp:86
int width
Dimension of the input / number of columns in the matrix.
Definition operator.hpp:28
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
int NumRows() const
Get the number of rows (size of output) of the Operator. Synonym with Height().
Definition operator.hpp:69
Class for parallel bilinear form.
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
ParFiniteElementSpace * ParFESpace() const
Return the parallel FE space associated with the ParBilinearForm.
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.
General product operator: x -> (A*B)(x) = A(B(x)).
Definition operator.hpp:797
Base class for solvers.
Definition operator.hpp:683
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition operator.hpp:686
General linear combination operator: x -> a A(x) + b B(x).
Definition operator.hpp:774
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition operator.hpp:786
The transpose of a given operator. Switches the roles of the methods Mult() and MultTranspose().
Definition operator.hpp:750
Vector data type.
Definition vector.hpp:80
Vector & Set(const real_t a, const Vector &x)
(*this) = a * x
Definition vector.cpp:262
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition vector.hpp:136
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
real_t Min() const
Returns the minimal element of the vector.
Definition vector.cpp:1212
Bramble-Pasciak Conjugate Gradient.
void UpdateVectors()
Bramble-Pasciak CG.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
static HypreParMatrix * ConstructMassPreconditioner(ParBilinearForm &mVarf, real_t alpha=0.5)
Assemble a preconditioner for the mass matrix.
BramblePasciakSolver(ParBilinearForm *mVarf, ParMixedBilinearForm *bVarf, const BPSParameters &param)
System and mass preconditioner are constructed from bilinear forms.
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Abstract solver class for Darcy's flow.
Vector beta
const real_t alpha
Definition ex15.cpp:369
real_t b
Definition lissajous.cpp:42
real_t delta
Definition lissajous.cpp:43
void SetOptions(IterativeSolver &solver, const IterSolveParameters &param)
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
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition vector.cpp:316
bool IsFinite(const real_t &val)
Definition vector.hpp:507
MemoryType GetMemoryType(MemoryClass mc)
Return a suitable MemoryType for a given MemoryClass.
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
Definition hypre.cpp:2949
void subtract(const Vector &x, const Vector &y, Vector &z)
Definition vector.cpp:472
float real_t
Definition config.hpp:43
MemoryType
Memory types supported by MFEM.
bool iterations
Detailed information about each iteration will be reported to mfem::out.
Definition solvers.hpp:88
bool warnings
If a non-fatal problem has been detected some context-specific information will be reported to mfem::...
Definition solvers.hpp:85
bool first_and_last
Information about the first and last iteration will be printed to mfem::out.
Definition solvers.hpp:94
bool summary
A summary of the solver process will be reported after the last iteration to mfem::out.
Definition solvers.hpp:91
Parameters for the BramblePasciakSolver method.