PaStiX Handbook  6.3.2
PaStiX

PaStiX users documentation. More...

Modules

 Parameters
 This group describes parameters and constants to control the solver.
 
 Examples
 This group describes how to use the solver with some examples.
 
 Analyze step
 This group describes all the analyse routines that prepare data structures for the factorization and solve steps.
 
 Factorization step
 This group describes all the routines of the internal block CSC and the solver matrix structure before applying the numerical factorization step.
 
 Solve step
 This group describes all the routines to solve the system.
 
 Refinement step
 This group describes all the routines to apply iterative refinement on the system.
 
 Schur complement routines
 

Data Structures

struct  pastix_data_s
 Main PaStiX data structure. More...
 
struct  pastix_rhs_s
 Main PaStiX RHS structure. More...
 

Functions

int pastix_task_analyze (pastix_data_t *pastix_data, const spmatrix_t *spm)
 Perform all the preprocessing steps: ordering, symbolic factorization, reordering, proportionnal mapping, ... More...
 
int pastix_task_numfact (pastix_data_t *pastix_data, spmatrix_t *spm)
 Perform all the numerical factorization steps: fill the internal block CSC and the solver matrix structures, then apply the factorization step. More...
 
int pastix_task_solve (pastix_data_t *pastix_data, pastix_int_t m, pastix_int_t nrhs, void *b, pastix_int_t ldb)
 Solve the given problem. More...
 
int pastix_task_refine (pastix_data_t *pastix_data, pastix_int_t m, pastix_int_t nrhs, void *b, pastix_int_t ldb, void *x, pastix_int_t ldx)
 Perform iterative refinement. More...
 
int pastix_task_solve_and_refine (pastix_data_t *pastix_data, pastix_int_t m, pastix_int_t nrhs, void *b, pastix_int_t ldb, void *x, pastix_int_t ldx)
 Performs solve and iterative refinement without unnecessary permutations. More...
 

Detailed Description

PaStiX users documentation.


Data Type Documentation

◆ pastix_data_s

struct pastix_data_s

Main PaStiX data structure.

This structure holds all informations related to the library and problem instance. It stores information from one step to another.

Warning
This structure should not be modified directly by the user.

Definition at line 67 of file pastixdata.h.

Data Fields
pastix_int_t id

Unique identifier of the pastix instance

pastix_int_t * iparm

Store integer parameters (input/output)

double * dparm

Store floating parameters (input/output)

pastix_int_t steps

Bitmask of the steps performed or not

pastix_scheduler_t sched

Indicates the scheduler used for the factorization step

PASTIX_Comm pastix_comm

PaStiX MPI communicator used for the ordering step

PASTIX_Comm intra_node_comm

PaStiX intra node MPI communicator used for synchronizations

PASTIX_Comm inter_node_comm

PaStiX inter node MPI communicator used for the factorization

int procnbr

Total number of MPI processes

int procnum

Local MPI rank

int intra_node_procnbr

Number of MPI tasks in intra node communicator

int intra_node_procnum

Local MPI rank in intra node communicator

int inter_node_procnbr

Number of MPI tasks in inter node communicator

int inter_node_procnum

Local MPI rank in inter node communicator

isched_t * isched

Internal scheduler structure that is always available

void * parsec

PaRSEC context if available

void * starpu

StarPU context if available

const spmatrix_t * csc

Pointer to the user csc structure used as input

pastix_graph_t * graph

Symmetrized graph of the problem used within ordering and symbolic factorization steps.

pastix_int_t schur_n

Number of entries for the Schur complement

pastix_int_t * schur_list

List of entries for the schur complement

pastix_int_t zeros_n

Number of diagonal entries considered as zeros

pastix_int_t * zeros_list

List of diagonal entries considered as zeros

pastix_order_t * ordemesh

Ordering structure

symbol_matrix_t * symbmtx

Symbol Matrix

pastix_bcsc_t * bcsc

Csc after reordering grouped by cblk

SolverMatrix * solvmatr

Solver informations associated to the matrix problem

SolverMatrix * solvloc

Solver informations associated to the matrix problem - Local

SolverMatrix * solvglob

Solver informations associated to the matrix problem - Global

pastix_model_t * cpu_models

CPU model coefficients for the kernels

pastix_model_t * gpu_models

GPU model coefficients for the kernels

char * dir_global

Unique directory name to store output files

char * dir_local

Unique directory name to store output specific to a MPI process

void * b
void * x0
int * bindtab

Former fields that are no longer used for now

