SpECTRE
v2025.01.30
|
An additive Schwarz subdomain solver for linear systems of equations \(Ax=b\). More...
#include <Schwarz.hpp>
Public Types | |
using | operand_tag = FieldsTag |
using | fields_tag = FieldsTag |
using | source_tag = SourceTag |
using | options_group = OptionsGroup |
using | component_list = tmpl::list<> |
using | observed_reduction_data_tags = observers::make_reduction_data_tags< tmpl::list< async_solvers::reduction_data, detail::reduction_data > > |
using | subdomain_solver = detail::subdomain_solver< FieldsTag, SubdomainOperator, SubdomainPreconditioners > |
using | initialize_element = tmpl::list< async_solvers::InitializeElement< FieldsTag, OptionsGroup, SourceTag >, detail::InitializeElement< FieldsTag, OptionsGroup, SubdomainOperator, SubdomainPreconditioners > > |
using | amr_projectors = initialize_element |
using | register_element = tmpl::list< async_solvers::RegisterElement< FieldsTag, OptionsGroup, SourceTag, ArraySectionIdTag >, detail::RegisterElement< FieldsTag, OptionsGroup, SourceTag, ArraySectionIdTag > > |
template<typename ApplyOperatorActions , typename Label = OptionsGroup> | |
using | solve = tmpl::list< async_solvers::PrepareSolve< FieldsTag, OptionsGroup, SourceTag, Label, ArraySectionIdTag >, detail::SendOverlapData< FieldsTag, OptionsGroup, SubdomainOperator >, detail::SolveSubdomain< FieldsTag, OptionsGroup, SubdomainOperator, ArraySectionIdTag >, detail::ReceiveOverlapSolution< FieldsTag, OptionsGroup, SubdomainOperator >, ApplyOperatorActions, async_solvers::CompleteStep< FieldsTag, OptionsGroup, SourceTag, Label, ArraySectionIdTag > > |
An additive Schwarz subdomain solver for linear systems of equations \(Ax=b\).
This Schwarz-type linear solver works by solving many sub-problems in parallel and combining their solutions as a weighted sum to converge towards the global solution. Each sub-problem is the restriction of the global problem to an element-centered subdomain, which consists of a central element and an overlap region with its neighbors. The decomposition into independent sub-problems makes this linear solver very parallelizable, avoiding global synchronization points altogether. It is commonly used as a preconditioner to Krylov-type solvers such as LinearSolver::gmres::Gmres
or LinearSolver::cg::ConjugateGradient
, or as the "smoother" for a Multigrid solver (which may in turn precondition a Krylov-type solver).
This linear solver relies on an implementation of the global linear operator \(A(x)\) and its restriction to a subdomain \(A_{ss}(x)\). Each step of the algorithm expects that \(A(x)\) is computed and stored in the DataBox as db::add_tag_prefix<LinearSolver::Tags::OperatorAppliedTo, operand_tag>
. To perform a solve, add the solve
action list to an array parallel component. Pass the actions that compute \(A(x)\), as well as any further actions you wish to run in each step of the algorithm, as the first template parameter to solve
. If you add the solve
action list multiple times, use the second template parameter to label each solve with a different type. Pass the subdomain operator as the SubdomainOperator
template parameter to this class (not to the solve
action list template). See the paragraph below for information on implementing a subdomain operator.
Fig. 1 shows part of a 2-dimensional computational domain. The domain is composed of elements (light gray) that each carry a Legendre-Gauss-Lobatto mesh of grid points (black dots). The diagonally shaded region to the left illustrates an external domain boundary. The lines between neighboring element faces illustrate mortar meshes, which are relevant when implementing a subdomain operator in a DG context but play no role in the Schwarz algorithm. The dashed line gives an example of an element-centered subdomain with a maximum of 2 grid points of overlap into neighboring elements. For details on how the number of overlap points is chosen see LinearSolver::Schwarz::overlap_extent
. Note that the subdomain does not extend into corner- or edge-neighbors, which is different to both [185] and [198]. The reason we don't include such diagonal couplings is that in a DG context information only propagates across faces, as noted in [185]. By eliminating the corner- and edge-neighbors we significantly reduce the complexity of the subdomain geometry and potentially also the communication costs. However, the advantages and disadvantages of this choice have yet to be carefully evaluated.
LinearSolver::Schwarz::SubdomainOperator
for details on how to implement a subdomain operator for your problem.LinearSolver::gmres::Gmres
that require global inner products in each iteration). The algorithm assumes that the action list passed to solve
applies the global linear operator, as described above.LinearSolver::Serial
are available to solve subdomains. They include Krylov-type methods that support nested preconditioning with any other linear solver. Additional problem-specific subdomain preconditioners can be added with the SubdomainPreconditioners
template parameter. Note that the choice of subdomain solver (and, by extension, the choice of subdomain preconditioner) affects only the performance of the Schwarz solver, not its convergence or parallelization properties (assuming the subdomain solutions it produces are sufficiently precise).LinearSolver::Schwarz::element_weight
and Section 3.1 of [185] for details. Note that we must account for the missing corner- and edge-neighbors when constructing the weights. See LinearSolver::Schwarz::intruding_weight
for a discussion.ArraySectionIdTag
template parameter if residual norms should be computed over a section. Pass void
(default) to compute residual norms over all elements in the array.