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 #include "Ifpack_ConfigDefs.h"
47 #ifdef HAVE_HYPRE
48 
49 #include "HYPRE_IJ_mv.h"
50 #include "HYPRE_parcsr_ls.h"
51 #include "krylov.h"
52 #include "_hypre_parcsr_mv.h"
53 #include "_hypre_IJ_mv.h"
54 #include "HYPRE_parcsr_mv.h"
55 #include "HYPRE.h"
56 #include "Ifpack_Preconditioner.h"
57 #include "Ifpack_Condest.h"
58 #include "Ifpack_ScalingType.h"
59 #include "Epetra_CompObject.h"
60 #include "Epetra_MultiVector.h"
61 #include "Epetra_Vector.h"
62 #include "Epetra_CrsGraph.h"
63 #include "Epetra_CrsMatrix.h"
64 #include "Epetra_BlockMap.h"
65 #include "Epetra_Map.h"
66 #include "Epetra_Object.h"
67 #include "Epetra_Comm.h"
68 #include "Epetra_CrsMatrix.h"
69 #include "Epetra_Time.h"
70 #include "Teuchos_RefCountPtr.hpp"
71 #include "Epetra_MpiComm.h"
72 
73 #ifndef HYPRE_ENUMS
74 #define HYPRE_ENUMS
75 enum Hypre_Solver{
77  BoomerAMG,
78  ParaSails,
79  Euclid,
80  AMS,
81  Hybrid,
82  PCG,
83  GMRES,
84  FlexGMRES,
85  LGMRES,
86  BiCGSTAB
87 };
88 
90 enum Hypre_Chooser{
91  Solver,
92  Preconditioner
93 };
94 #endif //HYPRE_ENUMS
95 
97 class FunctionParameter{
98  public:
100  FunctionParameter(Hypre_Chooser chooser, int (*funct_name)(HYPRE_Solver, int), int param1) :
101  chooser_(chooser),
102  option_(0),
103  int_func_(funct_name),
104  int_param1_(param1) {}
105 
107  FunctionParameter(Hypre_Chooser chooser, int (*funct_name)(HYPRE_Solver, double), double param1):
108  chooser_(chooser),
109  option_(1),
110  double_func_(funct_name),
111  double_param1_(param1) {}
112 
114  FunctionParameter(Hypre_Chooser chooser, int (*funct_name)(HYPRE_Solver, double, int), double param1, int param2):
115  chooser_(chooser),
116  option_(2),
117  double_int_func_(funct_name),
118  int_param1_(param2),
119  double_param1_(param1) {}
120 
122  FunctionParameter(Hypre_Chooser chooser, int (*funct_name)(HYPRE_Solver, int, int), int param1, int param2):
123  chooser_(chooser),
124  option_(3),
125  int_int_func_(funct_name),
126  int_param1_(param1),
127  int_param2_(param2) {}
128 
130  FunctionParameter(Hypre_Chooser chooser, int (*funct_name)(HYPRE_Solver, int*), int *param1):
131  chooser_(chooser),
132  option_(4),
133  int_star_func_(funct_name),
134  int_star_param_(param1) {}
135 
137  FunctionParameter(Hypre_Chooser chooser, int (*funct_name)(HYPRE_Solver, double*), double* param1):
138  chooser_(chooser),
139  option_(5),
140  double_star_func_(funct_name),
141  double_star_param_(param1) {}
142 
144  FunctionParameter(Hypre_Chooser chooser, int (*funct_name)(HYPRE_Solver, int**), int ** param1):
145  chooser_(chooser),
146  option_(6),
147  int_star_star_func_(funct_name),
148  int_star_star_param_(param1) {}
149 
151  int CallFunction(HYPRE_Solver solver, HYPRE_Solver precond){
152  if(chooser_ == Solver){
153  if(option_ == 0){
154  return int_func_(solver, int_param1_);
155  } else if(option_ == 1){
156  return double_func_(solver, double_param1_);
157  } else if(option_ == 2){
158  return double_int_func_(solver, double_param1_, int_param1_);
159  } else if (option_ == 3){
160  return int_int_func_(solver, int_param1_, int_param2_);
161  } else if (option_ == 4){
162  return int_star_func_(solver, int_star_param_);
163  } else if (option_ == 5){
164  return double_star_func_(solver, double_star_param_);
165  } else {
166  return int_star_star_func_(solver,int_star_star_param_);
167  }
168  } else {
169  if(option_ == 0){
170  return int_func_(precond, int_param1_);
171  } else if(option_ == 1){
172  return double_func_(precond, double_param1_);
173  } else if(option_ == 2){
174  return double_int_func_(precond, double_param1_, int_param1_);
175  } else if(option_ == 3) {
176  return int_int_func_(precond, int_param1_, int_param2_);
177  } else if(option_ == 4) {
178  return int_star_func_(precond, int_star_param_);
179  } else if (option_ == 5){
180  return double_star_func_(precond, double_star_param_);
181  } else {
182  return int_star_star_func_(precond,int_star_star_param_);
183  }
184  }
185  }
186 
187  private:
188  Hypre_Chooser chooser_;
189  int option_;
190  int (*int_func_)(HYPRE_Solver, int);
191  int (*double_func_)(HYPRE_Solver, double);
192  int (*double_int_func_)(HYPRE_Solver, double, int);
193  int (*int_int_func_)(HYPRE_Solver, int, int);
194  int (*int_star_func_)(HYPRE_Solver, int*);
195  int (*double_star_func_)(HYPRE_Solver, double*);
196  int (*int_star_star_func_)(HYPRE_Solver, int **);
197  int int_param1_;
198  int int_param2_;
199  double double_param1_;
200  int *int_star_param_;
201  double *double_star_param_;
202  int ** int_star_star_param_;
203 };
204 
205 namespace Teuchos {
206  class ParameterList;
207 }
208 
210 
215 class Ifpack_Hypre: public Ifpack_Preconditioner {
216 
217 public:
218  // @{ Constructors and destructors.
220  Ifpack_Hypre(Epetra_RowMatrix* A);
221 
223  ~Ifpack_Hypre(){ Destroy();}
224 
225  // @}
226  // @{ Construction methods
227 
229  int Initialize();
230 
232  bool IsInitialized() const{ return(IsInitialized_);}
233 
235 
237  int Compute();
238 
240  bool IsComputed() const{ return(IsComputed_);}
241 
242 
244  /* This method is only available if the Teuchos package is enabled.
245  This method recognizes six parameter names: Solver,
246  Preconditioner, SolveOrPrecondition, SetPreconditioner, NumFunctions and Functions. These names are
247  case sensitive. Solver requires an enumerated parameter of type Hypre_Solver. Preconditioner is similar
248  except requires the type be a preconditioner. The options are listed below:
249  Solvers Preconditioners
250  BoomerAMG BoomerAMG
251  AMS ParaSails
252  Hybrid AMS
253  PCG (Default) Euclid (Default)
254  GMRES
255  FlexGMRES
256  LGMRES
257  BiCGSTAB
258  SolveOrPrecondition takes enumerated type Hypre_Chooser, Solver will solve the system, Preconditioner will apply the preconditioner.
259  SetPreconditioner takes a boolean, true means the solver will use the preconditioner.
260  NumFunctions takes an int that describes how many parameters will be passed into Functions. (This needs to be correct.)
261  Functions takes an array of Ref Counted Pointers to an object called FunctionParameter. This class is implemented in Ifpack_Hypre.h.
262  The object takes whether it is Solver or Preconditioner that we are setting a parameter for.
263  The function in Hypre that sets the parameter, and the parameters for that function. An example is below:
264 
265  RCP<FunctionParameter> functs[2];
266  functs[0] = rcp(new FunctionParameter(Solver, &HYPRE_PCGSetMaxIter, 1000)); // max iterations
267  functs[1] = rcp(new FunctionParameter(Solver, &HYPRE_PCGSetTol, 1e-7)); // conv. tolerance
268  list.set("NumFunctions", 2);
269  list.set<RCP<FunctionParameter>*>("Functions", functs);
270  NOTE: SetParameters() must be called to use ApplyInverse(), the solvers will not be created otherwise. An empty list is acceptable to use defaults.
271  */
272  int SetParameters(Teuchos::ParameterList& parameterlist);
273 
275 
283  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, int), int parameter);
284 
286 
294  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, double), double parameter);
295 
297 
306  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, double, int), double parameter1, int parameter2);
307 
309 
318  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, int, int), int parameter1, int parameter2);
319 
321 
329  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, double*), double* parameter);
330 
332 
340  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, int*), int* parameter);
341 
343 
352  int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, int**), int** parameter);
353 
355 
364  int SetParameter(Hypre_Chooser chooser, Hypre_Solver Solver);
365 
367 
375  int SetParameter(bool UsePreconditioner){ UsePreconditioner_ = UsePreconditioner; return 0;}
376 
378 
384  int SetParameter(Hypre_Chooser chooser) { SolveOrPrec_ = chooser; return 0;}
385 
387  int CallFunctions() const;
388 
390 
399  int SetUseTranspose(bool UseTranspose_in) {UseTranspose_ = UseTranspose_in; return(0);};
400 
401  // @}
402 
403  // @{ Mathematical functions.
404  // Applies the matrix to X, returns the result in Y.
405  int Apply(const Epetra_MultiVector& X,
406  Epetra_MultiVector& Y) const{ return(Multiply(false,X,Y));}
407 
409 
420  int Multiply(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
421 
423 
436  int ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
437 
439  double Condest(const Ifpack_CondestType CT = Ifpack_Cheap,
440  const int MaxIters = 1550,
441  const double Tol = 1e-9,
442  Epetra_RowMatrix* Matrix_in = 0);
443 
445  double Condest() const{ return(Condest_);}
446 
447  // @}
448  // @{ Query methods
449 
451  const char* Label() const {return(Label_);}
452 
454  int SetLabel(const char* Label_in)
455  {
456  strcpy(Label_,Label_in);
457  return(0);
458  }
459 
461  const Epetra_Map& OperatorDomainMap() const{ return *GloballyContiguousRowMap_;}
462 
464  const Epetra_Map& OperatorRangeMap() const{ return *GloballyContiguousRowMap_;}
465 
467  double NormInf() const {return(0.0);};
468 
470  bool HasNormInf() const {return(false);};
471 
473  bool UseTranspose() const {return(UseTranspose_);};
474 
476  const Epetra_Comm & Comm() const{return(A_->Comm());};
477 
479  const Epetra_RowMatrix& Matrix() const{ return(*A_);}
480 
482  const HYPRE_IJMatrix& HypreMatrix()
483  {
484  if(IsInitialized() == false)
485  Initialize();
486  return(HypreA_);
487  }
488 
490  virtual std::ostream& Print(std::ostream& os) const;
491 
493  virtual int NumInitialize() const{ return(NumInitialize_);}
494 
496  virtual int NumCompute() const{ return(NumCompute_);}
497 
499  virtual int NumApplyInverse() const{ return(NumApplyInverse_);}
500 
502  virtual double InitializeTime() const{ return(InitializeTime_);}
503 
505  virtual double ComputeTime() const{ return(ComputeTime_);}
506 
508  virtual double ApplyInverseTime() const{ return(ApplyInverseTime_);}
509 
511  virtual double InitializeFlops() const{ return(0.0);}
512 
514  virtual double ComputeFlops() const{ return(ComputeFlops_);}
515 
517  virtual double ApplyInverseFlops() const{ return(ApplyInverseFlops_);}
518 
519 private:
520 
521  // @}
522  // @{ Private methods
523 
525  Ifpack_Hypre(const Ifpack_Hypre& RHS) : Time_(RHS.Comm()){}
526 
528  Ifpack_Hypre& operator=(const Ifpack_Hypre& /*RHS*/){ return(*this);}
529 
531  void Destroy();
532 
534  MPI_Comm GetMpiComm() const
535  { return (dynamic_cast<const Epetra_MpiComm*>(&A_->Comm()))->GetMpiComm();}
536 
538 
548  int Solve(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
549 
550 
552  int NumGlobalRows() const {return(A_->NumGlobalRows());};
553 
555  int NumGlobalCols() const {return(A_->NumGlobalCols());};
556 
558  int NumMyRows() const {return(A_->NumMyRows());};
559 
561  int NumMyCols() const {return(A_->NumMyCols());};
562 
564  int SetSolverType(Hypre_Solver solver);
565 
567  int SetPrecondType(Hypre_Solver precond);
568 
570  int CreateSolver();
571 
573  int CreatePrecond();
574 
576  int CopyEpetraToHypre();
577 
579  int AddFunToList(Teuchos::RCP<FunctionParameter> NewFun);
580 
582  int Hypre_BoomerAMGCreate(MPI_Comm /*comm*/, HYPRE_Solver *solver)
583  { return HYPRE_BoomerAMGCreate(solver);}
584 
586  int Hypre_ParaSailsCreate(MPI_Comm comm, HYPRE_Solver *solver)
587  { return HYPRE_ParaSailsCreate(comm, solver);}
588 
590  int Hypre_EuclidCreate(MPI_Comm comm, HYPRE_Solver *solver)
591  { return HYPRE_EuclidCreate(comm, solver);}
592 
594  int Hypre_AMSCreate(MPI_Comm /*comm*/, HYPRE_Solver *solver)
595  { return HYPRE_AMSCreate(solver);}
596 
598  int Hypre_ParCSRHybridCreate(MPI_Comm /*comm*/, HYPRE_Solver *solver)
599  { return HYPRE_ParCSRHybridCreate(solver);}
600 
602  int Hypre_ParCSRPCGCreate(MPI_Comm comm, HYPRE_Solver *solver)
603  { return HYPRE_ParCSRPCGCreate(comm, solver);}
604 
606  int Hypre_ParCSRGMRESCreate(MPI_Comm comm, HYPRE_Solver *solver)
607  { return HYPRE_ParCSRGMRESCreate(comm, solver);}
608 
610  int Hypre_ParCSRFlexGMRESCreate(MPI_Comm comm, HYPRE_Solver *solver)
611  { return HYPRE_ParCSRFlexGMRESCreate(comm, solver);}
612 
614  int Hypre_ParCSRLGMRESCreate(MPI_Comm comm, HYPRE_Solver *solver)
615  { return HYPRE_ParCSRLGMRESCreate(comm, solver);}
616 
618  int Hypre_ParCSRBiCGSTABCreate(MPI_Comm comm, HYPRE_Solver *solver)
619  { return HYPRE_ParCSRBiCGSTABCreate(comm, solver);}
620 
621  // @}
622  // @{ Internal data
623 
625  Teuchos::RefCountPtr<Epetra_RowMatrix> A_;
627  Teuchos::ParameterList List_;
629  bool UseTranspose_;
631  double Condest_;
633  bool IsInitialized_;
635  bool IsComputed_;
637  char Label_[160];
639  int NumInitialize_;
641  int NumCompute_;
643  mutable int NumApplyInverse_;
645  double InitializeTime_;
647  double ComputeTime_;
649  mutable double ApplyInverseTime_;
651  double ComputeFlops_;
653  mutable double ApplyInverseFlops_;
655  mutable Epetra_Time Time_;
656 
658  mutable HYPRE_IJMatrix HypreA_;
660  mutable HYPRE_ParCSRMatrix ParMatrix_;
662  mutable HYPRE_IJVector XHypre_;
664  mutable HYPRE_IJVector YHypre_;
665  mutable HYPRE_ParVector ParX_;
666  mutable HYPRE_ParVector ParY_;
667  mutable hypre_ParVector *XVec_;
668  mutable hypre_ParVector *YVec_;
669  mutable hypre_Vector *XLocal_;
670  mutable hypre_Vector *YLocal_;
672  mutable HYPRE_Solver Solver_;
674  mutable HYPRE_Solver Preconditioner_;
675  // The following are pointers to functions to use the solver and preconditioner.
676  int (Ifpack_Hypre::*SolverCreatePtr_)(MPI_Comm, HYPRE_Solver*);
677  int (*SolverDestroyPtr_)(HYPRE_Solver);
678  int (*SolverSetupPtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
679  int (*SolverSolvePtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
680  int (*SolverPrecondPtr_)(HYPRE_Solver, HYPRE_PtrToParSolverFcn, HYPRE_PtrToParSolverFcn, HYPRE_Solver);
681  int (Ifpack_Hypre::*PrecondCreatePtr_)(MPI_Comm, HYPRE_Solver*);
682  int (*PrecondDestroyPtr_)(HYPRE_Solver);
683  int (*PrecondSetupPtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
684  int (*PrecondSolvePtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
685 
686  bool *IsSolverSetup_;
687  bool *IsPrecondSetup_;
689  Hypre_Chooser SolveOrPrec_;
691  Teuchos::RCP<const Epetra_Map> GloballyContiguousRowMap_;
692  Teuchos::RCP<const Epetra_Map> GloballyContiguousColMap_;
694  int NumFunsToCall_;
696  Hypre_Solver SolverType_;
698  Hypre_Solver PrecondType_;
700  bool UsePreconditioner_;
702  std::vector<Teuchos::RCP<FunctionParameter> > FunsToCall_;
703 };
704 
705 #endif // HAVE_HYPRE
706 #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().