18#ifdef DEAL_II_WITH_TRILINOS
26# include <AztecOO_StatusTest.h>
27# include <AztecOO_StatusTestCombo.h>
28# include <AztecOO_StatusTestMaxIters.h>
29# include <AztecOO_StatusTestResNorm.h>
30# include <AztecOO_StatusType.h>
41 const bool output_solver_details,
42 const unsigned int gmres_restart_parameter)
43 : output_solver_details(output_solver_details)
44 , gmres_restart_parameter(gmres_restart_parameter)
60 : solver_name(solver_name)
62 , additional_data(data)
84 const_cast<Epetra_CrsMatrix *
>(&
A.trilinos_matrix()),
86 const_cast<Epetra_MultiVector *
>(&
b.trilinos_vector()));
106 const_cast<Epetra_MultiVector *
>(
107 &
b.trilinos_vector()));
127 const_cast<Epetra_MultiVector *
>(
128 &
b.trilinos_vector()));
139 Epetra_MultiVector & x,
140 const Epetra_MultiVector &
b,
148 const_cast<Epetra_MultiVector *
>(
160 Epetra_MultiVector & x,
161 const Epetra_MultiVector &
b,
169 const_cast<Epetra_MultiVector *
>(
180 const ::Vector<double> &
b,
187 Assert(
A.local_range().second ==
A.m(),
188 ExcMessage(
"Can only work in serial when using deal.II vectors."));
189 Assert(
A.trilinos_matrix().Filled(),
190 ExcMessage(
"Matrix is not compressed. Call compress() method."));
192 Epetra_Vector ep_x(
View,
A.trilinos_matrix().DomainMap(), x.
begin());
193 Epetra_Vector ep_b(
View,
194 A.trilinos_matrix().RangeMap(),
195 const_cast<double *
>(
b.begin()));
200 const_cast<Epetra_CrsMatrix *
>(&
A.trilinos_matrix()), &ep_x, &ep_b);
210 const ::Vector<double> &
b,
213 Epetra_Vector ep_x(
View,
A.OperatorDomainMap(), x.
begin());
214 Epetra_Vector ep_b(
View,
215 A.OperatorRangeMap(),
216 const_cast<double *
>(
b.begin()));
220 linear_problem = std::make_unique<Epetra_LinearProblem>(&
A, &ep_x, &ep_b);
230 const ::LinearAlgebra::distributed::Vector<double> &
b,
236 A.trilinos_matrix().DomainMap().NumMyElements());
238 A.trilinos_matrix().RangeMap().NumMyElements());
240 Epetra_Vector ep_x(
View,
A.trilinos_matrix().DomainMap(), x.
begin());
241 Epetra_Vector ep_b(
View,
242 A.trilinos_matrix().RangeMap(),
243 const_cast<double *
>(
b.begin()));
248 const_cast<Epetra_CrsMatrix *
>(&
A.trilinos_matrix()), &ep_x, &ep_b);
258 const ::LinearAlgebra::distributed::Vector<double> &
b,
259 const PreconditionBase &preconditioner)
264 Epetra_Vector ep_x(
View,
A.OperatorDomainMap(), x.
begin());
265 Epetra_Vector ep_b(
View,
266 A.OperatorRangeMap(),
267 const_cast<double *
>(
b.begin()));
271 linear_problem = std::make_unique<Epetra_LinearProblem>(&
A, &ep_x, &ep_b);
282 compute_residual(
const Epetra_MultiVector *
const residual_vector)
284 Assert(residual_vector->NumVectors() == 1,
285 ExcMessage(
"Residual multivector holds more than one vector"));
287 const int ierr = residual_vector->Norm2(&res_l2_norm);
292 class TrilinosReductionControl :
public AztecOO_StatusTest
295 TrilinosReductionControl(
const int max_steps,
296 const double tolerance,
298 const Epetra_LinearProblem &linear_problem);
300 virtual ~TrilinosReductionControl()
override =
default;
303 ResidualVectorRequired()
const override
308 virtual AztecOO_StatusType
309 CheckStatus(
int CurrentIter,
310 Epetra_MultiVector *CurrentResVector,
311 double CurrentResNormEst,
312 bool SolutionUpdated)
override
317 (CurrentResNormEst < 0.0 ? compute_residual(CurrentResVector) :
319 if (CurrentIter == 0)
328 virtual AztecOO_StatusType
329 GetStatus()
const override
334 virtual std::ostream &
335 Print(std::ostream &stream,
int indent = 0)
const override
341 get_initial_residual()
const
347 get_current_residual()
const
362 TrilinosReductionControl::TrilinosReductionControl(
364 const double tolerance,
366 const Epetra_LinearProblem &linear_problem)
372 AztecOO_StatusTestCombo::OR))
377 std::make_unique<AztecOO_StatusTestMaxIters>(max_steps);
380 Assert(linear_problem.GetRHS()->NumVectors() == 1,
381 ExcMessage(
"RHS multivector holds more than one vector"));
385 *linear_problem.GetOperator(),
386 *(linear_problem.GetLHS()->operator()(0)),
387 *(linear_problem.GetRHS()->operator()(0)),
390 AztecOO_StatusTestResNorm::TwoNorm);
392 AztecOO_StatusTestResNorm::None, AztecOO_StatusTestResNorm::TwoNorm);
397 *linear_problem.GetOperator(),
398 *(linear_problem.GetLHS()->operator()(0)),
399 *(linear_problem.GetRHS()->operator()(0)),
402 AztecOO_StatusTestResNorm::TwoNorm);
404 AztecOO_StatusTestResNorm::NormOfInitRes,
405 AztecOO_StatusTestResNorm::TwoNorm);
413 template <
typename Preconditioner>
415 SolverBase::do_solve(
const Preconditioner &preconditioner)
420 solver.SetProblem(*linear_problem);
426 solver.SetAztecOption(AZ_solver, AZ_cg);
429 solver.SetAztecOption(AZ_solver, AZ_cgs);
432 solver.SetAztecOption(AZ_solver, AZ_gmres);
433 solver.SetAztecOption(AZ_kspace,
434 additional_data.gmres_restart_parameter);
437 solver.SetAztecOption(AZ_solver, AZ_bicgstab);
440 solver.SetAztecOption(AZ_solver, AZ_tfqmr);
447 set_preconditioner(solver, preconditioner);
450 solver.SetAztecOption(AZ_output,
451 additional_data.output_solver_details ? AZ_all :
453 solver.SetAztecOption(AZ_conv, AZ_noscaled);
473 status_test = std::make_unique<internal::TrilinosReductionControl>(
474 reduction_control->max_steps(),
475 reduction_control->tolerance(),
476 reduction_control->reduction(),
478 solver.SetStatusTest(status_test.get());
484 solver.Iterate(solver_control.max_steps(), solver_control.tolerance());
494 "option not implemented"));
499 "numerical breakdown"));
504 "loss of precision"));
509 "GMRES Hessenberg ill-conditioned"));
522 if (
const internal::TrilinosReductionControl
523 *
const reduction_control_status =
524 dynamic_cast<const internal::TrilinosReductionControl *const
>(
533 if (solver.NumIters() > 0)
538 solver_control.check(
539 0, reduction_control_status->get_initial_residual());
540 solver_control.check(
542 reduction_control_status->get_current_residual());
545 solver_control.check(
547 reduction_control_status->get_current_residual());
552 solver_control.check(solver.NumIters(), solver.TrueResidual());
558 solver_control.last_value()));
565 SolverBase::set_preconditioner(AztecOO & solver,
572 const int ierr = solver.SetPrecOperator(
577 solver.SetAztecOption(AZ_precond, AZ_none);
583 SolverBase::set_preconditioner(AztecOO & solver,
587 solver.SetPrecOperator(
const_cast<Epetra_Operator *
>(&preconditioner));
594 SolverCG::AdditionalData::AdditionalData(
const bool output_solver_details)
611 const bool output_solver_details,
612 const unsigned int restart_parameter)
629 const bool output_solver_details)
679 const std::string &solver_type)
680 : output_solver_details(output_solver_details)
681 , solver_type(solver_type)
710 const_cast<Epetra_CrsMatrix *
>(&
A.trilinos_matrix()));
727 std::string(
"You tried to select the solver type <") +
729 "> but this solver is not supported by Trilinos either "
730 "because it does not exist, or because Trilinos was not "
731 "configured for its use."));
736 verbose_cout <<
"Starting symbolic factorization" << std::endl;
737 ierr =
solver->SymbolicFactorization();
740 verbose_cout <<
"Starting numeric factorization" << std::endl;
741 ierr =
solver->NumericFactorization();
754 const_cast<Epetra_MultiVector *
>(&
b.trilinos_vector()));
762 verbose_cout <<
"Starting solve" << std::endl;
765 int ierr =
solver->Solve();
777 const ::LinearAlgebra::distributed::Vector<double> &
b)
779 Epetra_Vector ep_x(
View,
782 Epetra_Vector ep_b(
View,
784 const_cast<double *
>(
b.begin()));
797 verbose_cout <<
"Starting solve" << std::endl;
800 int ierr =
solver->Solve();
827 std::string(
"You tried to select the solver type <") +
829 "> but this solver is not supported by Trilinos either "
830 "because it does not exist, or because Trilinos was not "
831 "configured for its use."));
836 verbose_cout <<
"Starting symbolic factorization" << std::endl;
837 ierr =
solver->SymbolicFactorization();
840 verbose_cout <<
"Starting numeric factorization" << std::endl;
841 ierr =
solver->NumericFactorization();
844 verbose_cout <<
"Starting solve" << std::endl;
868 const_cast<Epetra_CrsMatrix *
>(&
A.trilinos_matrix()),
870 const_cast<Epetra_MultiVector *
>(&
b.trilinos_vector()));
880 const ::Vector<double> &
b)
886 Assert(
A.local_range().second ==
A.m(),
887 ExcMessage(
"Can only work in serial when using deal.II vectors."));
888 Epetra_Vector ep_x(
View,
A.trilinos_matrix().DomainMap(), x.
begin());
889 Epetra_Vector ep_b(
View,
890 A.trilinos_matrix().RangeMap(),
891 const_cast<double *
>(
b.begin()));
896 const_cast<Epetra_CrsMatrix *
>(&
A.trilinos_matrix()), &ep_x, &ep_b);
907 const ::LinearAlgebra::distributed::Vector<double> &
b)
910 A.trilinos_matrix().DomainMap().NumMyElements());
912 A.trilinos_matrix().RangeMap().NumMyElements());
913 Epetra_Vector ep_x(
View,
A.trilinos_matrix().DomainMap(), x.
begin());
914 Epetra_Vector ep_b(
View,
915 A.trilinos_matrix().RangeMap(),
916 const_cast<double *
>(
b.begin()));
921 const_cast<Epetra_CrsMatrix *
>(&
A.trilinos_matrix()), &ep_x, &ep_b);
unsigned int last_step() const
double last_value() const
virtual State check(const unsigned int step, const double check_value)
@ success
Stop iteration, goal reached.
SolverControl & control() const
void solve(const SparseMatrix &A, MPI::Vector &x, const MPI::Vector &b, const PreconditionBase &preconditioner)
const AdditionalData additional_data
SolverControl & solver_control
void do_solve(const Preconditioner &preconditioner)
std::unique_ptr< Epetra_LinearProblem > linear_problem
SolverBase(SolverControl &cn, const AdditionalData &data=AdditionalData())
enum TrilinosWrappers::SolverBase::SolverName solver_name
const AdditionalData additional_data
SolverBicgstab(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
SolverCGS(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
SolverCG(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
SolverDirect(SolverControl &cn, const AdditionalData &data=AdditionalData())
void initialize(const SparseMatrix &A)
std::unique_ptr< Amesos_BaseSolver > solver
SolverControl & solver_control
void solve(MPI::Vector &x, const MPI::Vector &b)
std::unique_ptr< Epetra_LinearProblem > linear_problem
SolverControl & control() const
const AdditionalData additional_data
SolverGMRES(SolverControl &cn, const AdditionalData &data=AdditionalData())
SolverTFQMR(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
__global__ void reduction(Number *result, const Number *v, const size_type N)
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
const Epetra_MultiVector & trilinos_vector() const
Teuchos::RCP< Epetra_Operator > preconditioner
size_type local_size() const
typename ::internal::FEValuesViews::ViewType< dim, spacedim, Extractor >::type View
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
AdditionalData(const bool output_solver_details=false, const unsigned int gmres_restart_parameter=30)
AdditionalData(const bool output_solver_details=false)
AdditionalData(const bool output_solver_details=false)
AdditionalData(const bool output_solver_details=false, const std::string &solver_type="Amesos_Klu")
bool output_solver_details
AdditionalData(const bool output_solver_details=false, const unsigned int restart_parameter=30)
AdditionalData(const bool output_solver_details=false)
std::unique_ptr< AztecOO_StatusTestResNorm > status_test_abs_tol
std::unique_ptr< AztecOO_StatusTestResNorm > status_test_rel_tol
std::unique_ptr< AztecOO_StatusTestMaxIters > status_test_max_steps
std::unique_ptr< AztecOO_StatusTestCombo > status_test_collection