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> 24 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
26 # include <Ifpack_Chebyshev.h> 27 # include <Teuchos_ParameterList.hpp> 28 # include <Teuchos_RCP.hpp> 29 # include <Epetra_MultiVector.h> 30 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
32 DEAL_II_NAMESPACE_OPEN
38 #ifdef DEAL_II_WITH_MPI 40 communicator (MPI_COMM_SELF)
49 preconditioner (base.preconditioner),
50 #ifdef DEAL_II_WITH_MPI
51 communicator (base.communicator),
53 vector_distributor (new Epetra_Map(*base.vector_distributor))
66 #ifdef DEAL_II_WITH_MPI 76 #ifdef DEAL_II_WITH_MPI 109 const double min_diagonal,
110 const unsigned int n_sweeps)
113 min_diagonal (min_diagonal),
127 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
130 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 137 Teuchos::ParameterList parameter_list;
138 parameter_list.set (
"relaxation: sweeps", static_cast<int>(additional_data.
n_sweeps));
139 parameter_list.set (
"relaxation: type",
"Jacobi");
140 parameter_list.set (
"relaxation: damping factor", additional_data.
omega);
141 parameter_list.set (
"relaxation: min diagonal value",
144 ierr = ifpack->SetParameters(parameter_list);
147 ierr = ifpack->Initialize();
150 ierr = ifpack->Compute();
160 const double min_diagonal,
161 const unsigned int overlap,
162 const unsigned int n_sweeps)
165 min_diagonal (min_diagonal),
179 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
182 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 189 Teuchos::ParameterList parameter_list;
190 parameter_list.set (
"relaxation: sweeps", static_cast<int>(additional_data.
n_sweeps));
191 parameter_list.set (
"relaxation: type",
"symmetric Gauss-Seidel");
192 parameter_list.set (
"relaxation: damping factor", additional_data.
omega);
193 parameter_list.set (
"relaxation: min diagonal value",
195 parameter_list.set (
"schwarz: combine mode",
"Add");
197 ierr = ifpack->SetParameters(parameter_list);
200 ierr = ifpack->Initialize();
203 ierr = ifpack->Compute();
213 const double min_diagonal,
214 const unsigned int overlap,
215 const unsigned int n_sweeps)
218 min_diagonal (min_diagonal),
232 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
235 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 242 Teuchos::ParameterList parameter_list;
243 parameter_list.set (
"relaxation: sweeps", static_cast<int>(additional_data.
n_sweeps));
244 parameter_list.set (
"relaxation: type",
"Gauss-Seidel");
245 parameter_list.set (
"relaxation: damping factor", additional_data.
omega);
246 parameter_list.set (
"relaxation: min diagonal value",
248 parameter_list.set (
"schwarz: combine mode",
"Add");
250 ierr = ifpack->SetParameters(parameter_list);
253 ierr = ifpack->Initialize();
256 ierr = ifpack->Compute();
266 const std::string &block_creation_type,
268 const double min_diagonal,
269 const unsigned int n_sweeps)
271 block_size(block_size),
272 block_creation_type(block_creation_type),
274 min_diagonal (min_diagonal),
288 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
291 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 298 Teuchos::ParameterList parameter_list;
299 parameter_list.set (
"relaxation: sweeps", static_cast<int>(additional_data.
n_sweeps));
300 parameter_list.set (
"relaxation: type",
"Jacobi");
301 parameter_list.set (
"relaxation: damping factor", additional_data.
omega);
302 parameter_list.set (
"relaxation: min diagonal value",
305 int n_local_parts = (matrix.trilinos_matrix().NumMyRows()+additional_data.
307 parameter_list.set (
"partitioner: local parts", n_local_parts);
309 ierr = ifpack->SetParameters(parameter_list);
312 ierr = ifpack->Initialize();
315 ierr = ifpack->Compute();
325 const std::string &block_creation_type,
327 const double min_diagonal,
328 const unsigned int overlap,
329 const unsigned int n_sweeps)
331 block_size(block_size),
332 block_creation_type(block_creation_type),
334 min_diagonal (min_diagonal),
348 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
351 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 358 Teuchos::ParameterList parameter_list;
359 parameter_list.set (
"relaxation: sweeps", static_cast<int>(additional_data.
n_sweeps));
360 parameter_list.set (
"relaxation: type",
"symmetric Gauss-Seidel");
361 parameter_list.set (
"relaxation: damping factor", additional_data.
omega);
362 parameter_list.set (
"relaxation: min diagonal value",
364 parameter_list.set (
"schwarz: combine mode",
"Add");
366 int n_local_parts = (matrix.trilinos_matrix().NumMyRows()+additional_data.
368 parameter_list.set (
"partitioner: local parts", n_local_parts);
370 ierr = ifpack->SetParameters(parameter_list);
373 ierr = ifpack->Initialize();
376 ierr = ifpack->Compute();
386 const std::string &block_creation_type,
388 const double min_diagonal,
389 const unsigned int overlap,
390 const unsigned int n_sweeps)
392 block_size(block_size),
393 block_creation_type(block_creation_type),
395 min_diagonal (min_diagonal),
409 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
412 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 419 Teuchos::ParameterList parameter_list;
420 parameter_list.set (
"relaxation: sweeps", static_cast<int>(additional_data.
n_sweeps));
421 parameter_list.set (
"relaxation: type",
"Gauss-Seidel");
422 parameter_list.set (
"relaxation: damping factor", additional_data.
omega);
423 parameter_list.set (
"relaxation: min diagonal value",
425 parameter_list.set (
"schwarz: combine mode",
"Add");
427 int n_local_parts = (matrix.trilinos_matrix().NumMyRows()+additional_data.
429 parameter_list.set (
"partitioner: local parts", n_local_parts);
431 ierr = ifpack->SetParameters(parameter_list);
434 ierr = ifpack->Initialize();
437 ierr = ifpack->Compute();
447 const double ic_atol,
448 const double ic_rtol,
449 const unsigned int overlap)
466 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
469 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 476 Teuchos::ParameterList parameter_list;
477 parameter_list.set (
"fact: level-of-fill",(
int)additional_data.
ic_fill);
478 parameter_list.set (
"fact: absolute threshold",additional_data.
ic_atol);
479 parameter_list.set (
"fact: relative threshold",additional_data.
ic_rtol);
480 parameter_list.set (
"schwarz: combine mode",
"Add");
482 ierr = ifpack->SetParameters(parameter_list);
485 ierr = ifpack->Initialize();
488 ierr = ifpack->Compute();
498 const double ilu_atol,
499 const double ilu_rtol,
500 const unsigned int overlap)
517 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
520 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 527 Teuchos::ParameterList parameter_list;
528 parameter_list.set (
"fact: level-of-fill", static_cast<int>(additional_data.
ilu_fill));
529 parameter_list.set (
"fact: absolute threshold", additional_data.
ilu_atol);
530 parameter_list.set (
"fact: relative threshold", additional_data.
ilu_rtol);
531 parameter_list.set (
"schwarz: combine mode",
"Add");
533 ierr = ifpack->SetParameters(parameter_list);
536 ierr = ifpack->Initialize();
539 ierr = ifpack->Compute();
549 const unsigned int ilut_fill,
550 const double ilut_atol,
551 const double ilut_rtol,
552 const unsigned int overlap)
554 ilut_drop (ilut_drop),
555 ilut_fill (ilut_fill),
556 ilut_atol (ilut_atol),
557 ilut_rtol (ilut_rtol),
570 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
573 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 580 Teuchos::ParameterList parameter_list;
581 parameter_list.set (
"fact: drop value",additional_data.
ilut_drop);
582 parameter_list.set (
"fact: level-of-fill",(
int)additional_data.
ilut_fill);
583 parameter_list.set (
"fact: absolute threshold",additional_data.
ilut_atol);
584 parameter_list.set (
"fact: relative threshold",additional_data.
ilut_rtol);
585 parameter_list.set (
"schwarz: combine mode",
"Add");
587 ierr = ifpack->SetParameters(parameter_list);
590 ierr = ifpack->Initialize();
593 ierr = ifpack->Compute();
616 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
619 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 626 Teuchos::ParameterList parameter_list;
627 parameter_list.set (
"schwarz: combine mode",
"Add");
629 ierr = ifpack->SetParameters(parameter_list);
632 ierr = ifpack->Initialize();
635 ierr = ifpack->Compute();
645 const double max_eigenvalue,
646 const double eigenvalue_ratio,
647 const double min_eigenvalue,
648 const double min_diagonal,
649 const bool nonzero_starting)
652 max_eigenvalue (max_eigenvalue),
653 eigenvalue_ratio (eigenvalue_ratio),
654 min_eigenvalue (min_eigenvalue),
655 min_diagonal (min_diagonal),
656 nonzero_starting (nonzero_starting)
666 preconditioner.reset (
new Ifpack_Chebyshev (&matrix.trilinos_matrix()));
668 Ifpack_Chebyshev *ifpack =
static_cast<Ifpack_Chebyshev *
> 675 Teuchos::ParameterList parameter_list;
676 parameter_list.set (
"chebyshev: ratio eigenvalue",
678 parameter_list.set (
"chebyshev: min eigenvalue",
680 parameter_list.set (
"chebyshev: max eigenvalue",
682 parameter_list.set (
"chebyshev: degree",
683 (
int)additional_data.
degree);
684 parameter_list.set (
"chebyshev: min diagonal value",
686 parameter_list.set (
"chebyshev: zero starting solution",
689 ierr = ifpack->SetParameters(parameter_list);
692 ierr = ifpack->Initialize();
695 ierr = ifpack->Compute();
725 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
728 Ifpack_Preconditioner *ifpack =
static_cast<Ifpack_Preconditioner *
> 735 Teuchos::ParameterList parameter_list;
736 parameter_list.set (
"relaxation: sweeps", 1);
737 parameter_list.set (
"relaxation: type",
"Jacobi");
738 parameter_list.set (
"relaxation: damping factor", 1.0);
739 parameter_list.set (
"relaxation: min diagonal value", 0.0);
741 ierr = ifpack->SetParameters(parameter_list);
744 ierr = ifpack->Initialize();
747 ierr = ifpack->Compute();
767 const ::Vector<double> &src)
const 774 const ::Vector<double> &src)
const 794 DEAL_II_NAMESPACE_CLOSE
796 #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())
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_cxx11::shared_ptr< Epetra_Map > vector_distributor
std::string block_creation_type
static ::ExceptionBase & ExcMessage(std::string arg1)
std_cxx11::shared_ptr< Epetra_Operator > preconditioner
#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 vmult(VectorBase &dst, const VectorBase &src) const
void Tvmult(VectorBase &dst, const VectorBase &src) const
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())
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 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