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 EXADG_OPERATORS_MULTIGRID_OPERATOR_H_
23#define EXADG_OPERATORS_MULTIGRID_OPERATOR_H_
24
25// ExaDG
26#include <exadg/operators/multigrid_operator_base.h>
27
28namespace ExaDG
29{
30template<int dim, typename Number, typename Operator>
31class MultigridOperator : public MultigridOperatorBase<dim, Number>
32{
33public:
34 typedef MultigridOperatorBase<dim, Number> Base;
35 typedef typename Base::value_type value_type;
36 typedef typename Base::VectorType VectorType;
37
38 MultigridOperator(std::shared_ptr<Operator> op) : pde_operator(op)
39 {
40 }
41
42 virtual ~MultigridOperator()
43 {
44 }
45
46 std::shared_ptr<Operator>
47 get_pde_operator() const
48 {
49 AssertThrow(pde_operator.get() != 0, dealii::ExcMessage("Invalid pointer"));
50
51 return pde_operator;
52 }
53
54 dealii::AffineConstraints<typename Operator::value_type> const &
55 get_affine_constraints() const final
56 {
57 return pde_operator->get_affine_constraints();
58 }
59
60 dealii::MatrixFree<dim, Number> const &
61 get_matrix_free() const final
62 {
63 return pde_operator->get_matrix_free();
64 }
65
66 unsigned int
67 get_dof_index() const final
68 {
69 return pde_operator->get_dof_index();
70 }
71
72 dealii::types::global_dof_index
73 m() const final
74 {
75 return pde_operator->m();
76 }
77
78 dealii::types::global_dof_index
79 n() const final
80 {
81 return pde_operator->n();
82 }
83
84 Number
85 el(unsigned int const i, unsigned int const j) const final
86 {
87 return pde_operator->el(i, j);
88 }
89
90 void
91 initialize_dof_vector(VectorType & vector) const final
92 {
93 pde_operator->initialize_dof_vector(vector);
94 }
95
96 void
97 vmult(VectorType & dst, VectorType const & src) const final
98 {
99 pde_operator->vmult(dst, src);
100 }
101
102 void
103 vmult_add(VectorType & dst, VectorType const & src) const final
104 {
105 pde_operator->vmult_add(dst, src);
106 }
107
108 void
109 vmult_interface_down(VectorType & dst, VectorType const & src) const final
110 {
111 pde_operator->vmult_interface_down(dst, src);
112 }
113
114 void
115 vmult_add_interface_up(VectorType & dst, VectorType const & src) const final
116 {
117 pde_operator->vmult_add_interface_up(dst, src);
118 }
119
120 void
121 calculate_inverse_diagonal(VectorType & inverse_diagonal_entries) const final
122 {
123 pde_operator->calculate_inverse_diagonal(inverse_diagonal_entries);
124 }
125
126 void
127 initialize_block_diagonal_preconditioner(bool const initialize) const final
128 {
129 pde_operator->initialize_block_diagonal_preconditioner(initialize);
130 }
131
132 void
133 update_block_diagonal_preconditioner() const final
134 {
135 pde_operator->update_block_diagonal_preconditioner();
136 }
137
138 void
139 apply_inverse_block_diagonal(VectorType & dst, VectorType const & src) const final
140 {
141 pde_operator->apply_inverse_block_diagonal(dst, src);
142 }
143
144 virtual void
145 apply_inverse_additive_schwarz_matrices(VectorType & dst, VectorType const & src) const final
146 {
147 pde_operator->apply_inverse_additive_schwarz_matrices(dst, src);
148 }
149
150 virtual void
151 compute_factorized_additive_schwarz_matrices() const final
152 {
153 pde_operator->compute_factorized_additive_schwarz_matrices();
154 }
155
156#ifdef DEAL_II_WITH_TRILINOS
157 void
158 init_system_matrix(dealii::TrilinosWrappers::SparseMatrix & system_matrix,
159 MPI_Comm const & mpi_comm) const final
160 {
161 pde_operator->init_system_matrix(system_matrix, mpi_comm);
162 }
163
164 void
165 calculate_system_matrix(dealii::TrilinosWrappers::SparseMatrix & system_matrix) const final
166 {
167 pde_operator->calculate_system_matrix(system_matrix);
168 }
169#endif
170
171#ifdef DEAL_II_WITH_PETSC
172 void
173 init_system_matrix(dealii::PETScWrappers::MPI::SparseMatrix & system_matrix,
174 MPI_Comm const & mpi_comm) const final
175 {
176 pde_operator->init_system_matrix(system_matrix, mpi_comm);
177 }
178
179 void
180 calculate_system_matrix(dealii::PETScWrappers::MPI::SparseMatrix & system_matrix) const final
181 {
182 pde_operator->calculate_system_matrix(system_matrix);
183 }
184#endif
185
186 void
187 get_constant_modes(std::vector<std::vector<bool>> & constant_modes,
188 std::vector<std::vector<double>> & constant_modes_values) const final
189 {
190 pde_operator->get_constant_modes(constant_modes, constant_modes_values);
191 }
192
193private:
194 std::shared_ptr<Operator> pde_operator;
195};
196
197} // namespace ExaDG
198
199#endif /* EXADG_OPERATORS_MULTIGRID_OPERATOR_H_ */
Definition driver.cpp:33