15 #include "../general/forall.hpp"
23 :
Solver(fes_orig.GetTrueVSize()),
24 fec(fes_orig.GetMaxElementOrder(),
25 fes_orig.
GetMesh()->Dimension(),
27 fes_orig.GetFE(0)->GetMapType()),
33 const int btype_orig =
36 if (btype_orig == btype)
97 MFEM_ABORT(
"SetOperator not supported with DGMassInverse.")
118 template<
int DIM,
int D1D,
int Q1D>
121 using namespace internal;
127 const int ND =
static_cast<int>(pow(d1d,
DIM));
142 const bool CHANGE_BASIS = (
d2q !=
nullptr);
149 double *b2 =
nullptr;
150 const double *b_orig =
nullptr;
151 const double *d2q_B =
nullptr;
152 const double *q2d_B =
nullptr;
153 const double *q2d_Bt =
nullptr;
169 constexpr
int NB = Q1D ? Q1D : 1;
171 MFEM_FORALL_2D(e, NE, NB, NB, 1,
173 constexpr
int NB = Q1D ? Q1D : 1;
179 DGMassBasis<DIM,D1D,MAX_D1D>(e, NE, q2d_Bt, b_orig, b2, d1d);
183 DGMassBasis<DIM,D1D,MAX_D1D>(e, NE, d2q_B,
u,
u, d1d);
187 const int tid = MFEM_THREAD_ID(x) + NB*MFEM_THREAD_ID(y);
192 DGMassApply<DIM,D1D,Q1D>(e, NE, B, Bt, pa_data,
u, r, d1d, q1d);
193 DGMassAxpy(e, NE, ND, 1.0,
b, -1.0, r, r);
198 const int BX = MFEM_THREAD_SIZE(x);
199 const int BY = MFEM_THREAD_SIZE(y);
200 const int bxy = BX*BY;
204 for (
int i = tid; i < ND; i += bxy)
212 DGMassPreconditioner(e, NE, ND, dinv, r, z);
213 DGMassAxpy(e, NE, ND, 1.0, z, 0.0, z, d);
215 double nom = DGMassDot<NB>(e, NE, ND, d, r);
216 if (nom < 0.0) {
return; }
217 double r0 = fmax(nom*RELTOL*RELTOL, ABSTOL*ABSTOL);
218 if (nom <= r0) {
return; }
220 DGMassApply<DIM,D1D,Q1D>(e, NE, B, Bt, pa_data, d, z, d1d, q1d);
221 double den = DGMassDot<NB>(e, NE, ND, z, d);
224 DGMassDot<NB>(e, NE, ND, d, d);
226 if (den == 0.0) {
return; }
233 const double alpha = nom/den;
234 DGMassAxpy(e, NE, ND, 1.0,
u, alpha, d,
u);
235 DGMassAxpy(e, NE, ND, 1.0, r, -alpha, z, r);
237 DGMassPreconditioner(e, NE, ND, dinv, r, z);
239 double betanom = DGMassDot<NB>(e, NE, ND, r, z);
240 if (betanom < 0.0) {
return; }
241 if (betanom <= r0) {
break; }
243 if (++i > MAXIT) {
break; }
245 const double beta = betanom/nom;
246 DGMassAxpy(e, NE, ND, 1.0, z, beta, d, d);
247 DGMassApply<DIM,D1D,Q1D>(e, NE, B, Bt, pa_data, d, z, d1d, q1d);
248 den = DGMassDot<NB>(e, NE, ND, d, z);
251 DGMassDot<NB>(e, NE, ND, d, d);
253 if (den == 0.0) {
break; }
260 DGMassBasis<DIM,D1D,MAX_D1D>(e, NE, q2d_B,
u,
u, d1d);
272 const int id = (d1d << 4) | q1d;
278 case 0x11:
return DGMassCGIteration<2,1,1>(Mu,
u);
279 case 0x22:
return DGMassCGIteration<2,2,2>(Mu,
u);
280 case 0x33:
return DGMassCGIteration<2,3,3>(Mu,
u);
281 case 0x35:
return DGMassCGIteration<2,3,5>(Mu,
u);
282 case 0x44:
return DGMassCGIteration<2,4,4>(Mu,
u);
283 case 0x46:
return DGMassCGIteration<2,4,6>(Mu,
u);
284 case 0x55:
return DGMassCGIteration<2,5,5>(Mu,
u);
285 case 0x57:
return DGMassCGIteration<2,5,7>(Mu,
u);
286 case 0x66:
return DGMassCGIteration<2,6,6>(Mu,
u);
287 case 0x68:
return DGMassCGIteration<2,6,8>(Mu,
u);
288 default:
return DGMassCGIteration<2>(Mu,
u);
295 case 0x22:
return DGMassCGIteration<3,2,2>(Mu,
u);
296 case 0x23:
return DGMassCGIteration<3,2,3>(Mu,
u);
297 case 0x33:
return DGMassCGIteration<3,3,3>(Mu,
u);
298 case 0x34:
return DGMassCGIteration<3,3,4>(Mu,
u);
299 case 0x35:
return DGMassCGIteration<3,3,5>(Mu,
u);
300 case 0x44:
return DGMassCGIteration<3,4,4>(Mu,
u);
301 case 0x45:
return DGMassCGIteration<3,4,5>(Mu,
u);
302 case 0x46:
return DGMassCGIteration<3,4,6>(Mu,
u);
303 case 0x48:
return DGMassCGIteration<3,4,8>(Mu,
u);
304 case 0x55:
return DGMassCGIteration<3,5,5>(Mu,
u);
305 case 0x56:
return DGMassCGIteration<3,5,6>(Mu,
u);
306 case 0x57:
return DGMassCGIteration<3,5,7>(Mu,
u);
307 case 0x58:
return DGMassCGIteration<3,5,8>(Mu,
u);
308 case 0x66:
return DGMassCGIteration<3,6,6>(Mu,
u);
309 case 0x67:
return DGMassCGIteration<3,6,7>(Mu,
u);
310 default:
return DGMassCGIteration<3>(Mu,
u);
T * HostReadWrite()
Shortcut for mfem::ReadWrite(a.GetMemory(), a.Size(), false).
Class for an integration rule - an Array of IntegrationPoint.
DeviceTensor< 2, const double > ConstDeviceMatrix
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Tensor product representation using 1D matrices/tensors with dimensions using 1D number of quadrature...
void SetSize(int s)
Resize the vector to size s.
void SetOperator(const Operator &op)
Not implemented. Aborts.
Data type dense matrix using column-major storage.
int Size() const
Returns the size of the vector.
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
DGMassInverse(FiniteElementSpace &fes_, Coefficient *coeff, const IntegrationRule *ir, int btype)
Protected constructor, used internally.
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
void SetAbsTol(const double abs_tol_)
Set the absolute tolerance.
const DofToQuad * maps
Not owned.
DeviceTensor< 2, double > DeviceMatrix
Array< double > Bt_
Inverse of change of basis, transposed.
int GetNE() const
Returns number of elements in the mesh.
void SetMaxIter(const double max_iter_)
Set the maximum number of iterations.
double abs_tol
Absolute CG tolerance.
Mesh * GetMesh() const
Returns the mesh.
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
double rel_tol
Relative CG tolerance.
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
bool IsDGSpace() const
Return whether or not the space is discontinuous (L2)
Array< double > B_
Inverse of change of basis.
Array< double > Bt
Transpose of B.
void Transpose()
(*this) = (*this)^t
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void Mult(const Vector &b, Vector &u) const
Solve the system M b = u.
class BilinearForm * M
Mass bilinear form, owned.
int max_iter
Maximum number of CG iterations;.
int height
Dimension of the output / number of rows in the matrix.
Array< double > B
Basis functions evaluated at quadrature points.
void DGMassCGIteration(const Vector &b_, Vector &u_) const
Solve the system M b = u. Not part of the public interface.
T * HostWrite()
Shortcut for mfem::Write(a.GetMemory(), a.Size(), false).
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...
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
class MassIntegrator * m
Mass integrator, owned by the form M.
virtual bool Factor(int m, double TOL=0.0)
Compute the LU factorization of the current matrix.
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
const DofToQuad * d2q
Change of basis. Not owned.
Vector diag_inv
Jacobi preconditioner.
const FiniteElementCollection * FEColl() const
double u(const Vector &xvec)
FiniteElementSpace fes
FE space in requested basis.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
void SetRelTol(const double rel_tol_)
Set the relative tolerance.
Solver for the discontinuous Galerkin mass matrix.
void Update()
Recompute operator and preconditioner (when coefficient or mesh changes).
Arbitrary order "L2-conforming" discontinuous finite elements.