31class MultigridOperator :
public MultigridOperatorBase<dim, Number>
34 typedef MultigridOperatorBase<dim, Number> Base;
35 typedef typename Base::value_type value_type;
36 typedef typename Base::VectorType VectorType;
38 MultigridOperator(std::shared_ptr<Operator> op) : pde_operator(op)
42 virtual ~MultigridOperator()
46 std::shared_ptr<Operator>
47 get_pde_operator()
const
49 AssertThrow(pde_operator.get() != 0, dealii::ExcMessage(
"Invalid pointer"));
54 dealii::AffineConstraints<typename Operator::value_type>
const &
55 get_affine_constraints()
const final
57 return pde_operator->get_affine_constraints();
60 dealii::MatrixFree<dim, Number>
const &
61 get_matrix_free()
const final
63 return pde_operator->get_matrix_free();
67 get_dof_index()
const final
69 return pde_operator->get_dof_index();
72 dealii::types::global_dof_index
75 return pde_operator->m();
78 dealii::types::global_dof_index
81 return pde_operator->n();
85 el(
unsigned int const i,
unsigned int const j)
const final
87 return pde_operator->el(i, j);
91 initialize_dof_vector(VectorType & vector)
const final
93 pde_operator->initialize_dof_vector(vector);
97 vmult(VectorType & dst, VectorType
const & src)
const final
99 pde_operator->vmult(dst, src);
103 vmult_add(VectorType & dst, VectorType
const & src)
const final
105 pde_operator->vmult_add(dst, src);
109 vmult_interface_down(VectorType & dst, VectorType
const & src)
const final
111 pde_operator->vmult_interface_down(dst, src);
115 vmult_add_interface_up(VectorType & dst, VectorType
const & src)
const final
117 pde_operator->vmult_add_interface_up(dst, src);
121 calculate_inverse_diagonal(VectorType & inverse_diagonal_entries)
const final
123 pde_operator->calculate_inverse_diagonal(inverse_diagonal_entries);
127 initialize_block_diagonal_preconditioner(
bool const initialize)
const final
129 pde_operator->initialize_block_diagonal_preconditioner(initialize);
133 update_block_diagonal_preconditioner()
const final
135 pde_operator->update_block_diagonal_preconditioner();
139 apply_inverse_block_diagonal(VectorType & dst, VectorType
const & src)
const final
141 pde_operator->apply_inverse_block_diagonal(dst, src);
145 apply_inverse_additive_schwarz_matrices(VectorType & dst, VectorType
const & src)
const final
147 pde_operator->apply_inverse_additive_schwarz_matrices(dst, src);
151 compute_factorized_additive_schwarz_matrices()
const final
153 pde_operator->compute_factorized_additive_schwarz_matrices();
156#ifdef DEAL_II_WITH_TRILINOS
158 init_system_matrix(dealii::TrilinosWrappers::SparseMatrix & system_matrix,
159 MPI_Comm
const & mpi_comm)
const final
161 pde_operator->init_system_matrix(system_matrix, mpi_comm);
165 calculate_system_matrix(dealii::TrilinosWrappers::SparseMatrix & system_matrix)
const final
167 pde_operator->calculate_system_matrix(system_matrix);
171#ifdef DEAL_II_WITH_PETSC
173 init_system_matrix(dealii::PETScWrappers::MPI::SparseMatrix & system_matrix,
174 MPI_Comm
const & mpi_comm)
const final
176 pde_operator->init_system_matrix(system_matrix, mpi_comm);
180 calculate_system_matrix(dealii::PETScWrappers::MPI::SparseMatrix & system_matrix)
const final
182 pde_operator->calculate_system_matrix(system_matrix);
187 get_constant_modes(std::vector<std::vector<bool>> & constant_modes,
188 std::vector<std::vector<double>> & constant_modes_values)
const final
190 pde_operator->get_constant_modes(constant_modes, constant_modes_values);
194 std::shared_ptr<Operator> pde_operator;