programmer's 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-2015 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_sles.h"
36 #include "cs_sles_it.h"
37 #include "cs_sles_pc.h"
38 #include "cs_time_plot.h"
39 
40 /*----------------------------------------------------------------------------*/
41 
43 
44 /*============================================================================
45  * Macro definitions
46  *============================================================================*/
47 
48 /*============================================================================
49  * Type definitions
50  *============================================================================*/
51 
52 /* Multigrid linear solver context (opaque) */
53 
54 typedef struct _cs_multigrid_t cs_multigrid_t;
55 
56 /*============================================================================
57  * Global variables
58  *============================================================================*/
59 
60 /*=============================================================================
61  * Public function prototypes
62  *============================================================================*/
63 
64 /*----------------------------------------------------------------------------
65  * Initialize multigrid solver API.
66  *----------------------------------------------------------------------------*/
67 
68 void
70 
71 /*----------------------------------------------------------------------------
72  * Finalize multigrid solver API.
73  *----------------------------------------------------------------------------*/
74 
75 void
77 
78 /*----------------------------------------------------------------------------
79  * Indicate if multigrid solver API is used for at least one system.
80  *
81  * returns:
82  * true if at least one system uses a multigrid solver, false otherwise
83  *----------------------------------------------------------------------------*/
84 
85 bool
87 
88 /*----------------------------------------------------------------------------
89  * Define and associate a multigrid sparse linear system solver
90  * for a given field or equation name.
91  *
92  * If this system did not previously exist, it is added to the list of
93  * "known" systems. Otherwise, its definition is replaced by the one
94  * defined here.
95  *
96  * This is a utility function: if finer control is needed, see
97  * cs_sles_define() and cs_multigrid_create().
98  *
99  * Note that this function returns a pointer directly to the multigrid solver
100  * management structure. This may be used to set further options, for
101  * example calling cs_multigrid_set_coarsening_options() and
102  * cs_multigrid_set_solver_options().
103  * If needed, cs_sles_find() may be used to obtain a pointer to the
104  * matching cs_sles_t container.
105  *
106  * parameters:
107  * f_id <-- associated field id, or < 0
108  * name <-- associated name if f_id < 0, or NULL
109  *
110  * \return pointer to new multigrid info and context
111  */
112 /*----------------------------------------------------------------------------*/
113 
115 cs_multigrid_define(int f_id,
116  const char *name);
117 
118 /*----------------------------------------------------------------------------
119  * Create multigrid linear system solver info and context.
120  *
121  * The multigrid variant is an ACM (Additive Corrective Multigrid) method.
122  *
123  * returns:
124  * pointer to new multigrid info and context
125  *----------------------------------------------------------------------------*/
126 
128 cs_multigrid_create(void);
129 
130 /*----------------------------------------------------------------------------
131  * Destroy multigrid linear system solver info and context.
132  *
133  * parameters:
134  * context <-> pointer to multigrid linear solver info
135  * (actual type: cs_multigrid_t **)
136  *----------------------------------------------------------------------------*/
137 
138 void
139 cs_multigrid_destroy(void **context);
140 
141 /*----------------------------------------------------------------------------
142  * Create multigrid sparse linear system solver info and context
143  * based on existing info and context.
144  *
145  * parameters:
146  * context <-- pointer to reference info and context
147  * (actual type: cs_multigrid_t *)
148  *
149  * returns:
150  * pointer to newly created solver info object
151  * (actual type: cs_multigrid_t *)
152  *----------------------------------------------------------------------------*/
153 
154 void *
155 cs_multigrid_copy(const void *context);
156 
157 /*----------------------------------------------------------------------------
158  * Set multigrid coarsening parameters.
159  *
160  * parameters:
161  * mg <-> pointer to multigrid info and context
162  * aggregation_limit <-- maximum allowed fine cells per coarse cell
163  * coarsening_type <-- coarsening type:
164  * 0: algebraic, natural face traversal;
165  * 1: algebraic, face traveral by criteria;
166  * 2: algebraic, Hilbert face traversal;
167  * n_max_levels <-- maximum number of grid levels
168  * min_g_cells <-- global number of cells on coarse grids
169  * under which no coarsening occurs
170  * p0p1_relax <-- p0/p1 relaxation_parameter
171  * postprocess_block_size <-- if > 0, postprocess coarsening
172  * (using coarse cell numbers modulo this value)
173  *----------------------------------------------------------------------------*/
174 
175 void
177  int aggregation_limit,
178  int coarsening_type,
179  int n_max_levels,
180  cs_gnum_t min_g_cells,
181  double p0p1_relax,
182  int postprocess_block_size);
183 
184 /*----------------------------------------------------------------------------
185  * Set multigrid parameters for associated iterative solvers.
186  *
187  * parameters:
188  * mg <-> pointer to multigrid info and context
189  * descent_smoother_type <-- type of smoother for descent
190  * ascent_smoother_type <-- type of smoother for ascent
191  * coarse_solver_type <-- type of solver
192  * n_max_cycles <-- maximum number of cycles
193  * n_max_iter_descent <-- maximum iterations per descent phase
194  * n_max_iter_ascent <-- maximum iterations per descent phase
195  * n_max_iter_coarse <-- maximum iterations per coarsest solution
196  * poly_degree_descent <-- preconditioning polynomial degree
197  * for descent phases (0: diagonal)
198  * poly_degree_ascent <-- preconditioning polynomial degree
199  * for ascent phases (0: diagonal)
200  * poly_degree_coarse <-- preconditioning polynomial degree
201  * for coarse solver (0: diagonal)
202  * precision_mult_descent <-- precision multiplier for descent phases
203  * (levels >= 1)
204  * precision_mult_ascent <-- precision multiplier for ascent phases
205  * precision_mult_coarse <-- precision multiplier for coarsest grid
206  *----------------------------------------------------------------------------*/
207 
208 void
210  cs_sles_it_type_t descent_smoother_type,
211  cs_sles_it_type_t ascent_smoother_type,
212  cs_sles_it_type_t coarse_solver_type,
213  int n_max_cycles,
214  int n_max_iter_descent,
215  int n_max_iter_ascent,
216  int n_max_iter_coarse,
217  int poly_degree_descent,
218  int poly_degree_ascent,
219  int poly_degree_coarse,
220  double precision_mult_descent,
221  double precision_mult_ascent,
222  double precision_mult_coarse);
223 
224 /*----------------------------------------------------------------------------
225  * Return solver type used on fine mesh.
226  *
227  * parameters:
228  * mg <-- pointer to multigrid info and context
229  *
230  * returns:
231  * type of smoother for descent (used for fine mesh)
232  *----------------------------------------------------------------------------*/
233 
236 
237 /*----------------------------------------------------------------------------
238  * Setup multigrid sparse linear equation solver.
239  *
240  * parameters:
241  * context <-> pointer to multigrid info and context
242  * (actual type: cs_multigrid_t *)
243  * name <-- pointer to name of linear system
244  * a <-- associated matrix
245  * verbosity <-- associated verbosity
246  *----------------------------------------------------------------------------*/
247 
248 void
249 cs_multigrid_setup(void *context,
250  const char *name,
251  const cs_matrix_t *a,
252  int verbosity);
253 
254 /*----------------------------------------------------------------------------
255  * Call multigrid sparse linear equation solver.
256  *
257  * parameters:
258  * context <-> pointer to iterative sparse linear solver info
259  * (actual type: cs_multigrid_t *)
260  * name <-- pointer to name of linear system
261  * a <-- matrix
262  * verbosity <-- associated verbosity
263  * rotation_mode <-- halo update option for rotational periodicity
264  * precision <-- solver precision
265  * r_norm <-- residue normalization
266  * n_iter --> number of iterations
267  * residue --> residue
268  * rhs <-- right hand side
269  * vx <-> system solution
270  * aux_size <-- number of elements in aux_vectors
271  * aux_vectors --- optional working area (internal allocation if NULL)
272  *
273  * returns:
274  * convergence state
275  *----------------------------------------------------------------------------*/
276 
278 cs_multigrid_solve(void *context,
279  const char *name,
280  const cs_matrix_t *a,
281  int verbosity,
282  cs_halo_rotation_t rotation_mode,
283  double precision,
284  double r_norm,
285  int *n_iter,
286  double *residue,
287  const cs_real_t *rhs,
288  cs_real_t *vx,
289  size_t aux_size,
290  void *aux_vectors);
291 
292 /*----------------------------------------------------------------------------
293  * Free iterative sparse linear equation solver setup context.
294  *
295  * Note that this function should free resolution-related data, such as
296  * buffers and preconditioning but doesd not free the whole context,
297  * as info used for logging (especially performance data) is maintained.
298 
299  * parameters:
300  * context <-> pointer to iterative sparse linear solver info
301  * (actual type: cs_multigrid_t *)
302  *----------------------------------------------------------------------------*/
303 
304 void
305 cs_multigrid_free(void *context);
306 
307 /*----------------------------------------------------------------------------
308  * Log sparse linear equation solver info.
309  *
310  * parameters:
311  * context <-> pointer to iterative sparse linear solver info
312  * (actual type: cs_multigrid_t *)
313  * log_type <-- log type
314  *----------------------------------------------------------------------------*/
315 
316 void
317 cs_multigrid_log(const void *context,
318  cs_log_t log_type);
319 
320 /*----------------------------------------------------------------------------
321  * Create a multigrid preconditioner.
322  *
323  * returns:
324  * pointer to newly created preconditioner object.
325  *----------------------------------------------------------------------------*/
326 
327 cs_sles_pc_t *
329 
330 /*----------------------------------------------------------------------------
331  * Error handler for multigrid sparse linear equation solver.
332  *
333  * In case of divergence or breakdown, this error handler outputs
334  * postprocessing data to assist debugging, then aborts the run.
335  * It does nothing in case the maximum iteration count is reached.
336  *
337  * parameters:
338  * context <-> pointer to multigrid sparse linear system solver info
339  * (actual type: cs_multigrid_t *)
340  * name <-- pointer to name of linear system
341  * state <-- convergence status
342  * a <-- matrix
343  * rotation_mode <-- halo update option for rotational periodicity
344  * rhs <-- right hand side
345  * vx <-> system solution
346  *----------------------------------------------------------------------------*/
347 
348 void
351  const char *name,
352  const cs_matrix_t *a,
353  cs_halo_rotation_t rotation_mode,
354  const cs_real_t *rhs,
355  cs_real_t *vx);
356 
357 /*----------------------------------------------------------------------------
358  * Set plotting options for multigrid.
359  *
360  * parameters:
361  * mg <-> pointer to multigrid info and context
362  * base_name <-- base plot name to activate, NULL otherwise
363  * use_iteration <-- if true, use iteration as time stamp
364  * otherwise, use wall clock time
365  *----------------------------------------------------------------------------*/
366 
367 void
369  const char *base_name,
370  bool use_iteration);
371 
372 /*----------------------------------------------------------------------------*/
373 
375 
376 #endif /* __CS_MULTIGRID_H__ */
unsigned long cs_gnum_t
global mesh entity number
Definition: cs_defs.h:280
cs_halo_rotation_t
Definition: cs_halo.h:59
cs_sles_it_type_t
Iterative solver types.
Definition: cs_sles_it.h:57
cs_multigrid_t * cs_multigrid_define(int f_id, const char *name)
Define and associate a multigrid sparse linear system solver for a given field or equation name...
Definition: cs_multigrid.c:2143
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:2289
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:3172
cs_sles_pc_t * cs_multigrid_pc_create(void)
Create a multigrid preconditioner.
Definition: cs_multigrid.c:2962
#define BEGIN_C_DECLS
Definition: cs_defs.h:429
void cs_multigrid_destroy(void **context)
Destroy multigrid linear system solver info and context.
Definition: cs_multigrid.c:2233
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:2483
struct _cs_sles_pc_t cs_sles_pc_t
Definition: cs_sles_pc.h:66
bool cs_multigrid_needed(void)
Indicate if multigrid solver API is used for at least one system.
Definition: cs_multigrid.c:2111
void cs_multigrid_log(const void *context, cs_log_t log_type)
Log multigrid solver info.
Definition: cs_multigrid.c:2317
struct _cs_matrix_t cs_matrix_t
Definition: cs_matrix.h:89
void cs_multigrid_error_post_and_abort(void *context, cs_sles_convergence_state_t state, const char *name, const cs_matrix_t *a, cs_halo_rotation_t rotation_mode, const cs_real_t *rhs, cs_real_t *vx)
cs_multigrid_t * cs_multigrid_create(void)
Create multigrid linear system solver info and context.
Definition: cs_multigrid.c:2177
cs_sles_convergence_state_t
Convergence status indicator.
Definition: cs_sles.h:55
void cs_multigrid_finalize(void)
Finalize multigrid solver API.
Definition: cs_multigrid.c:2095
double precision, save a
Definition: cs_fuel_incl.f90:146
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:2763
struct _cs_multigrid_t cs_multigrid_t
Definition: cs_multigrid.h:54
void cs_multigrid_free(void *context)
Free multigrid sparse linear equation solver setup context.
Definition: cs_multigrid.c:2902
cs_log_t
Definition: cs_log.h:47
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:2460
void cs_multigrid_initialize(void)
Initialize multigrid solver API.
Definition: cs_multigrid.c:2084
#define END_C_DECLS
Definition: cs_defs.h:430
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:2351
double cs_real_t
Definition: cs_defs.h:296
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:2404