| Amesos2 - Direct Sparse Solver Interfaces Version of the Day
    | 
Amesos2 interface to the distributed memory version of SuperLU. More...
#include <Amesos2_Superludist_decl.hpp>


| Public Types | |
| typedef Superludist< Matrix, Vector > | type | 
| typedef SolverCore< Amesos2::Superludist, Matrix, Vector > | super_type | 
| typedef Matrix | matrix_type | 
| typedef Vector | vector_type | 
| typedef super_type::scalar_type | scalar_type | 
| typedef super_type::local_ordinal_type | local_ordinal_type | 
| typedef super_type::global_ordinal_type | global_ordinal_type | 
| typedef super_type::global_size_type | global_size_type | 
| typedef super_type::node_type | node_type | 
| typedef TypeMap< Amesos2::Superludist, scalar_type > | type_map | 
| typedef type_map::type | slu_type | 
| typedef type_map::magnitude_type | magnitude_type | 
| typedef FunctionMap< Amesos2::Superludist, slu_type > | function_map | 
| typedef Kokkos::DefaultHostExecutionSpace | HostExecSpaceType | 
| typedef Kokkos::View< SLUD::int_t *, HostExecSpaceType > | host_size_type_array | 
| typedef Kokkos::View< SLUD::int_t *, HostExecSpaceType > | host_ordinal_type_array | 
| typedef Kokkos::View< slu_type *, HostExecSpaceType > | host_value_type_array | 
| typedef Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > | map_type | 
| typedef Amesos2::Superludist< Matrix, Vector > | solver_type | 
| Public Member Functions | |
| Mathematical functions | |
| super_type & | preOrdering () override | 
| Pre-orders the matrix A for minimal fill-in. | |
| super_type & | symbolicFactorization () override | 
| Performs symbolic factorization on the matrix A. | |
| super_type & | numericFactorization () override | 
| Performs numeric factorization on the matrix A. | |
| void | solve () override | 
| Solves   | |
| void | solve (const Teuchos::Ptr< Vector > X, const Teuchos::Ptr< const Vector > B) const override | 
| Solve  | |
| void | solve (Vector *X, const Vector *B) const override | 
| Solve  | |
| int | solve_ir (const Teuchos::Ptr< Vector > X, const Teuchos::Ptr< const Vector > B, const int maxNumIters, const bool verbose) const | 
| int | solve_ir (Vector *X, const Vector *B, const int maxNumIters, const bool verbose) const | 
| int | solve_ir (const int maxNumIters, const bool verbose) | 
| bool | matrixShapeOK () override | 
| Returns trueif the solver can handle this matrix shape. | |
| void | setA (const Teuchos::RCP< const Matrix > a, EPhase keep_phase=CLEAN) override | 
| Sets the matrix A of this solver. | |
| void | setA (const Matrix *a, EPhase keep_phase=CLEAN) override | 
| Sets the matrix A of this solver. | |
| void | setX (const Teuchos::RCP< Vector > x) override | 
| Sets the LHS vector X. | |
| void | setX (Vector *x) override | 
| Sets the LHS vector X using a raw pointer. | |
| const Teuchos::RCP< Vector > | getX () override | 
| Returns the vector that is the LHS of the linear system. | |
| Vector * | getXRaw () override | 
| Returns a raw pointer to the LHS of the linear system. | |
| void | setB (const Teuchos::RCP< const Vector > b) override | 
| Sets the RHS vector B. | |
| void | setB (const Vector *b) override | 
| Sets the RHS vector B using a raw pointer. | |
| const Teuchos::RCP< const Vector > | getB () override | 
| Returns the vector that is the RHS of the linear system. | |
| const Vector * | getBRaw () override | 
| Returns a raw pointer to the RHS of the linear system. | |
| Parameter methods | |
| super_type & | setParameters (const Teuchos::RCP< Teuchos::ParameterList > ¶meterList) override | 
| Set/update internal variables and solver options. | |
| Teuchos::RCP< const Teuchos::ParameterList > | getValidParameters () const override | 
| Return a const parameter list of all of the valid parameters that this->setParameterList(...) will accept. | |
| void | setParameterList (const Teuchos::RCP< Teuchos::ParameterList > ¶meterList) | 
| Set or update internal variables and solver options. | |
| Teuchos::RCP< Teuchos::ParameterList > | getNonconstParameterList () | 
| This is a empty stub. | |
| Teuchos::RCP< Teuchos::ParameterList > | unsetParameterList () | 
| This is an empty stub. | |
| Static Public Attributes | |
| static const char * | name = "SuperLU_DIST" | 
| Name of this solver interface. | |
| Friends | |
| class | SolverCore< Amesos2::Superludist, Matrix, Vector > | 
| Constructor/Destructor methods | |
| struct Amesos2::Superludist::SLUData | data_ | 
| host_value_type_array | nzvals_view_ | 
| Stores the values of the nonzero entries for SuperLU_DIST. | |
| host_value_type_array | nzvals_temp_ | 
| host_ordinal_type_array | colind_view_ | 
| Stores the row indices of the nonzero entries. | |
| host_size_type_array | rowptr_view_ | 
| Stores the location in Ai_and Aval_ that starts row j. | |
| Teuchos::Array< slu_type > | bvals_ | 
| 1D store for B values | |
| Teuchos::Array< slu_type > | xvals_ | 
| 1D store for X values | |
| bool | in_grid_ | 
| trueif this processor is in SuperLU_DISTS's 2D process grid | |
| bool | same_symbolic_ | 
| bool | force_symbfact_ | 
| bool | same_solve_struct_ | 
| Teuchos::RCP< const map_type > | superlu_rowmap_ | 
| Maps rows of the matrix to processors in the SuperLU_DIST processor grid. | |
| Teuchos::RCP< const map_type > | superlu_contig_rowmap_ | 
| Teuchos::RCP< const map_type > | superlu_contig_colmap_ | 
| bool | is_contiguous_ | 
| Superludist (Teuchos::RCP< const Matrix > A, Teuchos::RCP< Vector > X, Teuchos::RCP< const Vector > B) | |
| Initialize from Teuchos::RCP. | |
| ~Superludist () | |
| Destructor. | |
| void | computeRowPermutationLargeDiagMC64 (SLUD::SuperMatrix &GA) | 
| Compute the row permutation for option LargeDiag-MC64. | |
| int | preOrdering_impl () | 
| Performs pre-ordering on the matrix to increase efficiency. | |
| int | symbolicFactorization_impl () | 
| Perform symbolic factorization of the matrix using SuperLU_DIST. | |
| int | numericFactorization_impl () | 
| SuperLU_DIST specific numeric factorization. | |
| int | solve_impl (const Teuchos::Ptr< MultiVecAdapter< Vector > > X, const Teuchos::Ptr< const MultiVecAdapter< Vector > > B) const | 
| SuperLU_DIST specific solve. | |
| bool | matrixShapeOK_impl () const | 
| Determines whether the shape of the matrix is OK for this solver. | |
| void | setParameters_impl (const Teuchos::RCP< Teuchos::ParameterList > ¶meterList) | 
| Teuchos::RCP< const Teuchos::ParameterList > | getValidParameters_impl () const | 
| void | get_default_grid_size (int nprocs, SLUD::int_t &nprow, SLUD::int_t &npcol) const | 
| bool | loadA_impl (EPhase current_phase) | 
| Reads matrix data into internal solver structures. | |
| Accessor methods | |
| Teuchos::RCP< const Teuchos::Comm< int > > | getComm () const override | 
| Returns a pointer to the Teuchos::Comm communicator with this operator. | |
| Status & | getStatus () const override | 
| Returns a reference to this solver's internal status object. | |
| std::string | description () const override | 
| Returns a short description of this Solver. | |
| void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override | 
| void | printTiming (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override | 
| Prints timing information about the current solver. | |
| void | getTiming (Teuchos::ParameterList &timingParameterList) const override | 
| Extracts timing information from the current solver. | |
| void | loadA (EPhase current_phase) | 
| Refresh this solver's internal data about A. | |
| typedef Kokkos::View< scalar_type *, HostExecSpaceType > | host_scalar_type_array | 
| Teuchos::RCP< const MatrixAdapter< Matrix > > | matrixA_ | 
| The LHS operator. | |
| bool | matrix_loaded_ | 
| Teuchos::RCP< Vector > | multiVecX_ | 
| The LHS vector/multi-vector. | |
| Teuchos::RCP< const Vector > | multiVecB_ | 
| The RHS vector/multi-vector. | |
| global_size_type | globalNumRows_ | 
| Number of global rows in matrixA_. | |
| global_size_type | globalNumCols_ | 
| Number of global columns in matrixA_. | |
| global_size_type | globalNumNonZeros_ | 
| Number of global non-zero values in matrixA_. | |
| global_size_type | rowIndexBase_ | 
| Index base of rowmap of matrixA_. | |
| global_size_type | columnIndexBase_ | 
| Index base of column map of matrixA_. | |
| Status | status_ | 
| Holds status information about a solver. | |
| Control | control_ | 
| Parameters for solving. | |
| Timers | timers_ | 
| Various timing statistics. | |
| int | rank_ | 
| The MPI rank of this image. | |
| bool | root_ | 
| If true, then this is the root processor. | |
| int | nprocs_ | 
| Number of process images in the matrix communicator. | |
| Teuchos::RCP< const map_type > | contig_rowmap_ | 
| Teuchos::RCP< const map_type > | contig_colmap_ | 
| host_ordinal_type_array | perm_g2l | 
| host_ordinal_type_array | recvCountRows | 
| host_ordinal_type_array | recvDisplRows | 
| host_ordinal_type_array | recvCounts | 
| host_ordinal_type_array | recvDispls | 
| host_ordinal_type_array | transpose_map | 
| host_scalar_type_array | nzvals_t | 
| void | setNnzLU (size_t nnz) | 
| Set the number of non-zero values in the   | |
Amesos2 interface to the distributed memory version of SuperLU.
The distributed memory version of SuperLU, SuperLU_DIST, is supported by this Amesos2 interface. Currently support is for the SuperLU_DIST 2.5 version.
This interface to SuperLU_DIST currently does not support row permutations due to a sequential bottleneck present in SuperLU_DIST when calculating such a row permutation (i.e. the matrix must be brought to the root processor, which then finds the row permutation and broadcasts this permutation to the other processors). In the future we may support row permutations through Zoltan. By not supporting this option, we make use of the entirely distributed interface to SuperLU_DIST. On the other hand, if you absolutely need row permutations and your matrix will fit on a single node, then you may consider using Amesos2's SuperLU_MT interface instead.
SuperLU_DIST does not provide a means to directly get the number of non-zeros in the L and U factors.
| Amesos2::Superludist< Matrix, Vector >::Superludist | ( | Teuchos::RCP< const Matrix > | A, | 
| Teuchos::RCP< Vector > | X, | ||
| Teuchos::RCP< const Vector > | B | ||
| ) | 
Initialize from Teuchos::RCP.
References Amesos2::Superludist< Matrix, Vector >::get_default_grid_size(), Amesos2::SolverCore< Amesos2::Superludist, Matrix, Vector >::getComm(), Amesos2::SolverCore< Amesos2::Superludist, Matrix, Vector >::getValidParameters(), Amesos2::SolverCore< Amesos2::Superludist, Matrix, Vector >::globalNumCols_, Amesos2::SolverCore< Amesos2::Superludist, Matrix, Vector >::globalNumRows_, Amesos2::Superludist< Matrix, Vector >::in_grid_, Amesos2::SolverCore< Amesos2::Superludist, Matrix, Vector >::matrixA_, Amesos2::SolverCore< Amesos2::Superludist, Matrix, Vector >::setParameters(), and Amesos2::Superludist< Matrix, Vector >::superlu_rowmap_.
| Amesos2::Superludist< Matrix, Vector >::~Superludist | ( | ) | 
Destructor.
These are created by superlu malloc and should be deallocated by superlu free
| 
 | private | 
Compute the row permutation for option LargeDiag-MC64.
SuperLU_DIST supports several forms of row permutations. Refer to slu_mt_options for the available RowPerm options. 
| 
 | private | 
Performs pre-ordering on the matrix to increase efficiency.
SuperLU_DIST supports several forms of column permutations. Refer to slu_mt_options for the available ColPerm options. 
| 
 | private | 
Perform symbolic factorization of the matrix using SuperLU_DIST.
Called second in the sequence before numericFactorization.
| std::runtime_error | SuperLU_DIST is not able to factor the matrix. | 
| 
 | private | 
SuperLU_DIST specific numeric factorization.
SuperLU_DIST factors the matrix in a shared memory environment using nprocs threads, where nprocs defaults to 1 if it is not changed through setParameters().
| std::runtime_error | SuperLU_DIST is not able to factor the matrix | 
| 
 | private | 
SuperLU_DIST specific solve.
Uses the symbolic and numeric factorizations, along with the RHS vector B to solve the sparse system of equations. The solution is placed in X.
| std::runtime_error | SuperLU_DIST is not able to solve the system. | 
| 
 | private | 
Determines whether the shape of the matrix is OK for this solver.
SuperLU_DIST supports square matrices.
| 
 | private | 
Currently, the following SuperLU_DIST parameters/options are recognized:
"npcol"(int) and "nprow"(int) : Specified together, these parameters set the size of the SuperLU_DIST processor grid to nprow rows by npcol columns. If these parameters are not set, the SuperLU_DIST interface uses a heuristic to pick the grid dimensions based on the number of processors in the matrix' communicator. "ColPerm" which takes one of the following: "NATURAL" : natural column ordering. "PARMETIS" : use the ParMETIS TPL to order the columns. (default) "ReplaceTinyPivot" : { true | false }. Specifies whether to replace tiny diagonals with 
true) | 
 | private | 
