My Project
programmer's documentation
Functions
cs_convection_diffusion.c File Reference

Convection-diffusion operators. More...

#include "cs_defs.h"
#include <assert.h>
#include <errno.h>
#include <stdio.h>
#include <stdarg.h>
#include <string.h>
#include <math.h>
#include <float.h>
#include "bft_error.h"
#include "bft_mem.h"
#include "bft_printf.h"
#include "cs_blas.h"
#include "cs_halo.h"
#include "cs_halo_perio.h"
#include "cs_log.h"
#include "cs_math.h"
#include "cs_mesh.h"
#include "cs_field.h"
#include "cs_field_operator.h"
#include "cs_field_pointer.h"
#include "cs_gradient.h"
#include "cs_gradient_perio.h"
#include "cs_ext_neighborhood.h"
#include "cs_mesh_quantities.h"
#include "cs_parall.h"
#include "cs_parameters.h"
#include "cs_prototypes.h"
#include "cs_timer.h"
#include "cs_stokes_model.h"
#include "cs_boundary_conditions.h"
#include "cs_internal_coupling.h"
#include "cs_bad_cells_regularisation.h"
#include "cs_convection_diffusion.h"
Include dependency graph for cs_convection_diffusion.c:

Functions

void CS_PROCF (itrmas, ITRMAS) const
 
void CS_PROCF (itrmav, ITRMAV) const
 
void CS_PROCF (itrgrp, ITRGRP) const
 
void CS_PROCF (itrgrv, ITRGRV) const
 
void cs_slope_test_gradient (int f_id, 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_real_t *coefap, const cs_real_t *coefbp, const cs_real_t *i_massflux)
 Compute the upwind gradient used in the slope tests. More...
 
void cs_upwind_gradient (const int f_id, const int inc, const cs_halo_type_t halo_type, const cs_real_t coefap[], const cs_real_t coefbp[], 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 in order to cope with SOLU schemes observed in the litterature. More...
 
void cs_slope_test_gradient_vector (const int inc, const cs_halo_type_t halo_type, const cs_real_33_t *grad, cs_real_33_t *grdpa, const cs_real_3_t *pvar, const cs_real_3_t *coefa, const cs_real_33_t *coefb, const cs_real_t *i_massflux)
 Compute the upwind gradient used in the slope tests. More...
 
void cs_slope_test_gradient_tensor (const int inc, const cs_halo_type_t halo_type, const cs_real_63_t *grad, cs_real_63_t *grdpa, const cs_real_6_t *pvar, const cs_real_6_t *coefa, const cs_real_66_t *coefb, const cs_real_t *i_massflux)
 Compute the upwind gradient used in the slope tests. More...
 
void cs_max_limiter_building (int f_id, int inc, const cs_real_t rovsdt[])
 Compute a coefficient for blending that ensures the positivity of the scalar. More...
 
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 imasac, cs_real_t *restrict pvar, const cs_real_t *restrict pvara, 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_face_convection_scalar (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int icvflb, int inc, int iccocg, int imasac, cs_real_t *restrict pvar, const cs_real_t *restrict pvara, const cs_int_t icvfli[], const cs_real_t coefap[], const cs_real_t coefbp[], 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[])
 Update face flux with convection contribution 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 ivisep, int imasac, cs_real_3_t *restrict pvar, const cs_real_3_t *restrict pvara, 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 i_secvis[], const cs_real_t b_secvis[], 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_tensor (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int icvflb, int inc, int imasac, cs_real_6_t *restrict pvar, const cs_real_6_t *restrict pvara, const cs_real_6_t coefa[], const cs_real_66_t coefb[], const cs_real_6_t cofaf[], const cs_real_66_t cofbf[], 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 *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 imasac, 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_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_left_diffusion_vector (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int inc, int ivisep, cs_real_3_t *restrict pvar, const cs_real_3_t *restrict pvara, 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 i_secvis[], cs_real_3_t *restrict rhs)
 Add explicit part of the terms of diffusion by a left-multiplying symmetric tensorial diffusivity for a transport equation of a vector field $ \vect{\varia} $. More...
 
void cs_anisotropic_right_diffusion_vector (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int inc, cs_real_3_t *restrict pvar, const cs_real_3_t *restrict pvara, 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_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_3_t *restrict rhs)
 Add explicit part of the terms of diffusion by a right-multiplying symmetric tensorial diffusivity for a transport equation of a vector field $ \vect{\varia} $. More...
 
void cs_anisotropic_diffusion_tensor (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int inc, cs_real_6_t *restrict pvar, const cs_real_6_t *restrict pvara, const cs_real_6_t coefa[], const cs_real_66_t coefb[], const cs_real_6_t cofaf[], const cs_real_66_t cofbf[], 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_6_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_face_diffusion_potential (const int f_id, 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 iwgrp, 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_t *restrict visel, 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 int f_id, 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 iwgrp, 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 int f_id, 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[], cs_real_t visel[], 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 int f_id, 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 diften). More...
 

Detailed Description

Convection-diffusion operators.

Please refer to the convection-diffusion section of the theory guide for more informations.

Function Documentation

◆ cs_anisotropic_diffusion_potential()

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  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 diften).

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]f_idfield id (or -1)
[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 squares 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]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

◆ cs_anisotropic_diffusion_scalar()

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} $

