22#ifndef DOXYGEN_SHOULD_SKIP_THIS
28#include "order/order_internal.h"
30#include "sopalin/sopalin_data.h"
31#include "pastix_papi.h"
40#ifndef DOXYGEN_SHOULD_SKIP_THIS
42#if defined(PASTIX_DEBUG_SOLVE)
70 static volatile int32_t
step = 0;
72 int32_t lstep = pastix_atomic_inc_32b( &
step );
77 rc = asprintf( &fullname,
"%02d_%s.%02d.rhs",
78 lstep, name, pastix_data->
solvmatr->clustnum );
87 z_spmDensePrint( f, B->
m, B->
n, B->
b, B->
ld );
90 c_spmDensePrint( f, B->
m, B->
n, B->
b, B->
ld );
93 d_spmDensePrint( f, B->
m, B->
n, B->
b, B->
ld );
96 s_spmDensePrint( f, B->
m, B->
n, B->
b, B->
ld );
182 if (pastix_data == NULL) {
183 pastix_print_error(
"pastix_subtask_applyorder: wrong pastix_data parameter" );
187 pastix_print_error(
"pastix_subtask_applyorder: wrong Bp parameter" );
193 pastix_print_error(
"pastix_subtask_applyorder: ordermesh must be 0-based" );
197 flttype = pastix_data->
csc->flttype;
198#if defined(PASTIX_DEBUG_SOLVE)
208 pastix_rhs_dump( pastix_data,
"applyorder_Forward_input", &rhsb );
211 pastix_rhs_dump( pastix_data,
"applyorder_Backward_input", Bp );
217 case PastixComplex64:
221 case PastixComplex32:
234 pastix_print_error(
"The floating type of the rhs is not defined\n" );
238#if defined(PASTIX_DEBUG_SOLVE)
240 pastix_rhs_dump( pastix_data,
"applyorder_Forward_output", Bp );
251 pastix_rhs_dump( pastix_data,
"applyorder_Backward_output", &rhsb );
304 sopalin_data_t sopalin_data;
311 if (pastix_data == NULL) {
312 pastix_print_error(
"pastix_subtask_trsm: wrong pastix_data parameter" );
316 pastix_print_error(
"pastix_subtask_trsm: wrong Bp parameter" );
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" );
327 if ( Bp->
cblkb == NULL ) {
331 if ( nbbuffers > 0 ) {
332 Bp->
cblkb = calloc( nbbuffers,
sizeof(
void*) );
340 pastix_check_and_correct_scheduler( pastix_data );
342 sopalin_data.solvmtx = solvmtx;
347 case PastixComplex64:
353 case PastixComplex32:
374 fprintf(stderr,
"Unknown floating point arithmetic\n" );
379#if defined(PASTIX_DEBUG_SOLVE)
382 asprintf( &fullname,
"solve_trsm_%c%c%c%c",
388 pastix_rhs_dump( pastix_data, fullname, Bp );
425 sopalin_data_t sopalin_data;
431 if (pastix_data == NULL) {
432 pastix_print_error(
"pastix_subtask_diag: wrong pastix_data parameter" );
436 pastix_print_error(
"pastix_subtask_diag: wrong Bp parameter" );
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" );
448 pastix_check_and_correct_scheduler( pastix_data );
450 sopalin_data.solvmtx = pastix_data->
solvmatr;
455 case PastixComplex64:
458 case PastixComplex32:
468 fprintf(stderr,
"Unknown floating point arithmetic\n" );
473 pastix_rhs_dump( pastix_data,
"solve_diag", Bp );
529 if (pastix_data == NULL) {
530 pastix_print_error(
"pastix_task_solve: wrong pastix_data parameter" );
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" );
538 bcsc = pastix_data->
bcsc;
560 if ( ((bcsc->flttype == PastixComplex32) || (bcsc->flttype == PastixComplex64)) &&
567 (transA != transfact) )
569 pastix_print_error(
"pastix_task_solve: transA incompatible with the factotype used (require extra conj(L) not handled)" );
574 double timer, energy;
584 ( ( Bp->
flttype == PastixComplex64 ) || ( Bp->
flttype == PastixDouble ) ) )
652 ( ( Bp->
flttype == PastixComplex64 ) || ( Bp->
flttype == PastixDouble ) ) )
663 energy = papiEnergyStop();
671#if defined(PASTIX_WITH_PAPI)
770 if (pastix_data == NULL) {
771 pastix_print_error(
"pastix_task_solve: wrong pastix_data parameter" );
775 if ( !(pastix_data->
steps & STEP_NUMFACT) ) {
776 pastix_print_error(
"pastix_task_solve: Numerical factorization hasn't been done." );
BEGIN_C_DECLS typedef int pastix_int_t
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.
int pastixRhsFinalize(pastix_rhs_t rhs)
Cleanup an RHS data structure.
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.
int pastixBlasSetNumThreadsOne(void)
Set the number of threads for blas calls (BLIS, MKL, OpenBLAS) to 1 and return the previous number of...
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.
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...
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.
@ PASTIX_ERR_BADPARAMETER
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.
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.
pastix_order_t * ordemesh
pastix_coeftype_t flttype
PASTIX_Comm inter_node_comm
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.
Main PaStiX RHS structure.
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.
Solver column block structure.