8.1
general documentation
cs_multigrid.h
Go to the documentation of this file.
1 #ifndef __CS_MULTIGRID_H__
2 #define __CS_MULTIGRID_H__
3 
4 /*============================================================================
5  * Multigrid solver.
6  *============================================================================*/
7 
8 /*
9  This file is part of code_saturne, a general-purpose CFD tool.
10 
11  Copyright (C) 1998-2023 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_grid.h"
36 #include "cs_sles.h"
37 #include "cs_sles_it.h"
38 #include "cs_sles_pc.h"
39 #include "cs_time_plot.h"
40 
41 /*----------------------------------------------------------------------------*/
42 
44 
45 /*============================================================================
46  * Macro definitions
47  *============================================================================*/
48 
49 /*============================================================================
50  * Type definitions
51  *============================================================================*/
52 
53 /*----------------------------------------------------------------------------
54  * Multigrid types
55  *----------------------------------------------------------------------------*/
56 
57 typedef enum {
58 
65 
66 /* Multigrid linear solver context (opaque) */
67 
68 typedef struct _cs_multigrid_t cs_multigrid_t;
69 
70 /*============================================================================
71  * Global variables
72  *============================================================================*/
73 
74 /* Names for multigrid types */
75 
76 extern const char *cs_multigrid_type_name[];
77 
78 /*=============================================================================
79  * Public function prototypes
80  *============================================================================*/
81 
82 /*----------------------------------------------------------------------------
83  * Initialize multigrid solver API.
84  *----------------------------------------------------------------------------*/
85 
86 void
88 
89 /*----------------------------------------------------------------------------
90  * Finalize multigrid solver API.
91  *----------------------------------------------------------------------------*/
92 
93 void
95 
96 /*----------------------------------------------------------------------------
97  * Define and associate a multigrid sparse linear system solver
98  * for a given field or equation name.
99  *
100  * If this system did not previously exist, it is added to the list of
101  * "known" systems. Otherwise, its definition is replaced by the one
102  * defined here.
103  *
104  * This is a utility function: if finer control is needed, see
105  * cs_sles_define() and cs_multigrid_create().
106  *
107  * Note that this function returns a pointer directly to the multigrid solver
108  * management structure. This may be used to set further options, for
109  * example calling cs_multigrid_set_coarsening_options() and
110  * cs_multigrid_set_solver_options().
111  * If needed, cs_sles_find() may be used to obtain a pointer to the
112  * matching cs_sles_t container.
113  *
114  * \param[in] f_id associated field id, or < 0
115  * \param[in] name associated name if f_id < 0, or NULL
116  * \param[in] mg_type type of multigrid algorithm to use
117  *
118  * \return pointer to new multigrid info and context
119  */
120 /*----------------------------------------------------------------------------*/
121 
123 cs_multigrid_define(int f_id,
124  const char *name,
125  cs_multigrid_type_t mg_type);
126 
127 /*----------------------------------------------------------------------------*/
137 /*----------------------------------------------------------------------------*/
138 
141 
142 /*----------------------------------------------------------------------------
143  * Destroy multigrid linear system solver info and context.
144  *
145  * parameters:
146  * context <-> pointer to multigrid linear solver info
147  * (actual type: cs_multigrid_t **)
148  *----------------------------------------------------------------------------*/
149 
150 void
151 cs_multigrid_destroy(void **context);
152 
153 /*----------------------------------------------------------------------------
154  * Create multigrid sparse linear system solver info and context
155  * based on existing info and context.
156  *
157  * parameters:
158  * context <-- pointer to reference info and context
159  * (actual type: cs_multigrid_t *)
160  *
161  * returns:
162  * pointer to newly created solver info object
163  * (actual type: cs_multigrid_t *)
164  *----------------------------------------------------------------------------*/
165 
166 void *
167 cs_multigrid_copy(const void *context);
168 
169 /*----------------------------------------------------------------------------
170  * Set multigrid coarsening parameters.
171  *
172  * parameters:
173  * mg <-> pointer to multigrid info and context
174  * aggregation_limit <-- maximum allowed fine rows per coarse cell
175  * coarsening_type <-- coarsening type; see cs_grid_coarsening_t
176  * n_max_levels <-- maximum number of grid levels
177  * min_g_rows <-- global number of rows on coarse grids
178  * under which no coarsening occurs
179  * p0p1_relax <-- p0/p1 relaxation_parameter
180  * postprocess_block_size <-- if > 0, postprocess coarsening
181  * (using coarse cell numbers modulo this value)
182  *----------------------------------------------------------------------------*/
183 
184 void
186  int aggregation_limit,
187  int coarsening_type,
188  int n_max_levels,
189  cs_gnum_t min_g_rows,
190  double p0p1_relax,
191  int postprocess_block_size);
192 
193 /*----------------------------------------------------------------------------
194  * Set multigrid parameters for associated iterative solvers.
195  *
196  * parameters:
197  * mg <-> pointer to multigrid info and context
198  * descent_smoother_type <-- type of smoother for descent
199  * ascent_smoother_type <-- type of smoother for ascent
200  * coarse_solver_type <-- type of solver
201  * n_max_cycles <-- maximum number of cycles
202  * n_max_iter_descent <-- maximum iterations per descent phase
203  * n_max_iter_ascent <-- maximum iterations per descent phase
204  * n_max_iter_coarse <-- maximum iterations per coarsest solution
205  * poly_degree_descent <-- preconditioning polynomial degree
206  * for descent phases (0: diagonal)
207  * poly_degree_ascent <-- preconditioning polynomial degree
208  * for ascent phases (0: diagonal)
209  * poly_degree_coarse <-- preconditioning polynomial degree
210  * for coarse solver (0: diagonal)
211  * precision_mult_descent <-- precision multiplier for descent phases
212  * (levels >= 1)
213  * precision_mult_ascent <-- precision multiplier for ascent phases
214  * precision_mult_coarse <-- precision multiplier for coarsest grid
215  *----------------------------------------------------------------------------*/
216 
217 void
219  cs_sles_it_type_t descent_smoother_type,
220  cs_sles_it_type_t ascent_smoother_type,
221  cs_sles_it_type_t coarse_solver_type,
222  int n_max_cycles,
223  int n_max_iter_descent,
224  int n_max_iter_ascent,
225  int n_max_iter_coarse,
226  int poly_degree_descent,
227  int poly_degree_ascent,
228  int poly_degree_coarse,
229  double precision_mult_descent,
230  double precision_mult_ascent,
231  double precision_mult_coarse);
232 
233 /*----------------------------------------------------------------------------*/
240 /*----------------------------------------------------------------------------*/
241 
242 void
244  int n_max_cycles);
245 
246 /*----------------------------------------------------------------------------
247  * Return solver type used on fine mesh.
248  *
249  * parameters:
250  * mg <-- pointer to multigrid info and context
251  *
252  * returns:
253  * type of smoother for descent (used for fine mesh)
254  *----------------------------------------------------------------------------*/
255 
258 
259 /*----------------------------------------------------------------------------
260  * Setup multigrid sparse linear equation solver.
261  *
262  * parameters:
263  * context <-> pointer to multigrid info and context
264  * (actual type: cs_multigrid_t *)
265  * name <-- pointer to name of linear system
266  * a <-- associated matrix
267  * verbosity <-- associated verbosity
268  *----------------------------------------------------------------------------*/
269 
270 void
271 cs_multigrid_setup(void *context,
272  const char *name,
273  const cs_matrix_t *a,
274  int verbosity);
275 
276 /*----------------------------------------------------------------------------
277  * Setup multigrid sparse linear equation solver with separate
278  * convection-diffusion matrixes
279  *
280  * parameters:
281  * context <-> pointer to multigrid info and context
282  * (actual type: cs_multigrid_t *)
283  * name <-- pointer to name of linear system
284  * a <-- associated matrix
285  * conv_diff <-- convection-diffusion mode
286  * verbosity <-- associated verbosity
287  *----------------------------------------------------------------------------*/
288 
289 void
290 cs_multigrid_setup_conv_diff(void *context,
291  const char *name,
292  const cs_matrix_t *a,
293  bool conv_diff,
294  int verbosity);
295 
296 /*----------------------------------------------------------------------------
297  * Call multigrid sparse linear equation solver.
298  *
299  * parameters:
300  * context <-> pointer to iterative sparse linear solver info
301  * (actual type: cs_multigrid_t *)
302  * name <-- pointer to name of linear system
303  * a <-- matrix
304  * verbosity <-- associated verbosity
305  * precision <-- solver precision
306  * r_norm <-- residual normalization
307  * n_iter --> number of iterations
308  * residual --> residual
309  * rhs <-- right hand side
310  * vx <-> system solution
311  * aux_size <-- number of elements in aux_vectors
312  * aux_vectors --- optional working area (internal allocation if NULL)
313  *
314  * returns:
315  * convergence state
316  *----------------------------------------------------------------------------*/
317 
319 cs_multigrid_solve(void *context,
320  const char *name,
321  const cs_matrix_t *a,
322  int verbosity,
323  double precision,
324  double r_norm,
325  int *n_iter,
326  double *residual,
327  const cs_real_t *rhs,
328  cs_real_t *vx,
329  size_t aux_size,
330  void *aux_vectors);
331 
332 /*----------------------------------------------------------------------------
333  * Free iterative sparse linear equation solver setup context.
334  *
335  * Note that this function should free resolution-related data, such as
336  * buffers and preconditioning but does not free the whole context,
337  * as info used for logging (especially performance data) is maintained.
338 
339  * parameters:
340  * context <-> pointer to iterative sparse linear solver info
341  * (actual type: cs_multigrid_t *)
342  *----------------------------------------------------------------------------*/
343 
344 void
345 cs_multigrid_free(void *context);
346 
347 /*----------------------------------------------------------------------------
348  * Log sparse linear equation solver info.
349  *
350  * parameters:
351  * context <-> pointer to iterative sparse linear solver info
352  * (actual type: cs_multigrid_t *)
353  * log_type <-- log type
354  *----------------------------------------------------------------------------*/
355 
356 void
357 cs_multigrid_log(const void *context,
358  cs_log_t log_type);
359 
360 /*----------------------------------------------------------------------------*/
368 /*----------------------------------------------------------------------------*/
369 
370 cs_sles_pc_t *
372 
373 /*----------------------------------------------------------------------------
374  * Error handler for multigrid sparse linear equation solver.
375  *
376  * In case of divergence or breakdown, this error handler outputs
377  * postprocessing data to assist debugging, then aborts the run.
378  * It does nothing in case the maximum iteration count is reached.
379  *
380  * parameters:
381  * sles <-> pointer to solver object
382  * state <-- convergence status
383  * a <-- matrix
384  * rhs <-- right hand side
385  * vx <-> system solution
386  *
387  * returns:
388  * false (do not attempt new solve)
389  *----------------------------------------------------------------------------*/
390 
391 bool
394  const cs_matrix_t *a,
395  const cs_real_t *rhs,
396  cs_real_t *vx);
397 
398 /*----------------------------------------------------------------------------
399  * Set plotting options for multigrid.
400  *
401  * parameters:
402  * mg <-> pointer to multigrid info and context
403  * base_name <-- base plot name to activate, NULL otherwise
404  * use_iteration <-- if true, use iteration as time stamp
405  * otherwise, use wall clock time
406  *----------------------------------------------------------------------------*/
407 
408 void
410  const char *base_name,
411  bool use_iteration);
412 
413 /*----------------------------------------------------------------------------*/
425 /*----------------------------------------------------------------------------*/
426 
427 void
429  int *rank_stride,
430  int *rows_mean_threshold,
431  cs_gnum_t *rows_glob_threshold);
432 
433 /*----------------------------------------------------------------------------*/
445 /*----------------------------------------------------------------------------*/
446 
447 void
449  int rank_stride,
450  int rows_mean_threshold,
451  cs_gnum_t rows_glob_threshold);
452 
453 /*----------------------------------------------------------------------------*/
466 /*----------------------------------------------------------------------------*/
467 
468 const cs_grid_t *
470  int level);
471 
472 /*----------------------------------------------------------------------------*/
473 
475 
476 #endif /* __CS_MULTIGRID_H__ */
#define BEGIN_C_DECLS
Definition: cs_defs.h:514
double cs_real_t
Floating-point value.
Definition: cs_defs.h:319
unsigned long cs_gnum_t
global mesh entity number
Definition: cs_defs.h:298
#define END_C_DECLS
Definition: cs_defs.h:515
struct _cs_grid_t cs_grid_t
Definition: cs_grid.h:67
cs_log_t
Definition: cs_log.h:48
struct _cs_matrix_t cs_matrix_t
Definition: cs_matrix.h:110
void cs_multigrid_set_coarsening_options(cs_multigrid_t *mg, int aggregation_limit, int coarsening_type, int n_max_levels, cs_gnum_t min_g_rows, double p0p1_relax, int postprocess_block_size)
Set multigrid coarsening parameters.
Definition: cs_multigrid.c:3993
void cs_multigrid_setup(void *context, const char *name, const cs_matrix_t *a, int verbosity)
Setup multigrid sparse linear equation solver.
Definition: cs_multigrid.c:4152
void cs_multigrid_setup_conv_diff(void *context, const char *name, const cs_matrix_t *a, bool conv_diff, int verbosity)
Setup multigrid sparse linear equation solver.
Definition: cs_multigrid.c:4178
cs_sles_it_type_t cs_multigrid_get_fine_solver_type(const cs_multigrid_t *mg)
Return solver type used on fine mesh.
Definition: cs_multigrid.c:4129
struct _cs_multigrid_t cs_multigrid_t
Definition: cs_multigrid.h:68
void * cs_multigrid_copy(const void *context)
Create multigrid sparse linear system solver info and context based on existing info and context.
Definition: cs_multigrid.c:3929
cs_multigrid_type_t
Definition: cs_multigrid.h:57
@ CS_MULTIGRID_K_CYCLE
Definition: cs_multigrid.h:60
@ CS_MULTIGRID_N_TYPES
Definition: cs_multigrid.h:62
@ CS_MULTIGRID_K_CYCLE_HPC
Definition: cs_multigrid.h:61
@ CS_MULTIGRID_V_CYCLE
Definition: cs_multigrid.h:59
void cs_multigrid_get_merge_options(const cs_multigrid_t *mg, int *rank_stride, int *rows_mean_threshold, cs_gnum_t *rows_glob_threshold)
Query the global multigrid parameters for parallel grid merging.
Definition: cs_multigrid.c:4783
cs_sles_pc_t * cs_multigrid_pc_create(cs_multigrid_type_t mg_type)
Create a multigrid preconditioner.
Definition: cs_multigrid.c:4495
cs_multigrid_t * cs_multigrid_define(int f_id, const char *name, cs_multigrid_type_t mg_type)
Define and associate a multigrid sparse linear system solver for a given field or equation name.
Definition: cs_multigrid.c:3724
cs_multigrid_t * cs_multigrid_create(cs_multigrid_type_t mg_type)
Create multigrid linear system solver info and context.
Definition: cs_multigrid.c:3761
const cs_grid_t * cs_multigrid_get_grid(const cs_multigrid_t *mg, int level)
Return a pointer to a grid associated with a given multigrid setup and level.
Definition: cs_multigrid.c:4850
void cs_multigrid_set_merge_options(cs_multigrid_t *mg, int rank_stride, int rows_mean_threshold, cs_gnum_t rows_glob_threshold)
Set global multigrid parameters for parallel grid merging behavior.
Definition: cs_multigrid.c:4820
cs_sles_convergence_state_t cs_multigrid_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, size_t aux_size, void *aux_vectors)
Call multigrid sparse linear equation solver.
Definition: cs_multigrid.c:4265
void cs_multigrid_set_plot_options(cs_multigrid_t *mg, const char *base_name, bool use_iteration)
Set plotting options for multigrid.
Definition: cs_multigrid.c:4737
bool cs_multigrid_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)
void cs_multigrid_set_solver_options(cs_multigrid_t *mg, cs_sles_it_type_t descent_smoother_type, cs_sles_it_type_t ascent_smoother_type, cs_sles_it_type_t coarse_solver_type, int n_max_cycles, int n_max_iter_descent, int n_max_iter_ascent, int n_max_iter_coarse, int poly_degree_descent, int poly_degree_ascent, int poly_degree_coarse, double precision_mult_descent, double precision_mult_ascent, double precision_mult_coarse)
Set multigrid parameters for associated iterative solvers.
Definition: cs_multigrid.c:4046
void cs_multigrid_free(void *context)
Free multigrid sparse linear equation solver setup context.
Definition: cs_multigrid.c:4427
void cs_multigrid_initialize(void)
Initialize multigrid solver API.
Definition: cs_multigrid.c:3680
void cs_multigrid_log(const void *context, cs_log_t log_type)
Log multigrid solver info.
Definition: cs_multigrid.c:3961
void cs_multigrid_set_max_cycles(cs_multigrid_t *mg, int n_max_cycles)
Set the max. number of cycles for a multigrid.
Definition: cs_multigrid.c:4107
void cs_multigrid_finalize(void)
Finalize multigrid solver API.
Definition: cs_multigrid.c:3691
const char * cs_multigrid_type_name[]
void cs_multigrid_destroy(void **context)
Destroy multigrid linear system solver info and context.
Definition: cs_multigrid.c:3876
cs_sles_convergence_state_t
Definition: cs_sles.h:56
struct _cs_sles_t cs_sles_t
Definition: cs_sles.h:68
cs_sles_it_type_t
Definition: cs_sles_it.h:55
struct _cs_sles_pc_t cs_sles_pc_t
Definition: cs_sles_pc.h:66