PaStiX Handbook  6.3.2

Functions to generate and manipulate the graph structure. More...

Graph basic subroutines

int graphPrepare (pastix_data_t *pastix_data, const spmatrix_t *spm, pastix_graph_t **graph)
 This routine initializes the graph. More...
 
void graphBase (pastix_graph_t *graph, pastix_int_t baseval)
 Rebase the graph to the given value. More...
 
void graphExit (pastix_graph_t *graph)
 Free the content of the graph structure. More...
 
void graphInitEmpty (pastix_graph_t *graph, PASTIX_Comm comm)
 Initialize an empty graph. More...
 

Graph IO subroutines

void graphLoad (const pastix_data_t *pastix_data, pastix_graph_t *graph)
 Load a graph from a file. More...
 
void graphSave (pastix_data_t *pastix_data, const pastix_graph_t *graph)
 Save a graph to file. More...
 

Graph manipulation subroutines

int graphCopy (pastix_graph_t *graphdst, const pastix_graph_t *graphsrc)
 This routine copy a given ordering in a new one. More...
 
int graphSpm2Graph (pastix_graph_t *graph, const spmatrix_t *spm)
 This routine build a graph thanks to an spm;. More...
 
void graphSort (pastix_graph_t *graph)
 This routine sortes the subarray of edges of each vertex. More...
 
void graphNoDiag (pastix_graph_t *graph)
 This routine removes the diagonal edges from a centralized graph. More...
 
int graphSymmetrize (pastix_graph_t *graph)
 Symmetrize a given graph. More...
 
int graphUpdateComputedFields (pastix_graph_t *graph)
 Update dofs, nnz, nnzexp, gnnz, n, nexp, gN of a given graph. More...
 
int graphScatterInPlace (pastix_graph_t *graph, PASTIX_Comm comm)
 This routine scatter a graph from node root to the other nodes. More...
 
int graphGatherInPlace (pastix_graph_t *graph)
 This routine gather a distributed graph on each node in place. More...
 
int graphIsolate (const pastix_graph_t *ingraph, pastix_graph_t *outgraph, pastix_int_t isolate_n, pastix_int_t *isolate_list, pastix_int_t **new_perm, pastix_int_t **new_invp)
 Isolate a subset of vertices from a given graph. More...
 
int graphIsolateRange (const pastix_graph_t *graph, const pastix_order_t *order, pastix_graph_t *out_graph, pastix_int_t fnode, pastix_int_t lnode, pastix_int_t distance)
 Isolate the subgraph associated to a range of unknowns in the permuted graph. More...
 
void graphComputeProjection (const pastix_graph_t *graph, const int *vertlvl, pastix_order_t *order, const pastix_graph_t *subgraph, pastix_order_t *suborder, pastix_int_t fnode, pastix_int_t lnode, pastix_int_t sn_level, pastix_int_t distance, pastix_int_t maxdepth, pastix_int_t maxwidth, pastix_int_t *depthsze)
 TODO. More...
 
pastix_int_t graphIsolateConnectedComponents (const pastix_graph_t *graph, pastix_int_t *comp_vtx, pastix_int_t *comp_sze)
 Isolate the connected components from the original graph. More...
 
int graphComputeKway (const pastix_graph_t *graph, pastix_order_t *order, pastix_int_t *peritab, pastix_int_t *comp_nbr, pastix_int_t *comp_sze, pastix_int_t *comp_vtx, pastix_int_t comp_id, pastix_int_t nbpart)
 Compute the K-way partition of a given set of unknowns. More...
 
pastix_int_tgraphGetWeights (const pastix_graph_t *graph)
 Build the vertex weight array out of the dof array. More...
 

Detailed Description

Functions to generate and manipulate the graph structure.

This module provides the set of function to prepare the graph structure associated to a given sparse matrix. It is possible to symmetrize a graph, to extract a subgraph and to apply a new permutation.

Function Documentation

◆ graphPrepare()

