IFPACK  Development
 All Classes Namespaces Files Functions Variables Enumerations Friends Pages
Ifpack_Hypre.h
1 /*@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack: Object-Oriented Algebraic Preconditioner Package
5 // Copyright (2002) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 */
42 
43 #ifndef IFPACK_HYPRE_H
44 #define IFPACK_HYPRE_H
45 
46 #if defined(Ifpack_SHOW_DEPRECATED_WARNINGS)
47 #ifdef __GNUC__
48 #warning "The Ifpack package is deprecated"
49 #endif
50 #endif
51 
52 #include "Ifpack_ConfigDefs.h"
53 #ifdef HAVE_HYPRE
54 
55 #include "Ifpack_Preconditioner.h"
56 #include "Ifpack_Condest.h"
57 #include "Ifpack_ScalingType.h"
58 #include "Epetra_CompObject.h"
59 #include "Epetra_MultiVector.h"
60 #include "Epetra_Vector.h"
61 #include "Epetra_CrsGraph.h"
62 #include "Epetra_CrsMatrix.h"
63 #include "Epetra_BlockMap.h"
64 #include "Epetra_Map.h"
65 #include "Epetra_Object.h"
66 #include "Epetra_Comm.h"
67 #include "Epetra_CrsMatrix.h"
68 #include "Epetra_Time.h"
69 #include "Teuchos_RefCountPtr.hpp"
70 #include "Teuchos_ArrayRCP.hpp"
71 #include "Epetra_MpiComm.h"
72 
73 
74 #include <map>
75 
76 // Hypre forward declarations (to avoid downstream header pollution)
77 
78 struct hypre_IJMatrix_struct;
79 typedef struct hypre_IJMatrix_struct *HYPRE_IJMatrix;
80 struct hypre_IJVector_struct;
81 typedef struct hypre_IJVector_struct *HYPRE_IJVector;
82 struct hypre_ParCSRMatrix_struct;
83 typedef struct hypre_ParCSRMatrix_struct* HYPRE_ParCSRMatrix;
84 struct hypre_ParVector_struct;
85 typedef struct hypre_ParVector_struct * HYPRE_ParVector;
86 struct hypre_Solver_struct;
87 typedef struct hypre_Solver_struct *HYPRE_Solver;
88 struct hypre_ParVector_struct;
89 typedef struct hypre_ParVector_struct hypre_ParVector;
90 //struct hypre_Vector;
91 
92 // Will only work if Hypre is built with HYPRE_BIGINT=OFF
93 typedef int HYPRE_Int;
94 
95 typedef HYPRE_Int (*HYPRE_PtrToParSolverFcn)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
96 
97 #ifndef HYPRE_ENUMS
98 #define HYPRE_ENUMS
99 enum Hypre_Solver{
101  BoomerAMG,
102  ParaSails,
103  Euclid,
104  AMS,
105  Hybrid,
106  PCG,
107  GMRES,
108  FlexGMRES,
109  LGMRES,
110  BiCGSTAB
111 };
112 
114 enum Hypre_Chooser{
115  Solver,
116  Preconditioner
117 };
118 #endif //HYPRE_ENUMS
119 
120 class FunctionParameter;
121 
122 namespace Teuchos {
123  class ParameterList;
124 }
125 
127 
132 class Ifpack_Hypre: public Ifpack_Preconditioner {
133 
134 public:
135  // @{ Constructors and destructors.
137  Ifpack_Hypre(Epetra_RowMatrix* A);
138 
140  ~Ifpack_Hypre(){ Destroy();}
141 
142  // @}
143  // @{ Construction methods
144 
146  int Initialize();
147 
149  bool IsInitialized() const{ return(IsInitialized_);}
150 
152 
154  int Compute();
155 
157  bool IsComputed() const{ return(IsComputed_);}
158 
159 
161  /* This method is only available if the Teuchos package is enabled.
162  This method recognizes six parameter names: Solver,
163  Preconditioner, SolveOrPrecondition, SetPreconditioner, NumFunctions and Functions. These names are
164  case sensitive. Solver requires an enumerated parameter of type Hypre_Solver. Preconditioner is similar
165  except requires the type be a preconditioner. The options are listed below:
166  Solvers Preconditioners
167  BoomerAMG BoomerAMG
168  AMS ParaSails
169  Hybrid AMS
170  PCG (Default) Euclid (Default)
171  GMRES
172  FlexGMRES
173  LGMRES
174  BiCGSTAB
175  SolveOrPrecondition takes enumerated type Hypre_Chooser, Solver will solve the system, Preconditioner will apply the preconditioner.
176  SetPreconditioner takes a boolean, true means the solver will use the preconditioner.
177  NumFunctions takes an int that describes how many parameters will be passed into Functions. (This needs to be correct.)
178  Functions takes an array of Ref Counted Pointers to an object called FunctionParameter. This class is implemented in Ifpack_Hypre.h.
179  The object takes whether it is Solver or Preconditioner that we are setting a parameter for.
180  The function in Hypre that sets the parameter, and the parameters for that function. An example is below:
181 
182  RCP<FunctionParameter> functs[2];
183  functs[0] = rcp(new FunctionParameter(Solver, &HYPRE_PCGSetMaxIter, 1000)); // max iterations
184  functs[1] = rcp(new FunctionParameter(Solver, &HYPRE_PCGSetTol, 1e-7)); // conv. tolerance
185  list.set("NumFunctions", 2);
186  list.set<RCP<FunctionParameter>*>("Functions", functs);
187  NOTE: SetParameters() must be called to use ApplyInverse(), the solvers will not be created otherwise. An empty list is acceptable to use defaults.
188  */
189  int SetParameters(Teuchos::ParameterList& parameterlist);
190 
192 
200  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, int), int parameter);
201 
203 
211  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, double), double parameter);
212 
214 
223  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, double, int), double parameter1, int parameter2);
224 
226 
235  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, int, double), int parameter1, double parameter2);
236 
238 
247  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, int, int), int parameter1, int parameter2);
248 
250 
258  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, double*), double* parameter);
259 
261 
269  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, int*), int* parameter);
270 
272 
281  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, int**), int** parameter);
282 
284 
293  int SetParameter(Hypre_Chooser chooser, Hypre_Solver Solver);
294 
296 
304  int SetParameter(bool UsePreconditioner){ UsePreconditioner_ = UsePreconditioner; return 0;}
305 
307 
313  int SetParameter(Hypre_Chooser chooser) { SolveOrPrec_ = chooser; return 0;}
314 
316  int SetCoordinates(Teuchos::RCP<Epetra_MultiVector> coords);
317 
319  int SetDiscreteGradient(Teuchos::RCP<const Epetra_CrsMatrix> G);
320 
322  int CallFunctions() const;
323 
325 
334  int SetUseTranspose(bool UseTranspose_in) {UseTranspose_ = UseTranspose_in; return(0);};
335 
336  // @}
337 
338  // @{ Mathematical functions.
339  // Applies the matrix to X, returns the result in Y.
340  int Apply(const Epetra_MultiVector& X,
341  Epetra_MultiVector& Y) const{ return(Multiply(false,X,Y));}
342 
344 
355  int Multiply(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
356 
358 
371  int ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
372 
374  double Condest(const Ifpack_CondestType CT = Ifpack_Cheap,
375  const int MaxIters = 1550,
376  const double Tol = 1e-9,
377  Epetra_RowMatrix* Matrix_in = 0);
378 
380  double Condest() const{ return(Condest_);}
381 
382  // @}
383  // @{ Query methods
384 
386  const char* Label() const {return(Label_);}
387 
389  int SetLabel(const char* Label_in)
390  {
391  strcpy(Label_,Label_in);
392  return(0);
393  }
394 
396  const Epetra_Map& OperatorDomainMap() const{ return *GloballyContiguousRowMap_;}
397 
399  const Epetra_Map& OperatorRangeMap() const{ return *GloballyContiguousRowMap_;}
400 
402  double NormInf() const {return(0.0);};
403 
405  bool HasNormInf() const {return(false);};
406 
408  bool UseTranspose() const {return(UseTranspose_);};
409 
411  const Epetra_Comm & Comm() const{return(A_->Comm());};
412 
414  const Epetra_RowMatrix& Matrix() const{ return(*A_);}
415 
417 #if 0
418  const HYPRE_IJMatrix& HypreMatrix()
419  {
420  if(IsInitialized() == false)
421  Initialize();
422  return(HypreA_);
423  }
424 #endif
425 
427  virtual std::ostream& Print(std::ostream& os) const;
428 
430  virtual int NumInitialize() const{ return(NumInitialize_);}
431 
433  virtual int NumCompute() const{ return(NumCompute_);}
434 
436  virtual int NumApplyInverse() const{ return(NumApplyInverse_);}
437 
439  virtual double InitializeTime() const{ return(InitializeTime_);}
440 
442  virtual double ComputeTime() const{ return(ComputeTime_);}
443 
445  virtual double ApplyInverseTime() const{ return(ApplyInverseTime_);}
446 
448  virtual double InitializeFlops() const{ return(0.0);}
449 
451  virtual double ComputeFlops() const{ return(ComputeFlops_);}
452 
454  virtual double ApplyInverseFlops() const{ return(ApplyInverseFlops_);}
455 
456 private:
457 
458  // @}
459  // @{ Private methods
460 
462  Ifpack_Hypre(const Ifpack_Hypre& RHS) : Time_(RHS.Comm()){}
463 
465  Ifpack_Hypre& operator=(const Ifpack_Hypre& /*RHS*/){ return(*this);}
466 
468  void Destroy();
469 
471  MPI_Comm GetMpiComm() const
472  { return (dynamic_cast<const Epetra_MpiComm*>(&A_->Comm()))->GetMpiComm();}
473 
475 
485  int Solve(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
486 
487 
489  int NumGlobalRows() const {return(A_->NumGlobalRows());};
490 
492  int NumGlobalCols() const {return(A_->NumGlobalCols());};
493 
495  int NumMyRows() const {return(A_->NumMyRows());};
496 
498  int NumMyCols() const {return(A_->NumMyCols());};
499 
501  int SetSolverType(Hypre_Solver solver);
502 
504  int SetPrecondType(Hypre_Solver precond);
505 
507  int CreateSolver();
508 
510  int CreatePrecond();
511 
513  int CopyEpetraToHypre();
514 
516  int AddFunToList(Teuchos::RCP<FunctionParameter> NewFun);
517 
519  int Hypre_BoomerAMGCreate(MPI_Comm comm, HYPRE_Solver *solver);
520 
522  int Hypre_ParaSailsCreate(MPI_Comm comm, HYPRE_Solver *solver);
523 
525  int Hypre_EuclidCreate(MPI_Comm comm, HYPRE_Solver *solver);
526 
528  int Hypre_AMSCreate(MPI_Comm comm, HYPRE_Solver *solver);
529 
531  int Hypre_ParCSRHybridCreate(MPI_Comm comm, HYPRE_Solver *solver);
532 
534  int Hypre_ParCSRPCGCreate(MPI_Comm comm, HYPRE_Solver *solver);
535 
537  int Hypre_ParCSRGMRESCreate(MPI_Comm comm, HYPRE_Solver *solver);
538 
540  int Hypre_ParCSRFlexGMRESCreate(MPI_Comm comm, HYPRE_Solver *solver);
541 
543  int Hypre_ParCSRLGMRESCreate(MPI_Comm comm, HYPRE_Solver *solver);
544 
546  int Hypre_ParCSRBiCGSTABCreate(MPI_Comm comm, HYPRE_Solver *solver);
547 
549  Teuchos::RCP<const Epetra_Map> MakeContiguousColumnMap(Teuchos::RCP<const Epetra_RowMatrix> &Matrix) const;
550  // @}
551  // @{ Internal data
552 
554  Teuchos::RCP<Epetra_RowMatrix> A_;
556  Teuchos::ParameterList List_;
558  bool UseTranspose_;
560  double Condest_;
562  bool IsInitialized_;
564  bool IsComputed_;
566  char Label_[160];
568  int NumInitialize_;
570  int NumCompute_;
572  mutable int NumApplyInverse_;
574  double InitializeTime_;
576  double ComputeTime_;
578  mutable double ApplyInverseTime_;
580  double ComputeFlops_;
582  mutable double ApplyInverseFlops_;
584  mutable Epetra_Time Time_;
585 
587  mutable HYPRE_IJMatrix HypreA_;
589  mutable HYPRE_ParCSRMatrix ParMatrix_;
590 
592  Teuchos::RCP<const Epetra_CrsMatrix> G_;
594  mutable HYPRE_IJMatrix HypreG_;
596  mutable HYPRE_ParCSRMatrix ParMatrixG_;
597 
599  mutable HYPRE_IJVector XHypre_;
601  mutable HYPRE_IJVector YHypre_;
602  mutable HYPRE_ParVector ParX_;
603  mutable HYPRE_ParVector ParY_;
604  mutable Teuchos::RCP<hypre_ParVector> XVec_;
605  mutable Teuchos::RCP<hypre_ParVector> YVec_;
606 
607  Teuchos::RCP<Epetra_MultiVector> Coords_;
608  mutable HYPRE_IJVector xHypre_;
609  mutable HYPRE_IJVector yHypre_;
610  mutable HYPRE_IJVector zHypre_;
611  mutable HYPRE_ParVector xPar_;
612  mutable HYPRE_ParVector yPar_;
613  mutable HYPRE_ParVector zPar_;
614 
616  mutable HYPRE_Solver Solver_;
618  mutable HYPRE_Solver Preconditioner_;
619  // The following are pointers to functions to use the solver and preconditioner.
620  int (Ifpack_Hypre::*SolverCreatePtr_)(MPI_Comm, HYPRE_Solver*);
621  int (*SolverDestroyPtr_)(HYPRE_Solver);
622  int (*SolverSetupPtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
623  int (*SolverSolvePtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
624  int (*SolverPrecondPtr_)(HYPRE_Solver, HYPRE_PtrToParSolverFcn, HYPRE_PtrToParSolverFcn, HYPRE_Solver);
625  int (Ifpack_Hypre::*PrecondCreatePtr_)(MPI_Comm, HYPRE_Solver*);
626  int (*PrecondDestroyPtr_)(HYPRE_Solver);
627  int (*PrecondSetupPtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
628  int (*PrecondSolvePtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
629 
630  bool IsSolverCreated_;
631  bool IsPrecondCreated_;
633  Hypre_Chooser SolveOrPrec_;
635  Teuchos::RCP<const Epetra_Map> GloballyContiguousRowMap_;
636  Teuchos::RCP<const Epetra_Map> GloballyContiguousColMap_;
637  Teuchos::RCP<const Epetra_Map> GloballyContiguousNodeRowMap_;
638  Teuchos::RCP<const Epetra_Map> GloballyContiguousNodeColMap_;
640  int NumFunsToCall_;
642  Hypre_Solver SolverType_;
644  Hypre_Solver PrecondType_;
646  bool UsePreconditioner_;
648  std::vector<Teuchos::RCP<FunctionParameter> > FunsToCall_;
650  bool Dump_;
652  mutable Teuchos::ArrayRCP<double> VectorCache_;
653 
654 };
655 
656 #endif // HAVE_HYPRE
657 #endif /* IFPACK_HYPRE_H */
virtual int NumInitialize() const =0
Returns the number of calls to Initialize().
virtual int SetUseTranspose(bool UseTranspose)=0
virtual double ComputeTime() const =0
Returns the time spent in Compute().
virtual double ComputeFlops() const =0
Returns the number of flops in the computation phase.
virtual double ApplyInverseTime() const =0
Returns the time spent in ApplyInverse().
virtual double ApplyInverseFlops() const =0
Returns the number of flops in the application of the preconditioner.
virtual const Epetra_RowMatrix & Matrix() const =0
Returns a pointer to the matrix to be preconditioned.
virtual bool IsInitialized() const =0
Returns true if the preconditioner has been successfully initialized, false otherwise.
virtual std::ostream & Print(std::ostream &os) const =0
Prints basic information on iostream. This function is used by operator&lt;&lt;.
virtual const Epetra_Map & OperatorDomainMap() const =0
virtual const char * Label() const =0
virtual int Initialize()=0
Computes all it is necessary to initialize the preconditioner.
virtual double InitializeTime() const =0
Returns the time spent in Initialize().
virtual int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0
virtual const Epetra_Map & OperatorRangeMap() const =0
virtual const Epetra_Comm & Comm() const =0
virtual int SetParameters(Teuchos::ParameterList &List)=0
Sets all parameters for the preconditioner.
virtual double Condest() const =0
Returns the computed condition number estimate, or -1.0 if not computed.
virtual bool UseTranspose() const =0
Ifpack_ScalingType enumerable type.
Ifpack_Preconditioner: basic class for preconditioning in Ifpack.
virtual int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0
Applies the preconditioner to vector X, returns the result in Y.
virtual double InitializeFlops() const =0
Returns the number of flops in the initialization phase.
virtual bool HasNormInf() const =0
virtual int NumCompute() const =0
Returns the number of calls to Compute().
virtual double NormInf() const =0
virtual bool IsComputed() const =0
Returns true if the preconditioner has been successfully computed, false otherwise.
virtual int Compute()=0
Computes all it is necessary to apply the preconditioner.
virtual int NumApplyInverse() const =0
Returns the number of calls to ApplyInverse().