PaStiX Handbook  6.3.2
solver.h
Go to the documentation of this file.
1 /**
2  * @file solver.h
3  *
4  * PaStiX solver structure header.
5  *
6  * @copyright 2004-2023 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
7  * Univ. Bordeaux. All rights reserved.
8  *
9  * @version 6.3.2
10  * @author David Goudin
11  * @author Pascal Henon
12  * @author Francois Pellegrini
13  * @author Pierre Ramet
14  * @author Mathieu Faverge
15  * @author Xavier Lacoste
16  * @author Esragul Korkmaz
17  * @author Gregoire Pichon
18  * @author Tony Delarue
19  * @date 2023-10-23
20  *
21  **/
22 #ifndef _solver_h_
23 #define _solver_h_
24 
25 #ifndef DOXYGEN_SHOULD_SKIP_THIS
26 struct blendctrl_s;
27 typedef struct blendctrl_s BlendCtrl;
28 
29 struct simuctrl_s;
30 typedef struct simuctrl_s SimuCtrl;
31 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
32 
33 #include "pastix_lowrank.h"
34 
35 /**
36  * @name Cblk properties
37  * @{
38  * The type and structure definitions.
39  * Define the mask for the cblks in the cblktype field:
40  * - 1st bit: The cblk is a local fanin accumulating the contributions for a remote cblk
41  * - 2nd bit: The cblk is stored in a 2D layout fashion as in a tiled matrix, otherwise the standard 1D lapack layout is used
42  * - 3rd bit: The cblk generates 2D granularity tasks, instead of a single 1D tasks that perform factorization, solves and updates
43  * - 4th bit: The cblk is compressed in Low-Rank (implies CBLK_LAYOUT_2D), otherwise it is stored in dense
44  * - 5th bit: The cblk is part of the Schur complement if set
45  * - 6th bit: The cblk is a local copy of a remote fanin to be received.
46  */
47 #ifndef DOXYGEN_SHOULD_SKIP_THIS
48 #define CBLK_FANIN (1 << 0)
49 #define CBLK_LAYOUT_2D (1 << 1)
50 #define CBLK_TASKS_2D (1 << 2)
51 #define CBLK_COMPRESSED (1 << 3)
52 #define CBLK_IN_SCHUR (1 << 4)
53 #define CBLK_IN_LAST (1 << 5)
54 #define CBLK_RECV (1 << 6)
55 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
56 
57 /**
58  *@}
59  */
60 
61 /*
62  * The type and structure definitions.
63  */
64 #ifndef DOXYGEN_SHOULD_SKIP_THIS
65 #define COMP_1D 0
66 #define DIAG 1
67 #define E1 2
68 #define E2 3
69 #define DRUNK 4
70 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
71 
72 /**
73  * @brief Tags used in MPI communications.
74  */
75 typedef enum solve_step_e {
76  PastixSolveForward,
77  PastixSolveBackward,
78  PastixFacto,
80 
81 /**
82  * @brief Arguments for the solve.
83  */
84 typedef struct args_solve_s
85 {
86  solve_step_t solve_step;
87  pastix_solv_mode_t mode;
88  pastix_side_t side;
89  pastix_uplo_t uplo;
90  pastix_trans_t trans;
91  pastix_diag_t diag;
93 
94 /**
95  * @brief Solver recv block structure.
96  */
97 typedef struct solver_blok_recv_s {
98  pastix_int_t frownum; /**< First row contributed in the block (lrownum+1 of the original block if none) */
99  pastix_int_t lrownum; /**< Last row contributed in the block (frownum-1 of the original block if none) */
101 
102 /**
103  * @brief Solver recv column block structure.
104  *
105  * The structure works as a list of cblk for each remote owner to compute the exact contributions.
106  */
107 typedef struct solver_cblk_recv_s {
108  struct solver_cblk_recv_s *next;
109  pastix_int_t ownerid;
110  pastix_int_t fcolnum; /**< First column contributed in the block (lcolnum+1 of the original cblk if none) */
111  pastix_int_t lcolnum; /**< Last column contributed in the block (fcolnum-1 of the original cblk if none) */
112  solver_blok_recv_t bloktab[1]; /**< Array of reception blocks */
114 
115 /**
116  * @brief The task structure for the numerical factorization
117  */
118 typedef struct task_s {
119  pastix_int_t taskid; /**< COMP_1D DIAG E1 E2 */
120  pastix_int_t prionum; /**< Priority value for the factorization */
121  pastix_int_t cblknum; /**< Attached column block */
122  pastix_int_t bloknum; /**< Attached block */
123  pastix_int_t volatile ctrbcnt; /**< Total number of contributions */
124 #if defined(PASTIX_DYNSCHED)
125  int threadid;/**< Index of the bubble which contains the task */
126 #endif
128 
129 #ifndef DOXYGEN_SHOULD_SKIP_THIS
130 #define GPUID_UNDEFINED -2 /**< GPU still undefined */
131 #define GPUID_NONE -1 /**< Block not computed on GPU */
132 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
133 
134 /**
135  * @brief Solver block structure.
136  */
137 typedef struct solver_blok_s {
138  void *handler[2]; /**< Runtime data handler */
139  pastix_int_t lcblknm; /**< Local column block */
140  pastix_int_t fcblknm; /**< Facing column block */
141  pastix_int_t gbloknm; /**< Index in global bloktab (UNUSED) */
142  pastix_int_t frownum; /**< First row index */
143  pastix_int_t lrownum; /**< Last row index (inclusive) */
144  pastix_int_t coefind; /**< Index in coeftab */
145  pastix_int_t browind; /**< Index in browtab */
146  int8_t inlast; /**< Index of the block among last separator (2), coupling with last separator (1) or other blocks (0) */
147  int iluklvl; /**< The block ILU(k) level */
148 
149  /* LR structures */
150  pastix_lrblock_t *LRblock[2]; /**< Store the blok (L/U) in LR format. Allocated for the cblk. */
152 
153 /**
154  * @brief Solver column block structure.
155  */
156 typedef struct solver_cblk_s {
157  pastix_atomic_lock_t lock; /**< Lock to protect computation on the cblk */
158  volatile int32_t ctrbcnt; /**< Number of contribution to receive */
159  int8_t cblktype; /**< Type of cblk */
160  int8_t partitioned; /**< Bitmask to know if one of the side has been partitioned or not */
161  pastix_int_t fcolnum; /**< First column index (Global numbering) */
162  pastix_int_t lcolnum; /**< Last column index (Global numbering, inclusive) */
163  SolverBlok *fblokptr; /**< First block in column (diagonal) */
164  pastix_int_t stride; /**< Column block stride */
165  pastix_int_t lcolidx; /**< First column index (Local numbering), used for the rhs vectors */
166  pastix_int_t brownum; /**< First block in row facing the diagonal block in browtab, 0-based */
167  pastix_int_t brown2d; /**< First 2D-block in row facing the diagonal block in browtab, 0-based */
168  pastix_int_t sndeidx; /**< Global index of the original supernode the cblk belongs to */
169  pastix_int_t gcblknum; /**< Global column block index */
170  pastix_int_t bcscnum; /**< Index in the bcsctab if local cblk, -1 otherwise (FANIN | RECV) */
171  void *lcoeftab; /**< Coefficients access vector, lower part */
172  void *ucoeftab; /**< Coefficients access vector, upper part */
173  void *handler[2]; /**< Runtime data handler */
174  pastix_int_t selevtx; /**< Index to identify selected cblk for which intra-separator contributions are not compressed */
175  int ownerid; /**< Rank of the owner */
176  int threadid; /**< Rank of the accessing thread */
177  pastix_int_t priority; /**< Priority of the cblk */
179 
180 
181 #ifndef DOXYGEN_SHOULD_SKIP_THIS
184 
187 
190 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
191 
192 /**
193  * @brief Solver column block structure.
194  *
195  * This structure stores all the numerical information about the factorization,
196  * as well as the structure of the problem. Only local information to each
197  * process is stored in this structure.
198  *
199  */
201  int restore; /**< Flag to indicate if it is require to restore data with
202  solverBackupRestore: 0: No need, 1:After solve,
203  2:After Factorization */
204  pastix_int_t baseval; /**< Base value for numberings */
205  pastix_int_t nodenbr; /**< Number of local columns (after dof extension) */
206  pastix_int_t coefnbr; /**< Number of locally stored values (after dof extension) */
207  pastix_int_t gcblknbr; /**< Global number of column blocks */
208  pastix_int_t cblknbr; /**< Local number of column blocks */
209  pastix_int_t faninnbr; /**< Local number of fanin cblk (included in cblknbr) */
210  pastix_int_t fanincnt; /**< Number of sends to realize */
211  pastix_int_t maxrecv; /**< Maximum blok size for a cblk_recv */
212  pastix_int_t recvnbr; /**< Local number of recv cblk (included in cblknbr) */
213  pastix_int_t recvcnt; /**< Number of receptions to realize */
214  pastix_int_t cblkmax1d; /**< Rank of the last cblk not beeing enabled for 2D computations */
215  pastix_int_t cblkmin2d; /**< Rank of the first cblk beeing enabled for 2D computations */
216  pastix_int_t cblkmaxblk; /**< Maximum number of blocks per cblk */
217  pastix_int_t cblkschur; /**< Index of the first local cblk in Schur */
218  pastix_int_t nb2dcblk; /**< Number of 2D cblks */
219  pastix_int_t nb2dblok; /**< Number of 2D blocks */
220  pastix_int_t bloknbr; /**< Number of blocks */
221  pastix_int_t brownbr; /**< Size of the browtab array */
222  SolverCblk * restrict cblktab; /**< Array of solver column blocks [+1] */
223  SolverBlok * restrict bloktab; /**< Array of solver blocks [+1] */
224  pastix_int_t * restrict browtab; /**< Array of blocks */
225  pastix_coeftype_t flttype; /**< valtab datatype: PastixFloat, PastixDouble, PastixComplex32 or PastixComplex64 */
226  int globalalloc; /**< Boolean for global allocation of coeftab */
227 
228  pastix_int_t *gcbl2loc; /**< Array of local cblknum corresponding to gcblknum */
229 
230  pastix_lr_t lowrank; /**< Low-rank parameters */
231  pastix_factotype_t factotype; /**< General or symmetric factorization? */
232  double diagthreshold; /**< Diagonal threshold for pivoting */
233  volatile int32_t nbpivots; /**< Number of pivots during the factorization */
234 
235 #if defined(PASTIX_WITH_PARSEC)
236  parsec_sparse_matrix_desc_t *parsec_desc;
237 #endif
238 #if defined(PASTIX_WITH_STARPU)
239  starpu_sparse_matrix_desc_t *starpu_desc;
240  starpu_dense_matrix_desc_t *starpu_desc_rhs;
241 #endif
242 
243  pastix_int_t offdmax; /*+ Maximum size of the off-diagonal blocks for hetrf/sytrf temporary buffers +*/
244  pastix_int_t gemmmax; /*+ Maximum size of the GEMM update for 1d GEMM computations +*/
245  pastix_int_t blokmax; /*+ Maximum size of 2D blocks +*/
246  pastix_int_t colmax; /*+ Maximum column width in solvmtx +*/
247 
248  int clustnum; /*+ current processor number +*/
249  int clustnbr; /*+ number of processors +*/
250  pastix_int_t procnbr; /*+ Number of physical processor used +*/
251  pastix_int_t thrdnbr; /*+ Number of local computation threads +*/
252  pastix_int_t bublnbr; /*+ Number of local computation threads +*/
253  /* BubbleTree * restrict btree; /\*+ Bubbles tree +*\/ */
254 
255  Task * restrict tasktab; /*+ Task access vector +*/
256  pastix_int_t tasknbr; /*+ Number of Tasks in 1d mode +*/
257  pastix_int_t tasknbr_1dp; /*+ Number of Tasks in 1d plus mode +*/
258  pastix_int_t ** ttsktab; /*+ Task access vector by thread +*/
259  pastix_int_t * ttsknbr; /*+ Number of tasks by thread +*/
260  pastix_queue_t ** computeQueue; /*+ Queue of task to compute by thread +*/
261 
262  pastix_int_t *selevtx; /*+ Array to identify which cblk are pre-selected +*/
263 
264  MPI_Request *reqtab; /**< Array of requests for MPI asynchronous messages */
265  pastix_int_t *reqidx; /**< Array of local cblknum index corresponding to reqtab */
266  pastix_int_t reqnbr; /**< Length of the reqtab/reqidx arrays */
267  volatile int32_t reqnum; /**< Current amount of active requests (packed) */
268  pastix_atomic_lock_t reqlock; /**< Lock to access the request arrays */
269  void *rcoeftab; /**< Reception buffer for the communication */
270 
271  size_t *com_vector; /**< Matrix of communications between nodes. */
272 
273  PASTIX_Comm solv_comm; /*+ Copy of the pastix_data->inter_node_comm +*/
274 };
275 
276 /**
277  *******************************************************************************
278  *
279  * @brief Computes the current solve step.
280  *
281  *******************************************************************************
282  *
283  * @param[in] side
284  * Specify whether the off-diagonal blocks appear on the left or right
285  * in the equation. It has to be either PastixLeft or PastixRight.
286  *
287  * @param[in] uplo
288  * Specify whether the off-diagonal blocks are upper or lower
289  * triangular. It has to be either PastixUpper or PastixLower.
290  *
291  * @param[in] trans
292  * Specify the transposition used for the off-diagonal blocks. It has
293  * to be either PastixTrans or PastixConjTrans.
294  *
295  *******************************************************************************
296  *
297  * @return PastixSolveForward or PastixSolveBackward
298  *
299  *******************************************************************************/
300 static inline solve_step_t
302  pastix_uplo_t uplo,
303  pastix_trans_t trans )
304 {
305  if ( ( (side == PastixLeft) && (uplo == PastixUpper) && (trans == PastixNoTrans) ) ||
306  ( (side == PastixLeft) && (uplo == PastixLower) && (trans != PastixNoTrans) ) ||
307  ( (side == PastixRight) && (uplo == PastixUpper) && (trans != PastixNoTrans) ) ||
308  ( (side == PastixRight) && (uplo == PastixLower) && (trans == PastixNoTrans) ) )
309  {
310  return PastixSolveBackward;
311  }
312  else {
313  return PastixSolveForward;
314  }
315 }
316 
317 /**
318  * @brief Compute the number of columns in a column block.
319  * @param[in] cblk
320  * The pointer to the column block.
321  * @return The number of columns in the cblk.
322  */
323 static inline pastix_int_t
324 cblk_colnbr( const SolverCblk *cblk )
325 {
326  return cblk->lcolnum - cblk->fcolnum + 1;
327 }
328 
329 /**
330  * @brief Get the pointer to the data associated to the lower part of the cblk.
331  * @param[in] cblk
332  * The pointer to the column block.
333  * @return lcoeftab if the cblk is not compressed.
334  * fblokptr->LRblock[0] if the cblk is compressed.
335  */
336 static inline void *
337 cblk_getdataL( const SolverCblk *cblk ) {
338  return (cblk->cblktype & CBLK_COMPRESSED) ? cblk->fblokptr->LRblock[0] : cblk->lcoeftab;
339 }
340 
341 /**
342  * @brief Get the pointer to the data associated to the upper part of the cblk.
343  * @param[in] cblk
344  * The pointer to the column block.
345  * @return lcoeftab if the cblk is not compressed.
346  * fblokptr->LRblock[1] if the cblk is compressed.
347  */
348 static inline void *
349 cblk_getdataU( const SolverCblk *cblk ) {
350  return (cblk->cblktype & CBLK_COMPRESSED) ? cblk->fblokptr->LRblock[1] : cblk->ucoeftab;
351 }
352 
353 /**
354  * @brief Get the pointer to the data associated to the side part of the cblk.
355  * @param[in] cblk
356  * The pointer to the column block.
357  * @param[in] side
358  * side of the data needed.
359  * PastixLCoef for the lower part PastixUCoef for the upper part.
360  * @return [lu]coeftab if the cblk is not compressed.
361  * fblokptr->LRblock[side] if the cblk is compressed.
362  */
363 static inline void *
365  return cblk->cblktype & CBLK_COMPRESSED
366  ? cblk->fblokptr->LRblock[side]
367  : side == PastixUCoef ? cblk->ucoeftab
368  : cblk->lcoeftab;
369 }
370 
371 /**
372  * @brief Compute the number of blocks in a column block.
373  * @param[in] cblk
374  * The pointer to the column block.
375  * @return The number of blocks in the cblk including the diagonal block.
376  */
377 static inline pastix_int_t
378 cblk_bloknbr( const SolverCblk *cblk )
379 {
380  return (cblk+1)->fblokptr - cblk->fblokptr + 1;
381 }
382 
383 /**
384  * @brief Compute the number of rows of a block.
385  * @param[in] blok
386  * The pointer to the block.
387  * @return The number of rows in the block.
388  */
389 static inline pastix_int_t
390 blok_rownbr( const SolverBlok *blok )
391 {
392  return blok->lrownum - blok->frownum + 1;
393 }
394 
395 /**
396  * @brief Compute the number of rows of a contiguous block in front of the same cblk.
397  * @param[in] blok
398  * The pointer to the block.
399  * @return The number of rows in the block.
400  */
401 static inline pastix_int_t
403 {
404  pastix_int_t rownbr = blok_rownbr( blok );
405 
406  while( (blok[0].lcblknm == blok[1].lcblknm) &&
407  (blok[0].fcblknm == blok[1].fcblknm) )
408  {
409  blok++;
410  rownbr += blok_rownbr( blok );
411  }
412  return rownbr;
413 }
414 
415 /**
416  * @brief Return if a block is preselected as either part of the projection, or
417  * as a sub-diagonal block.
418  * @param[in] cblk
419  * The pointer to the cblk to which belong the block.
420  * @param[in] blok
421  * The pointer to the block.
422  * @param[in] fcbk
423  * The pointer to the facing cblk of the block.
424  * @return True is the block is preselected, false otherwise.
425  */
426 static inline int
428  const SolverBlok *blok,
429  const SolverCblk *fcbk )
430 {
431  int is_preselected = ( fcbk->selevtx );
432  int is_firstoffd = ( blok == (cblk->fblokptr + 1) );
433 
434  return ( fcbk->sndeidx == cblk->sndeidx ) && ( is_preselected | is_firstoffd );
435 }
436 
437 /**
438  * @brief Compute the number of rows of a column block.
439  * @param[in] cblk
440  * The pointer to the column block.
441  * @return The number of rows in the column block.
442  */
443 static inline pastix_int_t
444 cblk_rownbr( const SolverCblk *cblk )
445 {
446  pastix_int_t rownbr = 0;
447  SolverBlok * blok;
448  for (blok = cblk->fblokptr; blok < cblk[1].fblokptr; blok++) {
449  rownbr += blok_rownbr(blok);
450  }
451  return rownbr;
452 }
453 
454 /**
455  * @brief Task stealing method.
456  *
457  * @param[inout] solvmtx
458  * The pointer to the solverMatrix.
459  * @param[in] rank
460  * Rank of the computeQueue.
461  * @param[in] nbthreads
462  * Total amount of threads in the node.
463  * @return The concerned cblk if it exists. -1 otherwhise.
464  */
465 static inline pastix_int_t
467  int rank,
468  int nbthreads )
469 {
470  int rk = (rank + 1)%nbthreads;
471  pastix_queue_t *stoleQueue;
472  pastix_int_t cblknum = -1;
473  while( rk != rank )
474  {
475  assert( solvmtx->computeQueue[ rk ] );
476  stoleQueue = solvmtx->computeQueue[ rk ];
477  if( (cblknum = pqueuePop(stoleQueue)) != -1 ){
478  return cblknum;
479  }
480  rk = (rk + 1)%nbthreads;
481  }
482  return cblknum;
483 }
484 
485 /**
486  * @brief Check if a block is included inside another one.
487  *
488  * Indicate if a blok is included inside an other block.
489  * i.e. indicate if the row range of the first block is included in the
490  * one of the second.
491  *
492  * @param[in] blok The block that is tested for inclusion.
493  * @param[in] fblok The block that is suppose to include the first one.
494  *
495  * @retval true if the first block is included in the second one.
496  * @retval false if the first block is not included in the second one.
497  */
498 static inline int
500  const SolverBlok *fblok )
501 {
502 # if defined(NAPA_SOPALIN)
503  return (((blok->frownum >= fblok->frownum) &&
504  (blok->lrownum <= fblok->lrownum)) ||
505  ((blok->frownum <= fblok->frownum) &&
506  (blok->lrownum >= fblok->lrownum)) ||
507  ((blok->frownum <= fblok->frownum) &&
508  (blok->lrownum >= fblok->frownum)) ||
509  ((blok->frownum <= fblok->lrownum) &&
510  (blok->lrownum >= fblok->lrownum)));
511 # else
512  return ((blok->frownum >= fblok->frownum) &&
513  (blok->lrownum <= fblok->lrownum));
514 # endif /* defined(NAPA_SOPALIN) */
515 }
516 
517 void solverInit( SolverMatrix *solvmtx );
518 void solverExit( SolverMatrix *solvmtx );
519 
520 int solverMatrixGen( SolverMatrix *solvmtx,
521  const symbol_matrix_t *symbmtx,
522  const pastix_order_t *ordeptr,
523  const SimuCtrl *simuctl,
524  const BlendCtrl *ctrl,
525  PASTIX_Comm comm,
526  isched_t *isched );
527 
528 int solverMatrixGenSeq( SolverMatrix *solvmtx,
529  const symbol_matrix_t *symbmtx,
530  const pastix_order_t *ordeptr,
531  const SimuCtrl *simuctl,
532  const BlendCtrl *ctrl,
533  PASTIX_Comm comm,
534  isched_t *isched,
535  pastix_int_t is_dbg );
536 
537 int solverLoad( SolverMatrix *solvptr,
538  FILE *stream );
539 int solverSave( const SolverMatrix *solvptr,
540  FILE *stream );
541 
542 void solverRealloc( SolverMatrix *solvptr);
543 SolverMatrix *solverCopy ( const SolverMatrix *solvptr,
544  pastix_coeftype_t flttype );
545 
546 int solverCheck ( const SolverMatrix *solvmtx );
547 int solverDraw ( const SolverMatrix *solvptr,
548  FILE *stream,
549  int verbose,
550  const char *directory );
551 void solverPrintStats( const SolverMatrix *solvptr );
552 
553 void solverRequestInit( solve_step_t solve_step,
554  SolverMatrix *solvmtx );
555 void solverRequestExit( SolverMatrix *solvmtx );
556 
558  SolverMatrix *solvmtx,
559  pastix_coeftype_t flttype );
560 void solverRecvExit( SolverMatrix *solvmtx );
561 
562 void solverRhsRecvInit( solve_step_t solve_step,
563  SolverMatrix *solvmtx,
564  pastix_coeftype_t flttype,
565  pastix_rhs_t rhsb );
566 void solverRhsRecvExit( SolverMatrix *solvmtx );
567 
568 /*
569  * Solver backup
570  */
571 #ifndef DOXYGEN_SHOULD_SKIP_THIS
572 struct SolverBackup_s;
573 typedef struct SolverBackup_s SolverBackup_t;
574 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
575 
576 SolverBackup_t *solverBackupInit ( const SolverMatrix *solvmtx );
577 int solverBackupRestore( SolverMatrix *solvmtx, const SolverBackup_t *b );
578 void solverBackupExit ( SolverBackup_t *b );
579 
580 #endif /* _solver_h_ */
BEGIN_C_DECLS typedef int pastix_int_t
Definition: datatypes.h:51
The type and structure definitions.
Definition: blendctrl.h:28
static pastix_int_t pqueuePop(pastix_queue_t *q)
Pop the head of the queue whithout returning the keys.
Definition: queue.h:75
Queue structure.
Definition: queue.h:38
Control structure for the simulation.
Definition: simu.h:116
void solverRealloc(SolverMatrix *solvptr)
Realloc in a contiguous way a given solver structure.
Definition: solver_copy.c:205
int solverBackupRestore(SolverMatrix *solvmtx, const SolverBackup_t *b)
Restore initial values.
SolverMatrix * solverCopy(const SolverMatrix *solvptr, pastix_coeftype_t flttype)
Generate a copy of a solver matrix structure.
Definition: solver_copy.c:171
void solverRecvExit(SolverMatrix *solvmtx)
Free the array linked to pending reception.
Definition: solver.c:563
int solverSave(const SolverMatrix *solvptr, FILE *stream)
Save a solver matrix structure into a file.
Definition: solver_io.c:261
void solverRequestExit(SolverMatrix *solvmtx)
Free the arrays related to the requests.
Definition: solver.c:481
void solverBackupExit(SolverBackup_t *b)
Free the solver backup data structure.
int solverLoad(SolverMatrix *solvptr, FILE *stream)
Load a solver matrix structure from a file.
Definition: solver_io.c:45
int solverDraw(const SolverMatrix *solvptr, FILE *stream, int verbose, const char *directory)
Writes a PostScript picture of the low-rank solver matrix.
Definition: solver_draw.c:56
SolverBackup_t * solverBackupInit(const SolverMatrix *solvmtx)
Initialize the backup structure.
Definition: solver_backup.c:59
void solverRhsRecvExit(SolverMatrix *solvmtx)
Frees the array linked to pending reception.
Definition: solver.c:676
void solverRecvInit(pastix_coefside_t side, SolverMatrix *solvmtx, pastix_coeftype_t flttype)
Allocate the reception buffer, and initiate the first persistant reception.
Definition: solver.c:522
void solverRhsRecvInit(solve_step_t solve_step, SolverMatrix *solvmtx, pastix_coeftype_t flttype, pastix_rhs_t rhsb)
Allocates the reception buffer, and initiate the first persistant reception.
Definition: solver.c:633
void solverPrintStats(const SolverMatrix *solvptr)
Print statistical information about the solver matrix structure.
Definition: solver.c:196
void solverRequestInit(solve_step_t solve_step, SolverMatrix *solvmtx)
Instanciate the arrays for the requests according to the scheduler.
Definition: solver.c:433
void solverInit(SolverMatrix *solvmtx)
Initialize the solver structure.
Definition: solver.c:118
int solverCheck(const SolverMatrix *solvmtx)
Checks the consistency of the given solver matrix structure.
Definition: solver_check.c:54
void solverExit(SolverMatrix *solvmtx)
Free the content of the solver matrix structure.
Definition: solver.c:143
Structure to define the type of function to use for the low-rank kernels and their parameters.
The block low-rank structure to hold a matrix in low-rank form.
spm_coeftype_t pastix_coeftype_t
Arithmetic types.
Definition: api.h:294
enum pastix_diag_e pastix_diag_t
Diagonal.
enum pastix_solv_mode_e pastix_solv_mode_t
Solve Schur modes.
enum pastix_uplo_e pastix_uplo_t
Upper/Lower part.
enum pastix_side_e pastix_side_t
Side of the operation.
enum pastix_factotype_e pastix_factotype_t
Factorization algorithms available for IPARM_FACTORIZATION parameter.
enum pastix_trans_e pastix_trans_t
Transpostion.
enum pastix_coefside_e pastix_coefside_t
Data blocks used in the kernel.
@ PastixUCoef
Definition: api.h:479
@ PastixUpper
Definition: api.h:466
@ PastixLower
Definition: api.h:467
@ PastixRight
Definition: api.h:496
@ PastixLeft
Definition: api.h:495
@ PastixNoTrans
Definition: api.h:445
int solverMatrixGen(SolverMatrix *solvmtx, const symbol_matrix_t *symbmtx, const pastix_order_t *ordeptr, const SimuCtrl *simuctl, 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 *simuctl, 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
PaRSEC descriptor stucture for the sparse matrix.
Definition: pastix_parsec.h:35
StarPU descriptor for the vectors linked to a given sparse matrix.
StarPU descriptor stucture for the sparse matrix.
Symbol matrix structure.
Definition: symbol.h:77
Main PaStiX RHS structure.
Definition: pastixdata.h:150
pastix_int_t nodenbr
Definition: solver.h:205
static pastix_int_t blok_rownbr_ext(const SolverBlok *blok)
Compute the number of rows of a contiguous block in front of the same cblk.
Definition: solver.h:402
pastix_int_t reqnbr
Definition: solver.h:266
pastix_int_t lrownum
Definition: solver.h:99
pastix_int_t browind
Definition: solver.h:145
pastix_int_t nb2dcblk
Definition: solver.h:218
pastix_atomic_lock_t reqlock
Definition: solver.h:268
pastix_int_t maxrecv
Definition: solver.h:211
pastix_int_t cblkmin2d
Definition: solver.h:215
pastix_coeftype_t flttype
Definition: solver.h:225
static pastix_int_t blok_rownbr(const SolverBlok *blok)
Compute the number of rows of a block.
Definition: solver.h:390
pastix_int_t baseval
Definition: solver.h:204
void * ucoeftab
Definition: solver.h:172
pastix_int_t taskid
Definition: solver.h:119
pastix_int_t cblknum
Definition: solver.h:121
static pastix_int_t cblk_bloknbr(const SolverCblk *cblk)
Compute the number of blocks in a column block.
Definition: solver.h:378
pastix_int_t lrownum
Definition: solver.h:143
pastix_int_t gcblknbr
Definition: solver.h:207
void * handler[2]
Definition: solver.h:138
pastix_lr_t lowrank
Definition: solver.h:230
pastix_int_t brownum
Definition: solver.h:166
solve_step_e
Tags used in MPI communications.
Definition: solver.h:75
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
struct task_s Task
The task structure for the numerical factorization.
pastix_int_t nb2dblok
Definition: solver.h:219
int iluklvl
Definition: solver.h:147
static pastix_int_t stealQueue(SolverMatrix *solvmtx, int rank, int nbthreads)
Task stealing method.
Definition: solver.h:466
pastix_int_t fcblknm
Definition: solver.h:140
pastix_int_t frownum
Definition: solver.h:98
static int is_block_inside_fblock(const SolverBlok *blok, const SolverBlok *fblok)
Check if a block is included inside another one.
Definition: solver.h:499
static solve_step_t compute_solve_step(pastix_side_t side, pastix_uplo_t uplo, pastix_trans_t trans)
Computes the current solve step.
Definition: solver.h:301
MPI_Request * reqtab
Definition: solver.h:264
pastix_int_t recvcnt
Definition: solver.h:213
void * rcoeftab
Definition: solver.h:269
pastix_int_t sndeidx
Definition: solver.h:168
pastix_int_t cblknbr
Definition: solver.h:208
pastix_int_t fcolnum
Definition: solver.h:110
pastix_int_t * gcbl2loc
Definition: solver.h:228
pastix_int_t cblkmaxblk
Definition: solver.h:216
SolverBlok *restrict bloktab
Definition: solver.h:223
int globalalloc
Definition: solver.h:226
pastix_int_t cblkmax1d
Definition: solver.h:214
pastix_int_t frownum
Definition: solver.h:142
static void * cblk_getdataU(const SolverCblk *cblk)
Get the pointer to the data associated to the upper part of the cblk.
Definition: solver.h:349
static pastix_int_t cblk_rownbr(const SolverCblk *cblk)
Compute the number of rows of a column block.
Definition: solver.h:444
pastix_int_t brownbr
Definition: solver.h:221
static void * cblk_getdata(const SolverCblk *cblk, pastix_coefside_t side)
Get the pointer to the data associated to the side part of the cblk.
Definition: solver.h:364
pastix_int_t gbloknm
Definition: solver.h:141
pastix_int_t selevtx
Definition: solver.h:174
pastix_factotype_t factotype
Definition: solver.h:231
pastix_int_t lcolnum
Definition: solver.h:111
static int blok_is_preselected(const SolverCblk *cblk, const SolverBlok *blok, const SolverCblk *fcbk)
Return if a block is preselected as either part of the projection, or as a sub-diagonal block.
Definition: solver.h:427
pastix_int_t faninnbr
Definition: solver.h:209
pastix_int_t fanincnt
Definition: solver.h:210
pastix_atomic_lock_t lock
Definition: solver.h:157
struct solver_cblk_s SolverCblk
Solver column block structure.
struct solver_blok_s SolverBlok
Solver block structure.
pastix_int_t prionum
Definition: solver.h:120
struct solver_blok_recv_s solver_blok_recv_t
Solver recv block structure.
pastix_int_t coefind
Definition: solver.h:144
volatile int32_t ctrbcnt
Definition: solver.h:158
SolverBlok * fblokptr
Definition: solver.h:163
pastix_int_t bloknbr
Definition: solver.h:220
static void * cblk_getdataL(const SolverCblk *cblk)
Get the pointer to the data associated to the lower part of the cblk.
Definition: solver.h:337
pastix_int_t recvnbr
Definition: solver.h:212
pastix_lrblock_t * LRblock[2]
Definition: solver.h:150
pastix_int_t *restrict browtab
Definition: solver.h:224
pastix_int_t volatile ctrbcnt
Definition: solver.h:123
volatile int32_t reqnum
Definition: solver.h:267
pastix_int_t lcblknm
Definition: solver.h:139
solver_blok_recv_t bloktab[1]
Definition: solver.h:112
struct args_solve_s args_solve_t
Arguments for the solve.
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
SolverCblk *restrict cblktab
Definition: solver.h:222
pastix_int_t stride
Definition: solver.h:164
void * handler[2]
Definition: solver.h:173
int8_t cblktype
Definition: solver.h:159
enum solve_step_e solve_step_t
Tags used in MPI communications.
pastix_int_t * reqidx
Definition: solver.h:265
int ownerid
Definition: solver.h:175
pastix_int_t cblkschur
Definition: solver.h:217
size_t * com_vector
Definition: solver.h:271
pastix_int_t lcolnum
Definition: solver.h:162
pastix_int_t bloknum
Definition: solver.h:122
struct solver_cblk_recv_s solver_cblk_recv_t
Solver recv column block structure.
void * lcoeftab
Definition: solver.h:171
double diagthreshold
Definition: solver.h:232
int8_t partitioned
Definition: solver.h:160
volatile int32_t nbpivots
Definition: solver.h:233
pastix_int_t fcolnum
Definition: solver.h:161
pastix_int_t coefnbr
Definition: solver.h:206
Arguments for the solve.
Definition: solver.h:85
Solver recv block structure.
Definition: solver.h:97
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
The task structure for the numerical factorization.
Definition: solver.h:118