Reference documentation for deal.II version 9.0.0
trilinos_solver.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/lac/trilinos_solver.h>
17 
18 #ifdef DEAL_II_WITH_TRILINOS
19 
20 # include <deal.II/base/conditional_ostream.h>
21 # include <deal.II/base/std_cxx14/memory.h>
22 # include <deal.II/lac/trilinos_sparse_matrix.h>
23 # include <deal.II/lac/trilinos_vector.h>
24 # include <deal.II/lac/trilinos_precondition.h>
25 
26 # include <AztecOO_StatusTest.h>
27 # include <AztecOO_StatusTestMaxIters.h>
28 # include <AztecOO_StatusTestResNorm.h>
29 # include <AztecOO_StatusTestCombo.h>
30 # include <AztecOO_StatusType.h>
31 
32 # include <cmath>
33 # include <limits>
34 
35 DEAL_II_NAMESPACE_OPEN
36 
37 namespace TrilinosWrappers
38 {
39 
40  SolverBase::AdditionalData::AdditionalData (const bool output_solver_details,
41  const unsigned int gmres_restart_parameter)
42  :
43  output_solver_details (output_solver_details),
44  gmres_restart_parameter (gmres_restart_parameter)
45  {}
46 
47 
48 
50  const AdditionalData &data)
51  :
52  solver_name (gmres),
53  solver_control (cn),
54  additional_data (data)
55  {}
56 
57 
58 
60  SolverControl &cn,
61  const AdditionalData &data)
62  :
63  solver_name (solver_name),
64  solver_control (cn),
65  additional_data (data)
66  {}
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_cxx14::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
99  SolverBase::solve (const Epetra_Operator &A,
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.
106  linear_problem = std_cxx14::make_unique<Epetra_LinearProblem>
107  (const_cast<Epetra_Operator *>(&A),
108  &x.trilinos_vector(),
109  const_cast<Epetra_MultiVector *>(&b.trilinos_vector()));
110 
111  do_solve(preconditioner);
112  }
113 
114 
115 
116  // Note: "A" is set as a constant reference so that all patterns for ::solve
117  // can be used by the inverse_operator of LinearOperator
118  void
119  SolverBase::solve (const Epetra_Operator &A,
120  MPI::Vector &x,
121  const MPI::Vector &b,
122  const Epetra_Operator &preconditioner)
123  {
124  // We need an Epetra_LinearProblem object to let the AztecOO solver know
125  // about the matrix and vectors.
126  linear_problem = std_cxx14::make_unique<Epetra_LinearProblem>
127  (const_cast<Epetra_Operator *>(&A),
128  &x.trilinos_vector(),
129  const_cast<Epetra_MultiVector *>(&b.trilinos_vector()));
130 
131  do_solve(preconditioner);
132  }
133 
134 
135 
136  // Note: "A" is set as a constant reference so that all patterns for ::solve
137  // can be used by the inverse_operator of LinearOperator
138  void
139  SolverBase::solve (const Epetra_Operator &A,
140  Epetra_MultiVector &x,
141  const Epetra_MultiVector &b,
142  const PreconditionBase &preconditioner)
143  {
144  // We need an Epetra_LinearProblem object to let the AztecOO solver know
145  // about the matrix and vectors.
146  linear_problem = std_cxx14::make_unique<Epetra_LinearProblem>
147  (const_cast<Epetra_Operator *>(&A),
148  &x,
149  const_cast<Epetra_MultiVector *>(&b));
150 
151  do_solve(preconditioner);
152  }
153 
154 
155 
156  // Note: "A" is set as a constant reference so that all patterns for ::solve
157  // can be used by the inverse_operator of LinearOperator
158  void
159  SolverBase::solve (const Epetra_Operator &A,
160  Epetra_MultiVector &x,
161  const Epetra_MultiVector &b,
162  const Epetra_Operator &preconditioner)
163  {
164  // We need an Epetra_LinearProblem object to let the AztecOO solver know
165  // about the matrix and vectors.
166  linear_problem = std_cxx14::make_unique<Epetra_LinearProblem>
167  (const_cast<Epetra_Operator *>(&A),
168  &x,
169  const_cast<Epetra_MultiVector *>(&b));
170 
171  do_solve(preconditioner);
172  }
173 
174 
175 
176  void
178  ::Vector<double> &x,
179  const ::Vector<double> &b,
180  const PreconditionBase &preconditioner)
181  {
182  // In case we call the solver with deal.II vectors, we create views of the
183  // vectors in Epetra format.
184  Assert (x.size() == A.n(),
185  ExcDimensionMismatch(x.size(), A.n()));
186  Assert (b.size() == A.m(),
187  ExcDimensionMismatch(b.size(), A.m()));
188  Assert (A.local_range ().second == A.m(),
189  ExcMessage ("Can only work in serial when using deal.II vectors."));
190  Assert (A.trilinos_matrix().Filled(),
191  ExcMessage ("Matrix is not compressed. Call compress() method."));
192 
193  Epetra_Vector ep_x (View, A.domain_partitioner(), x.begin());
194  Epetra_Vector ep_b (View, A.range_partitioner(), const_cast<double *>(b.begin()));
195 
196  // We need an Epetra_LinearProblem object to let the AztecOO solver know
197  // about the matrix and vectors.
198  linear_problem = std_cxx14::make_unique<Epetra_LinearProblem>
199  (const_cast<Epetra_CrsMatrix *>(&A.trilinos_matrix()),
200  &ep_x, &ep_b);
201 
202  do_solve(preconditioner);
203  }
204 
205 
206 
207  void
208  SolverBase::solve (Epetra_Operator &A,
209  ::Vector<double> &x,
210  const ::Vector<double> &b,
211  const PreconditionBase &preconditioner)
212  {
213  Epetra_Vector ep_x (View, A.OperatorDomainMap(), x.begin());
214  Epetra_Vector ep_b (View, A.OperatorRangeMap(), const_cast<double *>(b.begin()));
215 
216  // We need an Epetra_LinearProblem object to let the AztecOO solver know
217  // about the matrix and vectors.
218  linear_problem = std_cxx14::make_unique<Epetra_LinearProblem>(&A,&ep_x, &ep_b);
219 
220  do_solve(preconditioner);
221  }
222 
223 
224 
225  void
228  const ::LinearAlgebra::distributed::Vector<double> &b,
229  const PreconditionBase &preconditioner)
230  {
231  // In case we call the solver with deal.II vectors, we create views of the
232  // vectors in Epetra format.
233  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(x.local_size()),
234  A.domain_partitioner().NumMyElements());
235  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(b.local_size()),
236  A.range_partitioner().NumMyElements());
237 
238  Epetra_Vector ep_x (View, A.domain_partitioner(), x.begin());
239  Epetra_Vector ep_b (View, A.range_partitioner(), const_cast<double *>(b.begin()));
240 
241  // We need an Epetra_LinearProblem object to let the AztecOO solver know
242  // about the matrix and vectors.
243  linear_problem = std_cxx14::make_unique<Epetra_LinearProblem>
244  (const_cast<Epetra_CrsMatrix *>(&A.trilinos_matrix()),
245  &ep_x, &ep_b);
246 
247  do_solve(preconditioner);
248  }
249 
250 
251 
252  void
253  SolverBase::solve (Epetra_Operator &A,
255  const ::LinearAlgebra::distributed::Vector<double> &b,
256  const PreconditionBase &preconditioner)
257  {
258  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(x.local_size()),
259  A.OperatorDomainMap().NumMyElements());
260  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(b.local_size()),
261  A.OperatorRangeMap().NumMyElements());
262 
263  Epetra_Vector ep_x (View, A.OperatorDomainMap(), x.begin());
264  Epetra_Vector ep_b (View, A.OperatorRangeMap(), const_cast<double *>(b.begin()));
265 
266  // We need an Epetra_LinearProblem object to let the AztecOO solver know
267  // about the matrix and vectors.
268  linear_problem = std_cxx14::make_unique<Epetra_LinearProblem>(&A,&ep_x, &ep_b);
269 
270  do_solve(preconditioner);
271  }
272 
273 
274  namespace internal
275  {
276  namespace
277  {
278  double
279  compute_residual (const Epetra_MultiVector *const residual_vector)
280  {
281  Assert(residual_vector->NumVectors() == 1,
282  ExcMessage("Residual multivector holds more than one vector"));
283  TrilinosScalar res_l2_norm = 0.0;
284  const int ierr = residual_vector->Norm2 (&res_l2_norm);
285  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
286  return res_l2_norm;
287  }
288 
289  class TrilinosReductionControl : public AztecOO_StatusTest
290  {
291  public:
292  TrilinosReductionControl (const int &max_steps,
293  const double &tolerance,
294  const double &reduction,
295  const Epetra_LinearProblem &linear_problem);
296 
297  virtual ~TrilinosReductionControl() = default;
298 
299  virtual bool
300  ResidualVectorRequired () const
301  {
302  return status_test_collection->ResidualVectorRequired();
303  }
304 
305  virtual AztecOO_StatusType
306  CheckStatus (int CurrentIter,
307  Epetra_MultiVector *CurrentResVector,
308  double CurrentResNormEst,
309  bool SolutionUpdated)
310  {
311  // Note: CurrentResNormEst is set to -1.0 if no estimate of the
312  // residual value is available
313  current_residual = (CurrentResNormEst < 0.0 ?
314  compute_residual(CurrentResVector) :
315  CurrentResNormEst);
316  if (CurrentIter == 0)
317  initial_residual = current_residual;
318 
319  return status_test_collection->CheckStatus(CurrentIter,
320  CurrentResVector,
321  CurrentResNormEst,
322  SolutionUpdated);
323 
324  }
325 
326  virtual AztecOO_StatusType
327  GetStatus () const
328  {
329  return status_test_collection->GetStatus();
330  }
331 
332  virtual std::ostream &
333  Print (std::ostream &stream,
334  int indent = 0) const
335  {
336  return status_test_collection->Print(stream,indent);
337  }
338 
339  double
340  get_initial_residual() const
341  {
342  return initial_residual;
343  }
344 
345  double
346  get_current_residual() const
347  {
348  return current_residual;
349  }
350 
351  private:
352  double initial_residual;
353  double current_residual;
354  std::unique_ptr<AztecOO_StatusTestCombo> status_test_collection;
355  std::unique_ptr<AztecOO_StatusTestMaxIters> status_test_max_steps;
356  std::unique_ptr<AztecOO_StatusTestResNorm> status_test_abs_tol;
357  std::unique_ptr<AztecOO_StatusTestResNorm> status_test_rel_tol;
358  };
359 
360 
361  TrilinosReductionControl::TrilinosReductionControl(
362  const int &max_steps,
363  const double &tolerance,
364  const double &reduction,
365  const Epetra_LinearProblem &linear_problem )
366  :
367  initial_residual (std::numeric_limits<double>::max()),
368  current_residual(std::numeric_limits<double>::max())
369  {
370  // Consider linear problem converged if any of the collection
371  // of criterion are met
372  status_test_collection = std_cxx14::make_unique<AztecOO_StatusTestCombo>
373  (AztecOO_StatusTestCombo::OR);
374 
375  // Maximum number of iterations
376  Assert (max_steps >=0, ExcInternalError());
377  status_test_max_steps = std_cxx14::make_unique<AztecOO_StatusTestMaxIters>(max_steps);
378  status_test_collection->AddStatusTest(*status_test_max_steps);
379 
380  Assert(linear_problem.GetRHS()->NumVectors() == 1,
381  ExcMessage("RHS multivector holds more than one vector"));
382 
383  // Residual norm is below some absolute value
384  status_test_abs_tol = std_cxx14::make_unique<AztecOO_StatusTestResNorm>
385  (*linear_problem.GetOperator(),
386  *(linear_problem.GetLHS()->operator()(0)),
387  *(linear_problem.GetRHS()->operator()(0)),
388  tolerance);
389  status_test_abs_tol->DefineResForm(AztecOO_StatusTestResNorm::Explicit,
390  AztecOO_StatusTestResNorm::TwoNorm);
391  status_test_abs_tol->DefineScaleForm(AztecOO_StatusTestResNorm::None,
392  AztecOO_StatusTestResNorm::TwoNorm);
393  status_test_collection->AddStatusTest(*status_test_abs_tol);
394 
395  // Residual norm, scaled by some initial value, is below some threshold
396  status_test_rel_tol = std_cxx14::make_unique<AztecOO_StatusTestResNorm>
397  (*linear_problem.GetOperator(),
398  *(linear_problem.GetLHS()->operator()(0)),
399  *(linear_problem.GetRHS()->operator()(0)),
400  reduction);
401  status_test_rel_tol->DefineResForm(AztecOO_StatusTestResNorm::Explicit,
402  AztecOO_StatusTestResNorm::TwoNorm);
403  status_test_rel_tol->DefineScaleForm(AztecOO_StatusTestResNorm::NormOfInitRes,
404  AztecOO_StatusTestResNorm::TwoNorm);
405  status_test_collection->AddStatusTest(*status_test_rel_tol);
406  }
407 
408  }
409  }
410 
411 
412  template <typename Preconditioner>
413  void
414  SolverBase::do_solve(const Preconditioner &preconditioner)
415  {
416  int ierr;
417 
418  // Next we can allocate the AztecOO solver...
419  solver.SetProblem(*linear_problem);
420 
421  // ... and we can specify the solver to be used.
422  switch (solver_name)
423  {
424  case cg:
425  solver.SetAztecOption(AZ_solver, AZ_cg);
426  break;
427  case cgs:
428  solver.SetAztecOption(AZ_solver, AZ_cgs);
429  break;
430  case gmres:
431  solver.SetAztecOption(AZ_solver, AZ_gmres);
432  solver.SetAztecOption(AZ_kspace, additional_data.gmres_restart_parameter);
433  break;
434  case bicgstab:
435  solver.SetAztecOption(AZ_solver, AZ_bicgstab);
436  break;
437  case tfqmr:
438  solver.SetAztecOption(AZ_solver, AZ_tfqmr);
439  break;
440  default:
441  Assert (false, ExcNotImplemented());
442  }
443 
444  // Set the preconditioner
445  set_preconditioner(solver, preconditioner);
446 
447  // ... set some options, ...
448  solver.SetAztecOption (AZ_output, additional_data.output_solver_details ?
449  AZ_all : AZ_none);
450  solver.SetAztecOption (AZ_conv, AZ_noscaled);
451 
452  // By default, the Trilinos solver chooses convergence criterion based on
453  // the number of iterations made and an absolute tolerance.
454  // This implies that the use of the standard Trilinos convergence test
455  // actually coincides with ::IterationNumberControl because the
456  // solver, unless explicitly told otherwise, will Iterate() until a number
457  // of max_steps() are taken or an absolute tolerance() is attained.
458  // It is therefore suitable for use with both SolverControl or
459  // IterationNumberControl. The final check at the end will determine whether
460  // failure to converge to the defined residual norm constitutes failure
461  // (SolverControl) or is alright (IterationNumberControl).
462  // In the case that the SolverControl wants to perform ReductionControl,
463  // then we have to do a little extra something by prescribing a custom
464  // status test.
465  if (!status_test)
466  {
467  if (const ReductionControl* const reduction_control
468  = dynamic_cast<const ReductionControl *const>(&solver_control))
469  {
470  status_test = std_cxx14::make_unique<internal::TrilinosReductionControl>
471  (reduction_control->max_steps(),
472  reduction_control->tolerance(),
473  reduction_control->reduction(),
474  *linear_problem);
475  solver.SetStatusTest(status_test.get());
476  }
477  }
478 
479  // ... and then solve!
480  ierr = solver.Iterate (solver_control.max_steps(),
481  solver_control.tolerance());
482 
483  // report errors in more detail than just by checking whether the return
484  // status is zero or greater. the error strings are taken from the
485  // implementation of the AztecOO::Iterate function
486  switch (ierr)
487  {
488  case -1:
489  AssertThrow (false, ExcMessage("AztecOO::Iterate error code -1: "
490  "option not implemented"));
491  break;
492  case -2:
493  AssertThrow (false, ExcMessage("AztecOO::Iterate error code -2: "
494  "numerical breakdown"));
495  break;
496  case -3:
497  AssertThrow (false, ExcMessage("AztecOO::Iterate error code -3: "
498  "loss of precision"));
499  break;
500  case -4:
501  AssertThrow (false, ExcMessage("AztecOO::Iterate error code -4: "
502  "GMRES Hessenberg ill-conditioned"));
503  break;
504  default:
505  AssertThrow (ierr >= 0, ExcTrilinosError(ierr));
506  }
507 
508  // Finally, let the deal.II SolverControl object know what has
509  // happened. If the solve succeeded, the status of the solver control will
510  // turn into SolverControl::success.
511  // If the residual is not computed/stored by the solver, as can happen for
512  // certain choices of solver or if a custom status test is set, then the
513  // result returned by TrueResidual() is equal to -1. In this case we must
514  // compute it ourself.
515  if (const internal::TrilinosReductionControl* const reduction_control_status
516  = dynamic_cast<const internal::TrilinosReductionControl *const>(status_test.get()))
517  {
518  Assert(dynamic_cast<const ReductionControl *const>(&solver_control), ExcInternalError());
519 
520  // Check to see if solver converged in one step
521  // This can happen if the matrix is diagonal and a non-trivial
522  // preconditioner is used.
523  if (solver.NumIters() > 0)
524  {
525  // For ReductionControl, we must first register the initial residual
526  // value. This is the basis from which it will determine whether the
527  // current residual corresponds to a converged state.
528  solver_control.check (0, reduction_control_status->get_initial_residual());
529  solver_control.check (solver.NumIters(), reduction_control_status->get_current_residual());
530  }
531  else
532  solver_control.check (solver.NumIters(), reduction_control_status->get_current_residual());
533  }
534  else
535  {
536  Assert(solver.TrueResidual() >= 0.0, ExcInternalError());
537  solver_control.check (solver.NumIters(), solver.TrueResidual());
538  }
539 
540  if (solver_control.last_check() != SolverControl::success)
541  AssertThrow(false, SolverControl::NoConvergence (solver_control.last_step(),
542  solver_control.last_value()));
543  }
544 
545 
546 
547  template <>
548  void
549  SolverBase::set_preconditioner(AztecOO &solver,
550  const PreconditionBase &preconditioner)
551  {
552  // Introduce the preconditioner, if the identity preconditioner is used,
553  // the precondioner is set to none, ...
554  if (preconditioner.preconditioner.use_count()!=0)
555  {
556  const int ierr = solver.SetPrecOperator (const_cast<Epetra_Operator *>
557  (preconditioner.preconditioner.get()));
558  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
559  }
560  else
561  solver.SetAztecOption(AZ_precond,AZ_none);
562  }
563 
564 
565  template <>
566  void
567  SolverBase::set_preconditioner(AztecOO &solver,
568  const Epetra_Operator &preconditioner)
569  {
570  const int ierr = solver.SetPrecOperator (const_cast<Epetra_Operator *>(&preconditioner));
571  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
572  }
573 
574 
575  /* ---------------------- SolverCG ------------------------ */
576 
577  SolverCG::AdditionalData::
578  AdditionalData (const bool output_solver_details)
579  :
580  SolverBase::AdditionalData (output_solver_details)
581  {}
582 
583 
584 
586  const AdditionalData &data)
587  :
588  SolverBase (cn, data),
589  additional_data (data)
590  {
591  solver_name = cg;
592  }
593 
594 
595  /* ---------------------- SolverGMRES ------------------------ */
596 
598  AdditionalData (const bool output_solver_details,
599  const unsigned int restart_parameter)
600  :
601  SolverBase::AdditionalData (output_solver_details,
602  restart_parameter)
603  {}
604 
605 
606 
608  const AdditionalData &data)
609  :
610  SolverBase (cn, data),
611  additional_data (data)
612  {
613  solver_name = gmres;
614  }
615 
616 
617  /* ---------------------- SolverBicgstab ------------------------ */
618 
620  AdditionalData (const bool output_solver_details)
621  :
622  SolverBase::AdditionalData (output_solver_details)
623  {}
624 
625 
626 
627 
629  const AdditionalData &data)
630  :
631  SolverBase (cn, data),
632  additional_data (data)
633  {
634  solver_name = bicgstab;
635  }
636 
637 
638  /* ---------------------- SolverCGS ------------------------ */
639 
641  AdditionalData (const bool output_solver_details)
642  :
643  SolverBase::AdditionalData (output_solver_details)
644  {}
645 
646 
647 
649  const AdditionalData &data)
650  :
651  SolverBase (cn, data),
652  additional_data (data)
653  {
654  solver_name = cgs;
655  }
656 
657 
658  /* ---------------------- SolverTFQMR ------------------------ */
659 
661  AdditionalData (const bool output_solver_details)
662  :
663  SolverBase::AdditionalData (output_solver_details)
664  {}
665 
666 
667 
669  const AdditionalData &data)
670  :
671  SolverBase (cn, data),
672  additional_data (data)
673  {
674  solver_name = tfqmr;
675  }
676 
677 
678 
679  /* ---------------------- SolverDirect ------------------------ */
680 
682  AdditionalData (const bool output_solver_details,
683  const std::string &solver_type)
684  :
685  output_solver_details (output_solver_details),
686  solver_type(solver_type)
687  {}
688 
689 
690 
691 
693  const AdditionalData &data)
694  :
695  solver_control (cn),
696  additional_data (data.output_solver_details,data.solver_type)
697  {}
698 
699 
700 
701  SolverControl &
703  {
704  return solver_control;
705  }
706 
707 
708 
710  {
711  // We need an Epetra_LinearProblem object to let the Amesos solver know
712  // about the matrix and vectors.
713  linear_problem = std_cxx14::make_unique<Epetra_LinearProblem >();
714 
715  // Assign the matrix operator to the Epetra_LinearProblem object
716  linear_problem->SetOperator(const_cast<Epetra_CrsMatrix *>(&A.trilinos_matrix()));
717 
718  // Fetch return value of Amesos Solver functions
719  int ierr;
720 
721  // First set whether we want to print the solver information to screen or
722  // not.
723  ConditionalOStream verbose_cout (std::cout,
725 
726  // Next allocate the Amesos solver, this is done in two steps, first we
727  // create a solver Factory and and generate with that the concrete Amesos
728  // solver, if possible.
729  Amesos Factory;
730 
731  AssertThrow(
732  Factory.Query(additional_data.solver_type.c_str()),
733  ExcMessage (std::string ("You tried to select the solver type <") +
735  "> but this solver is not supported by Trilinos either "
736  "because it does not exist, or because Trilinos was not "
737  "configured for its use.")
738  );
739 
740  solver.reset (
741  Factory.Create(additional_data.solver_type.c_str(), *linear_problem)
742  );
743 
744  verbose_cout << "Starting symbolic factorization" << std::endl;
745  ierr = solver->SymbolicFactorization();
746  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
747 
748  verbose_cout << "Starting numeric factorization" << std::endl;
749  ierr = solver->NumericFactorization();
750  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
751  }
752 
753 
755  {
756  // Assign the empty LHS vector to the Epetra_LinearProblem object
757  linear_problem->SetLHS(&x.trilinos_vector());
758 
759  // Assign the RHS vector to the Epetra_LinearProblem object
760  linear_problem->SetRHS(const_cast<Epetra_MultiVector *>(&b.trilinos_vector()));
761 
762  // First set whether we want to print the solver information to screen or
763  // not.
764  ConditionalOStream verbose_cout (std::cout,
766 
767 
768  verbose_cout << "Starting solve" << std::endl;
769 
770  // Fetch return value of Amesos Solver functions
771  int ierr = solver->Solve ();
772  AssertThrow (ierr == 0, ExcTrilinosError (ierr));
773 
774  // Finally, force the SolverControl object to report convergence
775  solver_control.check (0, 0);
776  }
777 
778 
779 
781  const ::LinearAlgebra::distributed::Vector<double> &b)
782  {
783  Epetra_Vector ep_x (View, linear_problem->GetOperator()->OperatorDomainMap(), x.begin());
784  Epetra_Vector ep_b (View, linear_problem->GetOperator()->OperatorRangeMap(), const_cast<double *>(b.begin()));
785 
786  // Assign the empty LHS vector to the Epetra_LinearProblem object
787  linear_problem->SetLHS(&ep_x);
788 
789  // Assign the RHS vector to the Epetra_LinearProblem object
790  linear_problem->SetRHS(&ep_b);
791 
792  // First set whether we want to print the solver information to screen or
793  // not.
794  ConditionalOStream verbose_cout (std::cout,
796 
797  verbose_cout << "Starting solve" << std::endl;
798 
799  // Fetch return value of Amesos Solver functions
800  int ierr = solver->Solve ();
801  AssertThrow (ierr == 0, ExcTrilinosError (ierr));
802 
803  // Finally, force the SolverControl object to report convergence
804  solver_control.check (0, 0);
805  }
806 
807 
808 
809  void
811  {
812  // Fetch return value of Amesos Solver functions
813  int ierr;
814 
815  // First set whether we want to print the solver information to screen or
816  // not.
817  ConditionalOStream verbose_cout (std::cout,
819 
820  // Next allocate the Amesos solver, this is done in two steps, first we
821  // create a solver Factory and and generate with that the concrete Amesos
822  // solver, if possible.
823  Amesos Factory;
824 
825  AssertThrow (Factory.Query (additional_data.solver_type.c_str ()),
826  ExcMessage (std::
827  string ("You tried to select the solver type <")
829  "> but this solver is not supported by Trilinos either "
830  "because it does not exist, or because Trilinos was not "
831  "configured for its use."));
832 
833  solver.reset (Factory.
834  Create (additional_data.solver_type.c_str (),
835  *linear_problem));
836 
837  verbose_cout << "Starting symbolic factorization" << std::endl;
838  ierr = solver->SymbolicFactorization ();
839  AssertThrow (ierr == 0, ExcTrilinosError (ierr));
840 
841  verbose_cout << "Starting numeric factorization" << std::endl;
842  ierr = solver->NumericFactorization ();
843  AssertThrow (ierr == 0, ExcTrilinosError (ierr));
844 
845  verbose_cout << "Starting solve" << std::endl;
846  ierr = solver->Solve();
847  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
848 
849  // Finally, let the deal.II SolverControl object know what has
850  // happened. If the solve succeeded, the status of the solver control will
851  // turn into SolverControl::success.
852  solver_control.check (0, 0);
853 
857  }
858 
859 
860  void
862  MPI::Vector &x,
863  const MPI::Vector &b)
864  {
865  // We need an Epetra_LinearProblem object to let the Amesos solver know
866  // about the matrix and vectors.
867  linear_problem = std_cxx14::make_unique<Epetra_LinearProblem>
868  (const_cast<Epetra_CrsMatrix *>(&A.trilinos_matrix()),
869  &x.trilinos_vector(),
870  const_cast<Epetra_MultiVector *>(&b.trilinos_vector()));
871 
872  do_solve();
873  }
874 
875 
876 
877  void
879  ::Vector<double> &x,
880  const ::Vector<double> &b)
881  {
882 
883  // In case we call the solver with deal.II vectors, we create views of the
884  // vectors in Epetra format.
885  Assert (x.size() == A.n(),
886  ExcDimensionMismatch(x.size(), A.n()));
887  Assert (b.size() == A.m(),
888  ExcDimensionMismatch(b.size(), A.m()));
889  Assert (A.local_range ().second == A.m(),
890  ExcMessage ("Can only work in serial when using deal.II vectors."));
891  Epetra_Vector ep_x (View, A.domain_partitioner(), x.begin());
892  Epetra_Vector ep_b (View, A.range_partitioner(), const_cast<double *>(b.begin()));
893 
894  // We need an Epetra_LinearProblem object to let the Amesos solver know
895  // about the matrix and vectors.
896  linear_problem = std_cxx14::make_unique<Epetra_LinearProblem>
897  (const_cast<Epetra_CrsMatrix *>(&A.trilinos_matrix()),
898  &ep_x, &ep_b);
899 
900  do_solve();
901  }
902 
903 
904 
905  void
908  const ::LinearAlgebra::distributed::Vector<double> &b)
909  {
910  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(x.local_size()),
911  A.domain_partitioner().NumMyElements());
912  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(b.local_size()),
913  A.range_partitioner().NumMyElements());
914  Epetra_Vector ep_x (View, A.domain_partitioner(), x.begin());
915  Epetra_Vector ep_b (View, A.range_partitioner(), const_cast<double *>(b.begin()));
916 
917  // We need an Epetra_LinearProblem object to let the Amesos solver know
918  // about the matrix and vectors.
919  linear_problem = std_cxx14::make_unique<Epetra_LinearProblem>
920  (const_cast<Epetra_CrsMatrix *>(&A.trilinos_matrix()),
921  &ep_x, &ep_b);
922 
923  do_solve();
924  }
925 }
926 
927 
928 // explicit instantiations
929 // TODO: put these instantiations into generic file
930 namespace TrilinosWrappers
931 {
932  template void
933  SolverBase::do_solve(const PreconditionBase &preconditioner);
934 
935  template void
936  SolverBase::do_solve(const Epetra_Operator &preconditioner);
937 }
938 
939 DEAL_II_NAMESPACE_CLOSE
940 
941 #endif // DEAL_II_WITH_PETSC
static ::ExceptionBase & ExcTrilinosError(int arg1)
SolverCG(SolverControl &cn, const AdditionalData &data=AdditionalData())
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1248
virtual State check(const unsigned int step, const double check_value)
std::unique_ptr< Amesos_BaseSolver > solver
const Epetra_Map & range_partitioner() const
const AdditionalData additional_data
const AdditionalData additional_data
AdditionalData(const bool output_solver_details=false, const unsigned int restart_parameter=30)
const AdditionalData additional_data
STL namespace.
void initialize(const SparseMatrix &A)
void solve(MPI::Vector &x, const MPI::Vector &b)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
const AdditionalData additional_data
void do_solve(const Preconditioner &preconditioner)
static ::ExceptionBase & ExcMessage(std::string arg1)
SolverCGS(SolverControl &cn, const AdditionalData &data=AdditionalData())
Stop iteration, goal reached.
#define Assert(cond, exc)
Definition: exceptions.h:1142
AdditionalData(const bool output_solver_details=false)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
SolverBicgstab(SolverControl &cn, const AdditionalData &data=AdditionalData())
std::unique_ptr< Epetra_LinearProblem > linear_problem
SolverControl & control() const
void solve(const SparseMatrix &A, MPI::Vector &x, const MPI::Vector &b, const PreconditionBase &preconditioner)
const Epetra_CrsMatrix & trilinos_matrix() const
std::unique_ptr< Epetra_LinearProblem > linear_problem
iterator begin()
double last_value() const
SolverTFQMR(SolverControl &cn, const AdditionalData &data=AdditionalData())
AdditionalData(const bool output_solver_details=false)
SolverBase(SolverControl &cn, const AdditionalData &data=AdditionalData())
unsigned int last_step() const
const AdditionalData additional_data
AdditionalData(const bool output_solver_details=false, const std::string &solver_type="Amesos_Klu")
SolverDirect(SolverControl &cn, const AdditionalData &data=AdditionalData())
State last_check() const
const Epetra_Map & domain_partitioner() const
const AdditionalData additional_data
size_type size() const
std::pair< size_type, size_type > local_range() const
AdditionalData(const bool output_solver_details=false)
SolverControl & control() const
static ::ExceptionBase & ExcNotImplemented()
AdditionalData(const bool output_solver_details=false, const unsigned int gmres_restart_parameter=30)
std::shared_ptr< Epetra_Operator > preconditioner
const Epetra_MultiVector & trilinos_vector() const
SolverGMRES(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
T max(const T &t, const MPI_Comm &mpi_communicator)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcTrilinosError(int arg1)