linear_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-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 // The actual solve functions for dense LU linear solvers.
27 
28 // Config header generated by autoconfig
29 #ifdef HAVE_CONFIG_H
30 #include <oomph-lib-config.h>
31 #endif
32 
33 #ifdef OOMPH_HAS_MPI
34 #include "mpi.h"
35 #endif
36 
37 // oomph-lib includes
38 #include "Vector.h"
39 #include "linear_solver.h"
40 #include "matrices.h"
41 #include "problem.h"
42 
43 
44 namespace oomph
45 {
46  //=============================================================================
47  /// Solver: Takes pointer to problem and returns the results Vector
48  /// which contains the solution of the linear system defined by
49  /// the problem's fully assembled Jacobian and residual Vector.
50  //=============================================================================
51  void DenseLU::solve(Problem* const& problem_pt, DoubleVector& result)
52  {
53  // Initialise timer
54  double t_start = TimingHelpers::timer();
55 
56  // Find # of degrees of freedom (variables)
57  const unsigned n_dof = problem_pt->ndof();
58 
59  // Allocate storage for the residuals vector and the jacobian matrix
60  DoubleVector residuals;
61  DenseDoubleMatrix jacobian(n_dof);
62 
63  // initialise timer
64  double t_start_jacobian = TimingHelpers::timer();
65 
66  // Get the full jacobian and residuals of the problem
67  problem_pt->get_jacobian(residuals, jacobian);
68 
69  // compute jacobian setup time
70  double t_end_jacobian = TimingHelpers::timer();
71  Jacobian_setup_time = t_end_jacobian - t_start_jacobian;
72 
73  // Report the time
74  if (Doc_time)
75  {
76  oomph_info << std::endl
77  << "CPU for setup of Dense Jacobian: "
80  << std::endl;
81  }
82 
83  // Solve by dense LU decomposition VERY INEFFICIENT!
84  solve(&jacobian, residuals, result);
85 
86  // Set the sign of the determinant of the jacobian
88 
89  // Finalise/doc timings
90  double t_end = TimingHelpers::timer();
91  double total_time = t_end - t_start;
92  if (Doc_time)
93  {
94  oomph_info << "CPU for DenseLU LinearSolver: "
96  << std::endl
97  << std::endl;
98  }
99  }
100 
101 
102  //=============================================================================
103  /// Delete the storage that has been allocated for the LU factors, if
104  /// the matrix data is not itself being overwritten.
105  //=============================================================================
107  {
108  // delete the Distribution_pt
109  this->clear_distribution();
110 
111  // Clean up the LU factor storage, if it has been allocated
112  // N.B. we don't need to check the index storage as well.
113  if (LU_factors != 0)
114  {
115  // Delete the pointer to the LU factors
116  delete[] LU_factors;
117  // Null out the vector
118  LU_factors = 0;
119  // Delete the pointer to the Index
120  delete[] Index;
121  // Null out
122  Index = 0;
123  }
124  }
125 
126  //=============================================================================
127  /// LU decompose the matrix.
128  /// WARNING: this class does not perform any PARANOID checks on the vectors -
129  /// these are all performed in the solve(...) method.
130  //=============================================================================
131  void DenseLU::factorise(DoubleMatrixBase* const& matrix_pt)
132  {
133  // Set the number of unknowns
134  const unsigned long n = matrix_pt->nrow();
135 
136  // Small constant
137  const double small_number = 1.0e-20;
138 
139  // Vector scaling stores the implicit scaling of each row
140  Vector<double> scaling(n);
141 
142  // Integer to store the sign that must multiply the determinant as
143  // a consequence of the row/column interchanges
144  int signature = 1;
145 
146  // Loop over rows to get implicit scaling information
147  for (unsigned long i = 0; i < n; i++)
148  {
149  double largest_entry = 0.0;
150  for (unsigned long j = 0; j < n; j++)
151  {
152  double tmp = std::fabs((*matrix_pt)(i, j));
153  if (tmp > largest_entry) largest_entry = tmp;
154  }
155  if (largest_entry == 0.0)
156  {
157  throw OomphLibError(
158  "Singular Matrix", OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
159  }
160  // Save the scaling
161  scaling[i] = 1.0 / largest_entry;
162  }
163 
164  // Firsly, we shall delete any previous LU storage.
165  // If the user calls this function twice without changing the matrix
166  // then it is their own inefficiency, not ours (this time).
167  clean_up_memory();
168 
169  // Allocate storage for the LU factors, the index and store
170  // the number of unknowns
171  LU_factors = new double[n * n];
172  Index = new long[n];
173 
174  // Now we know that memory has been allocated, copy over
175  // the matrix values
176  unsigned count = 0;
177  for (unsigned long i = 0; i < n; i++)
178  {
179  for (unsigned long j = 0; j < n; j++)
180  {
181  LU_factors[count] = (*matrix_pt)(i, j);
182  ++count;
183  }
184  }
185 
186  // Loop over columns
187  for (unsigned long j = 0; j < n; j++)
188  {
189  // Initialise imax
190  unsigned long imax = 0;
191 
192  for (unsigned long i = 0; i < j; i++)
193  {
194  double sum = LU_factors[n * i + j];
195  for (unsigned long k = 0; k < i; k++)
196  {
197  sum -= LU_factors[n * i + k] * LU_factors[n * k + j];
198  }
199  LU_factors[n * i + j] = sum;
200  }
201 
202  // Initialise search for largest pivot element
203  double largest_entry = 0.0;
204  for (unsigned long i = j; i < n; i++)
205  {
206  double sum = LU_factors[n * i + j];
207  for (unsigned long k = 0; k < j; k++)
208  {
209  sum -= LU_factors[n * i + k] * LU_factors[n * k + j];
210  }
211  LU_factors[n * i + j] = sum;
212  // Set temporary
213  double tmp = scaling[i] * std::fabs(sum);
214  if (tmp >= largest_entry)
215  {
216  largest_entry = tmp;
217  imax = i;
218  }
219  }
220 
221  // Test to see if we need to interchange rows
222  if (j != imax)
223  {
224  for (unsigned long k = 0; k < n; k++)
225  {
226  double tmp = LU_factors[n * imax + k];
227  LU_factors[n * imax + k] = LU_factors[n * j + k];
228  LU_factors[n * j + k] = tmp;
229  }
230  // Change the parity of signature
231  signature = -signature;
232 
233  // Interchange scale factor
234  scaling[imax] = scaling[j];
235  }
236 
237  // Set the index
238  Index[j] = imax;
239  if (LU_factors[n * j + j] == 0.0)
240  {
241  LU_factors[n * j + j] = small_number;
242  }
243  // Divide by pivot element
244  if (j != n - 1)
245  {
246  double tmp = 1.0 / LU_factors[n * j + j];
247  for (unsigned long i = j + 1; i < n; i++)
248  {
249  LU_factors[n * i + j] *= tmp;
250  }
251  }
252 
253  } // End of loop over columns
254 
255 
256  // Now multiply all the diagonal terms together to get the determinant
257  // Note that we need to use the mantissa, exponent formulation to
258  // avoid underflow errors
259  double determinant_mantissa = 1.0;
260  int determinant_exponent = 0, iexp;
261  for (unsigned i = 0; i < n; i++)
262  {
263  // Multiply by the next diagonal entry's mantissa
264  // and return the exponent
265  determinant_mantissa *= frexp(LU_factors[n * i + i], &iexp);
266 
267  // Add the new exponent to the current exponent
268  determinant_exponent += iexp;
269 
270  // normalise
271  determinant_mantissa = frexp(determinant_mantissa, &iexp);
272  determinant_exponent += iexp;
273  }
274 
275  // If paranoid issue a warning that the matrix is near singular
276  // #ifdef PARANOID
277  // int tiny_exponent = -60;
278  // if(determinant_exponent < tiny_exponent)
279  // {
280  // std::ostringstream warning_stream;
281  // warning_stream << "The determinant of the matrix is very close to
282  // zero.\n"
283  // << "It is " << determinant_mantissa << " x 2^"
284  // << determinant_exponent << "\n";
285  // warning_stream << "The results will depend on the exact details of
286  // the\n"
287  // << "floating point implementation ... just to let you
288  // know\n";
289  // OomphLibWarning(warning_stream.str(),
290  // "DenseLU::factorise()",
291  // OOMPH_EXCEPTION_LOCATION);
292  // }
293  // #endif
294 
295  // Integer to store the sign of the determinant
296  int sign = 0;
297 
298  // Find the sign of the determinant
299  if (determinant_mantissa > 0.0)
300  {
301  sign = 1;
302  }
303  if (determinant_mantissa < 0.0)
304  {
305  sign = -1;
306  }
307 
308  // Multiply the sign by the signature
309  sign *= signature;
310 
311  // Return the sign of the determinant
313  }
314 
315  //=============================================================================
316  /// Do the backsubstitution for the DenseLU solver.
317  /// WARNING: this class does not perform any PARANOID checks on the vectors -
318  /// these are all performed in the solve(...) method.
319  //=============================================================================
320  void DenseLU::backsub(const DoubleVector& rhs, DoubleVector& result)
321  {
322  // Get pointers to first entries
323  const double* rhs_pt = rhs.values_pt();
324  double* result_pt = result.values_pt();
325 
326  // Copy the rhs vector into the result vector
327  const unsigned long n = rhs.nrow();
328  for (unsigned long i = 0; i < n; ++i)
329  {
330  result_pt[i] = rhs_pt[i];
331  }
332 
333  // Loop over all rows for forward substition
334  unsigned long k = 0;
335  for (unsigned long i = 0; i < n; i++)
336  {
337  unsigned long ip = Index[i];
338  double sum = result_pt[ip];
339  result_pt[ip] = result_pt[i];
340  if (k != 0)
341  {
342  for (unsigned long j = k - 1; j < i; j++)
343  {
344  sum -= LU_factors[n * i + j] * result_pt[j];
345  }
346  }
347  else if (sum != 0.0)
348  {
349  k = i + 1;
350  }
351  result_pt[i] = sum;
352  }
353 
354  // Now do the back substitution
355  for (long i = long(n) - 1; i >= 0; i--)
356  {
357  double sum = result_pt[i];
358  for (long j = i + 1; j < long(n); j++)
359  {
360  sum -= LU_factors[n * i + j] * result_pt[j];
361  }
362  result_pt[i] = sum / LU_factors[n * i + i];
363  }
364  }
365 
366  //=============================================================================
367  /// Do the backsubstitution for the DenseLU solver.
368  /// WARNING: this class does not perform any PARANOID checks on the vectors -
369  /// these are all performed in the solve(...) method. So, if you call backsub
370  /// directly, you have been warned...
371  //=============================================================================
373  {
374  // Copy the rhs vector into the result vector
375  const unsigned long n = rhs.size();
376  for (unsigned long i = 0; i < n; ++i)
377  {
378  result[i] = rhs[i];
379  }
380 
381  // Loop over all rows for forward substition
382  unsigned long k = 0;
383  for (unsigned long i = 0; i < n; i++)
384  {
385  unsigned long ip = Index[i];
386  double sum = result[ip];
387  result[ip] = result[i];
388  if (k != 0)
389  {
390  for (unsigned long j = k - 1; j < i; j++)
391  {
392  sum -= LU_factors[n * i + j] * result[j];
393  }
394  }
395  else if (sum != 0.0)
396  {
397  k = i + 1;
398  }
399  result[i] = sum;
400  }
401 
402  // Now do the back substitution
403  for (long i = long(n) - 1; i >= 0; i--)
404  {
405  double sum = result[i];
406  for (long j = i + 1; j < long(n); j++)
407  {
408  sum -= LU_factors[n * i + j] * result[j];
409  }
410  result[i] = sum / LU_factors[n * i + i];
411  }
412  }
413 
414 
415  //=============================================================================
416  /// Linear-algebra-type solver: Takes pointer to a matrix and rhs
417  /// vector and returns the solution of the linear system.
418  //============================================================================
419  void DenseLU::solve(DoubleMatrixBase* const& matrix_pt,
420  const DoubleVector& rhs,
421  DoubleVector& result)
422  {
423 #ifdef PARANOID
424  // check that the rhs vector is not distributed
425  if (rhs.distribution_pt()->distributed())
426  {
427  std::ostringstream error_message_stream;
428  error_message_stream
429  << "The vectors rhs and result must not be distributed";
430  throw OomphLibError(error_message_stream.str(),
431  OOMPH_CURRENT_FUNCTION,
432  OOMPH_EXCEPTION_LOCATION);
433  }
434 
435  // check that the matrix is square
436  if (matrix_pt->nrow() != matrix_pt->ncol())
437  {
438  std::ostringstream error_message_stream;
439  error_message_stream << "The matrix at matrix_pt must be square.";
440  throw OomphLibError(error_message_stream.str(),
441  OOMPH_CURRENT_FUNCTION,
442  OOMPH_EXCEPTION_LOCATION);
443  }
444  // check that the matrix and the rhs vector have the same nrow()
445  if (matrix_pt->nrow() != rhs.nrow())
446  {
447  std::ostringstream error_message_stream;
448  error_message_stream
449  << "The matrix and the rhs vector must have the same number of rows.";
450  throw OomphLibError(error_message_stream.str(),
451  OOMPH_CURRENT_FUNCTION,
452  OOMPH_EXCEPTION_LOCATION);
453  }
454 
455  // if the matrix is distributable then it too should have the same
456  // communicator as the rhs vector and should not be distributed
457  DistributableLinearAlgebraObject* dist_matrix_pt =
458  dynamic_cast<DistributableLinearAlgebraObject*>(matrix_pt);
459  if (dist_matrix_pt != 0)
460  {
461  if (dist_matrix_pt->distribution_pt()->communicator_pt()->nproc() > 1 &&
462  dist_matrix_pt->distribution_pt()->distributed() == true)
463  {
464  throw OomphLibError(
465  "Matrix must not be distributed or only one processor",
466  OOMPH_CURRENT_FUNCTION,
467  OOMPH_EXCEPTION_LOCATION);
468  }
470  if (!(temp_comm == *dist_matrix_pt->distribution_pt()->communicator_pt()))
471  {
472  std::ostringstream error_message_stream;
473  error_message_stream
474  << "The matrix matrix_pt must have the same communicator as the "
475  "vectors"
476  << " rhs and result must have the same communicator";
477  throw OomphLibError(error_message_stream.str(),
478  OOMPH_CURRENT_FUNCTION,
479  OOMPH_EXCEPTION_LOCATION);
480  }
481  }
482  // if the result vector is setup then check it is not distributed and has
483  // the same communicator as the rhs vector
484  if (result.distribution_built())
485  {
486  if (!(*result.distribution_pt() == *rhs.distribution_pt()))
487  {
488  std::ostringstream error_message_stream;
489  error_message_stream
490  << "The result vector distribution has been setup; it must have the "
491  << "same distribution as the rhs vector.";
492  throw OomphLibError(error_message_stream.str(),
493  OOMPH_CURRENT_FUNCTION,
494  OOMPH_EXCEPTION_LOCATION);
495  }
496  }
497 #endif
498 
499  if (!result.distribution_built())
500  {
501  result.build(rhs.distribution_pt(), 0.0);
502  }
503 
504  // set the distribution
505  this->build_distribution(rhs.distribution_pt());
506 
507  // Time the solver
508  double t_start = TimingHelpers::timer();
509 
510  // factorise
511  factorise(matrix_pt);
512 
513  // backsubstitute
514  backsub(rhs, result);
515 
516  // Doc time for solver
517  double t_end = TimingHelpers::timer();
518 
519  Solution_time = t_end - t_start;
520  if (Doc_time)
521  {
522  oomph_info << std::endl
523  << "CPU for solve with DenseLU : "
526  << std::endl
527  << std::endl;
528  }
529 
530  // If we are not resolving then delete storage
531  if (!Enable_resolve)
532  {
533  clean_up_memory();
534  }
535  }
536 
537  //=============================================================================
538  /// Linear-algebra-type solver: Takes pointer to a matrix and rhs
539  /// vector and returns the solution of the linear system.
540  //=============================================================================
541  void DenseLU::solve(DoubleMatrixBase* const& matrix_pt,
542  const Vector<double>& rhs,
543  Vector<double>& result)
544  {
545  // Time the solver
546  clock_t t_start = clock();
547 
548  factorise(matrix_pt);
549  backsub(rhs, result);
550 
551  // Doc time for solver
552  clock_t t_end = clock();
553 
554  Solution_time = double(t_end - t_start) / CLOCKS_PER_SEC;
555  if (Doc_time)
556  {
557  oomph_info << "CPU for solve with DenseLU : "
560  << std::endl;
561  }
562 
563  // If we are not resolving then delete storage
564  if (!Enable_resolve)
565  {
566  clean_up_memory();
567  }
568  }
569 
570  //==================================================================
571  /// Solver: Takes pointer to problem and returns the results Vector
572  /// which contains the solution of the linear system defined by
573  /// the problem's residual Vector. (Jacobian assembled by FD).
574  //==================================================================
575  void FD_LU::solve(Problem* const& problem_pt, DoubleVector& result)
576  {
577  // Initialise timer
578  clock_t t_start = clock();
579 
580 #ifdef PARANOID
581  // if the result vector is setup then check it is not distributed and has
582  // the same communicator as the rhs vector
583  if (result.built())
584  {
585  if (result.distributed())
586  {
587  std::ostringstream error_message_stream;
588  error_message_stream << "The result vector must not be distributed";
589  throw OomphLibError(error_message_stream.str(),
590  OOMPH_CURRENT_FUNCTION,
591  OOMPH_EXCEPTION_LOCATION);
592  }
593  }
594 #endif
595 
596  // Find # of degrees of freedom
597  unsigned long n_dof = problem_pt->ndof();
598 
599  // Allocate storage for the residuals vector and the jacobian matrix
600  DoubleVector residuals;
601  DenseDoubleMatrix jacobian(n_dof);
602 
603  {
604  // initialise timer
605  clock_t t_start = clock();
606 
607  // Get the full jacobian by finite differencing) VERY INEFFICIENT!
608  problem_pt->get_fd_jacobian(residuals, jacobian);
609 
610  // compute jacobian setup time
611  clock_t t_end = clock();
612  Jacobian_setup_time = double(t_end - t_start) / CLOCKS_PER_SEC;
613 
614  // Report the time
615  if (Doc_time)
616  {
617  oomph_info << std::endl
618  << "CPU for setup of Dense Jacobian: "
621  << std::endl
622  << std::endl;
623  }
624  }
625 
626  // Solve by dense LU decomposition (not efficient)
627  solve(&jacobian, residuals, result);
628 
629  // Set the sign of the determinant of the jacobian
631 
632  // Finalise/doc timings
633  clock_t t_end = clock();
634  double total_time = double(t_end - t_start) / CLOCKS_PER_SEC;
635  if (Doc_time)
636  {
637  oomph_info << "CPU for FD DenseLU LinearSolver: "
639  << std::endl
640  << std::endl;
641  }
642  }
643 
644 
645  //===================================================================
646  // Interface to SuperLU wrapper
647  //===================================================================
648  extern "C"
649  {
650  int superlu(int*,
651  int*,
652  int*,
653  int*,
654  double*,
655  int*,
656  int*,
657  double*,
658  int*,
659  int*,
660  int*,
661  void*,
662  int*);
663  }
664 
665 
666 #ifdef OOMPH_HAS_MPI
667  //===================================================================
668  // Interface to SuperLU_DIST wrapper
669  //===================================================================
670  extern "C"
671  {
672  // Interface to distributed SuperLU solver where each processor
673  // holds the entire matrix
674  void superlu_dist_global_matrix(int opt_flag,
675  int allow_permutations,
676  int n,
677  int nnz,
678  double* values,
679  int* row_index,
680  int* col_start,
681  double* b,
682  int nprow,
683  int npcol,
684  int doc,
685  void** data,
686  int* info,
687  MPI_Comm comm);
688 
689  // Interface to distributed SuperLU solver where each processor
690  // holds part of the matrix
692  int allow_permutations,
693  int n,
694  int nnz_local,
695  int nrow_local,
696  int first_row,
697  double* values,
698  int* col_index,
699  int* row_start,
700  double* b,
701  int nprow,
702  int npcol,
703  int doc,
704  void** data,
705  int* info,
706  MPI_Comm comm);
707 
708  // helper method - just calls the superlu method dCompRow_to_CompCol to
709  // convert the c-style vectors of a cr matrix to a cc matrix
710  void superlu_cr_to_cc(int nrow,
711  int ncol,
712  int nnz,
713  double* cr_values,
714  int* cr_index,
715  int* cr_start,
716  double** cc_values,
717  int** cc_index,
718  int** cc_start);
719  }
720 #endif
721 
722 
723  //===================================================================
724  // Interface to SuperLU wrapper extras
725  //===================================================================
726  extern "C"
727  {
728  /// Function to calculate the number of bytes used to store the
729  /// LU factors
731 
732  /// Function to calculate the number of bytes used in calculating
733  /// and storing the LU factors
735  }
736 
737 #ifdef OOMPH_HAS_MPI
738  //===================================================================
739  // Interface to SuperLU_DIST wrapper extras
740  //===================================================================
741  extern "C"
742  {
743  /// Function to calculate the number of bytes used to store the
744  /// LU factors
746 
747  /// Function to calculate the number of bytes used in calculating
748  /// and storing the LU factors
750  }
751 #endif
752 
753  //=============================================================================
754  /// How much memory do the LU factors take up? In bytes
755  /// NOTE: This has been scraped from dQuerySpace(...) in dmemory.c in
756  /// external_src/oomph_superlu_4.3
757  //=============================================================================
759  {
760  // If we're using the non-distributed version of SuperLU and the LU
761  // factors have also been computed
762  if ((Solver_type != Distributed) && (Serial_f_factors != 0))
763  {
765  }
766 #ifdef OOMPH_HAS_MPI
767  // If we're using SuperLU dist and the LU factors have been computed
768  if ((Solver_type == Distributed) && (Dist_solver_data_pt != 0))
769  {
771  }
772 #endif
773  // If the factors haven't been computed we can't do anything
774  else
775  {
776  return 0.0;
777  }
778  } // End of get_memory_usage_for_lu_factors
779 
780 
781  //=============================================================================
782  /// How much memory was used in total? In bytes
783  /// NOTE: This has been scraped from dQuerySpace(...) in dmemory.c in
784  /// external_src/oomph_superlu_4.3
785  //=============================================================================
787  {
788  // If we're using the non-distributed version of SuperLU and the LU
789  // factors have also been computed
790  if ((Solver_type != Distributed) && (Serial_f_factors != 0))
791  {
793  }
794 #ifdef OOMPH_HAS_MPI
795  // If we're using SuperLU dist and the LU factors have been computed
796  if ((Solver_type == Distributed) && (Dist_solver_data_pt != 0))
797  {
799  }
800 #endif
801  // If the factors haven't been computed we can't do anything
802  else
803  {
804  return 0.0;
805  }
806  } // End of get_total_needed_memory
807 
808 
809  //==========================================================================
810  /// Solver: Takes pointer to problem and returns the results Vector
811  /// which contains the solution of the linear system defined by
812  /// the problem's fully assembled Jacobian and residual Vector.
813  //==========================================================================
814  void SuperLUSolver::solve(Problem* const& problem_pt, DoubleVector& result)
815  {
816  // wipe memory
817  this->clean_up_memory();
818 
819 #ifdef OOMPH_HAS_MPI
820  // USING SUPERLU DIST
821  /// //////////////////
822  if (Solver_type == Distributed ||
823  (Solver_type == Default && problem_pt->communicator_pt()->nproc() > 1))
824  {
825  // init the timers
826  double t_start = TimingHelpers::timer();
827 
828  // number of dofs
829  unsigned n_dof = problem_pt->ndof();
830 
831  // set the distribution
833  problem_pt->communicator_pt(), n_dof, !Dist_use_global_solver);
834  this->build_distribution(dist);
835 
836  // Take a copy of Delete_matrix_data
837  bool copy_of_Delete_matrix_data = Dist_delete_matrix_data;
838 
839  // Set Delete_matrix to true
841 
842  // Use the distributed version of SuperLU_DIST?
844  {
845  // Initialise timer
846  double t_start = TimingHelpers::timer();
847 
848  // Storage for the residuals vector
849  DoubleVector residuals(this->distribution_pt(), 0.0);
850 
851  // Get the sparse jacobian and residuals of the problem
852  CRDoubleMatrix jacobian(this->distribution_pt());
853  problem_pt->get_jacobian(residuals, jacobian);
854 
855  // Doc time for setup
856  double t_end = TimingHelpers::timer();
857  Jacobian_setup_time = t_end - t_start;
858  if (Doc_time)
859  {
860  oomph_info << "Time to set up CRDoubleMatrix Jacobian : "
863  << std::endl;
864  }
865 
866  // Now call the linear algebra solve, if desired
867  if (!Suppress_solve)
868  {
869  // If the distribution of the result has been build and
870  // does not match that of
871  // the solver then redistribute before the solve and return
872  // to the incoming distribution afterwards.
873  if ((result.built()) &&
874  (!(*result.distribution_pt() == *this->distribution_pt())))
875  {
876  LinearAlgebraDistribution temp_global_dist(
877  result.distribution_pt());
878  result.build(this->distribution_pt(), 0.0);
879  solve(&jacobian, residuals, result);
880  result.redistribute(&temp_global_dist);
881  }
882  else
883  {
884  solve(&jacobian, residuals, result);
885  }
886  }
887  }
888  // Otherwise its the global solve version
889  else
890  {
891  // Storage for the residuals vector
892  // A non-distriubted residuals vector
894  problem_pt->communicator_pt(), problem_pt->ndof(), false);
895  DoubleVector residuals(&dist, 0.0);
896  CRDoubleMatrix jacobian(&dist);
897 
898  // Get the sparse jacobian and residuals of the problem
899  problem_pt->get_jacobian(residuals, jacobian);
900 
901  // Doc time for setup
902  double t_end = TimingHelpers::timer();
903  Jacobian_setup_time = t_end - t_start;
904  if (Doc_time)
905  {
906  oomph_info << "Time to set up CR Jacobian : "
909  << std::endl;
910  }
911 
912  // Now call the linear algebra solve, if desired
913  if (!Suppress_solve)
914  {
915  // If the result distribution has been built and
916  // does not match the global distribution
917  // the redistribute before the solve and then return to the
918  // distributed version afterwards
919  if ((result.built()) && (!(*result.distribution_pt() == dist)))
920  {
921  LinearAlgebraDistribution temp_global_dist(
922  result.distribution_pt());
923  result.build(&dist, 0.0);
924  solve(&jacobian, residuals, result);
925  result.redistribute(&temp_global_dist);
926  }
927  else
928  {
929  solve(&jacobian, residuals, result);
930  }
931  }
932  }
933  // Set Delete_matrix back to original value
934  Dist_delete_matrix_data = copy_of_Delete_matrix_data;
935  }
936 
937  // OTHERWISE WE ARE USING SUPERLU (SERIAL)
938  /// ///////////////////////////////////////
939  else
940 #endif
941  {
942  // set the solver distribution
944  problem_pt->communicator_pt(), problem_pt->ndof(), false);
945  this->build_distribution(dist);
946 
947  // Allocate storage for the residuals vector
948  DoubleVector residuals(dist, 0.0);
949 
950  // Use the compressed row version?
952  {
953  // Initialise timer
954  double t_start = TimingHelpers::timer();
955 
956  // Get the sparse jacobian and residuals of the problem
957  CRDoubleMatrix CR_jacobian(this->distribution_pt());
958  problem_pt->get_jacobian(residuals, CR_jacobian);
959 
960  // If we want to compute the gradient for the globally convergent
961  // Newton method, then do it here
962  if (Compute_gradient)
963  {
964  // Compute it
965  CR_jacobian.multiply_transpose(residuals,
967  // Set the flag
969  }
970 
971  // Doc time for setup
972  double t_end = TimingHelpers::timer();
973  Jacobian_setup_time = t_end - t_start;
974  if (Doc_time)
975  {
976  oomph_info << std::endl
977  << "Time to set up CRDoubleMatrix Jacobian : "
980  << std::endl;
981  }
982 
983  // Now call the linear algebra solve, if desired
984  if (!Suppress_solve)
985  {
986  // If the result vector is built and distributed
987  // then need to redistribute into the same form as the
988  // RHS (non-distributed)
989  if ((result.built()) &&
990  (!(*result.distribution_pt() == *this->distribution_pt())))
991  {
992  LinearAlgebraDistribution temp_global_dist(
993  result.distribution_pt());
994  result.build(this->distribution_pt(), 0.0);
995  solve(&CR_jacobian, residuals, result);
996  result.redistribute(&temp_global_dist);
997  }
998  // Otherwise just solve
999  else
1000  {
1001  solve(&CR_jacobian, residuals, result);
1002  }
1003  }
1004  }
1005  // Otherwise its the compressed column version
1006  else
1007  {
1008  // Initialise timer
1009  double t_start = TimingHelpers::timer();
1010 
1011  // Get the sparse jacobian and residuals of the problem
1012  CCDoubleMatrix CC_jacobian;
1013  problem_pt->get_jacobian(residuals, CC_jacobian);
1014 
1015  // If we want to compute the gradient for the globally convergent
1016  // Newton method, then do it here
1017  if (Compute_gradient)
1018  {
1019  // Compute it
1020  CC_jacobian.multiply_transpose(residuals,
1022  // Set the flag
1024  }
1025 
1026  // Doc time for setup
1027  double t_end = TimingHelpers::timer();
1028  Jacobian_setup_time = t_end - t_start;
1029  if (Doc_time)
1030  {
1031  oomph_info << "\nTime to set up CCDoubleMatrix Jacobian: "
1034  << std::endl;
1035  }
1036 
1037  // Now call the linear algebra solve, if desired
1038  if (!Suppress_solve)
1039  {
1040  // If the result vector is built and distributed
1041  // then need to redistribute into the same form as the
1042  // RHS
1043  if ((result.built()) &&
1044  (!(*result.distribution_pt() == *this->distribution_pt())))
1045  {
1046  LinearAlgebraDistribution temp_global_dist(
1047  result.distribution_pt());
1048  result.build(this->distribution_pt(), 0.0);
1049  solve(&CC_jacobian, residuals, result);
1050  result.redistribute(&temp_global_dist);
1051  }
1052  // Otherwise just solve
1053  else
1054  {
1055  solve(&CC_jacobian, residuals, result);
1056  }
1057  }
1058  }
1059 
1060  // Set the sign of the jacobian
1061  //(this is computed in the LU decomposition phase)
1063  }
1064  }
1065 
1066  //=========================================================================
1067  /// Linear-algebra-type solver: Takes pointer to a matrix and rhs
1068  /// vector and returns the solution of the linear system. Problem pointer
1069  /// defaults to NULL and can be omitted. The function returns the global
1070  /// result Vector.
1071  /// Note: if Delete_matrix_data is true the function
1072  /// matrix_pt->clean_up_memory() will be used to wipe the matrix data.
1073  //=========================================================================
1074  void SuperLUSolver::solve(DoubleMatrixBase* const& matrix_pt,
1075  const DoubleVector& rhs,
1076  DoubleVector& result)
1077  {
1078  // Initialise timer
1079  double t_start = TimingHelpers::timer();
1080 
1081  // Pointer used in various places
1082  CRDoubleMatrix* cr_pt = 0;
1083 
1084 
1085 #ifdef PARANOID
1086  // check that the rhs vector is setup
1087  if (!rhs.built())
1088  {
1089  std::ostringstream error_message_stream;
1090  error_message_stream << "The vectors rhs must be setup";
1091  throw OomphLibError(error_message_stream.str(),
1092  OOMPH_CURRENT_FUNCTION,
1093  OOMPH_EXCEPTION_LOCATION);
1094  }
1095 
1096  // check that the matrix is square
1097  if (matrix_pt->nrow() != matrix_pt->ncol())
1098  {
1099  std::ostringstream error_message_stream;
1100  error_message_stream << "The matrix at matrix_pt must be square.";
1101  throw OomphLibError(error_message_stream.str(),
1102  OOMPH_CURRENT_FUNCTION,
1103  OOMPH_EXCEPTION_LOCATION);
1104  }
1105 
1106  // check that the matrix has some entries, and so has a values_pt that
1107  // makes sense (only for CR because CC is never used I think dense
1108  // matrices will be safe since they don't use a values pointer).
1109  cr_pt = dynamic_cast<CRDoubleMatrix*>(matrix_pt);
1110  if (cr_pt != 0)
1111  {
1112  if (cr_pt->nnz() == 0)
1113  {
1114  std::ostringstream error_message_stream;
1115  error_message_stream
1116  << "Attempted to call SuperLu on a CRDoubleMatrix with no entries, "
1117  << "SuperLU would segfault (because the values array pt is "
1118  << "uninitialised or null).";
1119  throw OomphLibError(error_message_stream.str(),
1120  OOMPH_CURRENT_FUNCTION,
1121  OOMPH_EXCEPTION_LOCATION);
1122  }
1123  }
1124 
1125  // check that the matrix and the rhs vector have the same nrow()
1126  if (matrix_pt->nrow() != rhs.nrow())
1127  {
1128  std::ostringstream error_message_stream;
1129  error_message_stream
1130  << "The matrix and the rhs vector must have the same number of rows.";
1131  throw OomphLibError(error_message_stream.str(),
1132  OOMPH_CURRENT_FUNCTION,
1133  OOMPH_EXCEPTION_LOCATION);
1134  }
1135 
1136  // if the matrix is distributable then should have the same distribution
1137  // as the rhs vector
1138  DistributableLinearAlgebraObject* dist_matrix_pt =
1139  dynamic_cast<DistributableLinearAlgebraObject*>(matrix_pt);
1140  if (dist_matrix_pt != 0)
1141  {
1142  if (!(*dist_matrix_pt->distribution_pt() == *rhs.distribution_pt()))
1143  {
1144  std::ostringstream error_message_stream;
1145  error_message_stream
1146  << "The matrix matrix_pt must have the same distribution as the "
1147  << "rhs vector.";
1148  throw OomphLibError(error_message_stream.str(),
1149  OOMPH_CURRENT_FUNCTION,
1150  OOMPH_EXCEPTION_LOCATION);
1151  }
1152  }
1153  // if the matrix is not distributable then it the rhs vector should not be
1154  // distributed
1155  else
1156  {
1157  if (rhs.distribution_pt()->distributed())
1158  {
1159  std::ostringstream error_message_stream;
1160  error_message_stream
1161  << "The matrix (matrix_pt) is not distributable and therefore the rhs"
1162  << " vector must not be distributed";
1163  throw OomphLibError(error_message_stream.str(),
1164  OOMPH_CURRENT_FUNCTION,
1165  OOMPH_EXCEPTION_LOCATION);
1166  }
1167  }
1168  // if the result vector is setup then check it has the same distribution
1169  // as the rhs
1170  if (result.built())
1171  {
1172  if (!(*result.distribution_pt() == *rhs.distribution_pt()))
1173  {
1174  std::ostringstream error_message_stream;
1175  error_message_stream
1176  << "The result vector distribution has been setup; it must have the "
1177  << "same distribution as the rhs vector.";
1178  throw OomphLibError(error_message_stream.str(),
1179  OOMPH_CURRENT_FUNCTION,
1180  OOMPH_EXCEPTION_LOCATION);
1181  }
1182  }
1183 #endif
1184 
1185  // set the distribution
1186  if (dynamic_cast<DistributableLinearAlgebraObject*>(matrix_pt))
1187  {
1188  // the solver has the same distribution as the matrix if possible
1189  this->build_distribution(
1190  dynamic_cast<DistributableLinearAlgebraObject*>(matrix_pt)
1191  ->distribution_pt());
1192  }
1193  else
1194  {
1195  // the solver has the same distribution as the RHS
1196  this->build_distribution(rhs.distribution_pt());
1197  }
1198 
1199  // Doc time for solve
1200  double t_factorise_start = TimingHelpers::timer();
1201 
1202  // Factorise the matrix
1203  factorise(matrix_pt);
1204 
1205  // Doc the end time
1206  double t_factorise_end = TimingHelpers::timer();
1207 
1208  // How long did the factorisation take?
1209  double factorise_time = t_factorise_end - t_factorise_start;
1210 
1211  // Try and upcast the matrix to a CRDoubleMatrix
1212  // CRDoubleMatrix*
1213  cr_pt = dynamic_cast<CRDoubleMatrix*>(matrix_pt);
1214 
1215  // If the input matrix is a CRDoubleMatrix
1216  if (cr_pt != 0)
1217  {
1218  // ...and actually has an entry
1219  if (cr_pt->nnz() != 0)
1220  {
1221  // Find out how many rows there are in the global Jacobian
1222  unsigned n_row = cr_pt->nrow();
1223 
1224  // And how many non-zeros there are in the global Jacobian
1225  unsigned n_nnz = cr_pt->nnz();
1226 
1227  // Get the memory usage (in bytes) for the global Jacobian storage
1228  double memory_usage_for_jacobian =
1229  ((2 * ((n_row + 1) * sizeof(int))) +
1230  (n_nnz * (sizeof(int) + sizeof(double))));
1231 
1232  // Get the memory usage (in bytes) for storage of the LU factors in
1233  // SuperLU
1234  double memory_usage_for_lu_storage = get_total_needed_memory();
1235 
1236  // Get the memory usage (in bytes) for storage of the LU factors in
1237  // SuperLU
1238  double total_memory_usage =
1239  memory_usage_for_jacobian + memory_usage_for_lu_storage;
1240 
1241 
1242  // How much memory have we used?
1243  if (Doc_stats)
1244  {
1245  oomph_info << "\nMemory statistics:"
1246  << "\n - Memory used to store the Jacobian (MB): "
1247  << memory_usage_for_jacobian / 1.0e+06
1248  << "\n - Memory used to store the LU factors (MB): "
1249  << memory_usage_for_lu_storage / 1.0e+06
1250  << "\n - Total memory used for matrix storage (MB): "
1251  << total_memory_usage / 1.0e+06 << "\n"
1252  << std::endl;
1253  }
1254  }
1255  } // if (cr_pt!=0)
1256 
1257  // Doc the start time
1258  double t_backsub_start = TimingHelpers::timer();
1259 
1260  // Now do the back solve
1261  backsub(rhs, result);
1262 
1263  // Doc the end time
1264  double t_backsub_end = TimingHelpers::timer();
1265 
1266  // How long did the back substitution take?
1267  double backsub_time = t_backsub_end - t_backsub_start;
1268 
1269  // Doc time for solve
1270  double t_end = TimingHelpers::timer();
1271  Solution_time = t_end - t_start;
1272  if (Doc_time)
1273  {
1274  oomph_info
1275  << "Time for LU factorisation : "
1277  << "\nTime for back-substitution: "
1279  << "\nTime for SuperLUSolver solve (ndof=" << matrix_pt->nrow() << "): "
1281  << std::endl;
1282  }
1283 
1284  // If we are not storing the solver data for resolves, delete it
1285  if (!Enable_resolve)
1286  {
1287  clean_up_memory();
1288  }
1289  }
1290 
1291 
1292  //=============================================================================
1293  /// Solver: Takes pointer to problem and returns the results Vector
1294  /// which contains the solution of the linear system defined by
1295  /// the problem's fully assembled Jacobian and residual Vector.
1296  //=============================================================================
1297  void SuperLUSolver::solve_transpose(Problem* const& problem_pt,
1298  DoubleVector& result)
1299  {
1300  // wipe memory
1301  this->clean_up_memory();
1302 
1303 #ifdef OOMPH_HAS_MPI
1304  // USING SUPERLU DIST
1305  /// //////////////////
1306  if (Solver_type == Distributed ||
1307  (Solver_type == Default && problem_pt->communicator_pt()->nproc() > 1))
1308  {
1309  // init the timers
1310  double t_start = TimingHelpers::timer();
1311 
1312  // number of dofs
1313  unsigned n_dof = problem_pt->ndof();
1314 
1315  // set the distribution
1317  problem_pt->communicator_pt(), n_dof, !Dist_use_global_solver);
1318  this->build_distribution(dist);
1319 
1320  // Take a copy of Delete_matrix_data
1321  bool copy_of_Delete_matrix_data = Dist_delete_matrix_data;
1322 
1323  // Set Delete_matrix to true
1324  Dist_delete_matrix_data = true;
1325 
1326  // Use the distributed version of SuperLU_DIST?
1328  {
1329  // Initialise timer
1330  double t_start = TimingHelpers::timer();
1331 
1332  // Storage for the residuals vector
1333  DoubleVector residuals(this->distribution_pt(), 0.0);
1334 
1335  // Get the sparse jacobian and residuals of the problem
1336  CRDoubleMatrix jacobian(this->distribution_pt());
1337  problem_pt->get_jacobian(residuals, jacobian);
1338 
1339  // Doc time for setup
1340  double t_end = TimingHelpers::timer();
1341  Jacobian_setup_time = t_end - t_start;
1342  if (Doc_time)
1343  {
1344  oomph_info << "Time to set up CRDoubleMatrix Jacobian : "
1347  << std::endl;
1348  }
1349 
1350  // Now call the linear algebra solve, if desired
1351  if (!Suppress_solve)
1352  {
1353  // If the distribution of the result has been build and
1354  // does not match that of
1355  // the solver then redistribute before the solve and return
1356  // to the incoming distribution afterwards.
1357  if ((result.built()) &&
1358  (!(*result.distribution_pt() == *this->distribution_pt())))
1359  {
1360  LinearAlgebraDistribution temp_global_dist(
1361  result.distribution_pt());
1362  result.build(this->distribution_pt(), 0.0);
1363  solve_transpose(&jacobian, residuals, result);
1364  result.redistribute(&temp_global_dist);
1365  }
1366  else
1367  {
1368  solve_transpose(&jacobian, residuals, result);
1369  }
1370  }
1371  }
1372  // Otherwise its the global solve version
1373  else
1374  {
1375  // Storage for the residuals vector
1376  // A non-distriubted residuals vector
1378  problem_pt->communicator_pt(), problem_pt->ndof(), false);
1379  DoubleVector residuals(&dist, 0.0);
1380  CRDoubleMatrix jacobian(&dist);
1381 
1382  // Get the sparse jacobian and residuals of the problem
1383  problem_pt->get_jacobian(residuals, jacobian);
1384 
1385  // Doc time for setup
1386  double t_end = TimingHelpers::timer();
1387  Jacobian_setup_time = t_end - t_start;
1388  if (Doc_time)
1389  {
1390  oomph_info << "Time to set up CR Jacobian : "
1393  << std::endl;
1394  }
1395 
1396  // Now call the linear algebra solve, if desired
1397  if (!Suppress_solve)
1398  {
1399  // If the result distribution has been built and
1400  // does not match the global distribution
1401  // the redistribute before the solve and then return to the
1402  // distributed version afterwards
1403  if ((result.built()) && (!(*result.distribution_pt() == dist)))
1404  {
1405  LinearAlgebraDistribution temp_global_dist(
1406  result.distribution_pt());
1407  result.build(&dist, 0.0);
1408  solve_transpose(&jacobian, residuals, result);
1409  result.redistribute(&temp_global_dist);
1410  }
1411  else
1412  {
1413  solve_transpose(&jacobian, residuals, result);
1414  }
1415  }
1416  }
1417  // Set Delete_matrix back to original value
1418  Dist_delete_matrix_data = copy_of_Delete_matrix_data;
1419  }
1420 
1421  // OTHERWISE WE ARE USING SUPERLU (SERIAL)
1422  /// ///////////////////////////////////////
1423  else
1424 #endif
1425  {
1426  // set the solver distribution
1428  problem_pt->communicator_pt(), problem_pt->ndof(), false);
1429  this->build_distribution(dist);
1430 
1431  // Allocate storage for the residuals vector
1432  DoubleVector residuals(dist, 0.0);
1433 
1434  // Use the compressed row version?
1436  {
1437  // Initialise timer
1438  double t_start = TimingHelpers::timer();
1439 
1440  // Get the sparse jacobian and residuals of the problem
1441  CRDoubleMatrix CR_jacobian(this->distribution_pt());
1442  problem_pt->get_jacobian(residuals, CR_jacobian);
1443 
1444  // If we want to compute the gradient for the globally convergent
1445  // Newton method, then do it here
1446  if (Compute_gradient)
1447  {
1448  // Compute it
1449  CR_jacobian.multiply_transpose(residuals,
1451  // Set the flag
1453  }
1454 
1455  // Doc time for setup
1456  double t_end = TimingHelpers::timer();
1457  Jacobian_setup_time = t_end - t_start;
1458  if (Doc_time)
1459  {
1460  oomph_info << std::endl
1461  << "Time to set up CRDoubleMatrix Jacobian: "
1464  << std::endl;
1465  }
1466 
1467  // Now call the linear algebra solve, if desired
1468  if (!Suppress_solve)
1469  {
1470  // If the result vector is built and distributed
1471  // then need to redistribute into the same form as the
1472  // RHS (non-distributed)
1473  if ((result.built()) &&
1474  (!(*result.distribution_pt() == *this->distribution_pt())))
1475  {
1476  LinearAlgebraDistribution temp_global_dist(
1477  result.distribution_pt());
1478  result.build(this->distribution_pt(), 0.0);
1479  solve_transpose(&CR_jacobian, residuals, result);
1480  result.redistribute(&temp_global_dist);
1481  }
1482  // Otherwise just solve
1483  else
1484  {
1485  solve_transpose(&CR_jacobian, residuals, result);
1486  }
1487  }
1488  }
1489  // Otherwise its the compressed column version
1490  else
1491  {
1492  // Initialise timer
1493  double t_start = TimingHelpers::timer();
1494 
1495  // Get the sparse jacobian and residuals of the problem
1496  CCDoubleMatrix CC_jacobian;
1497  problem_pt->get_jacobian(residuals, CC_jacobian);
1498 
1499  // If we want to compute the gradient for the globally convergent
1500  // Newton method, then do it here
1501  if (Compute_gradient)
1502  {
1503  // Compute it
1504  CC_jacobian.multiply_transpose(residuals,
1506  // Set the flag
1508  }
1509 
1510  // Doc time for setup
1511  double t_end = TimingHelpers::timer();
1512  Jacobian_setup_time = t_end - t_start;
1513  if (Doc_time)
1514  {
1515  oomph_info << "\nTime to set up CCDoubleMatrix Jacobian: "
1518  << std::endl;
1519  }
1520 
1521  // Now call the linear algebra solve, if desired
1522  if (!Suppress_solve)
1523  {
1524  // If the result vector is built and distributed
1525  // then need to redistribute into the same form as the
1526  // RHS
1527  if ((result.built()) &&
1528  (!(*result.distribution_pt() == *this->distribution_pt())))
1529  {
1530  LinearAlgebraDistribution temp_global_dist(
1531  result.distribution_pt());
1532  result.build(this->distribution_pt(), 0.0);
1533  solve_transpose(&CC_jacobian, residuals, result);
1534  result.redistribute(&temp_global_dist);
1535  }
1536  // Otherwise just solve
1537  else
1538  {
1539  solve_transpose(&CC_jacobian, residuals, result);
1540  }
1541  }
1542  }
1543 
1544  // Set the sign of the jacobian
1545  //(this is computed in the LU decomposition phase)
1547  }
1548  }
1549 
1550  //=========================================================================
1551  /// Linear-algebra-type solver: Takes pointer to a matrix and rhs
1552  /// vector and returns the solution of the linear system. Problem pointer
1553  /// defaults to NULL and can be omitted. The function returns the global
1554  /// result Vector.
1555  /// Note: if Delete_matrix_data is true the function
1556  /// matrix_pt->clean_up_memory() will be used to wipe the matrix data.
1557  //=========================================================================
1559  const DoubleVector& rhs,
1560  DoubleVector& result)
1561  {
1562  // Initialise timer
1563  double t_start = TimingHelpers::timer();
1564 
1565  // Pointer used in various places
1566  CRDoubleMatrix* cr_pt = 0;
1567 
1568 #ifdef PARANOID
1569  // check that the rhs vector is setup
1570  if (!rhs.built())
1571  {
1572  std::ostringstream error_message_stream;
1573  error_message_stream << "The vectors rhs must be setup";
1574  throw OomphLibError(error_message_stream.str(),
1575  OOMPH_CURRENT_FUNCTION,
1576  OOMPH_EXCEPTION_LOCATION);
1577  }
1578 
1579  // check that the matrix is square
1580  if (matrix_pt->nrow() != matrix_pt->ncol())
1581  {
1582  std::ostringstream error_message_stream;
1583  error_message_stream << "The matrix at matrix_pt must be square.";
1584  throw OomphLibError(error_message_stream.str(),
1585  OOMPH_CURRENT_FUNCTION,
1586  OOMPH_EXCEPTION_LOCATION);
1587  }
1588 
1589  // check that the matrix has some entries, and so has a values_pt that
1590  // makes sense (only for CR because CC is never used I think dense
1591  // matrices will be safe since they don't use a values pointer).
1592  cr_pt = dynamic_cast<CRDoubleMatrix*>(matrix_pt);
1593  if (cr_pt != 0)
1594  {
1595  if (cr_pt->nnz() == 0)
1596  {
1597  std::ostringstream error_message_stream;
1598  error_message_stream
1599  << "Attempted to call SuperLu on a CRDoubleMatrix with no entries, "
1600  << "SuperLU would segfault (because the values array pt is "
1601  << "uninitialised or null).";
1602  throw OomphLibError(error_message_stream.str(),
1603  OOMPH_CURRENT_FUNCTION,
1604  OOMPH_EXCEPTION_LOCATION);
1605  }
1606  }
1607 
1608  // check that the matrix and the rhs vector have the same nrow()
1609  if (matrix_pt->nrow() != rhs.nrow())
1610  {
1611  std::ostringstream error_message_stream;
1612  error_message_stream
1613  << "The matrix and the rhs vector must have the same number of rows.";
1614  throw OomphLibError(error_message_stream.str(),
1615  OOMPH_CURRENT_FUNCTION,
1616  OOMPH_EXCEPTION_LOCATION);
1617  }
1618 
1619  // if the matrix is distributable then should have the same distribution
1620  // as the rhs vector
1621  DistributableLinearAlgebraObject* dist_matrix_pt =
1622  dynamic_cast<DistributableLinearAlgebraObject*>(matrix_pt);
1623  if (dist_matrix_pt != 0)
1624  {
1625  if (!(*dist_matrix_pt->distribution_pt() == *rhs.distribution_pt()))
1626  {
1627  std::ostringstream error_message_stream;
1628  error_message_stream
1629  << "The matrix matrix_pt must have the same distribution as the "
1630  << "rhs vector.";
1631  throw OomphLibError(error_message_stream.str(),
1632  OOMPH_CURRENT_FUNCTION,
1633  OOMPH_EXCEPTION_LOCATION);
1634  }
1635  }
1636  // if the matrix is not distributable then it the rhs vector should not be
1637  // distributed
1638  else
1639  {
1640  if (rhs.distribution_pt()->distributed())
1641  {
1642  std::ostringstream error_message_stream;
1643  error_message_stream
1644  << "The matrix (matrix_pt) is not distributable and therefore the rhs"
1645  << " vector must not be distributed";
1646  throw OomphLibError(error_message_stream.str(),
1647  OOMPH_CURRENT_FUNCTION,
1648  OOMPH_EXCEPTION_LOCATION);
1649  }
1650  }
1651  // if the result vector is setup then check it has the same distribution
1652  // as the rhs
1653  if (result.built())
1654  {
1655  if (!(*result.distribution_pt() == *rhs.distribution_pt()))
1656  {
1657  std::ostringstream error_message_stream;
1658  error_message_stream
1659  << "The result vector distribution has been setup; it must have the "
1660  << "same distribution as the rhs vector.";
1661  throw OomphLibError(error_message_stream.str(),
1662  OOMPH_CURRENT_FUNCTION,
1663  OOMPH_EXCEPTION_LOCATION);
1664  }
1665  }
1666 #endif
1667 
1668  // set the distribution
1669  if (dynamic_cast<DistributableLinearAlgebraObject*>(matrix_pt))
1670  {
1671  // the solver has the same distribution as the matrix if possible
1672  this->build_distribution(
1673  dynamic_cast<DistributableLinearAlgebraObject*>(matrix_pt)
1674  ->distribution_pt());
1675  }
1676  else
1677  {
1678  // the solver has the same distribution as the RHS
1679  this->build_distribution(rhs.distribution_pt());
1680  }
1681 
1682  // Doc time for solve
1683  double t_factorise_start = TimingHelpers::timer();
1684 
1685  // Factorise the matrix
1686  factorise(matrix_pt);
1687 
1688  // Doc the end time
1689  double t_factorise_end = TimingHelpers::timer();
1690 
1691  // How long did the factorisation take?
1692  double factorise_time = t_factorise_end - t_factorise_start;
1693 
1694  // Try and upcast the matrix to a CRDoubleMatrix
1695  // CRDoubleMatrix*
1696  cr_pt = dynamic_cast<CRDoubleMatrix*>(matrix_pt);
1697 
1698  // If the input matrix is a CRDoubleMatrix
1699  if (cr_pt != 0)
1700  {
1701  // ...and actually has an entry
1702  if (cr_pt->nnz() != 0)
1703  {
1704  // Find out how many rows there are in the global Jacobian
1705  unsigned n_row = cr_pt->nrow();
1706 
1707  // And how many non-zeros there are in the global Jacobian
1708  unsigned n_nnz = cr_pt->nnz();
1709 
1710  // Get the memory usage (in bytes) for the global Jacobian storage
1711  double memory_usage_for_jacobian =
1712  ((2 * ((n_row + 1) * sizeof(int))) +
1713  (n_nnz * (sizeof(int) + sizeof(double))));
1714 
1715  // Get the memory usage (in bytes) for storage of the LU factors in
1716  // SuperLU
1717  double memory_usage_for_lu_storage = get_total_needed_memory();
1718 
1719  // Get the memory usage (in bytes) for storage of the LU factors in
1720  // SuperLU
1721  double total_memory_usage =
1722  memory_usage_for_jacobian + memory_usage_for_lu_storage;
1723 
1724  // How much memory have we used?
1725  if (Doc_stats)
1726  {
1727  oomph_info << "\nMemory statistics:"
1728  << "\n - Memory used to store the Jacobian (MB): "
1729  << memory_usage_for_jacobian / 1.0e+06
1730  << "\n - Memory used to store the LU factors (MB): "
1731  << memory_usage_for_lu_storage / 1.0e+06
1732  << "\n - Total memory used for matrix storage (MB): "
1733  << total_memory_usage / 1.0e+06 << "\n"
1734  << std::endl;
1735  }
1736  }
1737  } // if (cr_pt!=0)
1738 
1739  // Doc the start time
1740  double t_backsub_start = TimingHelpers::timer();
1741 
1742  // Now do the back solve
1743  backsub_transpose(rhs, result);
1744 
1745  // Doc the end time
1746  double t_backsub_end = TimingHelpers::timer();
1747 
1748  // How long did the back substitution take?
1749  double backsub_time = t_backsub_end - t_backsub_start;
1750 
1751  // Doc time for solve
1752  double t_end = TimingHelpers::timer();
1753  Solution_time = t_end - t_start;
1754  if (Doc_time)
1755  {
1756  oomph_info
1757  << "Time for LU factorisation : "
1759  << "\nTime for back-substitution: "
1761  << "\nTime for SuperLUSolver solve (ndof=" << matrix_pt->nrow() << "): "
1763  << std::endl;
1764  }
1765 
1766  // If we are not storing the solver data for resolves, delete it
1767  if (!Enable_resolve)
1768  {
1769  clean_up_memory();
1770  }
1771  } // End of solve_transpose
1772 
1773  //===============================================================
1774  /// Resolve the system for a given RHS
1775  //===============================================================
1777  {
1778  // Store starting time for solve
1779  double t_start = TimingHelpers::timer();
1780 
1781  // backsub
1782  backsub(rhs, result);
1783 
1784  // Doc time for solve
1785  double t_end = TimingHelpers::timer();
1786  Solution_time = t_end - t_start;
1787  if (Doc_time)
1788  {
1789  oomph_info << "Time for SuperLUSolver solve (ndof=" << rhs.nrow() << "): "
1791  t_start)
1792  << std::endl;
1793  }
1794  }
1795 
1796 
1797  //===============================================================
1798  /// Resolve the (transposed) system for a given RHS
1799  //===============================================================
1801  DoubleVector& result)
1802  {
1803  // Store starting time for solve
1804  double t_start = TimingHelpers::timer();
1805 
1806  // Backsub (but solve the transposed system)
1807  backsub_transpose(rhs, result);
1808 
1809  // Doc time for solve
1810  double t_end = TimingHelpers::timer();
1811  Solution_time = t_end - t_start;
1812  if (Doc_time)
1813  {
1814  oomph_info << "Time for SuperLUSolver solve (ndof=" << rhs.nrow() << "): "
1816  t_start)
1817  << std::endl;
1818  }
1819  }
1820 
1821 
1822  //===================================================================
1823  /// LU decompose the matrix addressed by matrix_pt by using
1824  /// the SuperLU solver. The resulting matrix factors are stored
1825  /// internally.
1826  //===================================================================
1828  {
1829  // wipe memory
1830  this->clean_up_memory();
1831 
1832  // if we have mpi and the solver is distributed or default and nproc
1833  // gt 1
1834 #ifdef OOMPH_HAS_MPI
1835  DistributableLinearAlgebraObject* dist_matrix_pt =
1836  dynamic_cast<DistributableLinearAlgebraObject*>(matrix_pt);
1837  unsigned nproc = 1;
1838  if (dist_matrix_pt != 0)
1839  {
1840  nproc = dist_matrix_pt->distribution_pt()->communicator_pt()->nproc();
1841  }
1842  if (Solver_type == Distributed || (Solver_type == Default && nproc > 1 &&
1844  {
1845  // if the matrix is a distributed linear algebra object then use SuperLU
1846  // dist
1847  if (dist_matrix_pt != 0)
1848  {
1849  factorise_distributed(matrix_pt);
1850  Using_dist = true;
1851  }
1852  else
1853  {
1854  factorise_serial(matrix_pt);
1855  Using_dist = false;
1856  }
1857  }
1858  else
1859 #endif
1860  {
1861  factorise_serial(matrix_pt);
1862  Using_dist = false;
1863  }
1864  }
1865 
1866 #ifdef OOMPH_HAS_MPI
1867  //=============================================================================
1868  /// LU decompose the matrix addressed by matrix_pt using
1869  /// the SuperLU_DIST solver. The resulting matrix factors are stored
1870  /// internally.
1871  //=============================================================================
1873  {
1874  // Check that we have a square matrix
1875 #ifdef PARANOID
1876  int m = matrix_pt->ncol();
1877  int n = matrix_pt->nrow();
1878  if (n != m)
1879  {
1880  std::ostringstream error_message_stream;
1881  error_message_stream << "Can only solve for square matrices\n"
1882  << "N, M " << n << " " << m << std::endl;
1883 
1884  throw OomphLibError(error_message_stream.str(),
1885  OOMPH_CURRENT_FUNCTION,
1886  OOMPH_EXCEPTION_LOCATION);
1887  }
1888 #endif
1889 
1890  // number of processors
1891  unsigned nproc = MPI_Helpers::communicator_pt()->nproc();
1892  if (dynamic_cast<DistributableLinearAlgebraObject*>(matrix_pt) != 0)
1893  {
1894  nproc = dynamic_cast<DistributableLinearAlgebraObject*>(matrix_pt)
1895  ->distribution_pt()
1896  ->communicator_pt()
1897  ->nproc();
1898  }
1899 
1900  // Find number of rows and columns for the process grid
1901  // First guess at number of rows:
1902  int nprow = int(sqrt(double(nproc)));
1903 
1904  // Does this evenly divide the processor grid?
1905  while (nprow > 1)
1906  {
1907  if (nproc % nprow == 0) break;
1908  nprow -= 1;
1909  }
1910 
1911  // Store Number of rows/columns for process grid
1912  Dist_nprow = nprow;
1913  Dist_npcol = nproc / Dist_nprow;
1914 
1915  // Make sure any existing factors are deleted
1916  clean_up_memory();
1917 
1918  // Doc (0/1) = (true/false)
1919  int doc = !Doc_stats;
1920 
1921  // Rset Info
1922  Dist_info = 0;
1923 
1924  // Flag for row and column permutations
1925  int allow_permutations = Dist_allow_row_and_col_permutations;
1926 
1927  // Is it a DistributedCRDoubleMatrix?
1928  if (dynamic_cast<CRDoubleMatrix*>(matrix_pt) != 0)
1929  {
1930  // Get a cast pointer to the matrix
1931  CRDoubleMatrix* cr_matrix_pt = dynamic_cast<CRDoubleMatrix*>(matrix_pt);
1932 
1933  // Get the distribution from the matrix
1934  this->build_distribution(cr_matrix_pt->distribution_pt());
1935 
1936 #ifdef PARANOID
1937  // paranoid check that the matrix has been setup
1938  if (!cr_matrix_pt->built())
1939  {
1940  throw OomphLibError(
1941  "To apply SuperLUSolver to a CRDoubleMatrix - it must be built",
1942  OOMPH_CURRENT_FUNCTION,
1943  OOMPH_EXCEPTION_LOCATION);
1944  }
1945 #endif
1946 
1947  // if the matrix is distributed then setup setup superlu dist distributed
1948  if (cr_matrix_pt->distributed())
1949  {
1950  // Find the number of non-zero entries in the matrix
1951  const int nnz_local = int(cr_matrix_pt->nnz());
1952 
1953  // Set up the pointers to the matrix.
1954  // NOTE: these arrays (accessed via value_pt, index_pt and
1955  // start_pt) may be modified by the SuperLU_DIST routines, and so
1956  // a copy must be taken if the matrix is to be preserved.
1957 
1958  // Copy values
1959  Dist_value_pt = new double[nnz_local];
1960  double* matrix_value_pt = cr_matrix_pt->value();
1961  for (int i = 0; i < nnz_local; i++)
1962  {
1963  Dist_value_pt[i] = matrix_value_pt[i];
1964  }
1965 
1966  // Copy column indices
1967  Dist_index_pt = new int[nnz_local];
1968  int* matrix_index_pt = cr_matrix_pt->column_index();
1969  for (int i = 0; i < nnz_local; i++)
1970  {
1971  Dist_index_pt[i] = matrix_index_pt[i];
1972  }
1973 
1974  // Copy row starts
1975  int nrow_local = cr_matrix_pt->nrow_local();
1976  Dist_start_pt = new int[nrow_local + 1];
1977  int* matrix_start_pt = cr_matrix_pt->row_start();
1978  for (int i = 0; i <= nrow_local; i++)
1979  {
1980  Dist_start_pt[i] = matrix_start_pt[i];
1981  }
1982 
1983  // cache
1984  int ndof = cr_matrix_pt->distribution_pt()->nrow();
1985  int first_row = cr_matrix_pt->first_row();
1986 
1987  // Now delete the matrix if we are allowed
1988  if (Dist_delete_matrix_data == true)
1989  {
1990  cr_matrix_pt->clear();
1991  }
1992 
1993  // Factorize
1995  1,
1996  allow_permutations,
1997  ndof,
1998  nnz_local,
1999  nrow_local,
2000  first_row,
2001  Dist_value_pt,
2002  Dist_index_pt,
2003  Dist_start_pt,
2004  0,
2005  Dist_nprow,
2006  Dist_npcol,
2007  doc,
2009  &Dist_info,
2010  this->distribution_pt()->communicator_pt()->mpi_comm());
2011 
2012  // Record that data is stored
2014  }
2015  // else the CRDoubleMatrix is not distributed
2016  else
2017  {
2018  // Find the number of non-zero entries in the matrix
2019  const int nnz = int(cr_matrix_pt->nnz());
2020 
2021  // cache the number of rows
2022  int nrow = cr_matrix_pt->nrow();
2023 
2024  // Set up the pointers to the matrix.
2025  // NOTE: these arrays (accessed via value_pt, index_pt and
2026  // start_pt) may be modified by the SuperLU_DIST routines, and so
2027  // a copy must be taken if the matrix is to be preserved.
2028 
2029  // create the corresponing cc matrix
2031  nrow,
2032  nnz,
2033  cr_matrix_pt->value(),
2034  cr_matrix_pt->column_index(),
2035  cr_matrix_pt->row_start(),
2036  &Dist_value_pt,
2037  &Dist_index_pt,
2038  &Dist_start_pt);
2039 
2040  // Delete the matrix if we are allowed
2041  if (Dist_delete_matrix_data == true)
2042  {
2043  cr_matrix_pt->clear();
2044  }
2045 
2046  // do the factorization
2048  1,
2049  allow_permutations,
2050  nrow,
2051  nnz,
2052  Dist_value_pt,
2053  Dist_index_pt,
2054  Dist_start_pt,
2055  0,
2056  Dist_nprow,
2057  Dist_npcol,
2058  doc,
2060  &Dist_info,
2061  this->distribution_pt()->communicator_pt()->mpi_comm());
2062 
2063  // Record that data is stored
2065  }
2066  }
2067 
2068  // Or is it a CCDoubleMatrix?
2069  else if (dynamic_cast<CCDoubleMatrix*>(matrix_pt))
2070  {
2071  // Get a cast pointer to the matrix
2072  CCDoubleMatrix* serial_matrix_pt =
2073  dynamic_cast<CCDoubleMatrix*>(matrix_pt);
2074 
2075  // Find the number of non-zero entries in the matrix
2076  const int nnz = int(serial_matrix_pt->nnz());
2077 
2078  // Find # of degrees of freedom (variables)
2079  int ndof = int(serial_matrix_pt->nrow());
2080 
2081  // Find the local number of degrees of freedom in the linear system
2082  int ndof_local = ndof;
2083 
2084  // Set up the pointers to the matrix.
2085  // NOTE: these arrays (accessed via value_pt, index_pt and
2086  // start_pt) may be modified by the SuperLU_DIST routines, and so
2087  // a copy must be taken if the matrix is to be preserved.
2088 
2089  // Copy values
2090  Dist_value_pt = new double[nnz];
2091  double* matrix_value_pt = serial_matrix_pt->value();
2092  for (int i = 0; i < nnz; i++)
2093  {
2094  Dist_value_pt[i] = matrix_value_pt[i];
2095  }
2096 
2097  // copy row indices
2098  Dist_index_pt = new int[nnz];
2099  int* matrix_index_pt = serial_matrix_pt->row_index();
2100  for (int i = 0; i < nnz; i++)
2101  {
2102  Dist_index_pt[i] = matrix_index_pt[i];
2103  }
2104 
2105  // copy column starts
2106  Dist_start_pt = new int[ndof_local + 1];
2107  int* matrix_start_pt = serial_matrix_pt->column_start();
2108  for (int i = 0; i <= ndof_local; i++)
2109  {
2110  Dist_start_pt[i] = matrix_start_pt[i];
2111  }
2112 
2113  // Delete the matrix if we are allowed
2114  if (Dist_delete_matrix_data == true)
2115  {
2116  serial_matrix_pt->clean_up_memory();
2117  }
2118 
2119  // do the factorization
2121  1,
2122  allow_permutations,
2123  ndof,
2124  nnz,
2125  Dist_value_pt,
2126  Dist_index_pt,
2127  Dist_start_pt,
2128  0,
2129  Dist_nprow,
2130  Dist_npcol,
2131  doc,
2133  &Dist_info,
2134  this->distribution_pt()->communicator_pt()->mpi_comm());
2135 
2136  // Record that data is stored
2138  }
2139  // Otherwise throw an error
2140  else
2141  {
2142  std::ostringstream error_message_stream;
2143  error_message_stream << "SuperLUSolver implemented only for "
2144  << " CCDoubleMatrix, CRDoubleMatrix\n"
2145  << "and DistributedCRDoubleMatrix matrices\n";
2146  throw OomphLibError(error_message_stream.str(),
2147  OOMPH_CURRENT_FUNCTION,
2148  OOMPH_EXCEPTION_LOCATION);
2149  }
2150 
2151  // Throw an error if superLU returned an error status in info.
2152  if (Dist_info != 0)
2153  {
2154  std::ostringstream error_msg;
2155  error_msg << "SuperLU returned the error status code " << Dist_info
2156  << " . See the SuperLU documentation for what this means.";
2157  throw OomphLibError(
2158  error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2159  }
2160  }
2161 #endif
2162 
2163  //===================================================================
2164  /// LU decompose the matrix addressed by matrix_pt by using
2165  /// the SuperLU solver. The resulting matrix factors are stored
2166  /// internally.
2167  //===================================================================
2169  {
2170 #ifdef PARANOID
2171  // PARANOID check that if the matrix is distributable then it should not be
2172  // then it should not be distributed
2173  if (dynamic_cast<DistributableLinearAlgebraObject*>(matrix_pt) != 0)
2174  {
2175  if (dynamic_cast<DistributableLinearAlgebraObject*>(matrix_pt)
2176  ->distributed())
2177  {
2178  std::ostringstream error_message_stream;
2179  error_message_stream << "The matrix must not be distributed.";
2180  throw OomphLibError(error_message_stream.str(),
2181  OOMPH_CURRENT_FUNCTION,
2182  OOMPH_EXCEPTION_LOCATION);
2183  }
2184  }
2185 #endif
2186 
2187  // Find # of degrees of freedom (variables)
2188  int n = matrix_pt->nrow();
2189 
2190  // Check that we have a square matrix
2191 #ifdef PARANOID
2192  int m = matrix_pt->ncol();
2193  if (n != m)
2194  {
2195  std::ostringstream error_message_stream;
2196  error_message_stream << "Can only solve for square matrices\n"
2197  << "N, M " << n << " " << m << std::endl;
2198 
2199  throw OomphLibError(error_message_stream.str(),
2200  OOMPH_CURRENT_FUNCTION,
2201  OOMPH_EXCEPTION_LOCATION);
2202  }
2203 #endif
2204 
2205  // Storage for the values, rows and column indices
2206  // required by SuplerLU
2207  double* value = 0;
2208  int *index = 0, *start = 0;
2209 
2210  // Integer used to represent compressed row or column format
2211  // Default compressed row
2212  int transpose = 0;
2213 
2214  // Number of non-zero entries in the matrix
2215  int nnz = 0;
2216 
2217  // Doc flag (convert to int for SuperLU)
2218  int doc = Doc_stats;
2219 
2220  // Is it a CR matrix
2221  if (dynamic_cast<CRDoubleMatrix*>(matrix_pt))
2222  {
2223  // Set the appropriate row flags
2225  transpose = 1;
2226  // Get a cast pointer to the matrix
2227  CRDoubleMatrix* CR_matrix_pt = dynamic_cast<CRDoubleMatrix*>(matrix_pt);
2228 
2229  // Now set the pointers to the interanally stored values
2230  // and indices
2231  nnz = CR_matrix_pt->nnz();
2232  value = CR_matrix_pt->value();
2233  index = CR_matrix_pt->column_index();
2234  start = CR_matrix_pt->row_start();
2235  }
2236  // Otherwise is it the compressed column version?
2237  else if (dynamic_cast<CCDoubleMatrix*>(matrix_pt))
2238  {
2239  // Set the compressed row flag to false
2241  // Get a cast pointer to the matrix
2242  CCDoubleMatrix* CC_matrix_pt = dynamic_cast<CCDoubleMatrix*>(matrix_pt);
2243 
2244  // Now set the pointers to the interanally stored values
2245  // and indices
2246  nnz = CC_matrix_pt->nnz();
2247  value = CC_matrix_pt->value();
2248  index = CC_matrix_pt->row_index();
2249  start = CC_matrix_pt->column_start();
2250  }
2251  // Otherwise throw and error
2252  else
2253  {
2254  throw OomphLibError("SuperLU only works with CR or CC Double matrices",
2255  OOMPH_CURRENT_FUNCTION,
2256  OOMPH_EXCEPTION_LOCATION);
2257  }
2258 
2259  // Clean up any previous storage so that if this is called twice with
2260  // the same matrix, we don't get a memory leak
2261  clean_up_memory();
2262 
2263  // Perform the lu decompose phase (i=1)
2264  int i = 1;
2266  &n,
2267  &nnz,
2268  0,
2269  value,
2270  index,
2271  start,
2272  0,
2273  &n,
2274  &transpose,
2275  &doc,
2277  &Serial_info);
2278 
2279  // Throw an error if superLU returned an error status in info.
2280  if (Serial_info != 0)
2281  {
2282  std::ostringstream error_msg;
2283  error_msg << "SuperLU returned the error status code " << Serial_info
2284  << " . See the SuperLU documentation for what this means.";
2285  throw OomphLibError(
2286  error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2287  }
2288 
2289 
2290  // Set the number of degrees of freedom in the linear system
2291  Serial_n_dof = n;
2292  }
2293 
2294  //=============================================================================
2295  /// Do the backsubstitution for SuperLUSolver.
2296  /// Note - this method performs no paranoid checks - these are all performed
2297  /// in solve(...) and resolve(...)
2298  //=============================================================================
2300  {
2301 #ifdef OOMPH_HAS_MPI
2302  if (Using_dist)
2303  {
2304  backsub_distributed(rhs, result);
2305  }
2306  else
2307 #endif
2308  {
2309  backsub_serial(rhs, result);
2310  }
2311  }
2312 
2313 
2314  //=============================================================================
2315  /// Do the backsubstitution of the transposed system for SuperLUSolver.
2316  /// Note - this method performs no paranoid checks - these are all performed
2317  /// in solve(...) and resolve(...)
2318  //=============================================================================
2320  DoubleVector& result)
2321  {
2322 #ifdef OOMPH_HAS_MPI
2323  if (Using_dist)
2324  {
2325  backsub_transpose_distributed(rhs, result);
2326  }
2327  else
2328 #endif
2329  {
2330  backsub_transpose_serial(rhs, result);
2331  }
2332  }
2333 
2334 #ifdef OOMPH_HAS_MPI
2335  //=========================================================================
2336  /// Static warning to suppress warnings about incorrect distribution of
2337  /// RHS vector. Default is false
2338  //=========================================================================
2340  false;
2341 
2342  //=============================================================================
2343  /// Do the backsubstitution for SuperLU solver.
2344  /// Note - this method performs no paranoid checks - these are all performed
2345  /// in solve(...) and resolve(...)
2346  //=============================================================================
2348  DoubleVector& result)
2349  {
2350 #ifdef PARANOID
2351  // check that the rhs vector is setup
2352  if (!rhs.distribution_pt()->built())
2353  {
2354  std::ostringstream error_message_stream;
2355  error_message_stream << "The vectors rhs must be setup";
2356  throw OomphLibError(error_message_stream.str(),
2357  OOMPH_CURRENT_FUNCTION,
2358  OOMPH_EXCEPTION_LOCATION);
2359  }
2360 #endif
2361  // check that the rhs distribution is the same as the distribution as this
2362  // solver. If not redistribute and issue a warning
2363  LinearAlgebraDistribution rhs_distribution(rhs.distribution_pt());
2364  if (!(*rhs.distribution_pt() == *this->distribution_pt()))
2365  {
2367  {
2368  std::ostringstream warning_stream;
2369  warning_stream << "The distribution of rhs vector does not match that "
2370  "ofthe solver.\n";
2371  warning_stream << "The rhs will be redistributed, which is likely to "
2372  "be inefficient\n";
2373  warning_stream
2374  << "To remove this warning you can either:\n"
2375  << " i) Ensure that the rhs vector has the correct distribution\n"
2376  << " before calling the resolve() function\n"
2377  << "or ii) Set the flag \n"
2378  << " SuperLUSolver::Suppress_incorrect_rhs_distribution_warning_in_"
2379  "resolve\n"
2380  << " to be true\n\n";
2381 
2382  OomphLibWarning(warning_stream.str(),
2383  "SuperLUSolver::resolve()",
2384  OOMPH_EXCEPTION_LOCATION);
2385  }
2386 
2387  // Have to cast away const-ness (which tells us that we shouldn't really
2388  // be doing this!)
2389  const_cast<DoubleVector&>(rhs).redistribute(this->distribution_pt());
2390  }
2391 
2392 #ifdef PARANOID
2393  // if the result vector is setup then check it has the same distribution
2394  // as the rhs
2395  if (result.distribution_built())
2396  {
2397  if (!(*result.distribution_pt() == *rhs.distribution_pt()))
2398  {
2399  std::ostringstream error_message_stream;
2400  error_message_stream
2401  << "The result vector distribution has been setup; it must have the "
2402  << "same distribution as the rhs vector.";
2403  throw OomphLibError(error_message_stream.str(),
2404  OOMPH_CURRENT_FUNCTION,
2405  OOMPH_EXCEPTION_LOCATION);
2406  }
2407  }
2408 #endif
2409  // Doc (0/1) = (true/false)
2410  int doc = !Doc_stats;
2411 
2412  // Reset Info
2413  Dist_info = 0;
2414 
2415  // number of DOFs
2416  int ndof = this->distribution_pt()->nrow();
2417 
2418  // Copy the rhs values to result
2419  result = rhs;
2420 
2421  // Do the backsubsitition phase
2423  {
2424  // Call distributed solver
2426  2,
2427  -1,
2428  ndof,
2429  0,
2430  0,
2431  0,
2432  0,
2433  0,
2434  0,
2435  result.values_pt(),
2436  Dist_nprow,
2437  Dist_npcol,
2438  doc,
2440  &Dist_info,
2441  this->distribution_pt()->communicator_pt()->mpi_comm());
2442  }
2444  {
2445  // Call global solver
2447  2,
2448  -1,
2449  ndof,
2450  0,
2451  0,
2452  0,
2453  0,
2454  result.values_pt(),
2455  Dist_nprow,
2456  Dist_npcol,
2457  doc,
2459  &Dist_info,
2460  this->distribution_pt()->communicator_pt()->mpi_comm());
2461  }
2462  else
2463  {
2464  throw OomphLibError("The matrix factors have not been stored",
2465  OOMPH_CURRENT_FUNCTION,
2466  OOMPH_EXCEPTION_LOCATION);
2467  }
2468 
2469  // Throw an error if superLU returned an error status in info.
2470  if (Dist_info != 0)
2471  {
2472  std::ostringstream error_msg;
2473  error_msg << "SuperLU returned the error status code " << Dist_info
2474  << " . See the SuperLU documentation for what this means.";
2475  throw OomphLibError(
2476  error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2477  }
2478 
2479  // Redistribute to original distribution
2480  // Have to cast away const-ness (which tells us that we shouldn't really
2481  // be doing this!)
2482  const_cast<DoubleVector&>(rhs).redistribute(&rhs_distribution);
2483  }
2484 
2485  //=============================================================================
2486  /// Do the backsubstitution for SuperLU solver.
2487  /// Note - this method performs no paranoid checks - these are all performed
2488  /// in solve(...) and resolve(...)
2489  //=============================================================================
2491  DoubleVector& result)
2492  {
2493  // Create an output stream
2494  std::ostringstream error_message_stream;
2495 
2496  // Create the error message
2497  error_message_stream << "This function hasn't been implemented yet. If you "
2498  << "need it, implement it!" << std::endl;
2499 
2500  // Throw the error message
2501  throw OomphLibError(error_message_stream.str(),
2502  OOMPH_CURRENT_FUNCTION,
2503  OOMPH_EXCEPTION_LOCATION);
2504  }
2505 #endif
2506 
2507  //================================================================
2508  /// Do the backsubstitution for SuperLU
2509  //================================================================
2511  DoubleVector& result)
2512  {
2513  // Find the number of unknowns
2514  int n = rhs.nrow();
2515 
2516 #ifdef PARANOID
2517  // PARANOID check that this rhs distribution is setup
2518  if (!rhs.built())
2519  {
2520  std::ostringstream error_message_stream;
2521  error_message_stream << "The rhs vector distribution must be setup.";
2522  throw OomphLibError(error_message_stream.str(),
2523  OOMPH_CURRENT_FUNCTION,
2524  OOMPH_EXCEPTION_LOCATION);
2525  }
2526  // PARANOID check that the rhs has the right number of global rows
2527  if (static_cast<int>(Serial_n_dof) != n)
2528  {
2529  throw OomphLibError(
2530  "RHS does not have the same dimension as the linear system",
2531  OOMPH_CURRENT_FUNCTION,
2532  OOMPH_EXCEPTION_LOCATION);
2533  }
2534  // PARANOID check that the rhs is not distributed
2535  if (rhs.distribution_pt()->distributed())
2536  {
2537  std::ostringstream error_message_stream;
2538  error_message_stream << "The rhs vector must not be distributed.";
2539  throw OomphLibError(error_message_stream.str(),
2540  OOMPH_CURRENT_FUNCTION,
2541  OOMPH_EXCEPTION_LOCATION);
2542  }
2543  // PARANOID check that if the result is setup it matches the distribution
2544  // of the rhs
2545  if (result.built())
2546  {
2547  if (!(*rhs.distribution_pt() == *result.distribution_pt()))
2548  {
2549  std::ostringstream error_message_stream;
2550  error_message_stream << "If the result distribution is setup then it "
2551  "must be the same as the "
2552  << "rhs distribution";
2553  throw OomphLibError(error_message_stream.str(),
2554  OOMPH_CURRENT_FUNCTION,
2555  OOMPH_EXCEPTION_LOCATION);
2556  }
2557  }
2558 #endif
2559 
2560  // copy result to rhs
2561  result = rhs;
2562 
2563  // Number of RHSs
2564  int nrhs = 1;
2565 
2566  // Cast the boolean flags to ints for SuperLU
2567  int transpose = Serial_compressed_row_flag;
2568  int doc = Doc_stats;
2569 
2570  // Do the backsubsitition phase
2571  int i = 2;
2572  superlu(&i,
2573  &n,
2574  0,
2575  &nrhs,
2576  0,
2577  0,
2578  0,
2579  result.values_pt(),
2580  &n,
2581  &transpose,
2582  &doc,
2584  &Serial_info);
2585 
2586  // Throw an error if superLU returned an error status in info.
2587  if (Serial_info != 0)
2588  {
2589  std::ostringstream error_msg;
2590  error_msg << "SuperLU returned the error status code " << Serial_info
2591  << " . See the SuperLU documentation for what this means.";
2592  throw OomphLibError(
2593  error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2594  }
2595  }
2596 
2597  //================================================================
2598  /// Do the backsubstitution for SuperLU
2599  //================================================================
2601  DoubleVector& result)
2602  {
2603  // Find the number of unknowns
2604  int n = rhs.nrow();
2605 
2606 #ifdef PARANOID
2607  // PARANOID check that this rhs distribution is setup
2608  if (!rhs.built())
2609  {
2610  std::ostringstream error_message_stream;
2611  error_message_stream << "The rhs vector distribution must be setup.";
2612  throw OomphLibError(error_message_stream.str(),
2613  OOMPH_CURRENT_FUNCTION,
2614  OOMPH_EXCEPTION_LOCATION);
2615  }
2616  // PARANOID check that the rhs has the right number of global rows
2617  if (static_cast<int>(Serial_n_dof) != n)
2618  {
2619  throw OomphLibError(
2620  "RHS does not have the same dimension as the linear system",
2621  OOMPH_CURRENT_FUNCTION,
2622  OOMPH_EXCEPTION_LOCATION);
2623  }
2624  // PARANOID check that the rhs is not distributed
2625  if (rhs.distribution_pt()->distributed())
2626  {
2627  std::ostringstream error_message_stream;
2628  error_message_stream << "The rhs vector must not be distributed.";
2629  throw OomphLibError(error_message_stream.str(),
2630  OOMPH_CURRENT_FUNCTION,
2631  OOMPH_EXCEPTION_LOCATION);
2632  }
2633  // PARANOID check that if the result is setup it matches the distribution
2634  // of the rhs
2635  if (result.built())
2636  {
2637  if (!(*rhs.distribution_pt() == *result.distribution_pt()))
2638  {
2639  std::ostringstream error_message_stream;
2640  error_message_stream << "If the result distribution is setup then it "
2641  "must be the same as the "
2642  << "rhs distribution";
2643  throw OomphLibError(error_message_stream.str(),
2644  OOMPH_CURRENT_FUNCTION,
2645  OOMPH_EXCEPTION_LOCATION);
2646  }
2647  }
2648 #endif
2649 
2650  // copy result to rhs
2651  result = rhs;
2652 
2653  // Number of RHSs
2654  int nrhs = 1;
2655 
2656  // Cast the boolean flags to ints for SuperLU
2657  int transpose = (!Serial_compressed_row_flag);
2658  int doc = Doc_stats;
2659 
2660  // Do the backsubsitition phase
2661  int i = 2;
2662  superlu(&i,
2663  &n,
2664  0,
2665  &nrhs,
2666  0,
2667  0,
2668  0,
2669  result.values_pt(),
2670  &n,
2671  &transpose,
2672  &doc,
2674  &Serial_info);
2675 
2676  // Throw an error if superLU returned an error status in info.
2677  if (Serial_info != 0)
2678  {
2679  std::ostringstream error_msg;
2680  error_msg << "SuperLU returned the error status code " << Serial_info
2681  << " . See the SuperLU documentation for what this means.";
2682  throw OomphLibError(
2683  error_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
2684  }
2685  }
2686 
2687  //=============================================================================
2688  /// Clean up the memory
2689  //=============================================================================
2691  {
2692  // If we have non-zero LU factors stored
2693  if (Serial_f_factors != 0)
2694  {
2695  // Clean up those factors
2696  int i = 3;
2697  int transpose = Serial_compressed_row_flag;
2698  superlu(&i,
2699  0,
2700  0,
2701  0,
2702  0,
2703  0,
2704  0,
2705  0,
2706  0,
2707  &transpose,
2708  0,
2710  &Serial_info);
2711 
2712  // Set the F_factors to zero
2713  Serial_f_factors = 0;
2714  Serial_n_dof = 0;
2715  }
2716 
2717 #ifdef OOMPH_HAS_MPI
2718  // If we have non-zero LU factors stored
2719  if (Dist_solver_data_pt != 0)
2720  {
2721  // Clean up any stored solver data
2722 
2723  // Doc (0/1) = (true/false)
2724  int doc = !Doc_stats;
2725 
2726  // Reset Info flag
2727  Dist_info = 0;
2728 
2729  // number of DOFs
2730  int ndof = this->distribution_pt()->nrow();
2731 
2733  {
2735  3,
2736  -1,
2737  ndof,
2738  0,
2739  0,
2740  0,
2741  0,
2742  0,
2743  0,
2744  0,
2745  Dist_nprow,
2746  Dist_npcol,
2747  doc,
2749  &Dist_info,
2750  this->distribution_pt()->communicator_pt()->mpi_comm());
2752  }
2754  {
2756  3,
2757  -1,
2758  ndof,
2759  0,
2760  0,
2761  0,
2762  0,
2763  0,
2764  Dist_nprow,
2765  Dist_npcol,
2766  doc,
2768  &Dist_info,
2769  this->distribution_pt()->communicator_pt()->mpi_comm());
2771  }
2772 
2773  Dist_solver_data_pt = 0;
2774 
2775  // Delete internal copy of the matrix
2776  delete[] Dist_value_pt;
2777  delete[] Dist_index_pt;
2778  delete[] Dist_start_pt;
2779  Dist_value_pt = 0;
2780  Dist_index_pt = 0;
2781  Dist_start_pt = 0;
2782 
2783  // and the distribution
2784  this->clear_distribution();
2785  }
2786 #endif
2787  }
2788 
2789 } // namespace oomph
cstr elem_len * i
Definition: cfortran.h:603
//////////////////////////////////////////////////////////////// ////////////////////////////////////...
Definition: matrices.h:2791
void multiply_transpose(const DoubleVector &x, DoubleVector &soln) const
Multiply the transposed matrix by the vector x: soln=A^T x.
Definition: matrices.cc:715
unsigned long nrow() const
Return the number of rows of the matrix.
Definition: matrices.h:2817
int * column_start()
Access to C-style column_start array.
Definition: matrices.h:2692
int * row_index()
Access to C-style row index array.
Definition: matrices.h:2704
void clean_up_memory()
Wipe matrix data and set all values to 0.
Definition: matrices.h:3169
A class for compressed row matrices. This is a distributable object.
Definition: matrices.h:888
int * column_index()
Access to C-style column index array.
Definition: matrices.h:1072
int * row_start()
Access to C-style row_start array.
Definition: matrices.h:1060
void multiply_transpose(const DoubleVector &x, DoubleVector &soln) const
Multiply the transposed matrix by the vector x: soln=A^T x.
Definition: matrices.cc:1882
unsigned long nnz() const
Return the number of nonzero entries (the local nnz)
Definition: matrices.h:1096
bool built() const
access function to the Built flag - indicates whether the matrix has been build - i....
Definition: matrices.h:1210
double * value()
Access to C-style value array.
Definition: matrices.h:1084
unsigned long nrow() const
Return the number of rows of the matrix.
Definition: matrices.h:1002
Class of matrices containing doubles, and stored as a DenseMatrix<double>, but with solving functiona...
Definition: matrices.h:1271
double Jacobian_setup_time
Jacobian setup time.
void backsub(const DoubleVector &rhs, DoubleVector &result)
Do the backsubstitution step to solve the system LU result = rhs.
void clean_up_memory()
Clean up the stored LU factors.
long * Index
Pointer to storage for the index of permutations in the LU solve.
void factorise(DoubleMatrixBase *const &matrix_pt)
Perform the LU decomposition of the matrix.
double * LU_factors
Pointer to storage for the LU decomposition.
double Solution_time
Solution time.
void solve(Problem *const &problem_pt, DoubleVector &result)
Solver: Takes pointer to problem and returns the results Vector which contains the solution of the li...
int Sign_of_determinant_of_matrix
Sign of the determinant of the matrix (obtained during the LU decomposition)
Base class for any linear algebra object that is distributable. Just contains storage for the LinearA...
void clear_distribution()
clear the distribution of this distributable linear algebra object
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.
bool distribution_built() const
if the communicator_pt is null then the distribution is not setup then false is returned,...
unsigned nrow_local() const
access function for the num of local rows on this processor.
unsigned first_row() const
access function for the first row on this processor
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.
double * values_pt()
access function to the underlying values
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
void solve(Problem *const &problem_pt, DoubleVector &result)
Solver: Takes pointer to problem and returns the results Vector which contains the solution of the li...
Describes the distribution of a distributable linear algebra type object. Typically this is a contain...
bool distributed() const
access function to the distributed - indicates whether the distribution is serial or distributed
OomphCommunicator * communicator_pt() const
const access to the communicator pointer
bool built() const
if the communicator_pt is null then the distribution is not setup then false is returned,...
unsigned nrow() const
access function to the number of global rows.
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
bool Gradient_has_been_computed
flag that indicates whether the gradient was computed or not
Definition: linear_solver.h:84
DoubleVector Gradient_for_glob_conv_newton_solve
DoubleVector storing the gradient for the globally convergent Newton method.
Definition: linear_solver.h:88
bool Compute_gradient
flag that indicates whether the gradient required for the globally convergent Newton method should be...
Definition: linear_solver.h:81
static bool mpi_has_been_initialised()
return true if MPI has been initialised
static OomphCommunicator * communicator_pt()
access to global communicator. This is the oomph-lib equivalent of MPI_COMM_WORLD
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 OomphLibWarning object which should be created as a temporary object to issue a warning....
////////////////////////////////////////////////////////////////// //////////////////////////////////...
Definition: problem.h:151
unsigned long ndof() const
Return the number of dofs.
Definition: problem.h:1678
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:3921
OomphCommunicator * communicator_pt()
access function to the oomph-lib communicator
Definition: problem.h:1246
void get_fd_jacobian(DoubleVector &residuals, DenseMatrix< double > &jacobian)
Return the fully-assembled Jacobian and residuals, generated by finite differences.
Definition: problem.cc:7800
int & sign_of_jacobian()
Access function for the sign of the global jacobian matrix. This will be set by the linear solver,...
Definition: problem.h:2553
T * value()
Access to C-style value array.
Definition: matrices.h:616
unsigned long nnz() const
Return the number of nonzero entries.
Definition: matrices.h:640
bool Doc_stats
Set to true to output statistics (false by default).
bool Dist_delete_matrix_data
Delete_matrix_data flag. SuperLU_dist needs its own copy of the input matrix, therefore a copy must b...
void solve_transpose(Problem *const &problem_pt, DoubleVector &result)
Solver: Takes pointer to problem and returns the results Vector which contains the solution of the li...
static bool Suppress_incorrect_rhs_distribution_warning_in_resolve
Static flag that determines whether the warning about incorrect distribution of RHSs will be printed ...
int * Dist_index_pt
Pointer for storage of matrix rows or column indices required by SuperLU_DIST.
Type Solver_type
the solver type. see SuperLU_solver_type for details.
void backsub_serial(const DoubleVector &rhs, DoubleVector &result)
backsub method for SuperLU (serial)
void backsub_transpose_distributed(const DoubleVector &rhs, DoubleVector &result)
backsub method for SuperLU Dist
double Solution_time
Solution time.
double * Dist_value_pt
Pointer for storage of the matrix values required by SuperLU_DIST.
bool Serial_compressed_row_flag
Use compressed row version?
void solve(Problem *const &problem_pt, DoubleVector &result)
Solver: Takes pointer to problem and returns the results Vector which contains the solution of the li...
void factorise(DoubleMatrixBase *const &matrix_pt)
Do the factorisation stage Note: if Delete_matrix_data is true the function matrix_pt->clean_up_memor...
void factorise_serial(DoubleMatrixBase *const &matrix_pt)
factorise method for SuperLU (serial)
void backsub_distributed(const DoubleVector &rhs, DoubleVector &result)
backsub method for SuperLU Dist
bool Using_dist
boolean flag indicating whether superlu dist is being used
int * Dist_start_pt
Pointers for storage of matrix column or row starts.
unsigned long Serial_n_dof
The number of unknowns in the linear system.
void * Serial_f_factors
Storage for the LU factors as required by SuperLU.
void resolve_transpose(const DoubleVector &rhs, DoubleVector &result)
Resolve the (transposed) system defined by the last assembled Jacobian and the specified rhs vector i...
bool Suppress_solve
Suppress solve?
int Serial_sign_of_determinant_of_matrix
Sign of the determinant of the matrix.
bool Dist_use_global_solver
Flag that determines whether the MPIProblem based solve function uses the global or distributed versi...
double get_memory_usage_for_lu_factors()
How much memory do the LU factors take up? In bytes.
void factorise_distributed(DoubleMatrixBase *const &matrix_pt)
factorise method for SuperLU Dist
bool Dist_allow_row_and_col_permutations
If true then SuperLU_DIST is allowed to permute matrix rows and columns during factorisation....
int Dist_info
Info flag for the SuperLU solver.
bool Dist_global_solve_data_allocated
Flag is true if solve data has been generated for a global matrix.
double get_total_needed_memory()
How much memory was allocated by SuperLU? In bytes.
double Jacobian_setup_time
Jacobian setup time.
void backsub_transpose_serial(const DoubleVector &rhs, DoubleVector &result)
backsub method for SuperLU (serial)
void resolve(const DoubleVector &rhs, DoubleVector &result)
Resolve the system defined by the last assembled jacobian and the specified rhs vector if resolve has...
void * Dist_solver_data_pt
Storage for the LU factors and other data required by SuperLU.
bool Dist_distributed_solve_data_allocated
Flag is true if solve data has been generated for distributed matrix.
int Dist_npcol
Number of columns for the process grid.
void backsub(const DoubleVector &rhs, DoubleVector &result)
Do the backsubstitution for SuperLU solver Note: returns the global result Vector.
void clean_up_memory()
Clean up the memory allocated by the solver.
int Dist_nprow
Number of rows for the process grid.
void backsub_transpose(const DoubleVector &rhs, DoubleVector &result)
Do the back-substitution for transposed system of the SuperLU solver Note: Returns the global result ...
int Serial_info
Info flag for the SuperLU solver.
void start(const unsigned &i)
(Re-)start i-th timer
double timer()
returns the time in seconds after some point in past
std::string convert_secs_to_formatted_string(const double &time_in_sec)
Returns a nicely formatted string from an input time in seconds; the format depends on the size of ti...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
double get_total_memory_usage_in_bytes()
Function to calculate the number of bytes used in calculating and storing the LU factors.
double get_lu_factor_memory_usage_in_bytes_dist()
Function to calculate the number of bytes used to store the LU factors.
double get_total_memory_usage_in_bytes_dist()
Function to calculate the number of bytes used in calculating and storing the LU factors.
void superlu_dist_global_matrix(int opt_flag, int allow_permutations, int n, int nnz, double *values, int *row_index, int *col_start, double *b, int nprow, int npcol, int doc, void **data, int *info, MPI_Comm comm)
void superlu_cr_to_cc(int nrow, int ncol, int nnz, double *cr_values, int *cr_index, int *cr_start, double **cc_values, int **cc_index, int **cc_start)
int superlu(int *, int *, int *, int *, double *, int *, int *, double *, int *, int *, int *, void *, int *)
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...
double get_lu_factor_memory_usage_in_bytes()
Function to calculate the number of bytes used to store the LU factors.
void superlu_dist_distributed_matrix(int opt_flag, int allow_permutations, int n, int nnz_local, int nrow_local, int first_row, double *values, int *col_index, int *row_start, double *b, int nprow, int npcol, int doc, void **data, int *info, MPI_Comm comm)