12 #include "../general/forall.hpp"
29 if (mesh->
GetNE() == 0) {
return; }
53 GeometricFactors::JACOBIANS);
57 pa_data.SetSize(ne*nq, Device::GetDeviceMemoryType());
62 MFEM_VERIFY(cQ != NULL,
"Only ConstantCoefficient is supported.");
65 if (!(
dim == 2 ||
dim == 3))
67 MFEM_ABORT(
"Dimension not supported.");
71 const double constant = coeff;
75 auto J =
Reshape(geom->J.Read(), NQ,2,2,NE);
76 auto v =
Reshape(pa_data.Write(), NQ, NE);
79 for (
int q = 0; q < NQ; ++q)
81 const double J11 = J(q,0,0,e);
82 const double J12 = J(q,1,0,e);
83 const double J21 = J(q,0,1,e);
84 const double J22 = J(q,1,1,e);
85 const double detJ = (J11*J22)-(J21*J12);
86 v(q,e) = w[q] * constant * detJ;
92 const double constant = coeff;
96 auto J =
Reshape(geom->J.Read(), NQ,3,3,NE);
97 auto v =
Reshape(pa_data.Write(), NQ,NE);
100 for (
int q = 0; q < NQ; ++q)
102 const double J11 = J(q,0,0,e), J12 = J(q,0,1,e), J13 = J(q,0,2,e);
103 const double J21 = J(q,1,0,e), J22 = J(q,1,1,e), J23 = J(q,1,2,e);
104 const double J31 = J(q,2,0,e), J32 = J(q,2,1,e), J33 = J(q,2,2,e);
105 const double detJ = J11 * (J22 * J33 - J32 * J23) -
106 J21 * (J12 * J33 - J32 * J13) +
107 J31 * (J12 * J23 - J22 * J13);
108 v(q,e) = W[q] * constant * detJ;
114 template<
const int T_D1D = 0,
116 static void PAVectorMassApply2D(
const int NE,
125 const int D1D = T_D1D ? T_D1D : d1d;
126 const int Q1D = T_Q1D ? T_Q1D : q1d;
127 constexpr
int VDIM = 2;
128 MFEM_VERIFY(D1D <=
MAX_D1D,
"");
129 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
137 const int D1D = T_D1D ? T_D1D : d1d;
138 const int Q1D = T_Q1D ? T_Q1D : q1d;
140 constexpr
int max_D1D = T_D1D ? T_D1D :
MAX_D1D;
141 constexpr
int max_Q1D = T_Q1D ? T_Q1D :
MAX_Q1D;
142 double sol_xy[max_Q1D][max_Q1D];
143 for (
int c = 0; c < VDIM; ++c)
145 for (
int qy = 0; qy < Q1D; ++qy)
147 for (
int qx = 0; qx < Q1D; ++qx)
149 sol_xy[qy][qx] = 0.0;
152 for (
int dy = 0; dy < D1D; ++dy)
154 double sol_x[max_Q1D];
155 for (
int qy = 0; qy < Q1D; ++qy)
159 for (
int dx = 0; dx < D1D; ++dx)
161 const double s = x(dx,dy,c,e);
162 for (
int qx = 0; qx < Q1D; ++qx)
164 sol_x[qx] += B(qx,dx)*
s;
167 for (
int qy = 0; qy < Q1D; ++qy)
169 const double d2q = B(qy,dy);
170 for (
int qx = 0; qx < Q1D; ++qx)
172 sol_xy[qy][qx] += d2q * sol_x[qx];
176 for (
int qy = 0; qy < Q1D; ++qy)
178 for (
int qx = 0; qx < Q1D; ++qx)
180 sol_xy[qy][qx] *= op(qx,qy,e);
183 for (
int qy = 0; qy < Q1D; ++qy)
185 double sol_x[max_D1D];
186 for (
int dx = 0; dx < D1D; ++dx)
190 for (
int qx = 0; qx < Q1D; ++qx)
192 const double s = sol_xy[qy][qx];
193 for (
int dx = 0; dx < D1D; ++dx)
195 sol_x[dx] += Bt(dx,qx) *
s;
198 for (
int dy = 0; dy < D1D; ++dy)
200 const double q2d = Bt(dy,qy);
201 for (
int dx = 0; dx < D1D; ++dx)
203 y(dx,dy,c,e) += q2d * sol_x[dx];
211 template<
const int T_D1D = 0,
213 static void PAVectorMassApply3D(
const int NE,
214 const Array<double> &B_,
215 const Array<double> &Bt_,
222 const int D1D = T_D1D ? T_D1D : d1d;
223 const int Q1D = T_Q1D ? T_Q1D : q1d;
224 constexpr
int VDIM = 3;
225 MFEM_VERIFY(D1D <=
MAX_D1D,
"");
226 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
227 auto B =
Reshape(B_.Read(), Q1D, D1D);
228 auto Bt =
Reshape(Bt_.Read(), D1D, Q1D);
229 auto op =
Reshape(op_.Read(), Q1D, Q1D, Q1D, NE);
230 auto x =
Reshape(x_.Read(), D1D, D1D, D1D, VDIM, NE);
231 auto y =
Reshape(y_.ReadWrite(), D1D, D1D, D1D, VDIM, NE);
234 const int D1D = T_D1D ? T_D1D : d1d;
235 const int Q1D = T_Q1D ? T_Q1D : q1d;
236 constexpr
int max_D1D = T_D1D ? T_D1D :
MAX_D1D;
237 constexpr
int max_Q1D = T_Q1D ? T_Q1D :
MAX_Q1D;
238 double sol_xyz[max_Q1D][max_Q1D][max_Q1D];
239 for (
int c = 0; c < VDIM; ++ c)
241 for (
int qz = 0; qz < Q1D; ++qz)
243 for (
int qy = 0; qy < Q1D; ++qy)
245 for (
int qx = 0; qx < Q1D; ++qx)
247 sol_xyz[qz][qy][qx] = 0.0;
251 for (
int dz = 0; dz < D1D; ++dz)
253 double sol_xy[max_Q1D][max_Q1D];
254 for (
int qy = 0; qy < Q1D; ++qy)
256 for (
int qx = 0; qx < Q1D; ++qx)
258 sol_xy[qy][qx] = 0.0;
261 for (
int dy = 0; dy < D1D; ++dy)
263 double sol_x[max_Q1D];
264 for (
int qx = 0; qx < Q1D; ++qx)
268 for (
int dx = 0; dx < D1D; ++dx)
270 const double s = x(dx,dy,dz,c,e);
271 for (
int qx = 0; qx < Q1D; ++qx)
273 sol_x[qx] += B(qx,dx) *
s;
276 for (
int qy = 0; qy < Q1D; ++qy)
278 const double wy = B(qy,dy);
279 for (
int qx = 0; qx < Q1D; ++qx)
281 sol_xy[qy][qx] += wy * sol_x[qx];
285 for (
int qz = 0; qz < Q1D; ++qz)
287 const double wz = B(qz,dz);
288 for (
int qy = 0; qy < Q1D; ++qy)
290 for (
int qx = 0; qx < Q1D; ++qx)
292 sol_xyz[qz][qy][qx] += wz * sol_xy[qy][qx];
297 for (
int qz = 0; qz < Q1D; ++qz)
299 for (
int qy = 0; qy < Q1D; ++qy)
301 for (
int qx = 0; qx < Q1D; ++qx)
303 sol_xyz[qz][qy][qx] *= op(qx,qy,qz,e);
307 for (
int qz = 0; qz < Q1D; ++qz)
309 double sol_xy[max_D1D][max_D1D];
310 for (
int dy = 0; dy < D1D; ++dy)
312 for (
int dx = 0; dx < D1D; ++dx)
317 for (
int qy = 0; qy < Q1D; ++qy)
319 double sol_x[max_D1D];
320 for (
int dx = 0; dx < D1D; ++dx)
324 for (
int qx = 0; qx < Q1D; ++qx)
326 const double s = sol_xyz[qz][qy][qx];
327 for (
int dx = 0; dx < D1D; ++dx)
329 sol_x[dx] += Bt(dx,qx) *
s;
332 for (
int dy = 0; dy < D1D; ++dy)
334 const double wy = Bt(dy,qy);
335 for (
int dx = 0; dx < D1D; ++dx)
337 sol_xy[dy][dx] += wy * sol_x[dx];
341 for (
int dz = 0; dz < D1D; ++dz)
343 const double wz = Bt(dz,qz);
344 for (
int dy = 0; dy < D1D; ++dy)
346 for (
int dx = 0; dx < D1D; ++dx)
348 y(dx,dy,dz,c,e) += wz * sol_xy[dy][dx];
357 static void PAVectorMassApply(
const int dim,
361 const Array<double> &B,
362 const Array<double> &Bt,
369 return PAVectorMassApply2D(NE, B, Bt, op, x, y, D1D, Q1D);
373 return PAVectorMassApply3D(NE, B, Bt, op, x, y, D1D, Q1D);
375 MFEM_ABORT(
"Unknown kernel.");
382 ceedOp->AddMult(x, y);
386 PAVectorMassApply(dim, dofs1D, quad1D, ne, maps->B, maps->Bt, pa_data, x, y);
390 template<const
int T_D1D = 0, const
int T_Q1D = 0>
391 static void PAVectorMassAssembleDiagonal2D(
const int NE,
399 const int D1D = T_D1D ? T_D1D : d1d;
400 const int Q1D = T_Q1D ? T_Q1D : q1d;
401 constexpr
int VDIM = 2;
402 MFEM_VERIFY(D1D <=
MAX_D1D,
"");
403 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
409 const int D1D = T_D1D ? T_D1D : d1d;
410 const int Q1D = T_Q1D ? T_Q1D : q1d;
411 constexpr
int max_D1D = T_D1D ? T_D1D :
MAX_D1D;
412 constexpr
int max_Q1D = T_Q1D ? T_Q1D :
MAX_Q1D;
414 double temp[max_Q1D][max_D1D];
415 for (
int qx = 0; qx < Q1D; ++qx)
417 for (
int dy = 0; dy < D1D; ++dy)
420 for (
int qy = 0; qy < Q1D; ++qy)
422 temp[qx][dy] += B(qy, dy) * B(qy, dy) * op(qx, qy, e);
426 for (
int dy = 0; dy < D1D; ++dy)
428 for (
int dx = 0; dx < D1D; ++dx)
431 for (
int qx = 0; qx < Q1D; ++qx)
433 temp1 += B(qx, dx) * B(qx, dx) * temp[qx][dy];
435 y(dx, dy, 0, e) = temp1;
436 y(dx, dy, 1, e) = temp1;
442 template<const
int T_D1D = 0, const
int T_Q1D = 0>
443 static void PAVectorMassAssembleDiagonal3D(
const int NE,
444 const Array<double> &B_,
445 const Array<double> &Bt_,
451 const int D1D = T_D1D ? T_D1D : d1d;
452 const int Q1D = T_Q1D ? T_Q1D : q1d;
453 constexpr
int VDIM = 3;
454 MFEM_VERIFY(D1D <=
MAX_D1D,
"");
455 MFEM_VERIFY(Q1D <=
MAX_Q1D,
"");
456 auto B =
Reshape(B_.Read(), Q1D, D1D);
457 auto op =
Reshape(op_.Read(), Q1D, Q1D, Q1D, NE);
458 auto y =
Reshape(diag_.ReadWrite(), D1D, D1D, D1D, VDIM, NE);
461 const int D1D = T_D1D ? T_D1D : d1d;
462 const int Q1D = T_Q1D ? T_Q1D : q1d;
464 constexpr
int max_D1D = T_D1D ? T_D1D :
MAX_D1D;
465 constexpr
int max_Q1D = T_Q1D ? T_Q1D :
MAX_Q1D;
467 double temp[max_Q1D][max_Q1D][max_D1D];
468 for (
int qx = 0; qx < Q1D; ++qx)
470 for (
int qy = 0; qy < Q1D; ++qy)
472 for (
int dz = 0; dz < D1D; ++dz)
474 temp[qx][qy][dz] = 0.0;
475 for (
int qz = 0; qz < Q1D; ++qz)
477 temp[qx][qy][dz] += B(qz, dz) * B(qz, dz) * op(qx, qy, qz, e);
482 double temp2[max_Q1D][max_D1D][max_D1D];
483 for (
int qx = 0; qx < Q1D; ++qx)
485 for (
int dz = 0; dz < D1D; ++dz)
487 for (
int dy = 0; dy < D1D; ++dy)
489 temp2[qx][dy][dz] = 0.0;
490 for (
int qy = 0; qy < Q1D; ++qy)
492 temp2[qx][dy][dz] += B(qy, dy) * B(qy, dy) * temp[qx][qy][dz];
497 for (
int dz = 0; dz < D1D; ++dz)
499 for (
int dy = 0; dy < D1D; ++dy)
501 for (
int dx = 0; dx < D1D; ++dx)
504 for (
int qx = 0; qx < Q1D; ++qx)
506 temp3 += B(qx, dx) * B(qx, dx)
509 y(dx, dy, dz, 0, e) = temp3;
510 y(dx, dy, dz, 1, e) = temp3;
511 y(dx, dy, dz, 2, e) = temp3;
518 static void PAVectorMassAssembleDiagonal(
const int dim,
522 const Array<double> &B,
523 const Array<double> &Bt,
529 return PAVectorMassAssembleDiagonal2D(NE, B, Bt, op, y, D1D, Q1D);
533 return PAVectorMassAssembleDiagonal3D(NE, B, Bt, op, y, D1D, Q1D);
535 MFEM_ABORT(
"Dimension not implemented.");
538 void VectorMassIntegrator::AssembleDiagonalPA(
Vector &diag)
542 ceedOp->GetDiagonal(diag);
546 PAVectorMassAssembleDiagonal(dim,
int GetNPoints() const
Returns the number of the points in the integration rule.
Abstract class for all finite elements.
Class for an integration rule - an Array of IntegrationPoint.
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
A coefficient that is constant across space and time.
const GeometricFactors * GetGeometricFactors(const IntegrationRule &ir, const int flags, MemoryType d_mt=MemoryType::DEFAULT)
Return the mesh geometric factors corresponding to the given integration rule.
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
int GetNE() const
Returns number of elements.
const Array< double > & GetWeights() const
Return the quadrature weights in a contiguous array.
Mesh * GetMesh() const
Returns the mesh.
Represent a MassIntegrator with AssemblyLevel::Partial using libCEED.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
bool DeviceCanUseCeed()
Function that determines if a CEED kernel should be used, based on the current mfem::Device configura...
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...
const IntegrationRule & GetRule(const Integrator &integ, const FiniteElement &trial_fe, const FiniteElement &test_fe, ElementTransformation &Trans)
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.