int graphPrepare ( pastix_data_t pastix_data,
const spmatrix_t *  spm,
pastix_graph_t **  graph 
)

This routine initializes the graph.

This routine will also symmetrize the graph, remove duplicates, remove the diagonal elements, and keep only the lower part.

Parameters
[in,out]pastix_dataThe pointer to the solver instance. On exit, the fields n, cols, rows and loc2glob are initialized for future steps of the solver.
[in]spmThe initial user spm that needs to be transformed in a correct graph for future call in ordering and symbolic factorization routines.
[out]graphOn exit, the pointer to the allocated graph structure is returned. The graph can then be used with ordering and symbol factorization tools. The graph is symmetrized without diagonal elements and rows array is sorted.
Return values
PASTIX_SUCCESSon success,
!0on failure.

Definition at line 115 of file graph_prepare.c.

References graphLoad(), graphNoDiag(), graphSort(), graphSpm2Graph(), graphSymmetrize(), pastix_data_s::iparm, IPARM_IO_STRATEGY, IPARM_VERBOSE, pastix_int_t, PASTIX_SUCCESS, PastixIOLoadGraph, and PastixVerboseNo.

Referenced by pastix_subtask_order().

◆ graphBase()

void graphBase ( pastix_graph_t *  graph,
pastix_int_t  baseval 
)

Rebase the graph to the given value.

Parameters
[in,out]graphThe graph to rebase.
[in]basevalThe base value to use in the graph (0 or 1).

Definition at line 102 of file graph.c.

Referenced by graphIsolateRange(), ocpts_graph_init(), ocs_graph_init(), orderAmalgamate(), orderComputePTScotch(), pastix_subtask_order(), and pastix_subtask_symbfact().

◆ graphExit()

void graphExit ( pastix_graph_t *  graph)

Free the content of the graph structure.

Parameters
[in,out]graphThe pointer graph structure to free.

Definition at line 73 of file graph.c.

Referenced by graphComputeKway(), graphGatherInPlace(), graphScatterInPlace(), orderSupernodes(), pastix_subtask_blend(), pastix_subtask_order(), and pastixFinalize().

◆ graphInitEmpty()

void graphInitEmpty ( pastix_graph_t *  graph,
PASTIX_Comm  comm 
)

Initialize an empty graph.

Parameters
[in,out]graphThe empty graph to init.
[in]commThe MPI communicator used for the graph.

Definition at line 48 of file graph.c.

References graphUpdateComputedFields(), and pastix_int_t.

Referenced by graphIsolate().

◆ graphLoad()

void graphLoad ( const pastix_data_t pastix_data,
pastix_graph_t *  graph 
)

Load a graph from a file.

The file is named 'graphname' in the local directory.

Parameters
[in]pastix_dataThe pointer to the solver instance to get options as rank, communicators, ...
[in,out]graphThe graph structure to store the loaded graph. The graph is read from the file named by the environment variable PASTIX_FILE_GRAPH, and if PASTIX_FILE_GRAPH is not defined, the default filename "graphname" in the local directory is used.

Definition at line 45 of file graph_io.c.

Referenced by graphPrepare().

◆ graphSave()

void graphSave ( pastix_data_t pastix_data,
const pastix_graph_t *  graph 
)

Save a graph to file.

The file is named 'graphgen' in the local directory.

Parameters
[in]pastix_dataThe pointer to the solver instance to get options as rank, communicators, ...
[in]graphThe graph structure to store the loaded graph. The graph is written to the file named by the environment variable PASTIX_FILE_GRAPH, and if PASTIX_FILE_GRAPH is not defined, the default filename "graphname" in the local directory is used.

Definition at line 102 of file graph_io.c.

References pastix_data_s::dir_local, pastix_fname(), and pastix_gendirectories().

◆ graphCopy()

int graphCopy ( pastix_graph_t *  graphdst,
const pastix_graph_t *  graphsrc 
)

This routine copy a given ordering in a new one.

