MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
div_free_solver.hpp
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#ifndef MFEM_DIVFREE_SOLVER_HPP
13#define MFEM_DIVFREE_SOLVER_HPP
14
15#include "darcy_solver.hpp"
16
17namespace mfem
18{
19namespace blocksolvers
20{
21/// Parameters for the divergence free solver
23{
24 /** There are three components in the solver: a particular solution
25 satisfying the divergence constraint, the remaining div-free component of
26 the flux, and the pressure. When coupled_solve == false, the three
27 components will be solved one by one in the aforementioned order.
28 Otherwise, they will be solved at the same time. */
29 bool coupled_solve = false;
30 bool verbose = false;
33};
34
35/// Data for the divergence free solver
36struct DFSData
37{
38 std::vector<OperatorPtr> agg_hdivdof; // agglomerates to H(div) dofs table
39 std::vector<OperatorPtr> agg_l2dof; // agglomerates to L2 dofs table
40 std::vector<OperatorPtr> P_hdiv; // Interpolation matrix for H(div) space
41 std::vector<OperatorPtr> P_l2; // Interpolation matrix for L2 space
42 std::vector<OperatorPtr> P_hcurl; // Interpolation for kernel space of div
43 std::vector<OperatorPtr> Q_l2; // Q_l2[l] = (W_{l+1})^{-1} P_l2[l]^T W_l
44 Array<int> coarsest_ess_hdivdofs; // coarsest level essential H(div) dofs
45 std::vector<OperatorPtr> C; // discrete curl: ND -> RT, map to Null(B)
47};
48
49/// Finite element spaces concerning divergence free solvers
50/// The main usage of this class is to collect data needed for the solver.
52{
53 RT_FECollection hdiv_fec_;
54 L2_FECollection l2_fec_;
55 std::unique_ptr<FiniteElementCollection> hcurl_fec_;
56 L2_FECollection l2_0_fec_;
57
58 std::unique_ptr<ParFiniteElementSpace> coarse_hdiv_fes_;
59 std::unique_ptr<ParFiniteElementSpace> coarse_l2_fes_;
60 std::unique_ptr<ParFiniteElementSpace> coarse_hcurl_fes_;
61 std::unique_ptr<ParFiniteElementSpace> l2_0_fes_;
62
63 std::unique_ptr<ParFiniteElementSpace> hdiv_fes_;
64 std::unique_ptr<ParFiniteElementSpace> l2_fes_;
65 std::unique_ptr<ParFiniteElementSpace> hcurl_fes_;
66
67 std::vector<SparseMatrix> el_l2dof_;
68 const Array<int>& ess_bdr_attr_;
69 Array<int> all_bdr_attr_;
70
71 int level_;
72 DFSData data_;
73
74 void MakeDofRelationTables(int level);
75 void DataFinalize();
76public:
77 DFSSpaces(int order, int num_refine, ParMesh *mesh,
78 const Array<int>& ess_attr, const DFSParameters& param);
79
80 /** This should be called each time when the mesh (where the FE spaces are
81 defined) is refined. The spaces will be updated, and the prolongation for
82 the spaces and other data needed for the div-free solver are stored. */
83 void CollectDFSData();
84
85 const DFSData& GetDFSData() const { return data_; }
86 ParFiniteElementSpace* GetHdivFES() const { return hdiv_fes_.get(); }
87 ParFiniteElementSpace* GetL2FES() const { return l2_fes_.get(); }
88};
89
90/// Solvers for DFS
91/// Solver for B * B^T
92/// Compute the product B * B^T and solve it with CG preconditioned by BoomerAMG
93class BBTSolver : public Solver
94{
95 OperatorPtr BBT_;
96 OperatorPtr BBT_prec_;
97 CGSolver BBT_solver_;
98public:
100 virtual void Mult(const Vector &x, Vector &y) const { BBT_solver_.Mult(x, y); }
101 virtual void SetOperator(const Operator &op) { }
102};
103
104/// Block diagonal solver for symmetric A, each block is inverted by direct solver
106{
107public:
109 : DirectSubBlockSolver(A, block_dof) { }
110 virtual void MultTranspose(const Vector &x, Vector &y) const { Mult(x, y); }
111};
112
113/// non-overlapping additive Schwarz smoother for saddle point systems
114/// [ M B^T ]
115/// [ B 0 ]
117{
118 const SparseMatrix& agg_hdivdof_;
119 const SparseMatrix& agg_l2dof_;
120 OperatorPtr coarse_l2_projector_;
121
122 Array<int> offsets_;
123 mutable Array<int> offsets_loc_;
124 mutable Array<int> hdivdofs_loc_;
125 mutable Array<int> l2dofs_loc_;
126 std::vector<OperatorPtr> solvers_loc_;
127public:
128 /** SaddleSchwarzSmoother solves local saddle point problems defined on a
129 list of non-overlapping aggregates (of elements).
130 @param M the [1,1] block of the saddle point system
131 @param B the [2,1] block of the saddle point system
132 @param agg_hdivdof aggregate to H(div) dofs relation table (boolean matrix)
133 @param agg_l2dof aggregate to L2 dofs relation table (boolean matrix)
134 @param P_l2 prolongation matrix of the discrete L2 space
135 @param Q_l2 projection matrix of the discrete L2 space:
136 Q_l2 := (P_l2 W P_l2)^{-1} * P_l2 * W,
137 where W is the mass matrix of the discrete L2 space. */
139 const HypreParMatrix& B,
140 const SparseMatrix& agg_hdivdof,
141 const SparseMatrix& agg_l2dof,
142 const HypreParMatrix& P_l2,
143 const HypreParMatrix& Q_l2);
144 virtual void Mult(const Vector &x, Vector &y) const;
145 virtual void MultTranspose(const Vector &x, Vector &y) const { Mult(x, y); }
146 virtual void SetOperator(const Operator &op) { }
147};
148
149/// Solver for local problems in SaddleSchwarzSmoother
150class LocalSolver : public Solver
151{
152 DenseMatrix local_system_;
153 DenseMatrixInverse local_solver_;
154 const int offset_;
155public:
156 LocalSolver(const DenseMatrix &M, const DenseMatrix &B);
157 virtual void Mult(const Vector &x, Vector &y) const;
158 virtual void SetOperator(const Operator &op) { }
159};
160
161/// Divergence free solver.
162/** Divergence free solver.
163 The basic idea of the solver is to exploit a multilevel decomposition of
164 Raviart-Thomas space to find a particular solution satisfying the divergence
165 constraint, and then solve the remaining (divergence-free) component in the
166 kernel space of the discrete divergence operator.
167
168 For more details see
169 1. Vassilevski, Multilevel Block Factorization Preconditioners (Appendix F.3),
170 Springer, 2008.
171 2. Voronin, Lee, Neumuller, Sepulveda, Vassilevski, Space-time discretizations
172 using constrained first-order system least squares (CFOSLS).
173 J. Comput. Phys. 373: 863-876, 2018. */
175{
176 const DFSData& data_;
177 DFSParameters param_;
178 OperatorPtr BT_;
179 BBTSolver BBT_solver_;
180 std::vector<Array<int>> ops_offsets_;
182 Array<BlockOperator*> blk_Ps_;
183 Array<Solver*> smoothers_;
184 OperatorPtr prec_;
185 OperatorPtr solver_;
186
187 void SolveParticular(const Vector& rhs, Vector& sol) const;
188 void SolveDivFree(const Vector& rhs, Vector& sol) const;
189 void SolvePotential(const Vector &rhs, Vector& sol) const;
190public:
192 const DFSData& data);
194 virtual void Mult(const Vector &x, Vector &y) const;
195 virtual void SetOperator(const Operator &op) { }
196 virtual int GetNumIterations() const;
197};
198
199} // namespace blocksolvers
200
201} // namespace mfem
202
203#endif // MFEM_DIVFREE_SOLVER_HPP
Conjugate gradient method.
Definition solvers.hpp:513
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the Conjugate Gradient method.
Definition solvers.cpp:718
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
Block diagonal solver for A, each block is inverted by direct solver.
Definition solvers.hpp:1174
virtual void Mult(const Vector &x, Vector &y) const
Direct solution of the block diagonal linear system.
Definition solvers.cpp:3375
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
Pointer to an Operator of a specified type.
Definition handle.hpp:34
Abstract operator.
Definition operator.hpp:25
Abstract parallel finite element space.
Definition pfespace.hpp:29
Class for parallel meshes.
Definition pmesh.hpp:34
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition fe_coll.hpp:386
Base class for solvers.
Definition operator.hpp:683
Data type sparse matrix.
Definition sparsemat.hpp:51
Vector data type.
Definition vector.hpp:80
BBTSolver(const HypreParMatrix &B, IterSolveParameters param)
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
DFSSpaces(int order, int num_refine, ParMesh *mesh, const Array< int > &ess_attr, const DFSParameters &param)
const DFSData & GetDFSData() const
ParFiniteElementSpace * GetHdivFES() const
ParFiniteElementSpace * GetL2FES() const
Abstract solver class for Darcy's flow.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
DivFreeSolver(const HypreParMatrix &M, const HypreParMatrix &B, const DFSData &data)
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Solver for local problems in SaddleSchwarzSmoother.
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
LocalSolver(const DenseMatrix &M, const DenseMatrix &B)
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
SaddleSchwarzSmoother(const HypreParMatrix &M, const HypreParMatrix &B, const SparseMatrix &agg_hdivdof, const SparseMatrix &agg_l2dof, const HypreParMatrix &P_l2, const HypreParMatrix &Q_l2)
Block diagonal solver for symmetric A, each block is inverted by direct solver.
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
SymDirectSubBlockSolver(const SparseMatrix &A, const SparseMatrix &block_dof)
Data for the divergence free solver.
std::vector< OperatorPtr > agg_l2dof
std::vector< OperatorPtr > P_hcurl
std::vector< OperatorPtr > P_hdiv
std::vector< OperatorPtr > P_l2
std::vector< OperatorPtr > Q_l2
std::vector< OperatorPtr > agg_hdivdof
std::vector< OperatorPtr > C
Parameters for the divergence free solver.