PaStiX Handbook  6.3.2
Symbolic Factorization

Functions to generate and manipulate the symbolic factorization structure. More...

Modules

 Symbolic internal function documentation
 

Data Structures

struct  symbol_cblk_s
 Symbol column block structure. More...
 
struct  symbol_blok_s
 Symbol block structure. More...
 
struct  symbol_matrix_s
 Symbol matrix structure. More...
 

Typedefs

typedef struct symbol_cblk_s symbol_cblk_t
 Symbol column block structure.
 
typedef struct symbol_blok_s symbol_blok_t
 Symbol block structure.
 
typedef struct symbol_matrix_s symbol_matrix_t
 Symbol matrix structure. More...
 

Functions

static pastix_int_t symbol_cblk_get_colnum (const symbol_matrix_t *symbmtx, symbol_cblk_t *symbcblk, pastix_int_t *fcolnum, pastix_int_t *lcolnum)
 Get the expanded column indexes of a symbol_cblk. More...
 
static pastix_int_t symbol_blok_get_rownum (const symbol_matrix_t *symbmtx, symbol_blok_t *symbblok, pastix_int_t *frownum, pastix_int_t *lrownum)
 Get the expanded row index of a symbol_blok. More...
 
static int is_symbblock_inside_fblock (const symbol_blok_t *blok, const symbol_blok_t *fblok)
 Check if a block is included inside another one. More...
 
void pastixSymbolPrint (const symbol_matrix_t *symbptr, FILE *file)
 Print the given block matrix structure in human readable format. More...
 

Symbol basic subroutines

void pastixSymbolInit (const pastix_graph_t *graph, const pastix_order_t *order, symbol_matrix_t *symbptr)
 Initialize the symbol structure. More...
 
void pastixSymbolExit (symbol_matrix_t *symbptr)
 Free the content of symbolic matrix. More...
 
void pastixSymbolBase (symbol_matrix_t *const symbptr, const pastix_int_t baseval)
 Sets the base of the given symbol matrix structure to the given base value. More...
 
void pastixSymbolRealloc (symbol_matrix_t *symbptr)
 Reallocate the data structure to optimize the memory alignment. More...
 
int pastixSymbolCheck (const symbol_matrix_t *const symbptr)
 Checks the consistency of the given symbolic block matrix. More...
 
void pastixSymbolExpand (symbol_matrix_t *symbptr)
 Expand the symbol matrix structure based on the dof information (compressed -> expanded) More...
 

Symbol IO subroutines

int pastixSymbolSave (const symbol_matrix_t *const symbptr, FILE *const stream)
 Save the given block matrix structure to the given stream. More...
 
int pastixSymbolLoad (symbol_matrix_t *const symbptr, FILE *const stream)
 Load the given block matrix structure from the given stream. More...
 
int pastixSymbolDraw (const symbol_matrix_t *const symbptr, FILE *const stream)
 Export the symbol structure in a PostScript format. More...
 
void pastixSymbolDrawMap (pastix_data_t *pastix_data, const char *extname, pastix_int_t sndeidx)
 Dump a separator mapping into a map file. More...
 

Symbol statistical information subroutines

void pastixSymbolPrintStats (const symbol_matrix_t *symbptr)
 Print statistical information about the symbolic matrix structure. More...
 
size_t pastixSymbolGetNNZ (const symbol_matrix_t *symbptr)
 Computes the number of non-zero elements in L. More...
 
void pastixSymbolGetFlops (const symbol_matrix_t *symbmtx, pastix_coeftype_t flttype, pastix_factotype_t factotype, double *thflops, double *rlflops)
 Computes the number of theoretical and real flops. More...
 
void pastixSymbolGetTimes (const symbol_matrix_t *symbmtx, pastix_coeftype_t flttype, pastix_factotype_t factotype, double *cblkcost, double *blokcost)
 Computes the cost of structure for the costMatrixBuild() function. More...
 

Symbol reordering subroutines

void pastixSymbolReordering (pastix_data_t *pastix_data)
 Compute the reordering on the complete matrix. More...
 
