MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
restriction.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 "restriction.hpp"
13 #include "gridfunc.hpp"
14 #include "fespace.hpp"
15 #include "../general/forall.hpp"
16 
17 namespace mfem
18 {
19 
21  ElementDofOrdering e_ordering)
22  : fes(f),
23  ne(fes.GetNE()),
24  vdim(fes.GetVDim()),
25  byvdim(fes.GetOrdering() == Ordering::byVDIM),
26  ndofs(fes.GetNDofs()),
27  dof(ne > 0 ? fes.GetFE(0)->GetDof() : 0),
28  nedofs(ne*dof),
29  offsets(ndofs+1),
30  indices(ne*dof),
31  gatherMap(ne*dof)
32 {
33  // Assuming all finite elements are the same.
34  height = vdim*ne*dof;
35  width = fes.GetVSize();
36  const bool dof_reorder = (e_ordering == ElementDofOrdering::LEXICOGRAPHIC);
37  const int *dof_map = NULL;
38  if (dof_reorder && ne > 0)
39  {
40  for (int e = 0; e < ne; ++e)
41  {
42  const FiniteElement *fe = fes.GetFE(e);
43  const TensorBasisElement* el =
44  dynamic_cast<const TensorBasisElement*>(fe);
45  if (el) { continue; }
46  mfem_error("Finite element not suitable for lexicographic ordering");
47  }
48  const FiniteElement *fe = fes.GetFE(0);
49  const TensorBasisElement* el =
50  dynamic_cast<const TensorBasisElement*>(fe);
51  const Array<int> &fe_dof_map = el->GetDofMap();
52  MFEM_VERIFY(fe_dof_map.Size() > 0, "invalid dof map");
53  dof_map = fe_dof_map.GetData();
54  }
55  const Table& e2dTable = fes.GetElementToDofTable();
56  const int* elementMap = e2dTable.GetJ();
57  // We will be keeping a count of how many local nodes point to its global dof
58  for (int i = 0; i <= ndofs; ++i)
59  {
60  offsets[i] = 0;
61  }
62  for (int e = 0; e < ne; ++e)
63  {
64  for (int d = 0; d < dof; ++d)
65  {
66  const int sgid = elementMap[dof*e + d]; // signed
67  const int gid = (sgid >= 0) ? sgid : -1 - sgid;
68  ++offsets[gid + 1];
69  }
70  }
71  // Aggregate to find offsets for each global dof
72  for (int i = 1; i <= ndofs; ++i)
73  {
74  offsets[i] += offsets[i - 1];
75  }
76  // For each global dof, fill in all local nodes that point to it
77  for (int e = 0; e < ne; ++e)
78  {
79  for (int d = 0; d < dof; ++d)
80  {
81  const int sdid = dof_reorder ? dof_map[d] : 0; // signed
82  const int did = (!dof_reorder)?d:(sdid >= 0 ? sdid : -1-sdid);
83  const int sgid = elementMap[dof*e + did]; // signed
84  const int gid = (sgid >= 0) ? sgid : -1-sgid;
85  const int lid = dof*e + d;
86  const bool plus = (sgid >= 0 && sdid >= 0) || (sgid < 0 && sdid < 0);
87  gatherMap[lid] = plus ? gid : -1-gid;
88  indices[offsets[gid]++] = plus ? lid : -1-lid;
89  }
90  }
91  // We shifted the offsets vector by 1 by using it as a counter.
92  // Now we shift it back.
93  for (int i = ndofs; i > 0; --i)
94  {
95  offsets[i] = offsets[i - 1];
96  }
97  offsets[0] = 0;
98 }
99 
100 void ElementRestriction::Mult(const Vector& x, Vector& y) const
101 {
102  // Assumes all elements have the same number of dofs
103  const int nd = dof;
104  const int vd = vdim;
105  const bool t = byvdim;
106  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
107  auto d_y = Reshape(y.Write(), nd, vd, ne);
108  auto d_gatherMap = gatherMap.Read();
109  MFEM_FORALL(i, dof*ne,
110  {
111  const int gid = d_gatherMap[i];
112  const bool plus = gid >= 0;
113  const int j = plus ? gid : -1-gid;
114  for (int c = 0; c < vd; ++c)
115  {
116  const double dofValue = d_x(t?c:j, t?j:c);
117  d_y(i % nd, c, i / nd) = plus ? dofValue : -dofValue;
118  }
119  });
120 }
121 
123 {
124  // Assumes all elements have the same number of dofs
125  const int nd = dof;
126  const int vd = vdim;
127  const bool t = byvdim;
128  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
129  auto d_y = Reshape(y.Write(), nd, vd, ne);
130  auto d_gatherMap = gatherMap.Read();
131 
132  MFEM_FORALL(i, dof*ne,
133  {
134  const int gid = d_gatherMap[i];
135  const int j = gid >= 0 ? gid : -1-gid;
136  for (int c = 0; c < vd; ++c)
137  {
138  d_y(i % nd, c, i / nd) = d_x(t?c:j, t?j:c);
139  }
140  });
141 }
142 
144 {
145  // Assumes all elements have the same number of dofs
146  const int nd = dof;
147  const int vd = vdim;
148  const bool t = byvdim;
149  auto d_offsets = offsets.Read();
150  auto d_indices = indices.Read();
151  auto d_x = Reshape(x.Read(), nd, vd, ne);
152  auto d_y = Reshape(y.Write(), t?vd:ndofs, t?ndofs:vd);
153  MFEM_FORALL(i, ndofs,
154  {
155  const int offset = d_offsets[i];
156  const int nextOffset = d_offsets[i + 1];
157  for (int c = 0; c < vd; ++c)
158  {
159  double dofValue = 0;
160  for (int j = offset; j < nextOffset; ++j)
161  {
162  const int idx_j = (d_indices[j] >= 0) ? d_indices[j] : -1 - d_indices[j];
163  dofValue += (d_indices[j] >= 0) ? d_x(idx_j % nd, c,
164  idx_j / nd) : -d_x(idx_j % nd, c, idx_j / nd);
165  }
166  d_y(t?c:i,t?i:c) = dofValue;
167  }
168  });
169 }
170 
172 {
173  // Assumes all elements have the same number of dofs
174  const int nd = dof;
175  const int vd = vdim;
176  const bool t = byvdim;
177  auto d_offsets = offsets.Read();
178  auto d_indices = indices.Read();
179  auto d_x = Reshape(x.Read(), nd, vd, ne);
180  auto d_y = Reshape(y.Write(), t?vd:ndofs, t?ndofs:vd);
181  MFEM_FORALL(i, ndofs,
182  {
183  const int offset = d_offsets[i];
184  const int nextOffset = d_offsets[i + 1];
185  for (int c = 0; c < vd; ++c)
186  {
187  double dofValue = 0;
188  for (int j = offset; j < nextOffset; ++j)
189  {
190  const int idx_j = (d_indices[j] >= 0) ? d_indices[j] : -1 - d_indices[j];
191  dofValue += d_x(idx_j % nd, c, idx_j / nd);
192  }
193  d_y(t?c:i,t?i:c) = dofValue;
194  }
195  });
196 }
197 
199 {
200  // Assumes all elements have the same number of dofs
201  const int nd = dof;
202  const int vd = vdim;
203  const bool t = byvdim;
204 
205  Array<char> processed(vd * ndofs);
206  processed = 0;
207 
208  auto d_offsets = offsets.HostRead();
209  auto d_indices = indices.HostRead();
210  auto d_x = Reshape(processed.HostReadWrite(), t?vd:ndofs, t?ndofs:vd);
211  auto d_y = Reshape(y.HostWrite(), nd, vd, ne);
212  for (int i = 0; i < ndofs; ++i)
213  {
214  const int offset = d_offsets[i];
215  const int nextOffset = d_offsets[i+1];
216  for (int c = 0; c < vd; ++c)
217  {
218  for (int j = offset; j < nextOffset; ++j)
219  {
220  const int idx_j = d_indices[j];
221  if (d_x(t?c:i,t?i:c))
222  {
223  d_y(idx_j % nd, c, idx_j / nd) = 0.0;
224  }
225  else
226  {
227  d_y(idx_j % nd, c, idx_j / nd) = 1.0;
228  d_x(t?c:i,t?i:c) = 1;
229  }
230  }
231  }
232  }
233 }
234 
236  SparseMatrix &mat) const
237 {
238  mat.GetMemoryI().New(mat.Height()+1, mat.GetMemoryI().GetMemoryType());
239  const int nnz = FillI(mat);
240  mat.GetMemoryJ().New(nnz, mat.GetMemoryJ().GetMemoryType());
241  mat.GetMemoryData().New(nnz, mat.GetMemoryData().GetMemoryType());
242  FillJAndData(mat_ea, mat);
243 }
244 
245 template <int MaxNbNbr>
246 static MFEM_HOST_DEVICE int GetMinElt(const int *my_elts, const int nbElts,
247  const int *nbr_elts, const int nbrNbElts)
248 {
249  // Building the intersection
250  int inter[MaxNbNbr];
251  int cpt = 0;
252  for (int i = 0; i < nbElts; i++)
253  {
254  const int e_i = my_elts[i];
255  for (int j = 0; j < nbrNbElts; j++)
256  {
257  if (e_i==nbr_elts[j])
258  {
259  inter[cpt] = e_i;
260  cpt++;
261  }
262  }
263  }
264  // Finding the minimum
265  int min = inter[0];
266  for (int i = 1; i < cpt; i++)
267  {
268  if (inter[i] < min)
269  {
270  min = inter[i];
271  }
272  }
273  return min;
274 }
275 
276 /** Returns the index where a non-zero entry should be added and increment the
277  number of non-zeros for the row i_L. */
278 static MFEM_HOST_DEVICE int GetAndIncrementNnzIndex(const int i_L, int* I)
279 {
280  int ind = AtomicAdd(I[i_L],1);
281  return ind;
282 }
283 
285 {
286  static constexpr int Max = MaxNbNbr;
287  const int all_dofs = ndofs;
288  const int vd = vdim;
289  const int elt_dofs = dof;
290  auto I = mat.ReadWriteI();
291  auto d_offsets = offsets.Read();
292  auto d_indices = indices.Read();
293  auto d_gatherMap = gatherMap.Read();
294  MFEM_FORALL(i_L, vd*all_dofs+1,
295  {
296  I[i_L] = 0;
297  });
298  MFEM_FORALL(e, ne,
299  {
300  for (int i = 0; i < elt_dofs; i++)
301  {
302  int i_elts[Max];
303  const int i_E = e*elt_dofs + i;
304  const int i_L = d_gatherMap[i_E];
305  const int i_offset = d_offsets[i_L];
306  const int i_nextOffset = d_offsets[i_L+1];
307  const int i_nbElts = i_nextOffset - i_offset;
308  for (int e_i = 0; e_i < i_nbElts; ++e_i)
309  {
310  const int i_E = d_indices[i_offset+e_i];
311  i_elts[e_i] = i_E/elt_dofs;
312  }
313  for (int j = 0; j < elt_dofs; j++)
314  {
315  const int j_E = e*elt_dofs + j;
316  const int j_L = d_gatherMap[j_E];
317  const int j_offset = d_offsets[j_L];
318  const int j_nextOffset = d_offsets[j_L+1];
319  const int j_nbElts = j_nextOffset - j_offset;
320  if (i_nbElts == 1 || j_nbElts == 1) // no assembly required
321  {
322  GetAndIncrementNnzIndex(i_L, I);
323  }
324  else // assembly required
325  {
326  int j_elts[Max];
327  for (int e_j = 0; e_j < j_nbElts; ++e_j)
328  {
329  const int j_E = d_indices[j_offset+e_j];
330  const int elt = j_E/elt_dofs;
331  j_elts[e_j] = elt;
332  }
333  int min_e = GetMinElt<Max>(i_elts, i_nbElts, j_elts, j_nbElts);
334  if (e == min_e) // add the nnz only once
335  {
336  GetAndIncrementNnzIndex(i_L, I);
337  }
338  }
339  }
340  }
341  });
342  // We need to sum the entries of I, we do it on CPU as it is very sequential.
343  auto h_I = mat.HostReadWriteI();
344  const int nTdofs = vd*all_dofs;
345  int sum = 0;
346  for (int i = 0; i < nTdofs; i++)
347  {
348  const int nnz = h_I[i];
349  h_I[i] = sum;
350  sum+=nnz;
351  }
352  h_I[nTdofs] = sum;
353  // We return the number of nnz
354  return h_I[nTdofs];
355 }
356 
358  SparseMatrix &mat) const
359 {
360  static constexpr int Max = MaxNbNbr;
361  const int all_dofs = ndofs;
362  const int vd = vdim;
363  const int elt_dofs = dof;
364  auto I = mat.ReadWriteI();
365  auto J = mat.WriteJ();
366  auto Data = mat.WriteData();
367  auto d_offsets = offsets.Read();
368  auto d_indices = indices.Read();
369  auto d_gatherMap = gatherMap.Read();
370  auto mat_ea = Reshape(ea_data.Read(), elt_dofs, elt_dofs, ne);
371  MFEM_FORALL(e, ne,
372  {
373  for (int i = 0; i < elt_dofs; i++)
374  {
375  int i_elts[Max];
376  int i_B[Max];
377  const int i_E = e*elt_dofs + i;
378  const int i_L = d_gatherMap[i_E];
379  const int i_offset = d_offsets[i_L];
380  const int i_nextOffset = d_offsets[i_L+1];
381  const int i_nbElts = i_nextOffset - i_offset;
382  for (int e_i = 0; e_i < i_nbElts; ++e_i)
383  {
384  const int i_E = d_indices[i_offset+e_i];
385  i_elts[e_i] = i_E/elt_dofs;
386  i_B[e_i] = i_E%elt_dofs;
387  }
388  for (int j = 0; j < elt_dofs; j++)
389  {
390  const int j_E = e*elt_dofs + j;
391  const int j_L = d_gatherMap[j_E];
392  const int j_offset = d_offsets[j_L];
393  const int j_nextOffset = d_offsets[j_L+1];
394  const int j_nbElts = j_nextOffset - j_offset;
395  if (i_nbElts == 1 || j_nbElts == 1) // no assembly required
396  {
397  const int nnz = GetAndIncrementNnzIndex(i_L, I);
398  J[nnz] = j_L;
399  Data[nnz] = mat_ea(j,i,e);
400  }
401  else // assembly required
402  {
403  int j_elts[Max];
404  int j_B[Max];
405  for (int e_j = 0; e_j < j_nbElts; ++e_j)
406  {
407  const int j_E = d_indices[j_offset+e_j];
408  const int elt = j_E/elt_dofs;
409  j_elts[e_j] = elt;
410  j_B[e_j] = j_E%elt_dofs;
411  }
412  int min_e = GetMinElt<Max>(i_elts, i_nbElts, j_elts, j_nbElts);
413  if (e == min_e) // add the nnz only once
414  {
415  double val = 0.0;
416  for (int i = 0; i < i_nbElts; i++)
417  {
418  const int e_i = i_elts[i];
419  const int i_Bloc = i_B[i];
420  for (int j = 0; j < j_nbElts; j++)
421  {
422  const int e_j = j_elts[j];
423  const int j_Bloc = j_B[j];
424  if (e_i == e_j)
425  {
426  val += mat_ea(j_Bloc, i_Bloc, e_i);
427  }
428  }
429  }
430  const int nnz = GetAndIncrementNnzIndex(i_L, I);
431  J[nnz] = j_L;
432  Data[nnz] = val;
433  }
434  }
435  }
436  }
437  });
438  // We need to shift again the entries of I, we do it on CPU as it is very
439  // sequential.
440  auto h_I = mat.HostReadWriteI();
441  const int size = vd*all_dofs;
442  for (int i = 0; i < size; i++)
443  {
444  h_I[size-i] = h_I[size-(i+1)];
445  }
446  h_I[0] = 0;
447 }
448 
450  : ne(fes.GetNE()),
451  vdim(fes.GetVDim()),
452  byvdim(fes.GetOrdering() == Ordering::byVDIM),
453  ndof(ne > 0 ? fes.GetFE(0)->GetDof() : 0),
454  ndofs(fes.GetNDofs())
455 {
456  height = vdim*ne*ndof;
457  width = vdim*ne*ndof;
458 }
459 
460 void L2ElementRestriction::Mult(const Vector &x, Vector &y) const
461 {
462  const int nd = ndof;
463  const int vd = vdim;
464  const bool t = byvdim;
465  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
466  auto d_y = Reshape(y.Write(), nd, vd, ne);
467  MFEM_FORALL(i, ndofs,
468  {
469  const int idx = i;
470  const int dof = idx % nd;
471  const int e = idx / nd;
472  for (int c = 0; c < vd; ++c)
473  {
474  d_y(dof, c, e) = d_x(t?c:idx, t?idx:c);
475  }
476  });
477 }
478 
480 {
481  const int nd = ndof;
482  const int vd = vdim;
483  const bool t = byvdim;
484  auto d_x = Reshape(x.Read(), nd, vd, ne);
485  auto d_y = Reshape(y.Write(), t?vd:ndofs, t?ndofs:vd);
486  MFEM_FORALL(i, ndofs,
487  {
488  const int idx = i;
489  const int dof = idx % nd;
490  const int e = idx / nd;
491  for (int c = 0; c < vd; ++c)
492  {
493  d_y(t?c:idx,t?idx:c) = d_x(dof, c, e);
494  }
495  });
496 }
497 
499 {
500  const int elem_dofs = ndof;
501  const int vd = vdim;
502  auto I = mat.WriteI();
503  MFEM_FORALL(dof, ne*elem_dofs*vd,
504  {
505  I[dof] = elem_dofs;
506  });
507 }
508 
509 static MFEM_HOST_DEVICE int AddNnz(const int iE, int *I, const int dofs)
510 {
511  int val = AtomicAdd(I[iE],dofs);
512  return val;
513 }
514 
516  SparseMatrix &mat) const
517 {
518  const int elem_dofs = ndof;
519  const int vd = vdim;
520  auto I = mat.ReadWriteI();
521  auto J = mat.WriteJ();
522  auto Data = mat.WriteData();
523  auto mat_ea = Reshape(ea_data.Read(), elem_dofs, elem_dofs, ne);
524  MFEM_FORALL(iE, ne*elem_dofs*vd,
525  {
526  const int offset = AddNnz(iE,I,elem_dofs);
527  const int e = iE/elem_dofs;
528  const int i = iE%elem_dofs;
529  for (int j = 0; j < elem_dofs; j++)
530  {
531  J[offset+j] = e*elem_dofs+j;
532  Data[offset+j] = mat_ea(j,i,e);
533  }
534  });
535 }
536 
537 // Return the face degrees of freedom returned in Lexicographic order.
538 void GetFaceDofs(const int dim, const int face_id,
539  const int dof1d, Array<int> &faceMap)
540 {
541  switch (dim)
542  {
543  case 1:
544  switch (face_id)
545  {
546  case 0: // WEST
547  faceMap[0] = 0;
548  break;
549  case 1: // EAST
550  faceMap[0] = dof1d-1;
551  break;
552  }
553  break;
554  case 2:
555  switch (face_id)
556  {
557  case 0: // SOUTH
558  for (int i = 0; i < dof1d; ++i)
559  {
560  faceMap[i] = i;
561  }
562  break;
563  case 1: // EAST
564  for (int i = 0; i < dof1d; ++i)
565  {
566  faceMap[i] = dof1d-1 + i*dof1d;
567  }
568  break;
569  case 2: // NORTH
570  for (int i = 0; i < dof1d; ++i)
571  {
572  faceMap[i] = (dof1d-1)*dof1d + i;
573  }
574  break;
575  case 3: // WEST
576  for (int i = 0; i < dof1d; ++i)
577  {
578  faceMap[i] = i*dof1d;
579  }
580  break;
581  }
582  break;
583  case 3:
584  switch (face_id)
585  {
586  case 0: // BOTTOM
587  for (int i = 0; i < dof1d; ++i)
588  {
589  for (int j = 0; j < dof1d; ++j)
590  {
591  faceMap[i+j*dof1d] = i + j*dof1d;
592  }
593  }
594  break;
595  case 1: // SOUTH
596  for (int i = 0; i < dof1d; ++i)
597  {
598  for (int j = 0; j < dof1d; ++j)
599  {
600  faceMap[i+j*dof1d] = i + j*dof1d*dof1d;
601  }
602  }
603  break;
604  case 2: // EAST
605  for (int i = 0; i < dof1d; ++i)
606  {
607  for (int j = 0; j < dof1d; ++j)
608  {
609  faceMap[i+j*dof1d] = dof1d-1 + i*dof1d + j*dof1d*dof1d;
610  }
611  }
612  break;
613  case 3: // NORTH
614  for (int i = 0; i < dof1d; ++i)
615  {
616  for (int j = 0; j < dof1d; ++j)
617  {
618  faceMap[i+j*dof1d] = (dof1d-1)*dof1d + i + j*dof1d*dof1d;
619  }
620  }
621  break;
622  case 4: // WEST
623  for (int i = 0; i < dof1d; ++i)
624  {
625  for (int j = 0; j < dof1d; ++j)
626  {
627  faceMap[i+j*dof1d] = i*dof1d + j*dof1d*dof1d;
628  }
629  }
630  break;
631  case 5: // TOP
632  for (int i = 0; i < dof1d; ++i)
633  {
634  for (int j = 0; j < dof1d; ++j)
635  {
636  faceMap[i+j*dof1d] = (dof1d-1)*dof1d*dof1d + i + j*dof1d;
637  }
638  }
639  break;
640  }
641  break;
642  }
643 }
644 
646  const ElementDofOrdering e_ordering,
647  const FaceType type)
648  : fes(fes),
649  nf(fes.GetNFbyType(type)),
650  vdim(fes.GetVDim()),
651  byvdim(fes.GetOrdering() == Ordering::byVDIM),
652  ndofs(fes.GetNDofs()),
653  dof(nf > 0 ? fes.GetFaceElement(0)->GetDof() : 0),
654  nfdofs(nf*dof),
655  scatter_indices(nf*dof),
656  offsets(ndofs+1),
657  gather_indices(nf*dof)
658 {
659  if (nf==0) { return; }
660  // If fespace == H1
661  const FiniteElement *fe = fes.GetFE(0);
662  const TensorBasisElement *tfe = dynamic_cast<const TensorBasisElement*>(fe);
663  MFEM_VERIFY(tfe != NULL &&
666  "Only Gauss-Lobatto and Bernstein basis are supported in "
667  "H1FaceRestriction.");
668  MFEM_VERIFY(fes.GetMesh()->Conforming(),
669  "Non-conforming meshes not yet supported with partial assembly.");
670  // Assuming all finite elements are using Gauss-Lobatto.
671  height = vdim*nf*dof;
672  width = fes.GetVSize();
673  const bool dof_reorder = (e_ordering == ElementDofOrdering::LEXICOGRAPHIC);
674  if (dof_reorder && nf > 0)
675  {
676  for (int f = 0; f < fes.GetNF(); ++f)
677  {
678  const FiniteElement *fe = fes.GetFaceElement(f);
679  const TensorBasisElement* el =
680  dynamic_cast<const TensorBasisElement*>(fe);
681  if (el) { continue; }
682  mfem_error("Finite element not suitable for lexicographic ordering");
683  }
684  const FiniteElement *fe = fes.GetFaceElement(0);
685  const TensorBasisElement* el =
686  dynamic_cast<const TensorBasisElement*>(fe);
687  const Array<int> &fe_dof_map = el->GetDofMap();
688  MFEM_VERIFY(fe_dof_map.Size() > 0, "invalid dof map");
689  }
690  const TensorBasisElement* el =
691  dynamic_cast<const TensorBasisElement*>(fe);
692  const int *dof_map = el->GetDofMap().GetData();
693  const Table& e2dTable = fes.GetElementToDofTable();
694  const int* elementMap = e2dTable.GetJ();
695  Array<int> faceMap(dof);
696  int e1, e2;
697  int inf1, inf2;
698  int face_id;
699  int orientation;
700  const int dof1d = fes.GetFE(0)->GetOrder()+1;
701  const int elem_dofs = fes.GetFE(0)->GetDof();
702  const int dim = fes.GetMesh()->SpaceDimension();
703  // Computation of scatter_indices
704  int f_ind = 0;
705  for (int f = 0; f < fes.GetNF(); ++f)
706  {
707  fes.GetMesh()->GetFaceElements(f, &e1, &e2);
708  fes.GetMesh()->GetFaceInfos(f, &inf1, &inf2);
709  orientation = inf1 % 64;
710  face_id = inf1 / 64;
711  if ((type==FaceType::Interior && (e2>=0 || (e2<0 && inf2>=0))) ||
712  (type==FaceType::Boundary && e2<0 && inf2<0) )
713  {
714  // Assumes Gauss-Lobatto basis
715  if (dof_reorder)
716  {
717  if (orientation != 0)
718  {
719  mfem_error("FaceRestriction used on degenerated mesh.");
720  }
721  GetFaceDofs(dim, face_id, dof1d, faceMap); // Only for hex
722  }
723  else
724  {
725  mfem_error("FaceRestriction not yet implemented for this type of "
726  "element.");
727  // TODO Something with GetFaceDofs?
728  }
729  for (int d = 0; d < dof; ++d)
730  {
731  const int face_dof = faceMap[d];
732  const int did = (!dof_reorder)?face_dof:dof_map[face_dof];
733  const int gid = elementMap[e1*elem_dofs + did];
734  const int lid = dof*f_ind + d;
735  scatter_indices[lid] = gid;
736  }
737  f_ind++;
738  }
739  }
740  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
741  // Computation of gather_indices
742  for (int i = 0; i <= ndofs; ++i)
743  {
744  offsets[i] = 0;
745  }
746  f_ind = 0;
747  for (int f = 0; f < fes.GetNF(); ++f)
748  {
749  fes.GetMesh()->GetFaceElements(f, &e1, &e2);
750  fes.GetMesh()->GetFaceInfos(f, &inf1, &inf2);
751  orientation = inf1 % 64;
752  face_id = inf1 / 64;
753  if ((type==FaceType::Interior && (e2>=0 || (e2<0 && inf2>=0))) ||
754  (type==FaceType::Boundary && e2<0 && inf2<0) )
755  {
756  GetFaceDofs(dim, face_id, dof1d, faceMap);
757  for (int d = 0; d < dof; ++d)
758  {
759  const int face_dof = faceMap[d];
760  const int did = (!dof_reorder)?face_dof:dof_map[face_dof];
761  const int gid = elementMap[e1*elem_dofs + did];
762  ++offsets[gid + 1];
763  }
764  f_ind++;
765  }
766  }
767  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
768  for (int i = 1; i <= ndofs; ++i)
769  {
770  offsets[i] += offsets[i - 1];
771  }
772  f_ind = 0;
773  for (int f = 0; f < fes.GetNF(); ++f)
774  {
775  fes.GetMesh()->GetFaceElements(f, &e1, &e2);
776  fes.GetMesh()->GetFaceInfos(f, &inf1, &inf2);
777  orientation = inf1 % 64;
778  face_id = inf1 / 64;
779  if ((type==FaceType::Interior && (e2>=0 || (e2<0 && inf2>=0))) ||
780  (type==FaceType::Boundary && e2<0 && inf2<0) )
781  {
782  GetFaceDofs(dim, face_id, dof1d, faceMap);
783  for (int d = 0; d < dof; ++d)
784  {
785  const int face_dof = faceMap[d];
786  const int did = (!dof_reorder)?face_dof:dof_map[face_dof];
787  const int gid = elementMap[e1*elem_dofs + did];
788  const int lid = dof*f_ind + d;
789  gather_indices[offsets[gid]++] = lid;
790  }
791  f_ind++;
792  }
793  }
794  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
795  for (int i = ndofs; i > 0; --i)
796  {
797  offsets[i] = offsets[i - 1];
798  }
799  offsets[0] = 0;
800 }
801 
802 void H1FaceRestriction::Mult(const Vector& x, Vector& y) const
803 {
804  // Assumes all elements have the same number of dofs
805  const int nd = dof;
806  const int vd = vdim;
807  const bool t = byvdim;
808  auto d_indices = scatter_indices.Read();
809  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
810  auto d_y = Reshape(y.Write(), nd, vd, nf);
811  MFEM_FORALL(i, nfdofs,
812  {
813  const int idx = d_indices[i];
814  const int dof = i % nd;
815  const int face = i / nd;
816  for (int c = 0; c < vd; ++c)
817  {
818  d_y(dof, c, face) = d_x(t?c:idx, t?idx:c);
819  }
820  });
821 }
822 
824 {
825  // Assumes all elements have the same number of dofs
826  const int nd = dof;
827  const int vd = vdim;
828  const bool t = byvdim;
829  auto d_offsets = offsets.Read();
830  auto d_indices = gather_indices.Read();
831  auto d_x = Reshape(x.Read(), nd, vd, nf);
832  auto d_y = Reshape(y.Write(), t?vd:ndofs, t?ndofs:vd);
833  MFEM_FORALL(i, ndofs,
834  {
835  const int offset = d_offsets[i];
836  const int nextOffset = d_offsets[i + 1];
837  for (int c = 0; c < vd; ++c)
838  {
839  double dofValue = 0;
840  for (int j = offset; j < nextOffset; ++j)
841  {
842  const int idx_j = d_indices[j];
843  dofValue += d_x(idx_j % nd, c, idx_j / nd);
844  }
845  d_y(t?c:i,t?i:c) += dofValue;
846  }
847  });
848 }
849 
850 static int ToLexOrdering2D(const int face_id, const int size1d, const int i)
851 {
852  if (face_id==2 || face_id==3)
853  {
854  return size1d-1-i;
855  }
856  else
857  {
858  return i;
859  }
860 }
861 
862 static int PermuteFace2D(const int face_id1, const int face_id2,
863  const int orientation,
864  const int size1d, const int index)
865 {
866  int new_index;
867  // Convert from lex ordering
868  if (face_id1==2 || face_id1==3)
869  {
870  new_index = size1d-1-index;
871  }
872  else
873  {
874  new_index = index;
875  }
876  // Permute based on face orientations
877  if (orientation==1)
878  {
879  new_index = size1d-1-new_index;
880  }
881  return ToLexOrdering2D(face_id2, size1d, new_index);
882 }
883 
884 static int ToLexOrdering3D(const int face_id, const int size1d, const int i,
885  const int j)
886 {
887  if (face_id==2 || face_id==1 || face_id==5)
888  {
889  return i + j*size1d;
890  }
891  else if (face_id==3 || face_id==4)
892  {
893  return (size1d-1-i) + j*size1d;
894  }
895  else // face_id==0
896  {
897  return i + (size1d-1-j)*size1d;
898  }
899 }
900 
901 static int PermuteFace3D(const int face_id1, const int face_id2,
902  const int orientation,
903  const int size1d, const int index)
904 {
905  int i=0, j=0, new_i=0, new_j=0;
906  i = index%size1d;
907  j = index/size1d;
908  // Convert from lex ordering
909  if (face_id1==3 || face_id1==4)
910  {
911  i = size1d-1-i;
912  }
913  else if (face_id1==0)
914  {
915  j = size1d-1-j;
916  }
917  // Permute based on face orientations
918  switch (orientation)
919  {
920  case 0:
921  new_i = i;
922  new_j = j;
923  break;
924  case 1:
925  new_i = j;
926  new_j = i;
927  break;
928  case 2:
929  new_i = j;
930  new_j = (size1d-1-i);
931  break;
932  case 3:
933  new_i = (size1d-1-i);
934  new_j = j;
935  break;
936  case 4:
937  new_i = (size1d-1-i);
938  new_j = (size1d-1-j);
939  break;
940  case 5:
941  new_i = (size1d-1-j);
942  new_j = (size1d-1-i);
943  break;
944  case 6:
945  new_i = (size1d-1-j);
946  new_j = i;
947  break;
948  case 7:
949  new_i = i;
950  new_j = (size1d-1-j);
951  break;
952  }
953  return ToLexOrdering3D(face_id2, size1d, new_i, new_j);
954 }
955 
956 // Permute dofs or quads on a face for e2 to match with the ordering of e1
957 int PermuteFaceL2(const int dim, const int face_id1,
958  const int face_id2, const int orientation,
959  const int size1d, const int index)
960 {
961  switch (dim)
962  {
963  case 1:
964  return 0;
965  case 2:
966  return PermuteFace2D(face_id1, face_id2, orientation, size1d, index);
967  case 3:
968  return PermuteFace3D(face_id1, face_id2, orientation, size1d, index);
969  default:
970  mfem_error("Unsupported dimension.");
971  return 0;
972  }
973 }
974 
976  const FaceType type,
977  const L2FaceValues m)
978  : fes(fes),
979  nf(fes.GetNFbyType(type)),
980  ne(fes.GetNE()),
981  vdim(fes.GetVDim()),
982  byvdim(fes.GetOrdering() == Ordering::byVDIM),
983  ndofs(fes.GetNDofs()),
984  dof(nf > 0 ?
985  fes.GetTraceElement(0, fes.GetMesh()->GetFaceBaseGeometry(0))->GetDof()
986  : 0),
987  elemDofs(fes.GetFE(0)->GetDof()),
988  m(m),
989  nfdofs(nf*dof),
990  scatter_indices1(nf*dof),
991  scatter_indices2(m==L2FaceValues::DoubleValued?nf*dof:0),
992  offsets(ndofs+1),
993  gather_indices((m==L2FaceValues::DoubleValued? 2 : 1)*nf*dof)
994 {
995 }
996 
998  const ElementDofOrdering e_ordering,
999  const FaceType type,
1000  const L2FaceValues m)
1001  : L2FaceRestriction(fes, type, m)
1002 {
1003  // If fespace == L2
1004  const FiniteElement *fe = fes.GetFE(0);
1005  const TensorBasisElement *tfe = dynamic_cast<const TensorBasisElement*>(fe);
1006  MFEM_VERIFY(tfe != NULL &&
1009  "Only Gauss-Lobatto and Bernstein basis are supported in "
1010  "L2FaceRestriction.");
1011  MFEM_VERIFY(fes.GetMesh()->Conforming(),
1012  "Non-conforming meshes not yet supported with partial assembly.");
1013  if (nf==0) { return; }
1015  width = fes.GetVSize();
1016  const bool dof_reorder = (e_ordering == ElementDofOrdering::LEXICOGRAPHIC);
1017  if (!dof_reorder)
1018  {
1019  mfem_error("Non-Tensor L2FaceRestriction not yet implemented.");
1020  }
1021  if (dof_reorder && nf > 0)
1022  {
1023  for (int f = 0; f < fes.GetNF(); ++f)
1024  {
1025  const FiniteElement *fe =
1027  const TensorBasisElement* el =
1028  dynamic_cast<const TensorBasisElement*>(fe);
1029  if (el) { continue; }
1030  mfem_error("Finite element not suitable for lexicographic ordering");
1031  }
1032  }
1033  const Table& e2dTable = fes.GetElementToDofTable();
1034  const int* elementMap = e2dTable.GetJ();
1035  Array<int> faceMap1(dof), faceMap2(dof);
1036  int e1, e2;
1037  int inf1, inf2;
1038  int face_id1, face_id2;
1039  int orientation;
1040  const int dof1d = fes.GetFE(0)->GetOrder()+1;
1041  const int elem_dofs = fes.GetFE(0)->GetDof();
1042  const int dim = fes.GetMesh()->SpaceDimension();
1043  // Computation of scatter indices
1044  int f_ind=0;
1045  for (int f = 0; f < fes.GetNF(); ++f)
1046  {
1047  fes.GetMesh()->GetFaceElements(f, &e1, &e2);
1048  fes.GetMesh()->GetFaceInfos(f, &inf1, &inf2);
1049  if (dof_reorder)
1050  {
1051  orientation = inf1 % 64;
1052  face_id1 = inf1 / 64;
1053  GetFaceDofs(dim, face_id1, dof1d, faceMap1); // Only for hex
1054  orientation = inf2 % 64;
1055  face_id2 = inf2 / 64;
1056  GetFaceDofs(dim, face_id2, dof1d, faceMap2); // Only for hex
1057  }
1058  else
1059  {
1060  mfem_error("FaceRestriction not yet implemented for this type of "
1061  "element.");
1062  // TODO Something with GetFaceDofs?
1063  }
1064  if ((type==FaceType::Interior && e2>=0) ||
1065  (type==FaceType::Boundary && e2<0))
1066  {
1067  for (int d = 0; d < dof; ++d)
1068  {
1069  const int face_dof = faceMap1[d];
1070  const int did = face_dof;
1071  const int gid = elementMap[e1*elem_dofs + did];
1072  const int lid = dof*f_ind + d;
1073  scatter_indices1[lid] = gid;
1074  }
1076  {
1077  for (int d = 0; d < dof; ++d)
1078  {
1079  if (type==FaceType::Interior && e2>=0) // interior face
1080  {
1081  const int pd = PermuteFaceL2(dim, face_id1, face_id2,
1082  orientation, dof1d, d);
1083  const int face_dof = faceMap2[pd];
1084  const int did = face_dof;
1085  const int gid = elementMap[e2*elem_dofs + did];
1086  const int lid = dof*f_ind + d;
1087  scatter_indices2[lid] = gid;
1088  }
1089  else if (type==FaceType::Boundary && e2<0) // true boundary face
1090  {
1091  const int lid = dof*f_ind + d;
1092  scatter_indices2[lid] = -1;
1093  }
1094  }
1095  }
1096  f_ind++;
1097  }
1098  }
1099  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
1100  // Computation of gather_indices
1101  for (int i = 0; i <= ndofs; ++i)
1102  {
1103  offsets[i] = 0;
1104  }
1105  f_ind = 0;
1106  for (int f = 0; f < fes.GetNF(); ++f)
1107  {
1108  fes.GetMesh()->GetFaceElements(f, &e1, &e2);
1109  fes.GetMesh()->GetFaceInfos(f, &inf1, &inf2);
1110  if ((type==FaceType::Interior && (e2>=0 || (e2<0 && inf2>=0))) ||
1111  (type==FaceType::Boundary && e2<0 && inf2<0) )
1112  {
1113  orientation = inf1 % 64;
1114  face_id1 = inf1 / 64;
1115  GetFaceDofs(dim, face_id1, dof1d, faceMap1);
1116  orientation = inf2 % 64;
1117  face_id2 = inf2 / 64;
1118  GetFaceDofs(dim, face_id2, dof1d, faceMap2);
1119  for (int d = 0; d < dof; ++d)
1120  {
1121  const int did = faceMap1[d];
1122  const int gid = elementMap[e1*elem_dofs + did];
1123  ++offsets[gid + 1];
1124  }
1126  {
1127  for (int d = 0; d < dof; ++d)
1128  {
1129  if (type==FaceType::Interior && e2>=0) // interior face
1130  {
1131  const int pd = PermuteFaceL2(dim, face_id1, face_id2,
1132  orientation, dof1d, d);
1133  const int did = faceMap2[pd];
1134  const int gid = elementMap[e2*elem_dofs + did];
1135  ++offsets[gid + 1];
1136  }
1137  }
1138  }
1139  f_ind++;
1140  }
1141  }
1142  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
1143  for (int i = 1; i <= ndofs; ++i)
1144  {
1145  offsets[i] += offsets[i - 1];
1146  }
1147  f_ind = 0;
1148  for (int f = 0; f < fes.GetNF(); ++f)
1149  {
1150  fes.GetMesh()->GetFaceElements(f, &e1, &e2);
1151  fes.GetMesh()->GetFaceInfos(f, &inf1, &inf2);
1152  if ((type==FaceType::Interior && (e2>=0 || (e2<0 && inf2>=0))) ||
1153  (type==FaceType::Boundary && e2<0 && inf2<0) )
1154  {
1155  orientation = inf1 % 64;
1156  face_id1 = inf1 / 64;
1157  GetFaceDofs(dim, face_id1, dof1d, faceMap1);
1158  orientation = inf2 % 64;
1159  face_id2 = inf2 / 64;
1160  GetFaceDofs(dim, face_id2, dof1d, faceMap2);
1161  for (int d = 0; d < dof; ++d)
1162  {
1163  const int did = faceMap1[d];
1164  const int gid = elementMap[e1*elem_dofs + did];
1165  const int lid = dof*f_ind + d;
1166  // We don't shift lid to express that it's e1 of f
1167  gather_indices[offsets[gid]++] = lid;
1168  }
1170  {
1171  for (int d = 0; d < dof; ++d)
1172  {
1173  if (type==FaceType::Interior && e2>=0) // interior face
1174  {
1175  const int pd = PermuteFaceL2(dim, face_id1, face_id2,
1176  orientation, dof1d, d);
1177  const int did = faceMap2[pd];
1178  const int gid = elementMap[e2*elem_dofs + did];
1179  const int lid = dof*f_ind + d;
1180  // We shift lid to express that it's e2 of f
1181  gather_indices[offsets[gid]++] = nfdofs + lid;
1182  }
1183  }
1184  }
1185  f_ind++;
1186  }
1187  }
1188  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
1189  for (int i = ndofs; i > 0; --i)
1190  {
1191  offsets[i] = offsets[i - 1];
1192  }
1193  offsets[0] = 0;
1194 }
1195 
1196 void L2FaceRestriction::Mult(const Vector& x, Vector& y) const
1197 {
1198  // Assumes all elements have the same number of dofs
1199  const int nd = dof;
1200  const int vd = vdim;
1201  const bool t = byvdim;
1202 
1204  {
1205  auto d_indices1 = scatter_indices1.Read();
1206  auto d_indices2 = scatter_indices2.Read();
1207  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
1208  auto d_y = Reshape(y.Write(), nd, vd, 2, nf);
1209  MFEM_FORALL(i, nfdofs,
1210  {
1211  const int dof = i % nd;
1212  const int face = i / nd;
1213  const int idx1 = d_indices1[i];
1214  for (int c = 0; c < vd; ++c)
1215  {
1216  d_y(dof, c, 0, face) = d_x(t?c:idx1, t?idx1:c);
1217  }
1218  const int idx2 = d_indices2[i];
1219  for (int c = 0; c < vd; ++c)
1220  {
1221  d_y(dof, c, 1, face) = idx2==-1 ? 0.0 : d_x(t?c:idx2, t?idx2:c);
1222  }
1223  });
1224  }
1225  else
1226  {
1227  auto d_indices1 = scatter_indices1.Read();
1228  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
1229  auto d_y = Reshape(y.Write(), nd, vd, nf);
1230  MFEM_FORALL(i, nfdofs,
1231  {
1232  const int dof = i % nd;
1233  const int face = i / nd;
1234  const int idx1 = d_indices1[i];
1235  for (int c = 0; c < vd; ++c)
1236  {
1237  d_y(dof, c, face) = d_x(t?c:idx1, t?idx1:c);
1238  }
1239  });
1240  }
1241 }
1242 
1244 {
1245  // Assumes all elements have the same number of dofs
1246  const int nd = dof;
1247  const int vd = vdim;
1248  const bool t = byvdim;
1249  const int dofs = nfdofs;
1250  auto d_offsets = offsets.Read();
1251  auto d_indices = gather_indices.Read();
1252 
1254  {
1255  auto d_x = Reshape(x.Read(), nd, vd, 2, nf);
1256  auto d_y = Reshape(y.Write(), t?vd:ndofs, t?ndofs:vd);
1257  MFEM_FORALL(i, ndofs,
1258  {
1259  const int offset = d_offsets[i];
1260  const int nextOffset = d_offsets[i + 1];
1261  for (int c = 0; c < vd; ++c)
1262  {
1263  double dofValue = 0;
1264  for (int j = offset; j < nextOffset; ++j)
1265  {
1266  int idx_j = d_indices[j];
1267  bool isE1 = idx_j < dofs;
1268  idx_j = isE1 ? idx_j : idx_j - dofs;
1269  dofValue += isE1 ?
1270  d_x(idx_j % nd, c, 0, idx_j / nd)
1271  :d_x(idx_j % nd, c, 1, idx_j / nd);
1272  }
1273  d_y(t?c:i,t?i:c) += dofValue;
1274  }
1275  });
1276  }
1277  else
1278  {
1279  auto d_x = Reshape(x.Read(), nd, vd, nf);
1280  auto d_y = Reshape(y.Write(), t?vd:ndofs, t?ndofs:vd);
1281  MFEM_FORALL(i, ndofs,
1282  {
1283  const int offset = d_offsets[i];
1284  const int nextOffset = d_offsets[i + 1];
1285  for (int c = 0; c < vd; ++c)
1286  {
1287  double dofValue = 0;
1288  for (int j = offset; j < nextOffset; ++j)
1289  {
1290  int idx_j = d_indices[j];
1291  dofValue += d_x(idx_j % nd, c, idx_j / nd);
1292  }
1293  d_y(t?c:i,t?i:c) += dofValue;
1294  }
1295  });
1296  }
1297 }
1298 
1300  SparseMatrix &face_mat) const
1301 {
1302  const int face_dofs = dof;
1303  auto d_indices1 = scatter_indices1.Read();
1304  auto d_indices2 = scatter_indices2.Read();
1305  auto I = mat.ReadWriteI();
1306  MFEM_FORALL(fdof, nf*face_dofs,
1307  {
1308  const int iE1 = d_indices1[fdof];
1309  const int iE2 = d_indices2[fdof];
1310  AddNnz(iE1,I,face_dofs);
1311  AddNnz(iE2,I,face_dofs);
1312  });
1313 }
1314 
1316  SparseMatrix &mat,
1317  SparseMatrix &face_mat) const
1318 {
1319  const int face_dofs = dof;
1320  auto d_indices1 = scatter_indices1.Read();
1321  auto d_indices2 = scatter_indices2.Read();
1322  auto I = mat.ReadWriteI();
1323  auto mat_fea = Reshape(ea_data.Read(), face_dofs, face_dofs, 2, nf);
1324  auto J = mat.WriteJ();
1325  auto Data = mat.WriteData();
1326  MFEM_FORALL(fdof, nf*face_dofs,
1327  {
1328  const int f = fdof/face_dofs;
1329  const int iF = fdof%face_dofs;
1330  const int iE1 = d_indices1[f*face_dofs+iF];
1331  const int iE2 = d_indices2[f*face_dofs+iF];
1332  const int offset1 = AddNnz(iE1,I,face_dofs);
1333  const int offset2 = AddNnz(iE2,I,face_dofs);
1334  for (int jF = 0; jF < face_dofs; jF++)
1335  {
1336  const int jE1 = d_indices1[f*face_dofs+jF];
1337  const int jE2 = d_indices2[f*face_dofs+jF];
1338  J[offset2+jF] = jE1;
1339  J[offset1+jF] = jE2;
1340  Data[offset2+jF] = mat_fea(jF,iF,0,f);
1341  Data[offset1+jF] = mat_fea(jF,iF,1,f);
1342  }
1343  });
1344 }
1345 
1347  Vector &ea_data) const
1348 {
1349  const int face_dofs = dof;
1350  const int elem_dofs = elemDofs;
1351  const int NE = ne;
1353  {
1354  auto d_indices1 = scatter_indices1.Read();
1355  auto d_indices2 = scatter_indices2.Read();
1356  auto mat_fea = Reshape(fea_data.Read(), face_dofs, face_dofs, 2, nf);
1357  auto mat_ea = Reshape(ea_data.ReadWrite(), elem_dofs, elem_dofs, ne);
1358  MFEM_FORALL(f, nf,
1359  {
1360  const int e1 = d_indices1[f*face_dofs]/elem_dofs;
1361  const int e2 = d_indices2[f*face_dofs]/elem_dofs;
1362  for (int j = 0; j < face_dofs; j++)
1363  {
1364  const int jB1 = d_indices1[f*face_dofs+j]%elem_dofs;
1365  for (int i = 0; i < face_dofs; i++)
1366  {
1367  const int iB1 = d_indices1[f*face_dofs+i]%elem_dofs;
1368  AtomicAdd(mat_ea(iB1,jB1,e1), mat_fea(i,j,0,f));
1369  }
1370  }
1371  if (e2 < NE)
1372  {
1373  for (int j = 0; j < face_dofs; j++)
1374  {
1375  const int jB2 = d_indices2[f*face_dofs+j]%elem_dofs;
1376  for (int i = 0; i < face_dofs; i++)
1377  {
1378  const int iB2 = d_indices2[f*face_dofs+i]%elem_dofs;
1379  AtomicAdd(mat_ea(iB2,jB2,e2), mat_fea(i,j,1,f));
1380  }
1381  }
1382  }
1383  });
1384  }
1385  else
1386  {
1387  auto d_indices = scatter_indices1.Read();
1388  auto mat_fea = Reshape(fea_data.Read(), face_dofs, face_dofs, nf);
1389  auto mat_ea = Reshape(ea_data.ReadWrite(), elem_dofs, elem_dofs, ne);
1390  MFEM_FORALL(f, nf,
1391  {
1392  const int e = d_indices[f*face_dofs]/elem_dofs;
1393  for (int j = 0; j < face_dofs; j++)
1394  {
1395  const int jE = d_indices[f*face_dofs+j]%elem_dofs;
1396  for (int i = 0; i < face_dofs; i++)
1397  {
1398  const int iE = d_indices[f*face_dofs+i]%elem_dofs;
1399  AtomicAdd(mat_ea(iE,jE,e), mat_fea(i,j,f));
1400  }
1401  }
1402  });
1403  }
1404 }
1405 
1406 int ToLexOrdering(const int dim, const int face_id, const int size1d,
1407  const int index)
1408 {
1409  switch (dim)
1410  {
1411  case 1:
1412  return 0;
1413  case 2:
1414  return ToLexOrdering2D(face_id, size1d, index);
1415  case 3:
1416  return ToLexOrdering3D(face_id, size1d, index%size1d, index/size1d);
1417  default:
1418  mfem_error("Unsupported dimension.");
1419  return 0;
1420  }
1421 }
1422 
1423 } // namespace mfem
Abstract class for all finite elements.
Definition: fe.hpp:235
T * HostReadWrite()
Shortcut for mfem::ReadWrite(a.GetMemory(), a.Size(), false).
Definition: array.hpp:310
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:400
int * GetJ()
Definition: table.hpp:114
L2FaceValues
Definition: restriction.hpp:26
void MultTransposeUnsigned(const Vector &x, Vector &y) const
Compute MultTranspose without applying signs based on DOF orientations.
const FiniteElement * GetTraceElement(int i, Geometry::Type geom_type) const
Return the trace element from element &#39;i&#39; to the given &#39;geom_type&#39;.
Definition: fespace.cpp:2108
virtual void FillI(SparseMatrix &mat, SparseMatrix &face_mat) const
bool Conforming() const
Definition: mesh.hpp:1213
Array< int > gather_indices
const T * HostRead() const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), false).
Definition: array.hpp:294
MFEM_HOST_DEVICE T AtomicAdd(T &add, const T val)
Definition: backends.hpp:59
T * GetData()
Returns the data.
Definition: array.hpp:98
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe.hpp:319
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void AddFaceMatricesToElementMatrices(Vector &fea_data, Vector &ea_data) const
This methods adds the DG face matrices to the element matrices.
H1FaceRestriction(const FiniteElementSpace &, const ElementDofOrdering, const FaceType)
Geometry::Type GetFaceBaseGeometry(int i) const
Definition: mesh.hpp:829
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
const FiniteElement * GetFaceElement(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th face in the ...
Definition: fespace.cpp:2073
double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:380
Memory< double > & GetMemoryData()
Definition: sparsemat.hpp:220
const L2FaceValues m
Memory< int > & GetMemoryI()
Definition: sparsemat.hpp:188
MemoryType GetMemoryType() const
Return a MemoryType that is currently valid. If both the host and the device pointers are currently v...
Array< int > gather_indices
L2ElementRestriction(const FiniteElementSpace &)
Memory< int > & GetMemoryJ()
Definition: sparsemat.hpp:204
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 GetBasisType() const
Definition: fe.hpp:2021
FaceType
Definition: mesh.hpp:45
double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:388
Bernstein polynomials.
Definition: fe.hpp:35
L2FaceRestriction(const FiniteElementSpace &, const FaceType, const L2FaceValues m=L2FaceValues::DoubleValued)
Data type sparse matrix.
Definition: sparsemat.hpp:46
int FillI(SparseMatrix &mat) const
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:65
void FillI(SparseMatrix &mat) const
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:314
Operator that extracts Face degrees of freedom.
void BooleanMask(Vector &y) const
Fills the E-vector y with boolean values 0.0 and 1.0 such that each each entry of the L-vector is uni...
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:153
int * ReadWriteI(bool on_dev=true)
Definition: sparsemat.hpp:194
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:290
void MultUnsigned(const Vector &x, Vector &y) const
Compute Mult without applying signs based on DOF orientations.
int * WriteI(bool on_dev=true)
Definition: sparsemat.hpp:192
ElementRestriction(const FiniteElementSpace &, ElementDofOrdering)
Definition: restriction.cpp:20
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
Definition: fespace.hpp:433
int SpaceDimension() const
Definition: mesh.hpp:789
int ToLexOrdering(const int dim, const int face_id, const int size1d, const int index)
const FiniteElementSpace & fes
Definition: restriction.hpp:39
void GetFaceElements(int Face, int *Elem1, int *Elem2) const
Definition: mesh.cpp:1031
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void GetFaceInfos(int Face, int *Inf1, int *Inf2) const
Definition: mesh.cpp:1037
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:315
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void GetFaceDofs(const int dim, const int face_id, const int dof1d, Array< int > &faceMap)
The ordering method used when the number of unknowns per mesh node (vector dimension) is bigger than ...
Definition: fespace.hpp:28
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
int * WriteJ(bool on_dev=true)
Definition: sparsemat.hpp:208
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
const Array< int > & GetDofMap() const
Get an Array&lt;int&gt; that maps lexicographically ordered indices to the indices of the respective nodes/...
Definition: fe.hpp:2029
Array< int > scatter_indices1
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:372
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:65
int dim
Definition: ex24.cpp:53
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:1798
int index(int i, int j, int nx, int ny)
Definition: life.cpp:241
Lexicographic ordering for tensor-product FiniteElements.
virtual void FillJAndData(const Vector &ea_data, SparseMatrix &mat, SparseMatrix &face_mat) const
double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:384
Vector data type.
Definition: vector.hpp:51
int PermuteFaceL2(const int dim, const int face_id1, const int face_id2, const int orientation, const int size1d, const int index)
const Table & GetElementToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each mesh element...
Definition: fespace.hpp:543
int * HostReadWriteI()
Definition: sparsemat.hpp:200
Array< int > scatter_indices
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
Array< int > scatter_indices2
void FillSparseMatrix(const Vector &mat_ea, SparseMatrix &mat) const
Fill a Sparse Matrix with Element Matrices.
VDIM x NQPT x NE.
double f(const Vector &p)
double * WriteData(bool on_dev=true)
Definition: sparsemat.hpp:224
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...