MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
lor_nd.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 "lor_nd.hpp"
13 #include "lor_util.hpp"
14 #include "../../linalg/dtensor.hpp"
15 #include "../../general/forall.hpp"
16 
17 namespace mfem
18 {
19 
20 template <int ORDER>
22 {
23  const int nel_ho = fes_ho.GetNE();
24 
25  static constexpr int nv = 4;
26  static constexpr int ne = 4;
27  static constexpr int dim = 2;
28  static constexpr int ddm2 = (dim*(dim+1))/2;
29  static constexpr int ngeom = ddm2 + 1;
30  static constexpr int o = ORDER;
31  static constexpr int op1 = ORDER + 1;
32  static constexpr int ndof_per_el = dim*o*op1;
33  static constexpr int nnz_per_row = 7;
34  static constexpr int sz_local_mat = ne*ne;
35 
36  const bool const_mq = c1.Size() == 1;
37  const auto MQ = const_mq
38  ? Reshape(c1.Read(), 1, 1, 1)
39  : Reshape(c1.Read(), op1, op1, nel_ho);
40  const bool const_dq = c2.Size() == 1;
41  const auto DQ = const_dq
42  ? Reshape(c2.Read(), 1, 1, 1)
43  : Reshape(c2.Read(), op1, op1, nel_ho);
44 
45  sparse_ij.SetSize(nnz_per_row*ndof_per_el*nel_ho);
46  auto V = Reshape(sparse_ij.Write(), nnz_per_row, o*op1, dim, nel_ho);
47 
48  auto X = X_vert.Read();
49 
50  MFEM_FORALL_2D(iel_ho, nel_ho, ORDER, ORDER, 1,
51  {
52  // Assemble a sparse matrix over the macro-element by looping over each
53  // subelement.
54  // V(j,ix,iy) stores the jth nonzero in the row of the sparse matrix
55  // corresponding to local DOF (ix, iy).
56  MFEM_FOREACH_THREAD(iy,y,o)
57  {
58  MFEM_FOREACH_THREAD(ix,x,op1)
59  {
60  for (int c=0; c<2; ++c)
61  {
62  for (int j=0; j<nnz_per_row; ++j)
63  {
64  V(j,ix+iy*op1,c,iel_ho) = 0.0;
65  }
66  }
67  }
68  }
69  MFEM_SYNC_THREAD;
70 
71  // Loop over the sub-elements
72  MFEM_FOREACH_THREAD(ky,y,ORDER)
73  {
74  MFEM_FOREACH_THREAD(kx,x,ORDER)
75  {
76  // Compute geometric factors at quadrature points
77  double Q_[ngeom*nv];
78  double local_mat_[sz_local_mat];
79 
80  DeviceTensor<3> Q(Q_, ngeom, 2, 2);
81  DeviceTensor<2> local_mat(local_mat_, ne, ne);
82 
83  // local_mat is the local (dense) stiffness matrix
84  for (int i=0; i<sz_local_mat; ++i) { local_mat[i] = 0.0; }
85 
86  double vx[4], vy[4];
87  LORVertexCoordinates2D<ORDER>(X, iel_ho, kx, ky, vx, vy);
88 
89  for (int iqx=0; iqx<2; ++iqx)
90  {
91  for (int iqy=0; iqy<2; ++iqy)
92  {
93  const double x = iqx;
94  const double y = iqy;
95  const double w = 1.0/4.0;
96 
97  double J_[2*2];
98  DeviceTensor<2> J(J_, 2, 2);
99 
100  Jacobian2D(x, y, vx, vy, J);
101 
102  const double detJ = Det2D(J);
103  const double w_detJ = w/detJ;
104 
105  Q(0,iqy,iqx) = w_detJ * (J(0,1)*J(0,1) + J(1,1)*J(1,1)); // 1,1
106  Q(1,iqy,iqx) = -w_detJ * (J(0,1)*J(0,0) + J(1,1)*J(1,0)); // 1,2
107  Q(2,iqy,iqx) = w_detJ * (J(0,0)*J(0,0) + J(1,0)*J(1,0)); // 2,2
108  Q(3,iqy,iqx) = w_detJ;
109  }
110  }
111  for (int iqx=0; iqx<2; ++iqx)
112  {
113  for (int iqy=0; iqy<2; ++iqy)
114  {
115  const double mq = const_mq ? MQ(0,0,0) : MQ(kx+iqx, ky+iqy, iel_ho);
116  const double dq = const_dq ? DQ(0,0,0) : DQ(kx+iqx, ky+iqy, iel_ho);
117  // Loop over x,y components. c=0 => x, c=1 => y
118  for (int cj=0; cj<dim; ++cj)
119  {
120  for (int bj=0; bj<2; ++bj)
121  {
122  const double curl_j = ((cj == 0) ? 1 : -1)*((bj == 0) ? 1 : -1);
123  const double bxj = (cj == 0) ? ((bj == iqy) ? 1 : 0) : 0;
124  const double byj = (cj == 1) ? ((bj == iqx) ? 1 : 0) : 0;
125 
126  const double jj_loc = bj + 2*cj;
127 
128  for (int ci=0; ci<dim; ++ci)
129  {
130  for (int bi=0; bi<2; ++bi)
131  {
132  const double curl_i = ((ci == 0) ? 1 : -1)*((bi == 0) ? 1 : -1);
133  const double bxi = (ci == 0) ? ((bi == iqy) ? 1 : 0) : 0;
134  const double byi = (ci == 1) ? ((bi == iqx) ? 1 : 0) : 0;
135 
136  const double ii_loc = bi + 2*ci;
137 
138  // Only store the lower-triangular part of
139  // the matrix (by symmetry).
140  if (jj_loc > ii_loc) { continue; }
141 
142  double val = 0.0;
143  val += bxi*bxj*Q(0,iqy,iqx);
144  val += byi*bxj*Q(1,iqy,iqx);
145  val += bxi*byj*Q(1,iqy,iqx);
146  val += byi*byj*Q(2,iqy,iqx);
147  val *= mq;
148  val += dq*curl_i*curl_j*Q(3,iqy,iqx);
149 
150  local_mat(ii_loc, jj_loc) += val;
151  }
152  }
153  }
154  }
155  }
156  }
157  // Assemble the local matrix into the macro-element sparse matrix
158  // in a format similar to coordinate format. The (I,J) arrays
159  // are implicit (not stored explicitly).
160  for (int ii_loc=0; ii_loc<ne; ++ii_loc)
161  {
162  const int ci = ii_loc/2;
163  const int bi = ii_loc%2;
164  const int ix = (ci == 0) ? 0 : bi;
165  const int iy = (ci == 1) ? 0 : bi;
166 
167  int ii = (ci == 0) ? kx+ix + (ky+iy)*o : kx+ix + (ky+iy)*op1;
168 
169  for (int jj_loc=0; jj_loc<ne; ++jj_loc)
170  {
171  const int cj = jj_loc/2;
172  const int bj = jj_loc%2;
173 
174  const int jj_off = (ci == cj) ? (bj - bi + 1) : (3 + bj + (1-bi)*2);
175 
176  // Symmetry
177  const double val = (jj_loc <= ii_loc)
178  ? local_mat(ii_loc, jj_loc)
179  : local_mat(jj_loc, ii_loc);
180  AtomicAdd(V(jj_off, ii, ci, iel_ho), val);
181  }
182  }
183  }
184  }
185  });
186 
187  sparse_mapping.SetSize(nnz_per_row*ndof_per_el);
188  sparse_mapping = -1;
189  auto map = Reshape(sparse_mapping.HostReadWrite(), nnz_per_row, ndof_per_el);
190  for (int ci=0; ci<2; ++ci)
191  {
192  for (int i1=0; i1<o; ++i1)
193  {
194  for (int i2=0; i2<op1; ++i2)
195  {
196  const int ii_el = (ci == 0) ? i1 + i2*o : i2 + i1*op1 + o*op1;
197  for (int cj=0; cj<2; ++cj)
198  {
199  const int j1_begin = (ci == cj) ? i1 : ((i2 > 0) ? i2-1 : i2);
200  const int j1_end = (ci == cj) ? i1 : ((i2 < o) ? i2 : i2-1);
201  const int j2_begin = (ci == cj) ? ((i2 > 0) ? i2-1 : i2) : i1;
202  const int j2_end = (ci == cj) ? ((i2 < o) ? i2+1 : i2) : i1+1;
203 
204  for (int j1=j1_begin; j1<=j1_end; ++j1)
205  {
206  for (int j2=j2_begin; j2<=j2_end; ++j2)
207  {
208  const int jj_el = (cj == 0) ? j1 + j2*o : j2 + j1*op1 + o*op1;
209  int jj_off = (ci == cj) ? (j2-i2+1) : 3 + (j2-i1) + 2*(j1-i2+1);
210  map(jj_off, ii_el) = jj_el;
211  }
212  }
213  }
214  }
215  }
216  }
217 }
218 
219 template <int ORDER>
221 {
222  const int nel_ho = fes_ho.GetNE();
223 
224  static constexpr int nv = 8; // number of vertices in hexahedron
225  static constexpr int ne = 12; // number of edges in hexahedron
226  static constexpr int dim = 3;
227  static constexpr int ddm2 = (dim*(dim+1))/2;
228  static constexpr int ngeom = 2*ddm2; // number of geometric factors stored
229  static constexpr int o = ORDER;
230  static constexpr int op1 = ORDER + 1;
231  static constexpr int ndof_per_el = dim*o*op1*op1;
232  static constexpr int nnz_per_row = 33;
233  static constexpr int sz_local_mat = ne*ne;
234 
235  const bool const_mq = c1.Size() == 1;
236  const auto MQ = const_mq
237  ? Reshape(c1.Read(), 1, 1, 1, 1)
238  : Reshape(c1.Read(), op1, op1, op1, nel_ho);
239  const bool const_dq = c2.Size() == 1;
240  const auto DQ = const_dq
241  ? Reshape(c2.Read(), 1, 1, 1, 1)
242  : Reshape(c2.Read(), op1, op1, op1, nel_ho);
243 
244  sparse_ij.SetSize(nnz_per_row*ndof_per_el*nel_ho);
245  auto V = Reshape(sparse_ij.Write(), nnz_per_row, o*op1*op1, dim, nel_ho);
246 
247  auto X = X_vert.Read();
248 
249  // Last thread dimension is lowered to avoid "too many resources" error
250  MFEM_FORALL_3D(iel_ho, nel_ho, ORDER, ORDER, (ORDER>6)?4:ORDER,
251  {
252  MFEM_FOREACH_THREAD(iz,z,o)
253  {
254  MFEM_FOREACH_THREAD(iy,y,op1)
255  {
256  MFEM_FOREACH_THREAD(ix,x,op1)
257  {
258  for (int c=0; c<dim; ++c)
259  {
260  for (int j=0; j<nnz_per_row; ++j)
261  {
262  V(j,ix+iy*op1+iz*op1*op1,c,iel_ho) = 0.0;
263  }
264  }
265  }
266  }
267  }
268  MFEM_SYNC_THREAD;
269 
270  // Loop over the sub-elements
271  MFEM_FOREACH_THREAD(kz,z,ORDER)
272  {
273  MFEM_FOREACH_THREAD(ky,y,ORDER)
274  {
275  MFEM_FOREACH_THREAD(kx,x,ORDER)
276  {
277  // Geometric factors at quadrature points (element vertices)
278  double Q_[ngeom*nv];
279  DeviceTensor<4> Q(Q_, ngeom, 2, 2, 2);
280 
281  double local_mat_[sz_local_mat];
282  DeviceTensor<2> local_mat(local_mat_, ne, ne);
283  for (int i=0; i<sz_local_mat; ++i) { local_mat[i] = 0.0; }
284 
285  double vx[8], vy[8], vz[8];
286  LORVertexCoordinates3D<ORDER>(X, iel_ho, kx, ky, kz, vx, vy, vz);
287 
288  for (int iqz=0; iqz<2; ++iqz)
289  {
290  for (int iqy=0; iqy<2; ++iqy)
291  {
292  for (int iqx=0; iqx<2; ++iqx)
293  {
294  const double x = iqx;
295  const double y = iqy;
296  const double z = iqz;
297  const double w = 1.0/8.0;
298 
299  double J_[3*3];
300  DeviceTensor<2> J(J_, 3, 3);
301 
302  Jacobian3D(x, y, z, vx, vy, vz, J);
303 
304  const double detJ = Det3D(J);
305  const double w_detJ = w/detJ;
306 
307  // adj(J)
308  double A_[3*3];
309  DeviceTensor<2> A(A_, 3, 3);
310  Adjugate3D(J, A);
311 
312  Q(0,iqz,iqy,iqx) = w_detJ*(A(0,0)*A(0,0)+A(0,1)*A(0,1)+A(0,2)*A(0,2)); // 1,1
313  Q(1,iqz,iqy,iqx) = w_detJ*(A(0,0)*A(1,0)+A(0,1)*A(1,1)+A(0,2)*A(1,2)); // 2,1
314  Q(2,iqz,iqy,iqx) = w_detJ*(A(0,0)*A(2,0)+A(0,1)*A(2,1)+A(0,2)*A(2,2)); // 3,1
315  Q(3,iqz,iqy,iqx) = w_detJ*(A(1,0)*A(1,0)+A(1,1)*A(1,1)+A(1,2)*A(1,2)); // 2,2
316  Q(4,iqz,iqy,iqx) = w_detJ*(A(1,0)*A(2,0)+A(1,1)*A(2,1)+A(1,2)*A(2,2)); // 3,2
317  Q(5,iqz,iqy,iqx) = w_detJ*(A(2,0)*A(2,0)+A(2,1)*A(2,1)+A(2,2)*A(2,2)); // 3,3
318 
319  // w J^T J / det(J)
320  Q(6,iqz,iqy,iqx) = w_detJ*(J(0,0)*J(0,0)+J(1,0)*J(1,0)+J(2,0)*J(2,0)); // 1,1
321  Q(7,iqz,iqy,iqx) = w_detJ*(J(0,0)*J(0,1)+J(1,0)*J(1,1)+J(2,0)*J(2,1)); // 2,1
322  Q(8,iqz,iqy,iqx) = w_detJ*(J(0,0)*J(0,2)+J(1,0)*J(1,2)+J(2,0)*J(2,2)); // 3,1
323  Q(9,iqz,iqy,iqx) = w_detJ*(J(0,1)*J(0,1)+J(1,1)*J(1,1)+J(2,1)*J(2,1)); // 2,2
324  Q(10,iqz,iqy,iqx) = w_detJ*(J(0,1)*J(0,2)+J(1,1)*J(1,2)+J(2,1)*J(2,2)); // 3,2
325  Q(11,iqz,iqy,iqx) = w_detJ*(J(0,2)*J(0,2)+J(1,2)*J(1,2)+J(2,2)*J(2,2)); // 3,3
326  }
327  }
328  }
329  for (int iqz=0; iqz<2; ++iqz)
330  {
331  for (int iqy=0; iqy<2; ++iqy)
332  {
333  for (int iqx=0; iqx<2; ++iqx)
334  {
335  const double mq = const_mq ? MQ(0,0,0,0) : MQ(kx+iqx, ky+iqy, kz+iqz, iel_ho);
336  const double dq = const_dq ? DQ(0,0,0,0) : DQ(kx+iqx, ky+iqy, kz+iqz, iel_ho);
337  // Loop over x,y,z components. 0 => x, 1 => y, 2 => z
338  for (int cj=0; cj<dim; ++cj)
339  {
340  const double jq1 = (cj == 0) ? iqy : ((cj == 1) ? iqz : iqx);
341  const double jq2 = (cj == 0) ? iqz : ((cj == 1) ? iqx : iqy);
342 
343  const int jd_0 = cj;
344  const int jd_1 = (cj + 1)%3;
345  const int jd_2 = (cj + 2)%3;
346 
347  for (int bj=0; bj<4; ++bj) // 4 edges in each dim
348  {
349  const int bj1 = bj%2;
350  const int bj2 = bj/2;
351 
352  double curl_j[3];
353  curl_j[jd_0] = 0.0;
354  curl_j[jd_1] = ((bj1 == 0) ? jq1 - 1 : -jq1)*((bj2 == 0) ? 1 : -1);
355  curl_j[jd_2] = ((bj2 == 0) ? 1 - jq2 : jq2)*((bj1 == 0) ? 1 : -1);
356 
357  double basis_j[3];
358  basis_j[jd_0] = ((bj1 == 0) ? 1 - jq1 : jq1)*((bj2 == 0) ? 1 - jq2 : jq2);
359  basis_j[jd_1] = 0.0;
360  basis_j[jd_2] = 0.0;
361 
362  const int jj_loc = bj + 4*cj;
363 
364  for (int ci=0; ci<dim; ++ci)
365  {
366  const double iq1 = (ci == 0) ? iqy : ((ci == 1) ? iqz : iqx);
367  const double iq2 = (ci == 0) ? iqz : ((ci == 1) ? iqx : iqy);
368 
369  const int id_0 = ci;
370  const int id_1 = (ci + 1)%3;
371  const int id_2 = (ci + 2)%3;
372 
373  for (int bi=0; bi<4; ++bi)
374  {
375  const int bi1 = bi%2;
376  const int bi2 = bi/2;
377 
378  double curl_i[3];
379  curl_i[id_0] = 0.0;
380  curl_i[id_1] = ((bi1 == 0) ? iq1 - 1 : -iq1)*((bi2 == 0) ? 1 : -1);
381  curl_i[id_2] = ((bi2 == 0) ? 1 - iq2 : iq2)*((bi1 == 0) ? 1 : -1);
382 
383  double basis_i[3];
384  basis_i[id_0] = ((bi1 == 0) ? 1 - iq1 : iq1)*((bi2 == 0) ? 1 - iq2 : iq2);
385  basis_i[id_1] = 0.0;
386  basis_i[id_2] = 0.0;
387 
388  const int ii_loc = bi + 4*ci;
389 
390  // Only store the lower-triangular part of
391  // the matrix (by symmetry).
392  if (jj_loc > ii_loc) { continue; }
393 
394  double curl_curl = 0.0;
395  curl_curl += Q(6,iqz,iqy,iqx)*curl_i[0]*curl_j[0];
396  curl_curl += Q(7,iqz,iqy,iqx)*(curl_i[0]*curl_j[1] + curl_i[1]*curl_j[0]);
397  curl_curl += Q(8,iqz,iqy,iqx)*(curl_i[0]*curl_j[2] + curl_i[2]*curl_j[0]);
398  curl_curl += Q(9,iqz,iqy,iqx)*curl_i[1]*curl_j[1];
399  curl_curl += Q(10,iqz,iqy,iqx)*(curl_i[1]*curl_j[2] + curl_i[2]*curl_j[1]);
400  curl_curl += Q(11,iqz,iqy,iqx)*curl_i[2]*curl_j[2];
401 
402  double basis_basis = 0.0;
403  basis_basis += Q(0,iqz,iqy,iqx)*basis_i[0]*basis_j[0];
404  basis_basis += Q(1,iqz,iqy,iqx)*(basis_i[0]*basis_j[1] + basis_i[1]*basis_j[0]);
405  basis_basis += Q(2,iqz,iqy,iqx)*(basis_i[0]*basis_j[2] + basis_i[2]*basis_j[0]);
406  basis_basis += Q(3,iqz,iqy,iqx)*basis_i[1]*basis_j[1];
407  basis_basis += Q(4,iqz,iqy,iqx)*(basis_i[1]*basis_j[2] + basis_i[2]*basis_j[1]);
408  basis_basis += Q(5,iqz,iqy,iqx)*basis_i[2]*basis_j[2];
409 
410  const double val = dq*curl_curl + mq*basis_basis;
411 
412  local_mat(ii_loc, jj_loc) += val;
413  }
414  }
415  }
416  }
417  }
418  }
419  }
420  // Assemble the local matrix into the macro-element sparse matrix
421  // The nonzeros of the macro-element sparse matrix are ordered as
422  // follows:
423  //
424  // The axes are ordered relative to the direction of the basis
425  // vector, e.g. for x-vectors, the axes are (x,y,z), for
426  // y-vectors the axes are (y,z,x), and for z-vectors the axes are
427  // (z,x,y).
428  //
429  // The nonzeros are then given in "rotated lexicographic"
430  // ordering, according to these axes.
431  for (int ii_loc=0; ii_loc<ne; ++ii_loc)
432  {
433  const int ci = ii_loc/4;
434  const int bi = ii_loc%4;
435 
436  const int id0 = ci;
437  const int id1 = (ci+1)%3;
438  const int id2 = (ci+2)%3;
439 
440  const int i0 = 0;
441  const int i1 = bi%2;
442  const int i2 = bi/2;
443 
444  int ii_lex[3];
445  ii_lex[id0] = i0;
446  ii_lex[id1] = i1;
447  ii_lex[id2] = i2;
448 
449  const int nx = (ci == 0) ? o : op1;
450  const int ny = (ci == 1) ? o : op1;
451 
452  const int ii = kx+ii_lex[0] + (ky+ii_lex[1])*nx + (kz+ii_lex[2])*nx*ny;
453 
454  for (int jj_loc=0; jj_loc<ne; ++jj_loc)
455  {
456  const int cj = jj_loc/4;
457  // add 3 to take modulus (rather than remainder) when
458  // (cj - ci) is negative
459  const int cj_rel = (3 + cj - ci)%3;
460 
461  const int bj = jj_loc%4;
462 
463  const int jd0 = cj_rel;
464  const int jd1 = (cj_rel+1)%3;
465  const int jd2 = (cj_rel+2)%3;
466 
467  int jj_rel[3];
468  jj_rel[jd0] = 0;
469  jj_rel[jd1] = bj%2;
470  jj_rel[jd2] = bj/2;
471 
472  const int d0 = jj_rel[0] - i0;
473  const int d1 = 1 + jj_rel[1] - i1;
474  const int d2 = 1 + jj_rel[2] - i2;
475  int jj_off;
476  if (cj_rel == 0) { jj_off = d1 + 3*d2; }
477  else if (cj_rel == 1) { jj_off = 9 + d0 + 2*d1 + 4*d2; }
478  else /* if (cj_rel == 2) */ { jj_off = 21 + d0 + 2*d1 + 6*d2; }
479 
480  // Symmetry
481  const double val = (jj_loc <= ii_loc)
482  ? local_mat(ii_loc, jj_loc)
483  : local_mat(jj_loc, ii_loc);
484  AtomicAdd(V(jj_off, ii, ci, iel_ho), val);
485  }
486  }
487  }
488  }
489  }
490  });
491 
492  sparse_mapping.SetSize(nnz_per_row*ndof_per_el);
493  sparse_mapping = -1;
494  auto map = Reshape(sparse_mapping.HostReadWrite(), nnz_per_row, ndof_per_el);
495  for (int ci=0; ci<dim; ++ci)
496  {
497  const int i_off = ci*o*op1*op1;
498  const int id0 = ci;
499  const int id1 = (ci+1)%3;
500  const int id2 = (ci+2)%3;
501 
502  const int nxi = (ci == 0) ? o : op1;
503  const int nyi = (ci == 1) ? o : op1;
504 
505  for (int i0=0; i0<o; ++i0)
506  {
507  for (int i1=0; i1<op1; ++i1)
508  {
509  for (int i2=0; i2<op1; ++i2)
510  {
511  int ii_lex[3];
512  ii_lex[id0] = i0;
513  ii_lex[id1] = i1;
514  ii_lex[id2] = i2;
515  const int ii_el = i_off + ii_lex[0] + ii_lex[1]*nxi + ii_lex[2]*nxi*nyi;
516 
517  for (int cj_rel=0; cj_rel<dim; ++cj_rel)
518  {
519  const int cj = (ci + cj_rel) % 3;
520  const int j_off = cj*o*op1*op1;
521 
522  const int nxj = (cj == 0) ? o : op1;
523  const int nyj = (cj == 1) ? o : op1;
524 
525  const int j0_begin = i0;
526  const int j0_end = (cj_rel == 0) ? i0 : i0 + 1;
527  const int j1_begin = (i1 > 0) ? i1-1 : i1;
528  const int j1_end = (cj_rel == 1)
529  ? ((i1 < o) ? i1 : i1-1)
530  : ((i1 < o) ? i1+1 : i1);
531  const int j2_begin = (i2 > 0) ? i2-1 : i2;
532  const int j2_end = (cj_rel == 2)
533  ? ((i2 < o) ? i2 : i2-1)
534  : ((i2 < o) ? i2+1 : i2);
535 
536  for (int j0=j0_begin; j0<=j0_end; ++j0)
537  {
538  const int d0 = j0 - i0;
539  for (int j1=j1_begin; j1<=j1_end; ++j1)
540  {
541  const int d1 = j1 - i1 + 1;
542  for (int j2=j2_begin; j2<=j2_end; ++j2)
543  {
544  const int d2 = j2 - i2 + 1;
545  int jj_lex[3];
546  jj_lex[id0] = j0;
547  jj_lex[id1] = j1;
548  jj_lex[id2] = j2;
549  const int jj_el = j_off + jj_lex[0] + jj_lex[1]*nxj + jj_lex[2]*nxj*nyj;
550  int jj_off;
551  if (cj_rel == 0) { jj_off = d1 + 3*d2; }
552  else if (cj_rel == 1) { jj_off = 9 + d0 + 2*d1 + 4*d2; }
553  else /* if (cj_rel == 2) */ { jj_off = 21 + d0 + 2*d1 + 6*d2; }
554  map(jj_off, ii_el) = jj_el;
555  }
556  }
557  }
558  }
559  }
560  }
561  }
562  }
563 }
564 
565 // Explicit template instantiations
566 template void BatchedLOR_ND::Assemble2D<1>();
567 template void BatchedLOR_ND::Assemble2D<2>();
568 template void BatchedLOR_ND::Assemble2D<3>();
569 template void BatchedLOR_ND::Assemble2D<4>();
570 template void BatchedLOR_ND::Assemble2D<5>();
571 template void BatchedLOR_ND::Assemble2D<6>();
572 template void BatchedLOR_ND::Assemble2D<7>();
573 template void BatchedLOR_ND::Assemble2D<8>();
574 
575 template void BatchedLOR_ND::Assemble3D<1>();
576 template void BatchedLOR_ND::Assemble3D<2>();
577 template void BatchedLOR_ND::Assemble3D<3>();
578 template void BatchedLOR_ND::Assemble3D<4>();
579 template void BatchedLOR_ND::Assemble3D<5>();
580 template void BatchedLOR_ND::Assemble3D<6>();
581 template void BatchedLOR_ND::Assemble3D<7>();
582 template void BatchedLOR_ND::Assemble3D<8>();
583 
585  FiniteElementSpace &fes_ho_,
586  Vector &X_vert_,
587  Vector &sparse_ij_,
588  Array<int> &sparse_mapping_)
589  : BatchedLORKernel(fes_ho_, X_vert_, sparse_ij_, sparse_mapping_)
590 {
591  ProjectLORCoefficient<VectorFEMassIntegrator>(a, c1);
592  ProjectLORCoefficient<CurlCurlIntegrator>(a, c2);
593 }
594 
595 } // namespace mfem
T * HostReadWrite()
Shortcut for mfem::ReadWrite(a.GetMemory(), a.Size(), false).
Definition: array.hpp:324
Abstract base class for the batched LOR assembly kernels.
BatchedLOR_ND(BilinearForm &a, FiniteElementSpace &fes_ho_, Vector &X_vert_, Vector &sparse_ij_, Array< int > &sparse_mapping_)
Definition: lor_nd.cpp:584
Array< int > & sparse_mapping
Local element sparsity pattern.
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
MFEM_HOST_DEVICE T AtomicAdd(T &add, const T val)
Definition: backends.hpp:84
int Size() const
Returns the size of the vector.
Definition: vector.hpp:200
CoefficientVector c2
Coefficient of second integrator.
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:614
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:457
Vector & sparse_ij
Local element sparsity matrix data.
MFEM_HOST_DEVICE double Det3D(DeviceMatrix &J)
Definition: lor_util.hpp:188
MFEM_HOST_DEVICE double Det2D(DeviceMatrix &J)
Definition: lor_util.hpp:183
MFEM_HOST_DEVICE void Jacobian3D(const double x, const double y, const double z, const double vx[8], const double vy[8], const double vz[8], DeviceMatrix &J)
Definition: lor_util.hpp:126
FiniteElementSpace & fes_ho
The associated high-order space.
A basic generic Tensor class, appropriate for use on the GPU.
Definition: dtensor.hpp:81
MFEM_HOST_DEVICE void Jacobian2D(const double x, const double y, const double vx[4], const double vy[4], DeviceMatrix &J)
Definition: lor_util.hpp:115
Vector & X_vert
Mesh coordinate vector.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
double a
Definition: lissajous.cpp:41
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
int dim
Definition: ex24.cpp:53
Vector data type.
Definition: vector.hpp:60
MFEM_HOST_DEVICE void Adjugate3D(const DeviceMatrix &J, DeviceMatrix &A)
Definition: lor_util.hpp:170
CoefficientVector c1
Coefficient of first integrator.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:449
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition: dtensor.hpp:131