15 #include "../general/forall.hpp"
27 ndofs(fes.GetNDofs()),
28 dof(ne > 0 ? fes.GetFE(0)->GetDof() : 0),
38 const int *dof_map = NULL;
39 if (dof_reorder &&
ne > 0)
41 for (
int e = 0; e <
ne; ++e)
47 mfem_error(
"Finite element not suitable for lexicographic ordering");
53 MFEM_VERIFY(fe_dof_map.
Size() > 0,
"invalid dof map");
57 const int* elementMap = e2dTable.
GetJ();
59 for (
int i = 0; i <=
ndofs; ++i)
63 for (
int e = 0; e <
ne; ++e)
65 for (
int d = 0; d <
dof; ++d)
67 const int sgid = elementMap[dof*e + d];
68 const int gid = (sgid >= 0) ? sgid : -1 - sgid;
73 for (
int i = 1; i <=
ndofs; ++i)
78 for (
int e = 0; e <
ne; ++e)
80 for (
int d = 0; d <
dof; ++d)
82 const int sdid = dof_reorder ? dof_map[d] : 0;
83 const int did = (!dof_reorder)?d:(sdid >= 0 ? sdid : -1-sdid);
84 const int sgid = elementMap[dof*e + did];
85 const int gid = (sgid >= 0) ? sgid : -1-sgid;
86 const int lid = dof*e + d;
87 const bool plus = (sgid >= 0 && sdid >= 0) || (sgid < 0 && sdid < 0);
94 for (
int i =
ndofs; i > 0; --i)
110 MFEM_FORALL(i,
dof*
ne,
112 const int gid = d_gatherMap[i];
113 const bool plus = gid >= 0;
114 const int j = plus ? gid : -1-gid;
115 for (
int c = 0; c < vd; ++c)
117 const double dofValue = d_x(t?c:j, t?j:c);
118 d_y(i % nd, c, i / nd) = plus ? dofValue : -dofValue;
133 MFEM_FORALL(i,
dof*
ne,
135 const int gid = d_gatherMap[i];
136 const int j = gid >= 0 ? gid : -1-gid;
137 for (
int c = 0; c < vd; ++c)
139 d_y(i % nd, c, i / nd) = d_x(t?c:j, t?j:c);
154 MFEM_FORALL(i,
ndofs,
156 const int offset = d_offsets[i];
157 const int nextOffset = d_offsets[i + 1];
158 for (
int c = 0; c < vd; ++c)
161 for (
int j = offset; j < nextOffset; ++j)
163 const int idx_j = (d_indices[j] >= 0) ? d_indices[j] : -1 - d_indices[j];
164 dofValue += (d_indices[j] >= 0) ? d_x(idx_j % nd, c,
165 idx_j / nd) : -d_x(idx_j % nd, c, idx_j / nd);
167 d_y(t?c:i,t?i:c) = dofValue;
182 MFEM_FORALL(i,
ndofs,
184 const int offset = d_offsets[i];
185 const int nextOffset = d_offsets[i + 1];
186 for (
int c = 0; c < vd; ++c)
189 for (
int j = offset; j < nextOffset; ++j)
191 const int idx_j = (d_indices[j] >= 0) ? d_indices[j] : -1 - d_indices[j];
192 dofValue += d_x(idx_j % nd, c, idx_j / nd);
194 d_y(t?c:i,t?i:c) = dofValue;
209 MFEM_FORALL(i,
ndofs,
211 const int nextOffset = d_offsets[i + 1];
212 for (
int c = 0; c < vd; ++c)
215 const int j = nextOffset - 1;
216 const int idx_j = (d_indices[j] >= 0) ? d_indices[j] : -1 - d_indices[j];
217 dofValue = (d_indices[j] >= 0) ? d_x(idx_j % nd, c, idx_j / nd) :
218 -d_x(idx_j % nd, c, idx_j / nd);
219 d_y(t?c:i,t?i:c) = dofValue;
238 for (
int i = 0; i <
ndofs; ++i)
240 const int offset = d_offsets[i];
241 const int nextOffset = d_offsets[i+1];
242 for (
int c = 0; c < vd; ++c)
244 for (
int j = offset; j < nextOffset; ++j)
246 const int idx_j = d_indices[j];
247 if (d_x(t?c:i,t?i:c))
249 d_y(idx_j % nd, c, idx_j / nd) = 0.0;
253 d_y(idx_j % nd, c, idx_j / nd) = 1.0;
254 d_x(t?c:i,t?i:c) = 1;
265 const int nnz =
FillI(mat);
271 static MFEM_HOST_DEVICE
int GetMinElt(
const int *my_elts,
const int nbElts,
272 const int *nbr_elts,
const int nbrNbElts)
275 int min_el = INT_MAX;
276 for (
int i = 0; i < nbElts; i++)
278 const int e_i = my_elts[i];
279 if (e_i >= min_el) {
continue; }
280 for (
int j = 0; j < nbrNbElts; j++)
282 if (e_i==nbr_elts[j])
294 static MFEM_HOST_DEVICE
int GetAndIncrementNnzIndex(
const int i_L,
int* I)
302 static constexpr
int Max = MaxNbNbr;
303 const int all_dofs =
ndofs;
305 const int elt_dofs =
dof;
310 MFEM_FORALL(i_L, vd*all_dofs+1,
316 for (
int i = 0; i < elt_dofs; i++)
319 const int i_E = e*elt_dofs + i;
320 const int i_L = d_gatherMap[i_E];
321 const int i_offset = d_offsets[i_L];
322 const int i_nextOffset = d_offsets[i_L+1];
323 const int i_nbElts = i_nextOffset - i_offset;
324 for (
int e_i = 0; e_i < i_nbElts; ++e_i)
326 const int i_E = d_indices[i_offset+e_i];
327 i_elts[e_i] = i_E/elt_dofs;
329 for (
int j = 0; j < elt_dofs; j++)
331 const int j_E = e*elt_dofs + j;
332 const int j_L = d_gatherMap[j_E];
333 const int j_offset = d_offsets[j_L];
334 const int j_nextOffset = d_offsets[j_L+1];
335 const int j_nbElts = j_nextOffset - j_offset;
336 if (i_nbElts == 1 || j_nbElts == 1)
338 GetAndIncrementNnzIndex(i_L, I);
343 for (
int e_j = 0; e_j < j_nbElts; ++e_j)
345 const int j_E = d_indices[j_offset+e_j];
346 const int elt = j_E/elt_dofs;
349 int min_e = GetMinElt(i_elts, i_nbElts, j_elts, j_nbElts);
352 GetAndIncrementNnzIndex(i_L, I);
360 const int nTdofs = vd*all_dofs;
362 for (
int i = 0; i < nTdofs; i++)
364 const int nnz = h_I[i];
376 static constexpr
int Max = MaxNbNbr;
377 const int all_dofs =
ndofs;
379 const int elt_dofs =
dof;
386 auto mat_ea =
Reshape(ea_data.
Read(), elt_dofs, elt_dofs,
ne);
389 for (
int i = 0; i < elt_dofs; i++)
393 const int i_E = e*elt_dofs + i;
394 const int i_L = d_gatherMap[i_E];
395 const int i_offset = d_offsets[i_L];
396 const int i_nextOffset = d_offsets[i_L+1];
397 const int i_nbElts = i_nextOffset - i_offset;
398 for (
int e_i = 0; e_i < i_nbElts; ++e_i)
400 const int i_E = d_indices[i_offset+e_i];
401 i_elts[e_i] = i_E/elt_dofs;
402 i_B[e_i] = i_E%elt_dofs;
404 for (
int j = 0; j < elt_dofs; j++)
406 const int j_E = e*elt_dofs + j;
407 const int j_L = d_gatherMap[j_E];
408 const int j_offset = d_offsets[j_L];
409 const int j_nextOffset = d_offsets[j_L+1];
410 const int j_nbElts = j_nextOffset - j_offset;
411 if (i_nbElts == 1 || j_nbElts == 1)
413 const int nnz = GetAndIncrementNnzIndex(i_L, I);
415 Data[nnz] = mat_ea(j,i,e);
421 for (
int e_j = 0; e_j < j_nbElts; ++e_j)
423 const int j_E = d_indices[j_offset+e_j];
424 const int elt = j_E/elt_dofs;
426 j_B[e_j] = j_E%elt_dofs;
428 int min_e = GetMinElt(i_elts, i_nbElts, j_elts, j_nbElts);
432 for (
int i = 0; i < i_nbElts; i++)
434 const int e_i = i_elts[i];
435 const int i_Bloc = i_B[i];
436 for (
int j = 0; j < j_nbElts; j++)
438 const int e_j = j_elts[j];
439 const int j_Bloc = j_B[j];
442 val += mat_ea(j_Bloc, i_Bloc, e_i);
446 const int nnz = GetAndIncrementNnzIndex(i_L, I);
457 const int size = vd*all_dofs;
458 for (
int i = 0; i < size; i++)
460 h_I[size-i] = h_I[size-(i+1)];
469 ndof(ne > 0 ? fes.GetFE(0)->GetDof() : 0),
470 ndofs(fes.GetNDofs())
473 width = vdim*ne*ndof;
480 const bool t = byvdim;
481 auto d_x =
Reshape(x.
Read(), t?vd:ndofs, t?ndofs:vd);
483 MFEM_FORALL(i, ndofs,
486 const int dof = idx % nd;
487 const int e = idx / nd;
488 for (
int c = 0; c < vd; ++c)
490 d_y(dof, c, e) = d_x(t?c:idx, t?idx:c);
499 const bool t = byvdim;
502 MFEM_FORALL(i, ndofs,
505 const int dof = idx % nd;
506 const int e = idx / nd;
507 for (
int c = 0; c < vd; ++c)
509 d_y(t?c:idx,t?idx:c) = d_x(dof, c, e);
516 const int elem_dofs = ndof;
519 const int isize = mat.
Height() + 1;
520 const int interior_dofs = ne*elem_dofs*vd;
521 MFEM_FORALL(dof, isize,
523 I[dof] = dof<interior_dofs ? elem_dofs : 0;
527 static MFEM_HOST_DEVICE
int AddNnz(
const int iE,
int *I,
const int dofs)
536 const int elem_dofs = ndof;
541 auto mat_ea =
Reshape(ea_data.
Read(), elem_dofs, elem_dofs, ne);
542 MFEM_FORALL(iE, ne*elem_dofs*vd,
544 const int offset = AddNnz(iE,I,elem_dofs);
545 const int e = iE/elem_dofs;
546 const int i = iE%elem_dofs;
547 for (
int j = 0; j < elem_dofs; j++)
549 J[offset+j] = e*elem_dofs+j;
550 Data[offset+j] = mat_ea(j,i,e);
568 faceMap[0] = dof1d-1;
576 for (
int i = 0; i < dof1d; ++i)
582 for (
int i = 0; i < dof1d; ++i)
584 faceMap[i] = dof1d-1 + i*dof1d;
588 for (
int i = 0; i < dof1d; ++i)
590 faceMap[i] = (dof1d-1)*dof1d + i;
594 for (
int i = 0; i < dof1d; ++i)
596 faceMap[i] = i*dof1d;
605 for (
int i = 0; i < dof1d; ++i)
607 for (
int j = 0; j < dof1d; ++j)
609 faceMap[i+j*dof1d] = i + j*dof1d;
614 for (
int i = 0; i < dof1d; ++i)
616 for (
int j = 0; j < dof1d; ++j)
618 faceMap[i+j*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] = dof1d-1 + 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 + i + j*dof1d*dof1d;
641 for (
int i = 0; i < dof1d; ++i)
643 for (
int j = 0; j < dof1d; ++j)
645 faceMap[i+j*dof1d] = i*dof1d + j*dof1d*dof1d;
650 for (
int i = 0; i < dof1d; ++i)
652 for (
int j = 0; j < dof1d; ++j)
654 faceMap[i+j*dof1d] = (dof1d-1)*dof1d*dof1d + i + j*dof1d;
667 nf(fes.GetNFbyType(type)),
670 ndofs(fes.GetNDofs()),
671 dof(nf > 0 ? fes.GetFaceElement(0)->GetDof() : 0),
673 scatter_indices(nf*dof),
675 gather_indices(nf*dof)
677 if (
nf==0) {
return; }
681 MFEM_VERIFY(tfe != NULL &&
684 "Only Gauss-Lobatto and Bernstein basis are supported in "
685 "H1FaceRestriction.");
687 "Non-conforming meshes not yet supported with partial assembly.");
692 if (dof_reorder &&
nf > 0)
694 for (
int f = 0;
f < fes.
GetNF(); ++
f)
699 if (el) {
continue; }
700 mfem_error(
"Finite element not suitable for lexicographic ordering");
706 MFEM_VERIFY(fe_dof_map.
Size() > 0,
"invalid dof map");
712 const int* elementMap = e2dTable.
GetJ();
723 for (
int f = 0;
f < fes.
GetNF(); ++
f)
727 orientation = inf1 % 64;
735 if (orientation != 0)
737 mfem_error(
"FaceRestriction used on degenerated mesh.");
743 mfem_error(
"FaceRestriction not yet implemented for this type of "
747 for (
int d = 0; d <
dof; ++d)
749 const int face_dof = faceMap[d];
750 const int did = (!dof_reorder)?face_dof:dof_map[face_dof];
751 const int gid = elementMap[e1*elem_dofs + did];
752 const int lid = dof*f_ind + d;
758 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
760 for (
int i = 0; i <=
ndofs; ++i)
765 for (
int f = 0;
f < fes.
GetNF(); ++
f)
769 orientation = inf1 % 64;
775 for (
int d = 0; d <
dof; ++d)
777 const int face_dof = faceMap[d];
778 const int did = (!dof_reorder)?face_dof:dof_map[face_dof];
779 const int gid = elementMap[e1*elem_dofs + did];
785 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
786 for (
int i = 1; i <=
ndofs; ++i)
791 for (
int f = 0;
f < fes.
GetNF(); ++
f)
795 orientation = inf1 % 64;
801 for (
int d = 0; d <
dof; ++d)
803 const int face_dof = faceMap[d];
804 const int did = (!dof_reorder)?face_dof:dof_map[face_dof];
805 const int gid = elementMap[e1*elem_dofs + did];
806 const int lid = dof*f_ind + d;
812 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
813 for (
int i =
ndofs; i > 0; --i)
831 const int idx = d_indices[i];
832 const int dof = i % nd;
833 const int face = i / nd;
834 for (
int c = 0; c < vd; ++c)
836 d_y(dof, c, face) = d_x(t?c:idx, t?idx:c);
851 MFEM_FORALL(i,
ndofs,
853 const int offset = d_offsets[i];
854 const int nextOffset = d_offsets[i + 1];
855 for (
int c = 0; c < vd; ++c)
858 for (
int j = offset; j < nextOffset; ++j)
860 const int idx_j = d_indices[j];
861 dofValue += d_x(idx_j % nd, c, idx_j / nd);
863 d_y(t?c:i,t?i:c) += dofValue;
868 static int ToLexOrdering2D(
const int face_id,
const int size1d,
const int i)
870 if (face_id==2 || face_id==3)
880 static int PermuteFace2D(
const int face_id1,
const int face_id2,
881 const int orientation,
882 const int size1d,
const int index)
886 if (face_id1==2 || face_id1==3)
888 new_index = size1d-1-
index;
897 new_index = size1d-1-new_index;
899 return ToLexOrdering2D(face_id2, size1d, new_index);
902 static int ToLexOrdering3D(
const int face_id,
const int size1d,
const int i,
905 if (face_id==2 || face_id==1 || face_id==5)
909 else if (face_id==3 || face_id==4)
911 return (size1d-1-i) + j*size1d;
915 return i + (size1d-1-j)*size1d;
919 static int PermuteFace3D(
const int face_id1,
const int face_id2,
920 const int orientation,
921 const int size1d,
const int index)
923 int i=0, j=0, new_i=0, new_j=0;
927 if (face_id1==3 || face_id1==4)
931 else if (face_id1==0)
948 new_j = (size1d-1-i);
951 new_i = (size1d-1-i);
955 new_i = (size1d-1-i);
956 new_j = (size1d-1-j);
959 new_i = (size1d-1-j);
960 new_j = (size1d-1-i);
963 new_i = (size1d-1-j);
968 new_j = (size1d-1-j);
971 return ToLexOrdering3D(face_id2, size1d, new_i, new_j);
976 const int face_id2,
const int orientation,
977 const int size1d,
const int index)
984 return PermuteFace2D(face_id1, face_id2, orientation, size1d, index);
986 return PermuteFace3D(face_id1, face_id2, orientation, size1d, index);
997 nf(fes.GetNFbyType(type)),
1001 ndofs(fes.GetNDofs()),
1003 fes.GetTraceElement(0, fes.
GetMesh()->GetFaceBaseGeometry(0))->GetDof()
1005 elemDofs(fes.GetFE(0)->GetDof()),
1008 scatter_indices1(nf*dof),
1024 MFEM_VERIFY(tfe != NULL &&
1027 "Only Gauss-Lobatto and Bernstein basis are supported in "
1028 "L2FaceRestriction.");
1030 "Non-conforming meshes not yet supported with partial assembly.");
1031 if (
nf==0) {
return; }
1037 mfem_error(
"Non-Tensor L2FaceRestriction not yet implemented.");
1039 if (dof_reorder &&
nf > 0)
1041 for (
int f = 0;
f < fes.
GetNF(); ++
f)
1047 if (el) {
continue; }
1048 mfem_error(
"Finite element not suitable for lexicographic ordering");
1052 const int* elementMap = e2dTable.
GetJ();
1056 int face_id1 = -1, face_id2 = -1;
1057 int orientation = -1;
1063 for (
int f = 0;
f < fes.
GetNF(); ++
f)
1069 orientation = inf1 % 64;
1070 face_id1 = inf1 / 64;
1072 orientation = inf2 % 64;
1073 face_id2 = inf2 / 64;
1078 mfem_error(
"FaceRestriction not yet implemented for this type of "
1085 for (
int d = 0; d <
dof; ++d)
1087 const int face_dof = faceMap1[d];
1088 const int did = face_dof;
1089 const int gid = elementMap[e1*elem_dofs + did];
1090 const int lid = dof*f_ind + d;
1095 for (
int d = 0; d <
dof; ++d)
1100 orientation, dof1d, d);
1101 const int face_dof = faceMap2[pd];
1102 const int did = face_dof;
1103 const int gid = elementMap[e2*elem_dofs + did];
1104 const int lid = dof*f_ind + d;
1109 const int lid = dof*f_ind + d;
1117 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
1119 for (
int i = 0; i <=
ndofs; ++i)
1124 for (
int f = 0;
f < fes.
GetNF(); ++
f)
1131 orientation = inf1 % 64;
1132 face_id1 = inf1 / 64;
1134 orientation = inf2 % 64;
1135 face_id2 = inf2 / 64;
1137 for (
int d = 0; d <
dof; ++d)
1139 const int did = faceMap1[d];
1140 const int gid = elementMap[e1*elem_dofs + did];
1145 for (
int d = 0; d <
dof; ++d)
1150 orientation, dof1d, d);
1151 const int did = faceMap2[pd];
1152 const int gid = elementMap[e2*elem_dofs + did];
1160 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
1161 for (
int i = 1; i <=
ndofs; ++i)
1166 for (
int f = 0;
f < fes.
GetNF(); ++
f)
1173 orientation = inf1 % 64;
1174 face_id1 = inf1 / 64;
1176 orientation = inf2 % 64;
1177 face_id2 = inf2 / 64;
1179 for (
int d = 0; d <
dof; ++d)
1181 const int did = faceMap1[d];
1182 const int gid = elementMap[e1*elem_dofs + did];
1183 const int lid = dof*f_ind + d;
1189 for (
int d = 0; d <
dof; ++d)
1194 orientation, dof1d, d);
1195 const int did = faceMap2[pd];
1196 const int gid = elementMap[e2*elem_dofs + did];
1197 const int lid = dof*f_ind + d;
1206 MFEM_VERIFY(f_ind==
nf,
"Unexpected number of faces.");
1207 for (
int i =
ndofs; i > 0; --i)
1218 const int vd =
vdim;
1229 const int dof = i % nd;
1230 const int face = i / nd;
1231 const int idx1 = d_indices1[i];
1232 for (
int c = 0; c < vd; ++c)
1234 d_y(dof, c, 0, face) = d_x(t?c:idx1, t?idx1:c);
1236 const int idx2 = d_indices2[i];
1237 for (
int c = 0; c < vd; ++c)
1239 d_y(dof, c, 1, face) = idx2==-1 ? 0.0 : d_x(t?c:idx2, t?idx2:c);
1250 const int dof = i % nd;
1251 const int face = i / nd;
1252 const int idx1 = d_indices1[i];
1253 for (
int c = 0; c < vd; ++c)
1255 d_y(dof, c, face) = d_x(t?c:idx1, t?idx1:c);
1265 const int vd =
vdim;
1275 MFEM_FORALL(i,
ndofs,
1277 const int offset = d_offsets[i];
1278 const int nextOffset = d_offsets[i + 1];
1279 for (
int c = 0; c < vd; ++c)
1281 double dofValue = 0;
1282 for (
int j = offset; j < nextOffset; ++j)
1284 int idx_j = d_indices[j];
1285 bool isE1 = idx_j < dofs;
1286 idx_j = isE1 ? idx_j : idx_j - dofs;
1288 d_x(idx_j % nd, c, 0, idx_j / nd)
1289 :d_x(idx_j % nd, c, 1, idx_j / nd);
1291 d_y(t?c:i,t?i:c) += dofValue;
1299 MFEM_FORALL(i,
ndofs,
1301 const int offset = d_offsets[i];
1302 const int nextOffset = d_offsets[i + 1];
1303 for (
int c = 0; c < vd; ++c)
1305 double dofValue = 0;
1306 for (
int j = offset; j < nextOffset; ++j)
1308 int idx_j = d_indices[j];
1309 dofValue += d_x(idx_j % nd, c, idx_j / nd);
1311 d_y(t?c:i,t?i:c) += dofValue;
1318 const bool keep_nbr_block)
const
1320 const int face_dofs =
dof;
1324 MFEM_FORALL(fdof,
nf*face_dofs,
1326 const int iE1 = d_indices1[fdof];
1327 const int iE2 = d_indices2[fdof];
1328 AddNnz(iE1,I,face_dofs);
1329 AddNnz(iE2,I,face_dofs);
1335 const bool keep_nbr_block)
const
1337 const int face_dofs =
dof;
1341 auto mat_fea =
Reshape(ea_data.
Read(), face_dofs, face_dofs, 2,
nf);
1344 MFEM_FORALL(fdof,
nf*face_dofs,
1346 const int f = fdof/face_dofs;
1347 const int iF = fdof%face_dofs;
1348 const int iE1 = d_indices1[f*face_dofs+iF];
1349 const int iE2 = d_indices2[f*face_dofs+iF];
1350 const int offset1 = AddNnz(iE1,I,face_dofs);
1351 const int offset2 = AddNnz(iE2,I,face_dofs);
1352 for (
int jF = 0; jF < face_dofs; jF++)
1354 const int jE1 = d_indices1[f*face_dofs+jF];
1355 const int jE2 = d_indices2[f*face_dofs+jF];
1356 J[offset2+jF] = jE1;
1357 J[offset1+jF] = jE2;
1358 Data[offset2+jF] = mat_fea(jF,iF,0,f);
1359 Data[offset1+jF] = mat_fea(jF,iF,1,f);
1367 const int face_dofs =
dof;
1374 auto mat_fea =
Reshape(fea_data.
Read(), face_dofs, face_dofs, 2,
nf);
1378 const int e1 = d_indices1[
f*face_dofs]/elem_dofs;
1379 const int e2 = d_indices2[
f*face_dofs]/elem_dofs;
1380 for (
int j = 0; j < face_dofs; j++)
1382 const int jB1 = d_indices1[
f*face_dofs+j]%elem_dofs;
1383 for (
int i = 0; i < face_dofs; i++)
1385 const int iB1 = d_indices1[
f*face_dofs+i]%elem_dofs;
1386 AtomicAdd(mat_ea(iB1,jB1,e1), mat_fea(i,j,0,
f));
1391 for (
int j = 0; j < face_dofs; j++)
1393 const int jB2 = d_indices2[
f*face_dofs+j]%elem_dofs;
1394 for (
int i = 0; i < face_dofs; i++)
1396 const int iB2 = d_indices2[
f*face_dofs+i]%elem_dofs;
1397 AtomicAdd(mat_ea(iB2,jB2,e2), mat_fea(i,j,1,
f));
1406 auto mat_fea =
Reshape(fea_data.
Read(), face_dofs, face_dofs,
nf);
1410 const int e = d_indices[
f*face_dofs]/elem_dofs;
1411 for (
int j = 0; j < face_dofs; j++)
1413 const int jE = d_indices[
f*face_dofs+j]%elem_dofs;
1414 for (
int i = 0; i < face_dofs; i++)
1416 const int iE = d_indices[
f*face_dofs+i]%elem_dofs;
1432 return ToLexOrdering2D(face_id, size1d, index);
1434 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 Mult(const Vector &x, Vector &y) const override
Extract the face degrees of freedom from x into y.
void MultTransposeUnsigned(const Vector &x, Vector &y) const
Compute MultTranspose without applying signs based on DOF orientations.
virtual void FillJAndData(const Vector &ea_data, SparseMatrix &mat, const bool keep_nbr_block=false) const
const FiniteElement * GetTraceElement(int i, Geometry::Type geom_type) const
Return the trace element from element 'i' to the given 'geom_type'.
Array< int > gather_indices
void AddMultTranspose(const Vector &x, Vector &y) const override
Add the face degrees of freedom x to the element degrees of freedom y.
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 double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
void AddFaceMatricesToElementMatrices(Vector &fea_data, Vector &ea_data) const
This methods adds the DG face matrices to the element matrices.
void MultLeftInverse(const Vector &x, Vector &y) const
Geometry::Type GetFaceBaseGeometry(int i) const
const FiniteElement * GetFaceElement(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th face in the ...
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 f(const Vector &xvec)
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 on L2 FiniteElementSpaces.
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
void Mult(const Vector &x, Vector &y) const override
Extract the face degrees of freedom from x into y.
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)
H1FaceRestriction(const FiniteElementSpace &fes, const ElementDofOrdering ordering, const FaceType type)
Constructor for a H1FaceRestriction.
void AddMultTranspose(const Vector &x, Vector &y) const override
Add the face degrees of freedom x to the element degrees of freedom y.
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 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 ...
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
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(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.
virtual 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.
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...
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Array< int > scatter_indices
int width
Dimension of the input / number of columns in the matrix.
Array< int > scatter_indices2
virtual void FillI(SparseMatrix &mat, const bool keep_nbr_block=false) const
void FillSparseMatrix(const Vector &mat_ea, SparseMatrix &mat) const
Fill a Sparse Matrix with Element Matrices.
VDIM x NQPT x NE (values) / VDIM x DIM x NQPT x NE (grads)
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 ...