◆ cs_anisotropic_diffusion_tensor()

void cs_anisotropic_diffusion_tensor ( int  idtvar,
int  f_id,
const cs_var_cal_opt_t  var_cal_opt,
int  inc,
cs_real_6_t *restrict  pvar,
const cs_real_6_t *restrict  pvara,
const cs_real_6_t  coefa[],
const cs_real_66_t  coefb[],
const cs_real_6_t  cofaf[],
const cs_real_66_t  cofbf[],
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_6_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]pvarsolved variable (current time step)
[in]pvarasolved variable (previous time step)
[in]coefaboundary condition array for the variable (explicit part)
[in]coefbboundary condition array for the variable (implicit part)
[in]cofafboundary condition array for the diffusion of the variable (explicit part)
[in]cofbfboundary 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} $

◆ cs_anisotropic_left_diffusion_vector()

void cs_anisotropic_left_diffusion_vector ( int  idtvar,
int  f_id,
const cs_var_cal_opt_t  var_cal_opt,
int  inc,
int  ivisep,
cs_real_3_t *restrict  pvar,
const cs_real_3_t *restrict  pvara,
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  i_secvis[],
cs_real_3_t *restrict  rhs 
)

Add explicit part of the terms of diffusion by a left-multiplying 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( - \gradt_\fij \vect{\varia} \tens{\mu}_\fij \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 the present function
  • 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]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]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]i_secvissecondary viscosity at interior faces
[in,out]rhsright hand side $ \vect{Rhs} $

◆ cs_anisotropic_right_diffusion_vector()

void cs_anisotropic_right_diffusion_vector ( int  idtvar,
int  f_id,
const cs_var_cal_opt_t  var_cal_opt,
int  inc,
cs_real_3_t *restrict  pvar,
const cs_real_3_t *restrict  pvara,
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_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_3_t *restrict  rhs 
)

Add explicit part of the terms of diffusion by a right-multiplying 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( - \gradt_\fij \vect{\varia} \tens{\mu}_\fij \cdot \vect{S}_\ij \right) \]

Warning:

  • $ \vect{Rhs} $ has already been initialized before calling the present function
  • 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]pvarsolved variable (current time step)
[in]pvarasolved variable (previous time step)
[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]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} $

◆ cs_convection_diffusion_scalar()

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  imasac,
cs_real_t *restrict  pvar,
const cs_real_t *restrict  pvara,
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