Hooked in by Amesos2::Solver parent class.
| 
 | private | 
Calculates a SuperLU_DIST grid size of nprow by npcol processes which will try to utilize all nprocs available processes, but in case of failure, will return a square grid that may not use all nprocs processes.
If you're ever not pleased with how the algorithm's heuristics treat prime numbers, don't give a prime for nprocs.
nprocs , nprow and npcol parameters may be set together directly with setParameters() Referenced by Amesos2::Superludist< Matrix, Vector >::Superludist().
| 
 | private | 
Reads matrix data into internal solver structures.
Loads data from the matrix A into the internal SuperLU_DIST matrix structure. This function requires communication accross all processors as the matrix is redistributed as necessary to the processors in the SuperLU_DIST process grid.
true if the matrix was loaded, false if not References Amesos2::ARBITRARY, Amesos2::DISTRIBUTED_NO_OVERLAP, Amesos2::ROOTED, and Amesos2::SORTED_INDICES.
| 
 | overridevirtualinherited | 
Pre-orders the matrix A for minimal fill-in.
Rearranges the rows and columns of the matrix A to minimize the amount of fill-in of the non-zero entries of the matrix. Pre-ordering may or may not be supported by the underlying solver. If not supported, a call to this method simply does nothing.
this . Implements Amesos2::Solver< Matrix, Vector >.
| 
 | overridevirtualinherited | 
