MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
bilininteg_vectorfe.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 
15 namespace mfem
16 {
17 
18 void PADiffusionSetup3D(const int Q1D,
19  const int coeffDim,
20  const int NE,
21  const Array<double> &w,
22  const Vector &j,
23  const Vector &coeff_,
24  Vector &op);
25 
26 void PAHcurlMassAssembleDiagonal2D(const int D1D,
27  const int Q1D,
28  const int NE,
29  const bool symmetric,
30  const Array<double> &bo,
31  const Array<double> &bc,
32  const Vector &pa_data,
33  Vector &diag);
34 
35 void PAHcurlMassAssembleDiagonal3D(const int D1D,
36  const int Q1D,
37  const int NE,
38  const bool symmetric,
39  const Array<double> &bo,
40  const Array<double> &bc,
41  const Vector &pa_data,
42  Vector &diag);
43 
44 template<int T_D1D = 0, int T_Q1D = 0>
45 void SmemPAHcurlMassAssembleDiagonal3D(const int D1D,
46  const int Q1D,
47  const int NE,
48  const bool symmetric,
49  const Array<double> &bo,
50  const Array<double> &bc,
51  const Vector &pa_data,
52  Vector &diag);
53 
54 void PAHcurlMassApply2D(const int D1D,
55  const int Q1D,
56  const int NE,
57  const bool symmetric,
58  const Array<double> &bo,
59  const Array<double> &bc,
60  const Array<double> &bot,
61  const Array<double> &bct,
62  const Vector &pa_data,
63  const Vector &x,
64  Vector &y);
65 
66 void PAHcurlMassApply3D(const int D1D,
67  const int Q1D,
68  const int NE,
69  const bool symmetric,
70  const Array<double> &bo,
71  const Array<double> &bc,
72  const Array<double> &bot,
73  const Array<double> &bct,
74  const Vector &pa_data,
75  const Vector &x,
76  Vector &y);
77 
78 template<int T_D1D = 0, int T_Q1D = 0>
79 void SmemPAHcurlMassApply3D(const int D1D,
80  const int Q1D,
81  const int NE,
82  const bool symmetric,
83  const Array<double> &bo,
84  const Array<double> &bc,
85  const Array<double> &bot,
86  const Array<double> &bct,
87  const Vector &pa_data,
88  const Vector &x,
89  Vector &y);
90 
91 void PAHdivSetup2D(const int Q1D,
92  const int NE,
93  const Array<double> &w,
94  const Vector &j,
95  Vector &coeff_,
96  Vector &op);
97 
98 void PAHdivSetup3D(const int Q1D,
99  const int NE,
100  const Array<double> &w,
101  const Vector &j,
102  Vector &coeff_,
103  Vector &op);
104 
105 void PAHcurlH1Apply2D(const int D1D,
106  const int Q1D,
107  const int NE,
108  const Array<double> &bc,
109  const Array<double> &gc,
110  const Array<double> &bot,
111  const Array<double> &bct,
112  const Vector &pa_data,
113  const Vector &x,
114  Vector &y);
115 
116 void PAHcurlH1ApplyTranspose2D(const int D1D,
117  const int Q1D,
118  const int NE,
119  const Array<double> &bc,
120  const Array<double> &bo,
121  const Array<double> &bct,
122  const Array<double> &gct,
123  const Vector &pa_data,
124  const Vector &x,
125  Vector &y);
126 
127 void PAHcurlH1Apply3D(const int D1D,
128  const int Q1D,
129  const int NE,
130  const Array<double> &bc,
131  const Array<double> &gc,
132  const Array<double> &bot,
133  const Array<double> &bct,
134  const Vector &pa_data,
135  const Vector &x,
136  Vector &y);
137 
138 void PAHcurlH1ApplyTranspose3D(const int D1D,
139  const int Q1D,
140  const int NE,
141  const Array<double> &bc,
142  const Array<double> &bo,
143  const Array<double> &bct,
144  const Array<double> &gct,
145  const Vector &pa_data,
146  const Vector &x,
147  Vector &y);
148 
149 void PAHdivMassAssembleDiagonal2D(const int D1D,
150  const int Q1D,
151  const int NE,
152  const Array<double> &Bo_,
153  const Array<double> &Bc_,
154  const Vector &op_,
155  Vector &diag_);
156 
157 void PAHdivMassAssembleDiagonal3D(const int D1D,
158  const int Q1D,
159  const int NE,
160  const Array<double> &Bo_,
161  const Array<double> &Bc_,
162  const Vector &op_,
163  Vector &diag_);
164 
165 void PAHdivMassApply2D(const int D1D,
166  const int Q1D,
167  const int NE,
168  const Array<double> &Bo_,
169  const Array<double> &Bc_,
170  const Array<double> &Bot_,
171  const Array<double> &Bct_,
172  const Vector &op_,
173  const Vector &x_,
174  Vector &y_);
175 
176 void PAHdivMassApply3D(const int D1D,
177  const int Q1D,
178  const int NE,
179  const Array<double> &Bo_,
180  const Array<double> &Bc_,
181  const Array<double> &Bot_,
182  const Array<double> &Bct_,
183  const Vector &op_,
184  const Vector &x_,
185  Vector &y_);
186 
187 void PAHcurlL2Setup(const int NQ,
188  const int coeffDim,
189  const int NE,
190  const Array<double> &w,
191  Vector &coeff_,
192  Vector &op);
193 
194 // PA H(curl) x H(div) mass assemble 3D kernel, with factor
195 // dF^{-1} C dF for a vector or matrix coefficient C.
196 // If transpose, use dF^T C dF^{-T} for H(div) x H(curl).
197 void PAHcurlHdivSetup3D(const int Q1D,
198  const int coeffDim,
199  const int NE,
200  const bool transpose,
201  const Array<double> &w_,
202  const Vector &j,
203  Vector &coeff_,
204  Vector &op)
205 {
206  const bool symmetric = (coeffDim != 9);
207  auto W = Reshape(w_.Read(), Q1D, Q1D, Q1D);
208  auto J = Reshape(j.Read(), Q1D, Q1D, Q1D, 3, 3, NE);
209  auto coeff = Reshape(coeff_.Read(), coeffDim, Q1D, Q1D, Q1D, NE);
210  auto y = Reshape(op.Write(), 9, Q1D, Q1D, Q1D, NE);
211 
212  const int i11 = 0;
213  const int i12 = transpose ? 3 : 1;
214  const int i13 = transpose ? 6 : 2;
215  const int i21 = transpose ? 1 : 3;
216  const int i22 = 4;
217  const int i23 = transpose ? 7 : 5;
218  const int i31 = transpose ? 2 : 6;
219  const int i32 = transpose ? 5 : 7;
220  const int i33 = 8;
221 
222  MFEM_FORALL_3D(e, NE, Q1D, Q1D, Q1D,
223  {
224  MFEM_FOREACH_THREAD(qx,x,Q1D)
225  {
226  MFEM_FOREACH_THREAD(qy,y,Q1D)
227  {
228  MFEM_FOREACH_THREAD(qz,z,Q1D)
229  {
230  const double J11 = J(qx,qy,qz,0,0,e);
231  const double J21 = J(qx,qy,qz,1,0,e);
232  const double J31 = J(qx,qy,qz,2,0,e);
233  const double J12 = J(qx,qy,qz,0,1,e);
234  const double J22 = J(qx,qy,qz,1,1,e);
235  const double J32 = J(qx,qy,qz,2,1,e);
236  const double J13 = J(qx,qy,qz,0,2,e);
237  const double J23 = J(qx,qy,qz,1,2,e);
238  const double J33 = J(qx,qy,qz,2,2,e);
239  const double detJ = J11 * (J22 * J33 - J32 * J23) -
240  /* */ J21 * (J12 * J33 - J32 * J13) +
241  /* */ J31 * (J12 * J23 - J22 * J13);
242  const double w_detJ = W(qx,qy,qz) / detJ;
243  // adj(J)
244  const double A11 = (J22 * J33) - (J23 * J32);
245  const double A12 = (J32 * J13) - (J12 * J33);
246  const double A13 = (J12 * J23) - (J22 * J13);
247  const double A21 = (J31 * J23) - (J21 * J33);
248  const double A22 = (J11 * J33) - (J13 * J31);
249  const double A23 = (J21 * J13) - (J11 * J23);
250  const double A31 = (J21 * J32) - (J31 * J22);
251  const double A32 = (J31 * J12) - (J11 * J32);
252  const double A33 = (J11 * J22) - (J12 * J21);
253 
254  if (coeffDim == 6 || coeffDim == 9) // Matrix coefficient version
255  {
256  // First compute entries of R = M^T J
257  const double M11 = (!symmetric) ? coeff(i11,qx,qy,qz,e) : coeff(0,qx,qy,qz,e);
258  const double M12 = (!symmetric) ? coeff(i12,qx,qy,qz,e) : coeff(1,qx,qy,qz,e);
259  const double M13 = (!symmetric) ? coeff(i13,qx,qy,qz,e) : coeff(2,qx,qy,qz,e);
260  const double M21 = (!symmetric) ? coeff(i21,qx,qy,qz,e) : M12;
261  const double M22 = (!symmetric) ? coeff(i22,qx,qy,qz,e) : coeff(3,qx,qy,qz,e);
262  const double M23 = (!symmetric) ? coeff(i23,qx,qy,qz,e) : coeff(4,qx,qy,qz,e);
263  const double M31 = (!symmetric) ? coeff(i31,qx,qy,qz,e) : M13;
264  const double M32 = (!symmetric) ? coeff(i32,qx,qy,qz,e) : M23;
265  const double M33 = (!symmetric) ? coeff(i33,qx,qy,qz,e) : coeff(5,qx,qy,qz,e);
266 
267  const double R11 = M11*J11 + M21*J21 + M31*J31;
268  const double R12 = M11*J12 + M21*J22 + M31*J32;
269  const double R13 = M11*J13 + M21*J23 + M31*J33;
270  const double R21 = M12*J11 + M22*J21 + M32*J31;
271  const double R22 = M12*J12 + M22*J22 + M32*J32;
272  const double R23 = M12*J13 + M22*J23 + M32*J33;
273  const double R31 = M13*J11 + M23*J21 + M33*J31;
274  const double R32 = M13*J12 + M23*J22 + M33*J32;
275  const double R33 = M13*J13 + M23*J23 + M33*J33;
276 
277  // y = (J^{-1} M^T J)^T
278  y(i11,qx,qy,qz,e) = w_detJ * (A11*R11 + A12*R21 + A13*R31); // 1,1
279  y(i21,qx,qy,qz,e) = w_detJ * (A11*R12 + A12*R22 + A13*R32); // 1,2
280  y(i31,qx,qy,qz,e) = w_detJ * (A11*R13 + A12*R23 + A13*R33); // 1,3
281  y(i12,qx,qy,qz,e) = w_detJ * (A21*R11 + A22*R21 + A23*R31); // 2,1
282  y(i22,qx,qy,qz,e) = w_detJ * (A21*R12 + A22*R22 + A23*R32); // 2,2
283  y(i32,qx,qy,qz,e) = w_detJ * (A21*R13 + A22*R23 + A23*R33); // 2,3
284  y(i13,qx,qy,qz,e) = w_detJ * (A31*R11 + A32*R21 + A33*R31); // 3,1
285  y(i23,qx,qy,qz,e) = w_detJ * (A31*R12 + A32*R22 + A33*R32); // 3,2
286  y(i33,qx,qy,qz,e) = w_detJ * (A31*R13 + A32*R23 + A33*R33); // 3,3
287  }
288  else if (coeffDim == 3) // Vector coefficient version
289  {
290  const double D1 = coeff(0,qx,qy,qz,e);
291  const double D2 = coeff(1,qx,qy,qz,e);
292  const double D3 = coeff(2,qx,qy,qz,e);
293  // detJ J^{-1} DJ = adj(J) DJ
294  // transpose
295  y(i11,qx,qy,qz,e) = w_detJ * (D1*A11*J11 + D2*A12*J21 + D3*A13*J31); // 1,1
296  y(i21,qx,qy,qz,e) = w_detJ * (D1*A11*J12 + D2*A12*J22 + D3*A13*J32); // 1,2
297  y(i31,qx,qy,qz,e) = w_detJ * (D1*A11*J13 + D2*A12*J23 + D3*A13*J33); // 1,3
298  y(i12,qx,qy,qz,e) = w_detJ * (D1*A21*J11 + D2*A22*J21 + D3*A23*J31); // 2,1
299  y(i22,qx,qy,qz,e) = w_detJ * (D1*A21*J12 + D2*A22*J22 + D3*A23*J32); // 2,2
300  y(i32,qx,qy,qz,e) = w_detJ * (D1*A21*J13 + D2*A22*J23 + D3*A23*J33); // 2,3
301  y(i13,qx,qy,qz,e) = w_detJ * (D1*A31*J11 + D2*A32*J21 + D3*A33*J31); // 3,1
302  y(i23,qx,qy,qz,e) = w_detJ * (D1*A31*J12 + D2*A32*J22 + D3*A33*J32); // 3,2
303  y(i33,qx,qy,qz,e) = w_detJ * (D1*A31*J13 + D2*A32*J23 + D3*A33*J33); // 3,3
304  }
305  }
306  }
307  }
308  });
309 }
310 
311 // PA H(curl) x H(div) mass assemble 2D kernel, with factor
312 // dF^{-1} C dF for a vector or matrix coefficient C.
313 // If transpose, use dF^T C dF^{-T} for H(div) x H(curl).
314 void PAHcurlHdivSetup2D(const int Q1D,
315  const int coeffDim,
316  const int NE,
317  const bool transpose,
318  const Array<double> &w_,
319  const Vector &j,
320  Vector &coeff_,
321  Vector &op)
322 {
323  const bool symmetric = (coeffDim != 4);
324  auto W = Reshape(w_.Read(), Q1D, Q1D);
325  auto J = Reshape(j.Read(), Q1D, Q1D, 2, 2, NE);
326  auto coeff = Reshape(coeff_.Read(), coeffDim, Q1D, Q1D, NE);
327  auto y = Reshape(op.Write(), 4, Q1D, Q1D, NE);
328 
329  const int i11 = 0;
330  const int i12 = transpose ? 2 : 1;
331  const int i21 = transpose ? 1 : 2;
332  const int i22 = 3;
333 
334  MFEM_FORALL_2D(e, NE, Q1D, Q1D, 1,
335  {
336  MFEM_FOREACH_THREAD(qx,x,Q1D)
337  {
338  MFEM_FOREACH_THREAD(qy,y,Q1D)
339  {
340  const double J11 = J(qx,qy,0,0,e);
341  const double J21 = J(qx,qy,1,0,e);
342  const double J12 = J(qx,qy,0,1,e);
343  const double J22 = J(qx,qy,1,1,e);
344  const double w_detJ = W(qx,qy) / ((J11*J22) - (J21*J12));
345 
346  if (coeffDim == 3 || coeffDim == 4) // Matrix coefficient version
347  {
348  // First compute entries of R = MJ
349  const double M11 = coeff(i11,qx,qy,e);
350  const double M12 = (!symmetric) ? coeff(i12,qx,qy,e) : coeff(1,qx,qy,e);
351  const double M21 = (!symmetric) ? coeff(i21,qx,qy,e) : M12;
352  const double M22 = (!symmetric) ? coeff(i22,qx,qy,e) : coeff(2,qx,qy,e);
353 
354  // J^{-1} M^T
355  const double R11 = ( J22*M11 - J12*M12); // 1,1
356  const double R12 = ( J22*M21 - J12*M22); // 1,2
357  const double R21 = (-J21*M11 + J11*M12); // 2,1
358  const double R22 = (-J21*M21 + J11*M22); // 2,2
359 
360  // (RJ)^T
361  y(i11,qx,qy,e) = w_detJ * (R11*J11 + R12*J21); // 1,1
362  y(i21,qx,qy,e) = w_detJ * (R11*J12 + R12*J22); // 1,2 (transpose)
363  y(i12,qx,qy,e) = w_detJ * (R21*J11 + R22*J21); // 2,1 (transpose)
364  y(i22,qx,qy,e) = w_detJ * (R21*J12 + R22*J22); // 2,2
365  }
366  else if (coeffDim == 2) // Vector coefficient version
367  {
368  const double D1 = coeff(0,qx,qy,e);
369  const double D2 = coeff(1,qx,qy,e);
370  const double R11 = D1*J11;
371  const double R12 = D1*J12;
372  const double R21 = D2*J21;
373  const double R22 = D2*J22;
374  y(i11,qx,qy,e) = w_detJ * ( J22*R11 - J12*R21); // 1,1
375  y(i12,qx,qy,e) = w_detJ * ( J22*R12 - J12*R22); // 1,2
376  y(i21,qx,qy,e) = w_detJ * (-J21*R11 + J11*R21); // 2,1
377  y(i22,qx,qy,e) = w_detJ * (-J21*R12 + J11*R22); // 2,2
378  }
379  }
380  }
381  });
382 }
383 
384 // Mass operator for H(curl) and H(div) functions, using Piola transformations
385 // u = dF^{-T} \hat{u} in H(curl), v = (1 / det dF) dF \hat{v} in H(div).
386 void PAHcurlHdivMassApply3D(const int D1D,
387  const int D1Dtest,
388  const int Q1D,
389  const int NE,
390  const bool scalarCoeff,
391  const bool trialHcurl,
392  const Array<double> &Bo_,
393  const Array<double> &Bc_,
394  const Array<double> &Bot_,
395  const Array<double> &Bct_,
396  const Vector &op_,
397  const Vector &x_,
398  Vector &y_)
399 {
400  constexpr static int MAX_D1D = HCURL_MAX_D1D;
401  constexpr static int MAX_Q1D = HCURL_MAX_Q1D;
402 
403  MFEM_VERIFY(D1D <= MAX_D1D, "Error: D1D > MAX_D1D");
404  MFEM_VERIFY(Q1D <= MAX_Q1D, "Error: Q1D > MAX_Q1D");
405  constexpr static int VDIM = 3;
406 
407  auto Bo = Reshape(Bo_.Read(), Q1D, D1D-1);
408  auto Bc = Reshape(Bc_.Read(), Q1D, D1D);
409  auto Bot = Reshape(Bot_.Read(), D1Dtest-1, Q1D);
410  auto Bct = Reshape(Bct_.Read(), D1Dtest, Q1D);
411  auto op = Reshape(op_.Read(), scalarCoeff ? 1 : 9, Q1D, Q1D, Q1D, NE);
412  auto x = Reshape(x_.Read(), 3*(D1D-1)*D1D*(trialHcurl ? D1D : D1D-1), NE);
413  auto y = Reshape(y_.ReadWrite(), 3*(D1Dtest-1)*D1Dtest*
414  (trialHcurl ? D1Dtest-1 : D1Dtest), NE);
415 
416  MFEM_FORALL(e, NE,
417  {
418  double mass[MAX_Q1D][MAX_Q1D][MAX_Q1D][VDIM];
419 
420  for (int qz = 0; qz < Q1D; ++qz)
421  {
422  for (int qy = 0; qy < Q1D; ++qy)
423  {
424  for (int qx = 0; qx < Q1D; ++qx)
425  {
426  for (int c = 0; c < VDIM; ++c)
427  {
428  mass[qz][qy][qx][c] = 0.0;
429  }
430  }
431  }
432  }
433 
434  int osc = 0;
435  for (int c = 0; c < VDIM; ++c) // loop over x, y, z trial components
436  {
437  const int D1Dz = trialHcurl ? ((c == 2) ? D1D - 1 : D1D) :
438  ((c == 2) ? D1D : D1D - 1);
439  const int D1Dy = trialHcurl ? ((c == 1) ? D1D - 1 : D1D) :
440  ((c == 1) ? D1D : D1D - 1);
441  const int D1Dx = trialHcurl ? ((c == 0) ? D1D - 1 : D1D) :
442  ((c == 0) ? D1D : D1D - 1);
443 
444  for (int dz = 0; dz < D1Dz; ++dz)
445  {
446  double massXY[MAX_Q1D][MAX_Q1D];
447  for (int qy = 0; qy < Q1D; ++qy)
448  {
449  for (int qx = 0; qx < Q1D; ++qx)
450  {
451  massXY[qy][qx] = 0.0;
452  }
453  }
454 
455  for (int dy = 0; dy < D1Dy; ++dy)
456  {
457  double massX[MAX_Q1D];
458  for (int qx = 0; qx < Q1D; ++qx)
459  {
460  massX[qx] = 0.0;
461  }
462 
463  for (int dx = 0; dx < D1Dx; ++dx)
464  {
465  const double t = x(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e);
466  for (int qx = 0; qx < Q1D; ++qx)
467  {
468  massX[qx] += t * (trialHcurl ? ((c == 0) ? Bo(qx,dx) : Bc(qx,dx)) :
469  ((c == 0) ? Bc(qx,dx) : Bo(qx,dx)));
470  }
471  }
472 
473  for (int qy = 0; qy < Q1D; ++qy)
474  {
475  const double wy = trialHcurl ? ((c == 1) ? Bo(qy,dy) : Bc(qy,dy)) :
476  ((c == 1) ? Bc(qy,dy) : Bo(qy,dy));
477  for (int qx = 0; qx < Q1D; ++qx)
478  {
479  const double wx = massX[qx];
480  massXY[qy][qx] += wx * wy;
481  }
482  }
483  }
484 
485  for (int qz = 0; qz < Q1D; ++qz)
486  {
487  const double wz = trialHcurl ? ((c == 2) ? Bo(qz,dz) : Bc(qz,dz)) :
488  ((c == 2) ? Bc(qz,dz) : Bo(qz,dz));
489  for (int qy = 0; qy < Q1D; ++qy)
490  {
491  for (int qx = 0; qx < Q1D; ++qx)
492  {
493  mass[qz][qy][qx][c] += massXY[qy][qx] * wz;
494  }
495  }
496  }
497  }
498 
499  osc += D1Dx * D1Dy * D1Dz;
500  } // loop (c) over components
501 
502  // Apply D operator.
503  for (int qz = 0; qz < Q1D; ++qz)
504  {
505  for (int qy = 0; qy < Q1D; ++qy)
506  {
507  for (int qx = 0; qx < Q1D; ++qx)
508  {
509  const double O11 = op(0,qx,qy,qz,e);
510  const double O12 = scalarCoeff ? 0.0 : op(1,qx,qy,qz,e);
511  const double O13 = scalarCoeff ? 0.0 : op(2,qx,qy,qz,e);
512  const double O21 = scalarCoeff ? 0.0 : op(3,qx,qy,qz,e);
513  const double O22 = scalarCoeff ? O11 : op(4,qx,qy,qz,e);
514  const double O23 = scalarCoeff ? 0.0 : op(5,qx,qy,qz,e);
515  const double O31 = scalarCoeff ? 0.0 : op(6,qx,qy,qz,e);
516  const double O32 = scalarCoeff ? 0.0 : op(7,qx,qy,qz,e);
517  const double O33 = scalarCoeff ? O11 : op(8,qx,qy,qz,e);
518  const double massX = mass[qz][qy][qx][0];
519  const double massY = mass[qz][qy][qx][1];
520  const double massZ = mass[qz][qy][qx][2];
521  mass[qz][qy][qx][0] = (O11*massX)+(O12*massY)+(O13*massZ);
522  mass[qz][qy][qx][1] = (O21*massX)+(O22*massY)+(O23*massZ);
523  mass[qz][qy][qx][2] = (O31*massX)+(O32*massY)+(O33*massZ);
524  }
525  }
526  }
527 
528  for (int qz = 0; qz < Q1D; ++qz)
529  {
530  double massXY[HDIV_MAX_D1D][HDIV_MAX_D1D];
531 
532  osc = 0;
533  for (int c = 0; c < VDIM; ++c) // loop over x, y, z test components
534  {
535  const int D1Dz = trialHcurl ? ((c == 2) ? D1Dtest : D1Dtest - 1) :
536  ((c == 2) ? D1Dtest - 1 : D1Dtest);
537  const int D1Dy = trialHcurl ? ((c == 1) ? D1Dtest : D1Dtest - 1) :
538  ((c == 1) ? D1Dtest - 1 : D1Dtest);
539  const int D1Dx = trialHcurl ? ((c == 0) ? D1Dtest : D1Dtest - 1) :
540  ((c == 0) ? D1Dtest - 1 : D1Dtest);
541 
542  for (int dy = 0; dy < D1Dy; ++dy)
543  {
544  for (int dx = 0; dx < D1Dx; ++dx)
545  {
546  massXY[dy][dx] = 0.0;
547  }
548  }
549  for (int qy = 0; qy < Q1D; ++qy)
550  {
551  double massX[HDIV_MAX_D1D];
552  for (int dx = 0; dx < D1Dx; ++dx)
553  {
554  massX[dx] = 0.0;
555  }
556  for (int qx = 0; qx < Q1D; ++qx)
557  {
558  for (int dx = 0; dx < D1Dx; ++dx)
559  {
560  massX[dx] += mass[qz][qy][qx][c] * (trialHcurl ?
561  ((c == 0) ? Bct(dx,qx) : Bot(dx,qx)) :
562  ((c == 0) ? Bot(dx,qx) : Bct(dx,qx)));
563  }
564  }
565  for (int dy = 0; dy < D1Dy; ++dy)
566  {
567  const double wy = trialHcurl ? ((c == 1) ? Bct(dy,qy) : Bot(dy,qy)) :
568  ((c == 1) ? Bot(dy,qy) : Bct(dy,qy));
569  for (int dx = 0; dx < D1Dx; ++dx)
570  {
571  massXY[dy][dx] += massX[dx] * wy;
572  }
573  }
574  }
575 
576  for (int dz = 0; dz < D1Dz; ++dz)
577  {
578  const double wz = trialHcurl ? ((c == 2) ? Bct(dz,qz) : Bot(dz,qz)) :
579  ((c == 2) ? Bot(dz,qz) : Bct(dz,qz));
580  for (int dy = 0; dy < D1Dy; ++dy)
581  {
582  for (int dx = 0; dx < D1Dx; ++dx)
583  {
584  y(dx + ((dy + (dz * D1Dy)) * D1Dx) + osc, e) +=
585  massXY[dy][dx] * wz;
586  }
587  }
588  }
589 
590  osc += D1Dx * D1Dy * D1Dz;
591  } // loop c
592  } // loop qz
593  }); // end of element loop
594 }
595 
596 // Mass operator for H(curl) and H(div) functions, using Piola transformations
597 // u = dF^{-T} \hat{u} in H(curl), v = (1 / det dF) dF \hat{v} in H(div).
598 void PAHcurlHdivMassApply2D(const int D1D,
599  const int D1Dtest,
600  const int Q1D,
601  const int NE,
602  const bool scalarCoeff,
603  const bool trialHcurl,
604  const Array<double> &Bo_,
605  const Array<double> &Bc_,
606  const Array<double> &Bot_,
607  const Array<double> &Bct_,
608  const Vector &op_,
609  const Vector &x_,
610  Vector &y_)
611 {
612  constexpr static int MAX_D1D = HCURL_MAX_D1D;
613  constexpr static int MAX_Q1D = HCURL_MAX_Q1D;
614 
615  MFEM_VERIFY(D1D <= MAX_D1D, "Error: D1D > MAX_D1D");
616  MFEM_VERIFY(Q1D <= MAX_Q1D, "Error: Q1D > MAX_Q1D");
617  constexpr static int VDIM = 2;
618 
619  auto Bo = Reshape(Bo_.Read(), Q1D, D1D-1);
620  auto Bc = Reshape(Bc_.Read(), Q1D, D1D);
621  auto Bot = Reshape(Bot_.Read(), D1Dtest-1, Q1D);
622  auto Bct = Reshape(Bct_.Read(), D1Dtest, Q1D);
623  auto op = Reshape(op_.Read(), scalarCoeff ? 1 : 4, Q1D, Q1D, NE);
624  auto x = Reshape(x_.Read(), 2*(D1D-1)*D1D, NE);
625  auto y = Reshape(y_.ReadWrite(), 2*(D1Dtest-1)*D1Dtest, NE);
626 
627  MFEM_FORALL(e, NE,
628  {
629  double mass[MAX_Q1D][MAX_Q1D][VDIM];
630 
631  for (int qy = 0; qy < Q1D; ++qy)
632  {
633  for (int qx = 0; qx < Q1D; ++qx)
634  {
635  for (int c = 0; c < VDIM; ++c)
636  {
637  mass[qy][qx][c] = 0.0;
638  }
639  }
640  }
641 
642  int osc = 0;
643  for (int c = 0; c < VDIM; ++c) // loop over x, y trial components
644  {
645  const int D1Dy = trialHcurl ? ((c == 1) ? D1D - 1 : D1D) :
646  ((c == 1) ? D1D : D1D - 1);
647  const int D1Dx = trialHcurl ? ((c == 0) ? D1D - 1 : D1D) :
648  ((c == 0) ? D1D : D1D - 1);
649 
650  for (int dy = 0; dy < D1Dy; ++dy)
651  {
652  double massX[MAX_Q1D];
653  for (int qx = 0; qx < Q1D; ++qx)
654  {
655  massX[qx] = 0.0;
656  }
657 
658  for (int dx = 0; dx < D1Dx; ++dx)
659  {
660  const double t = x(dx + (dy * D1Dx) + osc, e);
661  for (int qx = 0; qx < Q1D; ++qx)
662  {
663  massX[qx] += t * (trialHcurl ? ((c == 0) ? Bo(qx,dx) : Bc(qx,dx)) :
664  ((c == 0) ? Bc(qx,dx) : Bo(qx,dx)));
665  }
666  }
667 
668  for (int qy = 0; qy < Q1D; ++qy)
669  {
670  const double wy = trialHcurl ? ((c == 1) ? Bo(qy,dy) : Bc(qy,dy)) :
671  ((c == 1) ? Bc(qy,dy) : Bo(qy,dy));
672  for (int qx = 0; qx < Q1D; ++qx)
673  {
674  mass[qy][qx][c] += massX[qx] * wy;
675  }
676  }
677  }
678 
679  osc += D1Dx * D1Dy;
680  } // loop (c) over components
681 
682  // Apply D operator.
683  for (int qy = 0; qy < Q1D; ++qy)
684  {
685  for (int qx = 0; qx < Q1D; ++qx)
686  {
687  const double O11 = op(0,qx,qy,e);
688  const double O12 = scalarCoeff ? 0.0 : op(1,qx,qy,e);
689  const double O21 = scalarCoeff ? 0.0 : op(2,qx,qy,e);
690  const double O22 = scalarCoeff ? O11 : op(3,qx,qy,e);
691  const double massX = mass[qy][qx][0];
692  const double massY = mass[qy][qx][1];
693  mass[qy][qx][0] = (O11*massX)+(O12*massY);
694  mass[qy][qx][1] = (O21*massX)+(O22*massY);
695  }
696  }
697 
698  osc = 0;
699  for (int c = 0; c < VDIM; ++c) // loop over x, y test components
700  {
701  const int D1Dy = trialHcurl ? ((c == 1) ? D1Dtest : D1Dtest - 1) :
702  ((c == 1) ? D1Dtest - 1 : D1Dtest);
703  const int D1Dx = trialHcurl ? ((c == 0) ? D1Dtest : D1Dtest - 1) :
704  ((c == 0) ? D1Dtest - 1 : D1Dtest);
705 
706  for (int qy = 0; qy < Q1D; ++qy)
707  {
708  double massX[HDIV_MAX_D1D];
709  for (int dx = 0; dx < D1Dx; ++dx)
710  {
711  massX[dx] = 0.0;
712  }
713  for (int qx = 0; qx < Q1D; ++qx)
714  {
715  for (int dx = 0; dx < D1Dx; ++dx)
716  {
717  massX[dx] += mass[qy][qx][c] * (trialHcurl ?
718  ((c == 0) ? Bct(dx,qx) : Bot(dx,qx)) :
719  ((c == 0) ? Bot(dx,qx) : Bct(dx,qx)));
720  }
721  }
722  for (int dy = 0; dy < D1Dy; ++dy)
723  {
724  const double wy = trialHcurl ? ((c == 1) ? Bct(dy,qy) : Bot(dy,qy)) :
725  ((c == 1) ? Bot(dy,qy) : Bct(dy,qy));
726  for (int dx = 0; dx < D1Dx; ++dx)
727  {
728  y(dx + (dy * D1Dx) + osc, e) += massX[dx] * wy;
729  }
730  }
731  }
732 
733  osc += D1Dx * D1Dy;
734  } // loop c
735  }); // end of element loop
736 }
737 
738 void VectorFEMassIntegrator::AssemblePA(const FiniteElementSpace &fes)
739 {
740  AssemblePA(fes, fes);
741 }
742 
743 void VectorFEMassIntegrator::AssemblePA(const FiniteElementSpace &trial_fes,
744  const FiniteElementSpace &test_fes)
745 {
746  // Assumes tensor-product elements
747  Mesh *mesh = trial_fes.GetMesh();
748 
749  const FiniteElement *trial_fel = trial_fes.GetFE(0);
750  const VectorTensorFiniteElement *trial_el =
751  dynamic_cast<const VectorTensorFiniteElement*>(trial_fel);
752  MFEM_VERIFY(trial_el != NULL, "Only VectorTensorFiniteElement is supported!");
753 
754  const FiniteElement *test_fel = test_fes.GetFE(0);
755  const VectorTensorFiniteElement *test_el =
756  dynamic_cast<const VectorTensorFiniteElement*>(test_fel);
757  MFEM_VERIFY(test_el != NULL, "Only VectorTensorFiniteElement is supported!");
758 
759  const IntegrationRule *ir
760  = IntRule ? IntRule : &MassIntegrator::GetRule(*trial_el, *trial_el,
761  *mesh->GetElementTransformation(0));
762  const int dims = trial_el->GetDim();
763  MFEM_VERIFY(dims == 2 || dims == 3, "");
764 
765  const int symmDims = (dims * (dims + 1)) / 2; // 1x1: 1, 2x2: 3, 3x3: 6
766  nq = ir->GetNPoints();
767  dim = mesh->Dimension();
768  MFEM_VERIFY(dim == 2 || dim == 3, "");
769 
770  ne = trial_fes.GetNE();
771  MFEM_VERIFY(ne == test_fes.GetNE(),
772  "Different meshes for test and trial spaces");
773  geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS);
774  mapsC = &trial_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
775  mapsO = &trial_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
776  dofs1D = mapsC->ndof;
777  quad1D = mapsC->nqpt;
778 
779  mapsCtest = &test_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
780  mapsOtest = &test_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
781  dofs1Dtest = mapsCtest->ndof;
782 
783  MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
784 
785  trial_fetype = trial_el->GetDerivType();
786  test_fetype = test_el->GetDerivType();
787 
788  const int MQsymmDim = SMQ ? (SMQ->GetSize() * (SMQ->GetSize() + 1)) / 2 : 0;
789  const int MQfullDim = MQ ? (MQ->GetHeight() * MQ->GetWidth()) : 0;
790  const int MQdim = MQ ? MQfullDim : MQsymmDim;
791  const int coeffDim = (MQ || SMQ) ? MQdim : (DQ ? DQ->GetVDim() : 1);
792 
793  symmetric = (MQ == NULL);
794 
795  const bool trial_curl = (trial_fetype == mfem::FiniteElement::CURL);
796  const bool trial_div = (trial_fetype == mfem::FiniteElement::DIV);
797  const bool test_curl = (test_fetype == mfem::FiniteElement::CURL);
798  const bool test_div = (test_fetype == mfem::FiniteElement::DIV);
799 
800  if ((trial_curl && test_div) || (trial_div && test_curl))
801  pa_data.SetSize((coeffDim == 1 ? 1 : dim*dim) * nq * ne,
802  Device::GetMemoryType());
803  else
804  pa_data.SetSize((symmetric ? symmDims : MQfullDim) * nq * ne,
805  Device::GetMemoryType());
806 
807  Vector coeff(coeffDim * ne * nq);
808  coeff = 1.0;
809  auto coeffh = Reshape(coeff.HostWrite(), coeffDim, nq, ne);
810  if (Q || DQ || MQ || SMQ)
811  {
812  Vector DM(DQ ? coeffDim : 0);
813  DenseMatrix M;
814  DenseSymmetricMatrix SM;
815 
816  if (DQ)
817  {
818  MFEM_VERIFY(coeffDim == dim, "");
819  }
820  if (MQ)
821  {
822  MFEM_VERIFY(coeffDim == MQdim, "");
823  MFEM_VERIFY(MQ->GetHeight() == dim && MQ->GetWidth() == dim, "");
824  M.SetSize(dim);
825  }
826  if (SMQ)
827  {
828  MFEM_VERIFY(SMQ->GetSize() == dim, "");
829  SM.SetSize(dim);
830  }
831 
832  for (int e=0; e<ne; ++e)
833  {
834  ElementTransformation *tr = mesh->GetElementTransformation(e);
835  for (int p=0; p<nq; ++p)
836  {
837  if (MQ)
838  {
839  MQ->Eval(M, *tr, ir->IntPoint(p));
840 
841  for (int i=0; i<dim; ++i)
842  for (int j=0; j<dim; ++j)
843  {
844  coeffh(j+(i*dim), p, e) = M(i,j);
845  }
846  }
847  else if (SMQ)
848  {
849  SMQ->Eval(SM, *tr, ir->IntPoint(p));
850  int cnt = 0;
851  for (int i=0; i<dim; ++i)
852  for (int j=i; j<dim; ++j, ++cnt)
853  {
854  coeffh(cnt, p, e) = SM(i,j);
855  }
856  }
857  else if (DQ)
858  {
859  DQ->Eval(DM, *tr, ir->IntPoint(p));
860  for (int i=0; i<coeffDim; ++i)
861  {
862  coeffh(i, p, e) = DM[i];
863  }
864  }
865  else
866  {
867  coeffh(0, p, e) = Q->Eval(*tr, ir->IntPoint(p));
868  }
869  }
870  }
871  }
872 
873  if (trial_curl && test_curl && dim == 3)
874  {
875  PADiffusionSetup3D(quad1D, coeffDim, ne, ir->GetWeights(), geom->J,
876  coeff, pa_data);
877  }
878  else if (trial_curl && test_curl && dim == 2)
879  {
880  PADiffusionSetup2D<2>(quad1D, coeffDim, ne, ir->GetWeights(), geom->J,
881  coeff, pa_data);
882  }
883  else if (trial_div && test_div && dim == 3)
884  {
885  PAHdivSetup3D(quad1D, ne, ir->GetWeights(), geom->J,
886  coeff, pa_data);
887  }
888  else if (trial_div && test_div && dim == 2)
889  {
890  PAHdivSetup2D(quad1D, ne, ir->GetWeights(), geom->J,
891  coeff, pa_data);
892  }
893  else if (((trial_curl && test_div) || (trial_div && test_curl)) &&
894  test_fel->GetOrder() == trial_fel->GetOrder())
895  {
896  if (coeffDim == 1)
897  {
898  PAHcurlL2Setup(nq, coeffDim, ne, ir->GetWeights(), coeff, pa_data);
899  }
900  else
901  {
902  const bool tr = (trial_div && test_curl);
903  if (dim == 3)
904  PAHcurlHdivSetup3D(quad1D, coeffDim, ne, tr, ir->GetWeights(),
905  geom->J, coeff, pa_data);
906  else
907  PAHcurlHdivSetup2D(quad1D, coeffDim, ne, tr, ir->GetWeights(),
908  geom->J, coeff, pa_data);
909  }
910  }
911  else
912  {
913  MFEM_ABORT("Unknown kernel.");
914  }
915 }
916 
917 void VectorFEMassIntegrator::AssembleDiagonalPA(Vector& diag)
918 {
919  if (dim == 3)
920  {
921  if (trial_fetype == mfem::FiniteElement::CURL && test_fetype == trial_fetype)
922  {
923  if (Device::Allows(Backend::DEVICE_MASK))
924  {
925  const int ID = (dofs1D << 4) | quad1D;
926  switch (ID)
927  {
928  case 0x23: return SmemPAHcurlMassAssembleDiagonal3D<2,3>(dofs1D, quad1D, ne,
929  symmetric,
930  mapsO->B, mapsC->B, pa_data, diag);
931  case 0x34: return SmemPAHcurlMassAssembleDiagonal3D<3,4>(dofs1D, quad1D, ne,
932  symmetric,
933  mapsO->B, mapsC->B, pa_data, diag);
934  case 0x45: return SmemPAHcurlMassAssembleDiagonal3D<4,5>(dofs1D, quad1D, ne,
935  symmetric,
936  mapsO->B, mapsC->B, pa_data, diag);
937  case 0x56: return SmemPAHcurlMassAssembleDiagonal3D<5,6>(dofs1D, quad1D, ne,
938  symmetric,
939  mapsO->B, mapsC->B, pa_data, diag);
940  default: return SmemPAHcurlMassAssembleDiagonal3D(dofs1D, quad1D, ne, symmetric,
941  mapsO->B, mapsC->B, pa_data, diag);
942  }
943  }
944  else
945  PAHcurlMassAssembleDiagonal3D(dofs1D, quad1D, ne, symmetric,
946  mapsO->B, mapsC->B, pa_data, diag);
947  }
948  else if (trial_fetype == mfem::FiniteElement::DIV &&
949  test_fetype == trial_fetype)
950  {
951  PAHdivMassAssembleDiagonal3D(dofs1D, quad1D, ne,
952  mapsO->B, mapsC->B, pa_data, diag);
953  }
954  else
955  {
956  MFEM_ABORT("Unknown kernel.");
957  }
958  }
959  else
960  {
961  if (trial_fetype == mfem::FiniteElement::CURL && test_fetype == trial_fetype)
962  {
963  PAHcurlMassAssembleDiagonal2D(dofs1D, quad1D, ne, symmetric,
964  mapsO->B, mapsC->B, pa_data, diag);
965  }
966  else if (trial_fetype == mfem::FiniteElement::DIV &&
967  test_fetype == trial_fetype)
968  {
969  PAHdivMassAssembleDiagonal2D(dofs1D, quad1D, ne,
970  mapsO->B, mapsC->B, pa_data, diag);
971  }
972  else
973  {
974  MFEM_ABORT("Unknown kernel.");
975  }
976  }
977 }
978 
979 void VectorFEMassIntegrator::AddMultPA(const Vector &x, Vector &y) const
980 {
981  const bool trial_curl = (trial_fetype == mfem::FiniteElement::CURL);
982  const bool trial_div = (trial_fetype == mfem::FiniteElement::DIV);
983  const bool test_curl = (test_fetype == mfem::FiniteElement::CURL);
984  const bool test_div = (test_fetype == mfem::FiniteElement::DIV);
985 
986  if (dim == 3)
987  {
988  if (trial_curl && test_curl)
989  {
990  if (Device::Allows(Backend::DEVICE_MASK))
991  {
992  const int ID = (dofs1D << 4) | quad1D;
993  switch (ID)
994  {
995  case 0x23: return SmemPAHcurlMassApply3D<2,3>(dofs1D, quad1D, ne, symmetric,
996  mapsO->B,
997  mapsC->B, mapsO->Bt,
998  mapsC->Bt, pa_data, x, y);
999  case 0x34: return SmemPAHcurlMassApply3D<3,4>(dofs1D, quad1D, ne, symmetric,
1000  mapsO->B,
1001  mapsC->B, mapsO->Bt,
1002  mapsC->Bt, pa_data, x, y);
1003  case 0x45: return SmemPAHcurlMassApply3D<4,5>(dofs1D, quad1D, ne, symmetric,
1004  mapsO->B,
1005  mapsC->B, mapsO->Bt,
1006  mapsC->Bt, pa_data, x, y);
1007  case 0x56: return SmemPAHcurlMassApply3D<5,6>(dofs1D, quad1D, ne, symmetric,
1008  mapsO->B,
1009  mapsC->B, mapsO->Bt,
1010  mapsC->Bt, pa_data, x, y);
1011  default: return SmemPAHcurlMassApply3D(dofs1D, quad1D, ne, symmetric, mapsO->B,
1012  mapsC->B,
1013  mapsO->Bt, mapsC->Bt, pa_data, x, y);
1014  }
1015  }
1016  else
1017  PAHcurlMassApply3D(dofs1D, quad1D, ne, symmetric, mapsO->B, mapsC->B, mapsO->Bt,
1018  mapsC->Bt, pa_data, x, y);
1019  }
1020  else if (trial_div && test_div)
1021  {
1022  PAHdivMassApply3D(dofs1D, quad1D, ne, mapsO->B, mapsC->B, mapsO->Bt,
1023  mapsC->Bt, pa_data, x, y);
1024  }
1025  else if (trial_curl && test_div)
1026  {
1027  const bool scalarCoeff = !(DQ || MQ || SMQ);
1028  PAHcurlHdivMassApply3D(dofs1D, dofs1Dtest, quad1D, ne, scalarCoeff,
1029  true, mapsO->B, mapsC->B, mapsOtest->Bt,
1030  mapsCtest->Bt, pa_data, x, y);
1031  }
1032  else if (trial_div && test_curl)
1033  {
1034  const bool scalarCoeff = !(DQ || MQ || SMQ);
1035  PAHcurlHdivMassApply3D(dofs1D, dofs1Dtest, quad1D, ne, scalarCoeff,
1036  false, mapsO->B, mapsC->B, mapsOtest->Bt,
1037  mapsCtest->Bt, pa_data, x, y);
1038  }
1039  else
1040  {
1041  MFEM_ABORT("Unknown kernel.");
1042  }
1043  }
1044  else
1045  {
1046  if (trial_curl && test_curl)
1047  {
1048  PAHcurlMassApply2D(dofs1D, quad1D, ne, symmetric, mapsO->B, mapsC->B,
1049  mapsO->Bt, mapsC->Bt, pa_data, x, y);
1050  }
1051  else if (trial_div && test_div)
1052  {
1053  PAHdivMassApply2D(dofs1D, quad1D, ne, mapsO->B, mapsC->B, mapsO->Bt,
1054  mapsC->Bt, pa_data, x, y);
1055  }
1056  else if ((trial_curl && test_div) || (trial_div && test_curl))
1057  {
1058  const bool scalarCoeff = !(DQ || MQ || SMQ);
1059  PAHcurlHdivMassApply2D(dofs1D, dofs1Dtest, quad1D, ne, scalarCoeff,
1060  trial_curl, mapsO->B, mapsC->B, mapsOtest->Bt,
1061  mapsCtest->Bt, pa_data, x, y);
1062  }
1063  else
1064  {
1065  MFEM_ABORT("Unknown kernel.");
1066  }
1067  }
1068 }
1069 
1070 void MixedVectorGradientIntegrator::AssemblePA(const FiniteElementSpace
1071  &trial_fes,
1072  const FiniteElementSpace &test_fes)
1073 {
1074  // Assumes tensor-product elements, with a vector test space and H^1 trial space.
1075  Mesh *mesh = trial_fes.GetMesh();
1076  const FiniteElement *trial_fel = trial_fes.GetFE(0);
1077  const FiniteElement *test_fel = test_fes.GetFE(0);
1078 
1079  const NodalTensorFiniteElement *trial_el =
1080  dynamic_cast<const NodalTensorFiniteElement*>(trial_fel);
1081  MFEM_VERIFY(trial_el != NULL, "Only NodalTensorFiniteElement is supported!");
1082 
1083  const VectorTensorFiniteElement *test_el =
1084  dynamic_cast<const VectorTensorFiniteElement*>(test_fel);
1085  MFEM_VERIFY(test_el != NULL, "Only VectorTensorFiniteElement is supported!");
1086 
1087  const IntegrationRule *ir
1088  = IntRule ? IntRule : &MassIntegrator::GetRule(*trial_el, *trial_el,
1089  *mesh->GetElementTransformation(0));
1090  const int dims = trial_el->GetDim();
1091  MFEM_VERIFY(dims == 2 || dims == 3, "");
1092 
1093  const int symmDims = (dims * (dims + 1)) / 2; // 1x1: 1, 2x2: 3, 3x3: 6
1094  const int nq = ir->GetNPoints();
1095  dim = mesh->Dimension();
1096  MFEM_VERIFY(dim == 2 || dim == 3, "");
1097 
1098  MFEM_VERIFY(trial_el->GetOrder() == test_el->GetOrder(), "");
1099 
1100  ne = trial_fes.GetNE();
1101  geom = mesh->GetGeometricFactors(*ir, GeometricFactors::JACOBIANS);
1102  mapsC = &test_el->GetDofToQuad(*ir, DofToQuad::TENSOR);
1103  mapsO = &test_el->GetDofToQuadOpen(*ir, DofToQuad::TENSOR);
1104  dofs1D = mapsC->ndof;
1105  quad1D = mapsC->nqpt;
1106 
1107  MFEM_VERIFY(dofs1D == mapsO->ndof + 1 && quad1D == mapsO->nqpt, "");
1108 
1109  pa_data.SetSize(symmDims * nq * ne, Device::GetMemoryType());
1110 
1111  Vector coeff(ne * nq);
1112  coeff = 1.0;
1113  if (Q)
1114  {
1115  for (int e=0; e<ne; ++e)
1116  {
1117  ElementTransformation *tr = mesh->GetElementTransformation(e);
1118  for (int p=0; p<nq; ++p)
1119  {
1120  coeff[p + (e * nq)] = Q->Eval(*tr, ir->IntPoint(p));
1121  }
1122  }
1123  }
1124 
1125  // Use the same setup functions as VectorFEMassIntegrator.
1126  if (test_el->GetDerivType() == mfem::FiniteElement::CURL && dim == 3)
1127  {
1128  PADiffusionSetup3D(quad1D, 1, ne, ir->GetWeights(), geom->J,
1129  coeff, pa_data);
1130  }
1131  else if (test_el->GetDerivType() == mfem::FiniteElement::CURL && dim == 2)
1132  {
1133  PADiffusionSetup2D<2>(quad1D, 1, ne, ir->GetWeights(), geom->J,
1134  coeff, pa_data);
1135  }
1136  else
1137  {
1138  MFEM_ABORT("Unknown kernel.");
1139  }
1140 }
1141 
1142 void MixedVectorGradientIntegrator::AddMultPA(const Vector &x, Vector &y) const
1143 {
1144  if (dim == 3)
1145  PAHcurlH1Apply3D(dofs1D, quad1D, ne, mapsC->B, mapsC->G,
1146  mapsO->Bt, mapsC->Bt, pa_data, x, y);
1147  else if (dim == 2)
1148  PAHcurlH1Apply2D(dofs1D, quad1D, ne, mapsC->B, mapsC->G,
1149  mapsO->Bt, mapsC->Bt, pa_data, x, y);
1150  else
1151  {
1152  MFEM_ABORT("Unsupported dimension!");
1153  }
1154 }
1155 
1156 void MixedVectorGradientIntegrator::AddMultTransposePA(const Vector &x,
1157  Vector &y) const
1158 {
1159  if (dim == 3)
1160  PAHcurlH1ApplyTranspose3D(dofs1D, quad1D, ne, mapsC->B, mapsO->B,
1161  mapsC->Bt, mapsC->Gt, pa_data, x, y);
1162  else if (dim == 2)
1163  PAHcurlH1ApplyTranspose2D(dofs1D, quad1D, ne, mapsC->B, mapsO->B,
1164  mapsC->Bt, mapsC->Gt, pa_data, x, y);
1165  else
1166  {
1167  MFEM_ABORT("Unsupported dimension!");
1168  }
1169 }
1170 
1171 } // namespace mfem
void PAHdivSetup3D(const int Q1D, const int NE, const Array< double > &w, const Vector &j, Vector &coeff_, Vector &op)
void PAHdivMassApply2D(const int D1D, const int Q1D, const int NE, const Array< double > &Bo_, const Array< double > &Bc_, const Array< double > &Bot_, const Array< double > &Bct_, const Vector &op_, const Vector &x_, Vector &y_)
void PAHcurlH1Apply3D(const int D1D, const int Q1D, const int NE, const Array< double > &bc, const Array< double > &gc, const Array< double > &bot, const Array< double > &bct, const Vector &pa_data, const Vector &x, Vector &y)
void PAHcurlMassApply3D(const int D1D, const int Q1D, const int NE, const bool symmetric, const Array< double > &bo, const Array< double > &bc, const Array< double > &bot, const Array< double > &bct, const Vector &pa_data, const Vector &x, Vector &y)
void PAHcurlL2Setup(const int NQ, const int coeffDim, const int NE, const Array< double > &w, Vector &coeff, Vector &op)
DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition: dtensor.hpp:131
void SmemPAHcurlMassApply3D(const int D1D, const int Q1D, const int NE, const bool symmetric, const Array< double > &bo, const Array< double > &bc, const Array< double > &bot, const Array< double > &bct, const Vector &pa_data, const Vector &x, Vector &y)
void PAHcurlHdivMassApply3D(const int D1D, const int D1Dtest, const int Q1D, const int NE, const bool scalarCoeff, const bool trialHcurl, const Array< double > &Bo_, const Array< double > &Bc_, const Array< double > &Bot_, const Array< double > &Bct_, const Vector &op_, const Vector &x_, Vector &y_)
void PAHcurlMassAssembleDiagonal2D(const int D1D, const int Q1D, const int NE, const bool symmetric, const Array< double > &bo, const Array< double > &bc, const Vector &pa_data, Vector &diag)
const int MAX_Q1D
Definition: forall.hpp:29
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:446
constexpr int HCURL_MAX_D1D
Definition: bilininteg.hpp:23
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:304
void PAHdivMassApply3D(const int D1D, const int Q1D, const int NE, const Array< double > &Bo_, const Array< double > &Bc_, const Array< double > &Bot_, const Array< double > &Bct_, const Vector &op_, const Vector &x_, Vector &y_)
void PAHcurlH1ApplyTranspose2D(const int D1D, const int Q1D, const int NE, const Array< double > &bc, const Array< double > &bo, const Array< double > &bct, const Array< double > &gct, const Vector &pa_data, const Vector &x, Vector &y)
void PAHdivSetup2D(const int Q1D, const int NE, const Array< double > &w, const Vector &j, Vector &coeff_, Vector &op)
void PAHcurlMassAssembleDiagonal3D(const int D1D, const int Q1D, const int NE, const bool symmetric, const Array< double > &bo, const Array< double > &bc, const Vector &pa_data, Vector &diag)
void PADiffusionSetup3D(const int Q1D, const int coeffDim, const int NE, const Array< double > &w, const Vector &j, const Vector &c, Vector &d)
const int MAX_D1D
Definition: forall.hpp:28
void PAHcurlMassApply2D(const int D1D, const int Q1D, const int NE, const bool symmetric, const Array< double > &bo, const Array< double > &bc, const Array< double > &bot, const Array< double > &bct, const Vector &pa_data, const Vector &x, Vector &y)
void PAHdivMassAssembleDiagonal2D(const int D1D, const int Q1D, const int NE, const Array< double > &Bo_, const Array< double > &Bc_, const Vector &op_, Vector &diag_)
void SmemPAHcurlMassAssembleDiagonal3D(const int D1D, const int Q1D, const int NE, const bool symmetric, const Array< double > &bo, const Array< double > &bc, const Vector &pa_data, Vector &diag)
void PAHcurlHdivSetup3D(const int Q1D, const int coeffDim, const int NE, const bool transpose, const Array< double > &w_, const Vector &j, Vector &coeff_, Vector &op)
Vector data type.
Definition: vector.hpp:60
void PAHdivMassAssembleDiagonal3D(const int D1D, const int Q1D, const int NE, const Array< double > &Bo_, const Array< double > &Bc_, const Vector &op_, Vector &diag_)
void PAHcurlHdivSetup2D(const int Q1D, const int coeffDim, const int NE, const bool transpose, const Array< double > &w_, const Vector &j, Vector &coeff_, Vector &op)
void PAHcurlH1Apply2D(const int D1D, const int Q1D, const int NE, const Array< double > &bc, const Array< double > &gc, const Array< double > &bot, const Array< double > &bct, const Vector &pa_data, const Vector &x, Vector &y)
constexpr int HCURL_MAX_Q1D
Definition: bilininteg.hpp:25
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:438
void PAHcurlH1ApplyTranspose3D(const int D1D, const int Q1D, const int NE, const Array< double > &bc, const Array< double > &bo, const Array< double > &bct, const Array< double > &gct, const Vector &pa_data, const Vector &x, Vector &y)