12 #include "../config/config.hpp"
21 #include "../general/forall.hpp"
31 nf(fes.GetNFbyType(type)),
34 ndofs(fes.GetNDofs()),
36 fes.GetTraceElement(0, fes.GetMesh()->GetFaceBaseGeometry(0))->GetDof()
40 scatter_indices1(nf*dof),
45 if (
nf==0) {
return; }
49 MFEM_VERIFY(tfe != NULL &&
52 "Only Gauss-Lobatto and Bernstein basis are supported in "
53 "ParL2FaceRestriction.");
55 "Non-conforming meshes not yet supported with partial assembly.");
62 mfem_error(
"Non-Tensor L2FaceRestriction not yet implemented.");
64 if (dof_reorder &&
nf > 0)
66 for (
int f = 0; f < fes.
GetNF(); ++f)
73 mfem_error(
"Finite element not suitable for lexicographic ordering");
77 const int* elementMap = e2dTable.
GetJ();
81 int face_id1, face_id2;
88 for (
int f = 0; f < fes.
GetNF(); ++f)
94 orientation = inf1 % 64;
97 orientation = inf2 % 64;
103 mfem_error(
"FaceRestriction not yet implemented for this type of "
107 face_id1 = face_id2 = 0;
110 (e2>=0 || (e2<0 && inf2>=0) ))
112 for (
int d = 0; d <
dof; ++d)
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;
124 for (
int d = 0; d <
dof; ++d)
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;
137 const int se2 = -1 - e2;
140 for (
int d = 0; d <
dof; ++d)
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;
157 for (
int d = 0; d <
dof; ++d)
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;
167 for (
int d = 0; d <
dof; ++d)
169 const int lid = dof*f_ind + d;
176 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
178 for (
int i = 0; i <=
ndofs; ++i)
183 for (
int f = 0; f < fes.
GetNF(); ++f)
190 orientation = inf1 % 64;
191 face_id1 = inf1 / 64;
193 orientation = inf2 % 64;
194 face_id2 = inf2 / 64;
196 for (
int d = 0; d <
dof; ++d)
198 const int did = faceMap1[d];
199 const int gid = elementMap[e1*elem_dofs + did];
204 for (
int d = 0; d <
dof; ++d)
209 orientation, dof1d, d);
210 const int did = faceMap2[pd];
211 const int gid = elementMap[e2*elem_dofs + did];
219 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
220 for (
int i = 1; i <=
ndofs; ++i)
225 for (
int f = 0; f < fes.
GetNF(); ++f)
232 orientation = inf1 % 64;
233 face_id1 = inf1 / 64;
235 orientation = inf2 % 64;
236 face_id2 = inf2 / 64;
238 for (
int d = 0; d <
dof; ++d)
240 const int did = faceMap1[d];
241 const int gid = elementMap[e1*elem_dofs + did];
242 const int lid = dof*f_ind + d;
248 for (
int d = 0; d <
dof; ++d)
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;
265 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
266 for (
int i =
ndofs; i > 0; --i)
276 x_gf.
MakeRef(const_cast<ParFiniteElementSpace*>(&
fes),
277 const_cast<Vector&>(x), 0);
284 const int threshold =
ndofs;
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)
301 d_y(dof, c, 0, face) = d_x(t?c:idx1, t?idx1:c);
303 const int idx2 = d_indices2[i];
304 for (
int c = 0; c < vd; ++c)
306 if (idx2>-1 && idx2<threshold)
308 d_y(dof, c, 1, face) = d_x(t?c:idx2, t?idx2:c);
310 else if (idx2>=threshold)
312 d_y(dof, c, 1, face) = d_x_shared(t?c:(idx2-threshold),
313 t?(idx2-threshold):c);
317 d_y(dof, c, 1, face) = 0.0;
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)
334 d_y(dof, c, face) = d_x(t?c:idx1, t?idx1:c);
351 MFEM_FORALL(i,
ndofs,
353 const int offset = d_offsets[i];
354 const int nextOffset = d_offsets[i + 1];
355 for (
int c = 0; c < vd; ++c)
358 for (
int j = offset; j < nextOffset; ++j)
360 int idx_j = d_indices[j];
361 bool isE1 = idx_j < dofs;
362 idx_j = isE1 ? idx_j : idx_j - dofs;
364 d_x(idx_j % nd, c, 0, idx_j / nd)
365 :d_x(idx_j % nd, c, 1, idx_j / nd);
367 d_y(t?c:i,t?i:c) += dofValue;
Abstract class for Finite Elements.
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
void ExchangeFaceNbrData()
const FiniteElement * GetTraceElement(int i, Geometry::Type geom_type) const
Return the trace element from element 'i' to the given 'geom_type'.
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Geometry::Type GetFaceBaseGeometry(int i) const
Abstract parallel finite element space.
double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Array< int > scatter_indices1
Array< int > gather_indices
DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
void GetFaceNbrElementVDofs(int i, Array< int > &vdofs) const
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Mesh * GetMesh() const
Returns the mesh.
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
int SpaceDimension() const
void GetFaceElements(int Face, int *Elem1, int *Elem2) const
void GetFaceInfos(int Face, int *Inf1, int *Inf2) const
int GetDof() const
Returns the number of degrees of freedom in the finite element.
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 ...
int height
Dimension of the output / number of rows in the matrix.
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement associated with i'th element.
Lexicographic ordering for tensor-product FiniteElements.
Array< int > scatter_indices2
ParL2FaceRestriction(const ParFiniteElementSpace &, ElementDofOrdering, FaceType type, L2FaceValues m=L2FaceValues::DoubleValued)
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 Table & GetElementToDofTable() const
Class for parallel grid function.
int width
Dimension of the input / number of columns in the matrix.