22 MFEM_ABORT(
"the assembly level has already been set!");
37 mfem_error(
"Unknown assembly level for this form.");
79 MFEM_ABORT(
"internal MFEM error");
90 MFEM_VERIFY(!
fnfi.Size(),
"Interior faces terms not yet implemented!");
91 MFEM_VERIFY(!
bfnfi.Size(),
"Boundary face terms not yet implemented!");
108 for (
int k = 0; k <
dnfi.Size(); k++)
116 MFEM_ASSERT(marker.
Size() == attr_marker.
Size(),
117 "invalid marker for domain integrator #"
118 << k <<
", counting from zero");
119 for (
int i = 0; i < attr_marker.
Size(); i++)
121 attr_marker[i] |= marker[i];
126 for (
int i = 0; i <
fes->
GetNE(); i++)
129 if (attr_marker[attr-1] == 0) {
continue; }
136 for (
int k = 0; k <
dnfi.Size(); k++)
141 energy +=
dnfi[k]->GetElementEnergy(*fe, *T, el_x);
152 for (
int k = 0; k <
bnfi.Size(); k++)
160 MFEM_ASSERT(bdr_marker.
Size() == bdr_attr_marker.
Size(),
161 "invalid boundary marker for boundary integrator #"
162 << k <<
", counting from zero");
163 for (
int i = 0; i < bdr_attr_marker.
Size(); i++)
165 bdr_attr_marker[i] |= bdr_marker[i];
173 if (bdr_attr_marker[bdr_attr-1] == 0) {
continue; }
180 for (
int k = 0; k <
bnfi.Size(); k++)
185 energy +=
bnfi[k]->GetElementEnergy(*fe, *T, el_x);
193 MFEM_ABORT(
"TODO: add energy contribution from interior face terms");
198 MFEM_ABORT(
"TODO: add energy contribution from boundary face terms");
206 MFEM_VERIFY(x.
Size() ==
Width(),
"invalid input Vector size");
234 mfem::forall(N, [=] MFEM_HOST_DEVICE (
int i) { Y[tdof[i]] = 0.0; });
254 for (
int k = 0; k <
dnfi.Size(); k++)
262 MFEM_ASSERT(marker.
Size() == attr_marker.
Size(),
263 "invalid marker for domain integrator #"
264 << k <<
", counting from zero");
265 for (
int i = 0; i < attr_marker.
Size(); i++)
267 attr_marker[i] |= marker[i];
272 for (
int i = 0; i <
fes->
GetNE(); i++)
275 if (attr_marker[attr-1] == 0) {
continue; }
282 for (
int k = 0; k <
dnfi.Size(); k++)
287 dnfi[k]->AssembleElementVector(*fe, *T, el_x, el_y);
300 for (
int k = 0; k <
bnfi.Size(); k++)
308 MFEM_ASSERT(bdr_marker.
Size() == bdr_attr_marker.
Size(),
309 "invalid boundary marker for boundary integrator #"
310 << k <<
", counting from zero");
311 for (
int i = 0; i < bdr_attr_marker.
Size(); i++)
313 bdr_attr_marker[i] |= bdr_marker[i];
321 if (bdr_attr_marker[bdr_attr-1] == 0) {
continue; }
329 for (
int k = 0; k <
bnfi.Size(); k++)
334 bnfi[k]->AssembleElementVector(*fe, *T, el_x, el_y);
361 for (
int k = 0; k <
fnfi.Size(); k++)
363 fnfi[k]->AssembleFaceVector(*fe1, *fe2, *tr, el_x, el_y);
379 for (
int k = 0; k <
bfnfi.Size(); k++)
387 MFEM_ASSERT(bdr_marker.
Size() == bdr_attr_marker.
Size(),
388 "invalid boundary marker for boundary face integrator #"
389 << k <<
", counting from zero");
390 for (
int i = 0; i < bdr_attr_marker.
Size(); i++)
392 bdr_attr_marker[i] |= bdr_marker[i];
396 for (
int i = 0; i <
fes -> GetNBE(); i++)
399 if (bdr_attr_marker[bdr_attr-1] == 0) {
continue; }
412 for (
int k = 0; k <
bfnfi.Size(); k++)
417 bfnfi[k]->AssembleFaceVector(*fe1, *fe2, *tr, el_x, el_y);
453 const int skip_zeros = 0;
477 for (
int k = 0; k <
dnfi.Size(); k++)
485 MFEM_ASSERT(marker.
Size() == attr_marker.
Size(),
486 "invalid marker for domain integrator #"
487 << k <<
", counting from zero");
488 for (
int i = 0; i < attr_marker.
Size(); i++)
490 attr_marker[i] |= marker[i];
495 for (
int i = 0; i <
fes->
GetNE(); i++)
498 if (attr_marker[attr-1] == 0) {
continue; }
505 for (
int k = 0; k <
dnfi.Size(); k++)
510 dnfi[k]->AssembleElementGrad(*fe, *T, el_x, elmat);
524 for (
int k = 0; k <
bnfi.Size(); k++)
532 MFEM_ASSERT(bdr_marker.
Size() == bdr_attr_marker.
Size(),
533 "invalid boundary marker for boundary integrator #"
534 << k <<
", counting from zero");
535 for (
int i = 0; i < bdr_attr_marker.
Size(); i++)
537 bdr_attr_marker[i] |= bdr_marker[i];
545 if (bdr_attr_marker[bdr_attr-1] == 0) {
continue; }
552 for (
int k = 0; k <
bnfi.Size(); k++)
557 bnfi[k]->AssembleElementGrad(*fe, *T, el_x, elmat);
584 for (
int k = 0; k <
fnfi.Size(); k++)
586 fnfi[k]->AssembleFaceGrad(*fe1, *fe2, *tr, el_x, elmat);
602 for (
int k = 0; k <
bfnfi.Size(); k++)
610 MFEM_ASSERT(bdr_marker.
Size() == bdr_attr_marker.
Size(),
611 "invalid boundary marker for boundary face integrator #"
612 << k <<
", counting from zero");
613 for (
int i = 0; i < bdr_attr_marker.
Size(); i++)
615 bdr_attr_marker[i] |= bdr_marker[i];
619 for (
int i = 0; i <
fes -> GetNBE(); i++)
622 if (bdr_attr_marker[bdr_attr-1] == 0) {
continue; }
635 for (
int k = 0; k <
bfnfi.Size(); k++)
640 bfnfi[k]->AssembleFaceGrad(*fe1, *fe2, *tr, el_x, elmat);
647 if (!finalize) {
return *
Grad; }
700 for (
int i = 0; i <
dnfi.Size(); i++) {
delete dnfi[i]; }
701 for (
int i = 0; i <
bnfi.Size(); i++) {
delete bnfi[i]; }
702 for (
int i = 0; i <
fnfi.Size(); i++) {
delete fnfi[i]; }
703 for (
int i = 0; i <
bfnfi.Size(); i++) {
delete bfnfi[i]; }
710 fes(0), BlockGrad(NULL)
720 for (
int i=0; i<
Grads.NumRows(); ++i)
722 for (
int j=0; j<
Grads.NumCols(); ++j)
728 for (
int i = 0; i <
ess_tdofs.Size(); ++i)
741 for (
int i=0; i<
fes.Size(); ++i)
759 P.SetSize(
fes.Size());
760 cP.SetSize(
fes.Size());
762 for (
int s = 0; s <
fes.Size(); ++s)
765 P[s] =
fes[s]->GetProlongationMatrix();
787 fes(0), BlockGrad(NULL)
795 for (
int s = 0; s <
fes.Size(); ++s)
797 fes[s]->GetEssentialTrueDofs(*bdr_attr_is_ess[s], *
ess_tdofs[s]);
801 rhs[s]->SetSubVector(*
ess_tdofs[s], 0.0);
809 for (
int s = 0; s <
fes.Size(); ++s)
814 rhs[s]->SetSubVector(*
ess_tdofs[s], 0.0);
826 Mesh *mesh =
fes[0]->GetMesh();
829 for (
int i=0; i<
fes.Size(); ++i)
831 el_x_const[i] = el_x[i] =
new Vector();
841 for (
int k = 0; k <
dnfi.Size(); k++)
849 MFEM_ASSERT(marker.
Size() == attr_marker.
Size(),
850 "invalid marker for domain integrator #"
851 << k <<
", counting from zero");
852 for (
int i = 0; i < attr_marker.
Size(); i++)
854 attr_marker[i] |= marker[i];
859 for (
int i = 0; i <
fes[0]->GetNE(); ++i)
862 if (attr_marker[attr-1] == 0) {
continue; }
864 T =
fes[0]->GetElementTransformation(i);
865 for (
int s=0; s<
fes.Size(); ++s)
867 fe[s] =
fes[s]->GetFE(i);
868 fes[s]->GetElementVDofs(i, *vdofs[s], doftrans);
873 for (
int k = 0; k <
dnfi.Size(); ++k)
878 energy +=
dnfi[k]->GetElementEnergy(fe, *T, el_x_const);
889 for (
int k = 0; k <
bnfi.Size(); k++)
897 MFEM_ASSERT(bdr_marker.
Size() == bdr_attr_marker.
Size(),
898 "invalid boundary marker for boundary integrator #"
899 << k <<
", counting from zero");
900 for (
int i = 0; i < bdr_attr_marker.
Size(); i++)
902 bdr_attr_marker[i] |= bdr_marker[i];
907 for (
int i = 0; i < mesh->
GetNBE(); i++)
910 if (bdr_attr_marker[bdr_attr-1] == 0) {
continue; }
912 T =
fes[0]->GetBdrElementTransformation(i);
913 for (
int s = 0; s <
fes.Size(); ++s)
915 fe[s] =
fes[s]->GetBE(i);
916 fes[s]->GetBdrElementVDofs(i, *(vdofs[s]), doftrans);
921 for (
int k = 0; k <
bnfi.Size(); k++)
926 energy +=
bnfi[k]->GetElementEnergy(fe, *T, el_x_const);
932 for (
int i = 0; i <
fes.Size(); ++i)
940 MFEM_ABORT(
"TODO: add energy contribution from interior face terms");
945 MFEM_ABORT(
"TODO: add energy contribution from boundary face terms");
968 std::vector<DofTransformation> doftrans(
fes.Size());
969 Mesh *mesh =
fes[0]->GetMesh();
974 for (
int s=0; s<
fes.Size(); ++s)
976 el_x_const[s] = el_x[s] =
new Vector();
988 for (
int k = 0; k <
dnfi.Size(); k++)
996 MFEM_ASSERT(marker.
Size() == attr_marker.
Size(),
997 "invalid marker for domain integrator #"
998 << k <<
", counting from zero");
999 for (
int i = 0; i < attr_marker.
Size(); i++)
1001 attr_marker[i] |= marker[i];
1005 for (
int i = 0; i <
fes[0]->GetNE(); ++i)
1008 if (attr_marker[attr-1] == 0) {
continue; }
1010 T =
fes[0]->GetElementTransformation(i);
1011 for (
int s = 0; s <
fes.Size(); ++s)
1013 fes[s]->GetElementVDofs(i, *(vdofs[s]), doftrans[s]);
1014 fe[s] =
fes[s]->GetFE(i);
1016 doftrans[s].InvTransformPrimal(*el_x[s]);
1019 for (
int k = 0; k <
dnfi.Size(); ++k)
1024 dnfi[k]->AssembleElementVector(fe, *T,
1027 for (
int s=0; s<
fes.Size(); ++s)
1029 if (el_y[s]->Size() == 0) {
continue; }
1030 doftrans[s].TransformDual(*el_y[s]);
1042 bdr_attr_marker = 0;
1043 for (
int k = 0; k <
bnfi.Size(); k++)
1047 bdr_attr_marker = 1;
1051 MFEM_ASSERT(bdr_marker.
Size() == bdr_attr_marker.
Size(),
1052 "invalid boundary marker for boundary integrator #"
1053 << k <<
", counting from zero");
1054 for (
int i = 0; i < bdr_attr_marker.
Size(); i++)
1056 bdr_attr_marker[i] |= bdr_marker[i];
1060 for (
int i = 0; i < mesh->
GetNBE(); i++)
1063 if (bdr_attr_marker[bdr_attr-1] == 0) {
continue; }
1065 T =
fes[0]->GetBdrElementTransformation(i);
1066 for (
int s = 0; s <
fes.Size(); ++s)
1068 fes[s]->GetBdrElementVDofs(i, *(vdofs[s]), doftrans[s]);
1069 fe[s] =
fes[s]->GetBE(i);
1071 doftrans[s].InvTransformPrimal(*el_x[s]);
1074 for (
int k = 0; k <
bnfi.Size(); k++)
1079 bnfi[k]->AssembleElementVector(fe, *T, el_x_const, el_y);
1081 for (
int s=0; s<
fes.Size(); ++s)
1083 if (el_y[s]->Size() == 0) {
continue; }
1084 doftrans[s].TransformDual(*el_y[s]);
1101 for (
int s=0; s<
fes.Size(); ++s)
1103 fe[s] =
fes[s]->GetFE(tr->Elem1No);
1104 fe2[s] =
fes[s]->GetFE(tr->Elem2No);
1106 fes[s]->GetElementVDofs(tr->Elem1No, *(vdofs[s]));
1107 fes[s]->GetElementVDofs(tr->Elem2No, *(vdofs2[s]));
1109 vdofs[s]->
Append(*(vdofs2[s]));
1114 for (
int k = 0; k <
fnfi.Size(); ++k)
1117 fnfi[k]->AssembleFaceVector(fe, fe2, *tr, el_x_const, el_y);
1119 for (
int s=0; s<
fes.Size(); ++s)
1121 if (el_y[s]->Size() == 0) {
continue; }
1136 bdr_attr_marker = 0;
1137 for (
int k = 0; k <
bfnfi.Size(); ++k)
1141 bdr_attr_marker = 1;
1145 MFEM_ASSERT(bdr_marker.
Size() == bdr_attr_marker.
Size(),
1146 "invalid boundary marker for boundary face integrator #"
1147 << k <<
", counting from zero");
1148 for (
int i = 0; i < bdr_attr_marker.
Size(); ++i)
1150 bdr_attr_marker[i] |= bdr_marker[i];
1154 for (
int i = 0; i < mesh->
GetNBE(); ++i)
1157 if (bdr_attr_marker[bdr_attr-1] == 0) {
continue; }
1162 for (
int s=0; s<
fes.Size(); ++s)
1164 fe[s] =
fes[s]->GetFE(tr->Elem1No);
1165 fe2[s] =
fes[s]->GetFE(tr->Elem1No);
1167 fes[s]->GetElementVDofs(tr->Elem1No, *(vdofs[s]));
1171 for (
int k = 0; k <
bfnfi.Size(); ++k)
1176 bfnfi[k]->AssembleFaceVector(fe, fe2, *tr, el_x_const, el_y);
1178 for (
int s=0; s<
fes.Size(); ++s)
1180 if (el_y[s]->Size() == 0) {
continue; }
1188 for (
int s=0; s<
fes.Size(); ++s)
1201 MFEM_VERIFY(bx.
Size() ==
Width(),
"invalid input BlockVector size");
1206 for (
int s = 0; s <
fes.Size(); s++)
1238 for (
int s = 0; s <
fes.Size(); s++)
1253 bool finalize)
const
1255 const int skip_zeros = 0;
1264 std::vector<DofTransformation> doftrans(
fes.Size());
1265 Mesh *mesh =
fes[0]->GetMesh();
1267 for (
int i=0; i<
fes.Size(); ++i)
1269 el_x_const[i] = el_x[i] =
new Vector();
1272 for (
int j=0; j<
fes.Size(); ++j)
1278 for (
int i=0; i<
fes.Size(); ++i)
1280 for (
int j=0; j<
fes.Size(); ++j)
1282 if (
Grads(i,j) != NULL)
1289 fes[j]->GetVSize());
1300 for (
int k = 0; k <
dnfi.Size(); k++)
1308 MFEM_ASSERT(marker.
Size() == attr_marker.
Size(),
1309 "invalid marker for domain integrator #"
1310 << k <<
", counting from zero");
1311 for (
int i = 0; i < attr_marker.
Size(); i++)
1313 attr_marker[i] |= marker[i];
1317 for (
int i = 0; i <
fes[0]->GetNE(); ++i)
1320 if (attr_marker[attr-1] == 0) {
continue; }
1322 T =
fes[0]->GetElementTransformation(i);
1323 for (
int s = 0; s <
fes.Size(); ++s)
1325 fe[s] =
fes[s]->GetFE(i);
1326 fes[s]->GetElementVDofs(i, *vdofs[s], doftrans[s]);
1328 doftrans[s].InvTransformPrimal(*el_x[s]);
1331 for (
int k = 0; k <
dnfi.Size(); ++k)
1336 dnfi[k]->AssembleElementGrad(fe, *T, el_x_const, elmats);
1338 for (
int j=0; j<
fes.Size(); ++j)
1340 for (
int l=0; l<
fes.Size(); ++l)
1342 if (elmats(j,l)->
Height() == 0) {
continue; }
1344 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1345 *elmats(j,l), skip_zeros);
1357 bdr_attr_marker = 0;
1358 for (
int k = 0; k <
bnfi.Size(); k++)
1362 bdr_attr_marker = 1;
1366 MFEM_ASSERT(bdr_marker.
Size() == bdr_attr_marker.
Size(),
1367 "invalid boundary marker for boundary integrator #"
1368 << k <<
", counting from zero");
1369 for (
int i = 0; i < bdr_attr_marker.
Size(); i++)
1371 bdr_attr_marker[i] |= bdr_marker[i];
1375 for (
int i = 0; i < mesh->
GetNBE(); i++)
1378 if (bdr_attr_marker[bdr_attr-1] == 0) {
continue; }
1380 T =
fes[0]->GetBdrElementTransformation(i);
1381 for (
int s = 0; s <
fes.Size(); ++s)
1383 fe[s] =
fes[s]->GetBE(i);
1384 fes[s]->GetBdrElementVDofs(i, *(vdofs[s]), doftrans[s]);
1386 doftrans[s].InvTransformPrimal(*el_x[s]);
1389 for (
int k = 0; k <
bnfi.Size(); k++)
1394 bnfi[k]->AssembleElementGrad(fe, *T, el_x_const, elmats);
1396 for (
int j=0; j<
fes.Size(); ++j)
1398 for (
int l=0; l<
fes.Size(); ++l)
1400 if (elmats(j,l)->
Height() == 0) {
continue; }
1402 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1403 *elmats(j,l), skip_zeros);
1418 for (
int s=0; s <
fes.Size(); ++s)
1420 fe[s] =
fes[s]->GetFE(tr->Elem1No);
1421 fe2[s] =
fes[s]->GetFE(tr->Elem2No);
1423 fes[s]->GetElementVDofs(tr->Elem1No, *vdofs[s]);
1424 fes[s]->GetElementVDofs(tr->Elem2No, *vdofs2[s]);
1425 vdofs[s]->
Append(*(vdofs2[s]));
1430 for (
int k = 0; k <
fnfi.Size(); ++k)
1432 fnfi[k]->AssembleFaceGrad(fe, fe2, *tr, el_x_const, elmats);
1433 for (
int j=0; j<
fes.Size(); ++j)
1435 for (
int l=0; l<
fes.Size(); ++l)
1437 if (elmats(j,l)->
Height() == 0) {
continue; }
1438 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1439 *elmats(j,l), skip_zeros);
1453 bdr_attr_marker = 0;
1454 for (
int k = 0; k <
bfnfi.Size(); ++k)
1458 bdr_attr_marker = 1;
1462 MFEM_ASSERT(bdr_marker.
Size() == bdr_attr_marker.
Size(),
1463 "invalid boundary marker for boundary face integrator #"
1464 << k <<
", counting from zero");
1465 for (
int i = 0; i < bdr_attr_marker.
Size(); ++i)
1467 bdr_attr_marker[i] |= bdr_marker[i];
1471 for (
int i = 0; i < mesh->
GetNBE(); ++i)
1474 if (bdr_attr_marker[bdr_attr-1] == 0) {
continue; }
1479 for (
int s = 0; s <
fes.Size(); ++s)
1481 fe[s] =
fes[s]->GetFE(tr->Elem1No);
1484 fes[s]->GetElementVDofs(tr->Elem1No, *vdofs[s]);
1488 for (
int k = 0; k <
bfnfi.Size(); ++k)
1492 bfnfi[k]->AssembleFaceGrad(fe, fe2, *tr, el_x_const, elmats);
1493 for (
int l=0; l<
fes.Size(); ++l)
1495 for (
int j=0; j<
fes.Size(); ++j)
1497 if (elmats(j,l)->
Height() == 0) {
continue; }
1498 Grads(j,l)->AddSubMatrix(*vdofs[j], *vdofs[l],
1499 *elmats(j,l), skip_zeros);
1507 if (finalize && !
Grads(0,0)->Finalized())
1509 for (
int i=0; i<
fes.Size(); ++i)
1511 for (
int j=0; j<
fes.Size(); ++j)
1513 Grads(i,j)->Finalize(skip_zeros);
1518 for (
int i=0; i<
fes.Size(); ++i)
1520 for (
int j=0; j<
fes.Size(); ++j)
1541 for (
int s1 = 0; s1 <
fes.Size(); ++s1)
1543 for (
int s2 = 0; s2 <
fes.Size(); ++s2)
1546 if (
cP[s1] &&
cP[s2])
1563 mGrads(s1, s2) =
cGrads(s1, s2);
1568 for (
int s = 0; s <
fes.Size(); ++s)
1570 for (
int i = 0; i <
ess_tdofs[s]->Size(); ++i)
1572 for (
int j = 0; j <
fes.Size(); ++j)
1576 mGrads(s, s)->EliminateRowCol((*
ess_tdofs[s])[i],
1581 mGrads(s, j)->EliminateRow((*
ess_tdofs[s])[i]);
1582 mGrads(j, s)->EliminateCol((*
ess_tdofs[s])[i]);
1590 for (
int i = 0; i <
fes.Size(); ++i)
1592 for (
int j = 0; j <
fes.Size(); ++j)
1603 for (
int i=0; i<
fes.Size(); ++i)
1605 for (
int j=0; j<
fes.Size(); ++j)
1613 for (
int i = 0; i <
dnfi.Size(); ++i)
1618 for (
int i = 0; i <
bnfi.Size(); ++i)
1623 for (
int i = 0; i <
fnfi.Size(); ++i)
1628 for (
int i = 0; i <
bfnfi.Size(); ++i)
Dynamic 2D array using row-major layout.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
const T * HostRead() const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), false).
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int Size() const
Return the logical size of the array.
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
int Append(const T &el)
Append element 'el' to array, resize if necessary.
void Copy(Array ©) const
Create a copy of the internal array to the provided copy.
const T * Read(bool on_dev=true) const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), on_dev).
A class to handle Block systems in a matrix-free implementation.
void SetBlock(int iRow, int iCol, Operator *op, real_t c=1.0)
Add a block op in the block-entry (iblock, jblock).
A class to handle Vectors in a block fashion.
void Update(real_t *data, const Array< int > &bOffsets)
Update method.
void SyncFromBlocks() const
Synchronize the memory location flags (i.e. the memory validity flags) of the big/monolithic block-ve...
void SyncToBlocks() const
Synchronize the memory location flags (i.e. the memory validity flags) of the big/monolithic block-ve...
Vector & GetBlock(int i)
Get the i-th vector in the block.
Data type dense matrix using column-major storage.
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th boundary fac...
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
ElementTransformation * GetElementTransformation(int i) const
static void ListToMarker(const Array< int > &list, int marker_size, Array< int > &marker, int mark_val=-1)
Convert an array of indices (list) to a Boolean marker array where all indices in the list are marked...
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
int GetNBE() const
Returns number of boundary elements in the mesh.
virtual const Operator * GetProlongationMatrix() const
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i-th boundary element.
int GetNE() const
Returns number of elements in the mesh.
static void MarkerToList(const Array< int > &marker, Array< int > &list)
Convert a Boolean marker array to a list containing all marked indices.
void ConvertToConformingVDofs(const Array< int > &dofs, Array< int > &cdofs)
For a partially conforming FE space, convert a marker array (nonzero entries are true) on the partial...
Mesh * GetMesh() const
Returns the mesh.
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
DofTransformation * GetBdrElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for i'th boundary element. The returned indices are offsets int...
Abstract class for all finite elements.
void BooleanMultTranspose(int alpha, const int *x, int beta, int *y)
The "Boolean" analog of y = alpha * A^T * x + beta * y, where elements in the sparsity pattern of the...
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
int GetAttribute(int i) const
Return the attribute of element i.
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Builds the transformation defining the given boundary face.
int GetNBE() const
Returns number of boundary elements.
FaceElementTransformations * GetInteriorFaceTransformations(int FaceNo)
See GetFaceElementTransformations().
Array< int > attributes
A list of all unique element attributes used by the Mesh.
void Clear()
Clear the OperatorHandle, deleting the held Operator (if owned), while leaving the type id unchanged.
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
int width
Dimension of the input / number of columns in the matrix.
void FormSystemOperator(const Array< int > &ess_tdof_list, Operator *&A)
Return in A a parallel (on truedofs) version of this square operator.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
int height
Dimension of the output / number of rows in the matrix.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
@ DIAG_ONE
Set the diagonal value to one.
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Abstract parallel finite element space.
int GetTrueVSize() const override
Return the number of local vector true dofs.
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
void EliminateRowCol(int rc, const real_t sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Eliminate row rc and column rc and modify the rhs using sol.
void MultTranspose(const Vector &x, Vector &y) const override
Multiply a vector with the transposed matrix. y = At * x.
bool Finalized() const
Returns whether or not CSR format has been finalized.
void AddSubMatrix(const Array< int > &rows, const Array< int > &cols, const DenseMatrix &subm, int skip_zeros=1)
void Finalize(int skip_zeros=1) override
Finalize the matrix initialization, switching the storage format from LIL to CSR.
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add elements of the elemvect Vector to the entries listed in dofs. Negative dof values cause the -dof...
int Size() const
Returns the size of the vector.
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
void SetSize(int s)
Resize the vector to size s.
virtual real_t * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
void mfem_error(const char *msg)
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
void TransformDual(const DofTransformation &ran_dof_trans, const DofTransformation &dom_dof_trans, DenseMatrix &elmat)
AssemblyLevel
Enumeration defining the assembly level for bilinear and nonlinear form classes derived from Operator...
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
std::function< real_t(const Vector &)> f(real_t mass_coeff)
SparseMatrix * TransposeMult(const SparseMatrix &A, const SparseMatrix &B)
C = A^T B.
void forall(int N, lambda &&body)