MFEM  v4.2.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  : L2FaceRestriction(fes, type, m)
31 {
32  if (nf==0) { return; }
33  // If fespace == L2
34  const ParFiniteElementSpace &pfes =
35  static_cast<const ParFiniteElementSpace&>(this->fes);
36  const FiniteElement *fe = pfes.GetFE(0);
37  const TensorBasisElement *tfe = dynamic_cast<const TensorBasisElement*>(fe);
38  MFEM_VERIFY(tfe != NULL &&
41  "Only Gauss-Lobatto and Bernstein basis are supported in "
42  "ParL2FaceRestriction.");
43  MFEM_VERIFY(pfes.GetMesh()->Conforming(),
44  "Non-conforming meshes not yet supported with partial assembly.");
45  // Assuming all finite elements are using Gauss-Lobatto dofs
47  width = pfes.GetVSize();
48  const bool dof_reorder = (e_ordering == ElementDofOrdering::LEXICOGRAPHIC);
49  if (!dof_reorder)
50  {
51  mfem_error("Non-Tensor L2FaceRestriction not yet implemented.");
52  }
53  if (dof_reorder && nf > 0)
54  {
55  for (int f = 0; f < pfes.GetNF(); ++f)
56  {
57  const FiniteElement *fe =
59  const TensorBasisElement* el =
60  dynamic_cast<const TensorBasisElement*>(fe);
61  if (el) { continue; }
62  mfem_error("Finite element not suitable for lexicographic ordering");
63  }
64  }
65  const Table& e2dTable = pfes.GetElementToDofTable();
66  const int* elementMap = e2dTable.GetJ();
67  Array<int> faceMap1(dof), faceMap2(dof);
68  int e1, e2;
69  int inf1, inf2;
70  int face_id1, face_id2;
71  int orientation;
72  const int dof1d = pfes.GetFE(0)->GetOrder()+1;
73  const int elem_dofs = pfes.GetFE(0)->GetDof();
74  const int dim = pfes.GetMesh()->SpaceDimension();
75  // Computation of scatter indices
76  int f_ind=0;
77  for (int f = 0; f < pfes.GetNF(); ++f)
78  {
79  pfes.GetMesh()->GetFaceElements(f, &e1, &e2);
80  pfes.GetMesh()->GetFaceInfos(f, &inf1, &inf2);
81  if (dof_reorder)
82  {
83  orientation = inf1 % 64;
84  face_id1 = inf1 / 64;
85  GetFaceDofs(dim, face_id1, dof1d, faceMap1); // only for hex
86  orientation = inf2 % 64;
87  face_id2 = inf2 / 64;
88  GetFaceDofs(dim, face_id2, dof1d, faceMap2); // only for hex
89  }
90  else
91  {
92  mfem_error("FaceRestriction not yet implemented for this type of "
93  "element.");
94  // TODO Something with GetFaceDofs?
95  orientation = 0; // suppress compiler warning
96  face_id1 = face_id2 = 0; // suppress compiler warning
97  }
98  if (type==FaceType::Interior &&
99  (e2>=0 || (e2<0 && inf2>=0) )) // interior/shared face
100  {
101  for (int d = 0; d < dof; ++d)
102  {
103  const int face_dof = faceMap1[d];
104  const int did = face_dof;
105  const int gid = elementMap[e1*elem_dofs + did];
106  const int lid = dof*f_ind + d;
107  scatter_indices1[lid] = gid;
108  }
110  {
111  if (e2>=0) // interior face
112  {
113  for (int d = 0; d < dof; ++d)
114  {
115  const int pd = PermuteFaceL2(dim, face_id1, face_id2,
116  orientation, dof1d, d);
117  const int face_dof = faceMap2[pd];
118  const int did = face_dof;
119  const int gid = elementMap[e2*elem_dofs + did];
120  const int lid = dof*f_ind + d;
121  scatter_indices2[lid] = gid;
122  }
123  }
124  else if (inf2>=0) // shared boundary
125  {
126  const int se2 = -1 - e2;
127  Array<int> sharedDofs;
128  pfes.GetFaceNbrElementVDofs(se2, sharedDofs);
129  for (int d = 0; d < dof; ++d)
130  {
131  const int pd = PermuteFaceL2(dim, face_id1, face_id2,
132  orientation, dof1d, d);
133  const int face_dof = faceMap2[pd];
134  const int did = face_dof;
135  const int gid = sharedDofs[did];
136  const int lid = dof*f_ind + d;
137  // Trick to differentiate dof location inter/shared
138  scatter_indices2[lid] = ndofs+gid;
139  }
140  }
141  }
142  f_ind++;
143  }
144  else if (type==FaceType::Boundary && e2<0 && inf2<0) // true boundary
145  {
146  for (int d = 0; d < dof; ++d)
147  {
148  const int face_dof = faceMap1[d];
149  const int did = face_dof;
150  const int gid = elementMap[e1*elem_dofs + did];
151  const int lid = dof*f_ind + d;
152  scatter_indices1[lid] = gid;
153  }
155  {
156  for (int d = 0; d < dof; ++d)
157  {
158  const int lid = dof*f_ind + d;
159  scatter_indices2[lid] = -1;
160  }
161  }
162  f_ind++;
163  }
164  }
165  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
166  // Computation of gather_indices
167  for (int i = 0; i <= ndofs; ++i)
168  {
169  offsets[i] = 0;
170  }
171  f_ind = 0;
172  for (int f = 0; f < pfes.GetNF(); ++f)
173  {
174  pfes.GetMesh()->GetFaceElements(f, &e1, &e2);
175  pfes.GetMesh()->GetFaceInfos(f, &inf1, &inf2);
176  if ((type==FaceType::Interior && (e2>=0 || (e2<0 && inf2>=0))) ||
177  (type==FaceType::Boundary && e2<0 && inf2<0) )
178  {
179  orientation = inf1 % 64;
180  face_id1 = inf1 / 64;
181  GetFaceDofs(dim, face_id1, dof1d, faceMap1);
182  orientation = inf2 % 64;
183  face_id2 = inf2 / 64;
184  GetFaceDofs(dim, face_id2, dof1d, faceMap2);
185  for (int d = 0; d < dof; ++d)
186  {
187  const int did = faceMap1[d];
188  const int gid = elementMap[e1*elem_dofs + did];
189  ++offsets[gid + 1];
190  }
192  {
193  for (int d = 0; d < dof; ++d)
194  {
195  if (type==FaceType::Interior && e2>=0) // interior face
196  {
197  const int pd = PermuteFaceL2(dim, face_id1, face_id2,
198  orientation, dof1d, d);
199  const int did = faceMap2[pd];
200  const int gid = elementMap[e2*elem_dofs + did];
201  ++offsets[gid + 1];
202  }
203  }
204  }
205  f_ind++;
206  }
207  }
208  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
209  for (int i = 1; i <= ndofs; ++i)
210  {
211  offsets[i] += offsets[i - 1];
212  }
213  f_ind = 0;
214  for (int f = 0; f < pfes.GetNF(); ++f)
215  {
216  pfes.GetMesh()->GetFaceElements(f, &e1, &e2);
217  pfes.GetMesh()->GetFaceInfos(f, &inf1, &inf2);
218  if ((type==FaceType::Interior && (e2>=0 || (e2<0 && inf2>=0))) ||
219  (type==FaceType::Boundary && e2<0 && inf2<0) )
220  {
221  orientation = inf1 % 64;
222  face_id1 = inf1 / 64;
223  GetFaceDofs(dim, face_id1, dof1d, faceMap1);
224  orientation = inf2 % 64;
225  face_id2 = inf2 / 64;
226  GetFaceDofs(dim, face_id2, dof1d, faceMap2);
227  for (int d = 0; d < dof; ++d)
228  {
229  const int did = faceMap1[d];
230  const int gid = elementMap[e1*elem_dofs + did];
231  const int lid = dof*f_ind + d;
232  // We don't shift lid to express that it's e1 of f
233  gather_indices[offsets[gid]++] = lid;
234  }
236  {
237  for (int d = 0; d < dof; ++d)
238  {
239  if (type==FaceType::Interior && e2>=0) // interior face
240  {
241  const int pd = PermuteFaceL2(dim, face_id1, face_id2,
242  orientation, dof1d, d);
243  const int did = faceMap2[pd];
244  const int gid = elementMap[e2*elem_dofs + did];
245  const int lid = dof*f_ind + d;
246  // We shift lid to express that it's e2 of f
247  gather_indices[offsets[gid]++] = nfdofs + lid;
248  }
249  }
250  }
251  f_ind++;
252  }
253  }
254  MFEM_VERIFY(f_ind==nf, "Unexpected number of faces.");
255  for (int i = ndofs; i > 0; --i)
256  {
257  offsets[i] = offsets[i - 1];
258  }
259  offsets[0] = 0;
260 }
261 
262 void ParL2FaceRestriction::Mult(const Vector& x, Vector& y) const
263 {
264  const ParFiniteElementSpace &pfes =
265  static_cast<const ParFiniteElementSpace&>(this->fes);
266  ParGridFunction x_gf;
267  x_gf.MakeRef(const_cast<ParFiniteElementSpace*>(&pfes),
268  const_cast<Vector&>(x), 0);
269  x_gf.ExchangeFaceNbrData();
270 
271  // Assumes all elements have the same number of dofs
272  const int nd = dof;
273  const int vd = vdim;
274  const bool t = byvdim;
275  const int threshold = ndofs;
276 
278  {
279  auto d_indices1 = scatter_indices1.Read();
280  auto d_indices2 = scatter_indices2.Read();
281  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
282  auto d_x_shared = Reshape(x_gf.FaceNbrData().Read(),
283  t?vd:ndofs, t?ndofs:vd);
284  auto d_y = Reshape(y.Write(), nd, vd, 2, nf);
285  MFEM_FORALL(i, nfdofs,
286  {
287  const int dof = i % nd;
288  const int face = i / nd;
289  const int idx1 = d_indices1[i];
290  for (int c = 0; c < vd; ++c)
291  {
292  d_y(dof, c, 0, face) = d_x(t?c:idx1, t?idx1:c);
293  }
294  const int idx2 = d_indices2[i];
295  for (int c = 0; c < vd; ++c)
296  {
297  if (idx2>-1 && idx2<threshold) // interior face
298  {
299  d_y(dof, c, 1, face) = d_x(t?c:idx2, t?idx2:c);
300  }
301  else if (idx2>=threshold) // shared boundary
302  {
303  d_y(dof, c, 1, face) = d_x_shared(t?c:(idx2-threshold),
304  t?(idx2-threshold):c);
305  }
306  else // true boundary
307  {
308  d_y(dof, c, 1, face) = 0.0;
309  }
310  }
311  });
312  }
313  else
314  {
315  auto d_indices1 = scatter_indices1.Read();
316  auto d_x = Reshape(x.Read(), t?vd:ndofs, t?ndofs:vd);
317  auto d_y = Reshape(y.Write(), nd, vd, nf);
318  MFEM_FORALL(i, nfdofs,
319  {
320  const int dof = i % nd;
321  const int face = i / nd;
322  const int idx1 = d_indices1[i];
323  for (int c = 0; c < vd; ++c)
324  {
325  d_y(dof, c, face) = d_x(t?c:idx1, t?idx1:c);
326  }
327  });
328  }
329 }
330 
331 static MFEM_HOST_DEVICE int AddNnz(const int iE, int *I, const int dofs)
332 {
333  int val = AtomicAdd(I[iE],dofs);
334  return val;
335 }
336 
338  SparseMatrix &face_mat) const
339 {
340  const int face_dofs = dof;
341  const int Ndofs = ndofs;
342  auto d_indices1 = scatter_indices1.Read();
343  auto d_indices2 = scatter_indices2.Read();
344  auto I = mat.ReadWriteI();
345  auto I_face = face_mat.ReadWriteI();
346  MFEM_FORALL(i, ne*elemDofs*vdim+1,
347  {
348  I_face[i] = 0;
349  });
350  MFEM_FORALL(fdof, nf*face_dofs,
351  {
352  const int f = fdof/face_dofs;
353  const int iF = fdof%face_dofs;
354  const int iE1 = d_indices1[f*face_dofs+iF];
355  if (iE1 < Ndofs)
356  {
357  for (int jF = 0; jF < face_dofs; jF++)
358  {
359  const int jE2 = d_indices2[f*face_dofs+jF];
360  if (jE2 < Ndofs)
361  {
362  AddNnz(iE1,I,1);
363  }
364  else
365  {
366  AddNnz(iE1,I_face,1);
367  }
368  }
369  }
370  const int iE2 = d_indices2[f*face_dofs+iF];
371  if (iE2 < Ndofs)
372  {
373  for (int jF = 0; jF < face_dofs; jF++)
374  {
375  const int jE1 = d_indices1[f*face_dofs+jF];
376  if (jE1 < Ndofs)
377  {
378  AddNnz(iE2,I,1);
379  }
380  else
381  {
382  AddNnz(iE2,I_face,1);
383  }
384  }
385  }
386  });
387 }
388 
390  SparseMatrix &mat,
391  SparseMatrix &face_mat) const
392 {
393  const int face_dofs = dof;
394  const int Ndofs = ndofs;
395  auto d_indices1 = scatter_indices1.Read();
396  auto d_indices2 = scatter_indices2.Read();
397  auto mat_fea = Reshape(ea_data.Read(), face_dofs, face_dofs, 2, nf);
398  auto I = mat.ReadWriteI();
399  auto I_face = face_mat.ReadWriteI();
400  auto J = mat.WriteJ();
401  auto J_face = face_mat.WriteJ();
402  auto Data = mat.WriteData();
403  auto Data_face = face_mat.WriteData();
404  MFEM_FORALL(fdof, nf*face_dofs,
405  {
406  const int f = fdof/face_dofs;
407  const int iF = fdof%face_dofs;
408  const int iE1 = d_indices1[f*face_dofs+iF];
409  if (iE1 < Ndofs)
410  {
411  for (int jF = 0; jF < face_dofs; jF++)
412  {
413  const int jE2 = d_indices2[f*face_dofs+jF];
414  if (jE2 < Ndofs)
415  {
416  const int offset = AddNnz(iE1,I,1);
417  J[offset] = jE2;
418  Data[offset] = mat_fea(jF,iF,1,f);
419  }
420  else
421  {
422  const int offset = AddNnz(iE1,I_face,1);
423  J_face[offset] = jE2-Ndofs;
424  Data_face[offset] = mat_fea(jF,iF,1,f);
425  }
426  }
427  }
428  const int iE2 = d_indices2[f*face_dofs+iF];
429  if (iE2 < Ndofs)
430  {
431  for (int jF = 0; jF < face_dofs; jF++)
432  {
433  const int jE1 = d_indices1[f*face_dofs+jF];
434  if (jE1 < Ndofs)
435  {
436  const int offset = AddNnz(iE2,I,1);
437  J[offset] = jE1;
438  Data[offset] = mat_fea(jF,iF,0,f);
439  }
440  else
441  {
442  const int offset = AddNnz(iE2,I_face,1);
443  J_face[offset] = jE1-Ndofs;
444  Data_face[offset] = mat_fea(jF,iF,0,f);
445  }
446  }
447  }
448  });
449 }
450 
451 } // namespace mfem
452 
453 #endif
Abstract class for all finite elements.
Definition: fe.hpp:235
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:400
int * GetJ()
Definition: table.hpp:114
L2FaceValues
Definition: restriction.hpp:26
const FiniteElement * GetTraceElement(int i, Geometry::Type geom_type) const
Return the trace element from element &#39;i&#39; to the given &#39;geom_type&#39;.
Definition: fespace.cpp:2108
bool Conforming() const
Definition: mesh.hpp:1213
MFEM_HOST_DEVICE T AtomicAdd(T &add, const T val)
Definition: backends.hpp:59
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe.hpp:319
Geometry::Type GetFaceBaseGeometry(int i) const
Definition: mesh.hpp:829
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:380
const L2FaceValues m
Array< int > gather_indices
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:1165
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:2021
FaceType
Definition: mesh.hpp:45
const FiniteElementSpace & fes
Bernstein polynomials.
Definition: fe.hpp:35
Data type sparse matrix.
Definition: sparsemat.hpp:46
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:314
Operator that extracts Face degrees of freedom.
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:153
int * ReadWriteI(bool on_dev=true)
Definition: sparsemat.hpp:194
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
Definition: array.hpp:290
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
Definition: fespace.hpp:433
Vector & FaceNbrData()
Definition: pgridfunc.hpp:203
int SpaceDimension() const
Definition: mesh.hpp:789
void GetFaceElements(int Face, int *Elem1, int *Elem2) const
Definition: mesh.cpp:1031
void GetFaceInfos(int Face, int *Inf1, int *Inf2) const
Definition: mesh.cpp:1037
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:315
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void GetFaceDofs(const int dim, const int face_id, const int dof1d, Array< int > &faceMap)
int * WriteJ(bool on_dev=true)
Definition: sparsemat.hpp:208
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
Array< int > scatter_indices1
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:372
void FillI(SparseMatrix &mat, SparseMatrix &face_mat) const
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
Definition: fespace.hpp:65
int dim
Definition: ex24.cpp:53
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:1798
Lexicographic ordering for tensor-product FiniteElements.
ParL2FaceRestriction(const ParFiniteElementSpace &, ElementDofOrdering, FaceType type, L2FaceValues m=L2FaceValues::DoubleValued)
Vector data type.
Definition: vector.hpp:51
void FillJAndData(const Vector &ea_data, SparseMatrix &mat, SparseMatrix &face_mat) const
int PermuteFaceL2(const int dim, const int face_id1, const int face_id2, const int orientation, const int size1d, const int index)
const Table & GetElementToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each mesh element...
Definition: fespace.hpp:543
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
Array< int > scatter_indices2
double f(const Vector &p)
double * WriteData(bool on_dev=true)
Definition: sparsemat.hpp:224