programmer's documentation
 All Data Structures Files Functions Variables Typedefs Enumerations Enumerator Macros Groups Pages
Functions
cs_convection_diffusion.h File Reference
#include "cs_base.h"
#include "cs_halo.h"
Include dependency graph for cs_convection_diffusion.h:

Go to the source code of this file.

Functions

void bilsc2 (const cs_int_t *const idtvar, const cs_int_t *const f_id, const cs_var_cal_opt_t *const var_cal_opt, const cs_int_t *const icvflb, const cs_int_t *const inc, const cs_int_t *const iccocg, const cs_int_t *const ifaccp, cs_real_t pvar[], const cs_real_t pvara[], const cs_int_t bc_type[], const cs_int_t icvfli[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], 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[])
 
void bilsc4 (const cs_int_t *const idtvar, const cs_int_t *const f_id, const cs_var_cal_opt_t *const var_cal_opt, const cs_int_t *const icvflb, const cs_int_t *const inc, const cs_int_t *const ifaccp, const cs_int_t *const ivisep, cs_real_3_t pvar[], const cs_real_3_t pvara[], const cs_int_t bc_type[], const cs_int_t icvfli[], const cs_real_3_t coefav[], const cs_real_33_t coefbv[], const cs_real_3_t cofafv[], const cs_real_33_t cofbfv[], 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 secvif[], cs_real_3_t rhs[])
 
void bilsct (const cs_int_t *const idtvar, const cs_int_t *const f_id, const cs_var_cal_opt_t *const var_cal_opt, const cs_int_t *const inc, const cs_int_t *const iccocg, const cs_int_t *const ifaccp, cs_real_t pvar[], const cs_real_t pvara[], const cs_int_t bc_type[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], 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[])
 
