Reference documentation for deal.II version 8.5.1
trilinos_precondition.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2017 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 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
25 # include <Ifpack.h>
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
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 namespace TrilinosWrappers
35 {
36 
38 #ifdef DEAL_II_WITH_MPI
39  :
40  communicator (MPI_COMM_SELF)
41 #endif
42  {}
43 
44 
45 
47  :
48  Subscriptor (),
49  preconditioner (base.preconditioner),
50 #ifdef DEAL_II_WITH_MPI
51  communicator (base.communicator),
52 #endif
53  vector_distributor (new Epetra_Map(*base.vector_distributor))
54  {}
55 
56 
57 
59  {}
60 
61 
62 
64  {
65  preconditioner.reset();
66 #ifdef DEAL_II_WITH_MPI
67  communicator = MPI_COMM_SELF;
68 #endif
69  vector_distributor.reset();
70  }
71 
72 
73  MPI_Comm
75  {
76 #ifdef DEAL_II_WITH_MPI
77  return communicator.Comm();
78 #else
79  return MPI_COMM_SELF;
80 #endif
81  }
82 
83 
84  Epetra_Operator &
86  {
87  AssertThrow (preconditioner, ExcMessage("Trying to dereference a null pointer."));
88  return (*preconditioner);
89  }
90 
91 
92  IndexSet
94  {
95  return IndexSet(preconditioner->OperatorDomainMap());
96  }
97 
98 
99  IndexSet
101  {
102  return IndexSet(preconditioner->OperatorRangeMap());
103  }
104 
105  /* -------------------------- PreconditionJacobi -------------------------- */
106 
108  AdditionalData (const double omega,
109  const double min_diagonal,
110  const unsigned int n_sweeps)
111  :
112  omega (omega),
113  min_diagonal (min_diagonal),
114  n_sweeps (n_sweeps)
115  {}
116 
117 
118 
119  void
121  const AdditionalData &additional_data)
122  {
123  // release memory before reallocation
124  preconditioner.reset ();
125  preconditioner.reset (Ifpack().Create
126  ("point relaxation",
127  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
128  0));
129 
130  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
131  (preconditioner.get());
132  Assert (ifpack != 0, ExcMessage ("Trilinos could not create this "
133  "preconditioner"));
134 
135  int ierr;
136 
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",
142  additional_data.min_diagonal);
143 
144  ierr = ifpack->SetParameters(parameter_list);
145  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
146 
147  ierr = ifpack->Initialize();
148  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
149 
150  ierr = ifpack->Compute();
151  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
152  }
153 
154 
155 
156  /* -------------------------- PreconditionSSOR -------------------------- */
157 
159  AdditionalData (const double omega,
160  const double min_diagonal,
161  const unsigned int overlap,
162  const unsigned int n_sweeps)
163  :
164  omega (omega),
165  min_diagonal (min_diagonal),
166  overlap (overlap),
167  n_sweeps (n_sweeps)
168  {}
169 
170 
171 
172  void
174  const AdditionalData &additional_data)
175  {
176  preconditioner.reset ();
177  preconditioner.reset (Ifpack().Create
178  ("point relaxation",
179  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
180  additional_data.overlap));
181 
182  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
183  (preconditioner.get());
184  Assert (ifpack != 0, ExcMessage ("Trilinos could not create this "
185  "preconditioner"));
186 
187  int ierr;
188 
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",
194  additional_data.min_diagonal);
195  parameter_list.set ("schwarz: combine mode", "Add");
196 
197  ierr = ifpack->SetParameters(parameter_list);
198  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
199 
200  ierr = ifpack->Initialize();
201  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
202 
203  ierr = ifpack->Compute();
204  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
205  }
206 
207 
208 
209  /* -------------------------- PreconditionSOR -------------------------- */
210 
212  AdditionalData (const double omega,
213  const double min_diagonal,
214  const unsigned int overlap,
215  const unsigned int n_sweeps)
216  :
217  omega (omega),
218  min_diagonal (min_diagonal),
219  overlap (overlap),
220  n_sweeps (n_sweeps)
221  {}
222 
223 
224 
225  void
227  const AdditionalData &additional_data)
228  {
229  preconditioner.reset ();
230  preconditioner.reset (Ifpack().Create
231  ("point relaxation",
232  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
233  additional_data.overlap));
234 
235  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
236  (preconditioner.get());
237  Assert (ifpack != 0, ExcMessage ("Trilinos could not create this "
238  "preconditioner"));
239 
240  int ierr;
241 
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",
247  additional_data.min_diagonal);
248  parameter_list.set ("schwarz: combine mode", "Add");
249 
250  ierr = ifpack->SetParameters(parameter_list);
251  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
252 
253  ierr = ifpack->Initialize();
254  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
255 
256  ierr = ifpack->Compute();
257  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
258  }
259 
260 
261 
262  /* ----------------------- PreconditionBlockJacobi ---------------------- */
263 
265  AdditionalData (const unsigned int block_size,
266  const std::string &block_creation_type,
267  const double omega,
268  const double min_diagonal,
269  const unsigned int n_sweeps)
270  :
271  block_size(block_size),
272  block_creation_type(block_creation_type),
273  omega (omega),
274  min_diagonal (min_diagonal),
275  n_sweeps (n_sweeps)
276  {}
277 
278 
279 
280  void
282  const AdditionalData &additional_data)
283  {
284  // release memory before reallocation
285  preconditioner.reset ();
286  preconditioner.reset (Ifpack().Create
287  ("block relaxation",
288  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
289  0));
290 
291  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
292  (preconditioner.get());
293  Assert (ifpack != 0, ExcMessage ("Trilinos could not create this "
294  "preconditioner"));
295 
296  int ierr;
297 
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",
303  additional_data.min_diagonal);
304  parameter_list.set ("partitioner: type", additional_data.block_creation_type);
305  int n_local_parts = (matrix.trilinos_matrix().NumMyRows()+additional_data.
306  block_size-1)/additional_data.block_size;
307  parameter_list.set ("partitioner: local parts", n_local_parts);
308 
309  ierr = ifpack->SetParameters(parameter_list);
310  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
311 
312  ierr = ifpack->Initialize();
313  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
314 
315  ierr = ifpack->Compute();
316  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
317  }
318 
319 
320 
321  /* ----------------------- PreconditionBlockSSOR ------------------------ */
322 
324  AdditionalData (const unsigned int block_size,
325  const std::string &block_creation_type,
326  const double omega,
327  const double min_diagonal,
328  const unsigned int overlap,
329  const unsigned int n_sweeps)
330  :
331  block_size(block_size),
332  block_creation_type(block_creation_type),
333  omega (omega),
334  min_diagonal (min_diagonal),
335  overlap (overlap),
336  n_sweeps (n_sweeps)
337  {}
338 
339 
340 
341  void
343  const AdditionalData &additional_data)
344  {
345  preconditioner.reset ();
346  preconditioner.reset (Ifpack().Create
347  ("block relaxation",
348  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
349  additional_data.overlap));
350 
351  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
352  (preconditioner.get());
353  Assert (ifpack != 0, ExcMessage ("Trilinos could not create this "
354  "preconditioner"));
355 
356  int ierr;
357 
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",
363  additional_data.min_diagonal);
364  parameter_list.set ("schwarz: combine mode", "Add");
365  parameter_list.set ("partitioner: type", additional_data.block_creation_type);
366  int n_local_parts = (matrix.trilinos_matrix().NumMyRows()+additional_data.
367  block_size-1)/additional_data.block_size;
368  parameter_list.set ("partitioner: local parts", n_local_parts);
369 
370  ierr = ifpack->SetParameters(parameter_list);
371  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
372 
373  ierr = ifpack->Initialize();
374  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
375 
376  ierr = ifpack->Compute();
377  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
378  }
379 
380 
381 
382  /* ------------------------ PreconditionBlockSOR ------------------------ */
383 
385  AdditionalData (const unsigned int block_size,
386  const std::string &block_creation_type,
387  const double omega,
388  const double min_diagonal,
389  const unsigned int overlap,
390  const unsigned int n_sweeps)
391  :
392  block_size(block_size),
393  block_creation_type(block_creation_type),
394  omega (omega),
395  min_diagonal (min_diagonal),
396  overlap (overlap),
397  n_sweeps (n_sweeps)
398  {}
399 
400 
401 
402  void
404  const AdditionalData &additional_data)
405  {
406  preconditioner.reset ();
407  preconditioner.reset (Ifpack().Create
408  ("block relaxation",
409  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
410  additional_data.overlap));
411 
412  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
413  (preconditioner.get());
414  Assert (ifpack != 0, ExcMessage ("Trilinos could not create this "
415  "preconditioner"));
416 
417  int ierr;
418 
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",
424  additional_data.min_diagonal);
425  parameter_list.set ("schwarz: combine mode", "Add");
426  parameter_list.set ("partitioner: type", additional_data.block_creation_type);
427  int n_local_parts = (matrix.trilinos_matrix().NumMyRows()+additional_data.
428  block_size-1)/additional_data.block_size;
429  parameter_list.set ("partitioner: local parts", n_local_parts);
430 
431  ierr = ifpack->SetParameters(parameter_list);
432  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
433 
434  ierr = ifpack->Initialize();
435  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
436 
437  ierr = ifpack->Compute();
438  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
439  }
440 
441 
442 
443  /* -------------------------- PreconditionIC -------------------------- */
444 
446  AdditionalData (const unsigned int ic_fill,
447  const double ic_atol,
448  const double ic_rtol,
449  const unsigned int overlap)
450  :
451  ic_fill (ic_fill),
452  ic_atol (ic_atol),
453  ic_rtol (ic_rtol),
454  overlap (overlap)
455  {}
456 
457 
458 
459  void
461  const AdditionalData &additional_data)
462  {
463  preconditioner.reset ();
464  preconditioner.reset (Ifpack().Create
465  ("IC",
466  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
467  additional_data.overlap));
468 
469  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
470  (preconditioner.get());
471  Assert (ifpack != 0, ExcMessage ("Trilinos could not create this "
472  "preconditioner"));
473 
474  int ierr;
475 
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");
481 
482  ierr = ifpack->SetParameters(parameter_list);
483  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
484 
485  ierr = ifpack->Initialize();
486  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
487 
488  ierr = ifpack->Compute();
489  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
490  }
491 
492 
493 
494  /* -------------------------- PreconditionILU -------------------------- */
495 
497  AdditionalData (const unsigned int ilu_fill,
498  const double ilu_atol,
499  const double ilu_rtol,
500  const unsigned int overlap)
501  :
502  ilu_fill (ilu_fill),
503  ilu_atol (ilu_atol),
504  ilu_rtol (ilu_rtol),
505  overlap (overlap)
506  {}
507 
508 
509 
510  void
512  const AdditionalData &additional_data)
513  {
514  preconditioner.reset ();
515  preconditioner.reset (Ifpack().Create
516  ("ILU",
517  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
518  additional_data.overlap));
519 
520  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
521  (preconditioner.get());
522  Assert (ifpack != 0, ExcMessage ("Trilinos could not create this "
523  "preconditioner"));
524 
525  int ierr;
526 
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");
532 
533  ierr = ifpack->SetParameters(parameter_list);
534  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
535 
536  ierr = ifpack->Initialize();
537  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
538 
539  ierr = ifpack->Compute();
540  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
541  }
542 
543 
544 
545  /* -------------------------- PreconditionILUT -------------------------- */
546 
548  AdditionalData (const double ilut_drop,
549  const unsigned int ilut_fill,
550  const double ilut_atol,
551  const double ilut_rtol,
552  const unsigned int overlap)
553  :
554  ilut_drop (ilut_drop),
555  ilut_fill (ilut_fill),
556  ilut_atol (ilut_atol),
557  ilut_rtol (ilut_rtol),
558  overlap (overlap)
559  {}
560 
561 
562 
563  void
565  const AdditionalData &additional_data)
566  {
567  preconditioner.reset ();
568  preconditioner.reset (Ifpack().Create
569  ("ILUT",
570  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
571  additional_data.overlap));
572 
573  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
574  (preconditioner.get());
575  Assert (ifpack != 0, ExcMessage ("Trilinos could not create this "
576  "preconditioner"));
577 
578  int ierr;
579 
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");
586 
587  ierr = ifpack->SetParameters(parameter_list);
588  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
589 
590  ierr = ifpack->Initialize();
591  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
592 
593  ierr = ifpack->Compute();
594  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
595  }
596 
597 
598 
599  /* ---------------------- PreconditionBlockDirect --------------------- */
600 
602  AdditionalData (const unsigned int overlap)
603  :
604  overlap (overlap)
605  {}
606 
607 
608 
609  void
611  const AdditionalData &additional_data)
612  {
613  preconditioner.reset ();
614  preconditioner.reset (Ifpack().Create
615  ("Amesos",
616  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
617  additional_data.overlap));
618 
619  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
620  (preconditioner.get());
621  Assert (ifpack != 0, ExcMessage ("Trilinos could not create this "
622  "preconditioner"));
623 
624  int ierr;
625 
626  Teuchos::ParameterList parameter_list;
627  parameter_list.set ("schwarz: combine mode", "Add");
628 
629  ierr = ifpack->SetParameters(parameter_list);
630  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
631 
632  ierr = ifpack->Initialize();
633  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
634 
635  ierr = ifpack->Compute();
636  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
637  }
638 
639 
640 
641  /* ---------------------- PreconditionBlockDirect --------------------- */
642 
644  AdditionalData (const unsigned int degree,
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)
650  :
651  degree (degree),
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)
657  {}
658 
659 
660 
661  void
663  const AdditionalData &additional_data)
664  {
665  preconditioner.reset ();
666  preconditioner.reset (new Ifpack_Chebyshev (&matrix.trilinos_matrix()));
667 
668  Ifpack_Chebyshev *ifpack = static_cast<Ifpack_Chebyshev *>
669  (preconditioner.get());
670  Assert (ifpack != 0, ExcMessage ("Trilinos could not create this "
671  "preconditioner"));
672 
673  int ierr;
674 
675  Teuchos::ParameterList parameter_list;
676  parameter_list.set ("chebyshev: ratio eigenvalue",
677  additional_data.eigenvalue_ratio);
678  parameter_list.set ("chebyshev: min eigenvalue",
679  additional_data.min_eigenvalue);
680  parameter_list.set ("chebyshev: max eigenvalue",
681  additional_data.max_eigenvalue);
682  parameter_list.set ("chebyshev: degree",
683  (int)additional_data.degree);
684  parameter_list.set ("chebyshev: min diagonal value",
685  additional_data.min_diagonal);
686  parameter_list.set ("chebyshev: zero starting solution",
687  !additional_data.nonzero_starting);
688 
689  ierr = ifpack->SetParameters(parameter_list);
690  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
691 
692  ierr = ifpack->Initialize();
693  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
694 
695  ierr = ifpack->Compute();
696  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
697  }
698 
699 
700 
701 
702 
703 
704 
705 
706  /* -------------------------- PreconditionIdentity --------------------- */
707 
708  void
710  const AdditionalData &)
711  {
712  // What follows just configures a dummy preconditioner that
713  // sets up the domain and range maps, as well as the communicator.
714  // It is never used as the vmult, Tvmult operations are
715  // given a custom defintion.
716  // Note: This is only required in order to wrap this
717  // preconditioner in a LinearOperator without an exemplar
718  // matrix.
719 
720  // From PreconditionJacobi:
721  // release memory before reallocation
722  preconditioner.reset ();
723  preconditioner.reset (Ifpack().Create
724  ("point relaxation",
725  const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
726  0));
727 
728  Ifpack_Preconditioner *ifpack = static_cast<Ifpack_Preconditioner *>
729  (preconditioner.get());
730  Assert (ifpack != 0, ExcMessage ("Trilinos could not create this "
731  "preconditioner"));
732 
733  int ierr;
734 
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);
740 
741  ierr = ifpack->SetParameters(parameter_list);
742  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
743 
744  ierr = ifpack->Initialize();
745  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
746 
747  ierr = ifpack->Compute();
748  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
749  }
750 
751  void
753  const VectorBase &src) const
754  {
755  dst = src;
756  }
757 
758  void
760  const VectorBase &src) const
761  {
762  dst = src;
763  }
764 
765  void
767  const ::Vector<double> &src) const
768  {
769  dst = src;
770  }
771 
772  void
774  const ::Vector<double> &src) const
775  {
776  dst = src;
777  }
778 
779  void
782  {
783  dst = src;
784  }
785 
786  void
789  {
790  dst = src;
791  }
792 }
793 
794 DEAL_II_NAMESPACE_CLOSE
795 
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)
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)
Definition: exceptions.h:369
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_cxx11::shared_ptr< Epetra_Map > vector_distributor
static ::ExceptionBase & ExcMessage(std::string arg1)
std_cxx11::shared_ptr< Epetra_Operator > preconditioner
#define Assert(cond, exc)
Definition: exceptions.h:313
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())
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())
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)