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