void diften (const cs_int_t *const idtvar, const cs_int_t *const f_id, const cs_var_cal_opt_t *const var_cal_opt, const cs_int_t *const inc, const cs_int_t *const iccocg, cs_real_t pvar[], const cs_real_t pvara[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], 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 diftnv (const cs_int_t *const idtvar, const cs_int_t *const f_id, const cs_var_cal_opt_t *const var_cal_opt, const cs_int_t *const inc, const cs_int_t *const ifaccp, const cs_int_t *const ivisep, cs_real_3_t pvar[], const cs_real_3_t pvara[], const cs_int_t bc_type[], const cs_real_3_t coefav[], const cs_real_33_t coefbv[], const cs_real_3_t cofafv[], const cs_real_33_t cofbfv[], const cs_real_33_t i_visc[], const cs_real_t b_visc[], const cs_real_t secvif[], cs_real_3_t rhs[])
 
void itrmas (const cs_int_t *const init, const cs_int_t *const inc, const cs_int_t *const imrgra, const cs_int_t *const iccocg, const cs_int_t *const nswrgp, const cs_int_t *const imligp, const cs_int_t *const iphydp, const cs_int_t *const iwarnp, const cs_real_t *const epsrgp, const cs_real_t *const climgp, const cs_real_t *const extrap, cs_real_3_t frcxt[], cs_real_t pvar[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t viselx[], const cs_real_t visely[], const cs_real_t viselz[], cs_real_t i_massflux[], cs_real_t b_massflux[])
 
void itrmav (const cs_int_t *const init, const cs_int_t *const inc, const cs_int_t *const imrgra, const cs_int_t *const iccocg, const cs_int_t *const nswrgp, const cs_int_t *const imligp, const cs_int_t *const ircflp, const cs_int_t *const iphydp, const cs_int_t *const iwarnp, const cs_real_t *const epsrgp, const cs_real_t *const climgp, const cs_real_t *const extrap, cs_real_3_t frcxt[], cs_real_t pvar[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], 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 itrgrp (const cs_int_t *const init, const cs_int_t *const inc, const cs_int_t *const imrgra, const cs_int_t *const iccocg, const cs_int_t *const nswrgp, const cs_int_t *const imligp, const cs_int_t *const iphydp, const cs_int_t *const iwarnp, const cs_real_t *const epsrgp, const cs_real_t *const climgp, const cs_real_t *const extrap, cs_real_3_t frcxt[], cs_real_t pvar[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t viselx[], const cs_real_t visely[], const cs_real_t viselz[], cs_real_t diverg[])
 
void itrgrv (const cs_int_t *const init, const cs_int_t *const inc, const cs_int_t *const imrgra, const cs_int_t *const iccocg, const cs_int_t *const nswrgp, const cs_int_t *const imligp, const cs_int_t *const ircflp, const cs_int_t *const iphydp, const cs_int_t *const iwarnp, const cs_real_t *const epsrgp, const cs_real_t *const climgp, const cs_real_t *const extrap, cs_real_3_t frcxt[], cs_real_t pvar[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], 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[])
 
void cs_convection_diffusion_scalar (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int icvflb, int inc, int iccocg, int ifaccp, cs_real_t *restrict pvar, const cs_real_t *restrict pvara, const cs_int_t bc_type[], const cs_int_t icvfli[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], 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 *restrict rhs)
 Add the explicit part of the convection/diffusion terms of a standard transport equation of a scalar field $ \varia $. More...
 
void cs_convection_diffusion_vector (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int icvflb, int inc, int ifaccp, int ivisep, cs_real_3_t *restrict pvar, const cs_real_3_t *restrict pvara, const cs_int_t bc_type[], const cs_int_t icvfli[], const cs_real_3_t coefav[], const cs_real_33_t coefbv[], const cs_real_3_t cofafv[], const cs_real_33_t cofbfv[], 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 secvif[], cs_real_3_t *restrict rhs)
 Add the explicit part of the convection/diffusion terms of a transport equation of a vector field $ \vect{\varia} $. More...
 
void cs_convection_diffusion_thermal (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int inc, int iccocg, int ifaccp, cs_real_t *restrict pvar, const cs_real_t *restrict pvara, const cs_int_t bc_type[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], 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 *restrict rhs)
 Add the explicit part of the convection/diffusion terms of a transport equation of a scalar field $ \varia $ such as the temperature. More...
 
void cs_anisotropic_diffusion_scalar (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int inc, int iccocg, cs_real_t *restrict pvar, const cs_real_t *restrict pvara, const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t *restrict viscel, const cs_real_2_t weighf[], const cs_real_t weighb[], cs_real_t *restrict rhs)
 Add the explicit part of the diffusion terms with a symmetric tensor diffusivity for a transport equation of a scalar field $ \varia $. More...
 
void cs_anisotropic_diffusion_vector (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int inc, int ifaccp, int ivisep, cs_real_3_t *restrict pvar, const cs_real_3_t *restrict pvara, const cs_int_t bc_type[], const cs_real_3_t coefav[], const cs_real_33_t coefbv[], const cs_real_3_t cofafv[], const cs_real_33_t cofbfv[], const cs_real_33_t i_visc[], const cs_real_t b_visc[], const cs_real_t secvif[], cs_real_3_t *restrict rhs)
 Add the explicit part of the diffusion terms with a symmetric tensorial diffusivity for a transport equation of a vector field $ \vect{\varia} $. More...
 
void cs_face_diffusion_potential (const cs_mesh_t *m, cs_mesh_quantities_t *fvq, int init, int inc, int imrgra, int iccocg, int nswrgp, int imligp, int iphydp, int iwarnp, double epsrgp, double climgp, double extrap, cs_real_3_t *restrict frcxt, cs_real_t *restrict pvar, const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t viselx[], const cs_real_t visely[], const cs_real_t viselz[], cs_real_t *restrict i_massflux, cs_real_t *restrict b_massflux)
 Update the face mass flux with the face pressure (or pressure increment, or pressure double increment) gradient. More...
 
void cs_face_anisotropic_diffusion_potential (const cs_mesh_t *m, cs_mesh_quantities_t *fvq, int init, int inc, int imrgra, int iccocg, int nswrgp, int imligp, int ircflp, int iphydp, int iwarnp, double epsrgp, double climgp, double extrap, cs_real_3_t *restrict frcxt, cs_real_t *restrict pvar, const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t *restrict viscel, const cs_real_2_t weighf[], const cs_real_t weighb[], cs_real_t *restrict i_massflux, cs_real_t *restrict b_massflux)
 Add the explicit part of the pressure gradient term to the mass flux in case of anisotropic diffusion of the pressure field $ P $. More...
 
void cs_diffusion_potential (const cs_mesh_t *m, cs_mesh_quantities_t *fvq, int init, int inc, int imrgra, int iccocg, int nswrgp, int imligp, int iphydp, int iwarnp, double epsrgp, double climgp, double extrap, cs_real_3_t *restrict frcxt, cs_real_t *restrict pvar, const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t viselx[], const cs_real_t visely[], const cs_real_t viselz[], cs_real_t *restrict diverg)
 Update the cell mass flux divergence with the face pressure (or pressure increment, or pressure double increment) gradient. More...
 
void cs_anisotropic_diffusion_potential (const cs_mesh_t *m, cs_mesh_quantities_t *fvq, int init, int inc, int imrgra, int iccocg, int nswrgp, int imligp, int ircflp, int iphydp, int iwarnp, double epsrgp, double climgp, double extrap, cs_real_3_t *restrict frcxt, cs_real_t *restrict pvar, const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t *restrict viscel, const cs_real_2_t weighf[], const cs_real_t weighb[], cs_real_t *restrict diverg)
 Add the explicit part of the divergence of the mass flux due to the pressure gradient (routine analog to cs_anisotropic_diffusion_scalar). More...
 

Function Documentation

void bilsc2 ( const cs_int_t *const  idtvar,
const cs_int_t *const  f_id,
const cs_var_cal_opt_t *const  var_cal_opt,
const cs_int_t *const  icvflb,
const cs_int_t *const  inc,
const cs_int_t *const  iccocg,
const cs_int_t *const  ifaccp,
cs_real_t  pvar[],
const cs_real_t  pvara[],
const cs_int_t  bc_type[],
const cs_int_t  icvfli[],
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
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[] 
)

(end ignore by Doxygen)

void bilsc4 ( const cs_int_t *const  idtvar,
const cs_int_t *const  f_id,
const cs_var_cal_opt_t *const  var_cal_opt,
const cs_int_t *const  icvflb,
const cs_int_t *const  inc,
const cs_int_t *const  ifaccp,
const cs_int_t *const  ivisep,
cs_real_3_t  pvar[],
const cs_real_3_t  pvara[],
const cs_int_t  bc_type[],
const cs_int_t  icvfli[],
const cs_real_3_t  coefav[],
const cs_real_33_t  coefbv[],
const cs_real_3_t  cofafv[],
const cs_real_33_t  cofbfv[],
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  secvif[],
cs_real_3_t  rhs[] 
)
void bilsct ( const cs_int_t *const  idtvar,
const cs_int_t *const  f_id,
const cs_var_cal_opt_t *const  var_cal_opt,
const cs_int_t *const  inc,
const cs_int_t *const  iccocg,
const cs_int_t *const  ifaccp,
cs_real_t  pvar[],
const cs_real_t  pvara[],
const cs_int_t  bc_type[],
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
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[] 
)
void cs_anisotropic_diffusion_potential ( const cs_mesh_t m,
cs_mesh_quantities_t fvq,
int  init,
int  inc,
int  imrgra,
int  iccocg,
int  nswrgp,
int  imligp,
int  ircflp,
int  iphydp,
int  iwarnp,
double  epsrgp,
double  climgp,
double  extrap,
cs_real_3_t *restrict  frcxt,
cs_real_t *restrict  pvar,
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
const cs_real_t  i_visc[],
const cs_real_t  b_visc[],
cs_real_6_t *restrict  viscel,
const cs_real_2_t  weighf[],
const cs_real_t  weighb[],
cs_real_t *restrict  diverg 
)

Add the explicit part of the divergence of the mass flux due to the pressure gradient (routine analog to cs_anisotropic_diffusion_scalar).

More precisely, the divergence of the mass flux side $ \sum_{\fij \in \Facei{\celli}} \dot{m}_\fij $ is updated as follows:

\[ \sum_{\fij \in \Facei{\celli}} \dot{m}_\fij = \sum_{\fij \in \Facei{\celli}} \dot{m}_\fij - \sum_{\fij \in \Facei{\celli}} \left( \tens{\mu}_\fij \gradv_\fij P \cdot \vect{S}_\ij \right) \]

Parameters
[in]mpointer to mesh
[in]fvqpointer to finite volume quantities
[in]initindicator
  • 1 initialize the mass flux to 0
  • 0 otherwise
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]imrgraindicator
  • 0 iterative gradient
  • 1 least square gradient
[in]iccocgindicator
  • 1 re-compute cocg matrix (for iterativ gradients)
  • 0 otherwise
[in]nswrgpnumber of reconstruction sweeps for the gradients
[in]imligpclipping gradient method
  • < 0 no clipping
  • = 0 thank to neighbooring gradients
  • = 1 thank to the mean gradient
[in]ircflpindicator
  • 1 flux reconstruction,
  • 0 otherwise
[in]iphydpindicator
  • 1 hydrostatic pressure taken into account
  • 0 otherwise
[in]iwarnpverbosity
[in]epsrgprelative precision for the gradient reconstruction
[in]climgpclipping coeffecient for the computation of the gradient
[in]extrapcoefficient for extrapolation of the gradient
[in]frcxtbody force creating the hydrostatic pressure
[in]pvarsolved variable (pressure)
[in]coefapboundary condition array for the variable (Explicit part)
[in]coefbpboundary condition array for the variable (Implicit part)
[in]cofafpboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbfpboundary condition array for the diffusion of the variable (Implicit part)
[in]i_visc$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the r.h.s.
[in]b_visc$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at border faces for the r.h.s.
[in]viscelsymmetric cell tensor $ \tens{\mu}_\celli $
[in]weighfinternal face weight between cells i j in case of tensor diffusion
[in]weighbboundary face weight for cells i in case of tensor diffusion
[in,out]divergdivergence of the mass flux

Add the explicit part of the divergence of the mass flux due to the pressure gradient (routine analog to cs_anisotropic_diffusion_scalar).

More precisely, the divergence of the mass flux side $ \sum_{\fij \in \Facei{\celli}} \dot{m}_\fij $ is updated as follows:

\[ \sum_{\fij \in \Facei{\celli}} \dot{m}_\fij = \sum_{\fij \in \Facei{\celli}} \dot{m}_\fij - \sum_{\fij \in \Facei{\celli}} \left( \tens{\mu}_\fij \gradv_\fij P \cdot \vect{S}_\ij \right) \]

Parameters
[in]mpointer to mesh
[in]fvqpointer to finite volume quantities
[in]initindicator
  • 1 initialize the mass flux to 0
  • 0 otherwise
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]imrgraindicator
  • 0 iterative gradient
  • 1 least square gradient
[in]iccocgindicator
  • 1 re-compute cocg matrix (for iterativ gradients)
  • 0 otherwise
[in]nswrgpnumber of reconstruction sweeps for the gradients
[in]imligpclipping gradient method
  • < 0 no clipping
  • = 0 thank to neighbooring gradients
  • = 1 thank to the mean gradient
[in]ircflpindicator
  • 1 flux reconstruction,
  • 0 otherwise
[in]iphydpindicator
  • 1 hydrostatic pressure taken into account
  • 0 otherwise
[in]iwarnpverbosity
[in]epsrgprelative precision for the gradient reconstruction
[in]climgpclipping coeffecient for the computation of the gradient
[in]extrapcoefficient for extrapolation of the gradient
[in]frcxtbody force creating the hydrostatic pressure
[in]pvarsolved variable (pressure)
[in]coefapboundary condition array for the variable (Explicit part)
[in]coefbpboundary condition array for the variable (Implicit part)
[in]cofafpboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbfpboundary condition array for the diffusion of the variable (Implicit part)
[in]i_visc$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the r.h.s.
[in]b_visc$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at border faces for the r.h.s.
[in]viscelsymmetric cell tensor $ \tens{\mu}_\celli $
[in]weighfinternal face weight between cells i j in case of tensor diffusion
[in]weighbboundary face weight for cells i in case of tensor diffusion
[in,out]divergdivergence of the mass flux
void cs_anisotropic_diffusion_scalar ( int  idtvar,
int  f_id,
const cs_var_cal_opt_t  var_cal_opt,
int  inc,
int  iccocg,
cs_real_t *restrict  pvar,
const cs_real_t *restrict  pvara,
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
const cs_real_t  i_visc[],
const cs_real_t  b_visc[],
cs_real_6_t *restrict  viscel,
const cs_real_2_t  weighf[],
const cs_real_t  weighb[],
cs_real_t *restrict  rhs 
)

Add the explicit part of the diffusion terms with a symmetric tensor diffusivity for a transport equation of a scalar field $ \varia $.

More precisely, the right hand side $ Rhs $ is updated as follows:

\[ Rhs = Rhs - \sum_{\fij \in \Facei{\celli}} \left( - \tens{\mu}_\fij \gradv_\fij \varia \cdot \vect{S}_\ij \right) \]

Warning:

  • $ Rhs $ has already been initialized before calling cs_anisotropic_diffusion_scalar!
  • mind the sign minus
Parameters
[in]idtvarindicator of the temporal scheme
[in]f_idindex of the current variable
[in]var_cal_optvariable calculation options
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]iccocgindicator
  • 1 re-compute cocg matrix (for iterativ gradients)
  • 0 otherwise
[in]pvarsolved variable (current time step)
[in]pvarasolved variable (previous time step)
[in]coefapboundary condition array for the variable (Explicit part)
[in]coefbpboundary condition array for the variable (Implicit part)
[in]cofafpboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbfpboundary condition array for the diffusion of the variable (Implicit part)
[in]i_visc$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the r.h.s.
[in]b_visc$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at border faces for the r.h.s.
[in]viscelsymmetric cell tensor $ \tens{\mu}_\celli $
[in]weighfinternal face weight between cells i j in case of tensor diffusion
[in]weighbboundary face weight for cells i in case of tensor diffusion
[in,out]rhsright hand side $ \vect{Rhs} $
void cs_anisotropic_diffusion_vector ( int  idtvar,
int  f_id,
const cs_var_cal_opt_t  var_cal_opt,
int  inc,
int  ifaccp,
int  ivisep,
cs_real_3_t *restrict  pvar,
const cs_real_3_t *restrict  pvara,
const cs_int_t  bc_type[],
const cs_real_3_t  coefav[],
const cs_real_33_t  coefbv[],
const cs_real_3_t  cofafv[],
const cs_real_33_t  cofbfv[],
const cs_real_33_t  i_visc[],
const cs_real_t  b_visc[],
const cs_real_t  secvif[],
cs_real_3_t *restrict  rhs 
)

Add the explicit part of the diffusion terms with a symmetric tensorial diffusivity for a transport equation of a vector field $ \vect{\varia} $.

More precisely, the right hand side $ \vect{Rhs} $ is updated as follows:

\[ \vect{Rhs} = \vect{Rhs} - \sum_{\fij \in \Facei{\celli}} \left( - \tens{\mu}_\fij \gradt_\fij \vect{\varia} \cdot \vect{S}_\ij \right) \]

Warning:

  • $ \vect{Rhs} $ has already been initialized before calling diftnv!
  • mind the sign minus
Parameters
[in]idtvarindicator of the temporal scheme
[in]f_idindex of the current variable
[in]var_cal_optvariable calculation options
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]ifaccpindicator
  • 1 coupling activated
  • 0 coupling not activated
