schwz  Generated automatically from develop
Public Member Functions | List of all members
schwz::SolverRAS< ValueType, IndexType, MixedValueType > Class Template Reference

An implementation of the solver interface using the RAS solver. More...

#include <restricted_schwarz.hpp>

Collaboration diagram for schwz::SolverRAS< ValueType, IndexType, MixedValueType >:
[legend]

Public Member Functions

 SolverRAS (Settings &settings, Metadata< ValueType, IndexType > &metadata)
 The constructor that takes in the user settings and a metadata struct containing the solver metadata. More...
 
void setup_local_matrices (Settings &settings, Metadata< ValueType, IndexType > &metadata, std::vector< unsigned int > &partition_indices, std::shared_ptr< gko::matrix::Csr< ValueType, IndexType >> &global_matrix, std::shared_ptr< gko::matrix::Csr< ValueType, IndexType >> &local_matrix, std::shared_ptr< gko::matrix::Csr< ValueType, IndexType >> &interface_matrix) override
 Sets up the local and the interface matrices from the global matrix and the partition indices. More...
 
void setup_comm_buffers () override
 Sets up the communication buffers needed for the boundary exchange.
 
void setup_windows (const Settings &settings, const Metadata< ValueType, IndexType > &metadata, std::shared_ptr< gko::matrix::Dense< ValueType >> &main_buffer) override
 Sets up the windows needed for the asynchronous communication. More...
 
void exchange_boundary (const Settings &settings, const Metadata< ValueType, IndexType > &metadata, std::shared_ptr< gko::matrix::Dense< ValueType >> &global_solution) override
 Exchanges the elements of the solution vector. More...
 
void update_boundary (const Settings &settings, const Metadata< ValueType, IndexType > &metadata, std::shared_ptr< gko::matrix::Dense< ValueType >> &local_solution, const std::shared_ptr< gko::matrix::Dense< ValueType >> &local_rhs, const std::shared_ptr< gko::matrix::Dense< ValueType >> &global_solution, const std::shared_ptr< gko::matrix::Csr< ValueType, IndexType >> &interface_matrix) override
 Update the values into local vector from obtained from the neighboring sub-domains using the interface matrix. More...
 
- Public Member Functions inherited from schwz::SchwarzBase< ValueType, IndexType, MixedValueType >
 SchwarzBase (Settings &settings, Metadata< ValueType, IndexType > &metadata)
 The constructor that takes in the user settings and a metadata struct containing the solver metadata. More...
 
void initialize ()
 Initialize the matrix and vectors.
 
void run (std::shared_ptr< gko::matrix::Dense< ValueType >> &solution)
 The function that runs the actual solver and obtains the final solution. More...
 
void print_vector (const std::shared_ptr< gko::matrix::Dense< ValueType >> &vector, int subd, std::string name)
 The auxiliary function that prints a passed in vector. More...
 
void print_matrix (const std::shared_ptr< gko::matrix::Csr< ValueType, IndexType >> &matrix, int rank, std::string name)
 The auxiliary function that prints a passed in CSR matrix. More...
 
- Public Member Functions inherited from schwz::Initialize< ValueType, IndexType >
 Initialize (Settings &settings, Metadata< ValueType, IndexType > &metadata)
 
void generate_rhs (std::vector< ValueType > &rhs)
 Generates the right hand side vector. More...
 
void setup_global_matrix (const std::string &filename, const gko::size_type &oned_laplacian_size, std::shared_ptr< gko::matrix::Csr< ValueType, IndexType >> &global_matrix)
 Generates the 2D global laplacian matrix. More...
 
void partition (const Settings &settings, const Metadata< ValueType, IndexType > &metadata, const std::shared_ptr< gko::matrix::Csr< ValueType, IndexType >> &global_matrix, std::vector< unsigned int > &partition_indices)
 The partitioning function. More...
 
void setup_vectors (const Settings &settings, const Metadata< ValueType, IndexType > &metadata, std::vector< ValueType > &rhs, std::shared_ptr< gko::matrix::Dense< ValueType >> &local_rhs, std::shared_ptr< gko::matrix::Dense< ValueType >> &global_rhs, std::shared_ptr< gko::matrix::Dense< ValueType >> &local_solution)
 Setup the vectors with default values and allocate mameory if not allocated. More...
 
- Public Member Functions inherited from schwz::Settings
 Settings (std::string executor_string="reference")
 
