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