[in]ivisepindicator to take $ \divv \left(\mu \gradt \transpose{\vect{a}} \right) -2/3 \grad\left( \mu \dive \vect{a} \right)$
  • 1 take into account,
[in]pvarsolved variable (current time step)
[in]pvarasolved variable (previous time step)
[in]bc_typeboundary condition type
[in]coefavboundary condition array for the variable (Explicit part)
[in]coefbvboundary condition array for the variable (Implicit part)
[in]cofafvboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbfvboundary condition array for the diffusion of the variable (Implicit part)
[in]i_visc$ \tens{\mu}_\fij \dfrac{S_\fij}{\ipf\jpf} $ at interior faces for the r.h.s.
[in]b_visc$ \dfrac{S_\fib}{\ipf \centf} $ at border faces for the r.h.s.
[in]secvifsecondary viscosity at interior faces
[in,out]rhsright hand side $ \vect{Rhs} $
void cs_convection_diffusion_scalar ( int  idtvar,
int  f_id,
const cs_var_cal_opt_t  var_cal_opt,
int  icvflb,
int  inc,
int  iccocg,
int  ifaccp,
cs_real_t *restrict  pvar,
const cs_real_t *restrict  pvara,
const cs_int_t  bc_type[],
const cs_int_t  icvfli[],
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
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 *restrict  rhs 
)

