1 #ifndef __CS_NAVSTO_PARAM_H__
2 #define __CS_NAVSTO_PARAM_H__
3 
4 /*============================================================================
5  * Functions to handle cs_navsto_param_t structure
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  *  Local headers
30  *----------------------------------------------------------------------------*/
31 
32 #include "cs_boundary.h"
33 #include "cs_cdo_turbulence.h"
34 #include "cs_equation_param.h"
35 #include "cs_iter_algo.h"
36 #include "cs_math.h"
37 #include "cs_param_sles.h"
38 #include "cs_physical_constants.h"
39 
40 /*----------------------------------------------------------------------------*/
41 
42 BEGIN_C_DECLS
43 
44 /*============================================================================
45  * Macro definitions
46  *============================================================================*/
47 
48 /* Manage the naming of properties, variables and equations related to the
49  * Navier-Stokes module
50  */
51 
52 #define CS_NAVSTO_STREAM_EQNAME      "streamfunction_eq"
53 
54 /*============================================================================
55  * Type definitions
56  *============================================================================*/
57 
58 typedef cs_flag_t  cs_navsto_param_model_flag_t;
59 typedef cs_flag_t  cs_navsto_param_post_flag_t;
60 
61 /*! \enum cs_navsto_param_model_t
62  *  \brief Describe the system of equations related to the Navier-Stokes to be
63  *  solved
64  *
65  * \var CS_NAVSTO_MODEL_STOKES
66  * Stokes equations (mass and momentum) with the classical choice of variables
67  * i.e. velocity and pressure. Mass density is assumed to be constant.
68  *
69  * \var CS_NAVSTO_MODEL_INCOMPRESSIBLE_NAVIER_STOKES
70  * Navier-Stokes equations: mass and momentum with a constant mass density
71  * Mass equation is equivalent to an incompressibility constraint
72  *
73  * \var CS_NAVSTO_MODEL_OSEEN
74  * Like the incompressible Navier-Stokes equations (mass and momentum) but with
75  * a velocity field which is given. Thus the advection term in the momentum
76  * equation is linear. Unknowns: velocity and pressure. Mass density is assumed
77  * to be constant. The advection field is set using
78  * \ref cs_navsto_add_oseen_field
79  *
80  */
81 
82 typedef enum {
83 
84   CS_NAVSTO_MODEL_STOKES,
85   CS_NAVSTO_MODEL_OSEEN,
86   CS_NAVSTO_MODEL_INCOMPRESSIBLE_NAVIER_STOKES,
87 
88   CS_NAVSTO_N_MODELS
89 
90 } cs_navsto_param_model_t;
91 
92 /*! \enum cs_navsto_param_model_bit_t
93  *  \brief Bit values for additional physical modelling related to the
94  *  Navier-Stokes system of equations
95  *
96  * \var CS_NAVSTO_MODEL_STEADY
97  * There is no time dependency in the model
98  *
99  * \var CS_NAVSTO_MODEL_GRAVITY_EFFECTS
100  * Take into account the gravity effects (add a constant source term equal to
101  * rho*vect(g))
102  *
103  * \var CS_NAVSTO_MODEL_CORIOLIS_EFFECTS
104  * Take into account the Coriolis effects (add a source term)
105  *
106  * \var CS_NAVSTO_MODEL_PASSIVE_THERMAL_TRACER
107  * An additional equation is created involving the thermal equation.
108  * The advection field is automatically set as the mass flux. By default,
109  * the temperature is solved in Celsius and homogeneous Neumann boundary
110  * conditions (no flux) are set. This option is not compatible with the option
111  * \ref CS_NAVSTO_MODEL_BOUSSINESQ
112  *
113  * \var CS_NAVSTO_MODEL_BOUSSINESQ
114  * An additional equation is created involving the thermal equation. The
115  * advection field is automatically set as the mass flux. By default, the
116  * temperature is solved in Celsius and homogeneous Neumann boundary conditions
117  * (no flux) are set.  The variation of mass density around a reference value is
118  * related to the variation of temperature w.r.t. a reference temperature and it
119  * plays a role in the gravity effects (rho.vect(g)) The gradient of temperature
120  * is assumed to have a small norm and the mass density variates in a small
121  * range. In this case, an additional momentum source term is added.
122  *
123  * \var CS_NAVSTO_MODEL_WITH_SOLIDIFICATION
124  * The Navier-Stokes is modified to take into account solidification process.
125  * A boussinesq term is added as well as a head loss term derived from a Darcy
126  * approximation. Please see the \ref cs_solidification_model_t structure and
127  * related structures for more details
128  *
129  */
130 
131 typedef enum {
132 
133   CS_NAVSTO_MODEL_STEADY                          = 1<<0, /* =   1 */
134   CS_NAVSTO_MODEL_GRAVITY_EFFECTS                 = 1<<1, /* =   2 */
135   CS_NAVSTO_MODEL_CORIOLIS_EFFECTS                = 1<<2, /* =   4 */
136   CS_NAVSTO_MODEL_PASSIVE_THERMAL_TRACER          = 1<<3, /* =   8 */
137   CS_NAVSTO_MODEL_BOUSSINESQ                      = 1<<4, /* =  16 */
138   CS_NAVSTO_MODEL_WITH_SOLIDIFICATION             = 1<<5  /* =  32 */
139 
140 } cs_navsto_param_model_bit_t;
141 
142 /*! \enum cs_navsto_param_post_bit_t
143  *  \brief Bit values for additional generic postprocessing related to the
144  *  Navier-Stokes module. In what follows, w denotes the vorticity vector, u
145  *  the velocity vector and k the kinetic energy defined by 1/2 * u \cdot u
146  *
147  * \var CS_NAVSTO_POST_VELOCITY_DIVERGENCE
148  * Compute div(u) and associate a field to this quantity
149  *
150  * \var CS_NAVSTO_POST_KINETIC_ENERGY
151  * Compute rho/2 u\cdot u and associate a field to this quantity
152  *
153  * \var CS_NAVSTO_POST_VORTICITY
154  * Compute w = curl(u) and associate a field to this quantity
155  *
156  * \var CS_NAVSTO_POST_VELOCITY_GRADIENT
157  * Compute the tensor grad(u) and associate a field to this quantity
158  *
159  * \var CS_NAVSTO_POST_STREAM_FUNCTION
160  * Add an equation to compute the stream function and associate a field to this
161  * quantity. This equation is a scalar valued diffusion equation:
162  * -Lap(Psi) = w_z
163  * This is only useful for a 2D computation (assuming that the z-azis is the
164  * extruded one, i.e. the flow is in the x-y plane
165  *
166  * \var CS_NAVSTO_POST_HELICITY
167  * Compute h = u \cdot w and associate a field to this quantity
168  *
169  * \var CS_NAVSTO_POST_ENSTROPHY
170  * Compute  w \cdot w and associate a field to this quantity
171  *
172  * \var CS_NAVSTO_POST_MASS_DENSITY
173  * Compute the mass density in each cell and monitor the evolution of the
174  * integral of the mass in the full computational domain. If a Boussinesq
175  * approximation is set, then the mass density used in the buoyancy term is
176  * considered (not the reference value assumed to be constant).
177  *
178  * \var CS_NAVSTO_POST_CELL_MASS_FLUX_BALANCE
179  * Compute the balance of the mass flux for each cell. When the mass density is
180  * constant, this quantity is only a scaling (by the mass density) of the
181  * quantity computed with the flag CS_NAVSTO_POST_VELOCITY_DIVERGENCE.
182  *
183  * \var CS_NAVSTO_POST_PRESSURE_GRADIENT
184  * Compute the presure gradient in each cell.
185  */
186 
187 typedef enum {
188 
189   CS_NAVSTO_POST_VELOCITY_DIVERGENCE      = 1<< 0, /* =   1 */
190   CS_NAVSTO_POST_KINETIC_ENERGY           = 1<< 1, /* =   2 */
191   CS_NAVSTO_POST_VORTICITY                = 1<< 2, /* =   4 */
192   CS_NAVSTO_POST_VELOCITY_GRADIENT        = 1<< 3, /* =   8 */
193   CS_NAVSTO_POST_STREAM_FUNCTION          = 1<< 4, /* =  16 */
194   CS_NAVSTO_POST_HELICITY                 = 1<< 5, /* =  32 */
195   CS_NAVSTO_POST_ENSTROPHY                = 1<< 6, /* =  64 */
196   CS_NAVSTO_POST_MASS_DENSITY             = 1<< 7, /* = 128 */
197   CS_NAVSTO_POST_CELL_MASS_FLUX_BALANCE   = 1<< 8, /* = 256 */
198   CS_NAVSTO_POST_PRESSURE_GRADIENT        = 1<< 9, /* = 512 */
199 
200 } cs_navsto_param_post_bit_t;
201 
202 /*! \enum cs_navsto_sles_t
203  *
204  *  \brief High-level information about the way of settings the SLES for solving
205  *  the Navier-Stokes system. When the system is treated as a saddle-point
206  *  problem (monolithic approach in what follows), then one uses these
207  *  notations: A_{00} is the upper-left block and A_{11} (should be 0 but the
208  *  preconditioner may have entries for the approximation of the inverse of the
209  *  Schur complement).
210  *
211  *
212  * \var CS_NAVSTO_SLES_ADDITIVE_GMRES_BY_BLOCK
213  * Associated keyword: "additive_gmres"
214  *
215  * Available choice when a monolithic approach is used (i.e. with the parameter
216  * CS_NAVSTO_COUPLING_MONOLITHIC is set as coupling algorithm) The
217  * Navier-Stokes system of equations is solved an additive preconditioner
218  * (block diagonal matrix where the block 00 is A_{00}) and the block 11 is set
219  * to the identity.  Preconditioner/solver for the block 00 is set using the
220  * momentum equation.  This option is only available with the support to the
221  * PETSc library up to now.
222  *
223  *
224  * \var CS_NAVSTO_SLES_BLOCK_MULTIGRID_CG
225  * Associated keyword: "block_amg_cg"
226  *
227  * The Navier-Stokes system of equations is solved using a multigrid on each
228  * diagonal block as a preconditioner and applying a conjugate gradient as
229  * solver. Use this strategy when the saddle-point problem has been reformulated
230  * into a "classical" linear system. For instance when a Uzawa or an Artificial
231  * Compressibility coupling algorithm is used. (i.e. with the parameter
232  * \ref CS_NAVSTO_COUPLING_ARTIFICIAL_COMPRESSIBILITY is set as coupling
233  * algorithm). This option is only available with the support to the PETSc
234  * library up to now.
235  *
236  *
237  * \var CS_NAVSTO_SLES_DIAG_SCHUR_GMRES
238  * Associated keyword: "diag_schur_gmres"
239  *
240  * Available choice when a monolithic approach is used (i.e. with the parameter
241  * CS_NAVSTO_COUPLING_MONOLITHIC is set as coupling algorithm). The
242  * Navier-Stokes system of equations is solved using a block diagonal
243  * preconditioner where the block 00 is A_{00} preconditioned with one multigrid
244  * iteration and the block 11 is an approximation of the Schur complement
245  * preconditioned with one multigrid iteration. The main iterative solver is a
246  * flexible GMRES. This option is only available with the support to the PETSc
247  * library up to now.
248  *
249  *
250  * \var CS_NAVSTO_SLES_DIAG_SCHUR_GCR
251  * Associated keyword: "diag_schur_gcr"
252  *
253  * The Stokes or Navier-Stokes system is solved using a GCR algorithm with a
254  * block diagonal preconditioner using a Schur approximation for the pressure
255  * block (the block 22). The system is stored using a hybrid
256  * assembled/unassembled blocks. The velocity block is assembled (with
257  * potentially sub-blocks for each component) and the velocity
258  * divergence/pressure gradient operators are unassembled.
259  *
260  *
261  * \var CS_NAVSTO_SLES_DIAG_SCHUR_MINRES
262  * Associated keyword: "diag_schur_minres"
263  *
264  * The Stokes or Navier-Stokes system with an explicit advection is solved
265  * using a MINRES algorithm with a block diagonal preconditioner using a Schur
266  * approximation for the pressure block (the block 22). The system is stored
267  * using a hybrid assembled/unassembled blocks. The velocity block is assembled
268  * (with potentially sub-blocks for each component) and the velocity
269  * divergence/pressure gradient operators are unassembled.
270  *
271  *
272  * \var CS_NAVSTO_SLES_EQ_WITHOUT_BLOCK
273  * Associated keyword: "no_block"
274  *
275  * Use the same mechanism as for a stand-alone equation. In this case, the
276  * setting relies on the function \ref cs_equation_set_sles and the different
277  * options for solving a linear system such as the choice of the iterative
278  * solver or the choice of the preconditioner or the type of residual
279  * normalization
280  *
281  *
282  * \var CS_NAVSTO_SLES_GCR
283  * Associated keyword: "gcr"
284  *
285  * The Stokes or Navier-Stokes system is solved using a GCR algorithm without
286  * preconditioning. The system is stored using a hybrid assembled/unassembled
287  * blocks. The velocity block is assembled (with potentially sub-blocks for
288  * each component) and the velocity divergence/pressure gradient operators are
289  * unassembled.
290  *
291  *
292  * \var CS_NAVSTO_SLES_GKB_PETSC
293  * Associated keyword: "gkb_petsc"
294  *
295  * Available choice when a monolithic approach is used (i.e. with the parameter
296  * CS_NAVSTO_COUPLING_MONOLITHIC is set as coupling algorithm). The
297  * Navier-Stokes system of equations is solved using a Golub-Kahan
298  * bi-diagonalization. One assumes that the saddle-point system is symmetric.
299  * By default, the block A_{00} may be augmented (this is not the default
300  * choice) and is solved with a conjugate gradient algorithm preconditioned
301  * with a multigrid. The residual is computed in the energy norm. This option is
302  * only available with the support to the PETSc library up to now.
303  *
304  * * \var CS_NAVSTO_SLES_GKB_GMRES
305  * Associated keyword: "gkb_gmres"
306  *
307  * Available choice when a monolithic approach is used (i.e. with the parameter
308  * CS_NAVSTO_COUPLING_MONOLITHIC is set as coupling algorithm). The
309  * Navier-Stokes system of equations is solved using a Golub-Kahan
310  * bi-diagonalization (GKB) as preconditioner of a flexible GMRES solver. The
311  * GKB algorithm is solved with a reduced tolerance as well as the CG+Multigrid
312  * used as an inner solver in the GKB algorithm. One assumes that the
313  * saddle-point system is symmetric. The residual for the GKB part is computed
314  * in the energy norm. This option is only available with the support to the
315  * PETSc library up to now.
316  *
317  * \var CS_NAVSTO_SLES_GKB_SATURNE
318  * Associated keyword: "gkb" or "gkb_saturne"
319  *
320  * Available choice when a monolithic approach is used (i.e. with the parameter
321  * CS_NAVSTO_COUPLING_MONOLITHIC is set as coupling algorithm). The
322  * Navier-Stokes system of equations is solved using a Golub-Kahan
323  * bi-diagonalization.
324  * By default, the block A_{00} may be augmented (this is not the default
325  * choice) and is solved with the SLES settings given to the momentum equation
326  * A conjugate gradient algorithm preconditioned
327  * with a multigrid for Stokes for instance.
328  * The residual is computed in the energy norm.
329  *
330  *
331  * \var CS_NAVSTO_SLES_LOWER_SCHUR_GCR
332  * Associated keyword: "lower_schur_gcr"
333  *
334  * The Stokes or Navier-Stokes system is solved using a GCR algorithm with an
335  * lower triangular block preconditioner using a Schur approximation for the
336  * pressure block (the block 22). The system is stored using a hybrid
337  * assembled/unassembled blocks. The velocity block is assembled (with
338  * potentially sub-blocks for each component) and the velocity
339  * divergence/pressure gradient operators are unassembled.
340  *
341  *
342  * \var CS_NAVSTO_SLES_MINRES
343  * Associated keyword: "minres"
344  *
345  * The Stokes or Navier-Stokes system with an explicit advection is solved
346  * using a MINRES algorithm without preconditioning. The system is stored using
347  * a hybrid assembled/unassembled blocks. The velocity block is assembled (with
348  * potentially sub-blocks for each component) and the velocity
349  * divergence/pressure gradient operators are unassembled.
350  *
351  *
352  * \var CS_NAVSTO_SLES_MULTIPLICATIVE_GMRES_BY_BLOCK
353  * Associated keyword: "multiplicative_gmres"
354  *
355  * Available choice when a monolithic approach is used (i.e. with the parameter
356  * CS_NAVSTO_COUPLING_MONOLITHIC is set as coupling algorithm) The Navier-Stokes
357  * system of equations is solved a multiplicative preconditioner (block diagonal
358  * matrix where the block 00 is A_{00}) and the block 11 is set to the identity.
359  * Block 01 is also considered in the block preconditioner.
360  * Preconditioner/solver for the block 00 is set using the momentum equation.
361  * This option is only available with the support to the PETSc library up to
362  * now.
363  *
364  * \var CS_NAVSTO_SLES_MUMPS
365  * Associated keyword: "mumps"
366  *
367  * Direct solver to solve the full (saddle-point) system arising from the
368  * discretization of the Navier-Stokes equations
369  *
370  * \var CS_NAVSTO_SLES_NOTAY_TRANSFORM
371  * Associated keyword: "notay"
372  *
373  * Transform the saddle-point problem into an equivalent system without a zero
374  * block for the (2,2) block. This transformation allows one to consider
375  * standard preconditionner and iterative Krylow solver.
376  *
377  * \var CS_NAVSTO_SLES_SGS_SCHUR_GCR
378  * Associated keyword: "sgs_schur_gcr"
379  *
380  * The Stokes or Navier-Stokes system is solved using a GCR algorithm with a
381  * symmetric Gauss-Seidel block preconditioner using a Schur approximation for
382  * the pressure block (the block 22). The system is stored using a hybrid
383  * assembled/unassembled blocks. The velocity block is assembled (with
384  * potentially sub-blocks for each component) and the velocity
385  * divergence/pressure gradient operators are unassembled.
386  *
387  *
388  * \var CS_NAVSTO_SLES_UPPER_SCHUR_GCR
389  * Associated keyword: "upper_schur_gcr"
390  *
391  * The Stokes or Navier-Stokes system is solved using a GCR algorithm with an
392  * upper triangular block preconditioner using a Schur approximation for the
393  * pressure block (the block 22). The system is stored using a hybrid
394  * assembled/unassembled blocks. The velocity block is assembled (with
395  * potentially sub-blocks for each component) and the velocity
396  * divergence/pressure gradient operators are unassembled.
397  *
398  *
399  * \var CS_NAVSTO_SLES_UPPER_SCHUR_GMRES
400  * Associated keyword: "upper_schur_gmres"
401  *
402  * Available choice when a monolithic approach is used (i.e. with the parameter
403  * CS_NAVSTO_COUPLING_MONOLITHIC is set as coupling algorithm). The
404  * Navier-Stokes system of equations is solved using a upper triangular block
405  * preconditioner where the block 00 is A_{00} preconditioned with one multigrid
406  * iteration and the block 11 is an approximation of the Schur complement
407  * preconditioned with a minres. The main iterative solver is a flexible
408  * GMRES. This option is only available with the support to the PETSc
409  * library up to now.
410  *
411  * \var CS_NAVSTO_SLES_USER
412  * Associated keyword: "user"
413  *
414  * Resolution using a user-defined function. This function fulfills a
415  * pre-defined prototype
416  *
417  * \var CS_NAVSTO_SLES_UZAWA_AL
418  * Associated keyword: "uzawa_al"
419  *
420  * Resolution using an uzawa algorithm with an Augmented Lagrangian approach
421  *
422  * \var CS_NAVSTO_SLES_UZAWA_CG
423  * Associated keyword: "uzawa_cg"
424  *
425  * Resolution using an uzawa algorithm optimized using a conjugate gradient
426  * reformulation. Two systems are solved at each iteration (one related to the
427  * velocity block, one related to the Schur complement approximation - size of
428  * the pressure space).
429  *
430  * \var CS_NAVSTO_SLES_UZAWA_SCHUR_GCR
431  * Associated keyword: "uza_schur_gcr"
432  *
433  * The Stokes or Navier-Stokes system is solved using a GCR algorithm with an
434  * Uzawa algorithm tuned for block preconditioning and using a Schur
435  * approximation for the pressure block (the block 22). The system is stored
436  * using a hybrid assembled/unassembled blocks. The velocity block is assembled
437  * (with potentially sub-blocks for each component) and the velocity
438  * divergence/pressure gradient operators are unassembled.
439  */
440 
441 typedef enum {
442 
443   CS_NAVSTO_SLES_ADDITIVE_GMRES_BY_BLOCK,
444   CS_NAVSTO_SLES_BLOCK_MULTIGRID_CG,
445   CS_NAVSTO_SLES_BY_BLOCKS,     /* deprecated */
446   CS_NAVSTO_SLES_DIAG_SCHUR_GCR,
447   CS_NAVSTO_SLES_DIAG_SCHUR_GMRES,
448   CS_NAVSTO_SLES_DIAG_SCHUR_MINRES,
449   CS_NAVSTO_SLES_EQ_WITHOUT_BLOCK,
450   CS_NAVSTO_SLES_GCR,
451   CS_NAVSTO_SLES_GKB_PETSC,
452   CS_NAVSTO_SLES_GKB_GMRES,
453   CS_NAVSTO_SLES_GKB_SATURNE,
454   CS_NAVSTO_SLES_LOWER_SCHUR_GCR,
455   CS_NAVSTO_SLES_MINRES,
456   CS_NAVSTO_SLES_MULTIPLICATIVE_GMRES_BY_BLOCK,
457   CS_NAVSTO_SLES_MUMPS,
458   CS_NAVSTO_SLES_NOTAY_TRANSFORM,
459   CS_NAVSTO_SLES_SGS_SCHUR_GCR,
460   CS_NAVSTO_SLES_UPPER_SCHUR_GCR,
461   CS_NAVSTO_SLES_UPPER_SCHUR_GMRES,
462   CS_NAVSTO_SLES_USER,
463   CS_NAVSTO_SLES_UZAWA_SCHUR_GCR,
464   CS_NAVSTO_SLES_UZAWA_AL,
465   CS_NAVSTO_SLES_UZAWA_CG,
466 
467   CS_NAVSTO_SLES_N_TYPES
468 
469 } cs_navsto_sles_t;
470 
471 /*! \struct cs_navsto_param_sles_t
472  *  \brief Structure storing the parameters for solving the Navier-Stokes system
473  */
474 
475 typedef struct {
476 
477   /*! \var strategy
478    *  Choice of strategy for solving the Navier--Stokes system
479    */
480 
481   cs_navsto_sles_t              strategy;
482 
483   /*! \var schur_strategy
484    *  Choice of the way of preconditioning the schur approximation
485    */
486 
487   cs_param_schur_approx_t       schur_approximation;
488 
489   /*!
490    * @name Inner and linear algorithm
491    * Set of parameters to drive the resolution of the (inner) linear system
492    * @{
493    */
494 
495   /*! \var il_algo_param
496    *  Structure storing several tolerances and metadata to drive the inner
497    *  (linear) iterative algorithm used to solve either the Oseen or the Stokes
498    *  system. This algorithm is for instance an Uzawa or GKB algorithm. This is
499    *  incorporated in a non-linear process in case of Navier--Stokes equations.
500    */
501 
502   cs_iter_algo_param_t          il_algo_param;
503 
504   /*! \var il_algo_restart
505    *  Number of iterations before restarting the iterative solver associated to
506    *  the inner linear system
507    */
508 
509   int                           il_algo_restart;
510 
511   /*!
512    * @}
513    * @name Non-linear algorithm
514    * Set of parameters to drive the resolution of the non-linearity arising from
515    * the Navier--Stokes system
516    * @{
517    */
518 
519   /*! \var nl_algo
520    *  Type of algorithm used to tackle the non-linearity arising from the
521    *  advection term
522    */
523 
524   cs_param_nl_algo_t            nl_algo_type;
525 
526   /*! \var nl_algo_param
527    *  Structure storing several tolerances and metadata to drive the non-linear
528    *  iterative algorithm used to solve te Navier-Stokes when the advection
529    *  term is implicit and non linearized.
530    */
531 
532   cs_iter_algo_param_t          nl_algo_param;
533 
534   /*! \var anderson_param
535    * Set of parameters to drive the Anderson acceleration (useful if the type
536    * of non-linear algorithm is set to Anderson acceleration).
537    */
538 
539   cs_iter_algo_param_aa_t       anderson_param;
540 
541   /*! @} */
542 
543   /*!
544    * @}
545    * @name Block preconditioning or Schur complement approximation
546    * Set of parameters to drive the resolution of the pressure-related
547    * block. This is often a Schur complement approximation to B.A^-1.Bt
548    * @{
549    */
550 
551   cs_param_sles_t              *schur_sles_param;
552 
553   /*! @} */
554 
555 } cs_navsto_param_sles_t;
556 
557 /*! \enum cs_navsto_param_coupling_t
558  *  \brief Choice of algorithm for solving the system
559  *
560  * \var CS_NAVSTO_COUPLING_ARTIFICIAL_COMPRESSIBILITY
561  * The system is solved using an artificial compressibility algorithm.
562  * One vectorial equation is solved followed by a pressure update.
563  *
564  * \var CS_NAVSTO_COUPLING_MONOLITHIC
565  * The system is treated as a "monolithic" matrix
566  *
567  * \var CS_NAVSTO_COUPLING_PROJECTION
568  * The system is solved using an incremental projection algorithm
569  */
570 
571 typedef enum {
572 
573   CS_NAVSTO_COUPLING_ARTIFICIAL_COMPRESSIBILITY,
574   CS_NAVSTO_COUPLING_MONOLITHIC,
575   CS_NAVSTO_COUPLING_PROJECTION,
576 
577   CS_NAVSTO_N_COUPLINGS
578 
579 } cs_navsto_param_coupling_t;
580 
581 /*! \struct cs_navsto_param_boussinesq_t
582  *  \brief Structure storing the parameters related to the Boussinesq source
583  *         term in the momentum equation
584  */
585 
586 typedef struct {
587 
588   cs_real_t    beta;      /* Dilatation coefficient */
589   cs_real_t    var0;      /* Reference value of the variable */
590 
591   /* Array of values of the variable (for instance the temperature). This is a
592    * shared pointer. The lifecycle of this array is not managed by this
593    * structure.
594    */
595 
596   const cs_real_t   *var;
597 
598 } cs_navsto_param_boussinesq_t;
599 
600 /*! \struct cs_navsto_param_t
601  *  \brief Structure storing the parameters related to the resolution of the
602  *         Navier-Stokes system
603  */
604 
605 typedef struct {
606 
607   /*!
608    * @name Physical modelling
609    * Which equations to solve ?  Properties and their related fields are
610    * allocated according to the choice of model for Navier-Stokes
611    * @{
612    */
613 
614   /*! \var model
615    * Modelling related to the Navier-Stokes system of equations
616    */
617 
618   cs_navsto_param_model_t        model;
619 
620   /*! \var model_flag
621    * Flag storing high-level option related to the Navier-Stokes system
622    */
623 
624   cs_navsto_param_model_flag_t   model_flag;
625 
626   /*! \var turbulence
627    *  Structure storing all information needed to set the turbulence modelling
628    */
629 
630   cs_turbulence_param_t         *turbulence;
631 
632   /*!
633    * @}
634    * @name Properties and fields related to the Navier-Stokes module
635    * @{
636    */
637 
638   /*! \var phys_constants
639    * Main physical constants (gravity vector and Coriolis source term). This
640    * structure is shared with the legacy part.
641    */
642 
643   cs_physical_constants_t       *phys_constants;
644 
645   /*! \var mass_density
646    * Mass_density of the fluid, pointer to \ref cs_property_t used in several
647    * terms in the Navier-Stokes equations
648    */
649 
650   cs_property_t                 *mass_density;
651 
652   /*! \var tot_viscosity
653    *  Laminar viscosity + if needed the turbulent viscosity
654    *  Pointer to \ref cs_property_t associated to the
655    *  diffusion term for the momentum equation
656    */
657 
658   cs_property_t                *tot_viscosity;
659 
660   /*! \var lami_viscosity
661    *  Laminar viscosity
662    */
663 
664   cs_property_t                *lam_viscosity;
665 
666   /*!
667    * @}
668    * @name Numerical options
669    * Set of numerical options to build the linear system and how to solve it
670    * @{
671    */
672 
673   /*! \var dof_reduction_mode
674    *  How are defined the Degrees of freedom
675    */
676 
677   cs_param_dof_reduction_t       dof_reduction_mode;
678 
679   /*! \var coupling
680    * Choice of algorithm for solving the system
681    */
682 
683   cs_navsto_param_coupling_t     coupling;
684 
685   /*! \var gd_scale_coef
686    *  Default value to set the scaling of the grad-div term when an
687    *  artificial compressibility algorithm or an Uzawa-Augmented Lagrangian
688    *  method is used
689    */
690 
691   cs_real_t                      gd_scale_coef;
692 
693   /*! \var time_scheme
694    * Discretization scheme for time
695    *
696    * \var theta
697    * Value of the parameter for the time scheme when a theta-scheme is used
698    *
699    */
700 
701   cs_param_time_scheme_t         time_scheme;
702   cs_real_t                      theta;
703 
704   /*! \var space_scheme
705    * Discretization scheme for space
706    */
707 
708   cs_param_space_scheme_t        space_scheme;
709 
710   /*! \var adv_form
711    *  Type of formulation for the advection term
712    *
713    *  \var adv_scheme
714    *  Type of scheme for the advection term
715    *
716    *  \var adv_strategy
717    *  Strategy to handle the advection term
718    *
719    *  \var adv_extrapol
720    *  Extrapolation used to estimate the advection field
721    */
722 
723   cs_param_advection_form_t      adv_form;
724   cs_param_advection_scheme_t    adv_scheme;
725   cs_param_advection_strategy_t  adv_strategy;
726   cs_param_advection_extrapol_t  adv_extrapol;
727 
728   /* Boussinesq approximation:
729    *
730    * Take into account buoyancy terms (variation of mass density w.r.t. the
731    * variation of a field (for instance the temperature but can be also a
732    * concentration as in segregation model in the solidification module)
733    *
734    * \var n_boussinesq_terms
735    * Number of contributions to the buoyancy source term in the Boussinesq
736    * approximation
737    *
738    * \var boussinesq_param
739    * Structure storing elements used to compute the Boussinesq approximation
740    */
741 
742   int                                 n_boussinesq_terms;
743   cs_navsto_param_boussinesq_t       *boussinesq_param;
744 
745   /*! \var qtype
746    *  A \ref cs_quadrature_type_t indicating the type of quadrature to use in
747    *  all functions involving quadratures
748    */
749 
750   cs_quadrature_type_t           qtype;
751 
752   /*! \var sles_param
753    * Set of choices to control the resolution of the Navier--Stokes system
754    */
755 
756   cs_navsto_param_sles_t        *sles_param;
757 
758   /*! \var delta_thermal_tolerance
759    * Value under which one considers that the thermal equation is converged
760    * max_{c \in Cells} |T_c - T_{ref}|/|T_{ref}| < eps => stop iteration
761    */
762 
763   cs_real_t                      delta_thermal_tolerance;
764 
765   /*! \var n_max_outer_iter
766    * Stopping crierion related to the maximum number of outer iterations
767    * allowed. This outer iteration encompasses the Navier-Stokes system,
768    * and (according to the case settings) the turbulence system and/or
769    * the thermal system.
770    */
771 
772   int                            n_max_outer_iter;
773 
774   /*!
775    * @}
776    * @name Output
777    * @{
778    *
779    * \var verbosity
780    * Level of display of the information related to the Navier-Stokes system
781    */
782 
783   int                            verbosity;
784 
785   /*! \var post_flag
786    * Flag storing which predefined post-processing has to be done
787    */
788 
789   cs_navsto_param_post_flag_t    post_flag;
790 
791   /*!
792    * @}
793    * @name Initial conditions (IC)
794    *
795    * Set of parameters used to take into account the initial condition on the
796    * pressure and/or the velocity.
797    * CAUTION: so far, there is no check if the different IC are compatible
798    * with the boundary conditions for instance
799    * @{
800    */
801 
802   /*! \var velocity_ic_is_owner
803    *  True if the definitions are stored inside this structure, otherwise
804    *  the definitions are stored inside a \ref cs_equation_param_t
805    *  structure dedicated to the momentum equation.
806    *
807    * \var n_velocity_ic_defs
808    *  Number of initial conditions associated to the velocity
809    *
810    * \var velocity_ic_defs
811    *  Pointers to the definitions of the initial conditions associated to the
812    *  velocity.
813    *  The code does not check if the resulting initial velocity satisfies the
814    *  divergence constraint.
815    */
816 
817   bool         velocity_ic_is_owner;
818   int          n_velocity_ic_defs;
819   cs_xdef_t  **velocity_ic_defs;
820 
821   /*! \var pressure_ic_is_owner
822    *  True if the definitions are stored inside this structure, otherwise
823    *  the definitions are stored inside a dedicated \ref cs_equation_param_t
824    *
825    * \var n_pressure_ic_defs
826    *  Number of initial conditions associated to the pressure
827    *
828    * \var pressure_ic_defs
829    *  Pointers to the definitions of the initial conditions associated to the
830    *  pressure.
831    *  In order to force a zero-mean pressure, the code can compute the average
832    *  of the resulting pressure and subtract it
833    */
834 
835   bool         pressure_ic_is_owner;
836   int          n_pressure_ic_defs;
837   cs_xdef_t  **pressure_ic_defs;
838 
839   /*! @}
840    * @name Boundary conditions (BC)
841    *
842    * Set of parameters used to take into account the boundary conditions on the
843    * pressure and/or the velocity.
844    * @{
845    */
846 
847   /* \var boundaries
848    * Pointer to a \ref cs_boundary_t structure shared with the domain
849    */
850   const cs_boundary_t   *boundaries;
851 
852   /*! \var velocity_bc_is_owner
853    *  True if the definitions are stored inside this structure, otherwise
854    *  the definitions are stored inside a dedicated \ref cs_equation_param_t
855    *  Most of the time this should be false since an equation is associated to
856    *  to the resolution of the velocity field (the momentum equation).
857    *
858    * \var n_velocity_bc_defs
859    * Number of definitions related to the settings of the boundary conditions
860    * for the velocity field.
861    *
862    * \var velocity_bc_defs
863    * Array of pointers to the definition of boundary conditions for the velocity
864    * field
865    */
866 
867   bool         velocity_bc_is_owner;
868   int          n_velocity_bc_defs;
869   cs_xdef_t  **velocity_bc_defs;
870 
871   /*! \var pressure_bc_is_owner
872    *  True if the definitions are stored inside this structure, otherwise
873    *  the definitions are stored inside a dedicated \ref cs_equation_param_t
874    *  if an equation solves the pressure field.
875    *
876    * \var n_pressure_bc_defs
877    *  Number of boundary conditions associated to the pressure field.
878    *
879    * \var pressure_bc_defs
880    *  Pointers to the definitions of the boundary conditions associated to the
881    *  pressure field.
882    */
883 
884   bool         pressure_bc_is_owner;
885   int          n_pressure_bc_defs;
886   cs_xdef_t  **pressure_bc_defs;
887 
888   /*! @}
889    * @name Other conditions
890    * @{
891    */
892 
893   /*! \var reference_pressure
894    *  Value of the reference pressure p0 (used for rescaling or during update
895    *  of physical quantities). By default: 0.
896    */
897 
898   cs_real_t    reference_pressure;
899 
900   /*! @} */
901 
902 } cs_navsto_param_t;
903 
904 /*! \enum cs_navsto_key_t
905  *  \brief List of available keys for setting the parameters of the
906  *         Navier-Stokes system
907  *
908  * \var CS_NSKEY_ADVECTION_EXTRAPOL
909  * Set the extrapolation to use for the estimation of the advection field
910  * (cf. \ref cs_param_advection_extrapol_t))
911  *
912  * \var CS_NSKEY_ADVECTION_FORMULATION
913  * Set the type of formulation for the advection term, for example in the Oseen
914  * problem. (cf. \ref cs_param_advection_form_t)
915  *
916  * \var CS_NSKEY_ADVECTION_SCHEME
917  * Set the type of scheme for the advection term, for example in the Oseen
918  * problem. (cf. \ref cs_param_advection_scheme_t)
919  *
920  * \var CS_NSKEY_ADVECTION_STRATEGY
921  * Set the strategy to handle the advection term
922  * (cf. \ref cs_param_advection_strategy_t)
923  *
924  * \var CS_NSKEY_DOF_REDUCTION
925  * Set how the DoFs are defined (similar to \ref CS_EQKEY_DOF_REDUCTION)
926  * Enable to set this type of DoFs definition for all related equations
927  *
928  * \var CS_NSKEY_GD_SCALE_COEF
929  * Set the scaling of the grad-div term when an artificial compressibility
930  * algorithm or an Uzawa-Augmented Lagrangian method is used
931  *
932  * \var CS_NSKEY_IL_ALGO_ATOL
933  * Absolute tolerance at which the Oseen or Stokes system is resolved. These
934  * systems corresponds to an inner linear system to solve when considering the
935  * Navier-Stokes system since one has to handle the non-linearity in addition as
936  * an outer process.
937  *
938  * \var CS_NSKEY_IL_ALGO_DTOL
939  * Divergence tolerance at which the Oseen or Stokes system is resolved. These
940  * systems corresponds to an inner linear system to solve when considering the
941  * Navier-Stokes system since one has to handle the non-linearity in addition as
942  * an outer process.
943  *
944  * \var CS_NSKEY_IL_ALGO_RTOL
945  * Relative tolerance at which the Oseen or Stokes system is resolved. These
946  * systems corresponds to an inner linear system to solve when considering the
947  * Navier-Stokes system since one has to handle the non-linearity in addition as
948  * an outer process.
949  *
950  * \var CS_NSKEY_IL_ALGO_RESTART
951  * Number of iterations before restarting a Krylov solver as the main solver
952  * (useful if the strategy implied a GMRES, flexible GMRES or GCR)
953  *
954  * \var CS_NSKEY_IL_ALGO_VERBOSITY
955  * Level of verbosity related to the inner linear algorithm (cf. \ref
956  * CS_NSKEY_SLES_STRATEGY)
957  *
958  * \var CS_NSKEY_MAX_IL_ALGO_ITER
959  * Set the maximal number of iteration for solving the inner linear system.
960  *
961  * \var CS_NSKEY_MAX_NL_ALGO_ITER
962  * Set the maximal number of iterations for solving the non-linearity
963  * arising from the advection form
964  *
965  * \var CS_NSKEY_MAX_OUTER_ITER
966  * Set the maximal number of outer iterations for solving the full system
967  * including the turbulence modelling or the thermal system for instance
968  *
969  * \var CS_NSKEY_NL_ALGO
970  * Type of algorithm to consider to solve the non-linearity arising from the
971  * Navier-Stokes system (Picard or Anderson)
972  *
973  * \var CS_NSKEY_NL_ALGO_ATOL
974  * Absolute tolerance at which the non-linearity arising from the advection
975  * term is resolved
976  *
977  * \var CS_NSKEY_NL_ALGO_DTOL
978  * Divergence tolerance at which the non-linearity arising from the advection
979  * term is detected
980  *
981  * \var CS_NSKEY_NL_ALGO_RTOL
982  * Relative tolerance at which the non-linearity arising from the advection
983  * term is resolved
984  *
985  * \var CS_NSKEY_NL_ALGO_VERBOSITY
986  * Level of verbosity related to the non-linear algorithm
987  *
988  * \var CS_NSKEY_QUADRATURE
989  * Set the type to use in all functions involving quadrature (similar to
990  * \ref CS_EQKEY_BC_QUADRATURE)
991  *
992  * \var CS_NSKEY_SCHUR_STRATEGY
993  * Set the way to define the Schur complement approximation
994  * (cf. \ref cs_param_schur_approx_t)
995  *
996  * \var CS_NSKEY_SLES_STRATEGY
997  * Strategy for solving the SLES arising from the discretization of the
998  * Navier-Stokes system
999  *
1000  * \var CS_NSKEY_SPACE_SCHEME
1001  * Numerical scheme for the space discretization. Available choices are:
1002  * - "cdo_fb" or "cdofb" for CDO face-based scheme
1003  *
1004  * \var CS_NSKEY_THERMAL_TOLERANCE
1005  * Value of the tolerance criterion under which one stops iterating on
1006  * the Navier-Stokes and thermal systems
1007  *
1008  * \var CS_NSKEY_TIME_SCHEME
1009  * Numerical scheme for the time discretization
1010  *
1011  * \var CS_NSKEY_TIME_THETA
1012  * Set the value of theta. Only useful if CS_NSKEY_TIME_SCHEME is set to
1013  * "theta_scheme"
1014  * - Example: "0.75" (keyval must be between 0 and 1)
1015  *
1016  * \var CS_NSKEY_VERBOSITY
1017  * Set the level of details for the specific part related to the Navier-Stokes
1018  * system
1019  */
1020 
1021 typedef enum {
1022 
1023   CS_NSKEY_ADVECTION_EXTRAPOL,
1024   CS_NSKEY_ADVECTION_FORMULATION,
1025   CS_NSKEY_ADVECTION_SCHEME,
1026   CS_NSKEY_ADVECTION_STRATEGY,
1027   CS_NSKEY_DOF_REDUCTION,
1028   CS_NSKEY_GD_SCALE_COEF,
1029   CS_NSKEY_IL_ALGO_ATOL,
1030   CS_NSKEY_IL_ALGO_DTOL,
1031   CS_NSKEY_IL_ALGO_RTOL,
1032   CS_NSKEY_IL_ALGO_RESTART,
1033   CS_NSKEY_IL_ALGO_VERBOSITY,
1034   CS_NSKEY_MAX_IL_ALGO_ITER,
1035   CS_NSKEY_MAX_NL_ALGO_ITER,
1036   CS_NSKEY_MAX_OUTER_ITER,
1037   CS_NSKEY_NL_ALGO,
1038   CS_NSKEY_NL_ALGO_ATOL,
1039   CS_NSKEY_NL_ALGO_DTOL,
1040   CS_NSKEY_NL_ALGO_RTOL,
1041   CS_NSKEY_NL_ALGO_VERBOSITY,
1042   CS_NSKEY_QUADRATURE,
1043   CS_NSKEY_SCHUR_STRATEGY,
1044   CS_NSKEY_SLES_STRATEGY,
1045   CS_NSKEY_SPACE_SCHEME,
1046   CS_NSKEY_THERMAL_TOLERANCE,
1047   CS_NSKEY_TIME_SCHEME,
1048   CS_NSKEY_TIME_THETA,
1049   CS_NSKEY_VERBOSITY,
1050 
1051   CS_NSKEY_N_KEYS
1052 
1053 } cs_navsto_key_t;
1054 
1055 /*============================================================================
1056  * Inline static public function prototypes
1057  *============================================================================*/
1058 
1059 /*----------------------------------------------------------------------------*/
1060 /*!
1061  * \brief  Ask \ref cs_navsto_param_t structure if the settings correspond to
1062  *         a steady computation
1063  *
1064  * \param[in]  nsp     pointer to a \ref cs_navsto_param_t structure
1065 *
1066  * \return true or false
1067  */
1068 /*----------------------------------------------------------------------------*/
1069 
1070 static inline bool
cs_navsto_param_is_steady(const cs_navsto_param_t * nsp)1071 cs_navsto_param_is_steady(const cs_navsto_param_t       *nsp)
1072 {
1073   if (nsp == NULL)
1074     return true;
1075 
1076   if (nsp->model_flag & CS_NAVSTO_MODEL_STEADY)
1077     return true;
1078   else
1079     return false;
1080 }
1081 
1082 /*============================================================================
1083  * Public function prototypes
1084  *============================================================================*/
1085 
1086 /*----------------------------------------------------------------------------*/
1087 /*!
1088  * \brief  Retrieve the scaling coefficient used in the Notay's transformation
1089  *         devised in "Algebraic multigrid for Stokes equations" SIAM
1090  *         J. Sci. Comput. Vol. 39 (5), 2017
1091  *         In this article, this scaling is denoted by alpha
1092  *
1093  * \return the value of the scaling coefficient
1094  */
1095 /*----------------------------------------------------------------------------*/
1096 
1097 double
1098 cs_navsto_param_get_notay_scaling(void);
1099 
1100  /*----------------------------------------------------------------------------*/
1101 /*!
1102  * \brief  Set the scaling coefficient used in the Notay's transformation
1103  *         devised in "Algebraic multigrid for Stokes equations" SIAM
1104  *         J. Sci. Comput. Vol. 39 (5), 2017
1105  *         In this article, this scaling is denoted by alpha
1106  *
1107  * \param[in]  scaling_coef   valued of the scaling coefficient
1108  */
1109 /*----------------------------------------------------------------------------*/
1110 
1111 void
1112 cs_navsto_param_set_notay_scaling(double  scaling_coef);
1113 
1114 /*----------------------------------------------------------------------------*/
1115 /*!
1116  * \brief  Create a new structure to store all numerical parameters related
1117  *         to the resolution of the Navier-Stokes (NS) system
1118  *
1119  * \param[in] boundaries      pointer to a cs_boundary_t structure
1120  * \param[in] model           type of model related to the NS system
1121  * \param[in] model_flag      additional high-level model options
1122  * \param[in] algo_coupling   algorithm used for solving the NS system
1123  * \param[in] post_flag       predefined post-processings options
1124  *
1125  * \return a pointer to a new allocated structure
1126  */
1127 /*----------------------------------------------------------------------------*/
1128 
1129 cs_navsto_param_t *
1130 cs_navsto_param_create(const cs_boundary_t            *boundaries,
1131                        cs_navsto_param_model_t         model,
1132                        cs_navsto_param_model_flag_t    model_flag,
1133                        cs_navsto_param_coupling_t      algo_coupling,
1134                        cs_navsto_param_post_flag_t     post_flag);
1135 
1136 /*----------------------------------------------------------------------------*/
1137 /*!
1138  * \brief  Free a \ref cs_navsto_param_t structure
1139  *
1140  * \param[in, out]  param    pointer to a \ref cs_navsto_param_t structure
1141  *
1142  * \return a NULL pointer
1143  */
1144 /*----------------------------------------------------------------------------*/
1145 
1146 cs_navsto_param_t *
1147 cs_navsto_param_free(cs_navsto_param_t    *param);
1148 
1149 /*----------------------------------------------------------------------------*/
1150 /*!
1151  * \brief  Set a parameter attached to a keyname in a \ref cs_navsto_param_t
1152  *         structure
1153  *
1154  * \param[in, out] nsp      pointer to a \ref cs_navsto_param_t structure to set
1155  * \param[in]      key      key related to the member of eq to set
1156  * \param[in]      keyval   accessor to the value to set
1157  */
1158 /*----------------------------------------------------------------------------*/
1159 
1160 void
1161 cs_navsto_param_set(cs_navsto_param_t    *nsp,
1162                     cs_navsto_key_t       key,
1163                     const char           *keyval);
1164 
1165 /*----------------------------------------------------------------------------*/
1166 /*!
1167  * \brief  Apply the numerical settings defined for the Navier-Stokes system
1168  *         to an equation related to this system.
1169  *
1170  * \param[in]       nsp    pointer to a \ref cs_navsto_param_t structure
1171  * \param[in, out]  eqp    pointer to a \ref cs_equation_param_t structure
1172  */
1173 /*----------------------------------------------------------------------------*/
1174 
1175 void
1176 cs_navsto_param_transfer(const cs_navsto_param_t    *nsp,
1177                          cs_equation_param_t        *eqp);
1178 
1179 /*----------------------------------------------------------------------------*/
1180 /*!
1181  * \brief  Summary of the main cs_navsto_param_t structure
1182  *
1183  * \param[in]  nsp    pointer to a cs_navsto_param_t structure
1184  */
1185 /*----------------------------------------------------------------------------*/
1186 
1187 void
1188 cs_navsto_param_log(const cs_navsto_param_t    *nsp);
1189 
1190 /*----------------------------------------------------------------------------*/
1191 /*!
1192  * \brief  Add a new Boussinesq term (source term for the momemtum equation)
1193  *
1194  * \param[in, out]  nsp    pointer to a cs_navsto_param_t structure
1195  *
1196  * \return a pointer to the newly added structure
1197  */
1198 /*----------------------------------------------------------------------------*/
1199 
1200 cs_navsto_param_boussinesq_t *
1201 cs_navsto_param_add_boussinesq_term(cs_navsto_param_t    *nsp,
1202                                     cs_real_t             dilatation_coef,
1203                                     cs_real_t             reference_value);
1204 
1205 /*----------------------------------------------------------------------------*/
1206 /*!
1207  * \brief  Set the array of values to consider in the Boussinesq term
1208  *
1209  * \param[in, out]  bp    pointer to a cs_navsto_param_boussinesq_t structure
1210  * \param[in]       var   shared pointer to the array of values to consider
1211  */
1212 /*----------------------------------------------------------------------------*/
1213 
1214 void
1215 cs_navsto_param_set_boussinesq_array(cs_navsto_param_boussinesq_t   *bp,
1216                                      const cs_real_t                *var);
1217 
1218 /*----------------------------------------------------------------------------*/
1219 /*!
1220  * \brief  Retrieve the \ref cs_equation_param_t structure related to the
1221  *         velocity equation (momentum equation in most of the cases)
1222  *
1223  * \param[in]  nsp    pointer to a cs_navsto_param_t structure
1224  *
1225  * \return a pointer to the set of SLES parameters
1226  */
1227 /*----------------------------------------------------------------------------*/
1228 
1229 cs_navsto_param_sles_t *
1230 cs_navsto_param_get_sles_param(const cs_navsto_param_t    *nsp);
1231 
1232 /*----------------------------------------------------------------------------*/
1233 /*!
1234  * \brief  Retrieve the \ref cs_equation_param_t structure related to the
1235  *         velocity equation (momentum equation in most of the cases)
1236  *
1237  * \param[in]  nsp    pointer to a cs_navsto_param_t structure
1238  *
1239  * \return a pointer to the set of parameters related to the momentum equation
1240  */
1241 /*----------------------------------------------------------------------------*/
1242 
1243 cs_equation_param_t *
1244 cs_navsto_param_get_velocity_param(const cs_navsto_param_t    *nsp);
1245 
1246 /*----------------------------------------------------------------------------*/
1247 /*!
1248  * \brief  Retrieve the name of the model system of equations
1249  *
1250  * \param[in]   model    a \ref cs_navsto_param_model_t
1251  *
1252  * \return the corresponding name
1253  */
1254 /*----------------------------------------------------------------------------*/
1255 
1256 const char *
1257 cs_navsto_param_get_model_name(cs_navsto_param_model_t   model);
1258 
1259 /*----------------------------------------------------------------------------*/
1260 /*!
1261  * \brief  Retrieve the name of the coupling algorithm
1262  *
1263  * \param[in]     coupling    a \ref cs_navsto_param_coupling_t
1264  *
1265  * \return the name
1266  */
1267 /*----------------------------------------------------------------------------*/
1268 
1269 const char *
1270 cs_navsto_param_get_coupling_name(cs_navsto_param_coupling_t  coupling);
1271 
1272 /*----------------------------------------------------------------------------*/
1273 /*!
1274  * \brief  Set the value to consider for the reference pressure
1275  *
1276  * \param[in]  nsp       pointer to a \ref cs_navsto_param_t structure
1277  * \param[in]  pref      value of the reference pressure
1278  */
1279 /*----------------------------------------------------------------------------*/
1280 
1281 void
1282 cs_navsto_set_reference_pressure(cs_navsto_param_t    *nsp,
1283                                  cs_real_t             pref);
1284 
1285 /*----------------------------------------------------------------------------*/
1286 /*!
1287  * \brief  Define the initial condition for the velocity unknowns.
1288  *         This definition can be done on a specified mesh location.
1289  *         By default, the unknown is set to zero everywhere.
1290  *         Here the initial value is set to a constant value
1291  *
1292  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1293  * \param[in]      z_name    name of the associated zone (if NULL or "" if
1294  *                           all cells are considered)
1295  * \param[in]      val       pointer to the value
1296  *
1297  * \return a pointer to the new \ref cs_xdef_t structure
1298  */
1299 /*----------------------------------------------------------------------------*/
1300 
1301 cs_xdef_t *
1302 cs_navsto_add_velocity_ic_by_value(cs_navsto_param_t    *nsp,
1303                                    const char           *z_name,
1304                                    cs_real_t            *val);
1305 
1306 /*----------------------------------------------------------------------------*/
1307 /*!
1308  * \brief  Define the initial condition for the velocity unknowns.
1309  *         This definition can be done on a specified mesh location.
1310  *         By default, the unknown is set to zero everywhere.
1311  *         Here the initial value is set according to an analytical function
1312  *
1313  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1314  * \param[in]      z_name    name of the associated zone (if NULL or "" if
1315  *                           all cells are considered)
1316  * \param[in]      analytic  pointer to an analytic function
1317  * \param[in]      input     NULL or pointer to a structure cast on-the-fly
1318  *
1319  * \return a pointer to the new \ref cs_xdef_t structure
1320  */
1321 /*----------------------------------------------------------------------------*/
1322 
1323 cs_xdef_t *
1324 cs_navsto_add_velocity_ic_by_analytic(cs_navsto_param_t      *nsp,
1325                                       const char             *z_name,
1326                                       cs_analytic_func_t     *analytic,
1327                                       void                   *input);
1328 
1329 /*----------------------------------------------------------------------------*/
1330 /*!
1331  * \brief  Define the initial condition for the pressure unknowns.
1332  *         This definition can be done on a specified mesh location.
1333  *         By default, the unknown is set to zero everywhere.
1334  *         Here the initial value is set to a constant value
1335  *
1336  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1337  * \param[in]      z_name    name of the associated zone (if NULL or "" if
1338  *                           all cells are considered)
1339  * \param[in]      val       pointer to the value
1340  *
1341  * \return a pointer to the new \ref cs_xdef_t structure
1342  */
1343 /*----------------------------------------------------------------------------*/
1344 
1345 cs_xdef_t *
1346 cs_navsto_add_pressure_ic_by_value(cs_navsto_param_t    *nsp,
1347                                    const char           *z_name,
1348                                    cs_real_t            *val);
1349 
1350 /*----------------------------------------------------------------------------*/
1351 /*!
1352  * \brief  Define the initial condition for the pressure unknowns.
1353  *         This definition can be done on a specified mesh location.
1354  *         By default, the unknown is set to zero everywhere.
1355  *         Here the initial value is set according to an analytical function
1356  *
1357  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1358  * \param[in]      z_name    name of the associated zone (if NULL or "" if
1359  *                           all cells are considered)
1360  * \param[in]      analytic  pointer to an analytic function
1361  * \param[in]      input     NULL or pointer to a structure cast on-the-fly
1362  *
1363  * \return a pointer to the new \ref cs_xdef_t structure
1364  */
1365 /*----------------------------------------------------------------------------*/
1366 
1367 cs_xdef_t *
1368 cs_navsto_add_pressure_ic_by_analytic(cs_navsto_param_t      *nsp,
1369                                       const char             *z_name,
1370                                       cs_analytic_func_t     *analytic,
1371                                       void                   *input);
1372 
1373 /*----------------------------------------------------------------------------*/
1374 /*!
1375  * \brief  Add the definition of boundary conditions related to a fixed wall
1376  *         into the set of parameters for the management of the Navier-Stokes
1377  *         system of equations
1378  *
1379  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1380  */
1381 /*----------------------------------------------------------------------------*/
1382 
1383 void
1384 cs_navsto_set_fixed_walls(cs_navsto_param_t    *nsp);
1385 
1386 /*----------------------------------------------------------------------------*/
1387 /*!
1388  * \brief  Add the definition of boundary conditions related to a symmetry
1389  *         into the set of parameters for the management of the Navier-Stokes
1390  *         system of equations
1391  *
1392  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1393  */
1394 /*----------------------------------------------------------------------------*/
1395 
1396 void
1397 cs_navsto_set_symmetries(cs_navsto_param_t    *nsp);
1398 
1399 /*----------------------------------------------------------------------------*/
1400 /*!
1401  * \brief  Add the definition of boundary conditions related to outlets
1402  *         into the set of parameters for the management of the Navier-Stokes
1403  *         system of equations
1404  *
1405  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1406  */
1407 /*----------------------------------------------------------------------------*/
1408 
1409 void
1410 cs_navsto_set_outlets(cs_navsto_param_t    *nsp);
1411 
1412 /*----------------------------------------------------------------------------*/
1413 /*!
1414  * \brief  Set the pressure field on a boundary using a uniform value.
1415  *
1416  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1417  * \param[in]      z_name    name of the associated zone (if NULL or "" all
1418  *                           boundary faces are considered)
1419  * \param[in]      value     value to set
1420  *
1421  * \return a pointer to the new \ref cs_xdef_t structure
1422  */
1423 /*----------------------------------------------------------------------------*/
1424 
1425 cs_xdef_t *
1426 cs_navsto_set_pressure_bc_by_value(cs_navsto_param_t    *nsp,
1427                                    const char           *z_name,
1428                                    cs_real_t            *values);
1429 
1430 /*----------------------------------------------------------------------------*/
1431 /*!
1432  * \brief  Define the velocity field for a sliding wall boundary using a
1433  *         uniform value
1434  *
1435  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1436  * \param[in]      z_name    name of the associated zone (if NULL or "" all
1437  *                           boundary faces are considered)
1438  * \param[in]      values    array of three real values
1439  *
1440  * \return a pointer to the new \ref cs_xdef_t structure
1441  */
1442 /*----------------------------------------------------------------------------*/
1443 
1444 cs_xdef_t *
1445 cs_navsto_set_velocity_wall_by_value(cs_navsto_param_t    *nsp,
1446                                      const char           *z_name,
1447                                      cs_real_t            *values);
1448 
1449 /*----------------------------------------------------------------------------*/
1450 /*!
1451  * \brief  Define the velocity field for an inlet boundary using a uniform
1452  *         value
1453  *
1454  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1455  * \param[in]      z_name    name of the associated zone (if NULL or "" all
1456  *                           boundary faces are considered)
1457  * \param[in]      values    array of three real values
1458  *
1459  * \return a pointer to the new \ref cs_xdef_t structure
1460  */
1461 /*----------------------------------------------------------------------------*/
1462 
1463 cs_xdef_t *
1464 cs_navsto_set_velocity_inlet_by_value(cs_navsto_param_t    *nsp,
1465                                       const char           *z_name,
1466                                       cs_real_t            *values);
1467 
1468 /*----------------------------------------------------------------------------*/
1469 /*!
1470  * \brief  Define the velocity field for an inlet boundary using an analytical
1471  *         function
1472  *
1473  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1474  * \param[in]      z_name    name of the associated zone (if NULL or "" all
1475  *                           boundary faces are considered)
1476  * \param[in]      ana       pointer to an analytical function
1477  * \param[in]      input     NULL or pointer to a structure cast on-the-fly
1478  *
1479  * \return a pointer to the new \ref cs_xdef_t structure
1480  */
1481 /*----------------------------------------------------------------------------*/
1482 
1483 cs_xdef_t *
1484 cs_navsto_set_velocity_inlet_by_analytic(cs_navsto_param_t    *nsp,
1485                                          const char           *z_name,
1486                                          cs_analytic_func_t   *ana,
1487                                          void                 *input);
1488 
1489 /*----------------------------------------------------------------------------*/
1490 /*!
1491  * \brief  Define the velocity field for an inlet boundary using an array
1492  *         of values
1493  *
1494  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1495  * \param[in]      z_name    name of the associated zone (if NULL or "" all
1496  *                           boundary faces are considered)
1497  * \param[in]      loc       information to know where are located values
1498  * \param[in]      array     pointer to an array
1499  * \param[in]      is_owner  transfer the lifecycle to the cs_xdef_t structure
1500  *                           (true or false)
1501  * \param[in]      index     optional pointer to the array index
1502  *
1503  * \return a pointer to the new \ref cs_xdef_t structure
1504  */
1505 /*----------------------------------------------------------------------------*/
1506 
1507 cs_xdef_t *
1508 cs_navsto_set_velocity_inlet_by_array(cs_navsto_param_t    *nsp,
1509                                       const char           *z_name,
1510                                       cs_flag_t             loc,
1511                                       cs_real_t            *array,
1512                                       bool                  is_owner,
1513                                       cs_lnum_t            *index);
1514 
1515 /*----------------------------------------------------------------------------*/
1516 /*!
1517  * \brief  Define the velocity field for an inlet boundary using a DoF function
1518  *
1519  * \param[in]  nsp         pointer to a \ref cs_navsto_param_t structure
1520  * \param[in]  z_name      name of the associated zone (if NULL or "" all
1521  *                         boundary faces are considered)
1522  * \param[in]  func        pointer to a \ref cs_dof_function_t
1523  * \param[in]  func_input  NULL or pointer to a structure cast on-the-fly
1524  *
1525  * \return a pointer to the new \ref cs_xdef_t structure
1526  */
1527 /*----------------------------------------------------------------------------*/
1528 
1529 cs_xdef_t *
1530 cs_navsto_set_velocity_inlet_by_dof_func(cs_navsto_param_t    *nsp,
1531                                          const char           *z_name,
1532                                          cs_dof_func_t        *func,
1533                                          void                 *func_input);
1534 
1535 /*----------------------------------------------------------------------------*/
1536 /*!
1537  * \brief  Define a new source term structure defined by an analytical function
1538  *
1539  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1540  * \param[in]      z_name    name of the associated zone (if NULL or "" all
1541  *                           cells are considered)
1542  * \param[in]      ana       pointer to an analytical function
1543  * \param[in]      input     NULL or pointer to a structure cast on-the-fly
1544  *
1545  * \return a pointer to the new \ref cs_xdef_t structure
1546  */
1547 /*----------------------------------------------------------------------------*/
1548 
1549 cs_xdef_t *
1550 cs_navsto_add_source_term_by_analytic(cs_navsto_param_t    *nsp,
1551                                       const char           *z_name,
1552                                       cs_analytic_func_t   *ana,
1553                                       void                 *input);
1554 
1555 /*----------------------------------------------------------------------------*/
1556 /*!
1557  * \brief  Define a new source term structure defined by a constant value
1558  *
1559  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1560  * \param[in]      z_name    name of the associated zone (if NULL or "" all
1561  *                           cells are considered)
1562  * \param[in]      val       pointer to the value to set
1563  *
1564  * \return a pointer to the new \ref cs_xdef_t structure
1565  */
1566 /*----------------------------------------------------------------------------*/
1567 
1568 cs_xdef_t *
1569 cs_navsto_add_source_term_by_val(cs_navsto_param_t    *nsp,
1570                                  const char           *z_name,
1571                                  cs_real_t            *val);
1572 
1573 /*----------------------------------------------------------------------------*/
1574 /*!
1575  * \brief  Define a new source term structure defined by an array
1576  *
1577  * \param[in]      nsp       pointer to a \ref cs_navsto_param_t structure
1578  * \param[in]      z_name    name of the associated zone (if NULL or "" all
1579  *                           cells are considered)
1580  * \param[in]      loc       information to know where are located values
1581  * \param[in]      array     pointer to an array
1582  * \param[in]      is_owner  transfer the lifecycle to the cs_xdef_t structure
1583  *                           (true or false)
1584  * \param[in]      index     optional pointer to the array index
1585  *
1586  * \return a pointer to the new \ref cs_xdef_t structure
1587  */
1588 /*----------------------------------------------------------------------------*/
1589 
1590 cs_xdef_t *
1591 cs_navsto_add_source_term_by_array(cs_navsto_param_t    *nsp,
1592                                    const char           *z_name,
1593                                    cs_flag_t             loc,
1594                                    cs_real_t            *array,
1595                                    bool                  is_owner,
1596                                    cs_lnum_t            *index);
1597 
1598 /*----------------------------------------------------------------------------*/
1599 /*!
1600  * \brief  Add a advection field for the Oseen problem
1601  *
1602  * \param[in, out]    nsp        pointer to a \ref cs_navsto_param_t
1603  * \param[in, out]    adv_fld    pointer to a \ref cs_adv_field_t
1604  */
1605 /*----------------------------------------------------------------------------*/
1606 
1607 void
1608 cs_navsto_add_oseen_field(cs_navsto_param_t   *nsp,
1609                           cs_adv_field_t      *adv_fld);
1610 
1611 /*----------------------------------------------------------------------------*/
1612 
1613 END_C_DECLS
1614 
1615 #endif /* __CS_NAVSTO_PARAM_H__ */
1616