Schur usage example.
/
#include <spm.h>
#include <lapacke.h>
void
void *S,
int **ipiv )
{
int info = 0;
assert( ipiv != NULL );
*ipiv = malloc( N * sizeof(int) );
}
switch (flttype) {
case PastixFloat:
switch (factotype) {
info = LAPACKE_spotrf_work( LAPACK_COL_MAJOR, 'L', N, S, lds );
break;
info = LAPACKE_sgetrf_work( LAPACK_COL_MAJOR, N, N, S, lds, *ipiv );
break;
default:
fprintf(stderr, "Factorization type not handled by Schur example\n");
}
break;
case PastixComplex32:
switch (factotype) {
info = LAPACKE_cpotrf_work( LAPACK_COL_MAJOR, 'L', N, S, lds );
break;
info = LAPACKE_cgetrf_work( LAPACK_COL_MAJOR, N, N, S, lds, *ipiv );
break;
default:
fprintf(stderr, "Factorization type not handled by Schur example\n");
}
break;
case PastixComplex64:
switch (factotype) {
info = LAPACKE_zpotrf_work( LAPACK_COL_MAJOR, 'L', N, S, lds );
break;
info = LAPACKE_zgetrf_work( LAPACK_COL_MAJOR, N, N, S, lds, *ipiv );
break;
default:
fprintf(stderr, "Factorization type not handled by Schur example\n");
}
break;
case PastixDouble:
switch (factotype) {
info = LAPACKE_dpotrf_work( LAPACK_COL_MAJOR, 'L', N, S, lds );
break;
info = LAPACKE_dgetrf_work( LAPACK_COL_MAJOR, N, N, S, lds, *ipiv );
break;
default:
fprintf(stderr, "Factorization type not handled by Schur example\n");
}
break;
default:
fprintf(stderr, "Incorrect arithmetic type\n");
}
if (info != 0) {
fprintf(stderr, "Error in schurFactorize with info =%d\n", info );
}
return;
}
void
void *S,
void *bptr,
int **ipiv )
{
int info = 0;
assert(ipiv != NULL);
switch (flttype) {
case PastixFloat:
{
float *b = (float *)bptr;
switch (factotype) {
info = LAPACKE_spotrs_work( LAPACK_COL_MAJOR, 'L', Nschur, NRHS, S, lds, b, ldb );
break;
info = LAPACKE_sgetrs_work( LAPACK_COL_MAJOR, 'N', Nschur, NRHS, S, lds, *ipiv, b, ldb );
break;
default:
fprintf(stderr, "Factorization type not handled by Schur example\n");
}
}
break;
case PastixComplex32:
{
switch (factotype) {
info = LAPACKE_cpotrs_work( LAPACK_COL_MAJOR, 'L', Nschur, NRHS, S, lds, b, ldb );
break;
info = LAPACKE_cgetrs_work( LAPACK_COL_MAJOR, 'N', Nschur, NRHS, S, lds, *ipiv, b, ldb );
break;
default:
fprintf(stderr, "Factorization type not handled by Schur example\n");
}
}
break;
case PastixComplex64:
{
pastix_complex64_t *b = (pastix_complex64_t *)bptr;
switch (factotype) {
info = LAPACKE_zpotrs_work( LAPACK_COL_MAJOR, 'L', Nschur, NRHS, S, lds, b, ldb );
break;
info = LAPACKE_zgetrs_work( LAPACK_COL_MAJOR, 'N', Nschur, NRHS, S, lds, *ipiv, b, ldb );
break;
default:
fprintf(stderr, "Factorization type not handled by Schur example\n");
}
}
break;
case PastixDouble:
{
double *b = (double *)bptr;
switch (factotype) {
info = LAPACKE_dpotrs_work( LAPACK_COL_MAJOR, 'L', Nschur, NRHS, S, lds, b, ldb );
break;
info = LAPACKE_dgetrs_work( LAPACK_COL_MAJOR, 'N', Nschur, NRHS, S, lds, *ipiv, b, ldb );
break;
default:
fprintf(stderr, "Factorization type not handled by Schur example\n");
}
}
break;
default:
fprintf(stderr, "Incorrect arithmetic type\n");
}
if (*ipiv != NULL) {
free( *ipiv );
*ipiv = NULL;
}
if (info != 0) {
fprintf(stderr, "Error in schurSolve with info =%d\n", info );
}
return;
}
int main (int argc, char **argv)
{
double dparm[DPARM_SIZE];
spm_driver_t driver;
char *filename = NULL;
spmatrix_t *spm, spm2;
void *x, *b, *S, *x0 = NULL;
size_t size;
int scatter = 0;
int check = 1;
int nrhs = 1;
int rc = 0;
int *ipiv = NULL;
nschur = (nschur > 5000) ? 5000 : nschur;
if ( nschur > 0 ) {
for (i=0; i<nschur; i++) {
list[i] = i+baseval;
}
free( list );
}
}
spmScalMat( 1./normA, spm, nrhs, b, spm->nexp );
}
spm->nexp, nrhs, x, ldb, Xp );
}
}
{
void *schur_x;
size = pastix_size_of( spm->flttype ) * nschur * nrhs;
schur_x = malloc( size );
nschur, nrhs, S, lds, schur_x, nschur, &ipiv );
free( schur_x );
}
Xp );
}
Xp );
}
spm->nexp, nrhs, x, ldb, Xp );
if ( check )
{
if ( x0 ) {
free( x0 );
}
}
spmExit( spm );
free( spm );
free( S );
free( x );
free( b );
return rc;
}
BEGIN_C_DECLS typedef int pastix_int_t
float _Complex pastix_complex32_t
int pastixRhsInit(pastix_rhs_t *rhs)
Initialize an RHS data structure.
int pastixRhsFinalize(pastix_rhs_t rhs)
Cleanup an RHS data structure.
spm_coeftype_t pastix_coeftype_t
Arithmetic types.
void pastixFinalize(pastix_data_t **pastix_data)
Finalize the solver instance.
enum pastix_diag_e pastix_diag_t
Diagonal.
enum pastix_factotype_e pastix_factotype_t
Factorization algorithms available for IPARM_FACTORIZATION parameter.
@ DPARM_EPSILON_REFINEMENT
void pastixSetSchurUnknownList(pastix_data_t *pastix_data, pastix_int_t n, const pastix_int_t *list)
Set the list of unknowns that belongs to the schur complement.
int pastixRhsSchurSet(const pastix_data_t *pastix_data, pastix_int_t m, pastix_int_t n, void *B, pastix_int_t ldb, pastix_rhs_t rhsB)
Set the vector in an RHS data structure.
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_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 b)
Apply a triangular solve on the right-and-side vectors.
int pastixRhsSchurGet(const pastix_data_t *pastix_data, pastix_int_t m, pastix_int_t n, pastix_rhs_t rhsB, void *B, pastix_int_t ldb)
Get the vector in an RHS data structure.
Main PaStiX data structure.
Main PaStiX RHS structure.