Add the explicit part of the convection/diffusion terms of a standard transport equation of a scalar field $ \varia $.

More precisely, the right hand side $ Rhs $ is updated as follows:

\[ Rhs = Rhs - \sum_{\fij \in \Facei{\celli}} \left( \dot{m}_\ij \left( \varia_\fij - \varia_\celli \right) - \mu_\fij \gradv_\fij \varia \cdot \vect{S}_\ij \right) \]

Warning:

  • $ Rhs $ has already been initialized before calling bilsc2!
  • mind the sign minus
Parameters
[in]idtvarindicator of the temporal scheme
[in]f_idfield id (or -1)
[in]var_cal_optvariable calculation options
[in]icvflbglobal indicator of boundary convection flux
  • 0 upwind scheme at all boundary faces
  • 1 imposed flux at some boundary faces
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]iccocgindicator
  • 1 re-compute cocg matrix (for iterative gradients)
  • 0 otherwise
[in]ifaccpindicator
  • 1 coupling activated
  • 0 coupling not activated
[in]pvarsolved variable (current time step)
[in]pvarasolved variable (previous time step)
[in]bc_typeboundary condition type
[in]icvfliboundary face indicator array of convection flux
  • 0 upwind scheme
  • 1 imposed flux
[in]coefapboundary condition array for the variable (Explicit part)
[in]coefbpboundary condition array for the variable (Implicit part)
[in]cofafpboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbfpboundary condition array for the diffusion of the variable (Implicit part)
[in]i_massfluxmass flux at interior faces
[in]b_massfluxmass flux at boundary faces
[in]i_visc$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the r.h.s.
[in]b_visc$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at border faces for the r.h.s.
[in,out]rhsright hand side $ \vect{Rhs} $
void cs_convection_diffusion_thermal ( int  idtvar,
int  f_id,
const cs_var_cal_opt_t  var_cal_opt,
int  inc,
int  iccocg,
int  ifaccp,
cs_real_t *restrict  pvar,
const cs_real_t *restrict  pvara,
const cs_int_t  bc_type[],
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
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 *restrict  rhs 
)