This function copies a graph structure into another one. If all subpointers are NULL, then they are all allocated and contains the original graphsrc values on exit. If one or more array pointers are not NULL, then, only those are copied to the graphdst structure.

Parameters
[out]graphdstThe destination graph. Must be allocated on entry.
[in]graphsrcThe source graph
Return values
PASTIX_SUCCESSon successful exit
PASTIX_ERR_BADPARAMETERif one parameter is incorrect.

Definition at line 151 of file graph.c.

References PASTIX_ERR_BADPARAMETER, and PASTIX_SUCCESS.

Referenced by graphIsolate(), and graphIsolateRange().

◆ graphSpm2Graph()

int graphSpm2Graph ( pastix_graph_t *  graph,
const spmatrix_t *  spm 
)

This routine build a graph thanks to an spm;.

This function copies an spm structure into a graph one. If all subpointers are NULL, then they are all allocated and contains the original spm values on exit. If one or more array pointers are not NULL, then, only those are copied to the graphdst structure. We will take care that our graph does not contain coefficients, therefore has SpmPattern floating type and is a SpmCSC format type.

Parameters
[in,out]graphThe destination graph.
[in]spmThe source Sparse Matrix.
Return values
PASTIX_SUCCESSon successful exit
PASTIX_ERR_BADPARAMETERif one parameter is incorrect.

Definition at line 379 of file graph.c.

References PASTIX_ERR_BADPARAMETER, and PASTIX_SUCCESS.

Referenced by graphPrepare().

◆ graphSort()

void graphSort ( pastix_graph_t *  graph)

This routine sortes the subarray of edges of each vertex.

WARNING: The sort is always performed, do not call this routine when it is not required.

Parameters
[in,out]graphOn entry, the pointer to the graph structure. On exit, the same graph with subarrays of edges sorted by ascending order.

Definition at line 282 of file graph.c.

Referenced by graphPrepare().

◆ graphNoDiag()

void graphNoDiag ( pastix_graph_t *  graph)

This routine removes the diagonal edges from a centralized graph.

Parameters
[in,out]graphOn entry, the pointer to the graph structure with possible diagonal edges (i,i). On exit, all entries of type (i,i) are removed from the graph.

Definition at line 37 of file graph_prepare.c.

References graphUpdateComputedFields(), and pastix_int_t.

Referenced by graphPrepare().

◆ graphSymmetrize()

int graphSymmetrize ( pastix_graph_t *  graph)

Symmetrize a given graph.

Parameters
[in,out]graphThe initialized graph structure that will be symmetrized.
Return values
PASTIX_SUCCESSon success,
PASTIX_ERR_BADPARAMETERif incorrect parameters are given.

Definition at line 307 of file graph.c.

References PASTIX_ERR_BADPARAMETER, and PASTIX_SUCCESS.

Referenced by graphPrepare().

◆ graphUpdateComputedFields()

int graphUpdateComputedFields ( pastix_graph_t *  graph)

Update dofs, nnz, nnzexp, gnnz, n, nexp, gN of a given graph.

Parameters
[in,out]graphThe initialized graph structure that will be updated.
Return values
PASTIX_SUCCESSon success,
PASTIX_ERR_BADPARAMETERif incorrect parameters are given.

Definition at line 338 of file graph.c.

References PASTIX_ERR_BADPARAMETER, and PASTIX_SUCCESS.

Referenced by graph_isolate_compress(), graphInitEmpty(), graphIsolateRange(), and graphNoDiag().

◆ graphScatterInPlace()

int graphScatterInPlace ( pastix_graph_t *  graph,
PASTIX_Comm  comm 
)

This routine scatter a graph from node root to the other nodes.

Parameters
[in,out]graphOn entry, the graph to scatter. On exit, the scattered graph
[in]commMPI communicator.
Return values
1if the graph has been scattered, 0 if untouched

Definition at line 194 of file graph.c.

References graphExit().

