ExaDG
Loading...
Searching...
No Matches
multigrid_operator.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef INCLUDE_EXADG_OPERATORS_MULTIGRID_OPERATOR_H_
23#define INCLUDE_EXADG_OPERATORS_MULTIGRID_OPERATOR_H_
24
25#include <exadg/operators/multigrid_operator_base.h>
26
27namespace ExaDG
28{
29template<int dim, typename Number, typename Operator>
30class MultigridOperator : public MultigridOperatorBase<dim, Number>
31{
32public:
34 typedef typename Base::value_type value_type;
35 typedef typename Base::VectorType VectorType;
36
37 MultigridOperator(std::shared_ptr<Operator> op) : pde_operator(op)
38 {
39 }
40
41 virtual ~MultigridOperator()
42 {
43 }
44
45 std::shared_ptr<Operator>
46 get_pde_operator() const
47 {
48 AssertThrow(pde_operator.get() != 0, dealii::ExcMessage("Invalid pointer"));
49
50 return pde_operator;
51 }
52
53 dealii::AffineConstraints<typename Operator::value_type> const &
54 get_affine_constraints() const final
55 {
56 return pde_operator->get_affine_constraints();
57 }
58
59 dealii::MatrixFree<dim, Number> const &
60 get_matrix_free() const final
61 {
62 return pde_operator->get_matrix_free();
63 }
64
65 unsigned int
66 get_dof_index() const final
67 {
68 return pde_operator->get_dof_index();
69 }
70
71 dealii::types::global_dof_index
72 m() const final
73 {
74 return pde_operator->m();
75 }
76
77 dealii::types::global_dof_index
78 n() const final
79 {
80 return pde_operator->n();
81 }
82
83 Number
84 el(unsigned int const i, unsigned int const j) const final
85 {
86 return pde_operator->el(i, j);
87 }
88
89 void
90 initialize_dof_vector(VectorType & vector) const final
91 {
92 pde_operator->initialize_dof_vector(vector);
93 }
94
95 void
96 vmult(VectorType & dst, VectorType const & src) const final
97 {
98 pde_operator->vmult(dst, src);
99 }
100
101 void
102 vmult_add(VectorType & dst, VectorType const & src) const final
103 {
104 pde_operator->vmult_add(dst, src);
105 }
106
107 void
108 vmult_interface_down(VectorType & dst, VectorType const & src) const final
109 {
110 pde_operator->vmult_interface_down(dst, src);
111 }
112
113 void
114 vmult_add_interface_up(VectorType & dst, VectorType const & src) const final
115 {
116 pde_operator->vmult_add_interface_up(dst, src);
117 }
118
119 void
120 calculate_inverse_diagonal(VectorType & inverse_diagonal_entries) const final
121 {
122 pde_operator->calculate_inverse_diagonal(inverse_diagonal_entries);
123 }
124
125 void
126 initialize_block_diagonal_preconditioner(bool const initialize) const final
127 {
128 pde_operator->initialize_block_diagonal_preconditioner(initialize);
129 }
130
131 void
132 update_block_diagonal_preconditioner() const final
133 {
134 pde_operator->update_block_diagonal_preconditioner();
135 }
136
137 void
138 apply_inverse_block_diagonal(VectorType & dst, VectorType const & src) const final
139 {
140 pde_operator->apply_inverse_block_diagonal(dst, src);
141 }
142
143 virtual void
144 apply_inverse_additive_schwarz_matrices(VectorType & dst, VectorType const & src) const final
145 {
146 pde_operator->apply_inverse_additive_schwarz_matrices(dst, src);
147 }
148
149 virtual void
150 compute_factorized_additive_schwarz_matrices() const final
151 {
152 pde_operator->compute_factorized_additive_schwarz_matrices();
153 }
154
155#ifdef DEAL_II_WITH_TRILINOS
156 void
157 init_system_matrix(dealii::TrilinosWrappers::SparseMatrix & system_matrix,
158 MPI_Comm const & mpi_comm) const final
159 {
160 pde_operator->init_system_matrix(system_matrix, mpi_comm);
161 }
162
163 void
164 calculate_system_matrix(dealii::TrilinosWrappers::SparseMatrix & system_matrix) const final
165 {
166 pde_operator->calculate_system_matrix(system_matrix);
167 }
168#endif
169
170#ifdef DEAL_II_WITH_PETSC
171 void
172 init_system_matrix(dealii::PETScWrappers::MPI::SparseMatrix & system_matrix,
173 MPI_Comm const & mpi_comm) const final
174 {
175 pde_operator->init_system_matrix(system_matrix, mpi_comm);
176 }
177
178 void
179 calculate_system_matrix(dealii::PETScWrappers::MPI::SparseMatrix & system_matrix) const final
180 {
181 pde_operator->calculate_system_matrix(system_matrix);
182 }
183#endif
184
185private:
186 std::shared_ptr<Operator> pde_operator;
187};
188
189} // namespace ExaDG
190
191#endif /* INCLUDE_EXADG_OPERATORS_MULTIGRID_OPERATOR_H_ */
Definition multigrid_operator_base.h:35
Definition multigrid_operator.h:31
Definition driver.cpp:33