- Public Member Functions inherited from schwz::Communicate< ValueType, IndexType, MixedValueType >
void local_to_global_vector (const Settings &settings, const Metadata< ValueType, IndexType > &metadata, const std::shared_ptr< gko::matrix::Dense< ValueType >> &local_vector, std::shared_ptr< gko::matrix::Dense< ValueType >> &global_vector)
 Transforms data from a local vector to a global vector. More...
 
void clear (Settings &settings)
 Clears the data.
 
- Public Member Functions inherited from schwz::Solve< ValueType, IndexType, MixedValueType >
 Solve (const Settings &settings)
 

Additional Inherited Members

- Public Types inherited from schwz::Settings
enum  partition_settings {
  partition_regular = 0x0, partition_regular2d = 0x4, partition_metis = 0x1, partition_zoltan = 0x2,
  partition_custom = 0x3
}
 The partition algorithm to be used for partitioning the matrix.
 
enum  local_solver_settings {
  direct_solver_cholmod = 0x0, direct_solver_umfpack = 0x5, direct_solver_ginkgo = 0x1, iterative_solver_ginkgo = 0x2,
  iterative_solver_dealii = 0x3, solver_custom = 0x4
}
 The local solver algorithm for the local subdomain solves.
 
- Public Types inherited from schwz::Solve< ValueType, IndexType, MixedValueType >
using ResidualCriterionFactory = typename gko::stop::ResidualNormReduction< ValueType >::Factory
 
using IterationCriterionFactory = typename gko::stop::Iteration::Factory
 
- Public Attributes inherited from schwz::SchwarzBase< ValueType, IndexType, MixedValueType >
std::shared_ptr< gko::matrix::Csr< ValueType, IndexType > > local_matrix
 The local subdomain matrix.
 
std::shared_ptr< gko::matrix::Permutation< IndexType > > local_perm
 The local subdomain permutation matrix/array.
 
std::shared_ptr< gko::matrix::Permutation< IndexType > > local_inv_perm
 The local subdomain inverse permutation matrix/array.
 
std::shared_ptr< gko::matrix::Csr< ValueType, IndexType > > triangular_factor_l
 The local lower triangular factor used for the triangular solves.
 
std::shared_ptr< gko::matrix::Csr< ValueType, IndexType > > triangular_factor_u
 The local upper triangular factor used for the triangular solves.
 
std::shared_ptr< gko::matrix::Csr< ValueType, IndexType > > interface_matrix
 The local interface matrix.
 
std::shared_ptr< gko::matrix::Csr< ValueType, IndexType > > global_matrix
 The global matrix.
 
std::shared_ptr< gko::matrix::Dense< ValueType > > local_rhs
 The local right hand side.
 
std::shared_ptr< gko::matrix::Dense< ValueType > > global_rhs
 The global right hand side.
 
std::shared_ptr< gko::matrix::Dense< ValueType > > local_solution
 The local solution vector.
 
std::shared_ptr< gko::matrix::Dense< ValueType > > global_solution
 The global solution vector.
 
std::vector< ValueType > local_residual_vector_out
 The global residual vector.
 
std::vector< std::vector< ValueType > > global_residual_vector_out
 The local residual vector.
 
- Public Attributes inherited from schwz::Initialize< ValueType, IndexType >
std::vector< unsigned int > partition_indices
 The partition indices containing the subdomains to which each row(vertex) of the matrix(graph) belongs to.
 
std::vector< unsigned int > cell_weights
 The cell weights for the partition algorithm.
 
- Public Attributes inherited from schwz::Settings
std::string executor_string
 The string that contains the ginkgo executor paradigm.
 
std::shared_ptr< gko::Executor > executor = gko::ReferenceExecutor::create()
 The ginkgo executor the code is to be executed on.
 
std::shared_ptr< device_guardcuda_device_guard
 The ginkgo executor the code is to be executed on.
 
partition_settings partition = partition_settings::partition_regular
 
gko::int32 overlap = 2
 The overlap between the subdomains.
 
std::string matrix_filename = "null"
 The string that contains the matrix file name to read from .
 
bool explicit_laplacian = true
 Flag if the laplacian matrix should be generated within the library. More...
 
bool use_mixed_precision = false
 Flag if mixed precision should be used.
 
bool enable_random_rhs = false
 Flag to enable a random rhs.
 
bool print_matrices = false
 Flag to enable printing of matrices.
 
bool debug_print = false
 Flag to enable some debug printing.
 
local_solver_settings local_solver
 
bool non_symmetric_matrix = false
 Is the matrix non-symmetric ? , Use GMRES for local solves.
 