Referenced by pastix_subtask_order().

◆ graphGatherInPlace()

int graphGatherInPlace ( pastix_graph_t *  graph)

This routine gather a distributed graph on each node in place.

Parameters
[in]graphOn entry, the distributed graph. On exit, the gathered graph.
Return values
1if the graph has been gathered, 0 if untouched

Definition at line 239 of file graph.c.

References graphExit().

Referenced by pastix_subtask_order(), and pastix_subtask_symbfact().

◆ graphIsolate()

int graphIsolate ( const pastix_graph_t *  ingraph,
pastix_graph_t *  outgraph,
pastix_int_t  isolate_n,
pastix_int_t isolate_list,
pastix_int_t **  new_perm,
pastix_int_t **  new_invp 
)

Isolate a subset of vertices from a given graph.

Return a new graph cleaned from those vertices.

Parameters
[in]ingraphTODO
[in]outgraphTODO
[in]isolate_nThe number of columns to isolate from the original graph.
[in,out]isolate_listArray of size isolate_n. List of columns to isolate. On exit, the list is sorted by ascending indexes. Must be based as the graph.
[out]new_permArray of size n-isolate_n. Contains permutation generated to isolate the columns at the end of the graph that is 0-based. If new_perm == NULL, nothing is returned, otherwise the pointer to the allocated structure.
[out]new_invpArray of size n-isolate_n. Contains the inverse permutation generated to isolate the columns at the end of the graph that is 0-based. If new_invp == NULL, nothing is returned, otherwise the pointer to the allocated structure.
Return values
PASTIX_SUCCESSon success,
PASTIX_ERR_ALLOCif allocation went wrong,
PASTIX_ERR_BADPARAMETERif incorrect parameters are given.

Definition at line 334 of file graph_isolate.c.

References graph_isolate_assign_ptr(), graph_isolate_compress(), graph_isolate_everything(), graph_isolate_permutations(), graphCopy(), graphInitEmpty(), PASTIX_ERR_BADPARAMETER, pastix_int_t, and PASTIX_SUCCESS.

Referenced by pastix_subtask_order().

◆ graphIsolateRange()

int graphIsolateRange ( const pastix_graph_t *  graph,
const pastix_order_t order,
pastix_graph_t *  out_graph,
pastix_int_t  fnode,
pastix_int_t  lnode,
pastix_int_t  distance 
)

Isolate the subgraph associated to a range of unknowns in the permuted graph.

This routine isolates a continuous subset of vertices from a given graph, and returns a new graph made of those vertices and internal connexions. Extra edges are created between vertices if they are connected through a halo at a distance d given in parameter.

Parameters
[in]graphThe original graph associated from which vertices and edges must be extracted.
[in]orderThe ordering structure associated to the graph.
[in,out]out_graphThe extracted graph. If the graph is allocated, it is freed before usage. On exit, contains the subgraph of the vertices invp[fnode] to invp[lnode-1].
[in]fnodeThe index of the first node to extract in the inverse permutation.
[in]lnodeThe index (+1) of the last node to extract in the inverse permutation.
[in]distanceDistance considered in number of edges to create an edge in isolated graph.
Return values
PASTIX_SUCCESSon success.
PASTIX_ERR_ALLOCif allocation went wrong.
PASTIX_ERR_BADPARAMETERif incorrect parameters are given.

Definition at line 540 of file graph_isolate.c.

References pastix_order_s::cblknbr, extendint_Exit(), extendint_Init(), graph_iRange_fill_outptr(), graphBase(), graphCopy(), graphUpdateComputedFields(), PASTIX_ERR_BADPARAMETER, pastix_int_t, and PASTIX_SUCCESS.

Referenced by graphComputeKway(), orderDraw(), and orderSupernodes().

◆ graphComputeProjection()

