PaStiX Handbook 6.4.0
Loading...
Searching...
No Matches
solver_matrix_gen.c
Go to the documentation of this file.
1/**
2 *
3 * @file solver_matrix_gen.c
4 *
5 * PaStiX solver structure generation function.
6 *
7 * @copyright 1998-2024 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
8 * Univ. Bordeaux. All rights reserved.
9 *
10 * @version 6.4.0
11 * @author Tony Delarue
12 * @author Pascal Henon
13 * @author Pierre Ramet
14 * @author Xavier Lacoste
15 * @author Mathieu Faverge
16 * @author Gregoire Pichon
17 * @author Nolan Bredel
18 * @author Alycia Lisito
19 * @date 2024-07-05
20 *
21 **/
22#include <stdio.h>
23#include <string.h>
24#include <strings.h>
25#include <math.h>
26#include <assert.h>
27#include <sys/stat.h>
28
29#include "common.h"
30#include "symbol/symbol.h"
31#include "blend/solver.h"
32#include "elimintree.h"
33#include "cost.h"
34#include "cand.h"
35#include "pastix/order.h"
36#include "extendVector.h"
37#include "simu.h"
38#include "blendctrl.h"
40
41/**
42 *******************************************************************************
43 *
44 * @ingroup pastix_blend
45 *
46 * @brief Initialize the solver matrix structure
47 *
48 * This function takes all the global preprocessing steps: the symbol matrix
49 * and the result of the simulation step to generate the solver matrix that holds
50 * only local information of each PaStiX process.
51 *
52 *******************************************************************************
53 *
54 * @param[inout] solvmtx
55 * On entry, the allocated pointer to a solver matrix structure.
56 * On exit, this structure holds alls the local information required to
57 * perform the numerical factorization.
58 *
59 * @param[in] symbmtx
60 * The global symbol matrix structure.
61 *
62 * @param[in] ordeptr
63 * The ordering structure.
64 *
65 * @param[in] simuctrl
66 * The information resulting from the simulation that will provide the
67 * data mapping, and the order of the task execution for the static
68 * scheduling.
69 *
70 * @param[in] ctrl
71 * The blend control structure that contains extra information
72 * computed during the analyze steps and the parameters of the analyze
73 * step.
74 *
75 * @param[in] comm
76 * TODO
77 *
78 * @param[in] isched
79 * TODO
80 *
81 *******************************************************************************
82 *
83 * @retval PASTIX_SUCCESS if success.
84 * @retval PASTIX_ERR_OUTOFMEMORY if one of the malloc failed.
85 *
86 *******************************************************************************/
87int
89 const symbol_matrix_t *symbmtx,
90 const pastix_order_t *ordeptr,
91 const SimuCtrl *simuctrl,
92 const BlendCtrl *ctrl,
93 PASTIX_Comm comm,
94 isched_t *isched )
95{
97 pastix_int_t *cblklocalnum;
98 pastix_int_t *bloklocalnum;
99 pastix_int_t *tasklocalnum;
100 pastix_int_t *browtmp;
101 solver_cblk_recv_t **ftgttab = NULL;
102 pastix_int_t *faninnbr_tab = NULL;
103 pastix_int_t clustnbr = ctrl->clustnbr;
104 (void)ordeptr;
105
106 assert( symbmtx->baseval == 0 );
107
108 solverInit( solvmtx );
109
110#ifdef PASTIX_DYNSCHED
111 solvmtx->btree = ctrl->btree;
112#endif
113 solvmtx->clustnum = ctrl->clustnum;
114 solvmtx->clustnbr = clustnbr;
115 solvmtx->procnbr = ctrl->total_nbcores;
116 solvmtx->thrdnbr = ctrl->local_nbthrds;
117 solvmtx->bublnbr = ctrl->local_nbctxts;
118 solvmtx->solv_comm = comm;
119
120 /* Allocate the different local numbering arrays */
121 MALLOC_INTERN( bloklocalnum, symbmtx->bloknbr, pastix_int_t );
122 MALLOC_INTERN( cblklocalnum, symbmtx->cblknbr, pastix_int_t );
123 MALLOC_INTERN( tasklocalnum, simuctrl->tasknbr, pastix_int_t );
124 MALLOC_INTERN( ftgttab, symbmtx->cblknbr, solver_cblk_recv_t* );
125 memset( ftgttab, 0, symbmtx->cblknbr * sizeof(solver_cblk_recv_t*) );
126
127 if ( clustnbr > 1 ) {
128 /*
129 * k is the proc number
130 * - faninnbr_tab[2*k] corresponds to the first index of clbk fanin for k.
131 * - faninnbr_tab[2*k+1] corresponds to the first index of blok fanin for k.
132 * - faninnbr_tab[2*clustnbr+2*k] corresponds to the first index of clbk recv for k.
133 * - faninnbr_tab[2*clustnbr+2*k+1] corresponds to the first index of blok recv for k.
134 */
135 MALLOC_INTERN( faninnbr_tab, clustnbr * 4, pastix_int_t );
136 memset( faninnbr_tab, 0, clustnbr * 4 * sizeof( pastix_int_t ) );
137 }
138
139 /* Compute local indexes to compress the symbol information into solver */
140 solvMatGen_fill_localnums( symbmtx, simuctrl, solvmtx,
141 cblklocalnum, bloklocalnum, tasklocalnum,
142 ftgttab, faninnbr_tab );
143
144 solvmtx->cblkmin2d = solvmtx->cblknbr;
145 solvmtx->cblkschur = solvmtx->cblknbr;
146 solvmtx->gcblknbr = symbmtx->cblknbr;
147 solvmtx->gbloknbr = symbmtx->bloknbr;
148
149 /***************************************************************************
150 * Fill in the local bloktab and cblktab
151 */
152 /* Allocate the cblktab and bloktab with the computed size */
153 MALLOC_INTERN( solvmtx->cblktab, solvmtx->cblknbr+1, SolverCblk );
154 MALLOC_INTERN( solvmtx->bloktab, solvmtx->bloknbr+1, SolverBlok );
155 MALLOC_INTERN( solvmtx->browtab, solvmtx->brownbr, pastix_int_t );
156 MALLOC_INTERN( browtmp, symbmtx->browmax, pastix_int_t );
157 MALLOC_INTERN( solvmtx->gcbl2loc, symbmtx->cblknbr, pastix_int_t );
158 memset( solvmtx->gcbl2loc, 0xff, symbmtx->cblknbr * sizeof(pastix_int_t) );
159 {
160 solver_cblk_recv_t *ftgtcblk;
161 SolverCblk *solvcblk = solvmtx->cblktab;
162 SolverBlok *solvblok = solvmtx->bloktab;
163 symbol_cblk_t *symbcblk = symbmtx->cblktab;
164 Cand *candcblk = ctrl->candtab;
165 pastix_int_t nbcblk2d = 0;
166 pastix_int_t nbblok2d = 0;
167 pastix_int_t bcscidx = 0; /* Index of the local classic cblk */
168 pastix_int_t sndeidx = 0; /* Index of the current supernode in the original elimination tree */
169 pastix_int_t cblknum = 0;
170 pastix_int_t brownum = 0;
171 pastix_int_t coefnbr = 0;
172 pastix_int_t nodenbr = 0;
173
174 for ( i = 0; i < symbmtx->cblknbr; i++, symbcblk++, candcblk++ ) {
175 SolverBlok *fblokptr;
176 pastix_int_t nbcols;
177 int recvcnt = 0;
178 int tasks2D, flaglocal;
179
180 flaglocal = ( simuctrl->cblktab[i].owned ) || ( ftgttab[i] != NULL );
181 if ( !flaglocal ) {
182 continue;
183 }
184
185 tasks2D = candcblk->cblktype & CBLK_TASKS_2D;
186
187 if ( simuctrl->cblktab[i].owned ) {
188 pastix_int_t cblksize, fcolnum, lcolnum;
189
190 symbol_cblk_get_colnum( symbmtx, symbcblk, &fcolnum, &lcolnum );
191
192 /*
193 * Register reception cblk
194 */
195 ftgtcblk = ftgttab[i];
196 while( ftgtcblk != NULL ) {
197 assert( (ftgtcblk->ownerid != -1) &&
198 (ftgtcblk->ownerid != ctrl->clustnum) );
199
200 /*
201 * solvcblk->cblktype set to CBLK_RECV temporarily for the blok
202 * faninnm initialisation.
203 */
204 solvcblk->cblktype = CBLK_RECV;
205 solvblok = solvMatGen_register_remote_cblk( solvmtx, symbmtx, ftgtcblk,
206 candcblk, cblklocalnum,
207 solvcblk, solvblok,
208 cblknum, brownum, i, faninnbr_tab );
209
210 /* Initialize missing fields and set to RECV */
211 solvcblk->brown2d = brownum;
212 solvcblk->cblktype |= CBLK_RECV;
213 solvcblk->priority = -1;
214
215 /* Update colmax is necessary */
216 solvmtx->colmax = pastix_imax( solvmtx->colmax, cblk_colnbr(solvcblk) );
217
218 /* Update the maximum reception buffer size */
219 cblksize = cblk_colnbr(solvcblk) * solvcblk->stride;
220
221 if ( solvcblk->cblktype & CBLK_COMPRESSED ) {
222 cblksize += (solvblok - solvcblk->fblokptr);
223 }
224 solvmtx->maxrecv = pastix_imax( solvmtx->maxrecv, cblksize );
225
226 /* Update information about 1d/2d tasks */
227 solvMatGen_cblkIs2D( solvmtx, &nbcblk2d, &nbblok2d,
228 (solvblok - solvcblk->fblokptr),
229 tasks2D, cblknum );
230
231 /*
232 * lcolidx of the RECV block must be shifted is the
233 * reception does not start at the firs column
234 */
235 solvcblk->lcolidx = nodenbr + solvcblk->fcolnum - fcolnum;
236
237 recvcnt++;
238 solvmtx->recvcnt++;
239 cblknum++;
240
241 solvcblk->bcscnum = -( solvmtx->fanincnt + solvmtx->recvcnt );
242 solvcblk->gfaninnum = faninnbr_tab[2 * (clustnbr + ftgtcblk->ownerid)] + solvmtx->gcblknbr;
243 faninnbr_tab[2 * (clustnbr + ftgtcblk->ownerid)]++;
244 solvcblk++;
245
246 {
247 solver_cblk_recv_t *current = ftgtcblk;
248 ftgtcblk = ftgtcblk->next;
249 free( current );
250 }
251 }
252 ftgttab[i] = NULL;
253
254 /*
255 * Register the local cblk
256 */
257 solvblok = solvMatGen_register_local_cblk( symbmtx, candcblk, cblklocalnum,
258 solvcblk, solvblok,
259 cblknum, brownum, i, ctrl->clustnum );
260
261 /* Store the information for the bcsc */
262 solvcblk->bcscnum = bcscidx;
263 bcscidx++;
264
265 /* Store index for the RHS */
266 solvcblk->lcolidx = nodenbr;
267 assert( ((clustnbr > 1) && (nodenbr <= solvcblk->fcolnum)) ||
268 ((clustnbr == 1) && (nodenbr == solvcblk->fcolnum)) );
269
270 /* Update local statistics */
271 nbcols = cblk_colnbr( solvcblk );
272 nodenbr += nbcols;
273 coefnbr += solvcblk->stride * nbcols;
274 }
275 /*
276 * Register a fanin
277 */
278 else {
279 ftgtcblk = ftgttab[i];
280
281 assert( ftgtcblk != NULL );
282 assert( (ftgtcblk->ownerid != -1) &&
283 (ftgtcblk->ownerid != ctrl->clustnum) );
284
285 /*
286 * solvcblk->cblktype set to CBLK_FANIN temporarily for the blok
287 * faninnm initialisation.
288 */
289 solvcblk->cblktype = CBLK_FANIN;
290 solvblok = solvMatGen_register_remote_cblk( solvmtx, symbmtx, ftgtcblk,
291 candcblk, cblklocalnum,
292 solvcblk, solvblok,
293 cblknum, brownum, i, faninnbr_tab );
294
295 /* Set to fan-in */
296 solvcblk->cblktype |= CBLK_FANIN;
297 solvcblk->priority = -1;
298 solvmtx->fanincnt++;
299 solvcblk->bcscnum = -( solvmtx->fanincnt + solvmtx->recvcnt );
300 solvcblk->gfaninnum = faninnbr_tab[2 * ftgtcblk->ownerid] + solvmtx->gcblknbr;
301 faninnbr_tab[2 * ftgtcblk->ownerid]++;
302
303 nbcols = cblk_colnbr( solvcblk );
304
305 /* Update colmax is necessary */
306 solvmtx->colmax = pastix_imax( solvmtx->colmax, nbcols );
307
308 free( ftgtcblk );
309 ftgttab[i] = NULL;
310 }
311
312 /* Update the 1D/2D infos of the solvmtx through a cblk. */
313 solvMatGen_cblkIs2D( solvmtx, &nbcblk2d, &nbblok2d,
314 (solvblok - solvcblk->fblokptr), tasks2D, cblknum );
315
316#if defined(PASTIX_WITH_MPI)
317 if ( clustnbr > 1 ) {
318#if defined(PASTIX_BLEND_FANIN_FR)
319 if ( ( solvcblk->cblktype & CBLK_COMPRESSED ) &&
320 ( solvcblk->cblktype & ( CBLK_FANIN | CBLK_RECV ) ) ) {
321 solvcblk->cblktype &= ( ~CBLK_COMPRESSED );
322 }
323#endif
324 /* No Schur complement in distributed for the moment */
325 if( solvcblk->cblktype & CBLK_IN_SCHUR ) {
326 static int warning_schur = 1;
327 if ( warning_schur && (solvmtx->clustnum == 0) ) {
328 fprintf( stderr,
329 "Warning: Schur complement support is not yet available with multiple MPI processes\n"
330 " It is thus disabled and the factorization will be fully performed\n" );
331 warning_schur = 0;
332 }
333 solvcblk->cblktype &= (~CBLK_IN_SCHUR);
334 }
335 }
336#endif
337
338 /* Store first local cblk in Schur */
339 if ( (cblknum < solvmtx->cblkschur) &&
340 (solvcblk->cblktype & CBLK_IN_SCHUR) )
341 {
342 solvmtx->cblkschur = cblknum;
343 }
344
345 solvmtx->gcbl2loc[i] = cblknum;
346 assert( cblknum == (solvcblk - solvmtx->cblktab) );
347
348 /* Compute the original supernode in the nested dissection */
349 sndeidx = solvMatGen_supernode_index( symbcblk, solvcblk, sndeidx, ordeptr );
350
351 /*
352 * Copy browtab information
353 * In case of 2D tasks, we reorder the browtab to first store
354 * the 1D contributions, and then the 2D updates.
355 * This might also be used for low rank compression, to first
356 * accumulate small dense contributions, and then, switch to a
357 * low rank - low rank update scheme.
358 */
359 {
360 pastix_int_t brownbr;
361 brownbr = solvMatGen_reorder_browtab( symbmtx, symbcblk, solvmtx, solvcblk,
362 browtmp, cblklocalnum, bloklocalnum, brownum );
363
364 /* Diagonal bloks of CBLK_RECV are added at the end of the browtab */
365 while( recvcnt ) {
366 fblokptr = solvcblk[-recvcnt].fblokptr;
367 solvmtx->browtab[brownum + brownbr] = fblokptr - solvmtx->bloktab;
368 fblokptr->browind = brownum + brownbr;
369 brownbr++;
370
371 /* Supernode index is copied in too */
372 solvcblk[-recvcnt].sndeidx = solvcblk->sndeidx;
373 recvcnt--;
374 }
375
376 brownum += brownbr;
377
378 assert( brownum <= solvmtx->brownbr );
379 assert( solvcblk->brown2d >= solvcblk->brownum );
380 assert( solvcblk->brown2d <= solvcblk->brownum + brownbr );
381 }
382
383 cblknum++;
384 solvcblk++;
385 }
386
387 /* Add a virtual cblk to avoid side effect in the loops on cblk bloks */
388 if ( cblknum > 0 ) {
389 solvMatGen_init_cblk( solvcblk, solvblok, candcblk, symbcblk,
390 solvcblk[-1].lcolnum + 1, solvcblk[-1].lcolnum + 1,
391 brownum, 0, -1, solvmtx->clustnum );
392 solvcblk->lcolidx = nodenbr;
393 }
394
395 /* Add a virtual blok to avoid side effect in the loops on cblk bloks */
396 if ( solvmtx->bloknbr > 0 ) {
397 solvMatGen_init_blok( solvblok, symbmtx->cblknbr + 1, symbmtx->cblknbr + 1,
398 solvcblk[-1].lcolnum + 1, solvcblk[-1].lcolnum + 1,
399 0, 0, 0 );
400 }
401
402 solvmtx->nodenbr = nodenbr;
403 solvmtx->coefnbr = coefnbr;
404
405 solvmtx->nb2dcblk = nbcblk2d;
406 solvmtx->nb2dblok = nbblok2d;
407
408 assert( solvmtx->cblkmax1d + 1 >= solvmtx->cblkmin2d );
409 assert( solvmtx->brownbr == brownum );
410 assert( solvmtx->cblknbr == cblknum );
411 assert( solvmtx->faninnbr == solvmtx->fanincnt );
412 assert( solvmtx->recvnbr == solvmtx->recvcnt );
413 assert( solvmtx->bloknbr == solvblok - solvmtx->bloktab );
414 }
415 memFree_null( browtmp );
416
417 /* Fill in tasktab */
418 solvMatGen_fill_tasktab( solvmtx, isched, simuctrl,
419 tasklocalnum, cblklocalnum,
420 bloklocalnum, ctrl->clustnum, 0 );
421
422 memFree_null(cblklocalnum);
423 memFree_null(bloklocalnum);
424 memFree_null(tasklocalnum);
425 memFree_null(ftgttab);
426
427 if ( clustnbr > 1 ) {
428 memFree_null( faninnbr_tab );
429 }
430
431 /* Compute the maximum area of the temporary buffer */
432 solvMatGen_max_buffers( solvmtx );
433 solvMatGen_stats_last( solvmtx );
434
435 return PASTIX_SUCCESS;
436}
437
438/**
439 *******************************************************************************
440 *
441 * @ingroup pastix_blend
442 *
443 * @brief Initialize the solver matrix structure in sequential
444 *
445 * This function takes all the global preprocessing steps: the symbol matrix,
446 * and the result of the simulation step to generate the solver matrix for one
447 * PaStiX process.
448 *
449 *******************************************************************************
450 *
451 * @param[inout] solvmtx
452 * On entry, the allocated pointer to a solver matrix structure.
453 * On exit, this structure holds alls the local information required to
454 * perform the numerical factorization.
455 *
456 * @param[in] symbmtx
457 * The global symbol matrix structure.
458 *
459 * @param[in] ordeptr
460 * The ordering structure.
461 *
462 * @param[in] simuctrl
463 * The information resulting from the simulation that will provide the
464 * data mapping, and the order of the task execution for the static
465 * scheduling.
466 *
467 * @param[in] ctrl
468 * The blend control structure that contains extra information
469 * computed during the analyze steps and the parameters of the analyze
470 * step.
471 *
472 * @param[in] comm
473 * TODO
474 *
475 * @param[in] isched
476 * TODO
477 *
478 * @param[in] is_dbg
479 * TODO
480 *
481 *******************************************************************************
482 *
483 * @retval PASTIX_SUCCESS if success.
484 * @retval PASTIX_ERR_OUTOFMEMORY if one of the malloc failed.
485 *
486 *******************************************************************************/
487int
489 const symbol_matrix_t *symbmtx,
490 const pastix_order_t *ordeptr,
491 const SimuCtrl *simuctrl,
492 const BlendCtrl *ctrl,
493 PASTIX_Comm comm,
494 isched_t *isched,
495 pastix_int_t is_dbg )
496{
497 pastix_int_t i;
498 pastix_int_t *browtmp = 0;
499 (void)ordeptr;
500
501 assert( symbmtx->baseval == 0 );
502
503 solverInit( solvmtx );
504
505 solvmtx->clustnum = ctrl->clustnum;
506 solvmtx->clustnbr = ctrl->clustnbr;
507 solvmtx->procnbr = ctrl->total_nbcores;
508 solvmtx->thrdnbr = ctrl->local_nbthrds;
509 solvmtx->bublnbr = ctrl->local_nbctxts;
510 solvmtx->solv_comm = comm;
511
512 /* Set values computed through solvMatGen_fill_localnum in distributed */
513 solvmtx->tasknbr = simuctrl->tasknbr;
514 solvmtx->cblknbr = symbmtx->cblknbr;
515 solvmtx->bloknbr = symbmtx->bloknbr;
516 solvmtx->brownbr = symbmtx->cblktab[ solvmtx->cblknbr ].brownum
517 - symbmtx->cblktab[0].brownum;
518
519 solvmtx->cblkmin2d = solvmtx->cblknbr;
520 solvmtx->cblkschur = solvmtx->cblknbr;
521 solvmtx->gcblknbr = symbmtx->cblknbr;
522
523 /***************************************************************************
524 * Fill in the local bloktab and cblktab
525 */
526 /* Allocate the cblktab and bloktab with the computed size */
527 MALLOC_INTERN(solvmtx->cblktab, solvmtx->cblknbr+1, SolverCblk );
528 MALLOC_INTERN(solvmtx->bloktab, solvmtx->bloknbr+1, SolverBlok );
529 MALLOC_INTERN(solvmtx->browtab, solvmtx->brownbr, pastix_int_t);
530 MALLOC_INTERN(browtmp, symbmtx->browmax, pastix_int_t);
531 {
532 SolverCblk *solvcblk = solvmtx->cblktab;
533 SolverBlok *solvblok = solvmtx->bloktab;
534 symbol_cblk_t *symbcblk = symbmtx->cblktab;
535 Cand *candcblk = ctrl->candtab;
536 pastix_int_t nbcblk2d = 0;
537 pastix_int_t nbblok2d = 0;
538 pastix_int_t sndeidx = 0; /* Index of the current supernode in the original elimination tree */
539 pastix_int_t cblknum = 0;
540 pastix_int_t brownum = 0;
541 pastix_int_t coefnbr = 0;
542 pastix_int_t nodenbr = 0;
543
544 for(i=0; i<symbmtx->cblknbr; i++, symbcblk++, candcblk++)
545 {
546 pastix_int_t nbcols;
547 int tasks2D = candcblk->cblktype & CBLK_TASKS_2D;
548
549 /*
550 * Register the local cblk
551 */
552 solvblok = solvMatGen_register_local_cblk( symbmtx, candcblk, NULL,
553 solvcblk, solvblok,
554 cblknum, brownum, i,
555 simuctrl->bloktab[ symbcblk->bloknum ].ownerclust );
556
557 /* Store the information for the bcsc */
558 solvcblk->bcscnum = i;
559 solvcblk->lcolidx = solvcblk->fcolnum;
560
561 /* Update local statistics */
562 assert( nodenbr == solvcblk->fcolnum );
563 nbcols = cblk_colnbr( solvcblk );
564 nodenbr += nbcols;
565 coefnbr += solvcblk->stride * nbcols;
566
567 solvMatGen_cblkIs2D( solvmtx, &nbcblk2d, &nbblok2d,
568 (solvblok - solvcblk->fblokptr), tasks2D, cblknum );
569
570#if defined(PASTIX_WITH_MPI)
571 if ( (solvcblk->cblktype & CBLK_COMPRESSED) &&
572 (solvcblk->cblktype & (CBLK_FANIN | CBLK_RECV)) )
573 {
574 solvcblk->cblktype &= (~CBLK_COMPRESSED);
575 }
576#endif
577
578 /* Store first local cblk in Schur */
579 if ( (cblknum < solvmtx->cblkschur) &&
580 (solvcblk->cblktype & CBLK_IN_SCHUR) )
581 {
582 solvmtx->cblkschur = cblknum;
583 }
584
585 /* Compute the original supernode in the nested dissection */
586 sndeidx = solvMatGen_supernode_index( symbcblk, solvcblk, sndeidx, ordeptr );
587
588 /*
589 * Copy browtab information
590 * In case of 2D tasks, we reorder the browtab to first store
591 * the 1D contributions, and then the 2D updates.
592 * This might also be used for low rank compression, to first
593 * accumulate small dense contributions, and then, switch to a
594 * low rank - low rank update scheme.
595 */
596 {
597 pastix_int_t brownbr;
598 brownbr = solvMatGen_reorder_browtab( symbmtx, symbcblk, solvmtx, solvcblk,
599 browtmp, NULL, NULL, brownum );
600
601 brownum += brownbr;
602
603 assert( brownum <= solvmtx->brownbr );
604 assert( solvcblk->brown2d >= solvcblk->brownum );
605 assert( solvcblk->brown2d <= solvcblk->brownum + brownbr );
606 }
607 cblknum++;
608 solvcblk++;
609 }
610
611 /* Add a virtual cblk to avoid side effect in the loops on cblk bloks */
612 if ( cblknum > 0 ) {
613 solvMatGen_init_cblk( solvcblk, solvblok, candcblk, symbcblk,
614 solvcblk[-1].lcolnum + 1, solvcblk[-1].lcolnum + 1,
615 symbcblk->brownum, 0, -1, ctrl->clustnum);
616 solvcblk->lcolidx = solvcblk->fcolnum;
617 }
618
619 /* Add a virtual blok to avoid side effect in the loops on cblk bloks */
620 if ( solvmtx->bloknbr > 0 ) {
621 solvMatGen_init_blok( solvblok, symbmtx->cblknbr + 1, symbmtx->cblknbr + 1,
622 solvcblk[-1].lcolnum + 1, solvcblk[-1].lcolnum + 1,
623 0, 0, 0 );
624 }
625
626 solvmtx->nodenbr = nodenbr;
627 solvmtx->coefnbr = coefnbr;
628
629 solvmtx->nb2dcblk = nbcblk2d;
630 solvmtx->nb2dblok = nbblok2d;
631
632 assert( solvmtx->cblkmax1d+1 >= solvmtx->cblkmin2d );
633 assert( solvmtx->cblknbr == cblknum );
634 assert( solvmtx->bloknbr == solvblok - solvmtx->bloktab );
635 }
636 memFree_null( browtmp );
637
638 /*
639 * Update browind fields
640 */
641 for(i=0; i<solvmtx->brownbr; i++)
642 {
643 pastix_int_t bloknum = solvmtx->browtab[i];
644 solvmtx->bloktab[ bloknum ].browind = i;
645 }
646
647 /* Fill in tasktab */
648 solvMatGen_fill_tasktab( solvmtx, isched, simuctrl,
649 NULL, NULL, NULL, ctrl->clustnum, is_dbg );
650
651 /* Compute the maximum area of the temporary buffer */
652 solvMatGen_max_buffers( solvmtx );
653 solvMatGen_stats_last( solvmtx );
654
655 return PASTIX_SUCCESS;
656}
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
pastix_int_t clustnbr
Definition blendctrl.h:71
pastix_int_t clustnum
Definition blendctrl.h:70
pastix_int_t local_nbctxts
Definition blendctrl.h:76
Cand * candtab
Definition blendctrl.h:98
pastix_int_t local_nbthrds
Definition blendctrl.h:75
pastix_int_t total_nbcores
Definition blendctrl.h:72
The type and structure definitions.
Definition blendctrl.h:28
int8_t owned
Definition simu.h:80
int ownerclust
Definition simu.h:93
SimuBlok * bloktab
Definition simu.h:126
SimuCblk * cblktab
Definition simu.h:125
pastix_int_t tasknbr
Definition simu.h:119
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)
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.
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.
void solverInit(SolverMatrix *solvmtx)
Initialize the solver structure.
Definition solver.c:118
@ PASTIX_SUCCESS
Definition api.h:367
int solverMatrixGen(SolverMatrix *solvmtx, const symbol_matrix_t *symbmtx, const pastix_order_t *ordeptr, const SimuCtrl *simuctrl, const BlendCtrl *ctrl, PASTIX_Comm comm, isched_t *isched)
Initialize the solver matrix structure.
int solverMatrixGenSeq(SolverMatrix *solvmtx, const symbol_matrix_t *symbmtx, const pastix_order_t *ordeptr, const SimuCtrl *simuctrl, const BlendCtrl *ctrl, PASTIX_Comm comm, isched_t *isched, pastix_int_t is_dbg)
Initialize the solver matrix structure in sequential.
Order structure.
Definition order.h:47
pastix_int_t bloknbr
Definition symbol.h:80
pastix_int_t baseval
Definition symbol.h:78
symbol_cblk_t * cblktab
Definition symbol.h:83
pastix_int_t browmax
Definition symbol.h:86
pastix_int_t bloknum
Definition symbol.h:48
pastix_int_t brownum
Definition symbol.h:49
pastix_int_t cblknbr
Definition symbol.h:79
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.
Definition symbol.h:116
Symbol column block structure.
Definition symbol.h:45
Symbol matrix structure.
Definition symbol.h:77
pastix_int_t nodenbr
Definition solver.h:208
pastix_int_t browind
Definition solver.h:150
pastix_int_t nb2dcblk
Definition solver.h:222
pastix_int_t maxrecv
Definition solver.h:215
pastix_int_t cblkmin2d
Definition solver.h:219
pastix_int_t gbloknbr
Definition solver.h:225
pastix_int_t gcblknbr
Definition solver.h:210
pastix_int_t brownum
Definition solver.h:171
pastix_int_t priority
Definition solver.h:183
static pastix_int_t cblk_colnbr(const SolverCblk *cblk)
Compute the number of columns in a column block.
Definition solver.h:329
pastix_int_t brown2d
Definition solver.h:172
pastix_int_t nb2dblok
Definition solver.h:223
pastix_int_t recvcnt
Definition solver.h:217
pastix_int_t gfaninnum
Definition solver.h:176
pastix_int_t sndeidx
Definition solver.h:173
pastix_int_t cblknbr
Definition solver.h:211
pastix_int_t * gcbl2loc
Definition solver.h:234
SolverBlok *restrict bloktab
Definition solver.h:229
pastix_int_t cblkmax1d
Definition solver.h:218
pastix_int_t brownbr
Definition solver.h:227
pastix_int_t faninnbr
Definition solver.h:213
pastix_int_t fanincnt
Definition solver.h:214
SolverBlok * fblokptr
Definition solver.h:168
pastix_int_t bloknbr
Definition solver.h:224
pastix_int_t recvnbr
Definition solver.h:216
pastix_int_t *restrict browtab
Definition solver.h:230
pastix_int_t lcolidx
Definition solver.h:170
pastix_int_t bcscnum
Definition solver.h:175
SolverCblk *restrict cblktab
Definition solver.h:228
pastix_int_t stride
Definition solver.h:169
int8_t cblktype
Definition solver.h:164
pastix_int_t cblkschur
Definition solver.h:221
pastix_int_t lcolnum
Definition solver.h:167
pastix_int_t fcolnum
Definition solver.h:166
pastix_int_t coefnbr
Definition solver.h:209
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