PaStiX Handbook  6.3.2
step-by-step.c
Go to the documentation of this file.
1 /**
2  * @file step-by-step.c
3  *
4  * @brief A step-by-step example that runs one full analyze (ordering, symbolic
5  * factorization, analyze), then loops over 2 factorizations that are both
6  * used for 2 solves each.
7  *
8  * @copyright 2015-2023 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
9  * Univ. Bordeaux. All rights reserved.
10  *
11  * @version 6.3.2
12  * @author Mathieu Faverge
13  * @author Pierre Ramet
14  * @author Tony Delarue
15  * @author Alycia Lisito
16  * @date 2023-11-07
17  *
18  * @ingroup pastix_examples
19  * @code
20  *
21  */
22 #include <pastix.h>
23 #include <spm.h>
24 
25 int main (int argc, char **argv)
26 {
27  pastix_data_t *pastix_data = NULL; /*< Pointer to the storage structure required by pastix */
28  pastix_int_t iparm[IPARM_SIZE]; /*< Integer in/out parameters for pastix */
29  double dparm[DPARM_SIZE]; /*< Floating in/out parameters for pastix */
30  spm_driver_t driver;
31  char *filename = NULL;
32  spmatrix_t *spm, spm2;
33  void *x, *b, *x0 = NULL;
34  size_t size;
35  int scatter = 0;
36  int check = 1;
37  int nrhs = 10;
38  int rc = PASTIX_SUCCESS;
39  int nfact = 2;
40  int nsolv = 2;
41  long i, j;
42 
43  /**
44  * Initialize parameters to default values
45  */
46  pastixInitParam( iparm, dparm );
47 
48  /**
49  * Get options from command line
50  */
51  pastixGetOptions( argc, argv,
52  iparm, dparm,
53  &check, &scatter, &driver, &filename );
54 
55  /**
56  * Startup PaStiX
57  */
58  pastixInit( &pastix_data, MPI_COMM_WORLD, iparm, dparm );
59 
60  /**
61  * Read the sparse matrix with the driver
62  */
63  spm = malloc( sizeof( spmatrix_t ) );
64  if ( scatter ) {
65  rc = spmReadDriverDist( driver, filename, spm, MPI_COMM_WORLD );
66  }
67  else {
68  rc = spmReadDriver( driver, filename, spm );
69  }
70  free( filename );
71  if ( rc != SPM_SUCCESS ) {
72  pastixFinalize( &pastix_data );
73  return rc;
74  }
75 
76  spmPrintInfo( spm, stdout );
77 
78  rc = spmCheckAndCorrect( spm, &spm2 );
79  if ( rc != 0 ) {
80  spmExit( spm );
81  *spm = spm2;
82  rc = 0;
83  }
84 
85  /**
86  * Generate a Fake values array if needed for the numerical part
87  */
88  if ( spm->flttype == SpmPattern ) {
89  spmGenFakeValues( spm );
90  }
91 
92  /**
93  * Perform ordering, symbolic factorization, and analyze steps
94  */
95  pastix_subtask_order( pastix_data, spm, NULL );
96  pastix_subtask_symbfact( pastix_data );
97  pastix_subtask_reordering( pastix_data );
98  pastix_subtask_blend( pastix_data );
99 
100  /**
101  * Normalize A matrix (optional, but recommended for low-rank functionality)
102  */
103  double normA = spmNorm( SpmFrobeniusNorm, spm );
104  spmScal( 1./normA, spm );
105 
106  size = pastix_size_of( spm->flttype ) * spm->nexp * nrhs;
107  x = malloc( size );
108  b = malloc( size );
109  if ( check > 1 ) {
110  x0 = malloc( size );
111  }
112 
113  /* Do nfact factorization */
114  for (i = 0; i < nfact; i++)
115  {
116  /**
117  * Perform the numerical factorization
118  */
119  pastix_subtask_spm2bcsc( pastix_data, spm );
120  pastix_subtask_bcsc2ctab( pastix_data );
121  pastix_subtask_sopalin( pastix_data );
122  for (j = 0; j < nsolv; j++)
123  {
124  /**
125  * Generates the b and x vector such that A * x = b
126  * Compute the norms of the initial vectors if checking purpose.
127  */
128  if ( check )
129  {
130  spmGenRHS( SpmRhsRndX, nrhs, spm, x0, spm->nexp, b, spm->nexp );
131  memcpy( x, b, size );
132  }
133  else {
134  spmGenRHS( SpmRhsRndB, nrhs, spm, NULL, spm->nexp, x, spm->nexp );
135 
136  /* Apply also normalization to b vectors */
137  spmScalMat( 1./normA, spm, nrhs, b, spm->nexp );
138 
139  /* Save b for refinement */
140  memcpy( b, x, size );
141  }
142 
143  /**
144  * Solve the linear system
145  */
146  pastix_task_solve( pastix_data, spm->nexp, nrhs, x, spm->nexp );
147  pastix_task_refine( pastix_data, spm->nexp, nrhs, b, spm->nexp, x, spm->nexp );
148 
149  if ( check ) {
150  rc |= spmCheckAxb( dparm[DPARM_EPSILON_REFINEMENT], nrhs, spm, x0, spm->nexp, b, spm->nexp, x, spm->nexp );
151  }
152  }
153  }
154 
155  spmExit( spm );
156  free( spm );
157  free( b );
158  free( x );
159  if ( x0 ) {
160  free( x0 );
161  }
162  pastixFinalize( &pastix_data );
163 
164  return rc;
165 }
166 
167 /**
168  * @endcode
169  */
BEGIN_C_DECLS typedef int pastix_int_t
Definition: datatypes.h:51
int pastix_subtask_symbfact(pastix_data_t *pastix_data)
Computes the symbolic factorization step.
int pastix_subtask_order(pastix_data_t *pastix_data, const spmatrix_t *spm, pastix_order_t *myorder)
Computes the ordering of the given graph in parameters.
int pastix_subtask_blend(pastix_data_t *pastix_data)
Compute the proportional mapping and the final solver structure.
int pastix_subtask_reordering(pastix_data_t *pastix_data)
Apply the reordering step to compact off-diagonal blocks.
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
@ PASTIX_SUCCESS
Definition: api.h:367
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_subtask_sopalin(pastix_data_t *pastix_data)
Factorize the given problem using Cholesky or LU decomposition.
int pastix_subtask_spm2bcsc(pastix_data_t *pastix_data, spmatrix_t *spm)
Fill the internal block CSC structure.
int pastix_subtask_bcsc2ctab(pastix_data_t *pastix_data)
Fill the internal solver matrix structure.
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_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:67