PaStiX Handbook
6.4.0
|
#include "pastix_lowrank.h"
Go to the source code of this file.
Data Structures | |
struct | args_solve_s |
Arguments for the solve. More... | |
struct | solver_blok_recv_s |
Solver recv block structure. More... | |
struct | solver_cblk_recv_s |
Solver recv column block structure. More... | |
struct | task_s |
The task structure for the numerical factorization. More... | |
struct | solver_blok_s |
Solver block structure. More... | |
struct | solver_cblk_s |
Solver column block structure. More... | |
struct | solver_matrix_s |
Solver column block structure. More... | |
Typedefs | |
typedef enum solve_step_e | solve_step_t |
Tags used in MPI communications. | |
typedef struct args_solve_s | args_solve_t |
Arguments for the solve. | |
typedef struct solver_blok_recv_s | solver_blok_recv_t |
Solver recv block structure. | |
typedef struct solver_cblk_recv_s | solver_cblk_recv_t |
Solver recv column block structure. More... | |
typedef struct task_s | Task |
The task structure for the numerical factorization. | |
typedef struct solver_blok_s | SolverBlok |
Solver block structure. | |
typedef struct solver_cblk_s | SolverCblk |
Solver column block structure. | |
Enumerations | |
enum | solve_step_e { PastixSolveForward , PastixSolveBackward , PastixFacto } |
Tags used in MPI communications. | |
Functions | |
static solve_step_t | compute_solve_step (pastix_side_t side, pastix_uplo_t uplo, pastix_trans_t trans) |
Computes the current solve step. More... | |
static pastix_int_t | cblk_colnbr (const SolverCblk *cblk) |
Compute the number of columns in a column block. More... | |
static void * | cblk_getdataL (const SolverCblk *cblk) |
Get the pointer to the data associated to the lower part of the cblk. More... | |
static void * | cblk_getdataU (const SolverCblk *cblk) |
Get the pointer to the data associated to the upper part of the cblk. More... | |
static void * | cblk_getdata (const SolverCblk *cblk, pastix_coefside_t side) |
Get the pointer to the data associated to the side part of the cblk. More... | |
static pastix_int_t | cblk_bloknbr (const SolverCblk *cblk) |
Compute the number of blocks in a column block. More... | |
static pastix_int_t | blok_rownbr (const SolverBlok *blok) |
Compute the number of rows of a block. More... | |
static pastix_int_t | blok_rownbr_ext (const SolverBlok *blok) |
Compute the number of rows of a contiguous block in front of the same cblk. More... | |
static int | blok_is_preselected (const SolverCblk *cblk, const SolverBlok *blok, const SolverCblk *fcbk) |
Return if a block is preselected as either part of the projection, or as a sub-diagonal block. More... | |
static pastix_int_t | cblk_rownbr (const SolverCblk *cblk) |
Compute the number of rows of a column block. More... | |
static pastix_int_t | stealQueue (SolverMatrix *solvmtx, int rank, int nbthreads) |
Task stealing method. More... | |
static int | is_block_inside_fblock (const SolverBlok *blok, const SolverBlok *fblok) |
Check if a block is included inside another one. More... | |
void | solverInit (SolverMatrix *solvmtx) |
Initialize the solver structure. More... | |
void | solverExit (SolverMatrix *solvmtx) |
Free the content of the solver matrix structure. More... | |
int | solverMatrixGen (SolverMatrix *solvmtx, const symbol_matrix_t *symbmtx, const pastix_order_t *ordeptr, const SimuCtrl *simuctl, const BlendCtrl *ctrl, PASTIX_Comm comm, isched_t *isched) |
Initialize the solver matrix structure. More... | |
int | solverMatrixGenSeq (SolverMatrix *solvmtx, const symbol_matrix_t *symbmtx, const pastix_order_t *ordeptr, const SimuCtrl *simuctl, const BlendCtrl *ctrl, PASTIX_Comm comm, isched_t *isched, pastix_int_t is_dbg) |
Initialize the solver matrix structure in sequential. More... | |
int | solverLoad (SolverMatrix *solvptr, FILE *stream) |
Load a solver matrix structure from a file. More... | |
int | solverSave (const SolverMatrix *solvptr, FILE *stream) |
Save a solver matrix structure into a file. More... | |
void | solverRealloc (SolverMatrix *solvptr) |
Realloc in a contiguous way a given solver structure. More... | |
SolverMatrix * | solverCopy (const SolverMatrix *solvptr, pastix_coeftype_t flttype) |
Generate a copy of a solver matrix structure. More... | |
int | solverCheck (const SolverMatrix *solvmtx) |
Checks the consistency of the given solver matrix structure. More... | |
int | solverDraw (const SolverMatrix *solvptr, FILE *stream, int verbose, const char *directory) |
Writes a PostScript picture of the low-rank solver matrix. More... | |
void | solverPrintStats (const SolverMatrix *solvptr) |
Print statistical information about the solver matrix structure. More... | |
void | solverRequestInit (solve_step_t solve_step, SolverMatrix *solvmtx) |
Instanciate the arrays for the requests according to the scheduler. More... | |
void | solverRequestExit (SolverMatrix *solvmtx) |
Free the arrays related to the requests. More... | |
void | solverRecvInit (pastix_coefside_t side, SolverMatrix *solvmtx, pastix_coeftype_t flttype) |
Allocate the reception buffer, and initiate the first persistant reception. More... | |
void | solverRecvExit (SolverMatrix *solvmtx) |
Free the array linked to pending reception. More... | |
void | solverRhsRecvInit (solve_step_t solve_step, SolverMatrix *solvmtx, pastix_coeftype_t flttype, pastix_rhs_t rhsb) |
Allocates the reception buffer, and initiate the first persistant reception. More... | |
void | solverRhsRecvExit (SolverMatrix *solvmtx) |
Frees the array linked to pending reception. More... | |
SolverBackup_t * | solverBackupInit (const SolverMatrix *solvmtx) |
Initialize the backup structure. More... | |
int | solverBackupRestore (SolverMatrix *solvmtx, const SolverBackup_t *b) |
Restore initial values. More... | |
void | solverBackupExit (SolverBackup_t *b) |
Free the solver backup data structure. More... | |
PaStiX solver structure header.
Definition in file solver.h.
struct args_solve_s |
Data Fields | ||
---|---|---|
pastix_scheduler_t | sched | |
solve_step_t | solve_step | |
pastix_solv_mode_t | mode | |
pastix_side_t | side | |
pastix_uplo_t | uplo | |
pastix_trans_t | trans | |
pastix_diag_t | diag |
struct solver_blok_recv_s |
Data Fields | ||
---|---|---|
pastix_int_t | frownum |
First row contributed in the block (lrownum+1 of the original block if none) |
pastix_int_t | lrownum |
Last row contributed in the block (frownum-1 of the original block if none) |
struct solver_cblk_recv_s |
Solver recv column block structure.
The structure works as a list of cblk for each remote owner to compute the exact contributions.
Data Fields | ||
---|---|---|
struct solver_cblk_recv_s * | next | |
pastix_int_t | ownerid | |
pastix_int_t | fcolnum |
First column contributed in the block (lcolnum+1 of the original cblk if none) |
pastix_int_t | lcolnum |
Last column contributed in the block (fcolnum-1 of the original cblk if none) |
solver_blok_recv_t | bloktab[1] |
Array of reception blocks |
struct task_s |
Data Fields | ||
---|---|---|
pastix_int_t | taskid |
COMP_1D DIAG E1 E2 |
pastix_int_t | prionum |
Priority value for the factorization |
pastix_int_t | cblknum |
Attached column block |
pastix_int_t | bloknum |
Attached block |
pastix_int_t volatile | ctrbcnt |
Total number of contributions |
struct solver_blok_s |
Data Fields | ||
---|---|---|
void * | handler[2] |
Runtime data handler |
pastix_int_t | lcblknm |
Local column block |
pastix_int_t | fcblknm |
Facing column block |
pastix_int_t | gbloknm |
Index in global bloktab (UNUSED) |
pastix_int_t | gfaninnm |
Global fanin block index, -1 if cblk not (FANIN | RECV) |
pastix_int_t | frownum |
First row index |
pastix_int_t | lrownum |
Last row index (inclusive) |
pastix_int_t | coefind |
Index in coeftab |
pastix_int_t | browind |
Index in browtab |
int8_t | inlast |
Index of the block among last separator (2), coupling with last separator (1) or other blocks (0) |
int | iluklvl |
The block ILU(k) level |
pastix_lrblock_t * | LRblock[2] |
Store the blok (L/U) in LR format. Allocated for the cblk. |
struct solver_cblk_s |
Data Fields | ||
---|---|---|
pastix_atomic_lock_t | lock |
Lock to protect computation on the cblk |
volatile int32_t | ctrbcnt |
Number of contribution to receive |
int8_t | cblktype |
Type of cblk |
int8_t | partitioned |
Bitmask to know if one of the side has been partitioned or not |
pastix_int_t | fcolnum |
First column index (Global numbering) |
pastix_int_t | lcolnum |
Last column index (Global numbering, inclusive) |
SolverBlok * | fblokptr |
First block in column (diagonal) |
pastix_int_t | stride |
Column block stride |
pastix_int_t | lcolidx |
First column index (Local numbering), used for the rhs vectors |
pastix_int_t | brownum |
First block in row facing the diagonal block in browtab, 0-based |
pastix_int_t | brown2d |
First 2D-block in row facing the diagonal block in browtab, 0-based |
pastix_int_t | sndeidx |
Global index of the original supernode the cblk belongs to |
pastix_int_t | gcblknum |
Global column block index |
pastix_int_t | bcscnum |
Index in the bcsctab if local cblk, -1 otherwise (FANIN | RECV) |
pastix_int_t | gfaninnum |
Global fanin column block index, -1 if not (FANIN | RECV) |
void * | lcoeftab |
Coefficients access vector, lower part |
void * | ucoeftab |
Coefficients access vector, upper part |
void * | handler[2] |
Runtime data handler |
pastix_int_t | selevtx |
Index to identify selected cblk for which intra-separator contributions are not compressed |
int | ownerid |
Rank of the owner |
int | threadid |
Rank of the accessing thread |
pastix_int_t | priority |
Priority of the cblk |
struct solver_matrix_s |
Solver column block structure.
This structure stores all the numerical information about the factorization, as well as the structure of the problem. Only local information to each process is stored in this structure.
Data Fields | ||
---|---|---|
int | restore |
Flag to indicate if it is require to restore data with solverBackupRestore: 0: No need, 1:After solve, 2:After Factorization |
pastix_int_t | baseval |
Base value for numberings |
pastix_int_t | nodenbr |
Number of local columns (after dof extension) |
pastix_int_t | coefnbr |
Number of locally stored values (after dof extension) |
pastix_int_t | gcblknbr |
Global number of column blocks |
pastix_int_t | cblknbr |
Local number of column blocks |
pastix_int_t | gfanincblknbr |
Global number of fanin column blocks |
pastix_int_t | faninnbr |
Local number of fanin cblk (included in cblknbr) |
pastix_int_t | fanincnt |
Number of sends to realize |
pastix_int_t | maxrecv |
Maximum blok size for a cblk_recv |
pastix_int_t | recvnbr |
Local number of recv cblk (included in cblknbr) |
pastix_int_t | recvcnt |
Number of receptions to realize |
pastix_int_t | cblkmax1d |
Rank of the last cblk not beeing enabled for 2D computations |
pastix_int_t | cblkmin2d |
Rank of the first cblk beeing enabled for 2D computations |
pastix_int_t | cblkmaxblk |
Maximum number of blocks per cblk |
pastix_int_t | cblkschur |
Index of the first local cblk in Schur |
pastix_int_t | nb2dcblk |
Number of 2D cblks |
pastix_int_t | nb2dblok |
Number of 2D blocks |
pastix_int_t | bloknbr |
Number of blocks |
pastix_int_t | gbloknbr |
Global number of blocks |
pastix_int_t | gfaninbloknbr |
Global number of fanin blocks |
pastix_int_t | brownbr |
Size of the browtab array |
SolverCblk *restrict | cblktab |
Array of solver column blocks [+1] |
SolverBlok *restrict | bloktab |
Array of solver blocks [+1] |
pastix_int_t *restrict | browtab |
Array of blocks |
pastix_coeftype_t | flttype |
valtab datatype: PastixFloat, PastixDouble, PastixComplex32 or PastixComplex64 |
int | globalalloc |
Boolean for global allocation of coeftab |
pastix_int_t * | gcbl2loc |
Array of local cblknum corresponding to gcblknum |
pastix_lr_t | lowrank |
Low-rank parameters |
pastix_factotype_t | factotype |
General or symmetric factorization? |
double | diagthreshold |
Diagonal threshold for pivoting |
volatile int32_t | nbpivots |
Number of pivots during the factorization |
pastix_int_t | offdmax | |
pastix_int_t | gemmmax | |
pastix_int_t | blokmax | |
pastix_int_t | colmax | |
int | clustnum | |
int | clustnbr | |
pastix_int_t | procnbr | |
pastix_int_t | thrdnbr | |
pastix_int_t | bublnbr | |
Task *restrict | tasktab | |
pastix_int_t | tasknbr | |
pastix_int_t | tasknbr_1dp | |
pastix_int_t ** | ttsktab | |
pastix_int_t * | ttsknbr | |
pastix_queue_t ** | computeQueue | |
pastix_int_t * | selevtx | |
MPI_Request * | reqtab |
Array of requests for MPI asynchronous messages |
pastix_int_t * | reqidx |
Array of local cblknum index corresponding to reqtab |
pastix_int_t | reqnbr |
Length of the reqtab/reqidx arrays |
volatile int32_t | reqnum |
Current amount of active requests (packed) |
pastix_atomic_lock_t | reqlock |
Lock to access the request arrays |
void * | rcoeftab |
Reception buffer for the communication |
size_t * | com_vector |
Matrix of communications between nodes. |
PASTIX_Comm | solv_comm |
typedef struct solver_cblk_recv_s solver_cblk_recv_t |
Solver recv column block structure.
The structure works as a list of cblk for each remote owner to compute the exact contributions.
|
inlinestatic |
Computes the current solve step.
[in] | side | Specify whether the off-diagonal blocks appear on the left or right in the equation. It has to be either PastixLeft or PastixRight. |
[in] | uplo | Specify whether the off-diagonal blocks are upper or lower triangular. It has to be either PastixUpper or PastixLower. |
[in] | trans | Specify the transposition used for the off-diagonal blocks. It has to be either PastixTrans or PastixConjTrans. |
Definition at line 306 of file solver.h.
References PastixLeft, PastixLower, PastixNoTrans, PastixRight, and PastixUpper.
|
inlinestatic |
Compute the number of columns in a column block.
[in] | cblk | The pointer to the column block. |
Definition at line 329 of file solver.h.
References solver_cblk_s::fcolnum, and solver_cblk_s::lcolnum.
Referenced by bcsc_init_coltab(), bvec_caxpy_seq(), bvec_ccopy_seq(), bvec_cdotu_seq(), bvec_cgemv_seq(), bvec_cnrm2_seq(), bvec_cscal_seq(), bvec_daxpy_seq(), bvec_dcopy_seq(), bvec_ddot_seq(), bvec_dgemv_seq(), bvec_dnrm2_seq(), bvec_dscal_seq(), bvec_saxpy_seq(), bvec_scopy_seq(), bvec_sdot_seq(), bvec_sgemv_seq(), bvec_snrm2_seq(), bvec_sscal_seq(), bvec_zaxpy_seq(), bvec_zcopy_seq(), bvec_zdotu_seq(), bvec_zgemv_seq(), bvec_znrm2_seq(), bvec_zscal_seq(), coeftab_cgetdiag(), coeftab_cmemory_fr(), coeftab_dgetdiag(), coeftab_dmemory_fr(), coeftab_sgetdiag(), coeftab_smemory_fr(), coeftab_zgetdiag(), coeftab_zmemory_fr(), core_cgemmsp_1d1d(), core_cgemmsp_1d2d(), core_cgemmsp_2d2d(), core_cgemmsp_block_frfr(), core_cgemmsp_block_frlr(), core_cgemmsp_block_lrlr(), core_cgemmsp_fulllr(), core_cgemmsp_lr(), core_cgemmsp_lrfr(), core_chetrfsp1d_gemm(), core_csytrfsp1d_gemm(), core_dgemmsp_1d1d(), core_dgemmsp_1d2d(), core_dgemmsp_2d2d(), core_dgemmsp_block_frfr(), core_dgemmsp_block_frlr(), core_dgemmsp_block_lrlr(), core_dgemmsp_fulllr(), core_dgemmsp_lr(), core_dgemmsp_lrfr(), core_dsytrfsp1d_gemm(), core_sgemmsp_1d1d(), core_sgemmsp_1d2d(), core_sgemmsp_2d2d(), core_sgemmsp_block_frfr(), core_sgemmsp_block_frlr(), core_sgemmsp_block_lrlr(), core_sgemmsp_fulllr(), core_sgemmsp_lr(), core_sgemmsp_lrfr(), core_ssytrfsp1d_gemm(), core_zgemmsp_1d1d(), core_zgemmsp_1d2d(), core_zgemmsp_2d2d(), core_zgemmsp_block_frfr(), core_zgemmsp_block_frlr(), core_zgemmsp_block_lrlr(), core_zgemmsp_fulllr(), core_zgemmsp_lr(), core_zgemmsp_lrfr(), core_zhetrfsp1d_gemm(), core_zsytrfsp1d_gemm(), cpublok_cadd(), cpublok_cadd_frfr(), cpublok_cadd_frlr(), cpublok_cadd_lrlr(), cpublok_calloc_lrws(), cpublok_cscalo(), cpublok_dadd(), cpublok_dadd_frfr(), cpublok_dadd_frlr(), cpublok_dadd_lrlr(), cpublok_dalloc_lrws(), cpublok_dscalo(), cpublok_sadd(), cpublok_sadd_frfr(), cpublok_sadd_frlr(), cpublok_sadd_lrlr(), cpublok_salloc_lrws(), cpublok_sscalo(), cpublok_zadd(), cpublok_zadd_frfr(), cpublok_zadd_frlr(), cpublok_zadd_lrlr(), cpublok_zalloc_lrws(), cpublok_zscalo(), cpucblk_cadd(), cpucblk_cadd_frfr(), cpucblk_cadd_frlr(), cpucblk_cadd_lrlr(), cpucblk_calloc_fr(), cpucblk_calloc_lr(), cpucblk_calloc_lrws(), cpucblk_ccompress(), cpucblk_ccompute_size(), cpucblk_ccompute_size_lr(), cpucblk_cdiff(), cpucblk_cgemmsp(), cpucblk_cgetschur_fr(), cpucblk_cgetschur_lr(), cpucblk_cmemory(), cpucblk_cpack_lr(), cpucblk_crecv_rhs_backward(), cpucblk_crecv_rhs_forward(), cpucblk_cscalo(), cpucblk_csend_rhs_backward(), cpucblk_csend_rhs_forward(), cpucblk_cuncompress(), cpucblk_cunpack_fr(), cpucblk_cunpack_lr(), cpucblk_dadd(), cpucblk_dadd_frfr(), cpucblk_dadd_frlr(), cpucblk_dadd_lrlr(), cpucblk_dalloc_fr(), cpucblk_dalloc_lr(), cpucblk_dalloc_lrws(), cpucblk_dcompress(), cpucblk_dcompute_size(), cpucblk_dcompute_size_lr(), cpucblk_ddiff(), cpucblk_dgemmsp(), cpucblk_dgetschur_fr(), cpucblk_dgetschur_lr(), cpucblk_dmemory(), cpucblk_dpack_lr(), cpucblk_drecv_rhs_backward(), cpucblk_drecv_rhs_forward(), cpucblk_dscalo(), cpucblk_dsend_rhs_backward(), cpucblk_dsend_rhs_forward(), cpucblk_duncompress(), cpucblk_dunpack_fr(), cpucblk_dunpack_lr(), cpucblk_sadd(), cpucblk_sadd_frfr(), cpucblk_sadd_frlr(), cpucblk_sadd_lrlr(), cpucblk_salloc_fr(), cpucblk_salloc_lr(), cpucblk_salloc_lrws(), cpucblk_scompress(), cpucblk_scompute_size(), cpucblk_scompute_size_lr(), cpucblk_sdiff(), cpucblk_sgemmsp(), cpucblk_sgetschur_fr(), cpucblk_sgetschur_lr(), cpucblk_smemory(), cpucblk_spack_lr(), cpucblk_srecv_rhs_backward(), cpucblk_srecv_rhs_forward(), cpucblk_sscalo(), cpucblk_ssend_rhs_backward(), cpucblk_ssend_rhs_forward(), cpucblk_suncompress(), cpucblk_sunpack_fr(), cpucblk_sunpack_lr(), cpucblk_zadd(), cpucblk_zadd_frfr(), cpucblk_zadd_frlr(), cpucblk_zadd_lrlr(), cpucblk_zalloc_fr(), cpucblk_zalloc_lr(), cpucblk_zalloc_lrws(), cpucblk_zcompress(), cpucblk_zcompute_size(), cpucblk_zcompute_size_lr(), cpucblk_zdiff(), cpucblk_zgemmsp(), cpucblk_zgetschur_fr(), cpucblk_zgetschur_lr(), cpucblk_zmemory(), cpucblk_zpack_lr(), cpucblk_zrecv_rhs_backward(), cpucblk_zrecv_rhs_forward(), cpucblk_zscalo(), cpucblk_zsend_rhs_backward(), cpucblk_zsend_rhs_forward(), cpucblk_zuncompress(), cpucblk_zunpack_fr(), cpucblk_zunpack_lr(), fct_blok_cadd_cost(), fct_blok_cgemmsp_cost(), fct_blok_cgetrfsp_cost(), fct_blok_chetrfsp_cost(), fct_blok_cpotrfsp_cost(), fct_blok_cscalo_cost(), fct_blok_ctrsmsp_cost(), fct_blok_dadd_cost(), fct_blok_dgemmsp_cost(), fct_blok_dgetrfsp_cost(), fct_blok_dpotrfsp_cost(), fct_blok_dscalo_cost(), fct_blok_dtrsmsp_cost(), fct_blok_sadd_cost(), fct_blok_sgemmsp_cost(), fct_blok_sgetrfsp_cost(), fct_blok_spotrfsp_cost(), fct_blok_sscalo_cost(), fct_blok_strsmsp_cost(), fct_blok_zadd_cost(), fct_blok_zgemmsp_cost(), fct_blok_zgetrfsp_cost(), fct_blok_zhetrfsp_cost(), fct_blok_zpotrfsp_cost(), fct_blok_zscalo_cost(), fct_blok_ztrsmsp_cost(), fct_cblk_cadd_cost(), fct_cblk_cgemmsp_cost(), fct_cblk_cgetrfsp_cost(), fct_cblk_chetrfsp_cost(), fct_cblk_cpotrfsp_cost(), fct_cblk_cpxtrfsp_cost(), fct_cblk_csytrfsp_cost(), fct_cblk_dadd_cost(), fct_cblk_dgemmsp_cost(), fct_cblk_dgetrfsp_cost(), fct_cblk_dpotrfsp_cost(), fct_cblk_dsytrfsp_cost(), fct_cblk_sadd_cost(), fct_cblk_sgemmsp_cost(), fct_cblk_sgetrfsp_cost(), fct_cblk_spotrfsp_cost(), fct_cblk_ssytrfsp_cost(), fct_cblk_zadd_cost(), fct_cblk_zgemmsp_cost(), fct_cblk_zgetrfsp_cost(), fct_cblk_zhetrfsp_cost(), fct_cblk_zpotrfsp_cost(), fct_cblk_zpxtrfsp_cost(), fct_cblk_zsytrfsp_cost(), parsec_sparse_matrix_init(), pastix_parsec_register_cblk_fr(), pastix_starpu_register(), pastix_starpu_rhs_data_register(), sequential_chetrf(), sequential_csytrf(), sequential_dsytrf(), sequential_ssytrf(), sequential_zhetrf(), sequential_zsytrf(), solve_blok_cgemm(), solve_blok_ctrsm(), solve_blok_dgemm(), solve_blok_dtrsm(), solve_blok_sgemm(), solve_blok_strsm(), solve_blok_zgemm(), solve_blok_ztrsm(), solver_copy(), solverCheck(), solverDraw(), solverMatrixGen(), solverMatrixGenSeq(), solverRhsRecvMax(), solvMatGen_max_buffers(), starpu_stask_blok_cadd_fwd_recv(), starpu_stask_blok_ccpy_bwd_recv(), starpu_stask_blok_cgemm(), starpu_stask_blok_dadd_fwd_recv(), starpu_stask_blok_dcpy_bwd_recv(), starpu_stask_blok_dgemm(), starpu_stask_blok_sadd_fwd_recv(), starpu_stask_blok_scpy_bwd_recv(), starpu_stask_blok_sgemm(), starpu_stask_blok_zadd_fwd_recv(), starpu_stask_blok_zcpy_bwd_recv(), starpu_stask_blok_zgemm(), thread_chetrf_dynamic(), thread_chetrf_static(), thread_csytrf_dynamic(), thread_csytrf_static(), thread_dsytrf_dynamic(), thread_dsytrf_static(), thread_ssytrf_dynamic(), thread_ssytrf_static(), thread_zhetrf_dynamic(), thread_zhetrf_static(), thread_zsytrf_dynamic(), and thread_zsytrf_static().
|
inlinestatic |
Get the pointer to the data associated to the lower part of the cblk.
[in] | cblk | The pointer to the column block. |
Definition at line 342 of file solver.h.
References solver_cblk_s::cblktype, solver_cblk_s::fblokptr, solver_cblk_s::lcoeftab, and solver_blok_s::LRblock.
Referenced by cpucblk_cgetrfsp1d(), cpucblk_cgetrfsp1dplus(), cpucblk_cgetrfsp1dplus_update(), cpucblk_chetrfsp1d(), cpucblk_chetrfsp1dplus(), cpucblk_chetrfsp1dplus_update(), cpucblk_cpotrfsp1d(), cpucblk_cpotrfsp1dplus(), cpucblk_cpotrfsp1dplus_update(), cpucblk_cpxtrfsp1d(), cpucblk_cpxtrfsp1dplus(), cpucblk_cpxtrfsp1dplus_update(), cpucblk_csytrfsp1d(), cpucblk_csytrfsp1dplus(), cpucblk_csytrfsp1dplus_update(), cpucblk_dgetrfsp1d(), cpucblk_dgetrfsp1dplus(), cpucblk_dgetrfsp1dplus_update(), cpucblk_dpotrfsp1d(), cpucblk_dpotrfsp1dplus(), cpucblk_dpotrfsp1dplus_update(), cpucblk_dsytrfsp1d(), cpucblk_dsytrfsp1dplus(), cpucblk_dsytrfsp1dplus_update(), cpucblk_sgetrfsp1d(), cpucblk_sgetrfsp1dplus(), cpucblk_sgetrfsp1dplus_update(), cpucblk_spotrfsp1d(), cpucblk_spotrfsp1dplus(), cpucblk_spotrfsp1dplus_update(), cpucblk_ssytrfsp1d(), cpucblk_ssytrfsp1dplus(), cpucblk_ssytrfsp1dplus_update(), cpucblk_zgetrfsp1d(), cpucblk_zgetrfsp1dplus(), cpucblk_zgetrfsp1dplus_update(), cpucblk_zhetrfsp1d(), cpucblk_zhetrfsp1dplus(), cpucblk_zhetrfsp1dplus_update(), cpucblk_zpotrfsp1d(), cpucblk_zpotrfsp1dplus(), cpucblk_zpotrfsp1dplus_update(), cpucblk_zpxtrfsp1d(), cpucblk_zpxtrfsp1dplus(), cpucblk_zpxtrfsp1dplus_update(), cpucblk_zsytrfsp1d(), cpucblk_zsytrfsp1dplus(), and cpucblk_zsytrfsp1dplus_update().
|
inlinestatic |
Get the pointer to the data associated to the upper part of the cblk.
[in] | cblk | The pointer to the column block. |
Definition at line 354 of file solver.h.
References solver_cblk_s::cblktype, solver_cblk_s::fblokptr, solver_blok_s::LRblock, and solver_cblk_s::ucoeftab.
Referenced by cpucblk_cgetrfsp1d(), cpucblk_cgetrfsp1dplus(), cpucblk_cgetrfsp1dplus_update(), cpucblk_chetrfsp1d(), cpucblk_csytrfsp1d(), cpucblk_dgetrfsp1d(), cpucblk_dgetrfsp1dplus(), cpucblk_dgetrfsp1dplus_update(), cpucblk_dsytrfsp1d(), cpucblk_sgetrfsp1d(), cpucblk_sgetrfsp1dplus(), cpucblk_sgetrfsp1dplus_update(), cpucblk_ssytrfsp1d(), cpucblk_zgetrfsp1d(), cpucblk_zgetrfsp1dplus(), cpucblk_zgetrfsp1dplus_update(), cpucblk_zhetrfsp1d(), and cpucblk_zsytrfsp1d().
|
inlinestatic |
Get the pointer to the data associated to the side part of the cblk.
[in] | cblk | The pointer to the column block. |
[in] | side | side of the data needed. PastixLCoef for the lower part PastixUCoef for the upper part. |
Definition at line 369 of file solver.h.
References solver_cblk_s::cblktype, solver_cblk_s::fblokptr, solver_cblk_s::lcoeftab, solver_blok_s::LRblock, PastixUCoef, and solver_cblk_s::ucoeftab.
|
inlinestatic |
Compute the number of blocks in a column block.
[in] | cblk | The pointer to the column block. |
Definition at line 383 of file solver.h.
References solver_cblk_s::fblokptr.
|
inlinestatic |
Compute the number of rows of a block.
[in] | blok | The pointer to the block. |
Definition at line 395 of file solver.h.
References solver_blok_s::frownum, and solver_blok_s::lrownum.
Referenced by blok_rownbr_ext(), cblk_rownbr(), core_cgemmsp_1d1d(), core_cgemmsp_1d2d(), core_cgemmsp_2d2d(), core_cgemmsp_block_frfr(), core_cgemmsp_block_frlr(), core_cgemmsp_block_lrlr(), core_cgemmsp_fulllr(), core_cgemmsp_lr(), core_cgemmsp_lrfr(), core_chetrfsp1d_gemm(), core_csytrfsp1d_gemm(), core_ctrsmsp_1d(), core_ctrsmsp_2d(), core_ctrsmsp_2dsub(), core_ctrsmsp_lr(), core_ctrsmsp_lrsub(), core_dgemmsp_1d1d(), core_dgemmsp_1d2d(), core_dgemmsp_2d2d(), core_dgemmsp_block_frfr(), core_dgemmsp_block_frlr(), core_dgemmsp_block_lrlr(), core_dgemmsp_fulllr(), core_dgemmsp_lr(), core_dgemmsp_lrfr(), core_dsytrfsp1d_gemm(), core_dtrsmsp_1d(), core_dtrsmsp_2d(), core_dtrsmsp_2dsub(), core_dtrsmsp_lr(), core_dtrsmsp_lrsub(), core_sgemmsp_1d1d(), core_sgemmsp_1d2d(), core_sgemmsp_2d2d(), core_sgemmsp_block_frfr(), core_sgemmsp_block_frlr(), core_sgemmsp_block_lrlr(), core_sgemmsp_fulllr(), core_sgemmsp_lr(), core_sgemmsp_lrfr(), core_ssytrfsp1d_gemm(), core_strsmsp_1d(), core_strsmsp_2d(), core_strsmsp_2dsub(), core_strsmsp_lr(), core_strsmsp_lrsub(), core_zgemmsp_1d1d(), core_zgemmsp_1d2d(), core_zgemmsp_2d2d(), core_zgemmsp_block_frfr(), core_zgemmsp_block_frlr(), core_zgemmsp_block_lrlr(), core_zgemmsp_fulllr(), core_zgemmsp_lr(), core_zgemmsp_lrfr(), core_zhetrfsp1d_gemm(), core_zsytrfsp1d_gemm(), core_ztrsmsp_1d(), core_ztrsmsp_2d(), core_ztrsmsp_2dsub(), core_ztrsmsp_lr(), core_ztrsmsp_lrsub(), cpublok_cadd_frfr(), cpublok_cadd_frlr(), cpublok_cadd_lrlr(), cpublok_calloc_lrws(), cpublok_ccompute_size_lr(), cpublok_cpack_lr(), cpublok_cscalo(), cpublok_cunpack_lr(), cpublok_dadd_frfr(), cpublok_dadd_frlr(), cpublok_dadd_lrlr(), cpublok_dalloc_lrws(), cpublok_dcompute_size_lr(), cpublok_dpack_lr(), cpublok_dscalo(), cpublok_dunpack_lr(), cpublok_sadd_frfr(), cpublok_sadd_frlr(), cpublok_sadd_lrlr(), cpublok_salloc_lrws(), cpublok_scompute_size_lr(), cpublok_spack_lr(), cpublok_sscalo(), cpublok_sunpack_lr(), cpublok_zadd_frfr(), cpublok_zadd_frlr(), cpublok_zadd_lrlr(), cpublok_zalloc_lrws(), cpublok_zcompute_size_lr(), cpublok_zpack_lr(), cpublok_zscalo(), cpublok_zunpack_lr(), cpucblk_cadd_frfr(), cpucblk_cadd_frlr(), cpucblk_cadd_lrlr(), cpucblk_calloc_lr(), cpucblk_calloc_lrws(), cpucblk_ccompress(), cpucblk_cdump(), cpucblk_cfillin_fr(), cpucblk_cfillin_lr(), cpucblk_cgemmsp(), cpucblk_cgetschur_fr(), cpucblk_cgetschur_lr(), cpucblk_cmemory(), cpucblk_cscalo(), cpucblk_cuncompress(), cpucblk_dadd_frfr(), cpucblk_dadd_frlr(), cpucblk_dadd_lrlr(), cpucblk_dalloc_lr(), cpucblk_dalloc_lrws(), cpucblk_dcompress(), cpucblk_ddump(), cpucblk_dfillin_fr(), cpucblk_dfillin_lr(), cpucblk_dgemmsp(), cpucblk_dgetschur_fr(), cpucblk_dgetschur_lr(), cpucblk_dmemory(), cpucblk_dscalo(), cpucblk_dsfillin_fr(), cpucblk_dsfillin_lr(), cpucblk_duncompress(), cpucblk_sadd_frfr(), cpucblk_sadd_frlr(), cpucblk_sadd_lrlr(), cpucblk_salloc_lr(), cpucblk_salloc_lrws(), cpucblk_scompress(), cpucblk_sdump(), cpucblk_sfillin_fr(), cpucblk_sfillin_lr(), cpucblk_sgemmsp(), cpucblk_sgetschur_fr(), cpucblk_sgetschur_lr(), cpucblk_smemory(), cpucblk_sscalo(), cpucblk_suncompress(), cpucblk_zadd_frfr(), cpucblk_zadd_frlr(), cpucblk_zadd_lrlr(), cpucblk_zalloc_lr(), cpucblk_zalloc_lrws(), cpucblk_zcfillin_fr(), cpucblk_zcfillin_lr(), cpucblk_zcompress(), cpucblk_zdump(), cpucblk_zfillin_fr(), cpucblk_zfillin_lr(), cpucblk_zgemmsp(), cpucblk_zgetschur_fr(), cpucblk_zgetschur_lr(), cpucblk_zmemory(), cpucblk_zscalo(), cpucblk_zuncompress(), fct_blok_cadd_cost(), fct_blok_dadd_cost(), fct_blok_sadd_cost(), fct_blok_zadd_cost(), fct_cblk_cgemmsp_cost(), fct_cblk_dgemmsp_cost(), fct_cblk_sgemmsp_cost(), fct_cblk_zgemmsp_cost(), parsec_sparse_matrix_init(), solve_blok_cgemm(), solve_blok_dgemm(), solve_blok_sgemm(), solve_blok_zgemm(), solve_cblk_cdiag(), solve_cblk_ddiag(), solve_cblk_sdiag(), solve_cblk_zdiag(), solverDraw(), and solvMatGen_max_buffers().
|
inlinestatic |
Compute the number of rows of a contiguous block in front of the same cblk.
[in] | blok | The pointer to the block. |
Definition at line 407 of file solver.h.
References blok_rownbr(), and pastix_int_t.
Referenced by fct_blok_cgemmsp_cost(), fct_blok_cscalo_cost(), fct_blok_ctrsmsp_cost(), fct_blok_dgemmsp_cost(), fct_blok_dscalo_cost(), fct_blok_dtrsmsp_cost(), fct_blok_sgemmsp_cost(), fct_blok_sscalo_cost(), fct_blok_strsmsp_cost(), fct_blok_zgemmsp_cost(), fct_blok_zscalo_cost(), and fct_blok_ztrsmsp_cost().
|
inlinestatic |
Return if a block is preselected as either part of the projection, or as a sub-diagonal block.
[in] | cblk | The pointer to the cblk to which belong the block. |
[in] | blok | The pointer to the block. |
[in] | fcbk | The pointer to the facing cblk of the block. |
Definition at line 432 of file solver.h.
References solver_cblk_s::fblokptr, solver_cblk_s::selevtx, and solver_cblk_s::sndeidx.
Referenced by coeftab_ccblkComputePreselect(), coeftab_dcblkComputePreselect(), coeftab_scblkComputePreselect(), coeftab_zcblkComputePreselect(), cpucblk_cmemory(), cpucblk_dmemory(), cpucblk_smemory(), and cpucblk_zmemory().
|
inlinestatic |
Compute the number of rows of a column block.
[in] | cblk | The pointer to the column block. |
Definition at line 449 of file solver.h.
References blok_rownbr(), solver_cblk_s::fblokptr, and pastix_int_t.
|
inlinestatic |
Task stealing method.
[in,out] | solvmtx | The pointer to the solverMatrix. |
[in] | rank | Rank of the computeQueue. |
[in] | nbthreads | Total amount of threads in the node. |
Definition at line 471 of file solver.h.
References pastix_int_t, and pqueuePop().
Referenced by thread_cgetrf_dynamic(), thread_chetrf_dynamic(), thread_cpotrf_dynamic(), thread_cpxtrf_dynamic(), thread_csytrf_dynamic(), thread_dgetrf_dynamic(), thread_dpotrf_dynamic(), thread_dsytrf_dynamic(), thread_sgetrf_dynamic(), thread_spotrf_dynamic(), thread_ssytrf_dynamic(), thread_zgetrf_dynamic(), thread_zhetrf_dynamic(), thread_zpotrf_dynamic(), thread_zpxtrf_dynamic(), and thread_zsytrf_dynamic().
|
inlinestatic |
Check if a block is included inside another one.
Indicate if a blok is included inside an other block. i.e. indicate if the row range of the first block is included in the one of the second.
[in] | blok | The block that is tested for inclusion. |
[in] | fblok | The block that is suppose to include the first one. |
true | if the first block is included in the second one. |
false | if the first block is not included in the second one. |
Definition at line 504 of file solver.h.
References solver_blok_s::frownum, and solver_blok_s::lrownum.
Referenced by coeftabComputeCblkILULevels(), core_cgemmsp_1d2d(), core_cgemmsp_2d2d(), core_cgemmsp_block_frfr(), core_cgemmsp_block_frlr(), core_cgemmsp_block_lrlr(), core_cgemmsp_fulllr(), core_cgemmsp_lr(), core_cgemmsp_lrfr(), core_chetrfsp1d_gemm(), core_csytrfsp1d_gemm(), core_dgemmsp_1d2d(), core_dgemmsp_2d2d(), core_dgemmsp_block_frfr(), core_dgemmsp_block_frlr(), core_dgemmsp_block_lrlr(), core_dgemmsp_fulllr(), core_dgemmsp_lr(), core_dgemmsp_lrfr(), core_dsytrfsp1d_gemm(), core_sgemmsp_1d2d(), core_sgemmsp_2d2d(), core_sgemmsp_block_frfr(), core_sgemmsp_block_frlr(), core_sgemmsp_block_lrlr(), core_sgemmsp_fulllr(), core_sgemmsp_lr(), core_sgemmsp_lrfr(), core_ssytrfsp1d_gemm(), core_zgemmsp_1d2d(), core_zgemmsp_2d2d(), core_zgemmsp_block_frfr(), core_zgemmsp_block_frlr(), core_zgemmsp_block_lrlr(), core_zgemmsp_fulllr(), core_zgemmsp_lr(), core_zgemmsp_lrfr(), core_zhetrfsp1d_gemm(), core_zsytrfsp1d_gemm(), cpublok_cadd_frfr(), cpublok_cadd_frlr(), cpublok_cadd_lrlr(), cpublok_dadd_frfr(), cpublok_dadd_frlr(), cpublok_dadd_lrlr(), cpublok_sadd_frfr(), cpublok_sadd_frlr(), cpublok_sadd_lrlr(), cpublok_zadd_frfr(), cpublok_zadd_frlr(), cpublok_zadd_lrlr(), cpucblk_cadd_frfr(), cpucblk_cadd_frlr(), cpucblk_cadd_lrlr(), cpucblk_cgeaddsp1d(), cpucblk_dadd_frfr(), cpucblk_dadd_frlr(), cpucblk_dadd_lrlr(), cpucblk_dgeaddsp1d(), cpucblk_sadd_frfr(), cpucblk_sadd_frlr(), cpucblk_sadd_lrlr(), cpucblk_sgeaddsp1d(), cpucblk_zadd_frfr(), cpucblk_zadd_frlr(), cpucblk_zadd_lrlr(), cpucblk_zgeaddsp1d(), and solverCheck().