Performs symbolic factorization on the matrix A.
In addition to performing symbolic factorization on the matrix A, the call to symbolicFactorization() implies that no change will be made to the non-zero structure of the underlying matrix without a subsequent call to symbolicFactorization().
this . Implements Amesos2::Solver< Matrix, Vector >.
| 
 | overridevirtualinherited | 
Performs numeric factorization on the matrix A.
In addition to performing numeric factorization on the matrix A, the call to numericFactorization() implies that no change will be made to the underlying matrix values without a subsequent call to numericFactorization().
this Implements Amesos2::Solver< Matrix, Vector >.
| 
 | overridevirtualinherited | 
Solves 



Implements Amesos2::Solver< Matrix, Vector >.
| 
 | overridevirtualinherited | 
Solve 
This overload of solve uses the given X and B vectors when solving. This X and B are used in place of any X and B that were given upon construction of the Amesos2 solver instance and are used only for this solve.
If a permanent change of X and B are required, see the setX() and setB() methods.
X contains the solution to the systemX and B given at construction time (if any) are unchanged. Implements Amesos2::Solver< Matrix, Vector >.
| 
 | overridevirtualinherited | 
Solve 
This overload of solve uses the given X and B vectors when solving. This X and B are used in place of any X and B that were given upon construction of the Amesos2 solver instance and are used only for this solve.
If a permanent change of X and B are required, see the setX() and setB() methods.
X contains the solution to the systemX and B given at construction time (if any) are unchanged. Implements Amesos2::Solver< Matrix, Vector >.
| 
 | overridevirtualinherited | 
