41 #include <settings.hpp> 43 #include <communicate.hpp> 44 #include <conv_tools.hpp> 45 #include <solver_tools.hpp> 61 template <
typename ValueType = gko::default_precision,
62 typename IndexType = gko::int32,
63 typename MixedValueType = gko::default_precision>
66 using ResidualCriterionFactory =
67 typename gko::stop::ResidualNormReduction<ValueType>::Factory;
68 using IterationCriterionFactory =
typename gko::stop::Iteration::Factory;
77 std::shared_ptr<gko::matrix::Dense<ValueType>> local_residual_vector;
79 std::shared_ptr<gko::matrix::Dense<ValueType>> residual_vector;
81 std::shared_ptr<gko::Array<IndexType>> convergence_vector;
83 std::shared_ptr<gko::Array<IndexType>> convergence_sent;
85 std::shared_ptr<gko::Array<IndexType>> convergence_local;
90 MPI_Win window_residual_vector;
95 MPI_Win window_convergence;
110 void setup_local_solver(
112 const std::shared_ptr<gko::matrix::Csr<ValueType, IndexType>>
114 std::shared_ptr<gko::matrix::Csr<ValueType, IndexType>>
115 &triangular_factor_l,
116 std::shared_ptr<gko::matrix::Csr<ValueType, IndexType>>
117 &triangular_factor_u,
118 std::shared_ptr<gko::matrix::Permutation<IndexType>> &local_perm,
119 std::shared_ptr<gko::matrix::Permutation<IndexType>> &local_inv_perm,
120 std::shared_ptr<gko::matrix::Dense<ValueType>> &local_rhs);
131 void compute_local_factors(
134 const std::shared_ptr<gko::matrix::Csr<ValueType, IndexType>>
152 const std::shared_ptr<gko::matrix::Csr<ValueType, IndexType>>
154 const std::shared_ptr<gko::matrix::Csr<ValueType, IndexType>>
155 &triangular_factor_l,
156 const std::shared_ptr<gko::matrix::Csr<ValueType, IndexType>>
157 &triangular_factor_u,
158 std::shared_ptr<gko::matrix::Permutation<IndexType>> &local_perm,
159 std::shared_ptr<gko::matrix::Permutation<IndexType>> &local_inv_perm,
160 std::shared_ptr<gko::matrix::Dense<ValueType>> &work_vector,
161 std::shared_ptr<gko::matrix::Dense<ValueType>> &init_guess,
162 std::shared_ptr<gko::matrix::Dense<ValueType>> &local_solution);
181 void check_convergence(
185 std::shared_ptr<gko::Array<IndexType>> &convergence_vector,
186 const std::shared_ptr<gko::matrix::Dense<ValueType>>
187 &global_old_solution,
188 const std::shared_ptr<gko::matrix::Dense<ValueType>> &global_solution,
189 const std::shared_ptr<gko::matrix::Csr<ValueType, IndexType>>
191 std::shared_ptr<gko::matrix::Dense<ValueType>> &work_vector,
192 ValueType &local_residual_norm, ValueType &local_residual_norm0,
193 ValueType &global_residual_norm, ValueType &global_residual_norm0,
194 int &num_converged_procs);
211 void check_global_convergence(
215 std::shared_ptr<gko::Array<IndexType>> &convergence_vector,
216 ValueType &local_resnorm, ValueType &local_resnorm0,
217 ValueType &global_resnorm, ValueType &global_resnorm0,
218 int &converged_all_local,
int &num_converged_procs);
234 bool check_local_convergence(
236 const std::shared_ptr<gko::matrix::Dense<ValueType>> &global_solution,
237 const std::shared_ptr<gko::matrix::Dense<ValueType>>
238 &global_old_solution,
239 const std::shared_ptr<gko::matrix::Csr<ValueType, IndexType>>
241 std::shared_ptr<gko::matrix::Dense<ValueType>> &work_vector,
242 ValueType &local_resnorm, ValueType &local_resnorm0);
253 void update_residual(
255 std::shared_ptr<gko::matrix::Dense<ValueType>> &global_solution,
256 const std::shared_ptr<gko::matrix::Csr<ValueType, IndexType>>
258 const std::shared_ptr<gko::matrix::Dense<ValueType>>
259 &global_old_solution);
274 void compute_residual_norm(
277 const std::shared_ptr<gko::matrix::Csr<ValueType, IndexType>>
279 const std::shared_ptr<gko::matrix::Dense<ValueType>> &global_solution,
280 const std::shared_ptr<gko::matrix::Dense<ValueType>> &global_rhs,
281 ValueType &mat_norm, ValueType &rhs_norm, ValueType &sol_norm,
282 ValueType &residual_norm);
293 std::shared_ptr<gko::LinOp> solver;
298 std::shared_ptr<gko::stop::Combined::Factory> combined_criterion;
308 std::shared_ptr<gko::log::Record> record_logger;
313 std::shared_ptr<gko::solver::LowerTrs<ValueType, IndexType>> L_solver;
318 std::shared_ptr<gko::solver::UpperTrs<ValueType, IndexType>> U_solver;
320 #if SCHW_HAVE_CHOLMOD 325 cholmod_common settings;
326 cholmod_sparse *system_matrix;
328 cholmod_factor *L_factor;
337 #if SCHW_HAVE_UMFPACK 346 std::shared_ptr<gko::matrix::Dense<ValueType>> row_scale;
348 double control[UMFPACK_CONTROL];
349 double info[UMFPACK_INFO];
The Solver class the provides the solver and the convergence checking methods.
Definition: solve.hpp:64
The communication struct used to store the communication data.
Definition: communicate.hpp:67
The struct that contains the solver settings and the parameters to be set by the user.
Definition: settings.hpp:77
The Schwarz wrappers namespace.
Definition: comm_helpers.hpp:49
The initialization class that provides methods for initialization of the solver.
Definition: initialization.hpp:72