Reference documentation for deal.II version 9.0.0
trilinos_precondition.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2018 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/lac/trilinos_precondition.h>
17 
18 #ifdef DEAL_II_WITH_TRILINOS
19 
20 # include <deal.II/lac/vector.h>
21 # include <deal.II/lac/sparse_matrix.h>
22 # include <deal.II/lac/trilinos_sparse_matrix.h>
23 
24 # include <Ifpack.h>
25 # include <Ifpack_Chebyshev.h>
26 # include <Teuchos_ParameterList.hpp>
27 # include <Teuchos_RCP.hpp>
28 # include <Epetra_MultiVector.h>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 namespace TrilinosWrappers
33 {
34 
36 #ifdef DEAL_II_WITH_MPI
37  :
38  communicator (MPI_COMM_SELF)
39 #endif
40  {}
41 
42 
43 
45  :
46  Subscriptor (),
47  preconditioner (base.preconditioner),
48 #ifdef DEAL_II_WITH_MPI
49  communicator (base.communicator),
50 #endif
51  vector_distributor (new Epetra_Map(*base.vector_distributor))
52  {}
53 
54 
55 
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, ExcMessage("Trying to dereference a null pointer."));
81  return (*preconditioner);
82  }
83 
84 
85  IndexSet
87  {
88  return IndexSet(preconditioner->OperatorDomainMap());
89  }
90 
91 
92  IndexSet
94  {
95  return IndexSet(preconditioner->OperatorRangeMap());
96  }
97 
98  /* -------------------------- PreconditionJacobi -------------------------- */
99 
101  AdditionalData (const double omega,
102  const double min_diagonal,
103  const unsigned int n_sweeps)
104  :
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 (Ifpack().Create
119  ("point relaxation",
120  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
121  0));
122 
123  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
124  (preconditioner.get());
125  Assert (ifpack != nullptr, ExcMessage ("Trilinos could not create this "
126  "preconditioner"));
127 
128  int ierr;
129 
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",
135  additional_data.min_diagonal);
136 
137  ierr = ifpack->SetParameters(parameter_list);
138  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
139 
140  ierr = ifpack->Initialize();
141  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
142 
143  ierr = ifpack->Compute();
144  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
145  }
146 
147 
148 
149  /* -------------------------- PreconditionSSOR -------------------------- */
150 
152  AdditionalData (const double omega,
153  const double min_diagonal,
154  const unsigned int overlap,
155  const unsigned int n_sweeps)
156  :
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 (Ifpack().Create
171  ("point relaxation",
172  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
173  additional_data.overlap));
174 
175  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
176  (preconditioner.get());
177  Assert (ifpack != nullptr, ExcMessage ("Trilinos could not create this "
178  "preconditioner"));
179 
180  int ierr;
181 
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",
187  additional_data.min_diagonal);
188  parameter_list.set ("schwarz: combine mode", "Add");
189 
190  ierr = ifpack->SetParameters(parameter_list);
191  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
192 
193  ierr = ifpack->Initialize();
194  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
195 
196  ierr = ifpack->Compute();
197  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
198  }
199 
200 
201 
202  /* -------------------------- PreconditionSOR -------------------------- */
203 
205  AdditionalData (const double omega,
206  const double min_diagonal,
207  const unsigned int overlap,
208  const unsigned int n_sweeps)
209  :
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 (Ifpack().Create
224  ("point relaxation",
225  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
226  additional_data.overlap));
227 
228  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
229  (preconditioner.get());
230  Assert (ifpack != nullptr, ExcMessage ("Trilinos could not create this "
231  "preconditioner"));
232 
233  int ierr;
234 
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",
240  additional_data.min_diagonal);
241  parameter_list.set ("schwarz: combine mode", "Add");
242 
243  ierr = ifpack->SetParameters(parameter_list);
244  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
245 
246  ierr = ifpack->Initialize();
247  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
248 
249  ierr = ifpack->Compute();
250  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
251  }
252 
253 
254 
255  /* ----------------------- PreconditionBlockJacobi ---------------------- */
256 
258  AdditionalData (const unsigned int block_size,
259  const std::string &block_creation_type,
260  const double omega,
261  const double min_diagonal,
262  const unsigned int n_sweeps)
263  :
264  block_size(block_size),
265  block_creation_type(block_creation_type),
266  omega (omega),
267  min_diagonal (min_diagonal),
268  n_sweeps (n_sweeps)
269  {}
270 
271 
272 
273  void
275  const AdditionalData &additional_data)
276  {
277  // release memory before reallocation
278  preconditioner.reset ();
279 
280  // Block relaxation setup fails if we have no locally owned rows. As a work-around
281  // we just pretend to use point relaxation on those processors:
282  preconditioner.reset (Ifpack().Create
283  (
284  (matrix.trilinos_matrix().NumMyRows()==0) ?
285  "point relaxation" : "block relaxation",
286  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
287  0));
288 
289  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
290  (preconditioner.get());
291  Assert (ifpack != nullptr, ExcMessage ("Trilinos could not create this "
292  "preconditioner"));
293 
294  int ierr;
295 
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",
301  additional_data.min_diagonal);
302  parameter_list.set ("partitioner: type", additional_data.block_creation_type);
303  int n_local_parts = (matrix.trilinos_matrix().NumMyRows()+additional_data.
304  block_size-1)/additional_data.block_size;
305  parameter_list.set ("partitioner: local parts", n_local_parts);
306 
307  ierr = ifpack->SetParameters(parameter_list);
308  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
309 
310  ierr = ifpack->Initialize();
311  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
312 
313  ierr = ifpack->Compute();
314  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
315  }
316 
317 
318 
319  /* ----------------------- PreconditionBlockSSOR ------------------------ */
320 
322  AdditionalData (const unsigned int block_size,
323  const std::string &block_creation_type,
324  const double omega,
325  const double min_diagonal,
326  const unsigned int overlap,
327  const unsigned int n_sweeps)
328  :
329  block_size(block_size),
330  block_creation_type(block_creation_type),
331  omega (omega),
332  min_diagonal (min_diagonal),
333  overlap (overlap),
334  n_sweeps (n_sweeps)
335  {}
336 
337 
338 
339  void
341  const AdditionalData &additional_data)
342  {
343  preconditioner.reset ();
344 
345  // Block relaxation setup fails if we have no locally owned rows. As a work-around
346  // we just pretend to use point relaxation on those processors:
347  preconditioner.reset (Ifpack().Create
348  (
349  (matrix.trilinos_matrix().NumMyRows()==0) ?
350  "point relaxation" : "block relaxation",
351  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
352  additional_data.overlap));
353 
354  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
355  (preconditioner.get());
356  Assert (ifpack != nullptr, ExcMessage ("Trilinos could not create this "
357  "preconditioner"));
358 
359  int ierr;
360 
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",
366  additional_data.min_diagonal);
367  parameter_list.set ("schwarz: combine mode", "Add");
368  parameter_list.set ("partitioner: type", additional_data.block_creation_type);
369  int n_local_parts = (matrix.trilinos_matrix().NumMyRows()+additional_data.
370  block_size-1)/additional_data.block_size;
371  parameter_list.set ("partitioner: local parts", n_local_parts);
372 
373  ierr = ifpack->SetParameters(parameter_list);
374  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
375 
376  ierr = ifpack->Initialize();
377  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
378 
379  ierr = ifpack->Compute();
380  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
381  }
382 
383 
384 
385  /* ------------------------ PreconditionBlockSOR ------------------------ */
386 
388  AdditionalData (const unsigned int block_size,
389  const std::string &block_creation_type,
390  const double omega,
391  const double min_diagonal,
392  const unsigned int overlap,
393  const unsigned int n_sweeps)
394  :
395  block_size(block_size),
396  block_creation_type(block_creation_type),
397  omega (omega),
398  min_diagonal (min_diagonal),
399  overlap (overlap),
400  n_sweeps (n_sweeps)
401  {}
402 
403 
404 
405  void
407  const AdditionalData &additional_data)
408  {
409  preconditioner.reset ();
410 
411  // Block relaxation setup fails if we have no locally owned rows. As a work-around
412  // we just pretend to use point relaxation on those processors:
413  preconditioner.reset (Ifpack().Create
414  (
415  (matrix.trilinos_matrix().NumMyRows()==0) ?
416  "point relaxation" : "block relaxation",
417  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
418  additional_data.overlap));
419 
420  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
421  (preconditioner.get());
422  Assert (ifpack != nullptr, ExcMessage ("Trilinos could not create this "
423  "preconditioner"));
424 
425  int ierr;
426 
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",
432  additional_data.min_diagonal);
433  parameter_list.set ("schwarz: combine mode", "Add");
434  parameter_list.set ("partitioner: type", additional_data.block_creation_type);
435  int n_local_parts = (matrix.trilinos_matrix().NumMyRows()+additional_data.
436  block_size-1)/additional_data.block_size;
437  parameter_list.set ("partitioner: local parts", n_local_parts);
438 
439  ierr = ifpack->SetParameters(parameter_list);
440  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
441 
442  ierr = ifpack->Initialize();
443  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
444 
445  ierr = ifpack->Compute();
446  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
447  }
448 
449 
450 
451  /* -------------------------- PreconditionIC -------------------------- */
452 
454  AdditionalData (const unsigned int ic_fill,
455  const double ic_atol,
456  const double ic_rtol,
457  const unsigned int overlap)
458  :
459  ic_fill (ic_fill),
460  ic_atol (ic_atol),
461  ic_rtol (ic_rtol),
462  overlap (overlap)
463  {}
464 
465 
466 
467  void
469  const AdditionalData &additional_data)
470  {
471  preconditioner.reset ();
472  preconditioner.reset (Ifpack().Create
473  ("IC",
474  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
475  additional_data.overlap));
476 
477  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
478  (preconditioner.get());
479  Assert (ifpack != nullptr, ExcMessage ("Trilinos could not create this "
480  "preconditioner"));
481 
482  int ierr;
483 
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");
489 
490  ierr = ifpack->SetParameters(parameter_list);
491  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
492 
493  ierr = ifpack->Initialize();
494  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
495 
496  ierr = ifpack->Compute();
497  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
498  }
499 
500 
501 
502  /* -------------------------- PreconditionILU -------------------------- */
503 
505  AdditionalData (const unsigned int ilu_fill,
506  const double ilu_atol,
507  const double ilu_rtol,
508  const unsigned int overlap)
509  :
510  ilu_fill (ilu_fill),
511  ilu_atol (ilu_atol),
512  ilu_rtol (ilu_rtol),
513  overlap (overlap)
514  {}
515 
516 
517 
518  void
520  const AdditionalData &additional_data)
521  {
522  preconditioner.reset ();
523  preconditioner.reset (Ifpack().Create
524  ("ILU",
525  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
526  additional_data.overlap));
527 
528  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
529  (preconditioner.get());
530  Assert (ifpack != nullptr, ExcMessage ("Trilinos could not create this "
531  "preconditioner"));
532 
533  int ierr;
534 
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");
540 
541  ierr = ifpack->SetParameters(parameter_list);
542  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
543 
544  ierr = ifpack->Initialize();
545  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
546 
547  ierr = ifpack->Compute();
548  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
549  }
550 
551 
552 
553  /* -------------------------- PreconditionILUT -------------------------- */
554 
556  AdditionalData (const double ilut_drop,
557  const unsigned int ilut_fill,
558  const double ilut_atol,
559  const double ilut_rtol,
560  const unsigned int overlap)
561  :
562  ilut_drop (ilut_drop),
563  ilut_fill (ilut_fill),
564  ilut_atol (ilut_atol),
565  ilut_rtol (ilut_rtol),
566  overlap (overlap)
567  {}
568 
569 
570 
571  void
573  const AdditionalData &additional_data)
574  {
575  preconditioner.reset ();
576  preconditioner.reset (Ifpack().Create
577  ("ILUT",
578  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
579  additional_data.overlap));
580 
581  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
582  (preconditioner.get());
583  Assert (ifpack != nullptr, ExcMessage ("Trilinos could not create this "
584  "preconditioner"));
585 
586  int ierr;
587 
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");
594 
595  ierr = ifpack->SetParameters(parameter_list);
596  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
597 
598  ierr = ifpack->Initialize();
599  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
600 
601  ierr = ifpack->Compute();
602  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
603  }
604 
605 
606 
607  /* ---------------------- PreconditionBlockDirect --------------------- */
608 
610  AdditionalData (const unsigned int overlap)
611  :
612  overlap (overlap)
613  {}
614 
615 
616 
617  void
619  const AdditionalData &additional_data)
620  {
621  preconditioner.reset ();
622  preconditioner.reset (Ifpack().Create
623  ("Amesos",
624  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
625  additional_data.overlap));
626 
627  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
628  (preconditioner.get());
629  Assert (ifpack != nullptr, ExcMessage ("Trilinos could not create this "
630  "preconditioner"));
631 
632  int ierr;
633 
634  Teuchos::ParameterList parameter_list;
635  parameter_list.set ("schwarz: combine mode", "Add");
636 
637  ierr = ifpack->SetParameters(parameter_list);
638  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
639 
640  ierr = ifpack->Initialize();
641  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
642 
643  ierr = ifpack->Compute();
644  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
645  }
646 
647 
648 
649  /* ---------------------- PreconditionBlockDirect --------------------- */
650 
652  AdditionalData (const unsigned int degree,
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)
658  :
659  degree (degree),
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)
665  {}
666 
667 
668 
669  void
671  const AdditionalData &additional_data)
672  {
673  preconditioner.reset ();
674  preconditioner = std::make_shared<Ifpack_Chebyshev> (&matrix.trilinos_matrix());
675 
676  Ifpack_Chebyshev *ifpack = static_cast<Ifpack_Chebyshev *>
677  (preconditioner.get());
678  Assert (ifpack != nullptr, ExcMessage ("Trilinos could not create this "
679  "preconditioner"));
680 
681  int ierr;
682 
683  Teuchos::ParameterList parameter_list;
684  parameter_list.set ("chebyshev: ratio eigenvalue",
685  additional_data.eigenvalue_ratio);
686  parameter_list.set ("chebyshev: min eigenvalue",
687  additional_data.min_eigenvalue);
688  parameter_list.set ("chebyshev: max eigenvalue",
689  additional_data.max_eigenvalue);
690  parameter_list.set ("chebyshev: degree",
691  (int)additional_data.degree);
692  parameter_list.set ("chebyshev: min diagonal value",
693  additional_data.min_diagonal);
694  parameter_list.set ("chebyshev: zero starting solution",
695  !additional_data.nonzero_starting);
696 
697  ierr = ifpack->SetParameters(parameter_list);
698  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
699 
700  ierr = ifpack->Initialize();
701  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
702 
703  ierr = ifpack->Compute();
704  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
705  }
706 
707 
708 
709 
710 
711 
712 
713 
714  /* -------------------------- PreconditionIdentity --------------------- */
715 
716  void
718  const AdditionalData &)
719  {
720  // What follows just configures a dummy preconditioner that
721  // sets up the domain and range maps, as well as the communicator.
722  // It is never used as the vmult, Tvmult operations are
723  // given a custom definition.
724  // Note: This is only required in order to wrap this
725  // preconditioner in a LinearOperator without an exemplar
726  // matrix.
727 
728  // From PreconditionJacobi:
729  // release memory before reallocation
730  preconditioner.reset ();
731  preconditioner.reset (Ifpack().Create
732  ("point relaxation",
733  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
734  0));
735 
736  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
737  (preconditioner.get());
738  Assert (ifpack != nullptr, ExcMessage ("Trilinos could not create this "
739  "preconditioner"));
740 
741  int ierr;
742 
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);
748 
749  ierr = ifpack->SetParameters(parameter_list);
750  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
751 
752  ierr = ifpack->Initialize();
753  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
754 
755  ierr = ifpack->Compute();
756  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
757  }
758 
759  void
761  const MPI::Vector &src) const
762  {
763  dst = src;
764  }
765 
766  void
768  const MPI::Vector &src) const
769  {
770  dst = src;
771  }
772 
773  void
775  const ::Vector<double> &src) const
776  {
777  dst = src;
778  }
779 
780  void
782  const ::Vector<double> &src) const
783  {
784  dst = src;
785  }
786 
787  void
790  {
791  dst = src;
792  }
793 
794  void
797  {
798  dst = src;
799  }
800 }
801 
802 DEAL_II_NAMESPACE_CLOSE
803 
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)
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:1221
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:1142
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())
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())
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)