PaStiX Handbook  6.3.2
solver_matrix_gen_utils.h
Go to the documentation of this file.
1 /**
2  *
3  * @file solver_matrix_gen_utils.h
4  *
5  * PaStiX solver structure generation functions to factorize
6  * solver_matric_gen.c .
7  *
8  * @copyright 1998-2023 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
9  * Univ. Bordeaux. All rights reserved.
10  *
11  * @version 6.3.2
12  * @author Tony Delarue
13  * @author Pascal Henon
14  * @author Pierre Ramet
15  * @author Xavier Lacoste
16  * @author Mathieu Faverge
17  * @date 2023-07-21
18  *
19  * @addtogroup blend_dev_solver
20  * @{
21  *
22  **/
23 #ifndef _solver_matrix_gen_utils_h_
24 #define _solver_matrix_gen_utils_h_
25 
26 #include "elimintree.h"
27 #include "cost.h"
28 #include "symbol/symbol.h"
29 #include "cand.h"
30 #include "blend/solver.h"
31 #include "pastix/order.h"
32 
33 /**
34  *******************************************************************************
35  *
36  * @brief Initialize a solver block.
37  *
38  *******************************************************************************
39  *
40  * @param[inout] solvblok
41  * The pointer to the solver block to initialize.
42  *
43  * @param[in] lcblknm
44  * Local column block index.
45  *
46  * @param[in] fcblknm
47  * Facing column block index.
48  *
49  * @param[in] frownum
50  * First row of the block.
51  *
52  * @param[in] lrownum
53  * Last row of the block.
54  *
55  * @param[in] stride
56  * Stride of the column block.
57  *
58  * @param[in] nbcols
59  * Number of columns of the cblk to which belong the current block.
60  *
61  * @param[in] layout2D
62  * Parameter which indicates if the cblk layout is 1D lapack or 2D tile
63  * layout.
64  *
65  ********************************************************************************/
66 static inline void
68  pastix_int_t lcblknm,
69  pastix_int_t fcblknm,
70  pastix_int_t frownum,
71  pastix_int_t lrownum,
72  pastix_int_t stride,
73  pastix_int_t nbcols,
74  pastix_int_t layout2D )
75 {
76  assert( fcblknm >= -1 );
77  assert( lcblknm >= 0 );
78  assert( (fcblknm == -1) || (lcblknm <= fcblknm) );
79  assert( frownum >= 0 );
80  assert( lrownum >= frownum );
81  assert( stride >= 0 );
82  assert( nbcols >= 0 );
83 
84  solvblok->handler[0] = NULL;
85  solvblok->handler[1] = NULL;
86  solvblok->fcblknm = fcblknm;
87  solvblok->lcblknm = lcblknm;
88  solvblok->frownum = frownum;
89  solvblok->lrownum = lrownum;
90  solvblok->coefind = layout2D ? stride * nbcols : stride;
91  solvblok->browind = -1;
92  solvblok->inlast = 0;
93  solvblok->LRblock[0] = NULL;
94  solvblok->LRblock[1] = NULL;
95 }
96 
97 /**
98  *******************************************************************************
99  *
100  * @brief Initialize a solver cblk.
101  *
102  *******************************************************************************
103  *
104  * @param[inout] solvcblk
105  * The pointer to the solver cblk to initialize.
106  *
107  * @param[in] fblokptr
108  * The pointer to the first block.
109  *
110  * @param[in] candcblk
111  * The associated cand structure to the cblk to know the type of the cblk.
112  *
113  * @param[in] symbcblk
114  * The associated symbol cblk structure to get symbol information.
115  *
116  * @param[in] fcolnum
117  * Index of the first column included in the cblk.
118  *
119  * @param[in] lcolnum
120  * Index of the last column included in the cblk.
121  *
122  * @param[in] brownum
123  * Index of the first contribution to this cblk in the browtab.
124  *
125  * @param[in] stride
126  * Stride of the cblk.
127  *
128  * @param[in] cblknum
129  * The global index of the cblk. -1 if virtual.
130  *
131  * @param[in] ownerid
132  * The owner if of the cluster that owns the main instance of the cblk.
133  *
134  *******************************************************************************/
135 static inline void
137  SolverBlok *fblokptr,
138  const Cand *candcblk,
139  const symbol_cblk_t *symbcblk,
140  pastix_int_t fcolnum,
141  pastix_int_t lcolnum,
142  pastix_int_t brownum,
143  pastix_int_t stride,
144  pastix_int_t cblknum,
145  int ownerid )
146 {
147  assert( fblokptr != NULL );
148  assert( fcolnum >= 0 );
149  assert( lcolnum >= fcolnum );
150  assert( stride >= 0 );
151  assert( brownum >= 0 );
152 
153  /* Init the cblk */
154  solvcblk->lock = PASTIX_ATOMIC_UNLOCKED;
155  solvcblk->ctrbcnt = -1;
156  solvcblk->cblktype = (cblknum == -1) ? 0 : candcblk->cblktype;
157  solvcblk->fcolnum = fcolnum;
158  solvcblk->lcolnum = lcolnum;
159  solvcblk->fblokptr = fblokptr;
160  solvcblk->stride = stride;
161  solvcblk->lcolidx = -1;
162  solvcblk->brownum = brownum;
163  solvcblk->gcblknum = cblknum;
164  solvcblk->bcscnum = -1;
165  solvcblk->selevtx = (symbcblk->selevtx == SYMBCBLK_PROJ) ? 1 : 0;
166  solvcblk->ownerid = ownerid;
167  solvcblk->lcoeftab = NULL;
168  solvcblk->ucoeftab = NULL;
169  solvcblk->handler[0] = NULL;
170  solvcblk->handler[1] = NULL;
171  solvcblk->threadid = -1;
172 }
173 
174 /**
175  *******************************************************************************
176  *
177  * @brief Register the original supernode index of the cblk.
178  *
179  * This computes the index in the original elimination tree before spliting the
180  * large supernodes.
181  *
182  *******************************************************************************
183  *
184  * @param[inout] symbcblk
185  * TODO
186  *
187  * @param[inout] solvcblk
188  * The pointer to the cblk.
189  *
190  * @param[in] sndeidx
191  * The index of the last visited supernode to reduce the complexity of
192  * the function.
193  *
194  * @param[in] ordeptr
195  * The ordering structure.
196  *
197  *******************************************************************************
198  *
199  * @return The supernode index of the given cblk. The value can be used in the
200  * future calls of the function to reduce its complexity cost.
201  *
202  *******************************************************************************/
203 static inline pastix_int_t
205  SolverCblk *solvcblk,
206  pastix_int_t sndeidx,
207  const pastix_order_t *ordeptr )
208 {
209  while ( (sndeidx < ordeptr->sndenbr ) &&
210  (ordeptr->sndetab[sndeidx+1] <= symbcblk->lcolnum) )
211  {
212  sndeidx++;
213  }
214  assert( (ordeptr->sndetab[sndeidx] <= symbcblk->fcolnum) &&
215  (ordeptr->sndetab[sndeidx+1] > symbcblk->lcolnum) );
216  solvcblk->sndeidx = sndeidx;
217 
218  /* Register the cblk as being part of the last supernode */
219  if ( solvcblk->sndeidx+1 == ordeptr->sndenbr ) {
220  solvcblk->cblktype = solvcblk->cblktype | CBLK_IN_LAST;
221  }
222 
223  return sndeidx;
224 }
225 
226 /**
227  *******************************************************************************
228  *
229  * @brief Update the 1D/2D infos of the solver matrix through a cblk.
230  *
231  *******************************************************************************
232  *
233  * @param[inout] solvmtx
234  * Pointer to the solver matrix.
235  *
236  * @param[inout] nbcblk2d
237  * Amount of 2D cblk.
238  * On exit, the number of cblk is updated if the given cblk is
239  * considered as 2D.
240  *
241  * @param[inout] nbblok2d
242  * Amount of 2D blocks.
243  * On exit, the number of blok is updated if the given cblk is
244  * considered as 2D.
245  *
246  * @param[in] nbbloks
247  * Amount blocks in the current cblk.
248  *
249  * @param[in] tasks2D
250  * Boolean which indicate if the task is 2D.
251  *
252  * @param[in] cblknum
253  * Current cblk index.
254  *
255  *******************************************************************************/
256 static inline void
258  pastix_int_t *nbcblk2d,
259  pastix_int_t *nbblok2d,
260  pastix_int_t nbbloks,
261  pastix_int_t tasks2D,
262  pastix_int_t cblknum )
263 {
264  /*
265  * 2D tasks: Compute the number of cblk split in 2D tasks, and
266  * the smallest id
267  */
268  if ( tasks2D ) {
269  solvmtx->cblkmin2d = pastix_imin( solvmtx->cblkmin2d, cblknum );
270  *nbcblk2d += 1;
271  *nbblok2d += nbbloks;
272  }
273  else {
274  solvmtx->cblkmax1d = pastix_imax( solvmtx->cblkmax1d, cblknum );
275  }
276 
277  /*
278  * Compute the maximum number of block per cblk for data
279  * structure in PaRSEC/StarPU
280  */
281  if ( cblknum >= solvmtx->cblkmin2d ) {
282  solvmtx->cblkmaxblk = pastix_imax( solvmtx->cblkmaxblk, nbbloks );
283  }
284 }
285 
286 void solvMatGen_fill_localnums( const symbol_matrix_t *symbmtx,
287  const SimuCtrl *simuctrl,
288  SolverMatrix *solvmtx,
289  pastix_int_t *cblklocalnum,
290  pastix_int_t *bloklocalnum,
291  pastix_int_t *tasklocalnum,
292  solver_cblk_recv_t **ftgttab );
293 
295  const Cand *candcblk,
296  const pastix_int_t *cblklocalnum,
297  SolverCblk *solvcblk,
298  SolverBlok *solvblok,
299  pastix_int_t lcblknm,
300  pastix_int_t brownum,
301  pastix_int_t gcblknm,
302  pastix_int_t ownerid );
303 
305  const solver_cblk_recv_t *recvcblk,
306  const Cand *candcblk,
307  const pastix_int_t *cblklocalnum,
308  SolverCblk *solvcblk,
309  SolverBlok *solvblok,
310  pastix_int_t lcblknm,
311  pastix_int_t brownum,
312  pastix_int_t gcblknm );
313 
315  const symbol_cblk_t *symbcblk,
316  SolverMatrix *solvmtx,
317  SolverCblk *solvcblk,
318  pastix_int_t *browtmp,
319  const pastix_int_t *cblklocalnum,
320  const pastix_int_t *bloklocalnum,
321  pastix_int_t brownum );
322 
324  isched_t *isched,
325  const SimuCtrl *simuctrl,
326  const pastix_int_t *tasklocalnum,
327  const pastix_int_t *cblklocalnum,
328  const pastix_int_t *bloklocalnum,
329  pastix_int_t clustnum,
330  int is_dbg );
331 
332 void solvMatGen_stats_last( SolverMatrix *solvmtx );
333 void solvMatGen_max_buffers( SolverMatrix *solvmtx );
334 
336  const symbol_matrix_t *symbmtx,
337  const symbol_cblk_t *cblk,
338  const symbol_blok_t *blok,
339  const symbol_cblk_t *fcblk,
340  int ownerid );
342  const symbol_matrix_t *symbmtx,
343  const symbol_cblk_t *cblk,
344  const symbol_blok_t *blok,
345  const symbol_cblk_t *fcblk,
346  int ownerid );
347 int solver_recv_get_bloknbr( const solver_cblk_recv_t *ftgtptr,
348  const symbol_cblk_t *symbcblk,
349  const symbol_blok_t *symbblok );
350 
351 #endif /* _solver_matrix_gen_utils_h_ */
352 
353 /**
354  *@}
355  */
BEGIN_C_DECLS typedef int pastix_int_t
Definition: datatypes.h:51
int8_t cblktype
Definition: cand.h:36
Processor candidate group to own a column blok.
Definition: cand.h:28
Control structure for the simulation.
Definition: simu.h:116
SolverBlok * solvMatGen_register_local_cblk(const symbol_matrix_t *symbmtx, const Cand *candcblk, const pastix_int_t *cblklocalnum, SolverCblk *solvcblk, SolverBlok *solvblok, pastix_int_t lcblknm, pastix_int_t brownum, pastix_int_t gcblknm, pastix_int_t ownerid)
Register a local cblk from a symbol_cblk_t structure !(Fanin|Recv)
void solvMatGen_fill_localnums(const symbol_matrix_t *symbmtx, const SimuCtrl *simuctrl, SolverMatrix *solvmtx, pastix_int_t *cblklocalnum, pastix_int_t *bloklocalnum, pastix_int_t *tasklocalnum, solver_cblk_recv_t **ftgttab)
Fill the local numbering arrays to compress the symbol information into solver.
int solver_recv_get_bloknbr(const solver_cblk_recv_t *ftgtptr, const symbol_cblk_t *symbcblk, const symbol_blok_t *symbblok)
Compute the number of valid blocks in fanin/recv cblk.
Definition: solver_recv.c:306
static pastix_int_t solvMatGen_supernode_index(const symbol_cblk_t *symbcblk, SolverCblk *solvcblk, pastix_int_t sndeidx, const pastix_order_t *ordeptr)
Register the original supernode index of the cblk.
void solver_recv_update_fanin(solver_cblk_recv_t **faninptr, const symbol_matrix_t *symbmtx, const symbol_cblk_t *cblk, const symbol_blok_t *blok, const symbol_cblk_t *fcblk, int ownerid)
Register a new contribution to a fanin cblk.
Definition: solver_recv.c:203
SolverBlok * solvMatGen_register_remote_cblk(const symbol_matrix_t *symbmtx, const solver_cblk_recv_t *recvcblk, const Cand *candcblk, const pastix_int_t *cblklocalnum, SolverCblk *solvcblk, SolverBlok *solvblok, pastix_int_t lcblknm, pastix_int_t brownum, pastix_int_t gcblknm)
Register a remote cblk from a solver_recv_cblk_t structure (Fanin|Recv)
void solver_recv_update_recv(solver_cblk_recv_t **recvptr, const symbol_matrix_t *symbmtx, const symbol_cblk_t *cblk, const symbol_blok_t *blok, const symbol_cblk_t *fcblk, int ownerid)
Register a new contribution to a recv cblk.
Definition: solver_recv.c:250
static void solvMatGen_init_cblk(SolverCblk *solvcblk, SolverBlok *fblokptr, const Cand *candcblk, const symbol_cblk_t *symbcblk, pastix_int_t fcolnum, pastix_int_t lcolnum, pastix_int_t brownum, pastix_int_t stride, pastix_int_t cblknum, int ownerid)
Initialize a solver cblk.
void solvMatGen_fill_tasktab(SolverMatrix *solvmtx, isched_t *isched, const SimuCtrl *simuctrl, const pastix_int_t *tasklocalnum, const pastix_int_t *cblklocalnum, const pastix_int_t *bloklocalnum, pastix_int_t clustnum, int is_dbg)
Fill the global tasktab array, as well as the thread ttsktab arrays.
void solvMatGen_stats_last(SolverMatrix *solvmtx)
Mark blocks if they belong to the last supernode, or if they are facing it for statistical purpose on...
static void solvMatGen_cblkIs2D(SolverMatrix *solvmtx, pastix_int_t *nbcblk2d, pastix_int_t *nbblok2d, pastix_int_t nbbloks, pastix_int_t tasks2D, pastix_int_t cblknum)
Update the 1D/2D infos of the solver matrix through a cblk.
pastix_int_t solvMatGen_reorder_browtab(const symbol_matrix_t *symbmtx, const symbol_cblk_t *symbcblk, SolverMatrix *solvmtx, SolverCblk *solvcblk, pastix_int_t *browtmp, const pastix_int_t *cblklocalnum, const pastix_int_t *bloklocalnum, pastix_int_t brownum)
Reorder the browtab from the symbol structure in a distributed way. First stock the 1D blocks and the...
void solvMatGen_max_buffers(SolverMatrix *solvmtx)
Compute the maximum area of the temporary buffers used during computation.
static void solvMatGen_init_blok(SolverBlok *solvblok, pastix_int_t lcblknm, pastix_int_t fcblknm, pastix_int_t frownum, pastix_int_t lrownum, pastix_int_t stride, pastix_int_t nbcols, pastix_int_t layout2D)
Initialize a solver block.
pastix_int_t * sndetab
Definition: order.h:57
pastix_int_t sndenbr
Definition: order.h:56
Order structure.
Definition: order.h:47
pastix_int_t fcolnum
Definition: symbol.h:46
pastix_int_t lcolnum
Definition: symbol.h:47
Symbol block structure.
Definition: symbol.h:59
Symbol column block structure.
Definition: symbol.h:45
Symbol matrix structure.
Definition: symbol.h:77
pastix_int_t browind
Definition: solver.h:145
pastix_int_t cblkmin2d
Definition: solver.h:215
void * ucoeftab
Definition: solver.h:172
pastix_int_t lrownum
Definition: solver.h:143
void * handler[2]
Definition: solver.h:138
pastix_int_t brownum
Definition: solver.h:166
pastix_int_t fcblknm
Definition: solver.h:140
pastix_int_t sndeidx
Definition: solver.h:168
pastix_int_t cblkmaxblk
Definition: solver.h:216
pastix_int_t cblkmax1d
Definition: solver.h:214
pastix_int_t frownum
Definition: solver.h:142
pastix_int_t selevtx
Definition: solver.h:174
pastix_atomic_lock_t lock
Definition: solver.h:157
pastix_int_t coefind
Definition: solver.h:144
volatile int32_t ctrbcnt
Definition: solver.h:158
SolverBlok * fblokptr
Definition: solver.h:163
pastix_lrblock_t * LRblock[2]
Definition: solver.h:150
pastix_int_t lcblknm
Definition: solver.h:139
pastix_int_t gcblknum
Definition: solver.h:169
int threadid
Definition: solver.h:176
pastix_int_t lcolidx
Definition: solver.h:165
int8_t inlast
Definition: solver.h:146
pastix_int_t bcscnum
Definition: solver.h:170
pastix_int_t stride
Definition: solver.h:164
void * handler[2]
Definition: solver.h:173
int8_t cblktype
Definition: solver.h:159
int ownerid
Definition: solver.h:175
pastix_int_t lcolnum
Definition: solver.h:162
void * lcoeftab
Definition: solver.h:171
pastix_int_t fcolnum
Definition: solver.h:161
Solver block structure.
Definition: solver.h:137
Solver recv column block structure.
Definition: solver.h:107
Solver column block structure.
Definition: solver.h:156
Solver column block structure.
Definition: solver.h:200