Returns true if the solver can handle this matrix shape. 
Returns true if the matrix shape is one that the underlying concrete sparse direct solver can handle. Classes that work only on square matrices should return false for rectangular matrices. Classes that work only on symmetric matrices would return false for non-symmetric matrices. etc.
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | overridevirtualinherited | 
Sets the matrix A of this solver.
| [in] | a | An RCP to a matrix will will be used for future computation steps | 
| [in] | keep_phase | This parameter tells the solver what state it should keep. For example, you may want to replace the matrix but keep the symbolic factorization because you know the structure of the new matrix is the same as the structure of the old matrix. In this case you would pass Amesos2::SYMBFACT as this parameter. | 
The default value for the second parameter is Amesos2::CLEAN, which means that the internal state of the solver will be completely reset. It will be as if no previous computational steps were performed.
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | inlineoverridevirtualinherited | 
Sets the matrix A of this solver.
| [in] | a | An raw C pointer to a matrix will will be used for future computation steps. | 
| [in] | keep_phase | This parameter tells the solver what state it should keep. For example, you may want to replace the matrix but keep the symbolic factorization because you know the structure of the new matrix is the same as the structure of the old matrix. In this case you would pass Amesos2::SYMBFACT as this parameter. | 
The default value for the second parameter is Amesos2::CLEAN, which means that the internal state of the solver will be completely reset. It will be as if no previous computational steps were performed.
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | inlineoverridevirtualinherited | 
Sets the LHS vector X.
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | inlineoverridevirtualinherited | 
Sets the LHS vector X using a raw pointer.
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | inlineoverridevirtualinherited | 
Returns the vector that is the LHS of the linear system.
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | inlineoverridevirtualinherited | 
Returns a raw pointer to the LHS of the linear system.
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | inlineoverridevirtualinherited | 
Sets the RHS vector B.
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | inlineoverridevirtualinherited | 
Sets the RHS vector B using a raw pointer.
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | inlineoverridevirtualinherited | 
Returns the vector that is the RHS of the linear system.
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | inlineoverridevirtualinherited | 
Returns a raw pointer to the RHS of the linear system.
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | overridevirtualinherited | 
Set/update internal variables and solver options.
The setParameters method is consistent over all concrete solvers. It accepts general status and control parameters, as well as parameters specific to a given solver. If the solver does not recognize the parameter, then it will simply be ignored
parameterList. Parameters not specified in parameterList revert to their default values.this Implements Amesos2::Solver< Matrix, Vector >.
| 
 | overridevirtualinherited | 
