deal.II version GIT relicensing-2250-g88cb8ba3cb 2024-12-13 12:20:00+00:00
\(\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// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2008 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
16
17#ifdef DEAL_II_WITH_TRILINOS
18
21# include <deal.II/lac/vector.h>
22
23# include <Epetra_MultiVector.h>
24# include <Ifpack.h>
25# include <Ifpack_Chebyshev.h>
26# include <Teuchos_ParameterList.hpp>
27# include <Teuchos_RCP.hpp>
28
30
31namespace TrilinosWrappers
32{
36
37
38 void
44
45
48 {
49 return communicator.Comm();
50 }
51
52
55 {
56 AssertThrow(!preconditioner.is_null(),
57 ExcMessage("Trying to dereference a null pointer."));
58 return (*preconditioner);
59 }
60
61
64 {
65 return IndexSet(preconditioner->OperatorDomainMap());
66 }
67
68
71 {
72 return IndexSet(preconditioner->OperatorRangeMap());
73 }
74
75 /* -------------------------- PreconditionJacobi -------------------------- */
76
78 const double omega,
79 const double min_diagonal,
80 const unsigned int n_sweeps)
81 : omega(omega)
82 , min_diagonal(min_diagonal)
83 , n_sweeps(n_sweeps)
84 {}
85
86
87
88 void
90 const AdditionalData &additional_data)
91 {
92 // release memory before reallocation
93 preconditioner.reset();
94 preconditioner.reset(
95 Ifpack().Create("point relaxation",
96 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
97 0));
98
100 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
101 Assert(ifpack != nullptr,
102 ExcMessage("Trilinos could not create this "
103 "preconditioner"));
104
105 int ierr;
106
107 Teuchos::ParameterList parameter_list;
108 parameter_list.set("relaxation: sweeps",
109 static_cast<int>(additional_data.n_sweeps));
110 parameter_list.set("relaxation: type", "Jacobi");
111 parameter_list.set("relaxation: damping factor", additional_data.omega);
112 parameter_list.set("relaxation: min diagonal value",
113 additional_data.min_diagonal);
114
115 ierr = ifpack->SetParameters(parameter_list);
117
118 ierr = ifpack->Initialize();
120
121 ierr = ifpack->Compute();
123 }
124
125
126
127 /* -------------------------- PreconditionSSOR -------------------------- */
128
130 const double min_diagonal,
131 const unsigned int overlap,
132 const unsigned int n_sweeps)
133 : omega(omega)
134 , min_diagonal(min_diagonal)
135 , overlap(overlap)
136 , n_sweeps(n_sweeps)
137 {}
138
139
140
141 void
143 const AdditionalData &additional_data)
144 {
145 preconditioner.reset();
146 preconditioner.reset(
147 Ifpack().Create("point relaxation",
148 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
149 additional_data.overlap));
150
152 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
153 Assert(ifpack != nullptr,
154 ExcMessage("Trilinos could not create this "
155 "preconditioner"));
156
157 int ierr;
158
159 Teuchos::ParameterList parameter_list;
160 parameter_list.set("relaxation: sweeps",
161 static_cast<int>(additional_data.n_sweeps));
162 parameter_list.set("relaxation: type", "symmetric Gauss-Seidel");
163 parameter_list.set("relaxation: damping factor", additional_data.omega);
164 parameter_list.set("relaxation: min diagonal value",
165 additional_data.min_diagonal);
166 parameter_list.set("schwarz: combine mode", "Add");
167
168 ierr = ifpack->SetParameters(parameter_list);
170
171 ierr = ifpack->Initialize();
173
174 ierr = ifpack->Compute();
176 }
177
178
179
180 /* -------------------------- PreconditionSOR -------------------------- */
181
183 const double min_diagonal,
184 const unsigned int overlap,
185 const unsigned int n_sweeps)
186 : omega(omega)
187 , min_diagonal(min_diagonal)
188 , overlap(overlap)
189 , n_sweeps(n_sweeps)
190 {}
191
192
193
194 void
196 const AdditionalData &additional_data)
197 {
198 preconditioner.reset();
199 preconditioner.reset(
200 Ifpack().Create("point relaxation",
201 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
202 additional_data.overlap));
203
205 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
206 Assert(ifpack != nullptr,
207 ExcMessage("Trilinos could not create this "
208 "preconditioner"));
209
210 int ierr;
211
212 Teuchos::ParameterList parameter_list;
213 parameter_list.set("relaxation: sweeps",
214 static_cast<int>(additional_data.n_sweeps));
215 parameter_list.set("relaxation: type", "Gauss-Seidel");
216 parameter_list.set("relaxation: damping factor", additional_data.omega);
217 parameter_list.set("relaxation: min diagonal value",
218 additional_data.min_diagonal);
219 parameter_list.set("schwarz: combine mode", "Add");
220
221 ierr = ifpack->SetParameters(parameter_list);
223
224 ierr = ifpack->Initialize();
226
227 ierr = ifpack->Compute();
229 }
230
231
232
233 /* ----------------------- PreconditionBlockJacobi ---------------------- */
234
236 const unsigned int block_size,
237 const std::string &block_creation_type,
238 const double omega,
239 const double min_diagonal,
240 const unsigned int n_sweeps)
241 : block_size(block_size)
242 , block_creation_type(block_creation_type)
243 , omega(omega)
244 , min_diagonal(min_diagonal)
245 , n_sweeps(n_sweeps)
246 {}
247
248
249
250 void
252 const AdditionalData &additional_data)
253 {
254 // release memory before reallocation
255 preconditioner.reset();
256
257 // Block relaxation setup fails if we have no locally owned rows. As a
258 // work-around we just pretend to use point relaxation on those processors:
259 preconditioner.reset(
260 Ifpack().Create((matrix.trilinos_matrix().NumMyRows() == 0) ?
261 "point relaxation" :
262 "block relaxation",
263 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
264 0));
265
267 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
268 Assert(ifpack != nullptr,
269 ExcMessage("Trilinos could not create this "
270 "preconditioner"));
271
272 int ierr;
273
274 Teuchos::ParameterList parameter_list;
275 parameter_list.set("relaxation: sweeps",
276 static_cast<int>(additional_data.n_sweeps));
277 parameter_list.set("relaxation: type", "Jacobi");
278 parameter_list.set("relaxation: damping factor", additional_data.omega);
279 parameter_list.set("relaxation: min diagonal value",
280 additional_data.min_diagonal);
281 parameter_list.set("partitioner: type",
282 additional_data.block_creation_type);
283 int n_local_parts =
284 (matrix.trilinos_matrix().NumMyRows() + additional_data.block_size - 1) /
285 additional_data.block_size;
286 parameter_list.set("partitioner: local parts", n_local_parts);
287
288 ierr = ifpack->SetParameters(parameter_list);
290
291 ierr = ifpack->Initialize();
293
294 ierr = ifpack->Compute();
296 }
297
298
299
300 /* ----------------------- PreconditionBlockSSOR ------------------------ */
301
303 const unsigned int block_size,
304 const std::string &block_creation_type,
305 const double omega,
306 const double min_diagonal,
307 const unsigned int overlap,
308 const unsigned int n_sweeps)
309 : block_size(block_size)
310 , block_creation_type(block_creation_type)
311 , omega(omega)
312 , min_diagonal(min_diagonal)
313 , overlap(overlap)
314 , n_sweeps(n_sweeps)
315 {}
316
317
318
319 void
321 const AdditionalData &additional_data)
322 {
323 preconditioner.reset();
324
325 // Block relaxation setup fails if we have no locally owned rows. As a
326 // work-around we just pretend to use point relaxation on those processors:
327 preconditioner.reset(
328 Ifpack().Create((matrix.trilinos_matrix().NumMyRows() == 0) ?
329 "point relaxation" :
330 "block relaxation",
331 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
332 additional_data.overlap));
333
335 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
336 Assert(ifpack != nullptr,
337 ExcMessage("Trilinos could not create this "
338 "preconditioner"));
339
340 int ierr;
341
342 Teuchos::ParameterList parameter_list;
343 parameter_list.set("relaxation: sweeps",
344 static_cast<int>(additional_data.n_sweeps));
345 parameter_list.set("relaxation: type", "symmetric Gauss-Seidel");
346 parameter_list.set("relaxation: damping factor", additional_data.omega);
347 parameter_list.set("relaxation: min diagonal value",
348 additional_data.min_diagonal);
349 parameter_list.set("schwarz: combine mode", "Add");
350 parameter_list.set("partitioner: type",
351 additional_data.block_creation_type);
352 int n_local_parts =
353 (matrix.trilinos_matrix().NumMyRows() + additional_data.block_size - 1) /
354 additional_data.block_size;
355 parameter_list.set("partitioner: local parts", n_local_parts);
356
357 ierr = ifpack->SetParameters(parameter_list);
359
360 ierr = ifpack->Initialize();
362
363 ierr = ifpack->Compute();
365 }
366
367
368
369 /* ------------------------ PreconditionBlockSOR ------------------------ */
370
372 const unsigned int block_size,
373 const std::string &block_creation_type,
374 const double omega,
375 const double min_diagonal,
376 const unsigned int overlap,
377 const unsigned int n_sweeps)
378 : block_size(block_size)
379 , block_creation_type(block_creation_type)
380 , omega(omega)
381 , min_diagonal(min_diagonal)
382 , overlap(overlap)
383 , n_sweeps(n_sweeps)
384 {}
385
386
387
388 void
390 const AdditionalData &additional_data)
391 {
392 preconditioner.reset();
393
394 // Block relaxation setup fails if we have no locally owned rows. As a
395 // work-around we just pretend to use point relaxation on those processors:
396 preconditioner.reset(
397 Ifpack().Create((matrix.trilinos_matrix().NumMyRows() == 0) ?
398 "point relaxation" :
399 "block relaxation",
400 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
401 additional_data.overlap));
402
404 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
405 Assert(ifpack != nullptr,
406 ExcMessage("Trilinos could not create this "
407 "preconditioner"));
408
409 int ierr;
410
411 Teuchos::ParameterList parameter_list;
412 parameter_list.set("relaxation: sweeps",
413 static_cast<int>(additional_data.n_sweeps));
414 parameter_list.set("relaxation: type", "Gauss-Seidel");
415 parameter_list.set("relaxation: damping factor", additional_data.omega);
416 parameter_list.set("relaxation: min diagonal value",
417 additional_data.min_diagonal);
418 parameter_list.set("schwarz: combine mode", "Add");
419 parameter_list.set("partitioner: type",
420 additional_data.block_creation_type);
421 int n_local_parts =
422 (matrix.trilinos_matrix().NumMyRows() + additional_data.block_size - 1) /
423 additional_data.block_size;
424 parameter_list.set("partitioner: local parts", n_local_parts);
425
426 ierr = ifpack->SetParameters(parameter_list);
428
429 ierr = ifpack->Initialize();
431
432 ierr = ifpack->Compute();
434 }
435
436
437
438 /* -------------------------- PreconditionIC -------------------------- */
439
441 const double ic_atol,
442 const double ic_rtol,
443 const unsigned int overlap)
444 : ic_fill(ic_fill)
445 , ic_atol(ic_atol)
446 , ic_rtol(ic_rtol)
447 , overlap(overlap)
448 {}
449
450
451
452 void
454 const AdditionalData &additional_data)
455 {
456 preconditioner.reset();
457 preconditioner.reset(
458 Ifpack().Create("IC",
459 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
460 additional_data.overlap));
461
463 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
464 Assert(ifpack != nullptr,
465 ExcMessage("Trilinos could not create this "
466 "preconditioner"));
467
468 int ierr;
469
470 Teuchos::ParameterList parameter_list;
471 parameter_list.set("fact: level-of-fill", additional_data.ic_fill);
472 parameter_list.set("fact: absolute threshold", additional_data.ic_atol);
473 parameter_list.set("fact: relative threshold", additional_data.ic_rtol);
474 parameter_list.set("schwarz: combine mode", "Add");
475
476 ierr = ifpack->SetParameters(parameter_list);
478
479 ierr = ifpack->Initialize();
481
482 ierr = ifpack->Compute();
484 }
485
486
487
488 /* -------------------------- PreconditionILU -------------------------- */
489
491 const double ilu_atol,
492 const double ilu_rtol,
493 const unsigned int overlap)
494 : ilu_fill(ilu_fill)
495 , ilu_atol(ilu_atol)
496 , ilu_rtol(ilu_rtol)
497 , overlap(overlap)
498 {}
499
500
501
502 void
504 const AdditionalData &additional_data)
505 {
506 preconditioner.reset();
507 preconditioner.reset(
508 Ifpack().Create("ILU",
509 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
510 additional_data.overlap));
511
513 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
514 Assert(ifpack != nullptr,
515 ExcMessage("Trilinos could not create this "
516 "preconditioner"));
517
518 int ierr;
519
520 Teuchos::ParameterList parameter_list;
521 parameter_list.set("fact: level-of-fill",
522 static_cast<int>(additional_data.ilu_fill));
523 parameter_list.set("fact: absolute threshold", additional_data.ilu_atol);
524 parameter_list.set("fact: relative threshold", additional_data.ilu_rtol);
525 parameter_list.set("schwarz: combine mode", "Add");
526
527 ierr = ifpack->SetParameters(parameter_list);
529
530 ierr = ifpack->Initialize();
532
533 ierr = ifpack->Compute();
535 }
536
537
538
539 /* -------------------------- PreconditionILUT -------------------------- */
540
542 const unsigned int ilut_fill,
543 const double ilut_atol,
544 const double ilut_rtol,
545 const unsigned int overlap)
546 : ilut_drop(ilut_drop)
547 , ilut_fill(ilut_fill)
548 , ilut_atol(ilut_atol)
549 , ilut_rtol(ilut_rtol)
550 , overlap(overlap)
551 {}
552
553
554
555 void
557 const AdditionalData &additional_data)
558 {
559 preconditioner.reset();
560 preconditioner.reset(
561 Ifpack().Create("ILUT",
562 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
563 additional_data.overlap));
564
566 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
567 Assert(ifpack != nullptr,
568 ExcMessage("Trilinos could not create this "
569 "preconditioner"));
570
571 int ierr;
572
573 Teuchos::ParameterList parameter_list;
574 parameter_list.set("fact: drop value", additional_data.ilut_drop);
575 parameter_list.set("fact: level-of-fill",
576 static_cast<int>(additional_data.ilut_fill));
577 parameter_list.set("fact: absolute threshold", additional_data.ilut_atol);
578 parameter_list.set("fact: relative threshold", additional_data.ilut_rtol);
579 parameter_list.set("schwarz: combine mode", "Add");
580
581 ierr = ifpack->SetParameters(parameter_list);
583
584 ierr = ifpack->Initialize();
586
587 ierr = ifpack->Compute();
589 }
590
591
592
593 /* ---------------------- PreconditionBlockDirect --------------------- */
594
596 const unsigned int overlap)
597 : overlap(overlap)
598 {}
599
600
601
602 void
604 const AdditionalData &additional_data)
605 {
606 preconditioner.reset();
607 preconditioner.reset(
608 Ifpack().Create("Amesos",
609 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
610 additional_data.overlap));
611
613 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
614 Assert(ifpack != nullptr,
615 ExcMessage("Trilinos could not create this "
616 "preconditioner"));
617
618 int ierr;
619
620 Teuchos::ParameterList parameter_list;
621 parameter_list.set("schwarz: combine mode", "Add");
622
623 ierr = ifpack->SetParameters(parameter_list);
625
626 ierr = ifpack->Initialize();
628
629 ierr = ifpack->Compute();
631 }
632
633
634
635 /* ---------------------- PreconditionBlockDirect --------------------- */
636
638 const unsigned int degree,
639 const double max_eigenvalue,
640 const double eigenvalue_ratio,
641 const double min_eigenvalue,
642 const double min_diagonal,
643 const bool nonzero_starting)
644 : degree(degree)
645 , max_eigenvalue(max_eigenvalue)
646 , eigenvalue_ratio(eigenvalue_ratio)
647 , min_eigenvalue(min_eigenvalue)
648 , min_diagonal(min_diagonal)
649 , nonzero_starting(nonzero_starting)
650 {}
651
652
653
654 void
656 const AdditionalData &additional_data)
657 {
659 Teuchos::rcp(new Ifpack_Chebyshev(&matrix.trilinos_matrix()));
660
662 dynamic_cast<Ifpack_Chebyshev *>(preconditioner.get());
663 Assert(ifpack != nullptr,
664 ExcMessage("Trilinos could not create this "
665 "preconditioner"));
666
667 int ierr;
668
669 Teuchos::ParameterList parameter_list;
670 parameter_list.set("chebyshev: ratio eigenvalue",
671 additional_data.eigenvalue_ratio);
672 parameter_list.set("chebyshev: min eigenvalue",
673 additional_data.min_eigenvalue);
674 parameter_list.set("chebyshev: max eigenvalue",
675 additional_data.max_eigenvalue);
676 parameter_list.set("chebyshev: degree",
677 static_cast<int>(additional_data.degree));
678 parameter_list.set("chebyshev: min diagonal value",
679 additional_data.min_diagonal);
680 parameter_list.set("chebyshev: zero starting solution",
681 !additional_data.nonzero_starting);
682
683 ierr = ifpack->SetParameters(parameter_list);
685
686 ierr = ifpack->Initialize();
688
689 ierr = ifpack->Compute();
691 }
692
693
694
695 /* -------------------------- PreconditionIdentity --------------------- */
696
697 void
699 const AdditionalData &)
700 {
701 // What follows just configures a dummy preconditioner that
702 // sets up the domain and range maps, as well as the communicator.
703 // It is never used as the vmult, Tvmult operations are
704 // given a custom definition.
705 // Note: This is only required in order to wrap this
706 // preconditioner in a LinearOperator without an exemplar
707 // matrix.
708
709 // From PreconditionJacobi:
710 // release memory before reallocation
711 preconditioner.reset();
712 preconditioner.reset(
713 Ifpack().Create("point relaxation",
714 const_cast<Epetra_CrsMatrix *>(&matrix.trilinos_matrix()),
715 0));
716
718 dynamic_cast<Ifpack_Preconditioner *>(preconditioner.get());
719 Assert(ifpack != nullptr,
720 ExcMessage("Trilinos could not create this "
721 "preconditioner"));
722
723 int ierr;
724
725 Teuchos::ParameterList parameter_list;
726 parameter_list.set("relaxation: sweeps", 1);
727 parameter_list.set("relaxation: type", "Jacobi");
728 parameter_list.set("relaxation: damping factor", 1.0);
729 parameter_list.set("relaxation: min diagonal value", 0.0);
730
731 ierr = ifpack->SetParameters(parameter_list);
733
734 ierr = ifpack->Initialize();
736
737 ierr = ifpack->Compute();
739 }
740
741 void
743 {
744 dst = src;
745 }
746
747 void
749 {
750 dst = src;
751 }
752
753 void
755 const ::Vector<double> &src) const
756 {
757 dst = src;
758 }
759
760 void
762 const ::Vector<double> &src) const
763 {
764 dst = src;
765 }
766
767# ifndef DOXYGEN
768 void
772 {
773 dst = src;
774 }
775
776 void
780 {
781 dst = src;
782 }
783# endif // DOXYGEN
784} // namespace TrilinosWrappers
785
787
788#endif // DEAL_II_WITH_TRILINOS
friend class Tensor
Definition tensor.h:865
Epetra_Operator & trilinos_operator() const
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:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
#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)