Add the explicit part of the convection/diffusion terms of a transport equation of a scalar field $ \varia $ such as the temperature.

More precisely, the right hand side $ Rhs $ is updated as follows:

\[ Rhs = Rhs + \sum_{\fij \in \Facei{\celli}} \left( C_p\dot{m}_\ij \varia_\fij - \lambda_\fij \gradv_\fij \varia \cdot \vect{S}_\ij \right) \]

Warning: $ Rhs $ has already been initialized before calling bilsct!

Parameters
[in]idtvarindicator of the temporal scheme
[in]f_idindex of the current variable
[in]var_cal_optvariable calculation options
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]iccocgindicator
  • 1 re-compute cocg matrix (for iterative gradients)
  • 0 otherwise
[in]ifaccpindicator
  • 1 coupling activated
  • 0 coupling not activated
[in]pvarsolved variable (current time step)
[in]pvarasolved variable (previous time step)
[in]bc_typeboundary condition type
[in]coefapboundary condition array for the variable (Explicit part)
[in]coefbpboundary condition array for the variable (Implicit part)
[in]cofafpboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbfpboundary condition array for the diffusion of the variable (Implicit part)
[in]i_massfluxmass flux at interior faces
[in]b_massfluxmass flux at boundary faces
[in]i_visc$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the r.h.s.
[in]b_visc$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at border faces for the r.h.s.
[in]xcpparray of specific heat ( $ C_p $)
[in,out]rhsright hand side $ \vect{Rhs} $
void cs_convection_diffusion_vector ( int  idtvar,
int  f_id,
const cs_var_cal_opt_t  var_cal_opt,
int  icvflb,
int  inc,
int  ifaccp,
int  ivisep,
cs_real_3_t *restrict  pvar,
const cs_real_3_t *restrict  pvara,
const cs_int_t  bc_type[],
const cs_int_t  icvfli[],
const cs_real_3_t  coefav[],
const cs_real_33_t  coefbv[],
const cs_real_3_t  cofafv[],
const cs_real_33_t  cofbfv[],
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  secvif[],
cs_real_3_t *restrict  rhs 
)

