MFEM  v4.4.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-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 "restriction.hpp"
13 #include "gridfunc.hpp"
14 #include "fespace.hpp"
15 #include "../general/forall.hpp"
16 #include <climits>
17 
18 #ifdef MFEM_USE_MPI
19 
20 #include "pfespace.hpp"
21 
22 #endif
23 
24 namespace mfem
25 {
26 
28  ElementDofOrdering e_ordering)
29  : fes(f),
30  ne(fes.GetNE()),
31  vdim(fes.GetVDim()),
32  byvdim(fes.GetOrdering() == Ordering::byVDIM),
33  ndofs(fes.GetNDofs()),
34  dof(ne > 0 ? fes.GetFE(0)->GetDof() : 0),
35  nedofs(ne*dof),
36  offsets(ndofs+1),
37  indices(ne*dof),
38  gather_map(ne*dof)
39 {
40  // Assuming all finite elements are the same.
41  height = vdim*ne*dof;
42  width = fes.GetVSize();
43  const bool dof_reorder = (e_ordering == ElementDofOrdering::LEXICOGRAPHIC);
44  const int *dof_map = NULL;
45  if (dof_reorder && ne > 0)
46  {
47  for (int e = 0; e < ne; ++e)
48  {
49  const FiniteElement *fe = fes.GetFE(e);
50  const TensorBasisElement* el =
51  dynamic_cast<const TensorBasisElement*>(fe);
52  if (el) { continue; }
53  MFEM_ABORT("Finite element not suitable for lexicographic ordering");
54  }
55  const FiniteElement *fe = fes.GetFE(0);
56  const TensorBasisElement* el =
57  dynamic_cast<const TensorBasisElement*>(fe);
58  const Array<int> &fe_dof_map = el->GetDofMap();
59  MFEM_VERIFY(fe_dof_map.Size() > 0, "invalid dof map");
60  dof_map = fe_dof_map.GetData();
61  }
62  const Table& e2dTable = fes.GetElementToDofTable();
63  const int* element_map = e2dTable.GetJ();
64  // We will be keeping a count of how many local nodes point to its global dof
65  for (int i = 0; i <= ndofs; ++i)
66  {
67  offsets[i] = 0;
68  }
69  for (int e = 0; e < ne; ++e)
70  {
71  for (int d = 0; d < dof; ++d)
72  {
73  const int sgid = element_map[dof*e + d]; // signed
74  const int gid = (sgid >= 0) ? sgid : -1 - sgid;
75  ++offsets[gid + 1];
76  }
77  }
78  // Aggregate to find offsets for each global dof
79  for (int i = 1; i <= ndofs; ++i)
80  {
81  offsets[i] += offsets[i - 1];
82  }
83  // For each global dof, fill in all local nodes that point to it
84  for (int e = 0; e < ne; ++e)
85  {
86  for (int d = 0; d < dof; ++d)
87  {
88  const int sdid = dof_reorder ? dof_map[d] : 0; // signed
89  const int did = (!dof_reorder)?d:(sdid >= 0 ? sdid : -1-sdid);
90  const int sgid = element_map[dof*e + did]; // signed
91  const int gid = (sgid >= 0) ? sgid : -1-sgid;
92  const int lid = dof*e + d;
93  const bool plus = (sgid >= 0 && sdid >= 0) || (sgid < 0 && sdid < 0);
94  gather_map[lid] = plus ? gid : -1-gid;
95  indices[offsets[gid]++] = plus ? lid : -1-lid;
96  }
97  }
98  // We shifted the offsets vector by 1 by using it as a counter.
99  // Now we shift it back.
100  for (int i = ndofs; i > 0; --i)
101  {
102  offsets[i] = offsets[i - 1];
103  }
104  offsets[0] = 0;
105 }
106 
107 void ElementRestriction::Mult(const Vector& x, Vector& y) const
108 {
109  // Assumes all elements have the same number of dofs
110  const int nd = dof;
111  const int vd = vdim;
112  const bool t = byvdim;
113  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
114  auto d_y = Reshape(y.Write(), nd, vd, ne);
115  auto d_gather_map = gather_map.Read();
116  MFEM_FORALL(i, dof*ne,
117  {
118  const int gid = d_gather_map[i];
119  const bool plus = gid >= 0;
120  const int j = plus ? gid : -1-gid;
121  for (int c = 0; c < vd; ++c)
122  {
123  const double dof_value = d_x(t?c:j, t?j:c);
124  d_y(i % nd, c, i / nd) = plus ? dof_value : -dof_value;
125  }
126  });
127 }
128 
130 {
131  // Assumes all elements have the same number of dofs
132  const int nd = dof;
133  const int vd = vdim;
134  const bool t = byvdim;
135  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
136  auto d_y = Reshape(y.Write(), nd, vd, ne);
137  auto d_gather_map = gather_map.Read();
138 
139  MFEM_FORALL(i, dof*ne,
140  {
141  const int gid = d_gather_map[i];
142  const int j = gid >= 0 ? gid : -1-gid;
143  for (int c = 0; c < vd; ++c)
144  {
145  d_y(i % nd, c, i / nd) = d_x(t?c:j, t?j:c);
146  }
147  });
148 }
149 
151 {
152  // Assumes all elements have the same number of dofs
153  const int nd = dof;
154  const int vd = vdim;
155  const bool t = byvdim;
156  auto d_offsets = offsets.Read();
157  auto d_indices = indices.Read();
158  auto d_x = Reshape(x.Read(), nd, vd, ne);
159  auto d_y = Reshape(y.Write(), t?vd:ndofs, t?ndofs:vd);
160  MFEM_FORALL(i, ndofs,
161  {
162  const int offset = d_offsets[i];
163  const int next_offset = d_offsets[i + 1];
164  for (int c = 0; c < vd; ++c)
165  {
166  double dof_value = 0;
167  for (int j = offset; j < next_offset; ++j)
168  {
169  const int idx_j = (d_indices[j] >= 0) ? d_indices[j] : -1 - d_indices[j];
170  dof_value += ((d_indices[j] >= 0) ? d_x(idx_j % nd, c, idx_j / nd) :
171  -d_x(idx_j % nd, c, idx_j / nd));
172  }
173  d_y(t?c:i,t?i:c) = dof_value;
174  }
175  });
176 }
177 
179 {
180  // Assumes all elements have the same number of dofs
181  const int nd = dof;
182  const int vd = vdim;
183  const bool t = byvdim;
184  auto d_offsets = offsets.Read();
185  auto d_indices = indices.Read();
186  auto d_x = Reshape(x.Read(), nd, vd, ne);
187  auto d_y = Reshape(y.Write(), t?vd:ndofs, t?ndofs:vd);
188  MFEM_FORALL(i, ndofs,
189  {
190  const int offset = d_offsets[i];
191  const int next_offset = d_offsets[i + 1];
192  for (int c = 0; c < vd; ++c)
193  {
194  double dof_value = 0;
195  for (int j = offset; j < next_offset; ++j)
196  {
197  const int idx_j = (d_indices[j] >= 0) ? d_indices[j] : -1 - d_indices[j];
198  dof_value += d_x(idx_j % nd, c, idx_j / nd);
199  }
200  d_y(t?c:i,t?i:c) = dof_value;
201  }
202  });
203 }
204 
206 {
207  // Assumes all elements have the same number of dofs
208  const int nd = dof;
209  const int vd = vdim;
210  const bool t = byvdim;
211  auto d_offsets = offsets.Read();
212  auto d_indices = indices.Read();
213  auto d_x = Reshape(x.Read(), nd, vd, ne);
214  auto d_y = Reshape(y.Write(), t?vd:ndofs, t?ndofs:vd);
215  MFEM_FORALL(i, ndofs,
216  {
217  const int next_offset = d_offsets[i + 1];
218  for (int c = 0; c < vd; ++c)
219  {
220  double dof_value = 0;
221  const int j = next_offset - 1;
222  const int idx_j = (d_indices[j] >= 0) ? d_indices[j] : -1 - d_indices[j];
223  dof_value = (d_indices[j] >= 0) ? d_x(idx_j % nd, c, idx_j / nd) :
224  -d_x(idx_j % nd, c, idx_j / nd);
225  d_y(t?c:i,t?i:c) = dof_value;
226  }
227  });
228 }
229 
231 {
232  // Assumes all elements have the same number of dofs
233  const int nd = dof;
234  const int vd = vdim;
235  const bool t = byvdim;
236 
237  Array<char> processed(vd * ndofs);
238  processed = 0;
239 
240  auto d_offsets = offsets.HostRead();
241  auto d_indices = indices.HostRead();
242  auto d_x = Reshape(processed.HostReadWrite(), t?vd:ndofs, t?ndofs:vd);
243  auto d_y = Reshape(y.HostWrite(), nd, vd, ne);
244  for (int i = 0; i < ndofs; ++i)
245  {
246  const int offset = d_offsets[i];
247  const int next_offset = d_offsets[i+1];
248  for (int c = 0; c < vd; ++c)
249  {
250  for (int j = offset; j < next_offset; ++j)
251  {
252  const int idx_j = d_indices[j];
253  if (d_x(t?c:i,t?i:c))
254  {
255  d_y(idx_j % nd, c, idx_j / nd) = 0.0;
256  }
257  else
258  {
259  d_y(idx_j % nd, c, idx_j / nd) = 1.0;
260  d_x(t?c:i,t?i:c) = 1;
261  }
262  }
263  }
264  }
265 }
266 
268  SparseMatrix &mat) const
269 {
270  mat.GetMemoryI().New(mat.Height()+1, mat.GetMemoryI().GetMemoryType());
271  const int nnz = FillI(mat);
272  mat.GetMemoryJ().New(nnz, mat.GetMemoryJ().GetMemoryType());
273  mat.GetMemoryData().New(nnz, mat.GetMemoryData().GetMemoryType());
274  FillJAndData(mat_ea, mat);
275 }
276 
277 static MFEM_HOST_DEVICE int GetMinElt(const int *my_elts, const int nbElts,
278  const int *nbr_elts, const int nbrNbElts)
279 {
280  // Find the minimal element index found in both my_elts[] and nbr_elts[]
281  int min_el = INT_MAX;
282  for (int i = 0; i < nbElts; i++)
283  {
284  const int e_i = my_elts[i];
285  if (e_i >= min_el) { continue; }
286  for (int j = 0; j < nbrNbElts; j++)
287  {
288  if (e_i==nbr_elts[j])
289  {
290  min_el = e_i; // we already know e_i < min_el
291  break;
292  }
293  }
294  }
295  return min_el;
296 }
297 
298 /** Returns the index where a non-zero entry should be added and increment the
299  number of non-zeros for the row i_L. */
300 static MFEM_HOST_DEVICE int GetAndIncrementNnzIndex(const int i_L, int* I)
301 {
302  int ind = AtomicAdd(I[i_L],1);
303  return ind;
304 }
305 
307 {
308  static constexpr int Max = MaxNbNbr;
309  const int all_dofs = ndofs;
310  const int vd = vdim;
311  const int elt_dofs = dof;
312  auto I = mat.ReadWriteI();
313  auto d_offsets = offsets.Read();
314  auto d_indices = indices.Read();
315  auto d_gather_map = gather_map.Read();
316  MFEM_FORALL(i_L, vd*all_dofs+1,
317  {
318  I[i_L] = 0;
319  });
320  MFEM_FORALL(e, ne,
321  {
322  for (int i = 0; i < elt_dofs; i++)
323  {
324  int i_elts[Max];
325  const int i_gm = e*elt_dofs + i;
326  const int i_L = d_gather_map[i_gm];
327  const int i_offset = d_offsets[i_L];
328  const int i_next_offset = d_offsets[i_L+1];
329  const int i_nbElts = i_next_offset - i_offset;
330  MFEM_ASSERT_KERNEL(
331  i_nbElts <= Max,
332  "The connectivity of this mesh is beyond the max, increase the "
333  "MaxNbNbr variable to comply with your mesh.");
334  for (int e_i = 0; e_i < i_nbElts; ++e_i)
335  {
336  const int i_E = d_indices[i_offset+e_i];
337  i_elts[e_i] = i_E/elt_dofs;
338  }
339  for (int j = 0; j < elt_dofs; j++)
340  {
341  const int j_gm = e*elt_dofs + j;
342  const int j_L = d_gather_map[j_gm];
343  const int j_offset = d_offsets[j_L];
344  const int j_next_offset = d_offsets[j_L+1];
345  const int j_nbElts = j_next_offset - j_offset;
346  if (i_nbElts == 1 || j_nbElts == 1) // no assembly required
347  {
348  GetAndIncrementNnzIndex(i_L, I);
349  }
350  else // assembly required
351  {
352  int j_elts[Max];
353  for (int e_j = 0; e_j < j_nbElts; ++e_j)
354  {
355  const int j_E = d_indices[j_offset+e_j];
356  const int elt = j_E/elt_dofs;
357  j_elts[e_j] = elt;
358  }
359  int min_e = GetMinElt(i_elts, i_nbElts, j_elts, j_nbElts);
360  if (e == min_e) // add the nnz only once
361  {
362  GetAndIncrementNnzIndex(i_L, I);
363  }
364  }
365  }
366  }
367  });
368  // We need to sum the entries of I, we do it on CPU as it is very sequential.
369  auto h_I = mat.HostReadWriteI();
370  const int nTdofs = vd*all_dofs;
371  int sum = 0;
372  for (int i = 0; i < nTdofs; i++)
373  {
374  const int nnz = h_I[i];
375  h_I[i] = sum;
376  sum+=nnz;
377  }
378  h_I[nTdofs] = sum;
379  // We return the number of nnz
380  return h_I[nTdofs];
381 }
382 
384  SparseMatrix &mat) const
385 {
386  static constexpr int Max = MaxNbNbr;
387  const int all_dofs = ndofs;
388  const int vd = vdim;
389  const int elt_dofs = dof;
390  auto I = mat.ReadWriteI();
391  auto J = mat.WriteJ();
392  auto Data = mat.WriteData();
393  auto d_offsets = offsets.Read();
394  auto d_indices = indices.Read();
395  auto d_gather_map = gather_map.Read();
396  auto mat_ea = Reshape(ea_data.Read(), elt_dofs, elt_dofs, ne);
397  MFEM_FORALL(e, ne,
398  {
399  for (int i = 0; i < elt_dofs; i++)
400  {
401  int i_elts[Max];
402  int i_B[Max];
403  const int i_gm = e*elt_dofs + i;
404  const int i_L = d_gather_map[i_gm];
405  const int i_offset = d_offsets[i_L];
406  const int i_next_offset = d_offsets[i_L+1];
407  const int i_nbElts = i_next_offset - i_offset;
408  MFEM_ASSERT_KERNEL(
409  i_nbElts <= Max,
410  "The connectivity of this mesh is beyond the max, increase the "
411  "MaxNbNbr variable to comply with your mesh.");
412  for (int e_i = 0; e_i < i_nbElts; ++e_i)
413  {
414  const int i_E = d_indices[i_offset+e_i];
415  i_elts[e_i] = i_E/elt_dofs;
416  i_B[e_i] = i_E%elt_dofs;
417  }
418  for (int j = 0; j < elt_dofs; j++)
419  {
420  const int j_gm = e*elt_dofs + j;
421  const int j_L = d_gather_map[j_gm];
422  const int j_offset = d_offsets[j_L];
423  const int j_next_offset = d_offsets[j_L+1];
424  const int j_nbElts = j_next_offset - j_offset;
425  if (i_nbElts == 1 || j_nbElts == 1) // no assembly required
426  {
427  const int nnz = GetAndIncrementNnzIndex(i_L, I);
428  J[nnz] = j_L;
429  Data[nnz] = mat_ea(j,i,e);
430  }
431  else // assembly required
432  {
433  int j_elts[Max];
434  int j_B[Max];
435  for (int e_j = 0; e_j < j_nbElts; ++e_j)
436  {
437  const int j_E = d_indices[j_offset+e_j];
438  const int elt = j_E/elt_dofs;
439  j_elts[e_j] = elt;
440  j_B[e_j] = j_E%elt_dofs;
441  }
442  int min_e = GetMinElt(i_elts, i_nbElts, j_elts, j_nbElts);
443  if (e == min_e) // add the nnz only once
444  {
445  double val = 0.0;
446  for (int k = 0; k < i_nbElts; k++)
447  {
448  const int e_i = i_elts[k];
449  const int i_Bloc = i_B[k];
450  for (int l = 0; l < j_nbElts; l++)
451  {
452  const int e_j = j_elts[l];
453  const int j_Bloc = j_B[l];
454  if (e_i == e_j)
455  {
456  val += mat_ea(j_Bloc, i_Bloc, e_i);
457  }
458  }
459  }
460  const int nnz = GetAndIncrementNnzIndex(i_L, I);
461  J[nnz] = j_L;
462  Data[nnz] = val;
463  }
464  }
465  }
466  }
467  });
468  // We need to shift again the entries of I, we do it on CPU as it is very
469  // sequential.
470  auto h_I = mat.HostReadWriteI();
471  const int size = vd*all_dofs;
472  for (int i = 0; i < size; i++)
473  {
474  h_I[size-i] = h_I[size-(i+1)];
475  }
476  h_I[0] = 0;
477 }
478 
480  : ne(fes.GetNE()),
481  vdim(fes.GetVDim()),
482  byvdim(fes.GetOrdering() == Ordering::byVDIM),
483  ndof(ne > 0 ? fes.GetFE(0)->GetDof() : 0),
484  ndofs(fes.GetNDofs())
485 {
486  height = vdim*ne*ndof;
487  width = vdim*ne*ndof;
488 }
489 
490 void L2ElementRestriction::Mult(const Vector &x, Vector &y) const
491 {
492  const int nd = ndof;
493  const int vd = vdim;
494  const bool t = byvdim;
495  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
496  auto d_y = Reshape(y.Write(), nd, vd, ne);
497  MFEM_FORALL(i, ndofs,
498  {
499  const int idx = i;
500  const int dof = idx % nd;
501  const int e = idx / nd;
502  for (int c = 0; c < vd; ++c)
503  {
504  d_y(dof, c, e) = d_x(t?c:idx, t?idx:c);
505  }
506  });
507 }
508 
510 {
511  const int nd = ndof;
512  const int vd = vdim;
513  const bool t = byvdim;
514  auto d_x = Reshape(x.Read(), nd, vd, ne);
515  auto d_y = Reshape(y.Write(), t?vd:ndofs, t?ndofs:vd);
516  MFEM_FORALL(i, ndofs,
517  {
518  const int idx = i;
519  const int dof = idx % nd;
520  const int e = idx / nd;
521  for (int c = 0; c < vd; ++c)
522  {
523  d_y(t?c:idx,t?idx:c) = d_x(dof, c, e);
524  }
525  });
526 }
527 
529 {
530  const int elem_dofs = ndof;
531  const int vd = vdim;
532  auto I = mat.WriteI();
533  const int isize = mat.Height() + 1;
534  const int interior_dofs = ne*elem_dofs*vd;
535  MFEM_FORALL(dof, isize,
536  {
537  I[dof] = dof<interior_dofs ? elem_dofs : 0;
538  });
539 }
540 
541 static MFEM_HOST_DEVICE int AddNnz(const int iE, int *I, const int dofs)
542 {
543  int val = AtomicAdd(I[iE],dofs);
544  return val;
545 }
546 
548  SparseMatrix &mat) const
549 {
550  const int elem_dofs = ndof;
551  const int vd = vdim;
552  auto I = mat.ReadWriteI();
553  auto J = mat.WriteJ();
554  auto Data = mat.WriteData();
555  auto mat_ea = Reshape(ea_data.Read(), elem_dofs, elem_dofs, ne);
556  MFEM_FORALL(iE, ne*elem_dofs*vd,
557  {
558  const int offset = AddNnz(iE,I,elem_dofs);
559  const int e = iE/elem_dofs;
560  const int i = iE%elem_dofs;
561  for (int j = 0; j < elem_dofs; j++)
562  {
563  J[offset+j] = e*elem_dofs+j;
564  Data[offset+j] = mat_ea(j,i,e);
565  }
566  });
567 }
568 
569 /** Return the face degrees of freedom returned in Lexicographic order.
570  Note: Only for quad and hex */
571 void GetFaceDofs(const int dim, const int face_id,
572  const int dof1d, Array<int> &face_map)
573 {
574  switch (dim)
575  {
576  case 1:
577  switch (face_id)
578  {
579  case 0: // WEST
580  face_map[0] = 0;
581  break;
582  case 1: // EAST
583  face_map[0] = dof1d-1;
584  break;
585  }
586  break;
587  case 2:
588  switch (face_id)
589  {
590  case 0: // SOUTH
591  for (int i = 0; i < dof1d; ++i)
592  {
593  face_map[i] = i;
594  }
595  break;
596  case 1: // EAST
597  for (int i = 0; i < dof1d; ++i)
598  {
599  face_map[i] = dof1d-1 + i*dof1d;
600  }
601  break;
602  case 2: // NORTH
603  for (int i = 0; i < dof1d; ++i)
604  {
605  face_map[i] = (dof1d-1)*dof1d + i;
606  }
607  break;
608  case 3: // WEST
609  for (int i = 0; i < dof1d; ++i)
610  {
611  face_map[i] = i*dof1d;
612  }
613  break;
614  }
615  break;
616  case 3:
617  switch (face_id)
618  {
619  case 0: // BOTTOM
620  for (int i = 0; i < dof1d; ++i)
621  {
622  for (int j = 0; j < dof1d; ++j)
623  {
624  face_map[i+j*dof1d] = i + j*dof1d;
625  }
626  }
627  break;
628  case 1: // SOUTH
629  for (int i = 0; i < dof1d; ++i)
630  {
631  for (int j = 0; j < dof1d; ++j)
632  {
633  face_map[i+j*dof1d] = i + j*dof1d*dof1d;
634  }
635  }
636  break;
637  case 2: // EAST
638  for (int i = 0; i < dof1d; ++i)
639  {
640  for (int j = 0; j < dof1d; ++j)
641  {
642  face_map[i+j*dof1d] = dof1d-1 + i*dof1d + j*dof1d*dof1d;
643  }
644  }
645  break;
646  case 3: // NORTH
647  for (int i = 0; i < dof1d; ++i)
648  {
649  for (int j = 0; j < dof1d; ++j)
650  {
651  face_map[i+j*dof1d] = (dof1d-1)*dof1d + i + j*dof1d*dof1d;
652  }
653  }
654  break;
655  case 4: // WEST
656  for (int i = 0; i < dof1d; ++i)
657  {
658  for (int j = 0; j < dof1d; ++j)
659  {
660  face_map[i+j*dof1d] = i*dof1d + j*dof1d*dof1d;
661  }
662  }
663  break;
664  case 5: // TOP
665  for (int i = 0; i < dof1d; ++i)
666  {
667  for (int j = 0; j < dof1d; ++j)
668  {
669  face_map[i+j*dof1d] = (dof1d-1)*dof1d*dof1d + i + j*dof1d;
670  }
671  }
672  break;
673  }
674  break;
675  }
676 }
677 
679  const ElementDofOrdering e_ordering,
680  const FaceType type,
681  bool build)
682  : fes(fes),
683  nf(fes.GetNFbyType(type)),
684  vdim(fes.GetVDim()),
685  byvdim(fes.GetOrdering() == Ordering::byVDIM),
686  face_dofs(nf > 0 ? fes.GetFaceElement(0)->GetDof() : 0),
687  elem_dofs(fes.GetFE(0)->GetDof()),
688  nfdofs(nf*face_dofs),
689  ndofs(fes.GetNDofs()),
690  scatter_indices(nf*face_dofs),
691  gather_offsets(ndofs+1),
692  gather_indices(nf*face_dofs),
693  face_map(face_dofs)
694 {
696  width = fes.GetVSize();
697  if (!build) { return; }
698  if (nf==0) { return; }
699 
700  CheckFESpace(e_ordering);
701 
702  ComputeScatterIndicesAndOffsets(e_ordering, type);
703 
704  ComputeGatherIndices(e_ordering,type);
705 }
706 
708  const ElementDofOrdering e_ordering,
709  const FaceType type)
710  : H1FaceRestriction(fes, e_ordering, type, true)
711 { }
712 
713 void H1FaceRestriction::Mult(const Vector& x, Vector& y) const
714 {
715  if (nf==0) { return; }
716  // Assumes all elements have the same number of dofs
717  const int nface_dofs = face_dofs;
718  const int vd = vdim;
719  const bool t = byvdim;
720  auto d_indices = scatter_indices.Read();
721  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
722  auto d_y = Reshape(y.Write(), nface_dofs, vd, nf);
723  MFEM_FORALL(i, nfdofs,
724  {
725  const int idx = d_indices[i];
726  const int dof = i % nface_dofs;
727  const int face = i / nface_dofs;
728  for (int c = 0; c < vd; ++c)
729  {
730  d_y(dof, c, face) = d_x(t?c:idx, t?idx:c);
731  }
732  });
733 }
734 
736 {
737  if (nf==0) { return; }
738  // Assumes all elements have the same number of dofs
739  const int nface_dofs = face_dofs;
740  const int vd = vdim;
741  const bool t = byvdim;
742  auto d_offsets = gather_offsets.Read();
743  auto d_indices = gather_indices.Read();
744  auto d_x = Reshape(x.Read(), nface_dofs, vd, nf);
745  auto d_y = Reshape(y.ReadWrite(), t?vd:ndofs, t?ndofs:vd);
746  MFEM_FORALL(i, ndofs,
747  {
748  const int offset = d_offsets[i];
749  const int next_offset = d_offsets[i + 1];
750  for (int c = 0; c < vd; ++c)
751  {
752  double dof_value = 0;
753  for (int j = offset; j < next_offset; ++j)
754  {
755  const int idx_j = d_indices[j];
756  dof_value += d_x(idx_j % nface_dofs, c, idx_j / nface_dofs);
757  }
758  d_y(t?c:i,t?i:c) += dof_value;
759  }
760  });
761 }
762 
764 {
765 #ifdef MFEM_USE_MPI
766 
767  // If the underlying finite element space is parallel, ensure the face
768  // neighbor information is generated.
769  if (const ParFiniteElementSpace *pfes
770  = dynamic_cast<const ParFiniteElementSpace*>(&fes))
771  {
772  pfes->GetParMesh()->ExchangeFaceNbrData();
773  }
774 
775 #endif
776 
777 #ifdef MFEM_DEBUG
778  // If fespace == H1
779  const FiniteElement *fe0 = fes.GetFE(0);
780  const TensorBasisElement *tfe = dynamic_cast<const TensorBasisElement*>(fe0);
781  MFEM_VERIFY(tfe != NULL &&
784  "Only Gauss-Lobatto and Bernstein basis are supported in "
785  "H1FaceRestriction.");
786 
787  // Assuming all finite elements are using Gauss-Lobatto.
788  const bool dof_reorder = (e_ordering == ElementDofOrdering::LEXICOGRAPHIC);
789  if (dof_reorder && nf > 0)
790  {
791  for (int f = 0; f < fes.GetNF(); ++f)
792  {
793  const FiniteElement *fe = fes.GetFaceElement(f);
794  const TensorBasisElement* el =
795  dynamic_cast<const TensorBasisElement*>(fe);
796  if (el) { continue; }
797  MFEM_ABORT("Finite element not suitable for lexicographic ordering");
798  }
799  const FiniteElement *fe = fes.GetFaceElement(0);
800  const TensorBasisElement* el =
801  dynamic_cast<const TensorBasisElement*>(fe);
802  const Array<int> &fe_dof_map = el->GetDofMap();
803  MFEM_VERIFY(fe_dof_map.Size() > 0, "invalid dof map");
804  }
805 #endif
806 }
807 
808 void H1FaceRestriction::ComputeScatterIndicesAndOffsets(
809  const ElementDofOrdering ordering,
810  const FaceType type)
811 {
812  Mesh &mesh = *fes.GetMesh();
813 
814  // Initialization of the offsets
815  for (int i = 0; i <= ndofs; ++i)
816  {
817  gather_offsets[i] = 0;
818  }
819 
820  // Computation of scatter indices and offsets
821  int f_ind = 0;
822  for (int f = 0; f < fes.GetNF(); ++f)
823  {
824  Mesh::FaceInformation face = mesh.GetFaceInformation(f);
825  if ( face.IsNonconformingCoarse() )
826  {
827  // We skip nonconforming coarse faces as they are treated
828  // by the corresponding nonconforming fine faces.
829  continue;
830  }
831  else if ( face.IsOfFaceType(type) )
832  {
833  SetFaceDofsScatterIndices(face, f_ind, ordering);
834  f_ind++;
835  }
836  }
837  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
838 
839  // Summation of the offsets
840  for (int i = 1; i <= ndofs; ++i)
841  {
842  gather_offsets[i] += gather_offsets[i - 1];
843  }
844 }
845 
846 void H1FaceRestriction::ComputeGatherIndices(
847  const ElementDofOrdering ordering,
848  const FaceType type)
849 {
850  Mesh &mesh = *fes.GetMesh();
851 
852  // Computation of gather_indices
853  int f_ind = 0;
854  for (int f = 0; f < fes.GetNF(); ++f)
855  {
856  Mesh::FaceInformation face = mesh.GetFaceInformation(f);
857  if ( face.IsNonconformingCoarse() )
858  {
859  // We skip nonconforming coarse faces as they are treated
860  // by the corresponding nonconforming fine faces.
861  continue;
862  }
863  else if ( face.IsOfFaceType(type) )
864  {
865  SetFaceDofsGatherIndices(face, f_ind, ordering);
866  f_ind++;
867  }
868  }
869  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
870 
871  // Reset offsets to their initial value
872  for (int i = ndofs; i > 0; --i)
873  {
874  gather_offsets[i] = gather_offsets[i - 1];
875  }
876  gather_offsets[0] = 0;
877 }
878 
880  const Mesh::FaceInformation &face,
881  const int face_index,
882  const ElementDofOrdering ordering)
883 {
884  MFEM_ASSERT(!(face.IsNonconformingCoarse()),
885  "This method should not be used on nonconforming coarse faces.");
886  MFEM_ASSERT(face.element[0].orientation==0,
887  "FaceRestriction used on degenerated mesh.");
888 
889  const TensorBasisElement* el =
890  dynamic_cast<const TensorBasisElement*>(fes.GetFE(0));
891  const int *dof_map = el->GetDofMap().GetData();
892  const Table& e2dTable = fes.GetElementToDofTable();
893  const int* elem_map = e2dTable.GetJ();
894  const int face_id = face.element[0].local_face_id;
895  const int dim = fes.GetMesh()->Dimension();
896  const int dof1d = fes.GetFE(0)->GetOrder()+1;
897  const int elem_index = face.element[0].index;
898  const bool dof_reorder = (ordering == ElementDofOrdering::LEXICOGRAPHIC);
899  GetFaceDofs(dim, face_id, dof1d, face_map); // Only for quad and hex
900 
901  for (int face_dof = 0; face_dof < face_dofs; ++face_dof)
902  {
903  const int nat_volume_dof = face_map[face_dof];
904  const int volume_dof = (!dof_reorder)?
905  nat_volume_dof:
906  dof_map[nat_volume_dof];
907  const int global_dof = elem_map[elem_index*elem_dofs + volume_dof];
908  const int restriction_dof = face_dofs*face_index + face_dof;
909  scatter_indices[restriction_dof] = global_dof;
910  ++gather_offsets[global_dof + 1];
911  }
912 }
913 
915  const Mesh::FaceInformation &face,
916  const int face_index,
917  const ElementDofOrdering ordering)
918 {
919  MFEM_ASSERT(!(face.IsNonconformingCoarse()),
920  "This method should not be used on nonconforming coarse faces.");
921 
922  const TensorBasisElement* el =
923  dynamic_cast<const TensorBasisElement*>(fes.GetFE(0));
924  const int *dof_map = el->GetDofMap().GetData();
925  const Table& e2dTable = fes.GetElementToDofTable();
926  const int* elem_map = e2dTable.GetJ();
927  const int face_id = face.element[0].local_face_id;
928  const int dim = fes.GetMesh()->Dimension();
929  const int dof1d = fes.GetFE(0)->GetOrder()+1;
930  const int elem_index = face.element[0].index;
931  const bool dof_reorder = (ordering == ElementDofOrdering::LEXICOGRAPHIC);
932  GetFaceDofs(dim, face_id, dof1d, face_map); // Only for quad and hex
933 
934  for (int face_dof = 0; face_dof < face_dofs; ++face_dof)
935  {
936  const int nat_volume_dof = face_map[face_dof];
937  const int volume_dof = (!dof_reorder)?nat_volume_dof:dof_map[nat_volume_dof];
938  const int global_dof = elem_map[elem_index*elem_dofs + volume_dof];
939  const int restriction_dof = face_dofs*face_index + face_dof;
940  gather_indices[gather_offsets[global_dof]++] = restriction_dof;
941  }
942 }
943 
944 static int ToLexOrdering2D(const int face_id, const int size1d, const int i)
945 {
946  if (face_id==2 || face_id==3)
947  {
948  return size1d-1-i;
949  }
950  else
951  {
952  return i;
953  }
954 }
955 
956 static int PermuteFace2D(const int face_id1, const int face_id2,
957  const int orientation,
958  const int size1d, const int index)
959 {
960  int new_index;
961  // Convert from lex ordering
962  if (face_id1==2 || face_id1==3)
963  {
964  new_index = size1d-1-index;
965  }
966  else
967  {
968  new_index = index;
969  }
970  // Permute based on face orientations
971  if (orientation==1)
972  {
973  new_index = size1d-1-new_index;
974  }
975  return ToLexOrdering2D(face_id2, size1d, new_index);
976 }
977 
978 static int ToLexOrdering3D(const int face_id, const int size1d, const int i,
979  const int j)
980 {
981  if (face_id==2 || face_id==1 || face_id==5)
982  {
983  return i + j*size1d;
984  }
985  else if (face_id==3 || face_id==4)
986  {
987  return (size1d-1-i) + j*size1d;
988  }
989  else // face_id==0
990  {
991  return i + (size1d-1-j)*size1d;
992  }
993 }
994 
995 static int PermuteFace3D(const int face_id1, const int face_id2,
996  const int orientation,
997  const int size1d, const int index)
998 {
999  int i=0, j=0, new_i=0, new_j=0;
1000  i = index%size1d;
1001  j = index/size1d;
1002  // Convert from lex ordering
1003  if (face_id1==3 || face_id1==4)
1004  {
1005  i = size1d-1-i;
1006  }
1007  else if (face_id1==0)
1008  {
1009  j = size1d-1-j;
1010  }
1011  // Permute based on face orientations
1012  switch (orientation)
1013  {
1014  case 0:
1015  new_i = i;
1016  new_j = j;
1017  break;
1018  case 1:
1019  new_i = j;
1020  new_j = i;
1021  break;
1022  case 2:
1023  new_i = j;
1024  new_j = (size1d-1-i);
1025  break;
1026  case 3:
1027  new_i = (size1d-1-i);
1028  new_j = j;
1029  break;
1030  case 4:
1031  new_i = (size1d-1-i);
1032  new_j = (size1d-1-j);
1033  break;
1034  case 5:
1035  new_i = (size1d-1-j);
1036  new_j = (size1d-1-i);
1037  break;
1038  case 6:
1039  new_i = (size1d-1-j);
1040  new_j = i;
1041  break;
1042  case 7:
1043  new_i = i;
1044  new_j = (size1d-1-j);
1045  break;
1046  }
1047  return ToLexOrdering3D(face_id2, size1d, new_i, new_j);
1048 }
1049 
1050 // Permute dofs or quads on a face for e2 to match with the ordering of e1
1051 int PermuteFaceL2(const int dim, const int face_id1,
1052  const int face_id2, const int orientation,
1053  const int size1d, const int index)
1054 {
1055  switch (dim)
1056  {
1057  case 1:
1058  return 0;
1059  case 2:
1060  return PermuteFace2D(face_id1, face_id2, orientation, size1d, index);
1061  case 3:
1062  return PermuteFace3D(face_id1, face_id2, orientation, size1d, index);
1063  default:
1064  MFEM_ABORT("Unsupported dimension.");
1065  return 0;
1066  }
1067 }
1068 
1070  const ElementDofOrdering e_ordering,
1071  const FaceType type,
1072  const L2FaceValues m,
1073  bool build)
1074  : fes(fes),
1075  nf(fes.GetNFbyType(type)),
1076  ne(fes.GetNE()),
1077  vdim(fes.GetVDim()),
1078  byvdim(fes.GetOrdering() == Ordering::byVDIM),
1079  face_dofs(nf > 0 ?
1080  fes.GetTraceElement(0, fes.GetMesh()->GetFaceBaseGeometry(0))->GetDof()
1081  : 0),
1082  elem_dofs(fes.GetFE(0)->GetDof()),
1083  nfdofs(nf*face_dofs),
1084  ndofs(fes.GetNDofs()),
1085  type(type),
1086  m(m),
1087  scatter_indices1(nf*face_dofs),
1088  scatter_indices2(m==L2FaceValues::DoubleValued?nf*face_dofs:0),
1089  gather_offsets(ndofs+1),
1090  gather_indices((m==L2FaceValues::DoubleValued? 2 : 1)*nf*face_dofs),
1091  face_map(face_dofs)
1092 {
1094  width = fes.GetVSize();
1095  if (!build) { return; }
1096 
1097  CheckFESpace(e_ordering);
1098 
1099  ComputeScatterIndicesAndOffsets(e_ordering,type);
1100 
1101  ComputeGatherIndices(e_ordering, type);
1102 }
1103 
1105  const ElementDofOrdering e_ordering,
1106  const FaceType type,
1107  const L2FaceValues m)
1108  : L2FaceRestriction(fes, e_ordering, type, m, true)
1109 { }
1110 
1112  Vector& y) const
1113 {
1114  MFEM_ASSERT(
1116  "This method should be called when m == L2FaceValues::SingleValued.");
1117  // Assumes all elements have the same number of dofs
1118  const int nface_dofs = face_dofs;
1119  const int vd = vdim;
1120  const bool t = byvdim;
1121  auto d_indices1 = scatter_indices1.Read();
1122  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
1123  auto d_y = Reshape(y.Write(), nface_dofs, vd, nf);
1124  MFEM_FORALL(i, nfdofs,
1125  {
1126  const int dof = i % nface_dofs;
1127  const int face = i / nface_dofs;
1128  const int idx1 = d_indices1[i];
1129  for (int c = 0; c < vd; ++c)
1130  {
1131  d_y(dof, c, face) = d_x(t?c:idx1, t?idx1:c);
1132  }
1133  });
1134 }
1135 
1137  Vector& y) const
1138 {
1139  MFEM_ASSERT(
1141  "This method should be called when m == L2FaceValues::DoubleValued.");
1142  // Assumes all elements have the same number of dofs
1143  const int nface_dofs = face_dofs;
1144  const int vd = vdim;
1145  const bool t = byvdim;
1146  auto d_indices1 = scatter_indices1.Read();
1147  auto d_indices2 = scatter_indices2.Read();
1148  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
1149  auto d_y = Reshape(y.Write(), nface_dofs, vd, 2, nf);
1150  MFEM_FORALL(i, nfdofs,
1151  {
1152  const int dof = i % nface_dofs;
1153  const int face = i / nface_dofs;
1154  const int idx1 = d_indices1[i];
1155  for (int c = 0; c < vd; ++c)
1156  {
1157  d_y(dof, c, 0, face) = d_x(t?c:idx1, t?idx1:c);
1158  }
1159  const int idx2 = d_indices2[i];
1160  for (int c = 0; c < vd; ++c)
1161  {
1162  d_y(dof, c, 1, face) = idx2==-1 ? 0.0 : d_x(t?c:idx2, t?idx2:c);
1163  }
1164  });
1165 }
1166 
1167 void L2FaceRestriction::Mult(const Vector& x, Vector& y) const
1168 {
1169  if (nf==0) { return; }
1171  {
1173  }
1174  else
1175  {
1177  }
1178 }
1179 
1181  const Vector& x, Vector& y) const
1182 {
1183  // Assumes all elements have the same number of dofs
1184  const int nface_dofs = face_dofs;
1185  const int vd = vdim;
1186  const bool t = byvdim;
1187  auto d_offsets = gather_offsets.Read();
1188  auto d_indices = gather_indices.Read();
1189  auto d_x = Reshape(x.Read(), nface_dofs, vd, nf);
1190  auto d_y = Reshape(y.ReadWrite(), t?vd:ndofs, t?ndofs:vd);
1191  MFEM_FORALL(i, ndofs,
1192  {
1193  const int offset = d_offsets[i];
1194  const int next_offset = d_offsets[i + 1];
1195  for (int c = 0; c < vd; ++c)
1196  {
1197  double dof_value = 0;
1198  for (int j = offset; j < next_offset; ++j)
1199  {
1200  int idx_j = d_indices[j];
1201  dof_value += d_x(idx_j % nface_dofs, c, idx_j / nface_dofs);
1202  }
1203  d_y(t?c:i,t?i:c) += dof_value;
1204  }
1205  });
1206 }
1207 
1209  const Vector& x, Vector& y) const
1210 {
1211  // Assumes all elements have the same number of dofs
1212  const int nface_dofs = face_dofs;
1213  const int vd = vdim;
1214  const bool t = byvdim;
1215  const int dofs = nfdofs;
1216  auto d_offsets = gather_offsets.Read();
1217  auto d_indices = gather_indices.Read();
1218  auto d_x = Reshape(x.Read(), nface_dofs, vd, 2, nf);
1219  auto d_y = Reshape(y.ReadWrite(), t?vd:ndofs, t?ndofs:vd);
1220  MFEM_FORALL(i, ndofs,
1221  {
1222  const int offset = d_offsets[i];
1223  const int next_offset = d_offsets[i + 1];
1224  for (int c = 0; c < vd; ++c)
1225  {
1226  double dof_value = 0;
1227  for (int j = offset; j < next_offset; ++j)
1228  {
1229  int idx_j = d_indices[j];
1230  bool isE1 = idx_j < dofs;
1231  idx_j = isE1 ? idx_j : idx_j - dofs;
1232  dof_value += isE1 ?
1233  d_x(idx_j % nface_dofs, c, 0, idx_j / nface_dofs)
1234  :d_x(idx_j % nface_dofs, c, 1, idx_j / nface_dofs);
1235  }
1236  d_y(t?c:i,t?i:c) += dof_value;
1237  }
1238  });
1239 }
1240 
1242 {
1243  if (nf==0) { return; }
1245  {
1247  }
1248  else
1249  {
1251  }
1252 }
1253 
1255  const bool keep_nbr_block) const
1256 {
1257  const int nface_dofs = face_dofs;
1258  auto d_indices1 = scatter_indices1.Read();
1259  auto d_indices2 = scatter_indices2.Read();
1260  auto I = mat.ReadWriteI();
1261  MFEM_FORALL(fdof, nf*nface_dofs,
1262  {
1263  const int iE1 = d_indices1[fdof];
1264  const int iE2 = d_indices2[fdof];
1265  AddNnz(iE1,I,nface_dofs);
1266  AddNnz(iE2,I,nface_dofs);
1267  });
1268 }
1269 
1271  SparseMatrix &mat,
1272  const bool keep_nbr_block) const
1273 {
1274  const int nface_dofs = face_dofs;
1275  auto d_indices1 = scatter_indices1.Read();
1276  auto d_indices2 = scatter_indices2.Read();
1277  auto I = mat.ReadWriteI();
1278  auto mat_fea = Reshape(fea_data.Read(), nface_dofs, nface_dofs, 2, nf);
1279  auto J = mat.WriteJ();
1280  auto Data = mat.WriteData();
1281  MFEM_FORALL(fdof, nf*nface_dofs,
1282  {
1283  const int f = fdof/nface_dofs;
1284  const int iF = fdof%nface_dofs;
1285  const int iE1 = d_indices1[f*nface_dofs+iF];
1286  const int iE2 = d_indices2[f*nface_dofs+iF];
1287  const int offset1 = AddNnz(iE1,I,nface_dofs);
1288  const int offset2 = AddNnz(iE2,I,nface_dofs);
1289  for (int jF = 0; jF < nface_dofs; jF++)
1290  {
1291  const int jE1 = d_indices1[f*nface_dofs+jF];
1292  const int jE2 = d_indices2[f*nface_dofs+jF];
1293  J[offset2+jF] = jE1;
1294  J[offset1+jF] = jE2;
1295  Data[offset2+jF] = mat_fea(jF,iF,0,f);
1296  Data[offset1+jF] = mat_fea(jF,iF,1,f);
1297  }
1298  });
1299 }
1300 
1302  Vector &ea_data) const
1303 {
1304  const int nface_dofs = face_dofs;
1305  const int nelem_dofs = elem_dofs;
1306  const int NE = ne;
1308  {
1309  auto d_indices1 = scatter_indices1.Read();
1310  auto d_indices2 = scatter_indices2.Read();
1311  auto mat_fea = Reshape(fea_data.Read(), nface_dofs, nface_dofs, 2, nf);
1312  auto mat_ea = Reshape(ea_data.ReadWrite(), nelem_dofs, nelem_dofs, ne);
1313  MFEM_FORALL(f, nf,
1314  {
1315  const int e1 = d_indices1[f*nface_dofs]/nelem_dofs;
1316  const int e2 = d_indices2[f*nface_dofs]/nelem_dofs;
1317  for (int j = 0; j < nface_dofs; j++)
1318  {
1319  const int jB1 = d_indices1[f*nface_dofs+j]%nelem_dofs;
1320  for (int i = 0; i < nface_dofs; i++)
1321  {
1322  const int iB1 = d_indices1[f*nface_dofs+i]%nelem_dofs;
1323  AtomicAdd(mat_ea(iB1,jB1,e1), mat_fea(i,j,0,f));
1324  }
1325  }
1326  if (e2 < NE)
1327  {
1328  for (int j = 0; j < nface_dofs; j++)
1329  {
1330  const int jB2 = d_indices2[f*nface_dofs+j]%nelem_dofs;
1331  for (int i = 0; i < nface_dofs; i++)
1332  {
1333  const int iB2 = d_indices2[f*nface_dofs+i]%nelem_dofs;
1334  AtomicAdd(mat_ea(iB2,jB2,e2), mat_fea(i,j,1,f));
1335  }
1336  }
1337  }
1338  });
1339  }
1340  else
1341  {
1342  auto d_indices = scatter_indices1.Read();
1343  auto mat_fea = Reshape(fea_data.Read(), nface_dofs, nface_dofs, nf);
1344  auto mat_ea = Reshape(ea_data.ReadWrite(), nelem_dofs, nelem_dofs, ne);
1345  MFEM_FORALL(f, nf,
1346  {
1347  const int e = d_indices[f*nface_dofs]/nelem_dofs;
1348  for (int j = 0; j < nface_dofs; j++)
1349  {
1350  const int jE = d_indices[f*nface_dofs+j]%nelem_dofs;
1351  for (int i = 0; i < nface_dofs; i++)
1352  {
1353  const int iE = d_indices[f*nface_dofs+i]%nelem_dofs;
1354  AtomicAdd(mat_ea(iE,jE,e), mat_fea(i,j,f));
1355  }
1356  }
1357  });
1358  }
1359 }
1360 
1362 {
1363 #ifdef MFEM_USE_MPI
1364 
1365  // If the underlying finite element space is parallel, ensure the face
1366  // neighbor information is generated.
1367  if (const ParFiniteElementSpace *pfes
1368  = dynamic_cast<const ParFiniteElementSpace*>(&fes))
1369  {
1370  pfes->GetParMesh()->ExchangeFaceNbrData();
1371  }
1372 
1373 #endif
1374 
1375 #ifdef MFEM_DEBUG
1376  // If fespace == L2
1377  const FiniteElement *fe0 = fes.GetFE(0);
1378  const TensorBasisElement *tfe = dynamic_cast<const TensorBasisElement*>(fe0);
1379  MFEM_VERIFY(tfe != NULL &&
1382  "Only Gauss-Lobatto and Bernstein basis are supported in "
1383  "L2FaceRestriction.");
1384  if (nf==0) { return; }
1385  const bool dof_reorder = (e_ordering == ElementDofOrdering::LEXICOGRAPHIC);
1386  if (!dof_reorder)
1387  {
1388  MFEM_ABORT("Non-Tensor L2FaceRestriction not yet implemented.");
1389  }
1390  if (dof_reorder && nf > 0)
1391  {
1392  for (int f = 0; f < fes.GetNF(); ++f)
1393  {
1394  const FiniteElement *fe =
1396  const TensorBasisElement* el =
1397  dynamic_cast<const TensorBasisElement*>(fe);
1398  if (el) { continue; }
1399  MFEM_ABORT("Finite element not suitable for lexicographic ordering");
1400  }
1401  }
1402 #endif
1403 }
1404 
1405 void L2FaceRestriction::ComputeScatterIndicesAndOffsets(
1406  const ElementDofOrdering ordering,
1407  const FaceType face_type)
1408 {
1409  Mesh &mesh = *fes.GetMesh();
1410  // Initialization of the offsets
1411  for (int i = 0; i <= ndofs; ++i)
1412  {
1413  gather_offsets[i] = 0;
1414  }
1415 
1416  // Computation of scatter indices and offsets
1417  int f_ind=0;
1418  for (int f = 0; f < fes.GetNF(); ++f)
1419  {
1420  Mesh::FaceInformation face = mesh.GetFaceInformation(f);
1421  MFEM_ASSERT(!face.IsShared(),
1422  "Unexpected shared face in L2FaceRestriction.");
1423  if ( face.IsOfFaceType(face_type) )
1424  {
1425  SetFaceDofsScatterIndices1(face,f_ind);
1427  {
1428  if ( face_type==FaceType::Interior && face.IsInterior() )
1429  {
1431  }
1432  else if ( face_type==FaceType::Boundary && face.IsBoundary() )
1433  {
1434  SetBoundaryDofsScatterIndices2(face,f_ind);
1435  }
1436  }
1437  f_ind++;
1438  }
1439  }
1440  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
1441 
1442  // Summation of the offsets
1443  for (int i = 1; i <= ndofs; ++i)
1444  {
1445  gather_offsets[i] += gather_offsets[i - 1];
1446  }
1447 }
1448 
1449 void L2FaceRestriction::ComputeGatherIndices(
1450  const ElementDofOrdering ordering,
1451  const FaceType face_type)
1452 {
1453  Mesh &mesh = *fes.GetMesh();
1454  // Computation of gather_indices
1455  int f_ind = 0;
1456  for (int f = 0; f < fes.GetNF(); ++f)
1457  {
1458  Mesh::FaceInformation face = mesh.GetFaceInformation(f);
1459  MFEM_ASSERT(!face.IsShared(),
1460  "Unexpected shared face in L2FaceRestriction.");
1461  if ( face.IsOfFaceType(face_type) )
1462  {
1463  SetFaceDofsGatherIndices1(face,f_ind);
1464  if ( m==L2FaceValues::DoubleValued &&
1465  face_type==FaceType::Interior &&
1466  face.IsLocal())
1467  {
1469  }
1470  f_ind++;
1471  }
1472  }
1473  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
1474 
1475  // Reset offsets to their correct value
1476  for (int i = ndofs; i > 0; --i)
1477  {
1478  gather_offsets[i] = gather_offsets[i - 1];
1479  }
1480  gather_offsets[0] = 0;
1481 }
1482 
1484  const Mesh::FaceInformation &face,
1485  const int face_index)
1486 {
1487  MFEM_ASSERT(!(face.IsNonconformingCoarse()),
1488  "This method should not be used on nonconforming coarse faces.");
1489  const Table& e2dTable = fes.GetElementToDofTable();
1490  const int* elem_map = e2dTable.GetJ();
1491  const int face_id1 = face.element[0].local_face_id;
1492  const int dim = fes.GetMesh()->Dimension();
1493  const int dof1d = fes.GetFE(0)->GetOrder()+1;
1494  const int elem_index = face.element[0].index;
1495  GetFaceDofs(dim, face_id1, dof1d, face_map); // Only for quad and hex
1496 
1497  for (int face_dof_elem1 = 0; face_dof_elem1 < face_dofs; ++face_dof_elem1)
1498  {
1499  const int volume_dof_elem1 = face_map[face_dof_elem1];
1500  const int global_dof_elem1 = elem_map[elem_index*elem_dofs + volume_dof_elem1];
1501  const int restriction_dof_elem1 = face_dofs*face_index + face_dof_elem1;
1502  scatter_indices1[restriction_dof_elem1] = global_dof_elem1;
1503  ++gather_offsets[global_dof_elem1 + 1];
1504  }
1505 }
1506 
1508  const Mesh::FaceInformation &face,
1509  const int face_index)
1510 {
1511  MFEM_ASSERT(face.IsLocal(),
1512  "This method should only be used on local faces.");
1513  const Table& e2dTable = fes.GetElementToDofTable();
1514  const int* elem_map = e2dTable.GetJ();
1515  const int elem_index = face.element[1].index;
1516  const int face_id1 = face.element[0].local_face_id;
1517  const int face_id2 = face.element[1].local_face_id;
1518  const int orientation = face.element[1].orientation;
1519  const int dim = fes.GetMesh()->Dimension();
1520  const int dof1d = fes.GetFE(0)->GetOrder()+1;
1521  GetFaceDofs(dim, face_id2, dof1d, face_map); // Only for quad and hex
1522 
1523  for (int face_dof_elem1 = 0; face_dof_elem1 < face_dofs; ++face_dof_elem1)
1524  {
1525  const int face_dof_elem2 = PermuteFaceL2(dim, face_id1, face_id2,
1526  orientation, dof1d,
1527  face_dof_elem1);
1528  const int volume_dof_elem2 = face_map[face_dof_elem2];
1529  const int global_dof_elem2 = elem_map[elem_index*elem_dofs + volume_dof_elem2];
1530  const int restriction_dof_elem2 = face_dofs*face_index + face_dof_elem1;
1531  scatter_indices2[restriction_dof_elem2] = global_dof_elem2;
1532  ++gather_offsets[global_dof_elem2 + 1];
1533  }
1534 }
1535 
1537  const Mesh::FaceInformation &face,
1538  const int face_index)
1539 {
1540 #ifdef MFEM_USE_MPI
1541  MFEM_ASSERT(face.IsShared(),
1542  "This method should only be used on shared faces.");
1543  const int elem_index = face.element[1].index;
1544  const int face_id1 = face.element[0].local_face_id;
1545  const int face_id2 = face.element[1].local_face_id;
1546  const int orientation = face.element[1].orientation;
1547  const int dim = fes.GetMesh()->Dimension();
1548  const int dof1d = fes.GetFE(0)->GetOrder()+1;
1549  GetFaceDofs(dim, face_id2, dof1d, face_map); // Only for quad and hex
1550  Array<int> face_nbr_dofs;
1551  const ParFiniteElementSpace &pfes =
1552  static_cast<const ParFiniteElementSpace&>(this->fes);
1553  pfes.GetFaceNbrElementVDofs(elem_index, face_nbr_dofs);
1554 
1555  for (int face_dof_elem1 = 0; face_dof_elem1 < face_dofs; ++face_dof_elem1)
1556  {
1557  const int face_dof_elem2 = PermuteFaceL2(dim, face_id1, face_id2,
1558  orientation, dof1d, face_dof_elem1);
1559  const int volume_dof_elem2 = face_map[face_dof_elem2];
1560  const int global_dof_elem2 = face_nbr_dofs[volume_dof_elem2];
1561  const int restriction_dof_elem2 = face_dofs*face_index + face_dof_elem1;
1562  // Trick to differentiate dof location inter/shared
1563  scatter_indices2[restriction_dof_elem2] = ndofs+global_dof_elem2;
1564  }
1565 #endif
1566 }
1567 
1569  const Mesh::FaceInformation &face,
1570  const int face_index)
1571 {
1572  MFEM_ASSERT(face.IsBoundary(),
1573  "This method should only be used on boundary faces.");
1574 
1575  for (int d = 0; d < face_dofs; ++d)
1576  {
1577  const int restriction_dof_elem2 = face_dofs*face_index + d;
1578  scatter_indices2[restriction_dof_elem2] = -1;
1579  }
1580 }
1581 
1583  const Mesh::FaceInformation &face,
1584  const int face_index)
1585 {
1586  MFEM_ASSERT(!(face.IsNonconformingCoarse()),
1587  "This method should not be used on nonconforming coarse faces.");
1588  const Table& e2dTable = fes.GetElementToDofTable();
1589  const int* elem_map = e2dTable.GetJ();
1590  const int face_id1 = face.element[0].local_face_id;
1591  const int dim = fes.GetMesh()->Dimension();
1592  const int dof1d = fes.GetFE(0)->GetOrder()+1;
1593  const int elem_index = face.element[0].index;
1594  GetFaceDofs(dim, face_id1, dof1d, face_map); // Only for quad and hex
1595 
1596  for (int face_dof_elem1 = 0; face_dof_elem1 < face_dofs; ++face_dof_elem1)
1597  {
1598  const int volume_dof_elem1 = face_map[face_dof_elem1];
1599  const int global_dof_elem1 = elem_map[elem_index*elem_dofs + volume_dof_elem1];
1600  const int restriction_dof_elem1 = face_dofs*face_index + face_dof_elem1;
1601  // We don't shift restriction_dof_elem1 to express that it's elem1 of the face
1602  gather_indices[gather_offsets[global_dof_elem1]++] = restriction_dof_elem1;
1603  }
1604 }
1605 
1607  const Mesh::FaceInformation &face,
1608  const int face_index)
1609 {
1610  MFEM_ASSERT(face.IsLocal(),
1611  "This method should only be used on local faces.");
1612  const Table& e2dTable = fes.GetElementToDofTable();
1613  const int* elem_map = e2dTable.GetJ();
1614  const int elem_index = face.element[1].index;
1615  const int face_id1 = face.element[0].local_face_id;
1616  const int face_id2 = face.element[1].local_face_id;
1617  const int orientation = face.element[1].orientation;
1618  const int dim = fes.GetMesh()->Dimension();
1619  const int dof1d = fes.GetFE(0)->GetOrder()+1;
1620  GetFaceDofs(dim, face_id2, dof1d, face_map); // Only for quad and hex
1621 
1622  for (int face_dof_elem1 = 0; face_dof_elem1 < face_dofs; ++face_dof_elem1)
1623  {
1624  const int face_dof_elem2 = PermuteFaceL2(dim, face_id1, face_id2,
1625  orientation, dof1d,
1626  face_dof_elem1);
1627  const int volume_dof_elem2 = face_map[face_dof_elem2];
1628  const int global_dof_elem2 = elem_map[elem_index*elem_dofs + volume_dof_elem2];
1629  const int restriction_dof_elem2 = face_dofs*face_index + face_dof_elem1;
1630  // We shift restriction_dof_elem2 to express that it's elem2 of the face
1631  gather_indices[gather_offsets[global_dof_elem2]++] = nfdofs +
1632  restriction_dof_elem2;
1633  }
1634 }
1635 
1637  ElementDofOrdering ordering,
1638  FaceType type)
1639  : fes(fes),
1640  ordering(ordering),
1641  interp_config(type==FaceType::Interior ? fes.GetNFbyType(type) : 0),
1642  nc_cpt(0)
1643 { }
1644 
1646  const Mesh::FaceInformation &face,
1647  int face_index)
1648 {
1649  interp_config[face_index] = InterpConfig();
1650 }
1651 
1653  const Mesh::FaceInformation &face,
1654  int face_index)
1655 {
1656  MFEM_ASSERT(!face.IsConforming(),
1657  "Registering face as nonconforming even though it is not.");
1658  const DenseMatrix* ptMat = face.point_matrix;
1659  // In the case of nonconforming slave shared face the master face is elem1.
1660  const int master_side =
1662  const int face_key = (master_side == 0 ? 1000 : 0) +
1663  face.element[0].local_face_id +
1664  6*face.element[1].local_face_id +
1665  36*face.element[1].orientation ;
1666  // Unfortunately we can't trust unicity of the ptMat to identify the transformation.
1667  Key key(ptMat, face_key);
1668  auto itr = interp_map.find(key);
1669  if ( itr == interp_map.end() )
1670  {
1671  const DenseMatrix* interpolator =
1672  GetCoarseToFineInterpolation(face,ptMat);
1673  interp_map[key] = {nc_cpt, interpolator};
1674  interp_config[face_index] = {master_side, nc_cpt};
1675  nc_cpt++;
1676  }
1677  else
1678  {
1679  interp_config[face_index] = {master_side, itr->second.first};
1680  }
1681 }
1682 
1683 const DenseMatrix* InterpolationManager::GetCoarseToFineInterpolation(
1684  const Mesh::FaceInformation &face,
1685  const DenseMatrix* ptMat)
1686 {
1687  MFEM_VERIFY(ordering == ElementDofOrdering::LEXICOGRAPHIC,
1688  "The following interpolation operator is only implemented for"
1689  "lexicographic ordering.");
1690  MFEM_VERIFY(!face.IsConforming(),
1691  "This method should not be called on conforming faces.")
1692  const int face_id1 = face.element[0].local_face_id;
1693  const int face_id2 = face.element[1].local_face_id;
1694 
1695  const bool is_ghost_slave =
1697  const int master_face_id = is_ghost_slave ? face_id1 : face_id2;
1698 
1699  // Computation of the interpolation matrix from master
1700  // (coarse) face to slave (fine) face.
1701  // Assumes all trace elements are the same.
1702  const FiniteElement *trace_fe =
1704  const int face_dofs = trace_fe->GetDof();
1705  const TensorBasisElement* el =
1706  dynamic_cast<const TensorBasisElement*>(trace_fe);
1707  const auto dof_map = el->GetDofMap();
1708  DenseMatrix* interpolator = new DenseMatrix(face_dofs,face_dofs);
1709  Vector shape(face_dofs);
1710 
1712  isotr.SetIdentityTransformation(trace_fe->GetGeomType());
1713  isotr.SetPointMat(*ptMat);
1714  DenseMatrix& trans_pt_mat = isotr.GetPointMat();
1715  // PointMatrix needs to be flipped in 2D
1716  if ( trace_fe->GetGeomType()==Geometry::SEGMENT && !is_ghost_slave )
1717  {
1718  std::swap(trans_pt_mat(0,0),trans_pt_mat(0,1));
1719  }
1720  DenseMatrix native_interpolator(face_dofs,face_dofs);
1721  trace_fe->GetLocalInterpolation(isotr, native_interpolator);
1722  const int dim = trace_fe->GetDim()+1;
1723  const int dof1d = trace_fe->GetOrder()+1;
1724  const int orientation = face.element[1].orientation;
1725  for (int i = 0; i < face_dofs; i++)
1726  {
1727  const int ni = (dof_map.Size()==0) ? i : dof_map[i];
1728  int li = ToLexOrdering(dim, master_face_id, dof1d, i);
1729  if ( !is_ghost_slave )
1730  {
1731  // master side is elem 2, so we permute to order dofs as elem 1.
1732  li = PermuteFaceL2(dim, face_id2, face_id1,
1733  orientation, dof1d, li);
1734  }
1735  for (int j = 0; j < face_dofs; j++)
1736  {
1737  int lj = ToLexOrdering(dim, master_face_id, dof1d, j);
1738  if ( !is_ghost_slave )
1739  {
1740  // master side is elem 2, so we permute to order dofs as elem 1.
1741  lj = PermuteFaceL2(dim, face_id2, face_id1,
1742  orientation, dof1d, lj);
1743  }
1744  const int nj = (dof_map.Size()==0) ? j : dof_map[j];
1745  (*interpolator)(li,lj) = native_interpolator(ni,nj);
1746  }
1747  }
1748  return interpolator;
1749 }
1750 
1752 {
1753  // Assumes all trace elements are the same.
1754  const FiniteElement *trace_fe =
1756  const int face_dofs = trace_fe->GetDof();
1757  const int nc_size = interp_map.size();
1758  MFEM_VERIFY(nc_cpt==nc_size, "Unexpected number of interpolators.");
1759  interpolators.SetSize(face_dofs*face_dofs*nc_size);
1760  auto d_interp = Reshape(interpolators.HostWrite(),face_dofs,face_dofs,nc_size);
1761  for (auto val : interp_map)
1762  {
1763  const int idx = val.second.first;
1764  const DenseMatrix &interpolator = *val.second.second;
1765  for (int i = 0; i < face_dofs; i++)
1766  {
1767  for (int j = 0; j < face_dofs; j++)
1768  {
1769  d_interp(i,j,idx) = interpolator(i,j);
1770  }
1771  }
1772  delete val.second.second;
1773  }
1774  interp_map.clear();
1775 }
1776 
1778  const ElementDofOrdering ordering,
1779  const FaceType type,
1780  const L2FaceValues m,
1781  bool build)
1782  : L2FaceRestriction(fes, ordering, type, m, false),
1783  interpolations(fes, ordering, type)
1784 {
1785  if (!build) { return; }
1786  x_interp.UseDevice(true);
1787 
1788  CheckFESpace(ordering);
1789 
1790  ComputeScatterIndicesAndOffsets(ordering, type);
1791 
1792  ComputeGatherIndices(ordering, type);
1793 }
1794 
1796  const ElementDofOrdering ordering,
1797  const FaceType type,
1798  const L2FaceValues m)
1799  : NCL2FaceRestriction(fes, ordering, type, m, true)
1800 { }
1801 
1803  const Vector& x, Vector& y) const
1804 {
1805  // Assumes all elements have the same number of dofs
1806  const int nface_dofs = face_dofs;
1807  const int vd = vdim;
1808  const bool t = byvdim;
1809  auto d_indices1 = scatter_indices1.Read();
1810  auto d_indices2 = scatter_indices2.Read();
1811  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
1812  auto d_y = Reshape(y.Write(), nface_dofs, vd, 2, nf);
1813  auto interp_config_ptr = interpolations.GetFaceInterpConfig().Read();
1814  const int nc_size = interpolations.GetNumInterpolators();
1815  auto d_interp = Reshape(interpolations.GetInterpolators().Read(),
1816  nface_dofs, nface_dofs, nc_size);
1817  static constexpr int max_nd = 16*16;
1818  MFEM_VERIFY(nface_dofs<=max_nd, "Too many degrees of freedom.");
1819  MFEM_FORALL_3D(face, nf, nface_dofs, 1, 1,
1820  {
1821  MFEM_SHARED double dof_values[max_nd];
1822  const InterpConfig conf = interp_config_ptr[face];
1823  const int master_side = conf.master_side;
1824  const int interp_index = conf.index;
1825  for (int side = 0; side < 2; side++)
1826  {
1827  if ( !conf.is_non_conforming || side!=master_side )
1828  {
1829  // No interpolation needed
1830  MFEM_FOREACH_THREAD(dof,x,nface_dofs)
1831  {
1832  const int i = face*nface_dofs + dof;
1833  const int idx = side==0 ? d_indices1[i] : d_indices2[i];
1834  for (int c = 0; c < vd; ++c)
1835  {
1836  d_y(dof, c, side, face) = d_x(t?c:idx, t?idx:c);
1837  }
1838  }
1839  }
1840  else // Interpolation from coarse to fine
1841  {
1842  for (int c = 0; c < vd; ++c)
1843  {
1844  MFEM_FOREACH_THREAD(dof,x,nface_dofs)
1845  {
1846  const int i = face*nface_dofs + dof;
1847  const int idx = side==0 ? d_indices1[i] : d_indices2[i];
1848  dof_values[dof] = d_x(t?c:idx, t?idx:c);
1849  }
1850  MFEM_SYNC_THREAD;
1851  MFEM_FOREACH_THREAD(dof_out,x,nface_dofs)
1852  {
1853  double res = 0.0;
1854  for (int dof_in = 0; dof_in<nface_dofs; dof_in++)
1855  {
1856  res += d_interp(dof_out, dof_in, interp_index)*dof_values[dof_in];
1857  }
1858  d_y(dof_out, c, side, face) = res;
1859  }
1860  MFEM_SYNC_THREAD;
1861  }
1862  }
1863  }
1864  });
1865 }
1866 
1867 void NCL2FaceRestriction::Mult(const Vector& x, Vector& y) const
1868 {
1869  if (nf==0) { return; }
1870  if ( type==FaceType::Interior && m==L2FaceValues::DoubleValued )
1871  {
1872  DoubleValuedNonconformingMult(x, y);
1873  }
1874  else if ( type==FaceType::Boundary && m==L2FaceValues::DoubleValued )
1875  {
1876  DoubleValuedConformingMult(x, y);
1877  }
1878  else // Single valued (assumes no nonconforming master on elem1)
1879  {
1880  SingleValuedConformingMult(x, y);
1881  }
1882 }
1883 
1884 void NCL2FaceRestriction::SingleValuedNonconformingTransposeInterpolation(
1885  const Vector& x) const
1886 {
1887  MFEM_ASSERT(
1888  m == L2FaceValues::SingleValued,
1889  "This method should be called when m == L2FaceValues::SingleValued.");
1890  if (x_interp.Size()==0)
1891  {
1892  x_interp.SetSize(x.Size());
1893  }
1894  x_interp = x;
1895  // Assumes all elements have the same number of dofs
1896  const int nface_dofs = face_dofs;
1897  const int vd = vdim;
1898  // Interpolation
1899  auto d_x = Reshape(x_interp.ReadWrite(), nface_dofs, vd, nf);
1900  auto interp_config_ptr = interpolations.GetFaceInterpConfig().Read();
1901  auto interpolators = interpolations.GetInterpolators().Read();
1902  const int nc_size = interpolations.GetNumInterpolators();
1903  auto d_interp = Reshape(interpolators, nface_dofs, nface_dofs, nc_size);
1904  static constexpr int max_nd = 16*16;
1905  MFEM_VERIFY(nface_dofs<=max_nd, "Too many degrees of freedom.");
1906  MFEM_FORALL_3D(face, nf, nface_dofs, 1, 1,
1907  {
1908  MFEM_SHARED double dof_values[max_nd];
1909  const InterpConfig conf = interp_config_ptr[face];
1910  const int master_side = conf.master_side;
1911  const int interp_index = conf.index;
1912  if ( conf.is_non_conforming && master_side==0 )
1913  {
1914  // Interpolation from fine to coarse
1915  for (int c = 0; c < vd; ++c)
1916  {
1917  MFEM_FOREACH_THREAD(dof,x,nface_dofs)
1918  {
1919  dof_values[dof] = d_x(dof, c, face);
1920  }
1921  MFEM_SYNC_THREAD;
1922  MFEM_FOREACH_THREAD(dof_out,x,nface_dofs)
1923  {
1924  double res = 0.0;
1925  for (int dof_in = 0; dof_in<nface_dofs; dof_in++)
1926  {
1927  res += d_interp(dof_in, dof_out, interp_index)*dof_values[dof_in];
1928  }
1929  d_x(dof_out, c, face) = res;
1930  }
1931  MFEM_SYNC_THREAD;
1932  }
1933  }
1934  });
1935 }
1936 
1937 void NCL2FaceRestriction::DoubleValuedNonconformingTransposeInterpolation(
1938  const Vector& x) const
1939 {
1940  MFEM_ASSERT(
1941  m == L2FaceValues::DoubleValued,
1942  "This method should be called when m == L2FaceValues::DoubleValued.");
1943  if (x_interp.Size()==0)
1944  {
1945  x_interp.SetSize(x.Size());
1946  }
1947  x_interp = x;
1948  // Assumes all elements have the same number of dofs
1949  const int nface_dofs = face_dofs;
1950  const int vd = vdim;
1951  // Interpolation
1952  auto d_x = Reshape(x_interp.ReadWrite(), nface_dofs, vd, 2, nf);
1953  auto interp_config_ptr = interpolations.GetFaceInterpConfig().Read();
1954  auto interpolators = interpolations.GetInterpolators().Read();
1955  const int nc_size = interpolations.GetNumInterpolators();
1956  auto d_interp = Reshape(interpolators, nface_dofs, nface_dofs, nc_size);
1957  static constexpr int max_nd = 16*16;
1958  MFEM_VERIFY(nface_dofs<=max_nd, "Too many degrees of freedom.");
1959  MFEM_FORALL_3D(face, nf, nface_dofs, 1, 1,
1960  {
1961  MFEM_SHARED double dof_values[max_nd];
1962  const InterpConfig conf = interp_config_ptr[face];
1963  const int master_side = conf.master_side;
1964  const int interp_index = conf.index;
1965  if ( conf.is_non_conforming )
1966  {
1967  // Interpolation from fine to coarse
1968  for (int c = 0; c < vd; ++c)
1969  {
1970  MFEM_FOREACH_THREAD(dof,x,nface_dofs)
1971  {
1972  dof_values[dof] = d_x(dof, c, master_side, face);
1973  }
1974  MFEM_SYNC_THREAD;
1975  MFEM_FOREACH_THREAD(dof_out,x,nface_dofs)
1976  {
1977  double res = 0.0;
1978  for (int dof_in = 0; dof_in<nface_dofs; dof_in++)
1979  {
1980  res += d_interp(dof_in, dof_out, interp_index)*dof_values[dof_in];
1981  }
1982  d_x(dof_out, c, master_side, face) = res;
1983  }
1984  MFEM_SYNC_THREAD;
1985  }
1986  }
1987  });
1988 }
1989 
1990 void NCL2FaceRestriction::AddMultTranspose(const Vector& x, Vector& y) const
1991 {
1992  if (nf==0) { return; }
1993  if (type==FaceType::Interior)
1994  {
1995  if ( m==L2FaceValues::DoubleValued )
1996  {
1997  DoubleValuedNonconformingTransposeInterpolation(x);
1998  DoubleValuedConformingAddMultTranspose(x_interp, y);
1999  }
2000  else if ( m==L2FaceValues::SingleValued )
2001  {
2002  SingleValuedNonconformingTransposeInterpolation(x);
2003  SingleValuedConformingAddMultTranspose(x_interp, y);
2004  }
2005  }
2006  else
2007  {
2008  if ( m==L2FaceValues::DoubleValued )
2009  {
2010  DoubleValuedConformingAddMultTranspose(x, y);
2011  }
2012  else if ( m==L2FaceValues::SingleValued )
2013  {
2014  SingleValuedConformingAddMultTranspose(x, y);
2015  }
2016  }
2017 }
2018 
2019 void NCL2FaceRestriction::FillI(SparseMatrix &mat,
2020  const bool keep_nbr_block) const
2021 {
2022  MFEM_ABORT("Not yet implemented.");
2023 }
2024 
2025 void NCL2FaceRestriction::FillJAndData(const Vector &ea_data,
2026  SparseMatrix &mat,
2027  const bool keep_nbr_block) const
2028 {
2029  MFEM_ABORT("Not yet implemented.");
2030 }
2031 
2032 void NCL2FaceRestriction::AddFaceMatricesToElementMatrices(
2033  const Vector &fea_data,
2034  Vector &ea_data)
2035 const
2036 {
2037  MFEM_ABORT("Not yet implemented.");
2038 }
2039 
2040 int ToLexOrdering(const int dim, const int face_id, const int size1d,
2041  const int index)
2042 {
2043  switch (dim)
2044  {
2045  case 1:
2046  return 0;
2047  case 2:
2048  return ToLexOrdering2D(face_id, size1d, index);
2049  case 3:
2050  return ToLexOrdering3D(face_id, size1d, index%size1d, index/size1d);
2051  default:
2052  MFEM_ABORT("Unsupported dimension.");
2053  return 0;
2054  }
2055 }
2056 
2057 void NCL2FaceRestriction::ComputeScatterIndicesAndOffsets(
2058  const ElementDofOrdering ordering,
2059  const FaceType type)
2060 {
2061  Mesh &mesh = *fes.GetMesh();
2062 
2063  // Initialization of the offsets
2064  for (int i = 0; i <= ndofs; ++i)
2065  {
2066  gather_offsets[i] = 0;
2067  }
2068 
2069  // Computation of scatter and offsets indices
2070  int f_ind=0;
2071  for (int f = 0; f < fes.GetNF(); ++f)
2072  {
2073  Mesh::FaceInformation face = mesh.GetFaceInformation(f);
2074  if ( face.IsNonconformingCoarse() )
2075  {
2076  // We skip nonconforming coarse faces as they are treated
2077  // by the corresponding nonconforming fine faces.
2078  continue;
2079  }
2080  else if ( type==FaceType::Interior && face.IsInterior() )
2081  {
2082  SetFaceDofsScatterIndices1(face,f_ind);
2083  if ( m==L2FaceValues::DoubleValued )
2084  {
2085  PermuteAndSetFaceDofsScatterIndices2(face,f_ind);
2086  if ( face.IsConforming() )
2087  {
2088  interpolations.RegisterFaceConformingInterpolation(face,f_ind);
2089  }
2090  else // Non-conforming face
2091  {
2092  interpolations.RegisterFaceCoarseToFineInterpolation(face,f_ind);
2093  }
2094  }
2095  f_ind++;
2096  }
2097  else if ( type==FaceType::Boundary && face.IsBoundary() )
2098  {
2099  SetFaceDofsScatterIndices1(face,f_ind);
2100  if ( m==L2FaceValues::DoubleValued )
2101  {
2102  SetBoundaryDofsScatterIndices2(face,f_ind);
2103  }
2104  f_ind++;
2105  }
2106  }
2107  MFEM_VERIFY(f_ind==nf, "Unexpected number of " <<
2108  (type==FaceType::Interior? "interior" : "boundary") <<
2109  " faces: " << f_ind << " vs " << nf );
2110 
2111  // Summation of the offsets
2112  for (int i = 1; i <= ndofs; ++i)
2113  {
2114  gather_offsets[i] += gather_offsets[i - 1];
2115  }
2116 
2117  // Transform the interpolation matrix map into a contiguous memory structure.
2118  interpolations.LinearizeInterpolatorMapIntoVector();
2119 }
2120 
2121 void NCL2FaceRestriction::ComputeGatherIndices(
2122  const ElementDofOrdering ordering,
2123  const FaceType type)
2124 {
2125  Mesh &mesh = *fes.GetMesh();
2126  // Computation of gather_indices
2127  int f_ind = 0;
2128  for (int f = 0; f < fes.GetNF(); ++f)
2129  {
2130  Mesh::FaceInformation face = mesh.GetFaceInformation(f);
2131  MFEM_ASSERT(!face.IsShared(),
2132  "Unexpected shared face in NCL2FaceRestriction.");
2133  if ( face.IsNonconformingCoarse() )
2134  {
2135  // We skip nonconforming coarse faces as they are treated
2136  // by the corresponding nonconforming fine faces.
2137  continue;
2138  }
2139  else if ( face.IsOfFaceType(type) )
2140  {
2141  SetFaceDofsGatherIndices1(face,f_ind);
2142  if ( m==L2FaceValues::DoubleValued &&
2143  type==FaceType::Interior &&
2144  face.IsInterior() )
2145  {
2146  PermuteAndSetFaceDofsGatherIndices2(face,f_ind);
2147  }
2148  f_ind++;
2149  }
2150  }
2151  MFEM_VERIFY(f_ind==nf, "Unexpected number of " <<
2152  (type==FaceType::Interior? "interior" : "boundary") <<
2153  " faces: " << f_ind << " vs " << nf );
2154 
2155  // Switch back offsets to their correct value
2156  for (int i = ndofs; i > 0; --i)
2157  {
2158  gather_offsets[i] = gather_offsets[i - 1];
2159  }
2160  gather_offsets[0] = 0;
2161 }
2162 
2163 } // namespace mfem
Abstract class for all finite elements.
Definition: fe_base.hpp:235
T * HostReadWrite()
Shortcut for mfem::ReadWrite(a.GetMemory(), a.Size(), false).
Definition: array.hpp:324
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:138
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:563
void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
int * GetJ()
Definition: table.hpp:114
FaceInformation GetFaceInformation(int f) const
Definition: mesh.cpp:1130
H1FaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering ordering, const FaceType type, bool build)
Construct an H1FaceRestriction.
void MultTransposeUnsigned(const Vector &x, Vector &y) const
Compute MultTranspose without applying signs based on DOF orientations.
void PermuteAndSetSharedFaceDofsScatterIndices2(const Mesh::FaceInformation &face, const int face_index)
Permute and set the scattering indices of elem2 for the shared face described by the face...
Operator that extracts Face degrees of freedom for H1 FiniteElementSpaces.
bool IsNonconformingCoarse() const
Return true if the face is a nonconforming coarse face.
Definition: mesh.hpp:1412
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:3175
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:521
Array< int > gather_indices
void AddMultTranspose(const Vector &x, Vector &y) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
const T * HostRead() const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), false).
Definition: array.hpp:308
void SetBoundaryDofsScatterIndices2(const Mesh::FaceInformation &face, const int face_index)
Set the scattering indices of elem2 for the boundary face described by the face.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:472
std::pair< const DenseMatrix *, int > Key
MFEM_HOST_DEVICE T AtomicAdd(T &add, const T val)
Definition: backends.hpp:84
T * GetData()
Returns the data.
Definition: array.hpp:112
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_base.hpp:327
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
virtual double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:450
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
Operator that extracts face degrees of freedom for L2 nonconforming spaces.
virtual void DoubleValuedNonconformingMult(const Vector &x, Vector &y) const
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with no...
void MultLeftInverse(const Vector &x, Vector &y) const
Geometry::Type GetFaceBaseGeometry(int i) const
Definition: mesh.hpp:1059
Abstract parallel finite element space.
Definition: pfespace.hpp:28
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:3135
void SetFaceDofsGatherIndices1(const Mesh::FaceInformation &face, const int face_index)
Set the gathering indices of elem1 for the interior face described by the face.
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
Definition: vector.hpp:117
Memory< double > & GetMemoryData()
Definition: sparsemat.hpp:256
const L2FaceValues m
Memory< int > & GetMemoryI()
Definition: sparsemat.hpp:224
Array< int > gather_offsets
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
void DoubleValuedConformingAddMultTranspose(const Vector &x, Vector &y) const
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector. Should only be used with con...
void PermuteAndSetFaceDofsGatherIndices2(const Mesh::FaceInformation &face, const int face_index)
Permute and set the gathering indices of elem2 for the interior face described by the face...
L2ElementRestriction(const FiniteElementSpace &)
Memory< int > & GetMemoryJ()
Definition: sparsemat.hpp:240
Array< InterpConfig > interp_config
DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition: dtensor.hpp:131
bool IsShared() const
Return true if the face is a shared interior face which is NOT a master nonconforming face...
Definition: mesh.hpp:1359
int GetBasisType() const
Definition: fe_base.hpp:1173
void SetPointMat(const DenseMatrix &pm)
Set the underlying point matrix describing the transformation.
Definition: eltrans.hpp:401
FaceType
Definition: mesh.hpp:45
void RegisterFaceConformingInterpolation(const Mesh::FaceInformation &face, int face_index)
Register the face with face and index face_index as a conforming face for the interpolation of the de...
const FiniteElementSpace & fes
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
void PermuteAndSetFaceDofsScatterIndices2(const Mesh::FaceInformation &face, const int face_index)
Permute and set the scattering indices of elem2, and increment the offsets for the face described by ...
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:66
void SingleValuedConformingMult(const Vector &x, Vector &y) const
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with co...
void FillI(SparseMatrix &mat) const
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:433
Operator that extracts Face degrees of freedom for L2 spaces.
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:446
void Mult(const Vector &x, Vector &y) const override
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector.
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 CheckFESpace(const ElementDofOrdering ordering)
Verify that H1FaceRestriction is build from an H1 FESpace.
const FiniteElementSpace & fes
void LinearizeInterpolatorMapIntoVector()
Transform the interpolation matrix map into a contiguous memory structure.
int * ReadWriteI(bool on_dev=true)
Definition: sparsemat.hpp:230
Geometry::Type GetFaceGeometry(int i) const
Definition: mesh.hpp:1043
void AddMultTranspose(const Vector &x, Vector &y) const override
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:304
bool IsBoundary() const
Return true if the face is a boundary face.
Definition: mesh.hpp:1375
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:228
ElementRestriction(const FiniteElementSpace &, ElementDofOrdering)
Definition: restriction.cpp:27
void GetFaceDofs(const int dim, const int face_id, const int dof1d, Array< int > &face_map)
Return the face map that extracts the degrees of freedom for the requested local face of a quad or he...
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
Definition: fespace.hpp:596
int Dimension() const
Definition: mesh.hpp:999
This structure is used as a human readable output format that decipheres the information contained in...
Definition: mesh.hpp:1333
int ToLexOrdering(const int dim, const int face_id, const int size1d, const int index)
Convert a dof face index from Native ordering to lexicographic ordering for quads and hexes...
const FiniteElementSpace & fes
Definition: restriction.hpp:35
void SetFaceDofsGatherIndices(const Mesh::FaceInformation &face, const int face_index, const ElementDofOrdering ordering)
Set the gathering indices of elem1 for the interior face described by the face.
Array< int > gather_offsets
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
struct mfem::Mesh::FaceInformation::@3 element[2]
void SetIdentityTransformation(Geometry::Type GeomType)
Set the FiniteElement Geometry for the reference elements being used.
Definition: eltrans.cpp:383
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
bool IsLocal() const
Return true if the face is a local interior face which is NOT a master nonconforming face...
Definition: mesh.hpp:1352
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
int GetNumInterpolators() const
Return the total number of interpolators.
NCL2FaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering ordering, const FaceType type, const L2FaceValues m, bool build)
Constructs an NCL2FaceRestriction, this is a specialization of a L2FaceRestriction for nonconforming ...
The ordering method used when the number of unknowns per mesh node (vector dimension) is bigger than ...
Definition: fespace.hpp:29
int * WriteJ(bool on_dev=true)
Definition: sparsemat.hpp:244
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_base.hpp:1181
void SingleValuedConformingAddMultTranspose(const Vector &x, Vector &y) const
Gather the degrees of freedom, i.e. goes from face E-Vector to L-Vector. Should only be used with con...
Array< int > scatter_indices1
void CheckFESpace(const ElementDofOrdering ordering)
Verify that L2FaceRestriction is build from an L2 FESpace.
L2FaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering ordering, const FaceType type, const L2FaceValues m, bool build)
Constructs an L2FaceRestriction.
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:454
A standard isoparametric element transformation.
Definition: eltrans.hpp:361
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
Mesh * GetMesh(int type)
Definition: ex29.cpp:218
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:66
int dim
Definition: ex24.cpp:53
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2783
int index(int i, int j, int nx, int ny)
Definition: life.cpp:237
const FiniteElementSpace & fes
Lexicographic ordering for tensor-product FiniteElements.
void SetFaceDofsScatterIndices(const Mesh::FaceInformation &face, const int face_index, const ElementDofOrdering ordering)
Set the scattering indices of elem1, and increment the offsets for the face described by the face...
Bernstein polynomials.
Definition: fe_base.hpp:32
RefCoord t[3]
ElementConformity conformity
Definition: mesh.hpp:1340
Vector data type.
Definition: vector.hpp:60
virtual void DoubleValuedConformingMult(const Vector &x, Vector &y) const
Scatter the degrees of freedom, i.e. goes from L-Vector to face E-Vector. Should only be used with co...
virtual void FillJAndData(const Vector &fea_data, SparseMatrix &mat, const bool keep_nbr_block=false) const
Fill the J and Data arrays of the SparseMatrix corresponding to the sparsity pattern given by this L2...
int PermuteFaceL2(const int dim, const int face_id1, const int face_id2, const int orientation, const int size1d, const int index)
Compute the dof face index of elem2 corresponding to the given dof face index.
virtual void AddFaceMatricesToElementMatrices(const Vector &fea_data, Vector &ea_data) const
This methods adds the DG face matrices to the element matrices.
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:726
int * HostReadWriteI()
Definition: sparsemat.hpp:236
void SetFaceDofsScatterIndices1(const Mesh::FaceInformation &face, const int face_index)
Set the scattering indices of elem1, and increment the offsets for the face described by the face...
void RegisterFaceCoarseToFineInterpolation(const Mesh::FaceInformation &face, int face_index)
Register the face with face and index face_index as a conforming face for the interpolation of the de...
const DenseMatrix & GetPointMat() const
Return the stored point matrix.
Definition: eltrans.hpp:404
uint32_t is_non_conforming
const Array< InterpConfig > & GetFaceInterpConfig() const
Return an array containing the interpolation configuration for each face registered with RegisterFace...
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:438
Array< int > scatter_indices
InterpolationManager interpolations
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
Array< int > scatter_indices2
virtual void FillI(SparseMatrix &mat, const bool keep_nbr_block=false) const
Fill the I array of SparseMatrix corresponding to the sparsity pattern given by this L2FaceRestrictio...
const DenseMatrix * point_matrix
Definition: mesh.hpp:1348
void FillSparseMatrix(const Vector &mat_ea, SparseMatrix &mat) const
Fill a Sparse Matrix with Element Matrices.
VDIM x NQPT x NE (values) / VDIM x DIM x NQPT x NE (grads)
DofTransformation * GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:1464
bool IsConforming() const
Return true if the face is a conforming face.
Definition: mesh.hpp:1395
double f(const Vector &p)
double * WriteData(bool on_dev=true)
Definition: sparsemat.hpp:260
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 Vector & GetInterpolators() const
Return an mfem::Vector containing the interpolators in the following format: face_dofs x face_dofs x ...