PaStiX Handbook 6.4.0
Loading...
Searching...
No Matches
bcsc_sspmv.c
Go to the documentation of this file.
1/**
2 *
3 * @file bcsc_sspmv.c
4 *
5 * Functions computing matrix-vector products for the BCSC
6 *
7 * @copyright 2004-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 Vincent Bridonneau
13 * @author Theophile Terraz
14 * @author Tony Delarue
15 * @date 2024-07-05
16 *
17 * @generated from /builds/2mk6rsew/0/solverstack/pastix/bcsc/bcsc_zspmv.c, normal z -> s, Tue Feb 25 14:36:02 2025
18 *
19 **/
20#include "common.h"
21#include <math.h>
22#include "bcsc/bcsc.h"
23#include "bcsc_s.h"
24#include "blend/solver.h"
25#include "pastix/datatypes.h"
26
27#ifndef DOXYGEN_SHOULD_SKIP_THIS
28
29typedef void ( *bcsc_sspmv_Ax_fct_t )( const pastix_bcsc_t *,
30 const bcsc_cblk_t *,
31 float,
32 const float *,
33 const float *,
34 float,
35 float * );
36
37static inline void
38__bcsc_sspmv_by( pastix_int_t n,
39 float beta,
40 float *y )
41{
42 if( beta != (float)0.0 )
43 {
45 for( j=0; j<n; j++, y++ )
46 {
47 (*y) *= beta;
48 }
49 }
50 else
51 {
52 memset( y, 0, n * sizeof(float) );
53 }
54}
55
56static inline void
57__bcsc_sspmv_Ax( const pastix_bcsc_t *bcsc,
58 const bcsc_cblk_t *cblk,
59 float alpha,
60 const float *A,
61 const float *x,
62 float beta,
63 float *y )
64{
65 pastix_int_t i, j;
66
67 __bcsc_sspmv_by( cblk->colnbr, beta, y );
68
69 for( j=0; j<cblk->colnbr; j++, y++ )
70 {
71 for( i=cblk->coltab[j]; i< cblk->coltab[j+1]; i++ )
72 {
73 *y += alpha * A[i] * x[ bcsc->rowtab[i] ];
74 }
75 }
76}
77
78static inline void
79__bcsc_sspmv_Ax_ind( const pastix_bcsc_t *bcsc,
80 float alpha,
81 const float *A,
82 const float *x,
83 float beta,
84 float *y )
85{
86 const float *xptr = x;
87 pastix_int_t bloc, i, j;
88
89 __bcsc_sspmv_by( bcsc->gN, beta, y );
90
91 for( bloc=0; bloc<bcsc->cscfnbr; bloc++ )
92 {
93 for( j=0; j < bcsc->cscftab[bloc].colnbr; j++, xptr++ )
94 {
95 for( i = bcsc->cscftab[bloc].coltab[j]; i < bcsc->cscftab[bloc].coltab[j+1]; i++ )
96 {
97 y[ bcsc->rowtab[i] ] += alpha * A[i] * (*xptr);
98 }
99 }
100 }
101}
102
103#if defined(PRECISION_z) || defined(PRECISION_c)
104static inline void
105__bcsc_sspmv_conjAx( const pastix_bcsc_t *bcsc,
106 const bcsc_cblk_t *cblk,
107 float alpha,
108 const float *A,
109 const float *x,
110 float beta,
111 float *y )
112{
113 pastix_int_t i, j;
114
115 __bcsc_sspmv_by( cblk->colnbr, beta, y );
116
117 for( j=0; j<cblk->colnbr; j++, y++ )
118 {
119 for( i=cblk->coltab[j]; i< cblk->coltab[j+1]; i++ )
120 {
121 *y += alpha * ( A[i] ) * x[ bcsc->rowtab[i] ];
122 }
123 }
124}
125#endif
126
127static inline void
128__bcsc_sspmv_loop( const SolverMatrix *solvmtx,
129 pastix_trans_t trans,
130 float alpha,
131 const pastix_bcsc_t *bcsc,
132 const float *x,
133 float beta,
134 float *y,
135 pastix_int_t rank,
136 pastix_int_t begin,
137 pastix_int_t end )
138{
139 bcsc_sspmv_Ax_fct_t zspmv_Ax = __bcsc_sspmv_Ax;
140 float *valptr = NULL;
141 pastix_int_t bloc;
142 bcsc_cblk_t *cblk;
143
144 /*
145 * There are three cases:
146 * We can use the Lvalues pointer directly:
147 * - The matrix is general and we use A^t
148 * - the matrix is symmetric or symmetric
149 * We can use the Uvalues pointer directly
150 * - The matrix is general and we use A
151 * We have to use Lvalues per row (instead of column)
152 * - The matrix A is general and Uvalues is unavailable
153 *
154 * To this, we have to add the call if Trans or Symmetric
155 *
156 * Mtxtype | trans asked | algo applied
157 * ++++++++++++++++++++++++++++++++++++
158 + General | NoTrans | U if possible, otherwise indirect L
159 + General | Trans | L
160 + General | Trans | (L)
161 + Symmetric | NoTrans | L
162 + Symmetric | Trans | L
163 + Symmetric | Trans | (L)
164 + Symmetric | NoTrans | (L)
165 + Symmetric | Trans | L
166 + Symmetric | Trans | (L)
167 */
168 cblk = bcsc->cscftab + begin;
169 valptr = (float*)bcsc->Lvalues;
170
171 if ( (bcsc->mtxtype == PastixGeneral) && (trans == PastixNoTrans) )
172 {
173 /* U */
174 if ( bcsc->Uvalues != NULL ) {
175 valptr = (float*)bcsc->Uvalues;
176 }
177 /* Indirect L */
178 else {
179 /* Execute in sequential */
180 if ( rank != 0 ) {
181 return;
182 }
183 __bcsc_sspmv_Ax_ind( bcsc, alpha, valptr, x, beta, y );
184 }
185 }
186#if defined(PRECISION_z) || defined(PRECISION_c)
187 /* Conj(L) */
188 else if ( ( (bcsc->mtxtype == PastixGeneral ) && (trans == PastixTrans) ) ||
189 ( (bcsc->mtxtype == PastixSymmetric) && (trans == PastixTrans) ) ||
190 ( (bcsc->mtxtype == PastixSymmetric) && (trans != PastixTrans ) ) )
191 {
192 zspmv_Ax = __bcsc_sspmv_conjAx;
193 }
194#endif /* defined(PRECISION_z) || defined(PRECISION_c) */
195
196 for( bloc=begin; bloc<end; bloc++, cblk++ )
197 {
198 const SolverCblk *solv_cblk = solvmtx->cblktab + cblk->cblknum;
199 float *yptr = y + solv_cblk->lcolidx;
200
201 assert( !(solv_cblk->cblktype & (CBLK_FANIN|CBLK_RECV)) );
202
203 zspmv_Ax( bcsc, cblk, alpha, valptr, x, beta, yptr );
204 }
205}
206
207#endif /* DOXYGEN_SHOULD_SKIP_THIS */
208
209/**
210 *******************************************************************************
211 *
212 * @ingroup bcsc
213 *
214 * @brief Compute the matrix-vector product y = alpha * A * x + beta * y
215 * (Sequential version)
216 *
217 * Where A is given in the bcsc format, x and y are two vectors of size n, and
218 * alpha and beta are two scalars.
219 *
220 *******************************************************************************
221 *
222 * @param[in] pastix_data
223 * Provide information about bcsc
224 *
225 * @param[in] trans
226 * Specifies whether the matrix A from the bcsc is transposed, not
227 * transposed or conjugate transposed:
228 * = PastixNoTrans: A is not transposed;
229 * = PastixTrans: A is transposed;
230 * = PastixTrans: A is conjugate transposed.
231 *
232 * @param[in] alpha
233 * alpha specifies the scalar alpha
234 *
235 * @param[in] x
236 * The vector x.
237 *
238 * @param[in] beta
239 * beta specifies the scalar beta
240 *
241 * @param[inout] y
242 * The vector y.
243 *
244 *******************************************************************************/
245void
246bcsc_sspmv_seq( const pastix_data_t *pastix_data,
247 pastix_trans_t trans,
248 float alpha,
249 const float *x,
250 float beta,
251 float *y )
252{
253 pastix_bcsc_t *bcsc = pastix_data->bcsc;
254 SolverMatrix *solvmtx = pastix_data->solvmatr;
255
256 if( (bcsc == NULL) || (y == NULL) || (x == NULL) ) {
257 return;
258 }
259
260 __bcsc_sspmv_loop( solvmtx,
261 trans, alpha, bcsc, x, beta, y,
262 0, 0, bcsc->cscfnbr );
263}
264
265/**
266 * @brief Data structure for parallel arguments of spmv functions
267 */
268struct s_argument_spmv_s {
269 pastix_trans_t trans;
270 float alpha;
271 const pastix_bcsc_t *bcsc;
272 const float *x;
273 float beta;
274 float *y;
275 SolverMatrix *mtx;
276 pastix_int_t *start_indexes; /* starting position for each thread*/
277 pastix_int_t *start_bloc;
278};
279
280/**
281 *******************************************************************************
282 *
283 * @ingroup bcsc_internal
284 *
285 * @brief Compute the matrix-vector product y = alpha * op(A) * x + beta * y
286 *
287 * Where A is given in the bcsc format, x and y are two vectors of size n, and
288 * alpha and beta are two scalars.
289 * The op function is specified by the trans parameter and performs the
290 * operation as follows:
291 * trans = PastixNoTrans y := alpha*A *x + beta*y
292 * trans = PastixTrans y := alpha*A' *x + beta*y
293 * trans = PastixTrans y := alpha*(A')*x + beta*y
294 *
295 *******************************************************************************
296 *
297 * @param[in] ctx
298 * the context of the current thread
299 *
300 * @param[inout] args
301 * The parameter as specified in bcsc_sspmv.
302 *
303 *******************************************************************************/
304void
305pthread_bcsc_sspmv( isched_thread_t *ctx,
306 void *args )
307{
308 struct s_argument_spmv_s *arg = (struct s_argument_spmv_s*)args;
309 const pastix_bcsc_t *bcsc = arg->bcsc;
310 pastix_int_t begin, end, size, rank;
311 pastix_int_t *start_indexes = arg->start_indexes;
312 pastix_int_t *start_bloc = arg->start_bloc;
313
314 rank = (pastix_int_t)ctx->rank;
315 size = (pastix_int_t)ctx->global_ctx->world_size;
316
317 begin = start_bloc[rank];
318 if ( rank == (size - 1) )
319 {
320 end = bcsc->cscfnbr;
321 }
322 else {
323 end = start_bloc[rank + 1];
324 }
325
326 __bcsc_sspmv_loop( arg->mtx,
327 arg->trans, arg->alpha, bcsc, arg->x,
328 arg->beta, arg->y + start_indexes[rank],
329 rank, begin, end );
330}
331
332/**
333 *******************************************************************************
334 *
335 * @ingroup bcsc_internal
336 *
337 * @brief Compute the matrix-vector product y = alpha * op(A) * x + beta * y
338 *
339 * Where A is given in the bcsc format, x and y are two vectors of size n, and
340 * alpha and beta are two scalars.
341 * The op function is specified by the trans parameter and performs the
342 * operation as follows:
343 * trans = PastixNoTrans y := alpha*A *x + beta*y
344 * trans = PastixTrans y := alpha*A' *x + beta*y
345 * trans = PastixTrans y := alpha*(A')*x + beta*y
346 *
347 *******************************************************************************
348 *
349 * @param[in] ctx
350 * the context of the current thread
351 *
352 * @param[inout] args
353 * The parameter as specified in bcsc_sspmv.
354 *
355 *******************************************************************************/
356void
357pthread_bcsc_sspmv_tasktab( isched_thread_t *ctx,
358 void *args )
359{
360 bcsc_sspmv_Ax_fct_t zspmv_Ax = __bcsc_sspmv_Ax;
361 struct s_argument_spmv_s *arg = (struct s_argument_spmv_s*)args;
362 pastix_trans_t trans = arg->trans;
363 float alpha = arg->alpha;
364 const pastix_bcsc_t *bcsc = arg->bcsc;
365 const float *x = arg->x;
366 float beta = arg->beta;
367 float *y = arg->y;
368 float *valptr = NULL;
369 float *yptr;
370 pastix_int_t rank;
371 SolverMatrix *mtx = arg->mtx;
372 pastix_int_t tasknbr, *tasktab;
373 pastix_int_t ii, task_id;
374 SolverCblk *solv_cblk;
375 bcsc_cblk_t *bcsc_cblk;
376 Task *t;
377
378 rank = (pastix_int_t)ctx->rank;
379
380 tasknbr = mtx->ttsknbr[rank];
381 tasktab = mtx->ttsktab[rank];
382
383 /*
384 * There are three cases:
385 * We can use the Lvalues pointer directly:
386 * - The matrix is general and we use A^t
387 * - The matrix is symmetric or symmetric
388 * We can use the Uvalues pointer directly
389 * - The matrix is general and we use A
390 * We have to use Lvalues per row (instead of column)
391 * - The matrix A is general and Uvalues is unavailable
392 *
393 * To this, we have to add the call if Trans or Symmetric
394 *
395 * Mtxtype | trans asked | algo applied
396 * ++++++++++++++++++++++++++++++++++++
397 + General | NoTrans | U if possible, otherwise indirect L
398 + General | Trans | L
399 + General | Trans | (L)
400 + Symmetric | NoTrans | L
401 + Symmetric | Trans | L
402 + Symmetric | Trans | (L)
403 + Symmetric | NoTrans | (L)
404 + Symmetric | Trans | L
405 + Symmetric | Trans | (L)
406 */
407 valptr = (float*)bcsc->Lvalues;
408
409 if ( (bcsc->mtxtype == PastixGeneral) && (trans == PastixNoTrans) )
410 {
411 /* U */
412 if ( bcsc->Uvalues != NULL ) {
413 valptr = (float*)bcsc->Uvalues;
414 }
415 /* Indirect L */
416 else {
417 /* Execute in sequential */
418 if ( rank != 0 ) {
419 return;
420 }
421 __bcsc_sspmv_Ax_ind( bcsc, alpha, valptr, x, beta, y );
422 return;
423 }
424 }
425#if defined(PRECISION_z) || defined(PRECISION_c)
426 /* Conj(L) */
427 else if ( ( (bcsc->mtxtype == PastixGeneral ) && (trans == PastixTrans) ) ||
428 ( (bcsc->mtxtype == PastixSymmetric) && (trans == PastixTrans) ) ||
429 ( (bcsc->mtxtype == PastixSymmetric) && (trans != PastixTrans ) ) )
430 {
431 zspmv_Ax = __bcsc_sspmv_conjAx;
432 }
433#endif /* defined(PRECISION_z) || defined(PRECISION_c) */
434
435 for (ii=0; ii<tasknbr; ii++)
436 {
437 task_id = tasktab[ii];
438 t = mtx->tasktab + task_id;
439
440 solv_cblk = mtx->cblktab + t->cblknum;
441 bcsc_cblk = bcsc->cscftab + solv_cblk->bcscnum;
442 yptr = y + solv_cblk->lcolidx;
443
444 zspmv_Ax( bcsc, bcsc_cblk, alpha, valptr, x, beta, yptr );
445 }
446}
447
448/**
449 *******************************************************************************
450 *
451 * @ingroup bcsc_internal
452 *
453 * @brief Initialize indexes for vector pointer and bloc indexes
454 * for parallel version of spmv.
455 *
456 * This function Initial indexes for each thread
457 * in order to computes it once instead of once per thread. This is a more
458 * sophisticated version trying to balance the load for each thread in terms
459 * of bloc size.
460 *
461 *******************************************************************************
462 *
463 * @param[in] pastix_data
464 * The pastix_data structure providing number of threads and holding
465 * the A matrix.
466 *
467 * @param[out] args
468 * The argument containing arrays to initialise (blocs and indexes).
469 *
470 *******************************************************************************/
471void
473 struct s_argument_spmv_s *args )
474{
475 pastix_int_t rank, bloc, size;
476 pastix_int_t ratio, total, load;
477 pastix_bcsc_t *bcsc = pastix_data->bcsc;
478 bcsc_cblk_t *cblk = bcsc->cscftab;
479
480 if ( bcsc->mtxtype != PastixGeneral ) {
481 total = 2 * pastix_data->csc->nnzexp - bcsc->gN;
482 } else {
483 total = pastix_data->csc->nnzexp;
484 }
485 size = pastix_data->isched->world_size;
486 ratio = pastix_iceil( total, size );
487 load = 0;
488
489 args->start_bloc[0] = 0;
490 args->start_indexes[0] = 0;
491
492 for ( bloc = 0, rank = 1; bloc < bcsc->cscfnbr; ++bloc, ++cblk )
493 {
494 if ( load >= ratio ) {
495 assert( rank < size );
496
497 args->start_bloc[rank] = bloc;
498 args->start_indexes[rank] = pastix_data->solvmatr->cblktab[bloc].fcolnum;
499
500 rank ++;
501 total -= load;
502 load = 0;
503 }
504 load += cblk->coltab[cblk->colnbr] - cblk->coltab[0];
505 }
506
507 total -= load;
508 assert( total == 0 );
509
510 for ( ; rank < size; rank ++ ) {
511 args->start_bloc[rank] = bcsc->cscfnbr;
512 args->start_indexes[rank] = bcsc->gN;
513 }
514}
515
516/**
517 *******************************************************************************
518 *
519 * @ingroup bcsc
520 *
521 * @brief Perform y = alpha A x + beta y (Parallel version)
522 *
523 * This functions is parallelized through the internal static scheduler.
524 *
525 *******************************************************************************
526 *
527 * @param[in] pastix_data
528 * The pastix_data structure that holds the A matrix.
529 *
530 * @param[in] trans
531 * Specifies whether the matrix A from the bcsc is transposed, not
532 * transposed or conjugate transposed:
533 * = PastixNoTrans: A is not transposed;
534 * = PastixTrans: A is transposed;
535 * = PastixTrans: A is conjugate transposed.
536 *
537 * @param[in] alpha
538 * The scalar alpha.
539 *
540 * @param[in] x
541 * The vector x
542 *
543 * @param[in] beta
544 * The scalar beta.
545 *
546 * @param[inout] y
547 * On entry, the vector y
548 * On exit, alpha A x + y
549 *
550 *******************************************************************************/
551void
552bcsc_sspmv_smp( const pastix_data_t *pastix_data,
553 pastix_trans_t trans,
554 float alpha,
555 const float *x,
556 float beta,
557 float *y )
558{
559 pastix_bcsc_t *bcsc = pastix_data->bcsc;
560 struct s_argument_spmv_s arg = { trans, alpha, bcsc, x, beta, y,
561 pastix_data->solvmatr, NULL, NULL };
562
563 if( (bcsc == NULL) || (y == NULL) || (x == NULL) ) {
564 return;
565 }
566
567 isched_parallel_call( pastix_data->isched, pthread_bcsc_sspmv_tasktab, &arg );
568
569#if 0
570 /*
571 * Version that balances the number of nnz per thread, instead of exploiting
572 * the tasktab array.
573 */
574 {
575 MALLOC_INTERN( arg.start_indexes, 2 * pastix_data->isched->world_size, pastix_int_t );
576 arg.start_bloc = arg.start_indexes + pastix_data->isched->world_size;
577
578 bcsc_sspmv_get_balanced_indexes( pastix_data, &arg );
579
580 isched_parallel_call( pastix_data->isched, pthread_bcsc_sspmv, &arg );
581
582 memFree_null( arg.start_indexes );
583 }
584#endif
585}
586
587/**
588 *******************************************************************************
589 *
590 * @brief Compute the matrix-vector product y = alpha * op(A) * x + beta * y
591 *
592 * Where A is given in the bcsc format, x and y are two vectors of size n, and
593 * alpha and beta are two scalars.
594 * The op function is specified by the trans parameter and performs the
595 * operation as follows:
596 * trans = PastixNoTrans y := alpha*A *x + beta*y
597 * trans = PastixTrans y := alpha*A' *x + beta*y
598 * trans = PastixTrans y := alpha*(A')*x + beta*y
599 *
600 * This function is used only in testings.
601 *
602 *******************************************************************************
603 *
604 * @param[in] pastix_data
605 * Provide information about bcsc, and select the scheduling version
606 * based on iparm[IPARM_SCHEDULER].
607 *
608 * @param[in] trans
609 * Specifies whether the matrix A from the bcsc is transposed, not
610 * transposed or conjugate transposed:
611 * = PastixNoTrans: A is not transposed;
612 * = PastixTrans: A is transposed;
613 * = PastixTrans: A is conjugate transposed.
614 *
615 * @param[in] alpha
616 * alpha specifies the scalar alpha
617 *
618 * @param[in] x
619 * The vector x.
620 *
621 * @param[in] beta
622 * beta specifies the scalar beta
623 *
624 * @param[inout] y
625 * The vector y.
626 *
627 *******************************************************************************/
628void
629bcsc_sspmv( const pastix_data_t *pastix_data,
630 pastix_trans_t trans,
631 float alpha,
632 const float *x,
633 float beta,
634 float *y )
635{
636 const float *xglobal;
637 const pastix_int_t *iparm = pastix_data->iparm;
639
640 /*
641 * trans | transA | Final
642 * ----------------+-----------------+-----------------
643 * PastixNoTrans | PastixNoTrans | PastixNoTrans
644 * PastixNoTrans | PastixTrans | PastixTrans
645 * PastixNoTrans | PastixTrans | PastixTrans
646 * PastixTrans | PastixNoTrans | PastixTrans
647 * PastixTrans | PastixTrans | PastixNoTrans
648 * PastixTrans | PastixTrans | INCORRECT
649 * PastixTrans | PastixNoTrans | PastixTrans
650 * PastixTrans | PastixTrans | INCORRECT
651 * PastixTrans | PastixTrans | PastixNoTrans
652 */
653 if ( trans == PastixNoTrans ) {
654 trans = transA;
655 }
656 else if ( trans == transA ) {
657 trans = PastixNoTrans;
658 }
659 else if ( transA != PastixNoTrans ) {
660 pastix_print_error( "bcsc_sspmv: incompatible trans and transA" );
661 return;
662 }
663
664 /* y is duplicated on all nodes. Set to 0 non local data */
665 xglobal = bvec_sgather_remote( pastix_data, x );
666
667 if ( (iparm[IPARM_SCHEDULER] == PastixSchedStatic) ||
668 (iparm[IPARM_SCHEDULER] == PastixSchedDynamic) ) {
669 bcsc_sspmv_smp( pastix_data, trans, alpha, xglobal, beta, y );
670 }
671 else {
672 bcsc_sspmv_seq( pastix_data, trans, alpha, xglobal, beta, y );
673 }
674
675 if ( x != xglobal ) {
676 free( (void*)xglobal );
677 }
678}
BEGIN_C_DECLS typedef int pastix_int_t
Definition datatypes.h:51
void bcsc_sspmv_get_balanced_indexes(const pastix_data_t *pastix_data, struct s_argument_spmv_s *args)
Initialize indexes for vector pointer and bloc indexes for parallel version of spmv.
Definition bcsc_sspmv.c:472
void pthread_bcsc_sspmv_tasktab(isched_thread_t *ctx, void *args)
Compute the matrix-vector product y = alpha * op(A) * x + beta * y.
Definition bcsc_sspmv.c:357
void pthread_bcsc_sspmv(isched_thread_t *ctx, void *args)
Compute the matrix-vector product y = alpha * op(A) * x + beta * y.
Definition bcsc_sspmv.c:305
pastix_int_t colnbr
Definition bcsc.h:111
pastix_int_t * coltab
Definition bcsc.h:113
pastix_int_t cblknum
Definition bcsc.h:112
const float * bvec_sgather_remote(const pastix_data_t *pastix_data, const float *y)
Gather a distributed right hand side (bvec storage) on all nodes.
void bcsc_sspmv_smp(const pastix_data_t *pastix_data, pastix_trans_t trans, float alpha, const float *x, float beta, float *y)
Perform y = alpha A x + beta y (Parallel version)
Definition bcsc_sspmv.c:552
void bcsc_sspmv(const pastix_data_t *pastix_data, pastix_trans_t trans, float alpha, const float *x, float beta, float *y)
Compute the matrix-vector product y = alpha * op(A) * x + beta * y.
Definition bcsc_sspmv.c:629
void bcsc_sspmv_seq(const pastix_data_t *pastix_data, pastix_trans_t trans, float alpha, const float *x, float beta, float *y)
Compute the matrix-vector product y = alpha * A * x + beta * y (Sequential version)
Definition bcsc_sspmv.c:246
Compressed colptr format for the bcsc.
Definition bcsc.h:110
#define PastixSymmetric
Definition api.h:459
#define PastixGeneral
Definition api.h:458
enum pastix_trans_e pastix_trans_t
Transpostion.
@ IPARM_SCHEDULER
Definition api.h:117
@ IPARM_TRANSPOSE_SOLVE
Definition api.h:106
@ PastixNoTrans
Definition api.h:445
@ PastixTrans
Definition api.h:446
@ PastixSchedStatic
Definition api.h:335
@ PastixSchedDynamic
Definition api.h:338
SolverMatrix * solvmatr
Definition pastixdata.h:103
pastix_int_t * iparm
Definition pastixdata.h:70
const spmatrix_t * csc
Definition pastixdata.h:90
isched_t * isched
Definition pastixdata.h:86
pastix_bcsc_t * bcsc
Definition pastixdata.h:102
Main PaStiX data structure.
Definition pastixdata.h:68
pastix_int_t cblknum
Definition solver.h:125
pastix_int_t lcolidx
Definition solver.h:170
pastix_int_t bcscnum
Definition solver.h:175
SolverCblk *restrict cblktab
Definition solver.h:228
int8_t cblktype
Definition solver.h:164
pastix_int_t fcolnum
Definition solver.h:166
Solver column block structure.
Definition solver.h:161
Solver column block structure.
Definition solver.h:203
The task structure for the numerical factorization.
Definition solver.h:122