8.0
general documentation
Loading...
Searching...
No Matches
cs_cdofb_navsto.h
Go to the documentation of this file.
1#ifndef __CS_CDOFB_NAVSTO_H__
2#define __CS_CDOFB_NAVSTO_H__
3
4/*============================================================================
5 * Functions shared among all face-based schemes for the discretization of the
6 * Navier--Stokes system
7 *============================================================================*/
8
9/*
10 This file is part of code_saturne, a general-purpose CFD tool.
11
12 Copyright (C) 1998-2023 EDF S.A.
13
14 This program is free software; you can redistribute it and/or modify it under
15 the terms of the GNU General Public License as published by the Free Software
16 Foundation; either version 2 of the License, or (at your option) any later
17 version.
18
19 This program is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
21 FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
22 details.
23
24 You should have received a copy of the GNU General Public License along with
25 this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
26 Street, Fifth Floor, Boston, MA 02110-1301, USA.
27*/
28
29/*----------------------------------------------------------------------------*/
30
31#include "cs_defs.h"
32
33/*----------------------------------------------------------------------------
34 * Standard C library headers
35 *----------------------------------------------------------------------------*/
36
37/*----------------------------------------------------------------------------
38 * Local headers
39 *----------------------------------------------------------------------------*/
40
41#include "cs_base.h"
42#include "cs_cdo_connect.h"
43#include "cs_cdo_quantities.h"
44#include "cs_field.h"
45#include "cs_iter_algo.h"
46#include "cs_math.h"
47#include "cs_matrix.h"
48#include "cs_mesh.h"
49#include "cs_navsto_param.h"
50#include "cs_time_plot.h"
51#include "cs_time_step.h"
52#include "cs_sdm.h"
53
54/*----------------------------------------------------------------------------*/
55
57
58/*============================================================================
59 * Macro definitions
60 *============================================================================*/
61
62/*============================================================================
63 * Type definitions
64 *============================================================================*/
65
70
96
97
106
107typedef struct {
108
114
116
127
129
138
140
158
159 cs_boundary_type_t *bf_type; /* Size: n_fc */
160 cs_real_t *pressure_bc_val; /* Size: n_fc */
161
165
167
168/*----------------------------------------------------------------------------*/
181/*----------------------------------------------------------------------------*/
182
183typedef void
185 const cs_cell_mesh_t *cm,
186 const cs_cdofb_navsto_builder_t *nsb,
187 cs_cell_sys_t *csys);
188
189/*============================================================================
190 * Static inline public function prototypes
191 *============================================================================*/
192
193/*----------------------------------------------------------------------------*/
202/*----------------------------------------------------------------------------*/
203
204static inline void
206 cs_real_t div[])
207{
208 /* D(\hat{u}) = \frac{1}{|c|} \sum_{f_c} \iota_{fc} u_f.f
209 * But, when integrating [[ p, q ]]_{P_c} = |c| p_c q_c
210 * Thus, the volume in the divergence drops
211 */
212
213 for (short int f = 0; f < cm->n_fc; f++) {
214
215 const cs_quant_t pfq = cm->face[f];
216 const cs_real_t i_f = cm->f_sgn[f] * pfq.meas;
217
218 cs_real_t *_div_f = div + 3*f;
219 _div_f[0] = i_f * pfq.unitv[0];
220 _div_f[1] = i_f * pfq.unitv[1];
221 _div_f[2] = i_f * pfq.unitv[2];
222
223 } /* Loop on cell faces */
224}
225
226/*============================================================================
227 * Public function prototypes
228 *============================================================================*/
229
230/*----------------------------------------------------------------------------*/
236/*----------------------------------------------------------------------------*/
237
238void
240
241/*----------------------------------------------------------------------------*/
250/*----------------------------------------------------------------------------*/
251
254 const cs_cdo_connect_t *connect);
255
256/*----------------------------------------------------------------------------*/
262/*----------------------------------------------------------------------------*/
263
264void
266
267/*----------------------------------------------------------------------------*/
279/*----------------------------------------------------------------------------*/
280
281void
283 const cs_navsto_param_t *nsp,
284 const cs_cell_mesh_t *cm,
285 const cs_cell_sys_t *csys,
286 const cs_cdo_bc_face_t *pr_bc,
287 const cs_boundary_type_t *bf_type,
289
290/*----------------------------------------------------------------------------*/
302/*----------------------------------------------------------------------------*/
303
304void
306 const cs_cdo_quantities_t *quant,
307 const cs_real_t *face_vel,
308 cs_real_t *mass_flux);
309
310/*----------------------------------------------------------------------------*/
323/*----------------------------------------------------------------------------*/
324
325double
327 const cs_cdo_quantities_t *quant,
328 const cs_adjacency_t *c2f,
329 const cs_real_t *f_vals);
330
331/*----------------------------------------------------------------------------*/
343/*----------------------------------------------------------------------------*/
344
345void
347 const cs_cdo_connect_t *connect,
348 const cs_cdo_quantities_t *quant,
349 const cs_time_step_t *ts,
350 const cs_navsto_param_t *nsp,
351 const cs_real_t *p_cell,
352 cs_real_t *p_face);
353
354/*----------------------------------------------------------------------------*/
364/*----------------------------------------------------------------------------*/
365
366void
367cs_cdofb_navsto_add_grad_div(short int n_fc,
368 const cs_real_t zeta,
369 const cs_real_t div[],
370 cs_sdm_t *mat);
371
372/*----------------------------------------------------------------------------*/
381/*----------------------------------------------------------------------------*/
382
383void
385 const cs_cdo_quantities_t *quant,
386 const cs_time_step_t *ts,
387 cs_field_t *pr);
388
389/*----------------------------------------------------------------------------*/
399/*----------------------------------------------------------------------------*/
400
401void
403 const cs_cdo_connect_t *connect,
404 const cs_time_step_t *ts,
405 cs_real_t *pr_f);
406
407/*----------------------------------------------------------------------------*/
416/*----------------------------------------------------------------------------*/
417
418void
420 const cs_cdo_quantities_t *quant,
421 cs_real_t values[]);
422
423/*----------------------------------------------------------------------------*/
431/*----------------------------------------------------------------------------*/
432
433void
435 cs_real_t values[]);
436
437/*----------------------------------------------------------------------------*/
466/*----------------------------------------------------------------------------*/
467
468void
470 const cs_mesh_t *mesh,
471 const cs_cdo_quantities_t *quant,
472 const cs_cdo_connect_t *connect,
473 const cs_time_step_t *ts,
474 cs_time_plot_t *time_plotter,
475 const cs_adv_field_t *adv_field,
476 const cs_real_t *mass_flux,
477 const cs_real_t *p_cell,
478 const cs_real_t *u_cell,
479 const cs_real_t *u_face);
480
481/*----------------------------------------------------------------------------*/
496/*----------------------------------------------------------------------------*/
497
498void
500 const cs_equation_param_t *eqp,
501 const cs_cell_mesh_t *cm,
502 const cs_property_data_t *pty,
504 cs_cell_sys_t *csys);
505
506/*----------------------------------------------------------------------------*/
523/*----------------------------------------------------------------------------*/
524
525void
527 const cs_equation_param_t *eqp,
528 const cs_cell_mesh_t *cm,
529 const cs_property_data_t *pty,
531 cs_cell_sys_t *csys);
532
533/*----------------------------------------------------------------------------*/
550/*----------------------------------------------------------------------------*/
551
552void
554 const cs_equation_param_t *eqp,
555 const cs_cell_mesh_t *cm,
556 const cs_property_data_t *pty,
558 cs_cell_sys_t *csys);
559
560/*----------------------------------------------------------------------------*/
577/*----------------------------------------------------------------------------*/
578
579void
581 const cs_equation_param_t *eqp,
582 const cs_cell_mesh_t *cm,
583 const cs_property_data_t *pty,
585 cs_cell_sys_t *csys);
586
587/*----------------------------------------------------------------------------*/
603/*----------------------------------------------------------------------------*/
604
605void
606cs_cdofb_symmetry(short int fb,
607 const cs_equation_param_t *eqp,
608 const cs_cell_mesh_t *cm,
609 const cs_property_data_t *pty,
611 cs_cell_sys_t *csys);
612
613/*----------------------------------------------------------------------------*/
626/*----------------------------------------------------------------------------*/
627
628void
629cs_cdofb_fixed_wall(short int fb,
630 const cs_equation_param_t *eqp,
631 const cs_cell_mesh_t *cm,
632 const cs_property_data_t *pty,
634 cs_cell_sys_t *csys);
635
636/*----------------------------------------------------------------------------*/
649/*----------------------------------------------------------------------------*/
650
653 const cs_real_t *pre_iterate,
654 cs_real_t *cur_iterate,
655 cs_iter_algo_t *algo);
656
657/*----------------------------------------------------------------------------*/
666/*----------------------------------------------------------------------------*/
667
668void
670 cs_cdofb_navsto_source_t **p_func);
671
672/*----------------------------------------------------------------------------*/
684/*----------------------------------------------------------------------------*/
685
686void
688 const cs_cell_mesh_t *cm,
689 const cs_cdofb_navsto_builder_t *nsb,
690 cs_cell_sys_t *csys);
691
692/*----------------------------------------------------------------------------*/
705/*----------------------------------------------------------------------------*/
706
707void
709 const cs_cell_mesh_t *cm,
710 const cs_cdofb_navsto_builder_t *nsb,
711 cs_cell_sys_t *csys);
712
713/*----------------------------------------------------------------------------*/
726/*----------------------------------------------------------------------------*/
727
728void
730 const cs_cell_mesh_t *cm,
731 const cs_cdofb_navsto_builder_t *nsb,
732 cs_cell_sys_t *csys);
733
734/*----------------------------------------------------------------------------*/
747/*----------------------------------------------------------------------------*/
748
749void
751 const cs_cell_mesh_t *cm,
752 const cs_cdofb_navsto_builder_t *nsb,
753 cs_cell_sys_t *csys);
754
755/*----------------------------------------------------------------------------*/
767/*----------------------------------------------------------------------------*/
768
769void
771 const cs_lnum_t *elt_ids,
772 bool dense_output,
773 void *input,
774 cs_real_t *retval);
775
776/*----------------------------------------------------------------------------*/
777
779
780#endif /* __CS_CDOFB_NAVSTO_H__ */
int cs_boundary_type_t
Definition cs_boundary.h:69
void cs_cdofb_navsto_gravity_term(const cs_navsto_param_t *nsp, const cs_cell_mesh_t *cm, const cs_cdofb_navsto_builder_t *nsb, cs_cell_sys_t *csys)
Take into account the gravity effects. Compute and add the source term to the local RHS....
Definition cs_cdofb_navsto.c:2050
void cs_cdofb_block_dirichlet_pena(short int f, const cs_equation_param_t *eqp, const cs_cell_mesh_t *cm, const cs_property_data_t *pty, cs_cell_builder_t *cb, cs_cell_sys_t *csys)
Take into account a Dirichlet BCs on the three velocity components. For instance for a velocity inlet...
Definition cs_cdofb_navsto.c:1505
void cs_cdofb_navsto_extra_op(const cs_navsto_param_t *nsp, const cs_mesh_t *mesh, const cs_cdo_quantities_t *quant, const cs_cdo_connect_t *connect, const cs_time_step_t *ts, cs_time_plot_t *time_plotter, const cs_adv_field_t *adv_field, const cs_real_t *mass_flux, const cs_real_t *p_cell, const cs_real_t *u_cell, const cs_real_t *u_face)
Perform extra-operation related to Fb schemes when solving Navier-Stokes. Computation of the followin...
Definition cs_cdofb_navsto.c:954
void cs_cdofb_navsto_add_grad_div(short int n_fc, const cs_real_t zeta, const cs_real_t div[], cs_sdm_t *mat)
Add the grad-div part to the local matrix (i.e. for the current cell)
Definition cs_cdofb_navsto.c:574
void cs_cdofb_navsto_boussinesq_at_face(const cs_navsto_param_t *nsp, const cs_cell_mesh_t *cm, const cs_cdofb_navsto_builder_t *nsb, cs_cell_sys_t *csys)
Take into account the buoyancy force with the Boussinesq approx. Compute and add the source term to t...
Definition cs_cdofb_navsto.c:2151
void cs_cdofb_navsto_compute_face_pressure(const cs_mesh_t *mesh, const cs_cdo_connect_t *connect, const cs_cdo_quantities_t *quant, const cs_time_step_t *ts, const cs_navsto_param_t *nsp, const cs_real_t *p_cell, cs_real_t *p_face)
Compute an estimation of the pressure at faces.
Definition cs_cdofb_navsto.c:502
void cs_cdofb_navsto_source_t(const cs_navsto_param_t *nsp, const cs_cell_mesh_t *cm, const cs_cdofb_navsto_builder_t *nsb, cs_cell_sys_t *csys)
Compute and add a source term to the local RHS. This is a special treatment to enable source involvin...
Definition cs_cdofb_navsto.h:184
void cs_cdofb_block_dirichlet_wsym(short int fb, const cs_equation_param_t *eqp, const cs_cell_mesh_t *cm, const cs_property_data_t *pty, cs_cell_builder_t *cb, cs_cell_sys_t *csys)
Take into account a Dirichlet BCs on the three velocity components. For instance for a velocity inlet...
Definition cs_cdofb_navsto.c:1660
void cs_cdofb_symmetry(short int fb, const cs_equation_param_t *eqp, const cs_cell_mesh_t *cm, const cs_property_data_t *pty, cs_cell_builder_t *cb, cs_cell_sys_t *csys)
Take into account a boundary defined as 'symmetry' (treated as a sliding BCs on the three velocity co...
Definition cs_cdofb_navsto.c:1768
cs_sles_convergence_state_t cs_cdofb_navsto_nl_algo_cvg(cs_param_nl_algo_t nl_algo_type, const cs_real_t *pre_iterate, cs_real_t *cur_iterate, cs_iter_algo_t *algo)
Test if one has to do one more non-linear iteration. Test if performed on the relative norm on the in...
Definition cs_cdofb_navsto.c:1925
void cs_cdofb_navsto_init_pressure(const cs_navsto_param_t *nsp, const cs_cdo_quantities_t *quant, const cs_time_step_t *ts, cs_field_t *pr)
Initialize the pressure values.
Definition cs_cdofb_navsto.c:648
void cs_cdofb_block_dirichlet_alge(short int f, const cs_equation_param_t *eqp, const cs_cell_mesh_t *cm, const cs_property_data_t *pty, cs_cell_builder_t *cb, cs_cell_sys_t *csys)
Take into account a Dirichlet BCs on the three velocity components. For instance for a velocity inlet...
Definition cs_cdofb_navsto.c:1400
void cs_cdofb_fixed_wall(short int fb, const cs_equation_param_t *eqp, const cs_cell_mesh_t *cm, const cs_property_data_t *pty, cs_cell_builder_t *cb, cs_cell_sys_t *csys)
Take into account a wall BCs by a weak enforcement using Nitsche technique plus a symmetric treatment...
Definition cs_cdofb_navsto.c:1881
cs_cdofb_navsto_boussinesq_type_t
Type of algorithm to compute the Boussinesq approximation.
Definition cs_cdofb_navsto.h:71
@ CS_CDOFB_NAVSTO_BOUSSINESQ_CELL_DOF
Boussinesq approximation relyong on a cell contribution.
Definition cs_cdofb_navsto.h:82
@ CS_CDOFB_NAVSTO_BOUSSINESQ_FACE_DOF
Boussinesq approximation relyong on face contributions.
Definition cs_cdofb_navsto.h:93
void cs_cdofb_navsto_rescale_pressure_to_ref(const cs_navsto_param_t *nsp, const cs_cdo_quantities_t *quant, cs_real_t values[])
Update the pressure field in order to get a field with a mean-value equal to the reference value.
Definition cs_cdofb_navsto.c:851
void cs_cdofb_navsto_free_builder(cs_cdofb_navsto_builder_t *nsb)
Destroy the given cs_cdofb_navsto_builder_t structure.
Definition cs_cdofb_navsto.c:245
double cs_cdofb_navsto_cell_divergence(const cs_lnum_t c_id, const cs_cdo_quantities_t *quant, const cs_adjacency_t *c2f, const cs_real_t *f_vals)
Compute the divergence in a cell of a vector-valued array defined at faces (values are defined both a...
Definition cs_cdofb_navsto.c:463
void cs_cdofb_navsto_set_boussinesq_algo(cs_cdofb_navsto_boussinesq_type_t type)
Set the way to compute the Boussinesq approximation.
Definition cs_cdofb_navsto.c:198
void cs_cdofb_navsto_boussinesq_at_cell(const cs_navsto_param_t *nsp, const cs_cell_mesh_t *cm, const cs_cdofb_navsto_builder_t *nsb, cs_cell_sys_t *csys)
Take into account the buoyancy force with the Boussinesq approx. Compute and add the source term to t...
Definition cs_cdofb_navsto.c:2091
void cs_cdofb_navsto_stream_source_term(cs_lnum_t n_elts, const cs_lnum_t *elt_ids, bool dense_output, void *input, cs_real_t *retval)
Get the source term for computing the stream function. This relies on the prototype associated to the...
Definition cs_cdofb_navsto.c:2207
void cs_cdofb_navsto_init_face_pressure(const cs_navsto_param_t *nsp, const cs_cdo_connect_t *connect, const cs_time_step_t *ts, cs_real_t *pr_f)
Initialize the pressure values when the pressure is defined at faces.
Definition cs_cdofb_navsto.c:746
void cs_cdofb_navsto_define_builder(cs_real_t t_eval, const cs_navsto_param_t *nsp, const cs_cell_mesh_t *cm, const cs_cell_sys_t *csys, const cs_cdo_bc_face_t *pr_bc, const cs_boundary_type_t *bf_type, cs_cdofb_navsto_builder_t *nsb)
Set the members of the cs_cdofb_navsto_builder_t structure.
Definition cs_cdofb_navsto.c:269
void cs_cdofb_navsto_set_gravity_func(const cs_navsto_param_t *nsp, cs_cdofb_navsto_source_t **p_func)
Set the function pointer computing the source term in the momentum equation related to the gravity ef...
Definition cs_cdofb_navsto.c:2006
cs_cdofb_navsto_builder_t cs_cdofb_navsto_create_builder(const cs_navsto_param_t *nsp, const cs_cdo_connect_t *connect)
Create and allocate a local NavSto builder when Fb schemes are used.
Definition cs_cdofb_navsto.c:215
static void cs_cdofb_navsto_divergence_vect(const cs_cell_mesh_t *cm, cs_real_t div[])
Compute the divergence vector associated to the current cell. WARNING: mind that, differently form th...
Definition cs_cdofb_navsto.h:205
void cs_cdofb_block_dirichlet_weak(short int fb, const cs_equation_param_t *eqp, const cs_cell_mesh_t *cm, const cs_property_data_t *pty, cs_cell_builder_t *cb, cs_cell_sys_t *csys)
Take into account a Dirichlet BCs on the three velocity components. For instance for a velocity inlet...
Definition cs_cdofb_navsto.c:1573
void cs_cdofb_navsto_set_zero_mean_pressure(const cs_cdo_quantities_t *quant, cs_real_t values[])
Update the pressure field in order to get a field with a zero-mean average.
Definition cs_cdofb_navsto.c:889
void cs_cdofb_navsto_mass_flux(const cs_navsto_param_t *nsp, const cs_cdo_quantities_t *quant, const cs_real_t *face_vel, cs_real_t *mass_flux)
Compute the mass flux playing the role of the advection field in the Navier-Stokes equations One cons...
Definition cs_cdofb_navsto.c:421
void cs_cdofb_navsto_boussinesq_by_part(const cs_navsto_param_t *nsp, const cs_cell_mesh_t *cm, const cs_cdofb_navsto_builder_t *nsb, cs_cell_sys_t *csys)
Take into account the buoyancy force with the Boussinesq approx. Compute and add the source term to t...
#define BEGIN_C_DECLS
Definition cs_defs.h:509
double cs_real_t
Floating-point value.
Definition cs_defs.h:319
#define END_C_DECLS
Definition cs_defs.h:510
int cs_lnum_t
local mesh entity id
Definition cs_defs.h:313
cs_param_nl_algo_t
Class of non-linear iterative algorithm.
Definition cs_param_types.h:551
cs_sles_convergence_state_t
Definition cs_sles.h:56
struct _cs_time_plot_t cs_time_plot_t
Definition cs_time_plot.h:48
Definition mesh.f90:26
Definition cs_mesh_adjacencies.h:68
Definition cs_advection_field.h:150
Definition cs_cdo_bc.h:109
Definition cs_cdo_connect.h:61
Definition cs_cdo_quantities.h:137
Structure storing additional arrays related to the building of the Navier-Stokes system.
Definition cs_cdofb_navsto.h:107
cs_real_t * pressure_bc_val
Definition cs_cdofb_navsto.h:160
cs_real_t mass_rhs
Definition cs_cdofb_navsto.h:139
cs_real_t rho_c
Definition cs_cdofb_navsto.h:115
cs_real_t * div_op
Definition cs_cdofb_navsto.h:128
cs_boundary_type_t * bf_type
Definition cs_cdofb_navsto.h:159
Set of local and temporary buffers.
Definition cs_cdo_local.h:60
Set of local quantities and connectivities related to a mesh cell.
Definition cs_cdo_local.h:203
short int n_fc
Definition cs_cdo_local.h:238
short int * f_sgn
Definition cs_cdo_local.h:241
cs_quant_t * face
Definition cs_cdo_local.h:244
Set of arrays and local (small) dense matrices related to a mesh cell This is a key structure for bui...
Definition cs_cdo_local.h:147
Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources.
Definition cs_equation_param.h:192
Field descriptor.
Definition cs_field.h:130
Structure to handle the convergence of an iterative algorithm.
Definition cs_iter_algo.h:60
Definition cs_mesh.h:85
Structure storing the parameters related to the resolution of the Navier-Stokes system.
Definition cs_navsto_param.h:611
Structure storing the evaluation of a property and its related data.
Definition cs_property.h:211
Definition cs_cdo_quantities.h:129
double meas
Definition cs_cdo_quantities.h:131
double unitv[3]
Definition cs_cdo_quantities.h:132
time step descriptor
Definition cs_time_step.h:64