MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
dgmassinv.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#include "dgmassinv.hpp"
13#include "bilinearform.hpp"
14#include "dgmassinv_kernels.hpp"
15#include "../general/forall.hpp"
16
17namespace mfem
18{
19
21 const IntegrationRule *ir,
22 int btype)
23 : Solver(fes_orig.GetTrueVSize()),
24 fec(fes_orig.GetMaxElementOrder(),
25 fes_orig.GetMesh()->Dimension(),
26 btype,
27 fes_orig.GetFE(0)->GetMapType()),
28 fes(fes_orig.GetMesh(), &fec)
29{
30 MFEM_VERIFY(fes.IsDGSpace(), "Space must be DG.");
31 MFEM_VERIFY(!fes.IsVariableOrder(), "Variable orders not supported.");
32
33 const int btype_orig =
34 static_cast<const L2_FECollection*>(fes_orig.FEColl())->GetBasisType();
35
36 if (btype_orig == btype)
37 {
38 // No change of basis required
39 d2q = nullptr;
40 }
41 else
42 {
43 // original basis to solver basis
44 const auto mode = DofToQuad::TENSOR;
45 d2q = &fes_orig.GetFE(0)->GetDofToQuad(fes.GetFE(0)->GetNodes(), mode);
46
47 int n = d2q->ndof;
48 Array<real_t> B_inv = d2q->B; // deep copy
49 Array<int> ipiv(n);
50 // solver basis to original
51 LUFactors lu(B_inv.HostReadWrite(), ipiv.HostWrite());
52 lu.Factor(n);
53 B_.SetSize(n*n);
55 Bt_.SetSize(n*n);
56 DenseMatrix B_matrix(B_.HostReadWrite(), n, n);
57 DenseMatrix Bt_matrix(Bt_.HostWrite(), n, n);
58 Bt_matrix.Transpose(B_matrix);
59 }
60
61 if (coeff) { m = new MassIntegrator(*coeff, ir); }
62 else { m = new MassIntegrator(ir); }
63
65 // Workspace vectors used for CG
69 // Only need transformed RHS if basis is different
70 if (btype_orig != btype) { b2_.SetSize(height); }
71
72 M = new BilinearForm(&fes);
73 M->AddDomainIntegrator(m); // M assumes ownership of m
75
76 // Assemble the bilinear form and its diagonal (for preconditioning).
77 Update();
78}
79
81 int btype)
82 : DGMassInverse(fes_, &coeff, nullptr, btype) { }
83
85 const IntegrationRule &ir, int btype)
86 : DGMassInverse(fes_, &coeff, &ir, btype) { }
87
89 const IntegrationRule &ir, int btype)
90 : DGMassInverse(fes_, nullptr, &ir, btype) { }
91
93 : DGMassInverse(fes_, nullptr, nullptr, btype) { }
94
96{
97 MFEM_ABORT("SetOperator not supported with DGMassInverse.")
98}
99
100void DGMassInverse::SetRelTol(const real_t rel_tol_) { rel_tol = rel_tol_; }
101
102void DGMassInverse::SetAbsTol(const real_t abs_tol_) { abs_tol = abs_tol_; }
103
104void DGMassInverse::SetMaxIter(const real_t max_iter_) { max_iter = max_iter_; }
105
112
114{
115 delete M;
116}
117
118template<int DIM, int D1D, int Q1D>
120{
121 using namespace internal; // host/device kernel functions
122
123 const int NE = fes.GetNE();
124 const int d1d = m->dofs1D;
125 const int q1d = m->quad1D;
126
127 const int ND = static_cast<int>(pow(d1d, DIM));
128
129 const auto B = m->maps->B.Read();
130 const auto Bt = m->maps->Bt.Read();
131 const auto pa_data = m->pa_data.Read();
132 const auto dinv = diag_inv.Read();
133 auto r = r_.Write();
134 auto d = d_.Write();
135 auto z = z_.Write();
136 auto u = u_.ReadWrite();
137
138 const real_t RELTOL = rel_tol;
139 const real_t ABSTOL = abs_tol;
140 const real_t MAXIT = max_iter;
141 const bool IT_MODE = iterative_mode;
142 const bool CHANGE_BASIS = (d2q != nullptr);
143
144 // b is the right-hand side (if no change of basis, this just points to the
145 // incoming RHS vector, if we have to change basis, this points to the
146 // internal b2 vector where we put the transformed RHS)
147 const real_t *b;
148 // the following are non-null if we have to change basis
149 real_t *b2 = nullptr; // non-const access to b2
150 const real_t *b_orig = nullptr; // RHS vector in "original" basis
151 const real_t *d2q_B = nullptr; // matrix to transform initial guess
152 const real_t *q2d_B = nullptr; // matrix to transform solution
153 const real_t *q2d_Bt = nullptr; // matrix to transform RHS
154 if (CHANGE_BASIS)
155 {
156 d2q_B = d2q->B.Read();
157 q2d_B = B_.Read();
158 q2d_Bt = Bt_.Read();
159
160 b2 = b2_.Write();
161 b_orig = b_.Read();
162 b = b2;
163 }
164 else
165 {
166 b = b_.Read();
167 }
168
169 static constexpr int NB = Q1D ? Q1D : 1; // block size
170
171 mfem::forall_2D(NE, NB, NB, [=] MFEM_HOST_DEVICE (int e)
172 {
173 // Perform change of basis if needed
174 if (CHANGE_BASIS)
175 {
176 // Transform RHS
177 DGMassBasis<DIM,D1D>(e, NE, q2d_Bt, b_orig, b2, d1d);
178 if (IT_MODE)
179 {
180 // Transform initial guess
181 DGMassBasis<DIM,D1D>(e, NE, d2q_B, u, u, d1d);
182 }
183 }
184
185 const int tid = MFEM_THREAD_ID(x) + NB*MFEM_THREAD_ID(y);
186
187 // Compute first residual
188 if (IT_MODE)
189 {
190 DGMassApply<DIM,D1D,Q1D>(e, NE, B, Bt, pa_data, u, r, d1d, q1d);
191 DGMassAxpy(e, NE, ND, 1.0, b, -1.0, r, r); // r = b - r
192 }
193 else
194 {
195 // if not in iterative mode, use zero initial guess
196 const int BX = MFEM_THREAD_SIZE(x);
197 const int BY = MFEM_THREAD_SIZE(y);
198 const int bxy = BX*BY;
199 const auto B = ConstDeviceMatrix(b, ND, NE);
200 auto U = DeviceMatrix(u, ND, NE);
201 auto R = DeviceMatrix(r, ND, NE);
202 for (int i = tid; i < ND; i += bxy)
203 {
204 U(i, e) = 0.0;
205 R(i, e) = B(i, e);
206 }
207 MFEM_SYNC_THREAD;
208 }
209
210 DGMassPreconditioner(e, NE, ND, dinv, r, z);
211 DGMassAxpy(e, NE, ND, 1.0, z, 0.0, z, d); // d = z
212
213 real_t nom = DGMassDot<NB>(e, NE, ND, d, r);
214 if (nom < 0.0) { return; /* Not positive definite */ }
215 real_t r0 = fmax(nom*RELTOL*RELTOL, ABSTOL*ABSTOL);
216 if (nom <= r0) { return; /* Converged */ }
217
218 DGMassApply<DIM,D1D,Q1D>(e, NE, B, Bt, pa_data, d, z, d1d, q1d);
219 real_t den = DGMassDot<NB>(e, NE, ND, z, d);
220 if (den <= 0.0)
221 {
222 DGMassDot<NB>(e, NE, ND, d, d);
223 // d2 > 0 => not positive definite
224 if (den == 0.0) { return; }
225 }
226
227 // start iteration
228 int i = 1;
229 while (true)
230 {
231 const real_t alpha = nom/den;
232 DGMassAxpy(e, NE, ND, 1.0, u, alpha, d, u); // u = u + alpha*d
233 DGMassAxpy(e, NE, ND, 1.0, r, -alpha, z, r); // r = r - alpha*A*d
234
235 DGMassPreconditioner(e, NE, ND, dinv, r, z);
236
237 real_t betanom = DGMassDot<NB>(e, NE, ND, r, z);
238 if (betanom < 0.0) { return; /* Not positive definite */ }
239 if (betanom <= r0) { break; /* Converged */ }
240
241 if (++i > MAXIT) { break; }
242
243 const real_t beta = betanom/nom;
244 DGMassAxpy(e, NE, ND, 1.0, z, beta, d, d); // d = z + beta*d
245 DGMassApply<DIM,D1D,Q1D>(e, NE, B, Bt, pa_data, d, z, d1d, q1d); // z = A d
246 den = DGMassDot<NB>(e, NE, ND, d, z);
247 if (den <= 0.0)
248 {
249 DGMassDot<NB>(e, NE, ND, d, d);
250 // d2 > 0 => not positive definite
251 if (den == 0.0) { break; }
252 }
253 nom = betanom;
254 }
255
256 if (CHANGE_BASIS)
257 {
258 DGMassBasis<DIM,D1D>(e, NE, q2d_B, u, u, d1d);
259 }
260 });
261}
262
263void DGMassInverse::Mult(const Vector &Mu, Vector &u) const
264{
265 // Dispatch to templated version based on dim, d1d, and q1d.
266 const int dim = fes.GetMesh()->Dimension();
267 const int d1d = m->dofs1D;
268 const int q1d = m->quad1D;
269
270 const int id = (d1d << 4) | q1d;
271
272 if (dim == 2)
273 {
274 switch (id)
275 {
276 case 0x11: return DGMassCGIteration<2,1,1>(Mu, u);
277 case 0x22: return DGMassCGIteration<2,2,2>(Mu, u);
278 case 0x33: return DGMassCGIteration<2,3,3>(Mu, u);
279 case 0x35: return DGMassCGIteration<2,3,5>(Mu, u);
280 case 0x44: return DGMassCGIteration<2,4,4>(Mu, u);
281 case 0x46: return DGMassCGIteration<2,4,6>(Mu, u);
282 case 0x55: return DGMassCGIteration<2,5,5>(Mu, u);
283 case 0x57: return DGMassCGIteration<2,5,7>(Mu, u);
284 case 0x66: return DGMassCGIteration<2,6,6>(Mu, u);
285 case 0x68: return DGMassCGIteration<2,6,8>(Mu, u);
286 default: return DGMassCGIteration<2>(Mu, u); // Fallback
287 }
288 }
289 else if (dim == 3)
290 {
291 switch (id)
292 {
293 case 0x22: return DGMassCGIteration<3,2,2>(Mu, u);
294 case 0x23: return DGMassCGIteration<3,2,3>(Mu, u);
295 case 0x33: return DGMassCGIteration<3,3,3>(Mu, u);
296 case 0x34: return DGMassCGIteration<3,3,4>(Mu, u);
297 case 0x35: return DGMassCGIteration<3,3,5>(Mu, u);
298 case 0x44: return DGMassCGIteration<3,4,4>(Mu, u);
299 case 0x45: return DGMassCGIteration<3,4,5>(Mu, u);
300 case 0x46: return DGMassCGIteration<3,4,6>(Mu, u);
301 case 0x48: return DGMassCGIteration<3,4,8>(Mu, u);
302 case 0x55: return DGMassCGIteration<3,5,5>(Mu, u);
303 case 0x56: return DGMassCGIteration<3,5,6>(Mu, u);
304 case 0x57: return DGMassCGIteration<3,5,7>(Mu, u);
305 case 0x58: return DGMassCGIteration<3,5,8>(Mu, u);
306 case 0x66: return DGMassCGIteration<3,6,6>(Mu, u);
307 case 0x67: return DGMassCGIteration<3,6,7>(Mu, u);
308 default: return DGMassCGIteration<3>(Mu, u); // Fallback
309 }
310 }
311}
312
313} // namespace mfem
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition array.hpp:317
T * HostReadWrite()
Shortcut for mfem::ReadWrite(a.GetMemory(), a.Size(), false).
Definition array.hpp:337
T * HostWrite()
Shortcut for mfem::Write(a.GetMemory(), a.Size(), false).
Definition array.hpp:329
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
virtual void AssembleDiagonal(Vector &diag) const
Assemble the diagonal of the bilinear form into diag. Note that diag is a tdof Vector.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Solver for the discontinuous Galerkin mass matrix.
Definition dgmassinv.hpp:28
Vector diag_inv
Jacobi preconditioner.
Definition dgmassinv.hpp:37
class BilinearForm * M
Mass bilinear form, owned.
Definition dgmassinv.hpp:35
class MassIntegrator * m
Mass integrator, owned by the form M.
Definition dgmassinv.hpp:36
real_t rel_tol
Relative CG tolerance.
Definition dgmassinv.hpp:38
int max_iter
Maximum number of CG iterations;.
Definition dgmassinv.hpp:40
void SetOperator(const Operator &op)
Not implemented. Aborts.
Definition dgmassinv.cpp:95
Array< real_t > B_
Inverse of change of basis.
Definition dgmassinv.hpp:33
real_t abs_tol
Absolute CG tolerance.
Definition dgmassinv.hpp:39
void Mult(const Vector &b, Vector &u) const
Solve the system M b = u.
DGMassInverse(FiniteElementSpace &fes_, Coefficient *coeff, const IntegrationRule *ir, int btype)
Protected constructor, used internally.
Definition dgmassinv.cpp:20
Array< real_t > Bt_
Inverse of change of basis, transposed.
Definition dgmassinv.hpp:34
void Update()
Recompute operator and preconditioner (when coefficient or mesh changes).
void DGMassCGIteration(const Vector &b_, Vector &u_) const
Solve the system M b = u. Not part of the public interface.
void SetRelTol(const real_t rel_tol_)
Set the relative tolerance.
void SetMaxIter(const real_t max_iter_)
Set the maximum number of iterations.
FiniteElementSpace fes
FE space in requested basis.
Definition dgmassinv.hpp:31
const DofToQuad * d2q
Change of basis. Not owned.
Definition dgmassinv.hpp:32
void SetAbsTol(const real_t abs_tol_)
Set the absolute tolerance.
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void Transpose()
(*this) = (*this)^t
@ TENSOR
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
Definition fe_base.hpp:161
Array< real_t > B
Basis functions evaluated at quadrature points.
Definition fe_base.hpp:189
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number.
Definition fe_base.hpp:174
Array< real_t > Bt
Transpose of B.
Definition fe_base.hpp:195
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition fespace.hpp:581
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition fespace.cpp:3168
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:740
const FiniteElementCollection * FEColl() const
Definition fespace.hpp:727
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:559
bool IsDGSpace() const
Return whether or not the space is discontinuous (L2)
Definition fespace.hpp:1313
virtual const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Return a DofToQuad structure corresponding to the given IntegrationRule using the given DofToQuad::Mo...
Definition fe_base.cpp:365
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition fe_base.hpp:395
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
virtual void GetInverseMatrix(int m, real_t *X) const
Assuming L.U = P.A factored data of size (m x m), compute X <- A^{-1}.
virtual bool Factor(int m, real_t TOL=0.0)
Compute the LU factorization of the current matrix.
const DofToQuad * maps
Not owned.
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
Abstract operator.
Definition operator.hpp:25
int height
Dimension of the output / number of rows in the matrix.
Definition operator.hpp:27
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
Vector data type.
Definition vector.hpp:80
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:474
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:490
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
void Reciprocal()
(*this)(i) = 1.0 / (*this)(i)
Definition vector.cpp:308
virtual real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:482
Vector beta
const real_t alpha
Definition ex15.cpp:369
int dim
Definition ex24.cpp:53
Mesh * GetMesh(int type)
Definition ex29.cpp:218
real_t b
Definition lissajous.cpp:42
constexpr int DIM
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
DeviceTensor< 2, real_t > DeviceMatrix
Definition dtensor.hpp:143
DeviceTensor< 2, const real_t > ConstDeviceMatrix
Definition dtensor.hpp:144
void forall_2D(int N, int X, int Y, lambda &&body)
Definition forall.hpp:763
float real_t
Definition config.hpp:43