16#ifndef dealii_mesh_worker_simple_h
17#define dealii_mesh_worker_simple_h
57 template <
typename VectorType>
84 template <
class DOFINFO>
94 template <
class DOFINFO>
101 template <
class DOFINFO>
103 assemble(
const DOFINFO &info1,
const DOFINFO &info2);
150 template <
typename MatrixType>
190 template <
class DOFINFO>
198 template <
class DOFINFO>
206 template <
class DOFINFO>
208 assemble(
const DOFINFO &info1,
const DOFINFO &info2);
214 std::vector<SmartPointer<MatrixType, MatrixSimple<MatrixType>>>
matrix;
230 const unsigned int index,
231 const std::vector<types::global_dof_index> &i1,
232 const std::vector<types::global_dof_index> &i2);
252 template <
typename MatrixType>
296 template <
class DOFINFO>
303 template <
class DOFINFO>
311 template <
class DOFINFO>
313 assemble(
const DOFINFO &info1,
const DOFINFO &info2);
322 const std::vector<types::global_dof_index> &i1,
323 const std::vector<types::global_dof_index> &i2);
331 const std::vector<types::global_dof_index> &i1,
332 const std::vector<types::global_dof_index> &i2,
333 const unsigned int level);
342 const std::vector<types::global_dof_index> &i1,
343 const std::vector<types::global_dof_index> &i2,
352 const std::vector<types::global_dof_index> &i1,
353 const std::vector<types::global_dof_index> &i2,
363 const std::vector<types::global_dof_index> &i1,
364 const std::vector<types::global_dof_index> &i2,
374 const std::vector<types::global_dof_index> &i1,
375 const std::vector<types::global_dof_index> &i2,
435 template <
typename MatrixType,
typename VectorType>
469 template <
class DOFINFO>
476 template <
class DOFINFO>
484 template <
class DOFINFO>
486 assemble(
const DOFINFO &info1,
const DOFINFO &info2);
496 const unsigned int index,
497 const std::vector<types::global_dof_index> &indices);
502 const unsigned int index,
503 const std::vector<types::global_dof_index> &i1,
504 const std::vector<types::global_dof_index> &i2);
510 template <
typename VectorType>
519 template <
typename VectorType>
529 template <
typename VectorType>
530 template <
class DOFINFO>
534 info.initialize_vectors(residuals.size());
539 template <
typename VectorType>
540 template <
class DOFINFO>
544 for (
unsigned int k = 0; k < residuals.size(); ++k)
546 VectorType *v = residuals.entry<VectorType *>(k);
547 for (
unsigned int i = 0; i != info.vector(k).n_blocks(); ++i)
549 const std::vector<types::global_dof_index> &ldi =
550 info.vector(k).n_blocks() == 1 ? info.indices :
551 info.indices_by_block[i];
553 if (constraints !=
nullptr)
554 constraints->distribute_local_to_global(info.vector(k).block(i),
558 v->add(ldi, info.vector(k).block(i));
563 template <
typename VectorType>
564 template <
class DOFINFO>
567 const DOFINFO &info2)
576 template <
typename MatrixType>
578 : threshold(threshold)
582 template <
typename MatrixType>
591 template <
typename MatrixType>
595 matrix.resize(m.size());
596 for (
unsigned int i = 0; i < m.size(); ++i)
601 template <
typename MatrixType>
610 template <
typename MatrixType>
611 template <
class DOFINFO>
617 const unsigned int n = info.indices_by_block.size();
620 info.initialize_matrices(matrix.size(), face);
623 info.initialize_matrices(matrix.size() * n * n, face);
625 for (
unsigned int m = 0; m < matrix.size(); ++m)
626 for (
unsigned int i = 0; i < n; ++i)
627 for (
unsigned int j = 0; j < n; ++j, ++k)
629 info.matrix(k,
false).row = i;
630 info.matrix(k,
false).column = j;
633 info.matrix(k,
true).row = i;
634 info.matrix(k,
true).column = j;
642 template <
typename MatrixType>
646 const unsigned int index,
647 const std::vector<types::global_dof_index> &i1,
648 const std::vector<types::global_dof_index> &i2)
653 if (constraints ==
nullptr)
655 for (
unsigned int j = 0; j < i1.size(); ++j)
656 for (
unsigned int k = 0; k < i2.size(); ++k)
657 if (std::fabs(M(j, k)) >= threshold)
658 matrix[index]->add(i1[j], i2[k], M(j, k));
661 constraints->distribute_local_to_global(M, i1, i2, *matrix[index]);
665 template <
typename MatrixType>
666 template <
class DOFINFO>
671 const unsigned int n = info.indices_by_block.size();
674 for (
unsigned int m = 0; m < matrix.size(); ++m)
675 assemble(info.matrix(m,
false).matrix, m, info.indices, info.indices);
678 for (
unsigned int m = 0; m < matrix.size(); ++m)
679 for (
unsigned int k = 0; k < n * n; ++k)
682 info.matrix(k + m * n * n,
false).matrix,
684 info.indices_by_block[info.matrix(k + m * n * n,
false).row],
685 info.indices_by_block[info.matrix(k + m * n * n,
false)
692 template <
typename MatrixType>
693 template <
class DOFINFO>
696 const DOFINFO &info2)
701 info2.indices_by_block.size());
703 const unsigned int n = info1.indices_by_block.size();
707 for (
unsigned int m = 0; m < matrix.size(); ++m)
709 assemble(info1.matrix(m,
false).matrix,
713 assemble(info1.matrix(m,
true).matrix,
717 assemble(info2.matrix(m,
false).matrix,
721 assemble(info2.matrix(m,
true).matrix,
729 for (
unsigned int m = 0; m < matrix.size(); ++m)
730 for (
unsigned int k = 0; k < n * n; ++k)
732 const unsigned int row = info1.matrix(k + m * n * n,
false).row;
733 const unsigned int column =
734 info1.matrix(k + m * n * n,
false).column;
736 assemble(info1.matrix(k + m * n * n,
false).matrix,
738 info1.indices_by_block[row],
739 info1.indices_by_block[column]);
740 assemble(info1.matrix(k + m * n * n,
true).matrix,
742 info1.indices_by_block[row],
743 info2.indices_by_block[column]);
744 assemble(info2.matrix(k + m * n * n,
false).matrix,
746 info2.indices_by_block[row],
747 info2.indices_by_block[column]);
748 assemble(info2.matrix(k + m * n * n,
true).matrix,
750 info2.indices_by_block[row],
751 info1.indices_by_block[column]);
759 template <
typename MatrixType>
761 : threshold(threshold)
765 template <
typename MatrixType>
772 template <
typename MatrixType>
776 mg_constrained_dofs = &c;
780 template <
typename MatrixType>
791 template <
typename MatrixType>
798 interface_out = &out;
802 template <
typename MatrixType>
803 template <
class DOFINFO>
807 const unsigned int n = info.indices_by_block.size();
810 info.initialize_matrices(1, face);
813 info.initialize_matrices(n * n, face);
815 for (
unsigned int i = 0; i < n; ++i)
816 for (
unsigned int j = 0; j < n; ++j, ++k)
818 info.matrix(k,
false).row = i;
819 info.matrix(k,
false).column = j;
822 info.matrix(k,
true).row = i;
823 info.matrix(k,
true).column = j;
830 template <
typename MatrixType>
835 const std::vector<types::global_dof_index> &i1,
836 const std::vector<types::global_dof_index> &i2)
843 for (
unsigned int j = 0; j < i1.size(); ++j)
844 for (
unsigned int k = 0; k < i2.size(); ++k)
845 if (std::fabs(M(j, k)) >= threshold)
846 G.add(i1[j], i2[k], M(j, k));
850 template <
typename MatrixType>
855 const std::vector<types::global_dof_index> &i1,
856 const std::vector<types::global_dof_index> &i2,
857 const unsigned int level)
862 if (mg_constrained_dofs ==
nullptr)
864 for (
unsigned int j = 0; j < i1.size(); ++j)
865 for (
unsigned int k = 0; k < i2.size(); ++k)
866 if (std::fabs(M(j, k)) >= threshold)
867 G.add(i1[j], i2[k], M(j, k));
871 for (
unsigned int j = 0; j < i1.size(); ++j)
872 for (
unsigned int k = 0; k < i2.size(); ++k)
876 if (std::fabs(M(j, k)) < threshold)
886 if (mg_constrained_dofs->at_refinement_edge(
level, i1[j]) ||
887 mg_constrained_dofs->at_refinement_edge(
level, i2[k]))
892 if ((mg_constrained_dofs->is_boundary_index(
level, i1[j]) ||
893 mg_constrained_dofs->is_boundary_index(
level, i2[k])) &&
897 G.add(i1[j], i2[k], M(j, k));
903 template <
typename MatrixType>
908 const std::vector<types::global_dof_index> &i1,
909 const std::vector<types::global_dof_index> &i2,
910 const unsigned int level)
915 if (mg_constrained_dofs ==
nullptr)
917 for (
unsigned int j = 0; j < i1.size(); ++j)
918 for (
unsigned int k = 0; k < i2.size(); ++k)
919 if (std::fabs(M(k, j)) >= threshold)
920 G.add(i1[j], i2[k], M(k, j));
924 for (
unsigned int j = 0; j < i1.size(); ++j)
925 for (
unsigned int k = 0; k < i2.size(); ++k)
926 if (std::fabs(M(k, j)) >= threshold)
927 if (!mg_constrained_dofs->at_refinement_edge(
level, i2[k]))
928 G.add(i1[j], i2[k], M(k, j));
932 template <
typename MatrixType>
937 const std::vector<types::global_dof_index> &i1,
938 const std::vector<types::global_dof_index> &i2,
939 const unsigned int level)
944 if (mg_constrained_dofs ==
nullptr)
946 for (
unsigned int j = 0; j < i1.size(); ++j)
947 for (
unsigned int k = 0; k < i2.size(); ++k)
948 if (std::fabs(M(j, k)) >= threshold)
949 G.add(i1[j], i2[k], M(j, k));
953 for (
unsigned int j = 0; j < i1.size(); ++j)
954 for (
unsigned int k = 0; k < i2.size(); ++k)
955 if (std::fabs(M(j, k)) >= threshold)
956 if (!mg_constrained_dofs->at_refinement_edge(
level, i2[k]))
957 G.add(i1[j], i2[k], M(j, k));
961 template <
typename MatrixType>
966 const std::vector<types::global_dof_index> &i1,
967 const std::vector<types::global_dof_index> &i2,
968 const unsigned int level)
974 for (
unsigned int j = 0; j < i1.size(); ++j)
975 for (
unsigned int k = 0; k < i2.size(); ++k)
976 if (std::fabs(M(j, k)) >= threshold)
988 if (mg_constrained_dofs->at_refinement_edge(
level, i1[j]) &&
989 !mg_constrained_dofs->at_refinement_edge(
level, i2[k]))
991 if ((!mg_constrained_dofs->is_boundary_index(
level, i1[j]) &&
992 !mg_constrained_dofs->is_boundary_index(
level, i2[k])) ||
993 (mg_constrained_dofs->is_boundary_index(
level, i1[j]) &&
994 mg_constrained_dofs->is_boundary_index(
level, i2[k]) &&
996 G.add(i1[j], i2[k], M(j, k));
1001 template <
typename MatrixType>
1006 const std::vector<types::global_dof_index> &i1,
1007 const std::vector<types::global_dof_index> &i2,
1008 const unsigned int level)
1014 for (
unsigned int j = 0; j < i1.size(); ++j)
1015 for (
unsigned int k = 0; k < i2.size(); ++k)
1016 if (std::fabs(M(k, j)) >= threshold)
1017 if (mg_constrained_dofs->at_refinement_edge(
level, i1[j]) &&
1018 !mg_constrained_dofs->at_refinement_edge(
level, i2[k]))
1020 if ((!mg_constrained_dofs->is_boundary_index(
level, i1[j]) &&
1021 !mg_constrained_dofs->is_boundary_index(
level, i2[k])) ||
1022 (mg_constrained_dofs->is_boundary_index(
level, i1[j]) &&
1023 mg_constrained_dofs->is_boundary_index(
level, i2[k]) &&
1025 G.add(i1[j], i2[k], M(k, j));
1030 template <
typename MatrixType>
1031 template <
class DOFINFO>
1036 const unsigned int level = info.cell->level();
1038 if (info.indices_by_block.empty())
1040 assemble((*matrix)[
level],
1041 info.matrix(0,
false).matrix,
1045 if (mg_constrained_dofs !=
nullptr)
1047 assemble_in((*interface_in)[
level],
1048 info.matrix(0,
false).matrix,
1052 assemble_out((*interface_out)[
level],
1053 info.matrix(0,
false).matrix,
1060 for (
unsigned int k = 0; k < info.n_matrices(); ++k)
1062 const unsigned int row = info.matrix(k,
false).row;
1063 const unsigned int column = info.matrix(k,
false).column;
1065 assemble((*matrix)[
level],
1066 info.matrix(k,
false).matrix,
1067 info.indices_by_block[row],
1068 info.indices_by_block[column],
1071 if (mg_constrained_dofs !=
nullptr)
1073 assemble_in((*interface_in)[
level],
1074 info.matrix(k,
false).matrix,
1075 info.indices_by_block[row],
1076 info.indices_by_block[column],
1078 assemble_out((*interface_out)[
level],
1079 info.matrix(k,
false).matrix,
1080 info.indices_by_block[column],
1081 info.indices_by_block[row],
1088 template <
typename MatrixType>
1089 template <
class DOFINFO>
1092 const DOFINFO &info2)
1096 const unsigned int level1 = info1.cell->level();
1097 const unsigned int level2 = info2.cell->level();
1099 if (info1.indices_by_block.empty())
1101 if (level1 == level2)
1103 assemble((*matrix)[level1],
1104 info1.matrix(0,
false).matrix,
1108 assemble((*matrix)[level1],
1109 info1.matrix(0,
true).matrix,
1113 assemble((*matrix)[level1],
1114 info2.matrix(0,
false).matrix,
1118 assemble((*matrix)[level1],
1119 info2.matrix(0,
true).matrix,
1130 assemble((*matrix)[level1],
1131 info1.matrix(0,
false).matrix,
1137 assemble_up((*flux_up)[level1],
1138 info1.matrix(0,
true).matrix,
1142 assemble_down((*flux_down)[level1],
1143 info2.matrix(0,
true).matrix,
1151 for (
unsigned int k = 0; k < info1.n_matrices(); ++k)
1153 const unsigned int row = info1.matrix(k,
false).row;
1154 const unsigned int column = info1.matrix(k,
false).column;
1156 if (level1 == level2)
1158 assemble((*matrix)[level1],
1159 info1.matrix(k,
false).matrix,
1160 info1.indices_by_block[row],
1161 info1.indices_by_block[column],
1163 assemble((*matrix)[level1],
1164 info1.matrix(k,
true).matrix,
1165 info1.indices_by_block[row],
1166 info2.indices_by_block[column],
1168 assemble((*matrix)[level1],
1169 info2.matrix(k,
false).matrix,
1170 info2.indices_by_block[row],
1171 info2.indices_by_block[column],
1173 assemble((*matrix)[level1],
1174 info2.matrix(k,
true).matrix,
1175 info2.indices_by_block[row],
1176 info1.indices_by_block[column],
1185 assemble((*matrix)[level1],
1186 info1.matrix(k,
false).matrix,
1187 info1.indices_by_block[row],
1188 info1.indices_by_block[column],
1192 assemble_up((*flux_up)[level1],
1193 info1.matrix(k,
true).matrix,
1194 info2.indices_by_block[column],
1195 info1.indices_by_block[row],
1197 assemble_down((*flux_down)[level1],
1198 info2.matrix(k,
true).matrix,
1199 info2.indices_by_block[row],
1200 info1.indices_by_block[column],
1209 template <
typename MatrixType,
typename VectorType>
1215 template <
typename MatrixType,
typename VectorType>
1221 VectorType *p = &rhs;
1222 data.
add(p,
"right hand side");
1228 template <
typename MatrixType,
typename VectorType>
1237 template <
typename MatrixType,
typename VectorType>
1238 template <
class DOFINFO>
1247 template <
typename MatrixType,
typename VectorType>
1252 const unsigned int index,
1253 const std::vector<types::global_dof_index> &indices)
1259 VectorType *v = residuals.
entry<VectorType *>(index);
1263 for (
unsigned int i = 0; i < indices.size(); ++i)
1264 (*v)(indices[i]) += vector(i);
1266 for (
unsigned int j = 0; j < indices.size(); ++j)
1267 for (
unsigned int k = 0; k < indices.size(); ++k)
1285 template <
typename MatrixType,
typename VectorType>
1290 const unsigned int index,
1291 const std::vector<types::global_dof_index> &i1,
1292 const std::vector<types::global_dof_index> &i2)
1298 VectorType *v = residuals.
entry<VectorType *>(index);
1302 for (
unsigned int j = 0; j < i1.size(); ++j)
1303 for (
unsigned int k = 0; k < i2.size(); ++k)
1312 vector, i1, i2, *v, M,
false);
1319 template <
typename MatrixType,
typename VectorType>
1320 template <
class DOFINFO>
1327 const unsigned int n = info.indices_by_block.size();
1331 for (
unsigned int m = 0; m < MatrixSimple<MatrixType>::matrix.size();
1333 assemble(info.matrix(m,
false).matrix,
1334 info.vector(m).block(0),
1340 for (
unsigned int m = 0; m < MatrixSimple<MatrixType>::matrix.size();
1342 for (
unsigned int k = 0; k < n * n; ++k)
1344 const unsigned int row = info.matrix(k + m * n * n,
false).row;
1345 const unsigned int column =
1346 info.matrix(k + m * n * n,
false).column;
1349 assemble(info.matrix(k + m * n * n,
false).matrix,
1350 info.vector(m).block(row),
1352 info.indices_by_block[row]);
1354 assemble(info.matrix(k + m * n * n,
false).matrix,
1355 info.vector(m).block(row),
1357 info.indices_by_block[row],
1358 info.indices_by_block[column]);
1364 template <
typename MatrixType,
typename VectorType>
1365 template <
class DOFINFO>
1368 const DOFINFO &info2)
1373 info2.indices_by_block.size());
1375 const unsigned int n = info1.indices_by_block.size();
1379 for (
unsigned int m = 0; m < MatrixSimple<MatrixType>::matrix.size();
1382 assemble(info1.matrix(m,
false).matrix,
1383 info1.vector(m).block(0),
1386 assemble(info1.matrix(m,
true).matrix,
1387 info1.vector(m).block(0),
1391 assemble(info2.matrix(m,
false).matrix,
1392 info2.vector(m).block(0),
1395 assemble(info2.matrix(m,
true).matrix,
1396 info2.vector(m).block(0),
1404 for (
unsigned int m = 0; m < MatrixSimple<MatrixType>::matrix.size();
1406 for (
unsigned int k = 0; k < n * n; ++k)
1408 const unsigned int row = info1.matrix(k + m * n * n,
false).row;
1409 const unsigned int column =
1410 info1.matrix(k + m * n * n,
false).column;
1414 assemble(info1.matrix(k + m * n * n,
false).matrix,
1415 info1.vector(m).block(row),
1417 info1.indices_by_block[row]);
1418 assemble(info2.matrix(k + m * n * n,
false).matrix,
1419 info2.vector(m).block(row),
1421 info2.indices_by_block[row]);
1425 assemble(info1.matrix(k + m * n * n,
false).matrix,
1426 info1.vector(m).block(row),
1428 info1.indices_by_block[row],
1429 info1.indices_by_block[column]);
1430 assemble(info2.matrix(k + m * n * n,
false).matrix,
1431 info2.vector(m).block(row),
1433 info2.indices_by_block[row],
1434 info2.indices_by_block[column]);
1436 assemble(info1.matrix(k + m * n * n,
true).matrix,
1437 info1.vector(m).block(row),
1439 info1.indices_by_block[row],
1440 info2.indices_by_block[column]);
1441 assemble(info2.matrix(k + m * n * n,
true).matrix,
1442 info2.vector(m).block(row),
1444 info2.indices_by_block[row],
1445 info1.indices_by_block[column]);
type entry(const std::string &name)
Access to stored data object by name.
void add(type entry, const std::string &name)
Add a new data object.
void assemble_out(MatrixType &G, const FullMatrix< double > &M, const std::vector< types::global_dof_index > &i1, const std::vector< types::global_dof_index > &i2, const unsigned int level=numbers::invalid_unsigned_int)
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > flux_down
void initialize_fluxes(MGLevelObject< MatrixType > &flux_up, MGLevelObject< MatrixType > &flux_down)
void initialize_interfaces(MGLevelObject< MatrixType > &interface_in, MGLevelObject< MatrixType > &interface_out)
void assemble(const DOFINFO &info)
void initialize_info(DOFINFO &info, bool face) const
void assemble_in(MatrixType &G, const FullMatrix< double > &M, const std::vector< types::global_dof_index > &i1, const std::vector< types::global_dof_index > &i2, const unsigned int level=numbers::invalid_unsigned_int)
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > interface_out
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > matrix
MGMatrixSimple(double threshold=1.e-12)
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > interface_in
void initialize(MGLevelObject< MatrixType > &m)
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > flux_up
void assemble_down(MatrixType &G, const FullMatrix< double > &M, const std::vector< types::global_dof_index > &i1, const std::vector< types::global_dof_index > &i2, const unsigned int level=numbers::invalid_unsigned_int)
void assemble_up(MatrixType &G, const FullMatrix< double > &M, const std::vector< types::global_dof_index > &i1, const std::vector< types::global_dof_index > &i2, const unsigned int level=numbers::invalid_unsigned_int)
SmartPointer< const MGConstrainedDoFs, MGMatrixSimple< MatrixType > > mg_constrained_dofs
SmartPointer< const AffineConstraints< typename MatrixType::value_type >, MatrixSimple< MatrixType > > constraints
void initialize(MatrixType &m)
MatrixSimple(double threshold=1.e-12)
void initialize_info(DOFINFO &info, bool face) const
void assemble(const DOFINFO &info)
std::vector< SmartPointer< MatrixType, MatrixSimple< MatrixType > > > matrix
void assemble(const DOFINFO &info)
void initialize_info(DOFINFO &info, bool face) const
SmartPointer< const AffineConstraints< typename VectorType::value_type >, ResidualSimple< VectorType > > constraints
void initialize(AnyData &results)
SystemSimple(double threshold=1.e-12)
void assemble(const DOFINFO &info)
void initialize(MatrixType &m, VectorType &rhs)
void initialize_info(DOFINFO &info, bool face) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcMessage(std::string arg1)
static const unsigned int invalid_unsigned_int