16 #include <deal.II/lac/trilinos_solver.h> 18 #ifdef DEAL_II_WITH_TRILINOS 20 # include <deal.II/base/conditional_ostream.h> 21 # include <deal.II/base/std_cxx14/memory.h> 22 # include <deal.II/lac/trilinos_sparse_matrix.h> 23 # include <deal.II/lac/trilinos_vector.h> 24 # include <deal.II/lac/trilinos_precondition.h> 26 # include <AztecOO_StatusTest.h> 27 # include <AztecOO_StatusTestMaxIters.h> 28 # include <AztecOO_StatusTestResNorm.h> 29 # include <AztecOO_StatusTestCombo.h> 30 # include <AztecOO_StatusType.h> 35 DEAL_II_NAMESPACE_OPEN
41 const unsigned int gmres_restart_parameter)
43 output_solver_details (output_solver_details),
44 gmres_restart_parameter (gmres_restart_parameter)
63 solver_name (solver_name),
65 additional_data (data)
89 const_cast<Epetra_MultiVector *
>(&b.trilinos_vector()));
107 (
const_cast<Epetra_Operator *
>(&A),
109 const_cast<Epetra_MultiVector *
>(&b.trilinos_vector()));
122 const Epetra_Operator &preconditioner)
127 (
const_cast<Epetra_Operator *
>(&A),
129 const_cast<Epetra_MultiVector *
>(&b.trilinos_vector()));
140 Epetra_MultiVector &x,
141 const Epetra_MultiVector &b,
147 (
const_cast<Epetra_Operator *
>(&A),
149 const_cast<Epetra_MultiVector *>(&b));
160 Epetra_MultiVector &x,
161 const Epetra_MultiVector &b,
162 const Epetra_Operator &preconditioner)
167 (
const_cast<Epetra_Operator *
>(&A),
169 const_cast<Epetra_MultiVector *>(&b));
179 const ::Vector<double> &b,
189 ExcMessage (
"Can only work in serial when using deal.II vectors."));
191 ExcMessage (
"Matrix is not compressed. Call compress() method."));
194 Epetra_Vector ep_b (View, A.
range_partitioner(),
const_cast<double *
>(b.begin()));
210 const ::Vector<double> &b,
213 Epetra_Vector ep_x (View, A.OperatorDomainMap(), x.
begin());
214 Epetra_Vector ep_b (View, A.OperatorRangeMap(),
const_cast<double *
>(b.begin()));
218 linear_problem = std_cxx14::make_unique<Epetra_LinearProblem>(&A,&ep_x, &ep_b);
228 const ::LinearAlgebra::distributed::Vector<double> &b,
235 AssertDimension (static_cast<TrilinosWrappers::types::int_type>(b.local_size()),
239 Epetra_Vector ep_b (View, A.
range_partitioner(),
const_cast<double *
>(b.begin()));
255 const ::LinearAlgebra::distributed::Vector<double> &b,
259 A.OperatorDomainMap().NumMyElements());
260 AssertDimension (static_cast<TrilinosWrappers::types::int_type>(b.local_size()),
261 A.OperatorRangeMap().NumMyElements());
263 Epetra_Vector ep_x (View, A.OperatorDomainMap(), x.
begin());
264 Epetra_Vector ep_b (View, A.OperatorRangeMap(),
const_cast<double *
>(b.begin()));
268 linear_problem = std_cxx14::make_unique<Epetra_LinearProblem>(&A,&ep_x, &ep_b);
279 compute_residual (
const Epetra_MultiVector *
const residual_vector)
281 Assert(residual_vector->NumVectors() == 1,
282 ExcMessage(
"Residual multivector holds more than one vector"));
283 TrilinosScalar res_l2_norm = 0.0;
284 const int ierr = residual_vector->Norm2 (&res_l2_norm);
289 class TrilinosReductionControl :
public AztecOO_StatusTest
292 TrilinosReductionControl (
const int &max_steps,
293 const double &tolerance,
294 const double &reduction,
295 const Epetra_LinearProblem &linear_problem);
297 virtual ~TrilinosReductionControl() =
default;
300 ResidualVectorRequired ()
const 302 return status_test_collection->ResidualVectorRequired();
305 virtual AztecOO_StatusType
306 CheckStatus (
int CurrentIter,
307 Epetra_MultiVector *CurrentResVector,
308 double CurrentResNormEst,
309 bool SolutionUpdated)
313 current_residual = (CurrentResNormEst < 0.0 ?
314 compute_residual(CurrentResVector) :
316 if (CurrentIter == 0)
317 initial_residual = current_residual;
319 return status_test_collection->CheckStatus(CurrentIter,
326 virtual AztecOO_StatusType
329 return status_test_collection->GetStatus();
332 virtual std::ostream &
333 Print (std::ostream &stream,
334 int indent = 0)
const 336 return status_test_collection->Print(stream,indent);
340 get_initial_residual()
const 342 return initial_residual;
346 get_current_residual()
const 348 return current_residual;
352 double initial_residual;
353 double current_residual;
354 std::unique_ptr<AztecOO_StatusTestCombo> status_test_collection;
355 std::unique_ptr<AztecOO_StatusTestMaxIters> status_test_max_steps;
356 std::unique_ptr<AztecOO_StatusTestResNorm> status_test_abs_tol;
357 std::unique_ptr<AztecOO_StatusTestResNorm> status_test_rel_tol;
361 TrilinosReductionControl::TrilinosReductionControl(
362 const int &max_steps,
363 const double &tolerance,
364 const double &reduction,
365 const Epetra_LinearProblem &linear_problem )
367 initial_residual (
std::numeric_limits<double>::
max()),
368 current_residual(
std::numeric_limits<double>::
max())
372 status_test_collection = std_cxx14::make_unique<AztecOO_StatusTestCombo>
373 (AztecOO_StatusTestCombo::OR);
377 status_test_max_steps = std_cxx14::make_unique<AztecOO_StatusTestMaxIters>(max_steps);
378 status_test_collection->AddStatusTest(*status_test_max_steps);
381 ExcMessage(
"RHS multivector holds more than one vector"));
384 status_test_abs_tol = std_cxx14::make_unique<AztecOO_StatusTestResNorm>
389 status_test_abs_tol->DefineResForm(AztecOO_StatusTestResNorm::Explicit,
390 AztecOO_StatusTestResNorm::TwoNorm);
391 status_test_abs_tol->DefineScaleForm(AztecOO_StatusTestResNorm::None,
392 AztecOO_StatusTestResNorm::TwoNorm);
393 status_test_collection->AddStatusTest(*status_test_abs_tol);
396 status_test_rel_tol = std_cxx14::make_unique<AztecOO_StatusTestResNorm>
401 status_test_rel_tol->DefineResForm(AztecOO_StatusTestResNorm::Explicit,
402 AztecOO_StatusTestResNorm::TwoNorm);
403 status_test_rel_tol->DefineScaleForm(AztecOO_StatusTestResNorm::NormOfInitRes,
404 AztecOO_StatusTestResNorm::TwoNorm);
405 status_test_collection->AddStatusTest(*status_test_rel_tol);
412 template <
typename Preconditioner>
414 SolverBase::do_solve(
const Preconditioner &preconditioner)
419 solver.SetProblem(*linear_problem);
425 solver.SetAztecOption(AZ_solver, AZ_cg);
428 solver.SetAztecOption(AZ_solver, AZ_cgs);
431 solver.SetAztecOption(AZ_solver, AZ_gmres);
432 solver.SetAztecOption(AZ_kspace, additional_data.gmres_restart_parameter);
435 solver.SetAztecOption(AZ_solver, AZ_bicgstab);
438 solver.SetAztecOption(AZ_solver, AZ_tfqmr);
445 set_preconditioner(solver, preconditioner);
448 solver.SetAztecOption (AZ_output, additional_data.output_solver_details ?
450 solver.SetAztecOption (AZ_conv, AZ_noscaled);
468 = dynamic_cast<const ReductionControl *const>(&solver_control))
470 status_test = std_cxx14::make_unique<internal::TrilinosReductionControl>
471 (reduction_control->max_steps(),
472 reduction_control->tolerance(),
473 reduction_control->reduction(),
475 solver.SetStatusTest(status_test.get());
480 ierr = solver.Iterate (solver_control.max_steps(),
481 solver_control.tolerance());
490 "option not implemented"));
494 "numerical breakdown"));
498 "loss of precision"));
502 "GMRES Hessenberg ill-conditioned"));
515 if (
const internal::TrilinosReductionControl*
const reduction_control_status
516 = dynamic_cast<const internal::TrilinosReductionControl *const>(status_test.get()))
523 if (solver.NumIters() > 0)
528 solver_control.check (0, reduction_control_status->get_initial_residual());
529 solver_control.check (solver.NumIters(), reduction_control_status->get_current_residual());
532 solver_control.check (solver.NumIters(), reduction_control_status->get_current_residual());
537 solver_control.check (solver.NumIters(), solver.TrueResidual());
542 solver_control.last_value()));
549 SolverBase::set_preconditioner(AztecOO &solver,
556 const int ierr = solver.SetPrecOperator (const_cast<Epetra_Operator *>
561 solver.SetAztecOption(AZ_precond,AZ_none);
567 SolverBase::set_preconditioner(AztecOO &solver,
568 const Epetra_Operator &preconditioner)
570 const int ierr = solver.SetPrecOperator (const_cast<Epetra_Operator *>(&preconditioner));
577 SolverCG::AdditionalData::
578 AdditionalData (
const bool output_solver_details)
599 const unsigned int restart_parameter)
683 const std::string &solver_type)
685 output_solver_details (output_solver_details),
686 solver_type(solver_type)
733 ExcMessage (std::string (
"You tried to select the solver type <") +
735 "> but this solver is not supported by Trilinos either " 736 "because it does not exist, or because Trilinos was not " 737 "configured for its use.")
744 verbose_cout <<
"Starting symbolic factorization" << std::endl;
745 ierr =
solver->SymbolicFactorization();
748 verbose_cout <<
"Starting numeric factorization" << std::endl;
749 ierr =
solver->NumericFactorization();
760 linear_problem->SetRHS(const_cast<Epetra_MultiVector *>(&b.trilinos_vector()));
768 verbose_cout <<
"Starting solve" << std::endl;
771 int ierr =
solver->Solve ();
781 const ::LinearAlgebra::distributed::Vector<double> &b)
784 Epetra_Vector ep_b (View,
linear_problem->GetOperator()->OperatorRangeMap(),
const_cast<double *
>(b.begin()));
797 verbose_cout <<
"Starting solve" << std::endl;
800 int ierr =
solver->Solve ();
827 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."));
837 verbose_cout <<
"Starting symbolic factorization" << std::endl;
838 ierr =
solver->SymbolicFactorization ();
841 verbose_cout <<
"Starting numeric factorization" << std::endl;
842 ierr =
solver->NumericFactorization ();
845 verbose_cout <<
"Starting solve" << std::endl;
870 const_cast<Epetra_MultiVector *
>(&b.trilinos_vector()));
880 const ::Vector<double> &b)
890 ExcMessage (
"Can only work in serial when using deal.II vectors."));
892 Epetra_Vector ep_b (View, A.
range_partitioner(),
const_cast<double *
>(b.begin()));
908 const ::LinearAlgebra::distributed::Vector<double> &b)
912 AssertDimension (static_cast<TrilinosWrappers::types::int_type>(b.local_size()),
915 Epetra_Vector ep_b (View, A.
range_partitioner(),
const_cast<double *
>(b.begin()));
939 DEAL_II_NAMESPACE_CLOSE
941 #endif // DEAL_II_WITH_PETSC static ::ExceptionBase & ExcTrilinosError(int arg1)
SolverCG(SolverControl &cn, const AdditionalData &data=AdditionalData())
#define AssertDimension(dim1, dim2)
virtual State check(const unsigned int step, const double check_value)
std::unique_ptr< Amesos_BaseSolver > solver
const Epetra_Map & range_partitioner() const
const AdditionalData additional_data
const AdditionalData additional_data
AdditionalData(const bool output_solver_details=false, const unsigned int restart_parameter=30)
const AdditionalData additional_data
void initialize(const SparseMatrix &A)
void solve(MPI::Vector &x, const MPI::Vector &b)
#define AssertThrow(cond, exc)
const AdditionalData additional_data
void do_solve(const Preconditioner &preconditioner)
SolverControl & solver_control
static ::ExceptionBase & ExcMessage(std::string arg1)
SolverCGS(SolverControl &cn, const AdditionalData &data=AdditionalData())
Stop iteration, goal reached.
#define Assert(cond, exc)
AdditionalData(const bool output_solver_details=false)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
SolverBicgstab(SolverControl &cn, const AdditionalData &data=AdditionalData())
SolverControl & solver_control
std::unique_ptr< Epetra_LinearProblem > linear_problem
SolverControl & control() const
void solve(const SparseMatrix &A, MPI::Vector &x, const MPI::Vector &b, const PreconditionBase &preconditioner)
size_type local_size() const
const Epetra_CrsMatrix & trilinos_matrix() const
std::unique_ptr< Epetra_LinearProblem > linear_problem
double last_value() const
SolverTFQMR(SolverControl &cn, const AdditionalData &data=AdditionalData())
AdditionalData(const bool output_solver_details=false)
SolverBase(SolverControl &cn, const AdditionalData &data=AdditionalData())
unsigned int last_step() const
const AdditionalData additional_data
AdditionalData(const bool output_solver_details=false, const std::string &solver_type="Amesos_Klu")
SolverDirect(SolverControl &cn, const AdditionalData &data=AdditionalData())
const Epetra_Map & domain_partitioner() const
const AdditionalData additional_data
bool output_solver_details
std::pair< size_type, size_type > local_range() const
AdditionalData(const bool output_solver_details=false)
SolverControl & control() const
static ::ExceptionBase & ExcNotImplemented()
AdditionalData(const bool output_solver_details=false, const unsigned int gmres_restart_parameter=30)
std::shared_ptr< Epetra_Operator > preconditioner
const Epetra_MultiVector & trilinos_vector() const
SolverGMRES(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
T max(const T &t, const MPI_Comm &mpi_communicator)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcTrilinosError(int arg1)