MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
prestriction.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "../config/config.hpp"
13 
14 #ifdef MFEM_USE_MPI
15 
16 #include "restriction.hpp"
17 #include "prestriction.hpp"
18 #include "pgridfunc.hpp"
19 #include "pfespace.hpp"
20 #include "fespace.hpp"
21 #include "../general/forall.hpp"
22 
23 namespace mfem
24 {
25 
27  ElementDofOrdering e_ordering,
28  FaceType type,
29  L2FaceValues m)
30  : fes(fes),
31  nf(fes.GetNFbyType(type)),
32  vdim(fes.GetVDim()),
33  byvdim(fes.GetOrdering() == Ordering::byVDIM),
34  ndofs(fes.GetNDofs()),
35  dof(nf>0 ?
36  fes.GetTraceElement(0, fes.GetMesh()->GetFaceBaseGeometry(0))->GetDof()
37  : 0),
38  m(m),
39  nfdofs(nf*dof),
40  scatter_indices1(nf*dof),
41  scatter_indices2(m==L2FaceValues::DoubleValued?nf*dof:0),
42  offsets(ndofs+1),
43  gather_indices((m==L2FaceValues::DoubleValued? 2 : 1)*nf*dof)
44 {
45  if (nf==0) { return; }
46  // If fespace == L2
47  const FiniteElement *fe = fes.GetFE(0);
48  const TensorBasisElement *tfe = dynamic_cast<const TensorBasisElement*>(fe);
49  MFEM_VERIFY(tfe != NULL &&
52  "Only Gauss-Lobatto and Bernstein basis are supported in "
53  "ParL2FaceRestriction.");
54  MFEM_VERIFY(fes.GetMesh()->Conforming(),
55  "Non-conforming meshes not yet supported with partial assembly.");
56  // Assuming all finite elements are using Gauss-Lobatto dofs
58  width = fes.GetVSize();
59  const bool dof_reorder = (e_ordering == ElementDofOrdering::LEXICOGRAPHIC);
60  if (!dof_reorder)
61  {
62  mfem_error("Non-Tensor L2FaceRestriction not yet implemented.");
63  }
64  if (dof_reorder && nf > 0)
65  {
66  for (int f = 0; f < fes.GetNF(); ++f)
67  {
68  const FiniteElement *fe =
70  const TensorBasisElement* el =
71  dynamic_cast<const TensorBasisElement*>(fe);
72  if (el) { continue; }
73  mfem_error("Finite element not suitable for lexicographic ordering");
74  }
75  }
76  const Table& e2dTable = fes.GetElementToDofTable();
77  const int* elementMap = e2dTable.GetJ();
78  Array<int> faceMap1(dof), faceMap2(dof);
79  int e1, e2;
80  int inf1, inf2;
81  int face_id1, face_id2;
82  int orientation;
83  const int dof1d = fes.GetFE(0)->GetOrder()+1;
84  const int elem_dofs = fes.GetFE(0)->GetDof();
85  const int dim = fes.GetMesh()->SpaceDimension();
86  // Computation of scatter indices
87  int f_ind=0;
88  for (int f = 0; f < fes.GetNF(); ++f)
89  {
90  fes.GetMesh()->GetFaceElements(f, &e1, &e2);
91  fes.GetMesh()->GetFaceInfos(f, &inf1, &inf2);
92  if (dof_reorder)
93  {
94  orientation = inf1 % 64;
95  face_id1 = inf1 / 64;
96  GetFaceDofs(dim, face_id1, dof1d, faceMap1); // only for hex
97  orientation = inf2 % 64;
98  face_id2 = inf2 / 64;
99  GetFaceDofs(dim, face_id2, dof1d, faceMap2); // only for hex
100  }
101  else
102  {
103  mfem_error("FaceRestriction not yet implemented for this type of "
104  "element.");
105  // TODO Something with GetFaceDofs?
106  orientation = 0; // suppress compiler warning
107  face_id1 = face_id2 = 0; // suppress compiler warning
108  }
109  if (type==FaceType::Interior &&
110  (e2>=0 || (e2<0 && inf2>=0) )) // interior/shared face
111  {
112  for (int d = 0; d < dof; ++d)
113  {
114  const int face_dof = faceMap1[d];
115  const int did = face_dof;
116  const int gid = elementMap[e1*elem_dofs + did];
117  const int lid = dof*f_ind + d;
118  scatter_indices1[lid] = gid;
119  }
121  {
122  if (e2>=0) // interior face
123  {
124  for (int d = 0; d < dof; ++d)
125  {
126  const int pd = PermuteFaceL2(dim, face_id1, face_id2,
127  orientation, dof1d, d);
128  const int face_dof = faceMap2[pd];
129  const int did = face_dof;
130  const int gid = elementMap[e2*elem_dofs + did];
131  const int lid = dof*f_ind + d;
132  scatter_indices2[lid] = gid;
133  }
134  }
135  else if (inf2>=0) // shared boundary
136  {
137  const int se2 = -1 - e2;
138  Array<int> sharedDofs;
139  fes.GetFaceNbrElementVDofs(se2, sharedDofs);
140  for (int d = 0; d < dof; ++d)
141  {
142  const int pd = PermuteFaceL2(dim, face_id1, face_id2,
143  orientation, dof1d, d);
144  const int face_dof = faceMap2[pd];
145  const int did = face_dof;
146  const int gid = sharedDofs[did];
147  const int lid = dof*f_ind + d;
148  // Trick to differentiate dof location inter/shared
149  scatter_indices2[lid] = ndofs+gid;
150  }
151  }
152  }
153  f_ind++;
154  }
155  else if (type==FaceType::Boundary && e2<0 && inf2<0) // true boundary
156  {
157  for (int d = 0; d < dof; ++d)
158  {
159  const int face_dof = faceMap1[d];
160  const int did = face_dof;
161  const int gid = elementMap[e1*elem_dofs + did];
162  const int lid = dof*f_ind + d;
163  scatter_indices1[lid] = gid;
164  }
166  {
167  for (int d = 0; d < dof; ++d)
168  {
169  const int lid = dof*f_ind + d;
170  scatter_indices2[lid] = -1;
171  }
172  }
173  f_ind++;
174  }
175  }
176  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
177  // Computation of gather_indices
178  for (int i = 0; i <= ndofs; ++i)
179  {
180  offsets[i] = 0;
181  }
182  f_ind = 0;
183  for (int f = 0; f < fes.GetNF(); ++f)
184  {
185  fes.GetMesh()->GetFaceElements(f, &e1, &e2);
186  fes.GetMesh()->GetFaceInfos(f, &inf1, &inf2);
187  if ((type==FaceType::Interior && (e2>=0 || (e2<0 && inf2>=0))) ||
188  (type==FaceType::Boundary && e2<0 && inf2<0) )
189  {
190  orientation = inf1 % 64;
191  face_id1 = inf1 / 64;
192  GetFaceDofs(dim, face_id1, dof1d, faceMap1);
193  orientation = inf2 % 64;
194  face_id2 = inf2 / 64;
195  GetFaceDofs(dim, face_id2, dof1d, faceMap2);
196  for (int d = 0; d < dof; ++d)
197  {
198  const int did = faceMap1[d];
199  const int gid = elementMap[e1*elem_dofs + did];
200  ++offsets[gid + 1];
201  }
203  {
204  for (int d = 0; d < dof; ++d)
205  {
206  if (type==FaceType::Interior && e2>=0) // interior face
207  {
208  const int pd = PermuteFaceL2(dim, face_id1, face_id2,
209  orientation, dof1d, d);
210  const int did = faceMap2[pd];
211  const int gid = elementMap[e2*elem_dofs + did];
212  ++offsets[gid + 1];
213  }
214  }
215  }
216  f_ind++;
217  }
218  }
219  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
220  for (int i = 1; i <= ndofs; ++i)
221  {
222  offsets[i] += offsets[i - 1];
223  }
224  f_ind = 0;
225  for (int f = 0; f < fes.GetNF(); ++f)
226  {
227  fes.GetMesh()->GetFaceElements(f, &e1, &e2);
228  fes.GetMesh()->GetFaceInfos(f, &inf1, &inf2);
229  if ((type==FaceType::Interior && (e2>=0 || (e2<0 && inf2>=0))) ||
230  (type==FaceType::Boundary && e2<0 && inf2<0) )
231  {
232  orientation = inf1 % 64;
233  face_id1 = inf1 / 64;
234  GetFaceDofs(dim, face_id1, dof1d, faceMap1);
235  orientation = inf2 % 64;
236  face_id2 = inf2 / 64;
237  GetFaceDofs(dim, face_id2, dof1d, faceMap2);
238  for (int d = 0; d < dof; ++d)
239  {
240  const int did = faceMap1[d];
241  const int gid = elementMap[e1*elem_dofs + did];
242  const int lid = dof*f_ind + d;
243  // We don't shift lid to express that it's e1 of f
244  gather_indices[offsets[gid]++] = lid;
245  }
247  {
248  for (int d = 0; d < dof; ++d)
249  {
250  if (type==FaceType::Interior && e2>=0) // interior face
251  {
252  const int pd = PermuteFaceL2(dim, face_id1, face_id2,
253  orientation, dof1d, d);
254  const int did = faceMap2[pd];
255  const int gid = elementMap[e2*elem_dofs + did];
256  const int lid = dof*f_ind + d;
257  // We shift lid to express that it's e2 of f
258  gather_indices[offsets[gid]++] = nfdofs + lid;
259  }
260  }
261  }
262  f_ind++;
263  }
264  }
265  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
266  for (int i = ndofs; i > 0; --i)
267  {
268  offsets[i] = offsets[i - 1];
269  }
270  offsets[0] = 0;
271 }
272 
273 void ParL2FaceRestriction::Mult(const Vector& x, Vector& y) const
274 {
275  ParGridFunction x_gf;
276  x_gf.MakeRef(const_cast<ParFiniteElementSpace*>(&fes),
277  const_cast<Vector&>(x), 0);
278  x_gf.ExchangeFaceNbrData();
279 
280  // Assumes all elements have the same number of dofs
281  const int nd = dof;
282  const int vd = vdim;
283  const bool t = byvdim;
284  const int threshold = ndofs;
285 
287  {
288  auto d_indices1 = scatter_indices1.Read();
289  auto d_indices2 = scatter_indices2.Read();
290  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
291  auto d_x_shared = Reshape(x_gf.FaceNbrData().Read(),
292  t?vd:ndofs, t?ndofs:vd);
293  auto d_y = Reshape(y.Write(), nd, vd, 2, nf);
294  MFEM_FORALL(i, nfdofs,
295  {
296  const int dof = i % nd;
297  const int face = i / nd;
298  const int idx1 = d_indices1[i];
299  for (int c = 0; c < vd; ++c)
300  {
301  d_y(dof, c, 0, face) = d_x(t?c:idx1, t?idx1:c);
302  }
303  const int idx2 = d_indices2[i];
304  for (int c = 0; c < vd; ++c)
305  {
306  if (idx2>-1 && idx2<threshold) // interior face
307  {
308  d_y(dof, c, 1, face) = d_x(t?c:idx2, t?idx2:c);
309  }
310  else if (idx2>=threshold) // shared boundary
311  {
312  d_y(dof, c, 1, face) = d_x_shared(t?c:(idx2-threshold),
313  t?(idx2-threshold):c);
314  }
315  else // true boundary
316  {
317  d_y(dof, c, 1, face) = 0.0;
318  }
319  }
320  });
321  }
322  else
323  {
324  auto d_indices1 = scatter_indices1.Read();
325  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
326  auto d_y = Reshape(y.Write(), nd, vd, nf);
327  MFEM_FORALL(i, nfdofs,
328  {
329  const int dof = i % nd;
330  const int face = i / nd;
331  const int idx1 = d_indices1[i];
332  for (int c = 0; c < vd; ++c)
333  {
334  d_y(dof, c, face) = d_x(t?c:idx1, t?idx1:c);
335  }
336  });
337  }
338 }
339 
341 {
342  // Assumes all elements have the same number of dofs
343  const int nd = dof;
344  const int vd = vdim;
345  const bool t = byvdim;
346  const int dofs = nfdofs;
347  auto d_offsets = offsets.Read();
348  auto d_indices = gather_indices.Read();
349  auto d_x = Reshape(x.Read(), nd, vd, 2, nf);
350  auto d_y = Reshape(y.Write(), t?vd:ndofs, t?ndofs:vd);
351  MFEM_FORALL(i, ndofs,
352  {
353  const int offset = d_offsets[i];
354  const int nextOffset = d_offsets[i + 1];
355  for (int c = 0; c < vd; ++c)
356  {
357  double dofValue = 0;
358  for (int j = offset; j < nextOffset; ++j)
359  {
360  int idx_j = d_indices[j];
361  bool isE1 = idx_j < dofs;
362  idx_j = isE1 ? idx_j : idx_j - dofs;
363  dofValue += isE1 ?
364  d_x(idx_j % nd, c, 0, idx_j / nd)
365  :d_x(idx_j % nd, c, 1, idx_j / nd);
366  }
367  d_y(t?c:i,t?i:c) += dofValue;
368  }
369  });
370 }
371 
372 } // namespace mfem
373 
374 #endif
Abstract class for Finite Elements.
Definition: fe.hpp:232
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:381
int * GetJ()
Definition: table.hpp:114
L2FaceValues
Definition: restriction.hpp:26
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:1936
bool Conforming() const
Definition: mesh.hpp:1167
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe.hpp:318
Geometry::Type GetFaceBaseGeometry(int i) const
Definition: mesh.hpp:785
Abstract parallel finite element space.
Definition: pfespace.hpp:28
double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:341
DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition: dtensor.hpp:134
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
Definition: pfespace.cpp:1160
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Definition: pgridfunc.cpp:102
int GetBasisType() const
Definition: fe.hpp:1842
FaceType
Definition: mesh.hpp:42
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:295
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:153
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:273
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
Definition: fespace.hpp:414
Vector & FaceNbrData()
Definition: pgridfunc.hpp:198
int SpaceDimension() const
Definition: mesh.hpp:745
void GetFaceElements(int Face, int *Elem1, int *Elem2) const
Definition: mesh.cpp:974
void GetFaceInfos(int Face, int *Inf1, int *Inf2) const
Definition: mesh.cpp:980
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:314
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
const ParFiniteElementSpace & fes
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 ...
void GetFaceDofs(const int dim, const int face_id, const int dof1d, Array< int > &faceMap)
Return the face degrees of freedom returned in Lexicographic order.
The ordering method used when the number of unknowns per mesh node (vector dimension) is bigger than ...
Definition: fespace.hpp:28
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:333
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:65
int dim
Definition: ex24.cpp:43
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i&#39;th element.
Definition: fespace.cpp:1642
Lexicographic ordering for tensor-product FiniteElements.
ParL2FaceRestriction(const ParFiniteElementSpace &, ElementDofOrdering, FaceType type, L2FaceValues m=L2FaceValues::DoubleValued)
Vector data type.
Definition: vector.hpp:48
int PermuteFaceL2(const int dim, const int face_id1, const int face_id2, const int orientation, const int size1d, const int index)
Permute dofs or quads on a face for e2 to match with the ordering of e1.
const L2FaceValues m
const Table & GetElementToDofTable() const
Definition: fespace.hpp:524
Bernstein polynomials.
Definition: fe.hpp:35
Class for parallel grid function.
Definition: pgridfunc.hpp:32
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
VDIM x NQPT x NE.