PaStiX Handbook  6.3.2
pastix_task_refine.c
Go to the documentation of this file.
1 /**
2  *
3  * @file pastix_task_refine.c
4  *
5  * PaStiX refinement functions implementations.
6  *
7  * @copyright 2015-2023 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
8  * Univ. Bordeaux. All rights reserved.
9  *
10  * @version 6.3.2
11  * @author Mathieu Faverge
12  * @author Pierre Ramet
13  * @author Xavier Lacoste
14  * @author Gregoire Pichon
15  * @author Tony Delarue
16  * @author Alycia Lisito
17  * @date 2023-07-21
18  *
19  **/
20 #include "common.h"
21 #include "bcsc/bcsc.h"
26 #include "pastix/order.h"
27 
28 /**
29  *******************************************************************************
30  *
31  * @ingroup pastix_dev_refine
32  *
33  * @brief Select the refinement function to call depending on the matrix type
34  * and the precision
35  *
36  *******************************************************************************/
38 
39 #ifndef DOXYGEN_SHOULD_SKIP_THIS
40 static refine_fct_t sopalinRefine[4][4] =
41 {
42  /* PastixRefineGMRES */
43  {
48  },
49  /* PastixRefineCG */
50  {
51  s_grad_smp,
52  d_grad_smp,
53  c_grad_smp,
55  },
56  /* PastixRefineSR */
57  {
62  },
63  /* PastixRefineBiCGSTAB */
64  {
69  }
70 };
71 #endif /* DOXYGEN_SHOULD_SKIP_THIS */
72 
73 /**
74  *******************************************************************************
75  *
76  * @ingroup pastix_refine
77  *
78  * @brief Perform the iterative refinement without apply the permutations.
79  *
80  * This routine is affected by the following parameters:
81  * IPARM_REFINEMENT, DPARM_EPSILON_REFINEMENT
82  *
83  *******************************************************************************
84  *
85  * @param[in] pastix_data
86  * The PaStiX data structure that describes the solver instance.
87  *
88  * @param[inout] Bp
89  * The right hand side matrix of size ldb-by-nrhs.
90  * B is noted as inout, as permutation might be performed on the
91  * matrix. On exit, the matrix is restored as it was on entry.
92  *
93  * @param[inout] Xp
94  * The matrix x of size ldx-by-nrhs.
95  * On entry, the initial guess x0 for the refinement step, that may be
96  * the solution returned by the solve step or any other initial guess.
97  * On exit, contains the final solution after the iterative refinement.
98  *
99  *******************************************************************************
100  *
101  * @retval PASTIX_SUCCESS on successful exit,
102  * @retval PASTIX_ERR_BADPARAMETER if one parameter is incorrect,
103  *
104  *******************************************************************************/
105 int
107  pastix_rhs_t Bp,
108  pastix_rhs_t Xp )
109 {
110  pastix_int_t *iparm = pastix_data->iparm;
111  pastix_bcsc_t *bcsc = pastix_data->bcsc;
112  double timer;
113  pastix_int_t n = Bp->m;
114  pastix_int_t nrhs = Bp->n;
115  pastix_int_t ldb = Bp->ld;
116  pastix_int_t ldx = Xp->ld;
117  const void *b = Bp->b;
118  void *x = Xp->b;
119 
120  if (nrhs > 1)
121  {
122  pastix_print_warning( "Refinement works only with 1 rhs, We will iterate on each RHS one by one\n" );
123  }
124 
125  if ( (pastix_data->schur_n > 0) && (iparm[IPARM_SCHUR_SOLV_MODE] != PastixSolvModeLocal))
126  {
127  fprintf(stderr, "Refinement is not available with Schur complement when non local solve is required\n");
129  }
130 
131  /* Prepare the refinement threshold, if not set by the user */
132  if ( pastix_data->dparm[DPARM_EPSILON_REFINEMENT] < 0. ) {
133  int isDouble = (bcsc->flttype == PastixDouble) || (bcsc->flttype == PastixComplex64);
134  if ( (!isDouble) ) {
135  pastix_data->dparm[DPARM_EPSILON_REFINEMENT] = 1e-6;
136  }
137  else {
138  pastix_data->dparm[DPARM_EPSILON_REFINEMENT] = 1e-12;
139  }
140  }
141 
142  clockSyncStart( timer, pastix_data->inter_node_comm );
143  {
144  refine_fct_t refinefct = sopalinRefine[iparm[IPARM_REFINEMENT]][pastix_data->bcsc->flttype -2];
145  char *xptr = (char *)x;
146  char *bptr = (char *)b;
147  size_t shiftx, shiftb;
148  int i;
149 
150  shiftx = ldx * pastix_size_of( Xp->flttype );
151  shiftb = ldb * pastix_size_of( Bp->flttype );
152  Bp->n = 1;
153  Xp->n = 1;
154 
155  for(i=0; i<nrhs; i++, xptr += shiftx, bptr += shiftb ) {
156  pastix_int_t it;
157  Bp->b = bptr;
158  Xp->b = xptr;
159  it = refinefct( pastix_data, Xp, Bp );
160  pastix_data->iparm[IPARM_NBITER] = pastix_imax( it, pastix_data->iparm[IPARM_NBITER] );
161  }
162 
163  Bp->n = nrhs;
164  Bp->b = (void*)b;
165  Xp->n = nrhs;
166  Xp->b = (void*)x;
167  }
168  clockSyncStop( timer, pastix_data->inter_node_comm );
169 
170  pastix_data->dparm[DPARM_REFINE_TIME] = clockVal(timer);
171  if ( iparm[IPARM_VERBOSE] > PastixVerboseNot ) {
172  pastix_print( pastix_data->inter_node_procnum,
173  0, OUT_TIME_REFINE,
174  pastix_data->dparm[DPARM_REFINE_TIME] );
175  }
176 
177  (void)n;
178  return PASTIX_SUCCESS;
179 }
180 
181 /**
182  *******************************************************************************
183  *
184  * @ingroup pastix_users
185  *
186  * @brief Perform iterative refinement.
187  *
188  * This routine performs the permutation of x, and b before and after the
189  * iterative refinement solution. To prevent extra permuation to happen, see
190  * pastix_subtask_refine().
191  * This routine is affected by the following parameters:
192  * IPARM_REFINEMENT, DPARM_EPSILON_REFINEMENT
193  *
194  *******************************************************************************
195  *
196  * @param[in] pastix_data
197  * The PaStiX data structure that describes the solver instance.
198  *
199  * @param[in] m
200  * The size of system to solve, and the number of rows of both
201  * matrices b and x.
202  *
203  * @param[in] nrhs
204  * The number of right hand side members, and the number of columns of
205  * b and x.
206  *
207  * @param[inout] b
208  * The right hand side matrix of size ldb-by-nrhs.
209  * B is noted as inout, as permutation might be performed on the
210  * matrix. On exit, the matrix is restored as it was on entry.
211  *
212  * @param[in] ldb
213  * The leading dimension of the matrix b. ldb >= n.
214  *
215  * @param[inout] x
216  * The matrix x of size ldx-by-nrhs.
217  * On entry, the initial guess x0 for the refinement step, that may be
218  * the solution returned by the solve step or any other initial guess.
219  * On exit, contains the final solution after the iterative refinement.
220  *
221  * @param[in] ldx
222  * The leading dimension of the matrix x. ldx >= n.
223  *
224  *******************************************************************************
225  *
226  * @retval PASTIX_SUCCESS on successful exit,
227  * @retval PASTIX_ERR_BADPARAMETER if one parameter is incorrect,
228  *
229  *******************************************************************************/
230 int
232  pastix_int_t m,
233  pastix_int_t nrhs,
234  void *b,
235  pastix_int_t ldb,
236  void *x,
237  pastix_int_t ldx )
238 {
239  pastix_int_t *iparm = pastix_data->iparm;
240  pastix_rhs_t Bp, Xp;
241  int rc;
242 
243  if ( (pastix_data->schur_n > 0) && (iparm[IPARM_SCHUR_SOLV_MODE] != PastixSolvModeLocal))
244  {
245  fprintf(stderr, "Refinement is not available with Schur complement when non local solve is required\n");
247  }
248 
249  /* Compute P * b */
250  rc = pastixRhsInit( &Bp );
251  if( rc != PASTIX_SUCCESS ) {
252  return rc;
253  }
254 
255  rc = pastix_subtask_applyorder( pastix_data, PastixDirForward, m, nrhs, b, ldb, Bp );
256  if( rc != PASTIX_SUCCESS ) {
257  return rc;
258  }
259 
260  /* Compute P * x */
261  rc = pastixRhsInit( &Xp );
262  if( rc != PASTIX_SUCCESS ) {
263  return rc;
264  }
265 
266  rc = pastix_subtask_applyorder( pastix_data, PastixDirForward, m, nrhs, x, ldx, Xp );
267  if( rc != PASTIX_SUCCESS ) {
268  return rc;
269  }
270 
271  /* Performe the iterative refinement */
272  rc = pastix_subtask_refine( pastix_data, Bp, Xp );
273  if( rc != PASTIX_SUCCESS ) {
274  return rc;
275  }
276 
277  /* Compute P * b */
278  rc = pastix_subtask_applyorder( pastix_data, PastixDirBackward, m, nrhs, b, ldb, Bp );
279  if( rc != PASTIX_SUCCESS ) {
280  return rc;
281  }
282  rc = pastixRhsFinalize( Bp );
283  if( rc != PASTIX_SUCCESS ) {
284  return rc;
285  }
286 
287  /* Compute P * x */
288  rc = pastix_subtask_applyorder( pastix_data, PastixDirBackward, m, nrhs, x, ldx, Xp );
289  if( rc != PASTIX_SUCCESS ) {
290  return rc;
291  }
292 
293  rc = pastixRhsFinalize( Xp );
294  if( rc != PASTIX_SUCCESS ) {
295  return rc;
296  }
297 
298  (void)m;
299  return rc;
300 }
301 
302 /**
303  *******************************************************************************
304  *
305  * @ingroup pastix_users
306  *
307  * @brief Performs solve and iterative refinement without unnecessary
308  * permutations.
309  *
310  * This routine performs the permutation of x and b before and after the solve
311  * and iterative refinement solution. This prevent the permutation of x at
312  * the end of the solve and the permutations of x and b at the beginnning of the
313  * refinement when the user wants to perfom both solve and refinement.
314  * This routine is affected by the following parameters:
315  * IPARM_REFINEMENT, DPARM_EPSILON_REFINEMENT, IPARM_VERBOSE,
316  * IPARM_FACTORIZATION, IPARM_TRANSPOSE_SOLVE
317  *
318  *******************************************************************************
319  *
320  * @param[in] pastix_data
321  * The PaStiX data structure that describes the solver instance.
322  *
323  * @param[in] m
324  * The size of system to solve, and the number of rows of both
325  * matrices b and x.
326  *
327  * @param[in] nrhs
328  * The number of right hand side members, and the number of columns of
329  * b and x.
330  *
331  * @param[inout] b
332  * The right hand side matrix of size ldb-by-nrhs.
333  * B is noted as inout, as permutation might be performed on the
334  * matrix. On exit, the matrix is restored as it was on entry.
335  *
336  * @param[in] ldb
337  * The leading dimension of the matrix b. ldb >= n.
338  *
339  * @param[inout] x
340  * The matrix x of size ldx-by-nrhs.
341  * On entry, a copy of b.
342  * On exit, contains the final solution after the iterative refinement.
343  *
344  * @param[in] ldx
345  * The leading dimension of the matrix x. ldx >= n.
346  *
347  *******************************************************************************
348  *
349  * @retval PASTIX_SUCCESS on successful exit,
350  * @retval PASTIX_ERR_BADPARAMETER if one parameter is incorrect,
351  *
352  *******************************************************************************/
353 int
355  pastix_int_t m,
356  pastix_int_t nrhs,
357  void *b,
358  pastix_int_t ldb,
359  void *x,
360  pastix_int_t ldx )
361 {
362  pastix_int_t *iparm;
363  pastix_rhs_t Bp, Xp;
364  int rc;
365 
366  /*
367  * Checks parameters
368  */
369  if (pastix_data == NULL) {
370  pastix_print_error( "pastix_task_solve_and_refine: wrong pastix_data parameter" );
372  }
373  iparm = pastix_data->iparm;
374 
375  if ( !(pastix_data->steps & STEP_NUMFACT) ) {
376  pastix_print_error( "pastix_task_solve_and_refine: Numerical factorization hasn't been done." );
378  }
379 
380  if ( (pastix_data->schur_n > 0) && (iparm[IPARM_SCHUR_SOLV_MODE] != PastixSolvModeLocal))
381  {
382  fprintf(stderr, "Refinement is not available with Schur complement when non local solve is required\n");
384  }
385 
386 
387  /* Computes P * b */
388  rc = pastixRhsInit( &Bp );
389  if( rc != PASTIX_SUCCESS ) {
390  return rc;
391  }
392 
393  rc = pastix_subtask_applyorder( pastix_data, PastixDirForward, m, nrhs, b, ldb, Bp );
394  if( rc != PASTIX_SUCCESS ) {
395  return rc;
396  }
397 
398  /* Computes P * x */
399  rc = pastixRhsInit( &Xp );
400  if( rc != PASTIX_SUCCESS ) {
401  return rc;
402  }
403 
404  rc = pastix_subtask_applyorder( pastix_data, PastixDirForward, m, nrhs, x, ldx, Xp );
405  if( rc != PASTIX_SUCCESS ) {
406  return rc;
407  }
408 
409  /* Solves A x = b */
410  rc = pastix_subtask_solve( pastix_data, Xp );
411  if( rc != PASTIX_SUCCESS ) {
412  return rc;
413  }
414 
415  /* Performs the iterative refinement */
416  rc = pastix_subtask_refine( pastix_data, Bp, Xp );
417  if( rc != PASTIX_SUCCESS ) {
418  return rc;
419  }
420 
421  /* Computes P^t * P * b */
422  rc = pastix_subtask_applyorder( pastix_data, PastixDirBackward, m, nrhs, b, ldb, Bp );
423  if( rc != PASTIX_SUCCESS ) {
424  return rc;
425  }
426  rc = pastixRhsFinalize( Bp );
427  if( rc != PASTIX_SUCCESS ) {
428  return rc;
429  }
430 
431  /* Computes P^t * P * x */
432  rc = pastix_subtask_applyorder( pastix_data, PastixDirBackward, m, nrhs, x, ldx, Xp );
433  if( rc != PASTIX_SUCCESS ) {
434  return rc;
435  }
436 
437  rc = pastixRhsFinalize( Xp );
438  if( rc != PASTIX_SUCCESS ) {
439  return rc;
440  }
441 
442  (void)m;
443  return rc;
444 }
BEGIN_C_DECLS typedef int pastix_int_t
Definition: datatypes.h:51
int pastixRhsInit(pastix_rhs_t *rhs)
Initialize an RHS data structure.
Definition: pastix_rhs.c:41
int pastixRhsFinalize(pastix_rhs_t rhs)
Cleanup an RHS data structure.
Definition: pastix_rhs.c:86
@ PastixDirForward
Definition: api.h:513
@ PastixDirBackward
Definition: api.h:514
@ DPARM_REFINE_TIME
Definition: api.h:182
@ DPARM_EPSILON_REFINEMENT
Definition: api.h:161
@ IPARM_REFINEMENT
Definition: api.h:111
@ IPARM_VERBOSE
Definition: api.h:36
@ IPARM_NBITER
Definition: api.h:112
@ IPARM_SCHUR_SOLV_MODE
Definition: api.h:107
@ PastixVerboseNot
Definition: api.h:220
@ PASTIX_SUCCESS
Definition: api.h:367
@ PASTIX_ERR_BADPARAMETER
Definition: api.h:374
pastix_int_t d_gmres_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
pastix_int_t c_gmres_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
pastix_int_t s_grad_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
Definition: s_refine_grad.c:49
pastix_int_t s_gmres_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
pastix_int_t d_pivot_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
pastix_int_t s_pivot_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
pastix_int_t z_grad_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
Definition: z_refine_grad.c:49
pastix_int_t c_bicgstab_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
pastix_int_t c_grad_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
Definition: c_refine_grad.c:49
pastix_int_t c_pivot_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
pastix_int_t z_bicgstab_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
pastix_int_t z_pivot_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
pastix_int_t z_gmres_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
pastix_int_t s_bicgstab_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
pastix_int_t(* refine_fct_t)(pastix_data_t *, pastix_rhs_t, pastix_rhs_t)
Select the refinement function to call depending on the matrix type and the precision.
pastix_int_t d_bicgstab_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
pastix_int_t d_grad_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
Definition: d_refine_grad.c:49
int pastix_subtask_refine(pastix_data_t *pastix_data, pastix_rhs_t Bp, pastix_rhs_t Xp)
Perform the iterative refinement without apply the permutations.
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_solve(pastix_data_t *pastix_data, pastix_rhs_t b)
Solve the given problem without applying the permutation.
int inter_node_procnum
Definition: pastixdata.h:83
pastix_int_t * iparm
Definition: pastixdata.h:69
pastix_int_t ld
Definition: pastixdata.h:155
double * dparm
Definition: pastixdata.h:70
pastix_coeftype_t flttype
Definition: pastixdata.h:152
pastix_int_t schur_n
Definition: pastixdata.h:93
PASTIX_Comm inter_node_comm
Definition: pastixdata.h:77
pastix_int_t m
Definition: pastixdata.h:153
pastix_bcsc_t * bcsc
Definition: pastixdata.h:101
pastix_int_t steps
Definition: pastixdata.h:72
pastix_int_t n
Definition: pastixdata.h:154
int pastix_task_refine(pastix_data_t *pastix_data, pastix_int_t m, pastix_int_t nrhs, void *b, pastix_int_t ldb, void *x, pastix_int_t ldx)
Perform iterative refinement.
int pastix_task_solve_and_refine(pastix_data_t *pastix_data, pastix_int_t m, pastix_int_t nrhs, void *b, pastix_int_t ldb, void *x, pastix_int_t ldx)
Performs solve and iterative refinement without unnecessary permutations.
Main PaStiX data structure.
Definition: pastixdata.h:67
Main PaStiX RHS structure.
Definition: pastixdata.h:150