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