Add the explicit part of the convection/diffusion terms of a transport equation of a vector field $ \vect{\varia} $.

More precisely, the right hand side $ \vect{Rhs} $ is updated as follows:

\[ \vect{Rhs} = \vect{Rhs} - \sum_{\fij \in \Facei{\celli}} \left( \dot{m}_\ij \left( \vect{\varia}_\fij - \vect{\varia}_\celli \right) - \mu_\fij \gradt_\fij \vect{\varia} \cdot \vect{S}_\ij \right) \]

Remark: if ivisep = 1, then we also take $ \mu \transpose{\gradt\vect{\varia}} + \lambda \trace{\gradt\vect{\varia}} $, where $ \lambda $ is the secondary viscosity, i.e. usually $ -\frac{2}{3} \mu $.

Warning:

  • $ \vect{Rhs} $ has already been initialized before calling bilsc!
  • mind the sign minus
Parameters
[in]idtvarindicator of the temporal scheme
[in]f_idindex of the current variable
[in]var_cal_optvariable calculation options
[in]icvflbglobal indicator of boundary convection flux
  • 0 upwind scheme at all boundary faces
  • 1 imposed flux at some boundary faces
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]ifaccpindicator
  • 1 coupling activated
  • 0 coupling not activated
[in]ivisepindicator to take $ \divv \left(\mu \gradt \transpose{\vect{a}} \right) -2/3 \grad\left( \mu \dive \vect{a} \right)$
  • 1 take into account,
  • 0 otherwise
[in]pvarsolved velocity (current time step)
[in]pvarasolved velocity (previous time step)
[in]bc_typeboundary condition type
[in]icvfliboundary face indicator array of convection flux
  • 0 upwind scheme
  • 1 imposed flux
[in]coefavboundary condition array for the variable (Explicit part)
[in]coefbvboundary condition array for the variable (Implicit part)
[in]cofafvboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbfvboundary condition array for the diffusion of the variable (Implicit part)
[in]i_massfluxmass flux at interior faces
[in]b_massfluxmass flux at boundary faces
[in]i_visc$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the r.h.s.
[in]b_visc$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at border faces for the r.h.s.
[in]secvifsecondary viscosity at interior faces
[in,out]rhsright hand side $ \vect{Rhs} $
void cs_diffusion_potential ( const cs_mesh_t m,
cs_mesh_quantities_t fvq,
int  init,
int  inc,
int  imrgra,
int  iccocg,
int  nswrgp,
int  imligp,
int  iphydp,
int  iwarnp,
double  epsrgp,
double  climgp,
double  extrap,
cs_real_3_t *restrict  frcxt,
cs_real_t *restrict  pvar,
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
const cs_real_t  i_visc[],
const cs_real_t  b_visc[],
const cs_real_t  viselx[],
const cs_real_t  visely[],
const cs_real_t  viselz[],
cs_real_t *restrict  diverg 
)

Update the cell mass flux divergence with the face pressure (or pressure increment, or pressure double increment) gradient.

\[ \dot{m}_\ij = \dot{m}_\ij - \sum_j \Delta t \grad_\fij p \cdot \vect{S}_\ij \]

Parameters
[in]mpointer to mesh
[in]fvqpointer to finite volume quantities
[in]initindicator
  • 1 initialize the mass flux to 0
  • 0 otherwise
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]imrgraindicator
  • 0 iterative gradient
  • 1 least square gradient
[in]iccocgindicator
  • 1 re-compute cocg matrix (for iterative gradients)
  • 0 otherwise
[in]nswrgpnumber of reconstruction sweeps for the gradients
[in]imligpclipping gradient method
  • < 0 no clipping
  • = 0 thank to neighbooring gradients
  • = 1 thank to the mean gradient