unsigned int restart_iter = 1u
 The restart iter for the GMRES solver.
 
int reset_local_crit_iter = -1
 The global iter at which to reset the local solver criterion.
 
bool naturally_ordered_factor = false
 Disables the re-ordering of the matrix before computing the triangular factors during the CHOLMOD factorization. More...
 
std::string metis_objtype
 This setting defines the objective type for the metis partitioning.
 
bool use_precond = false
 Enable the block jacobi local preconditioner for the local solver.
 
bool write_debug_out = false
 Enable the writing of debug out to file.
 
bool write_iters_and_residuals = false
 Enable writing the iters and residuals to a file.
 
bool enable_logging = false
 Flag to enable logging for local iterative solvers. More...
 
bool write_perm_data = false
 Enable the local permutations from CHOLMOD to a file.
 
int shifted_iter = 1
 Iteration shift for node local communication.
 
comm_settings comm_settings
 
convergence_settings convergence_settings
 
std::string factorization = "cholmod"
 The factorization for the local direct solver.
 
std::string reorder
 The reordering for the local solve.
 
- Public Attributes inherited from schwz::Communicate< ValueType, IndexType, MixedValueType >
comm_struct comm_struct
 

Detailed Description

template<typename ValueType = gko::default_precision, typename IndexType = gko::int32, typename MixedValueType = gko::default_precision>
class schwz::SolverRAS< ValueType, IndexType, MixedValueType >

An implementation of the solver interface using the RAS solver.

Template Parameters
ValueTypeThe type of the floating point values.
IndexTypeThe type of the index type values.

Constructor & Destructor Documentation

◆ SolverRAS()

template<typename ValueType , typename IndexType , typename MixedValueType >
schwz::SolverRAS< ValueType, IndexType, MixedValueType >::SolverRAS ( Settings settings,
Metadata< ValueType, IndexType > &  metadata 
)

The constructor that takes in the user settings and a metadata struct containing the solver metadata.

Parameters
settingsThe settings struct.
metadataThe metadata struct.
dataThe additional data struct.

Member Function Documentation

◆ exchange_boundary()

template<typename ValueType , typename IndexType , typename MixedValueType >
void schwz::SolverRAS< ValueType, IndexType, MixedValueType >::exchange_boundary ( const Settings settings,
const Metadata< ValueType, IndexType > &  metadata,
std::shared_ptr< gko::matrix::Dense< ValueType >> &  global_solution 
)
overridevirtual

Exchanges the elements of the solution vector.

Parameters
settingsThe settings struct.
metadataThe metadata struct.
global_solutionThe solution vector being exchanged between the subdomains.

Implements schwz::Communicate< ValueType, IndexType, MixedValueType >.

References schwz::Settings::comm_settings::enable_onesided, and schwz::SchwarzBase< ValueType, IndexType, MixedValueType >::global_solution.

◆ setup_local_matrices()

template<typename ValueType , typename IndexType , typename MixedValueType >
void schwz::SolverRAS< ValueType, IndexType, MixedValueType >::setup_local_matrices ( Settings settings,
Metadata< ValueType, IndexType > &  metadata,
std::vector< unsigned int > &  partition_indices,
std::shared_ptr< gko::matrix::Csr< ValueType, IndexType >> &  global_matrix,
std::shared_ptr< gko::matrix::Csr< ValueType, IndexType >> &  local_matrix,
std::shared_ptr< gko::matrix::Csr< ValueType, IndexType >> &  interface_matrix 
)
overridevirtual

Sets up the local and the interface matrices from the global matrix and the partition indices.

Parameters
settingsThe settings struct.
metadataThe metadata struct.
partition_indicesThe array containing the partition indices.
global_matrixThe global system matrix.
local_matrixThe local system matrix.
interface_matrixThe interface matrix containing the interface and the overlap data mainly used for exchanging values between different sub-domains.
local_permThe local permutation, obtained through RCM or METIS.

Implements schwz::Initialize< ValueType, IndexType >.

