MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
bilininteg_vecmass.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, 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 "../general/forall.hpp"
13 #include "bilininteg.hpp"
14 #include "gridfunc.hpp"
16 
17 using namespace std;
18 
19 namespace mfem
20 {
21 
22 // PA Mass Integrator
23 
24 // PA Mass Assemble kernel
25 void VectorMassIntegrator::AssemblePA(const FiniteElementSpace &fes)
26 {
27  // Assuming the same element type
28  Mesh *mesh = fes.GetMesh();
29  if (mesh->GetNE() == 0) { return; }
30  const FiniteElement &el = *fes.GetFE(0);
32  const IntegrationRule *ir
33  = IntRule ? IntRule : &MassIntegrator::GetRule(el, el, *T);
34  if (DeviceCanUseCeed())
35  {
36  delete ceedOp;
37  const bool mixed = mesh->GetNumGeometries(mesh->Dimension()) > 1 ||
38  fes.IsVariableOrder();
39  if (mixed)
40  {
41  ceedOp = new ceed::MixedPAMassIntegrator(*this, fes, Q);
42  }
43  else
44  {
45  ceedOp = new ceed::PAMassIntegrator(fes, *ir, Q);
46  }
47  return;
48  }
49  dim = mesh->Dimension();
50  ne = fes.GetMesh()->GetNE();
51  nq = ir->GetNPoints();
52  geom = mesh->GetGeometricFactors(*ir, GeometricFactors::COORDINATES |
53  GeometricFactors::JACOBIANS);
54  maps = &el.GetDofToQuad(*ir, DofToQuad::TENSOR);
55  dofs1D = maps->ndof;
56  quad1D = maps->nqpt;
57  pa_data.SetSize(ne*nq, Device::GetDeviceMemoryType());
58  double coeff = 1.0;
59  if (Q)
60  {
61  ConstantCoefficient *cQ = dynamic_cast<ConstantCoefficient*>(Q);
62  MFEM_VERIFY(cQ != NULL, "Only ConstantCoefficient is supported.");
63  coeff = cQ->constant;
64  }
65  if (!(dim == 2 || dim == 3))
66  {
67  MFEM_ABORT("Dimension not supported.");
68  }
69  if (dim == 2)
70  {
71  const double constant = coeff;
72  const int NE = ne;
73  const int NQ = nq;
74  auto w = ir->GetWeights().Read();
75  auto J = Reshape(geom->J.Read(), NQ,2,2,NE);
76  auto v = Reshape(pa_data.Write(), NQ, NE);
77  MFEM_FORALL(e, NE,
78  {
79  for (int q = 0; q < NQ; ++q)
80  {
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;
87  }
88  });
89  }
90  if (dim == 3)
91  {
92  const double constant = coeff;
93  const int NE = ne;
94  const int NQ = nq;
95  auto W = ir->GetWeights().Read();
96  auto J = Reshape(geom->J.Read(), NQ,3,3,NE);
97  auto v = Reshape(pa_data.Write(), NQ,NE);
98  MFEM_FORALL(e, NE,
99  {
100  for (int q = 0; q < NQ; ++q)
101  {
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;
109  }
110  });
111  }
112 }
113 
114 template<const int T_D1D = 0,
115  const int T_Q1D = 0>
116 static void PAVectorMassApply2D(const int NE,
117  const Array<double> &B_,
118  const Array<double> &Bt_,
119  const Vector &op_,
120  const Vector &x_,
121  Vector &y_,
122  const int d1d = 0,
123  const int q1d = 0)
124 {
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, "");
130  auto B = Reshape(B_.Read(), Q1D, D1D);
131  auto Bt = Reshape(Bt_.Read(), D1D, Q1D);
132  auto op = Reshape(op_.Read(), Q1D, Q1D, NE);
133  auto x = Reshape(x_.Read(), D1D, D1D, VDIM, NE);
134  auto y = Reshape(y_.ReadWrite(), D1D, D1D, VDIM, NE);
135  MFEM_FORALL(e, NE,
136  {
137  const int D1D = T_D1D ? T_D1D : d1d; // nvcc workaround
138  const int Q1D = T_Q1D ? T_Q1D : q1d;
139  // the following variables are evaluated at compile time
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)
144  {
145  for (int qy = 0; qy < Q1D; ++qy)
146  {
147  for (int qx = 0; qx < Q1D; ++qx)
148  {
149  sol_xy[qy][qx] = 0.0;
150  }
151  }
152  for (int dy = 0; dy < D1D; ++dy)
153  {
154  double sol_x[max_Q1D];
155  for (int qy = 0; qy < Q1D; ++qy)
156  {
157  sol_x[qy] = 0.0;
158  }
159  for (int dx = 0; dx < D1D; ++dx)
160  {
161  const double s = x(dx,dy,c,e);
162  for (int qx = 0; qx < Q1D; ++qx)
163  {
164  sol_x[qx] += B(qx,dx)* s;
165  }
166  }
167  for (int qy = 0; qy < Q1D; ++qy)
168  {
169  const double d2q = B(qy,dy);
170  for (int qx = 0; qx < Q1D; ++qx)
171  {
172  sol_xy[qy][qx] += d2q * sol_x[qx];
173  }
174  }
175  }
176  for (int qy = 0; qy < Q1D; ++qy)
177  {
178  for (int qx = 0; qx < Q1D; ++qx)
179  {
180  sol_xy[qy][qx] *= op(qx,qy,e);
181  }
182  }
183  for (int qy = 0; qy < Q1D; ++qy)
184  {
185  double sol_x[max_D1D];
186  for (int dx = 0; dx < D1D; ++dx)
187  {
188  sol_x[dx] = 0.0;
189  }
190  for (int qx = 0; qx < Q1D; ++qx)
191  {
192  const double s = sol_xy[qy][qx];
193  for (int dx = 0; dx < D1D; ++dx)
194  {
195  sol_x[dx] += Bt(dx,qx) * s;
196  }
197  }
198  for (int dy = 0; dy < D1D; ++dy)
199  {
200  const double q2d = Bt(dy,qy);
201  for (int dx = 0; dx < D1D; ++dx)
202  {
203  y(dx,dy,c,e) += q2d * sol_x[dx];
204  }
205  }
206  }
207  }
208  });
209 }
210 
211 template<const int T_D1D = 0,
212  const int T_Q1D = 0>
213 static void PAVectorMassApply3D(const int NE,
214  const Array<double> &B_,
215  const Array<double> &Bt_,
216  const Vector &op_,
217  const Vector &x_,
218  Vector &y_,
219  const int d1d = 0,
220  const int q1d = 0)
221 {
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);
232  MFEM_FORALL(e, NE,
233  {
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)
240  {
241  for (int qz = 0; qz < Q1D; ++qz)
242  {
243  for (int qy = 0; qy < Q1D; ++qy)
244  {
245  for (int qx = 0; qx < Q1D; ++qx)
246  {
247  sol_xyz[qz][qy][qx] = 0.0;
248  }
249  }
250  }
251  for (int dz = 0; dz < D1D; ++dz)
252  {
253  double sol_xy[max_Q1D][max_Q1D];
254  for (int qy = 0; qy < Q1D; ++qy)
255  {
256  for (int qx = 0; qx < Q1D; ++qx)
257  {
258  sol_xy[qy][qx] = 0.0;
259  }
260  }
261  for (int dy = 0; dy < D1D; ++dy)
262  {
263  double sol_x[max_Q1D];
264  for (int qx = 0; qx < Q1D; ++qx)
265  {
266  sol_x[qx] = 0;
267  }
268  for (int dx = 0; dx < D1D; ++dx)
269  {
270  const double s = x(dx,dy,dz,c,e);
271  for (int qx = 0; qx < Q1D; ++qx)
272  {
273  sol_x[qx] += B(qx,dx) * s;
274  }
275  }
276  for (int qy = 0; qy < Q1D; ++qy)
277  {
278  const double wy = B(qy,dy);
279  for (int qx = 0; qx < Q1D; ++qx)
280  {
281  sol_xy[qy][qx] += wy * sol_x[qx];
282  }
283  }
284  }
285  for (int qz = 0; qz < Q1D; ++qz)
286  {
287  const double wz = B(qz,dz);
288  for (int qy = 0; qy < Q1D; ++qy)
289  {
290  for (int qx = 0; qx < Q1D; ++qx)
291  {
292  sol_xyz[qz][qy][qx] += wz * sol_xy[qy][qx];
293  }
294  }
295  }
296  }
297  for (int qz = 0; qz < Q1D; ++qz)
298  {
299  for (int qy = 0; qy < Q1D; ++qy)
300  {
301  for (int qx = 0; qx < Q1D; ++qx)
302  {
303  sol_xyz[qz][qy][qx] *= op(qx,qy,qz,e);
304  }
305  }
306  }
307  for (int qz = 0; qz < Q1D; ++qz)
308  {
309  double sol_xy[max_D1D][max_D1D];
310  for (int dy = 0; dy < D1D; ++dy)
311  {
312  for (int dx = 0; dx < D1D; ++dx)
313  {
314  sol_xy[dy][dx] = 0;
315  }
316  }
317  for (int qy = 0; qy < Q1D; ++qy)
318  {
319  double sol_x[max_D1D];
320  for (int dx = 0; dx < D1D; ++dx)
321  {
322  sol_x[dx] = 0;
323  }
324  for (int qx = 0; qx < Q1D; ++qx)
325  {
326  const double s = sol_xyz[qz][qy][qx];
327  for (int dx = 0; dx < D1D; ++dx)
328  {
329  sol_x[dx] += Bt(dx,qx) * s;
330  }
331  }
332  for (int dy = 0; dy < D1D; ++dy)
333  {
334  const double wy = Bt(dy,qy);
335  for (int dx = 0; dx < D1D; ++dx)
336  {
337  sol_xy[dy][dx] += wy * sol_x[dx];
338  }
339  }
340  }
341  for (int dz = 0; dz < D1D; ++dz)
342  {
343  const double wz = Bt(dz,qz);
344  for (int dy = 0; dy < D1D; ++dy)
345  {
346  for (int dx = 0; dx < D1D; ++dx)
347  {
348  y(dx,dy,dz,c,e) += wz * sol_xy[dy][dx];
349  }
350  }
351  }
352  }
353  }
354  });
355 }
356 
357 static void PAVectorMassApply(const int dim,
358  const int D1D,
359  const int Q1D,
360  const int NE,
361  const Array<double> &B,
362  const Array<double> &Bt,
363  const Vector &op,
364  const Vector &x,
365  Vector &y)
366 {
367  if (dim == 2)
368  {
369  return PAVectorMassApply2D(NE, B, Bt, op, x, y, D1D, Q1D);
370  }
371  if (dim == 3)
372  {
373  return PAVectorMassApply3D(NE, B, Bt, op, x, y, D1D, Q1D);
374  }
375  MFEM_ABORT("Unknown kernel.");
376 }
377 
378 void VectorMassIntegrator::AddMultPA(const Vector &x, Vector &y) const
379 {
380  if (DeviceCanUseCeed())
381  {
382  ceedOp->AddMult(x, y);
383  }
384  else
385  {
386  PAVectorMassApply(dim, dofs1D, quad1D, ne, maps->B, maps->Bt, pa_data, x, y);
387  }
388 }
389 
390 template<const int T_D1D = 0, const int T_Q1D = 0>
391 static void PAVectorMassAssembleDiagonal2D(const int NE,
392  const Array<double> &B_,
393  const Array<double> &Bt_,
394  const Vector &op_,
395  Vector &diag_,
396  const int d1d = 0,
397  const int q1d = 0)
398 {
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, "");
404  auto B = Reshape(B_.Read(), Q1D, D1D);
405  auto op = Reshape(op_.Read(), Q1D, Q1D, NE);
406  auto y = Reshape(diag_.ReadWrite(), D1D, D1D, VDIM, NE);
407  MFEM_FORALL(e, NE,
408  {
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;
413 
414  double temp[max_Q1D][max_D1D];
415  for (int qx = 0; qx < Q1D; ++qx)
416  {
417  for (int dy = 0; dy < D1D; ++dy)
418  {
419  temp[qx][dy] = 0.0;
420  for (int qy = 0; qy < Q1D; ++qy)
421  {
422  temp[qx][dy] += B(qy, dy) * B(qy, dy) * op(qx, qy, e);
423  }
424  }
425  }
426  for (int dy = 0; dy < D1D; ++dy)
427  {
428  for (int dx = 0; dx < D1D; ++dx)
429  {
430  double temp1 = 0.0;
431  for (int qx = 0; qx < Q1D; ++qx)
432  {
433  temp1 += B(qx, dx) * B(qx, dx) * temp[qx][dy];
434  }
435  y(dx, dy, 0, e) = temp1;
436  y(dx, dy, 1, e) = temp1;
437  }
438  }
439  });
440 }
441 
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_,
446  const Vector &op_,
447  Vector &diag_,
448  const int d1d = 0,
449  const int q1d = 0)
450 {
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);
459  MFEM_FORALL(e, NE,
460  {
461  const int D1D = T_D1D ? T_D1D : d1d; // nvcc workaround
462  const int Q1D = T_Q1D ? T_Q1D : q1d;
463  // the following variables are evaluated at compile time
464  constexpr int max_D1D = T_D1D ? T_D1D : MAX_D1D;
465  constexpr int max_Q1D = T_Q1D ? T_Q1D : MAX_Q1D;
466 
467  double temp[max_Q1D][max_Q1D][max_D1D];
468  for (int qx = 0; qx < Q1D; ++qx)
469  {
470  for (int qy = 0; qy < Q1D; ++qy)
471  {
472  for (int dz = 0; dz < D1D; ++dz)
473  {
474  temp[qx][qy][dz] = 0.0;
475  for (int qz = 0; qz < Q1D; ++qz)
476  {
477  temp[qx][qy][dz] += B(qz, dz) * B(qz, dz) * op(qx, qy, qz, e);
478  }
479  }
480  }
481  }
482  double temp2[max_Q1D][max_D1D][max_D1D];
483  for (int qx = 0; qx < Q1D; ++qx)
484  {
485  for (int dz = 0; dz < D1D; ++dz)
486  {
487  for (int dy = 0; dy < D1D; ++dy)
488  {
489  temp2[qx][dy][dz] = 0.0;
490  for (int qy = 0; qy < Q1D; ++qy)
491  {
492  temp2[qx][dy][dz] += B(qy, dy) * B(qy, dy) * temp[qx][qy][dz];
493  }
494  }
495  }
496  }
497  for (int dz = 0; dz < D1D; ++dz)
498  {
499  for (int dy = 0; dy < D1D; ++dy)
500  {
501  for (int dx = 0; dx < D1D; ++dx)
502  {
503  double temp3 = 0.0;
504  for (int qx = 0; qx < Q1D; ++qx)
505  {
506  temp3 += B(qx, dx) * B(qx, dx)
507  * temp2[qx][dy][dz];
508  }
509  y(dx, dy, dz, 0, e) = temp3;
510  y(dx, dy, dz, 1, e) = temp3;
511  y(dx, dy, dz, 2, e) = temp3;
512  }
513  }
514  }
515  });
516 }
517 
518 static void PAVectorMassAssembleDiagonal(const int dim,
519  const int D1D,
520  const int Q1D,
521  const int NE,
522  const Array<double> &B,
523  const Array<double> &Bt,
524  const Vector &op,
525  Vector &y)
526 {
527  if (dim == 2)
528  {
529  return PAVectorMassAssembleDiagonal2D(NE, B, Bt, op, y, D1D, Q1D);
530  }
531  else if (dim == 3)
532  {
533  return PAVectorMassAssembleDiagonal3D(NE, B, Bt, op, y, D1D, Q1D);
534  }
535  MFEM_ABORT("Dimension not implemented.");
536 }
537 
538 void VectorMassIntegrator::AssembleDiagonalPA(Vector &diag)
539 {
540  if (DeviceCanUseCeed())
541  {
542  ceedOp->GetDiagonal(diag);
543  }
544  else
545  {
546  PAVectorMassAssembleDiagonal(dim,
547  dofs1D,
548  quad1D,
549  ne,
550  maps->B,
551  maps->Bt,
552  pa_data,
553  diag);
554  }
555 }
556 
557 } // namespace mfem
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
Abstract class for all finite elements.
Definition: fe_base.hpp:235
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
bool IsVariableOrder() const
Returns true if the space contains elements of varying polynomial orders.
Definition: fespace.hpp:463
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
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.
Definition: mesh.cpp:840
int GetNumGeometries(int dim) const
Return the number of geometries of the given dimension present in the mesh.
Definition: mesh.cpp:5908
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
Definition: fe_base.hpp:170
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:923
const Array< double > & GetWeights() const
Return the quadrature weights in a contiguous array.
Definition: intrules.cpp:83
const int MAX_Q1D
Definition: forall.hpp:29
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:441
Represent a MassIntegrator with AssemblyLevel::Partial using libCEED.
Definition: mass.hpp:26
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:304
int Dimension() const
Definition: mesh.hpp:1006
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
bool DeviceCanUseCeed()
Function that determines if a CEED kernel should be used, based on the current mfem::Device configura...
Definition: util.cpp:33
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:366
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).
Definition: vector.hpp:465
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:348
int dim
Definition: ex24.cpp:53
const int MAX_D1D
Definition: forall.hpp:28
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2781
Vector data type.
Definition: vector.hpp:60
RefCoord s[3]
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:449
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition: dtensor.hpp:131