trilinos_solver.cc
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2024 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
26 #include "trilinos_solver.h"
27 
28 
29 namespace oomph
30 {
31  /// ////////////////////////////////////////////////////////////////////////////
32  /// ////////////////////////////////////////////////////////////////////////////
33  // functions for TrilinosAztecOOSolver class
34  /// ////////////////////////////////////////////////////////////////////////////
35  /// ////////////////////////////////////////////////////////////////////////////
36 
37 
38  //=============================================================================
39  /// Function which uses problem_pt's get_jacobian(...) function to
40  /// generate a linear system which is then solved. This function deletes
41  /// any existing internal data and then generates a new AztecOO solver.
42  //=============================================================================
43  void TrilinosAztecOOSolver::solve(Problem* const& problem_pt,
44  DoubleVector& solution)
45  {
46  // MemoryUsage::doc_memory_usage("start of TrilinosAztecOOSolver::solve");
47  // MemoryUsage::insert_comment_to_continous_top(
48  // "start of TrilinosAztecOOSolver::solve");
49 
50  // clean up from previous solve
52 
53  // if we were previously using a problem based solve then we should delete
54  // the oomphlib matrix as it was created during the last solve
56  {
57  delete Oomph_matrix_pt;
58  Oomph_matrix_pt = NULL;
59  }
60 
61  // this is a problem based solve
63 
64  // store the problem_pt
65  Problem_pt = problem_pt;
66 
67  // Get oomph-lib Jacobian matrix and residual vector
68 
69  // record the start time
70  double start_t = TimingHelpers::timer();
71 
72  // MemoryUsage::doc_memory_usage("start of get_jacobian()");
73  // MemoryUsage::insert_comment_to_continous_top("start of
74  // get_jacobian()");
75 
76  // create the residual
77  DoubleVector residual;
78 
79  // create the jacobian
80  CRDoubleMatrix* cr_matrix_pt = new CRDoubleMatrix;
81  Oomph_matrix_pt = cr_matrix_pt;
82  problem_pt->get_jacobian(residual, *cr_matrix_pt);
83  this->build_distribution(residual.distribution_pt());
84 
85  // record the end time and compute the matrix setup time
86  double end_t = TimingHelpers::timer();
87  Jacobian_setup_time = end_t - start_t;
88  if (this->Doc_time)
89  {
90  oomph_info << "Time to generate Jacobian [sec] : "
91  << Jacobian_setup_time << std::endl;
92  }
93 
94 
95  // MemoryUsage::doc_memory_usage("after get_jacobian() in trilinos
96  // solver"); MemoryUsage::insert_comment_to_continous_top(
97  // "after get_jacobian() in trilinos solver");
98 
99  // store the distribution of the solution vector
100  if (!solution.built())
101  {
102  solution.build(this->distribution_pt(), 0.0);
103  }
104  LinearAlgebraDistribution solution_dist(solution.distribution_pt());
105 
106  // redistribute the distribution
107  solution.redistribute(this->distribution_pt());
108 
109  // MemoryUsage::doc_memory_usage("before trilinos solve");
110  // MemoryUsage::insert_comment_to_continous_top("before trilinos solve ");
111 
112  // continue solving using matrix based solve function
113  solve(Oomph_matrix_pt, residual, solution);
114 
115 
116  // MemoryUsage::doc_memory_usage("after trilinos solve");
117  // MemoryUsage::insert_comment_to_continous_top("after trilinos solve ");
118 
119  // return to the original distribution
120  solution.redistribute(&solution_dist);
121 
122  // MemoryUsage::doc_memory_usage("end of TrilinosAztecOOSolver::solve");
123  // MemoryUsage::insert_comment_to_continous_top(
124  // "end of TrilinosAztecOOSolver::solve");
125  }
126 
127 
128  //=============================================================================
129  /// Function to solve the linear system defined by matrix_pt and rhs.
130  /// \b NOTE 1. The matrix has to be of type CRDoubleMatrix or
131  /// DistributedCRDoubleMatrix.
132  /// \b NOTE 2. This function will delete any existing internal data and
133  /// generate a new AztecOO solver.
134  //=============================================================================
136  const DoubleVector& rhs,
137  DoubleVector& result)
138  {
139  // start the timer
140  double start_t = TimingHelpers::timer();
141 
142 #ifdef PARANOID
143  // check that the matrix is square
144  if (matrix_pt->nrow() != matrix_pt->ncol())
145  {
146  std::ostringstream error_message_stream;
147  error_message_stream << "The matrix at matrix_pt must be square.";
148  throw OomphLibError(error_message_stream.str(),
149  OOMPH_CURRENT_FUNCTION,
150  OOMPH_EXCEPTION_LOCATION);
151  }
152  // check that the matrix and the rhs vector have the same nrow()
153  if (matrix_pt->nrow() != rhs.nrow())
154  {
155  std::ostringstream error_message_stream;
156  error_message_stream
157  << "The matrix and the rhs vector must have the same number of rows.";
158  throw OomphLibError(error_message_stream.str(),
159  OOMPH_CURRENT_FUNCTION,
160  OOMPH_EXCEPTION_LOCATION);
161  }
162 
163  // if the matrix is distributable then it too should have the same
164  // communicator as the rhs vector and should not be distributed
165  CRDoubleMatrix* cr_matrix_pt = dynamic_cast<CRDoubleMatrix*>(matrix_pt);
166  if (cr_matrix_pt != 0)
167  {
169  if (!(temp_comm == *cr_matrix_pt->distribution_pt()->communicator_pt()))
170  {
171  std::ostringstream error_message_stream;
172  error_message_stream
173  << "The matrix matrix_pt must have the same communicator as the "
174  "vectors"
175  << " rhs and result must have the same communicator";
176  throw OomphLibError(error_message_stream.str(),
177  OOMPH_CURRENT_FUNCTION,
178  OOMPH_EXCEPTION_LOCATION);
179  }
180  }
181  else
182  {
183  throw OomphLibError("Matrix must be of type CRDoubleMatrix",
184  OOMPH_CURRENT_FUNCTION,
185  OOMPH_EXCEPTION_LOCATION);
186  }
187 
188  // if the result vector is setup then check it is not distributed and has
189  // the same communicator as the rhs vector
190  if (result.built())
191  {
192  if (!(*result.distribution_pt() == *rhs.distribution_pt()))
193  {
194  std::ostringstream error_message_stream;
195  error_message_stream
196  << "The result vector distribution has been setup; it must have the "
197  << "same distribution as the rhs vector.";
198  throw OomphLibError(error_message_stream.str(),
199  OOMPH_CURRENT_FUNCTION,
200  OOMPH_EXCEPTION_LOCATION);
201  }
202  }
203 #endif
204 
205  // build the result if not built
206  if (!result.built())
207  {
208  result.build(rhs.distribution_pt(), 0.0);
209  }
210 
212  {
213  // Only call the setup method if this is the first time we call the
214  // solve method
216  {
217  // setup the solver
218  solver_setup(matrix_pt);
219  // Do not call the setup again
221  // Enable resolve since we are not going to build the solver, the
222  // matrix and the wrapper to the preconditioner again
223  Enable_resolve = true;
224  }
225  }
226  else
227  {
228  // setup the solver
229  solver_setup(matrix_pt);
230  }
231 
232  // create Epetra version of r
233  Epetra_Vector* epetra_r_pt =
235 
236  // create an empty Epetra vector for z
237  Epetra_Vector* epetra_z_pt =
239 
240  double start_t_trilinos = TimingHelpers::timer();
241 
242  // solve the system
243  solve_using_AztecOO(epetra_r_pt, epetra_z_pt);
244 
245  double end_t_trilinos = TimingHelpers::timer();
246  if (this->Doc_time)
247  {
248  oomph_info << "Time for trilinos solve itself : "
249  << end_t_trilinos - start_t_trilinos << "s" << std::endl;
250  }
251  // Copy result to z
253 
254  // clean up memory
255  delete epetra_r_pt;
256  delete epetra_z_pt;
257 
258  // delete solver data if required
259  if (!Enable_resolve)
260  {
261  clean_up_memory();
262  }
263 
264  // stop timers and compute solve time
265  double end_t = TimingHelpers::timer();
266  Linear_solver_solution_time = end_t - start_t;
267 
268  // output timings and info
269  if (this->Doc_time)
270  {
271  oomph_info << "Time for complete trilinos solve : "
272  << Linear_solver_solution_time << "s" << std::endl;
273  }
274  }
275 
276 
277  //=============================================================================
278  /// Helper function for setting up the solver. Converts the oomph-lib
279  /// matrices to Epetra matrices, sets up the preconditioner, creates the
280  /// Trilinos Aztec00 solver and passes in the matrix, preconditioner and
281  /// parameters.
282  //=============================================================================
284  {
285  // clean up the memory
286  // - delete all except Oomph_matrix_pt, which may have been set in the
287  // problem based solve
288  clean_up_memory();
289 
290  // cast to CRDoubleMatrix
291  // note cast check performed in matrix based solve(...) method
292  CRDoubleMatrix* cast_matrix_pt = dynamic_cast<CRDoubleMatrix*>(matrix_pt);
293 
294  // store the distribution
295  // distribution of preconditioner is same as matrix
296  this->build_distribution(cast_matrix_pt->distribution_pt());
297 
298  // create the new solver
299  AztecOO_solver_pt = new AztecOO();
300 
301  // if the preconditioner is an oomph-lib preconditioner then we set it up
302  TrilinosPreconditionerBase* trilinos_prec_pt =
304  if (trilinos_prec_pt == 0)
305  {
307  {
308  // setup the preconditioner
309  // start of prec setup
310  double prec_setup_start_t = TimingHelpers::timer();
311  Preconditioner_pt->setup(matrix_pt);
312  // start of prec setup
313  double prec_setup_finish_t = TimingHelpers::timer();
314  if (Doc_time)
315  {
316  double t_prec_setup = prec_setup_finish_t - prec_setup_start_t;
317  oomph_info << "Time for preconditioner setup [sec]: " << t_prec_setup
318  << std::endl;
319  }
320 #ifdef PARANOID
321  if (*Preconditioner_pt->distribution_pt() != *this->distribution_pt())
322  {
323  std::ostringstream error_message;
324  error_message << "The oomph-lib preconditioner and the solver must "
325  << "have the same distribution";
326  throw OomphLibError(error_message.str(),
327  OOMPH_CURRENT_FUNCTION,
328  OOMPH_EXCEPTION_LOCATION);
329  }
330 #endif
331  }
332 
333  // wrap the oomphlib preconditioner in the Epetra_Operator derived
334  // OoomphLibPreconditionerEpetraOperator to allow it to be passed to the
335  // trilinos preconditioner
338  }
339 
340  // create the matrix
341  double start_t_matrix = TimingHelpers::timer();
343  {
346  cast_matrix_pt);
347  }
348  else
349  {
352  cast_matrix_pt, cast_matrix_pt->distribution_pt());
353  }
354 
355  // record the end time and compute the matrix setup time
356  double end_t_matrix = TimingHelpers::timer();
357  if (trilinos_prec_pt == 0)
358  {
360  {
361  dynamic_cast<CRDoubleMatrix*>(Oomph_matrix_pt)->clear();
362  delete Oomph_matrix_pt;
363  Oomph_matrix_pt = NULL;
364  }
365 
366  // delete Oomph-lib matrix if requested
367  else if (Delete_matrix)
368  {
369  dynamic_cast<CRDoubleMatrix*>(matrix_pt)->clear();
370  }
371  }
372 
373  // output times
374  if (Doc_time)
375  {
376  oomph_info << "Time to generate Trilinos matrix : "
377  << double(end_t_matrix - start_t_matrix) << "s" << std::endl;
378  }
379 
380  // set the matrix
381  AztecOO_solver_pt->SetUserMatrix(Epetra_matrix_pt);
382 
383  // set the preconditioner
384  if (trilinos_prec_pt == 0)
385  {
387  }
388 
389 #ifdef PARANOID
390  // paranoid check the preconditioner exists
391  if (Preconditioner_pt == 0)
392  {
393  std::ostringstream error_message;
394  error_message << "Preconditioner_pt == 0. (Remember default "
395  << "preconditioner is IdentityPreconditioner)";
396  throw OomphLibError(
397  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
398  }
399 #endif
400 
401  // if the preconditioner is a trilinos preconditioner then setup the
402  // preconditioner
403  if (trilinos_prec_pt != 0)
404  {
405  // only setup the preconditioner if required
407  {
408  // start of prec setup
409  double prec_setup_start_t = TimingHelpers::timer();
410 
411  // setup the preconditioner
412  trilinos_prec_pt->set_matrix_pt(matrix_pt);
413  trilinos_prec_pt->setup(Epetra_matrix_pt);
414 
415 
416  // set the preconditioner
417  AztecOO_solver_pt->SetPrecOperator(
418  trilinos_prec_pt->epetra_operator_pt());
419 
420  // start of prec setup
421  double prec_setup_finish_t = TimingHelpers::timer();
422  if (Doc_time)
423  {
424  double t_prec_setup = prec_setup_finish_t - prec_setup_start_t;
425  oomph_info << "Time for preconditioner setup [sec]: " << t_prec_setup
426  << std::endl;
427  }
428  }
429 
430  // delete the oomph-matrix if required
432  {
433  dynamic_cast<CRDoubleMatrix*>(Oomph_matrix_pt)->clear();
434  delete Oomph_matrix_pt;
435  Oomph_matrix_pt = NULL;
436  }
437 
438  // delete Oomph-lib matrix if requested
439  else if (Delete_matrix)
440  {
441  dynamic_cast<CRDoubleMatrix*>(matrix_pt)->clear();
442  }
443  }
444 
445  // set solver options
446  if (Doc_time)
447  {
448  AztecOO_solver_pt->SetAztecOption(AZ_output, AZ_warnings);
449  }
450  else
451  {
452  AztecOO_solver_pt->SetAztecOption(AZ_output, AZ_none);
453  }
454  AztecOO_solver_pt->SetAztecOption(AZ_kspace, Max_iter);
455 
456  // set solver type
457  switch (Solver_type)
458  {
459  case CG:
460  AztecOO_solver_pt->SetAztecOption(AZ_solver, AZ_cg);
461  break;
462 
463  case GMRES:
464  AztecOO_solver_pt->SetAztecOption(AZ_solver, AZ_gmres);
465  break;
466 
467  case BiCGStab:
468  AztecOO_solver_pt->SetAztecOption(AZ_solver, AZ_bicgstab);
469  break;
470 
471  default:
472  std::ostringstream error_message;
473  error_message << "Solver_type set to incorrect value. "
474  << "Acceptable values are " << CG << ", " << GMRES
475  << " and " << BiCGStab << ". Current value is "
476  << Solver_type << ".";
477  throw OomphLibError(error_message.str(),
478  OOMPH_CURRENT_FUNCTION,
479  OOMPH_EXCEPTION_LOCATION);
480  }
481  }
482 
483  //=============================================================================
484  /// Function to resolve a linear system using the existing solver
485  /// data, allowing a solve with a new right hand side vector. This
486  /// function must be used after a call to solve(...) with
487  /// enable_resolve set to true.
488  //=============================================================================
490  DoubleVector& solution)
491  {
492  // start the timer
493  double start_t = TimingHelpers::timer();
494 
495 #ifdef PARANOID
496  if (Epetra_matrix_pt->NumGlobalRows() != static_cast<int>(rhs.nrow()))
497  {
498  std::ostringstream error_message;
499  error_message
500  << "The rhs vector and the matrix must have the same number "
501  << "of rows.\n"
502  << "The rhs vector has " << rhs.nrow() << " rows.\n"
503  << "The matrix has " << Epetra_matrix_pt->NumGlobalRows() << " rows.\n";
504  throw OomphLibError(
505  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
506  }
507 
508  // if the result vector is setup then check it is not distributed and has
509  // the same communicator as the rhs vector
510  if (solution.built())
511  {
512  if (!(*solution.distribution_pt() == *rhs.distribution_pt()))
513  {
514  std::ostringstream error_message_stream;
515  error_message_stream
516  << "The result vector distribution has been setup; it must have the "
517  << "same distribution as the rhs vector.";
518  throw OomphLibError(error_message_stream.str(),
519  OOMPH_CURRENT_FUNCTION,
520  OOMPH_EXCEPTION_LOCATION);
521  }
522  }
523 #endif
524 
525  // build the result if not built
526  if (!solution.built())
527  {
528  solution.build(rhs.distribution_pt(), 0.0);
529  }
530 
531 
532  // create Epetra version of r
533  Epetra_Vector* epetra_r_pt =
535 
536  // create an empty Epetra vector for z
537  Epetra_Vector* epetra_z_pt =
539 
540  // solve the system
541  solve_using_AztecOO(epetra_r_pt, epetra_z_pt);
542 
543  // Copy result to z
544  if (!solution.distribution_built())
545  {
546  solution.build(rhs.distribution_pt(), 0.0);
547  }
548  TrilinosEpetraHelpers::copy_to_oomphlib_vector(epetra_z_pt, solution);
549 
550  // clean up memory
551  delete epetra_r_pt;
552  delete epetra_z_pt;
553 
554  double end_t = TimingHelpers::timer();
555  Linear_solver_solution_time = end_t - start_t;
556 
557  // output timings and info
558  if (this->Doc_time)
559  {
560  oomph_info << "Time for resolve : "
561  << Linear_solver_solution_time << "s" << std::endl;
562  }
563  }
564 
565 
566  //=============================================================================
567  /// Helper function performs the actual solve once the AztecOO
568  /// solver is set up (i.e. solver_setup() is called)
569  //=============================================================================
570  void TrilinosAztecOOSolver::solve_using_AztecOO(Epetra_Vector*& rhs_pt,
571  Epetra_Vector*& soln_pt)
572  {
573 #ifdef PARANOID
574  // check the matrix and rhs are of consistent sizes
575  if (AztecOO_solver_pt == 0)
576  {
577  std::ostringstream error_message;
578  error_message << "Solver must be called with solve(...) "
579  << "before resolve(...) to set it up.\n";
580  throw OomphLibError(
581  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
582  }
583 #endif
584 
585  // set the vectors
586  AztecOO_solver_pt->SetLHS(soln_pt);
587  AztecOO_solver_pt->SetRHS(rhs_pt);
588 
589  // perform solve
591 
592 
593  // output iterations and final norm
594  Iterations = AztecOO_solver_pt->NumIters();
595  if (Doc_time)
596  {
597  double norm = AztecOO_solver_pt->TrueResidual();
598  oomph_info << "Linear solver iterations : " << Iterations << std::endl;
599  oomph_info << "Final Relative Residual Norm: " << norm << std::endl;
600  }
601  }
602 
603 
604  /// ////////////////////////////////////////////////////////////////////////////
605  /// ////////////////////////////////////////////////////////////////////////////
606  /// ////////////////////////////////////////////////////////////////////////////
607 
608 
609 } // namespace oomph
//////////////////////////////////////////////////////////////////////////// ////////////////////////...
//////////////////////////////////////////////////////////////////////////// ////////////////////////...
A class for compressed row matrices. This is a distributable object.
Definition: matrices.h:888
LinearAlgebraDistribution * distribution_pt() const
access to the LinearAlgebraDistribution
unsigned nrow() const
access function to the number of global rows.
bool distribution_built() const
if the communicator_pt is null then the distribution is not setup then false is returned,...
void build_distribution(const LinearAlgebraDistribution *const dist_pt)
setup the distribution of this distributable linear algebra object
Abstract base class for matrices of doubles – adds abstract interfaces for solving,...
Definition: matrices.h:261
virtual unsigned long ncol() const =0
Return the number of columns of the matrix.
virtual unsigned long nrow() const =0
Return the number of rows of the matrix.
A vector in the mathematical sense, initially developed for linear algebra type applications....
Definition: double_vector.h:58
void build(const DoubleVector &old_vector)
Just copys the argument DoubleVector.
void redistribute(const LinearAlgebraDistribution *const &dist_pt)
The contents of the vector are redistributed to match the new distribution. In a non-MPI rebuild this...
bool built() const
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
bool Use_iterative_solver_as_preconditioner
Use the iterative solver as preconditioner.
bool First_time_solve_when_used_as_preconditioner
When the iterative solver is used a preconditioner then we call the setup of solver method only once ...
double Tolerance
Convergence tolerance.
Preconditioner * Preconditioner_pt
Pointer to the preconditioner.
unsigned Max_iter
Maximum number of iterations.
bool Setup_preconditioner_before_solve
indicates whether the preconditioner should be setup before solve. Default = true;
Describes the distribution of a distributable linear algebra type object. Typically this is a contain...
OomphCommunicator * communicator_pt() const
const access to the communicator pointer
bool Doc_time
Boolean flag that indicates whether the time taken.
Definition: linear_solver.h:77
bool Enable_resolve
Boolean that indicates whether the matrix (or its factors, in the case of direct solver) should be st...
Definition: linear_solver.h:73
An oomph-lib wrapper to the MPI_Comm communicator object. Just contains an MPI_Comm object (which is ...
Definition: communicator.h:54
An OomphLibError object which should be thrown when an run-time error is encountered....
An Epetra_Operator class for oomph-lib preconditioners. A helper class for TrilinosOomphLibPreconditi...
void setup(DoubleMatrixBase *matrix_pt)
Setup the preconditioner: store the matrix pointer and the communicator pointer then call preconditio...
virtual void set_matrix_pt(DoubleMatrixBase *matrix_pt)
Set the matrix pointer.
////////////////////////////////////////////////////////////////// //////////////////////////////////...
Definition: problem.h:154
virtual void get_jacobian(DoubleVector &residuals, DenseDoubleMatrix &jacobian)
Return the fully-assembled Jacobian and residuals for the problem Interface for the case when the Jac...
Definition: problem.cc:3922
bool Using_problem_based_solve
Helper flag keeping track of whether we called the linear algebra or problem-based solve function.
void resolve(const DoubleVector &rhs, DoubleVector &solution)
Function to resolve a linear system using the existing solver data, allowing a solve with a new right...
void solver_setup(DoubleMatrixBase *const &matrix_pt)
Helper function for setting up the solver. Converts the oomph-lib matrices to Epetra matrices,...
unsigned Iterations
Stores number of iterations used.
Epetra_CrsMatrix * Epetra_matrix_pt
A pointer for the linear system matrix in Epetra_CrsMatrix format.
void clean_up_memory()
Clean up method - deletes the solver, the matrices and the preconditioner.
void solve_using_AztecOO(Epetra_Vector *&rhs_pt, Epetra_Vector *&soln_pt)
Helper function performs the actual solve once the AztecOO solver is set up.
bool Delete_matrix
Trilinos copies matrix data from oomph-lib's own CRDoubleMatrix or DistributedCRDoubleMatrix to Trili...
DoubleMatrixBase * Oomph_matrix_pt
Oomph lib matrix pointer.
Epetra_Operator * Epetra_preconditioner_pt
A pointer to the Epetra_Operator for the preconditioner. This is only used if the preconditioner NOT ...
void solve(Problem *const &problem_pt, DoubleVector &solution)
Function which uses problem_pt's get_jacobian(...) function to generate a linear system which is then...
unsigned Solver_type
Defines which solver is set up - available types are defined in AztecOO_solver_types.
double Linear_solver_solution_time
Stores time for the solution (excludes time to set up preconditioner)
Problem * Problem_pt
A pointer to the underlying problem (NULL if MATRIX based solve) The problem_pt is stored here in a p...
AztecOO * AztecOO_solver_pt
Pointer to the AztecOO solver.
bool Use_aztecoo_workaround_for_epetra_matrix_setup
Use workaround for creating of epetra matrix that respects aztecoo's ordering requirements.
double Jacobian_setup_time
Stores set up time for Jacobian.
Base class for Trilinos preconditioners as oomph-lib preconditioner.
Epetra_Operator *& epetra_operator_pt()
Access function to Epetra_preconditioner_pt. For use with TrilinosAztecOOSolver.
void setup()
Broken assignment operator.
double timer()
returns the time in seconds after some point in past
Epetra_CrsMatrix * create_distributed_epetra_matrix(const CRDoubleMatrix *oomph_matrix_pt, const LinearAlgebraDistribution *dist_pt)
create an Epetra_CrsMatrix from an oomph-lib CRDoubleMatrix. If oomph_matrix_pt is NOT distributed (i...
Epetra_CrsMatrix * create_distributed_epetra_matrix_for_aztecoo(CRDoubleMatrix *oomph_matrix_pt)
create and Epetra_CrsMatrix from an oomph-lib CRDoubleMatrix. Specialisation for Trilinos AztecOO....
void copy_to_oomphlib_vector(const Epetra_Vector *epetra_vec_pt, DoubleVector &oomph_vec)
Helper function to copy the contents of a Trilinos vector to an oomph-lib distributed vector....
Epetra_Vector * create_distributed_epetra_vector(const DoubleVector &oomph_vec)
create an Epetra_Vector from an oomph-lib DoubleVector. If oomph_vec is NOT distributed (i....
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...