[in]iphydphydrostatic pressure indicator
[in]iwarnpverbosity
[in]epsrgprelative precision for the gradient reconstruction
[in]climgpclipping coeffecient for the computation of the gradient
[in]extrapcoefficient for extrapolation of the gradient
[in]frcxtbody force creating the hydrostatic pressure
[in]pvarsolved variable (current time step)
[in]coefapboundary condition array for the variable (Explicit part)
[in]coefbpboundary condition array for the variable (Implicit part)
[in]cofafpboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbfpboundary condition array for the diffusion of the variable (Implicit part)
[in]i_visc$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the r.h.s.
[in]b_visc$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at border faces for the r.h.s.
[in]viselxviscosity by cell, dir x
[in]viselyviscosity by cell, dir y
[in]viselzviscosity by cell, dir z
[in,out]divergmass flux divergence
void cs_face_anisotropic_diffusion_potential ( const cs_mesh_t m,
cs_mesh_quantities_t fvq,
int  init,
int  inc,
int  imrgra,
int  iccocg,
int  nswrgp,
int  imligp,
int  ircflp,
int  iphydp,
int  iwarnp,
double  epsrgp,
double  climgp,
double  extrap,
cs_real_3_t *restrict  frcxt,
cs_real_t *restrict  pvar,
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
const cs_real_t  i_visc[],
const cs_real_t  b_visc[],
cs_real_6_t *restrict  viscel,
const cs_real_2_t  weighf[],
const cs_real_t  weighb[],
cs_real_t *restrict  i_massflux,
cs_real_t *restrict  b_massflux 
)

Add the explicit part of the pressure gradient term to the mass flux in case of anisotropic diffusion of the pressure field $ P $.

More precisely, the mass flux side $ \dot{m}_\fij $ is updated as follows:

\[ \dot{m}_\fij = \dot{m}_\fij - \left( \tens{\mu}_\fij \gradv_\fij P \cdot \vect{S}_\ij \right) \]

Parameters
[in]mpointer to mesh
[in]fvqpointer to finite volume quantities
[in]initindicator
  • 1 initialize the mass flux to 0
  • 0 otherwise
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]imrgraindicator
  • 0 iterative gradient
  • 1 least square gradient
[in]iccocgindicator
  • 1 re-compute cocg matrix (for iterativ gradients)
  • 0 otherwise
[in]nswrgpnumber of reconstruction sweeps for the gradients
[in]imligpclipping gradient method
  • < 0 no clipping
  • = 0 thank to neighbooring gradients
  • = 1 thank to the mean gradient
[in]ircflpindicator
  • 1 flux reconstruction,
  • 0 otherwise
[in]iphydpindicator
  • 1 hydrostatic pressure taken into account
  • 0 otherwise
[in]iwarnpverbosity
[in]epsrgprelative precision for the gradient reconstruction
[in]climgpclipping coeffecient for the computation of the gradient
[in]extrapcoefficient for extrapolation of the gradient
[in]frcxtbody force creating the hydrostatic pressure
[in]pvarsolved variable (pressure)
[in]coefapboundary condition array for the variable (Explicit part)
[in]coefbpboundary condition array for the variable (Implicit part)
[in]cofafpboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbfpboundary condition array for the diffusion of the variable (Implicit part)
[in]i_visc$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the r.h.s.
[in]b_visc$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at border faces for the r.h.s.
[in]viscelsymmetric cell tensor $ \tens{\mu}_\celli $
[in]weighfinternal face weight between cells i j in case of tensor diffusion
[in]weighbboundary face weight for cells i in case of tensor diffusion
[in,out]i_massfluxmass flux at interior faces
[in,out]b_massfluxmass flux at boundary faces
void cs_face_diffusion_potential ( const cs_mesh_t m,
cs_mesh_quantities_t fvq,
int  init,
int  inc,
int  imrgra,
int  iccocg,
int  nswrgp,
int  imligp,
int  iphydp,
int  iwarnp,
double  epsrgp,
double  climgp,
double  extrap,
cs_real_3_t *restrict  frcxt,
cs_real_t *restrict  pvar,
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
const cs_real_t  i_visc[],
const cs_real_t  b_visc[],
const cs_real_t  viselx[],
const cs_real_t  visely[],
const cs_real_t  viselz[],
cs_real_t *restrict  i_massflux,
cs_real_t *restrict  b_massflux 
)

Update the face mass flux with the face pressure (or pressure increment, or pressure double increment) gradient.

\[ \dot{m}_\ij = \dot{m}_\ij - \Delta t \grad_\fij \delta p \cdot \vect{S}_\ij \]

Parameters
[in]mpointer to mesh
[in]fvqpointer to finite volume quantities
[in]initindicator
  • 1 initialize the mass flux to 0
  • 0 otherwise
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]imrgraindicator
  • 0 iterative gradient
  • 1 least square gradient
[in]iccocgindicator
  • 1 re-compute cocg matrix (for iterative gradients)
  • 0 otherwise
[in]nswrgpnumber of reconstruction sweeps for the gradients
[in]imligpclipping gradient method
  • < 0 no clipping
  • = 0 thank to neighbooring gradients
  • = 1 thank to the mean gradient
