8.2
general documentation
cs_convection_diffusion.cpp File Reference
#include "cs_defs.h"
#include <cmath>
#include <chrono>
#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_array.h"
#include "cs_blas.h"
#include "cs_bad_cells_regularisation.h"
#include "cs_boundary_conditions.h"
#include "cs_dispatch.h"
#include "cs_equation_param.h"
#include "cs_halo.h"
#include "cs_halo_perio.h"
#include "cs_log.h"
#include "cs_internal_coupling.h"
#include "cs_math.h"
#include "cs_mesh.h"
#include "cs_field.h"
#include "cs_field_default.h"
#include "cs_field_operator.h"
#include "cs_field_pointer.h"
#include "cs_gradient.h"
#include "cs_mesh_quantities.h"
#include "cs_parall.h"
#include "cs_parameters.h"
#include "cs_porous_model.h"
#include "cs_prototypes.h"
#include "cs_timer.h"
#include "cs_velocity_pressure.h"
#include "cs_convection_diffusion.h"
#include "cs_convection_diffusion_priv.h"
+ Include dependency graph for cs_convection_diffusion.cpp:

Functions

void cs_cell_courant_number (const int f_id, cs_real_t *courant)
 
cs_real_tcs_get_v_slope_test (int f_id, const cs_equation_param_t eqp)
 
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 pair of values). More...
 
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 *restrict pvar, const cs_real_t *restrict 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 *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_equation_param_t eqp, int icvflb, int inc, int imasac, cs_real_t *restrict pvar, const cs_real_t *restrict 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_t i_conv_flux[][2], 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_equation_param_t eqp, int icvflb, int inc, int ivisep, int imasac, cs_real_3_t *restrict pvar, const cs_real_3_t *restrict pvara, const int icvfli[], const cs_field_bc_coeffs_t *bc_coeffs_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 *restrict i_pvar, cs_real_3_t *restrict b_pvar, 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_equation_param_t eqp, int icvflb, int inc, int imasac, cs_real_6_t *restrict pvar, const cs_real_6_t *restrict pvara, const cs_field_bc_coeffs_t *bc_coeffs_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 *restrict rhs)
 Add the explicit part of the convection/diffusion terms of a transport equation of a tensor field \( \tens{\varia} \). More...
 
void cs_convection_diffusion_thermal (int idtvar, int f_id, const cs_equation_param_t eqp, int inc, int imasac, cs_real_t *restrict pvar, const cs_real_t *restrict 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 *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_equation_param_t eqp, int inc, cs_real_t *restrict pvar, const cs_real_t *restrict 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 *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_equation_param_t eqp, int inc, int ivisep, cs_real_3_t *restrict pvar, const cs_real_3_t *restrict pvara, const cs_field_bc_coeffs_t *bc_coeffs_v, 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_equation_param_t eqp, int inc, cs_real_3_t *restrict pvar, const cs_real_3_t *restrict pvara, const cs_field_bc_coeffs_t *bc_coeffs_v, 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_equation_param_t eqp, int inc, cs_real_6_t *restrict pvar, const cs_real_6_t *restrict pvara, const cs_field_bc_coeffs_t *bc_coeffs_ts, 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 nswrgp, int imligp, int iphydp, int iwgrp, int iwarnp, double epsrgp, double climgp, cs_real_3_t *restrict frcxt, cs_real_t *restrict pvar, const cs_field_bc_coeffs_t *bc_coeffs, 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 nswrgp, int imligp, int ircflp, int iphydp, int iwgrp, int iwarnp, double epsrgp, double climgp, cs_real_3_t *restrict frcxt, cs_real_t *restrict 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 *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 nswrgp, int imligp, int iphydp, int iwgrp, int iwarnp, double epsrgp, double climgp, cs_real_3_t *restrict frcxt, cs_real_t *restrict 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 *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 nswrgp, int imligp, int ircflp, int iphydp, int iwgrp, int iwarnp, double epsrgp, double climgp, cs_real_3_t *restrict frcxt, cs_real_t *restrict 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 *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...
 

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  nswrgp,
int  imligp,
int  ircflp,
int  iphydp,
int  iwgrp,
int  iwarnp,
double  epsrgp,
double  climgp,
cs_real_3_t *restrict  frcxt,
cs_real_t *restrict  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 *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 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]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]frcxtbody force creating the hydrostatic pressure
[in]pvarsolved variable (pressure)
[in]bc_coeffsboundary condition structure for the variable
[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_equation_param_t  eqp,
int  inc,
cs_real_t *restrict  pvar,
const cs_real_t *restrict  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 *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]eqpequation parameters
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]pvarsolved variable (current time step)
[in]pvarasolved variable (previous time step)
[in]bc_coeffsboundary condition structure for the variable
[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_equation_param_t  eqp,
int  inc,
cs_real_6_t *restrict  pvar,
const cs_real_6_t *restrict  pvara,
const cs_field_bc_coeffs_t bc_coeffs_ts,
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]eqpequation parameters
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]pvarsolved variable (current time step)
[in]pvarasolved variable (previous time step)
[in]bc_coeffs_tsboundary condition structure for the variable
[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_equation_param_t  eqp,
int  inc,
int  ivisep,
cs_real_3_t *restrict  pvar,
const cs_real_3_t *restrict  pvara,
const cs_field_bc_coeffs_t bc_coeffs_v,
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]eqpequation parameters
[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]bc_coeffs_vboundary condition structure for the variable
[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_equation_param_t  eqp,
int  inc,
cs_real_3_t *restrict  pvar,
const cs_real_3_t *restrict  pvara,
const cs_field_bc_coeffs_t bc_coeffs_v,
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]eqpequation parameters
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]pvarsolved variable (current time step)
[in]pvarasolved variable (previous time step)
[in]bc_coeffs_vboundary condition structure for the variable
[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_beta_limiter_building()

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 pair of values).

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

