Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Amesos_Scalapack.h
Go to the documentation of this file.
1 /*
2  TODO:
3  1) Eliminate mb_ and grid_mb_ DONE
4 
5  1) Allow user level control over whether to use a 1D or 2D data
6  distribution.
7  2) Figure out how to redistribute the vectors. (Allow no more than
8  nb vectors at a time for now). Done - I think
9  3) Move the code from software:MyExamples/TwodMap/TwoDMap.cpp in to
10  this routine (for the 2D case) Done - I think
11  4) Create the ScaLAPACK 2D dense matrix. Done - I think
12  */
13 
14 // @HEADER
15 // ***********************************************************************
16 //
17 // Amesos: Direct Sparse Solver Package
18 // Copyright (2004) Sandia Corporation
19 //
20 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
21 // license for use of this work by or on behalf of the U.S. Government.
22 //
23 // This library is free software; you can redistribute it and/or modify
24 // it under the terms of the GNU Lesser General Public License as
25 // published by the Free Software Foundation; either version 2.1 of the
26 // License, or (at your option) any later version.
27 //
28 // This library is distributed in the hope that it will be useful, but
29 // WITHOUT ANY WARRANTY; without even the implied warranty of
30 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
31 // Lesser General Public License for more details.
32 //
33 // You should have received a copy of the GNU Lesser General Public
34 // License along with this library; if not, write to the Free Software
35 // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301
36 // USA
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef _AMESOS_SCALAPACK_H_
43 #define _AMESOS_SCALAPACK_H_
44 
45 #if defined(Amesos_SHOW_DEPRECATED_WARNINGS)
46 #ifdef __GNUC__
47 #warning "The Amesos package is deprecated"
48 #endif
49 #endif
50 
51 #include "Amesos_ConfigDefs.h"
52 #include "Amesos_BaseSolver.h"
53 #include "Amesos_NoCopiable.h"
54 #include "Amesos_Utils.h"
55 #include "Amesos_Time.h"
56 #include "Amesos_Status.h"
57 #include "Amesos_Control.h"
58 #include "Epetra_LinearProblem.h"
59 #include "Epetra_Time.h"
60 #ifdef EPETRA_MPI
61 #include "Epetra_MpiComm.h"
62 #else
63 #include "Epetra_Comm.h"
64 #endif
65 #include "Epetra_CrsGraph.h"
66 
68 
153  private Amesos_Time,
154  private Amesos_NoCopiable,
155  private Amesos_Utils,
156  private Amesos_Control,
157  private Amesos_Status
158  {
159 
160 public:
161 
163 
171  Amesos_Scalapack( const Epetra_LinearProblem& LinearProblem );
172 
174 
176  ~Amesos_Scalapack(void);
178 
180 
182 
189  int SymbolicFactorization() ;
190 
192 
217  int NumericFactorization() ;
218 
220 
240  int Solve();
241 
243 
245 
246 #if 0
247  char * Label() const {return(Epetra_Object::Label());};
249 #endif
250 
252  const Epetra_LinearProblem *GetProblem() const { return(Problem_); };
253 
255 
258  bool MatrixShapeOK() const ;
259 
261 
274 
276  bool UseTranspose() const {return(UseTranspose_);};
277 
279  const Epetra_Comm & Comm() const {return(GetProblem()->GetOperator()->Comm());};
280 
282 
299 
302 
305 
307  int NumSolve() const { return( Amesos_Status::NumSolve_ ); }
308 
310  void PrintTiming() const;
311 
313  void PrintStatus() const;
314 
316  void GetTiming( Teuchos::ParameterList &TimingParameterList ) const { Amesos_Time::GetTiming(TimingParameterList); }
317 
319 
320  private:
321  /*
322  RedistributeA - Convert matrix to a dense ScaLAPACK matrix
323  Preconditions:
324  Problem_ must be set
325  SetParameters()
326  ScaLAPACK1DMap and ScaLAPACK1DMatrix must either be 0 or be pointers to
327  appropriatly allocate objects. If they are non-zero, those objects
328  will be deleted and recreated.
329 
330  Postconditions:
331  nprow_, npcol_, DescA_
332  RowMatrixA_
333  ScaLAPACK1DMap_
334  ScaLAPACK1DMatrix_
335  ImportToScaLAPACK1D_
336  ImportBackToOriginal_
337 
338  */
339  int RedistributeA();
340 
341  /*
342  ConvertToScalapack - Convert matirx to form expected by Scalapack: Ai, Ap, Aval
343  Preconditions:
344  Problem_
345  Postconditions:
346  nprow_, npcol_,
347  */
348  int ConvertToScalapack();
349 
350  /*
351  PerformNumericFactorization - Call Scalapack to perform numeric factorization
352  Preconditions:
353 
354  Postconditions:
355  DenseA_, DescA_
356  */
358 
359  protected:
360 
361  int iam_; // Process number (i.e. Comm().MyPID()
362 
363  int NumGlobalElements_; // Number of rows and columns in the Problem_->GetOperator()
365 
366  //
367  // The following variables are required for the ScaLAPACK interface:
368  //
369  int nprow_ ; // number of process rows: 1 for now
370  int npcol_ ; // number of process columns
371  int ictxt_ ; // BLACS context
372  int m_per_p_; // Number of columns per process
373  int DescA_[10]; // ScaLAPACK array descriptor
374 
375  Epetra_Map *ScaLAPACK1DMap_ ; // Points to a 1D Map which matches a ScaLAPACK 1D
376  // blocked (not block cyclic) distribution
377  Epetra_CrsMatrix *ScaLAPACK1DMatrix_ ; // Points to a ScaLAPACK 1D
378  // blocked (not block cyclic) distribution
379  Epetra_Map *VectorMap_ ; // Points to a Map for vectors X and B
380  std::vector<double> DenseA_; // The data in a ScaLAPACK 1D blocked format
381  std::vector<int> Ipiv_ ; // ScaLAPACK pivot information
384 
385 
388 
389 
390  // some timing internal to ScaLAPACK
391  double ConTime_; // time to convert to ScaLAPACKformat
392  double SymTime_; // time for symbolic factorization
393  double NumTime_; // time for numeric factorization
394  double SolTime_; // time for solution
395  double VecTime_; // time to redistribute vectors
396  double MatTime_; // time to redistribute matrix
397 
398  //
399  // Control of the data distribution
400  //
401  bool TwoD_distribution_; // True if 2D data distribution is used
402  int grid_nb_; // Row and Column blocking factor (only used in 2D distribution)
403  int mypcol_; // Process column in the ScaLAPACK2D grid
404  int myprow_; // Process row in the ScaLAPACK2D grid
406 
407  //
408  // Blocking factors (For both 1D and 2D data distributions)
409  //
410  int nb_;
411  int lda_;
412 
414 
415 }; // End of class Amesos_Scalapack
416 #endif /* _AMESOS_SCALAPACK_H_ */
int NumSymbolicFact_
Number of symbolic factorization phases.
Definition: Amesos_Status.h:73
int SetUseTranspose(bool UseTranspose)
SetUseTranpose(true) is more efficient in Amesos_Scalapack.
Amesos_Control: Container for some control variables.
void PrintTiming() const
Print timing information.
const Epetra_LinearProblem * GetProblem() const
Get a pointer to the Problem.
int SymbolicFactorization()
Performs SymbolicFactorization on the matrix A.
int Solve()
Solves A X = B (or AT X = B)
bool UseTranspose() const
Returns the current UseTranspose setting.
int NumNumericFact_
Number of numeric factorization phases.
Definition: Amesos_Status.h:75
Epetra_Map * VectorMap_
int NumSolve_
Number of solves.
Definition: Amesos_Status.h:77
int NumSolve() const
Returns the number of solves performed by this object.
void GetTiming(Teuchos::ParameterList &TimingParameterList) const
Extracts timing information from the current solver and places it in the parameter list...
void GetTiming(Teuchos::ParameterList &list) const
Load up the current timing information into the parameter list.
Definition: Amesos_Time.h:130
Amesos_Status: Container for some status variables.
Definition: Amesos_Status.h:26
bool MatrixShapeOK() const
Returns true if SCALAPACK can handle this matrix shape.
const Epetra_LinearProblem * Problem_
Amesos_Time: Container for timing information.
Definition: Amesos_Time.h:56
Epetra_Time * Time_
int NumNumericFact() const
Returns the number of numeric factorizations performed by this object.
int NumericFactorization()
Performs NumericFactorization on the matrix A.
Epetra_CrsMatrix * FatOut_
std::vector< int > Ipiv_
int SetParameters(Teuchos::ParameterList &ParameterList)
Updates internal variables.
int NumSymbolicFact() const
Returns the number of symbolic factorizations performed by this object.
std::vector< double > DenseA_
const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this matrix.
Amesos_Scalapack(const Epetra_LinearProblem &LinearProblem)
Amesos_Scalapack Constructor.
Epetra_Map * ScaLAPACK1DMap_
Amesos_BaseSolver: A pure virtual class for direct solution of real-valued double-precision operators...
Amesos_NoCopiable: Simple class to prevent the usage of copy constructor and operator =...
Epetra_CrsMatrix * ScaLAPACK1DMatrix_
~Amesos_Scalapack(void)
Amesos_Scalapack Destructor.
Amesos_Scalapack: A serial and parallel dense solver. For now, we implement only the unsymmetric ScaL...
Amesos_Utils: Collections of basic utilities.
Definition: Amesos_Utils.h:25
void PrintStatus() const
Print information about the factorization and solution phases.