[in]iphydphydrostatic pressure indicator
[in]iwarnpverbosity
[in]epsrgprelative precision for the gradient reconstruction
[in]climgpclipping coeffecient for the computation of the gradient
[in]extrapcoefficient for extrapolation of the gradient
[in]frcxtbody force creating the hydrostatic pressure
[in]pvarsolved variable (current time step)
[in]coefapboundary condition array for the variable (Explicit part)
[in]coefbpboundary condition array for the variable (Implicit part)
[in]cofafpboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbfpboundary condition array for the diffusion of the variable (Implicit part)
[in]i_visc$ \mu_\fij \dfrac{S_\fij}{\ipf \jpf} $ at interior faces for the r.h.s.
[in]b_visc$ \mu_\fib \dfrac{S_\fib}{\ipf \centf} $ at border faces for the r.h.s.
[in]viselxviscosity by cell, dir x
[in]viselyviscosity by cell, dir y
[in]viselzviscosity by cell, dir z
[in,out]i_massfluxmass flux at interior faces
[in,out]b_massfluxmass flux at boundary faces
void diften ( const cs_int_t *const  idtvar,
const cs_int_t *const  f_id,
const cs_var_cal_opt_t *const  var_cal_opt,
const cs_int_t *const  inc,
const cs_int_t *const  iccocg,
cs_real_t  pvar[],
const cs_real_t  pvara[],
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
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 diftnv ( const cs_int_t *const  idtvar,
const cs_int_t *const  f_id,
const cs_var_cal_opt_t *const  var_cal_opt,
const cs_int_t *const  inc,
const cs_int_t *const  ifaccp,
const cs_int_t *const  ivisep,
cs_real_3_t  pvar[],
const cs_real_3_t  pvara[],
const cs_int_t  bc_type[],
const cs_real_3_t  coefav[],
const cs_real_33_t  coefbv[],
const cs_real_3_t  cofafv[],
const cs_real_33_t  cofbfv[],
const cs_real_33_t  i_visc[],
const cs_real_t  b_visc[],
const cs_real_t  secvif[],
cs_real_3_t  rhs[] 
)
void itrgrp ( const cs_int_t *const  init,
const cs_int_t *const  inc,
const cs_int_t *const  imrgra,
const cs_int_t *const  iccocg,
const cs_int_t *const  nswrgp,
const cs_int_t *const  imligp,
const cs_int_t *const  iphydp,
const cs_int_t *const  iwarnp,
const cs_real_t *const  epsrgp,
const cs_real_t *const  climgp,
const cs_real_t *const  extrap,
cs_real_3_t  frcxt[],
cs_real_t  pvar[],
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
const cs_real_t  i_visc[],
const cs_real_t  b_visc[],
const cs_real_t  viselx[],
const cs_real_t  visely[],
const cs_real_t  viselz[],
cs_real_t  diverg[] 
)
void itrgrv ( const cs_int_t *const  init,
const cs_int_t *const  inc,
const cs_int_t *const  imrgra,
const cs_int_t *const  iccocg,
const cs_int_t *const  nswrgp,
const cs_int_t *const  imligp,
const cs_int_t *const  ircflp,
const cs_int_t *const  iphydp,
const cs_int_t *const  iwarnp,
const cs_real_t *const  epsrgp,
const cs_real_t *const  climgp,
const cs_real_t *const  extrap,
cs_real_3_t  frcxt[],
cs_real_t  pvar[],
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
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[] 
)
void itrmas ( const cs_int_t *const  init,
const cs_int_t *const  inc,
const cs_int_t *const  imrgra,
const cs_int_t *const  iccocg,
const cs_int_t *const  nswrgp,
const cs_int_t *const  imligp,
const cs_int_t *const  iphydp,
const cs_int_t *const  iwarnp,
const cs_real_t *const  epsrgp,
const cs_real_t *const  climgp,
const cs_real_t *const  extrap,
cs_real_3_t  frcxt[],
cs_real_t  pvar[],
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
const cs_real_t  i_visc[],
const cs_real_t  b_visc[],
const cs_real_t  viselx[],
const cs_real_t  visely[],
const cs_real_t  viselz[],
cs_real_t  i_massflux[],
cs_real_t  b_massflux[] 
)
void itrmav ( const cs_int_t *const  init,
const cs_int_t *const  inc,
const cs_int_t *const  imrgra,
const cs_int_t *const  iccocg,
const cs_int_t *const  nswrgp,
const cs_int_t *const  imligp,
const cs_int_t *const  ircflp,
const cs_int_t *const  iphydp,
const cs_int_t *const  iwarnp,
const cs_real_t *const  epsrgp,
const cs_real_t *const  climgp,
const cs_real_t *const  extrap,
cs_real_3_t  frcxt[],
cs_real_t  pvar[],
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
const cs_real_t  cofafp[],
const cs_real_t  cofbfp[],
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[] 
)