trilinos_eigen_solver.h
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-2023 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 #ifndef OOMPH_TRILINOS_EIGEN_SOLVER_HEADER
27 #define OOMPH_TRILINOS_EIGEN_SOLVER_HEADER
28 
29 // oomph-lib includes
30 #include "double_multi_vector.h"
31 #include "problem.h"
32 #include "linear_solver.h"
33 #include "eigen_solver.h"
34 
35 // Trilinos includes
36 #include "AnasaziOutputManager.hpp"
37 #include "AnasaziBasicOutputManager.hpp"
38 #include "AnasaziConfigDefs.hpp"
39 #include "AnasaziOperator.hpp"
40 #include "AnasaziMultiVec.hpp"
41 #include "AnasaziBasicEigenproblem.hpp"
42 #include "AnasaziBlockKrylovSchurSolMgr.hpp"
43 
44 namespace Anasazi
45 {
46  // hierher Andrew get rid of this!
48 
49  //========================================================================
50  /// Specialize the Anasazi traits class for the oomph-lib DoubleMultiVector.
51  /// This provides the interfaces required by the Anasazi eigensolvers.
52  //========================================================================
53  template<>
54  class MultiVecTraits<double, oomph::DoubleMultiVector>
55  {
56  public:
57  /// \brief Creates a new empty DoubleMultiVector containing numvecs columns.
58  /// Return a reference-counted pointer to the new multivector
59  static Teuchos::RCP<oomph::DoubleMultiVector> Clone(
60  const oomph::DoubleMultiVector& mv, const int numvecs)
61  {
62  return Teuchos::rcp(new oomph::DoubleMultiVector(numvecs, mv));
63  }
64 
65  /// \brief Creates a deep copy of the DoubleMultiVector mv
66  /// return Reference-counted pointer to the new DoubleMultiVector
67  static Teuchos::RCP<oomph::DoubleMultiVector> CloneCopy(
68  const oomph::DoubleMultiVector& mv)
69  {
70  return Teuchos::rcp(new oomph::DoubleMultiVector(mv));
71  }
72 
73  /// \brief Creates a new oomph::DoubleMultiVector and (deep)
74  /// copies the contents of
75  /// the vectors in index into the new vector
76  /// return Reference-counted pointer to the new oomph::DoubleMultiVector
77  static Teuchos::RCP<oomph::DoubleMultiVector> CloneCopy(
78  const oomph::DoubleMultiVector& mv, const std::vector<int>& index)
79  {
80  return Teuchos::rcp(new oomph::DoubleMultiVector(mv, index));
81  }
82 
83  /// \brief Deep copy of specified columns of oomph::DoubleMultiVector
84  /// return Reference-counted pointer to the new oomph::DoubleMultiVector
85  static Teuchos::RCP<oomph::DoubleMultiVector> CloneCopy(
86  const oomph::DoubleMultiVector& mv, const Teuchos::Range1D& index)
87  {
88  return Teuchos::rcp(new oomph::DoubleMultiVector(mv, index));
89  }
90 
91  /// \brief Creates a new oomph::DoubleMultiVector that contains
92  /// shallow copies
93  /// of selected entries of the oomph::DoubleMultiVector mv
94  /// return Reference-counted pointer to the new oomph::DoubleMultiVector
95  static Teuchos::RCP<oomph::DoubleMultiVector> CloneViewNonConst(
96  oomph::DoubleMultiVector& mv, const std::vector<int>& index)
97  {
98  return Teuchos::rcp(new oomph::DoubleMultiVector(mv, index, false));
99  }
100 
101  /// \brief Creates a new oomph::DoubleMultiVector that
102  /// contains shallow copies
103  /// of selected entries of the oomph::DoubleMultiVector mv
104  /// return Reference-counted pointer to the new oomph::DoubleMultiVector
105  static Teuchos::RCP<oomph::DoubleMultiVector> CloneViewNonConst(
106  oomph::DoubleMultiVector& mv, const Teuchos::Range1D& index)
107  {
108  return Teuchos::rcp(new oomph::DoubleMultiVector(mv, index, false));
109  }
110 
111  /// \brief Creates a new oomph::DoubleMultiVector that contains
112  /// shallow copies
113  /// of selected entries of the oomph::DoubleMultiVector mv
114  /// return Reference-counted pointer to the new oomph::DoubleMultiVector
115  /// (const version)
116  static Teuchos::RCP<const oomph::DoubleMultiVector> CloneView(
117  const oomph::DoubleMultiVector& mv, const std::vector<int>& index)
118  {
119  return Teuchos::rcp(new oomph::DoubleMultiVector(mv, index, false));
120  }
121 
122  /// \brief Creates a new oomph::DoubleMultiVector that contains
123  /// shallow copies
124  /// of selected entries of the oomph::DoubleMultiVector mv
125  /// return Reference-counted pointer to the new oomph::DoubleMultiVector
126  /// (Non-const version for Trilinos 9 interface)
127  static Teuchos::RCP<oomph::DoubleMultiVector> CloneView(
128  oomph::DoubleMultiVector& mv, const std::vector<int>& index)
129  {
130  return Teuchos::rcp(new oomph::DoubleMultiVector(mv, index, false));
131  }
132 
133  /// \brief Creates a new oomph::DoubleMultiVector that
134  /// contains shallow copies
135  /// of selected entries of the oomph::DoubleMultiVector mv
136  /// return Reference-counted pointer to the new oomph::DoubleMultiVector
137  /// (const version)
138  static Teuchos::RCP<oomph::DoubleMultiVector> CloneView(
139  oomph::DoubleMultiVector& mv, const Teuchos::Range1D& index)
140  {
141  return Teuchos::rcp(new oomph::DoubleMultiVector(mv, index, false));
142  }
143 
144  /// Obtain the global length of the vector
146  {
147  return static_cast<int>(mv.nrow());
148  }
149 
150  /// Obtain the number of vectors in the multivector
152  {
153  return static_cast<int>(mv.nvector());
154  }
155 
156  /// \brief Update \c mv with \f$ \alpha AB + \beta mv \f$.
157  static void MvTimesMatAddMv(
158  const double alpha,
159  const oomph::DoubleMultiVector& A,
160  const Teuchos::SerialDenseMatrix<int, double>& B,
161  const double beta,
163  {
164  // For safety let's (deep) copy A
166  // First pre-multiply by the scalars
167  mv *= beta;
168  C *= alpha;
169 
170  // Find the number of vectors in each multivector
171  int mv_n_vector = mv.nvector();
172  int A_n_vector = A.nvector();
173  // Find the number of rows
174  int n_row_local = A.nrow_local();
175  // Now need to multiply by the matrix and add the result
176  // to the appropriate entry in the multivector
177  for (int i = 0; i < n_row_local; ++i)
178  {
179  for (int j = 0; j < mv_n_vector; j++)
180  {
181  double res = 0.0;
182  for (int k = 0; k < A_n_vector; k++)
183  {
184  res += C(k, i) * B(k, j);
185  }
186  mv(j, i) += res;
187  }
188  }
189  }
190 
191  /// \brief Replace \c mv with \f$\alpha A + \beta B\f$.
192  static void MvAddMv(const double alpha,
193  const oomph::DoubleMultiVector& A,
194  const double beta,
197  {
198  const unsigned n_vector = A.nvector();
199  const unsigned n_row_local = A.nrow_local();
200  for (unsigned v = 0; v < n_vector; v++)
201  {
202  for (unsigned i = 0; i < n_row_local; i++)
203  {
204  mv(v, i) = alpha * A(v, i) + beta * B(v, i);
205  }
206  }
207  }
208 
209  /*! \brief Scale each element of the vectors in \c mv with \c alpha.
210  */
211  static void MvScale(oomph::DoubleMultiVector& mv, const double alpha)
212  {
213  mv *= alpha;
214  }
215 
216  /*! \brief Scale each element of the \c i-th vector in \c mv with \c
217  * alpha[i].
218  */
220  const std::vector<double>& alpha)
221  {
222  const unsigned n_vector = mv.nvector();
223  const unsigned n_row_local = mv.nrow_local();
224  for (unsigned v = 0; v < n_vector; v++)
225  {
226  for (unsigned i = 0; i < n_row_local; i++)
227  {
228  mv(v, i) *= alpha[v];
229  }
230  }
231  }
232 
233  /*! \brief Compute a dense matrix \c B through the matrix-matrix multiply
234  * \f$ \alpha A^Hmv \f$.
235  */
236  static void MvTransMv(const double alpha,
237  const oomph::DoubleMultiVector& A,
238  const oomph::DoubleMultiVector& mv,
239  Teuchos::SerialDenseMatrix<int, double>& B)
240  {
241  const unsigned A_nvec = A.nvector();
242  const unsigned A_nrow_local = A.nrow_local();
243  const unsigned mv_nvec = mv.nvector();
244  // const unsigned mv_nrow_local = mv.nrow_local();
245 
246  for (unsigned v1 = 0; v1 < A_nvec; v1++)
247  {
248  for (unsigned v2 = 0; v2 < mv_nvec; v2++)
249  {
250  B(v1, v2) = 0.0;
251  for (unsigned i = 0; i < A_nrow_local; i++)
252  {
253  B(v1, v2) += A(v1, i) * mv(v2, i);
254  }
255  B(v1, v2) *= alpha;
256  }
257  }
258 
259  // This will need to be communicated if we're in parallel and distributed
260 #ifdef OOMPH_HAS_MPI
261  if (A.distributed() &&
262  A.distribution_pt()->communicator_pt()->nproc() > 1)
263  {
264  const int n_total_val = A_nvec * mv_nvec;
265  // Pack entries into a vector
266  double b[n_total_val];
267  double b2[n_total_val];
268  unsigned count = 0;
269  for (unsigned v1 = 0; v1 < A_nvec; v1++)
270  {
271  for (unsigned v2 = 0; v2 < mv_nvec; v2++)
272  {
273  b[count] = B(v1, v2);
274  b2[count] = b[count];
275  ++count;
276  }
277  }
278 
279 
280  MPI_Allreduce(&b,
281  &b2,
282  n_total_val,
283  MPI_DOUBLE,
284  MPI_SUM,
285  A.distribution_pt()->communicator_pt()->mpi_comm());
286 
287  count = 0;
288  for (unsigned v1 = 0; v1 < A_nvec; v1++)
289  {
290  for (unsigned v2 = 0; v2 < mv_nvec; v2++)
291  {
292  B(v1, v2) = b2[count];
293  ++count;
294  }
295  }
296  }
297 #endif
298  }
299 
300 
301  /*! \brief Compute a vector \c b where the components are the individual
302  * dot-products of the \c i-th columns of \c A and \c mv, i.e.
303  * \f$ b[i] = * A[i]^Hmv[i] \f$.
304  */
305  static void MvDot(const oomph::DoubleMultiVector& mv,
306  const oomph::DoubleMultiVector& A,
307  std::vector<double>& b)
308  {
309  mv.dot(A, b);
310  }
311 
312 
313  /*! \brief Compute the 2-norm of each individual vector of \c mv.
314  Upon return, \c normvec[i] holds the value of \f$||mv_i||_2\f$, the \c
315  i-th column of \c mv.
316  */
317  static void MvNorm(const oomph::DoubleMultiVector& mv,
318  std::vector<double>& normvec)
319  {
320  mv.norm(normvec);
321  }
322 
323  //@}
324 
325  //! @name Initialization methods
326  //@{
327  /*! \brief Copy the vectors in \c A to a set of vectors in \c mv indicated
328  by the indices given in \c index.
329 
330  The \c numvecs vectors in \c A are copied to a subset of vectors in \c mv
331  indicated by the indices given in \c index, i.e.<tt> mv[index[i]] =
332  A[i]</tt>.
333  */
334  static void SetBlock(const oomph::DoubleMultiVector& A,
335  const std::vector<int>& index,
337  {
338  // Check some stuff
339  const unsigned n_index = index.size();
340  if (n_index == 0)
341  {
342  return;
343  }
344  const unsigned n_row_local = mv.nrow_local();
345  for (unsigned v = 0; v < n_index; v++)
346  {
347  for (unsigned i = 0; i < n_row_local; i++)
348  {
349  mv(index[v], i) = A(v, i);
350  }
351  }
352  }
353 
354  /// \brief Deep copy of A into specified columns of mv
355  ///
356  /// (Deeply) copy the first <tt>index.size()</tt> columns of \c A
357  /// into the columns of \c mv specified by the given index range.
358  ///
359  /// Postcondition: <tt>mv[i] = A[i - index.lbound()]</tt>
360  /// for all <tt>i</tt> in <tt>[index.lbound(), index.ubound()]</tt>
361  ///
362  /// \param A [in] Source multivector
363  /// \param index [in] Inclusive index range of columns of mv;
364  /// index set of the target
365  /// \param mv [out] Target multivector
366  static void SetBlock(const oomph::DoubleMultiVector& A,
367  const Teuchos::Range1D& index,
369  {
370  // Check some stuff
371  const unsigned n_index = index.size();
372  if (n_index == 0)
373  {
374  return;
375  }
376  const unsigned n_row_local = mv.nrow_local();
377  unsigned range_index = index.lbound();
378  for (unsigned v = 0; v < n_index; v++)
379  {
380  for (unsigned i = 0; i < n_row_local; i++)
381  {
382  mv(range_index, i) = A(v, i);
383  }
384  ++range_index;
385  }
386  }
387 
388  /// \brief mv := A
389  ///
390  /// Assign (deep copy) A into mv.
391  static void Assign(const oomph::DoubleMultiVector& A,
393  {
394  mv = A;
395  }
396 
397  /*! \brief Replace the vectors in \c mv with random vectors.
398  */
400  {
401  const unsigned n_vector = mv.nvector();
402  const unsigned n_row_local = mv.nrow_local();
403  for (unsigned v = 0; v < n_vector; v++)
404  {
405  for (unsigned i = 0; i < n_row_local; i++)
406  {
407  mv(v, i) = Teuchos::ScalarTraits<double>::random();
408  }
409  }
410  }
411 
412  /*! \brief Replace each element of the vectors in \c mv with \c alpha.
413  */
414  static void MvInit(
416  const double alpha = Teuchos::ScalarTraits<double>::zero())
417  {
418  mv.initialise(alpha);
419  }
420 
421  //@}
422 
423  //! @name Print method
424  //@{
425 
426  /*! \brief Print the \c mv multi-vector to the \c os output stream.
427  */
428  static void MvPrint(const oomph::DoubleMultiVector& mv, std::ostream& os)
429  {
430  mv.output(os);
431  }
432 
433  //@}
434  };
435 
436 
437 } // namespace Anasazi
438 
439 namespace oomph
440 {
441  //===================================================================
442  /// Base class for Oomph-lib's Vector Operator classes that will be
443  /// used with the DoubleMultiVector
444  //==================================================================
446  {
447  public:
448  /// Empty constructor
450 
451  /// virtual destructor
453 
454  /// The apply interface
455  virtual void apply(const DoubleMultiVector& x,
456  DoubleMultiVector& y) const = 0;
457  };
458 
459 } // namespace oomph
460 
461 namespace Anasazi
462 {
463  //======================================================================
464  /// Anasazi Traits class that specialises the apply operator based on
465  /// oomph-lib's DoubleVector class and double as the primitive type
466  //======================================================================
467  template<>
468  class OperatorTraits<double,
469  oomph::DoubleMultiVector,
471  {
472  public:
473  /// Matrix operator application method
475  const oomph::DoubleMultiVector& x,
477  {
478  Op.apply(x, y);
479  }
480 
481  }; // End of the specialised traits class
482 
483 } // namespace Anasazi
484 
485 
486 namespace oomph
487 {
488  //====================================================================
489  /// Class for the shift invert operation
490  //====================================================================
492  {
493  private:
494  // Pointer to the problem
496 
497  // Pointer to the linear solver used
499 
500  // Storage for the shift
501  double Sigma;
502 
503  // Storage for the matrices
505 
506  public:
508  LinearSolver* const& linear_solver_pt,
509  const double& sigma = 0.0)
510  : Problem_pt(problem_pt), Linear_solver_pt(linear_solver_pt), Sigma(sigma)
511  {
512  // Before we get into the Arnoldi loop solve the shifted matrix problem
513  // Allocated Row compressed matrices for the mass matrix and shifted main
514  // matrix
516  {
517  M_pt = new CRDoubleMatrix(problem_pt->dof_distribution_pt());
518  AsigmaM_pt = new CRDoubleMatrix(problem_pt->dof_distribution_pt());
519  }
520  else
521  {
522  // No need for a distribution; gets automatically set up by the Problem
523  M_pt = new CRDoubleMatrix;
525  }
526 
527 
528  // Assemble the matrices
530 
531  // Do not report the time taken
533  }
534 
535 
536  // Now specify how to apply the operator
537  void apply(const DoubleMultiVector& x, DoubleMultiVector& y) const
538  {
539  const unsigned n_vec = x.nvector();
540  const unsigned n_row_local = x.nrow_local();
541  if (n_vec > 1)
542  {
544  }
545 
546  // Solve for the first vector (no need for a distribution)
547  DoubleVector X;
548 
549  // Premultiply by mass matrix
550  M_pt->multiply(x.doublevector(0), X);
551 
552  // Create output vector
553  DoubleVector Y;
555 
556  // Need to synchronise
557  //#ifdef OOMPH_HAS_MPI
558  // Problem_pt->synchronise_all_dofs();
559  //#endif
560 
561  for (unsigned i = 0; i < n_row_local; i++)
562  {
563  y(0, i) = Y[i];
564  }
565 
566  // Now loop over the others and resolve
567  for (unsigned v = 1; v < n_vec; ++v)
568  {
569  M_pt->multiply(x.doublevector(v), X);
570  Linear_solver_pt->resolve(X, Y);
571 
572  //#ifdef OOMPH_HAS_MPI
573  // Problem_pt->synchronise_all_dofs();
574  //#endif
575 
576  for (unsigned i = 0; i < n_row_local; i++)
577  {
578  y(v, i) = Y[i];
579  }
580  }
581  }
582  };
583 
584 
585  //====================================================================
586  /// Class for the adjoing problem shift invert operation
587  //====================================================================
590  {
591  private:
592  /// Pointer to the problem
594 
595  /// Pointer to the linear solver used
597 
598  /// Storage for the shift
599  double Sigma;
600 
601  /// Storage for the matrices
603 
604  public:
606  Problem* const& problem_pt,
607  LinearSolver* const& linear_solver_pt,
608  const double& sigma = 0.0)
609  : Problem_pt(problem_pt), Linear_solver_pt(linear_solver_pt), Sigma(sigma)
610  {
611  /// Before we get into the Arnoldi loop solve the shifted matrix problem
612  /// Allocated Row compressed matrices for the mass matrix and shifted main
613  /// matrix
614  M_pt = new CRDoubleMatrix(problem_pt->dof_distribution_pt());
615  AsigmaM_pt = new CRDoubleMatrix(problem_pt->dof_distribution_pt());
616 
617  /// Assemble the matrices
619 
620  /// Get the transpose of the mass and jacobian
621  /// NB Only difference to ProblemBasedShiftInvertOperator
624 
625  /// Do not report the time taken
627  }
628 
629 
630  /// Now specify how to apply the operator
631  void apply(const DoubleMultiVector& x, DoubleMultiVector& y) const
632  {
633  const unsigned n_vec = x.nvector();
634  const unsigned n_row_local = x.nrow_local();
635  if (n_vec > 1)
636  {
638  }
639 
640  /// Solve the first vector
642 
643  /// Premultiply by mass matrix
644  M_pt->multiply(x.doublevector(0), X);
645 
646  /// Create output vector
649 
650  for (unsigned i = 0; i < n_row_local; i++)
651  {
652  y(0, i) = Y[i];
653  }
654 
655  /// Now loop over the others and resolve
656  for (unsigned v = 1; v < n_vec; ++v)
657  {
658  M_pt->multiply(x.doublevector(v), X);
659  Linear_solver_pt->resolve(X, Y);
660 
661  for (unsigned i = 0; i < n_row_local; i++)
662  {
663  y(v, i) = Y[i];
664  }
665  }
666  }
667  };
668 
669 
670  //=====================================================================
671  /// Class for the Anasazi eigensolver
672  //=====================================================================
673  class ANASAZI : public EigenSolver
674  {
675  private:
676  typedef double ST;
677  typedef Teuchos::ScalarTraits<ST> SCT;
678  typedef SCT::magnitudeType MT;
679  typedef Anasazi::MultiVec<ST> MV;
680  typedef Anasazi::Operator<ST> OP;
681  typedef Anasazi::MultiVecTraits<ST, MV> MVT;
682  typedef Anasazi::OperatorTraits<ST, MV, OP> OPT;
683 
684  /// Pointer to output manager
685  Anasazi::OutputManager<ST>* Output_manager_pt;
686 
687  /// Pointer to a linear solver
689 
690  /// Pointer to a default linear solver
692 
693  // /// Integer to set whether the real, imaginary or magnitude is
694  // /// required
695  // /// to be small or large.
696  // int Spectrum;
697 
698  // /// Number of Arnoldi vectors to compute
699  // int NArnoldi;
700 
701  // /// Boolean to set which part of the spectrum left (default) or right
702  // /// of the shifted value.
703  // bool Small;
704 
705  // /// Boolean to indicate whether or not to compute the eigenvectors
706  // bool Compute_eigenvectors;
707 
708 
709  public:
710  /// Constructor
712  // Spectrum(0),
713  // NArnoldi(10),
714  // Compute_eigenvectors(true)
715  {
716  Output_manager_pt = new Anasazi::BasicOutputManager<ST>();
717  // Set verbosity level
718  int verbosity = 0;
719  verbosity += Anasazi::Warnings;
720  verbosity += Anasazi::FinalSummary;
721  verbosity += Anasazi::TimingDetails;
722  Output_manager_pt->setVerbosity(verbosity);
723 
724  // print greeting
725  Output_manager_pt->stream(Anasazi::Warnings)
726  << Anasazi::Anasazi_Version() << std::endl
727  << std::endl;
728  }
729 
730  /// Empty copy constructor
731  ANASAZI(const ANASAZI&) {}
732 
733  /// Destructor, delete the linear solver
734  virtual ~ANASAZI() {}
735 
736  // /// Access function for the number of Arnoldi vectors
737  // int& narnoldi()
738  // {
739  // return NArnoldi;
740  // }
741 
742  // /// Access function for the number of Arnoldi vectors (const version)
743  // const int& narnoldi() const
744  // {
745  // return NArnoldi;
746  // }
747 
748  // /// Set to enable the computation of the eigenvectors (default)
749  // void enable_compute_eigenvectors()
750  // {
751  // Compute_eigenvectors = true;
752  // }
753 
754  // /// Set to disable the computation of the eigenvectors
755  // void disable_compute_eigenvectors()
756  // {
757  // Compute_eigenvectors = false;
758  // }
759 
760  /// Solve the real eigenproblem that is assembled by elements in
761  /// a mesh in a Problem object. Note that the assembled matrices include the
762  /// shift and are real. The eigenvalues and eigenvectors are,
763  /// in general, complex. Eigenvalues may be infinite and are therefore
764  /// returned as
765  /// \f$ \lambda_i = \alpha_i / \beta_i \f$ where \f$ \alpha_i \f$ is complex
766  /// while \f$ \beta_i \f$ is real. The actual eigenvalues may then be
767  /// computed by doing the division, checking for zero betas to avoid NaNs.
768  /// There's a convenience wrapper to this function that simply computes
769  /// these eigenvalues regardless. That version may die in NaN checking is
770  /// enabled (via the fenv.h header and the associated feenable function).
771  /// NOTE: While the above statement is true, the implementation of this
772  /// function is actually everse engineered -- trilinos actually computes the
773  /// eigenvalues directly (and being an Arnoldi method it wouldn't be able to
774  /// obtain any infinite/NaN eigenvalues anyway
775  void solve_eigenproblem(Problem* const& problem_pt,
776  const int& n_eval,
777  Vector<std::complex<double>>& alpha,
778  Vector<double>& beta,
779  Vector<DoubleVector>& eigenvector_real,
780  Vector<DoubleVector>& eigenvector_imag,
781  const bool& do_adjoint_problem)
782  {
783  // Reverse engineer the "safe" version of the eigenvalues
784  Vector<std::complex<double>> eigenvalue;
785  solve_eigenproblem(problem_pt,
786  n_eval,
787  eigenvalue,
788  eigenvector_real,
789  eigenvector_imag,
790  do_adjoint_problem);
791  unsigned n = eigenvalue.size();
792  alpha.resize(n);
793  beta.resize(n);
794  for (unsigned i = 0; i < n; i++)
795  {
796  alpha[i] = eigenvalue[i];
797  beta[i] = 1.0;
798  }
799  }
800 
801  /// Solve the eigen problem
802  void solve_eigenproblem(Problem* const& problem_pt,
803  const int& n_eval,
804  Vector<std::complex<double>>& eigenvalue,
805  Vector<DoubleVector>& eigenvector_real,
806  Vector<DoubleVector>& eigenvector_imag,
807  const bool& do_adjoint_problem)
808  {
809  // Initially be dumb here
810  Linear_solver_pt = problem_pt->linear_solver_pt();
811 
812  // Get a shiny new linear algebra distribution from the problem
813  LinearAlgebraDistribution* dist_pt = 0;
814  problem_pt->create_new_linear_algebra_distribution(dist_pt);
815 
816  // Let's make the initial vector
817  Teuchos::RCP<DoubleMultiVector> initial =
818  Teuchos::rcp(new DoubleMultiVector(1, dist_pt));
819  Anasazi::MultiVecTraits<double, DoubleMultiVector>::MvRandom(*initial);
820 
821  // Make the operator
822  Teuchos::RCP<DoubleMultiVectorOperator> Op_pt;
823  if (do_adjoint_problem)
824  {
825  Op_pt = Teuchos::rcp(new AdjointProblemBasedShiftInvertOperator(
826  problem_pt, this->linear_solver_pt(), Sigma_real));
827  }
828  else
829  {
830  Op_pt = Teuchos::rcp(new ProblemBasedShiftInvertOperator(
831  problem_pt, this->linear_solver_pt(), Sigma_real));
832  }
833 
834  // Create the basic problem
835  Teuchos::RCP<Anasazi::BasicEigenproblem<double,
838  anasazi_pt = Teuchos::rcp(
839  new Anasazi::BasicEigenproblem<double,
842  initial));
843 
844  // No assumptions about niceness of problem matrices
845  anasazi_pt->setHermitian(false);
846 
847  // set the number of eigenvalues requested
848  anasazi_pt->setNEV(n_eval);
849 
850  // Inform the eigenproblem that you are done passing it information
851  if (anasazi_pt->setProblem() != true)
852  {
853  oomph_info << "Anasazi::BasicEigenproblem::setProblem() had an error."
854  << std::endl
855  << "End Result: TEST FAILED" << std::endl;
856  }
857 
858  // Create the solver manager
859  // No need to have ncv specificed, Triliinos has a sensible default
860  // int ncv = 10;
861  MT tol = 1.0e-10;
862  int verbosity = 0;
863  verbosity += Anasazi::Warnings;
864  // MH has switched off the (overly) verbose output
865  // Could introduce handle to switch it back on.
866  // verbosity+=Anasazi::FinalSummary;
867  // verbosity+=Anasazi::TimingDetails;
868 
869 
870  Teuchos::ParameterList MyPL;
871  MyPL.set("Which", "LM");
872  MyPL.set("Block Size", 1);
873  // MyPL.set( "Num Blocks", ncv );
874  MyPL.set("Maximum Restarts", 500);
875  MyPL.set("Orthogonalization", "DGKS");
876  MyPL.set("Verbosity", verbosity);
877  MyPL.set("Convergence Tolerance", tol);
878  Anasazi::BlockKrylovSchurSolMgr<double,
881  BKS(anasazi_pt, MyPL);
882 
883  // Solve the problem to the specified tolerances or length
884  Anasazi::ReturnType ret = BKS.solve();
885  if (ret != Anasazi::Converged)
886  {
887  std::string error_message = "Eigensolver did not converge!\n";
888 
889  throw OomphLibError(
890  error_message, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
891  }
892 
893 
894  // Get the eigenvalues and eigenvectors from the eigenproblem
895  Anasazi::Eigensolution<double, DoubleMultiVector> sol =
896  anasazi_pt->getSolution();
897  std::vector<Anasazi::Value<double>> evals = sol.Evals;
898  Teuchos::RCP<DoubleMultiVector> evecs = sol.Evecs;
899 
900  // Here's the vector that contains the information
901  // about how to translate the vectors that are returned
902  // by anasazi into real and imag parts of the actual
903  // eigenvectors
904  // std::vector<int> Anasazi::Eigensolution< ScalarType, MV >::index
905  std::vector<int> index_vector = sol.index;
906 
907  // Here's what we're actually going to return.
908  unsigned n_eval_avail = evals.size();
909  eigenvalue.resize(n_eval_avail);
910  eigenvector_real.resize(n_eval_avail);
911  eigenvector_imag.resize(n_eval_avail);
912 
913  // Extract these from the raw data returned by trilinos
914  for (unsigned i = 0; i < n_eval_avail; i++)
915  {
916  // Undo shift and invert
917  double a = evals[i].realpart;
918  double b = evals[i].imagpart;
919  double det = a * a + b * b;
920  eigenvalue[i] = std::complex<double>(a / det + Sigma_real, -b / det);
921 
922  // Now set the eigenvectors
923  unsigned nrow_local = evecs->nrow_local();
924  eigenvector_real[i].build(evecs->distribution_pt(), 0.0);
925  eigenvector_imag[i].build(evecs->distribution_pt(), 0.0);
926 
927  // std::vector<int> Anasazi::Eigensolution< ScalarType, MV >::index
928  // to translate into proper complex vector; see
929  // https://docs.trilinos.org/dev/packages/anasazi/doc/html/structAnasazi_1_1Eigensolution.html#ac9d141d98adcba85fbad011a7b7bda6e
930 
931  // An index into Evecs to allow compressed storage of eigenvectors for
932  // real, non-Hermitian problems.
933 
934  // index has length numVecs, where each entry is 0, +1, or -1. These
935  // have the following interpretation:
936 
937  // index[i] == 0: signifies that the corresponding eigenvector is
938  // stored as the i column of Evecs. This will usually be the case
939  // when ScalarType is complex, an eigenproblem is Hermitian, or a
940  // real, non-Hermitian eigenproblem has a real eigenvector. index[i]
941  // == +1: signifies that the corresponding eigenvector is stored in
942  // two vectors: the real part in the i column of Evecs and the
943  // positive imaginary part in the i+1 column of Evecs. index[i] ==
944  // -1: signifies that the corresponding eigenvector is stored in two
945  // vectors: the real part in the i-1 column of Evecs and the
946  // negative imaginary part in the i column of Evecs
947 
948  // Real eigenvector
949  if (index_vector[i] == 0)
950  {
951  for (unsigned j = 0; j < nrow_local; j++)
952  {
953  eigenvector_real[i][j] = (*evecs)(i, j);
954  eigenvector_imag[i][j] = 0.0;
955  }
956  }
957  else if (index_vector[i] == 1)
958  {
959  for (unsigned j = 0; j < nrow_local; j++)
960  {
961  eigenvector_real[i][j] = (*evecs)(i, j);
962  eigenvector_imag[i][j] = (*evecs)(i + 1, j);
963  }
964  }
965  else if (index_vector[i] == -1)
966  {
967  for (unsigned j = 0; j < nrow_local; j++)
968  {
969  eigenvector_real[i][j] = (*evecs)(i - 1, j);
970  eigenvector_imag[i][j] = -(*evecs)(i, j);
971  }
972  }
973  else
974  {
975  oomph_info << "Never get here: index_vector[ " << i
976  << " ] = " << index_vector[i] << std::endl;
977  abort();
978  }
979  }
980 
981  // Tidy up
982  delete dist_pt;
983  }
984 
985 
986  /// Solve the eigen problem
987  void solve_eigenproblem_legacy(Problem* const& problem_pt,
988  const int& n_eval,
989  Vector<std::complex<double>>& eigenvalue,
990  Vector<DoubleVector>& eigenvector,
991  const bool& do_adjoint_problem)
992  {
993  // Initially be dumb here
994  Linear_solver_pt = problem_pt->linear_solver_pt();
995 
996  // Let's make the initial one-dimensional vector
997  unsigned multivector_dimension = 1;
998  Teuchos::RCP<DoubleMultiVector> initial =
999  Teuchos::rcp(new DoubleMultiVector(multivector_dimension,
1000  problem_pt->dof_distribution_pt()));
1001  Anasazi::MultiVecTraits<double, DoubleMultiVector>::MvRandom(*initial);
1002 
1003  // Make the operator
1004  Teuchos::RCP<DoubleMultiVectorOperator> Op_pt;
1005  if (do_adjoint_problem)
1006  {
1007  Op_pt = Teuchos::rcp(new AdjointProblemBasedShiftInvertOperator(
1008  problem_pt, this->linear_solver_pt(), Sigma_real));
1009  }
1010  else
1011  {
1012  Op_pt = Teuchos::rcp(new ProblemBasedShiftInvertOperator(
1013  problem_pt, this->linear_solver_pt(), Sigma_real));
1014  }
1015 
1016  // Create the basic problem
1017  Teuchos::RCP<Anasazi::BasicEigenproblem<double,
1020  anasazi_pt = Teuchos::rcp(
1021  new Anasazi::BasicEigenproblem<double,
1024  initial));
1025 
1026  // The problem is not Hermitian in general
1027  anasazi_pt->setHermitian(false);
1028 
1029  // set the number of eigenvalues requested
1030  anasazi_pt->setNEV(n_eval);
1031 
1032  // Inform the eigenproblem that you are done passing it information
1033  if (anasazi_pt->setProblem() != true)
1034  {
1035  oomph_info << "Anasazi::BasicEigenproblem::setProblem() had an error."
1036  << std::endl
1037  << "End Result: TEST FAILED" << std::endl;
1038  }
1039 
1040  // Create the solver manager
1041  // No need to have ncv specified, Trilinos has a sensible default
1042  // int ncv = 10;
1043  MT tol = 1.0e-10;
1044  int verbosity =
1045  Anasazi::Warnings + Anasazi::FinalSummary + Anasazi::TimingDetails;
1046 
1047  Teuchos::ParameterList MyPL;
1048  MyPL.set("Which", "LM");
1049  MyPL.set("Block Size", 1);
1050  // MyPL.set( "Num Blocks", ncv );
1051  MyPL.set("Maximum Restarts", 500);
1052  MyPL.set("Orthogonalization", "DGKS");
1053  MyPL.set("Verbosity", verbosity);
1054  MyPL.set("Convergence Tolerance", tol);
1055  Anasazi::BlockKrylovSchurSolMgr<double,
1058  BKS(anasazi_pt, MyPL);
1059 
1060  // Solve the problem to the specified tolerances or length
1061  Anasazi::ReturnType ret = BKS.solve();
1062 
1063  if (ret != Anasazi::Converged)
1064  {
1065  oomph_info << "Eigensolver not converged\n";
1066  }
1067 
1068  // Get the eigenvalues and eigenvectors from the eigenproblem
1069  Anasazi::Eigensolution<double, DoubleMultiVector> sol =
1070  anasazi_pt->getSolution();
1071  std::vector<Anasazi::Value<double>> evals = sol.Evals;
1072  Teuchos::RCP<DoubleMultiVector> evecs = sol.Evecs;
1073 
1074  eigenvalue.resize(evals.size());
1075  eigenvector.resize(evals.size());
1076 
1077  for (unsigned i = 0; i < evals.size(); i++)
1078  {
1079  // Undo shift and invert
1080  double a = evals[i].realpart;
1081  double b = evals[i].imagpart;
1082  double det = a * a + b * b;
1083  eigenvalue[i] = std::complex<double>(a / det + Sigma_real, -b / det);
1084 
1085  // Now set the eigenvectors
1086  eigenvector[i].build(evecs->distribution_pt());
1087  unsigned nrow_local = evecs->nrow_local();
1088 
1089  // This is the legacy stuff. Use
1090  // std::vector<int> Anasazi::Eigensolution< ScalarType, MV >::index
1091  // to translate into proper complex vector; see
1092  // https://docs.trilinos.org/dev/packages/anasazi/doc/html/structAnasazi_1_1Eigensolution.html#ac9d141d98adcba85fbad011a7b7bda6e
1093  for (unsigned n = 0; n < nrow_local; n++)
1094  {
1095  eigenvector[i][n] = (*evecs)(i, n);
1096  }
1097  }
1098  }
1099 
1100  // /// Set the desired eigenvalues to be left of the shift
1101  // void get_eigenvalues_left_of_shift()
1102  // {
1103  // Small = true;
1104  // }
1105 
1106  // /// Set the desired eigenvalues to be right of the shift
1107  // void get_eigenvalues_right_of_shift()
1108  // {
1109  // Small = false;
1110  // }
1111 
1112  // /// Set the real part to be the quantity of interest (default)
1113  // void track_eigenvalue_real_part()
1114  // {
1115  // Spectrum = 1;
1116  // }
1117 
1118  // /// Set the imaginary part fo the quantity of interest
1119  // void track_eigenvalue_imaginary_part()
1120  // {
1121  // Spectrum = -1;
1122  // }
1123 
1124  // /// Set the magnitude to be the quantity of interest
1125  // void track_eigenvalue_magnitude()
1126  // {
1127  // Spectrum = 0;
1128  // }
1129 
1130  /// Return a pointer to the linear solver object
1132  {
1133  return Linear_solver_pt;
1134  }
1135 
1136  /// Return a pointer to the linear solver object (const version)
1138  {
1139  return Linear_solver_pt;
1140  }
1141  };
1142 
1143 } // namespace oomph
1144 
1145 
1146 #endif
cstr elem_len * i
Definition: cfortran.h:603
static Teuchos::RCP< oomph::DoubleMultiVector > CloneCopy(const oomph::DoubleMultiVector &mv)
Creates a deep copy of the DoubleMultiVector mv return Reference-counted pointer to the new DoubleMul...
static int GetVecLength(const oomph::DoubleMultiVector &mv)
Obtain the global length of the vector.
static Teuchos::RCP< oomph::DoubleMultiVector > CloneCopy(const oomph::DoubleMultiVector &mv, const std::vector< int > &index)
Creates a new oomph::DoubleMultiVector and (deep) copies the contents of the vectors in index into th...
static void MvTransMv(const double alpha, const oomph::DoubleMultiVector &A, const oomph::DoubleMultiVector &mv, Teuchos::SerialDenseMatrix< int, double > &B)
Compute a dense matrix B through the matrix-matrix multiply .
static void MvTimesMatAddMv(const double alpha, const oomph::DoubleMultiVector &A, const Teuchos::SerialDenseMatrix< int, double > &B, const double beta, oomph::DoubleMultiVector &mv)
Update mv with .
static void MvPrint(const oomph::DoubleMultiVector &mv, std::ostream &os)
Print the mv multi-vector to the os output stream.
static Teuchos::RCP< oomph::DoubleMultiVector > CloneView(oomph::DoubleMultiVector &mv, const Teuchos::Range1D &index)
Creates a new oomph::DoubleMultiVector that contains shallow copies of selected entries of the oomph:...
static void SetBlock(const oomph::DoubleMultiVector &A, const Teuchos::Range1D &index, oomph::DoubleMultiVector &mv)
Deep copy of A into specified columns of mv.
static Teuchos::RCP< oomph::DoubleMultiVector > CloneViewNonConst(oomph::DoubleMultiVector &mv, const std::vector< int > &index)
Creates a new oomph::DoubleMultiVector that contains shallow copies of selected entries of the oomph:...
static void MvScale(oomph::DoubleMultiVector &mv, const std::vector< double > &alpha)
Scale each element of the i-th vector in mv with alpha[i].
static void MvScale(oomph::DoubleMultiVector &mv, const double alpha)
Scale each element of the vectors in mv with alpha.
static Teuchos::RCP< const oomph::DoubleMultiVector > CloneView(const oomph::DoubleMultiVector &mv, const std::vector< int > &index)
Creates a new oomph::DoubleMultiVector that contains shallow copies of selected entries of the oomph:...
static void MvAddMv(const double alpha, const oomph::DoubleMultiVector &A, const double beta, const oomph::DoubleMultiVector &B, oomph::DoubleMultiVector &mv)
Replace mv with .
static Teuchos::RCP< oomph::DoubleMultiVector > Clone(const oomph::DoubleMultiVector &mv, const int numvecs)
Creates a new empty DoubleMultiVector containing numvecs columns. Return a reference-counted pointer ...
static int GetNumberVecs(const oomph::DoubleMultiVector &mv)
Obtain the number of vectors in the multivector.
static void SetBlock(const oomph::DoubleMultiVector &A, const std::vector< int > &index, oomph::DoubleMultiVector &mv)
Copy the vectors in A to a set of vectors in mv indicated by the indices given in index.
static Teuchos::RCP< oomph::DoubleMultiVector > CloneViewNonConst(oomph::DoubleMultiVector &mv, const Teuchos::Range1D &index)
Creates a new oomph::DoubleMultiVector that contains shallow copies of selected entries of the oomph:...
static void MvNorm(const oomph::DoubleMultiVector &mv, std::vector< double > &normvec)
Compute the 2-norm of each individual vector of mv. Upon return, normvec[i] holds the value of ,...
static Teuchos::RCP< oomph::DoubleMultiVector > CloneView(oomph::DoubleMultiVector &mv, const std::vector< int > &index)
Creates a new oomph::DoubleMultiVector that contains shallow copies of selected entries of the oomph:...
static void MvDot(const oomph::DoubleMultiVector &mv, const oomph::DoubleMultiVector &A, std::vector< double > &b)
Compute a vector b where the components are the individual dot-products of the i-th columns of A and ...
static void MvRandom(oomph::DoubleMultiVector &mv)
Replace the vectors in mv with random vectors.
static Teuchos::RCP< oomph::DoubleMultiVector > CloneCopy(const oomph::DoubleMultiVector &mv, const Teuchos::Range1D &index)
Deep copy of specified columns of oomph::DoubleMultiVector return Reference-counted pointer to the ne...
static void MvInit(oomph::DoubleMultiVector &mv, const double alpha=Teuchos::ScalarTraits< double >::zero())
Replace each element of the vectors in mv with alpha.
static void Assign(const oomph::DoubleMultiVector &A, oomph::DoubleMultiVector &mv)
mv := A
static void Apply(const oomph::DoubleMultiVectorOperator &Op, const oomph::DoubleMultiVector &x, oomph::DoubleMultiVector &y)
Matrix operator application method.
Class for the Anasazi eigensolver.
void solve_eigenproblem(Problem *const &problem_pt, const int &n_eval, Vector< std::complex< double >> &alpha, Vector< double > &beta, Vector< DoubleVector > &eigenvector_real, Vector< DoubleVector > &eigenvector_imag, const bool &do_adjoint_problem)
Solve the real eigenproblem that is assembled by elements in a mesh in a Problem object....
LinearSolver * Default_linear_solver_pt
Pointer to a default linear solver.
Anasazi::OutputManager< ST > * Output_manager_pt
Pointer to output manager.
Teuchos::ScalarTraits< ST > SCT
Anasazi::OperatorTraits< ST, MV, OP > OPT
virtual ~ANASAZI()
Destructor, delete the linear solver.
LinearSolver * Linear_solver_pt
Pointer to a linear solver.
void solve_eigenproblem(Problem *const &problem_pt, const int &n_eval, Vector< std::complex< double >> &eigenvalue, Vector< DoubleVector > &eigenvector_real, Vector< DoubleVector > &eigenvector_imag, const bool &do_adjoint_problem)
Solve the eigen problem.
LinearSolver *const & linear_solver_pt() const
Return a pointer to the linear solver object (const version)
SCT::magnitudeType MT
void solve_eigenproblem_legacy(Problem *const &problem_pt, const int &n_eval, Vector< std::complex< double >> &eigenvalue, Vector< DoubleVector > &eigenvector, const bool &do_adjoint_problem)
Solve the eigen problem.
Anasazi::MultiVec< ST > MV
Anasazi::MultiVecTraits< ST, MV > MVT
ANASAZI(const ANASAZI &)
Empty copy constructor.
Anasazi::Operator< ST > OP
LinearSolver *& linear_solver_pt()
Return a pointer to the linear solver object.
Class for the adjoing problem shift invert operation.
AdjointProblemBasedShiftInvertOperator(Problem *const &problem_pt, LinearSolver *const &linear_solver_pt, const double &sigma=0.0)
void apply(const DoubleMultiVector &x, DoubleMultiVector &y) const
Now specify how to apply the operator.
LinearSolver * Linear_solver_pt
Pointer to the linear solver used.
CRDoubleMatrix * M_pt
Storage for the matrices.
A class for compressed row matrices. This is a distributable object.
Definition: matrices.h:888
void multiply(const DoubleVector &x, DoubleVector &soln) const
Multiply the matrix by the vector x: soln=Ax.
Definition: matrices.cc:1782
void get_matrix_transpose(CRDoubleMatrix *result) const
Returns the transpose of this matrix.
Definition: matrices.cc:3271
bool distributed() const
distribution is serial or distributed
LinearAlgebraDistribution * distribution_pt() const
access to the LinearAlgebraDistribution
unsigned nrow() const
access function to the number of global rows.
unsigned nrow_local() const
access function for the num of local rows on this processor.
Base class for Oomph-lib's Vector Operator classes that will be used with the DoubleMultiVector.
virtual ~DoubleMultiVectorOperator()
virtual destructor
virtual void apply(const DoubleMultiVector &x, DoubleMultiVector &y) const =0
The apply interface.
A multi vector in the mathematical sense, initially developed for linear algebra type applications....
void output(std::ostream &outfile) const
output the contents of the vector
void initialise(const double &initial_value)
initialise the whole vector with value v
void dot(const DoubleMultiVector &vec, std::vector< double > &result) const
compute the 2 norm of this vector
void norm(std::vector< double > &result) const
compute the 2 norm of this vector
DoubleVector & doublevector(const unsigned &i)
access to the DoubleVector representatoin
unsigned nvector() const
Return the number of vectors.
A vector in the mathematical sense, initially developed for linear algebra type applications....
Definition: double_vector.h:58
Base class for all EigenProblem solves. This simply defines standard interfaces so that different sol...
Definition: eigen_solver.h:61
double Sigma_real
Double value that represents the real part of the shift in shifted eigensolvers.
Definition: eigen_solver.h:65
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
Base class for all linear solvers. This merely defines standard interfaces for linear solvers,...
Definition: linear_solver.h:68
virtual void solve(Problem *const &problem_pt, DoubleVector &result)=0
Solver: Takes pointer to problem and returns the results vector which contains the solution of the li...
virtual void enable_resolve()
Enable resolve (i.e. store matrix and/or LU decomposition, say) Virtual so it can be overloaded to pe...
virtual void resolve(const DoubleVector &rhs, DoubleVector &result)
Resolve the system defined by the last assembled jacobian and the rhs vector. Solution is returned in...
void disable_doc_time()
Disable documentation of solve times.
An OomphLibError object which should be thrown when an run-time error is encountered....
Class for the shift invert operation.
void apply(const DoubleMultiVector &x, DoubleMultiVector &y) const
The apply interface.
ProblemBasedShiftInvertOperator(Problem *const &problem_pt, LinearSolver *const &linear_solver_pt, const double &sigma=0.0)
////////////////////////////////////////////////////////////////// //////////////////////////////////...
Definition: problem.h:151
LinearAlgebraDistribution *const & dof_distribution_pt() const
Return the pointer to the dof distribution (read-only)
Definition: problem.h:1672
virtual void get_eigenproblem_matrices(CRDoubleMatrix &mass_matrix, CRDoubleMatrix &main_matrix, const double &shift=0.0)
Get the matrices required by a eigensolver. If the shift parameter is non-zero the second matrix will...
Definition: problem.cc:8602
LinearSolver *& linear_solver_pt()
Return a pointer to the linear solver object.
Definition: problem.h:1466
void create_new_linear_algebra_distribution(LinearAlgebraDistribution *&dist_pt)
Get new linear algebra distribution (you're in charge of deleting it!)
Definition: problem.cc:299
A slight extension to the standard template vector class so that we can include "graceful" array rang...
Definition: Vector.h:58
bool Use_temporary_code_for_andrew_legacy_version
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...