libMesh
laspack_linear_solver.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #include "libmesh/libmesh_common.h"
21 
22 #if defined(LIBMESH_HAVE_LASPACK)
23 
24 
25 // Local Includes
26 #include "libmesh/laspack_linear_solver.h"
27 #include "libmesh/libmesh_logging.h"
28 #include "libmesh/enum_to_string.h"
29 #include "libmesh/enum_solver_type.h"
30 #include "libmesh/enum_preconditioner_type.h"
31 #include "libmesh/enum_convergence_flags.h"
32 
33 namespace libMesh
34 {
35 
36 // #ifndef LIBMESH_USE_COMPLEX_NUMBERS
37 // extern "C"
38 // {
39 // #endif
40 
41 // void print_iter_accuracy(int Iter,
42 // _LPReal rNorm,
43 // _LPReal bNorm,
44 // IterIdType IterId)
45 // /* put out accuracy reached after each solver iteration */
46 // {
47 
48 // //FILE * out = fopen("residual.hist", "a");
49 // static int icall=0;
50 
51 // if (!icall)
52 // {
53 // printf("Iter ||r||/||f||\n");
54 // printf("------------------\n");
55 // icall=1;
56 // }
57 
58 // if ( Iter%1==0 && (IterId == CGIterId ||
59 // IterId == CGNIterId ||
60 // IterId == GMRESIterId ||
61 // IterId == BiCGIterId ||
62 // IterId == QMRIterId ||
63 // IterId == CGSIterId ||
64 // IterId == BiCGSTABIterId) )
65 // {
66 // if (!_LPIsZeroReal(bNorm))
67 // printf("%d \t %g\n", Iter, rNorm/bNorm);
68 // else
69 // printf("%d (fnorm == 0)\n", Iter);
70 // }
71 
72 // //fclose(out);
73 // }
74 
75 // #ifndef LIBMESH_USE_COMPLEX_NUMBERS
76 // }
77 // #endif
78 
79 /*----------------------- functions ----------------------------------*/
80 template <typename T>
82 {
83  if (this->initialized())
84  {
85  this->_is_initialized = false;
86 
87  this->_solver_type = GMRES;
88  this->_preconditioner_type = ILU_PRECOND;
89  }
90 }
91 
92 
93 
94 template <typename T>
95 void LaspackLinearSolver<T>::init (const char * /* name */)
96 {
97  // Initialize the data structures if not done so already.
98  if (!this->initialized())
99  {
100  this->_is_initialized = true;
101  }
102 
103  // SetRTCAuxProc (print_iter_accuracy);
104 }
105 
106 
107 
108 template <typename T>
109 std::pair<unsigned int, Real>
111  NumericVector<T> & solution_in,
112  NumericVector<T> & rhs_in,
113  const std::optional<double> tol,
114  const std::optional<unsigned int> m_its)
115 {
116  LOG_SCOPE("solve()", "LaspackLinearSolver");
117  this->init ();
118 
119  // Make sure the data passed in are really in Laspack types
120  LaspackMatrix<T> * matrix = cast_ptr<LaspackMatrix<T> *>(&matrix_in);
121  LaspackVector<T> * solution = cast_ptr<LaspackVector<T> *>(&solution_in);
122  LaspackVector<T> * rhs = cast_ptr<LaspackVector<T> *>(&rhs_in);
123 
124  const int max_its = this->get_int_solver_setting("max_its", m_its);
125  const double abs_tol = this->get_real_solver_setting("abs_tol", tol);
126 
127  // Zero-out the solution to prevent the solver from exiting in 0
128  // iterations (?)
129  //TODO:[BSK] Why does Laspack do this? Comment out this and try ex13...
130  solution->zero();
131 
132  // Close the matrix and vectors in case this wasn't already done.
133  matrix->close ();
134  solution->close ();
135  rhs->close ();
136 
137  // Set the preconditioner type
138  this->set_laspack_preconditioner_type ();
139 
140  // Set the solver tolerance
141  SetRTCAccuracy (abs_tol);
142 
143  // Solve the linear system
144  switch (this->_solver_type)
145  {
146  // Conjugate-Gradient
147  case CG:
148  {
149  CGIter (&matrix->_QMat,
150  &solution->_vec,
151  &rhs->_vec,
152  max_its,
153  _precond_type,
154  1.);
155  break;
156  }
157 
158  // Conjugate-Gradient Normalized
159  case CGN:
160  {
161  CGNIter (&matrix->_QMat,
162  &solution->_vec,
163  &rhs->_vec,
164  max_its,
165  _precond_type,
166  1.);
167  break;
168  }
169 
170  // Conjugate-Gradient Squared
171  case CGS:
172  {
173  CGSIter (&matrix->_QMat,
174  &solution->_vec,
175  &rhs->_vec,
176  max_its,
177  _precond_type,
178  1.);
179  break;
180  }
181 
182  // Bi-Conjugate Gradient
183  case BICG:
184  {
185  BiCGIter (&matrix->_QMat,
186  &solution->_vec,
187  &rhs->_vec,
188  max_its,
189  _precond_type,
190  1.);
191  break;
192  }
193 
194  // Bi-Conjugate Gradient Stabilized
195  case BICGSTAB:
196  {
197  BiCGSTABIter (&matrix->_QMat,
198  &solution->_vec,
199  &rhs->_vec,
200  max_its,
201  _precond_type,
202  1.);
203  break;
204  }
205 
206  // Quasi-Minimum Residual
207  case QMR:
208  {
209  QMRIter (&matrix->_QMat,
210  &solution->_vec,
211  &rhs->_vec,
212  max_its,
213  _precond_type,
214  1.);
215  break;
216  }
217 
218  // Symmetric over-relaxation
219  case SSOR:
220  {
221  SSORIter (&matrix->_QMat,
222  &solution->_vec,
223  &rhs->_vec,
224  max_its,
225  _precond_type,
226  1.);
227  break;
228  }
229 
230  // Jacobi Relaxation
231  case JACOBI:
232  {
233  JacobiIter (&matrix->_QMat,
234  &solution->_vec,
235  &rhs->_vec,
236  max_its,
237  _precond_type,
238  1.);
239  break;
240  }
241 
242  // Generalized Minimum Residual
243  case GMRES:
244  {
245  SetGMRESRestart (30);
246  GMRESIter (&matrix->_QMat,
247  &solution->_vec,
248  &rhs->_vec,
249  max_its,
250  _precond_type,
251  1.);
252  break;
253  }
254 
255  // Unknown solver, use GMRES
256  default:
257  {
258  libMesh::err << "ERROR: Unsupported LASPACK Solver: "
259  << Utility::enum_to_string(this->_solver_type) << std::endl
260  << "Continuing with GMRES" << std::endl;
261 
262  this->_solver_type = GMRES;
263 
264  return this->solve (*matrix,
265  *solution,
266  *rhs,
267  tol,
268  m_its);
269  }
270  }
271 
272  // Check for an error
273  if (LASResult() != LASOK)
274  {
275  WriteLASErrDescr(stdout);
276  libmesh_error_msg("Exiting after LASPACK Error!");
277  }
278 
279  // Get the convergence step # and residual
280  return std::make_pair(GetLastNoIter(), GetLastAccuracy());
281 }
282 
283 
284 
285 template <typename T>
286 std::pair<unsigned int, Real>
288  NumericVector<T> & solution_in,
289  NumericVector<T> & rhs_in,
290  const std::optional<double> tol,
291  const std::optional<unsigned int> m_its)
292 {
293  LOG_SCOPE("adjoint_solve()", "LaspackLinearSolver");
294  this->init ();
295 
296  // Make sure the data passed in are really in Laspack types
297  LaspackMatrix<T> * matrix = cast_ptr<LaspackMatrix<T> *>(&matrix_in);
298  LaspackVector<T> * solution = cast_ptr<LaspackVector<T> *>(&solution_in);
299  LaspackVector<T> * rhs = cast_ptr<LaspackVector<T> *>(&rhs_in);
300 
301  const int max_its = this->get_int_solver_setting("max_its", m_its);
302  const double abs_tol = this->get_real_solver_setting("abs_tol", tol);
303 
304  // Zero-out the solution to prevent the solver from exiting in 0
305  // iterations (?)
306  //TODO:[BSK] Why does Laspack do this? Comment out this and try ex13...
307  solution->zero();
308 
309  // Close the matrix and vectors in case this wasn't already done.
310  matrix->close ();
311  solution->close ();
312  rhs->close ();
313 
314  // Set the preconditioner type
315  this->set_laspack_preconditioner_type ();
316 
317  // Set the solver tolerance
318  SetRTCAccuracy (abs_tol);
319 
320  // Solve the linear system
321  switch (this->_solver_type)
322  {
323  // Conjugate-Gradient
324  case CG:
325  {
326  CGIter (Transp_Q(&matrix->_QMat),
327  &solution->_vec,
328  &rhs->_vec,
329  max_its,
330  _precond_type,
331  1.);
332  break;
333  }
334 
335  // Conjugate-Gradient Normalized
336  case CGN:
337  {
338  CGNIter (Transp_Q(&matrix->_QMat),
339  &solution->_vec,
340  &rhs->_vec,
341  max_its,
342  _precond_type,
343  1.);
344  break;
345  }
346 
347  // Conjugate-Gradient Squared
348  case CGS:
349  {
350  CGSIter (Transp_Q(&matrix->_QMat),
351  &solution->_vec,
352  &rhs->_vec,
353  max_its,
354  _precond_type,
355  1.);
356  break;
357  }
358 
359  // Bi-Conjugate Gradient
360  case BICG:
361  {
362  BiCGIter (Transp_Q(&matrix->_QMat),
363  &solution->_vec,
364  &rhs->_vec,
365  max_its,
366  _precond_type,
367  1.);
368  break;
369  }
370 
371  // Bi-Conjugate Gradient Stabilized
372  case BICGSTAB:
373  {
374  BiCGSTABIter (Transp_Q(&matrix->_QMat),
375  &solution->_vec,
376  &rhs->_vec,
377  max_its,
378  _precond_type,
379  1.);
380  break;
381  }
382 
383  // Quasi-Minimum Residual
384  case QMR:
385  {
386  QMRIter (Transp_Q(&matrix->_QMat),
387  &solution->_vec,
388  &rhs->_vec,
389  max_its,
390  _precond_type,
391  1.);
392  break;
393  }
394 
395  // Symmetric over-relaxation
396  case SSOR:
397  {
398  SSORIter (Transp_Q(&matrix->_QMat),
399  &solution->_vec,
400  &rhs->_vec,
401  max_its,
402  _precond_type,
403  1.);
404  break;
405  }
406 
407  // Jacobi Relaxation
408  case JACOBI:
409  {
410  JacobiIter (Transp_Q(&matrix->_QMat),
411  &solution->_vec,
412  &rhs->_vec,
413  max_its,
414  _precond_type,
415  1.);
416  break;
417  }
418 
419  // Generalized Minimum Residual
420  case GMRES:
421  {
422  SetGMRESRestart (30);
423  GMRESIter (Transp_Q(&matrix->_QMat),
424  &solution->_vec,
425  &rhs->_vec,
426  max_its,
427  _precond_type,
428  1.);
429  break;
430  }
431 
432  // Unknown solver, use GMRES
433  default:
434  {
435  libMesh::err << "ERROR: Unsupported LASPACK Solver: "
436  << Utility::enum_to_string(this->_solver_type) << std::endl
437  << "Continuing with GMRES" << std::endl;
438 
439  this->_solver_type = GMRES;
440 
441  return this->solve (*matrix,
442  *solution,
443  *rhs,
444  tol,
445  m_its);
446  }
447  }
448 
449  // Check for an error
450  if (LASResult() != LASOK)
451  {
452  WriteLASErrDescr(stdout);
453  libmesh_error_msg("Exiting after LASPACK Error!");
454  }
455 
456  // Get the convergence step # and residual
457  return std::make_pair(GetLastNoIter(), GetLastAccuracy());
458 }
459 
460 
461 
462 
463 template <typename T>
464 std::pair<unsigned int, Real>
466  NumericVector<T> & /*solution_in*/,
467  NumericVector<T> & /*rhs_in*/,
468  const std::optional<double> /*tol*/,
469  const std::optional<unsigned int> /*m_its*/)
470 {
471  libmesh_not_implemented();
472  return std::make_pair(0,0.0);
473 }
474 
475 
476 
477 template <typename T>
478 std::pair<unsigned int, Real>
480  const SparseMatrix<T> & /*precond_matrix*/,
481  NumericVector<T> & /*solution_in*/,
482  NumericVector<T> & /*rhs_in*/,
483  const std::optional<double> /*tol*/,
484  const std::optional<unsigned int> /*m_its*/)
485 {
486  libmesh_not_implemented();
487  return std::make_pair(0,0.0);
488 }
489 
490 
491 
492 template <typename T>
494 {
495  switch (this->_preconditioner_type)
496  {
497  case IDENTITY_PRECOND:
498  _precond_type = nullptr; return;
499 
500  case ILU_PRECOND:
501  _precond_type = ILUPrecond; return;
502 
503  case JACOBI_PRECOND:
504  _precond_type = JacobiPrecond; return;
505 
506  case SSOR_PRECOND:
507  _precond_type = SSORPrecond; return;
508 
509 
510  default:
511  libMesh::err << "ERROR: Unsupported LASPACK Preconditioner: "
512  << this->_preconditioner_type << std::endl
513  << "Continuing with ILU" << std::endl;
514  this->_preconditioner_type = ILU_PRECOND;
515  this->set_laspack_preconditioner_type();
516  }
517 }
518 
519 
520 
521 template <typename T>
523 {
524  switch (LASResult())
525  {
526  case LASOK :
527  libMesh::out << "Laspack converged.\n";
528  break;
529  default :
530  libMesh::out << "Laspack diverged.\n";
531  }
532  libMesh::out << "Detailed reporting is currently only supported"
533  << "with Petsc 2.3.1 and later." << std::endl;
534 }
535 
536 
537 
538 template <typename T>
540 {
541  switch (LASResult())
542  {
543  case LASOK : return CONVERGED_RTOL_NORMAL;
544  default : return DIVERGED_NULL;
545  }
546 }
547 
548 
549 
550 //------------------------------------------------------------------
551 // Explicit instantiations
552 template class LIBMESH_EXPORT LaspackLinearSolver<Number>;
553 
554 } // namespace libMesh
555 
556 
557 #endif // #ifdef LIBMESH_HAVE_LASPACK
OStreamProxy err
virtual void print_converged_reason() const override
Prints a useful message about why the latest linear solve con(di)verged.
void set_laspack_preconditioner_type()
Tells LASPACK to use the user-specified preconditioner stored in _preconditioner_type.
LinearConvergenceReason
Linear solver convergence flags (taken from the PETSc flags).
virtual void init(const char *name=nullptr) override
Initialize data structures if not done so already.
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
Definition: vector_fe_ex5.C:44
QMatrix _QMat
The Laspack sparse matrix pointer.
The libMesh namespace provides an interface to certain functionality in the library.
Generic sparse matrix.
Definition: vector_fe_ex5.C:46
This class provides an interface to Laspack iterative solvers that is compatible with the libMesh Lin...
virtual void close() override
Calls the SparseMatrix&#39;s internal assembly routines, ensuring that the values are consistent across p...
bool _is_initialized
Flag that tells if init() has been called.
Definition: libmesh.C:257
virtual void clear() override
Release all memory and clear data structures.
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
virtual LinearConvergenceReason get_converged_reason() const override
std::string enum_to_string(const T e)
This class provides a nice interface to the Laspack C-based data structures for serial vectors...
OStreamProxy out
bool initialized()
Checks that library initialization has been done.
Definition: libmesh.C:276
virtual std::pair< unsigned int, Real > adjoint_solve(SparseMatrix< T > &matrix, NumericVector< T > &solution, NumericVector< T > &rhs, const std::optional< double > tol=std::nullopt, const std::optional< unsigned int > m_its=std::nullopt) override
Call the Laspack solver to solve A^T x = b.
virtual std::pair< unsigned int, Real > solve(SparseMatrix< T > &matrix, NumericVector< T > &solution, NumericVector< T > &rhs, const std::optional< double > tol=std::nullopt, const std::optional< unsigned int > m_its=std::nullopt) override
Call the Laspack solver.
The LaspackMatrix class wraps a QMatrix object from the Laspack library.
Generic shell matrix, i.e.