Please refer to the bilsc2 section of the theory guide for more informations.

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]imasactake mass accumulation into account?
[in]pvarsolved variable (current time step)
[in]pvarasolved variable (previous time step)
[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} $

◆ cs_convection_diffusion_tensor()

void cs_convection_diffusion_tensor ( int  idtvar,
int  f_id,
const cs_var_cal_opt_t  var_cal_opt,
int  icvflb,
int  inc,
int  imasac,
cs_real_6_t *restrict  pvar,
const cs_real_6_t *restrict  pvara,
const cs_real_6_t  coefa[],
const cs_real_66_t  coefb[],
const cs_real_6_t  cofaf[],
const cs_real_66_t  cofbf[],
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 *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) \]

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]imasactake mass accumulation into account?
[in]pvarsolved velocity (current time step)
[in]pvarasolved velocity (previous time step)
[in]coefaboundary condition array for the variable (Explicit part)
[in]coefbboundary condition array for the variable (Implicit part)
[in]cofafboundary condition array for the diffusion of the variable (Explicit part)
[in]cofbfboundary 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} $

◆ cs_convection_diffusion_thermal()

void cs_convection_diffusion_thermal ( int  idtvar,
int  f_id,
const cs_var_cal_opt_t  var_cal_opt,
int  inc,
int  iccocg,
int  imasac,
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_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]imasactake mass accumulation into account?
[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_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} $

◆ cs_convection_diffusion_vector()

void cs_convection_diffusion_vector ( int  idtvar,
int  f_id,
const cs_var_cal_opt_t  var_cal_opt,
int  icvflb,
int  inc,
int  ivisep,
int  imasac,
cs_real_3_t *restrict  pvar,
const cs_real_3_t *restrict  pvara,
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  i_secvis[],
const cs_real_t  b_secvis[],
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]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]imasactake mass accumulation into account?
[in]pvarsolved velocity (current time step)
[in]pvarasolved velocity (previous time step)
[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]i_secvissecondary viscosity at interior faces
[in]b_secvissecondary viscosity at boundary faces
[in,out]rhsright hand side $ \vect{Rhs} $

◆ cs_diffusion_potential()

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  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[],
cs_real_t  visel[],
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]f_idfield id (or -1)
[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 squares 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]viselviscosity by cell
[in,out]divergmass flux divergence

◆ cs_face_anisotropic_diffusion_potential()

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  iccocg,
int  nswrgp,
int  imligp,
int  ircflp,
int  iphydp,
int  iwgrp,
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]f_idfield id (or -1)
[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 squares 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]iwgrpindicator
  • 1 weight gradient by vicosity*porosity
  • weighting determined by field options
[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

◆ cs_face_convection_scalar()

void cs_face_convection_scalar ( int  idtvar,
int  f_id,
const cs_var_cal_opt_t  var_cal_opt,
int  icvflb,
int  inc,
int  iccocg,
int  imasac,
cs_real_t *restrict  pvar,
const cs_real_t *restrict  pvara,
const cs_int_t  icvfli[],
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
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[] 
)

Update face flux with convection contribution of a standard transport equation of a scalar field $ \varia $.

\[ C_\ij = \dot{m}_\ij \left( \varia_\fij - \varia_\celli \right) \]

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]imasactake mass accumulation into account?
[in]pvarsolved variable (current time step)
[in]pvarasolved variable (previous time step)
[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]i_massfluxmass flux at interior faces
[in]b_massfluxmass flux at boundary faces
[in,out]i_conv_fluxscalar convection flux at interior faces
[in,out]b_conv_fluxscalar convection flux at boundary faces

◆ cs_face_diffusion_potential()

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  iccocg,
int  nswrgp,
int  imligp,
int  iphydp,
int  iwgrp,
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_t *restrict  visel,
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 \]

Please refer to the itrmas/itrgrp section of the theory guide for more informations.

Parameters
[in]f_idfield id (or -1)
[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 squares 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]iwgrpindicator
  • 1 weight gradient by vicosity*porosity
  • weighting determined by field options
