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

Wrapper to the function which adds the explicit part of the convection/diffusion terms of a transport equation of a field. 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_convection_diffusion.h"
#include "cs_balance.h"
Include dependency graph for cs_balance.c:

Functions

void cs_balance_scalar (int idtvar, int f_id, int imucpp, int imasac, int inc, int iccocg, cs_var_cal_opt_t *var_cal_opt, 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_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t viscel[], const cs_real_t xcpp[], const cs_real_2_t weighf[], const cs_real_t weighb[], int icvflb, const int icvfli[], cs_real_t smbrp[])
 Wrapper to the function which adds the explicit part of the convection/diffusion terms of a transport equation of a scalar field $ \varia $. More...
 
void cs_balance_vector (int idtvar, int f_id, int imasac, int inc, int ivisep, cs_var_cal_opt_t *var_cal_opt, cs_real_3_t pvar[], const cs_real_3_t 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_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t secvif[], const cs_real_t secvib[], cs_real_6_t viscel[], const cs_real_2_t weighf[], const cs_real_t weighb[], int icvflb, const int icvfli[], cs_real_3_t smbr[])
 Wrapper to the function which adds the explicit part of the convection/diffusion terms of a transport equation of a vector field $ \vect{\varia} $. More...
 
void cs_balance_tensor (int idtvar, int f_id, int imasac, int inc, cs_var_cal_opt_t *var_cal_opt, cs_real_6_t pvar[], const cs_real_6_t 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 viscel[], const cs_real_2_t weighf[], const cs_real_t weighb[], int icvflb, const int icvfli[], cs_real_6_t smbrp[])
 Wrapper to the function which adds the explicit part of the convection/diffusion terms of a transport equation of a tensor field $ \tens{\varia} $. More...
 

Detailed Description

Wrapper to the function which adds the explicit part of the convection/diffusion terms of a transport equation of a field.

Function Documentation

◆ cs_balance_scalar()

void cs_balance_scalar ( int  idtvar,
int  f_id,
int  imucpp,
int  imasac,
int  inc,
int  iccocg,
cs_var_cal_opt_t var_cal_opt,
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_massflux[],
const cs_real_t  b_massflux[],
const cs_real_t  i_visc[],
const cs_real_t  b_visc[],
cs_real_6_t  viscel[],
const cs_real_t  xcpp[],
const cs_real_2_t  weighf[],
const cs_real_t  weighb[],
int  icvflb,
const int  icvfli[],
cs_real_t  smbrp[] 
)

Wrapper to the function which adds the explicit part of the convection/diffusion terms of 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( \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 bilsca!
  • mind the minus sign

Options for the convective scheme:

  • blencp = 0: upwind scheme for the advection
  • blencp = 1: no upwind scheme except in the slope test
  • ischcp = 0: second order
  • ischcp = 1: centered
  • imucpp = 0: do not multiply the convective part by $ C_p $
  • imucpp = 1: multiply the convective part by $ C_p $
Parameters
[in]idtvarindicator of the temporal scheme
[in]f_idfield id (or -1)
[in]imucppindicator
  • 0 do not multiply the convectiv term by Cp
  • 1 do multiply the convectiv term by Cp
[in]imasactake mass accumulation into account?
[in]incindicator
  • 0 when solving an increment
  • 1 otherwise
[in]iccocgindicator
  • 1 re-compute cocg matrix (for iterative gradients)
  • 0 otherwise
[in]var_cal_optpointer to a cs_var_cal_opt_t structure which contains variable calculation options
[in]pvarsolved variable (current time step) may be NULL if pvara != NULL
[in]pvarasolved variable (previous time step) may be NULL if pvar != NULL
[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 boundary faces for the r.h.s.
[in]viscelsymmetric cell tensor $ \tens{\mu}_\celli $
[in]xcpparray of specific heat (Cp)
[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]icvflbglobal indicator of boundary convection flux
  • 0 upwind scheme at all boundary faces
  • 1 imposed flux at some boundary faces
[in]icvfliboundary face indicator array of convection flux
  • 0 upwind scheme
  • 1 imposed flux
[in,out]smbrpright hand side $ \vect{Rhs} $

◆ cs_balance_tensor()

void cs_balance_tensor ( int  idtvar,
int  f_id,
int  imasac,
int  inc,
cs_var_cal_opt_t var_cal_opt,
cs_real_6_t  pvar[],
const cs_real_6_t  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  viscel[],
const cs_real_2_t  weighf[],
const cs_real_t  weighb[],
int  icvflb,
const int  icvfli[],
cs_real_6_t  smbrp[] 
)

Wrapper to the function which adds the explicit part of the convection/diffusion terms of a transport equation of a tensor field $ \tens{\varia} $.

More precisely, the right hand side $ \vect{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 bilscts!
  • mind the sign minus

Options for the convective scheme:

  • blencp = 0: upwind scheme for the advection
  • blencp = 1: no upwind scheme except in the slope test
  • ischcp = 0: second order
  • ischcp = 1: centered
Parameters
[in]idtvarindicator of the temporal scheme
[in]f_idfield id (or -1)
[in]imasactake mass accumulation into account?
[in]incindicator
[in]var_cal_optpointer to a cs_var_cal_opt_t structure which contains variable calculation options
[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 (Impplicit 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 boundary 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]icvflbglobal indicator of boundary convection flux
  • 0 upwind scheme at all boundary faces
  • 1 imposed flux at some boundary faces
[in]icvfliboundary face indicator array of convection flux
  • 0 upwind scheme
  • 1 imposed flux
[in,out]smbrpright hand side $ \vect{Rhs} $

◆ cs_balance_vector()

void cs_balance_vector ( int  idtvar,
int  f_id,
int  imasac,
int  inc,
int  ivisep,
cs_var_cal_opt_t var_cal_opt,
cs_real_3_t  pvar[],
const cs_real_3_t  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_massflux[],
const cs_real_t  b_massflux[],
const cs_real_t  i_visc[],
const cs_real_t  b_visc[],
const cs_real_t  secvif[],
const cs_real_t  secvib[],
cs_real_6_t  viscel[],
const cs_real_2_t  weighf[],
const cs_real_t  weighb[],
int  icvflb,
const int  icvfli[],
cs_real_3_t  smbr[] 
)

Wrapper to the function which adds 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 bilscv!
  • mind the sign minus

Options for the convective scheme:

  • blencp = 0: upwind scheme for the advection
  • blencp = 1: no upwind scheme except in the slope test
  • ischcp = 0: second order
  • ischcp = 1: centered
Parameters
[in]idtvarindicator of the temporal scheme
[in]f_idfield id (or -1)
[in]imasactake mass accumulation into account?
[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]var_cal_optpointer to a cs_var_cal_opt_t structure which contains variable calculation options
[in]pvarsolved velocity (current time step)
[in]pvarasolved velocity (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_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 boundary faces for the r.h.s.
[in]secvifsecondary viscosity at interior faces
[in]secvibsecondary viscosity at boundary faces
[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]icvflbglobal indicator of boundary convection flux
  • 0 upwind scheme at all boundary faces
  • 1 imposed flux at some boundary faces
[in]icvfliboundary face indicator array of convection flux
  • 0 upwind scheme
  • 1 imposed flux
[in,out]smbrright hand side $ \vect{Rhs} $