MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
bilininteg_hcurl.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, 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"
15 #include "libceed/mass.hpp"
16 
17 using namespace std;
18 
19 namespace mfem
20 {
21 
22 // Local maximum size of dofs and quads in 1D
23 constexpr int HCURL_MAX_D1D = 5;
24 constexpr int HCURL_MAX_Q1D = 6;
25 
26 // PA H(curl) Mass Assemble 2D kernel
27 static void PAHcurlSetup2D(const int Q1D,
28  const int NE,
29  const Array<double> &w,
30  const Vector &j,
31  Vector &_coeff,
32  Vector &op)
33 {
34  const int NQ = Q1D*Q1D;
35  auto W = w.Read();
36 
37  auto J = Reshape(j.Read(), NQ, 2, 2, NE);
38  auto coeff = Reshape(_coeff.Read(), NQ, NE);
39  auto y = Reshape(op.Write(), NQ, 3, NE);
40 
41  MFEM_FORALL(e, NE,
42  {
43  for (int q = 0; q < NQ; ++q)
44  {
45  const double J11 = J(q,0,0,e);
46  const double J21 = J(q,1,0,e);
47  const double J12 = J(q,0,1,e);
48  const double J22 = J(q,1,1,e);
49  const double c_detJ = W[q] * coeff(q, e) / ((J11*J22)-(J21*J12));
50  y(q,0,e) = c_detJ * (J12*J12 + J22*J22); // 1,1
51  y(q,1,e) = -c_detJ * (J12*J11 + J22*J21); // 1,2
52  y(q,2,e) = c_detJ * (J11*J11 + J21*J21); // 2,2
53  }
54  });
55 }
56 
57 // PA H(curl) Mass Assemble 3D kernel
58 static void PAHcurlSetup3D(const int Q1D,
59  const int NE,
60  const Array<double> &w,
61  const Vector &j,
62  Vector &_coeff,
63  Vector &op)
64 {
65  const int NQ = Q1D*Q1D*Q1D;
66  auto W = w.Read();
67  auto J = Reshape(j.Read(), NQ, 3, 3, NE);
68  auto coeff = Reshape(_coeff.Read(), NQ, NE);
69  auto y = Reshape(op.Write(), NQ, 6, NE);
70 
71  MFEM_FORALL(e, NE,
72  {
73  for (int q = 0; q < NQ; ++q)
74  {
75  const double J11 = J(q,0,0,e);
76  const double J21 = J(q,1,0,e);
77  const double J31 = J(q,2,0,e);
78  const double J12 = J(q,0,1,e);
79  const double J22 = J(q,1,1,e);
80  const double J32 = J(q,2,1,e);
81  const double J13 = J(q,0,2,e);
82  const double J23 = J(q,1,2,e);
83  const double J33 = J(q,2,2,e);
84  const double detJ = J11 * (J22 * J33 - J32 * J23) -
85  /* */ J21 * (J12 * J33 - J32 * J13) +
86  /* */ J31 * (J12 * J23 - J22 * J13);
87  const double c_detJ = W[q] * coeff(q, e) / detJ;
88  // adj(J)
89  const double A11 = (J22 * J33) - (J23 * J32);
90  const double A12 = (J32 * J13) - (J12 * J33);
91  const double A13 = (J12 * J23) - (J22 * J13);
92  const double A21 = (J31 * J23) - (J21 * J33);
93  const double A22 = (J11 * J33) - (J13 * J31);
94  const double A23 = (J21 * J13) - (J11 * J23);
95  const double A31 = (J21 * J32) - (J31 * J22);
96  const double A32 = (J31 * J12) - (J11 * J32);
97  const double A33 = (J11 * J22) - (J12 * J21);
98  // detJ J^{-1} J^{-T} = (1/detJ) adj(J) adj(J)^T
99  y(q,0,e) = c_detJ * (A11*A11 + A12*A12 + A13*A13); // 1,1
100  y(q,1,e) = c_detJ * (A11*A21 + A12*A22 + A13*A23); // 2,1
101  y(q,2,e) = c_detJ * (A11*A31 + A12*A32 + A13*A33); // 3,1
102  y(q,3,e) = c_detJ * (A21*A21 + A22*A22 + A23*A23); // 2,2
103  y(q,4,e) = c_detJ * (A21*A31 + A22*A32 + A23*A33); // 3,2
104  y(q,5,e) = c_detJ * (A31*A31 + A32*A32 + A33*A33); // 3,3
105  }
106  });
107 }
108 
109 void VectorFEMassIntegrator::AssemblePA(const FiniteElementSpace &fes)
110 {
111  // Assumes tensor-product elements
112  Mesh *mesh = fes.GetMesh();
113  const FiniteElement *fel = fes.GetFE(0);
114 
115  const VectorTensorFiniteElement *el =
116  dynamic_cast<const VectorTensorFiniteElement*>(fel);
117  MFEM_VERIFY(el != NULL, "Only VectorTensorFiniteElement is supported!");
118 
119  const IntegrationRule *ir
120  = IntRule ? IntRule : &MassIntegrator::GetRule(*el, *el,
121  *mesh->GetElementTransformation(0));
122  const int dims = el->GetDim();
123  MFEM_VERIFY(dims == 2 || dims == 3, "");
124 
125  const int symmDims = (dims * (dims + 1)) / 2; // 1x1: 1, 2x2: 3, 3x3: 6
126  const int nq = ir->GetNPoints();
127  dim = mesh->Dimension();
128  MFEM_VERIFY(dim == 2 || dim == 3, "");
129 
130  ne = fes.GetNE();
131  geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS);
132  mapsC = &el->GetDofToQuad(*ir, DofToQuad::TENSOR);
133  mapsO = &el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
134  dofs1D = mapsC->ndof;
135  quad1D = mapsC->nqpt;
136 
137  MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
138 
139  pa_data.SetSize(symmDims * nq * ne, Device::GetMemoryType());
140 
141  Vector coeff(ne * nq);
142  coeff = 1.0;
143  if (Q)
144  {
145  for (int e=0; e<ne; ++e)
146  {
148  for (int p=0; p<nq; ++p)
149  {
150  coeff[p + (e * nq)] = Q->Eval(*tr, ir->IntPoint(p));
151  }
152  }
153  }
154 
155  if (el->GetDerivType() == mfem::FiniteElement::CURL && dim == 3)
156  {
157  PAHcurlSetup3D(quad1D, ne, ir->GetWeights(), geom->J,
158  coeff, pa_data);
159  }
160  else if (el->GetDerivType() == mfem::FiniteElement::CURL && dim == 2)
161  {
162  PAHcurlSetup2D(quad1D, ne, ir->GetWeights(), geom->J,
163  coeff, pa_data);
164  }
165  else
166  {
167  MFEM_ABORT("Unknown kernel.");
168  }
169 }
170 
171 static void PAHcurlMassApply2D(const int D1D,
172  const int Q1D,
173  const int NE,
174  const Array<double> &_Bo,
175  const Array<double> &_Bc,
176  const Array<double> &_Bot,
177  const Array<double> &_Bct,
178  const Vector &_op,
179  const Vector &_x,
180  Vector &_y)
181 {
182  constexpr static int VDIM = 2;
183 
184  auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
185  auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
186  auto Bot = Reshape(_Bot.Read(), D1D-1, Q1D);
187  auto Bct = Reshape(_Bct.Read(), D1D, Q1D);
188  auto op = Reshape(_op.Read(), Q1D, Q1D, 3, NE);
189  auto x = Reshape(_x.Read(), 2*(D1D-1)*D1D, NE);
190  auto y = Reshape(_y.ReadWrite(), 2*(D1D-1)*D1D, NE);
191 
192  MFEM_FORALL(e, NE,
193  {
194  double mass[MAX_Q1D][MAX_Q1D][VDIM];
195 
196  for (int qy = 0; qy < Q1D; ++qy)
197  {
198  for (int qx = 0; qx < Q1D; ++qx)
199  {
200  for (int c = 0; c < VDIM; ++c)
201  {
202  mass[qy][qx][c] = 0.0;
203  }
204  }
205  }
206 
207  int osc = 0;
208 
209  for (int c = 0; c < VDIM; ++c) // loop over x, y components
210  {
211  const int D1Dy = (c == 1) ? D1D - 1 : D1D;
212  const int D1Dx = (c == 0) ? D1D - 1 : D1D;
213 
214  for (int dy = 0; dy < D1Dy; ++dy)
215  {
216  double massX[MAX_Q1D];
217  for (int qx = 0; qx < Q1D; ++qx)
218  {
219  massX[qx] = 0.0;
220  }
221 
222  for (int dx = 0; dx < D1Dx; ++dx)
223  {
224  const double t = x(dx + (dy * D1Dx) + osc, e);
225  for (int qx = 0; qx < Q1D; ++qx)
226  {
227  massX[qx] += t * ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
228  }
229  }
230 
231  for (int qy = 0; qy < Q1D; ++qy)
232  {
233  const double wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
234  for (int qx = 0; qx < Q1D; ++qx)
235  {
236  mass[qy][qx][c] += massX[qx] * wy;
237  }
238  }
239  }
240 
241  osc += D1Dx * D1Dy;
242  } // loop (c) over components
243 
244  // Apply D operator.
245  for (int qy = 0; qy < Q1D; ++qy)
246  {
247  for (int qx = 0; qx < Q1D; ++qx)
248  {
249  const double O11 = op(qx,qy,0,e);
250  const double O12 = op(qx,qy,1,e);
251  const double O22 = op(qx,qy,2,e);
252  const double massX = mass[qy][qx][0];
253  const double massY = mass[qy][qx][1];
254  mass[qy][qx][0] = (O11*massX)+(O12*massY);
255  mass[qy][qx][1] = (O12*massX)+(O22*massY);
256  }
257  }
258 
259  for (int qy = 0; qy < Q1D; ++qy)
260  {
261  osc = 0;
262 
263  for (int c = 0; c < VDIM; ++c) // loop over x, y components
264  {
265  const int D1Dy = (c == 1) ? D1D - 1 : D1D;
266  const int D1Dx = (c == 0) ? D1D - 1 : D1D;
267 
268  double massX[MAX_D1D];
269  for (int dx = 0; dx < D1Dx; ++dx)
270  {
271  massX[dx] = 0;
272  }
273  for (int qx = 0; qx < Q1D; ++qx)
274  {
275  for (int dx = 0; dx < D1Dx; ++dx)
276  {
277  massX[dx] += mass[qy][qx][c] * ((c == 0) ? Bot(dx,qx) : Bct(dx,qx));
278  }
279  }
280 
281  for (int dy = 0; dy < D1Dy; ++dy)
282  {
283  const double wy = (c == 1) ? Bot(dy,qy) : Bct(dy,qy);
284 
285  for (int dx = 0; dx < D1Dx; ++dx)
286  {
287  y(dx + (dy * D1Dx) + osc, e) += massX[dx] * wy;
288  }
289  }
290 
291  osc += D1Dx * D1Dy;
292  } // loop c
293  } // loop qy
294  }); // end of element loop
295 }
296 
297 static void PAHcurlMassAssembleDiagonal2D(const int D1D,
298  const int Q1D,
299  const int NE,
300  const Array<double> &_Bo,
301  const Array<double> &_Bc,
302  const Vector &_op,
303  Vector &_diag)
304 {
305  constexpr static int VDIM = 2;
306 
307  auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
308  auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
309  auto op = Reshape(_op.Read(), Q1D, Q1D, 3, NE);
310  auto diag = Reshape(_diag.ReadWrite(), 2*(D1D-1)*D1D, NE);
311 
312  MFEM_FORALL(e, NE,
313  {
314  int osc = 0;
315 
316  for (int c = 0; c < VDIM; ++c) // loop over x, y components
317  {
318  const int D1Dy = (c == 1) ? D1D - 1 : D1D;
319  const int D1Dx = (c == 0) ? D1D - 1 : D1D;
320 
321  double mass[MAX_Q1D];
322 
323  for (int dy = 0; dy < D1Dy; ++dy)
324  {
325  for (int qx = 0; qx < Q1D; ++qx)
326  {
327  mass[qx] = 0.0;
328  for (int qy = 0; qy < Q1D; ++qy)
329  {
330  const double wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
331 
332  mass[qx] += wy * wy * ((c == 0) ? op(qx,qy,0,e) : op(qx,qy,2,e));
333  }
334  }
335 
336  for (int dx = 0; dx < D1Dx; ++dx)
337  {
338  for (int qx = 0; qx < Q1D; ++qx)
339  {
340  const double wx = ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
341  diag(dx + (dy * D1Dx) + osc, e) += mass[qx] * wx * wx;
342  }
343  }
344  }
345 
346  osc += D1Dx * D1Dy;
347  } // loop c
348  }); // end of element loop
349 }
350 
351 template<int MAX_D1D = HCURL_MAX_D1D, int MAX_Q1D = HCURL_MAX_Q1D>
352 static void PAHcurlMassAssembleDiagonal3D(const int D1D,
353  const int Q1D,
354  const int NE,
355  const Array<double> &_Bo,
356  const Array<double> &_Bc,
357  const Vector &_op,
358  Vector &_diag)
359 {
360  MFEM_VERIFY(D1D <= MAX_D1D, "Error: D1D > MAX_D1D");
361  MFEM_VERIFY(Q1D <= MAX_Q1D, "Error: Q1D > MAX_Q1D");
362  constexpr static int VDIM = 3;
363 
364  auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
365  auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
366  auto op = Reshape(_op.Read(), Q1D, Q1D, Q1D, 6, NE);
367  auto diag = Reshape(_diag.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
368 
369  MFEM_FORALL(e, NE,
370  {
371  int osc = 0;
372 
373  for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
374  {
375  const int D1Dz = (c == 2) ? D1D - 1 : D1D;
376  const int D1Dy = (c == 1) ? D1D - 1 : D1D;
377  const int D1Dx = (c == 0) ? D1D - 1 : D1D;
378 
379  const int opc = (c == 0) ? 0 : ((c == 1) ? 3 : 5);
380 
381  double mass[MAX_Q1D];
382 
383  for (int dz = 0; dz < D1Dz; ++dz)
384  {
385  for (int dy = 0; dy < D1Dy; ++dy)
386  {
387  for (int qx = 0; qx < Q1D; ++qx)
388  {
389  mass[qx] = 0.0;
390  for (int qy = 0; qy < Q1D; ++qy)
391  {
392  const double wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
393 
394  for (int qz = 0; qz < Q1D; ++qz)
395  {
396  const double wz = (c == 2) ? Bo(qz,dz) : Bc(qz,dz);
397 
398  mass[qx] += wy * wy * wz * wz * op(qx,qy,qz,opc,e);
399  }
400  }
401  }
402 
403  for (int dx = 0; dx < D1Dx; ++dx)
404  {
405  for (int qx = 0; qx < Q1D; ++qx)
406  {
407  const double wx = ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
408  diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) += mass[qx] * wx * wx;
409  }
410  }
411  }
412  }
413 
414  osc += D1Dx * D1Dy * D1Dz;
415  } // loop c
416  }); // end of element loop
417 }
418 
419 void VectorFEMassIntegrator::AssembleDiagonalPA(Vector& diag)
420 {
421  if (dim == 3)
422  PAHcurlMassAssembleDiagonal3D(dofs1D, quad1D, ne,
423  mapsO->B, mapsC->B, pa_data, diag);
424  else
425  PAHcurlMassAssembleDiagonal2D(dofs1D, quad1D, ne,
426  mapsO->B, mapsC->B, pa_data, diag);
427 }
428 
429 template<int MAX_D1D = HCURL_MAX_D1D, int MAX_Q1D = HCURL_MAX_Q1D>
430 static void PAHcurlMassApply3D(const int D1D,
431  const int Q1D,
432  const int NE,
433  const Array<double> &_Bo,
434  const Array<double> &_Bc,
435  const Array<double> &_Bot,
436  const Array<double> &_Bct,
437  const Vector &_op,
438  const Vector &_x,
439  Vector &_y)
440 {
441  MFEM_VERIFY(D1D <= MAX_D1D, "Error: D1D > MAX_D1D");
442  MFEM_VERIFY(Q1D <= MAX_Q1D, "Error: Q1D > MAX_Q1D");
443  constexpr static int VDIM = 3;
444 
445  auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
446  auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
447  auto Bot = Reshape(_Bot.Read(), D1D-1, Q1D);
448  auto Bct = Reshape(_Bct.Read(), D1D, Q1D);
449  auto op = Reshape(_op.Read(), Q1D, Q1D, Q1D, 6, NE);
450  auto x = Reshape(_x.Read(), 3*(D1D-1)*D1D*D1D, NE);
451  auto y = Reshape(_y.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
452 
453  MFEM_FORALL(e, NE,
454  {
455  double mass[MAX_Q1D][MAX_Q1D][MAX_Q1D][VDIM];
456 
457  for (int qz = 0; qz < Q1D; ++qz)
458  {
459  for (int qy = 0; qy < Q1D; ++qy)
460  {
461  for (int qx = 0; qx < Q1D; ++qx)
462  {
463  for (int c = 0; c < VDIM; ++c)
464  {
465  mass[qz][qy][qx][c] = 0.0;
466  }
467  }
468  }
469  }
470 
471  int osc = 0;
472 
473  for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
474  {
475  const int D1Dz = (c == 2) ? D1D - 1 : D1D;
476  const int D1Dy = (c == 1) ? D1D - 1 : D1D;
477  const int D1Dx = (c == 0) ? D1D - 1 : D1D;
478 
479  for (int dz = 0; dz < D1Dz; ++dz)
480  {
481  double massXY[MAX_Q1D][MAX_Q1D];
482  for (int qy = 0; qy < Q1D; ++qy)
483  {
484  for (int qx = 0; qx < Q1D; ++qx)
485  {
486  massXY[qy][qx] = 0.0;
487  }
488  }
489 
490  for (int dy = 0; dy < D1Dy; ++dy)
491  {
492  double massX[MAX_Q1D];
493  for (int qx = 0; qx < Q1D; ++qx)
494  {
495  massX[qx] = 0.0;
496  }
497 
498  for (int dx = 0; dx < D1Dx; ++dx)
499  {
500  const double t = x(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
501  for (int qx = 0; qx < Q1D; ++qx)
502  {
503  massX[qx] += t * ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
504  }
505  }
506 
507  for (int qy = 0; qy < Q1D; ++qy)
508  {
509  const double wy = (c == 1) ? Bo(qy,dy) : Bc(qy,dy);
510  for (int qx = 0; qx < Q1D; ++qx)
511  {
512  const double wx = massX[qx];
513  massXY[qy][qx] += wx * wy;
514  }
515  }
516  }
517 
518  for (int qz = 0; qz < Q1D; ++qz)
519  {
520  const double wz = (c == 2) ? Bo(qz,dz) : Bc(qz,dz);
521  for (int qy = 0; qy < Q1D; ++qy)
522  {
523  for (int qx = 0; qx < Q1D; ++qx)
524  {
525  mass[qz][qy][qx][c] += massXY[qy][qx] * wz;
526  }
527  }
528  }
529  }
530 
531  osc += D1Dx * D1Dy * D1Dz;
532  } // loop (c) over components
533 
534  // Apply D operator.
535  for (int qz = 0; qz < Q1D; ++qz)
536  {
537  for (int qy = 0; qy < Q1D; ++qy)
538  {
539  for (int qx = 0; qx < Q1D; ++qx)
540  {
541  const double O11 = op(qx,qy,qz,0,e);
542  const double O12 = op(qx,qy,qz,1,e);
543  const double O13 = op(qx,qy,qz,2,e);
544  const double O22 = op(qx,qy,qz,3,e);
545  const double O23 = op(qx,qy,qz,4,e);
546  const double O33 = op(qx,qy,qz,5,e);
547  const double massX = mass[qz][qy][qx][0];
548  const double massY = mass[qz][qy][qx][1];
549  const double massZ = mass[qz][qy][qx][2];
550  mass[qz][qy][qx][0] = (O11*massX)+(O12*massY)+(O13*massZ);
551  mass[qz][qy][qx][1] = (O12*massX)+(O22*massY)+(O23*massZ);
552  mass[qz][qy][qx][2] = (O13*massX)+(O23*massY)+(O33*massZ);
553  }
554  }
555  }
556 
557  for (int qz = 0; qz < Q1D; ++qz)
558  {
559  double massXY[MAX_D1D][MAX_D1D];
560 
561  osc = 0;
562 
563  for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
564  {
565  const int D1Dz = (c == 2) ? D1D - 1 : D1D;
566  const int D1Dy = (c == 1) ? D1D - 1 : D1D;
567  const int D1Dx = (c == 0) ? D1D - 1 : D1D;
568 
569  for (int dy = 0; dy < D1Dy; ++dy)
570  {
571  for (int dx = 0; dx < D1Dx; ++dx)
572  {
573  massXY[dy][dx] = 0;
574  }
575  }
576  for (int qy = 0; qy < Q1D; ++qy)
577  {
578  double massX[MAX_D1D];
579  for (int dx = 0; dx < D1Dx; ++dx)
580  {
581  massX[dx] = 0;
582  }
583  for (int qx = 0; qx < Q1D; ++qx)
584  {
585  for (int dx = 0; dx < D1Dx; ++dx)
586  {
587  massX[dx] += mass[qz][qy][qx][c] * ((c == 0) ? Bot(dx,qx) : Bct(dx,qx));
588  }
589  }
590  for (int dy = 0; dy < D1Dy; ++dy)
591  {
592  const double wy = (c == 1) ? Bot(dy,qy) : Bct(dy,qy);
593  for (int dx = 0; dx < D1Dx; ++dx)
594  {
595  massXY[dy][dx] += massX[dx] * wy;
596  }
597  }
598  }
599 
600  for (int dz = 0; dz < D1Dz; ++dz)
601  {
602  const double wz = (c == 2) ? Bot(dz,qz) : Bct(dz,qz);
603  for (int dy = 0; dy < D1Dy; ++dy)
604  {
605  for (int dx = 0; dx < D1Dx; ++dx)
606  {
607  y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) += massXY[dy][dx] * wz;
608  }
609  }
610  }
611 
612  osc += D1Dx * D1Dy * D1Dz;
613  } // loop c
614  } // loop qz
615  }); // end of element loop
616 }
617 
618 void VectorFEMassIntegrator::AddMultPA(const Vector &x, Vector &y) const
619 {
620  if (dim == 3)
621  {
622  PAHcurlMassApply3D(dofs1D, quad1D, ne, mapsO->B, mapsC->B, mapsO->Bt,
623  mapsC->Bt, pa_data, x, y);
624  }
625  else
626  {
627  PAHcurlMassApply2D(dofs1D, quad1D, ne, mapsO->B, mapsC->B, mapsO->Bt,
628  mapsC->Bt, pa_data, x, y);
629  }
630 }
631 
632 // PA H(curl) curl-curl assemble 2D kernel
633 static void PACurlCurlSetup2D(const int Q1D,
634  const int NE,
635  const Array<double> &w,
636  const Vector &j,
637  Vector &_coeff,
638  Vector &op)
639 {
640  const int NQ = Q1D*Q1D;
641  auto W = w.Read();
642  auto J = Reshape(j.Read(), NQ, 2, 2, NE);
643  auto coeff = Reshape(_coeff.Read(), NQ, NE);
644  auto y = Reshape(op.Write(), NQ, NE);
645  MFEM_FORALL(e, NE,
646  {
647  for (int q = 0; q < NQ; ++q)
648  {
649  const double J11 = J(q,0,0,e);
650  const double J21 = J(q,1,0,e);
651  const double J12 = J(q,0,1,e);
652  const double J22 = J(q,1,1,e);
653  const double detJ = (J11*J22)-(J21*J12);
654  y(q,e) = W[q] * coeff(q,e) / detJ;
655  }
656  });
657 }
658 
659 // PA H(curl) curl-curl assemble 3D kernel
660 static void PACurlCurlSetup3D(const int Q1D,
661  const int NE,
662  const Array<double> &w,
663  const Vector &j,
664  Vector &_coeff,
665  Vector &op)
666 {
667  const int NQ = Q1D*Q1D*Q1D;
668  auto W = w.Read();
669  auto J = Reshape(j.Read(), NQ, 3, 3, NE);
670  auto coeff = Reshape(_coeff.Read(), NQ, NE);
671  auto y = Reshape(op.Write(), NQ, 6, NE);
672  MFEM_FORALL(e, NE,
673  {
674  for (int q = 0; q < NQ; ++q)
675  {
676  const double J11 = J(q,0,0,e);
677  const double J21 = J(q,1,0,e);
678  const double J31 = J(q,2,0,e);
679  const double J12 = J(q,0,1,e);
680  const double J22 = J(q,1,1,e);
681  const double J32 = J(q,2,1,e);
682  const double J13 = J(q,0,2,e);
683  const double J23 = J(q,1,2,e);
684  const double J33 = J(q,2,2,e);
685  const double detJ = J11 * (J22 * J33 - J32 * J23) -
686  /* */ J21 * (J12 * J33 - J32 * J13) +
687  /* */ J31 * (J12 * J23 - J22 * J13);
688 
689  // set y to the 6 entries of J^T J / det^2
690  const double c_detJ = W[q] * coeff(q,e) / detJ;
691 
692  y(q,0,e) = c_detJ * (J11*J11 + J21*J21 + J31*J31); // 1,1
693  y(q,1,e) = c_detJ * (J11*J12 + J21*J22 + J31*J32); // 1,2
694  y(q,2,e) = c_detJ * (J11*J13 + J21*J23 + J31*J33); // 1,3
695  y(q,3,e) = c_detJ * (J12*J12 + J22*J22 + J32*J32); // 2,2
696  y(q,4,e) = c_detJ * (J12*J13 + J22*J23 + J32*J33); // 2,3
697  y(q,5,e) = c_detJ * (J13*J13 + J23*J23 + J33*J33); // 3,3
698  }
699  });
700 }
701 
702 void CurlCurlIntegrator::AssemblePA(const FiniteElementSpace &fes)
703 {
704  // Assumes tensor-product elements
705  Mesh *mesh = fes.GetMesh();
706  const FiniteElement *fel = fes.GetFE(0);
707 
708  const VectorTensorFiniteElement *el =
709  dynamic_cast<const VectorTensorFiniteElement*>(fel);
710  MFEM_VERIFY(el != NULL, "Only VectorTensorFiniteElement is supported!");
711 
712  const IntegrationRule *ir
713  = IntRule ? IntRule : &MassIntegrator::GetRule(*el, *el,
714  *mesh->GetElementTransformation(0));
715  const int dims = el->GetDim();
716  MFEM_VERIFY(dims == 2 || dims == 3, "");
717 
718  const int nq = ir->GetNPoints();
719  dim = mesh->Dimension();
720  MFEM_VERIFY(dim == 2 || dim == 3, "");
721 
722  ne = fes.GetNE();
723  geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS);
724  mapsC = &el->GetDofToQuad(*ir, DofToQuad::TENSOR);
725  mapsO = &el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
726  dofs1D = mapsC->ndof;
727  quad1D = mapsC->nqpt;
728 
729  MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
730 
731  const int ndata = (dim == 2) ? 1 : 6;
732  pa_data.SetSize(ndata * nq * ne, Device::GetMemoryType());
733 
734  Vector coeff(ne * nq);
735  coeff = 1.0;
736  if (Q)
737  {
738  for (int e=0; e<ne; ++e)
739  {
740  ElementTransformation *tr = mesh->GetElementTransformation(e);
741  for (int p=0; p<nq; ++p)
742  {
743  coeff[p + (e * nq)] = Q->Eval(*tr, ir->IntPoint(p));
744  }
745  }
746  }
747 
748  if (el->GetDerivType() == mfem::FiniteElement::CURL && dim == 3)
749  {
750  // pa_data_2.SetSize(6 * nq * ne, Device::GetMemoryType());
751 
752  PACurlCurlSetup3D(quad1D, ne, ir->GetWeights(), geom->J,
753  coeff, pa_data);
754  }
755  else if (el->GetDerivType() == mfem::FiniteElement::CURL && dim == 2)
756  {
757  PACurlCurlSetup2D(quad1D, ne, ir->GetWeights(), geom->J,
758  coeff, pa_data);
759  }
760  else
761  {
762  MFEM_ABORT("Unknown kernel.");
763  }
764 }
765 
766 static void PACurlCurlApply2D(const int D1D,
767  const int Q1D,
768  const int NE,
769  const Array<double> &_Bo,
770  const Array<double> &_Bot,
771  const Array<double> &_Gc,
772  const Array<double> &_Gct,
773  const Vector &_op,
774  const Vector &_x,
775  Vector &_y)
776 {
777  constexpr static int VDIM = 2;
778 
779  auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
780  auto Bot = Reshape(_Bot.Read(), D1D-1, Q1D);
781  auto Gc = Reshape(_Gc.Read(), Q1D, D1D);
782  auto Gct = Reshape(_Gct.Read(), D1D, Q1D);
783  auto op = Reshape(_op.Read(), Q1D, Q1D, NE);
784  auto x = Reshape(_x.Read(), 2*(D1D-1)*D1D, NE);
785  auto y = Reshape(_y.ReadWrite(), 2*(D1D-1)*D1D, NE);
786 
787  MFEM_FORALL(e, NE,
788  {
789  double curl[MAX_Q1D][MAX_Q1D];
790 
791  // curl[qy][qx] will be computed as du_y/dx - du_x/dy
792 
793  for (int qy = 0; qy < Q1D; ++qy)
794  {
795  for (int qx = 0; qx < Q1D; ++qx)
796  {
797  curl[qy][qx] = 0;
798  }
799  }
800 
801  int osc = 0;
802 
803  for (int c = 0; c < VDIM; ++c) // loop over x, y components
804  {
805  const int D1Dy = (c == 1) ? D1D - 1 : D1D;
806  const int D1Dx = (c == 0) ? D1D - 1 : D1D;
807 
808  for (int dy = 0; dy < D1Dy; ++dy)
809  {
810  double gradX[MAX_Q1D];
811  for (int qx = 0; qx < Q1D; ++qx)
812  {
813  gradX[qx] = 0;
814  }
815 
816  for (int dx = 0; dx < D1Dx; ++dx)
817  {
818  const double t = x(dx + (dy * D1Dx) + osc, e);
819  for (int qx = 0; qx < Q1D; ++qx)
820  {
821  gradX[qx] += t * ((c == 0) ? Bo(qx,dx) : Gc(qx,dx));
822  }
823  }
824 
825  for (int qy = 0; qy < Q1D; ++qy)
826  {
827  const double wy = (c == 0) ? -Gc(qy,dy) : Bo(qy,dy);
828  for (int qx = 0; qx < Q1D; ++qx)
829  {
830  curl[qy][qx] += gradX[qx] * wy;
831  }
832  }
833  }
834 
835  osc += D1Dx * D1Dy;
836  } // loop (c) over components
837 
838  // Apply D operator.
839  for (int qy = 0; qy < Q1D; ++qy)
840  {
841  for (int qx = 0; qx < Q1D; ++qx)
842  {
843  curl[qy][qx] *= op(qx,qy,e);
844  }
845  }
846 
847  for (int qy = 0; qy < Q1D; ++qy)
848  {
849  osc = 0;
850 
851  for (int c = 0; c < VDIM; ++c) // loop over x, y components
852  {
853  const int D1Dy = (c == 1) ? D1D - 1 : D1D;
854  const int D1Dx = (c == 0) ? D1D - 1 : D1D;
855 
856  double gradX[MAX_D1D];
857  for (int dx = 0; dx < D1Dx; ++dx)
858  {
859  gradX[dx] = 0;
860  }
861  for (int qx = 0; qx < Q1D; ++qx)
862  {
863  for (int dx = 0; dx < D1Dx; ++dx)
864  {
865  gradX[dx] += curl[qy][qx] * ((c == 0) ? Bot(dx,qx) : Gct(dx,qx));
866  }
867  }
868  for (int dy = 0; dy < D1Dy; ++dy)
869  {
870  const double wy = (c == 0) ? -Gct(dy,qy) : Bot(dy,qy);
871 
872  for (int dx = 0; dx < D1Dx; ++dx)
873  {
874  y(dx + (dy * D1Dx) + osc, e) += gradX[dx] * wy;
875  }
876  }
877 
878  osc += D1Dx * D1Dy;
879  } // loop c
880  } // loop qy
881  }); // end of element loop
882 }
883 
884 template<int MAX_D1D = HCURL_MAX_D1D, int MAX_Q1D = HCURL_MAX_Q1D>
885 static void PACurlCurlApply3D(const int D1D,
886  const int Q1D,
887  const int NE,
888  const Array<double> &_Bo,
889  const Array<double> &_Bc,
890  const Array<double> &_Bot,
891  const Array<double> &_Bct,
892  const Array<double> &_Gc,
893  const Array<double> &_Gct,
894  const Vector &_op,
895  const Vector &_x,
896  Vector &_y)
897 {
898  MFEM_VERIFY(D1D <= MAX_D1D, "Error: D1D > MAX_D1D");
899  MFEM_VERIFY(Q1D <= MAX_Q1D, "Error: Q1D > MAX_Q1D");
900  // Using (\nabla\times u) F = 1/det(dF) dF \hat{\nabla}\times\hat{u} (p. 78 of Monk), we get
901  // (\nabla\times u) \cdot (\nabla\times v) = 1/det(dF)^2 \hat{\nabla}\times\hat{u}^T dF^T dF \hat{\nabla}\times\hat{v}
902  // If c = 0, \hat{\nabla}\times\hat{u} reduces to [0, (u_0)_{x_2}, -(u_0)_{x_1}]
903  // If c = 1, \hat{\nabla}\times\hat{u} reduces to [-(u_1)_{x_2}, 0, (u_1)_{x_0}]
904  // If c = 2, \hat{\nabla}\times\hat{u} reduces to [(u_2)_{x_1}, -(u_2)_{x_0}, 0]
905 
906  constexpr static int VDIM = 3;
907 
908  auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
909  auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
910  auto Bot = Reshape(_Bot.Read(), D1D-1, Q1D);
911  auto Bct = Reshape(_Bct.Read(), D1D, Q1D);
912  auto Gc = Reshape(_Gc.Read(), Q1D, D1D);
913  auto Gct = Reshape(_Gct.Read(), D1D, Q1D);
914  auto op = Reshape(_op.Read(), Q1D, Q1D, Q1D, 6, NE);
915  auto x = Reshape(_x.Read(), 3*(D1D-1)*D1D*D1D, NE);
916  auto y = Reshape(_y.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
917 
918  MFEM_FORALL(e, NE,
919  {
920  double curl[MAX_Q1D][MAX_Q1D][MAX_Q1D][VDIM];
921  // curl[qz][qy][qx] will be computed as the vector curl at each quadrature point.
922 
923  for (int qz = 0; qz < Q1D; ++qz)
924  {
925  for (int qy = 0; qy < Q1D; ++qy)
926  {
927  for (int qx = 0; qx < Q1D; ++qx)
928  {
929  for (int c = 0; c < VDIM; ++c)
930  {
931  curl[qz][qy][qx][c] = 0.0;
932  }
933  }
934  }
935  }
936 
937  // We treat x, y, z components separately for optimization specific to each.
938 
939  int osc = 0;
940 
941  {
942  // x component
943  const int D1Dz = D1D;
944  const int D1Dy = D1D;
945  const int D1Dx = D1D - 1;
946 
947  for (int dz = 0; dz < D1Dz; ++dz)
948  {
949  double gradXY[MAX_Q1D][MAX_Q1D][2];
950  for (int qy = 0; qy < Q1D; ++qy)
951  {
952  for (int qx = 0; qx < Q1D; ++qx)
953  {
954  for (int d = 0; d < 2; ++d)
955  {
956  gradXY[qy][qx][d] = 0.0;
957  }
958  }
959  }
960 
961  for (int dy = 0; dy < D1Dy; ++dy)
962  {
963  double massX[MAX_Q1D];
964  for (int qx = 0; qx < Q1D; ++qx)
965  {
966  massX[qx] = 0.0;
967  }
968 
969  for (int dx = 0; dx < D1Dx; ++dx)
970  {
971  const double t = x(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
972  for (int qx = 0; qx < Q1D; ++qx)
973  {
974  massX[qx] += t * Bo(qx,dx);
975  }
976  }
977 
978  for (int qy = 0; qy < Q1D; ++qy)
979  {
980  const double wy = Bc(qy,dy);
981  const double wDy = Gc(qy,dy);
982  for (int qx = 0; qx < Q1D; ++qx)
983  {
984  const double wx = massX[qx];
985  gradXY[qy][qx][0] += wx * wDy;
986  gradXY[qy][qx][1] += wx * wy;
987  }
988  }
989  }
990 
991  for (int qz = 0; qz < Q1D; ++qz)
992  {
993  const double wz = Bc(qz,dz);
994  const double wDz = Gc(qz,dz);
995  for (int qy = 0; qy < Q1D; ++qy)
996  {
997  for (int qx = 0; qx < Q1D; ++qx)
998  {
999  // \hat{\nabla}\times\hat{u} is [0, (u_0)_{x_2}, -(u_0)_{x_1}]
1000  curl[qz][qy][qx][1] += gradXY[qy][qx][1] * wDz; // (u_0)_{x_2}
1001  curl[qz][qy][qx][2] -= gradXY[qy][qx][0] * wz; // -(u_0)_{x_1}
1002  }
1003  }
1004  }
1005  }
1006 
1007  osc += D1Dx * D1Dy * D1Dz;
1008  }
1009 
1010  {
1011  // y component
1012  const int D1Dz = D1D;
1013  const int D1Dy = D1D - 1;
1014  const int D1Dx = D1D;
1015 
1016  for (int dz = 0; dz < D1Dz; ++dz)
1017  {
1018  double gradXY[MAX_Q1D][MAX_Q1D][2];
1019  for (int qy = 0; qy < Q1D; ++qy)
1020  {
1021  for (int qx = 0; qx < Q1D; ++qx)
1022  {
1023  for (int d = 0; d < 2; ++d)
1024  {
1025  gradXY[qy][qx][d] = 0.0;
1026  }
1027  }
1028  }
1029 
1030  for (int dx = 0; dx < D1Dx; ++dx)
1031  {
1032  double massY[MAX_Q1D];
1033  for (int qy = 0; qy < Q1D; ++qy)
1034  {
1035  massY[qy] = 0.0;
1036  }
1037 
1038  for (int dy = 0; dy < D1Dy; ++dy)
1039  {
1040  const double t = x(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
1041  for (int qy = 0; qy < Q1D; ++qy)
1042  {
1043  massY[qy] += t * Bo(qy,dy);
1044  }
1045  }
1046 
1047  for (int qx = 0; qx < Q1D; ++qx)
1048  {
1049  const double wx = Bc(qx,dx);
1050  const double wDx = Gc(qx,dx);
1051  for (int qy = 0; qy < Q1D; ++qy)
1052  {
1053  const double wy = massY[qy];
1054  gradXY[qy][qx][0] += wDx * wy;
1055  gradXY[qy][qx][1] += wx * wy;
1056  }
1057  }
1058  }
1059 
1060  for (int qz = 0; qz < Q1D; ++qz)
1061  {
1062  const double wz = Bc(qz,dz);
1063  const double wDz = Gc(qz,dz);
1064  for (int qy = 0; qy < Q1D; ++qy)
1065  {
1066  for (int qx = 0; qx < Q1D; ++qx)
1067  {
1068  // \hat{\nabla}\times\hat{u} is [-(u_1)_{x_2}, 0, (u_1)_{x_0}]
1069  curl[qz][qy][qx][0] -= gradXY[qy][qx][1] * wDz; // -(u_1)_{x_2}
1070  curl[qz][qy][qx][2] += gradXY[qy][qx][0] * wz; // (u_1)_{x_0}
1071  }
1072  }
1073  }
1074  }
1075 
1076  osc += D1Dx * D1Dy * D1Dz;
1077  }
1078 
1079  {
1080  // z component
1081  const int D1Dz = D1D - 1;
1082  const int D1Dy = D1D;
1083  const int D1Dx = D1D;
1084 
1085  for (int dx = 0; dx < D1Dx; ++dx)
1086  {
1087  double gradYZ[MAX_Q1D][MAX_Q1D][2];
1088  for (int qz = 0; qz < Q1D; ++qz)
1089  {
1090  for (int qy = 0; qy < Q1D; ++qy)
1091  {
1092  for (int d = 0; d < 2; ++d)
1093  {
1094  gradYZ[qz][qy][d] = 0.0;
1095  }
1096  }
1097  }
1098 
1099  for (int dy = 0; dy < D1Dy; ++dy)
1100  {
1101  double massZ[MAX_Q1D];
1102  for (int qz = 0; qz < Q1D; ++qz)
1103  {
1104  massZ[qz] = 0.0;
1105  }
1106 
1107  for (int dz = 0; dz < D1Dz; ++dz)
1108  {
1109  const double t = x(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
1110  for (int qz = 0; qz < Q1D; ++qz)
1111  {
1112  massZ[qz] += t * Bo(qz,dz);
1113  }
1114  }
1115 
1116  for (int qy = 0; qy < Q1D; ++qy)
1117  {
1118  const double wy = Bc(qy,dy);
1119  const double wDy = Gc(qy,dy);
1120  for (int qz = 0; qz < Q1D; ++qz)
1121  {
1122  const double wz = massZ[qz];
1123  gradYZ[qz][qy][0] += wz * wy;
1124  gradYZ[qz][qy][1] += wz * wDy;
1125  }
1126  }
1127  }
1128 
1129  for (int qx = 0; qx < Q1D; ++qx)
1130  {
1131  const double wx = Bc(qx,dx);
1132  const double wDx = Gc(qx,dx);
1133 
1134  for (int qy = 0; qy < Q1D; ++qy)
1135  {
1136  for (int qz = 0; qz < Q1D; ++qz)
1137  {
1138  // \hat{\nabla}\times\hat{u} is [(u_2)_{x_1}, -(u_2)_{x_0}, 0]
1139  curl[qz][qy][qx][0] += gradYZ[qz][qy][1] * wx; // (u_2)_{x_1}
1140  curl[qz][qy][qx][1] -= gradYZ[qz][qy][0] * wDx; // -(u_2)_{x_0}
1141  }
1142  }
1143  }
1144  }
1145  }
1146 
1147  // Apply D operator.
1148  for (int qz = 0; qz < Q1D; ++qz)
1149  {
1150  for (int qy = 0; qy < Q1D; ++qy)
1151  {
1152  for (int qx = 0; qx < Q1D; ++qx)
1153  {
1154  const double O11 = op(qx,qy,qz,0,e);
1155  const double O12 = op(qx,qy,qz,1,e);
1156  const double O13 = op(qx,qy,qz,2,e);
1157  const double O22 = op(qx,qy,qz,3,e);
1158  const double O23 = op(qx,qy,qz,4,e);
1159  const double O33 = op(qx,qy,qz,5,e);
1160 
1161  const double c1 = (O11 * curl[qz][qy][qx][0]) + (O12 * curl[qz][qy][qx][1]) +
1162  (O13 * curl[qz][qy][qx][2]);
1163  const double c2 = (O12 * curl[qz][qy][qx][0]) + (O22 * curl[qz][qy][qx][1]) +
1164  (O23 * curl[qz][qy][qx][2]);
1165  const double c3 = (O13 * curl[qz][qy][qx][0]) + (O23 * curl[qz][qy][qx][1]) +
1166  (O33 * curl[qz][qy][qx][2]);
1167 
1168  curl[qz][qy][qx][0] = c1;
1169  curl[qz][qy][qx][1] = c2;
1170  curl[qz][qy][qx][2] = c3;
1171  }
1172  }
1173  }
1174 
1175  // x component
1176  osc = 0;
1177  {
1178  const int D1Dz = D1D;
1179  const int D1Dy = D1D;
1180  const int D1Dx = D1D - 1;
1181 
1182  for (int qz = 0; qz < Q1D; ++qz)
1183  {
1184  double gradXY12[MAX_D1D][MAX_D1D];
1185  double gradXY21[MAX_D1D][MAX_D1D];
1186 
1187  for (int dy = 0; dy < D1Dy; ++dy)
1188  {
1189  for (int dx = 0; dx < D1Dx; ++dx)
1190  {
1191  gradXY12[dy][dx] = 0.0;
1192  gradXY21[dy][dx] = 0.0;
1193  }
1194  }
1195  for (int qy = 0; qy < Q1D; ++qy)
1196  {
1197  double massX[MAX_D1D][2];
1198  for (int dx = 0; dx < D1Dx; ++dx)
1199  {
1200  for (int n = 0; n < 2; ++n)
1201  {
1202  massX[dx][n] = 0.0;
1203  }
1204  }
1205  for (int qx = 0; qx < Q1D; ++qx)
1206  {
1207  for (int dx = 0; dx < D1Dx; ++dx)
1208  {
1209  const double wx = Bot(dx,qx);
1210 
1211  massX[dx][0] += wx * curl[qz][qy][qx][1];
1212  massX[dx][1] += wx * curl[qz][qy][qx][2];
1213  }
1214  }
1215  for (int dy = 0; dy < D1Dy; ++dy)
1216  {
1217  const double wy = Bct(dy,qy);
1218  const double wDy = Gct(dy,qy);
1219 
1220  for (int dx = 0; dx < D1Dx; ++dx)
1221  {
1222  gradXY21[dy][dx] += massX[dx][0] * wy;
1223  gradXY12[dy][dx] += massX[dx][1] * wDy;
1224  }
1225  }
1226  }
1227 
1228  for (int dz = 0; dz < D1Dz; ++dz)
1229  {
1230  const double wz = Bct(dz,qz);
1231  const double wDz = Gct(dz,qz);
1232  for (int dy = 0; dy < D1Dy; ++dy)
1233  {
1234  for (int dx = 0; dx < D1Dx; ++dx)
1235  {
1236  // \hat{\nabla}\times\hat{u} is [0, (u_0)_{x_2}, -(u_0)_{x_1}]
1237  // (u_0)_{x_2} * (op * curl)_1 - (u_0)_{x_1} * (op * curl)_2
1238  y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1239  e) += (gradXY21[dy][dx] * wDz) - (gradXY12[dy][dx] * wz);
1240  }
1241  }
1242  }
1243  } // loop qz
1244 
1245  osc += D1Dx * D1Dy * D1Dz;
1246  }
1247 
1248  // y component
1249  {
1250  const int D1Dz = D1D;
1251  const int D1Dy = D1D - 1;
1252  const int D1Dx = D1D;
1253 
1254  for (int qz = 0; qz < Q1D; ++qz)
1255  {
1256  double gradXY02[MAX_D1D][MAX_D1D];
1257  double gradXY20[MAX_D1D][MAX_D1D];
1258 
1259  for (int dy = 0; dy < D1Dy; ++dy)
1260  {
1261  for (int dx = 0; dx < D1Dx; ++dx)
1262  {
1263  gradXY02[dy][dx] = 0.0;
1264  gradXY20[dy][dx] = 0.0;
1265  }
1266  }
1267  for (int qx = 0; qx < Q1D; ++qx)
1268  {
1269  double massY[MAX_D1D][2];
1270  for (int dy = 0; dy < D1Dy; ++dy)
1271  {
1272  massY[dy][0] = 0.0;
1273  massY[dy][1] = 0.0;
1274  }
1275  for (int qy = 0; qy < Q1D; ++qy)
1276  {
1277  for (int dy = 0; dy < D1Dy; ++dy)
1278  {
1279  const double wy = Bot(dy,qy);
1280 
1281  massY[dy][0] += wy * curl[qz][qy][qx][2];
1282  massY[dy][1] += wy * curl[qz][qy][qx][0];
1283  }
1284  }
1285  for (int dx = 0; dx < D1Dx; ++dx)
1286  {
1287  const double wx = Bct(dx,qx);
1288  const double wDx = Gct(dx,qx);
1289 
1290  for (int dy = 0; dy < D1Dy; ++dy)
1291  {
1292  gradXY02[dy][dx] += massY[dy][0] * wDx;
1293  gradXY20[dy][dx] += massY[dy][1] * wx;
1294  }
1295  }
1296  }
1297 
1298  for (int dz = 0; dz < D1Dz; ++dz)
1299  {
1300  const double wz = Bct(dz,qz);
1301  const double wDz = Gct(dz,qz);
1302  for (int dy = 0; dy < D1Dy; ++dy)
1303  {
1304  for (int dx = 0; dx < D1Dx; ++dx)
1305  {
1306  // \hat{\nabla}\times\hat{u} is [-(u_1)_{x_2}, 0, (u_1)_{x_0}]
1307  // -(u_1)_{x_2} * (op * curl)_0 + (u_1)_{x_0} * (op * curl)_2
1308  y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1309  e) += (-gradXY20[dy][dx] * wDz) + (gradXY02[dy][dx] * wz);
1310  }
1311  }
1312  }
1313  } // loop qz
1314 
1315  osc += D1Dx * D1Dy * D1Dz;
1316  }
1317 
1318  // z component
1319  {
1320  const int D1Dz = D1D - 1;
1321  const int D1Dy = D1D;
1322  const int D1Dx = D1D;
1323 
1324  for (int qx = 0; qx < Q1D; ++qx)
1325  {
1326  double gradYZ01[MAX_D1D][MAX_D1D];
1327  double gradYZ10[MAX_D1D][MAX_D1D];
1328 
1329  for (int dy = 0; dy < D1Dy; ++dy)
1330  {
1331  for (int dz = 0; dz < D1Dz; ++dz)
1332  {
1333  gradYZ01[dz][dy] = 0.0;
1334  gradYZ10[dz][dy] = 0.0;
1335  }
1336  }
1337  for (int qy = 0; qy < Q1D; ++qy)
1338  {
1339  double massZ[MAX_D1D][2];
1340  for (int dz = 0; dz < D1Dz; ++dz)
1341  {
1342  for (int n = 0; n < 2; ++n)
1343  {
1344  massZ[dz][n] = 0.0;
1345  }
1346  }
1347  for (int qz = 0; qz < Q1D; ++qz)
1348  {
1349  for (int dz = 0; dz < D1Dz; ++dz)
1350  {
1351  const double wz = Bot(dz,qz);
1352 
1353  massZ[dz][0] += wz * curl[qz][qy][qx][0];
1354  massZ[dz][1] += wz * curl[qz][qy][qx][1];
1355  }
1356  }
1357  for (int dy = 0; dy < D1Dy; ++dy)
1358  {
1359  const double wy = Bct(dy,qy);
1360  const double wDy = Gct(dy,qy);
1361 
1362  for (int dz = 0; dz < D1Dz; ++dz)
1363  {
1364  gradYZ01[dz][dy] += wy * massZ[dz][1];
1365  gradYZ10[dz][dy] += wDy * massZ[dz][0];
1366  }
1367  }
1368  }
1369 
1370  for (int dx = 0; dx < D1Dx; ++dx)
1371  {
1372  const double wx = Bct(dx,qx);
1373  const double wDx = Gct(dx,qx);
1374 
1375  for (int dy = 0; dy < D1Dy; ++dy)
1376  {
1377  for (int dz = 0; dz < D1Dz; ++dz)
1378  {
1379  // \hat{\nabla}\times\hat{u} is [(u_2)_{x_1}, -(u_2)_{x_0}, 0]
1380  // (u_2)_{x_1} * (op * curl)_0 - (u_2)_{x_0} * (op * curl)_1
1381  y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1382  e) += (gradYZ10[dz][dy] * wx) - (gradYZ01[dz][dy] * wDx);
1383  }
1384  }
1385  }
1386  } // loop qx
1387  }
1388 
1389  }); // end of element loop
1390 }
1391 
1392 void CurlCurlIntegrator::AddMultPA(const Vector &x, Vector &y) const
1393 {
1394  if (dim == 3)
1395  {
1396  PACurlCurlApply3D(dofs1D, quad1D, ne, mapsO->B, mapsC->B, mapsO->Bt,
1397  mapsC->Bt, mapsC->G, mapsC->Gt, pa_data, x, y);
1398  }
1399  else if (dim == 2)
1400  {
1401  PACurlCurlApply2D(dofs1D, quad1D, ne, mapsO->B, mapsO->Bt,
1402  mapsC->G, mapsC->Gt, pa_data, x, y);
1403  }
1404  else
1405  {
1406  MFEM_ABORT("Unsupported dimension!");
1407  }
1408 }
1409 
1410 static void PACurlCurlAssembleDiagonal2D(const int D1D,
1411  const int Q1D,
1412  const int NE,
1413  const Array<double> &_Bo,
1414  const Array<double> &_Gc,
1415  const Vector &_op,
1416  Vector &_diag)
1417 {
1418  constexpr static int VDIM = 2;
1419 
1420  auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
1421  auto Gc = Reshape(_Gc.Read(), Q1D, D1D);
1422  auto op = Reshape(_op.Read(), Q1D, Q1D, NE);
1423  auto diag = Reshape(_diag.ReadWrite(), 2*(D1D-1)*D1D, NE);
1424 
1425  MFEM_FORALL(e, NE,
1426  {
1427  int osc = 0;
1428 
1429  for (int c = 0; c < VDIM; ++c) // loop over x, y components
1430  {
1431  const int D1Dy = (c == 1) ? D1D - 1 : D1D;
1432  const int D1Dx = (c == 0) ? D1D - 1 : D1D;
1433 
1434  double t[MAX_Q1D];
1435 
1436  for (int dy = 0; dy < D1Dy; ++dy)
1437  {
1438  for (int qx = 0; qx < Q1D; ++qx)
1439  {
1440  t[qx] = 0.0;
1441  for (int qy = 0; qy < Q1D; ++qy)
1442  {
1443  const double wy = (c == 1) ? Bo(qy,dy) : -Gc(qy,dy);
1444  t[qx] += wy * wy * op(qx,qy,e);
1445  }
1446  }
1447 
1448  for (int dx = 0; dx < D1Dx; ++dx)
1449  {
1450  for (int qx = 0; qx < Q1D; ++qx)
1451  {
1452  const double wx = ((c == 0) ? Bo(qx,dx) : Gc(qx,dx));
1453  diag(dx + (dy * D1Dx) + osc, e) += t[qx] * wx * wx;
1454  }
1455  }
1456  }
1457 
1458  osc += D1Dx * D1Dy;
1459  } // loop c
1460  }); // end of element loop
1461 }
1462 
1463 template<int MAX_D1D = HCURL_MAX_D1D, int MAX_Q1D = HCURL_MAX_Q1D>
1464 static void PACurlCurlAssembleDiagonal3D(const int D1D,
1465  const int Q1D,
1466  const int NE,
1467  const Array<double> &_Bo,
1468  const Array<double> &_Bc,
1469  const Array<double> &_Go,
1470  const Array<double> &_Gc,
1471  const Vector &_op,
1472  Vector &_diag)
1473 {
1474  constexpr static int VDIM = 3;
1475  MFEM_VERIFY(D1D <= MAX_D1D, "Error: D1D > MAX_D1D");
1476  MFEM_VERIFY(Q1D <= MAX_Q1D, "Error: Q1D > MAX_Q1D");
1477 
1478  auto Bo = Reshape(_Bo.Read(), Q1D, D1D-1);
1479  auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
1480  auto Go = Reshape(_Go.Read(), Q1D, D1D-1);
1481  auto Gc = Reshape(_Gc.Read(), Q1D, D1D);
1482  auto op = Reshape(_op.Read(), Q1D, Q1D, Q1D, 6, NE);
1483  auto diag = Reshape(_diag.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
1484 
1485  MFEM_FORALL(e, NE,
1486  {
1487  // Using (\nabla\times u) F = 1/det(dF) dF \hat{\nabla}\times\hat{u} (p. 78 of Monk), we get
1488  // (\nabla\times u) \cdot (\nabla\times u) = 1/det(dF)^2 \hat{\nabla}\times\hat{u}^T dF^T dF \hat{\nabla}\times\hat{u}
1489  // If c = 0, \hat{\nabla}\times\hat{u} reduces to [0, (u_0)_{x_2}, -(u_0)_{x_1}]
1490  // If c = 1, \hat{\nabla}\times\hat{u} reduces to [-(u_1)_{x_2}, 0, (u_1)_{x_0}]
1491  // If c = 2, \hat{\nabla}\times\hat{u} reduces to [(u_2)_{x_1}, -(u_2)_{x_0}, 0]
1492 
1493  // For each c, we will keep 6 arrays for derivatives multiplied by the 6 entries of the symmetric 3x3 matrix (dF^T dF).
1494 
1495  int osc = 0;
1496 
1497  for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
1498  {
1499  const int D1Dz = (c == 2) ? D1D - 1 : D1D;
1500  const int D1Dy = (c == 1) ? D1D - 1 : D1D;
1501  const int D1Dx = (c == 0) ? D1D - 1 : D1D;
1502 
1503  double zt[MAX_Q1D][MAX_Q1D][MAX_D1D][6][3];
1504 
1505  // z contraction
1506  for (int qx = 0; qx < Q1D; ++qx)
1507  {
1508  for (int qy = 0; qy < Q1D; ++qy)
1509  {
1510  for (int dz = 0; dz < D1Dz; ++dz)
1511  {
1512  for (int i=0; i<6; ++i)
1513  {
1514  for (int d=0; d<3; ++d)
1515  {
1516  zt[qx][qy][dz][i][d] = 0.0;
1517  }
1518  }
1519 
1520  for (int qz = 0; qz < Q1D; ++qz)
1521  {
1522  const double wz = ((c == 2) ? Bo(qz,dz) : Bc(qz,dz));
1523  const double wDz = ((c == 2) ? Go(qz,dz) : Gc(qz,dz));
1524 
1525  for (int i=0; i<6; ++i)
1526  {
1527  zt[qx][qy][dz][i][0] += wz * wz * op(qx,qy,qz,i,e);
1528  zt[qx][qy][dz][i][1] += wDz * wz * op(qx,qy,qz,i,e);
1529  zt[qx][qy][dz][i][2] += wDz * wDz * op(qx,qy,qz,i,e);
1530  }
1531  }
1532  }
1533  }
1534  } // end of z contraction
1535 
1536  double yt[MAX_Q1D][MAX_D1D][MAX_D1D][6][3][3];
1537 
1538  // y contraction
1539  for (int qx = 0; qx < Q1D; ++qx)
1540  {
1541  for (int dz = 0; dz < D1Dz; ++dz)
1542  {
1543  for (int dy = 0; dy < D1Dy; ++dy)
1544  {
1545  for (int i=0; i<6; ++i)
1546  {
1547  for (int d=0; d<3; ++d)
1548  for (int j=0; j<3; ++j)
1549  {
1550  yt[qx][dy][dz][i][d][j] = 0.0;
1551  }
1552  }
1553 
1554  for (int qy = 0; qy < Q1D; ++qy)
1555  {
1556  const double wy = ((c == 1) ? Bo(qy,dy) : Bc(qy,dy));
1557  const double wDy = ((c == 1) ? Go(qy,dy) : Gc(qy,dy));
1558 
1559  for (int i=0; i<6; ++i)
1560  {
1561  for (int d=0; d<3; ++d)
1562  {
1563  yt[qx][dy][dz][i][d][0] += wy * wy * zt[qx][qy][dz][i][d];
1564  yt[qx][dy][dz][i][d][1] += wDy * wy * zt[qx][qy][dz][i][d];
1565  yt[qx][dy][dz][i][d][2] += wDy * wDy * zt[qx][qy][dz][i][d];
1566  }
1567  }
1568  }
1569  }
1570  }
1571  } // end of y contraction
1572 
1573  // x contraction
1574  for (int dz = 0; dz < D1Dz; ++dz)
1575  {
1576  for (int dy = 0; dy < D1Dy; ++dy)
1577  {
1578  for (int dx = 0; dx < D1Dx; ++dx)
1579  {
1580  for (int qx = 0; qx < Q1D; ++qx)
1581  {
1582  const double wx = ((c == 0) ? Bo(qx,dx) : Bc(qx,dx));
1583  const double wDx = ((c == 0) ? Go(qx,dx) : Gc(qx,dx));
1584 
1585  // Using (\nabla\times u) F = 1/det(dF) dF \hat{\nabla}\times\hat{u} (p. 78 of Monk), we get
1586  // (\nabla\times u) \cdot (\nabla\times u) = 1/det(dF)^2 \hat{\nabla}\times\hat{u}^T dF^T dF \hat{\nabla}\times\hat{u}
1587  // If c = 0, \hat{\nabla}\times\hat{u} reduces to [0, (u_0)_{x_2}, -(u_0)_{x_1}]
1588  // If c = 1, \hat{\nabla}\times\hat{u} reduces to [-(u_1)_{x_2}, 0, (u_1)_{x_0}]
1589  // If c = 2, \hat{\nabla}\times\hat{u} reduces to [(u_2)_{x_1}, -(u_2)_{x_0}, 0]
1590 
1591  /*
1592  const double O11 = op(q,0,e);
1593  const double O12 = op(q,1,e);
1594  const double O13 = op(q,2,e);
1595  const double O22 = op(q,3,e);
1596  const double O23 = op(q,4,e);
1597  const double O33 = op(q,5,e);
1598  */
1599 
1600  if (c == 0)
1601  {
1602  // (u_0)_{x_2} (O22 (u_0)_{x_2} - O23 (u_0)_{x_1}) - (u_0)_{x_1} (O32 (u_0)_{x_2} - O33 (u_0)_{x_1})
1603 
1604  // (u_0)_{x_2} O22 (u_0)_{x_2}
1605  diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1606  e) += yt[qx][dy][dz][3][2][0] * wx * wx;
1607 
1608  // -(u_0)_{x_2} O23 (u_0)_{x_1} - (u_0)_{x_1} O32 (u_0)_{x_2}
1609  diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1610  e) += -2.0 * yt[qx][dy][dz][4][1][1] * wx * wx;
1611 
1612  // (u_0)_{x_1} O33 (u_0)_{x_1}
1613  diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1614  e) += yt[qx][dy][dz][5][0][2] * wx * wx;
1615  }
1616  else if (c == 1)
1617  {
1618  // (u_1)_{x_2} (O11 (u_1)_{x_2} - O13 (u_1)_{x_0}) + (u_1)_{x_0} (-O31 (u_1)_{x_2} + O33 (u_1)_{x_0})
1619 
1620  // (u_1)_{x_2} O11 (u_1)_{x_2}
1621  diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1622  e) += yt[qx][dy][dz][0][2][0] * wx * wx;
1623 
1624  // -(u_1)_{x_2} O13 (u_1)_{x_0} - (u_1)_{x_0} O31 (u_1)_{x_2}
1625  diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1626  e) += -2.0 * yt[qx][dy][dz][2][1][0] * wDx * wx;
1627 
1628  // (u_1)_{x_0} O33 (u_1)_{x_0})
1629  diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1630  e) += yt[qx][dy][dz][5][0][0] * wDx * wDx;
1631  }
1632  else
1633  {
1634  // (u_2)_{x_1} (O11 (u_2)_{x_1} - O12 (u_2)_{x_0}) - (u_2)_{x_0} (O21 (u_2)_{x_1} - O22 (u_2)_{x_0})
1635 
1636  // (u_2)_{x_1} O11 (u_2)_{x_1}
1637  diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1638  e) += yt[qx][dy][dz][0][0][2] * wx * wx;
1639 
1640  // -(u_2)_{x_1} O12 (u_2)_{x_0} - (u_2)_{x_0} O21 (u_2)_{x_1}
1641  diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1642  e) += -2.0 * yt[qx][dy][dz][1][0][1] * wDx * wx;
1643 
1644  // (u_2)_{x_0} O22 (u_2)_{x_0}
1645  diag(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc,
1646  e) += yt[qx][dy][dz][3][0][0] * wDx * wDx;
1647  }
1648  }
1649  }
1650  }
1651  } // end of x contraction
1652 
1653  osc += D1Dx * D1Dy * D1Dz;
1654  } // loop c
1655  }); // end of element loop
1656 }
1657 
1658 void CurlCurlIntegrator::AssembleDiagonalPA(Vector& diag)
1659 {
1660  if (dim == 3)
1661  {
1662  // Reduce HCURL_MAX_D1D/Q1D to avoid using too much memory
1663  constexpr int MAX_D1D = 4;
1664  constexpr int MAX_Q1D = 5;
1665  PACurlCurlAssembleDiagonal3D<MAX_D1D,MAX_Q1D>(dofs1D, quad1D, ne,
1666  mapsO->B, mapsC->B,
1667  mapsO->G, mapsC->G,
1668  pa_data, diag);
1669  }
1670  else if (dim == 2)
1671  {
1672  PACurlCurlAssembleDiagonal2D(dofs1D, quad1D, ne,
1673  mapsO->B, mapsC->G, pa_data, diag);
1674  }
1675  else
1676  {
1677  MFEM_ABORT("Unsupported dimension!");
1678  }
1679 }
1680 
1681 void MixedVectorGradientIntegrator::AssemblePA(const FiniteElementSpace
1682  &trial_fes,
1683  const FiniteElementSpace &test_fes)
1684 {
1685  // Assumes tensor-product elements, with a vector test space and H^1 trial space.
1686  Mesh *mesh = trial_fes.GetMesh();
1687  const FiniteElement *trial_fel = trial_fes.GetFE(0);
1688  const FiniteElement *test_fel = test_fes.GetFE(0);
1689 
1690  const NodalTensorFiniteElement *trial_el =
1691  dynamic_cast<const NodalTensorFiniteElement*>(trial_fel);
1692  MFEM_VERIFY(trial_el != NULL, "Only NodalTensorFiniteElement is supported!");
1693 
1694  const VectorTensorFiniteElement *test_el =
1695  dynamic_cast<const VectorTensorFiniteElement*>(test_fel);
1696  MFEM_VERIFY(test_el != NULL, "Only VectorTensorFiniteElement is supported!");
1697 
1698  const IntegrationRule *ir
1699  = IntRule ? IntRule : &MassIntegrator::GetRule(*trial_el, *trial_el,
1700  *mesh->GetElementTransformation(0));
1701  const int dims = trial_el->GetDim();
1702  MFEM_VERIFY(dims == 2 || dims == 3, "");
1703 
1704  const int symmDims = (dims * (dims + 1)) / 2; // 1x1: 1, 2x2: 3, 3x3: 6
1705  const int nq = ir->GetNPoints();
1706  dim = mesh->Dimension();
1707  MFEM_VERIFY(dim == 2 || dim == 3, "");
1708 
1709  MFEM_VERIFY(trial_el->GetOrder() == test_el->GetOrder(), "");
1710 
1711  ne = trial_fes.GetNE();
1712  geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS);
1713  mapsC = &test_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
1714  mapsO = &test_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
1715  dofs1D = mapsC->ndof;
1716  quad1D = mapsC->nqpt;
1717 
1718  MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
1719 
1720  pa_data.SetSize(symmDims * nq * ne, Device::GetMemoryType());
1721 
1722  Vector coeff(ne * nq);
1723  coeff = 1.0;
1724  if (Q)
1725  {
1726  for (int e=0; e<ne; ++e)
1727  {
1728  ElementTransformation *tr = mesh->GetElementTransformation(e);
1729  for (int p=0; p<nq; ++p)
1730  {
1731  coeff[p + (e * nq)] = Q->Eval(*tr, ir->IntPoint(p));
1732  }
1733  }
1734  }
1735 
1736  // Use the same setup functions as VectorFEMassIntegrator.
1737  if (test_el->GetDerivType() == mfem::FiniteElement::CURL && dim == 3)
1738  {
1739  PAHcurlSetup3D(quad1D, ne, ir->GetWeights(), geom->J,
1740  coeff, pa_data);
1741  }
1742  else if (test_el->GetDerivType() == mfem::FiniteElement::CURL && dim == 2)
1743  {
1744  PAHcurlSetup2D(quad1D, ne, ir->GetWeights(), geom->J,
1745  coeff, pa_data);
1746  }
1747  else
1748  {
1749  MFEM_ABORT("Unknown kernel.");
1750  }
1751 }
1752 
1753 // Apply to x corresponding to DOF's in H^1 (trial), whose gradients are integrated
1754 // against H(curl) test functions corresponding to y.
1755 template<int MAX_D1D = HCURL_MAX_D1D, int MAX_Q1D = HCURL_MAX_Q1D>
1756 static void PAHcurlH1Apply3D(const int D1D,
1757  const int Q1D,
1758  const int NE,
1759  const Array<double> &_Bc,
1760  const Array<double> &_Gc,
1761  const Array<double> &_Bot,
1762  const Array<double> &_Bct,
1763  const Vector &_op,
1764  const Vector &_x,
1765  Vector &_y)
1766 {
1767  constexpr static int VDIM = 3;
1768 
1769  auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
1770  auto Gc = Reshape(_Gc.Read(), Q1D, D1D);
1771  auto Bot = Reshape(_Bot.Read(), D1D-1, Q1D);
1772  auto Bct = Reshape(_Bct.Read(), D1D, Q1D);
1773  auto op = Reshape(_op.Read(), Q1D, Q1D, Q1D, 6, NE);
1774  auto x = Reshape(_x.Read(), D1D, D1D, D1D, NE);
1775  auto y = Reshape(_y.ReadWrite(), 3*(D1D-1)*D1D*D1D, NE);
1776 
1777  MFEM_FORALL(e, NE,
1778  {
1779  double mass[MAX_Q1D][MAX_Q1D][MAX_Q1D][VDIM];
1780 
1781  for (int qz = 0; qz < Q1D; ++qz)
1782  {
1783  for (int qy = 0; qy < Q1D; ++qy)
1784  {
1785  for (int qx = 0; qx < Q1D; ++qx)
1786  {
1787  for (int c = 0; c < VDIM; ++c)
1788  {
1789  mass[qz][qy][qx][c] = 0.0;
1790  }
1791  }
1792  }
1793  }
1794 
1795  for (int dz = 0; dz < D1D; ++dz)
1796  {
1797  double gradXY[MAX_Q1D][MAX_Q1D][3];
1798  for (int qy = 0; qy < Q1D; ++qy)
1799  {
1800  for (int qx = 0; qx < Q1D; ++qx)
1801  {
1802  gradXY[qy][qx][0] = 0.0;
1803  gradXY[qy][qx][1] = 0.0;
1804  gradXY[qy][qx][2] = 0.0;
1805  }
1806  }
1807  for (int dy = 0; dy < D1D; ++dy)
1808  {
1809  double gradX[MAX_Q1D][2];
1810  for (int qx = 0; qx < Q1D; ++qx)
1811  {
1812  gradX[qx][0] = 0.0;
1813  gradX[qx][1] = 0.0;
1814  }
1815  for (int dx = 0; dx < D1D; ++dx)
1816  {
1817  const double s = x(dx,dy,dz,e);
1818  for (int qx = 0; qx < Q1D; ++qx)
1819  {
1820  gradX[qx][0] += s * Bc(qx,dx);
1821  gradX[qx][1] += s * Gc(qx,dx);
1822  }
1823  }
1824  for (int qy = 0; qy < Q1D; ++qy)
1825  {
1826  const double wy = Bc(qy,dy);
1827  const double wDy = Gc(qy,dy);
1828  for (int qx = 0; qx < Q1D; ++qx)
1829  {
1830  const double wx = gradX[qx][0];
1831  const double wDx = gradX[qx][1];
1832  gradXY[qy][qx][0] += wDx * wy;
1833  gradXY[qy][qx][1] += wx * wDy;
1834  gradXY[qy][qx][2] += wx * wy;
1835  }
1836  }
1837  }
1838  for (int qz = 0; qz < Q1D; ++qz)
1839  {
1840  const double wz = Bc(qz,dz);
1841  const double wDz = Gc(qz,dz);
1842  for (int qy = 0; qy < Q1D; ++qy)
1843  {
1844  for (int qx = 0; qx < Q1D; ++qx)
1845  {
1846  mass[qz][qy][qx][0] += gradXY[qy][qx][0] * wz;
1847  mass[qz][qy][qx][1] += gradXY[qy][qx][1] * wz;
1848  mass[qz][qy][qx][2] += gradXY[qy][qx][2] * wDz;
1849  }
1850  }
1851  }
1852  }
1853 
1854  // Apply D operator.
1855  for (int qz = 0; qz < Q1D; ++qz)
1856  {
1857  for (int qy = 0; qy < Q1D; ++qy)
1858  {
1859  for (int qx = 0; qx < Q1D; ++qx)
1860  {
1861  const double O11 = op(qx,qy,qz,0,e);
1862  const double O12 = op(qx,qy,qz,1,e);
1863  const double O13 = op(qx,qy,qz,2,e);
1864  const double O22 = op(qx,qy,qz,3,e);
1865  const double O23 = op(qx,qy,qz,4,e);
1866  const double O33 = op(qx,qy,qz,5,e);
1867  const double massX = mass[qz][qy][qx][0];
1868  const double massY = mass[qz][qy][qx][1];
1869  const double massZ = mass[qz][qy][qx][2];
1870  mass[qz][qy][qx][0] = (O11*massX)+(O12*massY)+(O13*massZ);
1871  mass[qz][qy][qx][1] = (O12*massX)+(O22*massY)+(O23*massZ);
1872  mass[qz][qy][qx][2] = (O13*massX)+(O23*massY)+(O33*massZ);
1873  }
1874  }
1875  }
1876 
1877  for (int qz = 0; qz < Q1D; ++qz)
1878  {
1879  double massXY[MAX_D1D][MAX_D1D];
1880 
1881  int osc = 0;
1882 
1883  for (int c = 0; c < VDIM; ++c) // loop over x, y, z components
1884  {
1885  const int D1Dz = (c == 2) ? D1D - 1 : D1D;
1886  const int D1Dy = (c == 1) ? D1D - 1 : D1D;
1887  const int D1Dx = (c == 0) ? D1D - 1 : D1D;
1888 
1889  for (int dy = 0; dy < D1Dy; ++dy)
1890  {
1891  for (int dx = 0; dx < D1Dx; ++dx)
1892  {
1893  massXY[dy][dx] = 0;
1894  }
1895  }
1896  for (int qy = 0; qy < Q1D; ++qy)
1897  {
1898  double massX[MAX_D1D];
1899  for (int dx = 0; dx < D1Dx; ++dx)
1900  {
1901  massX[dx] = 0;
1902  }
1903  for (int qx = 0; qx < Q1D; ++qx)
1904  {
1905  for (int dx = 0; dx < D1Dx; ++dx)
1906  {
1907  massX[dx] += mass[qz][qy][qx][c] * ((c == 0) ? Bot(dx,qx) : Bct(dx,qx));
1908  }
1909  }
1910  for (int dy = 0; dy < D1Dy; ++dy)
1911  {
1912  const double wy = (c == 1) ? Bot(dy,qy) : Bct(dy,qy);
1913  for (int dx = 0; dx < D1Dx; ++dx)
1914  {
1915  massXY[dy][dx] += massX[dx] * wy;
1916  }
1917  }
1918  }
1919 
1920  for (int dz = 0; dz < D1Dz; ++dz)
1921  {
1922  const double wz = (c == 2) ? Bot(dz,qz) : Bct(dz,qz);
1923  for (int dy = 0; dy < D1Dy; ++dy)
1924  {
1925  for (int dx = 0; dx < D1Dx; ++dx)
1926  {
1927  y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) += massXY[dy][dx] * wz;
1928  }
1929  }
1930  }
1931 
1932  osc += D1Dx * D1Dy * D1Dz;
1933  } // loop c
1934  } // loop qz
1935  }); // end of element loop
1936 }
1937 
1938 // Apply to x corresponding to DOF's in H^1 (trial), whose gradients are integrated
1939 // against H(curl) test functions corresponding to y.
1940 static void PAHcurlH1Apply2D(const int D1D,
1941  const int Q1D,
1942  const int NE,
1943  const Array<double> &_Bc,
1944  const Array<double> &_Gc,
1945  const Array<double> &_Bot,
1946  const Array<double> &_Bct,
1947  const Vector &_op,
1948  const Vector &_x,
1949  Vector &_y)
1950 {
1951  constexpr static int VDIM = 2;
1952 
1953  auto Bc = Reshape(_Bc.Read(), Q1D, D1D);
1954  auto Gc = Reshape(_Gc.Read(), Q1D, D1D);
1955  auto Bot = Reshape(_Bot.Read(), D1D-1, Q1D);
1956  auto Bct = Reshape(_Bct.Read(), D1D, Q1D);
1957  auto op = Reshape(_op.Read(), Q1D, Q1D, 3, NE);
1958  auto x = Reshape(_x.Read(), D1D, D1D, NE);
1959  auto y = Reshape(_y.ReadWrite(), 2*(D1D-1)*D1D, NE);
1960 
1961  MFEM_FORALL(e, NE,
1962  {
1963  double mass[MAX_Q1D][MAX_Q1D][VDIM];
1964 
1965  for (int qy = 0; qy < Q1D; ++qy)
1966  {
1967  for (int qx = 0; qx < Q1D; ++qx)
1968  {
1969  for (int c = 0; c < VDIM; ++c)
1970  {
1971  mass[qy][qx][c] = 0.0;
1972  }
1973  }
1974  }
1975 
1976  for (int dy = 0; dy < D1D; ++dy)
1977  {
1978  double gradX[MAX_Q1D][2];
1979  for (int qx = 0; qx < Q1D; ++qx)
1980  {
1981  gradX[qx][0] = 0.0;
1982  gradX[qx][1] = 0.0;
1983  }
1984  for (int dx = 0; dx < D1D; ++dx)
1985  {
1986  const double s = x(dx,dy,e);
1987  for (int qx = 0; qx < Q1D; ++qx)
1988  {
1989  gradX[qx][0] += s * Bc(qx,dx);
1990  gradX[qx][1] += s * Gc(qx,dx);
1991  }
1992  }
1993  for (int qy = 0; qy < Q1D; ++qy)
1994  {
1995  const double wy = Bc(qy,dy);
1996  const double wDy = Gc(qy,dy);
1997  for (int qx = 0; qx < Q1D; ++qx)
1998  {
1999  const double wx = gradX[qx][0];
2000  const double wDx = gradX[qx][1];
2001  mass[qy][qx][0] += wDx * wy;
2002  mass[qy][qx][1] += wx * wDy;
2003  }
2004  }
2005  }
2006 
2007  // Apply D operator.
2008  for (int qy = 0; qy < Q1D; ++qy)
2009  {
2010  for (int qx = 0; qx < Q1D; ++qx)
2011  {
2012  const double O11 = op(qx,qy,0,e);
2013  const double O12 = op(qx,qy,1,e);
2014  const double O22 = op(qx,qy,2,e);
2015  const double massX = mass[qy][qx][0];
2016  const double massY = mass[qy][qx][1];
2017  mass[qy][qx][0] = (O11*massX)+(O12*massY);
2018  mass[qy][qx][1] = (O12*massX)+(O22*massY);
2019  }
2020  }
2021 
2022  for (int qy = 0; qy < Q1D; ++qy)
2023  {
2024  int osc = 0;
2025 
2026  for (int c = 0; c < VDIM; ++c) // loop over x, y components
2027  {
2028  const int D1Dy = (c == 1) ? D1D - 1 : D1D;
2029  const int D1Dx = (c == 0) ? D1D - 1 : D1D;
2030 
2031  double massX[MAX_D1D];
2032  for (int dx = 0; dx < D1Dx; ++dx)
2033  {
2034  massX[dx] = 0;
2035  }
2036  for (int qx = 0; qx < Q1D; ++qx)
2037  {
2038  for (int dx = 0; dx < D1Dx; ++dx)
2039  {
2040  massX[dx] += mass[qy][qx][c] * ((c == 0) ? Bot(dx,qx) : Bct(dx,qx));
2041  }
2042  }
2043 
2044  for (int dy = 0; dy < D1Dy; ++dy)
2045  {
2046  const double wy = (c == 1) ? Bot(dy,qy) : Bct(dy,qy);
2047 
2048  for (int dx = 0; dx < D1Dx; ++dx)
2049  {
2050  y(dx + (dy * D1Dx) + osc, e) += massX[dx] * wy;
2051  }
2052  }
2053 
2054  osc += D1Dx * D1Dy;
2055  } // loop c
2056  }
2057  }); // end of element loop
2058 }
2059 
2060 void MixedVectorGradientIntegrator::AddMultPA(const Vector &x, Vector &y) const
2061 {
2062  if (dim == 3)
2063  PAHcurlH1Apply3D(dofs1D, quad1D, ne, mapsC->B, mapsC->G,
2064  mapsO->Bt, mapsC->Bt, pa_data, x, y);
2065  else if (dim == 2)
2066  PAHcurlH1Apply2D(dofs1D, quad1D, ne, mapsC->B, mapsC->G,
2067  mapsO->Bt, mapsC->Bt, pa_data, x, y);
2068  else
2069  {
2070  MFEM_ABORT("Unsupported dimension!");
2071  }
2072 }
2073 
2074 } // namespace mfem
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:245
Abstract class for Finite Elements.
Definition: fe.hpp:232
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe.hpp:308
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const DofToQuad & GetDofToQuad(const IntegrationRule &ir, DofToQuad::Mode mode) const
Definition: fe.cpp:11515
const Geometry::Type geom
Definition: ex1.cpp:40
int ndof
Number of degrees of freedom = number of basis functions. When mode is TENSOR, this is the 1D number...
Definition: fe.hpp:162
const Array< double > & GetWeights() const
Return the quadrature weights in a contiguous array.
Definition: intrules.cpp:85
double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:341
DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition: dtensor.hpp:134
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:408
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:248
Implements CalcCurlShape methods.
Definition: fe.hpp:294
const GeometricFactors * GetGeometricFactors(const IntegrationRule &ir, const int flags)
Return the mesh geometric factors corresponding to the given integration rule.
Definition: mesh.cpp:763
double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:349
const int MAX_Q1D
Definition: forall.hpp:36
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:295
constexpr int HCURL_MAX_D1D
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:273
int Dimension() const
Definition: mesh.hpp:744
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
const DofToQuad & GetDofToQuadOpen(const IntegrationRule &ir, DofToQuad::Mode mode) const
Definition: fe.cpp:11524
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:333
void GetElementTransformation(int i, IsoparametricTransformation *ElTr)
Definition: mesh.cpp:338
int dim
Definition: ex24.cpp:43
MemoryType GetMemoryType(MemoryClass mc)
Return a suitable MemoryType for a given MemoryClass.
Definition: mem_manager.cpp:49
const int MAX_D1D
Definition: forall.hpp:35
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1642
Vector data type.
Definition: vector.hpp:48
int GetDerivType() const
Definition: fe.hpp:336
constexpr int HCURL_MAX_Q1D