SpECTRE  v2024.04.12
LinearSolver::Schwarz::SubdomainOperator< Dim > Class Template Reference

Abstract base class for the subdomain operator, i.e. the linear operator restricted to an element-centered Schwarz subdomain. More...

#include <SubdomainOperator.hpp>

Static Public Attributes

static constexpr size_t volume_dim = Dim
 

Detailed Description

template<size_t Dim>
class LinearSolver::Schwarz::SubdomainOperator< Dim >

Abstract base class for the subdomain operator, i.e. the linear operator restricted to an element-centered Schwarz subdomain.

A subdomain operator must implement these member function templates

Since the subdomain operator has access to the element's DataBox, it can retrieve any background data, i.e. data that does not depend on the variables it operates on. Background data on overlap regions with other elements can be either initialized in advance if possible or communicated with the LinearSolver::Schwarz::Actions::SendOverlapFields and LinearSolver::Schwarz::Actions::ReceiveOverlapFields actions. Note that such communication should not be necessary between iterations of the Schwarz solve, but only between successive solves, because background data should not change during the solve. The Schwarz algorithm takes care of communicating all variable data that the subdomain operator operates on, whenever necessary. This variable data is passed to the operator as the operand argument (see above). It includes the data on the central element of the subdomain, as well as the data on overlap regions with neighbors. Since data on the entire subdomain is available, applying the subdomain operator requires no communication between elements. This is the strength of the Schwarz algorithm: all subdomain solves are independent of each other (see LinearSolver::Schwarz::Schwarz for details).

Here's an example of a subdomain operator that is the restriction of an explicit global matrix:

// Applies the linear operator given explicitly in the input file to data on an
// element-centered subdomain, using standard matrix multiplication.
//
// We assume all elements have the same extents so we can use the element's
// intruding overlap extents instead of initializing extruding overlap extents.
struct SubdomainOperator : LinearSolver::Schwarz::SubdomainOperator<1> {
template <typename ResultTags, typename OperandTags, typename DbTagsList>
void operator()(
result,
operand,
const db::DataBox<DbTagsList>& box) const {
// Retrieve arguments from the DataBox
const auto& matrix_slices =
db::get<helpers_distributed::LinearOperator>(box);
const auto& element = db::get<domain::Tags::Element<1>>(box);
const auto& overlap_extents = db::get<
const size_t element_index = helpers_distributed::get_index(element.id());
const size_t num_elements = matrix_slices.size();
const size_t num_points_per_element = matrix_slices.begin()->columns();
// Re-size the result buffer if necessary
result->destructive_resize(operand);
// Assemble full operator matrix
const auto operator_matrix = combine_matrix_slices(matrix_slices);
// Extend subdomain data with zeros outside the subdomain
const auto extended_operand =
extend_subdomain_data(operand, element_index, num_elements,
num_points_per_element, overlap_extents[0]);
// Apply matrix to extended subdomain data
const blaze::DynamicVector<double> extended_result =
operator_matrix * extended_operand;
// Restrict the result back to the subdomain
restrict_to_subdomain(result, extended_result, element_index,
num_points_per_element, overlap_extents[0]);
}
};
Abstract base class for the subdomain operator, i.e. the linear operator restricted to an element-cen...
Definition: SubdomainOperator.hpp:50
Require a pointer to not be a nullptr
Definition: Gsl.hpp:183
const auto & get(const Access &box)
Retrieve a tag from a db::Access
Definition: Access.hpp:60
Data on an element-centered subdomain.
Definition: ElementCenteredSubdomainData.hpp:59
The number of points a neighbor's subdomain extends into the element.
Definition: Tags.hpp:146

The documentation for this class was generated from the following file: