15 #include "../general/forall.hpp"
26 ndofs(fes.GetNDofs()),
27 dof(ne > 0 ? fes.GetFE(0)->GetDof() : 0),
37 const int *dof_map = NULL;
38 if (dof_reorder &&
ne > 0)
40 for (
int e = 0; e <
ne; ++e)
46 mfem_error(
"Finite element not suitable for lexicographic ordering");
52 MFEM_VERIFY(fe_dof_map.
Size() > 0,
"invalid dof map");
56 const int* elementMap = e2dTable.
GetJ();
58 for (
int i = 0; i <=
ndofs; ++i)
62 for (
int e = 0; e <
ne; ++e)
64 for (
int d = 0; d <
dof; ++d)
66 const int sgid = elementMap[dof*e + d];
67 const int gid = (sgid >= 0) ? sgid : -1 - sgid;
72 for (
int i = 1; i <=
ndofs; ++i)
77 for (
int e = 0; e <
ne; ++e)
79 for (
int d = 0; d <
dof; ++d)
81 const int sdid = dof_reorder ? dof_map[d] : 0;
82 const int did = (!dof_reorder)?d:(sdid >= 0 ? sdid : -1-sdid);
83 const int sgid = elementMap[dof*e + did];
84 const int gid = (sgid >= 0) ? sgid : -1-sgid;
85 const int lid = dof*e + d;
86 const bool plus = (sgid >= 0 && sdid >= 0) || (sgid < 0 && sdid < 0);
93 for (
int i =
ndofs; i > 0; --i)
109 MFEM_FORALL(i,
dof*
ne,
111 const int gid = d_gatherMap[i];
112 const bool plus = gid >= 0;
113 const int j = plus ? gid : -1-gid;
114 for (
int c = 0; c < vd; ++c)
116 const double dofValue = d_x(t?c:j, t?j:c);
117 d_y(i % nd, c, i / nd) = plus ? dofValue : -dofValue;
132 MFEM_FORALL(i,
dof*
ne,
134 const int gid = d_gatherMap[i];
135 const int j = gid >= 0 ? gid : -1-gid;
136 for (
int c = 0; c < vd; ++c)
138 d_y(i % nd, c, i / nd) = d_x(t?c:j, t?j:c);
153 MFEM_FORALL(i,
ndofs,
155 const int offset = d_offsets[i];
156 const int nextOffset = d_offsets[i + 1];
157 for (
int c = 0; c < vd; ++c)
160 for (
int j = offset; j < nextOffset; ++j)
162 const int idx_j = (d_indices[j] >= 0) ? d_indices[j] : -1 - d_indices[j];
163 dofValue += (d_indices[j] >= 0) ? d_x(idx_j % nd, c,
164 idx_j / nd) : -d_x(idx_j % nd, c, idx_j / nd);
166 d_y(t?c:i,t?i:c) = dofValue;
181 MFEM_FORALL(i,
ndofs,
183 const int offset = d_offsets[i];
184 const int nextOffset = d_offsets[i + 1];
185 for (
int c = 0; c < vd; ++c)
188 for (
int j = offset; j < nextOffset; ++j)
190 const int idx_j = (d_indices[j] >= 0) ? d_indices[j] : -1 - d_indices[j];
191 dofValue += d_x(idx_j % nd, c, idx_j / nd);
193 d_y(t?c:i,t?i:c) = dofValue;
212 for (
int i = 0; i <
ndofs; ++i)
214 const int offset = d_offsets[i];
215 const int nextOffset = d_offsets[i+1];
216 for (
int c = 0; c < vd; ++c)
218 for (
int j = offset; j < nextOffset; ++j)
220 const int idx_j = d_indices[j];
221 if (d_x(t?c:i,t?i:c))
223 d_y(idx_j % nd, c, idx_j / nd) = 0.0;
227 d_y(idx_j % nd, c, idx_j / nd) = 1.0;
228 d_x(t?c:i,t?i:c) = 1;
239 const int nnz =
FillI(mat);
245 template <
int MaxNbNbr>
246 static MFEM_HOST_DEVICE
int GetMinElt(
const int *my_elts,
const int nbElts,
247 const int *nbr_elts,
const int nbrNbElts)
252 for (
int i = 0; i < nbElts; i++)
254 const int e_i = my_elts[i];
255 for (
int j = 0; j < nbrNbElts; j++)
257 if (e_i==nbr_elts[j])
266 for (
int i = 1; i < cpt; i++)
278 static MFEM_HOST_DEVICE
int GetAndIncrementNnzIndex(
const int i_L,
int* I)
286 static constexpr
int Max = MaxNbNbr;
287 const int all_dofs =
ndofs;
289 const int elt_dofs =
dof;
294 MFEM_FORALL(i_L, vd*all_dofs+1,
300 for (
int i = 0; i < elt_dofs; i++)
303 const int i_E = e*elt_dofs + i;
304 const int i_L = d_gatherMap[i_E];
305 const int i_offset = d_offsets[i_L];
306 const int i_nextOffset = d_offsets[i_L+1];
307 const int i_nbElts = i_nextOffset - i_offset;
308 for (
int e_i = 0; e_i < i_nbElts; ++e_i)
310 const int i_E = d_indices[i_offset+e_i];
311 i_elts[e_i] = i_E/elt_dofs;
313 for (
int j = 0; j < elt_dofs; j++)
315 const int j_E = e*elt_dofs + j;
316 const int j_L = d_gatherMap[j_E];
317 const int j_offset = d_offsets[j_L];
318 const int j_nextOffset = d_offsets[j_L+1];
319 const int j_nbElts = j_nextOffset - j_offset;
320 if (i_nbElts == 1 || j_nbElts == 1)
322 GetAndIncrementNnzIndex(i_L, I);
327 for (
int e_j = 0; e_j < j_nbElts; ++e_j)
329 const int j_E = d_indices[j_offset+e_j];
330 const int elt = j_E/elt_dofs;
333 int min_e = GetMinElt<Max>(i_elts, i_nbElts, j_elts, j_nbElts);
336 GetAndIncrementNnzIndex(i_L, I);
344 const int nTdofs = vd*all_dofs;
346 for (
int i = 0; i < nTdofs; i++)
348 const int nnz = h_I[i];
360 static constexpr
int Max = MaxNbNbr;
361 const int all_dofs =
ndofs;
363 const int elt_dofs =
dof;
370 auto mat_ea =
Reshape(ea_data.
Read(), elt_dofs, elt_dofs,
ne);
373 for (
int i = 0; i < elt_dofs; i++)
377 const int i_E = e*elt_dofs + i;
378 const int i_L = d_gatherMap[i_E];
379 const int i_offset = d_offsets[i_L];
380 const int i_nextOffset = d_offsets[i_L+1];
381 const int i_nbElts = i_nextOffset - i_offset;
382 for (
int e_i = 0; e_i < i_nbElts; ++e_i)
384 const int i_E = d_indices[i_offset+e_i];
385 i_elts[e_i] = i_E/elt_dofs;
386 i_B[e_i] = i_E%elt_dofs;
388 for (
int j = 0; j < elt_dofs; j++)
390 const int j_E = e*elt_dofs + j;
391 const int j_L = d_gatherMap[j_E];
392 const int j_offset = d_offsets[j_L];
393 const int j_nextOffset = d_offsets[j_L+1];
394 const int j_nbElts = j_nextOffset - j_offset;
395 if (i_nbElts == 1 || j_nbElts == 1)
397 const int nnz = GetAndIncrementNnzIndex(i_L, I);
399 Data[nnz] = mat_ea(j,i,e);
405 for (
int e_j = 0; e_j < j_nbElts; ++e_j)
407 const int j_E = d_indices[j_offset+e_j];
408 const int elt = j_E/elt_dofs;
410 j_B[e_j] = j_E%elt_dofs;
412 int min_e = GetMinElt<Max>(i_elts, i_nbElts, j_elts, j_nbElts);
416 for (
int i = 0; i < i_nbElts; i++)
418 const int e_i = i_elts[i];
419 const int i_Bloc = i_B[i];
420 for (
int j = 0; j < j_nbElts; j++)
422 const int e_j = j_elts[j];
423 const int j_Bloc = j_B[j];
426 val += mat_ea(j_Bloc, i_Bloc, e_i);
430 const int nnz = GetAndIncrementNnzIndex(i_L, I);
441 const int size = vd*all_dofs;
442 for (
int i = 0; i < size; i++)
444 h_I[size-i] = h_I[size-(i+1)];
453 ndof(ne > 0 ? fes.GetFE(0)->GetDof() : 0),
454 ndofs(fes.GetNDofs())
457 width = vdim*ne*ndof;
464 const bool t = byvdim;
465 auto d_x =
Reshape(x.
Read(), t?vd:ndofs, t?ndofs:vd);
467 MFEM_FORALL(i, ndofs,
470 const int dof = idx % nd;
471 const int e = idx / nd;
472 for (
int c = 0; c < vd; ++c)
474 d_y(dof, c, e) = d_x(t?c:idx, t?idx:c);
483 const bool t = byvdim;
486 MFEM_FORALL(i, ndofs,
489 const int dof = idx % nd;
490 const int e = idx / nd;
491 for (
int c = 0; c < vd; ++c)
493 d_y(t?c:idx,t?idx:c) = d_x(dof, c, e);
500 const int elem_dofs = ndof;
503 MFEM_FORALL(dof, ne*elem_dofs*vd,
509 static MFEM_HOST_DEVICE
int AddNnz(
const int iE,
int *I,
const int dofs)
518 const int elem_dofs = ndof;
523 auto mat_ea =
Reshape(ea_data.
Read(), elem_dofs, elem_dofs, ne);
524 MFEM_FORALL(iE, ne*elem_dofs*vd,
526 const int offset = AddNnz(iE,I,elem_dofs);
527 const int e = iE/elem_dofs;
528 const int i = iE%elem_dofs;
529 for (
int j = 0; j < elem_dofs; j++)
531 J[offset+j] = e*elem_dofs+j;
532 Data[offset+j] = mat_ea(j,i,e);
550 faceMap[0] = dof1d-1;
558 for (
int i = 0; i < dof1d; ++i)
564 for (
int i = 0; i < dof1d; ++i)
566 faceMap[i] = dof1d-1 + i*dof1d;
570 for (
int i = 0; i < dof1d; ++i)
572 faceMap[i] = (dof1d-1)*dof1d + i;
576 for (
int i = 0; i < dof1d; ++i)
578 faceMap[i] = i*dof1d;
587 for (
int i = 0; i < dof1d; ++i)
589 for (
int j = 0; j < dof1d; ++j)
591 faceMap[i+j*dof1d] = i + j*dof1d;
596 for (
int i = 0; i < dof1d; ++i)
598 for (
int j = 0; j < dof1d; ++j)
600 faceMap[i+j*dof1d] = i + j*dof1d*dof1d;
605 for (
int i = 0; i < dof1d; ++i)
607 for (
int j = 0; j < dof1d; ++j)
609 faceMap[i+j*dof1d] = dof1d-1 + i*dof1d + j*dof1d*dof1d;
614 for (
int i = 0; i < dof1d; ++i)
616 for (
int j = 0; j < dof1d; ++j)
618 faceMap[i+j*dof1d] = (dof1d-1)*dof1d + i + j*dof1d*dof1d;
623 for (
int i = 0; i < dof1d; ++i)
625 for (
int j = 0; j < dof1d; ++j)
627 faceMap[i+j*dof1d] = i*dof1d + j*dof1d*dof1d;
632 for (
int i = 0; i < dof1d; ++i)
634 for (
int j = 0; j < dof1d; ++j)
636 faceMap[i+j*dof1d] = (dof1d-1)*dof1d*dof1d + i + j*dof1d;
649 nf(fes.GetNFbyType(type)),
652 ndofs(fes.GetNDofs()),
653 dof(nf > 0 ? fes.GetFaceElement(0)->GetDof() : 0),
655 scatter_indices(nf*dof),
657 gather_indices(nf*dof)
659 if (
nf==0) {
return; }
663 MFEM_VERIFY(tfe != NULL &&
666 "Only Gauss-Lobatto and Bernstein basis are supported in "
667 "H1FaceRestriction.");
669 "Non-conforming meshes not yet supported with partial assembly.");
674 if (dof_reorder &&
nf > 0)
676 for (
int f = 0;
f < fes.
GetNF(); ++
f)
681 if (el) {
continue; }
682 mfem_error(
"Finite element not suitable for lexicographic ordering");
688 MFEM_VERIFY(fe_dof_map.
Size() > 0,
"invalid dof map");
694 const int* elementMap = e2dTable.
GetJ();
705 for (
int f = 0;
f < fes.
GetNF(); ++
f)
709 orientation = inf1 % 64;
717 if (orientation != 0)
719 mfem_error(
"FaceRestriction used on degenerated mesh.");
725 mfem_error(
"FaceRestriction not yet implemented for this type of "
729 for (
int d = 0; d <
dof; ++d)
731 const int face_dof = faceMap[d];
732 const int did = (!dof_reorder)?face_dof:dof_map[face_dof];
733 const int gid = elementMap[e1*elem_dofs + did];
734 const int lid = dof*f_ind + d;
740 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
742 for (
int i = 0; i <=
ndofs; ++i)
747 for (
int f = 0;
f < fes.
GetNF(); ++
f)
751 orientation = inf1 % 64;
757 for (
int d = 0; d <
dof; ++d)
759 const int face_dof = faceMap[d];
760 const int did = (!dof_reorder)?face_dof:dof_map[face_dof];
761 const int gid = elementMap[e1*elem_dofs + did];
767 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
768 for (
int i = 1; i <=
ndofs; ++i)
773 for (
int f = 0;
f < fes.
GetNF(); ++
f)
777 orientation = inf1 % 64;
783 for (
int d = 0; d <
dof; ++d)
785 const int face_dof = faceMap[d];
786 const int did = (!dof_reorder)?face_dof:dof_map[face_dof];
787 const int gid = elementMap[e1*elem_dofs + did];
788 const int lid = dof*f_ind + d;
794 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
795 for (
int i =
ndofs; i > 0; --i)
813 const int idx = d_indices[i];
814 const int dof = i % nd;
815 const int face = i / nd;
816 for (
int c = 0; c < vd; ++c)
818 d_y(dof, c, face) = d_x(t?c:idx, t?idx:c);
833 MFEM_FORALL(i,
ndofs,
835 const int offset = d_offsets[i];
836 const int nextOffset = d_offsets[i + 1];
837 for (
int c = 0; c < vd; ++c)
840 for (
int j = offset; j < nextOffset; ++j)
842 const int idx_j = d_indices[j];
843 dofValue += d_x(idx_j % nd, c, idx_j / nd);
845 d_y(t?c:i,t?i:c) += dofValue;
850 static int ToLexOrdering2D(
const int face_id,
const int size1d,
const int i)
852 if (face_id==2 || face_id==3)
862 static int PermuteFace2D(
const int face_id1,
const int face_id2,
863 const int orientation,
864 const int size1d,
const int index)
868 if (face_id1==2 || face_id1==3)
870 new_index = size1d-1-
index;
879 new_index = size1d-1-new_index;
881 return ToLexOrdering2D(face_id2, size1d, new_index);
884 static int ToLexOrdering3D(
const int face_id,
const int size1d,
const int i,
887 if (face_id==2 || face_id==1 || face_id==5)
891 else if (face_id==3 || face_id==4)
893 return (size1d-1-i) + j*size1d;
897 return i + (size1d-1-j)*size1d;
901 static int PermuteFace3D(
const int face_id1,
const int face_id2,
902 const int orientation,
903 const int size1d,
const int index)
905 int i=0, j=0, new_i=0, new_j=0;
909 if (face_id1==3 || face_id1==4)
913 else if (face_id1==0)
930 new_j = (size1d-1-i);
933 new_i = (size1d-1-i);
937 new_i = (size1d-1-i);
938 new_j = (size1d-1-j);
941 new_i = (size1d-1-j);
942 new_j = (size1d-1-i);
945 new_i = (size1d-1-j);
950 new_j = (size1d-1-j);
953 return ToLexOrdering3D(face_id2, size1d, new_i, new_j);
958 const int face_id2,
const int orientation,
959 const int size1d,
const int index)
966 return PermuteFace2D(face_id1, face_id2, orientation, size1d, index);
968 return PermuteFace3D(face_id1, face_id2, orientation, size1d, index);
979 nf(fes.GetNFbyType(type)),
983 ndofs(fes.GetNDofs()),
985 fes.GetTraceElement(0, fes.GetMesh()->GetFaceBaseGeometry(0))->GetDof()
987 elemDofs(fes.GetFE(0)->GetDof()),
990 scatter_indices1(nf*dof),
1006 MFEM_VERIFY(tfe != NULL &&
1009 "Only Gauss-Lobatto and Bernstein basis are supported in "
1010 "L2FaceRestriction.");
1012 "Non-conforming meshes not yet supported with partial assembly.");
1013 if (
nf==0) {
return; }
1019 mfem_error(
"Non-Tensor L2FaceRestriction not yet implemented.");
1021 if (dof_reorder &&
nf > 0)
1023 for (
int f = 0;
f < fes.
GetNF(); ++
f)
1029 if (el) {
continue; }
1030 mfem_error(
"Finite element not suitable for lexicographic ordering");
1034 const int* elementMap = e2dTable.
GetJ();
1038 int face_id1, face_id2;
1045 for (
int f = 0;
f < fes.
GetNF(); ++
f)
1051 orientation = inf1 % 64;
1052 face_id1 = inf1 / 64;
1054 orientation = inf2 % 64;
1055 face_id2 = inf2 / 64;
1060 mfem_error(
"FaceRestriction not yet implemented for this type of "
1067 for (
int d = 0; d <
dof; ++d)
1069 const int face_dof = faceMap1[d];
1070 const int did = face_dof;
1071 const int gid = elementMap[e1*elem_dofs + did];
1072 const int lid = dof*f_ind + d;
1077 for (
int d = 0; d <
dof; ++d)
1082 orientation, dof1d, d);
1083 const int face_dof = faceMap2[pd];
1084 const int did = face_dof;
1085 const int gid = elementMap[e2*elem_dofs + did];
1086 const int lid = dof*f_ind + d;
1091 const int lid = dof*f_ind + d;
1099 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
1101 for (
int i = 0; i <=
ndofs; ++i)
1106 for (
int f = 0;
f < fes.
GetNF(); ++
f)
1113 orientation = inf1 % 64;
1114 face_id1 = inf1 / 64;
1116 orientation = inf2 % 64;
1117 face_id2 = inf2 / 64;
1119 for (
int d = 0; d <
dof; ++d)
1121 const int did = faceMap1[d];
1122 const int gid = elementMap[e1*elem_dofs + did];
1127 for (
int d = 0; d <
dof; ++d)
1132 orientation, dof1d, d);
1133 const int did = faceMap2[pd];
1134 const int gid = elementMap[e2*elem_dofs + did];
1142 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
1143 for (
int i = 1; i <=
ndofs; ++i)
1148 for (
int f = 0;
f < fes.
GetNF(); ++
f)
1155 orientation = inf1 % 64;
1156 face_id1 = inf1 / 64;
1158 orientation = inf2 % 64;
1159 face_id2 = inf2 / 64;
1161 for (
int d = 0; d <
dof; ++d)
1163 const int did = faceMap1[d];
1164 const int gid = elementMap[e1*elem_dofs + did];
1165 const int lid = dof*f_ind + d;
1171 for (
int d = 0; d <
dof; ++d)
1176 orientation, dof1d, d);
1177 const int did = faceMap2[pd];
1178 const int gid = elementMap[e2*elem_dofs + did];
1179 const int lid = dof*f_ind + d;
1188 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
1189 for (
int i =
ndofs; i > 0; --i)
1200 const int vd =
vdim;
1211 const int dof = i % nd;
1212 const int face = i / nd;
1213 const int idx1 = d_indices1[i];
1214 for (
int c = 0; c < vd; ++c)
1216 d_y(dof, c, 0, face) = d_x(t?c:idx1, t?idx1:c);
1218 const int idx2 = d_indices2[i];
1219 for (
int c = 0; c < vd; ++c)
1221 d_y(dof, c, 1, face) = idx2==-1 ? 0.0 : d_x(t?c:idx2, t?idx2:c);
1232 const int dof = i % nd;
1233 const int face = i / nd;
1234 const int idx1 = d_indices1[i];
1235 for (
int c = 0; c < vd; ++c)
1237 d_y(dof, c, face) = d_x(t?c:idx1, t?idx1:c);
1247 const int vd =
vdim;
1257 MFEM_FORALL(i,
ndofs,
1259 const int offset = d_offsets[i];
1260 const int nextOffset = d_offsets[i + 1];
1261 for (
int c = 0; c < vd; ++c)
1263 double dofValue = 0;
1264 for (
int j = offset; j < nextOffset; ++j)
1266 int idx_j = d_indices[j];
1267 bool isE1 = idx_j < dofs;
1268 idx_j = isE1 ? idx_j : idx_j - dofs;
1270 d_x(idx_j % nd, c, 0, idx_j / nd)
1271 :d_x(idx_j % nd, c, 1, idx_j / nd);
1273 d_y(t?c:i,t?i:c) += dofValue;
1281 MFEM_FORALL(i,
ndofs,
1283 const int offset = d_offsets[i];
1284 const int nextOffset = d_offsets[i + 1];
1285 for (
int c = 0; c < vd; ++c)
1287 double dofValue = 0;
1288 for (
int j = offset; j < nextOffset; ++j)
1290 int idx_j = d_indices[j];
1291 dofValue += d_x(idx_j % nd, c, idx_j / nd);
1293 d_y(t?c:i,t?i:c) += dofValue;
1302 const int face_dofs =
dof;
1306 MFEM_FORALL(fdof,
nf*face_dofs,
1308 const int iE1 = d_indices1[fdof];
1309 const int iE2 = d_indices2[fdof];
1310 AddNnz(iE1,I,face_dofs);
1311 AddNnz(iE2,I,face_dofs);
1319 const int face_dofs =
dof;
1323 auto mat_fea =
Reshape(ea_data.
Read(), face_dofs, face_dofs, 2,
nf);
1326 MFEM_FORALL(fdof,
nf*face_dofs,
1328 const int f = fdof/face_dofs;
1329 const int iF = fdof%face_dofs;
1330 const int iE1 = d_indices1[f*face_dofs+iF];
1331 const int iE2 = d_indices2[f*face_dofs+iF];
1332 const int offset1 = AddNnz(iE1,I,face_dofs);
1333 const int offset2 = AddNnz(iE2,I,face_dofs);
1334 for (
int jF = 0; jF < face_dofs; jF++)
1336 const int jE1 = d_indices1[f*face_dofs+jF];
1337 const int jE2 = d_indices2[f*face_dofs+jF];
1338 J[offset2+jF] = jE1;
1339 J[offset1+jF] = jE2;
1340 Data[offset2+jF] = mat_fea(jF,iF,0,f);
1341 Data[offset1+jF] = mat_fea(jF,iF,1,f);
1349 const int face_dofs =
dof;
1356 auto mat_fea =
Reshape(fea_data.
Read(), face_dofs, face_dofs, 2,
nf);
1360 const int e1 = d_indices1[
f*face_dofs]/elem_dofs;
1361 const int e2 = d_indices2[
f*face_dofs]/elem_dofs;
1362 for (
int j = 0; j < face_dofs; j++)
1364 const int jB1 = d_indices1[
f*face_dofs+j]%elem_dofs;
1365 for (
int i = 0; i < face_dofs; i++)
1367 const int iB1 = d_indices1[
f*face_dofs+i]%elem_dofs;
1368 AtomicAdd(mat_ea(iB1,jB1,e1), mat_fea(i,j,0,
f));
1373 for (
int j = 0; j < face_dofs; j++)
1375 const int jB2 = d_indices2[
f*face_dofs+j]%elem_dofs;
1376 for (
int i = 0; i < face_dofs; i++)
1378 const int iB2 = d_indices2[
f*face_dofs+i]%elem_dofs;
1379 AtomicAdd(mat_ea(iB2,jB2,e2), mat_fea(i,j,1,
f));
1388 auto mat_fea =
Reshape(fea_data.
Read(), face_dofs, face_dofs,
nf);
1392 const int e = d_indices[
f*face_dofs]/elem_dofs;
1393 for (
int j = 0; j < face_dofs; j++)
1395 const int jE = d_indices[
f*face_dofs+j]%elem_dofs;
1396 for (
int i = 0; i < face_dofs; i++)
1398 const int iE = d_indices[
f*face_dofs+i]%elem_dofs;
1414 return ToLexOrdering2D(face_id, size1d, index);
1416 return ToLexOrdering3D(face_id, size1d, index%size1d, index/size1d);
Abstract class for all finite elements.
T * HostReadWrite()
Shortcut for mfem::ReadWrite(a.GetMemory(), a.Size(), false).
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
int Size() const
Return the logical size of the array.
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
void MultTransposeUnsigned(const Vector &x, Vector &y) const
Compute MultTranspose without applying signs based on DOF orientations.
const FiniteElement * GetTraceElement(int i, Geometry::Type geom_type) const
Return the trace element from element 'i' to the given 'geom_type'.
virtual void FillI(SparseMatrix &mat, SparseMatrix &face_mat) const
Array< int > gather_indices
const T * HostRead() const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), false).
MFEM_HOST_DEVICE T AtomicAdd(T &add, const T val)
T * GetData()
Returns the data.
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...
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void AddFaceMatricesToElementMatrices(Vector &fea_data, Vector &ea_data) const
This methods adds the DG face matrices to the element matrices.
H1FaceRestriction(const FiniteElementSpace &, const ElementDofOrdering, const FaceType)
Geometry::Type GetFaceBaseGeometry(int i) 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 ...
const FiniteElement * GetFaceElement(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th face in the ...
double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Memory< double > & GetMemoryData()
Memory< int > & GetMemoryI()
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
L2ElementRestriction(const FiniteElementSpace &)
Memory< int > & GetMemoryJ()
DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims...dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
L2FaceRestriction(const FiniteElementSpace &, const FaceType, const L2FaceValues m=L2FaceValues::DoubleValued)
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().
void FillI(SparseMatrix &mat) const
Mesh * GetMesh() const
Returns the mesh.
Operator that extracts Face degrees of freedom.
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 mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
int * ReadWriteI(bool on_dev=true)
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
void MultUnsigned(const Vector &x, Vector &y) const
Compute Mult without applying signs based on DOF orientations.
int * WriteI(bool on_dev=true)
ElementRestriction(const FiniteElementSpace &, ElementDofOrdering)
int GetNF() const
Returns number of faces (i.e. co-dimension 1 entities) in the mesh.
int SpaceDimension() const
int ToLexOrdering(const int dim, const int face_id, const int size1d, const int index)
const FiniteElementSpace & fes
void GetFaceElements(int Face, int *Elem1, int *Elem2) const
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void GetFaceInfos(int Face, int *Inf1, int *Inf2) const
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
void FillJAndData(const Vector &ea_data, SparseMatrix &mat) const
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)
The ordering method used when the number of unknowns per mesh node (vector dimension) is bigger than ...
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 * WriteJ(bool on_dev=true)
int height
Dimension of the output / number of rows in the matrix.
const Array< int > & GetDofMap() const
Get an Array<int> that maps lexicographically ordered indices to the indices of the respective nodes/...
Array< int > scatter_indices1
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
void New(int size)
Allocate host memory for size entries with the current host memory type returned by MemoryManager::Ge...
ElementDofOrdering
Constants describing the possible orderings of the DOFs in one element.
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
int index(int i, int j, int nx, int ny)
Lexicographic ordering for tensor-product FiniteElements.
virtual void FillJAndData(const Vector &ea_data, SparseMatrix &mat, SparseMatrix &face_mat) const
double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
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...
Array< int > scatter_indices
int width
Dimension of the input / number of columns in the matrix.
Array< int > scatter_indices2
void FillSparseMatrix(const Vector &mat_ea, SparseMatrix &mat) const
Fill a Sparse Matrix with Element Matrices.
double f(const Vector &p)
double * WriteData(bool on_dev=true)
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 ...