Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
trilinos_precondition.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2008 - 2021 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
17
18#ifdef DEAL_II_WITH_TRILINOS
19
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
31
32namespace TrilinosWrappers
33{
35 : communicator(MPI_COMM_SELF)
36 {}
37
38
39
41 : Subscriptor()
42 , preconditioner(base.preconditioner)
43 , communicator(base.communicator)
44 , vector_distributor(new Epetra_Map(*base.vector_distributor))
45 {}
46
47
48
49 void
51 {
52 preconditioner.reset();
53 communicator = MPI_COMM_SELF;
54 vector_distributor.reset();
55 }
56
57
60 {
61 return communicator.Comm();
62 }
63
64
67 {
68 AssertThrow(!preconditioner.is_null(),
69 ExcMessage("Trying to dereference a null pointer."));
70 return (*preconditioner);
71 }
72
73
76 {
77 return IndexSet(preconditioner->OperatorDomainMap());
78 }
79
80
83 {
84 return IndexSet(preconditioner->OperatorRangeMap());
85 }
86
87 /* -------------------------- PreconditionJacobi -------------------------- */
88
90 const double omega,
91 const double min_diagonal,
92 const unsigned int n_sweeps)
93 : omega(omega)
94 , min_diagonal(min_diagonal)
95 , n_sweeps(n_sweeps)
96 {}
97
98
99
100 void
102 const AdditionalData &additional_data)
103 {
104 // release memory before reallocation
105 preconditioner.reset();
106 preconditioner.reset(
107 Ifpack().Create("point relaxation",
108 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
109 0));
110
111 Ifpack_Preconditioner *ifpack =
112 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
113 Assert(ifpack != nullptr,
114 ExcMessage("Trilinos could not create this "
115 "preconditioner"));
116
117 int ierr;
118
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",
125 additional_data.min_diagonal);
126
127 ierr = ifpack->SetParameters(parameter_list);
128 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
129
130 ierr = ifpack->Initialize();
131 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
132
133 ierr = ifpack->Compute();
134 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
135 }
136
137
138
139 /* -------------------------- PreconditionSSOR -------------------------- */
140
142 const double min_diagonal,
143 const unsigned int overlap,
144 const unsigned int n_sweeps)
145 : omega(omega)
146 , min_diagonal(min_diagonal)
147 , overlap(overlap)
148 , n_sweeps(n_sweeps)
149 {}
150
151
152
153 void
155 const AdditionalData &additional_data)
156 {
157 preconditioner.reset();
158 preconditioner.reset(
159 Ifpack().Create("point relaxation",
160 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
161 additional_data.overlap));
162
163 Ifpack_Preconditioner *ifpack =
164 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
165 Assert(ifpack != nullptr,
166 ExcMessage("Trilinos could not create this "
167 "preconditioner"));
168
169 int ierr;
170
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",
177 additional_data.min_diagonal);
178 parameter_list.set("schwarz: combine mode", "Add");
179
180 ierr = ifpack->SetParameters(parameter_list);
181 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
182
183 ierr = ifpack->Initialize();
184 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
185
186 ierr = ifpack->Compute();
187 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
188 }
189
190
191
192 /* -------------------------- PreconditionSOR -------------------------- */
193
195 const double min_diagonal,
196 const unsigned int overlap,
197 const unsigned int n_sweeps)
198 : omega(omega)
199 , min_diagonal(min_diagonal)
200 , overlap(overlap)
201 , n_sweeps(n_sweeps)
202 {}
203
204
205
206 void
208 const AdditionalData &additional_data)
209 {
210 preconditioner.reset();
211 preconditioner.reset(
212 Ifpack().Create("point relaxation",
213 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
214 additional_data.overlap));
215
216 Ifpack_Preconditioner *ifpack =
217 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
218 Assert(ifpack != nullptr,
219 ExcMessage("Trilinos could not create this "
220 "preconditioner"));
221
222 int ierr;
223
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",
230 additional_data.min_diagonal);
231 parameter_list.set("schwarz: combine mode", "Add");
232
233 ierr = ifpack->SetParameters(parameter_list);
234 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
235
236 ierr = ifpack->Initialize();
237 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
238
239 ierr = ifpack->Compute();
240 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
241 }
242
243
244
245 /* ----------------------- PreconditionBlockJacobi ---------------------- */
246
248 const unsigned int block_size,
249 const std::string &block_creation_type,
250 const double omega,
251 const double min_diagonal,
252 const unsigned int n_sweeps)
253 : block_size(block_size)
254 , block_creation_type(block_creation_type)
255 , omega(omega)
256 , min_diagonal(min_diagonal)
257 , n_sweeps(n_sweeps)
258 {}
259
260
261
262 void
264 const AdditionalData &additional_data)
265 {
266 // release memory before reallocation
267 preconditioner.reset();
268
269 // Block relaxation setup fails if we have no locally owned rows. As a
270 // work-around we just pretend to use point relaxation on those processors:
271 preconditioner.reset(
272 Ifpack().Create((matrix.trilinos_matrix().NumMyRows() == 0) ?
273 "point relaxation" :
274 "block relaxation",
275 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
276 0));
277
278 Ifpack_Preconditioner *ifpack =
279 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
280 Assert(ifpack != nullptr,
281 ExcMessage("Trilinos could not create this "
282 "preconditioner"));
283
284 int ierr;
285
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",
292 additional_data.min_diagonal);
293 parameter_list.set("partitioner: type",
294 additional_data.block_creation_type);
295 int n_local_parts =
296 (matrix.trilinos_matrix().NumMyRows() + additional_data.block_size - 1) /
297 additional_data.block_size;
298 parameter_list.set("partitioner: local parts", n_local_parts);
299
300 ierr = ifpack->SetParameters(parameter_list);
301 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
302
303 ierr = ifpack->Initialize();
304 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
305
306 ierr = ifpack->Compute();
307 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
308 }
309
310
311
312 /* ----------------------- PreconditionBlockSSOR ------------------------ */
313
315 const unsigned int block_size,
316 const std::string &block_creation_type,
317 const double omega,
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)
323 , omega(omega)
324 , min_diagonal(min_diagonal)
325 , overlap(overlap)
326 , n_sweeps(n_sweeps)
327 {}
328
329
330
331 void
333 const AdditionalData &additional_data)
334 {
335 preconditioner.reset();
336
337 // Block relaxation setup fails if we have no locally owned rows. As a
338 // work-around we just pretend to use point relaxation on those processors:
339 preconditioner.reset(
340 Ifpack().Create((matrix.trilinos_matrix().NumMyRows() == 0) ?
341 "point relaxation" :
342 "block relaxation",
343 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
344 additional_data.overlap));
345
346 Ifpack_Preconditioner *ifpack =
347 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
348 Assert(ifpack != nullptr,
349 ExcMessage("Trilinos could not create this "
350 "preconditioner"));
351
352 int ierr;
353
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",
360 additional_data.min_diagonal);
361 parameter_list.set("schwarz: combine mode", "Add");
362 parameter_list.set("partitioner: type",
363 additional_data.block_creation_type);
364 int n_local_parts =
365 (matrix.trilinos_matrix().NumMyRows() + additional_data.block_size - 1) /
366 additional_data.block_size;
367 parameter_list.set("partitioner: local parts", n_local_parts);
368
369 ierr = ifpack->SetParameters(parameter_list);
370 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
371
372 ierr = ifpack->Initialize();
373 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
374
375 ierr = ifpack->Compute();
376 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
377 }
378
379
380
381 /* ------------------------ PreconditionBlockSOR ------------------------ */
382
384 const unsigned int block_size,
385 const std::string &block_creation_type,
386 const double omega,
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)
392 , omega(omega)
393 , min_diagonal(min_diagonal)
394 , overlap(overlap)
395 , n_sweeps(n_sweeps)
396 {}
397
398
399
400 void
402 const AdditionalData &additional_data)
403 {
404 preconditioner.reset();
405
406 // Block relaxation setup fails if we have no locally owned rows. As a
407 // work-around we just pretend to use point relaxation on those processors:
408 preconditioner.reset(
409 Ifpack().Create((matrix.trilinos_matrix().NumMyRows() == 0) ?
410 "point relaxation" :
411 "block relaxation",
412 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
413 additional_data.overlap));
414
415 Ifpack_Preconditioner *ifpack =
416 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
417 Assert(ifpack != nullptr,
418 ExcMessage("Trilinos could not create this "
419 "preconditioner"));
420
421 int ierr;
422
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",
429 additional_data.min_diagonal);
430 parameter_list.set("schwarz: combine mode", "Add");
431 parameter_list.set("partitioner: type",
432 additional_data.block_creation_type);
433 int n_local_parts =
434 (matrix.trilinos_matrix().NumMyRows() + additional_data.block_size - 1) /
435 additional_data.block_size;
436 parameter_list.set("partitioner: local parts", n_local_parts);
437
438 ierr = ifpack->SetParameters(parameter_list);
439 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
440
441 ierr = ifpack->Initialize();
442 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
443
444 ierr = ifpack->Compute();
445 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
446 }
447
448
449
450 /* -------------------------- PreconditionIC -------------------------- */
451
453 const double ic_atol,
454 const double ic_rtol,
455 const unsigned int overlap)
456 : ic_fill(ic_fill)
457 , ic_atol(ic_atol)
458 , ic_rtol(ic_rtol)
459 , overlap(overlap)
460 {}
461
462
463
464 void
466 const AdditionalData &additional_data)
467 {
468 preconditioner.reset();
469 preconditioner.reset(
470 Ifpack().Create("IC",
471 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
472 additional_data.overlap));
473
474 Ifpack_Preconditioner *ifpack =
475 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
476 Assert(ifpack != nullptr,
477 ExcMessage("Trilinos could not create this "
478 "preconditioner"));
479
480 int ierr;
481
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");
487
488 ierr = ifpack->SetParameters(parameter_list);
489 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
490
491 ierr = ifpack->Initialize();
492 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
493
494 ierr = ifpack->Compute();
495 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
496 }
497
498
499
500 /* -------------------------- PreconditionILU -------------------------- */
501
503 const double ilu_atol,
504 const double ilu_rtol,
505 const unsigned int overlap)
506 : ilu_fill(ilu_fill)
507 , ilu_atol(ilu_atol)
508 , ilu_rtol(ilu_rtol)
509 , overlap(overlap)
510 {}
511
512
513
514 void
516 const AdditionalData &additional_data)
517 {
518 preconditioner.reset();
519 preconditioner.reset(
520 Ifpack().Create("ILU",
521 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
522 additional_data.overlap));
523
524 Ifpack_Preconditioner *ifpack =
525 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
526 Assert(ifpack != nullptr,
527 ExcMessage("Trilinos could not create this "
528 "preconditioner"));
529
530 int ierr;
531
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");
538
539 ierr = ifpack->SetParameters(parameter_list);
540 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
541
542 ierr = ifpack->Initialize();
543 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
544
545 ierr = ifpack->Compute();
546 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
547 }
548
549
550
551 /* -------------------------- PreconditionILUT -------------------------- */
552
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)
562 , overlap(overlap)
563 {}
564
565
566
567 void
569 const AdditionalData &additional_data)
570 {
571 preconditioner.reset();
572 preconditioner.reset(
573 Ifpack().Create("ILUT",
574 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
575 additional_data.overlap));
576
577 Ifpack_Preconditioner *ifpack =
578 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
579 Assert(ifpack != nullptr,
580 ExcMessage("Trilinos could not create this "
581 "preconditioner"));
582
583 int ierr;
584
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");
592
593 ierr = ifpack->SetParameters(parameter_list);
594 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
595
596 ierr = ifpack->Initialize();
597 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
598
599 ierr = ifpack->Compute();
600 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
601 }
602
603
604
605 /* ---------------------- PreconditionBlockDirect --------------------- */
606
608 const unsigned int overlap)
609 : overlap(overlap)
610 {}
611
612
613
614 void
616 const AdditionalData &additional_data)
617 {
618 preconditioner.reset();
619 preconditioner.reset(
620 Ifpack().Create("Amesos",
621 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
622 additional_data.overlap));
623
624 Ifpack_Preconditioner *ifpack =
625 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
626 Assert(ifpack != nullptr,
627 ExcMessage("Trilinos could not create this "
628 "preconditioner"));
629
630 int ierr;
631
632 Teuchos::ParameterList parameter_list;
633 parameter_list.set("schwarz: combine mode", "Add");
634
635 ierr = ifpack->SetParameters(parameter_list);
636 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
637
638 ierr = ifpack->Initialize();
639 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
640
641 ierr = ifpack->Compute();
642 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
643 }
644
645
646
647 /* ---------------------- PreconditionBlockDirect --------------------- */
648
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)
656 : degree(degree)
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)
662 {}
663
664
665
666 void
668 const AdditionalData &additional_data)
669 {
671 Teuchos::rcp(new Ifpack_Chebyshev(&matrix.trilinos_matrix()));
672
673 Ifpack_Chebyshev *ifpack =
674 dynamic_cast<Ifpack_Chebyshev *>(preconditioner.get());
675 Assert(ifpack != nullptr,
676 ExcMessage("Trilinos could not create this "
677 "preconditioner"));
678
679 int ierr;
680
681 Teuchos::ParameterList parameter_list;
682 parameter_list.set("chebyshev: ratio eigenvalue",
683 additional_data.eigenvalue_ratio);
684 parameter_list.set("chebyshev: min eigenvalue",
685 additional_data.min_eigenvalue);
686 parameter_list.set("chebyshev: max eigenvalue",
687 additional_data.max_eigenvalue);
688 parameter_list.set("chebyshev: degree",
689 static_cast<int>(additional_data.degree));
690 parameter_list.set("chebyshev: min diagonal value",
691 additional_data.min_diagonal);
692 parameter_list.set("chebyshev: zero starting solution",
693 !additional_data.nonzero_starting);
694
695 ierr = ifpack->SetParameters(parameter_list);
696 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
697
698 ierr = ifpack->Initialize();
699 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
700
701 ierr = ifpack->Compute();
702 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
703 }
704
705
706
707 /* -------------------------- PreconditionIdentity --------------------- */
708
709 void
711 const AdditionalData &)
712 {
713 // What follows just configures a dummy preconditioner that
714 // sets up the domain and range maps, as well as the communicator.
715 // It is never used as the vmult, Tvmult operations are
716 // given a custom definition.
717 // Note: This is only required in order to wrap this
718 // preconditioner in a LinearOperator without an exemplar
719 // matrix.
720
721 // From PreconditionJacobi:
722 // release memory before reallocation
723 preconditioner.reset();
724 preconditioner.reset(
725 Ifpack().Create("point relaxation",
726 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
727 0));
728
729 Ifpack_Preconditioner *ifpack =
730 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
731 Assert(ifpack != nullptr,
732 ExcMessage("Trilinos could not create this "
733 "preconditioner"));
734
735 int ierr;
736
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);
742
743 ierr = ifpack->SetParameters(parameter_list);
744 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
745
746 ierr = ifpack->Initialize();
747 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
748
749 ierr = ifpack->Compute();
750 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
751 }
752
753 void
755 {
756 dst = src;
757 }
758
759 void
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# ifndef DOXYGEN
780 void
784 {
785 dst = src;
786 }
787
788 void
792 {
793 dst = src;
794 }
795# endif // DOXYGEN
796} // namespace TrilinosWrappers
797
799
800#endif // DEAL_II_WITH_TRILINOS
Epetra_Operator & trilinos_operator() const
std::shared_ptr< Epetra_Map > vector_distributor
Teuchos::RCP< Epetra_Operator > preconditioner
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
Definition config.h:472
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
#define Assert(cond, exc)
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
AdditionalData(const unsigned int block_size=1, const std::string &block_creation_type="linear", const double omega=1, const double min_diagonal=0, const unsigned int n_sweeps=1)
AdditionalData(const unsigned int 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)
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)