SpECTRE
v2025.03.17
|
Items related to the Schwarz linear solver. More...
Namespaces | |
namespace | Actions |
Actions related to the Schwarz solver. | |
namespace | OptionTags |
Option tags related to the Schwarz solver. | |
namespace | Tags |
Tags related to the Schwarz solver. | |
Classes | |
struct | ElementCenteredSubdomainData |
Data on an element-centered subdomain. More... | |
struct | ElementCenteredSubdomainDataIterator |
Iterate over LinearSolver::Schwarz::ElementCenteredSubdomainData More... | |
class | OverlapIterator |
Iterate over grid points in a region that extends partially into the volume. More... | |
struct | Schwarz |
An additive Schwarz subdomain solver for linear systems of equations | |
class | SubdomainOperator |
Abstract base class for the subdomain operator, i.e. the linear operator restricted to an element-centered Schwarz subdomain. More... | |
Typedefs | |
template<size_t Dim> | |
using | OverlapId = DirectionalId< Dim > |
Identifies a subdomain region that overlaps with another element. | |
template<size_t Dim, typename ValueType > | |
using | OverlapMap = DirectionalIdMap< Dim, ValueType > |
Data structure that can store the ValueType on each possible overlap of an element-centered subdomain with its neighbors. Overlaps are identified by their OverlapId . | |
Functions | |
template<size_t Dim, typename LhsTagsList , typename RhsTagsList > | |
decltype(auto) | operator+ (ElementCenteredSubdomainData< Dim, LhsTagsList > lhs, const ElementCenteredSubdomainData< Dim, RhsTagsList > &rhs) |
template<size_t Dim, typename LhsTagsList , typename RhsTagsList > | |
decltype(auto) | operator- (ElementCenteredSubdomainData< Dim, LhsTagsList > lhs, const ElementCenteredSubdomainData< Dim, RhsTagsList > &rhs) |
template<size_t Dim, typename TagsList > | |
decltype(auto) | operator* (const double scalar, ElementCenteredSubdomainData< Dim, TagsList > data) |
template<size_t Dim, typename TagsList > | |
decltype(auto) | operator* (ElementCenteredSubdomainData< Dim, TagsList > data, const double scalar) |
template<size_t Dim, typename TagsList > | |
decltype(auto) | operator/ (ElementCenteredSubdomainData< Dim, TagsList > data, const double scalar) |
template<size_t Dim, typename TagsList > | |
std::ostream & | operator<< (std::ostream &os, const ElementCenteredSubdomainData< Dim, TagsList > &subdomain_data) |
template<size_t Dim, typename TagsList > | |
bool | operator== (const ElementCenteredSubdomainData< Dim, TagsList > &lhs, const ElementCenteredSubdomainData< Dim, TagsList > &rhs) |
template<size_t Dim, typename TagsList > | |
bool | operator!= (const ElementCenteredSubdomainData< Dim, TagsList > &lhs, const ElementCenteredSubdomainData< Dim, TagsList > &rhs) |
size_t | overlap_extent (size_t volume_extent, size_t max_overlap) |
The number of points that an overlap extends into the volume_extent More... | |
template<size_t Dim> | |
size_t | overlap_num_points (const Index< Dim > &volume_extents, size_t overlap_extent, size_t overlap_dimension) |
Total number of grid points in an overlap region that extends overlap_extent points into the volume_extents from either side in the overlap_dimension More... | |
double | overlap_width (size_t overlap_extent, const DataVector &collocation_points) |
Width of an overlap extending overlap_extent points into the collocation_points from either side. More... | |
template<size_t Dim, typename VolumeTagsList , typename OverlapTagsList > | |
void | add_overlap_data (const gsl::not_null< Variables< VolumeTagsList > * > volume_data, const Variables< OverlapTagsList > &overlap_data, const Index< Dim > &volume_extents, const size_t overlap_extent, const Direction< Dim > &direction) |
Add the overlap_data to the volume_data | |
DataVector | extruding_weight (const DataVector &logical_coords, double width, const Side &side) |
Weights for the solution on an element-centered subdomain, decreasing from 1 to 0.5 towards the side over the logical distance width , and further to 0 over the same distance outside the element. More... | |
DataVector | intruding_weight (const DataVector &logical_coords, double width, const Side &side) |
Weights for the intruding solution of a neighboring element-centered subdomain, increasing from 0 to 0.5 towards the side over the logical distance width , and further to 1 over the same distance outside the element. More... | |
template<size_t Dim, typename DataType , typename... TensorStructure> | |
void | data_on_overlap (const gsl::not_null< Tensor< DataType, TensorStructure... > * > restricted_tensor, const Tensor< DataType, TensorStructure... > &tensor, const Index< Dim > &volume_extents, const size_t overlap_extent, const Direction< Dim > &direction) |
The part of the tensor data that lies within the overlap region. | |
template<size_t Dim, typename DataType , typename... TensorStructure> | |
Tensor< DataType, TensorStructure... > | data_on_overlap (const Tensor< DataType, TensorStructure... > &tensor, const Index< Dim > &volume_extents, const size_t overlap_extent, const Direction< Dim > &direction) |
The part of the tensor data that lies within the overlap region. | |
template<size_t Dim, typename OverlapTagsList , typename VolumeTagsList > | |
void | data_on_overlap (const gsl::not_null< Variables< OverlapTagsList > * > overlap_data, const Variables< VolumeTagsList > &volume_data, const Index< Dim > &volume_extents, const size_t overlap_extent, const Direction< Dim > &direction) |
The part of the tensor data that lies within the overlap region. | |
template<size_t Dim, typename TagsList > | |
Variables< TagsList > | data_on_overlap (const Variables< TagsList > &volume_data, const Index< Dim > &volume_extents, const size_t overlap_extent, const Direction< Dim > &direction) |
The part of the tensor data that lies within the overlap region. | |
template<size_t Dim, typename ExtendedTagsList , typename OverlapTagsList > | |
void | extended_overlap_data (const gsl::not_null< Variables< ExtendedTagsList > * > extended_data, const Variables< OverlapTagsList > &overlap_data, const Index< Dim > &volume_extents, const size_t overlap_extent, const Direction< Dim > &direction) |
Extend the overlap data to the full mesh by filling it with zeros outside the overlap region. | |
template<size_t Dim, typename TagsList > | |
Variables< TagsList > | extended_overlap_data (const Variables< TagsList > &overlap_data, const Index< Dim > &volume_extents, const size_t overlap_extent, const Direction< Dim > &direction) |
Extend the overlap data to the full mesh by filling it with zeros outside the overlap region. | |
template<size_t Dim> | |
void | element_weight (gsl::not_null< Scalar< DataVector > * > element_weight, const tnsr::I< DataVector, Dim, Frame::ElementLogical > &logical_coords, const std::array< double, Dim > &overlap_widths, const std::unordered_set< Direction< Dim > > &external_boundaries) |
Weights for data on the central element of an element-centered subdomain. More... | |
template<size_t Dim> | |
Scalar< DataVector > | element_weight (const tnsr::I< DataVector, Dim, Frame::ElementLogical > &logical_coords, const std::array< double, Dim > &overlap_widths, const std::unordered_set< Direction< Dim > > &external_boundaries) |
Weights for data on the central element of an element-centered subdomain. More... | |
template<size_t Dim> | |
void | intruding_weight (gsl::not_null< Scalar< DataVector > * > weight, const tnsr::I< DataVector, Dim, Frame::ElementLogical > &logical_coords, const Direction< Dim > &direction, const std::array< double, Dim > &overlap_widths, size_t num_intruding_overlaps, const std::unordered_set< Direction< Dim > > &external_boundaries) |
Weights for data on overlap regions intruding into an element-centered subdomain. More... | |
template<size_t Dim> | |
Scalar< DataVector > | intruding_weight (const tnsr::I< DataVector, Dim, Frame::ElementLogical > &logical_coords, const Direction< Dim > &direction, const std::array< double, Dim > &overlap_widths, size_t num_intruding_overlaps, const std::unordered_set< Direction< Dim > > &external_boundaries) |
Weights for data on overlap regions intruding into an element-centered subdomain. More... | |
Items related to the Schwarz linear solver.
LinearSolver::Schwarz::Schwarz
Scalar< DataVector > LinearSolver::Schwarz::element_weight | ( | const tnsr::I< DataVector, Dim, Frame::ElementLogical > & | logical_coords, |
const std::array< double, Dim > & | overlap_widths, | ||
const std::unordered_set< Direction< Dim > > & | external_boundaries | ||
) |
Weights for data on the central element of an element-centered subdomain.
Constructs the weighting field
where LinearSolver::Schwarz::extruding_weight
and
void LinearSolver::Schwarz::element_weight | ( | gsl::not_null< Scalar< DataVector > * > | element_weight, |
const tnsr::I< DataVector, Dim, Frame::ElementLogical > & | logical_coords, | ||
const std::array< double, Dim > & | overlap_widths, | ||
const std::unordered_set< Direction< Dim > > & | external_boundaries | ||
) |
Weights for data on the central element of an element-centered subdomain.
Constructs the weighting field
where LinearSolver::Schwarz::extruding_weight
and
DataVector LinearSolver::Schwarz::extruding_weight | ( | const DataVector & | logical_coords, |
double | width, | ||
const Side & | side | ||
) |
Weights for the solution on an element-centered subdomain, decreasing from 1 to 0.5 towards the side
over the logical distance width
, and further to 0 over the same distance outside the element.
The weighting function over a full element-centered subdomain is
where smoothstep
, i.e. the quintic polynomial
(see Eq. (39) in [200]).
The LinearSolver::Schwarz::extruding_weight
and LinearSolver::Schwarz::intruding_weight
functions each compute one of the two terms in A
that overlaps with a neighboring element-centered subdomain B
. To combine solutions on A
and B
to a weighted solution on A
, multiply the solution on A
with the extruding_weight
and the solution on B
with the intruding_weight
, both evaluated at the logical coordinates in A
and at the side
of A
that faces B
.
DataVector LinearSolver::Schwarz::intruding_weight | ( | const DataVector & | logical_coords, |
double | width, | ||
const Side & | side | ||
) |
Weights for the intruding solution of a neighboring element-centered subdomain, increasing from 0 to 0.5 towards the side
over the logical distance width
, and further to 1 over the same distance outside the element.
Scalar< DataVector > LinearSolver::Schwarz::intruding_weight | ( | const tnsr::I< DataVector, Dim, Frame::ElementLogical > & | logical_coords, |
const Direction< Dim > & | direction, | ||
const std::array< double, Dim > & | overlap_widths, | ||
size_t | num_intruding_overlaps, | ||
const std::unordered_set< Direction< Dim > > & | external_boundaries | ||
) |
Weights for data on overlap regions intruding into an element-centered subdomain.
Constructs the weighting field LinearSolver::Schwarz::element_weight
for the data that overlaps with the central element of an element-centered subdomain. The weights are constructed in such a way that all weights at a grid point sum to one, i.e. the weight is conserved. The logical_coords
are the element-logical coordinates of the central element.
This function assumes that corner- and edge-neighbors of the central element are not part of the subdomain, which means that no contributions from those neighbors are expected although the weighting field is non-zero in overlap regions with those neighbors. Therefore, to retain conservation we must account for this missing weight by adding it to the central element, to the intruding overlaps from face-neighbors, or split it between the two. We choose to add the weight to the intruding overlaps, since that's where information from the corner- and edge-regions propagates through in a DG context.
void LinearSolver::Schwarz::intruding_weight | ( | gsl::not_null< Scalar< DataVector > * > | weight, |
const tnsr::I< DataVector, Dim, Frame::ElementLogical > & | logical_coords, | ||
const Direction< Dim > & | direction, | ||
const std::array< double, Dim > & | overlap_widths, | ||
size_t | num_intruding_overlaps, | ||
const std::unordered_set< Direction< Dim > > & | external_boundaries | ||
) |
Weights for data on overlap regions intruding into an element-centered subdomain.
Constructs the weighting field LinearSolver::Schwarz::element_weight
for the data that overlaps with the central element of an element-centered subdomain. The weights are constructed in such a way that all weights at a grid point sum to one, i.e. the weight is conserved. The logical_coords
are the element-logical coordinates of the central element.
This function assumes that corner- and edge-neighbors of the central element are not part of the subdomain, which means that no contributions from those neighbors are expected although the weighting field is non-zero in overlap regions with those neighbors. Therefore, to retain conservation we must account for this missing weight by adding it to the central element, to the intruding overlaps from face-neighbors, or split it between the two. We choose to add the weight to the intruding overlaps, since that's where information from the corner- and edge-regions propagates through in a DG context.
size_t LinearSolver::Schwarz::overlap_extent | ( | size_t | volume_extent, |
size_t | max_overlap | ||
) |
The number of points that an overlap extends into the volume_extent
In a dimension where an element has volume_extent
points, the overlap extent is the largest number under these constraints:
max_overlap
.volume_extent
.This means the overlap extent is always smaller than the volume_extent
. The reason for this constraint is that we define the width of the overlap as the element-logical coordinate distance from the face of the element to the first collocation point outside the overlap extent. Therefore, even an overlap region that covers the full element in width does not include the collocation point on the opposite side of the element.
Here's a few notes on the definition of the overlap extent and width:
size_t LinearSolver::Schwarz::overlap_num_points | ( | const Index< Dim > & | volume_extents, |
size_t | overlap_extent, | ||
size_t | overlap_dimension | ||
) |
Total number of grid points in an overlap region that extends overlap_extent
points into the volume_extents
from either side in the overlap_dimension
The overlap region has overlap_extent
points in the overlap_dimension
, and volume_extents
points in the other dimensions. The number of grid points returned by this function is the product of these extents.
double LinearSolver::Schwarz::overlap_width | ( | size_t | overlap_extent, |
const DataVector & | collocation_points | ||
) |
Width of an overlap extending overlap_extent
points into the collocation_points
from either side.
The "width" of an overlap is the element-logical coordinate distance from the element boundary to the first collocation point outside the overlap region in the overlap dimension, i.e. the dimension perpendicular to the element face. See LinearSolver::Schwarz::overlap_extent
for details.
This function assumes the collocation_points
are mirrored around 0.