hypre/babel/Interfaces.idl

1999 lines
73 KiB
Plaintext

/*BHEADER**********************************************************************
* Copyright (c) 2008, Lawrence Livermore National Security, LLC.
* Produced at the Lawrence Livermore National Laboratory.
* This file is part of HYPRE. See file COPYRIGHT for details.
*
* HYPRE is free software; you can redistribute it and/or modify it under the
* terms of the GNU Lesser General Public License (as published by the Free
* Software Foundation) version 2.1 dated February 1999.
*
* $Revision$
***********************************************************************EHEADER*/
/*
* This is the HYPRE SIDL interface description file, used with the
* Babel tool to provide multiple language support in hypre. The
* prototype of this work was presented at the 2001 SIAM Parallel
* Processing conference. The original authors were Andy Cleary, Jeff
* Painter, and Cal Ribbens (from HYPRE) and Scott Kohn and Gary
* Kumfert (from Components).
*
* For questions about SIDL and Babel, see the Components project web
* page at http://www.llnl.gov/CASC/components/.
**/
/**
* The bHYPRE package defines interfaces for the HYPRE software package.
**/
package bHYPRE version 1.0.0
{
/*----------- Interface -------------------------------------------------*/
/**
* The purpose of a ProblemDefinition is to:
*
* \begin{itemize}
* \item provide a particular view of how to define a problem
* \item construct and return a {\it problem object}
* \end{itemize}
*
* A {\it problem object} is an intentionally vague term that
* corresponds to any useful object used to define a problem.
* Prime examples are:
*
* \begin{itemize}
* \item a LinearOperator object, i.e., something with a matvec
* \item a MatrixAccess object, i.e., something with a getrow
* \item a Vector, i.e., something with a dot, axpy, ...
* \end{itemize}
*
* Note that {\tt Initialize} and {\tt Assemble} are reserved here
* for defining problem objects through a particular interface.
**/
interface ProblemDefinition
{
/**
* Set the MPI Communicator. DEPRECATED, Use Create()
**/
int SetCommunicator(in MPICommunicator mpi_comm);
/**
* The Destroy function doesn't necessarily destroy anything.
* It is just another name for deleteRef. Thus it decrements the
* object's reference count. The Babel memory management system will
* destroy the object if the reference count goes to zero.
**/
void Destroy();
/**
* Prepare an object for setting coefficient values, whether for
* the first time or subsequently.
**/
int Initialize();
/**
* Finalize the construction of an object before using, either
* for the first time or on subsequent uses. {\tt Initialize}
* and {\tt Assemble} always appear in a matched set, with
* Initialize preceding Assemble. Values can only be set in
* between a call to Initialize and Assemble.
**/
int Assemble();
/**
* Create function, like the builtin __create() except you must specify
* the MPI communicator in an argument; and most classes have more arguments.
* This function is documented here because it will be defined for all
* classes. However, it is declared in the classes, not in this interface.
* static int Create( in MPICommunicator mpi_comm,... );
**/
}
/*----------- Interface -------------------------------------------------*/
/**
* This interface is defined to express the conceptual structure of the object
* system. Derived interfaces and classes have similar functions such as
* SetValues and Print, but the functions are not declared here because the
* function argument lists vary
**/
interface MatrixVectorView extends ProblemDefinition
{}
/*----------- Interface -------------------------------------------------*/
/**
* This interface represents a linear-algebraic conceptual view of a
* linear system. The 'I' and 'J' in the name are meant to be
* mnemonic for the traditional matrix notation A(I,J).
**/
interface IJMatrixView extends MatrixVectorView
{
/**
* Set the local range for a matrix object. Each process owns
* some unique consecutive range of rows, indicated by the
* global row indices {\tt ilower} and {\tt iupper}. The row
* data is required to be such that the value of {\tt ilower} on
* any process $p$ be exactly one more than the value of {\tt
* iupper} on process $p-1$. Note that the first row of the
* global matrix may start with any integer value. In
* particular, one may use zero- or one-based indexing.
*
* For square matrices, {\tt jlower} and {\tt jupper} typically
* should match {\tt ilower} and {\tt iupper}, respectively.
* For rectangular matrices, {\tt jlower} and {\tt jupper}
* should define a partitioning of the columns. This
* partitioning must be used for any vector $v$ that will be
* used in matrix-vector products with the rectangular matrix.
* The matrix data structure may use {\tt jlower} and {\tt
* jupper} to store the diagonal blocks (rectangular in general)
* of the matrix separately from the rest of the matrix.
*
* Collective.
**/
int SetLocalRange(in int ilower,
in int iupper,
in int jlower,
in int jupper);
/**
* Sets values for {\tt nrows} of the matrix. The arrays {\tt
* ncols} and {\tt rows} are of dimension {\tt nrows} and
* contain the number of columns in each row and the row
* indices, respectively. The array {\tt cols} contains the
* column indices for each of the {\tt rows}, and is ordered by
* rows. The data in the {\tt values} array corresponds
* directly to the column entries in {\tt cols}. The last argument
* is the size of the cols and values arrays, i.e. the total number
* of nonzeros being provided, i.e. the sum of all values in ncols.
* This functin erases any previous values at the specified locations and
* replaces them with new ones, or, if there was no value there before,
* inserts a new one.
*
* Not collective.
**/
int SetValues( in int nrows,
in rarray<int,1> ncols(nrows),
in rarray<int,1> rows(nrows),
in rarray<int,1> cols(nnonzeros),
in rarray<double,1> values(nnonzeros),
in int nnonzeros
);
/**
* Adds to values for {\tt nrows} of the matrix. Usage details
* are analogous to {\tt SetValues}. Adds to any previous
* values at the specified locations, or, if there was no value
* there before, inserts a new one.
*
* Not collective.
**/
int AddToValues( in int nrows,
in rarray<int,1> ncols(nrows),
in rarray<int,1> rows(nrows),
in rarray<int,1> cols(nnonzeros),
in rarray<double,1> values(nnonzeros),
in int nnonzeros
);
/**
* Gets range of rows owned by this processor and range of
* column partitioning for this processor.
**/
int GetLocalRange(out int ilower,
out int iupper,
out int jlower,
out int jupper);
/**
* Gets number of nonzeros elements for {\tt nrows} rows
* specified in {\tt rows} and returns them in {\tt ncols},
* which needs to be allocated by the user.
**/
int GetRowCounts( in int nrows,
in rarray<int,1> rows(nrows),
inout rarray<int,1> ncols(nrows) );
/**
* Gets values for {\tt nrows} rows or partial rows of the
* matrix. Usage details are analogous to {\tt SetValues}.
**/
int GetValues( in int nrows,
in rarray<int,1> ncols(nrows),
in rarray<int,1> rows(nrows),
in rarray<int,1> cols(nnonzeros),
inout rarray<double,1> values(nnonzeros),
in int nnonzeros
);
/**
* (Optional) Set the max number of nonzeros to expect in each
* row. The array {\tt sizes} contains estimated sizes for each
* row on this process. The integer nrows is the number of rows in
* the local matrix. This call can significantly improve the
* efficiency of matrix construction, and should always be
* utilized if possible.
*
* Not collective.
**/
int SetRowSizes( in rarray<int,1> sizes(nrows),
in int nrows );
/**
* Print the matrix to file. This is mainly for debugging
* purposes.
**/
int Print(in string filename);
/**
* Read the matrix from file. This is mainly for debugging
* purposes.
**/
int Read(in string filename,
in MPICommunicator comm);
}
/*----------- Interface -------------------------------------------------*/
interface IJVectorView extends MatrixVectorView
{
/**
* Set the local range for a vector object. Each process owns
* some unique consecutive range of vector unknowns, indicated
* by the global indices {\tt jlower} and {\tt jupper}. The
* data is required to be such that the value of {\tt jlower} on
* any process $p$ be exactly one more than the value of {\tt
* jupper} on process $p-1$. Note that the first index of the
* global vector may start with any integer value. In
* particular, one may use zero- or one-based indexing.
*
* Collective.
**/
int SetLocalRange(in int jlower,
in int jupper);
/**
* Sets values in vector. The arrays {\tt values} and {\tt
* indices} are of dimension {\tt nvalues} and contain the
* vector values to be set and the corresponding global vector
* indices, respectively. Erases any previous values at the
* specified locations and replaces them with new ones.
*
* Not collective.
**/
int SetValues( in int nvalues,
in rarray<int,1> indices(nvalues),
in rarray<double,1> values(nvalues) );
/**
* Adds to values in vector. Usage details are analogous to
* {\tt SetValues}.
*
* Not collective.
**/
int AddToValues( in int nvalues,
in rarray<int,1> indices(nvalues),
in rarray<double,1> values(nvalues) );
/**
* Returns range of the part of the vector owned by this
* processor.
**/
int GetLocalRange(out int jlower,
out int jupper);
/**
* Gets values in vector. Usage details are analogous to {\tt
* SetValues}.
*
* Not collective.
**/
int GetValues( in int nvalues,
in rarray<int,1> indices(nvalues),
inout rarray<double,1> values(nvalues) );
/**
* Print the vector to file. This is mainly for debugging
* purposes.
**/
int Print(in string filename);
/**
* Read the vector from file. This is mainly for debugging
* purposes.
**/
int Read(in string filename,
in MPICommunicator comm);
}
/*----------- Interface -------------------------------------------------*/
interface SStructMatrixVectorView extends MatrixVectorView
{
/**
* A semi-structured matrix or vector contains a Struct or IJ matrix
* or vector. GetObject returns it.
* The returned type is a sidl.BaseInterface.
* A cast must be used on the returned object to convert it into a known type.
**/
int GetObject(out sidl.BaseInterface A);
}
interface SStructMatrixView extends SStructMatrixVectorView
{
/**
* Set the matrix graph.
* DEPRECATED Use Create
**/
int SetGraph(in SStructGraph graph);
/**
* Set matrix coefficients index by index.
*
* NOTE: Users are required to set values on all processes that
* own the associated variables. This means that some data will
* be multiply defined.
*
* NOTE: The entries in this routine must all be of the same
* type: either stencil or non-stencil, but not both. Also, if
* they are stencil entries, they must all represent couplings
* to the same variable type (there are no such restrictions for
* non-stencil entries).
*
* If the matrix is complex, then {\tt values} consists of pairs
* of doubles representing the real and imaginary parts of each
* complex value.
*
**/
int SetValues( in int part,
in rarray<int,1> index(dim),
in int dim,
in int var,
in int nentries,
in rarray<int,1> entries(nentries),
in rarray<double,1> values(nentries)
);
/**
* Set matrix coefficients a box at a time.
*
* NOTE: Users are required to set values on all processes that
* own the associated variables. This means that some data will
* be multiply defined.
*
* NOTE: The entries in this routine must all be of the same
* type: either stencil or non-stencil, but not both. Also, if
* they are stencil entries, they must all represent couplings
* to the same variable type (there are no such restrictions for
* non-stencil entries).
*
* If the matrix is complex, then {\tt values} consists of pairs
* of doubles representing the real and imaginary parts of each
* complex value.
**/
int SetBoxValues( in int part,
in rarray<int,1> ilower(dim),
in rarray<int,1> iupper(dim),
in int dim,
in int var,
in int nentries,
in rarray<int,1> entries(nentries),
in rarray<double,1> values(nvalues),
in int nvalues /* = nentries*box_size */
);
/**
* Add to matrix coefficients index by index.
*
* NOTE: Users are required to set values on all processes that
* own the associated variables. This means that some data will
* be multiply defined.
*
* NOTE: The entries in this routine must all be of the same
* type: either stencil or non-stencil, but not both. Also, if
* they are stencil entries, they must all represent couplings
* to the same variable type.
*
* If the matrix is complex, then {\tt values} consists of pairs
* of doubles representing the real and imaginary parts of each
* complex value.
**/
int AddToValues( in int part,
in rarray<int,1> index(dim),
in int dim,
in int var,
in int nentries,
in rarray<int,1> entries(nentries),
in rarray<double,1> values(nentries)
);
/**
* Add to matrix coefficients a box at a time.
*
* NOTE: Users are required to set values on all processes that
* own the associated variables. This means that some data will
* be multiply defined.
*
* NOTE: The entries in this routine must all be of stencil
* type. Also, they must all represent couplings to the same
* variable type.
*
* If the matrix is complex, then {\tt values} consists of pairs
* of doubles representing the real and imaginary parts of each
* complex value.
**/
int AddToBoxValues( in int part,
in rarray<int,1> ilower(dim),
in rarray<int,1> iupper(dim),
in int dim,
in int var,
in int nentries,
in rarray<int,1> entries(nentries),
in rarray<double,1> values(nvalues),
in int nvalues /* = nentries*box_size */
);
/**
* Define symmetry properties for the stencil entries in the
* matrix. The boolean argument {\tt symmetric} is applied to
* stencil entries on part {\tt part} that couple variable {\tt
* var} to variable {\tt to\_var}. A value of -1 may be used
* for {\tt part}, {\tt var}, or {\tt to\_var} to specify
* ``all''. For example, if {\tt part} and {\tt to\_var} are
* set to -1, then the boolean is applied to stencil entries on
* all parts that couple variable {\tt var} to all other
* variables.
*
* By default, matrices are assumed to be nonsymmetric.
* Significant storage savings can be made if the matrix is
* symmetric.
**/
int SetSymmetric(in int part,
in int var,
in int to_var,
in int symmetric);
/**
* Define symmetry properties for all non-stencil matrix
* entries.
**/
int SetNSSymmetric(in int symmetric);
/**
* Set the matrix to be complex.
**/
int SetComplex();
/**
* Print the matrix to file. This is mainly for debugging
* purposes.
**/
int Print(in string filename,
in int all);
}
/*----------- Interface -------------------------------------------------*/
interface SStructVectorView extends SStructMatrixVectorView
{
/**
* Set the vector grid.
**/
int SetGrid(in SStructGrid grid);
/**
* Set vector coefficients index by index.
*
* NOTE: Users are required to set values on all processes that
* own the associated variables. This means that some data will
* be multiply defined.
*
* If the vector is complex, then {\tt value} consists of a pair
* of doubles representing the real and imaginary parts of the
* complex value.
**/
int SetValues( in int part,
in rarray<int,1> index(dim),
in int dim,
in int var,
in double value );
/**
* Set vector coefficients a box at a time.
*
* NOTE: Users are required to set values on all processes that
* own the associated variables. This means that some data will
* be multiply defined.
*
* If the vector is complex, then {\tt values} consists of pairs
* of doubles representing the real and imaginary parts of each
* complex value.
**/
int SetBoxValues( in int part,
in rarray<int,1> ilower(dim),
in rarray<int,1> iupper(dim),
in int dim,
in int var,
in rarray<double,1> values(nvalues),
in int nvalues /* = box size */
);
/**
* Set vector coefficients index by index.
*
* NOTE: Users are required to set values on all processes that
* own the associated variables. This means that some data will
* be multiply defined.
*
* If the vector is complex, then {\tt value} consists of a pair
* of doubles representing the real and imaginary parts of the
* complex value.
**/
int AddToValues( in int part,
in rarray<int,1> index(dim),
in int dim,
in int var,
in double value );
/**
* Set vector coefficients a box at a time.
*
* NOTE: Users are required to set values on all processes that
* own the associated variables. This means that some data will
* be multiply defined.
*
* If the vector is complex, then {\tt values} consists of pairs
* of doubles representing the real and imaginary parts of each
* complex value.
**/
int AddToBoxValues( in int part,
in rarray<int,1> ilower(dim),
in rarray<int,1> iupper(dim),
in int dim,
in int var,
in rarray<double,1> values(nvalues),
in int nvalues /* = box size */
);
/**
* Gather vector data before calling {\tt GetValues}.
**/
int Gather();
/**
* Get vector coefficients index by index.
*
* NOTE: Users may only get values on processes that own the
* associated variables.
*
* If the vector is complex, then {\tt value} consists of a pair
* of doubles representing the real and imaginary parts of the
* complex value.
**/
int GetValues( in int part,
in rarray<int,1> index(dim),
in int dim,
in int var,
out double value );
/**
* Get vector coefficients a box at a time.
*
* NOTE: Users may only get values on processes that own the
* associated variables.
*
* If the vector is complex, then {\tt values} consists of pairs
* of doubles representing the real and imaginary parts of each
* complex value.
**/
int GetBoxValues( in int part,
in rarray<int,1> ilower(dim),
in rarray<int,1> iupper(dim),
in int dim,
in int var,
inout rarray<double,1> values(nvalues),
in int nvalues /* = box size */
);
/**
* Set the vector to be complex.
**/
int SetComplex();
/**
* Print the vector to file. This is mainly for debugging
* purposes.
**/
int Print(in string filename,
in int all);
}
/*----------- Interface -------------------------------------------------*/
interface StructMatrixView extends MatrixVectorView
{
/** Set the grid on which vectors are defined. This and the stencil
determine the matrix structure. */
int SetGrid(in StructGrid grid);
/** Set the stencil. This and the grid determine the matrix structure. */
int SetStencil(in StructStencil stencil);
/** Set matrix values at grid point, given by "index".
You can supply values for one or more positions in the stencil.
"index" is an array of size "dim"; and "stencil_indices" and "values"
are arrays of size "num_stencil_indices".
*/
int SetValues( in rarray<int,1> index(dim),
in int dim,
in int num_stencil_indices,
in rarray<int,1> stencil_indices(num_stencil_indices),
in rarray<double,1> values(num_stencil_indices)
);
/** Set matrix values throughout a box in the grid, specified by its lower
and upper corners. You can supply these values for one or more positions
in the stencil. Thus the total number of matrix values you supply,
"nvalues", is num_stencil_indices x box_size, where box_size is the
number of grid points in the box. The values array should be organized
so all values for a given box point are together (i.e., the stencil
index is the most rapidly varying).
"ilower" and "iupper" are arrays of size "dim", "stencil_indices" is an
array of size "num_stencil_indices", and "values" is an array of size
"nvalues". */
int SetBoxValues( in rarray<int,1> ilower(dim),
in rarray<int,1> iupper(dim),
in int dim,
in int num_stencil_indices,
in rarray<int,1> stencil_indices(num_stencil_indices),
in rarray<double,1> values(nvalues),
in int nvalues /* = num_stencil_indices * box_size */
);
/** Set the number of ghost zones, separately on the lower and upper sides
for each dimension.
"num_ghost" is an array of size "dim2", twice the number of dimensions*/
int SetNumGhost( in rarray<int,1> num_ghost(dim2),
in int dim2 /* = dim * 2 */
);
/** Call SetSymmetric with symmetric=1 to turn on symmetric matrix storage if
available. */
int SetSymmetric(in int symmetric);
/** State which stencil entries are constant over the grid.
Supported options are: (i) none (the default),
(ii) all (stencil_constant_points should include all stencil points)
(iii) all entries but the diagonal. */
int SetConstantEntries(
in int num_stencil_constant_points,
in rarray<int,1> stencil_constant_points(num_stencil_constant_points) );
/** Provide values for matrix coefficients which are constant throughout
the grid, one value for each stencil point.
"stencil_indices" and "values" is each an array of length
"num_stencil_indices" */
int SetConstantValues( in int num_stencil_indices,
in rarray<int,1> stencil_indices(num_stencil_indices),
in rarray<double,1> values(num_stencil_indices) );
}
/*----------- Interface -------------------------------------------------*/
interface StructVectorView extends MatrixVectorView
{
/** Set the grid on which vectors are defined. */
int SetGrid(in StructGrid grid);
/** Set the number of ghost zones, separately on the lower and upper sides
for each dimension.
"num_ghost" is an array of size "dim2", twice the number of dimensions. */
int SetNumGhost( in rarray<int,1> num_ghost(dim2),
in int dim2 /* = dim*2 */
); /* call before Initialize */
/** Set the value of a single vector coefficient, given by "grid_index".
"grid_index" is an array of size "dim", where dim is the number
of dimensions. */
int SetValue( in rarray<int,1> grid_index(dim),
in int dim,
in double value );
/** Set the values of all vector coefficient for grid points in a box.
The box is defined by its lower and upper corners in the grid.
"ilower" and "iupper" are arrays of size "dim", where dim is the
number of dimensions. The "values" array has size "nvalues", which
is the number of grid points in the box. */
int SetBoxValues( in rarray<int,1> ilower(dim),
in rarray<int,1> iupper(dim),
in int dim,
in rarray<double,1> values(nvalues),
in int nvalues /* = box size */
);
}
/*----------- Interface -------------------------------------------------*/
/**
* An Operator is anything that maps one Vector to another. The
* terms {\tt Setup} and {\tt Apply} are reserved for Operators.
* The implementation is allowed to assume that supplied parameter
* arrays will not be destroyed.
**/
interface Operator
{
/**
* Set the MPI Communicator.
* DEPRECATED, use Create:
**/
int SetCommunicator(in MPICommunicator mpi_comm);
/**
* The Destroy function doesn't necessarily destroy anything.
* It is just another name for deleteRef. Thus it decrements the
* object's reference count. The Babel memory management system will
* destroy the object if the reference count goes to zero.
**/
void Destroy();
/**
* Set the int parameter associated with {\tt name}.
**/
int SetIntParameter(in string name,
in int value);
/**
* Set the double parameter associated with {\tt name}.
**/
int SetDoubleParameter(in string name,
in double value);
/**
* Set the string parameter associated with {\tt name}.
**/
int SetStringParameter(in string name,
in string value);
/**
* Set the int 1-D array parameter associated with {\tt name}.
**/
int SetIntArray1Parameter( in string name,
in rarray<int,1> value(nvalues),
in int nvalues );
/**
* Set the int 2-D array parameter associated with {\tt name}.
**/
int SetIntArray2Parameter( in string name,
in array<int,2,column-major> value );
/* rarrays don't work well for 2-D: 2-D Fortran and 2-D C arrays are too different. */
/**
* Set the double 1-D array parameter associated with {\tt name}.
**/
int SetDoubleArray1Parameter( in string name,
in rarray<double,1> value(nvalues),
in int nvalues );
/**
* Set the double 2-D array parameter associated with {\tt name}.
**/
int SetDoubleArray2Parameter( in string name,
in array<double,2,column-major> value );
/* rarrays don't work well for 2-D: 2-D Fortran and 2-D C arrays are too different. */
/**
* Set the int parameter associated with {\tt name}.
**/
int GetIntValue(in string name,
out int value);
/**
* Get the double parameter associated with {\tt name}.
**/
int GetDoubleValue(in string name,
out double value);
/**
* (Optional) Do any preprocessing that may be necessary in
* order to execute {\tt Apply}.
**/
int Setup(in Vector b,
in Vector x);
/**
* Apply the operator to {\tt b}, returning {\tt x}.
**/
int Apply(in Vector b,
inout Vector x);
/**
* Apply the adjoint of the operator to {\tt b}, returning {\tt x}.
**/
int ApplyAdjoint(in Vector b,
inout Vector x);
}
/*----------- Interface -------------------------------------------------*/
interface Vector
{
/**
* Set {\tt self} to 0.
**/
int Clear();
/**
* Copy data from x into {\tt self}.
**/
int Copy(in Vector x);
/**
* Create an {\tt x} compatible with {\tt self}.
* The new vector's data is not specified.
*
* NOTE: When this method is used in an inherited class, the
* cloned {\tt Vector} object can be cast to an object with the
* inherited class type.
**/
int Clone(out Vector x);
/**
* Scale {\tt self} by {\tt a}.
**/
int Scale(in double a);
/**
* Compute {\tt d}, the inner-product of {\tt self} and {\tt x}.
**/
int Dot(in Vector x,
out double d);
/**
* Add {\tt a}*{\tt x} to {\tt self}.
**/
int Axpy(in double a,
in Vector x);
}
/*----------- Interface -------------------------------------------------*/
interface Solver extends Operator
{
/**
* Set the operator for the linear system being solved.
* DEPRECATED. use Create
**/
int SetOperator(in Operator A);
/**
* (Optional) Set the convergence tolerance.
* DEPRECATED. use SetDoubleParameter
**/
int SetTolerance(in double tolerance);
/**
* (Optional) Set maximum number of iterations.
* DEPRECATED use SetIntParameter
**/
int SetMaxIterations(in int max_iterations);
/**
* (Optional) Set the {\it logging level}, specifying the degree
* of additional informational data to be accumulated. Does
* nothing by default (level = 0). Other levels (if any) are
* implementation-specific. Must be called before {\tt Setup}
* and {\tt Apply}.
* DEPRECATED use SetIntParameter
**/
int SetLogging(in int level);
/**
* (Optional) Set the {\it print level}, specifying the degree
* of informational data to be printed either to the screen or
* to a file. Does nothing by default (level=0). Other levels
* (if any) are implementation-specific. Must be called before
* {\tt Setup} and {\tt Apply}.
* DEPRECATED use SetIntParameter
**/
int SetPrintLevel(in int level);
/**
* (Optional) Return the number of iterations taken.
**/
int GetNumIterations(out int num_iterations);
/**
* (Optional) Return the norm of the relative residual.
**/
int GetRelResidualNorm(out double norm);
}
/*----------- Interface -------------------------------------------------*/
interface PreconditionedSolver extends Solver
{
/**
* Set the preconditioner.
**/
int SetPreconditioner(in Solver s);
int GetPreconditioner(out Solver s);
int Clone(out PreconditionedSolver x); /* used in Hybrid solver */
}
/*----------- Interface -------------------------------------------------*/
interface CoefficientAccess
{
/**
* The GetRow method will allocate space for its two output
* arrays on the first call. The space will be reused on
* subsequent calls. Thus the user must not delete them, yet
* must not depend on the data from GetRow to persist beyond the
* next GetRow call.
**/
int GetRow( in int row,
out int size,
out array<int,1,column-major> col_ind,
out array<double,1,column-major> values );
/* rarrays can't be used here unless the user knows the sizes - they can
be inout but not out; and the dimensions must be specified in args */
}
/*----------- Classes ---------------------------------------------------*/
/**
* MPICommunicator class
* - two general Create functions: use CreateC if called from C code,
* CreateF if called from Fortran code.
* - Create_MPICommWorld will create a MPICommunicator to represent
* MPI_Comm_World, and can be called from any language.
*
**/
class MPICommunicator
{
/** Create an MPICommunicator object from C code. */
static MPICommunicator CreateC( in opaque mpi_comm );
/** Create an MPICommunicator object from Fortran code. */
static MPICommunicator CreateF( in opaque mpi_comm );
/** Create an MPICommunicator object which represents MPI_Comm_World. */
static MPICommunicator Create_MPICommWorld();
/**
* The Destroy function doesn't necessarily destroy anything.
* It is just another name for deleteRef. Thus it decrements the
* object's reference count. The Babel memory management system will
* destroy the object if the reference count goes to zero.
**/
void Destroy();
/** Init and Finalize are to help debug MPI interfaces;
* you should normally use the MPI library more directly:
**/
static void Init();
static void Finalize();
}
/**
* The IJParCSR matrix class.
*
* Objects of this type can be cast to IJMatrixView, Operator, or
* CoefficientAccess objects using the {\tt \_\_cast} methods.
**/
class IJParCSRMatrix implements-all IJMatrixView, Operator,
CoefficientAccess
{
/** This function is the preferred way to create an IJParCSR Matrix. */
static IJParCSRMatrix Create( in MPICommunicator mpi_comm,
in int ilower,
in int iupper,
in int jlower,
in int jupper);
/* The following function is for convenience in writing test drivers,
and subject to change without notice. nvalues should be 4 */
static IJParCSRMatrix GenerateLaplacian(
in MPICommunicator mpi_comm, in int nx, in int ny, in int nz,
in int Px, in int Py, in int Pz, in int p, in int q, in int r,
in rarray<double, 1> values(nvalues), in int nvalues,
in int discretization );
/**
* (Optional) Set the max number of nonzeros to expect in each
* row of the diagonal and off-diagonal blocks. The diagonal
* block is the submatrix whose column numbers correspond to
* rows owned by this process, and the off-diagonal block is
* everything else. The arrays {\tt diag\_sizes} and {\tt
* offdiag\_sizes} contain estimated sizes for each row of the
* diagonal and off-diagonal blocks, respectively. This routine
* can significantly improve the efficiency of matrix
* construction, and should always be utilized if possible.
*
* Not collective.
**/
int SetDiagOffdSizes( in rarray<int,1> diag_sizes(local_nrows),
in rarray<int,1> offdiag_sizes(local_nrows),
in int local_nrows );
}
/**
* The IJParCSR vector class.
*
* Objects of this type can be cast to IJVectorView or Vector
* objects using the {\tt \_\_cast} methods.
**/
class IJParCSRVector implements-all IJVectorView, Vector
{
/** This function is the preferred way to create an IJParCSR Vector. */
static IJParCSRVector Create( in MPICommunicator mpi_comm,
in int jlower, in int jupper );
}
/*----------- Classes ---------------------------------------------------*/
/**
* The SStructVariable enumerated type.
*
* An enumerated type that supports cell centered, node
* centered, face centered, and edge centered variables. Face
* centered variables are split into x-face, y-face, and z-face
* variables, and edge centered variables are split into x-edge,
* y-edge, and z-edge variables. The edge centered variable
* types are only used in 3D. In 2D, edge centered variables
* are handled by the face centered types.
*
* Variables are referenced relative to an abstract (cell centered)
* index in the following way:
* \begin{itemize}
* \item cell centered variables are aligned with the index;
* \item node centered variables are aligned with the cell corner
* at relative index (1/2, 1/2, 1/2);
* \item x-face, y-face, and z-face centered variables are aligned
* with the faces at relative indexes (1/2, 0, 0), (0, 1/2, 0),
* and (0, 0, 1/2), respectively;
* \item x-edge, y-edge, and z-edge centered variables are aligned
* with the edges at relative indexes (0, 1/2, 1/2), (1/2, 0, 1/2),
* and (1/2, 1/2, 0), respectively.
* \end{itemize}
*
* The supported identifiers are:
* \begin{itemize}
* \item {\tt HYPRE\_SSTRUCT\_VARIABLE\_CELL}
* \item {\tt HYPRE\_SSTRUCT\_VARIABLE\_NODE}
* \item {\tt HYPRE\_SSTRUCT\_VARIABLE\_XFACE}
* \item {\tt HYPRE\_SSTRUCT\_VARIABLE\_YFACE}
* \item {\tt HYPRE\_SSTRUCT\_VARIABLE\_ZFACE}
* \item {\tt HYPRE\_SSTRUCT\_VARIABLE\_XEDGE}
* \item {\tt HYPRE\_SSTRUCT\_VARIABLE\_YEDGE}
* \item {\tt HYPRE\_SSTRUCT\_VARIABLE\_ZEDGE}
* \end{itemize}
*
* NOTE: Although variables are referenced relative to a unique
* abstract cell-centered index, some variables are associated
* with multiple grid cells. For example, node centered
* variables in 3D are associated with 8 cells (away from
* boundaries). Although grid cells are distributed uniquely to
* different processes, variables may be owned by multiple
* processes because they may be associated with multiple cells.
**/
enum SStructVariable
{
UNDEFINED = -1,
CELL = 0,
NODE = 1,
XFACE = 2,
YFACE = 3,
ZFACE = 4,
XEDGE = 5,
YEDGE = 6,
ZEDGE = 7
};
/**
* The semi-structured grid class.
**/
class SStructGrid
{
/**
* Set the number of dimensions {\tt ndim} and the number of
* structured parts {\tt nparts}.
**/
/** This function is the preferred way to create a SStruct Grid. */
static SStructGrid Create( in MPICommunicator mpi_comm,
in int ndim,
in int nparts );
/* DEPRECATED, use Create: */
int SetNumDimParts(in int ndim,
in int nparts);
/* DEPRECATED, use Create: */
int SetCommunicator(in MPICommunicator mpi_comm);
/**
* The Destroy function doesn't necessarily destroy anything.
* It is just another name for deleteRef. Thus it decrements the
* object's reference count. The Babel memory management system will
* destroy the object if the reference count goes to zero.
**/
void Destroy();
/**
* Set the extents for a box on a structured part of the grid.
**/
int SetExtents( in int part,
in rarray<int,1> ilower(dim),
in rarray<int,1> iupper(dim),
in int dim );
/**
* Describe the variables that live on a structured part of the
* grid. Input: part number, variable number, total number of
* variables on that part (needed for memory allocation),
* variable type.
**/
int SetVariable(in int part,
in int var, in int nvars,
in SStructVariable vartype);
/**
* Describe additional variables that live at a particular
* index. These variables are appended to the array of
* variables set in {\tt SetVariables}, and are referenced as
* such.
**/
int AddVariable( in int part,
in rarray<int,1> index(dim),
in int dim,
in int var,
in SStructVariable vartype);
/**
* Describe how regions just outside of a part relate to other
* parts. This is done a box at a time.
*
* The indexes {\tt ilower} and {\tt iupper} map directly to the
* indexes {\tt nbor\_ilower} and {\tt nbor\_iupper}. Although,
* it is required that indexes increase from {\tt ilower} to
* {\tt iupper}, indexes may increase and/or decrease from {\tt
* nbor\_ilower} to {\tt nbor\_iupper}.
*
* The {\tt index\_map} describes the mapping of indexes 0, 1,
* and 2 on part {\tt part} to the corresponding indexes on part
* {\tt nbor\_part}. For example, triple (1, 2, 0) means that
* indexes 0, 1, and 2 on part {\tt part} map to indexes 1, 2,
* and 0 on part {\tt nbor\_part}, respectively.
*
* NOTE: All parts related to each other via this routine must
* have an identical list of variables and variable types. For
* example, if part 0 has only two variables on it, a cell
* centered variable and a node centered variable, and we
* declare part 1 to be a neighbor of part 0, then part 1 must
* also have only two variables on it, and they must be of type
* cell and node.
**/
int SetNeighborBox( in int part,
in rarray<int,1> ilower(dim),
in rarray<int,1> iupper(dim),
in int nbor_part,
in rarray<int,1> nbor_ilower(dim),
in rarray<int,1> nbor_iupper(dim),
in rarray<int,1> index_map(dim),
in int dim );
/**
* Add an unstructured part to the grid. The variables in the
* unstructured part of the grid are referenced by a global rank
* between 0 and the total number of unstructured variables
* minus one. Each process owns some unique consecutive range
* of variables, defined by {\tt ilower} and {\tt iupper}.
*
* NOTE: This is just a placeholder. This part of the interface
* is not finished.
**/
int AddUnstructuredPart(in int ilower,
in int iupper);
/**
* (Optional) Set periodic for a particular part.
**/
int SetPeriodic( in int part,
in rarray<int,1> periodic(dim),
in int dim );
/**
* Setting ghost in the sgrids.
**/
int SetNumGhost( in rarray<int,1> num_ghost(dim2),
in int dim2 /* expected to be dim*2 */ );
/** final construction of the object before its use */
int Assemble();
}
/*-----------------------*/
/**
* The semi-structured grid stencil class.
**/
class SStructStencil
{
/** This function is the preferred way to create a SStruct Stencil. */
static SStructStencil Create( in int ndim,
in int size );
/**
* The Destroy function doesn't necessarily destroy anything.
* It is just another name for deleteRef. Thus it decrements the
* object's reference count. The Babel memory management system will
* destroy the object if the reference count goes to zero.
**/
void Destroy();
/**
* Set the number of spatial dimensions and stencil entries.
* DEPRECATED, use Create:
**/
int SetNumDimSize(in int ndim,
in int size);
/**
* Set a stencil entry.
**/
int SetEntry( in int entry,
in rarray<int,1> offset(dim),
in int dim,
in int var );
}
/*-----------------------*/
/**
* The semi-structured grid graph class.
**/
class SStructGraph implements-all ProblemDefinition
{
/** This function is the preferred way to create a SStruct Graph. */
static SStructGraph Create(in MPICommunicator mpi_comm,
in SStructGrid grid);
/**
* Set the grid and communicator.
* DEPRECATED, use Create:
**/
int SetCommGrid(in MPICommunicator mpi_comm, in SStructGrid grid);
/**
* Set the stencil for a variable on a structured part of the
* grid.
**/
int SetStencil(in int part,
in int var,
in SStructStencil stencil);
/**
* Add a non-stencil graph entry at a particular index. This
* graph entry is appended to the existing graph entries, and is
* referenced as such.
*
* NOTE: Users are required to set graph entries on all
* processes that own the associated variables. This means that
* some data will be multiply defined.
**/
int AddEntries( in int part,
in rarray<int,1> index(dim),
in int dim,
in int var,
in int to_part,
in rarray<int,1> to_index(dim),
in int to_var );
int SetObjectType( in int type );
}
/*-----------------------*/
/**
* The semi-structured grid matrix class.
*
* Objects of this type can be cast to SStructMatrixView or
* Operator objects using the {\tt \_\_cast} methods.
**/
class SStructMatrix implements-all SStructMatrixView, Operator
{
/** This function is the preferred way to create a SStruct Matrix. */
static SStructMatrix Create( in MPICommunicator mpi_comm,
in SStructGraph graph );
int SetObjectType( in int type ); /* object type has 2 possible values */
}
/*-----------------------*/
/**
* The semi-structured grid vector class.
*
* Objects of this type can be cast to SStructVectorView or Vector
* objects using the {\tt \_\_cast} methods.
**/
class SStructVector implements-all SStructVectorView, Vector
{
/** This function is the preferred way to create a SStruct Vector. */
static SStructVector Create( in MPICommunicator mpi_comm,
in SStructGrid grid );
int SetObjectType( in int type ); /* object type has 2 possible values */
}
/*----------- Classes ---------------------------------------------------*/
/**
* The SStructParCSR matrix class.
*
* Objects of this type can be cast to SStructMatrixView or
* Operator objects using the {\tt \_\_cast} methods.
**/
class SStructParCSRMatrix implements-all SStructMatrixView, Operator
{
/** This function is the preferred way to create a SStruct ParCSR Matrix. */
static SStructParCSRMatrix Create( in MPICommunicator mpi_comm,
in SStructGraph graph );
/* no SetObjectType, as object type has only 1 possible value */
}
/**
* The SStructParCSR vector class.
*
* Objects of this type can be cast to SStructVectorView or Vector
* objects using the {\tt \_\_cast} methods.
**/
class SStructParCSRVector implements-all SStructVectorView, Vector
{
/** This function is the preferred way to create a SStruct ParCSR Vector. */
static SStructParCSRVector Create( in MPICommunicator mpi_comm,
in SStructGrid grid );
/* no SetObjectType, as object type has only 1 possible value */
}
/*----------- Classes ---------------------------------------------------*/
/**
* Define a structured stencil for a structured problem
* description. More than one implementation is not envisioned,
* thus the decision has been made to make this a class rather than
* an interface.
**/
class StructStencil
{
/** This function is the preferred way to create a Struct Stencil.
You provide the number of spatial dimensions and the number of
stencil entries. */
static StructStencil Create( in int ndim,
in int size );
/**
* The Destroy function doesn't necessarily destroy anything.
* It is just another name for deleteRef. Thus it decrements the
* object's reference count. The Babel memory management system will
* destroy the object if the reference count goes to zero.
**/
void Destroy();
/** Set the number of dimensions. DEPRECATED, use StructStencilCreate */
int SetDimension(in int dim);
/** Set the number of stencil entries.
DEPRECATED, use StructStencilCreate */
int SetSize(in int size);
/** Set a stencil element. Specify the stencil index, and an array of
offsets. "offset" is an array of length "dim", the number of spatial
dimensions. */
int SetElement( in int index,
in rarray<int,1> offset(dim),
in int dim );
}
/*-----------------------*/
/**
* Define a structured grid class.
**/
class StructGrid
{
/** This function is the preferred way to create a Struct Grid. */
static StructGrid Create( in MPICommunicator mpi_comm,
in int dim );
/**
* Set the MPI Communicator.
* DEPRECATED, use Create:
**/
int SetCommunicator(in MPICommunicator mpi_comm);
/**
* The Destroy function doesn't necessarily destroy anything.
* It is just another name for deleteRef. Thus it decrements the
* object's reference count. The Babel memory management system will
* destroy the object if the reference count goes to zero.
**/
void Destroy();
/* DEPRECATED, use Create:*/
int SetDimension(in int dim);
/* For comparison, here is some former code using sidl arrays:
int SetExtents(in array<int,1,column-major> ilower,
in array<int,1,column-major> iupper); */
/** Define the lower and upper corners of a box of the grid.
"ilower" and "iupper" are arrays of size "dim", the number of spatial
dimensions. */
int SetExtents( in rarray<int,1> ilower(dim),
in rarray<int,1> iupper(dim),
in int dim );
/** Set the periodicity for the grid. Default is no periodicity.
*
* The argument {\tt periodic} is an {\tt dim}-dimensional integer array that
* contains the periodicity for each dimension. A zero value for a dimension
* means non-periodic, while a nonzero value means periodic and contains the
* actual period. For example, periodicity in the first and third dimensions
* for a 10x11x12 grid is indicated by the array [10,0,12].
*
* NOTE: Some of the solvers in hypre have power-of-two restrictions on the size
* of the periodic dimensions.
*/
int SetPeriodic( in rarray<int,1> periodic(dim),
in int dim );
/** Set the number of ghost zones, separately on the lower and upper sides
for each dimension.
"num_ghost" is an array of size "dim2", twice the number of dimensions. */
int SetNumGhost( in rarray<int,1> num_ghost(dim2),
in int dim2 );
/** final construction of the object before its use */
int Assemble();
}
/*-----------------------*/
/**
* A single class that implements both a view interface and an
* operator interface.
* A StructMatrix is a matrix on a structured grid.
* One function unique to a StructMatrix is SetConstantEntries.
* This declares that matrix entries corresponding to certain stencil points
* (supplied as stencil element indices) will be constant throughout the grid.
**/
class StructMatrix implements-all StructMatrixView, Operator
{
/** This function is the preferred way to create a Struct Matrix. */
static StructMatrix Create( in MPICommunicator mpi_comm,
in StructGrid grid,
in StructStencil stencil );
}
/*-----------------------*/
class StructVector implements-all StructVectorView, Vector
{
/** This function is the preferred way to create a Struct Vector. */
static StructVector Create( in MPICommunicator mpi_comm,
in StructGrid grid );
}
/*----------- Classes ---------------------------------------------------*/
/**
* Identity solver, just solves an identity matrix, for when you don't really
* want a preconditioner
*
* Objects of this type can be cast to Solver objects using the
* {\tt \_\_cast} methods.
**/
class IdentitySolver implements-all Solver
{
/** This function is the preferred way to create an Identity (null) solver. */
static IdentitySolver Create( in MPICommunicator mpi_comm );
}
/**
* Diagonal scaling preconditioner for ParCSR matrix class.
*
* Objects of this type can be cast to Solver objects using the
* {\tt \_\_cast} methods.
**/
class ParCSRDiagScale implements-all Solver
{
/** This function is the preferred way to create a ParCSR DiagScale solver. */
static ParCSRDiagScale Create( in MPICommunicator mpi_comm,
in IJParCSRMatrix A );
}
/**
* Diagonal scaling preconditioner for STruct matrix class.
*
* Objects of this type can be cast to Solver objects using the
* {\tt \_\_cast} methods.
**/
class StructDiagScale implements-all Solver
{
/** This function is the preferred way to create a Struct DiagScale solver. */
static StructDiagScale Create( in MPICommunicator mpi_comm,
in StructMatrix A );
}
class SStructDiagScale implements-all Solver
{
/** This function is the preferred way to create a SStruct DiagScale solver. */
static SStructDiagScale Create( in MPICommunicator mpi_comm,
in Operator A );
/* ... A should also be a SStructMatrixView */
}
/* note: the fine grid options Cycle0NumSweeps, Cycle0RelaxType have been
removed from BoomerAMG as of October 2006 */
/**
* Algebraic multigrid solver, based on classical Ruge-Stueben.
*
* BoomerAMG requires an IJParCSR matrix
*
* The following optional parameters are available and may be set
* using the appropriate {\tt Parameter} function (as indicated in
* parentheses):
*
* \begin{description}
*
* \item[MaxLevels] ({\tt Int}) - maximum number of multigrid
* levels.
*
* \item[StrongThreshold] ({\tt Double}) - AMG strength threshold.
*
* \item[MaxRowSum] ({\tt Double}) -
*
* \item[CoarsenType] ({\tt Int}) - type of parallel coarsening
* algorithm used.
*
* \item[MeasureType] ({\tt Int}) - type of measure used; local or
* global.
*
* \item[CycleType] ({\tt Int}) - type of cycle used; a V-cycle
* (default) or a W-cycle.
*
* \item[NumGridSweeps] ({\tt IntArray 1D}) - number of sweeps for
* fine and coarse grid, up and down cycle. DEPRECATED:
* Use NumSweeps or Cycle?NumSweeps instead.
*
* \item[NumSweeps] ({\tt Int}) - number of sweeps for fine grid, up and
* down cycle.
*
* \item[Cycle1NumSweeps] ({\tt Int}) - number of sweeps for down cycle
*
* \item[Cycle2NumSweeps] ({\tt Int}) - number of sweeps for up cycle
*
* \item[Cycle3NumSweeps] ({\tt Int}) - number of sweeps for coarse grid
*
* \item[GridRelaxType] ({\tt IntArray 1D}) - type of smoother used on
* fine and coarse grid, up and down cycle. DEPRECATED:
* Use RelaxType or Cycle?RelaxType instead.
*
* \item[RelaxType] ({\tt Int}) - type of smoother for fine grid, up and
* down cycle.
*
* \item[Cycle1RelaxType] ({\tt Int}) - type of smoother for down cycle
*
* \item[Cycle2RelaxType] ({\tt Int}) - type of smoother for up cycle
*
* \item[Cycle3RelaxType] ({\tt Int}) - type of smoother for coarse grid
*
* \item[GridRelaxPoints] ({\tt IntArray 2D}) - point ordering used in
* relaxation. DEPRECATED.
*
* \item[RelaxWeight] ({\tt DoubleArray 1D}) - relaxation weight for
* smoothed Jacobi and hybrid SOR. DEPRECATED:
* Instead, use the RelaxWt parameter and the SetLevelRelaxWt function.
*
* \item[RelaxWt] ({\tt Int}) - relaxation weight for all levels for
* smoothed Jacobi and hybrid SOR.
*
* \item[TruncFactor] ({\tt Double}) - truncation factor for
* interpolation.
*
* \item[JacobiTruncThreshold] ({\tt Double}) - threshold for truncation
* of Jacobi interpolation.
*
* \item[SmoothType] ({\tt Int}) - more complex smoothers.
*
* \item[SmoothNumLevels] ({\tt Int}) - number of levels for more
* complex smoothers.
*
* \item[SmoothNumSweeps] ({\tt Int}) - number of sweeps for more
* complex smoothers.
*
* \item[PrintFileName] ({\tt String}) - name of file printed to in
* association with {\tt SetPrintLevel}.
*
* \item[NumFunctions] ({\tt Int}) - size of the system of PDEs
* (when using the systems version).
*
* \item[DOFFunc] ({\tt IntArray 1D}) - mapping that assigns the
* function to each variable (when using the systems version).
*
* \item[Variant] ({\tt Int}) - variant of Schwarz used.
*
* \item[Overlap] ({\tt Int}) - overlap for Schwarz.
*
* \item[DomainType] ({\tt Int}) - type of domain used for Schwarz.
*
* \item[SchwarzRlxWeight] ({\tt Double}) - the smoothing parameter
* for additive Schwarz.
*
* \item[SchwarzUseNonSymm] ({\tt Int}) - defines whether to use a nonsymmetric
* Schwarz smoother. Default:0 (symmetric smoother)
*
* \item[EuLevel] ({\tt Int}) - defines number of levels for ILU(k) smoother.
* To be used with SmoothType 9 and SmoothNumLevels > 0. Default:0
*
* \item[EuSparseA] ({\tt Double}) - defines drop tolerance for ILU(k) smoother
* To be used with SmoothType 9 and SmoothNumLevels > 0. Default:0
*
* \item[EuBJ] ({\tt Int}) - defines use of block Jacobi ILUT smoother.
* To be used with SmoothType 9 and SmoothNumLevels > 0. Default:0
*
* \item[Tolerance] ({\tt Double}) - convergence tolerance, if this
* is used as a solver; ignored if this is used as a preconditioner
*
* \item[DebugFlag] ({\tt Int}) -
*
* \item[InterpType] ({\tt Int}) - Defines which parallel interpolation
* operator is used. There are the following options for interp\_type:
*
* \begin{tabular}{|c|l|} \hline
* 0 & classical modified interpolation \\
* 1 & LS interpolation (for use with GSMG) \\
* 2 & classical modified interpolation for hyperbolic PDEs \\
* 3 & direct interpolation (with separation of weights) \\
* 4 & multipass interpolation \\
* 5 & multipass interpolation (with separation of weights) \\
* 6 & extended classical modified interpolation \\
* 7 & extended (if no common C neighbor) classical modified interpolation \\
* 8 & standard interpolation \\
* 9 & standard interpolation (with separation of weights) \\
* 10 & classical block interpolation (for use with nodal systems version only) \\
* 11 & classical block interpolation (for use with nodal systems version only) \\
* & with diagonalized diagonal blocks \\
* 12 & FF interpolation \\
* 13 & FF1 interpolation \\
* \hline
* \end{tabular}
*
* The default is 0.
*
* \item[PMaxElmts] ({\tt Int}) - Defines the maximal number of nonzero entries
* allowed in interpolation. Default: 0 (i.e. no limit)
*
* \item[AggNumLevels] ({\tt Int}) - Defines the number of levels
* of aggressive coarsening performed. Default:0
*
* \item[NumPaths] ({\tt Int}) - Defines the number of paths used for
* aggressive coarsening. Default:1
*
* \item[NumSamples] ({\tt Int}) - Defines the number of sample vectors used
* in GSMG or LS interpolation.
*
* \item[MaxIterations] ({\tt Int}) - maximum number of iterations
*
* \item[Logging] ({\tt Int}) - Set the {\it logging level}, specifying the
* degree of additional informational data to be accumulated. Does
* nothing by default (level = 0). Other levels (if any) are
* implementation-specific. Must be called before {\tt Setup}
* and {\tt Apply}.
*
* \item[PrintLevel] ({\tt Int}) - Set the {\it print level}, specifying the
* degree of informational data to be printed either to the screen or
* to a file. Does nothing by default (level=0). Other levels
* (if any) are implementation-specific. Must be called before
* {\tt Setup} and {\tt Apply}.
*
* \end{description}
*
* The following function is specific to this class:
*
* \begin{description}
*
* \item[SetLevelRelxWeight] ({\tt Double , \tt Int}) -
* relaxation weight for one specified level of smoothed Jacobi and hybrid SOR.
*
* \end{description}
*
* Objects of this type can be cast to Solver objects using the
* {\tt \_\_cast} methods.
**/
class BoomerAMG implements-all Solver
{
/** This function is the preferred way to create a BoomerAMG solver. */
static BoomerAMG Create( in MPICommunicator mpi_comm,
in IJParCSRMatrix A );
int SetLevelRelaxWt( in double relax_wt, in int level );
/* The following function is just a tool for test drivers,
subject to change without notice... */
int InitGridRelaxation(
out array<int,1,column-major> num_grid_sweeps,
out array<int,1,column-major> grid_relax_type,
out array<int,2,column-major> grid_relax_points,
in int coarsen_type,
out array<double,1,column-major> relax_weights,
in int max_levels );
}
/* RDF: Documentation goes here...*/
/**
* Objects of this type can be cast to Solver objects using the
* {\tt \_\_cast} methods.
*
* Although the usual Solver Set*Parameter functions are available,
* a Euclid-stype parameter-setting function is also available, SetParameters.
*
*
**/
class Euclid implements-all Solver
{
/** This function is the preferred way to create a Euclid solver. */
static Euclid Create( in MPICommunicator mpi_comm, in IJParCSRMatrix A );
int SetParameters( in int argc, inout string argv );
/* ... argv is really in, but inout makes implementation simpler */
}
/* RDF: Documentation goes here...*/
/**
* Objects of this type can be cast to Solver objects using the
* {\tt \_\_cast} methods.
*
* Pilut has not been implemented yet.
*
**/
class Pilut implements-all Solver
{
/** This function is the preferred way to create a Pilut solver. */
static Pilut Create( in MPICommunicator mpi_comm,
in Operator A );
}
/* RDF: Documentation goes here...*/
/**
* Objects of this type can be cast to Solver objects using the
* {\tt \_\_cast} methods.
*
* ParaSails requires an IJParCSR matrix
*
**/
class ParaSails implements-all Solver
{
/** This function is the preferred way to create a ParaSails solver. */
static ParaSails Create( in MPICommunicator mpi_comm,
in IJParCSRMatrix A );
}
/* RDF: Documentation goes here...*/
/**
* Objects of this type can be cast to Solver objects using the
* {\tt \_\_cast} methods.
*
* Schwarz requires an IJParCSR matrix
*
**/
class Schwarz implements-all Solver
{
/** This function is the preferred way to create a Schwarz solver. */
static Schwarz Create( in IJParCSRMatrix A );
}
/**
* Objects of this type can be cast to PreconditionedSolver objects
* using the {\tt \_\_cast} methods.
*
* The regular PCG solver calls Babel-interface matrix and vector functions.
* The HPCG solver calls HYPRE interface functions.
* The regular solver will work with any consistent matrix, vector, and
* preconditioner classes. The HPCG solver will work with the more common ones.
*
* The HPCG solver checks whether the matrix, vectors, and preconditioner
* are of known types, and will not work with any other types.
* Presently, the recognized data types are:
* matrix, vector: IJParCSRMatrix, IJParCSRVector
* matrix, vector: StructMatrix, StructVector
* preconditioner: BoomerAMG, ParaSails, ParCSRDiagScale, IdentitySolver
* preconditioner: StructSMG, StructPFMG
*
*
**/
class HPCG implements-all PreconditionedSolver
{
/** This function is the preferred way to create a HPCG solver. */
static HPCG Create( in MPICommunicator mpi_comm );
}
/**
* PCG solver.
* This calls Babel-interface matrix and vector functions, so it will work
* with any consistent matrix, vector, and preconditioner classes.
**/
class PCG implements-all PreconditionedSolver
{
/** This function is the preferred way to create a PCG solver. */
static PCG Create( in MPICommunicator mpi_comm, in Operator A );
}
/**
* GMRES solver.
* This calls Babel-interface matrix and vector functions, so it will work
* with any consistent matrix, vector, and preconditioner classes.
**/
class GMRES implements-all PreconditionedSolver
{
/** This function is the preferred way to create a GMRES solver. */
static GMRES Create( in MPICommunicator mpi_comm, in Operator A );
}
/**
* Objects of this type can be cast to PreconditionedSolver objects
* using the {\tt \_\_cast} methods.
*
* The regular GMRES solver calls Babel-interface matrix and vector functions.
* The HGMRES solver calls HYPRE interface functions.
* The regular solver will work with any consistent matrix, vector, and
* preconditioner classes. The HGMRES solver will work with the more common ones.
*
* The HGMRES solver checks whether the matrix, vectors, and preconditioner
* are of known types, and will not work with any other types.
* Presently, the recognized data types are:
* matrix, vector: IJParCSRMatrix, IJParCSRVector
* preconditioner: BoomerAMG, ParCSRDiagScale
*
**/
class HGMRES implements-all PreconditionedSolver
{
/** This function is the preferred way to create a HGMRES solver. */
static HGMRES Create( in MPICommunicator mpi_comm );
}
/**
* CGNR solver.
* This calls Babel-interface matrix and vector functions, so it will work
* with any consistent matrix, vector, and preconditioner classes.
**/
class CGNR implements-all PreconditionedSolver
{
/** This function is the preferred way to create a CGNR solver. */
static CGNR Create( in MPICommunicator mpi_comm, in Operator A );
}
/**
* BiCGSTAB solver.
* This calls Babel-interface matrix and vector functions, so it will work
* with any consistent matrix, vector, and preconditioner classes.
**/
class BiCGSTAB implements-all PreconditionedSolver
{
/** This function is the preferred way to create a BiCGSTAB solver. */
static BiCGSTAB Create( in MPICommunicator mpi_comm, in Operator A );
}
/* RDF: Documentation goes here.*/
/**
* Objects of this type can be cast to Solver objects
* using the {\tt \_\_cast} methods.
*
* The StructSMG solver requires a Struct matrix.
*
**/
class StructSMG implements-all Solver
{
/** This function is the preferred way to create a Struct SMG solver. */
static StructSMG Create( in MPICommunicator mpi_comm, in StructMatrix A );
}
/* RDF: Documentation goes here.*/
/**
* Objects of this type can be cast to Solver objects
* using the {\tt \_\_cast} methods.
*
* The StructPFMG solver requires a Struct matrix.
*
**/
class StructPFMG implements-all Solver
{
/** This function is the preferred way to create a Struct PFMG solver. */
static StructPFMG Create( in MPICommunicator mpi_comm,
in StructMatrix A );
}
/* RDF: Documentation goes here. */
/**
* Objects of this type can be cast to Solver objects
* using the {\tt \_\_cast} methods.
*
* The StructJacobi solver requires a Struct matrix.
*
**/
class StructJacobi implements-all Solver
{
/** This function is the preferred way to create a Struct Jacobi solver. */
static StructJacobi Create( in MPICommunicator mpi_comm,
in StructMatrix A );
}
/* Documentation goes here. */
/**
*
* The SStructSplit solver requires a SStruct matrix.
*
**/
class SStructSplit implements-all Solver
{
/** This function is the preferred way to create a SStruct Split solver. */
static SStructSplit Create( in MPICommunicator mpi_comm,
in Operator A );
}
/**
* Hybrid solver
* first tries to solve with the specified Krylov solver, preconditioned by
* If that fails to converge, it will try again with the user-specified
*
* Specify the preconditioner by calling SecondSolver's SetPreconditioner
* method. If no preconditioner is specified (equivalently, if the
* preconditioner for SecondSolver is IdentitySolver), the preconditioner for
* the second try will be one of the following defaults.
* StructMatrix: SMG. other matrix types: not implemented
*
* The Hybrid solver's Setup method will call Setup on KrylovSolver, so the
* user should not call Setup on KrylovSolver.
*
**/
class Hybrid implements-all Solver
{
/** This function is the preferred way to create a Hybrid solver. */
static Hybrid Create( in MPICommunicator mpi_comm,
in PreconditionedSolver SecondSolver,
in Operator A );
int GetFirstSolver( out PreconditionedSolver FirstSolver );
int GetSecondSolver( out PreconditionedSolver SecondSolver );
}
/**
* The ErrorCode enumerated type is used with methods of the ErrorHandler class.
**/
enum ErrorCode
{
HYPRE_ERROR_GENERIC = 1, /* generic error */
HYPRE_ERROR_MEMORY = 2, /* unable to allocate memory */
HYPRE_ERROR_ARG = 4, /* argument error */
/* bits 4-8 are reserved for the index of the argument error */
HYPRE_ERROR_CONV = 256 /* method did not converge as expected */
};
/**
* ErrorHandler class is an interface to the hypre error handling system.
* Its methods help interpret the error flag ierr returned by hypre functions.
**/
class ErrorHandler
{
/**
* The Check method will return nonzero when the error flag ierr
* includes an error of type error\_code; and zero otherwise.
**/
static int Check( in int ierr, in ErrorCode error_code );
/**
* included in the error flag ierr.
**/
static void Describe(in int ierr, out string message );
}
}