My Project
programmer's documentation
Loading...
Searching...
No Matches
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-2019 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
57typedef enum {
58
59 CS_SLES_PCG, /* Preconditionned conjugate gradient */
60 CS_SLES_FCG, /* Preconditionned flexible conjugate gradient */
61 CS_SLES_IPCG, /* Preconditionned inexact conjugate gradient */
62 CS_SLES_JACOBI, /* Jacobi */
63 CS_SLES_BICGSTAB, /* Bi-conjugate gradient stabilized */
64 CS_SLES_BICGSTAB2, /* Bi-conjugate gradient stabilized - 2*/
65 CS_SLES_GMRES, /* Generalized minimal residual */
66 CS_SLES_P_GAUSS_SEIDEL, /* Process-local Gauss-Seidel */
67 CS_SLES_P_SYM_GAUSS_SEIDEL, /* Process-local symmetric Gauss-Seidel */
68 CS_SLES_TS_F_GAUSS_SEIDEL, /* Truncated forward Gauss-Seidel
69 smoothing */
70 CS_SLES_TS_B_GAUSS_SEIDEL, /* Truncated backward Gauss-Seidel
71 smoothing */
72 CS_SLES_PCR3, /* 3-layer conjugate residual */
73 CS_SLES_N_IT_TYPES /* Number of resolution algorithms */
74
76
77/* Iterative linear solver context (opaque) */
78
79typedef struct _cs_sles_it_t cs_sles_it_t;
80
81/*============================================================================
82 * Global variables
83 *============================================================================*/
84
85/* Names for iterative solver types */
86
87extern const char *cs_sles_it_type_name[];
88
89/*=============================================================================
90 * Public function prototypes
91 *============================================================================*/
92
93/*----------------------------------------------------------------------------
94 * Define and associate an iterative sparse linear system solver
95 * for a given field or equation name.
96 *
97 * If this system did not previously exist, it is added to the list of
98 * "known" systems. Otherwise, its definition is replaced by the one
99 * defined here.
100 *
101 * This is a utility function: if finer control is needed, see
102 * cs_sles_define() and cs_sles_it_create().
103 *
104 * Note that this function returns a pointer directly to the iterative solver
105 * management structure. This may be used to set further options,
106 * for example using cs_sles_set_plot_options(). If needed, cs_sles_find()
107 * may be used to obtain a pointer to the matching cs_sles_t container.
108 *
109 * parameters:
110 * f_id <-- associated field id, or < 0
111 * name <-- associated name if f_id < 0, or NULL
112 * solver_type <-- type of solver (PCG, Jacobi, ...)
113 * poly_degree <-- preconditioning polynomial degree
114 * (0: diagonal; -1: non-preconditioned)
115 * n_max_iter <-- maximum number of iterations
116 *
117 * returns:
118 * pointer to newly created iterative solver info object.
119 *----------------------------------------------------------------------------*/
120
122cs_sles_it_define(int f_id,
123 const char *name,
124 cs_sles_it_type_t solver_type,
125 int poly_degree,
126 int n_max_iter);
127
128/*----------------------------------------------------------------------------
129 * Create iterative sparse linear system solver info and context.
130 *
131 * parameters:
132 * solver_type <-- type of solver (PCG, Jacobi, ...)
133 * poly_degree <-- preconditioning polynomial degree
134 * (0: diagonal; -1: non-preconditioned)
135 * n_max_iter <-- maximum number of iterations
136 * update_stats <-- automatic solver statistics indicator
137 *
138 * returns:
139 * pointer to newly created solver info object.
140 *----------------------------------------------------------------------------*/
141
144 int poly_degree,
145 int n_max_iter,
146 bool update_stats);
147
148/*----------------------------------------------------------------------------
149 * Destroy iterative sparse linear system solver info and context.
150 *
151 * parameters:
152 * context <-> pointer to iterative sparse linear solver info
153 * (actual type: cs_sles_it_t **)
154 *----------------------------------------------------------------------------*/
155
156void
157cs_sles_it_destroy(void **context);
158
159/*----------------------------------------------------------------------------
160 * Create iterative sparse linear system solver info and context
161 * based on existing info and context.
162 *
163 * parameters:
164 * context <-- pointer to reference info and context
165 * (actual type: cs_sles_it_t *)
166 *
167 * returns:
168 * pointer to newly created solver info object
169 * (actual type: cs_sles_it_t *)
170 *----------------------------------------------------------------------------*/
171
172void *
173cs_sles_it_copy(const void *context);
174
175/*----------------------------------------------------------------------------
176 * Setup iterative sparse linear equation solver.
177 *
178 * parameters:
179 * context <-> pointer to iterative sparse linear solver info
180 * (actual type: cs_sles_it_t *)
181 * name <-- pointer to system name
182 * a <-- associated matrix
183 * verbosity <-- verbosity level
184 *----------------------------------------------------------------------------*/
185
186void
187cs_sles_it_setup(void *context,
188 const char *name,
189 const cs_matrix_t *a,
190 int verbosity);
191
192/*----------------------------------------------------------------------------
193 * Call iterative sparse linear equation solver.
194 *
195 * parameters:
196 * context <-> pointer to iterative sparse linear solver info
197 * (actual type: cs_sles_it_t *)
198 * name <-- pointer to system name
199 * a <-- matrix
200 * verbosity <-- verbosity level
201 * rotation_mode <-- halo update option for rotational periodicity
202 * precision <-- solver precision
203 * r_norm <-- residue normalization
204 * n_iter --> number of iterations
205 * residue --> residue
206 * rhs <-- right hand side
207 * vx <-> system solution
208 * aux_size <-- number of elements in aux_vectors (in bytes)
209 * aux_vectors --- optional working area (internal allocation if NULL)
210 *
211 * returns:
212 * convergence state
213 *----------------------------------------------------------------------------*/
214
216cs_sles_it_solve(void *context,
217 const char *name,
218 const cs_matrix_t *a,
219 int verbosity,
220 cs_halo_rotation_t rotation_mode,
221 double precision,
222 double r_norm,
223 int *n_iter,
224 double *residue,
225 const cs_real_t *rhs,
226 cs_real_t *vx,
227 size_t aux_size,
228 void *aux_vectors);
229
230/*----------------------------------------------------------------------------
231 * Free iterative sparse linear equation solver setup context.
232 *
233 * This function frees resolution-related data, such as
234 * buffers and preconditioning but does not free the whole context,
235 * as info used for logging (especially performance data) is maintained.
236
237 * parameters:
238 * context <-> pointer to iterative sparse linear solver info
239 * (actual type: cs_sles_it_t *)
240 *----------------------------------------------------------------------------*/
241
242void
243cs_sles_it_free(void *context);
244
245/*----------------------------------------------------------------------------
246 * Log sparse linear equation solver info.
247 *
248 * parameters:
249 * context <-> pointer to iterative sparse linear solver info
250 * (actual type: cs_sles_it_t *)
251 * log_type <-- log type
252 *----------------------------------------------------------------------------*/
253
254void
255cs_sles_it_log(const void *context,
256 cs_log_t log_type);
257
258/*----------------------------------------------------------------------------
259 * Return iterative solver type.
260 *
261 * parameters:
262 * context <-- pointer to iterative solver info and context
263 *
264 * returns:
265 * selected solver type
266 *----------------------------------------------------------------------------*/
267
269cs_sles_it_get_type(const cs_sles_it_t *context);
270
271/*----------------------------------------------------------------------------
272 * Return the initial residue for the previous solve operation with a solver.
273 *
274 * This is useful for convergence tests when this solver is used as
275 * a preconditioning smoother.
276 *
277 * This operation is only valid between calls to cs_sles_it_setup()
278 * (or cs_sles_it_solve()) and cs_sles_it_free().
279 * It returns -1 otherwise.
280 *
281 * parameters:
282 * context <-- pointer to iterative solver info and context
283 *
284 * returns:
285 * initial residue from last call to \ref cs_sles_solve with this solver
286 *----------------------------------------------------------------------------*/
287
288double
290
291/*----------------------------------------------------------------------------
292 * Return a preconditioner context for an iterative sparse linear
293 * equation solver.
294 *
295 * This allows modifying parameters of a non default (Jacobi or polynomial)
296 * preconditioner.
297 *
298 * parameters:
299 * context <-- pointer to iterative solver info and context
300 *
301 * returns:
302 * pointer to preconditoner context
303 *----------------------------------------------------------------------------*/
304
307
308/*----------------------------------------------------------------------------
309 * Assign a preconditioner to an iterative sparse linear equation
310 * solver, transfering its ownership to to solver context.
311 *
312 * This allows assigning a non default (Jacobi or polynomial) preconditioner.
313 *
314 * The input pointer is set to NULL to make it clear the caller does not
315 * own the preconditioner anymore, though the context can be accessed using
316 * cs_sles_it_get_cp().
317 *
318 * parameters:
319 * context <-> pointer to iterative solver info and context
320 * pc <-> pointer to preconditoner context
321 *----------------------------------------------------------------------------*/
322
323void
325 cs_sles_pc_t **pc);
326
327/*----------------------------------------------------------------------------
328 * Copy options from one iterative sparse linear system solver info
329 * and context to another.
330 *
331 * Optional plotting contexts are shared between the source and destination
332 * contexts.
333 *
334 * Preconditioner settings are to be handled separately.
335 *
336 * parameters:
337 * src <-- pointer to source info and context
338 * dest <-> pointer to destination info and context
339 *----------------------------------------------------------------------------*/
340
341void
343 cs_sles_it_t *dest);
344
345/*----------------------------------------------------------------------------
346 * Associate a similar info and context object with which some setup
347 * data may be shared.
348 *
349 * This is especially useful for sharing preconditioning data between
350 * similar solver contexts (for example ascending and descending multigrid
351 * smoothers based on the same matrix).
352 *
353 * For preconditioning data to be effectively shared, cs_sles_it_setup()
354 * (or cs_sles_it_solve()) must be called on "shareable" before being
355 * called on "context" (without cs_sles_it_free() being called in between,
356 * of course).
357 *
358 * It is the caller's responsibility to ensure the context is not used
359 * for a cs_sles_it_setup() or cs_sles_it_solve() operation after the
360 * shareable object has been destroyed (normally by cs_sles_it_destroy()).
361 *
362 * parameters:
363 * context <-> pointer to iterative sparse linear system solver info
364 * shareable <-- pointer to iterative solver info and context
365 * whose context may be shared
366 *----------------------------------------------------------------------------*/
367
368void
370 const cs_sles_it_t *shareable);
371
372#if defined(HAVE_MPI)
373
374/*----------------------------------------------------------------------------*/
386/*----------------------------------------------------------------------------*/
387
388void
389cs_sles_it_set_mpi_reduce_comm(cs_sles_it_t *context,
390 MPI_Comm comm,
391 MPI_Comm caller_comm);
392
393#endif /* defined(HAVE_MPI) */
394
395/*----------------------------------------------------------------------------
396 * Assign ordering to iterative solver.
397 *
398 * The solver context takes ownership of the order array (i.e. it will
399 * handle its later deallocation).
400 *
401 * This is useful only for Block Gauss-Seidel.
402 *
403 * parameters:
404 * context <-> pointer to iterative solver info and context
405 * order <-> pointer to ordering array
406 *----------------------------------------------------------------------------*/
407
408void
410 cs_lnum_t **order);
411
412/*----------------------------------------------------------------------------*/
426/*----------------------------------------------------------------------------*/
427
428void
431
432/*----------------------------------------------------------------------------
433 * Query mean number of rows under which Conjugate Gradient algorithm
434 * uses the single-reduction variant.
435 *
436 * The single-reduction variant requires only one parallel sum per
437 * iteration (instead of 2), at the cost of additional vector operations,
438 * so it tends to be more expensive when the number of matrix rows per
439 * MPI rank is high, then becomes cheaper when the MPI latency cost becomes
440 * more significant.
441 *
442 * This option is ignored for non-parallel runs, so 0 is returned.
443 *
444 * return:
445 * mean number of rows per active rank under which the
446 * single-reduction variant will be used
447 *----------------------------------------------------------------------------*/
448
451
452/*----------------------------------------------------------------------------
453 * Set mean number of rows under which Conjugate Gradient algorithm
454 * should use the single-reduction variant.
455 *
456 * The single-reduction variant requires only one parallel sum per
457 * iteration (instead of 2), at the cost of additional vector operations,
458 * so it tends to be more expensive when the number of matrix rows per
459 * MPI rank is high, then becomes cheaper when the MPI latency cost becomes
460 * more significant.
461 *
462 * This option is ignored for non-parallel runs.
463 *
464 * parameters:
465 * threshold <-- mean number of rows per active rank under which the
466 * single-reduction variant will be used
467 *----------------------------------------------------------------------------*/
468
469void
471
472/*----------------------------------------------------------------------------
473 * Log the current global settings relative to parallelism.
474 *----------------------------------------------------------------------------*/
475
476void
478
479/*----------------------------------------------------------------------------
480 * Error handler for iterative sparse linear equation solver.
481 *
482 * In case of divergence or breakdown, this error handler outputs
483 * postprocessing data to assist debugging, then aborts the run.
484 * It does nothing in case the maximum iteration count is reached.
485 *
486 * parameters:
487 * sles <-> pointer to solver object
488 * state <-- convergence state
489 * a <-- matrix
490 * rotation_mode <-- halo update option for rotational periodicity
491 * rhs <-- right hand side
492 * vx <-> system solution
493 *
494 * returns:
495 * false (do not attempt new solve)
496 *----------------------------------------------------------------------------*/
497
498bool
501 const cs_matrix_t *a,
502 cs_halo_rotation_t rotation_mode,
503 const cs_real_t *rhs,
504 cs_real_t *vx);
505
506/*----------------------------------------------------------------------------
507 * Set plotting options for an iterative sparse linear equation solver.
508 *
509 * parameters:
510 * context <-> pointer to iterative solver info and context
511 * base_name <-- base plot name to activate, NULL otherwise
512 * use_iteration <-- if true, use iteration as time stamp
513 * otherwise, use wall clock time
514 *----------------------------------------------------------------------------*/
515
516void
518 const char *base_name,
519 bool use_iteration);
520
521/*----------------------------------------------------------------------------
522 * Assign existing time plot to iterative sparse linear equation solver.
523 *
524 * This is useful mainly when a time plot has a longer lifecycle than
525 * the linear solver context, such as inside a multigrid solver.
526 *
527 * parameters:
528 * context <-> pointer to iterative solver info and context
529 * time_plot <-- pointer to time plot structure
530 * time_stamp <-- associated time stamp
531 *----------------------------------------------------------------------------*/
532
533void
535 cs_time_plot_t *time_plot,
536 int time_stamp);
537
538/*----------------------------------------------------------------------------*/
539
541
542#endif /* __CS_SLES_IT_H__ */
#define BEGIN_C_DECLS
Definition cs_defs.h:467
double cs_real_t
Floating-point value.
Definition cs_defs.h:302
#define END_C_DECLS
Definition cs_defs.h:468
int cs_lnum_t
local mesh entity id
Definition cs_defs.h:298
cs_halo_rotation_t
Definition cs_halo.h:60
cs_log_t
Definition cs_log.h:48
struct _cs_matrix_t cs_matrix_t
Definition cs_matrix.h:90
cs_sles_convergence_state_t
Definition cs_sles.h:56
struct _cs_sles_t cs_sles_t
Definition cs_sles.h:68
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:5354
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:5407
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:5378
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:5439
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:5787
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:5129
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:4687
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:4967
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:5484
void cs_sles_it_log(const void *context, cs_log_t log_type)
Log sparse linear equation solver info.
Definition cs_sles_it.c:4863
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:4632
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:4825
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:5743
struct _cs_sles_it_t cs_sles_it_t
Definition cs_sles_it.h:79
void cs_sles_it_log_parallel_options(void)
Log the current global settings relative to parallelism.
Definition cs_sles_it.c:5663
void cs_sles_it_free(void *context)
Free iterative sparse linear equation solver setup context.
Definition cs_sles_it.c:5293
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:5557
const char * cs_sles_it_type_name[]
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:5621
void cs_sles_it_set_fallback_threshold(cs_sles_it_t *context, cs_sles_convergence_state_t threshold)
Define convergence level under which the fallback to another solver may be used if applicable.
Definition cs_sles_it.c:5596
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:5695
cs_sles_it_type_t
Definition cs_sles_it.h:57
@ CS_SLES_PCR3
Definition cs_sles_it.h:72
@ CS_SLES_FCG
Definition cs_sles_it.h:60
@ CS_SLES_P_SYM_GAUSS_SEIDEL
Definition cs_sles_it.h:67
@ CS_SLES_P_GAUSS_SEIDEL
Definition cs_sles_it.h:66
@ CS_SLES_TS_F_GAUSS_SEIDEL
Definition cs_sles_it.h:68
@ CS_SLES_TS_B_GAUSS_SEIDEL
Definition cs_sles_it.h:70
@ CS_SLES_BICGSTAB2
Definition cs_sles_it.h:64
@ CS_SLES_GMRES
Definition cs_sles_it.h:65
@ CS_SLES_N_IT_TYPES
Definition cs_sles_it.h:73
@ CS_SLES_PCG
Definition cs_sles_it.h:59
@ CS_SLES_IPCG
Definition cs_sles_it.h:61
@ CS_SLES_JACOBI
Definition cs_sles_it.h:62
@ CS_SLES_BICGSTAB
Definition cs_sles_it.h:63
void cs_sles_it_destroy(void **context)
Destroy iterative sparse linear system solver info and context.
Definition cs_sles_it.c:4787
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:5649
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:5329
struct _cs_sles_pc_t cs_sles_pc_t
Definition cs_sles_pc.h:66
struct _cs_time_plot_t cs_time_plot_t
Definition cs_time_plot.h:48