◆ cs_cell_courant_number()

void cs_cell_courant_number ( const int  f_id,
cs_real_t courant 
)

◆ cs_convection_diffusion_scalar()

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 *restrict  pvar,
const cs_real_t *restrict  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 *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]eqpequation parameters
[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 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]bc_coeffsboundary condition structure for the variable
[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_equation_param_t  eqp,
int  icvflb,
int  inc,
int  imasac,
cs_real_6_t *restrict  pvar,
const cs_real_6_t *restrict  pvara,
const cs_field_bc_coeffs_t bc_coeffs_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 *restrict  rhs 
)

Add the explicit part of the convection/diffusion terms of a transport equation of a tensor field \( \tens{\varia} \).

More precisely, the right hand side \( \tens{Rhs} \) is updated as follows:

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

Warning:

  • \( \tens{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]eqpequation parameters
[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]bc_coeffs_tsboundary condition structure for the variable
[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 \( \tens{Rhs} \)

◆ cs_convection_diffusion_thermal()

void cs_convection_diffusion_thermal ( int  idtvar,
int  f_id,
const cs_equation_param_t  eqp,
int  inc,
int  imasac,
cs_real_t *restrict  pvar,
const cs_real_t *restrict  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 *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]eqpequation parameters)
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]imasactake mass accumulation into account?
[in]pvarsolved variable (current time step)
[in]pvarasolved variable (previous time step)
[in]bc_coeffsboundary condition structure for the variable
[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_equation_param_t  eqp,
int  icvflb,
int  inc,
int  ivisep,
int  imasac,
cs_real_3_t *restrict  pvar,
const cs_real_3_t *restrict  pvara,
const int  icvfli[],
const cs_field_bc_coeffs_t bc_coeffs_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 *restrict  i_pvar,
cs_real_3_t *restrict  b_pvar,
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]eqpequation parameters
[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]bc_coeffs_vboundary conditions structure for the variable
[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]i_pvarvelocity at interior faces
[in]b_pvarvelocity 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  nswrgp,
int  imligp,
int  iphydp,
int  iwgrp,
int  iwarnp,
double  epsrgp,
double  climgp,
cs_real_3_t *restrict  frcxt,
cs_real_t *restrict  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 *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]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]iwgrpindicator
  • 1 weight gradient by vicosity*porosity
  • weighting determined by field options
[in]epsrgprelative precision for the gradient reconstruction
[in]climgpclipping coeffecient for the computation of the gradient
[in]frcxtbody force creating the hydrostatic pressure
[in]pvarsolved variable (current time step)
[in]bc_coeffsboundary condition structure for the variable
[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  nswrgp,
int  imligp,
int  ircflp,
int  iphydp,
int  iwgrp,
int  iwarnp,
double  epsrgp,
double  climgp,
cs_real_3_t *restrict  frcxt,
cs_real_t *restrict  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 *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]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]frcxtbody force creating the hydrostatic pressure
[in]pvarsolved variable (pressure)
[in]bc_coeffsboundary condition structure for the variable
[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_equation_param_t  eqp,
int  icvflb,
int  inc,
int  imasac,
cs_real_t *restrict  pvar,
const cs_real_t *restrict  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_t  i_conv_flux[][2],
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]f_idpointer to field
[in]eqpequation parameters
[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 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]bc_coeffsboundary condition structure for the variable
[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  nswrgp,
int  imligp,
int  iphydp,
int  iwgrp,
int  iwarnp,
double  epsrgp,
double  climgp,
cs_real_3_t *restrict  frcxt,
cs_real_t *restrict  pvar,
const cs_field_bc_coeffs_t bc_coeffs,
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 cs_face_diffusion_potential/cs_diffusion_potential section of the theory guide for more information.

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]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]frcxtbody force creating the hydrostatic pressure
[in]pvarsolved variable (current time step)
[in]bc_coeffsboundary condition structure for the variable
[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_get_v_slope_test()

cs_real_t* cs_get_v_slope_test ( int  f_id,
const cs_equation_param_t  eqp 
)