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-2016 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 
38 /*----------------------------------------------------------------------------*/
39 
41 
42 /*============================================================================
43  * Macro definitions
44  *============================================================================*/
45 
46 /*============================================================================
47  * Type definitions
48  *============================================================================*/
49 
50 /*----------------------------------------------------------------------------
51  * Solver types
52  *----------------------------------------------------------------------------*/
53 
54 typedef enum {
55 
56  CS_SLES_PCG, /* Preconditionned conjugate gradient */
57  CS_SLES_JACOBI, /* Jacobi */
58  CS_SLES_BICGSTAB, /* Bi-conjugate gradient stabilized */
59  CS_SLES_BICGSTAB2, /* Bi-conjugate gradient stabilized - 2*/
60  CS_SLES_GMRES, /* Generalized minimal residual */
61  CS_SLES_B_GAUSS_SEIDEL, /* Block Gauss-Seidel */
62  CS_SLES_N_IT_TYPES /* Number of resolution algorithms */
63 
65 
66 /* Iterative linear solver context (opaque) */
67 
68 typedef struct _cs_sles_it_t cs_sles_it_t;
69 
70 /*============================================================================
71  * Global variables
72  *============================================================================*/
73 
74 /* Short names for matrix types */
75 
76 extern const char *cs_sles_it_type_name[];
77 
78 /*=============================================================================
79  * Public function prototypes
80  *============================================================================*/
81 
82 /*----------------------------------------------------------------------------
83  * Define and associate an iterative sparse linear system solver
84  * for a given field or equation name.
85  *
86  * If this system did not previously exist, it is added to the list of
87  * "known" systems. Otherwise, its definition is replaced by the one
88  * defined here.
89  *
90  * This is a utility function: if finer control is needed, see
91  * cs_sles_define() and cs_sles_it_create().
92  *
93  * Note that this function returns a pointer directly to the iterative solver
94  * management structure. This may be used to set further options,
95  * for example using cs_sles_it_set_pcg_single_reduction(). If needed,
96  * cs_sles_find() may be used to obtain a pointer to the matching
97  * cs_sles_t container.
98  *
99  * parameters:
100  * f_id <-- associated field id, or < 0
101  * name <-- associated name if f_id < 0, or NULL
102  * solver_type <-- type of solver (PCG, Jacobi, ...)
103  * poly_degree <-- preconditioning polynomial degree
104  * (0: diagonal; -1: non-preconditioned)
105  * n_max_iter <-- maximum number of iterations
106  *
107  * returns:
108  * pointer to newly created iterative solver info object.
109  *----------------------------------------------------------------------------*/
110 
111 cs_sles_it_t *
112 cs_sles_it_define(int f_id,
113  const char *name,
114  cs_sles_it_type_t solver_type,
115  int poly_degree,
116  int n_max_iter);
117 
118 /*----------------------------------------------------------------------------
119  * Create iterative sparse linear system solver info and context.
120  *
121  * parameters:
122  * solver_type <-- type of solver (PCG, Jacobi, ...)
123  * poly_degree <-- preconditioning polynomial degree
124  * (0: diagonal; -1: non-preconditioned)
125  * n_max_iter <-- maximum number of iterations
126  * update_stats <-- automatic solver statistics indicator
127  *
128  * returns:
129  * pointer to newly created solver info object.
130  *----------------------------------------------------------------------------*/
131 
132 cs_sles_it_t *
134  int poly_degree,
135  int n_max_iter,
136  bool update_stats);
137 
138 /*----------------------------------------------------------------------------
139  * Destroy iterative sparse linear system solver info and context.
140  *
141  * parameters:
142  * context <-> pointer to iterative sparse linear solver info
143  * (actual type: cs_sles_it_t **)
144  *----------------------------------------------------------------------------*/
145 
146 void
147 cs_sles_it_destroy(void **context);
148 
149 /*----------------------------------------------------------------------------
150  * Create iterative sparse linear system solver info and context
151  * based on existing info and context.
152  *
153  * parameters:
154  * context <-- pointer to reference info and context
155  * (actual type: cs_sles_it_t *)
156  *
157  * returns:
158  * pointer to newly created solver info object
159  * (actual type: cs_sles_it_t *)
160  *----------------------------------------------------------------------------*/
161 
162 void *
163 cs_sles_it_copy(const void *context);
164 
165 /*----------------------------------------------------------------------------
166  * Setup iterative sparse linear equation solver.
167  *
168  * parameters:
169  * context <-> pointer to iterative sparse linear solver info
170  * (actual type: cs_sles_it_t *)
171  * name <-- pointer to system name
172  * a <-- associated matrix
173  * verbosity <-- verbosity level
174  *----------------------------------------------------------------------------*/
175 
176 void
177 cs_sles_it_setup(void *context,
178  const char *name,
179  const cs_matrix_t *a,
180  int verbosity);
181 
182 /*----------------------------------------------------------------------------
183  * Call iterative sparse linear equation solver.
184  *
185  * parameters:
186  * context <-> pointer to iterative sparse linear solver info
187  * (actual type: cs_sles_it_t *)
188  * name <-- pointer to system name
189  * a <-- matrix
190  * verbosity <-- verbosity level
191  * rotation_mode <-- halo update option for rotational periodicity
192  * precision <-- solver precision
193  * r_norm <-- residue normalization
194  * n_iter --> number of iterations
195  * residue --> residue
196  * rhs <-- right hand side
197  * vx <-> system solution
198  * aux_size <-- number of elements in aux_vectors (in bytes)
199  * aux_vectors --- optional working area (internal allocation if NULL)
200  *
201  * returns:
202  * convergence state
203  *----------------------------------------------------------------------------*/
204 
206 cs_sles_it_solve(void *context,
207  const char *name,
208  const cs_matrix_t *a,
209  int verbosity,
210  cs_halo_rotation_t rotation_mode,
211  double precision,
212  double r_norm,
213  int *n_iter,
214  double *residue,
215  const cs_real_t *rhs,
216  cs_real_t *vx,
217  size_t aux_size,
218  void *aux_vectors);
219 
220 /*----------------------------------------------------------------------------
221  * Free iterative sparse linear equation solver setup context.
222  *
223  * This function frees resolution-related data, such as
224  * buffers and preconditioning but does not free the whole context,
225  * as info used for logging (especially performance data) is maintained.
226 
227  * parameters:
228  * context <-> pointer to iterative sparse linear solver info
229  * (actual type: cs_sles_it_t *)
230  *----------------------------------------------------------------------------*/
231 
232 void
233 cs_sles_it_free(void *context);
234 
235 /*----------------------------------------------------------------------------
236  * Log sparse linear equation solver info.
237  *
238  * parameters:
239  * context <-> pointer to iterative sparse linear solver info
240  * (actual type: cs_sles_it_t *)
241  * log_type <-- log type
242  *----------------------------------------------------------------------------*/
243 
244 void
245 cs_sles_it_log(const void *context,
246  cs_log_t log_type);
247 
248 /*----------------------------------------------------------------------------
249  * Return iterative solver type.
250  *
251  * parameters:
252  * context <-- pointer to iterative solver info and context
253  *
254  * returns:
255  * selected solver type
256  *----------------------------------------------------------------------------*/
257 
259 cs_sles_it_get_type(const cs_sles_it_t *context);
260 
261 /*----------------------------------------------------------------------------
262  * Return the initial residue for the previous solve operation with a solver.
263  *
264  * This is useful for convergence tests when this solver is used as
265  * a preconditioning smoother.
266  *
267  * This operation is only valid between calls to cs_sles_it_setup()
268  * (or cs_sles_it_solve()) and cs_sles_it_free().
269  * It returns -1 otherwise.
270  *
271  * parameters:
272  * context <-- pointer to iterative solver info and context
273  *
274  * returns:
275  * initial residue from last call to \ref cs_sles_solve with this solver
276  *----------------------------------------------------------------------------*/
277 
278 double
280 
281 /*----------------------------------------------------------------------------
282  * Associate a similar info and context object with which some setup
283  * data may be shared.
284  *
285  * This is especially useful for sharing preconditioning data between
286  * similar solver contexts (for example ascending and descending multigrid
287  * smoothers based on the same matrix).
288  *
289  * For preconditioning data to be effectively shared, cs_sles_it_setup()
290  * (or cs_sles_it_solve()) must be called on "shareable" before being
291  * called on "context" (without cs_sles_it_free() being called in between,
292  * of course).
293  *
294  * It is the caller's responsibility to ensure the context is not used
295  * for a cs_sles_it_setup() or cs_sles_it_solve() operation after the
296  * shareable object has been destroyed (normally by cs_sles_it_destroy()).
297  *
298  * parameters:
299  * context <-> pointer to iterative sparse linear system solver info
300  * shareable <-- pointer to iterative solver info and context
301  * whose context may be shared
302  *----------------------------------------------------------------------------*/
303 
304 void
306  const cs_sles_it_t *shareable);
307 
308 #if defined(HAVE_MPI)
309 
310 /*----------------------------------------------------------------------------
311  * Set MPI communicator for dot products.
312  *
313  * parameters:
314  * context <-> pointer to iterative sparse linear system solver info
315  * comm <-- MPI communicator
316  *----------------------------------------------------------------------------*/
317 
318 void
319 cs_sles_it_set_mpi_reduce_comm(cs_sles_it_t *context,
320  MPI_Comm comm);
321 
322 #endif /* defined(HAVE_MPI) */
323 
324 /*----------------------------------------------------------------------------
325  * Assign ordering to iterative solver.
326  *
327  * The solver context takes ownership of the order array (i.e. it will
328  * handle its later deallocation).
329  *
330  * This is useful only for Block Gauss-Seidel.
331  *
332  * parameters:
333  * context <-> pointer to iterative solver info and context
334  * order <-> pointer to ordering array
335  *----------------------------------------------------------------------------*/
336 
337 void
339  cs_lnum_t **order);
340 
341 /*----------------------------------------------------------------------------
342  * Query mean number of rows under which Conjugate Gradient algorithm
343  * uses the single-reduction variant.
344  *
345  * The single-reduction variant requires only one parallel sum per
346  * iteration (instead of 2), at the cost of additional vector operations,
347  * so it tends to be more expensive when the number of matrix rows per
348  * MPI rank is high, then becomes cheaper when the MPI latency cost becomes
349  * more significant.
350  *
351  * This option is ignored for non-parallel runs, so 0 is returned.
352  *
353  * return:
354  * mean number of rows per active rank under which the
355  * single-reduction variant will be used
356  *----------------------------------------------------------------------------*/
357 
358 cs_lnum_t
360 
361 /*----------------------------------------------------------------------------
362  * Set mean number of rows under which Conjugate Gradient algorithm
363  * should use the single-reduction variant.
364  *
365  * The single-reduction variant requires only one parallel sum per
366  * iteration (instead of 2), at the cost of additional vector operations,
367  * so it tends to be more expensive when the number of matrix rows per
368  * MPI rank is high, then becomes cheaper when the MPI latency cost becomes
369  * more significant.
370  *
371  * This option is ignored for non-parallel runs.
372  *
373  * parameters:
374  * threshold <-- mean number of rows per active rank under which the
375  * single-reduction variant will be used
376  *----------------------------------------------------------------------------*/
377 
378 void
380 
381 /*----------------------------------------------------------------------------
382  * Log the current global settings relative to parallelism.
383  *----------------------------------------------------------------------------*/
384 
385 void
387 
388 /*----------------------------------------------------------------------------
389  * Error handler for iterative sparse linear equation solver.
390  *
391  * In case of divergence or breakdown, this error handler outputs
392  * postprocessing data to assist debugging, then aborts the run.
393  * It does nothing in case the maximum iteration count is reached.
394  *
395  * parameters:
396  * context <-> pointer to iterative sparse linear system solver info
397  * (actual type: cs_sles_it_t *)
398  * state <-- convergence state
399  * name <-- pointer to name of linear system
400  * a <-- matrix
401  * rotation_mode <-- halo update option for rotational periodicity
402  * rhs <-- right hand side
403  * vx <-> system solution
404  */
405 /*----------------------------------------------------------------------------*/
406 
407 void
408 cs_sles_it_error_post_and_abort(void *context,
410  const char *name,
411  const cs_matrix_t *a,
412  cs_halo_rotation_t rotation_mode,
413  const cs_real_t *rhs,
414  cs_real_t *vx);
415 
416 /*----------------------------------------------------------------------------*/
417 
419 
420 #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:3343
cs_halo_rotation_t
Definition: cs_halo.h:59
cs_sles_it_type_t
Iterative solver types.
Definition: cs_sles_it.h:54
Definition: cs_sles_it.h:60
void cs_sles_it_error_post_and_abort(void *context, cs_sles_convergence_state_t state, const char *name, 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:3907
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:3700
#define BEGIN_C_DECLS
Definition: cs_defs.h:419
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:3833
struct _cs_sles_it_t cs_sles_it_t
Definition: cs_sles_it.h:68
struct _cs_matrix_t cs_matrix_t
Definition: cs_matrix.h:86
void cs_sles_it_destroy(void **context)
Destroy iterative sparse linear system solver info and context.
Definition: cs_sles_it.c:3220
void cs_sles_it_free(void *context)
Free iterative sparse linear equation solver setup context.
Definition: cs_sles_it.c:3645
void cs_sles_it_log(const void *context, cs_log_t log_type)
Log sparse linear equation solver info.
Definition: cs_sles_it.c:3277
cs_sles_convergence_state_t
Convergence status indicator.
Definition: cs_sles.h:55
Definition: cs_sles_it.h:62
Definition: cs_sles_it.h:59
Definition: cs_sles_it.h:56
double precision, save a
Definition: cs_fuel_incl.f90:146
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:3384
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:3675
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:3875
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:3861
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:292
Definition: cs_sles_it.h:58
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:3734
#define END_C_DECLS
Definition: cs_defs.h:420
double cs_real_t
Definition: cs_defs.h:296
Definition: cs_sles_it.h:61
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:3792
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:3169
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:3124
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:3248
Definition: cs_sles_it.h:57