void pastixSymbolReorderingPrintComplexity (const symbol_matrix_t *symbptr)
 Compute the number of operations required to compute the reordering on the complete matrix. More...
 

Symbol construction subroutines

int pastixSymbolFaxDirect (symbol_matrix_t *symbptr, const pastix_graph_t *graphA, const pastix_order_t *ordeptr)
 Compute the block symbolic factorization of the given matrix graph according to the given vertex ordering. More...
 
int pastixSymbolFaxILUk (symbol_matrix_t *symbptr, pastix_int_t levelk, const pastix_graph_t *graphA, const pastix_order_t *ordeptr)
 Create the symbol matrix from the graph of the non zero pattern of the factorized matrix and the supernode partition. More...
 
void pastixSymbolBuildRowtab (symbol_matrix_t *symbptr)
 Construct the browtab array that stores the blocks in a CSR way. More...
 
pastix_int_t pastixSymbolGetFacingBloknum (const symbol_matrix_t *symbptr, pastix_int_t bloksrc, pastix_int_t bloknum, pastix_int_t startsearch, int ricar)
 Search the targeted block C for a couple of blocks A and B. More...
 

Detailed Description

Functions to generate and manipulate the symbolic factorization structure.

This module provides the set of function to generate the symbolic factorization structure based on a given graph, and an associated ordering. The symbolic structure is described in the symbol_matrix_t structure, and it can be generated through two different algorithms respectively for direct and incomplete factorization.


Data Type Documentation

◆ symbol_cblk_s

struct symbol_cblk_s

Symbol column block structure.

Definition at line 45 of file symbol.h.

Data Fields
pastix_int_t fcolnum

First column index

pastix_int_t lcolnum

Last column index (inclusive)

pastix_int_t bloknum

First block in column (diagonal)

pastix_int_t brownum

First block in row facing the diagonal block in browtab, 0-based

int8_t selevtx

◆ symbol_blok_s

struct symbol_blok_s

Symbol block structure.

Definition at line 59 of file symbol.h.

Data Fields
pastix_int_t frownum

First row index

pastix_int_t lrownum

Last row index (inclusive)

pastix_int_t lcblknm

Local column block

pastix_int_t fcblknm

Facing column block

◆ symbol_matrix_s

struct symbol_matrix_s

Symbol matrix structure.

This structure describes the symbolic block structure of the factorized matrix L, U is never stored as it is a symmetry of L. This structure is global and replicated on all processes. The default way to number the block is the CSC format where block are continuously number per column, the browtab array stores the CSR representation of the L structure to have a faster access to the list of blocks updating a column block.

Definition at line 77 of file symbol.h.

Data Fields
pastix_int_t baseval

Base value for numbering

pastix_int_t cblknbr

Number of column blocks

pastix_int_t bloknbr

Number of blocks

pastix_int_t nodenbr

Number of nodes (Equal to gN in spm)

pastix_int_t schurfcol

First column of the schur complement

symbol_cblk_t * cblktab

Array of column blocks [+1,based]

symbol_blok_t * bloktab

Array of blocks in CSC format [based]

pastix_int_t * browtab

Array of blocks in CSR format [based]

pastix_int_t browmax

Maximum number of input edges per node

pastix_int_t dof

