9.0
general documentation
Loading...
Searching...
No Matches
cs_convection_diffusion.h
Go to the documentation of this file.
1#ifndef __CS_CONVECTION_DIFFUSION_H__
2#define __CS_CONVECTION_DIFFUSION_H__
3
4/*============================================================================
5 * Convection-diffusion operators.
6 *============================================================================*/
7
8/*
9 This file is part of code_saturne, a general-purpose CFD tool.
10
11 Copyright (C) 1998-2025 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#include "base/cs_defs.h"
31
32/*----------------------------------------------------------------------------
33 * Local headers
34 *----------------------------------------------------------------------------*/
35
36#include "bft/bft_mem.h"
37
38#include "base/cs_base.h"
39#include "base/cs_dispatch.h"
40#include "base/cs_halo.h"
41#include "base/cs_math.h"
44
45/*----------------------------------------------------------------------------*/
46
48
49/*=============================================================================
50 * Macro definitions
51 *============================================================================*/
52
53/*============================================================================
54 * Type definition
55 *============================================================================*/
56
57/*----------------------------------------------------------------------------
58 * NVD/TVD Advection Scheme
59 *----------------------------------------------------------------------------*/
60
61typedef enum {
62
63 CS_NVD_GAMMA = 0, /* GAMMA */
64 CS_NVD_SMART = 1, /* SMART */
65 CS_NVD_CUBISTA = 2, /* CUBISTA */
66 CS_NVD_SUPERBEE = 3, /* SUPERBEE */
67 CS_NVD_MUSCL = 4, /* MUSCL */
68 CS_NVD_MINMOD = 5, /* MINMOD */
69 CS_NVD_CLAM = 6, /* CLAM */
70 CS_NVD_STOIC = 7, /* STOIC */
71 CS_NVD_OSHER = 8, /* OSHER */
72 CS_NVD_WASEB = 9, /* WASEB */
73 CS_NVD_VOF_HRIC = 10, /* M-HRIC for VOF */
74 CS_NVD_VOF_CICSAM = 11, /* M-CICSAM for VOF */
75 CS_NVD_VOF_STACS = 12, /* STACS for VOF */
76 CS_NVD_N_TYPES = 13 /* number of NVD schemes */
77
79
80/*============================================================================
81 * Global variables
82 *============================================================================*/
83
84/*============================================================================
85 * Public function prototypes for Fortran API
86 *============================================================================*/
87
88/*=============================================================================
89 * Public function prototypes
90 *============================================================================*/
91
92/*----------------------------------------------------------------------------
93 * Return pointer to slope test indicator field values if active.
94 *
95 * parameters:
96 * f_id <-- field id (or -1)
97 * eqp <-- equation parameters
98 *
99 * return:
100 * pointer to local values array, or NULL;
101 *----------------------------------------------------------------------------*/
102
103cs_real_t *
104cs_get_v_slope_test(int f_id,
105 const cs_equation_param_t eqp);
106
107/*----------------------------------------------------------------------------*/
116/*----------------------------------------------------------------------------*/
117
118void
120 int inc,
121 const cs_real_t rovsdt[]);
122
123/*----------------------------------------------------------------------------*/
170/*----------------------------------------------------------------------------*/
171
172void
174 int f_id,
175 const cs_equation_param_t eqp,
176 int icvflb,
177 int inc,
178 int imasac,
179 cs_real_t *pvar,
180 const cs_real_t *pvara,
181 const int icvfli[],
182 const cs_field_bc_coeffs_t *bc_coeffs,
183 const cs_real_t i_massflux[],
184 const cs_real_t b_massflux[],
185 const cs_real_t i_visc[],
186 const cs_real_t b_visc[],
187 cs_real_t *rhs,
188 cs_real_2_t i_flux[],
189 cs_real_t b_flux[]);
190
191/*----------------------------------------------------------------------------*/
192/*
193 * \brief Update face flux with convection contribution of a standard transport
194 * equation of a scalar field \f$ \varia \f$.
195 *
196 * \f[
197 * C_\ij = \dot{m}_\ij \left( \varia_\fij - \varia_\celli \right)
198 * \f]
199 *
200 * \param[in] idtvar indicator of the temporal scheme
201 * \param[in] f_id field id (or -1)
202 * \param[in] eqp equation parameters
203 * \param[in] icvflb global indicator of boundary convection flux
204 * - 0 upwind scheme at all boundary faces
205 * - 1 imposed flux at some boundary faces
206 * \param[in] inc indicator
207 * - 0 when solving an increment
208 * - 1 otherwise
209 * \param[in] imasac take mass accumulation into account?
210 * \param[in] pvar solved variable (current time step)
211 * \param[in] pvara solved variable (previous time step)
212 * \param[in] icvfli boundary face indicator array of convection flux
213 * - 0 upwind scheme
214 * - 1 imposed flux
215 * \param[in] bc_coeffs boundary condition structure for the variable
216 * \param[in] i_massflux mass flux at interior faces
217 * \param[in] b_massflux mass flux at boundary faces
218 * \param[in,out] i_conv_flux scalar convection flux at interior faces
219 * \param[in,out] b_conv_flux scalar convection flux at boundary faces
220 */
221/*----------------------------------------------------------------------------*/
222
223void
225 int f_id,
226 const cs_equation_param_t eqp,
227 int icvflb,
228 int inc,
229 int imasac,
230 cs_real_t *pvar,
231 const cs_real_t *pvara,
232 const int icvfli[],
233 const cs_field_bc_coeffs_t *bc_coeffs,
234 const cs_real_t i_massflux[],
235 const cs_real_t b_massflux[],
236 cs_real_2_t i_conv_flux[],
237 cs_real_t b_conv_flux[]);
238
239/*----------------------------------------------------------------------------*/
240/*
241 * \brief Add the explicit part of the convection/diffusion terms of a transport
242 * equation of a vector field \f$ \vect{\varia} \f$.
243 *
244 * More precisely, the right hand side \f$ \vect{Rhs} \f$ is updated as
245 * follows:
246 * \f[
247 * \vect{Rhs} = \vect{Rhs} - \sum_{\fij \in \Facei{\celli}} \left(
248 * \dot{m}_\ij \left( \vect{\varia}_\fij - \vect{\varia}_\celli \right)
249 * - \mu_\fij \gradt_\fij \vect{\varia} \cdot \vect{S}_\ij \right)
250 * \f]
251 *
252 * Remark:
253 * if ivisep = 1, then we also take \f$ \mu \transpose{\gradt\vect{\varia}}
254 * + \lambda \trace{\gradt\vect{\varia}} \f$, where \f$ \lambda \f$ is
255 * the secondary viscosity, i.e. usually \f$ -\frac{2}{3} \mu \f$.
256 *
257 * Warning:
258 * - \f$ \vect{Rhs} \f$ has already been initialized before calling bilsc!
259 * - mind the sign minus
260 *
261 * \param[in] idtvar indicator of the temporal scheme
262 * \param[in] f_id index of the current variable
263 * \param[in] eqp equation parameters
264 * \param[in] icvflb global indicator of boundary convection flux
265 * - 0 upwind scheme at all boundary faces
266 * - 1 imposed flux at some boundary faces
267 * \param[in] inc indicator
268 * - 0 when solving an increment
269 * - 1 otherwise
270 * \param[in] ivisep indicator to take \f$ \divv
271 * \left(\mu \gradt \transpose{\vect{a}} \right)
272 * -2/3 \grad\left( \mu \dive \vect{a} \right)\f$
273 * - 1 take into account,
274 * - 0 otherwise
275 * \param[in] imasac take mass accumulation into account?
276 * \param[in] pvar solved velocity (current time step)
277 * \param[in] pvara solved velocity (previous time step)
278 * \param[in] icvfli boundary face indicator array of convection flux
279 * - 0 upwind scheme
280 * - 1 imposed flux
281 * \param[in] bc_coeffs_v boundary conditions structure for the variable
282 * \param[in] bc_coeffs_solve_v sweep loop boundary conditions structure
283 * \param[in] i_massflux mass flux at interior faces
284 * \param[in] b_massflux mass flux at boundary faces
285 * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
286 * at interior faces for the r.h.s.
287 * \param[in] b_visc \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
288 * at border faces for the r.h.s.
289 * \param[in] i_secvis secondary viscosity at interior faces
290 * \param[in] b_secvis secondary viscosity at boundary faces
291 * \param[in] i_pvar velocity at interior faces
292 * \param[in] b_pvar velocity at boundary faces
293 * \param[in,out] rhs right hand side \f$ \vect{Rhs} \f$
294 */
295/*----------------------------------------------------------------------------*/
296
297void
299 int f_id,
300 const cs_equation_param_t eqp,
301 int icvflb,
302 int inc,
303 int ivisep,
304 int imasac,
305 cs_real_3_t *pvar,
306 const cs_real_3_t *pvara,
307 const int icvfli[],
308 const cs_field_bc_coeffs_t *bc_coeffs_v,
309 const cs_bc_coeffs_solve_t *bc_coeffs_solve_v,
310 const cs_real_t i_massflux[],
311 const cs_real_t b_massflux[],
312 const cs_real_t i_visc[],
313 const cs_real_t b_visc[],
314 const cs_real_t i_secvis[],
315 const cs_real_t b_secvis[],
316 cs_real_3_t *i_pvar,
317 cs_real_3_t *b_pvar,
318 cs_real_3_t *rhs);
319
320/*----------------------------------------------------------------------------*/
359/*----------------------------------------------------------------------------*/
360
361void
363 int f_id,
364 const cs_equation_param_t eqp,
365 int icvflb,
366 int inc,
367 int imasac,
368 cs_real_6_t *pvar,
369 const cs_real_6_t *pvara,
370 const cs_field_bc_coeffs_t *bc_coeffs_ts,
371 const cs_bc_coeffs_solve_t *bc_coeffs_solve_ts,
372 const cs_real_t i_massflux[],
373 const cs_real_t b_massflux[],
374 const cs_real_t i_visc[],
375 const cs_real_t b_visc[],
376 cs_real_6_t *rhs);
377
378/*----------------------------------------------------------------------------*/
413/*----------------------------------------------------------------------------*/
414
415void
417 int f_id,
418 const cs_equation_param_t eqp,
419 int inc,
420 int imasac,
421 cs_real_t * pvar,
422 const cs_real_t * pvara,
423 const cs_field_bc_coeffs_t *bc_coeffs,
424 const cs_real_t i_massflux[],
425 const cs_real_t b_massflux[],
426 const cs_real_t i_visc[],
427 const cs_real_t b_visc[],
428 const cs_real_t xcpp[],
429 cs_real_t * rhs);
430
431/*----------------------------------------------------------------------------*/
432/*
433 * \brief Add the explicit part of the diffusion terms with a symmetric tensor
434 * diffusivity for a transport equation of a scalar field \f$ \varia \f$.
435 *
436 * More precisely, the right hand side \f$ Rhs \f$ is updated as
437 * follows:
438 * \f[
439 * Rhs = Rhs - \sum_{\fij \in \Facei{\celli}} \left(
440 * - \tens{\mu}_\fij \gradv_\fij \varia \cdot \vect{S}_\ij \right)
441 * \f]
442 *
443 * Warning:
444 * - \f$ Rhs \f$ has already been initialized before
445 * calling cs_anisotropic_diffusion_scalar!
446 * - mind the sign minus
447 *
448 * \param[in] idtvar indicator of the temporal scheme
449 * \param[in] f_id index of the current variable
450 * \param[in] eqp equation parameters
451 * \param[in] inc indicator
452 * - 0 when solving an increment
453 * - 1 otherwise
454 * \param[in] pvar solved variable (current time step)
455 * \param[in] pvara solved variable (previous time step)
456 * \param[in] bc_coeffs boundary condition structure for the variable
457 * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
458 * at interior faces for the r.h.s.
459 * \param[in] b_visc \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
460 * at border faces for the r.h.s.
461 * \param[in] viscel symmetric cell tensor \f$ \tens{\mu}_\celli \f$
462 * \param[in] weighf internal face weight between cells i j in case
463 * of tensor diffusion
464 * \param[in] weighb boundary face weight for cells i in case
465 * of tensor diffusion
466 * \param[in,out] rhs right hand side \f$ \vect{Rhs} \f$
467 */
468/*----------------------------------------------------------------------------*/
469
470void
472 int f_id,
473 const cs_equation_param_t eqp,
474 int inc,
475 cs_real_t * pvar,
476 const cs_real_t * pvara,
477 const cs_field_bc_coeffs_t *bc_coeffs,
478 const cs_real_t i_visc[],
479 const cs_real_t b_visc[],
480 cs_real_6_t * viscel,
481 const cs_real_2_t weighf[],
482 const cs_real_t weighb[],
483 cs_real_t * rhs);
484
485/*-----------------------------------------------------------------------------*/
486/*
487 * \brief Add explicit part of the terms of diffusion by a left-multiplying
488 * symmetric tensorial diffusivity for a transport equation of a vector field
489 * \f$ \vect{\varia} \f$.
490 *
491 * More precisely, the right hand side \f$ \vect{Rhs} \f$ is updated as
492 * follows:
493 * \f[
494 * \vect{Rhs} = \vect{Rhs} - \sum_{\fij \in \Facei{\celli}} \left(
495 * - \gradt_\fij \vect{\varia} \tens{\mu}_\fij \cdot \vect{S}_\ij \right)
496 * \f]
497 *
498 * Remark:
499 * if ivisep = 1, then we also take \f$ \mu \transpose{\gradt\vect{\varia}}
500 * + \lambda \trace{\gradt\vect{\varia}} \f$, where \f$ \lambda \f$ is
501 * the secondary viscosity, i.e. usually \f$ -\frac{2}{3} \mu \f$.
502 *
503 * Warning:
504 * - \f$ \vect{Rhs} \f$ has already been initialized before calling the present
505 * function
506 * - mind the sign minus
507 *
508 * \param[in] idtvar indicator of the temporal scheme
509 * \param[in] f_id index of the current variable
510 * \param[in] eqp equation parameters
511 * \param[in] inc indicator
512 * - 0 when solving an increment
513 * - 1 otherwise
514 * \param[in] ivisep indicator to take \f$ \divv
515 * \left(\mu \gradt \transpose{\vect{a}} \right)
516 * -2/3 \grad\left( \mu \dive \vect{a} \right)\f$
517 * - 1 take into account,
518 * \param[in] pvar solved variable (current time step)
519 * \param[in] pvara solved variable (previous time step)
520 * \param[in] bc_coeffs_v boundary conditions structure for the variable
521 * \param[in] bc_coeffs_solve_v sweep loop boundary conditions structure
522 * \param[in] i_visc \f$ \tens{\mu}_\fij \dfrac{S_\fij}{\ipf\jpf} \f$
523 * at interior faces for the r.h.s.
524 * \param[in] b_visc \f$ \dfrac{S_\fib}{\ipf \centf} \f$
525 * at border faces for the r.h.s.
526 * \param[in] i_secvis secondary viscosity at interior faces
527 * \param[in,out] rhs right hand side \f$ \vect{Rhs} \f$
528 */
529/*----------------------------------------------------------------------------*/
530
531void
533 int f_id,
534 const cs_equation_param_t eqp,
535 int inc,
536 int ivisep,
537 cs_real_3_t *pvar,
538 const cs_real_3_t *pvara,
539 const cs_field_bc_coeffs_t *bc_coeffs_v,
540 const cs_bc_coeffs_solve_t *bc_coeffs_solve_v,
541 const cs_real_33_t i_visc[],
542 const cs_real_t b_visc[],
543 const cs_real_t i_secvis[],
544 cs_real_3_t *rhs);
545
546/*-----------------------------------------------------------------------------*/
547/*
548 * \brief Add explicit part of the terms of diffusion by a right-multiplying
549 * symmetric tensorial diffusivity for a transport equation of a vector field
550 * \f$ \vect{\varia} \f$.
551 *
552 * More precisely, the right hand side \f$ \vect{Rhs} \f$ is updated as
553 * follows:
554 * \f[
555 * \vect{Rhs} = \vect{Rhs} - \sum_{\fij \in \Facei{\celli}} \left(
556 * - \gradt_\fij \vect{\varia} \tens{\mu}_\fij \cdot \vect{S}_\ij \right)
557 * \f]
558 *
559 * Warning:
560 * - \f$ \vect{Rhs} \f$ has already been initialized before calling the present
561 * function
562 * - mind the sign minus
563 *
564 * \param[in] idtvar indicator of the temporal scheme
565 * \param[in] f_id index of the current variable
566 * \param[in] eqp equation parameters
567 * \param[in] inc indicator
568 * - 0 when solving an increment
569 * - 1 otherwise
570 * \param[in] pvar solved variable (current time step)
571 * \param[in] pvara solved variable (previous time step)
572 * \param[in] bc_coeffs_v boundary condition structure for the variable
573 * \param[in] bc_coeffs_solve_v sweep loop boundary conditions structure
574 * \param[in] i_visc \f$ \tens{\mu}_\fij \dfrac{S_\fij}{\ipf\jpf} \f$
575 * at interior faces for the r.h.s.
576 * \param[in] b_visc \f$ \dfrac{S_\fib}{\ipf \centf} \f$
577 * at border faces for the r.h.s.
578 * \param[in] viscel symmetric cell tensor \f$ \tens{\mu}_\celli \f$
579 * \param[in] weighf internal face weight between cells i j in case
580 * of tensor diffusion
581 * \param[in] weighb boundary face weight for cells i in case
582 * of tensor diffusion
583 * \param[in,out] rhs right hand side \f$ \vect{Rhs} \f$
584 */
585/*----------------------------------------------------------------------------*/
586
587void
589 int f_id,
590 const cs_equation_param_t eqp,
591 int inc,
592 cs_real_3_t *pvar,
593 const cs_real_3_t *pvara,
594 const cs_field_bc_coeffs_t *bc_coeffs_v,
595 const cs_bc_coeffs_solve_t *bc_coeffs_solve_v,
596 const cs_real_t i_visc[],
597 const cs_real_t b_visc[],
598 cs_real_6_t *viscel,
599 const cs_real_2_t weighf[],
600 const cs_real_t weighb[],
601 cs_real_3_t *rhs);
602
603/*----------------------------------------------------------------------------*/
604/*
605 * \brief Add the explicit part of the diffusion terms with a symmetric tensor
606 * diffusivity for a transport equation of a scalar field \f$ \varia \f$.
607 *
608 * More precisely, the right hand side \f$ Rhs \f$ is updated as
609 * follows:
610 * \f[
611 * Rhs = Rhs - \sum_{\fij \in \Facei{\celli}} \left(
612 * - \tens{\mu}_\fij \gradv_\fij \varia \cdot \vect{S}_\ij \right)
613 * \f]
614 *
615 * Warning:
616 * - \f$ Rhs \f$ has already been initialized before
617 * calling cs_anisotropic_diffusion_scalar!
618 * - mind the sign minus
619 *
620 * \param[in] idtvar indicator of the temporal scheme
621 * \param[in] f_id index of the current variable
622 * \param[in] eqp equation parameters
623 * \param[in] inc indicator
624 * - 0 when solving an increment
625 * - 1 otherwise
626 * \param[in] pvar solved variable (current time step)
627 * \param[in] pvara solved variable (previous time step)
628 * \param[in] bc_coeffs_ts boundary condition structure for the variable
629 * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
630 * at interior faces for the r.h.s.
631 * \param[in] b_visc \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
632 * at border faces for the r.h.s.
633 * \param[in] viscel symmetric cell tensor \f$ \tens{\mu}_\celli \f$
634 * \param[in] weighf internal face weight between cells i j in case
635 * of tensor diffusion
636 * \param[in] weighb boundary face weight for cells i in case
637 * of tensor diffusion
638 * \param[in,out] rhs right hand side \f$ \vect{Rhs} \f$
639 */
640/*----------------------------------------------------------------------------*/
641
642void
644 int f_id,
645 const cs_equation_param_t eqp,
646 int inc,
647 cs_real_6_t *pvar,
648 const cs_real_6_t *pvara,
649 const cs_field_bc_coeffs_t *bc_coeffs_ts,
650 const cs_bc_coeffs_solve_t *bc_coeffs_solve_ts,
651 const cs_real_t i_visc[],
652 const cs_real_t b_visc[],
653 cs_real_6_t *viscel,
654 const cs_real_2_t weighf[],
655 const cs_real_t weighb[],
656 cs_real_6_t *rhs);
657
658/*----------------------------------------------------------------------------*/
659/*
660 * \brief Update the face mass flux with the face pressure (or pressure
661 * increment, or pressure double increment) gradient.
662 *
663 * \f[
664 * \dot{m}_\ij = \dot{m}_\ij
665 * - \Delta t \grad_\fij \delta p \cdot \vect{S}_\ij
666 * \f]
667 *
668 * Please refer to the
669 * <a href="../../theory.pdf#cs_face_diffusion_potential">
670 <b>cs_face_diffusion_potential/cs_diffusion_potential</b></a>
671 * section of the theory guide for more information.
672 *
673 * \param[in] f_id field id (or -1)
674 * \param[in] m pointer to mesh
675 * \param[in] fvq pointer to finite volume quantities
676 * \param[in] init indicator
677 * - 1 initialize the mass flux to 0
678 * - 0 otherwise
679 * \param[in] inc indicator
680 * - 0 when solving an increment
681 * - 1 otherwise
682 * \param[in] imrgra indicator
683 * - 0 iterative gradient
684 * - 1 least squares gradient
685 * \param[in] nswrgp number of reconstruction sweeps for the
686 * gradients
687 * \param[in] imligp clipping gradient method
688 * - < 0 no clipping
689 * - = 0 thank to neighbooring gradients
690 * - = 1 thank to the mean gradient
691 * \param[in] iphydp hydrostatic pressure indicator
692 * \param[in] iwgrp indicator
693 * - 1 weight gradient by vicosity*porosity
694 * - weighting determined by field options
695 * \param[in] iwarnp verbosity
696 * \param[in] epsrgp relative precision for the gradient
697 * reconstruction
698 * \param[in] climgp clipping coeffecient for the computation of
699 * the gradient
700 * \param[in] frcxt body force creating the hydrostatic pressure
701 * \param[in] pvar solved variable (current time step)
702 * \param[in] bc_coeffs boundary condition structure for the variable
703 * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
704 * at interior faces for the r.h.s.
705 * \param[in] b_visc \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
706 * at border faces for the r.h.s.
707 * \param[in] visel viscosity by cell
708 * \param[in,out] i_massflux mass flux at interior faces
709 * \param[in,out] b_massflux mass flux at boundary faces
710 */
711/*----------------------------------------------------------------------------*/
712
713void
715 const cs_mesh_t *m,
717 int init,
718 int inc,
719 int imrgra,
720 int nswrgp,
721 int imligp,
722 int iphydp,
723 int iwgrp,
724 int iwarnp,
725 double epsrgp,
726 double climgp,
727 cs_real_3_t *frcxt,
728 cs_real_t *pvar,
729 const cs_field_bc_coeffs_t *bc_coeffs,
730 const cs_real_t i_visc[],
731 const cs_real_t b_visc[],
732 cs_real_t *visel,
733 cs_real_t *i_massflux,
734 cs_real_t *b_massflux);
735
736/*----------------------------------------------------------------------------*/
737/*
738 * \brief Add the explicit part of the pressure gradient term to the mass flux
739 * in case of anisotropic diffusion of the pressure field \f$ P \f$.
740 *
741 * More precisely, the mass flux side \f$ \dot{m}_\fij \f$ is updated as
742 * follows:
743 * \f[
744 * \dot{m}_\fij = \dot{m}_\fij -
745 * \left( \tens{\mu}_\fij \gradv_\fij P \cdot \vect{S}_\ij \right)
746 * \f]
747 *
748 * \param[in] f_id field id (or -1)
749 * \param[in] m pointer to mesh
750 * \param[in] fvq pointer to finite volume quantities
751 * \param[in] init indicator
752 * - 1 initialize the mass flux to 0
753 * - 0 otherwise
754 * \param[in] inc indicator
755 * - 0 when solving an increment
756 * - 1 otherwise
757 * \param[in] imrgra indicator
758 * - 0 iterative gradient
759 * - 1 least squares gradient
760 * \param[in] nswrgp number of reconstruction sweeps for the
761 * gradients
762 * \param[in] imligp clipping gradient method
763 * - < 0 no clipping
764 * - = 0 thank to neighbooring gradients
765 * - = 1 thank to the mean gradient
766 * \param[in] ircflp indicator
767 * - 1 flux reconstruction,
768 * - 0 otherwise
769 * \param[in] ircflb indicator
770 * - 1 boundary flux reconstruction,
771 * - 0 otherwise
772 * \param[in] iphydp indicator
773 * - 1 hydrostatic pressure taken into account
774 * - 0 otherwise
775 * \param[in] iwgrp indicator
776 * - 1 weight gradient by vicosity*porosity
777 * - weighting determined by field options
778 * \param[in] iwarnp verbosity
779 * \param[in] epsrgp relative precision for the gradient
780 * reconstruction
781 * \param[in] climgp clipping coeffecient for the computation of
782 * the gradient
783 * \param[in] frcxt body force creating the hydrostatic pressure
784 * \param[in] pvar solved variable (pressure)
785 * \param[in] bc_coeffs boundary condition structure for the variable
786 * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
787 * at interior faces for the r.h.s.
788 * \param[in] b_visc \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
789 * at border faces for the r.h.s.
790 * \param[in] viscel symmetric cell tensor \f$ \tens{\mu}_\celli \f$
791 * \param[in] weighf internal face weight between cells i j in case
792 * of tensor diffusion
793 * \param[in] weighb boundary face weight for cells i in case
794 * of tensor diffusion
795 * \param[in,out] i_massflux mass flux at interior faces
796 * \param[in,out] b_massflux mass flux at boundary faces
797 */
798/*----------------------------------------------------------------------------*/
799
800void
802 const cs_mesh_t *m,
804 int init,
805 int inc,
806 int imrgra,
807 int nswrgp,
808 int imligp,
809 int ircflp,
810 int ircflb,
811 int iphydp,
812 int iwgrp,
813 int iwarnp,
814 double epsrgp,
815 double climgp,
816 cs_real_3_t *frcxt,
817 cs_real_t *pvar,
818 const cs_field_bc_coeffs_t *bc_coeffs,
819 const cs_real_t i_visc[],
820 const cs_real_t b_visc[],
821 cs_real_6_t *viscel,
822 const cs_real_2_t weighf[],
823 const cs_real_t weighb[],
824 cs_real_t *i_massflux,
825 cs_real_t *b_massflux);
826
827/*----------------------------------------------------------------------------*/
874/*----------------------------------------------------------------------------*/
875
876void
878 const cs_mesh_t *m,
880 int init,
881 int inc,
882 int imrgra,
883 int nswrgp,
884 int imligp,
885 int iphydp,
886 int iwgrp,
887 int iwarnp,
888 double epsrgp,
889 double climgp,
890 cs_real_3_t *frcxt,
891 cs_real_t *pvar,
892 const cs_field_bc_coeffs_t *bc_coeffs,
893 const cs_real_t i_visc[],
894 const cs_real_t b_visc[],
895 cs_real_t visel[],
896 cs_real_t *diverg);
897
898/*----------------------------------------------------------------------------*/
899/*
900 * \brief Add the explicit part of the divergence of the mass flux due to the
901 * pressure gradient (analog to cs_anisotropic_diffusion_scalar).
902 *
903 * More precisely, the divergence of the mass flux side
904 * \f$ \sum_{\fij \in \Facei{\celli}} \dot{m}_\fij \f$ is updated as follows:
905 * \f[
906 * \sum_{\fij \in \Facei{\celli}} \dot{m}_\fij
907 * = \sum_{\fij \in \Facei{\celli}} \dot{m}_\fij
908 * - \sum_{\fij \in \Facei{\celli}}
909 * \left( \tens{\mu}_\fij \gradv_\fij P \cdot \vect{S}_\ij \right)
910 * \f]
911 *
912 * \param[in] f_id field id (or -1)
913 * \param[in] m pointer to mesh
914 * \param[in] fvq pointer to finite volume quantities
915 * \param[in] init indicator
916 * - 1 initialize the mass flux to 0
917 * - 0 otherwise
918 * \param[in] inc indicator
919 * - 0 when solving an increment
920 * - 1 otherwise
921 * \param[in] imrgra indicator
922 * - 0 iterative gradient
923 * - 1 least squares gradient
924 * \param[in] nswrgp number of reconstruction sweeps for the
925 * gradients
926 * \param[in] imligp clipping gradient method
927 * - < 0 no clipping
928 * - = 0 thank to neighbooring gradients
929 * - = 1 thank to the mean gradient
930 * \param[in] ircflp indicator
931 * - 1 flux reconstruction,
932 * - 0 otherwise
933 * \param[in] ircflb indicator
934 * - 1 boundary flux reconstruction,
935 * - 0 otherwise
936 * \param[in] iphydp indicator
937 * - 1 hydrostatic pressure taken into account
938 * - 0 otherwise
939 * \param[in] iwgrp indicator
940 * - 1 weight gradient by vicosity*porosity
941 * - weighting determined by field options
942 * \param[in] iwarnp verbosity
943 * \param[in] epsrgp relative precision for the gradient
944 * reconstruction
945 * \param[in] climgp clipping coeffecient for the computation of
946 * the gradient
947 * \param[in] frcxt body force creating the hydrostatic pressure
948 * \param[in] pvar solved variable (pressure)
949 * \param[in] bc_coeffs boundary condition structure for the variable
950 * \param[in] i_visc \f$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} \f$
951 * at interior faces for the r.h.s.
952 * \param[in] b_visc \f$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} \f$
953 * at border faces for the r.h.s.
954 * \param[in] viscel symmetric cell tensor \f$ \tens{\mu}_\celli \f$
955 * \param[in] weighf internal face weight between cells i j in case
956 * of tensor diffusion
957 * \param[in] weighb boundary face weight for cells i in case
958 * of tensor diffusion
959 * \param[in,out] diverg divergence of the mass flux
960 */
961/*----------------------------------------------------------------------------*/
962
963void
965 const cs_mesh_t *m,
967 int init,
968 int inc,
969 int imrgra,
970 int nswrgp,
971 int imligp,
972 int ircflp,
973 int ircflb,
974 int iphydp,
975 int iwgrp,
976 int iwarnp,
977 double epsrgp,
978 double climgp,
979 cs_real_3_t *frcxt,
980 cs_real_t *pvar,
981 const cs_field_bc_coeffs_t *bc_coeffs,
982 const cs_real_t i_visc[],
983 const cs_real_t b_visc[],
984 cs_real_6_t *viscel,
985 const cs_real_2_t weighf[],
986 const cs_real_t weighb[],
987 cs_real_t *diverg);
988
989/*----------------------------------------------------------------------------*/
990
992
993#ifdef __cplusplus
994
995/*----------------------------------------------------------------------------*/
1012/*----------------------------------------------------------------------------*/
1013
1014void
1015cs_slope_test_gradient(int f_id,
1016 cs_dispatch_context &ctx,
1017 int inc,
1018 cs_halo_type_t halo_type,
1019 const cs_real_3_t *grad,
1020 cs_real_3_t *grdpa,
1021 const cs_real_t *pvar,
1022 const cs_field_bc_coeffs_t *bc_coeffs,
1023 const cs_real_t *i_massflux);
1024
1025/*----------------------------------------------------------------------------*/
1026/*
1027 * \brief Compute the upwind gradient used in the pure SOLU schemes
1028 * (observed in the litterature).
1029 *
1030 * \param[in] f_id field index
1031 * \param[in] ctx Reference to dispatch context
1032 * \param[in] inc Not an increment flag
1033 * \param[in] halo_type halo type
1034 * \param[in] bc_coeffs boundary condition structure for the variable
1035 * \param[in] i_massflux mass flux at interior faces
1036 * \param[in] b_massflux mass flux at boundary faces
1037 * \param[in] pvar values
1038 * \param[out] grdpa upwind gradient
1039 */
1040/*----------------------------------------------------------------------------*/
1041
1042void
1043cs_upwind_gradient(const int f_id,
1044 cs_dispatch_context &ctx,
1045 const int inc,
1046 const cs_halo_type_t halo_type,
1047 const cs_field_bc_coeffs_t *bc_coeffs,
1048 const cs_real_t i_massflux[],
1049 const cs_real_t b_massflux[],
1050 const cs_real_t *pvar,
1051 cs_real_3_t *grdpa);
1052
1053/*----------------------------------------------------------------------------
1054 * Compute the local cell Courant number as the maximum of all cell face based
1055 * Courant number at each cell.
1056 *
1057 * parameters:
1058 * f <-- pointer to field
1059 * ctx <-- reference to dispatch context
1060 * courant --> cell Courant number
1061 */
1062/*----------------------------------------------------------------------------*/
1063
1064void
1066 cs_dispatch_context &ctx,
1067 cs_real_t *courant);
1068
1069#endif /* cplusplus */
1070
1071#endif /* __CS_CONVECTION_DIFFUSION_H__ */
void cs_upwind_gradient(const int f_id, cs_dispatch_context &ctx, const int inc, const cs_halo_type_t halo_type, const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t *restrict pvar, cs_real_3_t *restrict grdpa)
Compute the upwind gradient used in the pure SOLU schemes (observed in the litterature).
Definition cs_convection_diffusion.cpp:12558
void cs_slope_test_gradient(int f_id, cs_dispatch_context &ctx, int inc, cs_halo_type_t halo_type, const cs_real_3_t *grad, cs_real_3_t *grdpa, const cs_real_t *pvar, const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t *i_massflux)
Compute the upwind gradient used in the slope tests.
Definition cs_convection_diffusion.cpp:12470
void cs_cell_courant_number(const cs_field_t *f, cs_dispatch_context &ctx, cs_real_t *courant)
Definition cs_convection_diffusion.cpp:12383
void cs_convection_diffusion_vector(int idtvar, int f_id, const cs_equation_param_t eqp, int icvflb, int inc, int ivisep, int imasac, cs_real_3_t *pvar, const cs_real_3_t *pvara, const int icvfli[], const cs_field_bc_coeffs_t *bc_coeffs_v, const cs_bc_coeffs_solve_t *bc_coeffs_solve_v, const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t i_secvis[], const cs_real_t b_secvis[], cs_real_3_t *i_pvar, cs_real_3_t *b_pvar, cs_real_3_t *rhs)
void cs_anisotropic_diffusion_scalar(int idtvar, int f_id, const cs_equation_param_t eqp, int inc, cs_real_t *pvar, const cs_real_t *pvara, const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t *viscel, const cs_real_2_t weighf[], const cs_real_t weighb[], cs_real_t *rhs)
void cs_anisotropic_diffusion_tensor(int idtvar, int f_id, const cs_equation_param_t eqp, int inc, cs_real_6_t *pvar, const cs_real_6_t *pvara, const cs_field_bc_coeffs_t *bc_coeffs_ts, const cs_bc_coeffs_solve_t *bc_coeffs_solve_ts, const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t *viscel, const cs_real_2_t weighf[], const cs_real_t weighb[], cs_real_6_t *rhs)
void cs_anisotropic_right_diffusion_vector(int idtvar, int f_id, const cs_equation_param_t eqp, int inc, cs_real_3_t *pvar, const cs_real_3_t *pvara, const cs_field_bc_coeffs_t *bc_coeffs_v, const cs_bc_coeffs_solve_t *bc_coeffs_solve_v, const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t *viscel, const cs_real_2_t weighf[], const cs_real_t weighb[], cs_real_3_t *rhs)
cs_real_t * cs_get_v_slope_test(int f_id, const cs_equation_param_t eqp)
Definition cs_convection_diffusion.cpp:7216
void cs_convection_diffusion_thermal(int idtvar, int f_id, const cs_equation_param_t eqp, int inc, int imasac, cs_real_t *pvar, const cs_real_t *pvara, const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t xcpp[], cs_real_t *rhs)
Add the explicit part of the convection/diffusion terms of a transport equation of a scalar field su...
void cs_anisotropic_diffusion_potential(const int f_id, const cs_mesh_t *m, cs_mesh_quantities_t *fvq, int init, int inc, int imrgra, int nswrgp, int imligp, int ircflp, int ircflb, int iphydp, int iwgrp, int iwarnp, double epsrgp, double climgp, cs_real_3_t *frcxt, cs_real_t *pvar, const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t *viscel, const cs_real_2_t weighf[], const cs_real_t weighb[], cs_real_t *diverg)
cs_nvd_type_t
Definition cs_convection_diffusion.h:61
@ CS_NVD_SUPERBEE
Definition cs_convection_diffusion.h:66
@ CS_NVD_SMART
Definition cs_convection_diffusion.h:64
@ CS_NVD_CUBISTA
Definition cs_convection_diffusion.h:65
@ CS_NVD_STOIC
Definition cs_convection_diffusion.h:70
@ CS_NVD_WASEB
Definition cs_convection_diffusion.h:72
@ CS_NVD_CLAM
Definition cs_convection_diffusion.h:69
@ CS_NVD_OSHER
Definition cs_convection_diffusion.h:71
@ CS_NVD_VOF_CICSAM
Definition cs_convection_diffusion.h:74
@ CS_NVD_VOF_STACS
Definition cs_convection_diffusion.h:75
@ CS_NVD_GAMMA
Definition cs_convection_diffusion.h:63
@ CS_NVD_MINMOD
Definition cs_convection_diffusion.h:68
@ CS_NVD_VOF_HRIC
Definition cs_convection_diffusion.h:73
@ CS_NVD_MUSCL
Definition cs_convection_diffusion.h:67
@ CS_NVD_N_TYPES
Definition cs_convection_diffusion.h:76
void cs_face_anisotropic_diffusion_potential(const int f_id, const cs_mesh_t *m, cs_mesh_quantities_t *fvq, int init, int inc, int imrgra, int nswrgp, int imligp, int ircflp, int ircflb, int iphydp, int iwgrp, int iwarnp, double epsrgp, double climgp, cs_real_3_t *frcxt, cs_real_t *pvar, const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t *viscel, const cs_real_2_t weighf[], const cs_real_t weighb[], cs_real_t *i_massflux, cs_real_t *b_massflux)
void cs_face_diffusion_potential(const int f_id, const cs_mesh_t *m, cs_mesh_quantities_t *fvq, int init, int inc, int imrgra, int nswrgp, int imligp, int iphydp, int iwgrp, int iwarnp, double epsrgp, double climgp, cs_real_3_t *frcxt, cs_real_t *pvar, const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_t *visel, cs_real_t *i_massflux, cs_real_t *b_massflux)
void cs_diffusion_potential(const int f_id, const cs_mesh_t *m, cs_mesh_quantities_t *fvq, int init, int inc, int imrgra, int nswrgp, int imligp, int iphydp, int iwgrp, int iwarnp, double epsrgp, double climgp, cs_real_3_t *frcxt, cs_real_t *pvar, const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_t visel[], cs_real_t *diverg)
Update the cell mass flux divergence with the face pressure (or pressure increment,...
void cs_convection_diffusion_tensor(int idtvar, int f_id, const cs_equation_param_t eqp, int icvflb, int inc, int imasac, cs_real_6_t *pvar, const cs_real_6_t *pvara, const cs_field_bc_coeffs_t *bc_coeffs_ts, const cs_bc_coeffs_solve_t *bc_coeffs_solve_ts, const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t *rhs)
Add the explicit part of the convection/diffusion terms of a transport equation of a tensor field .
void cs_anisotropic_left_diffusion_vector(int idtvar, int f_id, const cs_equation_param_t eqp, int inc, int ivisep, cs_real_3_t *pvar, const cs_real_3_t *pvara, const cs_field_bc_coeffs_t *bc_coeffs_v, const cs_bc_coeffs_solve_t *bc_coeffs_solve_v, const cs_real_33_t i_visc[], const cs_real_t b_visc[], const cs_real_t i_secvis[], cs_real_3_t *rhs)
void cs_convection_diffusion_scalar(int idtvar, int f_id, const cs_equation_param_t eqp, int icvflb, int inc, int imasac, cs_real_t *pvar, const cs_real_t *pvara, const int icvfli[], const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_t *rhs, cs_real_2_t i_flux[], cs_real_t b_flux[])
Add the explicit part of the convection/diffusion terms of a standard transport equation of a scalar ...
void cs_face_convection_scalar(int idtvar, int f_id, const cs_equation_param_t eqp, int icvflb, int inc, int imasac, cs_real_t *pvar, const cs_real_t *pvara, const int icvfli[], const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t i_massflux[], const cs_real_t b_massflux[], cs_real_2_t i_conv_flux[], cs_real_t b_conv_flux[])
void cs_beta_limiter_building(int f_id, int inc, const cs_real_t rovsdt[])
Compute the beta blending coefficient of the beta limiter (ensuring preservation of a given min/max p...
Definition cs_convection_diffusion.cpp:7265
cs_real_t cs_real_33_t[3][3]
3x3 matrix of floating-point values
Definition cs_defs.h:368
#define BEGIN_C_DECLS
Definition cs_defs.h:542
double cs_real_t
Floating-point value.
Definition cs_defs.h:342
cs_real_t cs_real_2_t[2]
vector of 2 floating-point values
Definition cs_defs.h:358
cs_real_t cs_real_6_t[6]
vector of 6 floating-point values
Definition cs_defs.h:361
#define END_C_DECLS
Definition cs_defs.h:543
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition cs_defs.h:359
cs_halo_type_t
Definition cs_halo.h:56
Definition cs_field.h:104
Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources.
Definition cs_equation_param.h:193
Field boundary condition descriptor (for variables)
Definition cs_field.h:121
Field descriptor.
Definition cs_field.h:158
Definition cs_mesh_quantities.h:92
Definition cs_mesh.h:85