References schwz::Metadata< ValueType, IndexType >::comm_size, schwz::Settings::executor, schwz::Metadata< ValueType, IndexType >::first_row, schwz::SchwarzBase< ValueType, IndexType, MixedValueType >::global_matrix, schwz::Metadata< ValueType, IndexType >::global_size, schwz::Metadata< ValueType, IndexType >::global_to_local, schwz::Metadata< ValueType, IndexType >::i_permutation, schwz::SchwarzBase< ValueType, IndexType, MixedValueType >::interface_matrix, schwz::SchwarzBase< ValueType, IndexType, MixedValueType >::local_matrix, schwz::Metadata< ValueType, IndexType >::local_size, schwz::Metadata< ValueType, IndexType >::local_size_o, schwz::Metadata< ValueType, IndexType >::local_size_x, schwz::Metadata< ValueType, IndexType >::local_to_global, schwz::Metadata< ValueType, IndexType >::my_rank, schwz::Metadata< ValueType, IndexType >::num_subdomains, schwz::Settings::overlap, schwz::Metadata< ValueType, IndexType >::overlap_row, schwz::Metadata< ValueType, IndexType >::overlap_size, and schwz::Metadata< ValueType, IndexType >::permutation.

◆ setup_windows()

template<typename ValueType , typename IndexType , typename MixedValueType >
void schwz::SolverRAS< ValueType, IndexType, MixedValueType >::setup_windows ( const Settings settings,
const Metadata< ValueType, IndexType > &  metadata,
std::shared_ptr< gko::matrix::Dense< ValueType >> &  main_buffer 
)
overridevirtual

Sets up the windows needed for the asynchronous communication.

Parameters
settingsThe settings struct.
metadataThe metadata struct.
main_bufferThe main buffer being exchanged between the subdomains.

Implements schwz::Communicate< ValueType, IndexType, MixedValueType >.

References schwz::Settings::comm_settings::enable_get, schwz::Settings::comm_settings::enable_lock_all, schwz::Settings::comm_settings::enable_one_by_one, schwz::Settings::comm_settings::enable_onesided, schwz::Settings::comm_settings::enable_overlap, schwz::Settings::comm_settings::enable_put, schwz::Settings::executor, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::get_displacements, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::get_request, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::global_get, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::global_put, schwz::SchwarzBase< ValueType, IndexType, MixedValueType >::global_solution, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::is_local_neighbor, schwz::Metadata< ValueType, IndexType >::iter_count, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::local_get, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::local_neighbors_in, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::local_neighbors_out, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::local_num_neighbors_in, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::local_num_neighbors_out, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::local_put, schwz::Metadata< ValueType, IndexType >::local_size_o, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::mixedt_recv_buffer, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::mixedt_send_buffer, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::neighbors_in, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::neighbors_out, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::num_neighbors_in, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::num_neighbors_out, schwz::Metadata< ValueType, IndexType >::num_subdomains, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::put_displacements, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::put_request, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::recv_buffer, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::send_buffer, schwz::Settings::comm_settings::stage_through_host, schwz::Settings::use_mixed_precision, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::window_recv_buffer, schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::window_send_buffer, and schwz::Communicate< ValueType, IndexType, MixedValueType >::comm_struct::window_x.

◆ update_boundary()

template<typename ValueType , typename IndexType , typename MixedValueType >
void schwz::SolverRAS< ValueType, IndexType, MixedValueType >::update_boundary ( const Settings settings,
const Metadata< ValueType, IndexType > &  metadata,
std::shared_ptr< gko::matrix::Dense< ValueType >> &  local_solution,
const std::shared_ptr< gko::matrix::Dense< ValueType >> &  local_rhs,
const std::shared_ptr< gko::matrix::Dense< ValueType >> &  global_solution,
const std::shared_ptr< gko::matrix::Csr< ValueType, IndexType >> &  interface_matrix 
)
overridevirtual

Update the values into local vector from obtained from the neighboring sub-domains using the interface matrix.

Parameters
settingsThe settings struct.
metadataThe metadata struct.
local_solutionThe local solution vector in the subdomain.
local_rhsThe local right hand side vector in the subdomain.
global_solutionThe workspace solution vector.
global_old_solutionThe global solution vector of the previous iteration.
interface_matrixThe interface matrix containing the interface and the overlap data mainly used for exchanging values between different sub-domains.

Implements schwz::Communicate< ValueType, IndexType, MixedValueType >.

References schwz::Settings::executor, schwz::SchwarzBase< ValueType, IndexType, MixedValueType >::global_solution, schwz::SchwarzBase< ValueType, IndexType, MixedValueType >::interface_matrix, schwz::SchwarzBase< ValueType, IndexType, MixedValueType >::local_rhs, schwz::Metadata< ValueType, IndexType >::local_size_x, schwz::SchwarzBase< ValueType, IndexType, MixedValueType >::local_solution, schwz::Metadata< ValueType, IndexType >::num_subdomains, and schwz::Settings::overlap.


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