PaStiX Handbook  6.3.2
simple.c
Go to the documentation of this file.
1 /**
2  * @file simple.c
3  *
4  * @brief A simple example that reads the matrix and then runs pastix in one call.
5  *
6  * @copyright 2015-2023 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
7  * Univ. Bordeaux. All rights reserved.
8  *
9  * @version 6.3.2
10  * @author Gregoire Pichon
11  * @author Mathieu Faverge
12  * @author Matias Hastaran
13  * @author Pierre Ramet
14  * @author Theophile Terraz
15  * @author Tony Delarue
16  * @author Alycia Lisito
17  * @date 2023-07-21
18  *
19  * @ingroup pastix_examples
20  * @code
21  *
22  */
23 #include <pastix.h>
24 #include <spm.h>
25 
26 int main (int argc, char **argv)
27 {
28  pastix_data_t *pastix_data = NULL; /*< Pointer to the storage structure required by pastix */
29  pastix_int_t iparm[IPARM_SIZE]; /*< Integer in/out parameters for pastix */
30  double dparm[DPARM_SIZE]; /*< Floating in/out parameters for pastix */
31  spm_driver_t driver;
32  char *filename = NULL;
33  spmatrix_t *spm, spm2;
34  void *x, *b, *x0 = NULL;
35  size_t size;
36  int check = 1;
37  int scatter = 0;
38  int nrhs = 1;
39  int rc = 0;
40 
41  /**
42  * Initialize parameters to default values
43  */
44  pastixInitParam( iparm, dparm );
45 
46  /**
47  * Get options from command line
48  */
49  pastixGetOptions( argc, argv,
50  iparm, dparm,
51  &check, &scatter, &driver, &filename );
52 
53  /**
54  * Startup PaStiX
55  */
56  pastixInit( &pastix_data, MPI_COMM_WORLD, iparm, dparm );
57 
58  /**
59  * Read the sparse matrix with the driver
60  */
61  spm = malloc( sizeof( spmatrix_t ) );
62  if ( scatter ) {
63  rc = spmReadDriverDist( driver, filename, spm, MPI_COMM_WORLD );
64  }
65  else {
66  rc = spmReadDriver( driver, filename, spm );
67  }
68  free( filename );
69  if ( rc != SPM_SUCCESS ) {
70  pastixFinalize( &pastix_data );
71  return rc;
72  }
73 
74  spmPrintInfo( spm, stdout );
75 
76  rc = spmCheckAndCorrect( spm, &spm2 );
77  if ( rc != 0 ) {
78  spmExit( spm );
79  *spm = spm2;
80  rc = 0;
81  }
82 
83  /**
84  * Generate a Fake values array if needed for the numerical part
85  */
86  if ( spm->flttype == SpmPattern ) {
87  spmGenFakeValues( spm );
88  }
89 
90  /**
91  * Perform ordering, symbolic factorization, and analyze steps
92  */
93  pastix_task_analyze( pastix_data, spm );
94 
95  /**
96  * Normalize A matrix (optional, but recommended for low-rank functionality)
97  */
98  double normA = spmNorm( SpmFrobeniusNorm, spm );
99  spmScal( 1./normA, spm );
100 
101  /**
102  * Perform the numerical factorization
103  */
104  pastix_task_numfact( pastix_data, spm );
105 
106  /**
107  * Generates the b and x vector such that A * x = b
108  * Compute the norms of the initial vectors if checking purpose.
109  */
110  size = pastix_size_of( spm->flttype ) * spm->nexp * nrhs;
111  x = malloc( size );
112  b = malloc( size );
113 
114  if ( check )
115  {
116  if ( check > 1 ) {
117  x0 = malloc( size );
118  }
119  spmGenRHS( SpmRhsRndX, nrhs, spm, x0, spm->nexp, b, spm->nexp );
120  memcpy( x, b, size );
121  }
122  else {
123  spmGenRHS( SpmRhsRndB, nrhs, spm, NULL, spm->nexp, x, spm->nexp );
124 
125  /* Apply also normalization to b vectors */
126  spmScalMat( 1./normA, spm, nrhs, b, spm->nexp );
127 
128  /* Save b for refinement */
129  memcpy( b, x, size );
130  }
131 
132  /**
133  * Solve the linear system (and perform the optional refinement)
134  */
135  pastix_task_solve( pastix_data, spm->nexp, nrhs, x, spm->nexp );
136  pastix_task_refine( pastix_data, spm->nexp, nrhs, b, spm->nexp, x, spm->nexp );
137 
138  if ( check )
139  {
140  rc = spmCheckAxb( dparm[DPARM_EPSILON_REFINEMENT], nrhs, spm, x0, spm->nexp, b, spm->nexp, x, spm->nexp );
141 
142  if ( x0 ) {
143  free( x0 );
144  }
145  }
146 
147  spmExit( spm );
148  free( spm );
149  free( x );
150  free( b );
151  pastixFinalize( &pastix_data );
152 
153  return rc;
154 }
155 
156 /**
157  * @endcode
158  */
BEGIN_C_DECLS typedef int pastix_int_t
Definition: datatypes.h:51
void pastixFinalize(pastix_data_t **pastix_data)
Finalize the solver instance.
Definition: api.c:919
void pastixInitParam(pastix_int_t *iparm, double *dparm)
Initialize the iparm and dparm arrays to their default values.
Definition: api.c:411
void pastixInit(pastix_data_t **pastix_data, PASTIX_Comm pastix_comm, pastix_int_t *iparm, double *dparm)
Initialize the solver instance.
Definition: api.c:896
@ DPARM_EPSILON_REFINEMENT
Definition: api.h:161
void pastixGetOptions(int argc, char **argv, pastix_int_t *iparm, double *dparm, int *check, int *scatter, spm_driver_t *driver, char **filename)
PaStiX helper function to read command line options in examples.
Definition: get_options.c:149
int pastix_task_refine(pastix_data_t *pastix_data, pastix_int_t n, pastix_int_t nrhs, void *B, pastix_int_t ldb, void *X, pastix_int_t ldx)
Perform iterative refinement.
int pastix_task_analyze(pastix_data_t *pastix_data, const spmatrix_t *spm)
Perform all the preprocessing steps: ordering, symbolic factorization, reordering,...
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.
int pastix_task_numfact(pastix_data_t *pastix_data, spmatrix_t *spm)
Perform all the numerical factorization steps: fill the internal block CSC and the solver matrix stru...
Main PaStiX data structure.
Definition: pastixdata.h:67