PaStiX Handbook  6.4.0
starpu_sparse_matrix.c
Go to the documentation of this file.
1 /**
2  *
3  * @file starpu_sparse_matrix.c
4  *
5  * PaStiX sparse matrix descriptor for StarPU.
6  *
7  * @copyright 2016-2024 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
8  * Univ. Bordeaux. All rights reserved.
9  *
10  * @version 6.4.0
11  * @author Mathieu Faverge
12  * @author Pierre Ramet
13  * @author Tony Delarue
14  * @author Alycia Lisito
15  * @author Nolan Bredel
16  * @date 2024-07-05
17  *
18  * @ingroup pastix_starpu
19  * @{
20  *
21  **/
22 #include "common.h"
23 #include "blend/solver.h"
24 #include "pastix_starpu.h"
25 #include "pastix_zcores.h"
26 #include <starpu_data.h>
27 
28 #ifndef DOXYGEN_SHOULD_SKIP_THIS
29 static inline void
30 pastix_starpu_filter_interface( void *father_interface,
31  void *child_interface,
32  struct starpu_data_filter *f,
33  unsigned id,
34  unsigned nchunks )
35 {
36  pastix_starpu_interface_t *father = (pastix_starpu_interface_t *)father_interface;
37  pastix_starpu_interface_t *child = (pastix_starpu_interface_t *)child_interface;
38  size_t *sizetab = (size_t *)f->filter_arg_ptr;
39  size_t childoff = 0;
40 
41  assert( father->id == PASTIX_STARPU_INTERFACE_ID );
42  assert( father->offset == -1 );
43  assert( father->cblk->cblktype & CBLK_LAYOUT_2D );
44 
45  child->id = father->id;
46  child->flttype = father->flttype;
47  child->offset = sizetab[id];
48  child->nbblok = sizetab[id+1] - sizetab[id];
49  child->allocsize = 0;
50  child->cblk = father->cblk;
51  child->dataptr = NULL;
52 
53  assert( child->offset >= 0 );
54 
55  if ( father->cblk->cblktype & CBLK_COMPRESSED ) {
56  childoff = sizetab[id] * sizeof( pastix_lrblock_t );
57  child->allocsize = child->nbblok * sizeof( pastix_lrblock_t );
58  }
59  else {
60  SolverBlok *blok = father->cblk->fblokptr + sizetab[id];
61  SolverBlok *lblk = father->cblk->fblokptr + sizetab[id+1];
62  childoff = pastix_size_of( father->flttype ) * blok->coefind;
63 
64  if ( lblk->coefind == 0 ) {
65  int i;
66  int nbrow = 0;
67  for ( i=0; i<child->nbblok; i++, blok++) {
68  nbrow += blok_rownbr( blok );
69  }
70  child->allocsize = pastix_size_of( father->flttype ) * nbrow * cblk_colnbr( father->cblk );
71  }
72  else {
73  child->allocsize = pastix_size_of( father->flttype ) * (lblk->coefind - blok->coefind);
74  }
75  }
76 
77 #if defined(PASTIX_STARPU_INTERFACE_DEBUG)
78  fprintf( stderr,
79  "blok (%9s, size=%8zu, nbblok=%2ld )\n",
80  child->cblk->cblktype & CBLK_COMPRESSED ? "Low-rank" : "Full-rank",
81  child->allocsize, (long)(child->nbblok) );
82 #endif
83 
84  assert( child->allocsize > 0 );
85 
86  if ( father->dataptr != NULL ) {
87  child->dataptr = ((char*)father->dataptr) + childoff;
88  }
89 
90  (void)nchunks;
91 }
92 
93 static inline void
94 pastix_starpu_register_interface( const starpu_sparse_matrix_desc_t *spmtx,
95  SolverCblk *cblk,
96  int side,
97  pastix_coeftype_t flttype )
98 {
99  starpu_data_handle_t *handler = ( (starpu_data_handle_t *)( cblk->handler ) ) + side;
100 
101  pastix_starpu_register( handler, cblk, side, flttype );
102 #if defined( PASTIX_WITH_MPI )
103  int64_t tag_cblk;
104  if ( cblk->cblktype & CBLK_FANIN ) {
105  tag_cblk = 2 * cblk->gfaninnum + side;
106  starpu_mpi_data_register( *handler, spmtx->mpitag + tag_cblk, spmtx->solvmtx->clustnum );
107  starpu_data_set_reduction_methods( *handler, NULL, &cl_fanin_init_cpu );
108  }
109  else {
110  if ( cblk->cblktype & CBLK_RECV ) {
111  tag_cblk = 2 * cblk->gfaninnum + side;
112  }
113  else {
114  tag_cblk = 2 * cblk->gcblknum + side;
115  }
116  starpu_mpi_data_register( *handler, spmtx->mpitag + tag_cblk, cblk->ownerid );
117  }
118 #if defined(PASTIX_DEBUG_STARPU)
119  fprintf( stderr, "[%2d][pastix][%s] Matrix cblk=%d, owner=%d, tag=%ld, size=%ld\n",
120  spmtx->solvmtx->clustnum, __func__, cblk->gcblknum, cblk->ownerid,
121  spmtx->mpitag + tag_cblk,
122  cblk->stride * cblk_colnbr( cblk ) * pastix_size_of( spmtx->solvmtx->flttype ) );
123 #endif
124 #endif
125  (void)spmtx;
126 }
127 
128 static inline void
129 pastix_starpu_register_cblk( const starpu_sparse_matrix_desc_t *spmtx,
130  SolverCblk *cblk,
131  pastix_coeftype_t flttype )
132 {
133  /* Temporary fix for the factorisation with 2D cblk. */
134  if ( ( cblk->cblktype & (CBLK_FANIN|CBLK_RECV) ) &&
135  ( cblk->cblktype & CBLK_TASKS_2D ) ) {
136  if ( spmtx->mtxtype == PastixGeneral ) {
137  cpucblk_zalloc( PastixLUCoef, cblk );
138  }
139  else {
140  cpucblk_zalloc( PastixLCoef, cblk );
141  }
142  }
143 
144  pastix_starpu_register_interface( spmtx, cblk, PastixLCoef, flttype );
145  if ( spmtx->mtxtype == PastixGeneral ) {
146  pastix_starpu_register_interface( spmtx, cblk, PastixUCoef, flttype );
147  }
148 }
149 
150 #if defined( PASTIX_WITH_MPI )
151 static inline void
152 pastix_starpu_mpi_register_blok( const starpu_sparse_matrix_desc_t *spmtx,
153  const SolverCblk *cblk,
154  SolverBlok *blok,
155  int64_t tag_desc )
156 {
157  int64_t tag;
158  int ownerid = cblk->ownerid;
159 
160  if ( cblk->cblktype & CBLK_FANIN ) {
161  tag = tag_desc + 2 * blok->gfaninnm;
162  ownerid = spmtx->solvmtx->clustnum;
163  }
164  else if ( cblk->cblktype & CBLK_RECV ) {
165  tag = tag_desc + 2 * blok->gfaninnm;
166  }
167  else {
168  tag = tag_desc + 2 * blok->gbloknm;
169  }
170 
171  starpu_mpi_data_register( blok->handler[0], tag, ownerid );
172  if ( spmtx->mtxtype == PastixGeneral ) {
173  tag = tag + 1;
174  starpu_mpi_data_register( blok->handler[1], tag, ownerid );
175  }
176 }
177 #else
178 static inline void
179 pastix_starpu_mpi_register_blok( const starpu_sparse_matrix_desc_t *spmtx,
180  const SolverCblk *cblk,
181  SolverBlok *blok,
182  int64_t tag_desc )
183 {
184  (void)spmtx;
185  (void)cblk;
186  (void)blok;
187  (void)tag_desc;
188 }
189 #endif
190 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
191 
192 /**
193  *******************************************************************************
194  *
195  * @brief Generate the StarPU descriptor of the sparse matrix.
196  *
197  * This function creates the StarPU descriptor that will provide tha data
198  * mapping and memory location to StarPU for the computation.
199  *
200  *******************************************************************************
201  *
202  * @param[inout] solvmtx
203  * The solver matrix structure that describes the sparse matrix for
204  * PaStiX.
205  *
206  * @param[in] mtxtype
207  * The type of sparse matrix to describe.
208  * @arg PastixGeneral: The sparse matrix is general.
209  * @arg PastixSymmetric: The sparse matrix is lower triangular symmetric.
210  * @arg PastixHermitian: The sparse matrix is lower triangular hermitian.
211  *
212  * @param[in] nodes
213  * The number of processes used to solve the problem.
214  *
215  * @param[in] myrank
216  * The rank of the calling process.
217  *
218  * @param[in] flttype
219  * The memory size of the arithmetic used to store the matrix
220  * coefficients.
221  *
222  ******************************************************************************/
223 void
225  pastix_mtxtype_t mtxtype,
226  int nodes,
227  int myrank,
228  pastix_coeftype_t flttype )
229 {
230  pastix_int_t cblknbr, cblkmin2d;
231  size_t key1, key2;
232  SolverCblk *cblk;
233  SolverBlok *blok, *lblok;
234  pastix_int_t n = 0, cblknum;
235  pastix_int_t nbrow;
236  size_t size;
237  int64_t tag_desc;
238 
239  starpu_sparse_matrix_desc_t *spmtx = solvmtx->starpu_desc;
240  if ( spmtx != NULL ) {
242  }
243  else {
244  spmtx = (starpu_sparse_matrix_desc_t *)malloc( sizeof( starpu_sparse_matrix_desc_t ) );
245  }
246 
247  tag_desc = ( (int64_t) (solvmtx->gcblknbr + solvmtx->gbloknbr + solvmtx->gfanincblknbr + solvmtx->gfaninbloknbr) ) * 2;
248  spmtx->mpitag = pastix_starpu_tag_book( tag_desc );
249  tag_desc = spmtx->mpitag + 2 * ( solvmtx->gcblknbr + solvmtx->gfanincblknbr );
250  spmtx->typesze = pastix_size_of( flttype );
251  spmtx->mtxtype = mtxtype;
252  spmtx->solvmtx = solvmtx;
253  spmtx->cblktab_handle = NULL;
254  spmtx->gpu_blocktab = NULL;
255 
256  cblknbr = solvmtx->cblknbr;
257  cblkmin2d = solvmtx->cblkmin2d;
258  key1 = 2 * cblknbr;
259 
260  /* Initialize 1D cblk handlers */
261  cblk = spmtx->solvmtx->cblktab;
262  for ( cblknum = 0; cblknum < cblkmin2d; cblknum++, n++, cblk++ ) {
263  pastix_starpu_register_cblk( spmtx, cblk, flttype );
264  }
265 
266  /* Initialize 2D cblk handlers */
267  if ( cblkmin2d < cblknbr ) {
268  struct starpu_data_filter filter = { .filter_func = pastix_starpu_filter_interface };
269  starpu_cblk_t *cblkhandle;
270  size_t *sizetab = NULL;
271  pastix_int_t nchildren, sizenbr = 0;
272 
273  spmtx->cblktab_handle =
274  (starpu_cblk_t *)malloc( ( cblknbr - cblkmin2d ) * sizeof( starpu_cblk_t ) );
275 
276  cblk = spmtx->solvmtx->cblktab + cblkmin2d;
277  cblkhandle = spmtx->cblktab_handle;
278 
279  sizenbr = ( cblk[1].fblokptr - cblk[0].fblokptr ) + 1;
280  sizetab = malloc( sizenbr * sizeof( size_t ) );
281  assert( sizenbr >= 1 );
282 
283  for ( cblknum = cblkmin2d, n = 0; cblknum < cblknbr;
284  cblknum++, n++, cblk++, cblkhandle++ ) {
285  pastix_starpu_register_cblk( spmtx, cblk, flttype );
286 
287  if ( !(cblk->cblktype & CBLK_TASKS_2D) ) {
288  continue;
289  }
290 
291  /* Let's build the sizetab array */
292  blok = cblk[0].fblokptr;
293  lblok = cblk[1].fblokptr;
294 
295  if ( ( lblok - blok ) >= sizenbr ) {
296  sizenbr = ( lblok - blok ) + 1;
297  free( sizetab );
298  sizetab = malloc( sizenbr * sizeof( size_t ) );
299  }
300  nchildren = 0;
301  sizetab[0] = 0;
302 
303  /*
304  * Diagonal block
305  */
306  sizetab[nchildren + 1] = 1;
307  nchildren++;
308 
309  /*
310  * Off-diagonal blocks
311  */
312  blok++;
313  for ( ; blok < lblok; blok++ ) {
314  nbrow = 1;
315 
316  while ( ( blok + 1 < lblok ) &&
317  ( blok[0].fcblknm == blok[1].fcblknm ) &&
318  ( blok[0].lcblknm == blok[1].lcblknm ) )
319  {
320  blok++;
321  nbrow++;
322  }
323  size = nbrow;
324 
325  sizetab[nchildren + 1] = sizetab[nchildren] + size;
326  nchildren++;
327  }
328  filter.nchildren = nchildren;
329  filter.filter_arg_ptr = sizetab;
330 
331  cblkhandle->handlenbr = nchildren;
332  if ( mtxtype == PastixGeneral ) {
333  cblkhandle->handletab = (starpu_data_handle_t *)malloc(
334  2 * nchildren * sizeof( starpu_data_handle_t ) );
335 
336  starpu_data_partition_plan( cblk->handler[0], &filter, cblkhandle->handletab );
337 
338  starpu_data_partition_plan(
339  cblk->handler[1], &filter, cblkhandle->handletab + nchildren );
340  }
341  else {
342  cblkhandle->handletab =
343  (starpu_data_handle_t *)malloc( nchildren * sizeof( starpu_data_handle_t ) );
344 
345  starpu_data_partition_plan( cblk->handler[0], &filter, cblkhandle->handletab );
346  }
347 
348  nchildren = 0;
349  blok = cblk[0].fblokptr;
350  lblok = cblk[1].fblokptr;
351 
352  /*
353  * Diagonal block
354  */
355  blok->handler[0] = cblkhandle->handletab[nchildren];
356  if ( mtxtype == PastixGeneral ) {
357  blok->handler[1] = cblkhandle->handletab[cblkhandle->handlenbr + nchildren];
358  }
359  else {
360  blok->handler[1] = NULL;
361  }
362  pastix_starpu_mpi_register_blok( spmtx, cblk, blok, tag_desc );
363  nchildren++;
364 
365  /*
366  * Off-diagonal blocks
367  */
368  blok++;
369  for ( ; blok < lblok; blok++ ) {
370  blok->handler[0] = cblkhandle->handletab[nchildren];
371  if ( mtxtype == PastixGeneral ) {
372  blok->handler[1] = cblkhandle->handletab[cblkhandle->handlenbr + nchildren];
373  }
374  else {
375  blok->handler[1] = NULL;
376  }
377  pastix_starpu_mpi_register_blok( spmtx, cblk, blok, tag_desc );
378  nchildren++;
379 
380  while ( ( blok < lblok ) && ( blok[0].fcblknm == blok[1].fcblknm ) &&
381  ( blok[0].lcblknm == blok[1].lcblknm ) ) {
382  blok++;
383  blok->handler[0] = NULL;
384  blok->handler[1] = NULL;
385  }
386  }
387  }
388 
389  if ( sizetab != NULL ) {
390  free( sizetab );
391  }
392  }
393  solvmtx->starpu_desc = spmtx;
394 
395  (void)key1;
396  (void)key2;
397  (void)nodes;
398  (void)myrank;
399  (void)tag_desc;
400 }
401 
402 /**
403  *******************************************************************************
404  *
405  * @brief Submit asynchronous calls to retrieve the data on main memory.
406  *
407  *******************************************************************************
408  *
409  * @param[in] side
410  * The cblk for which the data needs to be flushed.
411  *
412  * @param[inout] cblk
413  * The cblk for which the data needs to be flushed.
414  *
415  ******************************************************************************/
416 void
418 {
419  if ( side != PastixUCoef ) {
420  assert( cblk->handler[0] );
421 #if defined(PASTIX_WITH_MPI)
422  starpu_mpi_cache_flush( MPI_COMM_WORLD, cblk->handler[0] );
423 #endif
424  starpu_data_wont_use( cblk->handler[0] );
425  }
426  if ( side != PastixLCoef ) {
427  assert( cblk->handler[1] );
428 #if defined(PASTIX_WITH_MPI)
429  starpu_mpi_cache_flush( MPI_COMM_WORLD, cblk->handler[1] );
430 #endif
431  starpu_data_wont_use( cblk->handler[1] );
432  }
433 }
434 
435 /**
436  *******************************************************************************
437  *
438  * @brief Submit asynchronous calls to retrieve the data on main memory.
439  *
440  *******************************************************************************
441  *
442  * @param[inout] spmtx
443  * The sparse matrix descriptor to retrieve on main memory.
444  *
445  ******************************************************************************/
446 void
448 {
449  SolverCblk *cblk;
450  pastix_int_t i;
451 
452  cblk = spmtx->solvmtx->cblktab;
453  for ( i = 0; i < spmtx->solvmtx->cblknbr; i++, cblk++ ) {
454  assert( cblk->handler[0] );
455 
456  if ( cblk->ownerid != spmtx->solvmtx->clustnum ) {
457  continue;
458  }
459 
460  starpu_data_acquire_cb( cblk->handler[0],
461  STARPU_R,
462  ( void( * )( void * ) ) & starpu_data_release,
463  cblk->handler[0] );
464 
465  if ( cblk->ucoeftab ) {
466  starpu_data_acquire_cb( cblk->handler[1],
467  STARPU_R,
468  ( void( * )( void * ) ) & starpu_data_release,
469  cblk->handler[1] );
470  }
471  }
472 }
473 
474 /**
475  *******************************************************************************
476  *
477  * @brief Destroy a single cblk StarPU data structure of the sparse matrix
478  *
479  *******************************************************************************
480  *
481  * @param[in] is_owner
482  * Boolean to specify if the calling process owns the cblk or not
483  *
484  * @param[in,out] cblk
485  * Pointer to the cblk data structure
486  *
487  * @param[in,out] starpu_cblk
488  * Pointer to the starpu data structure associated to the cblk
489  *
490  ******************************************************************************/
491 void
493  SolverCblk *cblk,
494  starpu_cblk_t *cblkhandle )
495 {
496  if ( cblk->cblktype & CBLK_TASKS_2D ) {
497  int gather_node = -1;
498 
499  if ( is_owner && !(cblk->cblktype & CBLK_FANIN) ) {
500  gather_node = STARPU_MAIN_RAM;
501  }
502 
503  /*
504  * First, let's unpartition ourself as long as StarPU does not use the
505  * correct gather_node in starpu_data_partition_clean()
506  */
507  if ( cblk->handler[0] ) {
508  starpu_data_partition_clean_node( cblk->handler[0],
509  cblkhandle->handlenbr,
510  cblkhandle->handletab,
511  gather_node );
512  }
513 
514  if ( cblk->handler[1] ) {
515  starpu_data_partition_clean_node( cblk->handler[1],
516  cblkhandle->handlenbr,
517  cblkhandle->handletab + cblkhandle->handlenbr,
518  gather_node );
519  }
520 
521  free( cblkhandle->handletab );
522  cblkhandle->handletab = NULL;
523  }
524 
525  if ( cblk->handler[0] ) {
526  starpu_data_unregister( cblk->handler[0] );
527  }
528  if ( cblk->handler[1] ) {
529  starpu_data_unregister( cblk->handler[1] );
530  }
531  cblk->handler[0] = NULL;
532  cblk->handler[1] = NULL;
533 }
534 
535 /**
536  *******************************************************************************
537  *
538  * @brief Free the StarPU descriptor of the sparse matrix.
539  *
540  * This function destroys the StarPU descriptor, but do not free the matrix data
541  * that are managed by PaStiX.
542  *
543  *******************************************************************************
544  *
545  * @param[inout] spmtx
546  * The descriptor to free.
547  *
548  ******************************************************************************/
549 void
551 {
552  starpu_cblk_t *cblkhandle;
553  SolverCblk *cblk;
554  pastix_int_t i, cblkmin2d;
555  pastix_int_t rank = spmtx->solvmtx->clustnum;
556 
557  cblkmin2d = spmtx->solvmtx->cblkmin2d;
558  cblk = spmtx->solvmtx->cblktab;
559  for ( i = 0; i < cblkmin2d; i++, cblk++ ) {
560  pastix_starpu_cblk_destroy( cblk->ownerid == rank, cblk, NULL );
561  }
562 
563  cblkhandle = spmtx->cblktab_handle;
564  for ( i = cblkmin2d; i < spmtx->solvmtx->cblknbr; i++, cblk++, cblkhandle++ ) {
565  pastix_starpu_cblk_destroy( cblk->ownerid == rank, cblk, cblkhandle );
566  }
567 
568  if ( spmtx->cblktab_handle != NULL ) {
569  free( spmtx->cblktab_handle );
570  }
571 
573 }
574 
575 /**
576  * @}
577  */
BEGIN_C_DECLS typedef int pastix_int_t
Definition: datatypes.h:51
void cpucblk_zalloc(pastix_coefside_t side, SolverCblk *cblk)
Allocate the cblk structure to store the coefficient.
struct pastix_lrblock_s pastix_lrblock_t
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
#define PastixGeneral
Definition: api.h:458
spm_mtxtype_t pastix_mtxtype_t
Matrix symmetry type property.
Definition: api.h:457
enum pastix_coefside_e pastix_coefside_t
Data blocks used in the kernel.
@ PastixLCoef
Definition: api.h:478
@ PastixLUCoef
Definition: api.h:480
@ PastixUCoef
Definition: api.h:479
starpu_data_handle_t * handletab
starpu_cblk_t * cblktab_handle
enum starpu_data_interface_id id
const SolverCblk * cblk
pastix_coeftype_t flttype
pastix_int_t handlenbr
void pastix_starpu_cblk_destroy(int is_owner, SolverCblk *cblk, starpu_cblk_t *cblkhandle)
Destroy a single cblk StarPU data structure of the sparse matrix.
void starpu_sparse_matrix_getoncpu(starpu_sparse_matrix_desc_t *spmtx)
Submit asynchronous calls to retrieve the data on main memory.
struct starpu_cblk_s starpu_cblk_t
Additional StarPU handlers for a column-block when using 2D kernels.
int64_t pastix_starpu_tag_book(int64_t nbtags)
Book a range of StarPU unique tags of size nbtags.
Definition: starpu_tags.c:246
struct starpu_codelet cl_fanin_init_cpu
Main structure for all tasks of fanin_init type.
void starpu_sparse_matrix_init(SolverMatrix *solvmtx, pastix_mtxtype_t mtxtype, int nodes, int myrank, pastix_coeftype_t flttype)
Generate the StarPU descriptor of the sparse matrix.
void pastix_starpu_tag_release(int64_t min)
Release the set of tags starting by min.
Definition: starpu_tags.c:264
#define PASTIX_STARPU_INTERFACE_ID
Alias to get the Interface id.
void starpu_sparse_cblk_wont_use(pastix_coefside_t side, SolverCblk *cblk)
Submit asynchronous calls to retrieve the data on main memory.
void pastix_starpu_register(starpu_data_handle_t *handleptr, const SolverCblk *cblk, pastix_coefside_t side, pastix_coeftype_t flttype)
Register a cblk at the StarPU level.
void starpu_sparse_matrix_destroy(starpu_sparse_matrix_desc_t *spmtx)
Free the StarPU descriptor of the sparse matrix.
Interface data structure to register the pieces of data in StarPU.
Additional StarPU handlers for a column-block when using 2D kernels.
Definition: pastix_starpu.h:99
StarPU descriptor stucture for the sparse matrix.
pastix_int_t cblkmin2d
Definition: solver.h:219
pastix_int_t gfanincblknbr
Definition: solver.h:212
pastix_coeftype_t flttype
Definition: solver.h:231
static pastix_int_t blok_rownbr(const SolverBlok *blok)
Compute the number of rows of a block.
Definition: solver.h:395
void * ucoeftab
Definition: solver.h:178
pastix_int_t gbloknbr
Definition: solver.h:225
pastix_int_t gcblknbr
Definition: solver.h:210
void * handler[2]
Definition: solver.h:142
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 gfaninnm
Definition: solver.h:146
pastix_int_t gfaninnum
Definition: solver.h:176
pastix_int_t cblknbr
Definition: solver.h:211
pastix_int_t gbloknm
Definition: solver.h:145
pastix_int_t gfaninbloknbr
Definition: solver.h:226
pastix_int_t coefind
Definition: solver.h:149
SolverBlok * fblokptr
Definition: solver.h:168
pastix_int_t gcblknum
Definition: solver.h:174
SolverCblk *restrict cblktab
Definition: solver.h:228
pastix_int_t stride
Definition: solver.h:169
void * handler[2]
Definition: solver.h:179
int8_t cblktype
Definition: solver.h:164
int ownerid
Definition: solver.h:181
Solver block structure.
Definition: solver.h:141
Solver column block structure.
Definition: solver.h:161
Solver column block structure.
Definition: solver.h:203