[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]viselviscosity by cell
[in,out]i_massfluxmass flux at interior faces
[in,out]b_massfluxmass flux at boundary faces

◆ cs_max_limiter_building()

void cs_max_limiter_building ( int  f_id,
int  inc,
const cs_real_t  rovsdt[] 
)

Compute a coefficient for blending that ensures the positivity of the scalar.

Parameters
[in]f_idfield id
[in]inc"not an increment" flag
[in]rovsdtrho * volume / dt

◆ CS_PROCF() [1/4]

void CS_PROCF ( itrgrp  ,
ITRGRP   
) const

◆ CS_PROCF() [2/4]

void CS_PROCF ( itrgrv  ,
ITRGRV   
) const

◆ CS_PROCF() [3/4]

void CS_PROCF ( itrmas  ,
ITRMAS   
) const

◆ CS_PROCF() [4/4]

void CS_PROCF ( itrmav  ,
ITRMAV   
) const

◆ cs_slope_test_gradient()

void cs_slope_test_gradient ( int  f_id,
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_real_t coefap,
const cs_real_t coefbp,
const cs_real_t i_massflux 
)

Compute the upwind gradient used in the slope tests.

This function assumes the input gradient and pvar values have already been synchronized.

Parameters
[in]f_idfield id
[in]incNot an increment flag
[in]halo_typehalo type
[in]gradstandard gradient
[out]grdpaupwind gradient
[in]pvarvalues
[in]coefapboundary condition array for the variable (explicit part)
[in]coefbpboundary condition array for the variable (implicit part)
[in]i_massfluxmass flux at interior faces

◆ cs_slope_test_gradient_tensor()

void cs_slope_test_gradient_tensor ( const int  inc,
const cs_halo_type_t  halo_type,
const cs_real_63_t grad,
cs_real_63_t grdpa,
const cs_real_6_t pvar,
const cs_real_6_t coefa,
const cs_real_66_t coefb,
const cs_real_t i_massflux 
)

Compute the upwind gradient used in the slope tests.

This function assumes the input gradient and pvar values have already been synchronized.

Parameters
[in]incNot an increment flag
[in]halo_typehalo type
[in]gradstandard gradient
[out]grdpaupwind gradient
[in]pvarvalues
[in]coefaboundary condition array for the variable (Explicit part)
[in]coefbboundary condition array for the variable (Implicit part)
[in]i_massfluxmass flux at interior faces

◆ cs_slope_test_gradient_vector()

void cs_slope_test_gradient_vector ( const int  inc,
const cs_halo_type_t  halo_type,
const cs_real_33_t grad,
cs_real_33_t grdpa,
const cs_real_3_t pvar,
const cs_real_3_t coefa,
const cs_real_33_t coefb,
const cs_real_t i_massflux 
)

Compute the upwind gradient used in the slope tests.

This function assumes the input gradient and pvar values have already been synchronized.

Parameters
[in]incNot an increment flag
[in]halo_typehalo type
[in]gradstandard gradient
[out]grdpaupwind gradient
[in]pvarvalues
[in]coefaboundary condition array for the variable (Explicit part)
[in]coefbboundary condition array for the variable (Implicit part)
[in]i_massfluxmass flux at interior faces

◆ cs_upwind_gradient()

void cs_upwind_gradient ( const int  f_id,
const int  inc,
const cs_halo_type_t  halo_type,
const cs_real_t  coefap[],
const cs_real_t  coefbp[],
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 in order to cope with SOLU schemes observed in the litterature.

Parameters
[in]f_idfield index
[in]incNot an increment flag
[in]halo_typehalo type
[in]coefapboundary condition array for the variable (explicit part)
[in]coefbpboundary condition array for the variable (implicit part)
[in]i_massfluxmass flux at interior faces
[in]b_massfluxmass flux at boundary faces
[in]pvarvalues
[out]grdpaupwind gradient