My Project
programmer's documentation
cs_sles.h
Go to the documentation of this file.
1 #ifndef __CS_SLES_H__
2 #define __CS_SLES_H__
3 
4 /*============================================================================
5  * Sparse Linear Equation Solver driver
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_log.h"
36 #include "cs_halo_perio.h"
37 #include "cs_matrix.h"
38 #include "cs_matrix.h"
39 
40 /*----------------------------------------------------------------------------*/
41 
43 
44 /*============================================================================
45  * Macro definitions
46  *============================================================================*/
47 
48 /*============================================================================
49  * Type definitions
50  *============================================================================*/
51 
52 /*----------------------------------------------------------------------------
53  * Convergence status
54  *----------------------------------------------------------------------------*/
55 
56 typedef enum {
57 
63 
65 
66 /* General linear solver context (opaque) */
67 
68 typedef struct _cs_sles_t cs_sles_t;
69 
70 /*----------------------------------------------------------------------------
71  * Function pointer for pre-resolution setup of a linear system solvers's
72  * context.
73  *
74  * This setup may include building a multigrid hierarchy, or a preconditioner.
75  *
76  * Use of this type of function is optional: the context is expected to
77  * maintain state, so that if a cs_sles_solve_t function is called before a
78  * cs_sles_setup_t function, the latter will be called automatically.
79  *
80  * parameters:
81  * context <-> pointer to solver context
82  * name <-- pointer to name of linear system
83  * a <-- matrix
84  * verbosity <-- associated verbosity
85  *----------------------------------------------------------------------------*/
86 
87 typedef void
88 (cs_sles_setup_t) (void *context,
89  const char *name,
90  const cs_matrix_t *a,
91  int verbosity);
92 
93 /*----------------------------------------------------------------------------
94  * Function pointer for resolution of a linear system.
95  *
96  * If the associated cs_sles_setup_t function has not been called before
97  * this function, it will be called automatically.
98  *
99  * The solution context setup by this call (or that of the matching setup
100  * function) will be maintained until the matching cs_sles_free_t function
101  * is called.
102  *
103  * The matrix is not expected to change between successive calls, although
104  * the right hand side may. If the matrix changes, the associated
105  * cs_sles_setup_t or cs_sles_free_t function must be called between
106  * solves.
107  *
108  * The system is considered to have converged when
109  * residue/r_norm <= precision, residue being the L2 norm of a.vx-rhs.
110  *
111  * parameters:
112  * context <-> pointer to solver context
113  * name <-- pointer to name of linear system
114  * a <-- matrix
115  * verbosity <-- associated verbosity
116  * rotation_mode <-- halo update option for rotational periodicity
117  * precision <-- solver precision
118  * r_norm <-- residue normalization
119  * n_iter --> number of "equivalent" iterations
120  * residue --> residue
121  * rhs <-- right hand side
122  * vx <-- system solution
123  * aux_size <-- number of elements in aux_vectors
124  * aux_vectors <-- optional working area (internal allocation if NULL)
125  *
126  * returns:
127  * convergence status
128  *----------------------------------------------------------------------------*/
129 
131 (cs_sles_solve_t) (void *context,
132  const char *name,
133  const cs_matrix_t *a,
134  int verbosity,
135  cs_halo_rotation_t rotation_mode,
136  double precision,
137  double r_norm,
138  int *n_iter,
139  double *residue,
140  const cs_real_t *rhs,
141  cs_real_t *vx,
142  size_t aux_size,
143  void *aux_vectors);
144 
145 /*----------------------------------------------------------------------------
146  * Function pointer for freeing of a linear system's context data.
147  *
148  * Note that this function should free resolution-related data, such as
149  * multigrid hierarchy, preconditioning, and any other temporary arrays or
150  * objects required for resolution, but should not free the whole context,
151  * as info used for logging (especially performance data) should be
152  * maintained.
153  *
154  * parameters:
155  * context <-> pointer to solver context
156  *----------------------------------------------------------------------------*/
157 
158 typedef void
159 (cs_sles_free_t) (void *context);
160 
161 /*----------------------------------------------------------------------------
162  * Function pointer for logging of linear solver setup,
163  * history and performance data.
164  *
165  * This function will be called for each solver when cs_sles_finalize()
166  * is called.
167  *
168  * parameters:
169  * context <-- pointer to solver context
170  * log_type <-- log type
171  *----------------------------------------------------------------------------*/
172 
173 typedef void
174 (cs_sles_log_t) (const void *context,
175  cs_log_t log_type);
176 
177 /*----------------------------------------------------------------------------
178  * Function pointer for creation of a solver context based on the copy
179  * of another.
180  *
181  * The new context copies the settings of the copied context, but not
182  * its setup data and logged info, such as performance data.
183  *
184  * This type of function is optional, but enables associating different
185  * solvers to related systems (to differentiate logging) while using
186  * the same settings by default.
187  *
188  * parameters:
189  * context <-- source context
190  *
191  * returns:
192  * pointer to newly created context
193  *----------------------------------------------------------------------------*/
194 
195 typedef void *
196 (cs_sles_copy_t) (const void *context);
197 
198 /*----------------------------------------------------------------------------
199  * Function pointer for destruction of a linear system solver context.
200  *
201  * This function should free all context data, and will be called for each
202  * system when cs_sles_finalize() is called.
203  *
204  * parameters:
205  * context <-> pointer to solver context
206  *----------------------------------------------------------------------------*/
207 
208 typedef void
209 (cs_sles_destroy_t) (void **context);
210 
211 /*----------------------------------------------------------------------------
212  * Function pointer for handling of non-convegence when solving
213  * a linear system.
214  *
215  * Such a function is optional, and may be used for a variety of purposes,
216  * such as logging, postprocessing, re-trying with different parameters,
217  * aborting the run, or any combination thereof.
218  *
219  * An error handler may be associated with a given solver using
220  * cs_sles_set_error_handler(), in which case it will be called whenever
221  * convergence fails.
222  *
223  * parameters:
224  * sles <-> pointer to solver object
225  * state <-- convergence status
226  * a <-- matrix
227  * rotation_mode <-- Halo update option for rotational periodicity
228  * rhs <-- Right hand side
229  * vx <-- System solution
230  *
231  * returns:
232  * true if solve should be re-executed, false otherwise
233  *----------------------------------------------------------------------------*/
234 
235 typedef bool
238  const cs_matrix_t *a,
239  cs_halo_rotation_t rotation_mode,
240  const cs_real_t *rhs,
241  cs_real_t *vx);
242 
243 /*----------------------------------------------------------------------------
244  * Function pointer for the default definition of a sparse
245  * linear equation solver
246  *
247  * The function may be associated using cs_sles_set_default_define(), so
248  * that it may provide a definition that will be used when
249  * cs_sles_setup() or cs_sles_solve() is used for a system for which
250  * no matching call to cs_sles_define() has been done.
251  *
252  * The function should call cs_sles_define() with arguments f_id
253  * and name, and appropriately chosen function pointers.
254  *
255  * A pointer to the matrix of the system to be solved is also provided,
256  * so that the corresponding information may be used to better choose
257  * defaults.
258  *
259  * parameters:
260  * f_id <-- associated field id, or < 0
261  * name <-- associated name if f_id < 0, or NULL
262  * a <-- Matrix
263  *----------------------------------------------------------------------------*/
264 
265 typedef void
267  const char *name,
268  const cs_matrix_t *a);
269 
270 /*----------------------------------------------------------------------------
271  * Function pointer for the default definition of a sparse
272  * linear equation solver's verbosity
273  *
274  * The function may be associated using cs_sles_set_default_verbosity(), so
275  * that it may provide a definition that will be used when
276  * cs_sles_default_verbosity() is called.
277  *
278  * parameters:
279  * f_id <-- associated field id, or < 0
280  * name <-- associated name if f_id < 0, or NULL
281  *
282  * returns:
283  * default verbosity value
284  *----------------------------------------------------------------------------*/
285 
286 typedef int
288  const char *name);
289 
290 /*============================================================================
291  * Global variables
292  *============================================================================*/
293 
294 /*=============================================================================
295  * Public function prototypes for Fortran API
296  *============================================================================*/
297 
298 /*=============================================================================
299  * Public function prototypes
300  *============================================================================*/
301 
302 /*----------------------------------------------------------------------------
303  * \brief Initialize sparse linear equation solver API.
304  *----------------------------------------------------------------------------*/
305 
306 void
307 cs_sles_initialize(void);
308 
309 /*----------------------------------------------------------------------------
310  * \brief Finalize sparse linear equation solver API.
311  *----------------------------------------------------------------------------*/
312 
313 void
314 cs_sles_finalize(void);
315 
316 /*----------------------------------------------------------------------------*/
322 /*----------------------------------------------------------------------------*/
323 
324 void
325 cs_sles_log(cs_log_t log_type);
326 
327 /*----------------------------------------------------------------------------*/
339 /*----------------------------------------------------------------------------*/
340 
341 cs_sles_t *
342 cs_sles_find(int f_id,
343  const char *name);
344 
345 /*----------------------------------------------------------------------------*/
362 /*----------------------------------------------------------------------------*/
363 
364 cs_sles_t *
366  const char *name);
367 
368 /*----------------------------------------------------------------------------*/
384 /*----------------------------------------------------------------------------*/
385 
386 void
387 cs_sles_push(int f_id,
388  const char *name);
389 
390 /*----------------------------------------------------------------------------*/
398 /*----------------------------------------------------------------------------*/
399 
400 void
401 cs_sles_pop(int f_id);
402 
403 /*----------------------------------------------------------------------------*/
438 /*----------------------------------------------------------------------------*/
439 
440 cs_sles_t *
441 cs_sles_define(int f_id,
442  const char *name,
443  void *context,
444  const char *type_name,
445  cs_sles_setup_t *setup_func,
446  cs_sles_solve_t *solve_func,
447  cs_sles_free_t *free_func,
448  cs_sles_log_t *log_func,
449  cs_sles_copy_t *copy_func,
450  cs_sles_destroy_t *destroy_func);
451 
452 /*----------------------------------------------------------------------------*/
464 /*----------------------------------------------------------------------------*/
465 
466 void
468  int verbosity);
469 
470 /*----------------------------------------------------------------------------*/
480 /*----------------------------------------------------------------------------*/
481 
482 int
484 
485 /*----------------------------------------------------------------------------*/
496 /*----------------------------------------------------------------------------*/
497 
498 void
500  int writer_id);
501 
502 /*----------------------------------------------------------------------------*/
511 /*----------------------------------------------------------------------------*/
512 
513 int
515 
516 /*----------------------------------------------------------------------------*/
535 /*----------------------------------------------------------------------------*/
536 
537 const char *
539 
540 /*----------------------------------------------------------------------------*/
552 /*----------------------------------------------------------------------------*/
553 
554 void *
556 
557 /*----------------------------------------------------------------------------*/
565 /*----------------------------------------------------------------------------*/
566 
567 int
568 cs_sles_get_f_id(const cs_sles_t *sles);
569 
570 /*----------------------------------------------------------------------------*/
581 /*----------------------------------------------------------------------------*/
582 
583 const char *
584 cs_sles_get_name(const cs_sles_t *sles);
585 
586 /*----------------------------------------------------------------------------*/
600 /*----------------------------------------------------------------------------*/
601 
602 void
604  const cs_matrix_t *a);
605 
606 /*----------------------------------------------------------------------------*/
637 /*----------------------------------------------------------------------------*/
638 
641  const cs_matrix_t *a,
642  cs_halo_rotation_t rotation_mode,
643  double precision,
644  double r_norm,
645  int *n_iter,
646  double *residue,
647  const cs_real_t *rhs,
648  cs_real_t *vx,
649  size_t aux_size,
650  void *aux_vectors);
651 
652 /*----------------------------------------------------------------------------*/
663 /*----------------------------------------------------------------------------*/
664 
665 void
666 cs_sles_free(cs_sles_t *sles);
667 
668 /*----------------------------------------------------------------------------*/
688 /*----------------------------------------------------------------------------*/
689 
690 int
691 cs_sles_copy(cs_sles_t *dest,
692  const cs_sles_t *src);
693 
694 /*----------------------------------------------------------------------------*/
709 /*----------------------------------------------------------------------------*/
710 
711 void
713  cs_sles_error_handler_t *error_handler_func);
714 
715 /*----------------------------------------------------------------------------*/
725 /*----------------------------------------------------------------------------*/
726 
729 
730 /*----------------------------------------------------------------------------*/
740 /*----------------------------------------------------------------------------*/
741 
742 void
744 
745 /*----------------------------------------------------------------------------*/
754 /*----------------------------------------------------------------------------*/
755 
756 void
758 
759 /*----------------------------------------------------------------------------*/
770 /*----------------------------------------------------------------------------*/
771 
772 void
773 cs_sles_post_error_output_def(const char *name,
774  int mesh_id,
775  cs_halo_rotation_t rotation_mode,
776  const cs_matrix_t *a,
777  const cs_real_t *rhs,
778  cs_real_t *vx);
779 
780 /*----------------------------------------------------------------------------*/
792 /*----------------------------------------------------------------------------*/
793 
794 void
795 cs_sles_post_output_var(const char *name,
796  int mesh_id,
797  int location_id,
798  int writer_id,
799  int diag_block_size,
800  cs_real_t var[]);
801 
802 /*----------------------------------------------------------------------------*/
814 /*----------------------------------------------------------------------------*/
815 
816 const char *
818  const char *name);
819 
820 /*----------------------------------------------------------------------------*/
829 /*----------------------------------------------------------------------------*/
830 
831 const char *
832 cs_sles_name(int f_id,
833  const char *name);
834 
835 /*----------------------------------------------------------------------------*/
836 
838 
839 #endif /* __CS_SLES_H__ */
f_id
void const int * f_id
Definition: cs_gui.h:292
cs_sles_log_t
void() cs_sles_log_t(const void *context, cs_log_t log_type)
Function pointer for logging of linear solver history and performance data.
Definition: cs_sles.h:174
cs_sles_set_default_verbosity
void cs_sles_set_default_verbosity(cs_sles_verbosity_t *verbosity_func)
Set default verbosity definition function.
Definition: cs_sles.c:1843
cs_sles_set_verbosity
void cs_sles_set_verbosity(cs_sles_t *sles, int verbosity)
Set the verbosity for a given linear equation solver.
Definition: cs_sles.c:1345
cs_sles_get_type
const char * cs_sles_get_type(cs_sles_t *sles)
Return type name of solver context.
Definition: cs_sles.c:1443
var
void cs_int_t cs_real_t var[]
Definition: cs_parall.h:63
CS_SLES_CONVERGED
Definition: cs_sles.h:62
cs_fuel_incl::a
double precision, save a
Definition: cs_fuel_incl.f90:146
cs_sles_find_or_add
cs_sles_t * cs_sles_find_or_add(int f_id, const char *name)
Return pointer to linear system object, based on matching field id or system name.
Definition: cs_sles.c:1167
cs_log.h
cs_sles_post_output_var
void cs_sles_post_output_var(const char *name, int mesh_id, int location_id, int writer_id, int diag_block_size, cs_real_t var[])
Output post-processing variable related to system convergence.
Definition: cs_sles.c:1976
cs_sles_get_context
void * cs_sles_get_context(cs_sles_t *sles)
Return pointer to solver context structure pointer.
Definition: cs_sles.c:1463
cs_sles_verbosity_t
int() cs_sles_verbosity_t(int f_id, const char *name)
Function pointer for the default definition of a sparse linear equation solver's verbosity.
Definition: cs_sles.h:287
cs_sles_initialize
void cs_sles_initialize(void)
Initialize sparse linear equation solver API.
Definition: cs_sles.c:891
cs_sles_setup
void cs_sles_setup(cs_sles_t *sles, const cs_matrix_t *a)
Setup sparse linear equation solver.
Definition: cs_sles.c:1520
cs_sles_push
void cs_sles_push(int f_id, const char *name)
Temporarily replace field id with name for matching calls to cs_sles_setup, cs_sles_solve,...
Definition: cs_sles.c:1204
END_C_DECLS
#define END_C_DECLS
Definition: cs_defs.h:468
cs_sles_post_error_output_def
void cs_sles_post_error_output_def(const char *name, int mesh_id, cs_halo_rotation_t rotation_mode, const cs_matrix_t *a, const cs_real_t *rhs, cs_real_t *vx)
Output default post-processing data for failed system convergence.
Definition: cs_sles.c:1862
cs_sles_finalize
void cs_sles_finalize(void)
Finalize sparse linear equation solver API.
Definition: cs_sles.c:911
cs_sles_solve
cs_sles_convergence_state_t cs_sles_solve(cs_sles_t *sles, const cs_matrix_t *a, 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)
General sparse linear system resolution.
Definition: cs_sles.c:1589
cs_sles_get_name
const char * cs_sles_get_name(const cs_sles_t *sles)
Return name associated with a given sparse linear equation solver.
Definition: cs_sles.c:1498
cs_sles_destroy_t
void() cs_sles_destroy_t(void **context)
Definition: cs_sles.h:209
cs_real_t
double cs_real_t
Floating-point value.
Definition: cs_defs.h:302
cs_sles_find
cs_sles_t * cs_sles_find(int f_id, const char *name)
Return pointer to linear system object, based on matching field id or system name.
Definition: cs_sles.c:1102
CS_SLES_ITERATING
Definition: cs_sles.h:61
BEGIN_C_DECLS
#define BEGIN_C_DECLS
Definition: cs_defs.h:467
cs_matrix_t
struct _cs_matrix_t cs_matrix_t
Definition: cs_matrix.h:90
cs_sles_convergence_state_t
cs_sles_convergence_state_t
Convergence status indicator.
Definition: cs_sles.h:56
cs_sles_copy
int cs_sles_copy(cs_sles_t *dest, const cs_sles_t *src)
Copy the definition of a sparse linear equation solver to another.
Definition: cs_sles.c:1731
cs_sles_pop
void cs_sles_pop(int f_id)
Restore behavior temporarily modified by cs_sles_push.
Definition: cs_sles.c:1240
CS_SLES_DIVERGED
Definition: cs_sles.h:58
cs_sles_free
void cs_sles_free(cs_sles_t *sles)
Free sparse linear equation solver setup.
Definition: cs_sles.c:1700
cs_sles_get_post_output
int cs_sles_get_post_output(cs_sles_t *sles)
Return the id of the associated writer if postprocessing output is active for a given linear equation...
Definition: cs_sles.c:1411
cs_sles_get_verbosity
int cs_sles_get_verbosity(cs_sles_t *sles)
Get the verbosity for a given linear equation solver.
Definition: cs_sles.c:1364
cs_sles_name
const char * cs_sles_name(int f_id, const char *name)
Return name associated to a field id, name couple.
Definition: cs_sles.c:2124
cs_sles_setup_t
void() cs_sles_setup_t(void *context, const char *name, const cs_matrix_t *a, int verbosity)
Function pointer for pre-resolution setup of a linear system's context.
Definition: cs_sles.h:88
cs_sles_define_t
void() cs_sles_define_t(int f_id, const char *name, const cs_matrix_t *a)
Function pointer for the default definition of a sparse linear equation solver.
Definition: cs_sles.h:266
cs_sles_set_default_define
void cs_sles_set_default_define(cs_sles_define_t *define_func)
Set default sparse linear solver definition function.
Definition: cs_sles.c:1826
cs_sles_copy_t
void *() cs_sles_copy_t(const void *context)
Function pointer for creation of a solver context based on the copy of another.
Definition: cs_sles.h:196
cs_sles_get_f_id
int cs_sles_get_f_id(const cs_sles_t *sles)
Return field id associated with a given sparse linear equation solver.
Definition: cs_sles.c:1479
cs_halo_perio.h
cs_sles_error_handler_t
bool() cs_sles_error_handler_t(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)
Function pointer for handling of non-convergence when solving a linear system.
Definition: cs_sles.h:236
cs_sles_free_t
void() cs_sles_free_t(void *context)
Function pointer for freeing of a linear system's context data.
Definition: cs_sles.h:159
cs_log_t
cs_log_t
Definition: cs_log.h:48
cs_sles_t
struct _cs_sles_t cs_sles_t
Definition: cs_sles.h:68
cs_sles_define
cs_sles_t * cs_sles_define(int f_id, const char *name, void *context, const char *type_name, cs_sles_setup_t *setup_func, cs_sles_solve_t *solve_func, cs_sles_free_t *free_func, cs_sles_log_t *log_func, cs_sles_copy_t *copy_func, cs_sles_destroy_t *destroy_func)
Define sparse linear equation solver for a given field or equation name.
Definition: cs_sles.c:1292
cs_sles_set_error_handler
void cs_sles_set_error_handler(cs_sles_t *sles, cs_sles_error_handler_t *error_handler_func)
Associate a convergence error handler to a given sparse linear equation solver.
Definition: cs_sles.c:1788
cs_matrix.h
cs_sles_set_post_output
void cs_sles_set_post_output(cs_sles_t *sles, int writer_id)
Activate postprocessing output for a given linear equation solver.
Definition: cs_sles.c:1383
CS_SLES_BREAKDOWN
Definition: cs_sles.h:59
cs_sles_solve_t
cs_sles_convergence_state_t() cs_sles_solve_t(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)
Function pointer for resolution of a linear system.
Definition: cs_sles.h:131
cs_sles_log
void cs_sles_log(cs_log_t log_type)
Log sparse linear equation solver info.
Definition: cs_sles.c:951
cs_halo_rotation_t
cs_halo_rotation_t
Definition: cs_halo.h:60
CS_SLES_MAX_ITERATION
Definition: cs_sles.h:60
cs_sles_base_name
const char * cs_sles_base_name(int f_id, const char *name)
Return base name associated to a field id, name couple.
Definition: cs_sles.c:2099
cs_base.h
cs_sles_get_default_define
cs_sles_define_t * cs_sles_get_default_define(void)
Return pointer to default sparse linear solver definition function.
Definition: cs_sles.c:1808