void graphComputeProjection ( const pastix_graph_t *  graph,
const int *  vertlvl,
pastix_order_t order,
const pastix_graph_t *  subgraph,
pastix_order_t suborder,
pastix_int_t  fnode,
pastix_int_t  lnode,
pastix_int_t  sn_level,
pastix_int_t  distance,
pastix_int_t  maxdepth,
pastix_int_t  maxwidth,
pastix_int_t depthsze 
)

TODO.

Parameters
[in]graphTODO
[in]vertlvlTODO
[in]orderTODO
[in]subgraphTODO
[in]suborderTODO
[in]fnodeTODO
[in]lnodeTODO
[in]sn_levelTODO
[in]distanceTODO
[in]maxdepthTODO
[in]maxwidthTODO
[in]depthszeTODO

Definition at line 68 of file graph_compute_projection.c.

References extendint_Add(), extendint_Clear(), extendint_Exit(), extendint_Init(), extendint_Read(), extendint_Size(), pastix_int_t, pastix_order_s::peritab, and pastix_order_s::permtab.

Referenced by orderSupernodes().

◆ graphIsolateConnectedComponents()

pastix_int_t graphIsolateConnectedComponents ( const pastix_graph_t *  graph,
pastix_int_t comp_vtx,
pastix_int_t comp_sze 
)

Isolate the connected components from the original graph.

Parameters
[in]graphThe original graph associated from which vertices and edges must be extracted.
[in,out]comp_vtxArray of size n that hold the index of the component for each vertex of the graph.
[in,out]comp_szeThe size of each components in the graph.
Return values
Theamount of connected components in the graph

Definition at line 114 of file graph_connected_components.c.

References pastix_int_t.

Referenced by orderSupernodes().

◆ graphComputeKway()

int graphComputeKway ( const pastix_graph_t *  graph,
pastix_order_t order,
pastix_int_t peritab,
pastix_int_t comp_nbr,
pastix_int_t comp_sze,
pastix_int_t comp_vtx,
pastix_int_t  comp_id,
pastix_int_t  nbpart 
)

Compute the K-way partition of a given set of unknowns.

Take a graph with its own associated ordering, and generates a k-way partition of the graph.

Parameters
[in]graphThe pointer to the graph structure on wich to apply the K-way partitioning on the subcomponent comp_id.
[in,out]orderThe order structure associated to the given graph. On entry, the structure is initialized through previous operation on the graph, and on exit vertices are reordered by components in the K-way partitioning.
[in,out]peritabArray of size n. Pointer to the peritab array of the global graph in which the subproblem belongs to. This array is updated according to the eprmutation performed for the K-way partitioning.
[in,out]comp_nbrOn entry the number of components in the graph. On exit, the new number of components after K-way partitioning.
[in,out]comp_szeThe size of each components in the graph. On entry, only the first comp_nbr sizes are initialized, others are 0. On exit, the array is extended to store the size of the new components generated by the k-way partitionning.
[in,out]comp_vtxArray of size n that hold the index of the component for each vertex of the graph.
[in]comp_idThe index of the component that must be split by K-way paritioning.
[in]nbpartThe number of part wanted in the k-way partitioning.
Return values
PASTIX_SUCCESSon success,
PASTIX_ERR_UNKNOWNif something went wrong with Scotch.

Definition at line 80 of file graph_compute_kway.c.

References pastix_order_s::baseval, graphExit(), graphIsolateRange(), PASTIX_ERR_BADPARAMETER, PASTIX_ERR_UNKNOWN, pastix_int_t, PASTIX_SUCCESS, pastix_order_s::peritab, and pastix_order_s::permtab.

Referenced by orderSupernodes().

◆ graphGetWeights()

pastix_int_t * graphGetWeights ( const pastix_graph_t *  graph)

Build the vertex weight array out of the dof array.

Parameters
[in]graphPointer to the graph structure.
Return values
Thevertex weight array if graph->dof != 1, NULL otherwise.

Definition at line 433 of file graph.c.

References pastix_int_t.

Referenced by ocpts_graph_init(), and ocs_graph_init().