void * schur_tab
pastix_int_t schur_tab_set
int cscInternFilled
int scaling
void * scalerowtab
void * iscalerowtab
void * scalecoltab
void * iscalecoltab
pastix_int_t pastix_id

◆ pastix_rhs_s

struct pastix_rhs_s

Main PaStiX RHS structure.

This structure holds all informations about a set of permuted right-hand-sides to be used in the solve/refinement steps.

Warning
This structure should not be modified directly by the user.

Definition at line 150 of file pastixdata.h.

Data Fields
int8_t allocated

Flag to know if the vector b is allocated internally or not.

pastix_coeftype_t flttype

Floating type of the vector.

pastix_int_t m

Local number of rows in the right hand sides.

pastix_int_t n

Number of columns in the right hand sides.

pastix_int_t ld

Leading dimension of the right hand side matrix.

void * b

Right hand sides of size ldb-by-n.

void ** cblkb

Array to store the temporary buffers associated to fanin/recv.

bvec_handle_comm_t * rhs_comm

Structure which handles the MPI communication (= NULL if PASTIX_WITH_MPI=OFF).

pastix_int_t * Ploc2Pglob

Array containing the local permuted index corresponding to the global permuted index.

Function Documentation

◆ pastix_task_analyze()

int pastix_task_analyze ( pastix_data_t pastix_data,
const spmatrix_t *  spm 
)

Perform all the preprocessing steps: ordering, symbolic factorization, reordering, proportionnal mapping, ...

This function combines the calls to all the preprocessing steps. See pastix_subtask_order() for further details on the ordering step, pastix_subtask_symbfact() for further details on the symbolic factotrization step, pastix_subtask_reordering() for further details on the reordering process, and pastix_subtask_blend() for the final analyse step that performs proportionnal mapin and generates the numerical factorization data structures.

Parameters
[in,out]pastix_dataThe pastix_data structure that describes the solver instance.
[in]spmThe sparse matrix descriptor that describes problem instance.
Return values
PASTIX_SUCCESSon successful exit
PASTIX_ERR_BADPARAMETERif one parameter is incorrect.
PASTIX_ERR_OUTOFMEMORYif one allocation failed.

Definition at line 53 of file pastix_task_analyze.c.

References pastix_data_s::dparm, DPARM_ANALYZE_TIME, pastix_data_s::iparm, IPARM_VERBOSE, PASTIX_ERR_BADPARAMETER, pastix_subtask_blend(), pastix_subtask_order(), pastix_subtask_reordering(), pastix_subtask_symbfact(), PASTIX_SUCCESS, PastixVerboseNot, pastix_data_s::procnum, and pastix_data_s::steps.

◆ pastix_task_numfact()

int pastix_task_numfact ( pastix_data_t pastix_data,
spmatrix_t *  spm 
)

Perform all the numerical factorization steps: fill the internal block CSC and the solver matrix structures, then apply the factorization step.

The user can call the pastix_task_solve() to obtain the solution.

This routine is affected by the following parameters: IPARM_VERBOSE.

Parameters
[in,out]pastix_dataThe pastix_data structure that describes the solver instance. On exit, the internal block CSC is filled with entries from the spm matrix and the solver matrix structure stores the factorization of the given problem.
[in,out]spmThe sparse matrix descriptor that describes problem instance. On exit, the spm structure is freed if IPARM_FREE_CSCUSER is true.
Return values
PASTIX_SUCCESSon successful exit,
PASTIX_ERR_BADPARAMETERif one parameter is incorrect,
PASTIX_ERR_OUTOFMEMORYif one allocation failed.

Definition at line 644 of file pastix_task_sopalin.c.

References pastix_data_s::inter_node_procnum, pastix_data_s::iparm, IPARM_FACTORIZATION, IPARM_VERBOSE, PASTIX_ERR_BADPARAMETER, pastix_int_t, pastix_subtask_bcsc2ctab(), pastix_subtask_sopalin(), pastix_subtask_spm2bcsc(), PASTIX_SUCCESS, PastixVerboseNot, and pastix_data_s::steps.

Referenced by pastix().

◆ pastix_task_solve()

int pastix_task_solve ( pastix_data_t pastix_data,
pastix_int_t  m,
pastix_int_t  nrhs,
void *  b,
pastix_int_t  ldb 
)

Solve the given problem.

This routine is affected by the following parameters: IPARM_VERBOSE, IPARM_FACTORIZATION, IPARM_TRANSPOSE_SOLVE.

