EpetraExt Development
Loading...
Searching...
No Matches
EpetraExt_HypreIJMatrix.h
Go to the documentation of this file.
1//@HEADER
2// ***********************************************************************
3//
4// EpetraExt: Epetra Extended - Linear Algebra Services Package
5// Copyright (2011) Sandia Corporation
6//
7// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8// the U.S. Government retains certain rights in this software.
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#ifndef EPETRAEXT_HYPREIJMATRIX_H_
43#define EPETRAEXT_HYPREIJMATRIX_H_
44
45#if defined(EpetraExt_SHOW_DEPRECATED_WARNINGS)
46#ifdef __GNUC__
47#warning "The EpetraExt package is deprecated"
48#endif
49#endif
50
51// Trilinos source files
52#include "Epetra_Object.h"
53#include "Epetra_CompObject.h"
54#include "Epetra_BasicRowMatrix.h"
55#include "Epetra_Map.h"
56#include "Epetra_Import.h"
57#include "Epetra_MpiComm.h"
58
59//Hypre source files
60#include "krylov.h"
61#include "HYPRE_parcsr_ls.h"
62#include "_hypre_parcsr_mv.h"
63#include "HYPRE_parcsr_mv.h"
64#include "HYPRE_IJ_mv.h"
65#include "_hypre_IJ_mv.h"
66#include "HYPRE.h"
67
68class Epetra_Vector;
70class Epetra_Import;
71
73
80#ifndef HYPRE_ENUMS
81#define HYPRE_ENUMS
83
97
103#endif //HYPRE_ENUMS
104
106
107public:
108
110
111
116 EpetraExt_HypreIJMatrix(HYPRE_IJMatrix matrix);
117
121
123
124
126
135 int ExtractMyRowCopy(int MyRow, int Length, int & NumEntries, double *Values, int * Indices) const;
136
138
147 int ExtractMyEntryView(int CurEntry, double *&Value, int &RowIndex, int &ColIndex);
148
150
159 int ExtractMyEntryView(int CurEntry, const double *&Value, int &RowIndex, int &ColIndex) const;
161
163
164
166
174 int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, int), int parameter);
175
177
185 int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, double), double parameter);
186
188
197 int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, double, int), double parameter1, int parameter2);
198
200
209 int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, int, int), int parameter1, int parameter2);
210
212
220 int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, double*), double* parameter);
221
223
231 int SetParameter(Hypre_Chooser chooser, int (*pt2Func)(HYPRE_Solver, int*), int* parameter);
232
234
244 int SetParameter(Hypre_Chooser chooser, Hypre_Solver Solver, bool transpose=false);
245
247
255 int SetParameter(bool UsePreconditioner);
256
258
266
268
270
277 int Multiply(bool TransA, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
278
280
290 int Solve(bool Upper, bool Trans, bool UnitDiagonal, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
291
293
300
301
303
311
313
314
316
324
326
339
341 virtual bool UseTranspose() const {return(false);}
342
344
346
348
355 int NumMyRowEntries(int MyRow, int & NumEntries) const;
356
358 HYPRE_IJMatrix& GetMatrix(){ return Matrix_;}
359
361protected:
362
365
366 // These methods are needed only because the create methods in Hypre sometimes take an MPI_Comm but not always.
367 // They simply call the create solver in the correct way.
369 int Hypre_BoomerAMGCreate(MPI_Comm comm, HYPRE_Solver *solver)
370 { return HYPRE_BoomerAMGCreate(solver);}
371
373 int Hypre_ParaSailsCreate(MPI_Comm comm, HYPRE_Solver *solver)
374 { return HYPRE_ParaSailsCreate(comm, solver);}
375
377 int Hypre_EuclidCreate(MPI_Comm comm, HYPRE_Solver *solver)
378 { return HYPRE_EuclidCreate(comm, solver);}
379
381 int Hypre_AMSCreate(MPI_Comm comm, HYPRE_Solver *solver)
382 { return HYPRE_AMSCreate(solver);}
383
385 int Hypre_ParCSRHybridCreate(MPI_Comm comm, HYPRE_Solver *solver)
386 { return HYPRE_ParCSRHybridCreate(solver);}
387
389 int Hypre_ParCSRPCGCreate(MPI_Comm comm, HYPRE_Solver *solver)
390 { return HYPRE_ParCSRPCGCreate(comm, solver);}
391
393 int Hypre_ParCSRGMRESCreate(MPI_Comm comm, HYPRE_Solver *solver)
394 { return HYPRE_ParCSRGMRESCreate(comm, solver);}
395
397 int Hypre_ParCSRFlexGMRESCreate(MPI_Comm comm, HYPRE_Solver *solver)
398 { return HYPRE_ParCSRFlexGMRESCreate(comm, solver);}
399
401 int Hypre_ParCSRLGMRESCreate(MPI_Comm comm, HYPRE_Solver *solver)
402 { return HYPRE_ParCSRLGMRESCreate(comm, solver);}
403
405 int Hypre_ParCSRBiCGSTABCreate(MPI_Comm comm, HYPRE_Solver *solver)
406 { return HYPRE_ParCSRBiCGSTABCreate(comm, solver);}
407
412 // These two methods setup the solver or preconditioner by calling the pointer.
413 // They are const because they are called from the Solve() routine.
414 // They really aren't const because they change the value of IsSolverSetup_. This is because it should only be called if it isn't setup.
415 int SetupSolver() const;
416 int SetupPrecond() const;
417
418 // Hypre variables
419 mutable HYPRE_IJMatrix Matrix_;
420 mutable HYPRE_ParCSRMatrix ParMatrix_;
421 mutable HYPRE_IJVector X_hypre;
422 mutable HYPRE_IJVector Y_hypre;
423 mutable HYPRE_ParVector par_x;
424 mutable HYPRE_ParVector par_y;
425 mutable hypre_ParVector *x_vec;
426 mutable hypre_ParVector *y_vec;
427 mutable hypre_Vector *x_local;
428 mutable hypre_Vector *y_local;
429 mutable HYPRE_Solver Solver_;
430 mutable HYPRE_Solver Preconditioner_;
431 // The following are pointers to functions to use the solver and preconditioner.
432 int (EpetraExt_HypreIJMatrix::*SolverCreatePtr_)(MPI_Comm, HYPRE_Solver*);
433 int (*SolverDestroyPtr_)(HYPRE_Solver);
434 int (*SolverSetupPtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
435 int (*SolverSolvePtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
436 int (*SolverPrecondPtr_)(HYPRE_Solver, HYPRE_PtrToParSolverFcn, HYPRE_PtrToParSolverFcn, HYPRE_Solver);
437 int (EpetraExt_HypreIJMatrix::*PrecondCreatePtr_)(MPI_Comm, HYPRE_Solver*);
438 int (*PrecondDestroyPtr_)(HYPRE_Solver);
439 int (*PrecondSetupPtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
440 int (*PrecondSolvePtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector);
441
453 mutable int MatType_;
462};
463#endif /* EPETRAEXT_HYPREIJMATRIX_H_ */
Hypre_Chooser
Enumerated type to choose to solve or precondition.
Hypre_Solver
Enumerated type for Hypre solvers.
int SetParameter(Hypre_Chooser chooser, int(*pt2Func)(HYPRE_Solver, int, int), int parameter1, int parameter2)
Set a parameter that takes two int parameters.
virtual bool UseTranspose() const
Returns the current UseTranspose setting.
int Multiply(bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a EpetraExt_HypreIJMatrix multiplied by a Epetra_MultiVector X in Y.
int(* SolverPrecondPtr_)(HYPRE_Solver, HYPRE_PtrToParSolverFcn, HYPRE_PtrToParSolverFcn, HYPRE_Solver)
int Hypre_ParCSRBiCGSTABCreate(MPI_Comm comm, HYPRE_Solver *solver)
BiCGSTAB Create passing function.
bool * IsSolverSetup_
Flag to know if solver needs to be destoyed.
int(* SolverSolvePtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector)
int(* SolverSetupPtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector)
int MyRowStart_
First row on local processor.
int SetParameter(bool UsePreconditioner)
Sets the solver to use the selected preconditioner.
int InitializeDefaults()
Set global variables to default values.
int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y.
int CreatePrecond()
Create the preconditioner with selected type.
int Hypre_EuclidCreate(MPI_Comm comm, HYPRE_Solver *solver)
Euclid Create passing function.
int SetParameter(Hypre_Chooser chooser)
Choose to solve the problem or apply the preconditioner.
bool TransposeSolve_
Do a transpose solve, only BoomerAMG.
int NumGlobalCols_
Number of columns across all processors.
int(* PrecondSolvePtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector)
int(EpetraExt_HypreIJMatrix::* PrecondCreatePtr_)(MPI_Comm, HYPRE_Solver *)
int ExtractMyEntryView(int CurEntry, double *&Value, int &RowIndex, int &ColIndex)
Returns a reference to the ith entry in the matrix, along with its row and column index.
int CreateSolver()
Create the solver with selected type.
int Hypre_BoomerAMGCreate(MPI_Comm comm, HYPRE_Solver *solver)
AMG Create passing function.
int Hypre_ParCSRHybridCreate(MPI_Comm comm, HYPRE_Solver *solver)
Hybrid Create passing function.
int LeftScale(const Epetra_Vector &X)
Scales the EpetraExt_HypreIJMatrix on the left with a Epetra_Vector x.
int SetParameter(Hypre_Chooser chooser, int(*pt2Func)(HYPRE_Solver, int), int parameter)
Set a parameter that takes a single int.
int NumGlobalRows_
Number of rows across all processors.
int Hypre_AMSCreate(MPI_Comm comm, HYPRE_Solver *solver)
AMS Create passing function.
int Hypre_ParCSRFlexGMRESCreate(MPI_Comm comm, HYPRE_Solver *solver)
FlexGMRES Create passing function.
int(* SolverDestroyPtr_)(HYPRE_Solver)
EpetraExt_HypreIJMatrix(HYPRE_IJMatrix matrix)
Epetra_HypreIJMatrix constructor.
int SetParameter(Hypre_Chooser chooser, int(*pt2Func)(HYPRE_Solver, double), double parameter)
Set a parameter that takes a single double.
int(EpetraExt_HypreIJMatrix::* SolverCreatePtr_)(MPI_Comm, HYPRE_Solver *)
int(* PrecondDestroyPtr_)(HYPRE_Solver)
int SetParameter(Hypre_Chooser chooser, Hypre_Solver Solver, bool transpose=false)
Sets the solver that is used by the Solve() and ApplyInverse() methods. Until this is called,...
Hypre_Chooser SolveOrPrec_
Choose to solve or apply preconditioner.
int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y.
int Solve(bool Upper, bool Trans, bool UnitDiagonal, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a EpetraExt_HypreIJMatrix solving a Epetra_MultiVector X in Y.
int SetParameter(Hypre_Chooser chooser, int(*pt2Func)(HYPRE_Solver, double *), double *parameter)
Set a parameter that takes a double*.
virtual ~EpetraExt_HypreIJMatrix()
EpetraExt_HypreIJMatrix Destructor.
int NumMyRows_
Number of rows on local processor.
int Hypre_ParCSRPCGCreate(MPI_Comm comm, HYPRE_Solver *solver)
PCG Create passing function.
int Hypre_ParCSRLGMRESCreate(MPI_Comm comm, HYPRE_Solver *solver)
LGMRES Create passing function.
int SetParameter(Hypre_Chooser chooser, int(*pt2Func)(HYPRE_Solver, int *), int *parameter)
Set a parameter that takes an int*.
int SetParameter(Hypre_Chooser chooser, int(*pt2Func)(HYPRE_Solver, double, int), double parameter1, int parameter2)
Set a parameter that takes a double then an int.
int NumMyRowEntries(int MyRow, int &NumEntries) const
Return the current number of values stored for the specified local row.
int MyRowEnd_
Last row on local processor.
int Hypre_ParaSailsCreate(MPI_Comm comm, HYPRE_Solver *solver)
ParaSails Create passing function.
int Hypre_ParCSRGMRESCreate(MPI_Comm comm, HYPRE_Solver *solver)
GMRES Create passing function.
HYPRE_IJMatrix & GetMatrix()
Return a reference to the Hypre matrix.
int ExtractMyEntryView(int CurEntry, const double *&Value, int &RowIndex, int &ColIndex) const
Returns a const reference to the ith entry in the matrix, along with its row and column index.
bool * IsPrecondSetup_
Flag to know if preconditioner needs to be destroyed.
int ExtractMyRowCopy(int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const
Returns a copy of the specified local row in user-provided arrays.
int(* PrecondSetupPtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector)
int MatType_
Hypre matrix type (parCSR).
int RightScale(const Epetra_Vector &X)
Scales the EpetraExt_HypreIJMatrix on the right with a Epetra_Vector x.
virtual bool UpperTriangular() const