17 #ifndef dealii_mesh_worker_simple_h 18 #define dealii_mesh_worker_simple_h 20 #include <deal.II/algorithms/any_data.h> 21 #include <deal.II/base/smartpointer.h> 22 #include <deal.II/base/mg_level_object.h> 23 #include <deal.II/lac/block_vector.h> 24 #include <deal.II/meshworker/dof_info.h> 25 #include <deal.II/meshworker/functional.h> 26 #include <deal.II/multigrid/mg_constrained_dofs.h> 34 DEAL_II_NAMESPACE_OPEN
53 template <
typename VectorType>
92 template <
class DOFINFO>
101 template <
class DOFINFO>
107 template <
class DOFINFO>
109 const DOFINFO &info2);
154 template <
typename MatrixType>
190 template <
class DOFINFO>
197 template <
class DOFINFO>
204 template <
class DOFINFO>
206 const DOFINFO &info2);
211 std::vector<SmartPointer<MatrixType,MatrixSimple<MatrixType> > >
matrix;
226 const unsigned int index,
227 const std::vector<types::global_dof_index> &i1,
228 const std::vector<types::global_dof_index> &i2);
248 template <
typename MatrixType>
288 template <
class DOFINFO>
294 template <
class DOFINFO>
301 template <
class DOFINFO>
303 const DOFINFO &info2);
310 const std::vector<types::global_dof_index> &i1,
311 const std::vector<types::global_dof_index> &i2);
318 const std::vector<types::global_dof_index> &i1,
319 const std::vector<types::global_dof_index> &i2,
320 const unsigned int level);
328 const std::vector<types::global_dof_index> &i1,
329 const std::vector<types::global_dof_index> &i2,
337 const std::vector<types::global_dof_index> &i1,
338 const std::vector<types::global_dof_index> &i2,
347 const std::vector<types::global_dof_index> &i1,
348 const std::vector<types::global_dof_index> &i2,
357 const std::vector<types::global_dof_index> &i1,
358 const std::vector<types::global_dof_index> &i2,
414 template <
typename MatrixType,
typename VectorType>
428 void initialize(MatrixType &m, VectorType &rhs);
446 template <
class DOFINFO>
452 template <
class DOFINFO>
459 template <
class DOFINFO>
461 const DOFINFO &info2);
470 const unsigned int index,
471 const std::vector<types::global_dof_index> &indices);
475 const unsigned int index,
476 const std::vector<types::global_dof_index> &i1,
477 const std::vector<types::global_dof_index> &i2);
483 template <
typename VectorType>
490 template <
typename VectorType>
498 template <
typename MatrixType>
504 template <
typename VectorType>
505 template <
class DOFINFO>
509 info.initialize_vectors(residuals.size());
513 template <
typename VectorType>
514 template <
class DOFINFO>
518 for (
unsigned int k=0; k<residuals.size(); ++k)
520 VectorType *v = residuals.entry<VectorType *>(k);
521 for (
unsigned int i=0; i != info.vector(k).n_blocks(); ++i)
523 const std::vector<types::global_dof_index> &ldi = info.vector(k).n_blocks()==1?
525 info.indices_by_block[i];
527 if (constraints !=
nullptr)
528 constraints->distribute_local_to_global(info.vector(k).block(i), ldi, *v);
530 v->add(ldi, info.vector(k).block(i));
535 template <
typename VectorType>
536 template <
class DOFINFO>
539 const DOFINFO &info2)
548 template <
typename MatrixType>
556 template <
typename MatrixType>
565 template <
typename MatrixType>
569 matrix.resize(m.size());
570 for (
unsigned int i=0; i<m.size(); ++i)
575 template <
typename MatrixType>
583 template <
typename MatrixType >
584 template <
class DOFINFO>
590 const unsigned int n = info.indices_by_block.size();
593 info.initialize_matrices(matrix.size(), face);
596 info.initialize_matrices(matrix.size()*n*n, face);
598 for (
unsigned int m=0; m<matrix.size(); ++m)
599 for (
unsigned int i=0; i<n; ++i)
600 for (
unsigned int j=0; j<n; ++j,++k)
602 info.matrix(k,
false).row = i;
603 info.matrix(k,
false).column = j;
606 info.matrix(k,
true).row = i;
607 info.matrix(k,
true).column = j;
615 template <
typename MatrixType>
618 const unsigned int index,
619 const std::vector<types::global_dof_index> &i1,
620 const std::vector<types::global_dof_index> &i2)
625 if (constraints ==
nullptr)
627 for (
unsigned int j=0; j<i1.size(); ++j)
628 for (
unsigned int k=0; k<i2.size(); ++k)
629 if (std::fabs(M(j,k)) >= threshold)
630 matrix[index]->add(i1[j], i2[k], M(j,k));
633 constraints->distribute_local_to_global(M, i1, i2, *matrix[index]);
637 template <
typename MatrixType>
638 template <
class DOFINFO>
643 const unsigned int n = info.indices_by_block.size();
646 for (
unsigned int m=0; m<matrix.size(); ++m)
647 assemble(info.matrix(m,
false).matrix, m, info.indices, info.indices);
650 for (
unsigned int m=0; m<matrix.size(); ++m)
651 for (
unsigned int k=0; k<n*n; ++k)
653 assemble(info.matrix(k+m*n*n,
false).matrix, m,
654 info.indices_by_block[info.matrix(k+m*n*n,
false).row],
655 info.indices_by_block[info.matrix(k+m*n*n,
false).column]);
661 template <
typename MatrixType>
662 template <
class DOFINFO>
668 AssertDimension(info1.indices_by_block.size(),info2.indices_by_block.size());
670 const unsigned int n = info1.indices_by_block.size();
674 for (
unsigned int m=0; m<matrix.size(); ++m)
676 assemble(info1.matrix(m,
false).matrix, m, info1.indices, info1.indices);
677 assemble(info1.matrix(m,
true).matrix, m, info1.indices, info2.indices);
678 assemble(info2.matrix(m,
false).matrix, m, info2.indices, info2.indices);
679 assemble(info2.matrix(m,
true).matrix, m, info2.indices, info1.indices);
684 for (
unsigned int m=0; m<matrix.size(); ++m)
685 for (
unsigned int k=0; k<n*n; ++k)
687 const unsigned int row = info1.matrix(k+m*n*n,
false).row;
688 const unsigned int column = info1.matrix(k+m*n*n,
false).column;
690 assemble(info1.matrix(k+m*n*n,
false).matrix, m,
691 info1.indices_by_block[row], info1.indices_by_block[column]);
692 assemble(info1.matrix(k+m*n*n,
true).matrix, m,
693 info1.indices_by_block[row], info2.indices_by_block[column]);
694 assemble(info2.matrix(k+m*n*n,
false).matrix, m,
695 info2.indices_by_block[row], info2.indices_by_block[column]);
696 assemble(info2.matrix(k+m*n*n,
true).matrix, m,
697 info2.indices_by_block[row], info1.indices_by_block[column]);
705 template <
typename MatrixType>
713 template <
typename MatrixType>
720 template <
typename MatrixType>
724 mg_constrained_dofs = &c;
728 template <
typename MatrixType>
738 template <
typename MatrixType>
744 interface_out = &out;
748 template <
typename MatrixType >
749 template <
class DOFINFO>
753 const unsigned int n = info.indices_by_block.size();
756 info.initialize_matrices(1, face);
759 info.initialize_matrices(n*n, face);
761 for (
unsigned int i=0; i<n; ++i)
762 for (
unsigned int j=0; j<n; ++j,++k)
764 info.matrix(k,
false).row = i;
765 info.matrix(k,
false).column = j;
768 info.matrix(k,
true).row = i;
769 info.matrix(k,
true).column = j;
776 template <
typename MatrixType>
781 const std::vector<types::global_dof_index> &i1,
782 const std::vector<types::global_dof_index> &i2)
789 for (
unsigned int j=0; j<i1.size(); ++j)
790 for (
unsigned int k=0; k<i2.size(); ++k)
791 if (std::fabs(M(j,k)) >= threshold)
792 G.add(i1[j], i2[k], M(j,k));
796 template <
typename MatrixType>
801 const std::vector<types::global_dof_index> &i1,
802 const std::vector<types::global_dof_index> &i2,
803 const unsigned int level)
808 if (mg_constrained_dofs ==
nullptr)
810 for (
unsigned int j=0; j<i1.size(); ++j)
811 for (
unsigned int k=0; k<i2.size(); ++k)
812 if (std::fabs(M(j,k)) >= threshold)
813 G.add(i1[j], i2[k], M(j,k));
817 for (
unsigned int j=0; j<i1.size(); ++j)
818 for (
unsigned int k=0; k<i2.size(); ++k)
822 if (std::fabs(M(j,k)) < threshold)
832 if (mg_constrained_dofs->at_refinement_edge(level, i1[j]) ||
833 mg_constrained_dofs->at_refinement_edge(level, i2[k]))
838 if ((mg_constrained_dofs->is_boundary_index(level, i1[j]) ||
839 mg_constrained_dofs->is_boundary_index(level, i2[k])) &&
843 G.add(i1[j], i2[k], M(j,k));
849 template <
typename MatrixType>
854 const std::vector<types::global_dof_index> &i1,
855 const std::vector<types::global_dof_index> &i2,
856 const unsigned int level)
861 if (mg_constrained_dofs ==
nullptr)
863 for (
unsigned int j=0; j<i1.size(); ++j)
864 for (
unsigned int k=0; k<i2.size(); ++k)
865 if (std::fabs(M(k,j)) >= threshold)
866 G.add(i1[j], i2[k], M(k,j));
870 for (
unsigned int j=0; j<i1.size(); ++j)
871 for (
unsigned int k=0; k<i2.size(); ++k)
872 if (std::fabs(M(k,j)) >= threshold)
873 if (!mg_constrained_dofs->at_refinement_edge(level, i2[k]))
874 G.add(i1[j], i2[k], M(k,j));
878 template <
typename MatrixType>
883 const std::vector<types::global_dof_index> &i1,
884 const std::vector<types::global_dof_index> &i2,
885 const unsigned int level)
890 if (mg_constrained_dofs ==
nullptr)
892 for (
unsigned int j=0; j<i1.size(); ++j)
893 for (
unsigned int k=0; k<i2.size(); ++k)
894 if (std::fabs(M(j,k)) >= threshold)
895 G.add(i1[j], i2[k], M(j,k));
899 for (
unsigned int j=0; j<i1.size(); ++j)
900 for (
unsigned int k=0; k<i2.size(); ++k)
901 if (std::fabs(M(j,k)) >= threshold)
902 if (!mg_constrained_dofs->at_refinement_edge(level, i2[k]))
903 G.add(i1[j], i2[k], M(j,k));
907 template <
typename MatrixType>
912 const std::vector<types::global_dof_index> &i1,
913 const std::vector<types::global_dof_index> &i2,
914 const unsigned int level)
920 for (
unsigned int j=0; j<i1.size(); ++j)
921 for (
unsigned int k=0; k<i2.size(); ++k)
922 if (std::fabs(M(j,k)) >= threshold)
934 if (mg_constrained_dofs->at_refinement_edge(level, i1[j]) &&
935 !mg_constrained_dofs->at_refinement_edge(level, i2[k]))
937 if ((!mg_constrained_dofs->is_boundary_index(level, i1[j]) &&
938 !mg_constrained_dofs->is_boundary_index(level, i2[k]))
940 (mg_constrained_dofs->is_boundary_index(level, i1[j]) &&
941 mg_constrained_dofs->is_boundary_index(level, i2[k]) &&
943 G.add(i1[j], i2[k], M(j,k));
948 template <
typename MatrixType>
953 const std::vector<types::global_dof_index> &i1,
954 const std::vector<types::global_dof_index> &i2,
955 const unsigned int level)
961 for (
unsigned int j=0; j<i1.size(); ++j)
962 for (
unsigned int k=0; k<i2.size(); ++k)
963 if (std::fabs(M(k,j)) >= threshold)
964 if (mg_constrained_dofs->at_refinement_edge(level, i1[j]) &&
965 !mg_constrained_dofs->at_refinement_edge(level, i2[k]))
967 if ((!mg_constrained_dofs->is_boundary_index(level, i1[j]) &&
968 !mg_constrained_dofs->is_boundary_index(level, i2[k]))
970 (mg_constrained_dofs->is_boundary_index(level, i1[j]) &&
971 mg_constrained_dofs->is_boundary_index(level, i2[k]) &&
973 G.add(i1[j], i2[k], M(k,j));
978 template <
typename MatrixType>
979 template <
class DOFINFO>
984 const unsigned int level = info.cell->level();
986 if (info.indices_by_block.size() == 0)
988 assemble((*matrix)[level], info.matrix(0,
false).matrix,
989 info.indices, info.indices, level);
990 if (mg_constrained_dofs !=
nullptr)
992 assemble_in((*interface_in)[level], info.matrix(0,
false).matrix,
993 info.indices, info.indices, level);
994 assemble_out((*interface_out)[level],info.matrix(0,
false).matrix,
995 info.indices, info.indices, level);
999 for (
unsigned int k=0; k<info.n_matrices(); ++k)
1001 const unsigned int row = info.matrix(k,
false).row;
1002 const unsigned int column = info.matrix(k,
false).column;
1004 assemble((*matrix)[level], info.matrix(k,
false).matrix,
1005 info.indices_by_block[row], info.indices_by_block[column], level);
1007 if (mg_constrained_dofs !=
nullptr)
1009 assemble_in((*interface_in)[level], info.matrix(k,
false).matrix,
1010 info.indices_by_block[row], info.indices_by_block[column], level);
1011 assemble_out((*interface_out)[level],info.matrix(k,
false).matrix,
1012 info.indices_by_block[column], info.indices_by_block[row], level);
1018 template <
typename MatrixType>
1019 template <
class DOFINFO>
1022 const DOFINFO &info2)
1026 const unsigned int level1 = info1.cell->level();
1027 const unsigned int level2 = info2.cell->level();
1029 if (info1.indices_by_block.size() == 0)
1031 if (level1 == level2)
1033 assemble((*matrix)[level1], info1.matrix(0,
false).matrix, info1.indices, info1.indices, level1);
1034 assemble((*matrix)[level1], info1.matrix(0,
true).matrix, info1.indices, info2.indices, level1);
1035 assemble((*matrix)[level1], info2.matrix(0,
false).matrix, info2.indices, info2.indices, level1);
1036 assemble((*matrix)[level1], info2.matrix(0,
true).matrix, info2.indices, info1.indices, level1);
1044 assemble((*matrix)[level1], info1.matrix(0,
false).matrix, info1.indices, info1.indices, level1);
1047 assemble_up((*flux_up)[level1],info1.matrix(0,
true).matrix, info2.indices, info1.indices, level1);
1048 assemble_down((*flux_down)[level1], info2.matrix(0,
true).matrix, info2.indices, info1.indices, level1);
1053 for (
unsigned int k=0; k<info1.n_matrices(); ++k)
1055 const unsigned int row = info1.matrix(k,
false).row;
1056 const unsigned int column = info1.matrix(k,
false).column;
1058 if (level1 == level2)
1060 assemble((*matrix)[level1], info1.matrix(k,
false).matrix, info1.indices_by_block[row], info1.indices_by_block[column], level1);
1061 assemble((*matrix)[level1], info1.matrix(k,
true).matrix, info1.indices_by_block[row], info2.indices_by_block[column], level1);
1062 assemble((*matrix)[level1], info2.matrix(k,
false).matrix, info2.indices_by_block[row], info2.indices_by_block[column], level1);
1063 assemble((*matrix)[level1], info2.matrix(k,
true).matrix, info2.indices_by_block[row], info1.indices_by_block[column], level1);
1071 assemble((*matrix)[level1], info1.matrix(k,
false).matrix, info1.indices_by_block[row], info1.indices_by_block[column], level1);
1074 assemble_up((*flux_up)[level1],info1.matrix(k,
true).matrix, info2.indices_by_block[column], info1.indices_by_block[row], level1);
1075 assemble_down((*flux_down)[level1], info2.matrix(k,
true).matrix, info2.indices_by_block[row], info1.indices_by_block[column], level1);
1083 template <
typename MatrixType,
typename VectorType>
1090 template <
typename MatrixType,
typename VectorType>
1095 VectorType *p = &rhs;
1096 data.
add(p,
"right hand side");
1102 template <
typename MatrixType,
typename VectorType>
1110 template <
typename MatrixType,
typename VectorType>
1111 template <
class DOFINFO>
1120 template <
typename MatrixType,
typename VectorType>
1124 const unsigned int index,
1125 const std::vector<types::global_dof_index> &indices)
1131 VectorType *v = residuals.entry<VectorType *>(index);
1135 for (
unsigned int i=0; i<indices.size(); ++i)
1136 (*v)(indices[i]) += vector(i);
1138 for (
unsigned int j=0; j<indices.size(); ++j)
1139 for (
unsigned int k=0; k<indices.size(); ++k)
1149 template <
typename MatrixType,
typename VectorType>
1153 const unsigned int index,
1154 const std::vector<types::global_dof_index> &i1,
1155 const std::vector<types::global_dof_index> &i2)
1161 VectorType *v = residuals.entry<VectorType *>(index);
1165 for (
unsigned int j=0; j<i1.size(); ++j)
1166 for (
unsigned int k=0; k<i2.size(); ++k)
1178 template <
typename MatrixType,
typename VectorType>
1179 template <
class DOFINFO>
1185 const unsigned int n = info.indices_by_block.size();
1189 for (
unsigned int m=0; m<MatrixSimple<MatrixType>::matrix.size(); ++m)
1190 assemble(info.matrix(m,
false).matrix,info.vector(m).block(0), m, info.indices);
1194 for (
unsigned int m=0; m<MatrixSimple<MatrixType>::matrix.size(); ++m)
1195 for (
unsigned int k=0; k<n*n; ++k)
1197 const unsigned int row = info.matrix(k+m*n*n,
false).row;
1198 const unsigned int column = info.matrix(k+m*n*n,
false).column;
1201 assemble(info.matrix(k+m*n*n,
false).matrix,
1202 info.vector(m).block(row), m,
1203 info.indices_by_block[row]);
1205 assemble(info.matrix(k+m*n*n,
false).matrix,
1206 info.vector(m).block(row), m,
1207 info.indices_by_block[row],
1208 info.indices_by_block[column]);
1215 template <
typename MatrixType,
typename VectorType>
1216 template <
class DOFINFO>
1219 const DOFINFO &info2)
1223 AssertDimension(info1.indices_by_block.size(),info2.indices_by_block.size());
1225 const unsigned int n = info1.indices_by_block.size();
1229 for (
unsigned int m=0; m<MatrixSimple<MatrixType>::matrix.size(); ++m)
1231 assemble(info1.matrix(m,
false).matrix, info1.vector(m).block(0), m, info1.indices);
1232 assemble(info1.matrix(m,
true).matrix, info1.vector(m).block(0), m, info1.indices, info2.indices);
1233 assemble(info2.matrix(m,
false).matrix, info2.vector(m).block(0), m, info2.indices);
1234 assemble(info2.matrix(m,
true).matrix, info2.vector(m).block(0), m, info2.indices, info1.indices);
1239 for (
unsigned int m=0; m<MatrixSimple<MatrixType>::matrix.size(); ++m)
1240 for (
unsigned int k=0; k<n*n; ++k)
1242 const unsigned int row = info1.matrix(k+m*n*n,
false).row;
1243 const unsigned int column = info1.matrix(k+m*n*n,
false).column;
1247 assemble(info1.matrix(k+m*n*n,
false).matrix, info1.vector(m).block(row), m,
1248 info1.indices_by_block[row]);
1249 assemble(info2.matrix(k+m*n*n,
false).matrix, info2.vector(m).block(row), m,
1250 info2.indices_by_block[row]);
1254 assemble(info1.matrix(k+m*n*n,
false).matrix, info1.vector(m).block(row), m,
1255 info1.indices_by_block[row], info1.indices_by_block[column]);
1256 assemble(info2.matrix(k+m*n*n,
false).matrix, info2.vector(m).block(row), m,
1257 info2.indices_by_block[row], info2.indices_by_block[column]);
1259 assemble(info1.matrix(k+m*n*n,
true).matrix, info1.vector(m).block(row), m,
1260 info1.indices_by_block[row], info2.indices_by_block[column]);
1261 assemble(info2.matrix(k+m*n*n,
true).matrix, info2.vector(m).block(row), m,
1262 info2.indices_by_block[row], info1.indices_by_block[column]);
1269 DEAL_II_NAMESPACE_CLOSE
static const unsigned int invalid_unsigned_int
void initialize_interfaces(MGLevelObject< MatrixType > &interface_in, MGLevelObject< MatrixType > &interface_out)
SystemSimple(double threshold=1.e-12)
#define AssertDimension(dim1, dim2)
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< const MGConstrainedDoFs, MGMatrixSimple< MatrixType > > mg_constrained_dofs
void initialize(MGLevelObject< MatrixType > &m)
static ::ExceptionBase & ExcNotInitialized()
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)
void initialize_info(DOFINFO &info, bool face) const
void initialize(MatrixType &m, VectorType &rhs)
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)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > interface_in
void assemble(const DOFINFO &info)
void initialize_fluxes(MGLevelObject< MatrixType > &flux_up, MGLevelObject< MatrixType > &flux_down)
std::vector< SmartPointer< MatrixType, MatrixSimple< MatrixType > > > matrix
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > flux_up
void initialize(AnyData &results)
MGMatrixSimple(double threshold=1.e-12)
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > matrix
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > interface_out
void initialize_local_blocks(const BlockIndices &)
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)
SmartPointer< const ConstraintMatrix, ResidualSimple< VectorType > > constraints
void add(type entry, const std::string &name)
Add a new data object.
void assemble(const DOFINFO &info)
void initialize_info(DOFINFO &info, bool face) const
void initialize_info(DOFINFO &info, bool face) const
void assemble(const DOFINFO &info)
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > flux_down
void initialize(MatrixType &m)
MatrixSimple(double threshold=1.e-12)
SmartPointer< const ConstraintMatrix, MatrixSimple< MatrixType > > constraints
void assemble(const DOFINFO &info)
void initialize_info(DOFINFO &info, bool face) const
static ::ExceptionBase & ExcInternalError()