programmer's documentation
cs_sles_it.h
Go to the documentation of this file.
1 #ifndef __CS_SLES_IT_H__
2 #define __CS_SLES_IT_H__
3 
4 /*============================================================================
5  * Sparse Linear Equation Solvers
6  *============================================================================*/
7 
8 /*
9  This file is part of Code_Saturne, a general-purpose CFD tool.
10 
11  Copyright (C) 1998-2018 EDF S.A.
12 
13  This program is free software; you can redistribute it and/or modify it under
14  the terms of the GNU General Public License as published by the Free Software
15  Foundation; either version 2 of the License, or (at your option) any later
16  version.
17 
18  This program is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
20  FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
21  details.
22 
23  You should have received a copy of the GNU General Public License along with
24  this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
25  Street, Fifth Floor, Boston, MA 02110-1301, USA.
26 */
27 
28 /*----------------------------------------------------------------------------*/
29 
30 /*----------------------------------------------------------------------------
31  * Local headers
32  *----------------------------------------------------------------------------*/
33 
34 #include "cs_base.h"
35 #include "cs_halo_perio.h"
36 #include "cs_matrix.h"
37 #include "cs_time_plot.h"
38 #include "cs_sles.h"
39 #include "cs_sles_pc.h"
40 
41 /*----------------------------------------------------------------------------*/
42 
44 
45 /*============================================================================
46  * Macro definitions
47  *============================================================================*/
48 
49 /*============================================================================
50  * Type definitions
51  *============================================================================*/
52 
53 /*----------------------------------------------------------------------------
54  * Solver types
55  *----------------------------------------------------------------------------*/
56 
57 typedef enum {
58 
59  CS_SLES_PCG, /* Preconditionned conjugate gradient */
60  CS_SLES_JACOBI, /* Jacobi */
61  CS_SLES_BICGSTAB, /* Bi-conjugate gradient stabilized */
62  CS_SLES_BICGSTAB2, /* Bi-conjugate gradient stabilized - 2*/
63  CS_SLES_GMRES, /* Generalized minimal residual */
64  CS_SLES_P_GAUSS_SEIDEL, /* Process-local Gauss-Seidel */
65  CS_SLES_PCR3, /* 3-layer conjugate residual */
66  CS_SLES_IPCG, /* Inexact preconditionned conjugate gradient */
67  CS_SLES_N_IT_TYPES /* Number of resolution algorithms */
68 
70 
71 /* Iterative linear solver context (opaque) */
72 
73 typedef struct _cs_sles_it_t cs_sles_it_t;
74 
75 /*============================================================================
76  * Global variables
77  *============================================================================*/
78 
79 /* Short names for matrix types */
80 
81 extern const char *cs_sles_it_type_name[];
82 
83 /*=============================================================================
84  * Public function prototypes
85  *============================================================================*/
86 
87 /*----------------------------------------------------------------------------
88  * Define and associate an iterative sparse linear system solver
89  * for a given field or equation name.
90  *
91  * If this system did not previously exist, it is added to the list of
92  * "known" systems. Otherwise, its definition is replaced by the one
93  * defined here.
94  *
95  * This is a utility function: if finer control is needed, see
96  * cs_sles_define() and cs_sles_it_create().
97  *
98  * Note that this function returns a pointer directly to the iterative solver
99  * management structure. This may be used to set further options,
100  * for example using cs_sles_set_plot_options(). If needed, cs_sles_find()
101  * may be used to obtain a pointer to the matching cs_sles_t container.
102  *
103  * parameters:
104  * f_id <-- associated field id, or < 0
105  * name <-- associated name if f_id < 0, or NULL
106  * solver_type <-- type of solver (PCG, Jacobi, ...)
107  * poly_degree <-- preconditioning polynomial degree
108  * (0: diagonal; -1: non-preconditioned)
109  * n_max_iter <-- maximum number of iterations
110  *
111  * returns:
112  * pointer to newly created iterative solver info object.
113  *----------------------------------------------------------------------------*/
114 
115 cs_sles_it_t *
116 cs_sles_it_define(int f_id,
117  const char *name,
118  cs_sles_it_type_t solver_type,
119  int poly_degree,
120  int n_max_iter);
121 
122 /*----------------------------------------------------------------------------
123  * Create iterative sparse linear system solver info and context.
124  *
125  * parameters:
126  * solver_type <-- type of solver (PCG, Jacobi, ...)
127  * poly_degree <-- preconditioning polynomial degree
128  * (0: diagonal; -1: non-preconditioned)
129  * n_max_iter <-- maximum number of iterations
130  * update_stats <-- automatic solver statistics indicator
131  *
132  * returns:
133  * pointer to newly created solver info object.
134  *----------------------------------------------------------------------------*/
135 
136 cs_sles_it_t *
138  int poly_degree,
139  int n_max_iter,
140  bool update_stats);
141 
142 /*----------------------------------------------------------------------------
143  * Destroy iterative sparse linear system solver info and context.
144  *
145  * parameters:
146  * context <-> pointer to iterative sparse linear solver info
147  * (actual type: cs_sles_it_t **)
148  *----------------------------------------------------------------------------*/
149 
150 void
151 cs_sles_it_destroy(void **context);
152 
153 /*----------------------------------------------------------------------------
154  * Create iterative sparse linear system solver info and context
155  * based on existing info and context.
156  *
157  * parameters:
158  * context <-- pointer to reference info and context
159  * (actual type: cs_sles_it_t *)
160  *
161  * returns:
162  * pointer to newly created solver info object
163  * (actual type: cs_sles_it_t *)
164  *----------------------------------------------------------------------------*/
165 
166 void *
167 cs_sles_it_copy(const void *context);
168 
169 /*----------------------------------------------------------------------------
170  * Setup iterative sparse linear equation solver.
171  *
172  * parameters:
173  * context <-> pointer to iterative sparse linear solver info
174  * (actual type: cs_sles_it_t *)
175  * name <-- pointer to system name
176  * a <-- associated matrix
177  * verbosity <-- verbosity level
178  *----------------------------------------------------------------------------*/
179 
180 void
181 cs_sles_it_setup(void *context,
182  const char *name,
183  const cs_matrix_t *a,
184  int verbosity);
185 
186 /*----------------------------------------------------------------------------
187  * Call iterative sparse linear equation solver.
188  *
189  * parameters:
190  * context <-> pointer to iterative sparse linear solver info
191  * (actual type: cs_sles_it_t *)
192  * name <-- pointer to system name
193  * a <-- matrix
194  * verbosity <-- verbosity level
195  * rotation_mode <-- halo update option for rotational periodicity
196  * precision <-- solver precision
197  * r_norm <-- residue normalization
198  * n_iter --> number of iterations
199  * residue --> residue
200  * rhs <-- right hand side
201  * vx <-> system solution
202  * aux_size <-- number of elements in aux_vectors (in bytes)
203  * aux_vectors --- optional working area (internal allocation if NULL)
204  *
205  * returns:
206  * convergence state
207  *----------------------------------------------------------------------------*/
208 
210 cs_sles_it_solve(void *context,
211  const char *name,
212  const cs_matrix_t *a,
213  int verbosity,
214  cs_halo_rotation_t rotation_mode,
215  double precision,
216  double r_norm,
217  int *n_iter,
218  double *residue,
219  const cs_real_t *rhs,
220  cs_real_t *vx,
221  size_t aux_size,
222  void *aux_vectors);
223 
224 /*----------------------------------------------------------------------------
225  * Free iterative sparse linear equation solver setup context.
226  *
227  * This function frees resolution-related data, such as
228  * buffers and preconditioning but does not free the whole context,
229  * as info used for logging (especially performance data) is maintained.
230 
231  * parameters:
232  * context <-> pointer to iterative sparse linear solver info
233  * (actual type: cs_sles_it_t *)
234  *----------------------------------------------------------------------------*/
235 
236 void
237 cs_sles_it_free(void *context);
238 
239 /*----------------------------------------------------------------------------
240  * Log sparse linear equation solver info.
241  *
242  * parameters:
243  * context <-> pointer to iterative sparse linear solver info
244  * (actual type: cs_sles_it_t *)
245  * log_type <-- log type
246  *----------------------------------------------------------------------------*/
247 
248 void
249 cs_sles_it_log(const void *context,
250  cs_log_t log_type);
251 
252 /*----------------------------------------------------------------------------
253  * Return iterative solver type.
254  *
255  * parameters:
256  * context <-- pointer to iterative solver info and context
257  *
258  * returns:
259  * selected solver type
260  *----------------------------------------------------------------------------*/
261 
263 cs_sles_it_get_type(const cs_sles_it_t *context);
264 
265 /*----------------------------------------------------------------------------
266  * Return the initial residue for the previous solve operation with a solver.
267  *
268  * This is useful for convergence tests when this solver is used as
269  * a preconditioning smoother.
270  *
271  * This operation is only valid between calls to cs_sles_it_setup()
272  * (or cs_sles_it_solve()) and cs_sles_it_free().
273  * It returns -1 otherwise.
274  *
275  * parameters:
276  * context <-- pointer to iterative solver info and context
277  *
278  * returns:
279  * initial residue from last call to \ref cs_sles_solve with this solver
280  *----------------------------------------------------------------------------*/
281 
282 double
284 
285 /*----------------------------------------------------------------------------
286  * Return a preconditioner context for an iterative sparse linear
287  * equation solver.
288  *
289  * This allows modifying parameters of a non default (Jacobi or polynomial)
290  * preconditioner.
291  *
292  * parameters:
293  * context <-- pointer to iterative solver info and context
294  *
295  * returns:
296  * pointer to preconditoner context
297  *----------------------------------------------------------------------------*/
298 
299 cs_sles_pc_t *
301 
302 /*----------------------------------------------------------------------------
303  * Assign a preconditioner to an iterative sparse linear equation
304  * solver, transfering its ownership to to solver context.
305  *
306  * This allows assigning a non default (Jacobi or polynomial) preconditioner.
307  *
308  * The input pointer is set to NULL to make it clear the caller does not
309  * own the preconditioner anymore, though the context can be accessed using
310  * cs_sles_it_get_cp().
311  *
312  * parameters:
313  * context <-> pointer to iterative solver info and context
314  * pc <-> pointer to preconditoner context
315  *----------------------------------------------------------------------------*/
316 
317 void
319  cs_sles_pc_t **pc);
320 
321 /*----------------------------------------------------------------------------
322  * Copy options from one iterative sparse linear system solver info
323  * and context to another.
324  *
325  * Optional plotting contexts are shared between the source and destination
326  * contexts.
327  *
328  * Preconditioner settings are to be handled separately.
329  *
330  * parameters:
331  * src <-- pointer to source info and context
332  * dest <-> pointer to destination info and context
333  *----------------------------------------------------------------------------*/
334 
335 void
337  cs_sles_it_t *dest);
338 
339 /*----------------------------------------------------------------------------
340  * Associate a similar info and context object with which some setup
341  * data may be shared.
342  *
343  * This is especially useful for sharing preconditioning data between
344  * similar solver contexts (for example ascending and descending multigrid
345  * smoothers based on the same matrix).
346  *
347  * For preconditioning data to be effectively shared, cs_sles_it_setup()
348  * (or cs_sles_it_solve()) must be called on "shareable" before being
349  * called on "context" (without cs_sles_it_free() being called in between,
350  * of course).
351  *
352  * It is the caller's responsibility to ensure the context is not used
353  * for a cs_sles_it_setup() or cs_sles_it_solve() operation after the
354  * shareable object has been destroyed (normally by cs_sles_it_destroy()).
355  *
356  * parameters:
357  * context <-> pointer to iterative sparse linear system solver info
358  * shareable <-- pointer to iterative solver info and context
359  * whose context may be shared
360  *----------------------------------------------------------------------------*/
361 
362 void
364  const cs_sles_it_t *shareable);
365 
366 #if defined(HAVE_MPI)
367 
368 /*----------------------------------------------------------------------------
369  * Set MPI communicator for dot products.
370  *
371  * parameters:
372  * context <-> pointer to iterative sparse linear system solver info
373  * comm <-- MPI communicator
374  *----------------------------------------------------------------------------*/
375 
376 void
377 cs_sles_it_set_mpi_reduce_comm(cs_sles_it_t *context,
378  MPI_Comm comm);
379 
380 #endif /* defined(HAVE_MPI) */
381 
382 /*----------------------------------------------------------------------------
383  * Assign ordering to iterative solver.
384  *
385  * The solver context takes ownership of the order array (i.e. it will
386  * handle its later deallocation).
387  *
388  * This is useful only for Block Gauss-Seidel.
389  *
390  * parameters:
391  * context <-> pointer to iterative solver info and context
392  * order <-> pointer to ordering array
393  *----------------------------------------------------------------------------*/
394 
395 void
397  cs_lnum_t **order);
398 
399 /*----------------------------------------------------------------------------
400  * Assign symmetric option to iterative solver.
401  *
402  * This is useful only for Process-local Gauss-Seidel.
403  *
404  * parameters:
405  * context <-> pointer to iterative solver info and context
406  * symmetric <-- symmetric if true
407  *----------------------------------------------------------------------------*/
408 
409 void
411  bool symmetric);
412 
413 /*----------------------------------------------------------------------------
414  * Query mean number of rows under which Conjugate Gradient algorithm
415  * uses the single-reduction variant.
416  *
417  * The single-reduction variant requires only one parallel sum per
418  * iteration (instead of 2), at the cost of additional vector operations,
419  * so it tends to be more expensive when the number of matrix rows per
420  * MPI rank is high, then becomes cheaper when the MPI latency cost becomes
421  * more significant.
422  *
423  * This option is ignored for non-parallel runs, so 0 is returned.
424  *
425  * return:
426  * mean number of rows per active rank under which the
427  * single-reduction variant will be used
428  *----------------------------------------------------------------------------*/
429 
430 cs_lnum_t
432 
433 /*----------------------------------------------------------------------------
434  * Set mean number of rows under which Conjugate Gradient algorithm
435  * should use the single-reduction variant.
436  *
437  * The single-reduction variant requires only one parallel sum per
438  * iteration (instead of 2), at the cost of additional vector operations,
439  * so it tends to be more expensive when the number of matrix rows per
440  * MPI rank is high, then becomes cheaper when the MPI latency cost becomes
441  * more significant.
442  *
443  * This option is ignored for non-parallel runs.
444  *
445  * parameters:
446  * threshold <-- mean number of rows per active rank under which the
447  * single-reduction variant will be used
448  *----------------------------------------------------------------------------*/
449 
450 void
452 
453 /*----------------------------------------------------------------------------
454  * Log the current global settings relative to parallelism.
455  *----------------------------------------------------------------------------*/
456 
457 void
459 
460 /*----------------------------------------------------------------------------
461  * Error handler for iterative sparse linear equation solver.
462  *
463  * In case of divergence or breakdown, this error handler outputs
464  * postprocessing data to assist debugging, then aborts the run.
465  * It does nothing in case the maximum iteration count is reached.
466  *
467  * parameters:
468  * sles <-> pointer to solver object
469  * state <-- convergence state
470  * a <-- matrix
471  * rotation_mode <-- halo update option for rotational periodicity
472  * rhs <-- right hand side
473  * vx <-> system solution
474  *
475  * returns:
476  * false (do not attempt new solve)
477  *----------------------------------------------------------------------------*/
478 
479 bool
482  const cs_matrix_t *a,
483  cs_halo_rotation_t rotation_mode,
484  const cs_real_t *rhs,
485  cs_real_t *vx);
486 
487 /*----------------------------------------------------------------------------
488  * Set plotting options for an iterative sparse linear equation solver.
489  *
490  * parameters:
491  * context <-> pointer to iterative solver info and context
492  * base_name <-- base plot name to activate, NULL otherwise
493  * use_iteration <-- if true, use iteration as time stamp
494  * otherwise, use wall clock time
495  *----------------------------------------------------------------------------*/
496 
497 void
499  const char *base_name,
500  bool use_iteration);
501 
502 /*----------------------------------------------------------------------------
503  * Assign existing time plot to iterative sparse linear equation solver.
504  *
505  * This is useful mainly when a time plot has a longer lifecycle than
506  * the linear solver context, such as inside a multigrid solver.
507  *
508  * parameters:
509  * context <-> pointer to iterative solver info and context
510  * time_plot <-- pointer to time plot structure
511  * time_stamp <-- associated time stamp
512  *----------------------------------------------------------------------------*/
513 
514 void
516  cs_time_plot_t *time_plot,
517  int time_stamp);
518 
519 /*----------------------------------------------------------------------------*/
520 
522 
523 #endif /* __CS_SLES_IT_H__ */
void cs_sles_it_setup(void *context, const char *name, const cs_matrix_t *a, int verbosity)
Setup iterative sparse linear equation solver.
Definition: cs_sles_it.c:4396
void cs_sles_it_transfer_parameters(const cs_sles_it_t *src, cs_sles_it_t *dest)
Copy options from one iterative sparse linear system solver info and context to another.
Definition: cs_sles_it.c:4883
cs_halo_rotation_t
Definition: cs_halo.h:60
cs_sles_it_type_t
Iterative solver types.
Definition: cs_sles_it.h:57
Definition: cs_sles_it.h:63
double cs_sles_it_get_last_initial_residue(const cs_sles_it_t *context)
Return the initial residue for the previous solve operation with a solver.
Definition: cs_sles_it.c:4798
#define BEGIN_C_DECLS
Definition: cs_defs.h:451
Definition: cs_sles_it.h:66
void cs_sles_it_transfer_pc(cs_sles_it_t *context, cs_sles_pc_t **pc)
Assign a preconditioner to an iterative sparse linear equation solver, transfering its ownership to t...
Definition: cs_sles_it.c:4851
struct _cs_sles_pc_t cs_sles_pc_t
Definition: cs_sles_pc.h:66
void cs_sles_it_set_symmetric(cs_sles_it_t *context, bool symmetric)
Assign symmetric option to iterative solver.
Definition: cs_sles_it.c:5025
cs_lnum_t cs_sles_it_get_pcg_single_reduction(void)
Query mean number of rows under which Conjugate Gradient algorithm uses the single-reduction variant...
Definition: cs_sles_it.c:5051
struct _cs_sles_it_t cs_sles_it_t
Definition: cs_sles_it.h:73
double cs_real_t
Floating-point value.
Definition: cs_defs.h:297
struct _cs_matrix_t cs_matrix_t
Definition: cs_matrix.h:90
void cs_sles_it_destroy(void **context)
Destroy iterative sparse linear system solver info and context.
Definition: cs_sles_it.c:4254
struct _cs_time_plot_t cs_time_plot_t
Definition: cs_time_plot.h:48
void cs_sles_it_free(void *context)
Free iterative sparse linear equation solver setup context.
Definition: cs_sles_it.c:4740
void cs_sles_it_log(const void *context, cs_log_t log_type)
Log sparse linear equation solver info.
Definition: cs_sles_it.c:4325
cs_sles_convergence_state_t
Convergence status indicator.
Definition: cs_sles.h:55
struct _cs_sles_t cs_sles_t
Definition: cs_sles.h:67
Definition: cs_sles_it.h:67
void cs_sles_it_assign_plot(cs_sles_it_t *context, cs_time_plot_t *time_plot, int time_stamp)
Assign existing time plot to iterative sparse linear equation solver.
Definition: cs_sles_it.c:5217
Definition: cs_sles_it.h:62
Definition: cs_sles_it.h:59
double precision, save a
Definition: cs_fuel_incl.f90:146
Definition: cs_sles_it.h:65
cs_sles_pc_t * cs_sles_it_get_pc(cs_sles_it_t *context)
Return a preconditioner context for an iterative sparse linear equation solver.
Definition: cs_sles_it.c:4822
const char * cs_sles_it_type_name[]
cs_sles_convergence_state_t cs_sles_it_solve(void *context, const char *name, const cs_matrix_t *a, int verbosity, cs_halo_rotation_t rotation_mode, double precision, double r_norm, int *n_iter, double *residue, const cs_real_t *rhs, cs_real_t *vx, size_t aux_size, void *aux_vectors)
Call iterative sparse linear equation solver.
Definition: cs_sles_it.c:4443
bool cs_sles_it_error_post_and_abort(cs_sles_t *sles, cs_sles_convergence_state_t state, const cs_matrix_t *a, cs_halo_rotation_t rotation_mode, const cs_real_t *rhs, cs_real_t *vx)
Error handler for iterative sparse linear equation solver.
Definition: cs_sles_it.c:5125
cs_sles_it_type_t cs_sles_it_get_type(const cs_sles_it_t *context)
Return iterative solver type.
Definition: cs_sles_it.c:4773
cs_log_t
Definition: cs_log.h:47
void cs_sles_it_log_parallel_options(void)
Log the current global settings relative to parallelism.
Definition: cs_sles_it.c:5093
void cs_sles_it_set_plot_options(cs_sles_it_t *context, const char *base_name, bool use_iteration)
Set plotting options for an iterative sparse linear equation solver.
Definition: cs_sles_it.c:5173
void cs_sles_it_set_pcg_single_reduction(cs_lnum_t threshold)
Set mean number of rows under which Conjugate Gradient algorithm should use the single-reduction vari...
Definition: cs_sles_it.c:5079
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:293
Definition: cs_sles_it.h:61
void cs_sles_it_set_shareable(cs_sles_it_t *context, const cs_sles_it_t *shareable)
Associate a similar info and context object with which some setup data may be shared.
Definition: cs_sles_it.c:4928
#define END_C_DECLS
Definition: cs_defs.h:452
void cs_sles_it_assign_order(cs_sles_it_t *context, cs_lnum_t **order)
Assign ordering to iterative solver.
Definition: cs_sles_it.c:4991
cs_sles_it_t * cs_sles_it_create(cs_sles_it_type_t solver_type, int poly_degree, int n_max_iter, bool update_stats)
Create iterative sparse linear system solver info and context.
Definition: cs_sles_it.c:4179
cs_sles_it_t * cs_sles_it_define(int f_id, const char *name, cs_sles_it_type_t solver_type, int poly_degree, int n_max_iter)
Define and associate an iterative sparse linear system solver for a given field or equation name...
Definition: cs_sles_it.c:4134
Definition: cs_sles_it.h:64
void * cs_sles_it_copy(const void *context)
Create iterative sparse linear system solver info and context based on existing info and context...
Definition: cs_sles_it.c:4287
Definition: cs_sles_it.h:60