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);
153 template <
typename MatrixType>
189 template <
class DOFINFO>
196 template<
class DOFINFO>
203 template<
class DOFINFO>
205 const DOFINFO &info2);
212 const unsigned int index,
213 const std::vector<types::global_dof_index> &i1,
214 const std::vector<types::global_dof_index> &i2);
219 std::vector<SmartPointer<MatrixType,MatrixSimple<MatrixType> > >
matrix;
245 template <
typename MatrixType>
286 template <
class DOFINFO>
292 template<
class DOFINFO>
299 template<
class DOFINFO>
301 const DOFINFO &info2);
308 const std::vector<types::global_dof_index> &i1,
309 const std::vector<types::global_dof_index> &i2);
316 const std::vector<types::global_dof_index> &i1,
317 const std::vector<types::global_dof_index> &i2,
318 const unsigned int level);
324 void assemble_up(MatrixType &G,
326 const std::vector<types::global_dof_index> &i1,
327 const std::vector<types::global_dof_index> &i2,
333 void assemble_down(MatrixType &G,
335 const std::vector<types::global_dof_index> &i1,
336 const std::vector<types::global_dof_index> &i2,
343 void assemble_in(MatrixType &G,
345 const std::vector<types::global_dof_index> &i1,
346 const std::vector<types::global_dof_index> &i2,
353 void assemble_out(MatrixType &G,
355 const std::vector<types::global_dof_index> &i1,
356 const std::vector<types::global_dof_index> &i2,
412 template <
typename MatrixType,
typename VectorType>
426 void initialize(MatrixType &m, VectorType &rhs);
444 template <
class DOFINFO>
450 template<
class DOFINFO>
457 template<
class DOFINFO>
459 const DOFINFO &info2);
465 template <
typename VectorType>
472 template <
typename VectorType>
480 template <
typename MatrixType>
486 template <
typename VectorType>
487 template <
class DOFINFO>
495 template <
typename VectorType>
496 template <
class DOFINFO>
505 for (
unsigned int i=0; i<info.vector(k).block(0).size(); ++i)
506 (*v)(info.indices[i]) += info.vector(k).block(0)(i);
510 if (info.indices_by_block.size() == 0)
511 constraints->distribute_local_to_global(info.vector(k).block(0), info.indices, *v);
513 for (
unsigned int i=0; i != info.vector(k).n_blocks(); ++i)
514 constraints->distribute_local_to_global(info.vector(k).block(i), info.indices_by_block[i], *v);
519 template <
typename VectorType>
520 template <
class DOFINFO>
523 const DOFINFO &info2)
530 for (
unsigned int i=0; i<info1.vector(k).block(0).size(); ++i)
531 (*v)(info1.indices[i]) += info1.vector(k).block(0)(i);
532 for (
unsigned int i=0; i<info2.vector(k).block(0).size(); ++i)
533 (*v)(info2.indices[i]) += info2.vector(k).block(0)(i);
537 if (info1.indices_by_block.size() == 0 && info2.indices_by_block.size() == 0)
540 (info1.vector(k).block(0), info1.indices, *v);
542 (info2.vector(k).block(0), info2.indices, *v);
544 else if (info1.indices_by_block.size() != 0 && info2.indices_by_block.size() != 0)
546 for (
unsigned int i=0; i<info1.vector(k).n_blocks(); ++i)
549 (info1.vector(k).block(i), info1.indices_by_block[i], *v);
551 (info2.vector(k).block(i), info2.indices_by_block[i], *v);
556 Assert(
false, ExcNotImplemented());
565 template <
typename MatrixType>
573 template <
typename MatrixType>
582 template <
typename MatrixType>
587 for (
unsigned int i=0; i<m.size(); ++i)
592 template <
typename MatrixType>
600 template <
typename MatrixType >
601 template <
class DOFINFO>
607 const unsigned int n = info.indices_by_block.size();
610 info.initialize_matrices(
matrix.size(), face);
613 info.initialize_matrices(
matrix.size()*n*n, face);
615 for (
unsigned int m=0; m<
matrix.size(); ++m)
616 for (
unsigned int i=0; i<n; ++i)
617 for (
unsigned int j=0; j<n; ++j,++k)
619 info.matrix(k,
false).row = i;
620 info.matrix(k,
false).column = j;
623 info.matrix(k,
true).row = i;
624 info.matrix(k,
true).column = j;
632 template <
typename MatrixType>
635 const unsigned int index,
636 const std::vector<types::global_dof_index> &i1,
637 const std::vector<types::global_dof_index> &i2)
644 for (
unsigned int j=0; j<i1.size(); ++j)
645 for (
unsigned int k=0; k<i2.size(); ++k)
647 matrix[index]->add(i1[j], i2[k], M(j,k));
654 template <
typename MatrixType>
655 template <
class DOFINFO>
660 const unsigned int n = info.indices_by_block.size();
663 for (
unsigned int m=0; m<
matrix.size(); ++m)
664 assemble(info.matrix(m,
false).matrix, m, info.indices, info.indices);
667 for (
unsigned int m=0; m<
matrix.size(); ++m)
668 for (
unsigned int k=0; k<n*n; ++k)
670 assemble(info.matrix(k+m*n*n,
false).matrix, m,
671 info.indices_by_block[info.matrix(k+m*n*n,
false).row],
672 info.indices_by_block[info.matrix(k+m*n*n,
false).column]);
678 template <
typename MatrixType>
679 template <
class DOFINFO>
685 AssertDimension(info1.indices_by_block.size(),info2.indices_by_block.size());
687 const unsigned int n = info1.indices_by_block.size();
691 for (
unsigned int m=0; m<
matrix.size(); ++m)
693 assemble(info1.matrix(m,
false).matrix, m, info1.indices, info1.indices);
694 assemble(info1.matrix(m,
true).matrix, m, info1.indices, info2.indices);
695 assemble(info2.matrix(m,
false).matrix, m, info2.indices, info2.indices);
696 assemble(info2.matrix(m,
true).matrix, m, info2.indices, info1.indices);
701 for (
unsigned int m=0; m<
matrix.size(); ++m)
702 for (
unsigned int k=0; k<n*n; ++k)
704 const unsigned int row = info1.matrix(k+m*n*n,
false).row;
705 const unsigned int column = info1.matrix(k+m*n*n,
false).column;
707 assemble(info1.matrix(k+m*n*n,
false).matrix, m,
708 info1.indices_by_block[row], info1.indices_by_block[column]);
709 assemble(info1.matrix(k+m*n*n,
true).matrix, m,
710 info1.indices_by_block[row], info2.indices_by_block[column]);
711 assemble(info2.matrix(k+m*n*n,
false).matrix, m,
712 info2.indices_by_block[row], info2.indices_by_block[column]);
713 assemble(info2.matrix(k+m*n*n,
true).matrix, m,
714 info2.indices_by_block[row], info1.indices_by_block[column]);
722 template <
typename MatrixType>
730 template <
typename MatrixType>
737 template <
typename MatrixType>
745 template <
typename MatrixType>
755 template <
typename MatrixType>
765 template <
typename MatrixType >
766 template <
class DOFINFO>
770 const unsigned int n = info.indices_by_block.size();
773 info.initialize_matrices(1, face);
776 info.initialize_matrices(n*n, face);
778 for (
unsigned int i=0; i<n; ++i)
779 for (
unsigned int j=0; j<n; ++j,++k)
781 info.matrix(k,
false).row = i;
782 info.matrix(k,
false).column = j;
785 info.matrix(k,
true).row = i;
786 info.matrix(k,
true).column = j;
793 template <
typename MatrixType>
798 const std::vector<types::global_dof_index> &i1,
799 const std::vector<types::global_dof_index> &i2)
806 for (
unsigned int j=0; j<i1.size(); ++j)
807 for (
unsigned int k=0; k<i2.size(); ++k)
809 G.add(i1[j], i2[k], M(j,k));
813 template <
typename MatrixType>
818 const std::vector<types::global_dof_index> &i1,
819 const std::vector<types::global_dof_index> &i2,
820 const unsigned int level)
827 for (
unsigned int j=0; j<i1.size(); ++j)
828 for (
unsigned int k=0; k<i2.size(); ++k)
830 G.add(i1[j], i2[k], M(j,k));
834 for (
unsigned int j=0; j<i1.size(); ++j)
835 for (
unsigned int k=0; k<i2.size(); ++k)
860 G.add(i1[j], i2[k], M(j,k));
866 template <
typename MatrixType>
871 const std::vector<types::global_dof_index> &i1,
872 const std::vector<types::global_dof_index> &i2,
873 const unsigned int level)
880 for (
unsigned int j=0; j<i1.size(); ++j)
881 for (
unsigned int k=0; k<i2.size(); ++k)
883 G.add(i1[j], i2[k], M(k,j));
887 for (
unsigned int j=0; j<i1.size(); ++j)
888 for (
unsigned int k=0; k<i2.size(); ++k)
891 G.add(i1[j], i2[k], M(k,j));
895 template <
typename MatrixType>
900 const std::vector<types::global_dof_index> &i1,
901 const std::vector<types::global_dof_index> &i2,
902 const unsigned int level)
909 for (
unsigned int j=0; j<i1.size(); ++j)
910 for (
unsigned int k=0; k<i2.size(); ++k)
912 G.add(i1[j], i2[k], M(j,k));
916 for (
unsigned int j=0; j<i1.size(); ++j)
917 for (
unsigned int k=0; k<i2.size(); ++k)
920 G.add(i1[j], i2[k], M(j,k));
924 template <
typename MatrixType>
929 const std::vector<types::global_dof_index> &i1,
930 const std::vector<types::global_dof_index> &i2,
931 const unsigned int level)
937 for (
unsigned int j=0; j<i1.size(); ++j)
938 for (
unsigned int k=0; k<i2.size(); ++k)
960 G.add(i1[j], i2[k], M(j,k));
965 template <
typename MatrixType>
970 const std::vector<types::global_dof_index> &i1,
971 const std::vector<types::global_dof_index> &i2,
972 const unsigned int level)
978 for (
unsigned int j=0; j<i1.size(); ++j)
979 for (
unsigned int k=0; k<i2.size(); ++k)
990 G.add(i1[j], i2[k], M(k,j));
995 template <
typename MatrixType>
996 template <
class DOFINFO>
1001 const unsigned int level = info.cell->level();
1003 if (info.indices_by_block.size() == 0)
1006 info.indices, info.indices, level);
1010 info.indices, info.indices, level);
1012 info.indices, info.indices, level);
1016 for (
unsigned int k=0; k<info.n_matrices(); ++k)
1018 const unsigned int row = info.matrix(k,
false).row;
1019 const unsigned int column = info.matrix(k,
false).column;
1022 info.indices_by_block[row], info.indices_by_block[column], level);
1027 info.indices_by_block[row], info.indices_by_block[column], level);
1029 info.indices_by_block[column], info.indices_by_block[row], level);
1035 template <
typename MatrixType>
1036 template <
class DOFINFO>
1039 const DOFINFO &info2)
1043 const unsigned int level1 = info1.cell->level();
1044 const unsigned int level2 = info2.cell->level();
1046 if (info1.indices_by_block.size() == 0)
1048 if (level1 == level2)
1050 assemble((*
matrix)[level1], info1.matrix(0,
false).matrix, info1.indices, info1.indices, level1);
1051 assemble((*
matrix)[level1], info1.matrix(0,
true).matrix, info1.indices, info2.indices, level1);
1052 assemble((*
matrix)[level1], info2.matrix(0,
false).matrix, info2.indices, info2.indices, level1);
1053 assemble((*
matrix)[level1], info2.matrix(0,
true).matrix, info2.indices, info1.indices, level1);
1057 Assert(level1 > level2, ExcInternalError());
1061 assemble((*
matrix)[level1], info1.matrix(0,
false).matrix, info1.indices, info1.indices, level1);
1064 assemble_up((*
flux_up)[level1],info1.matrix(0,
true).matrix, info2.indices, info1.indices, level1);
1070 for (
unsigned int k=0; k<info1.n_matrices(); ++k)
1072 const unsigned int row = info1.matrix(k,
false).row;
1073 const unsigned int column = info1.matrix(k,
false).column;
1075 if (level1 == level2)
1077 assemble((*
matrix)[level1], info1.matrix(k,
false).matrix, info1.indices_by_block[row], info1.indices_by_block[column], level1);
1078 assemble((*
matrix)[level1], info1.matrix(k,
true).matrix, info1.indices_by_block[row], info2.indices_by_block[column], level1);
1079 assemble((*
matrix)[level1], info2.matrix(k,
false).matrix, info2.indices_by_block[row], info2.indices_by_block[column], level1);
1080 assemble((*
matrix)[level1], info2.matrix(k,
true).matrix, info2.indices_by_block[row], info1.indices_by_block[column], level1);
1084 Assert(level1 > level2, ExcInternalError());
1088 assemble((*
matrix)[level1], info1.matrix(k,
false).matrix, info1.indices_by_block[row], info1.indices_by_block[column], level1);
1091 assemble_up((*
flux_up)[level1],info1.matrix(k,
true).matrix, info2.indices_by_block[row], info1.indices_by_block[column], level1);
1092 assemble_down((*
flux_down)[level1], info2.matrix(k,
true).matrix, info2.indices_by_block[row], info1.indices_by_block[column], level1);
1100 template <
typename MatrixType,
typename VectorType>
1107 template <
typename MatrixType,
typename VectorType>
1112 VectorType *p = &rhs;
1113 data.
add(p,
"right hand side");
1119 template <
typename MatrixType,
typename VectorType>
1128 template <
typename MatrixType,
typename VectorType>
1129 template <
class DOFINFO>
1139 template <
typename MatrixType,
typename VectorType>
1140 template <
class DOFINFO>
1149 template <
typename MatrixType,
typename VectorType>
1150 template <
class DOFINFO>
1153 const DOFINFO &info2)
1161 DEAL_II_NAMESPACE_CLOSE
static const unsigned int invalid_unsigned_int
type entry(const std::string &name)
Access to stored data object by name.
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)
::ExceptionBase & ExcMessage(std::string arg1)
SmartPointer< const MGConstrainedDoFs, MGMatrixSimple< MatrixType > > mg_constrained_dofs
Auxiliary class aiding in the handling of block structures like in BlockVector or FESystem...
void initialize(MGLevelObject< MatrixType > &m)
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)
#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.
unsigned int size() const
Number of stored data objects.
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