Parameters
[in,out]pastix_dataThe pastix_data structure that describes the solver instance.
[in]mThe local number of rows of the right-and-side vectors. If the spm is centralized or replicated, m must be equal to spm->gNexp, otherwise m must be equal to spm->nexp.
[in]nrhsThe number of right-and-side vectors.
[in,out]bThe right-and-side vectors (can be multiple rhs). On exit, the solution is stored in place of the right-hand-side vector.
[in]ldbThe leading dimension of the right-and-side vectors.
Return values
PASTIX_SUCCESSon successful exit,
PASTIX_ERR_BADPARAMETERif one parameter is incorrect.

Definition at line 757 of file pastix_task_solve.c.

References pastix_rhs_s::b, pastix_rhs_s::m, PASTIX_ERR_BADPARAMETER, pastix_int_t, pastix_subtask_applyorder(), pastix_subtask_solve(), PASTIX_SUCCESS, PastixDirBackward, PastixDirForward, pastixRhsFinalize(), pastixRhsInit(), and pastix_data_s::steps.

◆ pastix_task_refine()

int pastix_task_refine ( pastix_data_t pastix_data,
pastix_int_t  m,
pastix_int_t  nrhs,
void *  b,
pastix_int_t  ldb,
void *  x,
pastix_int_t  ldx 
)

Perform iterative refinement.

This routine performs the permutation of x, and b before and after the iterative refinement solution. To prevent extra permuation to happen, see pastix_subtask_refine(). This routine is affected by the following parameters: IPARM_REFINEMENT, DPARM_EPSILON_REFINEMENT

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[in]mThe size of system to solve, and the number of rows of both matrices b and x.
[in]nrhsThe number of right hand side members, and the number of columns of b and x.
[in,out]bThe right hand side matrix of size ldb-by-nrhs. B is noted as inout, as permutation might be performed on the matrix. On exit, the matrix is restored as it was on entry.
[in]ldbThe leading dimension of the matrix b. ldb >= n.
[in,out]xThe matrix x of size ldx-by-nrhs. On entry, the initial guess x0 for the refinement step, that may be the solution returned by the solve step or any other initial guess. On exit, contains the final solution after the iterative refinement.
[in]ldxThe leading dimension of the matrix x. ldx >= n.
Return values
PASTIX_SUCCESSon successful exit,
PASTIX_ERR_BADPARAMETERif one parameter is incorrect,

Definition at line 231 of file pastix_task_refine.c.

References pastix_data_s::iparm, IPARM_SCHUR_SOLV_MODE, pastix_int_t, and pastix_data_s::schur_n.

◆ pastix_task_solve_and_refine()

int pastix_task_solve_and_refine ( pastix_data_t pastix_data,
pastix_int_t  m,
pastix_int_t  nrhs,
void *  b,
pastix_int_t  ldb,
void *  x,
pastix_int_t  ldx 
)

Performs solve and iterative refinement without unnecessary permutations.

This routine performs the permutation of x and b before and after the solve and iterative refinement solution. This prevent the permutation of x at the end of the solve and the permutations of x and b at the beginnning of the refinement when the user wants to perfom both solve and refinement. This routine is affected by the following parameters: IPARM_REFINEMENT, DPARM_EPSILON_REFINEMENT, IPARM_VERBOSE, IPARM_FACTORIZATION, IPARM_TRANSPOSE_SOLVE

Parameters
[in]pastix_dataThe PaStiX data structure that describes the solver instance.
[in]mThe size of system to solve, and the number of rows of both matrices b and x.
[in]nrhsThe number of right hand side members, and the number of columns of b and x.
[in,out]bThe right hand side matrix of size ldb-by-nrhs. B is noted as inout, as permutation might be performed on the matrix. On exit, the matrix is restored as it was on entry.
[in]ldbThe leading dimension of the matrix b. ldb >= n.
[in,out]xThe matrix x of size ldx-by-nrhs. On entry, a copy of b. On exit, contains the final solution after the iterative refinement.
[in]ldxThe leading dimension of the matrix x. ldx >= n.
Return values
PASTIX_SUCCESSon successful exit,
PASTIX_ERR_BADPARAMETERif one parameter is incorrect,

Definition at line 354 of file pastix_task_refine.c.

References pastix_data_s::iparm, IPARM_SCHUR_SOLV_MODE, PASTIX_ERR_BADPARAMETER, pastix_int_t, pastix_data_s::schur_n, and pastix_data_s::steps.