Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
trilinos_precondition.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2019 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/lac/trilinos_precondition.h>
17 
18 #ifdef DEAL_II_WITH_TRILINOS
19 
20 # include <deal.II/lac/sparse_matrix.h>
21 # include <deal.II/lac/trilinos_sparse_matrix.h>
22 # include <deal.II/lac/vector.h>
23 
24 # include <Epetra_MultiVector.h>
25 # include <Ifpack.h>
26 # include <Ifpack_Chebyshev.h>
27 # include <Teuchos_ParameterList.hpp>
28 # include <Teuchos_RCP.hpp>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 namespace TrilinosWrappers
33 {
35 # ifdef DEAL_II_WITH_MPI
36  : communicator(MPI_COMM_SELF)
37 # endif
38  {}
39 
40 
41 
43  : Subscriptor()
44  , preconditioner(base.preconditioner)
45  ,
46 # ifdef DEAL_II_WITH_MPI
47  communicator(base.communicator)
48  ,
49 # endif
50  vector_distributor(new Epetra_Map(*base.vector_distributor))
51  {}
52 
53 
54 
55  void
57  {
58  preconditioner.reset();
59 # ifdef DEAL_II_WITH_MPI
60  communicator = MPI_COMM_SELF;
61 # endif
62  vector_distributor.reset();
63  }
64 
65 
66  MPI_Comm
68  {
69 # ifdef DEAL_II_WITH_MPI
70  return communicator.Comm();
71 # else
72  return MPI_COMM_SELF;
73 # endif
74  }
75 
76 
77  Epetra_Operator &
79  {
80  AssertThrow(!preconditioner.is_null(),
81  ExcMessage("Trying to dereference a null pointer."));
82  return (*preconditioner);
83  }
84 
85 
86  IndexSet
88  {
89  return IndexSet(preconditioner->OperatorDomainMap());
90  }
91 
92 
93  IndexSet
95  {
96  return IndexSet(preconditioner->OperatorRangeMap());
97  }
98 
99  /* -------------------------- PreconditionJacobi -------------------------- */
100 
102  const double omega,
103  const double min_diagonal,
104  const unsigned int n_sweeps)
105  : omega(omega)
106  , min_diagonal(min_diagonal)
107  , n_sweeps(n_sweeps)
108  {}
109 
110 
111 
112  void
114  const AdditionalData &additional_data)
115  {
116  // release memory before reallocation
117  preconditioner.reset();
118  preconditioner.reset(
119  Ifpack().Create("point relaxation",
120  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
121  0));
122 
123  Ifpack_Preconditioner *ifpack =
124  dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
125  Assert(ifpack != nullptr,
126  ExcMessage("Trilinos could not create this "
127  "preconditioner"));
128 
129  int ierr;
130 
131  Teuchos::ParameterList parameter_list;
132  parameter_list.set("relaxation: sweeps",
133  static_cast<int>(additional_data.n_sweeps));
134  parameter_list.set("relaxation: type", "Jacobi");
135  parameter_list.set("relaxation: damping factor", additional_data.omega);
136  parameter_list.set("relaxation: min diagonal value",
137  additional_data.min_diagonal);
138 
139  ierr = ifpack->SetParameters(parameter_list);
140  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
141 
142  ierr = ifpack->Initialize();
143  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
144 
145  ierr = ifpack->Compute();
146  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
147  }
148 
149 
150 
151  /* -------------------------- PreconditionSSOR -------------------------- */
152 
154  const double min_diagonal,
155  const unsigned int overlap,
156  const unsigned int n_sweeps)
157  : omega(omega)
158  , min_diagonal(min_diagonal)
159  , overlap(overlap)
160  , n_sweeps(n_sweeps)
161  {}
162 
163 
164 
165  void
167  const AdditionalData &additional_data)
168  {
169  preconditioner.reset();
170  preconditioner.reset(
171  Ifpack().Create("point relaxation",
172  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
173  additional_data.overlap));
174 
175  Ifpack_Preconditioner *ifpack =
176  dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
177  Assert(ifpack != nullptr,
178  ExcMessage("Trilinos could not create this "
179  "preconditioner"));
180 
181  int ierr;
182 
183  Teuchos::ParameterList parameter_list;
184  parameter_list.set("relaxation: sweeps",
185  static_cast<int>(additional_data.n_sweeps));
186  parameter_list.set("relaxation: type", "symmetric Gauss-Seidel");
187  parameter_list.set("relaxation: damping factor", additional_data.omega);
188  parameter_list.set("relaxation: min diagonal value",
189  additional_data.min_diagonal);
190  parameter_list.set("schwarz: combine mode", "Add");
191 
192  ierr = ifpack->SetParameters(parameter_list);
193  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
194 
195  ierr = ifpack->Initialize();
196  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
197 
198  ierr = ifpack->Compute();
199  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
200  }
201 
202 
203 
204  /* -------------------------- PreconditionSOR -------------------------- */
205 
207  const double min_diagonal,
208  const unsigned int overlap,
209  const unsigned int n_sweeps)
210  : omega(omega)
211  , min_diagonal(min_diagonal)
212  , overlap(overlap)
213  , n_sweeps(n_sweeps)
214  {}
215 
216 
217 
218  void
220  const AdditionalData &additional_data)
221  {
222  preconditioner.reset();
223  preconditioner.reset(
224  Ifpack().Create("point relaxation",
225  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
226  additional_data.overlap));
227 
228  Ifpack_Preconditioner *ifpack =
229  dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
230  Assert(ifpack != nullptr,
231  ExcMessage("Trilinos could not create this "
232  "preconditioner"));
233 
234  int ierr;
235 
236  Teuchos::ParameterList parameter_list;
237  parameter_list.set("relaxation: sweeps",
238  static_cast<int>(additional_data.n_sweeps));
239  parameter_list.set("relaxation: type", "Gauss-Seidel");
240  parameter_list.set("relaxation: damping factor", additional_data.omega);
241  parameter_list.set("relaxation: min diagonal value",
242  additional_data.min_diagonal);
243  parameter_list.set("schwarz: combine mode", "Add");
244 
245  ierr = ifpack->SetParameters(parameter_list);
246  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
247 
248  ierr = ifpack->Initialize();
249  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
250 
251  ierr = ifpack->Compute();
252  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
253  }
254 
255 
256 
257  /* ----------------------- PreconditionBlockJacobi ---------------------- */
258 
260  const unsigned int block_size,
261  const std::string &block_creation_type,
262  const double omega,
263  const double min_diagonal,
264  const unsigned int n_sweeps)
265  : block_size(block_size)
266  , block_creation_type(block_creation_type)
267  , omega(omega)
268  , min_diagonal(min_diagonal)
269  , n_sweeps(n_sweeps)
270  {}
271 
272 
273 
274  void
276  const AdditionalData &additional_data)
277  {
278  // release memory before reallocation
279  preconditioner.reset();
280 
281  // Block relaxation setup fails if we have no locally owned rows. As a
282  // work-around we just pretend to use point relaxation on those processors:
283  preconditioner.reset(
284  Ifpack().Create((matrix.trilinos_matrix().NumMyRows() == 0) ?
285  "point relaxation" :
286  "block relaxation",
287  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
288  0));
289 
290  Ifpack_Preconditioner *ifpack =
291  dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
292  Assert(ifpack != nullptr,
293  ExcMessage("Trilinos could not create this "
294  "preconditioner"));
295 
296  int ierr;
297 
298  Teuchos::ParameterList parameter_list;
299  parameter_list.set("relaxation: sweeps",
300  static_cast<int>(additional_data.n_sweeps));
301  parameter_list.set("relaxation: type", "Jacobi");
302  parameter_list.set("relaxation: damping factor", additional_data.omega);
303  parameter_list.set("relaxation: min diagonal value",
304  additional_data.min_diagonal);
305  parameter_list.set("partitioner: type",
306  additional_data.block_creation_type);
307  int n_local_parts =
308  (matrix.trilinos_matrix().NumMyRows() + additional_data.block_size - 1) /
309  additional_data.block_size;
310  parameter_list.set("partitioner: local parts", n_local_parts);
311 
312  ierr = ifpack->SetParameters(parameter_list);
313  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
314 
315  ierr = ifpack->Initialize();
316  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
317 
318  ierr = ifpack->Compute();
319  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
320  }
321 
322 
323 
324  /* ----------------------- PreconditionBlockSSOR ------------------------ */
325 
327  const unsigned int block_size,
328  const std::string &block_creation_type,
329  const double omega,
330  const double min_diagonal,
331  const unsigned int overlap,
332  const unsigned int n_sweeps)
333  : block_size(block_size)
334  , block_creation_type(block_creation_type)
335  , omega(omega)
336  , min_diagonal(min_diagonal)
337  , overlap(overlap)
338  , n_sweeps(n_sweeps)
339  {}
340 
341 
342 
343  void
345  const AdditionalData &additional_data)
346  {
347  preconditioner.reset();
348 
349  // Block relaxation setup fails if we have no locally owned rows. As a
350  // work-around we just pretend to use point relaxation on those processors:
351  preconditioner.reset(
352  Ifpack().Create((matrix.trilinos_matrix().NumMyRows() == 0) ?
353  "point relaxation" :
354  "block relaxation",
355  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
356  additional_data.overlap));
357 
358  Ifpack_Preconditioner *ifpack =
359  dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
360  Assert(ifpack != nullptr,
361  ExcMessage("Trilinos could not create this "
362  "preconditioner"));
363 
364  int ierr;
365 
366  Teuchos::ParameterList parameter_list;
367  parameter_list.set("relaxation: sweeps",
368  static_cast<int>(additional_data.n_sweeps));
369  parameter_list.set("relaxation: type", "symmetric Gauss-Seidel");
370  parameter_list.set("relaxation: damping factor", additional_data.omega);
371  parameter_list.set("relaxation: min diagonal value",
372  additional_data.min_diagonal);
373  parameter_list.set("schwarz: combine mode", "Add");
374  parameter_list.set("partitioner: type",
375  additional_data.block_creation_type);
376  int n_local_parts =
377  (matrix.trilinos_matrix().NumMyRows() + additional_data.block_size - 1) /
378  additional_data.block_size;
379  parameter_list.set("partitioner: local parts", n_local_parts);
380 
381  ierr = ifpack->SetParameters(parameter_list);
382  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
383 
384  ierr = ifpack->Initialize();
385  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
386 
387  ierr = ifpack->Compute();
388  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
389  }
390 
391 
392 
393  /* ------------------------ PreconditionBlockSOR ------------------------ */
394 
396  const unsigned int block_size,
397  const std::string &block_creation_type,
398  const double omega,
399  const double min_diagonal,
400  const unsigned int overlap,
401  const unsigned int n_sweeps)
402  : block_size(block_size)
403  , block_creation_type(block_creation_type)
404  , omega(omega)
405  , min_diagonal(min_diagonal)
406  , overlap(overlap)
407  , n_sweeps(n_sweeps)
408  {}
409 
410 
411 
412  void
414  const AdditionalData &additional_data)
415  {
416  preconditioner.reset();
417 
418  // Block relaxation setup fails if we have no locally owned rows. As a
419  // work-around we just pretend to use point relaxation on those processors:
420  preconditioner.reset(
421  Ifpack().Create((matrix.trilinos_matrix().NumMyRows() == 0) ?
422  "point relaxation" :
423  "block relaxation",
424  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
425  additional_data.overlap));
426 
427  Ifpack_Preconditioner *ifpack =
428  dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
429  Assert(ifpack != nullptr,
430  ExcMessage("Trilinos could not create this "
431  "preconditioner"));
432 
433  int ierr;
434 
435  Teuchos::ParameterList parameter_list;
436  parameter_list.set("relaxation: sweeps",
437  static_cast<int>(additional_data.n_sweeps));
438  parameter_list.set("relaxation: type", "Gauss-Seidel");
439  parameter_list.set("relaxation: damping factor", additional_data.omega);
440  parameter_list.set("relaxation: min diagonal value",
441  additional_data.min_diagonal);
442  parameter_list.set("schwarz: combine mode", "Add");
443  parameter_list.set("partitioner: type",
444  additional_data.block_creation_type);
445  int n_local_parts =
446  (matrix.trilinos_matrix().NumMyRows() + additional_data.block_size - 1) /
447  additional_data.block_size;
448  parameter_list.set("partitioner: local parts", n_local_parts);
449 
450  ierr = ifpack->SetParameters(parameter_list);
451  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
452 
453  ierr = ifpack->Initialize();
454  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
455 
456  ierr = ifpack->Compute();
457  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
458  }
459 
460 
461 
462  /* -------------------------- PreconditionIC -------------------------- */
463 
465  const double ic_atol,
466  const double ic_rtol,
467  const unsigned int overlap)
468  : ic_fill(ic_fill)
469  , ic_atol(ic_atol)
470  , ic_rtol(ic_rtol)
471  , overlap(overlap)
472  {}
473 
474 
475 
476  void
478  const AdditionalData &additional_data)
479  {
480  preconditioner.reset();
481  preconditioner.reset(
482  Ifpack().Create("IC",
483  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
484  additional_data.overlap));
485 
486  Ifpack_Preconditioner *ifpack =
487  dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
488  Assert(ifpack != nullptr,
489  ExcMessage("Trilinos could not create this "
490  "preconditioner"));
491 
492  int ierr;
493 
494  Teuchos::ParameterList parameter_list;
495  parameter_list.set("fact: level-of-fill", additional_data.ic_fill);
496  parameter_list.set("fact: absolute threshold", additional_data.ic_atol);
497  parameter_list.set("fact: relative threshold", additional_data.ic_rtol);
498  parameter_list.set("schwarz: combine mode", "Add");
499 
500  ierr = ifpack->SetParameters(parameter_list);
501  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
502 
503  ierr = ifpack->Initialize();
504  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
505 
506  ierr = ifpack->Compute();
507  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
508  }
509 
510 
511 
512  /* -------------------------- PreconditionILU -------------------------- */
513 
515  const double ilu_atol,
516  const double ilu_rtol,
517  const unsigned int overlap)
518  : ilu_fill(ilu_fill)
519  , ilu_atol(ilu_atol)
520  , ilu_rtol(ilu_rtol)
521  , overlap(overlap)
522  {}
523 
524 
525 
526  void
528  const AdditionalData &additional_data)
529  {
530  preconditioner.reset();
531  preconditioner.reset(
532  Ifpack().Create("ILU",
533  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
534  additional_data.overlap));
535 
536  Ifpack_Preconditioner *ifpack =
537  dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
538  Assert(ifpack != nullptr,
539  ExcMessage("Trilinos could not create this "
540  "preconditioner"));
541 
542  int ierr;
543 
544  Teuchos::ParameterList parameter_list;
545  parameter_list.set("fact: level-of-fill",
546  static_cast<int>(additional_data.ilu_fill));
547  parameter_list.set("fact: absolute threshold", additional_data.ilu_atol);
548  parameter_list.set("fact: relative threshold", additional_data.ilu_rtol);
549  parameter_list.set("schwarz: combine mode", "Add");
550 
551  ierr = ifpack->SetParameters(parameter_list);
552  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
553 
554  ierr = ifpack->Initialize();
555  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
556 
557  ierr = ifpack->Compute();
558  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
559  }
560 
561 
562 
563  /* -------------------------- PreconditionILUT -------------------------- */
564 
566  const unsigned int ilut_fill,
567  const double ilut_atol,
568  const double ilut_rtol,
569  const unsigned int overlap)
570  : ilut_drop(ilut_drop)
571  , ilut_fill(ilut_fill)
572  , ilut_atol(ilut_atol)
573  , ilut_rtol(ilut_rtol)
574  , overlap(overlap)
575  {}
576 
577 
578 
579  void
581  const AdditionalData &additional_data)
582  {
583  preconditioner.reset();
584  preconditioner.reset(
585  Ifpack().Create("ILUT",
586  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
587  additional_data.overlap));
588 
589  Ifpack_Preconditioner *ifpack =
590  dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
591  Assert(ifpack != nullptr,
592  ExcMessage("Trilinos could not create this "
593  "preconditioner"));
594 
595  int ierr;
596 
597  Teuchos::ParameterList parameter_list;
598  parameter_list.set("fact: drop value", additional_data.ilut_drop);
599  parameter_list.set("fact: level-of-fill",
600  static_cast<int>(additional_data.ilut_fill));
601  parameter_list.set("fact: absolute threshold", additional_data.ilut_atol);
602  parameter_list.set("fact: relative threshold", additional_data.ilut_rtol);
603  parameter_list.set("schwarz: combine mode", "Add");
604 
605  ierr = ifpack->SetParameters(parameter_list);
606  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
607 
608  ierr = ifpack->Initialize();
609  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
610 
611  ierr = ifpack->Compute();
612  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
613  }
614 
615 
616 
617  /* ---------------------- PreconditionBlockDirect --------------------- */
618 
620  const unsigned int overlap)
621  : overlap(overlap)
622  {}
623 
624 
625 
626  void
628  const AdditionalData &additional_data)
629  {
630  preconditioner.reset();
631  preconditioner.reset(
632  Ifpack().Create("Amesos",
633  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
634  additional_data.overlap));
635 
636  Ifpack_Preconditioner *ifpack =
637  dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
638  Assert(ifpack != nullptr,
639  ExcMessage("Trilinos could not create this "
640  "preconditioner"));
641 
642  int ierr;
643 
644  Teuchos::ParameterList parameter_list;
645  parameter_list.set("schwarz: combine mode", "Add");
646 
647  ierr = ifpack->SetParameters(parameter_list);
648  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
649 
650  ierr = ifpack->Initialize();
651  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
652 
653  ierr = ifpack->Compute();
654  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
655  }
656 
657 
658 
659  /* ---------------------- PreconditionBlockDirect --------------------- */
660 
662  const unsigned int degree,
663  const double max_eigenvalue,
664  const double eigenvalue_ratio,
665  const double min_eigenvalue,
666  const double min_diagonal,
667  const bool nonzero_starting)
668  : degree(degree)
669  , max_eigenvalue(max_eigenvalue)
670  , eigenvalue_ratio(eigenvalue_ratio)
671  , min_eigenvalue(min_eigenvalue)
672  , min_diagonal(min_diagonal)
673  , nonzero_starting(nonzero_starting)
674  {}
675 
676 
677 
678  void
680  const AdditionalData &additional_data)
681  {
683  Teuchos::rcp(new Ifpack_Chebyshev(&matrix.trilinos_matrix()));
684 
685  Ifpack_Chebyshev *ifpack =
686  dynamic_cast<Ifpack_Chebyshev *>(preconditioner.get());
687  Assert(ifpack != nullptr,
688  ExcMessage("Trilinos could not create this "
689  "preconditioner"));
690 
691  int ierr;
692 
693  Teuchos::ParameterList parameter_list;
694  parameter_list.set("chebyshev: ratio eigenvalue",
695  additional_data.eigenvalue_ratio);
696  parameter_list.set("chebyshev: min eigenvalue",
697  additional_data.min_eigenvalue);
698  parameter_list.set("chebyshev: max eigenvalue",
699  additional_data.max_eigenvalue);
700  parameter_list.set("chebyshev: degree",
701  static_cast<int>(additional_data.degree));
702  parameter_list.set("chebyshev: min diagonal value",
703  additional_data.min_diagonal);
704  parameter_list.set("chebyshev: zero starting solution",
705  !additional_data.nonzero_starting);
706 
707  ierr = ifpack->SetParameters(parameter_list);
708  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
709 
710  ierr = ifpack->Initialize();
711  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
712 
713  ierr = ifpack->Compute();
714  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
715  }
716 
717 
718 
719  /* -------------------------- PreconditionIdentity --------------------- */
720 
721  void
723  const AdditionalData &)
724  {
725  // What follows just configures a dummy preconditioner that
726  // sets up the domain and range maps, as well as the communicator.
727  // It is never used as the vmult, Tvmult operations are
728  // given a custom definition.
729  // Note: This is only required in order to wrap this
730  // preconditioner in a LinearOperator without an exemplar
731  // matrix.
732 
733  // From PreconditionJacobi:
734  // release memory before reallocation
735  preconditioner.reset();
736  preconditioner.reset(
737  Ifpack().Create("point relaxation",
738  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
739  0));
740 
741  Ifpack_Preconditioner *ifpack =
742  dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
743  Assert(ifpack != nullptr,
744  ExcMessage("Trilinos could not create this "
745  "preconditioner"));
746 
747  int ierr;
748 
749  Teuchos::ParameterList parameter_list;
750  parameter_list.set("relaxation: sweeps", 1);
751  parameter_list.set("relaxation: type", "Jacobi");
752  parameter_list.set("relaxation: damping factor", 1.0);
753  parameter_list.set("relaxation: min diagonal value", 0.0);
754 
755  ierr = ifpack->SetParameters(parameter_list);
756  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
757 
758  ierr = ifpack->Initialize();
759  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
760 
761  ierr = ifpack->Compute();
762  AssertThrow(ierr == 0, ExcTrilinosError(ierr));
763  }
764 
765  void
767  {
768  dst = src;
769  }
770 
771  void
773  {
774  dst = src;
775  }
776 
777  void
779  const ::Vector<double> &src) const
780  {
781  dst = src;
782  }
783 
784  void
786  const ::Vector<double> &src) const
787  {
788  dst = src;
789  }
790 
791  void
795  {
796  dst = src;
797  }
798 
799  void
803  {
804  dst = src;
805  }
806 } // namespace TrilinosWrappers
807 
808 DEAL_II_NAMESPACE_CLOSE
809 
810 #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)
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)
Definition: exceptions.h:1519
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)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1407
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 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())
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())
void vmult(MPI::Vector &dst, const MPI::Vector &src) const override
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())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData(const double omega=1, const double min_diagonal=0, const unsigned int n_sweeps=1)
Teuchos::RCP< Epetra_Operator > preconditioner