Degrees of freedom per node (constant if > 0, variadic if < 1

pastix_int_t * dofs

Array of the first column of each element in the expanded matrix [+1,based]

Typedef Documentation

◆ symbol_matrix_t

Symbol matrix structure.

This structure describes the symbolic block structure of the factorized matrix L, U is never stored as it is a symmetry of L. This structure is global and replicated on all processes. The default way to number the block is the CSC format where block are continuously number per column, the browtab array stores the CSR representation of the L structure to have a faster access to the list of blocks updating a column block.

Function Documentation

◆ symbol_cblk_get_colnum()

static pastix_int_t symbol_cblk_get_colnum ( const symbol_matrix_t symbmtx,
symbol_cblk_t symbcblk,
pastix_int_t fcolnum,
pastix_int_t lcolnum 
)
inlinestatic

Get the expanded column indexes of a symbol_cblk.

Parameters
[in]symbmtxPointer to the symbol matrix.
[in]symbcblkThe pointer to the current symbol_cblk.
[in,out]fcolnumFirst column index of the current cblk.
[in,out]lcolnumLast column index of the current cblk.
Returns
The number of columns of the expanded cblk.

Definition at line 116 of file symbol.h.

References symbol_matrix_s::dof, symbol_matrix_s::dofs, symbol_cblk_s::fcolnum, and symbol_cblk_s::lcolnum.

Referenced by solverMatrixGen(), and solvMatGen_register_local_cblk().

◆ symbol_blok_get_rownum()

static pastix_int_t symbol_blok_get_rownum ( const symbol_matrix_t symbmtx,
symbol_blok_t symbblok,
pastix_int_t frownum,
pastix_int_t lrownum 
)
inlinestatic

Get the expanded row index of a symbol_blok.

Parameters
[in]symbmtxPointer to the symbol matrix.
[in]symbblokThe pointer to the current symbol_blok.
[in,out]frownumFirst row index of the current blok.
[in,out]lrownumLast row index of the current blok.
Returns
The number of rows of the expanded blok.

Definition at line 155 of file symbol.h.

References symbol_matrix_s::dof, symbol_matrix_s::dofs, symbol_blok_s::frownum, and symbol_blok_s::lrownum.

Referenced by solvMatGen_register_local_cblk(), and solvMatGen_reorder_browtab().

◆ is_symbblock_inside_fblock()

static int is_symbblock_inside_fblock ( const symbol_blok_t blok,
const symbol_blok_t fblok 
)
inlinestatic

Check if a block is included inside another one.

Indicate if a blok is included inside another block. i.e. indicate if the row range of the first block is included in the one of the second.

Parameters
[in]blokThe block that is tested for inclusion.
[in]fblokThe block that is suppose to include the first one.
Return values
trueif the first block is included in the second one.
falseif the first block is not included in the second one.

Definition at line 185 of file symbol.h.

References symbol_blok_s::frownum, and symbol_blok_s::lrownum.

Referenced by solver_recv_add_contrib().

◆ pastixSymbolInit()

void pastixSymbolInit ( const pastix_graph_t *  graph,
const pastix_order_t order,
symbol_matrix_t symbptr 
)

Initialize the symbol structure.

Initialized the permuted dof array if graph and order are provided. The symbol is considered as dof = 1 otherwise.

Parameters
[in]graphThe original graph of the matrix
[in]orderThe ordering structure describing the permutation of the unknowns in the compressed form.
[in,out]symbptrThe symbol structure to initialize.

Definition at line 105 of file symbol.c.

References symbol_matrix_s::dof, symbol_matrix_s::schurfcol, and symbol_init_adddofs().

Referenced by pastix_subtask_symbfact(), and pastixSymbolExit().

◆ pastixSymbolExit()

void pastixSymbolExit ( symbol_matrix_t symbptr)

Free the content of symbolic matrix.

All the arrays from the structure are freed and the structure is memset to 0 at exit, but the symbol itself is not freed. It will require a new call to pastixSymbolInit() if the memory space area needs to be reused for a new symbol matrix.

Parameters
[in,out]symbptrThe pointer to the structure to free.

Definition at line 137 of file symbol.c.

References symbol_matrix_s::bloktab, symbol_matrix_s::browtab, symbol_matrix_s::cblktab, symbol_matrix_s::dofs, and pastixSymbolInit().

Referenced by extraCblkMerge(), pastix_subtask_blend(), pastix_subtask_symbfact(), pastixFinalize(), and pastixSymbolLoad().

◆ pastixSymbolBase()

void pastixSymbolBase ( symbol_matrix_t *const  symbptr,
const pastix_int_t  baseval 
)

◆ pastixSymbolRealloc()

void pastixSymbolRealloc ( symbol_matrix_t symbptr)

Reallocate the data structure to optimize the memory alignment.

This function is used when the symbol need to be shrinked to a smaller or larger set of blocks and column blocks. The original data is copied in the new arrays.

Parameters
[in,out]symbptrThe pointer to the structure that needs to be reallocated.

Definition at line 170 of file symbol.c.

References symbol_matrix_s::bloknbr, symbol_matrix_s::bloktab, symbol_matrix_s::cblknbr, and symbol_matrix_s::cblktab.

Referenced by pastix_subtask_symbfact().

◆ pastixSymbolCheck()

int pastixSymbolCheck ( const symbol_matrix_t *const  symbptr)

Checks the consistency of the given symbolic block matrix.

Because of incomplete factorization, from version 1.0, no check is performed regarding the existence of facing blocks in facing columns.

Todo:
Complete test set to check the brow information
Parameters
[in]symbptrThe symbol structure to check.
Return values
0if the symbol matrix is correct
1if incorrect

Definition at line 47 of file symbol_check.c.

References symbol_matrix_s::baseval, symbol_matrix_s::bloknbr, symbol_cblk_s::bloknum, symbol_matrix_s::bloktab, symbol_cblk_s::brownum, symbol_matrix_s::cblknbr, symbol_matrix_s::cblktab, symbol_cblk_s::fcolnum, symbol_cblk_s::lcolnum, symbol_matrix_s::nodenbr, pastix_int_t, and symbol_matrix_s::schurfcol.

Referenced by pastix_subtask_blend(), pastix_subtask_symbfact(), pastixSymbolExpand(), and splitSymbol().

◆ pastixSymbolExpand()

void pastixSymbolExpand ( symbol_matrix_t symbptr)

Expand the symbol matrix structure based on the dof information (compressed -> expanded)

Parameters
[in,out]symbptrThe symbol structure to expand. On entry, the information is based on the compressed graph. On exit, the information are expanded to the size of the expanded matrix.

Definition at line 181 of file symbol_expand.c.

References symbol_matrix_s::dof, symbol_matrix_s::dofs, pastixSymbolBase(), pastixSymbolCheck(), symbol_expand_fix(), and symbol_expand_var().

Referenced by pastixExpand().

◆ pastixSymbolSave()

int pastixSymbolSave ( const symbol_matrix_t *const  symbptr,
FILE *const  stream 
)

Save the given block matrix structure to the given stream.

Parameters
[in,out]symbptrThe symbolic matrix structure to write.
[in,out]streamThe stream to which to write the structure.
Return values
0on success.
!0on failure.

Definition at line 147 of file symbol_io.c.

References symbol_matrix_s::baseval, symbol_matrix_s::bloknbr, symbol_cblk_s::bloknum, symbol_matrix_s::bloktab, symbol_matrix_s::cblknbr, symbol_matrix_s::cblktab, symbol_blok_s::fcblknm, symbol_cblk_s::fcolnum, symbol_blok_s::frownum, symbol_cblk_s::lcolnum, symbol_blok_s::lrownum, and symbol_matrix_s::nodenbr.

Referenced by pastix_subtask_blend(), and pastix_subtask_symbfact().

◆ pastixSymbolLoad()

int pastixSymbolLoad ( symbol_matrix_t *const  symbptr,
FILE *const  stream 
)

Load the given block matrix structure from the given stream.

Parameters
[in,out]symbptrThe symbolic matrix structure to fill in.
[in,out]streamThe stream from which to read the structure.
Return values
0on success.
!0on failure.

Definition at line 51 of file symbol_io.c.

References symbol_matrix_s::baseval, symbol_matrix_s::bloknbr, symbol_cblk_s::bloknum, symbol_matrix_s::bloktab, symbol_matrix_s::cblknbr, symbol_matrix_s::cblktab, symbol_matrix_s::dof, symbol_blok_s::fcblknm, symbol_cblk_s::fcolnum, symbol_blok_s::frownum, symbol_cblk_s::lcolnum, symbol_blok_s::lrownum, symbol_matrix_s::nodenbr, pastix_int_t, and pastixSymbolExit().

Referenced by pastix_subtask_symbfact().

◆ pastixSymbolDraw()

int pastixSymbolDraw ( const symbol_matrix_t *const  symbptr,
FILE *const  stream 
)

Export the symbol structure in a PostScript format.

This routine writes to the given stream a PostScript (tm) picture of the symbolic block matrix, with diagonal blocks in black and off-diagonal blocks in dark gray.

Parameters
[in]symbptrThe pointer to the symbolic structure to draw.
[in,out]streamThe stream where to write the output.
Return values
0on success
-1on error

Definition at line 248 of file symbol_draw.c.

References pastixSymbolDrawFunc().

Referenced by pastix_subtask_blend(), and pastix_subtask_symbfact().

◆ pastixSymbolDrawMap()

void pastixSymbolDrawMap ( pastix_data_t pastix_data,
const char *  extname,
pastix_int_t  sndeidx 
)

Dump a separator mapping into a map file.

Parameters
[in,out]pastix_dataThe pastix data structure that holds the symbol and the ordering. On exit the output directories may be initialized, if not previously.
[in]extnameFilename extension to specify the .map file if multiple.
[in]sndeidxThe index of the supernode to dump into file.

Definition at line 44 of file symbol_draw_map.c.

References symbol_matrix_s::cblknbr, symbol_matrix_s::cblktab, pastix_data_s::dir_global, symbol_cblk_s::fcolnum, symbol_cblk_s::lcolnum, pastix_data_s::ordemesh, pastix_fopenw(), pastix_gendirectories(), pastix_int_t, pastix_order_s::sndetab, and pastix_data_s::symbmtx.

Referenced by pastix_subtask_blend().

◆ pastixSymbolPrintStats()

void pastixSymbolPrintStats ( const symbol_matrix_t symbptr)

Print statistical information about the symbolic matrix structure.

Parameters
[in]symbptrThe pointer to the symbolic structure.

Definition at line 389 of file symbol.c.

References symbol_matrix_s::bloknbr, symbol_matrix_s::cblknbr, and pastix_int_t.

Referenced by pastix_subtask_blend(), pastix_subtask_symbfact(), and splitSymbol().

◆ pastixSymbolGetNNZ()

size_t pastixSymbolGetNNZ ( const symbol_matrix_t symbptr)

Computes the number of non-zero elements in L.

This computes the number of non-zero elements stored in the symbol matrix in order to compute the fill-in. This routines returns the number of non-zero of the strictly lower part of the matrix.

Parameters
[in]symbptrThe symbol structure to study.
Returns
The number of non zero elements in the strictly lower part of the full symbol matrix.

Definition at line 325 of file symbol_cost.c.

References symbol_matrix_s::baseval, symbol_cblk_s::bloknum, symbol_matrix_s::bloktab, symbol_matrix_s::cblknbr, symbol_matrix_s::cblktab, symbol_matrix_s::dof, symbol_matrix_s::dofs, symbol_cblk_s::fcolnum, symbol_blok_s::frownum, symbol_cblk_s::lcolnum, symbol_blok_s::lrownum, and pastix_int_t.

Referenced by pastix_subtask_blend(), and pastix_subtask_symbfact().

◆ pastixSymbolGetFlops()

void pastixSymbolGetFlops ( const symbol_matrix_t symbmtx,
pastix_coeftype_t  flttype,
pastix_factotype_t  factotype,
double *  thflops,
double *  rlflops 
)

Computes the number of theoretical and real flops.

Given a symbolic factorization structure, the type of factorization: Cholesky or LU, and the arithmetic, this function will return the total number of floating point operation that will be performed during the numerical factorization.

Parameters
[in]symbmtxThe symbol structure to study.
[in]flttypeThe floating type of the elements in the matrix. PastixPattern, PastixFloat, PastixDouble, PastixComplex32 or PastixComplex64. In case of PastixPattern, values for PastixDouble are returned.
[in]factotypeThe factorization algorithm to perform: PastixFactLLT, PastixFactLDLT, PastixFactLLH, PastixFactLDLH or PastixFactLU.
[out]thflopsReturns the number of theoretical flops to perform. NULL if not asked.
[out]rlflopsReturns the number of real flops to perform, taking into account copies and scatter operations. NULL if not asked.

Definition at line 422 of file symbol_cost.c.

References symbol_matrix_s::cblknbr, flopstable, recursive_sum(), sum1d(), and sum2d().

Referenced by pastix_subtask_blend(), and pastix_subtask_symbfact().

◆ pastixSymbolGetTimes()

void pastixSymbolGetTimes ( const symbol_matrix_t symbmtx,
pastix_coeftype_t  flttype,
pastix_factotype_t  factotype,
double *  cblkcost,
double *  blokcost 
)

Computes the cost of structure for the costMatrixBuild() function.

This function iterates on the column-blocks and blocks to compute the cost of the operation performed on each of those elements for the costMatrixBuild() function that is used in the simulation for the data mapping. It returns an array of cost for each type of element.

Parameters
[in]symbmtxThe symbol structure to study.
[in]flttypeThe floating type of the elements in the matrix. PastixPattern, PastixFloat, PastixDouble, PastixComplex32 or PastixComplex64. In case of PastixPattern, values for PastixDouble are returned.
[in]factotypeThe factorization algorithm to perform: PastixFactLLT, PastixFactLDLT, PastixFactLLH, PastixFactLDLH or PastixFactLU.
[in,out]cblkcostAn allocated array of size cblknbr that will holds the cost per cblk on exit.
[in,out]blokcostAn allocated array of size bloknbr that will holds the cost per blok on exit.

Definition at line 480 of file symbol_cost.c.

References symbol_matrix_s::bloknbr, symbol_cblk_s::bloknum, symbol_matrix_s::cblknbr, symbol_matrix_s::cblktab, pastix_int_t, perfstable, and sum2dext().

Referenced by costMatrixBuild().

◆ pastixSymbolReordering()

void pastixSymbolReordering ( pastix_data_t pastix_data)

Compute the reordering on the complete matrix.

Parameters
[in]pastix_dataThe pastix_data structure that describes the solver instance. On exit, the field symbmtx is updated with the new symbol matrix, and the field ordemesh is updated with the new ordering. The split_level field activates the split level heuristic, dividing distances computations into two stages: for upper and for lower contruibuting supernodes. If a resulting distance for upper supernodes is large enough, the computation is stopped, as long as it will be large taking into account lower contributing supernodes. The stop_criterion field disregards rows that are far away.

Compute the reordering using either sequential or parallel method

Definition at line 619 of file symbol_reordering.c.

References symbol_matrix_s::cblknbr, compute_cblklevel(), symbol_matrix_s::nodenbr, pastix_data_s::ordemesh, pastix_int_t, pastix_order_s::peritab, pastix_order_s::permtab, pastix_data_s::symbmtx, symbol_reorder(), and pastix_order_s::treetab.

◆ pastixSymbolReorderingPrintComplexity()

void pastixSymbolReorderingPrintComplexity ( const symbol_matrix_t symbptr)

Compute the number of operations required to compute the reordering on the complete matrix.

The number of operation is compuyted and then printed on the standard output.

Parameters
[in]symbptrThe pointer to the symbolic structure.

Definition at line 668 of file symbol_reordering.c.

References symbol_matrix_s::bloktab, symbol_cblk_s::brownum, symbol_matrix_s::browtab, symbol_matrix_s::cblknbr, symbol_matrix_s::cblktab, symbol_blok_s::fcblknm, symbol_cblk_s::fcolnum, symbol_blok_s::frownum, symbol_cblk_s::lcolnum, symbol_blok_s::lrownum, pastix_int_t, and symbol_matrix_s::schurfcol.

◆ pastixSymbolFaxDirect()

int pastixSymbolFaxDirect ( symbol_matrix_t symbptr,
const pastix_graph_t *  graphA,
const pastix_order_t ordeptr 
)

Compute the block symbolic factorization of the given matrix graph according to the given vertex ordering.

pastixSymbolFaxGraph() could have called pastixSymbolFax() in the regular way, as do all of the grid-like factorization routines. However, for efficiency reasons, we have decided to inline pastixSymbolFax(), to avoid a function call for every arc.

Parameters
[in,out]symbptrThe symbolic matrix structure to fill in.
[in]graphATODO
[in]ordeptrThe ordering structure that contains the permutation and inverse permutation vector, as well as the list of supernodes and the associated tree.
Return values
0on success.
!0on failure.

Definition at line 66 of file symbol_fax_direct.c.

References pastix_int_t.

Referenced by pastix_subtask_symbfact().

◆ pastixSymbolFaxILUk()

int pastixSymbolFaxILUk ( symbol_matrix_t symbptr,
pastix_int_t  levelk,
const pastix_graph_t *  graphA,
const pastix_order_t ordeptr 
)

◆ pastixSymbolBuildRowtab()

void pastixSymbolBuildRowtab ( symbol_matrix_t symbptr)

Construct the browtab array that stores the blocks in a CSR way.

The browtab is an equivalent of the columns array in a CSR for the symbolic structure in terms of block indexes.

Parameters
[in,out]symbptrThe pointer to the symbolic structure to update.

Definition at line 302 of file symbol.c.

References symbol_matrix_s::bloknbr, symbol_cblk_s::bloknum, symbol_matrix_s::bloktab, symbol_matrix_s::browmax, symbol_cblk_s::brownum, symbol_matrix_s::browtab, symbol_matrix_s::cblknbr, symbol_matrix_s::cblktab, symbol_blok_s::fcblknm, and pastix_int_t.

Referenced by extraCblkMerge(), and pastix_subtask_symbfact().

◆ pastixSymbolGetFacingBloknum()

pastix_int_t pastixSymbolGetFacingBloknum ( const symbol_matrix_t symbptr,
pastix_int_t  bloksrc,
pastix_int_t  bloknum,
pastix_int_t  startsearch,
int  ricar 
)

Search the targeted block C for a couple of blocks A and B.

When executing the simulation run to map the data on the cores, it requires to compute dependencies of each block. In that case for each couple of blocks A (defined by bloknum), and B (defined by bloksrc), we need to find the block that will receive the contribution for the matrix-matrix product. To speedup the search, the startsearch parameter can be given to specify that the index of the block searched is larger than this parameter. It returns the index of the C block when found, or -1 if no block is facing the update. This happens, only if ILU(k) factorization is applied and ricar is set to true. Otherwise, returning -1 is an error.

Parameters
[in]symbptrThe pointer to the structure that needs to be reallocated.
[in]bloksrcIndex of the block used as A in teh matrix-matrix product.
[in]bloknumIndex of the block used as B in teh matrix-matrix product.
[in]startsearchIndex of a block belonging to the facing cblk of B used as a starting point to find the index of C.
[in]ricarBooleen to enable ILU(k) factorization or not
Return values
-1No block where found. This is an error if ricar is disabled, otherwise it means that the path is longer than the k parameter in the ILU(k) factorization.
iThe index of the block C that will receive the contribution from A * B^t

We found the first block that matches

Definition at line 232 of file symbol.c.

References symbol_cblk_s::bloknum, symbol_matrix_s::bloktab, symbol_matrix_s::cblktab, symbol_blok_s::fcblknm, symbol_blok_s::frownum, symbol_blok_s::lrownum, and pastix_int_t.

Referenced by simu_computeBlockCtrbNbr(), and simu_computeTask().

◆ pastixSymbolPrint()

void pastixSymbolPrint ( const symbol_matrix_t symbptr,
FILE *  file 
)

Print the given block matrix structure in human readable format.

Parameters
[in,out]symbptrThe symbolic matrix structure to write.
[in,out]fileThe file stream to which to write the structure.

Definition at line 194 of file symbol_io.c.

References symbol_cblk_s::bloknum, symbol_matrix_s::bloktab, symbol_matrix_s::cblknbr, symbol_matrix_s::cblktab, symbol_cblk_s::fcolnum, symbol_blok_s::frownum, symbol_cblk_s::lcolnum, symbol_blok_s::lrownum, and pastix_int_t.