ASPECT
stokes_matrix_free_global_coarsening.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2011 - 2025 by the authors of the ASPECT code.
3 
4  This file is part of ASPECT.
5 
6  ASPECT is free software; you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation; either version 2, or (at your option)
9  any later version.
10 
11  ASPECT is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with ASPECT; see the file LICENSE. If not see
18  <http://www.gnu.org/licenses/>.
19 */
20 
21 #ifndef _aspect_simulator_solver_stokes_matrix_free_global_coarsening_h
22 #define _aspect_simulator_solver_stokes_matrix_free_global_coarsening_h
23 
24 #include <aspect/global.h>
27 
28 #include <deal.II/matrix_free/matrix_free.h>
29 #include <deal.II/matrix_free/operators.h>
30 #include <deal.II/matrix_free/fe_evaluation.h>
31 
32 #include <deal.II/multigrid/mg_constrained_dofs.h>
33 #include <deal.II/multigrid/multigrid.h>
34 #include <deal.II/multigrid/mg_transfer_matrix_free.h>
35 #include <deal.II/multigrid/mg_transfer_global_coarsening.templates.h>
36 #include <deal.II/multigrid/mg_tools.h>
37 #include <deal.II/multigrid/mg_coarse.h>
38 #include <deal.II/multigrid/mg_smoother.h>
39 #include <deal.II/multigrid/mg_matrix.h>
40 
41 namespace aspect
42 {
53  template <int dim, int velocity_degree>
55  {
56  public:
63  const Parameters<dim> &parameters);
64 
69 
73  void initialize() override;
74 
78  std::string name() const override;
79 
103  solve(const LinearAlgebra::BlockSparseMatrix &system_matrix,
104  const LinearAlgebra::BlockVector &system_rhs,
105  const bool solve_newton_system,
106  const double last_pressure_normalization_adjustment,
107  LinearAlgebra::BlockVector &solution_vector) override;
108 
113  void setup_dofs() override;
114 
122  void assemble() override;
123 
128  void build_preconditioner() override;
129 
133  static
134  void declare_parameters (ParameterHandler &prm);
135 
139  void parse_parameters (ParameterHandler &prm) override;
140 
144  std::size_t get_dof_handler_memory_consumption() const override;
145 
149  std::size_t get_mg_transfer_memory_consumption() const override;
150 
154  std::size_t get_constraint_memory_consumption() const override;
155 
160  std::size_t get_cell_data_memory_consumption() const override;
161 
162  private:
169 
174  void correct_stokes_rhs();
175 
176 
178 
180 
186 
192 
193  FESystem<dim> fe_v;
194  FESystem<dim> fe_p;
195  FESystem<dim> fe_projection;
196 
201 
205  MGLevelObject<MatrixFreeStokesOperators::OperatorCellData<dim, GMGNumberType>> level_cell_data;
206 
213 
218 
219  MGLevelObject<GMGABlockMatrixType> mg_matrices_A_block;
220  MGLevelObject<GMGSchurComplementMatrixType> mg_matrices_Schur_complement;
221 
222  MGConstrainedDoFs mg_constrained_dofs_A_block;
225 
226  unsigned int min_level;
227  unsigned int max_level;
228 
229  std::vector<std::shared_ptr<MatrixFree<dim,double>>> matrix_free_objects;
230 
231  std::vector<std::shared_ptr<const Triangulation<dim, dim>>> trias;
232 
233  MGLevelObject<DoFHandler<dim>> dofhandlers_v;
234  MGLevelObject<DoFHandler<dim>> dofhandlers_p;
235  MGLevelObject<DoFHandler<dim>> dofhandlers_projection;
236 
237  MGLevelObject<AffineConstraints<double>> constraints_v;
238  MGLevelObject<AffineConstraints<double>> constraints_p;
239 
240 #if DEAL_II_VERSION_GTE(9,6,0)
242 #else
243  using transfer_t = MGTransferGlobalCoarsening<dim, ::LinearAlgebra::distributed::Vector<GMGNumberType>>;
244 #endif
245  MGLevelObject<MGTwoLevelTransfer<dim, ::LinearAlgebra::distributed::Vector<GMGNumberType>>> transfers_v;
246  MGLevelObject<MGTwoLevelTransfer<dim, ::LinearAlgebra::distributed::Vector<GMGNumberType>>> transfers_p;
247  std::unique_ptr<transfer_t> mg_transfer_A_block;
248  std::unique_ptr<transfer_t> mg_transfer_Schur_complement;
249  };
250 }
251 
252 #endif
void parse_parameters(ParameterHandler &prm) override
const Simulator< dim > * simulator
::TrilinosWrappers::MPI::BlockVector BlockVector
Definition: global.h:269
std::vector< std::shared_ptr< const Triangulation< dim, dim > > > trias
MGTransferGlobalCoarsening< dim, ::LinearAlgebra::distributed::Vector< GMGNumberType >> transfer_t
MGLevelObject< MatrixFreeStokesOperators::OperatorCellData< dim, GMGNumberType > > level_cell_data
StokesSolver::SolverOutputs solve(const LinearAlgebra::BlockSparseMatrix &system_matrix, const LinearAlgebra::BlockVector &system_rhs, const bool solve_newton_system, const double last_pressure_normalization_adjustment, LinearAlgebra::BlockVector &solution_vector) override
MGLevelObject< MGTwoLevelTransfer< dim, ::LinearAlgebra::distributed::Vector< GMGNumberType > > > transfers_p
::TrilinosWrappers::BlockSparseMatrix BlockSparseMatrix
Definition: global.h:280
MatrixFreeStokesOperators::OperatorCellData< dim, GMGNumberType > active_cell_data
double GMGNumberType
::MGTransferMatrixFree< dim, NumberType > MGTransferMF
Definition: compat.h:64
MGLevelObject< MGTwoLevelTransfer< dim, ::LinearAlgebra::distributed::Vector< GMGNumberType > > > transfers_v
std::vector< std::shared_ptr< MatrixFree< dim, double > > > matrix_free_objects
StokesMatrixFreeHandlerGlobalCoarseningImplementation(Simulator< dim > &simulator, const Parameters< dim > &parameters)