|
| 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> |
| 25 | +#include <aspect/simulator/solver/stokes_matrix_free.h> |
| 26 | +#include <aspect/simulator/solver/matrix_free_operators.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 | +{ |
| 43 | + /** |
| 44 | + * Main class of the Matrix-free method. Here are all the functions for |
| 45 | + * setup, assembly and solving the Stokes system. |
| 46 | + * |
| 47 | + * We need to derive from StokesMatrixFreeHandler to be able to introduce a |
| 48 | + * second template argument for the degree of the Stokes finite |
| 49 | + * element. This way, the main simulator does not need to know about the |
| 50 | + * degree by using a pointer to the base class and we can pick the desired |
| 51 | + * velocity degree at runtime. |
| 52 | + */ |
| 53 | + template <int dim, int velocity_degree> |
| 54 | + class StokesMatrixFreeHandlerGlobalCoarseningImplementation: public StokesMatrixFreeHandler<dim> |
| 55 | + { |
| 56 | + public: |
| 57 | + /** |
| 58 | + * Constructor. Give it a reference to the |
| 59 | + * Simulator that owns it, since it needs to make fairly extensive |
| 60 | + * changes to the internals of the simulator. |
| 61 | + */ |
| 62 | + StokesMatrixFreeHandlerGlobalCoarseningImplementation(Simulator<dim> &simulator, |
| 63 | + const Parameters<dim> ¶meters); |
| 64 | + |
| 65 | + /** |
| 66 | + * Destructor. |
| 67 | + */ |
| 68 | + ~StokesMatrixFreeHandlerGlobalCoarseningImplementation() override = default; |
| 69 | + |
| 70 | + /** |
| 71 | + * Initialize the matrix-free solver. |
| 72 | + */ |
| 73 | + void initialize() override; |
| 74 | + |
| 75 | + /** |
| 76 | + * Return the name of the solver for screen output. |
| 77 | + */ |
| 78 | + std::string name() const override; |
| 79 | + |
| 80 | + /** |
| 81 | + * Solves the Stokes linear system using the matrix-free |
| 82 | + * solver. |
| 83 | + * |
| 84 | + * @param system_matrix The system matrix. Note that we do not actually |
| 85 | + * use this matrix for this matrix free solver. |
| 86 | + * @param system_rhs The right hand side vector of the system. |
| 87 | + * @param solve_newton_system A flag indicating whether the system to be |
| 88 | + * solved is the normal linear system or the Newton system. If the Newton |
| 89 | + * system is solved, some operations have to change, e.g. the residual |
| 90 | + * is computed differently. |
| 91 | + * @param solution_vector The solution vector that will be |
| 92 | + * updated with the new solution. This vector is expected to have the |
| 93 | + * block structure of the full solution vector, and its velocity and |
| 94 | + * pressure blocks will be updated with the new solution. |
| 95 | + * |
| 96 | + * @return A structure that contains information about the solver, like |
| 97 | + * the initial and final residual. |
| 98 | + */ |
| 99 | + StokesSolver::SolverOutputs |
| 100 | + solve(const LinearAlgebra::BlockSparseMatrix &system_matrix, |
| 101 | + const LinearAlgebra::BlockVector &system_rhs, |
| 102 | + const bool solve_newton_system, |
| 103 | + const double last_pressure_normalization_adjustment, |
| 104 | + LinearAlgebra::BlockVector &solution_vector) override; |
| 105 | + |
| 106 | + /** |
| 107 | + * Allocates and sets up the members of the StokesMatrixFreeHandler. This |
| 108 | + * is called by Simulator<dim>::setup_dofs() |
| 109 | + */ |
| 110 | + void setup_dofs() override; |
| 111 | + |
| 112 | + /** |
| 113 | + * Perform various tasks to update the linear system to solve |
| 114 | + * for. Note that we are not assembling a matrix (as this is a |
| 115 | + * matrix-free algorithm), but we are evaluating the material |
| 116 | + * model and storing the information necessary for a later call |
| 117 | + * to solve(). |
| 118 | + */ |
| 119 | + void assemble() override; |
| 120 | + |
| 121 | + /** |
| 122 | + * Computes and sets the diagonal for both the mass matrix operator and the A-block |
| 123 | + * operators on each level for the purpose of smoothing inside the multigrid v-cycle. |
| 124 | + */ |
| 125 | + void build_preconditioner() override; |
| 126 | + |
| 127 | + /** |
| 128 | + * Declare parameters. |
| 129 | + */ |
| 130 | + static |
| 131 | + void declare_parameters (ParameterHandler &prm); |
| 132 | + |
| 133 | + /** |
| 134 | + * Parse parameters. |
| 135 | + */ |
| 136 | + void parse_parameters (ParameterHandler &prm) override; |
| 137 | + |
| 138 | + /** |
| 139 | + * Return memory consumption in bytes for all DoFHandler objects. |
| 140 | + */ |
| 141 | + std::size_t get_dof_handler_memory_consumption() const override; |
| 142 | + |
| 143 | + /** |
| 144 | + * Return memory consumption in bytes for all transfer objects. |
| 145 | + */ |
| 146 | + std::size_t get_mg_transfer_memory_consumption() const override; |
| 147 | + |
| 148 | + /** |
| 149 | + * Return memory consumption in bytes for all transfer objects. |
| 150 | + */ |
| 151 | + std::size_t get_constraint_memory_consumption() const override; |
| 152 | + |
| 153 | + /** |
| 154 | + * Return the memory consumption in bytes that are used to store |
| 155 | + * equation data like viscosity to be able to apply the operators. |
| 156 | + */ |
| 157 | + std::size_t get_cell_data_memory_consumption() const override; |
| 158 | + |
| 159 | + private: |
| 160 | + /** |
| 161 | + * Evaluate the MaterialModel to query information like the viscosity and |
| 162 | + * project this viscosity to the multigrid hierarchy. Also queries |
| 163 | + * other parameters like pressure scaling. |
| 164 | + */ |
| 165 | + void evaluate_material_model(); |
| 166 | + |
| 167 | + /** |
| 168 | + * Add correction to system RHS for non-zero boundary condition. See description in |
| 169 | + * StokesMatrixFreeHandler::correct_stokes_rhs() for more information. |
| 170 | + */ |
| 171 | + void correct_stokes_rhs(); |
| 172 | + |
| 173 | + |
| 174 | + Simulator<dim> ∼ |
| 175 | + |
| 176 | + bool print_details; |
| 177 | + |
| 178 | + /** |
| 179 | + * If true, it will time the key components of this matrix-free implementation, such as |
| 180 | + * vmult of different matrices, solver IDR with the cheap preconditioner, etc. |
| 181 | + */ |
| 182 | + bool do_timings; |
| 183 | + |
| 184 | + /** |
| 185 | + * The max/min of the evaluated viscosities. |
| 186 | + */ |
| 187 | + double minimum_viscosity; |
| 188 | + double maximum_viscosity; |
| 189 | + |
| 190 | + FESystem<dim> fe_v; |
| 191 | + FESystem<dim> fe_p; |
| 192 | + FESystem<dim> fe_projection; |
| 193 | + |
| 194 | + /** |
| 195 | + * Store the data for the Stokes operator (viscosity, etc.) for the active cells. |
| 196 | + */ |
| 197 | + MatrixFreeStokesOperators::OperatorCellData<dim, GMGNumberType> active_cell_data; |
| 198 | + |
| 199 | + /** |
| 200 | + * Store the data for the Stokes operator (viscosity, etc.) for each multigrid level. |
| 201 | + */ |
| 202 | + MGLevelObject<MatrixFreeStokesOperators::OperatorCellData<dim, GMGNumberType>> level_cell_data; |
| 203 | + |
| 204 | + using StokesMatrixType = MatrixFreeStokesOperators::StokesOperator<dim,velocity_degree,double>; |
| 205 | + using SchurComplementMatrixType = MatrixFreeStokesOperators::MassMatrixOperator<dim,velocity_degree-1,double>; |
| 206 | + using ABlockMatrixType = MatrixFreeStokesOperators::ABlockOperator<dim,velocity_degree,double>; |
| 207 | + using BTBlockOperatorType = MatrixFreeStokesOperators::BTBlockOperator<dim,velocity_degree,double>; |
| 208 | + using GMGSchurComplementMatrixType = MatrixFreeStokesOperators::MassMatrixOperator<dim,velocity_degree-1,GMGNumberType>; |
| 209 | + using GMGABlockMatrixType = MatrixFreeStokesOperators::ABlockOperator<dim,velocity_degree,GMGNumberType>; |
| 210 | + |
| 211 | + StokesMatrixType stokes_matrix; |
| 212 | + ABlockMatrixType A_block_matrix; |
| 213 | + BTBlockOperatorType BT_block; |
| 214 | + SchurComplementMatrixType Schur_complement_block_matrix; |
| 215 | + |
| 216 | + MGLevelObject<GMGABlockMatrixType> mg_matrices_A_block; |
| 217 | + MGLevelObject<GMGSchurComplementMatrixType> mg_matrices_Schur_complement; |
| 218 | + |
| 219 | + MGConstrainedDoFs mg_constrained_dofs_A_block; |
| 220 | + MGConstrainedDoFs mg_constrained_dofs_Schur_complement; |
| 221 | + MGConstrainedDoFs mg_constrained_dofs_projection; |
| 222 | + |
| 223 | + unsigned int min_level; |
| 224 | + unsigned int max_level; |
| 225 | + |
| 226 | + std::vector<std::shared_ptr<MatrixFree<dim,double>>> matrix_free_objects; |
| 227 | + |
| 228 | + std::vector<std::shared_ptr<const Triangulation<dim, dim>>> trias; |
| 229 | + |
| 230 | + MGLevelObject<DoFHandler<dim>> dofhandlers_v; |
| 231 | + MGLevelObject<DoFHandler<dim>> dofhandlers_p; |
| 232 | + MGLevelObject<DoFHandler<dim>> dofhandlers_projection; |
| 233 | + |
| 234 | + MGLevelObject<AffineConstraints<double>> constraints_v; |
| 235 | + MGLevelObject<AffineConstraints<double>> constraints_p; |
| 236 | + |
| 237 | +#if DEAL_II_VERSION_GTE(9,6,0) |
| 238 | + using transfer_t = MGTransferMF<dim, GMGNumberType>; |
| 239 | +#else |
| 240 | + using transfer_t = MGTransferGlobalCoarsening<dim, dealii::LinearAlgebra::distributed::Vector<GMGNumberType>>; |
| 241 | +#endif |
| 242 | + MGLevelObject<MGTwoLevelTransfer<dim, dealii::LinearAlgebra::distributed::Vector<GMGNumberType>>> transfers_v; |
| 243 | + MGLevelObject<MGTwoLevelTransfer<dim, dealii::LinearAlgebra::distributed::Vector<GMGNumberType>>> transfers_p; |
| 244 | + std::unique_ptr<transfer_t> mg_transfer_A_block; |
| 245 | + std::unique_ptr<transfer_t> mg_transfer_Schur_complement; |
| 246 | + }; |
| 247 | +} |
| 248 | + |
| 249 | +#endif |
0 commit comments