PaStiX Handbook 6.4.0
Loading...
Searching...
No Matches
pastix_task_solve.c
Go to the documentation of this file.
1/**
2 *
3 * @file pastix_task_solve.c
4 *
5 * PaStiX solve routines
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 Pascal Henon
12 * @author Xavier Lacoste
13 * @author Pierre Ramet
14 * @author Mathieu Faverge
15 * @author Theophile Terraz
16 * @author Alycia Lisito
17 * @author Brieuc Nicolas
18 * @author Tony Delarue
19 * @date 2024-07-05
20 *
21 **/
22#ifndef DOXYGEN_SHOULD_SKIP_THIS
23#define _GNU_SOURCE 1
24#endif /* DOXYGEN_SHOULD_SKIP_THIS */
25#include "common.h"
26#include "bcsc/bcsc.h"
27#include "pastix/order.h"
28#include "order/order_internal.h"
29#include "blend/solver.h"
30#include "sopalin/sopalin_data.h"
31#include "pastix_papi.h"
32
33#include "bcsc/bcsc_z.h"
34#include "bcsc/bcsc_c.h"
35#include "bcsc/bcsc_d.h"
36#include "bcsc/bcsc_s.h"
37
38#include <lapacke.h>
39
40#ifndef DOXYGEN_SHOULD_SKIP_THIS
41
42#if defined(PASTIX_DEBUG_SOLVE)
43#include <spm/z_spm.h>
44#include <spm/c_spm.h>
45#include <spm/d_spm.h>
46#include <spm/s_spm.h>
47
48/**
49 *******************************************************************************
50 *
51 * @brief TODO
52 *
53 *******************************************************************************
54 *
55 * @param[in] pastix_data
56 * The pastix_data structure that describes the solver instance.
57 *
58 * @param[in] name
59 * TODO
60 *
61 * @param[in] B
62 * TODO
63 *
64 *******************************************************************************/
65static inline void
66pastix_rhs_dump( pastix_data_t *pastix_data,
67 const char *name,
68 pastix_rhs_t B )
69{
70 static volatile int32_t step = 0;
71
72 int32_t lstep = pastix_atomic_inc_32b( &step );
73 char *fullname;
74 FILE *f;
75 int rc;
76
77 rc = asprintf( &fullname, "%02d_%s.%02d.rhs",
78 lstep, name, pastix_data->solvmatr->clustnum );
79 assert( rc != -1 );
80
81 f = pastix_fopenw( (pastix_data->dir_local == NULL) ? "." : pastix_data->dir_local,
82 fullname, "w" );
83 free( fullname );
84
85 switch( B->flttype ) {
86 case SpmComplex64:
87 z_spmDensePrint( f, B->m, B->n, B->b, B->ld );
88 break;
89 case SpmComplex32:
90 c_spmDensePrint( f, B->m, B->n, B->b, B->ld );
91 break;
92 case SpmDouble:
93 d_spmDensePrint( f, B->m, B->n, B->b, B->ld );
94 break;
95 case SpmFloat:
96 s_spmDensePrint( f, B->m, B->n, B->b, B->ld );
97 break;
98 case SpmPattern:
99 ;
100 }
101
102 fclose(f);
103 (void)rc;
104}
105#else
106static inline void
107pastix_rhs_dump( pastix_data_t *pastix_data,
108 const char *name,
109 pastix_rhs_t B )
110{
111 (void)pastix_data;
112 (void)name;
113 (void)B;
114}
115#endif
116#endif /* DOXYGEN_SHOULD_SKIP_THIS */
117
118/**
119 *******************************************************************************
120 *
121 * @ingroup pastix_solve
122 *
123 * @brief Apply a permutation on the right-and-side vector before the solve step.
124 *
125 * This routine is affected by the following parameters:
126 * IPARM_VERBOSE, IPARM_FACTORIZATION, IPARM_APPLYPERM_WS.
127 *
128 *******************************************************************************
129 *
130 * @param[inout] pastix_data
131 * The pastix_data structure that describes the solver instance.
132 *
133 * @param[in] dir
134 * Forward or backword application of the permutation.
135 *
136 * @param[in] m
137 * Size of the right-and-side vectors.
138 *
139 * @param[in] n
140 * Number of right-and-side vectors.
141 *
142 * @param[inout] b
143 * The right-and-side vectors of size ldb-by-n.
144 * If dir == PastixDirForward, b is used as input to initialize the
145 * permuted b. If the matrix is replicated on all nodes or in shared
146 * memory, the permuted vector has the same size as b, and b is
147 * modified in-place.
148 * If dir == PastixDirBackward, b must be allocated on entry and is
149 * filled by the reverse permutation of Bp. Note that if the matrix is
150 * replicated on all nodes or in shared memory, b is modified in-place.
151 *
152 * @param[in] ldb
153 * The leading dimension of the right-and-side vectors.
154 *
155 * @param[inout] Bp
156 * The right-and-side vectors of size ldb-by-n. On entry, must be
157 * initialized through pastixRhsInit(). On exit, contains the
158 * information about the permuted vector when dir ==
159 * PastixDirForward. When dir == PastixDirBackward, the data is used as
160 * input to update b.
161 *
162 *******************************************************************************
163 *
164 * @retval PASTIX_SUCCESS on successful exit,
165 * @retval PASTIX_ERR_BADPARAMETER if one parameter is incorrect.
166 *
167 *******************************************************************************/
168int
170 pastix_dir_t dir,
171 pastix_int_t m,
172 pastix_int_t n,
173 void *b,
174 pastix_int_t ldb,
175 pastix_rhs_t Bp )
176{
177 pastix_coeftype_t flttype;
178
179 /*
180 * Checks parameters.
181 */
182 if (pastix_data == NULL) {
183 pastix_print_error( "pastix_subtask_applyorder: wrong pastix_data parameter" );
185 }
186 if (Bp == NULL) {
187 pastix_print_error( "pastix_subtask_applyorder: wrong Bp parameter" );
189 }
190
191 /* Makes sure ordering is 0 based. */
192 if ( pastix_data->ordemesh->baseval != 0 ) {
193 pastix_print_error( "pastix_subtask_applyorder: ordermesh must be 0-based" );
195 }
196
197 flttype = pastix_data->csc->flttype;
198#if defined(PASTIX_DEBUG_SOLVE)
199 if ( dir == PastixDirForward ) {
200 struct pastix_rhs_s rhsb = {
201 .allocated = 0,
202 .flttype = flttype,
203 .m = m,
204 .n = n,
205 .b = b,
206 .ld = ldb,
207 };
208 pastix_rhs_dump( pastix_data, "applyorder_Forward_input", &rhsb );
209 }
210 else {
211 pastix_rhs_dump( pastix_data, "applyorder_Backward_input", Bp );
212 }
213#endif
214
215 /* See also xlapmr and xlapmt */
216 switch( flttype ) {
217 case PastixComplex64:
218 bvec_zlapmr( pastix_data, dir, m, n, b, ldb, Bp );
219 break;
220
221 case PastixComplex32:
222 bvec_clapmr( pastix_data, dir, m, n, b, ldb, Bp );
223 break;
224
225 case PastixDouble:
226 bvec_dlapmr( pastix_data, dir, m, n, b, ldb, Bp );
227 break;
228
229 case PastixFloat:
230 bvec_slapmr( pastix_data, dir, m, n, b, ldb, Bp );
231 break;
232
233 default:
234 pastix_print_error( "The floating type of the rhs is not defined\n" );
236 }
237
238#if defined(PASTIX_DEBUG_SOLVE)
239 if ( dir == PastixDirForward ) {
240 pastix_rhs_dump( pastix_data, "applyorder_Forward_output", Bp );
241 }
242 else {
243 struct pastix_rhs_s rhsb = {
244 .allocated = 0,
245 .flttype = flttype,
246 .m = m,
247 .n = n,
248 .b = b,
249 .ld = ldb,
250 };
251 pastix_rhs_dump( pastix_data, "applyorder_Backward_output", &rhsb );
252 }
253#endif
254
255 return PASTIX_SUCCESS;
256}
257
258/**
259 *******************************************************************************
260 *
261 * @ingroup pastix_solve
262 *
263 * @brief Apply a triangular solve on the right-and-side vectors.
264 *
265 * This routine is affected by the following parameters:
266 * IPARM_VERBOSE, IPARM_FACTORIZATION.
267 *
268 *******************************************************************************
269 *
270 * @param[inout] pastix_data
271 * The pastix_data structure that describes the solver instance.
272 *
273 * @param[in] side
274 * Left or right application.
275 *
276 * @param[in] uplo
277 * Upper or Lower part.
278 *
279 * @param[in] trans
280 * With or without transposition (or conjugate transposition).
281 *
282 * @param[in] diag
283 * Diagonal terms are unit or not.
284 *
285 * @param[inout] Bp
286 * The right-and-side vector (can be multiple rhs).
287 * On exit, the solution is stored in place of the right-hand-side vector.
288 *
289 *******************************************************************************
290 *
291 * @retval PASTIX_SUCCESS on successful exit,
292 * @retval PASTIX_ERR_BADPARAMETER if one parameter is incorrect.
293 *
294 *******************************************************************************/
295int
297 pastix_side_t side,
298 pastix_uplo_t uplo,
299 pastix_trans_t trans,
300 pastix_diag_t diag,
301 pastix_rhs_t Bp )
302{
303 SolverMatrix *solvmtx;
304 sopalin_data_t sopalin_data;
306 int prevnt;
307
308 /*
309 * Check parameters
310 */
311 if (pastix_data == NULL) {
312 pastix_print_error( "pastix_subtask_trsm: wrong pastix_data parameter" );
314 }
315 if (Bp == NULL) {
316 pastix_print_error( "pastix_subtask_trsm: wrong Bp parameter" );
318 }
319 if ( !(pastix_data->steps & STEP_NUMFACT) ) {
320 pastix_print_error( "pastix_subtask_trsm: All steps from pastix_task_init() to pastix_task_numfact() have to be called before calling this function" );
322 }
323
324 flttype = Bp->flttype;
325 solvmtx = pastix_data->solvmatr;
326
327 if ( Bp->cblkb == NULL ) {
328 pastix_int_t nbbuffers;
329
330 nbbuffers = solvmtx->faninnbr + solvmtx->recvnbr;
331 if ( nbbuffers > 0 ) {
332 Bp->cblkb = calloc( nbbuffers, sizeof(void*) );
333 }
334 }
335
336 /*
337 * Ensure that the scheduler is correct and is in the same
338 * family that the one used for the factorization.
339 */
340 pastix_check_and_correct_scheduler( pastix_data );
341
342 sopalin_data.solvmtx = solvmtx;
343
345
346 switch (flttype) {
347 case PastixComplex64:
348 {
349 sopalin_ztrsm( pastix_data, side, uplo, trans, diag,
350 &sopalin_data, Bp );
351 }
352 break;
353 case PastixComplex32:
354 {
355 sopalin_ctrsm( pastix_data, side, uplo, trans, diag,
356 &sopalin_data, Bp );
357 }
358 break;
359 case PastixDouble:
360 {
361 trans = (trans == PastixConjTrans) ? PastixTrans : trans;
362 sopalin_dtrsm( pastix_data, side, uplo, trans, diag,
363 &sopalin_data, Bp );
364 }
365 break;
366 case PastixFloat:
367 {
368 trans = (trans == PastixConjTrans) ? PastixTrans : trans;
369 sopalin_strsm( pastix_data, side, uplo, trans, diag,
370 &sopalin_data, Bp );
371 }
372 break;
373 default:
374 fprintf(stderr, "Unknown floating point arithmetic\n" );
375 }
376
377 pastixBlasSetNumThreads( prevnt );
378
379#if defined(PASTIX_DEBUG_SOLVE)
380 {
381 char *fullname;
382 asprintf( &fullname, "solve_trsm_%c%c%c%c",
383 (side == PastixLeft) ? 'L' : 'R',
384 (uplo == PastixUpper) ? 'U' : 'L',
385 (trans == PastixConjTrans) ? 'C' : (trans == PastixTrans ? 'T' : 'N'),
386 (diag == PastixUnit) ? 'U' : 'N' );
387
388 pastix_rhs_dump( pastix_data, fullname, Bp );
389 free(fullname);
390 }
391#endif
392
393 return PASTIX_SUCCESS;
394}
395
396/**
397 *******************************************************************************
398 *
399 * @ingroup pastix_solve
400 *
401 * @brief Apply a diagonal operation on the right-and-side vectors.
402 *
403 * This routine is affected by the following parameters:
404 * IPARM_VERBOSE, IPARM_FACTORIZATION.
405 *
406 *******************************************************************************
407 *
408 * @param[inout] pastix_data
409 * The pastix_data structure that describes the solver instance.
410 *
411 * @param[inout] Bp
412 * The right-and-side vector (can be multiple rhs).
413 * On exit, the solution is stored in place of the right-hand-side vector.
414 *
415 *******************************************************************************
416 *
417 * @retval PASTIX_SUCCESS on successful exit,
418 * @retval PASTIX_ERR_BADPARAMETER if one parameter is incorrect.
419 *
420 *******************************************************************************/
421int
423 pastix_rhs_t Bp )
424{
425 sopalin_data_t sopalin_data;
426 int prevnt;
427
428 /*
429 * Check parameters
430 */
431 if (pastix_data == NULL) {
432 pastix_print_error( "pastix_subtask_diag: wrong pastix_data parameter" );
434 }
435 if (Bp == NULL) {
436 pastix_print_error( "pastix_subtask_diag: wrong Bp parameter" );
438 }
439 if ( !(pastix_data->steps & STEP_NUMFACT) ) {
440 pastix_print_error( "pastix_subtask_trsm: All steps from pastix_task_init() to pastix_task_numfact() have to be called before calling this function" );
442 }
443
444 /*
445 * Ensure that the scheduler is correct and is in the same
446 * family that the one used for the factorization.
447 */
448 pastix_check_and_correct_scheduler( pastix_data );
449
450 sopalin_data.solvmtx = pastix_data->solvmatr;
451
453
454 switch ( Bp->flttype ) {
455 case PastixComplex64:
456 sopalin_zdiag( pastix_data, &sopalin_data, Bp );
457 break;
458 case PastixComplex32:
459 sopalin_cdiag( pastix_data, &sopalin_data, Bp );
460 break;
461 case PastixDouble:
462 sopalin_ddiag( pastix_data, &sopalin_data, Bp );
463 break;
464 case PastixFloat:
465 sopalin_sdiag( pastix_data, &sopalin_data, Bp );
466 break;
467 default:
468 fprintf(stderr, "Unknown floating point arithmetic\n" );
469 }
470
471 pastixBlasSetNumThreads( prevnt );
472
473 pastix_rhs_dump( pastix_data, "solve_diag", Bp );
474
475 return PASTIX_SUCCESS;
476}
477
478/**
479 *******************************************************************************
480 *
481 * @ingroup pastix_solve
482 *
483 * @brief Solve the given problem without applying the permutation.
484 *
485 * @warning The input vector is considered already permuted. For a solve step
486 * with permutation, see pastix_task_solve()
487 *
488 * This routine is affected by the following parameters:
489 * IPARM_VERBOSE, IPARM_FACTORIZATION.
490 *
491 *******************************************************************************
492 *
493 * @param[inout] pastix_data
494 * The pastix_data structure that describes the solver instance.
495 *
496 * @param[in] transA
497 * PastixNoTrans: A is not transposed (CSC matrix)
498 * PastixTrans: A is transposed (CSR of symmetric/general matrix)
499 * PastixConjTrans: A is conjugate transposed (CSR of hermitian matrix)
500 *
501 * @param[inout] Bp
502 * The right-and-side vectors (can be multiple rhs).
503 * On exit, the solution is stored in place of the right-hand-side vector.
504 *
505 *******************************************************************************
506 *
507 * @retval PASTIX_SUCCESS on successful exit,
508 * @retval PASTIX_ERR_BADPARAMETER if one parameter is incorrect.
509 *
510 *******************************************************************************/
511int
513 pastix_trans_t transA,
514 pastix_rhs_t Bp )
515{
516 pastix_bcsc_t *bcsc;
517 pastix_factotype_t factotype;
518
519 pastix_uplo_t uplo;
520 pastix_diag_t diag;
521 pastix_trans_t trans;
522 pastix_trans_t transfact = PastixTrans;
523 pastix_rhs_t sBp;
524 pastix_rhs_t B;
525
526 /*
527 * Check parameters
528 */
529 if (pastix_data == NULL) {
530 pastix_print_error( "pastix_task_solve: wrong pastix_data parameter" );
532 }
533 if ( !(pastix_data->steps & STEP_NUMFACT) ) {
534 pastix_print_error( "pastix_task_solve: All steps from pastix_task_init() to pastix_task_numfact() have to be called before calling this function" );
536 }
537
538 bcsc = pastix_data->bcsc;
539 factotype = pastix_data->iparm[IPARM_FACTORIZATION];
540
541 /**
542 *
543 * Summary of the operations to perform the full solve if the original A was
544 * either transposed or conjugate transposed.
545 *
546 * op(A) | Factorization | Step 1 | Step 2
547 * ----------+---------------------+-----------+-----------
548 * NoTrans | L U | L y = b | U x = y
549 * NoTrans | L L^t | L y = b | L^t x = y
550 * NoTrans | L L^h | L y = b | L^h x = y
551 * Trans |(L U )^t = U^t L^t | U^t y = b | L^t x = y
552 * Trans |(L L^t)^t = L L^t | L y = b | L^t x = y
553 * Trans |(L L^h)^t = c(L) L^t | Not handled (c(L))
554 * ConjTrans |(L U )^h = U^h L^h | Not handled (U^h)
555 * ConjTrans |(L L^t)^h = c(L) L^h | Not handled (c(L))
556 * ConjTrans |(L L^h)^h = L L^h | L y = b | L^h x = y
557 *
558 */
559 /* Value of the transpose case */
560 if ( ((bcsc->flttype == PastixComplex32) || (bcsc->flttype == PastixComplex64)) &&
561 ((factotype == PastixFactLLH ) || ( factotype == PastixFactLDLH )) )
562 {
563 transfact = PastixConjTrans;
564 }
565
566 if ( (transA != PastixNoTrans) &&
567 (transA != transfact) )
568 {
569 pastix_print_error( "pastix_task_solve: transA incompatible with the factotype used (require extra conj(L) not handled)" );
571 }
572
573 {
574 double timer, energy;
575
576 papiEnergyStart();
577 if ( pastix_data->iparm[IPARM_TRACE] & PastixTraceSolve ) {
579 }
580 /* Start timer */
581 clockSyncStart( timer, pastix_data->inter_node_comm );
582
583 if ( pastix_data->iparm[IPARM_MIXED] &&
584 ( ( Bp->flttype == PastixComplex64 ) || ( Bp->flttype == PastixDouble ) ) )
585 {
586 pastixRhsInit( &sBp );
587 pastixRhsDoubletoSingle( Bp, sBp );
588
589 B = sBp;
590 }
591 else {
592 B = Bp;
593 }
594
595 /*
596 * Solve the first step
597 */
598 uplo = PastixLower;
599 trans = PastixNoTrans;
600 diag = PastixNonUnit;
601
602 if (( transA != PastixNoTrans ) && ( factotype == PastixFactLU )) {
603 uplo = PastixUpper;
604 trans = transA;
605 }
606 if (( transA == PastixNoTrans ) && ( factotype == PastixFactLU )) {
607 diag = PastixUnit;
608 }
609
610 if (( factotype == PastixFactLDLT ) ||
611 ( factotype == PastixFactLDLH ) )
612 {
613 diag = PastixUnit;
614 }
615
616 pastix_subtask_trsm( pastix_data, PastixLeft, uplo, trans, diag, B );
617
618 /*
619 * Solve the diagonal step
620 */
621 if( (factotype == PastixFactLDLT) ||
622 (factotype == PastixFactLDLH) )
623 {
624 /* Solve y = D z with z = ([L^t | L^h] P x) */
625 pastix_subtask_diag( pastix_data, B );
626 }
627
628 /*
629 * Solve the second step
630 */
631 uplo = PastixLower;
632 trans = transfact;
633 diag = PastixNonUnit;
634
635 if (( transA == PastixNoTrans ) && ( factotype == PastixFactLU )) {
636 uplo = PastixUpper;
637 trans = PastixNoTrans;
638 }
639 if (( transA != PastixNoTrans ) && ( factotype == PastixFactLU )) {
640 diag = PastixUnit;
641 }
642
643 if (( factotype == PastixFactLDLT ) ||
644 ( factotype == PastixFactLDLH ) )
645 {
646 diag = PastixUnit;
647 }
648
649 pastix_subtask_trsm( pastix_data, PastixLeft, uplo, trans, diag, B );
650
651 if ( pastix_data->iparm[IPARM_MIXED] &&
652 ( ( Bp->flttype == PastixComplex64 ) || ( Bp->flttype == PastixDouble ) ) )
653 {
654 pastixRhsSingleToDouble( sBp, Bp );
655 pastixRhsFinalize( sBp );
656 }
657
658 /* Stop Timer */
659 clockSyncStop( timer, pastix_data->inter_node_comm );
660 if ( pastix_data->iparm[IPARM_TRACE] & PastixTraceSolve ) {
662 }
663 energy = papiEnergyStop();
664
665 pastix_data->dparm[DPARM_SOLV_TIME] = clockVal(timer);
666 pastix_data->dparm[DPARM_SOLV_ENERGY] = energy;
667
668 if ( pastix_data->iparm[IPARM_VERBOSE] > PastixVerboseNot ) {
669 pastix_print( pastix_data->inter_node_procnum, 0, OUT_TIME_SOLV,
670 pastix_data->dparm[DPARM_SOLV_TIME] );
671#if defined(PASTIX_WITH_PAPI)
672 pastix_print( pastix_data->inter_node_procnum, 0, OUT_SOLVE_ENERGY,
673 pastix_print_value_deci( pastix_data->dparm[DPARM_SOLV_ENERGY] ),
674 pastix_print_unit_deci( pastix_data->dparm[DPARM_SOLV_ENERGY] ),
675 pastix_print_value_deci( pastix_data->dparm[DPARM_SOLV_ENERGY] / pastix_data->dparm[DPARM_SOLV_TIME] ),
676 pastix_print_unit_deci( pastix_data->dparm[DPARM_SOLV_ENERGY] / pastix_data->dparm[DPARM_SOLV_TIME] ) );
677#endif
678 }
679 }
680
681 return PASTIX_SUCCESS;
682}
683
684/**
685 *******************************************************************************
686 *
687 * @ingroup pastix_solve
688 *
689 * @brief Solve the given problem without applying the permutation.
690 *
691 * @warning The input vector is considered already permuted. For a solve step
692 * with permutation, see pastix_task_solve()
693 *
694 * This routine is affected by the following parameters:
695 * IPARM_VERBOSE, IPARM_FACTORIZATION, IPARM_TRANSPOSE_SOLVE.
696 *
697 *******************************************************************************
698 *
699 * @param[inout] pastix_data
700 * The pastix_data structure that describes the solver instance.
701 *
702 * @param[inout] Bp
703 * The right-and-side vectors (can be multiple rhs).
704 * On exit, the solution is stored in place of the right-hand-side vector.
705 *
706 *******************************************************************************
707 *
708 * @retval PASTIX_SUCCESS on successful exit,
709 * @retval PASTIX_ERR_BADPARAMETER if one parameter is incorrect.
710 *
711 *******************************************************************************/
712int
714 pastix_rhs_t Bp )
715{
716 return pastix_subtask_solve_adv( pastix_data,
717 pastix_data->iparm[IPARM_TRANSPOSE_SOLVE],
718 Bp );
719}
720
721/**
722 *******************************************************************************
723 *
724 * @ingroup pastix_users
725 *
726 * @brief Solve the given problem.
727 *
728 * This routine is affected by the following parameters:
729 * IPARM_VERBOSE, IPARM_FACTORIZATION, IPARM_TRANSPOSE_SOLVE.
730 *
731 *******************************************************************************
732 *
733 * @param[inout] pastix_data
734 * The pastix_data structure that describes the solver instance.
735 *
736 * @param[in] m
737 * The local number of rows of the right-and-side vectors.
738 * If the spm is centralized or replicated, m must be equal to
739 * spm->gNexp, otherwise m must be equal to spm->nexp.
740 *
741 * @param[in] nrhs
742 * The number of right-and-side vectors.
743 *
744 * @param[inout] b
745 * The right-and-side vectors (can be multiple rhs).
746 * On exit, the solution is stored in place of the right-hand-side vector.
747 *
748 * @param[in] ldb
749 * The leading dimension of the right-and-side vectors.
750 *
751 *******************************************************************************
752 *
753 * @retval PASTIX_SUCCESS on successful exit,
754 * @retval PASTIX_ERR_BADPARAMETER if one parameter is incorrect.
755 *
756 *******************************************************************************/
757int
760 pastix_int_t nrhs,
761 void *b,
762 pastix_int_t ldb )
763{
764 pastix_rhs_t Bp;
765 pastix_int_t rc;
766
767 /*
768 * Check parameters
769 */
770 if (pastix_data == NULL) {
771 pastix_print_error( "pastix_task_solve: wrong pastix_data parameter" );
773 }
774
775 if ( !(pastix_data->steps & STEP_NUMFACT) ) {
776 pastix_print_error( "pastix_task_solve: Numerical factorization hasn't been done." );
778 }
779
780 /* Compute P * b */
781 rc = pastixRhsInit( &Bp );
782 if( rc != PASTIX_SUCCESS ) {
783 return rc;
784 }
785
786 rc = pastix_subtask_applyorder( pastix_data, PastixDirForward, m, nrhs, b, ldb, Bp );
787 if( rc != PASTIX_SUCCESS ) {
788 return rc;
789 }
790
791 /* Solve A x = b */
792 rc = pastix_subtask_solve( pastix_data, Bp );
793 if( rc != PASTIX_SUCCESS ) {
794 return rc;
795 }
796
797 /* Compute P^t * b */
798 rc = pastix_subtask_applyorder( pastix_data, PastixDirBackward, m, nrhs, b, ldb, Bp );
799 if( rc != PASTIX_SUCCESS ) {
800 return rc;
801 }
802
803 rc = pastixRhsFinalize( Bp );
804 if( rc != PASTIX_SUCCESS ) {
805 return rc;
806 }
807
808 return rc;
809}
BEGIN_C_DECLS typedef int pastix_int_t
Definition datatypes.h:51
int bvec_slapmr(pastix_data_t *pastix_data, pastix_dir_t dir, pastix_int_t m, pastix_int_t n, float *A, pastix_int_t lda, pastix_rhs_t PA)
Apply a row permutation to a right hand side A (LAPACK xlatmr)
int pastixRhsInit(pastix_rhs_t *rhs)
Initialize an RHS data structure.
Definition pastix_rhs.c:44
int pastixRhsFinalize(pastix_rhs_t rhs)
Cleanup an RHS data structure.
Definition pastix_rhs.c:92
int bvec_dlapmr(pastix_data_t *pastix_data, pastix_dir_t dir, pastix_int_t m, pastix_int_t n, double *A, pastix_int_t lda, pastix_rhs_t PA)
Apply a row permutation to a right hand side A (LAPACK xlatmr)
int bvec_zlapmr(pastix_data_t *pastix_data, pastix_dir_t dir, pastix_int_t m, pastix_int_t n, pastix_complex64_t *A, pastix_int_t lda, pastix_rhs_t PA)
Apply a row permutation to a right hand side A (LAPACK xlatmr)
int bvec_clapmr(pastix_data_t *pastix_data, pastix_dir_t dir, pastix_int_t m, pastix_int_t n, pastix_complex32_t *A, pastix_int_t lda, pastix_rhs_t PA)
Apply a row permutation to a right hand side A (LAPACK xlatmr)
void kernelsTraceStart(void)
Resumes the trace module.
void kernelsTraceStop(void)
Pauses the trace module.
spm_coeftype_t pastix_coeftype_t
Arithmetic types.
Definition api.h:294
int pastixBlasSetNumThreadsOne(void)
Set the number of threads for blas calls (BLIS, MKL, OpenBLAS) to 1 and return the previous number of...
Definition api.c:1179
enum pastix_diag_e pastix_diag_t
Diagonal.
FILE * pastix_fopenw(const char *dirname, const char *filename, const char *mode)
Open a file in the unique directory of the pastix instance.
Definition api.c:251
enum pastix_dir_e pastix_dir_t
Direction.
enum pastix_uplo_e pastix_uplo_t
Upper/Lower part.
int pastixBlasSetNumThreads(int nt)
Set the number of threads for blas calls (BLIS, MKL, OpenBLAS) and return the previous number of blas...
Definition api.c:1147
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.
@ PastixDirForward
Definition api.h:513
@ PastixDirBackward
Definition api.h:514
@ PastixFactLDLH
Definition api.h:319
@ PastixFactLDLT
Definition api.h:316
@ PastixFactLLH
Definition api.h:315
@ PastixFactLU
Definition api.h:317
@ DPARM_SOLV_TIME
Definition api.h:177
@ DPARM_SOLV_ENERGY
Definition api.h:181
@ IPARM_FACTORIZATION
Definition api.h:99
@ IPARM_MIXED
Definition api.h:139
@ IPARM_VERBOSE
Definition api.h:36
@ IPARM_TRACE
Definition api.h:44
@ IPARM_TRANSPOSE_SOLVE
Definition api.h:106
@ PastixUpper
Definition api.h:466
@ PastixLower
Definition api.h:467
@ PastixLeft
Definition api.h:495
@ PastixUnit
Definition api.h:488
@ PastixNonUnit
Definition api.h:487
@ PastixConjTrans
Definition api.h:447
@ PastixNoTrans
Definition api.h:445
@ PastixTrans
Definition api.h:446
@ PastixVerboseNot
Definition api.h:220
@ PASTIX_SUCCESS
Definition api.h:367
@ PASTIX_ERR_BADPARAMETER
Definition api.h:374
@ PastixTraceSolve
Definition api.h:212
pastix_int_t baseval
Definition order.h:48
int pastix_subtask_solve_adv(pastix_data_t *pastix_data, pastix_trans_t transA, pastix_rhs_t Bp)
Solve the given problem without applying the permutation.
int pastixRhsDoubletoSingle(const pastix_rhs_t dB, pastix_rhs_t sB)
Reduces the precision of an RHS.
Definition pastix_rhs.c:162
int pastix_subtask_applyorder(pastix_data_t *pastix_data, pastix_dir_t dir, pastix_int_t m, pastix_int_t n, void *b, pastix_int_t ldb, pastix_rhs_t Bp)
Apply a permutation on the right-and-side vector before the solve step.
int pastix_subtask_diag(pastix_data_t *pastix_data, pastix_rhs_t Bp)
Apply a diagonal operation on the right-and-side vectors.
int pastix_subtask_solve(pastix_data_t *pastix_data, pastix_rhs_t Bp)
Solve the given problem without applying the permutation.
int pastix_subtask_trsm(pastix_data_t *pastix_data, pastix_side_t side, pastix_uplo_t uplo, pastix_trans_t trans, pastix_diag_t diag, pastix_rhs_t Bp)
Apply a triangular solve on the right-and-side vectors.
int pastixRhsSingleToDouble(const pastix_rhs_t sB, pastix_rhs_t dB)
Increases the precision of an RHS.
Definition pastix_rhs.c:242
int inter_node_procnum
Definition pastixdata.h:84
void ** cblkb
Definition pastixdata.h:162
SolverMatrix * solvmatr
Definition pastixdata.h:103
pastix_order_t * ordemesh
Definition pastixdata.h:98
pastix_int_t * iparm
Definition pastixdata.h:70
pastix_int_t ld
Definition pastixdata.h:160
char * dir_local
Definition pastixdata.h:111
double * dparm
Definition pastixdata.h:71
const spmatrix_t * csc
Definition pastixdata.h:90
pastix_coeftype_t flttype
Definition pastixdata.h:157
PASTIX_Comm inter_node_comm
Definition pastixdata.h:78
pastix_int_t m
Definition pastixdata.h:158
pastix_bcsc_t * bcsc
Definition pastixdata.h:102
pastix_int_t steps
Definition pastixdata.h:73
pastix_int_t n
Definition pastixdata.h:159
int8_t allocated
Definition pastixdata.h:156
int pastix_task_solve(pastix_data_t *pastix_data, pastix_int_t m, pastix_int_t nrhs, void *b, pastix_int_t ldb)
Solve the given problem.
Main PaStiX data structure.
Definition pastixdata.h:68
Main PaStiX RHS structure.
Definition pastixdata.h:155
-by-step
void sopalin_cdiag(pastix_data_t *pastix_data, sopalin_data_t *sopalin_data, pastix_rhs_t rhsb)
Apply the diagonal solve.
void sopalin_ctrsm(pastix_data_t *pastix_data, pastix_side_t side, pastix_uplo_t uplo, pastix_trans_t trans, pastix_diag_t diag, sopalin_data_t *sopalin_data, pastix_rhs_t rhsb)
Calls the sequential, static, dynamic or runtime solve according to scheduler.
void sopalin_ddiag(pastix_data_t *pastix_data, sopalin_data_t *sopalin_data, pastix_rhs_t rhsb)
Apply the diagonal solve.
void sopalin_dtrsm(pastix_data_t *pastix_data, pastix_side_t side, pastix_uplo_t uplo, pastix_trans_t trans, pastix_diag_t diag, sopalin_data_t *sopalin_data, pastix_rhs_t rhsb)
Calls the sequential, static, dynamic or runtime solve according to scheduler.
void sopalin_sdiag(pastix_data_t *pastix_data, sopalin_data_t *sopalin_data, pastix_rhs_t rhsb)
Apply the diagonal solve.
void sopalin_strsm(pastix_data_t *pastix_data, pastix_side_t side, pastix_uplo_t uplo, pastix_trans_t trans, pastix_diag_t diag, sopalin_data_t *sopalin_data, pastix_rhs_t rhsb)
Calls the sequential, static, dynamic or runtime solve according to scheduler.
void sopalin_zdiag(pastix_data_t *pastix_data, sopalin_data_t *sopalin_data, pastix_rhs_t rhsb)
Apply the diagonal solve.
void sopalin_ztrsm(pastix_data_t *pastix_data, pastix_side_t side, pastix_uplo_t uplo, pastix_trans_t trans, pastix_diag_t diag, sopalin_data_t *sopalin_data, pastix_rhs_t rhsb)
Calls the sequential, static, dynamic or runtime solve according to scheduler.
pastix_int_t faninnbr
Definition solver.h:213
pastix_int_t recvnbr
Definition solver.h:216
Solver column block structure.
Definition solver.h:203