17#ifndef dealii_mesh_worker_simple_h
18#define dealii_mesh_worker_simple_h
58 template <
typename VectorType>
85 template <
class DOFINFO>
95 template <
class DOFINFO>
102 template <
class DOFINFO>
104 assemble(
const DOFINFO &info1,
const DOFINFO &info2);
151 template <
typename MatrixType>
191 template <
class DOFINFO>
199 template <
class DOFINFO>
207 template <
class DOFINFO>
209 assemble(
const DOFINFO &info1,
const DOFINFO &info2);
215 std::vector<SmartPointer<MatrixType, MatrixSimple<MatrixType>>>
matrix;
231 const unsigned int index,
232 const std::vector<types::global_dof_index> &i1,
233 const std::vector<types::global_dof_index> &i2);
253 template <
typename MatrixType>
297 template <
class DOFINFO>
304 template <
class DOFINFO>
312 template <
class DOFINFO>
314 assemble(
const DOFINFO &info1,
const DOFINFO &info2);
323 const std::vector<types::global_dof_index> &i1,
324 const std::vector<types::global_dof_index> &i2);
332 const std::vector<types::global_dof_index> &i1,
333 const std::vector<types::global_dof_index> &i2,
334 const unsigned int level);
343 const std::vector<types::global_dof_index> &i1,
344 const std::vector<types::global_dof_index> &i2,
353 const std::vector<types::global_dof_index> &i1,
354 const std::vector<types::global_dof_index> &i2,
364 const std::vector<types::global_dof_index> &i1,
365 const std::vector<types::global_dof_index> &i2,
375 const std::vector<types::global_dof_index> &i1,
376 const std::vector<types::global_dof_index> &i2,
436 template <
typename MatrixType,
typename VectorType>
470 template <
class DOFINFO>
477 template <
class DOFINFO>
485 template <
class DOFINFO>
487 assemble(
const DOFINFO &info1,
const DOFINFO &info2);
497 const unsigned int index,
498 const std::vector<types::global_dof_index> &indices);
503 const unsigned int index,
504 const std::vector<types::global_dof_index> &i1,
505 const std::vector<types::global_dof_index> &i2);
511 template <
typename VectorType>
520 template <
typename VectorType>
530 template <
typename VectorType>
531 template <
class DOFINFO>
535 info.initialize_vectors(residuals.size());
540 template <
typename VectorType>
541 template <
class DOFINFO>
545 for (
unsigned int k = 0; k < residuals.size(); ++k)
547 VectorType *v = residuals.entry<VectorType *>(k);
548 for (
unsigned int i = 0; i != info.vector(k).n_blocks(); ++i)
550 const std::vector<types::global_dof_index> &ldi =
551 info.vector(k).n_blocks() == 1 ? info.indices :
552 info.indices_by_block[i];
554 if (constraints !=
nullptr)
555 constraints->distribute_local_to_global(info.vector(k).block(i),
559 v->add(ldi, info.vector(k).block(i));
564 template <
typename VectorType>
565 template <
class DOFINFO>
568 const DOFINFO &info2)
577 template <
typename MatrixType>
579 : threshold(threshold)
583 template <
typename MatrixType>
592 template <
typename MatrixType>
597 for (
unsigned int i = 0; i < m.size(); ++i)
602 template <
typename MatrixType>
611 template <
typename MatrixType>
612 template <
class DOFINFO>
618 const unsigned int n = info.indices_by_block.size();
621 info.initialize_matrices(
matrix.size(), face);
624 info.initialize_matrices(
matrix.size() * n * n, face);
626 for (
unsigned int m = 0; m <
matrix.size(); ++m)
627 for (
unsigned int i = 0; i < n; ++i)
628 for (
unsigned int j = 0; j < n; ++j, ++k)
630 info.matrix(k,
false).row = i;
631 info.matrix(k,
false).column = j;
634 info.matrix(k,
true).row = i;
635 info.matrix(k,
true).column = j;
643 template <
typename MatrixType>
647 const unsigned int index,
648 const std::vector<types::global_dof_index> &i1,
649 const std::vector<types::global_dof_index> &i2)
654 if (constraints ==
nullptr)
656 for (
unsigned int j = 0; j < i1.size(); ++j)
657 for (
unsigned int k = 0; k < i2.size(); ++k)
659 matrix[index]->add(i1[j], i2[k], M(j, k));
662 constraints->distribute_local_to_global(M, i1, i2, *
matrix[index]);
666 template <
typename MatrixType>
667 template <
class DOFINFO>
672 const unsigned int n = info.indices_by_block.size();
675 for (
unsigned int m = 0; m <
matrix.size(); ++m)
676 assemble(info.matrix(m,
false).matrix, m, info.indices, info.indices);
679 for (
unsigned int m = 0; m <
matrix.size(); ++m)
680 for (
unsigned int k = 0; k < n * n; ++k)
683 info.matrix(k + m * n * n,
false).matrix,
685 info.indices_by_block[info.matrix(k + m * n * n,
false).row],
686 info.indices_by_block[info.matrix(k + m * n * n,
false)
693 template <
typename MatrixType>
694 template <
class DOFINFO>
697 const DOFINFO &info2)
702 info2.indices_by_block.size());
704 const unsigned int n = info1.indices_by_block.size();
708 for (
unsigned int m = 0; m <
matrix.size(); ++m)
710 assemble(info1.matrix(m,
false).matrix,
714 assemble(info1.matrix(m,
true).matrix,
718 assemble(info2.matrix(m,
false).matrix,
722 assemble(info2.matrix(m,
true).matrix,
730 for (
unsigned int m = 0; m <
matrix.size(); ++m)
731 for (
unsigned int k = 0; k < n * n; ++k)
733 const unsigned int row = info1.matrix(k + m * n * n,
false).row;
734 const unsigned int column =
735 info1.matrix(k + m * n * n,
false).column;
737 assemble(info1.matrix(k + m * n * n,
false).matrix,
739 info1.indices_by_block[row],
740 info1.indices_by_block[column]);
741 assemble(info1.matrix(k + m * n * n,
true).matrix,
743 info1.indices_by_block[row],
744 info2.indices_by_block[column]);
745 assemble(info2.matrix(k + m * n * n,
false).matrix,
747 info2.indices_by_block[row],
748 info2.indices_by_block[column]);
749 assemble(info2.matrix(k + m * n * n,
true).matrix,
751 info2.indices_by_block[row],
752 info1.indices_by_block[column]);
760 template <
typename MatrixType>
762 : threshold(threshold)
766 template <
typename MatrixType>
773 template <
typename MatrixType>
777 mg_constrained_dofs = &c;
781 template <
typename MatrixType>
792 template <
typename MatrixType>
799 interface_out = &out;
803 template <
typename MatrixType>
804 template <
class DOFINFO>
808 const unsigned int n = info.indices_by_block.size();
811 info.initialize_matrices(1, face);
814 info.initialize_matrices(n * n, face);
816 for (
unsigned int i = 0; i < n; ++i)
817 for (
unsigned int j = 0; j < n; ++j, ++k)
819 info.matrix(k,
false).row = i;
820 info.matrix(k,
false).column = j;
823 info.matrix(k,
true).row = i;
824 info.matrix(k,
true).column = j;
831 template <
typename MatrixType>
836 const std::vector<types::global_dof_index> &i1,
837 const std::vector<types::global_dof_index> &i2)
844 for (
unsigned int j = 0; j < i1.size(); ++j)
845 for (
unsigned int k = 0; k < i2.size(); ++k)
847 G.add(i1[j], i2[k], M(j, k));
851 template <
typename MatrixType>
856 const std::vector<types::global_dof_index> &i1,
857 const std::vector<types::global_dof_index> &i2,
858 const unsigned int level)
863 if (mg_constrained_dofs ==
nullptr)
865 for (
unsigned int j = 0; j < i1.size(); ++j)
866 for (
unsigned int k = 0; k < i2.size(); ++k)
868 G.add(i1[j], i2[k], M(j, k));
872 for (
unsigned int j = 0; j < i1.size(); ++j)
873 for (
unsigned int k = 0; k < i2.size(); ++k)
887 if (mg_constrained_dofs->at_refinement_edge(
level, i1[j]) ||
888 mg_constrained_dofs->at_refinement_edge(
level, i2[k]))
893 if ((mg_constrained_dofs->is_boundary_index(
level, i1[j]) ||
894 mg_constrained_dofs->is_boundary_index(
level, i2[k])) &&
898 G.add(i1[j], i2[k], M(j, k));
904 template <
typename MatrixType>
909 const std::vector<types::global_dof_index> &i1,
910 const std::vector<types::global_dof_index> &i2,
911 const unsigned int level)
916 if (mg_constrained_dofs ==
nullptr)
918 for (
unsigned int j = 0; j < i1.size(); ++j)
919 for (
unsigned int k = 0; k < i2.size(); ++k)
921 G.add(i1[j], i2[k], M(k, j));
925 for (
unsigned int j = 0; j < i1.size(); ++j)
926 for (
unsigned int k = 0; k < i2.size(); ++k)
928 if (!mg_constrained_dofs->at_refinement_edge(
level, i2[k]))
929 G.add(i1[j], i2[k], M(k, j));
933 template <
typename MatrixType>
938 const std::vector<types::global_dof_index> &i1,
939 const std::vector<types::global_dof_index> &i2,
940 const unsigned int level)
945 if (mg_constrained_dofs ==
nullptr)
947 for (
unsigned int j = 0; j < i1.size(); ++j)
948 for (
unsigned int k = 0; k < i2.size(); ++k)
950 G.add(i1[j], i2[k], M(j, k));
954 for (
unsigned int j = 0; j < i1.size(); ++j)
955 for (
unsigned int k = 0; k < i2.size(); ++k)
957 if (!mg_constrained_dofs->at_refinement_edge(
level, i2[k]))
958 G.add(i1[j], i2[k], M(j, k));
962 template <
typename MatrixType>
967 const std::vector<types::global_dof_index> &i1,
968 const std::vector<types::global_dof_index> &i2,
969 const unsigned int level)
975 for (
unsigned int j = 0; j < i1.size(); ++j)
976 for (
unsigned int k = 0; k < i2.size(); ++k)
989 if (mg_constrained_dofs->at_refinement_edge(
level, i1[j]) &&
990 !mg_constrained_dofs->at_refinement_edge(
level, i2[k]))
992 if ((!mg_constrained_dofs->is_boundary_index(
level, i1[j]) &&
993 !mg_constrained_dofs->is_boundary_index(
level, i2[k])) ||
994 (mg_constrained_dofs->is_boundary_index(
level, i1[j]) &&
995 mg_constrained_dofs->is_boundary_index(
level, i2[k]) &&
997 G.add(i1[j], i2[k], M(j, k));
1002 template <
typename MatrixType>
1007 const std::vector<types::global_dof_index> &i1,
1008 const std::vector<types::global_dof_index> &i2,
1009 const unsigned int level)
1015 for (
unsigned int j = 0; j < i1.size(); ++j)
1016 for (
unsigned int k = 0; k < i2.size(); ++k)
1018 if (mg_constrained_dofs->at_refinement_edge(
level, i1[j]) &&
1019 !mg_constrained_dofs->at_refinement_edge(
level, i2[k]))
1021 if ((!mg_constrained_dofs->is_boundary_index(
level, i1[j]) &&
1022 !mg_constrained_dofs->is_boundary_index(
level, i2[k])) ||
1023 (mg_constrained_dofs->is_boundary_index(
level, i1[j]) &&
1024 mg_constrained_dofs->is_boundary_index(
level, i2[k]) &&
1026 G.add(i1[j], i2[k], M(k, j));
1031 template <
typename MatrixType>
1032 template <
class DOFINFO>
1037 const unsigned int level = info.cell->level();
1039 if (info.indices_by_block.size() == 0)
1042 info.matrix(0,
false).matrix,
1046 if (mg_constrained_dofs !=
nullptr)
1048 assemble_in((*interface_in)[
level],
1049 info.matrix(0,
false).matrix,
1053 assemble_out((*interface_out)[
level],
1054 info.matrix(0,
false).matrix,
1061 for (
unsigned int k = 0; k < info.n_matrices(); ++k)
1063 const unsigned int row = info.matrix(k,
false).row;
1064 const unsigned int column = info.matrix(k,
false).column;
1067 info.matrix(k,
false).matrix,
1068 info.indices_by_block[row],
1069 info.indices_by_block[column],
1072 if (mg_constrained_dofs !=
nullptr)
1074 assemble_in((*interface_in)[
level],
1075 info.matrix(k,
false).matrix,
1076 info.indices_by_block[row],
1077 info.indices_by_block[column],
1079 assemble_out((*interface_out)[
level],
1080 info.matrix(k,
false).matrix,
1081 info.indices_by_block[column],
1082 info.indices_by_block[row],
1089 template <
typename MatrixType>
1090 template <
class DOFINFO>
1093 const DOFINFO &info2)
1097 const unsigned int level1 = info1.cell->level();
1098 const unsigned int level2 = info2.cell->level();
1100 if (info1.indices_by_block.size() == 0)
1102 if (level1 == level2)
1105 info1.matrix(0,
false).matrix,
1110 info1.matrix(0,
true).matrix,
1115 info2.matrix(0,
false).matrix,
1120 info2.matrix(0,
true).matrix,
1132 info1.matrix(0,
false).matrix,
1138 assemble_up((*flux_up)[level1],
1139 info1.matrix(0,
true).matrix,
1143 assemble_down((*flux_down)[level1],
1144 info2.matrix(0,
true).matrix,
1152 for (
unsigned int k = 0; k < info1.n_matrices(); ++k)
1154 const unsigned int row = info1.matrix(k,
false).row;
1155 const unsigned int column = info1.matrix(k,
false).column;
1157 if (level1 == level2)
1160 info1.matrix(k,
false).matrix,
1161 info1.indices_by_block[row],
1162 info1.indices_by_block[column],
1165 info1.matrix(k,
true).matrix,
1166 info1.indices_by_block[row],
1167 info2.indices_by_block[column],
1170 info2.matrix(k,
false).matrix,
1171 info2.indices_by_block[row],
1172 info2.indices_by_block[column],
1175 info2.matrix(k,
true).matrix,
1176 info2.indices_by_block[row],
1177 info1.indices_by_block[column],
1187 info1.matrix(k,
false).matrix,
1188 info1.indices_by_block[row],
1189 info1.indices_by_block[column],
1193 assemble_up((*flux_up)[level1],
1194 info1.matrix(k,
true).matrix,
1195 info2.indices_by_block[column],
1196 info1.indices_by_block[row],
1198 assemble_down((*flux_down)[level1],
1199 info2.matrix(k,
true).matrix,
1200 info2.indices_by_block[row],
1201 info1.indices_by_block[column],
1210 template <
typename MatrixType,
typename VectorType>
1216 template <
typename MatrixType,
typename VectorType>
1222 VectorType *p = &rhs;
1223 data.
add(p,
"right hand side");
1229 template <
typename MatrixType,
typename VectorType>
1238 template <
typename MatrixType,
typename VectorType>
1239 template <
class DOFINFO>
1248 template <
typename MatrixType,
typename VectorType>
1253 const unsigned int index,
1254 const std::vector<types::global_dof_index> &indices)
1260 VectorType *v = residuals.
entry<VectorType *>(index);
1264 for (
unsigned int i = 0; i < indices.size(); ++i)
1265 (*v)(indices[i]) += vector(i);
1267 for (
unsigned int j = 0; j < indices.size(); ++j)
1268 for (
unsigned int k = 0; k < indices.size(); ++k)
1286 template <
typename MatrixType,
typename VectorType>
1291 const unsigned int index,
1292 const std::vector<types::global_dof_index> &i1,
1293 const std::vector<types::global_dof_index> &i2)
1299 VectorType *v = residuals.
entry<VectorType *>(index);
1303 for (
unsigned int j = 0; j < i1.size(); ++j)
1304 for (
unsigned int k = 0; k < i2.size(); ++k)
1313 vector, i1, i2, *v, M,
false);
1320 template <
typename MatrixType,
typename VectorType>
1321 template <
class DOFINFO>
1328 const unsigned int n = info.indices_by_block.size();
1334 assemble(info.matrix(m,
false).matrix,
1335 info.vector(m).block(0),
1343 for (
unsigned int k = 0; k < n * n; ++k)
1345 const unsigned int row = info.matrix(k + m * n * n,
false).row;
1346 const unsigned int column =
1347 info.matrix(k + m * n * n,
false).column;
1350 assemble(info.matrix(k + m * n * n,
false).matrix,
1351 info.vector(m).block(row),
1353 info.indices_by_block[row]);
1355 assemble(info.matrix(k + m * n * n,
false).matrix,
1356 info.vector(m).block(row),
1358 info.indices_by_block[row],
1359 info.indices_by_block[column]);
1365 template <
typename MatrixType,
typename VectorType>
1366 template <
class DOFINFO>
1369 const DOFINFO &info2)
1374 info2.indices_by_block.size());
1376 const unsigned int n = info1.indices_by_block.size();
1383 assemble(info1.matrix(m,
false).matrix,
1384 info1.vector(m).block(0),
1387 assemble(info1.matrix(m,
true).matrix,
1388 info1.vector(m).block(0),
1392 assemble(info2.matrix(m,
false).matrix,
1393 info2.vector(m).block(0),
1396 assemble(info2.matrix(m,
true).matrix,
1397 info2.vector(m).block(0),
1407 for (
unsigned int k = 0; k < n * n; ++k)
1409 const unsigned int row = info1.matrix(k + m * n * n,
false).row;
1410 const unsigned int column =
1411 info1.matrix(k + m * n * n,
false).column;
1415 assemble(info1.matrix(k + m * n * n,
false).matrix,
1416 info1.vector(m).block(row),
1418 info1.indices_by_block[row]);
1419 assemble(info2.matrix(k + m * n * n,
false).matrix,
1420 info2.vector(m).block(row),
1422 info2.indices_by_block[row]);
1426 assemble(info1.matrix(k + m * n * n,
false).matrix,
1427 info1.vector(m).block(row),
1429 info1.indices_by_block[row],
1430 info1.indices_by_block[column]);
1431 assemble(info2.matrix(k + m * n * n,
false).matrix,
1432 info2.vector(m).block(row),
1434 info2.indices_by_block[row],
1435 info2.indices_by_block[column]);
1437 assemble(info1.matrix(k + m * n * n,
true).matrix,
1438 info1.vector(m).block(row),
1440 info1.indices_by_block[row],
1441 info2.indices_by_block[column]);
1442 assemble(info2.matrix(k + m * n * n,
true).matrix,
1443 info2.vector(m).block(row),
1445 info2.indices_by_block[row],
1446 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)
Expression fabs(const Expression &x)
@ matrix
Contents is actually a matrix.
void assemble(const MeshWorker::DoFInfoBox< dim, DOFINFO > &dinfo, A *assembler)
static const unsigned int invalid_unsigned_int