PaStiX Handbook 6.4.0
Loading...
Searching...
No Matches
refinement.c
Go to the documentation of this file.
1/**
2 * @file refinement.c
3 *
4 * @brief A refinement example that runs iterative methods without applying
5 * the preconditioner (factorization and solve steps are not called).
6 * Based on the step-by-step example.
7 *
8 * @copyright 2015-2024 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
9 * Univ. Bordeaux. All rights reserved.
10 *
11 * @version 6.4.0
12 * @author Pierre Ramet
13 * @author Mathieu Faverge
14 * @author Tony Delarue
15 * @author Alycia Lisito
16 * @date 2024-07-05
17 *
18 * @ingroup pastix_examples
19 * @code
20 *
21 */
22#include <pastix.h>
23#include <spm.h>
24
25int 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 pastix_spm_t *spm, spm2;
33 void *x, *b, *x0 = NULL;
34 size_t size;
35 int scatter = 0;
36 int check = 1;
37 int nrhs = 1;
38 int rc = 0;
39
40 /**
41 * Initialize parameters to default values
42 */
43 pastixInitParam( iparm, dparm );
44
45 /**
46 * Get options from command line
47 */
48 pastixGetOptions( argc, argv,
49 iparm, dparm,
50 &check, &scatter, &driver, &filename );
51
52 /**
53 * Startup PaStiX
54 */
55 pastixInit( &pastix_data, MPI_COMM_WORLD, iparm, dparm );
56
57 /**
58 * Read the sparse matrix with the driver
59 */
60 spm = malloc( sizeof( pastix_spm_t ) );
61 if ( scatter ) {
62 rc = spmReadDriverDist( driver, filename, spm, MPI_COMM_WORLD );
63 }
64 else {
65 rc = spmReadDriver( driver, filename, spm );
66 }
67 free( filename );
68 if ( rc != SPM_SUCCESS ) {
69 pastixFinalize( &pastix_data );
70 return rc;
71 }
72
73 spmPrintInfo( spm, stdout );
74
75 rc = spmCheckAndCorrect( spm, &spm2 );
76 if ( rc != 0 ) {
77 spmExit( spm );
78 *spm = spm2;
79 rc = 0;
80 }
81
82 /**
83 * Generate a Fake values array if needed for the numerical part
84 */
85 if ( spm->flttype == PastixPattern ) {
86 spmGenFakeValues( spm );
87 }
88
89 /**
90 * Perform ordering, symbolic factorization, and analyze steps
91 */
92 pastix_subtask_order( pastix_data, spm, NULL );
93 pastix_subtask_symbfact( pastix_data );
94 pastix_subtask_reordering( pastix_data );
95 pastix_subtask_blend( pastix_data );
96
97 /**
98 * Normalize A matrix (optional, but recommended for low-rank functionality)
99 */
100 double normA = spmNorm( SpmFrobeniusNorm, spm );
101 spmScal( 1./normA, spm );
102
103 size = pastix_size_of( spm->flttype ) * spm->nexp * nrhs;
104 x = malloc( size );
105 b = malloc( size );
106 if ( check > 1 ) {
107 x0 = malloc( size );
108 }
109
110 pastix_subtask_spm2bcsc( pastix_data, spm );
111
112 /**
113 * Generates the b and x vector such that A * x = b
114 * Compute the norms of the initial vectors if checking purpose.
115 */
116 if ( check )
117 {
118 spmGenRHS( SpmRhsRndX, nrhs, spm, x0, spm->nexp, b, spm->nexp );
119 memcpy( x, b, size );
120 }
121 else
122 {
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 * Refine the linear system
134 */
135 pastix_task_refine( pastix_data, spm->nexp, nrhs, b, spm->nexp, x, spm->nexp );
136
137 if ( check ) {
138 rc |= spmCheckAxb( dparm[DPARM_EPSILON_REFINEMENT], nrhs, spm, x0, spm->nexp, b, spm->nexp, x, spm->nexp );
139 }
140 if (iparm[IPARM_NBITER]>=iparm[IPARM_ITERMAX]) {
141 rc |= 1;
142 }
143
144 spmExit( spm );
145 free( spm );
146 free( b );
147 free( x );
148 if ( x0 ) {
149 free( x0 );
150 }
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
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:928
void pastixInitParam(pastix_int_t *iparm, double *dparm)
Initialize the iparm and dparm arrays to their default values.
Definition api.c:420
void pastixInit(pastix_data_t **pastix_data, PASTIX_Comm pastix_comm, pastix_int_t *iparm, double *dparm)
Initialize the solver instance.
Definition api.c:905
@ DPARM_EPSILON_REFINEMENT
Definition api.h:161
@ IPARM_ITERMAX
Definition api.h:113
@ IPARM_NBITER
Definition api.h:112
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.
int pastix_subtask_spm2bcsc(pastix_data_t *pastix_data, spmatrix_t *spm)
Fill the internal block CSC 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.
Main PaStiX data structure.
Definition pastixdata.h:68