16 #include <deal.II/lac/trilinos_precondition.h> 18 #ifdef DEAL_II_WITH_TRILINOS 20 # include <deal.II/lac/vector.h> 21 # include <deal.II/lac/sparse_matrix.h> 22 # include <deal.II/lac/trilinos_sparse_matrix.h> 25 # include <Ifpack_Chebyshev.h> 26 # include <Teuchos_ParameterList.hpp> 27 # include <Teuchos_RCP.hpp> 28 # include <Epetra_MultiVector.h> 30 DEAL_II_NAMESPACE_OPEN
36 #ifdef DEAL_II_WITH_MPI 38 communicator (MPI_COMM_SELF)
47 preconditioner (base.preconditioner),
48 #ifdef DEAL_II_WITH_MPI
49 communicator (base.communicator),
51 vector_distributor (new Epetra_Map(*base.vector_distributor))
59 #ifdef DEAL_II_WITH_MPI 69 #ifdef DEAL_II_WITH_MPI 102 const double min_diagonal,
103 const unsigned int n_sweeps)
106 min_diagonal (min_diagonal),
120 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
123 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 130 Teuchos::ParameterList parameter_list;
131 parameter_list.set (
"relaxation: sweeps", static_cast<int>(additional_data.
n_sweeps));
132 parameter_list.set (
"relaxation: type",
"Jacobi");
133 parameter_list.set (
"relaxation: damping factor", additional_data.
omega);
134 parameter_list.set (
"relaxation: min diagonal value",
137 ierr = ifpack->SetParameters(parameter_list);
140 ierr = ifpack->Initialize();
143 ierr = ifpack->Compute();
153 const double min_diagonal,
154 const unsigned int overlap,
155 const unsigned int n_sweeps)
158 min_diagonal (min_diagonal),
172 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
175 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 182 Teuchos::ParameterList parameter_list;
183 parameter_list.set (
"relaxation: sweeps", static_cast<int>(additional_data.
n_sweeps));
184 parameter_list.set (
"relaxation: type",
"symmetric Gauss-Seidel");
185 parameter_list.set (
"relaxation: damping factor", additional_data.
omega);
186 parameter_list.set (
"relaxation: min diagonal value",
188 parameter_list.set (
"schwarz: combine mode",
"Add");
190 ierr = ifpack->SetParameters(parameter_list);
193 ierr = ifpack->Initialize();
196 ierr = ifpack->Compute();
206 const double min_diagonal,
207 const unsigned int overlap,
208 const unsigned int n_sweeps)
211 min_diagonal (min_diagonal),
225 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
228 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 235 Teuchos::ParameterList parameter_list;
236 parameter_list.set (
"relaxation: sweeps", static_cast<int>(additional_data.
n_sweeps));
237 parameter_list.set (
"relaxation: type",
"Gauss-Seidel");
238 parameter_list.set (
"relaxation: damping factor", additional_data.
omega);
239 parameter_list.set (
"relaxation: min diagonal value",
241 parameter_list.set (
"schwarz: combine mode",
"Add");
243 ierr = ifpack->SetParameters(parameter_list);
246 ierr = ifpack->Initialize();
249 ierr = ifpack->Compute();
259 const std::string &block_creation_type,
261 const double min_diagonal,
262 const unsigned int n_sweeps)
264 block_size(block_size),
265 block_creation_type(block_creation_type),
267 min_diagonal (min_diagonal),
284 (matrix.trilinos_matrix().NumMyRows()==0) ?
285 "point relaxation" :
"block relaxation",
286 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
289 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 296 Teuchos::ParameterList parameter_list;
297 parameter_list.set (
"relaxation: sweeps", static_cast<int>(additional_data.
n_sweeps));
298 parameter_list.set (
"relaxation: type",
"Jacobi");
299 parameter_list.set (
"relaxation: damping factor", additional_data.
omega);
300 parameter_list.set (
"relaxation: min diagonal value",
303 int n_local_parts = (matrix.trilinos_matrix().NumMyRows()+additional_data.
305 parameter_list.set (
"partitioner: local parts", n_local_parts);
307 ierr = ifpack->SetParameters(parameter_list);
310 ierr = ifpack->Initialize();
313 ierr = ifpack->Compute();
323 const std::string &block_creation_type,
325 const double min_diagonal,
326 const unsigned int overlap,
327 const unsigned int n_sweeps)
329 block_size(block_size),
330 block_creation_type(block_creation_type),
332 min_diagonal (min_diagonal),
349 (matrix.trilinos_matrix().NumMyRows()==0) ?
350 "point relaxation" :
"block relaxation",
351 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
354 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 361 Teuchos::ParameterList parameter_list;
362 parameter_list.set (
"relaxation: sweeps", static_cast<int>(additional_data.
n_sweeps));
363 parameter_list.set (
"relaxation: type",
"symmetric Gauss-Seidel");
364 parameter_list.set (
"relaxation: damping factor", additional_data.
omega);
365 parameter_list.set (
"relaxation: min diagonal value",
367 parameter_list.set (
"schwarz: combine mode",
"Add");
369 int n_local_parts = (matrix.trilinos_matrix().NumMyRows()+additional_data.
371 parameter_list.set (
"partitioner: local parts", n_local_parts);
373 ierr = ifpack->SetParameters(parameter_list);
376 ierr = ifpack->Initialize();
379 ierr = ifpack->Compute();
389 const std::string &block_creation_type,
391 const double min_diagonal,
392 const unsigned int overlap,
393 const unsigned int n_sweeps)
395 block_size(block_size),
396 block_creation_type(block_creation_type),
398 min_diagonal (min_diagonal),
415 (matrix.trilinos_matrix().NumMyRows()==0) ?
416 "point relaxation" :
"block relaxation",
417 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
420 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 427 Teuchos::ParameterList parameter_list;
428 parameter_list.set (
"relaxation: sweeps", static_cast<int>(additional_data.
n_sweeps));
429 parameter_list.set (
"relaxation: type",
"Gauss-Seidel");
430 parameter_list.set (
"relaxation: damping factor", additional_data.
omega);
431 parameter_list.set (
"relaxation: min diagonal value",
433 parameter_list.set (
"schwarz: combine mode",
"Add");
435 int n_local_parts = (matrix.trilinos_matrix().NumMyRows()+additional_data.
437 parameter_list.set (
"partitioner: local parts", n_local_parts);
439 ierr = ifpack->SetParameters(parameter_list);
442 ierr = ifpack->Initialize();
445 ierr = ifpack->Compute();
455 const double ic_atol,
456 const double ic_rtol,
457 const unsigned int overlap)
474 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
477 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 484 Teuchos::ParameterList parameter_list;
485 parameter_list.set (
"fact: level-of-fill",(
int)additional_data.
ic_fill);
486 parameter_list.set (
"fact: absolute threshold",additional_data.
ic_atol);
487 parameter_list.set (
"fact: relative threshold",additional_data.
ic_rtol);
488 parameter_list.set (
"schwarz: combine mode",
"Add");
490 ierr = ifpack->SetParameters(parameter_list);
493 ierr = ifpack->Initialize();
496 ierr = ifpack->Compute();
506 const double ilu_atol,
507 const double ilu_rtol,
508 const unsigned int overlap)
525 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
528 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 535 Teuchos::ParameterList parameter_list;
536 parameter_list.set (
"fact: level-of-fill", static_cast<int>(additional_data.
ilu_fill));
537 parameter_list.set (
"fact: absolute threshold", additional_data.
ilu_atol);
538 parameter_list.set (
"fact: relative threshold", additional_data.
ilu_rtol);
539 parameter_list.set (
"schwarz: combine mode",
"Add");
541 ierr = ifpack->SetParameters(parameter_list);
544 ierr = ifpack->Initialize();
547 ierr = ifpack->Compute();
557 const unsigned int ilut_fill,
558 const double ilut_atol,
559 const double ilut_rtol,
560 const unsigned int overlap)
562 ilut_drop (ilut_drop),
563 ilut_fill (ilut_fill),
564 ilut_atol (ilut_atol),
565 ilut_rtol (ilut_rtol),
578 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
581 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 588 Teuchos::ParameterList parameter_list;
589 parameter_list.set (
"fact: drop value",additional_data.
ilut_drop);
590 parameter_list.set (
"fact: level-of-fill",(
int)additional_data.
ilut_fill);
591 parameter_list.set (
"fact: absolute threshold",additional_data.
ilut_atol);
592 parameter_list.set (
"fact: relative threshold",additional_data.
ilut_rtol);
593 parameter_list.set (
"schwarz: combine mode",
"Add");
595 ierr = ifpack->SetParameters(parameter_list);
598 ierr = ifpack->Initialize();
601 ierr = ifpack->Compute();
624 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
627 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 634 Teuchos::ParameterList parameter_list;
635 parameter_list.set (
"schwarz: combine mode",
"Add");
637 ierr = ifpack->SetParameters(parameter_list);
640 ierr = ifpack->Initialize();
643 ierr = ifpack->Compute();
653 const double max_eigenvalue,
654 const double eigenvalue_ratio,
655 const double min_eigenvalue,
656 const double min_diagonal,
657 const bool nonzero_starting)
660 max_eigenvalue (max_eigenvalue),
661 eigenvalue_ratio (eigenvalue_ratio),
662 min_eigenvalue (min_eigenvalue),
663 min_diagonal (min_diagonal),
664 nonzero_starting (nonzero_starting)
674 preconditioner = std::make_shared<Ifpack_Chebyshev> (&matrix.trilinos_matrix());
676 Ifpack_Chebyshev *ifpack =
static_cast<Ifpack_Chebyshev *
> 683 Teuchos::ParameterList parameter_list;
684 parameter_list.set (
"chebyshev: ratio eigenvalue",
686 parameter_list.set (
"chebyshev: min eigenvalue",
688 parameter_list.set (
"chebyshev: max eigenvalue",
690 parameter_list.set (
"chebyshev: degree",
691 (
int)additional_data.
degree);
692 parameter_list.set (
"chebyshev: min diagonal value",
694 parameter_list.set (
"chebyshev: zero starting solution",
697 ierr = ifpack->SetParameters(parameter_list);
700 ierr = ifpack->Initialize();
703 ierr = ifpack->Compute();
733 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
736 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 743 Teuchos::ParameterList parameter_list;
744 parameter_list.set (
"relaxation: sweeps", 1);
745 parameter_list.set (
"relaxation: type",
"Jacobi");
746 parameter_list.set (
"relaxation: damping factor", 1.0);
747 parameter_list.set (
"relaxation: min diagonal value", 0.0);
749 ierr = ifpack->SetParameters(parameter_list);
752 ierr = ifpack->Initialize();
755 ierr = ifpack->Compute();
775 const ::Vector<double> &src)
const 782 const ::Vector<double> &src)
const 802 DEAL_II_NAMESPACE_CLOSE
804 #endif // DEAL_II_WITH_TRILINOS
static ::ExceptionBase & ExcTrilinosError(int arg1)
AdditionalData(const unsigned int block_size=1, const std::string &block_creation_type="linear", const double omega=1, const double min_diagonal=0, const unsigned int overlap=0, const unsigned int n_sweeps=1)
IndexSet locally_owned_range_indices() const
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData(const double omega=1, const double min_diagonal=0, const unsigned int overlap=0, const unsigned int n_sweeps=1)
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
std::shared_ptr< Epetra_Map > vector_distributor
AdditionalData(const unsigned int ilu_fill=0, const double ilu_atol=0., const double ilu_rtol=1., const unsigned int overlap=0)
#define AssertThrow(cond, exc)
AdditionalData(const unsigned int block_size=1, const std::string &block_creation_type="linear", const double omega=1, const double min_diagonal=0, const unsigned int n_sweeps=1)
AdditionalData(const unsigned int overlap=0)
std::string block_creation_type
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData(const unsigned int degree=1, const double max_eigenvalue=10., const double eigenvalue_ratio=30., const double min_eigenvalue=1., const double min_diagonal=1e-12, const bool nonzero_starting=false)
AdditionalData(const double ilut_drop=0., const unsigned int ilut_fill=0, const double ilut_atol=0., const double ilut_rtol=1., const unsigned int overlap=0)
AdditionalData(const double omega=1, const double min_diagonal=0, const unsigned int overlap=0, const unsigned int n_sweeps=1)
Epetra_Operator & trilinos_operator() const
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
IndexSet locally_owned_domain_indices() const
std::string block_creation_type
AdditionalData(const unsigned int ic_fill=0, const double ic_atol=0., const double ic_rtol=1., const unsigned int overlap=0)
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
std::shared_ptr< Epetra_Operator > preconditioner
AdditionalData(const unsigned int block_size=1, const std::string &block_creation_type="linear", const double omega=1, const double min_diagonal=0, const unsigned int overlap=0, const unsigned int n_sweeps=1)
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void vmult(MPI::Vector &dst, const MPI::Vector &src) const
void Tvmult(MPI::Vector &dst, const MPI::Vector &src) const
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
Epetra_MpiComm communicator
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
MPI_Comm get_mpi_communicator() const
AdditionalData(const double omega=1, const double min_diagonal=0, const unsigned int n_sweeps=1)
std::string block_creation_type