9.1
general 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-2025 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 "base/cs_base.h"
35#include "alge/cs_matrix.h"
36#include "alge/cs_sles.h"
37#include "alge/cs_sles_pc.h"
38
39/*----------------------------------------------------------------------------*/
40
42
43/*============================================================================
44 * Macro definitions
45 *============================================================================*/
46
47/*============================================================================
48 * Type definitions
49 *============================================================================*/
50
51/*----------------------------------------------------------------------------
52 * Solver types
53 *----------------------------------------------------------------------------*/
54
55typedef enum {
56
57 /* Iterative solvers */
58
78 /* Multigrid smoothers */
79
92
93/* Iterative linear solver context (opaque) */
94
95typedef struct _cs_sles_it_t cs_sles_it_t;
96
97/* Forward type declarations */
98
99typedef struct _cs_sles_it_convergence_t cs_sles_it_convergence_t;
100
101/*============================================================================
102 * Global variables
103 *============================================================================*/
104
105/* Names for iterative solver types */
106
107extern const char *cs_sles_it_type_name[];
108
109/*=============================================================================
110 * User function prototypes
111 *============================================================================*/
112
113/*----------------------------------------------------------------------------
114 * Solution of A.vx = Rhs using a user-defined iterative solver
115 *
116 * On entry, vx is considered initialized.
117 *
118 * parameters:
119 * c <-- pointer to solver context info
120 * a <-- matrix
121 * diag_block_size <-- diagonal block size (unused here)
122 * convergence <-- convergence information structure
123 * rhs <-- right hand side
124 * vx_ini <-- initial system solution
125 * (vx if nonzero, nullptr if zero)
126 * vx <-> system solution
127 * aux_size <-- number of elements in aux_vectors (in bytes)
128 * aux_vectors --- optional working area (allocation otherwise)
129 *
130 * returns:
131 * convergence state
132 *----------------------------------------------------------------------------*/
133
136 const cs_matrix_t *a,
137 cs_lnum_t diag_block_size,
138 cs_sles_it_convergence_t *convergence,
139 const cs_real_t *rhs,
140 cs_real_t *vx_ini,
141 cs_real_t *vx,
142 size_t aux_size,
143 void *aux_vectors);
144
145/*=============================================================================
146 * Public function prototypes
147 *============================================================================*/
148
149/*----------------------------------------------------------------------------
150 * Define and associate an iterative sparse linear system solver
151 * for a given field or equation name.
152 *
153 * If this system did not previously exist, it is added to the list of
154 * "known" systems. Otherwise, its definition is replaced by the one
155 * defined here.
156 *
157 * This is a utility function: if finer control is needed, see
158 * cs_sles_define() and cs_sles_it_create().
159 *
160 * Note that this function returns a pointer directly to the iterative solver
161 * management structure. This may be used to set further options,
162 * for example using cs_sles_set_plot_options(). If needed, cs_sles_find()
163 * may be used to obtain a pointer to the matching cs_sles_t container.
164 *
165 * parameters:
166 * f_id <-- associated field id, or < 0
167 * name <-- associated name if f_id < 0, or NULL
168 * solver_type <-- type of solver (PCG, Jacobi, ...)
169 * poly_degree <-- preconditioning polynomial degree
170 * (0: diagonal; -1: non-preconditioned)
171 * n_max_iter <-- maximum number of iterations
172 *
173 * returns:
174 * pointer to newly created iterative solver info object.
175 *----------------------------------------------------------------------------*/
176
178cs_sles_it_define(int f_id,
179 const char *name,
180 cs_sles_it_type_t solver_type,
181 int poly_degree,
182 int n_max_iter);
183
184/*----------------------------------------------------------------------------
185 * Create iterative sparse linear system solver info and context.
186 *
187 * parameters:
188 * solver_type <-- type of solver (PCG, Jacobi, ...)
189 * poly_degree <-- preconditioning polynomial degree
190 * (0: diagonal; -1: non-preconditioned)
191 * n_max_iter <-- maximum number of iterations
192 * update_stats <-- automatic solver statistics indicator
193 *
194 * returns:
195 * pointer to newly created solver info object.
196 *----------------------------------------------------------------------------*/
197
200 int poly_degree,
201 int n_max_iter,
202 bool update_stats);
203
204/*----------------------------------------------------------------------------
205 * Destroy iterative sparse linear system solver info and context.
206 *
207 * parameters:
208 * context <-> pointer to iterative sparse linear solver info
209 * (actual type: cs_sles_it_t **)
210 *----------------------------------------------------------------------------*/
211
212void
213cs_sles_it_destroy(void **context);
214
215/*----------------------------------------------------------------------------
216 * Create iterative sparse linear system solver info and context
217 * based on existing info and context.
218 *
219 * parameters:
220 * context <-- pointer to reference info and context
221 * (actual type: cs_sles_it_t *)
222 *
223 * returns:
224 * pointer to newly created solver info object
225 * (actual type: cs_sles_it_t *)
226 *----------------------------------------------------------------------------*/
227
228void *
229cs_sles_it_copy(const void *context);
230
231/*----------------------------------------------------------------------------
232 * Setup iterative sparse linear equation solver.
233 *
234 * parameters:
235 * context <-> pointer to iterative sparse linear solver info
236 * (actual type: cs_sles_it_t *)
237 * name <-- pointer to system name
238 * a <-- associated matrix
239 * verbosity <-- verbosity level
240 *----------------------------------------------------------------------------*/
241
242void
243cs_sles_it_setup(void *context,
244 const char *name,
245 const cs_matrix_t *a,
246 int verbosity);
247
248/*----------------------------------------------------------------------------
249 * Call iterative sparse linear equation solver.
250 *
251 * parameters:
252 * context <-> pointer to iterative sparse linear solver info
253 * (actual type: cs_sles_it_t *)
254 * name <-- pointer to system name
255 * a <-- matrix
256 * verbosity <-- verbosity level
257 * precision <-- solver precision
258 * r_norm <-- residual normalization
259 * n_iter --> number of iterations
260 * residual --> residual
261 * rhs <-- right hand side
262 * vx_ini <-- initial solution (vx if nonzero, nullptr if zero)
263 * vx <-> system solution
264 * aux_size <-- number of elements in aux_vectors (in bytes)
265 * aux_vectors --- optional working area (internal allocation if NULL)
266 *
267 * returns:
268 * convergence state
269 *----------------------------------------------------------------------------*/
270
272cs_sles_it_solve(void *context,
273 const char *name,
274 const cs_matrix_t *a,
275 int verbosity,
276 double precision,
277 double r_norm,
278 int *n_iter,
279 double *residual,
280 const cs_real_t *rhs,
281 cs_real_t *vx,
282 cs_real_t *vx_ini,
283 size_t aux_size,
284 void *aux_vectors);
285
286/*----------------------------------------------------------------------------
287 * Free iterative sparse linear equation solver setup context.
288 *
289 * This function frees resolution-related data, such as
290 * buffers and preconditioning but does not free the whole context,
291 * as info used for logging (especially performance data) is maintained.
292
293 * parameters:
294 * context <-> pointer to iterative sparse linear solver info
295 * (actual type: cs_sles_it_t *)
296 *----------------------------------------------------------------------------*/
297
298void
299cs_sles_it_free(void *context);
300
301/*----------------------------------------------------------------------------
302 * Log sparse linear equation solver info.
303 *
304 * parameters:
305 * context <-> pointer to iterative sparse linear solver info
306 * (actual type: cs_sles_it_t *)
307 * log_type <-- log type
308 *----------------------------------------------------------------------------*/
309
310void
311cs_sles_it_log(const void *context,
312 cs_log_t log_type);
313
314/*----------------------------------------------------------------------------
315 * Return iterative solver type.
316 *
317 * parameters:
318 * context <-- pointer to iterative solver info and context
319 *
320 * returns:
321 * selected solver type
322 *----------------------------------------------------------------------------*/
323
325cs_sles_it_get_type(const cs_sles_it_t *context);
326
327/*----------------------------------------------------------------------------
328 * Return the initial residual for the previous solve operation with a solver.
329 *
330 * This is useful for convergence tests when this solver is used as
331 * a preconditioning smoother.
332 *
333 * This operation is only valid between calls to cs_sles_it_setup()
334 * (or cs_sles_it_solve()) and cs_sles_it_free().
335 * It returns -1 otherwise.
336 *
337 * parameters:
338 * context <-- pointer to iterative solver info and context
339 *
340 * returns:
341 * initial residual from last call to \ref cs_sles_solve with this solver
342 *----------------------------------------------------------------------------*/
343
344double
346
347/*----------------------------------------------------------------------------
348 * Return a preconditioner context for an iterative sparse linear
349 * equation solver.
350 *
351 * This allows modifying parameters of a non default (Jacobi or polynomial)
352 * preconditioner.
353 *
354 * parameters:
355 * context <-- pointer to iterative solver info and context
356 *
357 * returns:
358 * pointer to preconditoner context
359 *----------------------------------------------------------------------------*/
360
363
364/*----------------------------------------------------------------------------
365 * \brief Assign a preconditioner to an iterative sparse linear equation
366 * solver, transfering its ownership to the solver context.
367 *
368 * This allows assigning a non default (Jacobi or polynomial) preconditioner.
369 *
370 * The input pointer is set to NULL to make it clear the caller does not
371 * own the preconditioner anymore, though the context can be accessed using
372 * \ref cs_sles_it_get_pc.
373 *
374 * \param[in, out] context pointer to iterative solver info and context
375 * \param[in, out] pc pointer to preconditoner context
376 *----------------------------------------------------------------------------*/
377
378void
380 cs_sles_pc_t **pc);
381
382/*----------------------------------------------------------------------------
383 * Copy options from one iterative sparse linear system solver info
384 * and context to another.
385 *
386 * Optional plotting contexts are shared between the source and destination
387 * contexts.
388 *
389 * Preconditioner settings are to be handled separately.
390 *
391 * parameters:
392 * src <-- pointer to source info and context
393 * dest <-> pointer to destination info and context
394 *----------------------------------------------------------------------------*/
395
396void
398 cs_sles_it_t *dest);
399
400/*----------------------------------------------------------------------------
401 * Associate a similar info and context object with which some setup
402 * data may be shared.
403 *
404 * This is especially useful for sharing preconditioning data between
405 * similar solver contexts (for example ascending and descending multigrid
406 * smoothers based on the same matrix).
407 *
408 * For preconditioning data to be effectively shared, cs_sles_it_setup()
409 * (or cs_sles_it_solve()) must be called on "shareable" before being
410 * called on "context" (without cs_sles_it_free() being called in between,
411 * of course).
412 *
413 * It is the caller's responsibility to ensure the context is not used
414 * for a cs_sles_it_setup() or cs_sles_it_solve() operation after the
415 * shareable object has been destroyed (normally by cs_sles_it_destroy()).
416 *
417 * parameters:
418 * context <-> pointer to iterative sparse linear system solver info
419 * shareable <-- pointer to iterative solver info and context
420 * whose context may be shared
421 *----------------------------------------------------------------------------*/
422
423void
425 const cs_sles_it_t *shareable);
426
427#if defined(HAVE_MPI)
428
429/*----------------------------------------------------------------------------*/
441/*----------------------------------------------------------------------------*/
442
443void
444cs_sles_it_set_mpi_reduce_comm(cs_sles_it_t *context,
445 MPI_Comm comm,
446 MPI_Comm caller_comm);
447
448#endif /* defined(HAVE_MPI) */
449
450/*----------------------------------------------------------------------------
451 * Assign ordering to iterative solver.
452 *
453 * The solver context takes ownership of the order array (i.e. it will
454 * handle its later deallocation).
455 *
456 * This is useful only for Block Gauss-Seidel.
457 *
458 * parameters:
459 * context <-> pointer to iterative solver info and context
460 * order <-> pointer to ordering array
461 *----------------------------------------------------------------------------*/
462
463void
465 cs_lnum_t **order);
466
467/*----------------------------------------------------------------------------*/
474/*----------------------------------------------------------------------------*/
475
476double
478
479/*----------------------------------------------------------------------------*/
486/*----------------------------------------------------------------------------*/
487
488void
489cs_sles_it_set_breakdown_threshold(double threshold);
490
491/*----------------------------------------------------------------------------*/
508/*----------------------------------------------------------------------------*/
509
510void
513 int n_iter_max);
514
515/*----------------------------------------------------------------------------*/
523/*----------------------------------------------------------------------------*/
524
525void
527 int interval);
528
529/*----------------------------------------------------------------------------*/
536/*----------------------------------------------------------------------------*/
537
538void
540 int n_max_iter);
541
542/*----------------------------------------------------------------------------
543 * Error handler for iterative sparse linear equation solver.
544 *
545 * In case of divergence or breakdown, this error handler outputs
546 * postprocessing data to assist debugging, then aborts the run.
547 * It does nothing in case the maximum iteration count is reached.
548 *
549 * parameters:
550 * sles <-> pointer to solver object
551 * state <-- convergence state
552 * a <-- matrix
553 * rhs <-- right hand side
554 * vx <-> system solution
555 *
556 * returns:
557 * false (do not attempt new solve)
558 *----------------------------------------------------------------------------*/
559
560bool
563 const cs_matrix_t *a,
564 const cs_real_t *rhs,
565 cs_real_t *vx);
566
567/*----------------------------------------------------------------------------
568 * Set plotting options for an iterative sparse linear equation solver.
569 *
570 * parameters:
571 * context <-> pointer to iterative solver info and context
572 * base_name <-- base plot name to activate, NULL otherwise
573 * use_iteration <-- if true, use iteration as time stamp
574 * otherwise, use wall clock time
575 *----------------------------------------------------------------------------*/
576
577void
579 const char *base_name,
580 bool use_iteration);
581
582/*----------------------------------------------------------------------------
583 * Convergence test.
584 *
585 * parameters:
586 * c <-- pointer to solver context info
587 * n_iter <-- Number of iterations done
588 * residual <-- Non normalized residual
589 * convergence <-> Convergence information structure
590 *
591 * returns:
592 * convergence status.
593 *----------------------------------------------------------------------------*/
594
597 unsigned n_iter,
598 double residual,
599 cs_sles_it_convergence_t *convergence);
600
601/*----------------------------------------------------------------------------*/
602
604
605#endif /* __CS_SLES_IT_H__ */
#define BEGIN_C_DECLS
Definition: cs_defs.h:554
double cs_real_t
Floating-point value.
Definition: cs_defs.h:357
#define END_C_DECLS
Definition: cs_defs.h:555
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:350
cs_log_t
Definition: cs_log.h:48
struct _cs_matrix_t cs_matrix_t
Definition: cs_matrix.h:114
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_residual(const cs_sles_it_t *context)
Return the initial residual for the previous solve operation with a solver.
Definition: cs_sles_it.cpp:5286
cs_sles_convergence_state_t cs_sles_it_solve(void *context, const char *name, const cs_matrix_t *a, int verbosity, double precision, double r_norm, int *n_iter, double *residual, const cs_real_t *rhs, cs_real_t *vx, cs_real_t *vx_ini, size_t aux_size, void *aux_vectors)
Call iterative sparse linear equation solver.
Definition: cs_sles_it.cpp:4975
double cs_sles_it_get_breakdown_threshold(void)
Retrieve the threshold value under which a breakdown happens in solvers like BiCGStab or BiCGStab2.
Definition: cs_sles_it.cpp:5532
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.cpp:5339
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.cpp:5310
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.cpp:5371
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.cpp:4463
bool cs_sles_it_error_post_and_abort(cs_sles_t *sles, cs_sles_convergence_state_t state, const cs_matrix_t *a, const cs_real_t *rhs, cs_real_t *vx)
Error handler for iterative sparse linear equation solver.
Definition: cs_sles_it.cpp:5638
void cs_sles_it_set_n_max_iter(cs_sles_it_t *context, int n_max_iter)
Define the max. number of iterations before stopping the algorithm.
Definition: cs_sles_it.cpp:5610
void cs_sles_it_set_fallback_threshold(cs_sles_it_t *context, cs_sles_convergence_state_t threshold, int n_iter_max)
Define convergence level under which the fallback to another solver may be used if applicable.
Definition: cs_sles_it.cpp:5572
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.cpp:4768
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.cpp:5420
void cs_sles_it_set_restart_interval(cs_sles_it_t *context, int interval)
Define the number of iterations to be done before restarting the solver. Useful only for GCR or GMRES...
Definition: cs_sles_it.cpp:5591
void cs_sles_it_log(const void *context, cs_log_t log_type)
Log sparse linear equation solver info.
Definition: cs_sles_it.cpp:4652
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.cpp:4410
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.cpp:4606
cs_sles_convergence_state_t cs_sles_it_convergence_test(cs_sles_it_t *c, unsigned n_iter, double residual, cs_sles_it_convergence_t *convergence)
Definition: cs_sles_it.cpp:5729
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.cpp:5685
struct _cs_sles_it_t cs_sles_it_t
Definition: cs_sles_it.h:95
void cs_sles_it_free(void *context)
Free iterative sparse linear equation solver setup context.
Definition: cs_sles_it.cpp:5217
void cs_sles_it_assign_order(cs_sles_it_t *context, cs_lnum_t **order)
Assign ordering to iterative solver.
Definition: cs_sles_it.cpp:5500
const char * cs_sles_it_type_name[]
struct _cs_sles_it_convergence_t cs_sles_it_convergence_t
Definition: cs_sles_it.h:99
cs_sles_it_type_t
Definition: cs_sles_it.h:55
@ CS_SLES_P_SYM_GAUSS_SEIDEL
Definition: cs_sles_it.h:71
@ CS_SLES_JACOBI
Definition: cs_sles_it.h:63
@ CS_SLES_RJ2
Definition: cs_sles_it.h:82
@ CS_SLES_PCR3
Definition: cs_sles_it.h:72
@ CS_SLES_GMRES
Definition: cs_sles_it.h:68
@ CS_SLES_GCR
Definition: cs_sles_it.h:67
@ CS_SLES_P_GAUSS_SEIDEL
Definition: cs_sles_it.h:70
@ CS_SLES_BICGSTAB
Definition: cs_sles_it.h:64
@ CS_SLES_N_IT_TYPES
Definition: cs_sles_it.h:75
@ CS_SLES_BICGSTAB2
Definition: cs_sles_it.h:66
@ CS_SLES_L1_JACOBI
Definition: cs_sles_it.h:80
@ CS_SLES_USER_DEFINED
Definition: cs_sles_it.h:73
@ CS_SLES_RJ3
Definition: cs_sles_it.h:83
@ CS_SLES_PCG
Definition: cs_sles_it.h:59
@ CS_SLES_TS_B_GAUSS_SEIDEL
Definition: cs_sles_it.h:86
@ CS_SLES_IPCG
Definition: cs_sles_it.h:62
@ CS_SLES_R_JACOBI
Definition: cs_sles_it.h:81
@ CS_SLES_N_SMOOTHER_TYPES
Definition: cs_sles_it.h:88
@ CS_SLES_FCG
Definition: cs_sles_it.h:60
@ CS_SLES_TS_F_GAUSS_SEIDEL
Definition: cs_sles_it.h:85
void cs_sles_it_set_breakdown_threshold(double threshold)
Define the threshold value under which a breakdown happens in solvers like BiCGStab or BiCGStab2.
Definition: cs_sles_it.cpp:5547
cs_sles_convergence_state_t cs_user_sles_it_solver(cs_sles_it_t *c, const cs_matrix_t *a, cs_lnum_t diag_block_size, cs_sles_it_convergence_t *convergence, const cs_real_t *rhs, cs_real_t *vx_ini, cs_real_t *vx, size_t aux_size, void *aux_vectors)
void cs_sles_it_destroy(void **context)
Destroy iterative sparse linear system solver info and context.
Definition: cs_sles_it.cpp:4568
cs_sles_it_type_t cs_sles_it_get_type(const cs_sles_it_t *context)
Return iterative solver type.
Definition: cs_sles_it.cpp:5258
struct _cs_sles_pc_t cs_sles_pc_t
Definition: cs_sles_pc.h:66