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-2021 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 
43 BEGIN_C_DECLS
44 
45 /*============================================================================
46  * Macro definitions
47  *============================================================================*/
48 
49 /*============================================================================
50  * Type definitions
51  *============================================================================*/
52 
53 /*----------------------------------------------------------------------------
54  * Multigrid types
55  *----------------------------------------------------------------------------*/
56 
57 typedef enum {
58 
59   CS_MULTIGRID_V_CYCLE,          /*!< Use a V-cycle */
60   CS_MULTIGRID_K_CYCLE,          /*!< Use a K-cycle */
61   CS_MULTIGRID_K_CYCLE_HPC,      /*!< Use a K-cycle tuned for HPC systems */
62   CS_MULTIGRID_N_TYPES           /*!< Number of multigrid types */
63 
64 } cs_multigrid_type_t;
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
87 cs_multigrid_initialize(void);
88 
89 /*----------------------------------------------------------------------------
90  * Finalize multigrid solver API.
91  *----------------------------------------------------------------------------*/
92 
93 void
94 cs_multigrid_finalize(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 
122 cs_multigrid_t *
123 cs_multigrid_define(int                   f_id,
124                     const char           *name,
125                     cs_multigrid_type_t   mg_type);
126 
127 /*----------------------------------------------------------------------------*/
128 /*!
129  * \brief Create multigrid linear system solver info and context.
130  *
131  * The multigrid variant is an ACM (Additive Corrective Multigrid) method.
132  *
133  * \param[in]  mg_type  type of multigrid algorithm to use
134  *
135  * \return  pointer to new multigrid info and context
136  */
137 /*----------------------------------------------------------------------------*/
138 
139 cs_multigrid_t *
140 cs_multigrid_create(cs_multigrid_type_t  mg_type);
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
185 cs_multigrid_set_coarsening_options(cs_multigrid_t  *mg,
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
218 cs_multigrid_set_solver_options(cs_multigrid_t     *mg,
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 /*----------------------------------------------------------------------------
234  * Return solver type used on fine mesh.
235  *
236  * parameters:
237  *   mg <-- pointer to multigrid info and context
238  *
239  * returns:
240  *   type of smoother for descent (used for fine mesh)
241  *----------------------------------------------------------------------------*/
242 
243 cs_sles_it_type_t
244 cs_multigrid_get_fine_solver_type(const cs_multigrid_t  *mg);
245 
246 /*----------------------------------------------------------------------------
247  * Setup multigrid sparse linear equation solver.
248  *
249  * parameters:
250  *   context   <-> pointer to multigrid info and context
251  *                 (actual type: cs_multigrid_t  *)
252  *   name      <-- pointer to name of linear system
253  *   a         <-- associated matrix
254  *   verbosity <-- associated verbosity
255  *----------------------------------------------------------------------------*/
256 
257 void
258 cs_multigrid_setup(void               *context,
259                    const char         *name,
260                    const cs_matrix_t  *a,
261                    int                 verbosity);
262 
263 /*----------------------------------------------------------------------------
264  * Setup multigrid sparse linear equation solver with separate
265  * convection-diffusion matrixes
266  *
267  * parameters:
268  *   context   <-> pointer to multigrid info and context
269  *                 (actual type: cs_multigrid_t  *)
270  *   name      <-- pointer to name of linear system
271  *   a         <-- associated matrix
272  *   a_conv    <-- associated matrix (convection)
273  *   a_diff    <-- associated matrix (diffusion)
274  *   verbosity <-- associated verbosity
275  *----------------------------------------------------------------------------*/
276 
277 void
278 cs_multigrid_setup_conv_diff(void               *context,
279                              const char         *name,
280                              const cs_matrix_t  *a,
281                              const cs_matrix_t  *a_conv,
282                              const cs_matrix_t  *a_diff,
283                              int                 verbosity);
284 
285 /*----------------------------------------------------------------------------
286  * Call multigrid sparse linear equation solver.
287  *
288  * parameters:
289  *   context       <-> pointer to iterative sparse linear solver info
290  *                     (actual type: cs_multigrid_t  *)
291  *   name          <-- pointer to name of linear system
292  *   a             <-- matrix
293  *   verbosity     <-- associated verbosity
294  *   precision     <-- solver precision
295  *   r_norm        <-- residue normalization
296  *   n_iter        --> number of iterations
297  *   residue       --> residue
298  *   rhs           <-- right hand side
299  *   vx            <-> system solution
300  *   aux_size      <-- number of elements in aux_vectors
301  *   aux_vectors   --- optional working area (internal allocation if NULL)
302  *
303  * returns:
304  *   convergence state
305  *----------------------------------------------------------------------------*/
306 
307 cs_sles_convergence_state_t
308 cs_multigrid_solve(void                *context,
309                    const char          *name,
310                    const cs_matrix_t   *a,
311                    int                  verbosity,
312                    double               precision,
313                    double               r_norm,
314                    int                 *n_iter,
315                    double              *residue,
316                    const cs_real_t     *rhs,
317                    cs_real_t           *vx,
318                    size_t               aux_size,
319                    void                *aux_vectors);
320 
321 /*----------------------------------------------------------------------------
322  * Free iterative sparse linear equation solver setup context.
323  *
324  * Note that this function should free resolution-related data, such as
325  * buffers and preconditioning but does not free the whole context,
326  * as info used for logging (especially performance data) is maintained.
327 
328  * parameters:
329  *   context <-> pointer to iterative sparse linear solver info
330  *               (actual type: cs_multigrid_t  *)
331  *----------------------------------------------------------------------------*/
332 
333 void
334 cs_multigrid_free(void  *context);
335 
336 /*----------------------------------------------------------------------------
337  * Log sparse linear equation solver info.
338  *
339  * parameters:
340  *   context  <-> pointer to iterative sparse linear solver info
341  *                (actual type: cs_multigrid_t  *)
342  *   log_type <-- log type
343  *----------------------------------------------------------------------------*/
344 
345 void
346 cs_multigrid_log(const void  *context,
347                  cs_log_t     log_type);
348 
349 /*----------------------------------------------------------------------------*/
350 /*!
351  * \brief Create a multigrid preconditioner.
352  *
353  * \param[in]  mg_type  type of multigrid algorithm to use
354  *
355  * \return  pointer to newly created preconditioner object.
356  */
357 /*----------------------------------------------------------------------------*/
358 
359 cs_sles_pc_t *
360 cs_multigrid_pc_create(cs_multigrid_type_t  mg_type);
361 
362 /*----------------------------------------------------------------------------
363  * Error handler for multigrid sparse linear equation solver.
364  *
365  * In case of divergence or breakdown, this error handler outputs
366  * postprocessing data to assist debugging, then aborts the run.
367  * It does nothing in case the maximum iteration count is reached.
368  *
369  * parameters:
370  *   sles          <-> pointer to solver object
371  *   state         <-- convergence status
372  *   a             <-- matrix
373  *   rhs           <-- right hand side
374  *   vx            <-> system solution
375  *
376  * returns:
377  *   false (do not attempt new solve)
378  *----------------------------------------------------------------------------*/
379 
380 bool
381 cs_multigrid_error_post_and_abort(cs_sles_t                    *sles,
382                                   cs_sles_convergence_state_t   state,
383                                   const cs_matrix_t            *a,
384                                   const cs_real_t              *rhs,
385                                   cs_real_t                    *vx);
386 
387 /*----------------------------------------------------------------------------
388  * Set plotting options for multigrid.
389  *
390  * parameters:
391  *   mg            <-> pointer to multigrid info and context
392  *   base_name     <-- base plot name to activate, NULL otherwise
393  *   use_iteration <-- if true, use iteration as time stamp
394  *                     otherwise, use wall clock time
395  *----------------------------------------------------------------------------*/
396 
397 void
398 cs_multigrid_set_plot_options(cs_multigrid_t  *mg,
399                               const char      *base_name,
400                               bool             use_iteration);
401 
402 /*----------------------------------------------------------------------------*/
403 /*!
404  * \brief Query the global multigrid parameters for parallel grid merging.
405  *
406  * \param[in]   mg                   pointer to multigrid info and context
407  * \param[out]  rank_stride          number of ranks over which merging
408  *                                   takes place, or NULL
409  * \param[out]  rows_mean_threshold  mean number of rows under which merging
410  *                                   should be applied, or NULL
411  * \param[out]  rows_glob_threshold  global number of rows under which
412  *                                   merging should be applied, or NULL
413  */
414 /*----------------------------------------------------------------------------*/
415 
416 void
417 cs_multigrid_get_merge_options(const cs_multigrid_t  *mg,
418                                int                   *rank_stride,
419                                int                   *rows_mean_threshold,
420                                cs_gnum_t             *rows_glob_threshold);
421 
422 /*----------------------------------------------------------------------------*/
423 /*!
424  * \brief Set global multigrid parameters for parallel grid merging behavior.
425  *
426  * \param[in, out]  mg                   pointer to multigrid info and context
427  * \param[in]       rank_stride          number of ranks over which merging
428  *                                       takes place
429  * \param[in]       rows_mean_threshold  mean number of rows under which
430  *                                       merging should be applied
431  * \param[in]       rows_glob_threshold  global number of rows under which
432  *                                       merging should be applied
433  */
434 /*----------------------------------------------------------------------------*/
435 
436 void
437 cs_multigrid_set_merge_options(cs_multigrid_t  *mg,
438                                int              rank_stride,
439                                int              rows_mean_threshold,
440                                cs_gnum_t        rows_glob_threshold);
441 
442 /*----------------------------------------------------------------------------*/
443 
444 END_C_DECLS
445 
446 #endif /* __CS_MULTIGRID_H__ */
447