deal.II version 9.7.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_solver.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2008 - 2025 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
20
24
26
27# include <AztecOO_StatusTest.h>
28# include <AztecOO_StatusTestCombo.h>
29# include <AztecOO_StatusTestMaxIters.h>
30# include <AztecOO_StatusTestResNorm.h>
31# include <AztecOO_StatusType.h>
32
34
35# include <cmath>
36# include <limits>
37# include <memory>
38
40
41namespace TrilinosWrappers
42{
49
50
51
57
58
59
67
68
69
72 {
73 return solver_control;
74 }
75
76
77
78 void
80 MPI::Vector &x,
81 const MPI::Vector &b,
82 const PreconditionBase &preconditioner)
83 {
84 // We need an Epetra_LinearProblem object to let the AztecOO solver know
85 // about the matrix and vectors.
86 linear_problem = std::make_unique<Epetra_LinearProblem>(
87 const_cast<Epetra_CrsMatrix *>(&A.trilinos_matrix()),
88 &x.trilinos_vector(),
89 const_cast<Epetra_MultiVector *>(&b.trilinos_vector()));
90
91 do_solve(preconditioner);
92 }
93
94
95
96 // Note: "A" is set as a constant reference so that all patterns for ::solve
97 // can be used by the inverse_operator of LinearOperator
98 void
100 MPI::Vector &x,
101 const MPI::Vector &b,
102 const PreconditionBase &preconditioner)
103 {
104 // We need an Epetra_LinearProblem object to let the AztecOO solver know
105 // about the matrix and vectors.
107 std::make_unique<Epetra_LinearProblem>(const_cast<Epetra_Operator *>(&A),
108 &x.trilinos_vector(),
109 const_cast<Epetra_MultiVector *>(
110 &b.trilinos_vector()));
111
112 do_solve(preconditioner);
113 }
114
115
116
117 // Note: "A" is set as a constant reference so that all patterns for ::solve
118 // can be used by the inverse_operator of LinearOperator
119 void
121 MPI::Vector &x,
122 const MPI::Vector &b,
123 const Epetra_Operator &preconditioner)
124 {
125 // We need an Epetra_LinearProblem object to let the AztecOO solver know
126 // about the matrix and vectors.
128 std::make_unique<Epetra_LinearProblem>(const_cast<Epetra_Operator *>(&A),
129 &x.trilinos_vector(),
130 const_cast<Epetra_MultiVector *>(
131 &b.trilinos_vector()));
132
133 do_solve(preconditioner);
134 }
135
136
137
138 // Note: "A" is set as a constant reference so that all patterns for ::solve
139 // can be used by the inverse_operator of LinearOperator
140 void
142 Epetra_MultiVector &x,
143 const Epetra_MultiVector &b,
144 const PreconditionBase &preconditioner)
145 {
146 // We need an Epetra_LinearProblem object to let the AztecOO solver know
147 // about the matrix and vectors.
149 std::make_unique<Epetra_LinearProblem>(const_cast<Epetra_Operator *>(&A),
150 &x,
151 const_cast<Epetra_MultiVector *>(
152 &b));
153
154 do_solve(preconditioner);
155 }
156
157
158
159 // Note: "A" is set as a constant reference so that all patterns for ::solve
160 // can be used by the inverse_operator of LinearOperator
161 void
163 Epetra_MultiVector &x,
164 const Epetra_MultiVector &b,
165 const Epetra_Operator &preconditioner)
166 {
167 // We need an Epetra_LinearProblem object to let the AztecOO solver know
168 // about the matrix and vectors.
170 std::make_unique<Epetra_LinearProblem>(const_cast<Epetra_Operator *>(&A),
171 &x,
172 const_cast<Epetra_MultiVector *>(
173 &b));
174
175 do_solve(preconditioner);
176 }
177
178
179
180 void
183 const ::Vector<double> &b,
184 const PreconditionBase &preconditioner)
185 {
186 // In case we call the solver with deal.II vectors, we create views of the
187 // vectors in Epetra format.
188 Assert(x.size() == A.n(), ExcDimensionMismatch(x.size(), A.n()));
189 Assert(b.size() == A.m(), ExcDimensionMismatch(b.size(), A.m()));
190 Assert(A.local_range().second == A.m(),
191 ExcMessage("Can only work in serial when using deal.II vectors."));
192 Assert(A.trilinos_matrix().Filled(),
193 ExcMessage("Matrix is not compressed. Call compress() method."));
194
195 Epetra_Vector ep_x(View, A.trilinos_matrix().DomainMap(), x.begin());
196 Epetra_Vector ep_b(View,
197 A.trilinos_matrix().RangeMap(),
198 const_cast<double *>(b.begin()));
199
200 // We need an Epetra_LinearProblem object to let the AztecOO solver know
201 // about the matrix and vectors.
202 linear_problem = std::make_unique<Epetra_LinearProblem>(
203 const_cast<Epetra_CrsMatrix *>(&A.trilinos_matrix()), &ep_x, &ep_b);
204
205 do_solve(preconditioner);
206 }
207
208
209
210 void
213 const ::Vector<double> &b,
214 const PreconditionBase &preconditioner)
215 {
216 Epetra_Vector ep_x(View, A.OperatorDomainMap(), x.begin());
217 Epetra_Vector ep_b(View,
218 A.OperatorRangeMap(),
219 const_cast<double *>(b.begin()));
220
221 // We need an Epetra_LinearProblem object to let the AztecOO solver know
222 // about the matrix and vectors.
223 linear_problem = std::make_unique<Epetra_LinearProblem>(&A, &ep_x, &ep_b);
224
225 do_solve(preconditioner);
226 }
227
228
229
230 void
233 const ::LinearAlgebra::distributed::Vector<double> &b,
234 const PreconditionBase &preconditioner)
235 {
236 // In case we call the solver with deal.II vectors, we create views of the
237 // vectors in Epetra format.
239 A.trilinos_matrix().DomainMap().NumMyElements());
240 AssertDimension(b.locally_owned_size(),
241 A.trilinos_matrix().RangeMap().NumMyElements());
242
243 Epetra_Vector ep_x(View, A.trilinos_matrix().DomainMap(), x.begin());
244 Epetra_Vector ep_b(View,
245 A.trilinos_matrix().RangeMap(),
246 const_cast<double *>(b.begin()));
247
248 // We need an Epetra_LinearProblem object to let the AztecOO solver know
249 // about the matrix and vectors.
250 linear_problem = std::make_unique<Epetra_LinearProblem>(
251 const_cast<Epetra_CrsMatrix *>(&A.trilinos_matrix()), &ep_x, &ep_b);
252
253 do_solve(preconditioner);
254 }
255
256
257
258 void
261 const ::LinearAlgebra::distributed::Vector<double> &b,
262 const PreconditionBase &preconditioner)
263 {
265 A.OperatorDomainMap().NumMyElements());
266 AssertDimension(b.locally_owned_size(),
267 A.OperatorRangeMap().NumMyElements());
268
269 Epetra_Vector ep_x(View, A.OperatorDomainMap(), x.begin());
270 Epetra_Vector ep_b(View,
271 A.OperatorRangeMap(),
272 const_cast<double *>(b.begin()));
273
274 // We need an Epetra_LinearProblem object to let the AztecOO solver know
275 // about the matrix and vectors.
276 linear_problem = std::make_unique<Epetra_LinearProblem>(&A, &ep_x, &ep_b);
277
278 do_solve(preconditioner);
279 }
280
281
282 namespace internal
283 {
284 namespace
285 {
286 double
287 compute_residual(const Epetra_MultiVector *const residual_vector)
288 {
289 Assert(residual_vector->NumVectors() == 1,
290 ExcMessage("Residual multivector holds more than one vector"));
291 TrilinosScalar res_l2_norm = 0.0;
292 const int ierr = residual_vector->Norm2(&res_l2_norm);
293 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
294 return res_l2_norm;
295 }
296
297 class TrilinosReductionControl : public AztecOO_StatusTest
298 {
299 public:
300 TrilinosReductionControl(const int max_steps,
301 const double tolerance,
302 const double reduction,
303 const Epetra_LinearProblem &linear_problem);
304
305 virtual ~TrilinosReductionControl() override = default;
306
307 virtual bool
308 ResidualVectorRequired() const override
309 {
310 return status_test_collection->ResidualVectorRequired();
311 }
312
313 virtual AztecOO_StatusType
314 CheckStatus(int CurrentIter,
315 Epetra_MultiVector *CurrentResVector,
316 double CurrentResNormEst,
317 bool SolutionUpdated) override
318 {
319 // Note: CurrentResNormEst is set to -1.0 if no estimate of the
320 // residual value is available
321 current_residual =
322 (CurrentResNormEst < 0.0 ? compute_residual(CurrentResVector) :
323 CurrentResNormEst);
324 if (CurrentIter == 0)
325 initial_residual = current_residual;
326
327 return status_test_collection->CheckStatus(CurrentIter,
328 CurrentResVector,
329 CurrentResNormEst,
330 SolutionUpdated);
331 }
332
333 virtual AztecOO_StatusType
334 GetStatus() const override
335 {
336 return status_test_collection->GetStatus();
337 }
338
339 virtual std::ostream &
340 Print(std::ostream &stream, int indent = 0) const override
341 {
342 return status_test_collection->Print(stream, indent);
343 }
344
345 double
346 get_initial_residual() const
347 {
348 return initial_residual;
349 }
350
351 double
352 get_current_residual() const
353 {
354 return current_residual;
355 }
356
357 private:
358 double initial_residual;
359 double current_residual;
360 std::unique_ptr<AztecOO_StatusTestCombo> status_test_collection;
361 std::unique_ptr<AztecOO_StatusTestMaxIters> status_test_max_steps;
362 std::unique_ptr<AztecOO_StatusTestResNorm> status_test_abs_tol;
363 std::unique_ptr<AztecOO_StatusTestResNorm> status_test_rel_tol;
364 };
365
366
367 TrilinosReductionControl::TrilinosReductionControl(
368 const int max_steps,
369 const double tolerance,
370 const double reduction,
371 const Epetra_LinearProblem &linear_problem)
372 : initial_residual(std::numeric_limits<double>::max())
373 , current_residual(std::numeric_limits<double>::max())
374 // Consider linear problem converged if any of the collection of
375 // criterion are met
376 , status_test_collection(std::make_unique<AztecOO_StatusTestCombo>(
377 AztecOO_StatusTestCombo::OR))
378 {
379 // Maximum number of iterations
380 Assert(max_steps >= 0, ExcInternalError());
381 status_test_max_steps =
382 std::make_unique<AztecOO_StatusTestMaxIters>(max_steps);
383 status_test_collection->AddStatusTest(*status_test_max_steps);
384
385 Assert(linear_problem.GetRHS()->NumVectors() == 1,
386 ExcMessage("RHS multivector holds more than one vector"));
387
388 // Residual norm is below some absolute value
389 status_test_abs_tol = std::make_unique<AztecOO_StatusTestResNorm>(
390 *linear_problem.GetOperator(),
391 *(linear_problem.GetLHS()->operator()(0)),
392 *(linear_problem.GetRHS()->operator()(0)),
393 tolerance);
394 status_test_abs_tol->DefineResForm(AztecOO_StatusTestResNorm::Explicit,
395 AztecOO_StatusTestResNorm::TwoNorm);
396 status_test_abs_tol->DefineScaleForm(
397 AztecOO_StatusTestResNorm::None, AztecOO_StatusTestResNorm::TwoNorm);
398 status_test_collection->AddStatusTest(*status_test_abs_tol);
399
400 // Residual norm, scaled by some initial value, is below some threshold
401 status_test_rel_tol = std::make_unique<AztecOO_StatusTestResNorm>(
402 *linear_problem.GetOperator(),
403 *(linear_problem.GetLHS()->operator()(0)),
404 *(linear_problem.GetRHS()->operator()(0)),
405 reduction);
406 status_test_rel_tol->DefineResForm(AztecOO_StatusTestResNorm::Explicit,
407 AztecOO_StatusTestResNorm::TwoNorm);
408 status_test_rel_tol->DefineScaleForm(
409 AztecOO_StatusTestResNorm::NormOfInitRes,
410 AztecOO_StatusTestResNorm::TwoNorm);
411 status_test_collection->AddStatusTest(*status_test_rel_tol);
412 }
413
414 } // namespace
415 } // namespace internal
416
417
418 template <typename Preconditioner>
419 void
420 SolverBase::do_solve(const Preconditioner &preconditioner)
421 {
422 int ierr;
423
424 // Next we can allocate the AztecOO solver...
425 solver.SetProblem(*linear_problem);
426
427 // ... and we can specify the solver to be used.
428 switch (solver_name)
429 {
430 case cg:
431 solver.SetAztecOption(AZ_solver, AZ_cg);
432 break;
433 case cgs:
434 solver.SetAztecOption(AZ_solver, AZ_cgs);
435 break;
436 case gmres:
437 solver.SetAztecOption(AZ_solver, AZ_gmres);
438 solver.SetAztecOption(AZ_kspace,
439 additional_data.gmres_restart_parameter);
440 break;
441 case bicgstab:
442 solver.SetAztecOption(AZ_solver, AZ_bicgstab);
443 break;
444 case tfqmr:
445 solver.SetAztecOption(AZ_solver, AZ_tfqmr);
446 break;
447 default:
449 }
450
451 // Set the preconditioner
452 set_preconditioner(solver, preconditioner);
453
454 // ... set some options, ...
455 solver.SetAztecOption(AZ_output,
456 additional_data.output_solver_details ? AZ_all :
457 AZ_none);
458 solver.SetAztecOption(AZ_conv, AZ_noscaled);
459
460 // By default, the Trilinos solver chooses convergence criterion based on
461 // the number of iterations made and an absolute tolerance.
462 // This implies that the use of the standard Trilinos convergence test
463 // actually coincides with ::IterationNumberControl because the
464 // solver, unless explicitly told otherwise, will Iterate() until a number
465 // of max_steps() are taken or an absolute tolerance() is attained.
466 // It is therefore suitable for use with both SolverControl or
467 // IterationNumberControl. The final check at the end will determine whether
468 // failure to converge to the defined residual norm constitutes failure
469 // (SolverControl) or is alright (IterationNumberControl).
470 // In the case that the SolverControl wants to perform ReductionControl,
471 // then we have to do a little extra something by prescribing a custom
472 // status test.
473 if (!status_test)
474 {
475 if (const ReductionControl *const reduction_control =
476 dynamic_cast<const ReductionControl *>(&solver_control))
477 {
478 status_test = std::make_unique<internal::TrilinosReductionControl>(
479 reduction_control->max_steps(),
480 reduction_control->tolerance(),
481 reduction_control->reduction(),
483 solver.SetStatusTest(status_test.get());
484 }
485 }
486
487 // ... and then solve!
488 ierr =
489 solver.Iterate(solver_control.max_steps(), solver_control.tolerance());
490
491 // report errors in more detail than just by checking whether the return
492 // status is zero or greater. the error strings are taken from the
493 // implementation of the AztecOO::Iterate function
494 switch (ierr)
495 {
496 case -1:
497 AssertThrow(false,
498 ExcMessage("AztecOO::Iterate error code -1: "
499 "option not implemented"));
500 break;
501 case -2:
502 AssertThrow(false,
503 ExcMessage("AztecOO::Iterate error code -2: "
504 "numerical breakdown"));
505 break;
506 case -3:
507 AssertThrow(false,
508 ExcMessage("AztecOO::Iterate error code -3: "
509 "loss of precision"));
510 break;
511 case -4:
512 AssertThrow(false,
513 ExcMessage("AztecOO::Iterate error code -4: "
514 "GMRES Hessenberg ill-conditioned"));
515 break;
516 default:
517 AssertThrow(ierr >= 0, ExcTrilinosError(ierr));
518 }
519
520 // Finally, let the deal.II SolverControl object know what has
521 // happened. If the solve succeeded, the status of the solver control will
522 // turn into SolverControl::success.
523 // If the residual is not computed/stored by the solver, as can happen for
524 // certain choices of solver or if a custom status test is set, then the
525 // result returned by TrueResidual() is equal to -1. In this case we must
526 // compute it ourself.
527 if (const internal::TrilinosReductionControl
528 *const reduction_control_status =
529 dynamic_cast<const internal::TrilinosReductionControl *>(
530 status_test.get()))
531 {
532 Assert(dynamic_cast<const ReductionControl *>(&solver_control),
534
535 // Check to see if solver converged in one step
536 // This can happen if the matrix is diagonal and a non-trivial
537 // preconditioner is used.
538 if (solver.NumIters() > 0)
539 {
540 // For ReductionControl, we must first register the initial residual
541 // value. This is the basis from which it will determine whether the
542 // current residual corresponds to a converged state.
543 solver_control.check(
544 0, reduction_control_status->get_initial_residual());
545 solver_control.check(
546 solver.NumIters(),
547 reduction_control_status->get_current_residual());
548 }
549 else
550 solver_control.check(
551 solver.NumIters(),
552 reduction_control_status->get_current_residual());
553 }
554 else
555 {
556 Assert(solver.TrueResidual() >= 0.0, ExcInternalError());
557 solver_control.check(solver.NumIters(), solver.TrueResidual());
558 }
559
560 if (solver_control.last_check() != SolverControl::success)
561 AssertThrow(false,
563 solver_control.last_value()));
564 }
565
566
567
568 template <>
569 void
571 const PreconditionBase &preconditioner)
572 {
573 // Introduce the preconditioner, if the identity preconditioner is used,
574 // the precondioner is set to none, ...
575 if (preconditioner.preconditioner.strong_count() != 0)
576 {
577 const int ierr = solver.SetPrecOperator(
578 const_cast<Epetra_Operator *>(preconditioner.preconditioner.get()));
579 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
580 }
581 else
582 solver.SetAztecOption(AZ_precond, AZ_none);
583 }
584
585
586 template <>
587 void
589 const Epetra_Operator &preconditioner)
590 {
591 const int ierr =
592 solver.SetPrecOperator(const_cast<Epetra_Operator *>(&preconditioner));
593 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
594 }
595
596
597 /* ---------------------- SolverCG ------------------------ */
598
600 : SolverBase(cg, cn, data)
601 {}
602
603
604 /* ---------------------- SolverGMRES ------------------------ */
605
607 : SolverBase(gmres, cn, data)
608 {}
609
610
611 /* ---------------------- SolverBicgstab ------------------------ */
612
616
617
618 /* ---------------------- SolverCGS ------------------------ */
619
621 : SolverBase(cgs, cn, data)
622 {}
623
624
625 /* ---------------------- SolverTFQMR ------------------------ */
626
628 : SolverBase(tfqmr, cn, data)
629 {}
630
631
632
633 /* ---------------------- SolverDirect ------------------------ */
634
640
641
642
645 , additional_data(data.output_solver_details, data.solver_type)
646 {}
647
648
649
651 : solver_control(cn)
652 , additional_data(data.output_solver_details, data.solver_type)
653 {}
654
655
656
659 {
660 return solver_control;
661 }
662
663
664
665 void
667 {
668 // We need an Epetra_LinearProblem object to let the Amesos solver know
669 // about the matrix and vectors.
670 linear_problem = std::make_unique<Epetra_LinearProblem>();
671
672 // Assign the matrix operator to the Epetra_LinearProblem object
673 linear_problem->SetOperator(
674 const_cast<Epetra_CrsMatrix *>(&A.trilinos_matrix()));
675
676 // Fetch return value of Amesos Solver functions
677 int ierr;
678
679 // First set whether we want to print the solver information to screen or
680 // not.
681 ConditionalOStream verbose_cout(std::cout,
682 additional_data.output_solver_details);
683
684 // Next allocate the Amesos solver, this is done in two steps, first we
685 // create a solver Factory and generate with that the concrete Amesos
686 // solver, if possible.
687 Amesos Factory;
688
689 AssertThrow(Factory.Query(additional_data.solver_type.c_str()),
691 "You tried to select the solver type <" +
692 additional_data.solver_type +
693 "> but this solver is not supported by Trilinos either "
694 "because it does not exist, or because Trilinos was not "
695 "configured for its use."));
696
697 solver.reset(
698 Factory.Create(additional_data.solver_type.c_str(), *linear_problem));
699
700 verbose_cout << "Starting symbolic factorization" << std::endl;
701 ierr = solver->SymbolicFactorization();
702 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
703
704 verbose_cout << "Starting numeric factorization" << std::endl;
705 ierr = solver->NumericFactorization();
706 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
707 }
708
709
710
711 void
713 {
714 this->additional_data = data;
715
716 this->initialize(A);
717 }
718
719
720 void
722 {
723 this->vmult(x, b);
724 }
725
726
727
728 void
731 const ::LinearAlgebra::distributed::Vector<double> &b)
732 {
733 this->vmult(x, b);
734 }
735
736
737 void
739 {
740 // Assign the empty LHS vector to the Epetra_LinearProblem object
741 linear_problem->SetLHS(&x.trilinos_vector());
742
743 // Assign the RHS vector to the Epetra_LinearProblem object
744 linear_problem->SetRHS(
745 const_cast<Epetra_MultiVector *>(&b.trilinos_vector()));
746
747 // First set whether we want to print the solver information to screen or
748 // not.
749 ConditionalOStream verbose_cout(std::cout,
750 additional_data.output_solver_details);
751
752
753 verbose_cout << "Starting solve" << std::endl;
754
755 // Fetch return value of Amesos Solver functions
756 int ierr = solver->Solve();
757 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
758
759 // Finally, force the SolverControl object to report convergence
760 solver_control.check(0, 0);
761 }
762
763
764
765 void
768 const ::LinearAlgebra::distributed::Vector<double> &b) const
769 {
770 Epetra_Vector ep_x(View,
771 linear_problem->GetOperator()->OperatorDomainMap(),
772 x.begin());
773 Epetra_Vector ep_b(View,
774 linear_problem->GetOperator()->OperatorRangeMap(),
775 const_cast<double *>(b.begin()));
776
777 // Assign the empty LHS vector to the Epetra_LinearProblem object
778 linear_problem->SetLHS(&ep_x);
779
780 // Assign the RHS vector to the Epetra_LinearProblem object
781 linear_problem->SetRHS(&ep_b);
782
783 // First set whether we want to print the solver information to screen or
784 // not.
785 ConditionalOStream verbose_cout(std::cout,
786 additional_data.output_solver_details);
787
788 verbose_cout << "Starting solve" << std::endl;
789
790 // Fetch return value of Amesos Solver functions
791 int ierr = solver->Solve();
792 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
793
794 // Finally, force the SolverControl object to report convergence
795 solver_control.check(0, 0);
796 }
797
798
799
800 void
802 {
803 // Fetch return value of Amesos Solver functions
804 int ierr;
805
806 // First set whether we want to print the solver information to screen or
807 // not.
808 ConditionalOStream verbose_cout(std::cout,
809 additional_data.output_solver_details);
810
811 // Next allocate the Amesos solver, this is done in two steps, first we
812 // create a solver Factory and generate with that the concrete Amesos
813 // solver, if possible.
814 Amesos Factory;
815
816 AssertThrow(Factory.Query(additional_data.solver_type.c_str()),
818 "You tried to select the solver type <" +
819 additional_data.solver_type +
820 "> but this solver is not supported by Trilinos either "
821 "because it does not exist, or because Trilinos was not "
822 "configured for its use."));
823
824 solver.reset(
825 Factory.Create(additional_data.solver_type.c_str(), *linear_problem));
826
827 verbose_cout << "Starting symbolic factorization" << std::endl;
828 ierr = solver->SymbolicFactorization();
829 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
830
831 verbose_cout << "Starting numeric factorization" << std::endl;
832 ierr = solver->NumericFactorization();
833 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
834
835 verbose_cout << "Starting solve" << std::endl;
836 ierr = solver->Solve();
837 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
838
839 // Finally, let the deal.II SolverControl object know what has
840 // happened. If the solve succeeded, the status of the solver control will
841 // turn into SolverControl::success.
842 solver_control.check(0, 0);
843
844 if (solver_control.last_check() != SolverControl::success)
845 AssertThrow(false,
847 solver_control.last_value()));
848 }
849
850
851 void
853 Epetra_MultiVector &x,
854 const Epetra_MultiVector &b)
855 {
856 // We need an Epetra_LinearProblem object to let the Amesos solver know
857 // about the matrix and vectors.
859 std::make_unique<Epetra_LinearProblem>(const_cast<Epetra_Operator *>(&A),
860 &x,
861 const_cast<Epetra_MultiVector *>(
862 &b));
863
864 do_solve();
865 }
866
867
868 void
869 SolverDirect::solve(const SparseMatrix &sparse_matrix,
870 FullMatrix<double> &solution,
871 const FullMatrix<double> &rhs)
872 {
873 Assert(sparse_matrix.m() == sparse_matrix.n(), ExcInternalError());
874 Assert(rhs.m() == sparse_matrix.m(), ExcInternalError());
875 Assert(rhs.m() == solution.m(), ExcInternalError());
876 Assert(rhs.n() == solution.n(), ExcInternalError());
877
878 const unsigned int m = rhs.m();
879 const unsigned int n = rhs.n();
880
881 FullMatrix<double> rhs_t(n, m);
882 FullMatrix<double> solution_t(n, m);
883
884 rhs_t.copy_transposed(rhs);
885 solution_t.copy_transposed(solution);
886
887 std::vector<double *> rhs_ptrs(n);
888 std::vector<double *> sultion_ptrs(n);
889
890 for (unsigned int i = 0; i < n; ++i)
891 {
892 rhs_ptrs[i] = &rhs_t[i][0];
893 sultion_ptrs[i] = &solution_t[i][0];
894 }
895
896 const Epetra_CrsMatrix &mat = sparse_matrix.trilinos_matrix();
897
898 Epetra_MultiVector trilinos_dst(View,
899 mat.OperatorRangeMap(),
900 sultion_ptrs.data(),
901 sultion_ptrs.size());
902 Epetra_MultiVector trilinos_src(View,
903 mat.OperatorDomainMap(),
904 rhs_ptrs.data(),
905 rhs_ptrs.size());
906
907 this->initialize(sparse_matrix);
908 this->solve(mat, trilinos_dst, trilinos_src);
909
910 solution.copy_transposed(solution_t);
911 }
912
913
914 void
916 MPI::Vector &x,
917 const MPI::Vector &b)
918 {
919 // We need an Epetra_LinearProblem object to let the Amesos solver know
920 // about the matrix and vectors.
921 linear_problem = std::make_unique<Epetra_LinearProblem>(
922 const_cast<Epetra_CrsMatrix *>(&A.trilinos_matrix()),
923 &x.trilinos_vector(),
924 const_cast<Epetra_MultiVector *>(&b.trilinos_vector()));
925
926 do_solve();
927 }
928
929
930
931 void
934 const ::Vector<double> &b)
935 {
936 // In case we call the solver with deal.II vectors, we create views of the
937 // vectors in Epetra format.
938 Assert(x.size() == A.n(), ExcDimensionMismatch(x.size(), A.n()));
939 Assert(b.size() == A.m(), ExcDimensionMismatch(b.size(), A.m()));
940 Assert(A.local_range().second == A.m(),
941 ExcMessage("Can only work in serial when using deal.II vectors."));
942 Epetra_Vector ep_x(View, A.trilinos_matrix().DomainMap(), x.begin());
943 Epetra_Vector ep_b(View,
944 A.trilinos_matrix().RangeMap(),
945 const_cast<double *>(b.begin()));
946
947 // We need an Epetra_LinearProblem object to let the Amesos solver know
948 // about the matrix and vectors.
949 linear_problem = std::make_unique<Epetra_LinearProblem>(
950 const_cast<Epetra_CrsMatrix *>(&A.trilinos_matrix()), &ep_x, &ep_b);
951
952 do_solve();
953 }
954
955
956
957 void
959 const SparseMatrix &A,
961 const ::LinearAlgebra::distributed::Vector<double> &b)
962 {
964 A.trilinos_matrix().DomainMap().NumMyElements());
965 AssertDimension(b.locally_owned_size(),
966 A.trilinos_matrix().RangeMap().NumMyElements());
967 Epetra_Vector ep_x(View, A.trilinos_matrix().DomainMap(), x.begin());
968 Epetra_Vector ep_b(View,
969 A.trilinos_matrix().RangeMap(),
970 const_cast<double *>(b.begin()));
971
972 // We need an Epetra_LinearProblem object to let the Amesos solver know
973 // about the matrix and vectors.
974 linear_problem = std::make_unique<Epetra_LinearProblem>(
975 const_cast<Epetra_CrsMatrix *>(&A.trilinos_matrix()), &ep_x, &ep_b);
976
977 do_solve();
978 }
979} // namespace TrilinosWrappers
980
981
982// explicit instantiations
983// TODO: put these instantiations into generic file
984namespace TrilinosWrappers
985{
986 template void
987 SolverBase::do_solve(const PreconditionBase &preconditioner);
988
989 template void
990 SolverBase::do_solve(const Epetra_Operator &preconditioner);
991} // namespace TrilinosWrappers
992
994
995#endif // DEAL_II_WITH_PETSC
size_type n() const
void copy_transposed(const MatrixType &)
size_type m() const
size_type locally_owned_size() const
SolverCG(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
@ success
Stop iteration, goal reached.
const Epetra_MultiVector & trilinos_vector() const
Teuchos::RCP< Epetra_Operator > preconditioner
SolverControl & control() const
void solve(const SparseMatrix &A, MPI::Vector &x, const MPI::Vector &b, const PreconditionBase &preconditioner)
const AdditionalData additional_data
std::unique_ptr< AztecOO_StatusTest > status_test
void set_preconditioner(AztecOO &solver, const Preconditioner &preconditioner)
void do_solve(const Preconditioner &preconditioner)
std::unique_ptr< Epetra_LinearProblem > linear_problem
SolverBase(SolverControl &cn, const AdditionalData &data=AdditionalData())
enum TrilinosWrappers::SolverBase::SolverName solver_name
SolverBicgstab(SolverControl &cn, const AdditionalData &data=AdditionalData())
SolverCGS(SolverControl &cn, const AdditionalData &data=AdditionalData())
void initialize(const SparseMatrix &A)
void vmult(MPI::Vector &x, const MPI::Vector &b) const
std::unique_ptr< Amesos_BaseSolver > solver
void solve(MPI::Vector &x, const MPI::Vector &b)
std::unique_ptr< Epetra_LinearProblem > linear_problem
SolverControl & control() const
SolverDirect(const AdditionalData &data=AdditionalData())
SolverGMRES(SolverControl &cn, const AdditionalData &data=AdditionalData())
SolverTFQMR(SolverControl &cn, const AdditionalData &data=AdditionalData())
const Epetra_CrsMatrix & trilinos_matrix() const
virtual size_type size() const override
iterator begin()
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition config.h:603
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition config.h:647
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcTrilinosError(int arg1)
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
AdditionalData(const bool output_solver_details=false, const unsigned int gmres_restart_parameter=30)
AdditionalData(const bool output_solver_details=false, const std::string &solver_type="Amesos_Klu")
double TrilinosScalar
Definition types.h:198