Return a const parameter list of all of the valid parameters that this->setParameterList(...) will accept.
"YES" and "NO" as well as true and false ). Implements Amesos2::Solver< Matrix, Vector >.
| 
 | inlineinherited | 
Set or update internal variables and solver options.
Redefined from Teuchos::ParameterListAcceptor
setParameters() | [in] | parameterList | 
| 
 | inlineinherited | 
This is a empty stub.
| 
 | inlineinherited | 
This is an empty stub.
| 
 | inlineoverridevirtualinherited | 
Returns a pointer to the Teuchos::Comm communicator with this operator.
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | inlineoverridevirtualinherited | 
Returns a reference to this solver's internal status object.
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | overridevirtualinherited | 
Returns a short description of this Solver.
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | overridevirtualinherited | 
Prints the status information about the current solver with some level of verbosity
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | overridevirtualinherited | 
Prints timing information about the current solver.
The Amesos2::SolverCore base class takes care of tracking total time spent in the Amesos2 interface. Concrete solver interface class are responsible for reporting other timing statistics, which include time spent in:
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | overridevirtualinherited | 
Extracts timing information from the current solver.
Results are placed into the parameter list timingParameterList.
| [out] | timingParameterList | Accepts timing information from the current solver | 
Implements Amesos2::Solver< Matrix, Vector >.
| 
 | privateinherited | 
Refresh this solver's internal data about A.
Called whenever it would be necessary to refresh a solver's internal storage of the matrix A, which is whenever a phase is called that is equal to or below the current call.
For example, say a user has just previously called solve(), then calls numericFactorization(). Since the solve phase is greater than the numericFactorization phase, this is an indication that the internal store of A needs refreshing, since the user (assuming the user know what she's doing) wouldn't otherwise need to call numericFactorization following a solve.
| 
 | inlineprotectedinherited | 
Set the number of non-zero values in the 

Concrete solver classes may call this method if they wish to (or are able to) report the number of conbined non-zero count for the 

| 
 | protectedinherited | 
If true indicates that the current matrix A has been loaded into internal solver structures. 
| 
 | protectedinherited | 
The RHS vector/multi-vector.
We point to a const Vector because Amesos2 should never directly modify B.