My Project
programmer's 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-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_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
64
65/* Multigrid linear solver context (opaque) */
66
67typedef struct _cs_multigrid_t cs_multigrid_t;
68
69/*============================================================================
70 * Global variables
71 *============================================================================*/
72
73/* Names for multigrid types */
74
75extern const char *cs_multigrid_type_name[];
76
77/*=============================================================================
78 * Public function prototypes
79 *============================================================================*/
80
81/*----------------------------------------------------------------------------
82 * Initialize multigrid solver API.
83 *----------------------------------------------------------------------------*/
84
85void
87
88/*----------------------------------------------------------------------------
89 * Finalize multigrid solver API.
90 *----------------------------------------------------------------------------*/
91
92void
94
95/*----------------------------------------------------------------------------
96 * Indicate if multigrid solver API is used for at least one system.
97 *
98 * returns:
99 * true if at least one system uses a multigrid solver, false otherwise
100 *----------------------------------------------------------------------------*/
101
102bool
104
105/*----------------------------------------------------------------------------
106 * Define and associate a multigrid sparse linear system solver
107 * for a given field or equation name.
108 *
109 * If this system did not previously exist, it is added to the list of
110 * "known" systems. Otherwise, its definition is replaced by the one
111 * defined here.
112 *
113 * This is a utility function: if finer control is needed, see
114 * cs_sles_define() and cs_multigrid_create().
115 *
116 * Note that this function returns a pointer directly to the multigrid solver
117 * management structure. This may be used to set further options, for
118 * example calling cs_multigrid_set_coarsening_options() and
119 * cs_multigrid_set_solver_options().
120 * If needed, cs_sles_find() may be used to obtain a pointer to the
121 * matching cs_sles_t container.
122 *
123 * \param[in] f_id associated field id, or < 0
124 * \param[in] name associated name if f_id < 0, or NULL
125 * \param[in] mg_type type of multigrid algorithm to use
126 *
127 * \return pointer to new multigrid info and context
128 */
129/*----------------------------------------------------------------------------*/
130
132cs_multigrid_define(int f_id,
133 const char *name,
134 cs_multigrid_type_t mg_type);
135
136/*----------------------------------------------------------------------------*/
146/*----------------------------------------------------------------------------*/
147
150
151/*----------------------------------------------------------------------------
152 * Destroy multigrid linear system solver info and context.
153 *
154 * parameters:
155 * context <-> pointer to multigrid linear solver info
156 * (actual type: cs_multigrid_t **)
157 *----------------------------------------------------------------------------*/
158
159void
160cs_multigrid_destroy(void **context);
161
162/*----------------------------------------------------------------------------
163 * Create multigrid sparse linear system solver info and context
164 * based on existing info and context.
165 *
166 * parameters:
167 * context <-- pointer to reference info and context
168 * (actual type: cs_multigrid_t *)
169 *
170 * returns:
171 * pointer to newly created solver info object
172 * (actual type: cs_multigrid_t *)
173 *----------------------------------------------------------------------------*/
174
175void *
176cs_multigrid_copy(const void *context);
177
178/*----------------------------------------------------------------------------
179 * Set multigrid coarsening parameters.
180 *
181 * parameters:
182 * mg <-> pointer to multigrid info and context
183 * aggregation_limit <-- maximum allowed fine cells per coarse cell
184 * coarsening_type <-- coarsening type:
185 * 0: algebraic, natural face traversal;
186 * 1: algebraic, face traveral by criteria;
187 * 2: algebraic, Hilbert face traversal;
188 * n_max_levels <-- maximum number of grid levels
189 * min_g_cells <-- global number of cells on coarse grids
190 * under which no coarsening occurs
191 * p0p1_relax <-- p0/p1 relaxation_parameter
192 * postprocess_block_size <-- if > 0, postprocess coarsening
193 * (using coarse cell numbers modulo this value)
194 *----------------------------------------------------------------------------*/
195
196void
198 int aggregation_limit,
199 int coarsening_type,
200 int n_max_levels,
201 cs_gnum_t min_g_cells,
202 double p0p1_relax,
203 int postprocess_block_size);
204
205/*----------------------------------------------------------------------------
206 * Set multigrid parameters for associated iterative solvers.
207 *
208 * parameters:
209 * mg <-> pointer to multigrid info and context
210 * descent_smoother_type <-- type of smoother for descent
211 * ascent_smoother_type <-- type of smoother for ascent
212 * coarse_solver_type <-- type of solver
213 * n_max_cycles <-- maximum number of cycles
214 * n_max_iter_descent <-- maximum iterations per descent phase
215 * n_max_iter_ascent <-- maximum iterations per descent phase
216 * n_max_iter_coarse <-- maximum iterations per coarsest solution
217 * poly_degree_descent <-- preconditioning polynomial degree
218 * for descent phases (0: diagonal)
219 * poly_degree_ascent <-- preconditioning polynomial degree
220 * for ascent phases (0: diagonal)
221 * poly_degree_coarse <-- preconditioning polynomial degree
222 * for coarse solver (0: diagonal)
223 * precision_mult_descent <-- precision multiplier for descent phases
224 * (levels >= 1)
225 * precision_mult_ascent <-- precision multiplier for ascent phases
226 * precision_mult_coarse <-- precision multiplier for coarsest grid
227 *----------------------------------------------------------------------------*/
228
229void
231 cs_sles_it_type_t descent_smoother_type,
232 cs_sles_it_type_t ascent_smoother_type,
233 cs_sles_it_type_t coarse_solver_type,
234 int n_max_cycles,
235 int n_max_iter_descent,
236 int n_max_iter_ascent,
237 int n_max_iter_coarse,
238 int poly_degree_descent,
239 int poly_degree_ascent,
240 int poly_degree_coarse,
241 double precision_mult_descent,
242 double precision_mult_ascent,
243 double precision_mult_coarse);
244
245/*----------------------------------------------------------------------------
246 * Return solver type used on fine mesh.
247 *
248 * parameters:
249 * mg <-- pointer to multigrid info and context
250 *
251 * returns:
252 * type of smoother for descent (used for fine mesh)
253 *----------------------------------------------------------------------------*/
254
257
258/*----------------------------------------------------------------------------
259 * Setup multigrid sparse linear equation solver.
260 *
261 * parameters:
262 * context <-> pointer to multigrid info and context
263 * (actual type: cs_multigrid_t *)
264 * name <-- pointer to name of linear system
265 * a <-- associated matrix
266 * verbosity <-- associated verbosity
267 *----------------------------------------------------------------------------*/
268
269void
270cs_multigrid_setup(void *context,
271 const char *name,
272 const cs_matrix_t *a,
273 int verbosity);
274
275/*----------------------------------------------------------------------------
276 * Setup multigrid sparse linear equation solver with seperate
277 * convection-diffusion matrixes
278 *
279 * parameters:
280 * context <-> pointer to multigrid info and context
281 * (actual type: cs_multigrid_t *)
282 * name <-- pointer to name of linear system
283 * a <-- associated matrix
284 * a_conv <-- associated matrix (convection)
285 * a_diff <-- associated matrix (diffusion)
286 * verbosity <-- associated verbosity
287 *----------------------------------------------------------------------------*/
288
289void
290cs_multigrid_setup_conv_diff(void *context,
291 const char *name,
292 const cs_matrix_t *a,
293 const cs_matrix_t *a_conv,
294 const cs_matrix_t *a_diff,
295 int verbosity);
296
297/*----------------------------------------------------------------------------
298 * Call multigrid sparse linear equation solver.
299 *
300 * parameters:
301 * context <-> pointer to iterative sparse linear solver info
302 * (actual type: cs_multigrid_t *)
303 * name <-- pointer to name of linear system
304 * a <-- matrix
305 * verbosity <-- associated verbosity
306 * rotation_mode <-- halo update option for rotational periodicity
307 * precision <-- solver precision
308 * r_norm <-- residue normalization
309 * n_iter --> number of iterations
310 * residue --> residue
311 * rhs <-- right hand side
312 * vx <-> system solution
313 * aux_size <-- number of elements in aux_vectors
314 * aux_vectors --- optional working area (internal allocation if NULL)
315 *
316 * returns:
317 * convergence state
318 *----------------------------------------------------------------------------*/
319
321cs_multigrid_solve(void *context,
322 const char *name,
323 const cs_matrix_t *a,
324 int verbosity,
325 cs_halo_rotation_t rotation_mode,
326 double precision,
327 double r_norm,
328 int *n_iter,
329 double *residue,
330 const cs_real_t *rhs,
331 cs_real_t *vx,
332 size_t aux_size,
333 void *aux_vectors);
334
335/*----------------------------------------------------------------------------
336 * Free iterative sparse linear equation solver setup context.
337 *
338 * Note that this function should free resolution-related data, such as
339 * buffers and preconditioning but doesd not free the whole context,
340 * as info used for logging (especially performance data) is maintained.
341
342 * parameters:
343 * context <-> pointer to iterative sparse linear solver info
344 * (actual type: cs_multigrid_t *)
345 *----------------------------------------------------------------------------*/
346
347void
348cs_multigrid_free(void *context);
349
350/*----------------------------------------------------------------------------
351 * Log sparse linear equation solver info.
352 *
353 * parameters:
354 * context <-> pointer to iterative sparse linear solver info
355 * (actual type: cs_multigrid_t *)
356 * log_type <-- log type
357 *----------------------------------------------------------------------------*/
358
359void
360cs_multigrid_log(const void *context,
361 cs_log_t log_type);
362
363/*----------------------------------------------------------------------------*/
371/*----------------------------------------------------------------------------*/
372
375
376/*----------------------------------------------------------------------------
377 * Error handler for multigrid sparse linear equation solver.
378 *
379 * In case of divergence or breakdown, this error handler outputs
380 * postprocessing data to assist debugging, then aborts the run.
381 * It does nothing in case the maximum iteration count is reached.
382 *
383 * parameters:
384 * sles <-> pointer to solver object
385 * state <-- convergence status
386 * a <-- matrix
387 * rotation_mode <-- halo update option for rotational periodicity
388 * rhs <-- right hand side
389 * vx <-> system solution
390 *
391 * returns:
392 * false (do not attempt new solve)
393 *----------------------------------------------------------------------------*/
394
395bool
398 const cs_matrix_t *a,
399 cs_halo_rotation_t rotation_mode,
400 const cs_real_t *rhs,
401 cs_real_t *vx);
402
403/*----------------------------------------------------------------------------
404 * Set plotting options for multigrid.
405 *
406 * parameters:
407 * mg <-> pointer to multigrid info and context
408 * base_name <-- base plot name to activate, NULL otherwise
409 * use_iteration <-- if true, use iteration as time stamp
410 * otherwise, use wall clock time
411 *----------------------------------------------------------------------------*/
412
413void
415 const char *base_name,
416 bool use_iteration);
417
418/*----------------------------------------------------------------------------*/
419
421
422#endif /* __CS_MULTIGRID_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
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
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:3235
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:3212
struct _cs_multigrid_t cs_multigrid_t
Definition cs_multigrid.h:67
cs_sles_convergence_state_t cs_multigrid_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 multigrid sparse linear equation solver.
Definition cs_multigrid.c:3580
void cs_multigrid_setup_conv_diff(void *context, const char *name, const cs_matrix_t *a, const cs_matrix_t *a_conv, const cs_matrix_t *a_diff, int verbosity)
Setup multigrid sparse linear equation solver.
Definition cs_multigrid.c:3263
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:61
@ CS_MULTIGRID_V_CYCLE
Definition cs_multigrid.h:59
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:2841
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:3037
cs_sles_pc_t * cs_multigrid_pc_create(cs_multigrid_type_t mg_type)
Create a multigrid preconditioner.
Definition cs_multigrid.c:3804
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:4051
cs_multigrid_t * cs_multigrid_create(cs_multigrid_type_t mg_type)
Create multigrid linear system solver info and context.
Definition cs_multigrid.c:2878
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:3150
void cs_multigrid_free(void *context)
Free multigrid sparse linear equation solver setup context.
Definition cs_multigrid.c:3742
bool cs_multigrid_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)
void cs_multigrid_initialize(void)
Initialize multigrid solver API.
Definition cs_multigrid.c:2781
void cs_multigrid_log(const void *context, cs_log_t log_type)
Log multigrid solver info.
Definition cs_multigrid.c:3065
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_cells, double p0p1_relax, int postprocess_block_size)
Set multigrid coarsening parameters.
Definition cs_multigrid.c:3097
bool cs_multigrid_needed(void)
Indicate if multigrid solver API is used for at least one system.
Definition cs_multigrid.c:2808
void cs_multigrid_finalize(void)
Finalize multigrid solver API.
Definition cs_multigrid.c:2792
const char * cs_multigrid_type_name[]
void cs_multigrid_destroy(void **context)
Destroy multigrid linear system solver info and context.
Definition cs_multigrid.c:2981
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:57
struct _cs_sles_pc_t cs_sles_pc_t
Definition cs_sles_pc.h:66