PaStiX Handbook  6.3.2
order_amalgamate.c
Go to the documentation of this file.
1 /**
2  *
3  * @file order_amalgamate.c
4  *
5  * PaStiX amalgamation main routine to apply after ordering strategies that do
6  * not provide supernodes.
7  *
8  * @copyright 2004-2023 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
9  * Univ. Bordeaux. All rights reserved.
10  *
11  * @version 6.3.2
12  * @author Pascal Henon
13  * @author Mathieu Faverge
14  * @date 2023-07-21
15  *
16  **/
17 #include "common.h"
18 #include "graph/graph.h"
19 #include "order/order_internal.h"
20 #include "symbol/fax_csr.h"
21 
22 /**
23  *******************************************************************************
24  *
25  * @ingroup pastix_order
26 
27  * @brief Update the order structure with an amalgamation algorithm.
28  *
29  * This algorithm almagamates small blocks together in larger ones, and generates
30  * an updated order structure. See after for the different parameters.
31  *
32  *******************************************************************************
33  *
34  * @param[in] verbose
35  * Adjust the level of verbosity of the function
36  *
37  * @param[in] ilu
38  * - 1: incomplete factorization will be performed.
39  * - 0: direct factorization will be performed.
40  *
41  * @param[in] levelk
42  * Unused if ilu == 0.
43  * - k >= 0: ordering for ILU(k) factorization will be generated.
44  * - < 0: ordering for direct factorization will be generated.
45  *
46  * @param[in] rat_cblk
47  * Must be >= 0. Fill ratio that limits the amalgamation process based
48  * on the graph structure.
49  *
50  * @param[in] rat_blas
51  * Must be >= rat_cblk. Fill ratio that limits the amalgamation process
52  * that merges blocks in order to reduce the BLAS computational time
53  * (see amalgamate() for further informations).
54  *
55  * @param[inout] csc
56  * The original csc for which the symbol matrix needs to be generated.
57  * Rebase to C numbering on exit.
58  *
59  * @param[inout] orderptr
60  * The oder structure that contains the perm and invp array generated
61  * by the ordering step. The supernode partition might be initialized
62  * or not.
63  * On exit, it is rebased to c numbering and contains the updated
64  * perm/invp arrays as well as the supernode partition.
65  *
66  * @param[in] pastix_comm
67  * The PaStiX instance communicator.
68  *
69  *******************************************************************************
70  *
71  * @retval PASTIX_SUCCESS on success.
72  * @retval PASTIX_ERR_ALLOC if allocation went wrong.
73  * @retval PASTIX_ERR_BADPARAMETER if incorrect parameters are given.
74  *
75  *******************************************************************************/
76 int
77 orderAmalgamate( int verbose,
78  int ilu,
79  int levelk,
80  int rat_cblk,
81  int rat_blas,
82  pastix_graph_t *csc,
83  pastix_order_t *orderptr,
84  PASTIX_Comm pastix_comm )
85 {
86  fax_csr_t graphPA, graphL;
87  pastix_int_t n;
88  pastix_int_t nnzA, nnzL;
89  Clock timer;
90  int procnum;
91 
92  MPI_Comm_rank( pastix_comm, &procnum );
93 
94  /* Check parameters correctness */
95  if ( ( ilu == 0 ) || ( levelk < 0 ) ) {
96  /* Forces levelk to -1 */
97  levelk = -1;
98  }
99  if ( csc == NULL ) {
100  pastix_print_warning( "orderAmalgamate: wrong parameter csc" );
102  }
103  if ( orderptr == NULL ) {
104  pastix_print_warning( "orderAmalgamate: wrong parameter orderptr" );
106  }
107 
108  /* Convert Fortran to C numbering if not already done */
109  graphBase( csc, 0 );
110  pastixOrderBase( orderptr, 0 );
111 
112  n = csc->n;
113 
114  /* Create the graph of P A */
115  faxCSRGenPA( csc, orderptr->permtab, &graphPA );
116 
117  if ( verbose > PastixVerboseYes ) {
118  pastix_print( procnum,
119  0,
120  "Level of fill = %ld\n"
121  "Amalgamation ratio: cblk = %d, blas = %d\n",
122  (long)levelk,
123  rat_cblk,
124  rat_blas );
125  }
126 
127  /*
128  * Compute the graph of the factorized matrix L before amalgamation
129  * and the associated treetab and nodetab
130  */
131  memset( &graphL, 0, sizeof( fax_csr_t ) );
132  /* Direct Factorization */
133  if ( ( ilu == 0 ) || ( levelk == -1 ) )
134  {
135  clockStart( timer );
136  nnzL = faxCSRFactDirect( &graphPA, orderptr, &graphL );
137  clockStop( timer );
138 
139  if ( verbose > PastixVerboseYes ) {
140  pastix_print( procnum,
141  0,
142  "Time to compute scalar symbolic direct factorization %.3g s\n",
143  clockVal( timer ) );
144  }
145  }
146  /* ILU(k) Factorization */
147  else
148  {
149  clockStart( timer );
150  nnzL = faxCSRFactILUk( &graphPA, orderptr, levelk, &graphL );
151  clockStop( timer );
152 
153  if ( verbose > PastixVerboseYes ) {
154  pastix_print( procnum,
155  0,
156  "Time to compute scalar symbolic factorization of ILU(%ld) %.3g s\n",
157  (long)levelk,
158  clockVal( timer ) );
159  }
160  }
161 
162  nnzA = ( faxCSRGetNNZ( &graphPA ) + n ) / 2;
163  faxCSRClean( &graphPA );
164 
165  if ( verbose > PastixVerboseYes ) {
166  pastix_print( procnum,
167  0,
168  "Scalar nnza = %ld nnzlk = %ld, fillrate0 = %.3g \n",
169  (long)nnzA,
170  (long)nnzL,
171  (double)nnzL / (double)nnzA );
172  }
173 
174  /*
175  * Amalgamate the blocks
176  */
177  clockStart( timer );
178  faxCSRAmalgamate( ilu,
179  (double)rat_cblk / 100.,
180  (double)rat_blas / 100.,
181  &graphL,
182  orderptr,
183  pastix_comm );
184  clockStop( timer );
185 
186  faxCSRClean( &graphL );
187 
188  if ( verbose > PastixVerboseYes ) {
189  pastix_print( procnum,
190  0,
191  "Time to compute the amalgamation of supernodes %.3g s\n",
192  clockVal( timer ) );
193  pastix_print( procnum,
194  0,
195  "Number of cblk in the amalgamated symbol matrix = %ld \n",
196  (long)orderptr->cblknbr );
197  }
198 
199  return PASTIX_SUCCESS;
200 }
BEGIN_C_DECLS typedef int pastix_int_t
Definition: datatypes.h:51
@ PastixVerboseYes
Definition: api.h:222
@ PASTIX_SUCCESS
Definition: api.h:367
@ PASTIX_ERR_BADPARAMETER
Definition: api.h:374
void graphBase(pastix_graph_t *graph, pastix_int_t baseval)
Rebase the graph to the given value.
Definition: graph.c:102
pastix_int_t * permtab
Definition: order.h:51
pastix_int_t cblknbr
Definition: order.h:50
void pastixOrderBase(pastix_order_t *ordeptr, pastix_int_t baseval)
This routine sets the base of the given ordering structure to the given base value.
Definition: order.c:322
int orderAmalgamate(int verbose, int ilu, int levelk, int rat_cblk, int rat_blas, pastix_graph_t *csc, pastix_order_t *orderptr, PASTIX_Comm pastix_comm)
Update the order structure with an amalgamation algorithm.
Order structure.
Definition: order.h:47
pastix_int_t faxCSRFactILUk(const fax_csr_t *graphA, const pastix_order_t *order, pastix_int_t level, fax_csr_t *graphL)
Compute the non zero pattern of the levelized incomplete factor for a sparse lower triangular matrix ...
Definition: fax_csr_iluk.c:57
int faxCSRGenPA(const pastix_graph_t *graphA, const pastix_int_t *perm, fax_csr_t *graphPA)
Generate the graph of P*A from the graph of A and the permutation vector.
Definition: fax_csr.c:183
pastix_int_t faxCSRGetNNZ(const fax_csr_t *csr)
Computes the number of non zero entries in the graph.
Definition: fax_csr.c:99
void faxCSRClean(fax_csr_t *csr)
Free the data store in the structure.
Definition: fax_csr.c:65
void faxCSRAmalgamate(int ilu, double rat_cblk, double rat_blas, fax_csr_t *graphL, pastix_order_t *order, PASTIX_Comm pastix_comm)
Amalgamate the given graph up to a given tolerance.
pastix_int_t faxCSRFactDirect(const fax_csr_t *graphA, const pastix_order_t *order, fax_csr_t *graphL)
Compute the non zero pattern of the direct factorization of a matrix A, given the supernode partition...
Fax blocked csr structure.
Definition: fax_csr.h:25