8.2
general documentation
Loading...
Searching...
No Matches
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-2024 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
65
66/* Multigrid linear solver context (opaque) */
67
68typedef struct _cs_multigrid_t cs_multigrid_t;
69
70/*============================================================================
71 * Global variables
72 *============================================================================*/
73
74/* Names for multigrid types */
75
76extern const char *cs_multigrid_type_name[];
77
78/*=============================================================================
79 * Public function prototypes
80 *============================================================================*/
81
82/*----------------------------------------------------------------------------
83 * Initialize multigrid solver API.
84 *----------------------------------------------------------------------------*/
85
86void
88
89/*----------------------------------------------------------------------------
90 * Finalize multigrid solver API.
91 *----------------------------------------------------------------------------*/
92
93void
95
96/*----------------------------------------------------------------------------*
97 * Set maximum grid level which should run on device (i.e. GPU).
98 *----------------------------------------------------------------------------*/
99
100void
102
103/*----------------------------------------------------------------------------
104 * Define and associate a multigrid sparse linear system solver
105 * for a given field or equation name.
106 *
107 * If this system did not previously exist, it is added to the list of
108 * "known" systems. Otherwise, its definition is replaced by the one
109 * defined here.
110 *
111 * This is a utility function: if finer control is needed, see
112 * cs_sles_define() and cs_multigrid_create().
113 *
114 * Note that this function returns a pointer directly to the multigrid solver
115 * management structure. This may be used to set further options, for
116 * example calling cs_multigrid_set_coarsening_options() and
117 * cs_multigrid_set_solver_options().
118 * If needed, cs_sles_find() may be used to obtain a pointer to the
119 * matching cs_sles_t container.
120 *
121 * \param[in] f_id associated field id, or < 0
122 * \param[in] name associated name if f_id < 0, or NULL
123 * \param[in] mg_type type of multigrid algorithm to use
124 *
125 * \return pointer to new multigrid info and context
126 */
127/*----------------------------------------------------------------------------*/
128
130cs_multigrid_define(int f_id,
131 const char *name,
132 cs_multigrid_type_t mg_type);
133
134/*----------------------------------------------------------------------------*/
144/*----------------------------------------------------------------------------*/
145
148
149/*----------------------------------------------------------------------------
150 * Destroy multigrid linear system solver info and context.
151 *
152 * parameters:
153 * context <-> pointer to multigrid linear solver info
154 * (actual type: cs_multigrid_t **)
155 *----------------------------------------------------------------------------*/
156
157void
158cs_multigrid_destroy(void **context);
159
160/*----------------------------------------------------------------------------
161 * Create multigrid sparse linear system solver info and context
162 * based on existing info and context.
163 *
164 * parameters:
165 * context <-- pointer to reference info and context
166 * (actual type: cs_multigrid_t *)
167 *
168 * returns:
169 * pointer to newly created solver info object
170 * (actual type: cs_multigrid_t *)
171 *----------------------------------------------------------------------------*/
172
173void *
174cs_multigrid_copy(const void *context);
175
176/*----------------------------------------------------------------------------
177 * Set multigrid coarsening parameters.
178 *
179 * parameters:
180 * mg <-> pointer to multigrid info and context
181 * aggregation_limit <-- maximum allowed fine rows per coarse cell
182 * coarsening_type <-- coarsening type; see cs_grid_coarsening_t
183 * n_max_levels <-- maximum number of grid levels
184 * min_g_rows <-- global number of rows on coarse grids
185 * under which no coarsening occurs
186 * p0p1_relax <-- p0/p1 relaxation_parameter
187 * postprocess_block_size <-- if > 0, postprocess coarsening
188 * (using coarse cell numbers modulo this value)
189 *----------------------------------------------------------------------------*/
190
191void
193 (cs_multigrid_t *mg,
194 int aggregation_limit,
195 cs_grid_coarsening_t coarsening_type,
196 int n_max_levels,
197 cs_gnum_t min_g_rows,
198 double p0p1_relax,
199 int postprocess_block_size);
200
201/*----------------------------------------------------------------------------
202 * Set multigrid parameters for associated iterative solvers.
203 *
204 * parameters:
205 * mg <-> pointer to multigrid info and context
206 * descent_smoother_type <-- type of smoother for descent
207 * ascent_smoother_type <-- type of smoother for ascent
208 * coarse_solver_type <-- type of solver
209 * n_max_cycles <-- maximum number of cycles
210 * n_max_iter_descent <-- maximum iterations per descent phase
211 * n_max_iter_ascent <-- maximum iterations per descent phase
212 * n_max_iter_coarse <-- maximum iterations per coarsest solution
213 * poly_degree_descent <-- preconditioning polynomial degree
214 * for descent phases (0: diagonal)
215 * poly_degree_ascent <-- preconditioning polynomial degree
216 * for ascent phases (0: diagonal)
217 * poly_degree_coarse <-- preconditioning polynomial degree
218 * for coarse solver (0: diagonal)
219 * precision_mult_descent <-- precision multiplier for descent phases
220 * (levels >= 1)
221 * precision_mult_ascent <-- precision multiplier for ascent phases
222 * precision_mult_coarse <-- precision multiplier for coarsest grid
223 *----------------------------------------------------------------------------*/
224
225void
227 cs_sles_it_type_t descent_smoother_type,
228 cs_sles_it_type_t ascent_smoother_type,
229 cs_sles_it_type_t coarse_solver_type,
230 int n_max_cycles,
231 int n_max_iter_descent,
232 int n_max_iter_ascent,
233 int n_max_iter_coarse,
234 int poly_degree_descent,
235 int poly_degree_ascent,
236 int poly_degree_coarse,
237 double precision_mult_descent,
238 double precision_mult_ascent,
239 double precision_mult_coarse);
240
241/*----------------------------------------------------------------------------*/
248/*----------------------------------------------------------------------------*/
249
250void
252 int n_max_cycles);
253
254/*----------------------------------------------------------------------------
255 * Return solver type used on fine mesh.
256 *
257 * parameters:
258 * mg <-- pointer to multigrid info and context
259 *
260 * returns:
261 * type of smoother for descent (used for fine mesh)
262 *----------------------------------------------------------------------------*/
263
266
267/*----------------------------------------------------------------------------
268 * Setup multigrid sparse linear equation solver.
269 *
270 * parameters:
271 * context <-> pointer to multigrid info and context
272 * (actual type: cs_multigrid_t *)
273 * name <-- pointer to name of linear system
274 * a <-- associated matrix
275 * verbosity <-- associated verbosity
276 *----------------------------------------------------------------------------*/
277
278void
279cs_multigrid_setup(void *context,
280 const char *name,
281 const cs_matrix_t *a,
282 int verbosity);
283
284/*----------------------------------------------------------------------------
285 * Setup multigrid sparse linear equation solver with separate
286 * convection-diffusion matrixes
287 *
288 * parameters:
289 * context <-> pointer to multigrid info and context
290 * (actual type: cs_multigrid_t *)
291 * name <-- pointer to name of linear system
292 * a <-- associated matrix
293 * conv_diff <-- convection-diffusion mode
294 * verbosity <-- associated verbosity
295 *----------------------------------------------------------------------------*/
296
297void
298cs_multigrid_setup_conv_diff(void *context,
299 const char *name,
300 const cs_matrix_t *a,
301 bool conv_diff,
302 int verbosity);
303
304/*----------------------------------------------------------------------------
305 * Call multigrid sparse linear equation solver.
306 *
307 * parameters:
308 * context <-> pointer to iterative sparse linear solver info
309 * (actual type: cs_multigrid_t *)
310 * name <-- pointer to name of linear system
311 * a <-- matrix
312 * verbosity <-- associated verbosity
313 * precision <-- solver precision
314 * r_norm <-- residual normalization
315 * n_iter --> number of iterations
316 * residual --> residual
317 * rhs <-- right hand side
318 * vx_ini <-- initial solution (vx if nonzero, nullptr if zero)
319 * vx <-> system solution
320 * aux_size <-- number of elements in aux_vectors
321 * aux_vectors --- optional working area (internal allocation if NULL)
322 *
323 * returns:
324 * convergence state
325 *----------------------------------------------------------------------------*/
326
328cs_multigrid_solve(void *context,
329 const char *name,
330 const cs_matrix_t *a,
331 int verbosity,
332 double precision,
333 double r_norm,
334 int *n_iter,
335 double *residual,
336 const cs_real_t *rhs,
337 cs_real_t *vx_ini,
338 cs_real_t *vx,
339 size_t aux_size,
340 void *aux_vectors);
341
342/*----------------------------------------------------------------------------
343 * Free iterative sparse linear equation solver setup context.
344 *
345 * Note that this function should free resolution-related data, such as
346 * buffers and preconditioning but does not free the whole context,
347 * as info used for logging (especially performance data) is maintained.
348
349 * parameters:
350 * context <-> pointer to iterative sparse linear solver info
351 * (actual type: cs_multigrid_t *)
352 *----------------------------------------------------------------------------*/
353
354void
355cs_multigrid_free(void *context);
356
357/*----------------------------------------------------------------------------
358 * Log sparse linear equation solver info.
359 *
360 * parameters:
361 * context <-> pointer to iterative sparse linear solver info
362 * (actual type: cs_multigrid_t *)
363 * log_type <-- log type
364 *----------------------------------------------------------------------------*/
365
366void
367cs_multigrid_log(const void *context,
368 cs_log_t log_type);
369
370/*----------------------------------------------------------------------------*/
378/*----------------------------------------------------------------------------*/
379
382
383/*----------------------------------------------------------------------------
384 * Error handler for multigrid sparse linear equation solver.
385 *
386 * In case of divergence or breakdown, this error handler outputs
387 * postprocessing data to assist debugging, then aborts the run.
388 * It does nothing in case the maximum iteration count is reached.
389 *
390 * parameters:
391 * sles <-> pointer to solver object
392 * state <-- convergence status
393 * a <-- matrix
394 * rhs <-- right hand side
395 * vx <-> system solution
396 *
397 * returns:
398 * false (do not attempt new solve)
399 *----------------------------------------------------------------------------*/
400
401bool
404 const cs_matrix_t *a,
405 const cs_real_t *rhs,
406 cs_real_t *vx);
407
408/*----------------------------------------------------------------------------
409 * Set plotting options for multigrid.
410 *
411 * parameters:
412 * mg <-> pointer to multigrid info and context
413 * base_name <-- base plot name to activate, NULL otherwise
414 * use_iteration <-- if true, use iteration as time stamp
415 * otherwise, use wall clock time
416 *----------------------------------------------------------------------------*/
417
418void
420 const char *base_name,
421 bool use_iteration);
422
423/*----------------------------------------------------------------------------*/
435/*----------------------------------------------------------------------------*/
436
437void
439 int *rank_stride,
440 int *rows_mean_threshold,
441 cs_gnum_t *rows_glob_threshold);
442
443/*----------------------------------------------------------------------------*/
455/*----------------------------------------------------------------------------*/
456
457void
459 int rank_stride,
460 int rows_mean_threshold,
461 cs_gnum_t rows_glob_threshold);
462
463/*----------------------------------------------------------------------------*/
476/*----------------------------------------------------------------------------*/
477
478const cs_grid_t *
480 int level);
481
482/*----------------------------------------------------------------------------*/
483
485
486#endif /* __CS_MULTIGRID_H__ */
#define BEGIN_C_DECLS
Definition cs_defs.h:528
double cs_real_t
Floating-point value.
Definition cs_defs.h:332
#define END_C_DECLS
Definition cs_defs.h:529
cs_grid_coarsening_t
Definition cs_grid.h:55
struct _cs_grid_t cs_grid_t
Definition cs_grid.h:68
cs_log_t
Definition cs_log.h:48
struct _cs_matrix_t cs_matrix_t
Definition cs_matrix.h:110
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.cpp:4775
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.cpp:4801
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.cpp:4752
struct _cs_multigrid_t cs_multigrid_t
Definition cs_multigrid.h:68
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.cpp:5483
void cs_multigrid_set_max_grid_level_for_device(int level)
Set maximum grid level which should run on device (i.e. GPU).
Definition cs_multigrid.cpp:4313
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.cpp:4346
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.cpp:4551
void cs_multigrid_set_coarsening_options(cs_multigrid_t *mg, int aggregation_limit, cs_grid_coarsening_t 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.cpp:4616
cs_sles_pc_t * cs_multigrid_pc_create(cs_multigrid_type_t mg_type)
Create a multigrid preconditioner.
Definition cs_multigrid.cpp:5195
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_ini, cs_real_t *vx, size_t aux_size, void *aux_vectors)
Call multigrid sparse linear equation solver.
Definition cs_multigrid.cpp:4886
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.cpp:5520
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.cpp:5437
cs_multigrid_t * cs_multigrid_create(cs_multigrid_type_t mg_type)
Create multigrid linear system solver info and context.
Definition cs_multigrid.cpp:4383
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)
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.cpp:5550
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.cpp:4669
void cs_multigrid_free(void *context)
Free multigrid sparse linear equation solver setup context.
Definition cs_multigrid.cpp:5125
void cs_multigrid_initialize(void)
Initialize multigrid solver API.
Definition cs_multigrid.cpp:4290
void cs_multigrid_log(const void *context, cs_log_t log_type)
Log multigrid solver info.
Definition cs_multigrid.cpp:4583
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.cpp:4730
void cs_multigrid_finalize(void)
Finalize multigrid solver API.
Definition cs_multigrid.cpp:4301
const char * cs_multigrid_type_name[]
void cs_multigrid_destroy(void **context)
Destroy multigrid linear system solver info and context.
Definition cs_multigrid.cpp:4498
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