18#ifdef DEAL_II_WITH_TRILINOS
24# include <Epetra_MultiVector.h>
26# include <Ifpack_Chebyshev.h>
27# include <Teuchos_ParameterList.hpp>
28# include <Teuchos_RCP.hpp>
35 : communicator(MPI_COMM_SELF)
42 , preconditioner(base.preconditioner)
43 , communicator(base.communicator)
44 , vector_distributor(new Epetra_Map(*base.vector_distributor))
69 ExcMessage(
"Trying to dereference a null pointer."));
91 const double min_diagonal,
92 const unsigned int n_sweeps)
94 , min_diagonal(min_diagonal)
107 Ifpack().Create(
"point relaxation",
108 const_cast<Epetra_CrsMatrix *
>(&matrix.trilinos_matrix()),
111 Ifpack_Preconditioner *ifpack =
119 Teuchos::ParameterList parameter_list;
120 parameter_list.set(
"relaxation: sweeps",
121 static_cast<int>(additional_data.
n_sweeps));
122 parameter_list.set(
"relaxation: type",
"Jacobi");
123 parameter_list.set(
"relaxation: damping factor", additional_data.
omega);
124 parameter_list.set(
"relaxation: min diagonal value",
127 ierr = ifpack->SetParameters(parameter_list);
130 ierr = ifpack->Initialize();
133 ierr = ifpack->Compute();
142 const double min_diagonal,
143 const unsigned int overlap,
144 const unsigned int n_sweeps)
146 , min_diagonal(min_diagonal)
159 Ifpack().Create(
"point relaxation",
160 const_cast<Epetra_CrsMatrix *
>(&matrix.trilinos_matrix()),
163 Ifpack_Preconditioner *ifpack =
171 Teuchos::ParameterList parameter_list;
172 parameter_list.set(
"relaxation: sweeps",
173 static_cast<int>(additional_data.
n_sweeps));
174 parameter_list.set(
"relaxation: type",
"symmetric Gauss-Seidel");
175 parameter_list.set(
"relaxation: damping factor", additional_data.
omega);
176 parameter_list.set(
"relaxation: min diagonal value",
178 parameter_list.set(
"schwarz: combine mode",
"Add");
180 ierr = ifpack->SetParameters(parameter_list);
183 ierr = ifpack->Initialize();
186 ierr = ifpack->Compute();
195 const double min_diagonal,
196 const unsigned int overlap,
197 const unsigned int n_sweeps)
199 , min_diagonal(min_diagonal)
212 Ifpack().Create(
"point relaxation",
213 const_cast<Epetra_CrsMatrix *
>(&matrix.trilinos_matrix()),
216 Ifpack_Preconditioner *ifpack =
224 Teuchos::ParameterList parameter_list;
225 parameter_list.set(
"relaxation: sweeps",
226 static_cast<int>(additional_data.
n_sweeps));
227 parameter_list.set(
"relaxation: type",
"Gauss-Seidel");
228 parameter_list.set(
"relaxation: damping factor", additional_data.
omega);
229 parameter_list.set(
"relaxation: min diagonal value",
231 parameter_list.set(
"schwarz: combine mode",
"Add");
233 ierr = ifpack->SetParameters(parameter_list);
236 ierr = ifpack->Initialize();
239 ierr = ifpack->Compute();
248 const unsigned int block_size,
249 const std::string &block_creation_type,
251 const double min_diagonal,
252 const unsigned int n_sweeps)
253 : block_size(block_size)
254 , block_creation_type(block_creation_type)
256 , min_diagonal(min_diagonal)
272 Ifpack().Create((matrix.trilinos_matrix().NumMyRows() == 0) ?
275 const_cast<Epetra_CrsMatrix *
>(&matrix.trilinos_matrix()),
278 Ifpack_Preconditioner *ifpack =
286 Teuchos::ParameterList parameter_list;
287 parameter_list.set(
"relaxation: sweeps",
288 static_cast<int>(additional_data.
n_sweeps));
289 parameter_list.set(
"relaxation: type",
"Jacobi");
290 parameter_list.set(
"relaxation: damping factor", additional_data.
omega);
291 parameter_list.set(
"relaxation: min diagonal value",
293 parameter_list.set(
"partitioner: type",
296 (matrix.trilinos_matrix().NumMyRows() + additional_data.
block_size - 1) /
298 parameter_list.set(
"partitioner: local parts", n_local_parts);
300 ierr = ifpack->SetParameters(parameter_list);
303 ierr = ifpack->Initialize();
306 ierr = ifpack->Compute();
315 const unsigned int block_size,
316 const std::string &block_creation_type,
318 const double min_diagonal,
319 const unsigned int overlap,
320 const unsigned int n_sweeps)
321 : block_size(block_size)
322 , block_creation_type(block_creation_type)
324 , min_diagonal(min_diagonal)
340 Ifpack().Create((matrix.trilinos_matrix().NumMyRows() == 0) ?
343 const_cast<Epetra_CrsMatrix *
>(&matrix.trilinos_matrix()),
346 Ifpack_Preconditioner *ifpack =
354 Teuchos::ParameterList parameter_list;
355 parameter_list.set(
"relaxation: sweeps",
356 static_cast<int>(additional_data.
n_sweeps));
357 parameter_list.set(
"relaxation: type",
"symmetric Gauss-Seidel");
358 parameter_list.set(
"relaxation: damping factor", additional_data.
omega);
359 parameter_list.set(
"relaxation: min diagonal value",
361 parameter_list.set(
"schwarz: combine mode",
"Add");
362 parameter_list.set(
"partitioner: type",
365 (matrix.trilinos_matrix().NumMyRows() + additional_data.
block_size - 1) /
367 parameter_list.set(
"partitioner: local parts", n_local_parts);
369 ierr = ifpack->SetParameters(parameter_list);
372 ierr = ifpack->Initialize();
375 ierr = ifpack->Compute();
384 const unsigned int block_size,
385 const std::string &block_creation_type,
387 const double min_diagonal,
388 const unsigned int overlap,
389 const unsigned int n_sweeps)
390 : block_size(block_size)
391 , block_creation_type(block_creation_type)
393 , min_diagonal(min_diagonal)
409 Ifpack().Create((matrix.trilinos_matrix().NumMyRows() == 0) ?
412 const_cast<Epetra_CrsMatrix *
>(&matrix.trilinos_matrix()),
415 Ifpack_Preconditioner *ifpack =
423 Teuchos::ParameterList parameter_list;
424 parameter_list.set(
"relaxation: sweeps",
425 static_cast<int>(additional_data.
n_sweeps));
426 parameter_list.set(
"relaxation: type",
"Gauss-Seidel");
427 parameter_list.set(
"relaxation: damping factor", additional_data.
omega);
428 parameter_list.set(
"relaxation: min diagonal value",
430 parameter_list.set(
"schwarz: combine mode",
"Add");
431 parameter_list.set(
"partitioner: type",
434 (matrix.trilinos_matrix().NumMyRows() + additional_data.
block_size - 1) /
436 parameter_list.set(
"partitioner: local parts", n_local_parts);
438 ierr = ifpack->SetParameters(parameter_list);
441 ierr = ifpack->Initialize();
444 ierr = ifpack->Compute();
453 const double ic_atol,
454 const double ic_rtol,
455 const unsigned int overlap)
470 Ifpack().Create(
"IC",
471 const_cast<Epetra_CrsMatrix *
>(&matrix.trilinos_matrix()),
474 Ifpack_Preconditioner *ifpack =
482 Teuchos::ParameterList parameter_list;
483 parameter_list.set(
"fact: level-of-fill", additional_data.
ic_fill);
484 parameter_list.set(
"fact: absolute threshold", additional_data.
ic_atol);
485 parameter_list.set(
"fact: relative threshold", additional_data.
ic_rtol);
486 parameter_list.set(
"schwarz: combine mode",
"Add");
488 ierr = ifpack->SetParameters(parameter_list);
491 ierr = ifpack->Initialize();
494 ierr = ifpack->Compute();
503 const double ilu_atol,
504 const double ilu_rtol,
505 const unsigned int overlap)
520 Ifpack().Create(
"ILU",
521 const_cast<Epetra_CrsMatrix *
>(&matrix.trilinos_matrix()),
524 Ifpack_Preconditioner *ifpack =
532 Teuchos::ParameterList parameter_list;
533 parameter_list.set(
"fact: level-of-fill",
534 static_cast<int>(additional_data.
ilu_fill));
535 parameter_list.set(
"fact: absolute threshold", additional_data.
ilu_atol);
536 parameter_list.set(
"fact: relative threshold", additional_data.
ilu_rtol);
537 parameter_list.set(
"schwarz: combine mode",
"Add");
539 ierr = ifpack->SetParameters(parameter_list);
542 ierr = ifpack->Initialize();
545 ierr = ifpack->Compute();
554 const unsigned int ilut_fill,
555 const double ilut_atol,
556 const double ilut_rtol,
557 const unsigned int overlap)
558 : ilut_drop(ilut_drop)
559 , ilut_fill(ilut_fill)
560 , ilut_atol(ilut_atol)
561 , ilut_rtol(ilut_rtol)
573 Ifpack().Create(
"ILUT",
574 const_cast<Epetra_CrsMatrix *
>(&matrix.trilinos_matrix()),
577 Ifpack_Preconditioner *ifpack =
585 Teuchos::ParameterList parameter_list;
586 parameter_list.set(
"fact: drop value", additional_data.
ilut_drop);
587 parameter_list.set(
"fact: level-of-fill",
588 static_cast<int>(additional_data.
ilut_fill));
589 parameter_list.set(
"fact: absolute threshold", additional_data.
ilut_atol);
590 parameter_list.set(
"fact: relative threshold", additional_data.
ilut_rtol);
591 parameter_list.set(
"schwarz: combine mode",
"Add");
593 ierr = ifpack->SetParameters(parameter_list);
596 ierr = ifpack->Initialize();
599 ierr = ifpack->Compute();
608 const unsigned int overlap)
620 Ifpack().Create(
"Amesos",
621 const_cast<Epetra_CrsMatrix *
>(&matrix.trilinos_matrix()),
624 Ifpack_Preconditioner *ifpack =
632 Teuchos::ParameterList parameter_list;
633 parameter_list.set(
"schwarz: combine mode",
"Add");
635 ierr = ifpack->SetParameters(parameter_list);
638 ierr = ifpack->Initialize();
641 ierr = ifpack->Compute();
650 const unsigned int degree,
651 const double max_eigenvalue,
652 const double eigenvalue_ratio,
653 const double min_eigenvalue,
654 const double min_diagonal,
655 const bool nonzero_starting)
657 , max_eigenvalue(max_eigenvalue)
658 , eigenvalue_ratio(eigenvalue_ratio)
659 , min_eigenvalue(min_eigenvalue)
660 , min_diagonal(min_diagonal)
661 , nonzero_starting(nonzero_starting)
671 Teuchos::rcp(
new Ifpack_Chebyshev(&matrix.trilinos_matrix()));
673 Ifpack_Chebyshev *ifpack =
681 Teuchos::ParameterList parameter_list;
682 parameter_list.set(
"chebyshev: ratio eigenvalue",
684 parameter_list.set(
"chebyshev: min eigenvalue",
686 parameter_list.set(
"chebyshev: max eigenvalue",
688 parameter_list.set(
"chebyshev: degree",
689 static_cast<int>(additional_data.
degree));
690 parameter_list.set(
"chebyshev: min diagonal value",
692 parameter_list.set(
"chebyshev: zero starting solution",
695 ierr = ifpack->SetParameters(parameter_list);
698 ierr = ifpack->Initialize();
701 ierr = ifpack->Compute();
725 Ifpack().Create(
"point relaxation",
726 const_cast<Epetra_CrsMatrix *
>(&matrix.trilinos_matrix()),
729 Ifpack_Preconditioner *ifpack =
737 Teuchos::ParameterList parameter_list;
738 parameter_list.set(
"relaxation: sweeps", 1);
739 parameter_list.set(
"relaxation: type",
"Jacobi");
740 parameter_list.set(
"relaxation: damping factor", 1.0);
741 parameter_list.set(
"relaxation: min diagonal value", 0.0);
743 ierr = ifpack->SetParameters(parameter_list);
746 ierr = ifpack->Initialize();
749 ierr = ifpack->Compute();
767 const ::Vector<double> &src)
const
774 const ::Vector<double> &src)
const
IndexSet locally_owned_range_indices() const
Epetra_MpiComm communicator
Epetra_Operator & trilinos_operator() const
std::shared_ptr< Epetra_Map > vector_distributor
Teuchos::RCP< Epetra_Operator > preconditioner
MPI_Comm get_mpi_communicator() const
IndexSet locally_owned_domain_indices() 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())
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())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void vmult(MPI::Vector &dst, const MPI::Vector &src) const override
void Tvmult(MPI::Vector &dst, const MPI::Vector &src) const override
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())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
std::string block_creation_type
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)
std::string block_creation_type
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)
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)
std::string block_creation_type
AdditionalData(const unsigned int overlap=0)
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 unsigned int ic_fill=0, const double ic_atol=0., const double ic_rtol=1., const unsigned int overlap=0)
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 unsigned int ilu_fill=0, const double ilu_atol=0., const double ilu_rtol=1., const unsigned int overlap=0)
AdditionalData(const double omega=1, const double min_diagonal=0, const unsigned int n_sweeps=1)
AdditionalData(const double omega=1, const double min_diagonal=0, const unsigned int overlap=0, const unsigned int n_sweeps=1)
AdditionalData(const double omega=1, const double min_diagonal=0, const unsigned int overlap=0, const unsigned int n_sweeps=1)