#include "cs_defs.h"
#include "cs_base.h"
#include "cs_halo.h"
#include "cs_math.h"
#include "cs_mesh_quantities.h"
#include "cs_parameters.h"
#include "cs_field_pointer.h"
Go to the source code of this file.
Enumerations | |
enum | cs_nvd_type_t { CS_NVD_GAMMA = 0, CS_NVD_SMART = 1, CS_NVD_CUBISTA = 2, CS_NVD_SUPERBEE = 3, CS_NVD_MUSCL = 4, CS_NVD_MINMOD = 5, CS_NVD_CLAM = 6, CS_NVD_STOIC = 7, CS_NVD_OSHER = 8, CS_NVD_WASEB = 9, CS_NVD_VOF_HRIC = 10, CS_NVD_VOF_CICSAM = 11, CS_NVD_VOF_STACS = 12, CS_NVD_N_TYPES = 13 } |
Functions | |
void | itrmas (const int *const f_id, const int *const init, const int *const inc, const int *const imrgra, const int *const iccocg, const int *const nswrgp, const int *const imligp, const int *const iphydp, const int *const iwgrp, const int *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_t visel[], cs_real_t i_massflux[], cs_real_t b_massflux[]) |
void | itrmav (const int *const f_id, const int *const init, const int *const inc, const int *const imrgra, const int *const iccocg, const int *const nswrgp, const int *const imligp, const int *const ircflp, const int *const iphydp, const int *const iwgrp, const int *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 int *const f_id, const int *const init, const int *const inc, const int *const imrgra, const int *const iccocg, const int *const nswrgp, const int *const imligp, const int *const iphydp, const int *const iwgrp, const int *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_t visel[], cs_real_t diverg[]) |
void | itrgrv (const int *const f_id, const int *const init, const int *const inc, const int *const imrgra, const int *const iccocg, const int *const nswrgp, const int *const imligp, const int *const ircflp, const int *const iphydp, const int *const iwgrp, const int *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_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_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_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 int 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 ![]() | |
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 int 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 ![]() | |
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 int 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 ![]() | |
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 ![]() | |
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 ![]() | |
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 ![]() | |
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 ![]() | |
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 ![]() | |
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 ![]() | |
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, 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, 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 ![]() | |
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 iwgrp, int iwarnp, double epsrgp, double climgp, 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 iwgrp, int iwarnp, double epsrgp, double climgp, 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... | |
enum cs_nvd_type_t |
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 | iwgrp, | ||
int | iwarnp, | ||
double | epsrgp, | ||
double | climgp, | ||
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 is updated as follows:
[in] | f_id | field id (or -1) |
[in] | m | pointer to mesh |
[in] | fvq | pointer to finite volume quantities |
[in] | init | indicator
|
[in] | inc | indicator
|
[in] | imrgra | indicator
|
[in] | iccocg | indicator
|
[in] | nswrgp | number of reconstruction sweeps for the gradients |
[in] | imligp | clipping gradient method
|
[in] | ircflp | indicator
|
[in] | iphydp | indicator
|
[in] | iwgrp | indicator
|
[in] | iwarnp | verbosity |
[in] | epsrgp | relative precision for the gradient reconstruction |
[in] | climgp | clipping coeffecient for the computation of the gradient |
[in] | frcxt | body force creating the hydrostatic pressure |
[in] | pvar | solved variable (pressure) |
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | cofafp | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfp | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in] | viscel | symmetric cell tensor ![]() |
[in] | weighf | internal face weight between cells i j in case of tensor diffusion |
[in] | weighb | boundary face weight for cells i in case of tensor diffusion |
[in,out] | diverg | divergence 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 is updated as follows:
[in] | f_id | field id (or -1) |
[in] | m | pointer to mesh |
[in] | fvq | pointer to finite volume quantities |
[in] | init | indicator
|
[in] | inc | indicator
|
[in] | imrgra | indicator
|
[in] | iccocg | indicator
|
[in] | nswrgp | number of reconstruction sweeps for the gradients |
[in] | imligp | clipping gradient method
|
[in] | ircflp | indicator
|
[in] | iphydp | indicator
|
[in] | iwgrp | indicator
|
[in] | iwarnp | verbosity |
[in] | epsrgp | relative precision for the gradient reconstruction |
[in] | climgp | clipping coeffecient for the computation of the gradient |
[in] | frcxt | body force creating the hydrostatic pressure |
[in] | pvar | solved variable (pressure) |
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | cofafp | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfp | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in] | viscel | symmetric cell tensor ![]() |
[in] | weighf | internal face weight between cells i j in case of tensor diffusion |
[in] | weighb | boundary face weight for cells i in case of tensor diffusion |
[in,out] | diverg | divergence 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 .
More precisely, the right hand side is updated as follows:
Warning:
[in] | idtvar | indicator of the temporal scheme |
[in] | f_id | index of the current variable |
[in] | var_cal_opt | variable calculation options |
[in] | inc | indicator
|
[in] | iccocg | indicator
|
[in] | pvar | solved variable (current time step) |
[in] | pvara | solved variable (previous time step) |
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | cofafp | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfp | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in] | viscel | symmetric cell tensor ![]() |
[in] | weighf | internal face weight between cells i j in case of tensor diffusion |
[in] | weighb | boundary face weight for cells i in case of tensor diffusion |
[in,out] | rhs | right hand side ![]() |
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 .
More precisely, the right hand side is updated as follows:
Warning:
[in] | idtvar | indicator of the temporal scheme |
[in] | f_id | index of the current variable |
[in] | var_cal_opt | variable calculation options |
[in] | inc | indicator
|
[in] | pvar | solved variable (current time step) |
[in] | pvara | solved variable (previous time step) |
[in] | coefa | boundary condition array for the variable (explicit part) |
[in] | coefb | boundary condition array for the variable (implicit part) |
[in] | cofaf | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbf | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in] | viscel | symmetric cell tensor ![]() |
[in] | weighf | internal face weight between cells i j in case of tensor diffusion |
[in] | weighb | boundary face weight for cells i in case of tensor diffusion |
[in,out] | rhs | right hand side ![]() |
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 .
More precisely, the right hand side is updated as follows:
Remark: if ivisep = 1, then we also take , where
is the secondary viscosity, i.e. usually
.
Warning:
[in] | idtvar | indicator of the temporal scheme |
[in] | f_id | index of the current variable |
[in] | var_cal_opt | variable calculation options |
[in] | inc | indicator
|
[in] | ivisep | indicator to take ![]()
|
[in] | pvar | solved variable (current time step) |
[in] | pvara | solved variable (previous time step) |
[in] | coefav | boundary condition array for the variable (explicit part) |
[in] | coefbv | boundary condition array for the variable (implicit part) |
[in] | cofafv | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfv | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in] | i_secvis | secondary viscosity at interior faces |
[in,out] | rhs | right hand side ![]() |
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 .
More precisely, the right hand side is updated as follows:
Warning:
[in] | idtvar | indicator of the temporal scheme |
[in] | f_id | index of the current variable |
[in] | var_cal_opt | variable calculation options |
[in] | inc | indicator
|
[in] | pvar | solved variable (current time step) |
[in] | pvara | solved variable (previous time step) |
[in] | coefav | boundary condition array for the variable (explicit part) |
[in] | coefbv | boundary condition array for the variable (implicit part) |
[in] | cofafv | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfv | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in] | viscel | symmetric cell tensor ![]() |
[in] | weighf | internal face weight between cells i j in case of tensor diffusion |
[in] | weighb | boundary face weight for cells i in case of tensor diffusion |
[in,out] | rhs | right hand side ![]() |
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).
[in] | f_id | field id |
[in] | inc | "not an increment" flag |
[in] | rovsdt | rho * volume / dt |
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 int | 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 .
More precisely, the right hand side is updated as follows:
Warning:
[in] | idtvar | indicator of the temporal scheme |
[in] | f_id | field id (or -1) |
[in] | var_cal_opt | variable calculation options |
[in] | icvflb | global indicator of boundary convection flux
|
[in] | inc | indicator
|
[in] | iccocg | indicator
|
[in] | imasac | take mass accumulation into account? |
[in] | pvar | solved variable (current time step) |
[in] | pvara | solved variable (previous time step) |
[in] | icvfli | boundary face indicator array of convection flux
|
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | cofafp | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfp | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_massflux | mass flux at interior faces |
[in] | b_massflux | mass flux at boundary faces |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in,out] | rhs | right hand side ![]() |
More precisely, the right hand side is updated as follows:
Warning:
Please refer to the bilsc2 section of the theory guide for more informations.
[in] | idtvar | indicator of the temporal scheme |
[in] | f_id | field id (or -1) |
[in] | var_cal_opt | variable calculation options |
[in] | icvflb | global indicator of boundary convection flux
|
[in] | inc | indicator
|
[in] | iccocg | indicator
|
[in] | imasac | take mass accumulation into account? |
[in] | pvar | solved variable (current time step) |
[in] | pvara | solved variable (previous time step) |
[in] | icvfli | boundary face indicator array of convection flux
|
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | cofafp | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfp | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_massflux | mass flux at interior faces |
[in] | b_massflux | mass flux at boundary faces |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in,out] | rhs | right hand side ![]() |
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 .
More precisely, the right hand side is updated as follows:
Warning:
[in] | idtvar | indicator of the temporal scheme |
[in] | f_id | index of the current variable |
[in] | var_cal_opt | variable calculation options |
[in] | icvflb | global indicator of boundary convection flux
|
[in] | inc | indicator
|
[in] | imasac | take mass accumulation into account? |
[in] | pvar | solved velocity (current time step) |
[in] | pvara | solved velocity (previous time step) |
[in] | coefa | boundary condition array for the variable (Explicit part) |
[in] | coefb | boundary condition array for the variable (Implicit part) |
[in] | cofaf | boundary condition array for the diffusion of the variable (Explicit part) |
[in] | cofbf | boundary condition array for the diffusion of the variable (Implicit part) |
[in] | i_massflux | mass flux at interior faces |
[in] | b_massflux | mass flux at boundary faces |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in,out] | rhs | right hand side ![]() |
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 such as the temperature.
More precisely, the right hand side is updated as follows:
Warning: has already been initialized before calling bilsct!
[in] | idtvar | indicator of the temporal scheme |
[in] | f_id | index of the current variable |
[in] | var_cal_opt | variable calculation options |
[in] | inc | indicator
|
[in] | iccocg | indicator
|
[in] | imasac | take mass accumulation into account? |
[in] | pvar | solved variable (current time step) |
[in] | pvara | solved variable (previous time step) |
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | cofafp | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfp | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_massflux | mass flux at interior faces |
[in] | b_massflux | mass flux at boundary faces |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in] | xcpp | array of specific heat ( ![]() |
[in,out] | rhs | right hand side ![]() |
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 int | 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 .
More precisely, the right hand side is updated as follows:
Remark: if ivisep = 1, then we also take , where
is the secondary viscosity, i.e. usually
.
Warning:
[in] | idtvar | indicator of the temporal scheme |
[in] | f_id | index of the current variable |
[in] | var_cal_opt | variable calculation options |
[in] | icvflb | global indicator of boundary convection flux
|
[in] | inc | indicator
|
[in] | ivisep | indicator to take ![]()
|
[in] | imasac | take mass accumulation into account? |
[in] | pvar | solved velocity (current time step) |
[in] | pvara | solved velocity (previous time step) |
[in] | icvfli | boundary face indicator array of convection flux
|
[in] | coefav | boundary condition array for the variable (explicit part) |
[in] | coefbv | boundary condition array for the variable (implicit part) |
[in] | cofafv | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfv | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_massflux | mass flux at interior faces |
[in] | b_massflux | mass flux at boundary faces |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in] | i_secvis | secondary viscosity at interior faces |
[in] | b_secvis | secondary viscosity at boundary faces |
[in,out] | rhs | right hand side ![]() |
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 | iwgrp, | ||
int | iwarnp, | ||
double | epsrgp, | ||
double | climgp, | ||
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.
[in] | f_id | field id (or -1) |
[in] | m | pointer to mesh |
[in] | fvq | pointer to finite volume quantities |
[in] | init | indicator
|
[in] | inc | indicator
|
[in] | imrgra | indicator
|
[in] | iccocg | indicator
|
[in] | nswrgp | number of reconstruction sweeps for the gradients |
[in] | imligp | clipping gradient method
|
[in] | iphydp | hydrostatic pressure indicator |
[in] | iwgrp | indicator
|
[in] | iwarnp | verbosity |
[in] | epsrgp | relative precision for the gradient reconstruction |
[in] | climgp | clipping coeffecient for the computation of the gradient |
[in] | frcxt | body force creating the hydrostatic pressure |
[in] | pvar | solved variable (current time step) |
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | cofafp | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfp | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in] | visel | viscosity by cell |
[in,out] | diverg | mass flux divergence |
[in] | f_id | field id (or -1) |
[in] | m | pointer to mesh |
[in] | fvq | pointer to finite volume quantities |
[in] | init | indicator
|
[in] | inc | indicator
|
[in] | imrgra | indicator
|
[in] | iccocg | indicator
|
[in] | nswrgp | number of reconstruction sweeps for the gradients |
[in] | imligp | clipping gradient method
|
[in] | iphydp | hydrostatic pressure indicator |
[in] | iwarnp | verbosity |
[in] | iwgrp | indicator
|
[in] | epsrgp | relative precision for the gradient reconstruction |
[in] | climgp | clipping coeffecient for the computation of the gradient |
[in] | frcxt | body force creating the hydrostatic pressure |
[in] | pvar | solved variable (current time step) |
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | cofafp | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfp | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in] | visel | viscosity by cell |
[in,out] | diverg | mass flux divergence |
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, | ||
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 .
More precisely, the mass flux side is updated as follows:
[in] | f_id | field id (or -1) |
[in] | m | pointer to mesh |
[in] | fvq | pointer to finite volume quantities |
[in] | init | indicator
|
[in] | inc | indicator
|
[in] | imrgra | indicator
|
[in] | iccocg | indicator
|
[in] | nswrgp | number of reconstruction sweeps for the gradients |
[in] | imligp | clipping gradient method
|
[in] | ircflp | indicator
|
[in] | iphydp | indicator
|
[in] | iwgrp | indicator
|
[in] | iwarnp | verbosity |
[in] | epsrgp | relative precision for the gradient reconstruction |
[in] | climgp | clipping coeffecient for the computation of the gradient |
[in] | frcxt | body force creating the hydrostatic pressure |
[in] | pvar | solved variable (pressure) |
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | cofafp | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfp | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in] | viscel | symmetric cell tensor ![]() |
[in] | weighf | internal face weight between cells i j in case of tensor diffusion |
[in] | weighb | boundary face weight for cells i in case of tensor diffusion |
[in,out] | i_massflux | mass flux at interior faces |
[in,out] | b_massflux | mass flux at boundary faces |
More precisely, the mass flux side is updated as follows:
[in] | f_id | field id (or -1) |
[in] | m | pointer to mesh |
[in] | fvq | pointer to finite volume quantities |
[in] | init | indicator
|
[in] | inc | indicator
|
[in] | imrgra | indicator
|
[in] | iccocg | indicator
|
[in] | nswrgp | number of reconstruction sweeps for the gradients |
[in] | imligp | clipping gradient method
|
[in] | ircflp | indicator
|
[in] | iphydp | indicator
|
[in] | iwgrp | indicator
|
[in] | iwarnp | verbosity |
[in] | epsrgp | relative precision for the gradient reconstruction |
[in] | climgp | clipping coeffecient for the computation of the gradient |
[in] | frcxt | body force creating the hydrostatic pressure |
[in] | pvar | solved variable (pressure) |
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | cofafp | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfp | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in] | viscel | symmetric cell tensor ![]() |
[in] | weighf | internal face weight between cells i j in case of tensor diffusion |
[in] | weighb | boundary face weight for cells i in case of tensor diffusion |
[in,out] | i_massflux | mass flux at interior faces |
[in,out] | b_massflux | mass flux at boundary faces |
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 int | 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 .
[in] | idtvar | indicator of the temporal scheme |
[in] | f_id | field id (or -1) |
[in] | var_cal_opt | variable calculation options |
[in] | icvflb | global indicator of boundary convection flux
|
[in] | inc | indicator
|
[in] | iccocg | indicator
|
[in] | imasac | take mass accumulation into account? |
[in] | pvar | solved variable (current time step) |
[in] | pvara | solved variable (previous time step) |
[in] | icvfli | boundary face indicator array of convection flux
|
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | i_massflux | mass flux at interior faces |
[in] | b_massflux | mass flux at boundary faces |
[in,out] | i_conv_flux | scalar convection flux at interior faces |
[in,out] | b_conv_flux | scalar convection flux at boundary faces |
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, | ||
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.
[in] | f_id | field id (or -1) |
[in] | m | pointer to mesh |
[in] | fvq | pointer to finite volume quantities |
[in] | init | indicator
|
[in] | inc | indicator
|
[in] | imrgra | indicator
|
[in] | iccocg | indicator
|
[in] | nswrgp | number of reconstruction sweeps for the gradients |
[in] | imligp | clipping gradient method
|
[in] | iphydp | hydrostatic pressure indicator |
[in] | iwgrp | indicator
|
[in] | iwarnp | verbosity |
[in] | epsrgp | relative precision for the gradient reconstruction |
[in] | climgp | clipping coeffecient for the computation of the gradient |
[in] | frcxt | body force creating the hydrostatic pressure |
[in] | pvar | solved variable (current time step) |
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | cofafp | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfp | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in] | visel | viscosity by cell |
[in,out] | i_massflux | mass flux at interior faces |
[in,out] | b_massflux | mass flux at boundary faces |
Please refer to the itrmas/itrgrp section of the theory guide for more information.
[in] | f_id | field id (or -1) |
[in] | m | pointer to mesh |
[in] | fvq | pointer to finite volume quantities |
[in] | init | indicator
|
[in] | inc | indicator
|
[in] | imrgra | indicator
|
[in] | iccocg | indicator
|
[in] | nswrgp | number of reconstruction sweeps for the gradients |
[in] | imligp | clipping gradient method
|
[in] | iphydp | hydrostatic pressure indicator |
[in] | iwgrp | indicator
|
[in] | iwarnp | verbosity |
[in] | epsrgp | relative precision for the gradient reconstruction |
[in] | climgp | clipping coeffecient for the computation of the gradient |
[in] | frcxt | body force creating the hydrostatic pressure |
[in] | pvar | solved variable (current time step) |
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | cofafp | boundary condition array for the diffusion of the variable (explicit part) |
[in] | cofbfp | boundary condition array for the diffusion of the variable (implicit part) |
[in] | i_visc | ![]() |
[in] | b_visc | ![]() |
[in] | visel | viscosity by cell |
[in,out] | i_massflux | mass flux at interior faces |
[in,out] | b_massflux | mass flux at boundary faces |
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.
[in] | f_id | field id |
[in] | inc | Not an increment flag |
[in] | halo_type | halo type |
[in] | grad | standard gradient |
[out] | grdpa | upwind gradient |
[in] | pvar | values |
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | i_massflux | mass flux at interior faces |
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.
[in] | inc | Not an increment flag |
[in] | halo_type | halo type |
[in] | grad | standard gradient |
[out] | grdpa | upwind gradient |
[in] | pvar | values |
[in] | coefa | boundary condition array for the variable (explicit part) |
[in] | coefb | boundary condition array for the variable (implicit part) |
[in] | i_massflux | mass flux at interior faces |
This function assumes the input gradient and pvar values have already been synchronized.
[in] | inc | Not an increment flag |
[in] | halo_type | halo type |
[in] | grad | standard gradient |
[out] | grdpa | upwind gradient |
[in] | pvar | values |
[in] | coefa | boundary condition array for the variable (Explicit part) |
[in] | coefb | boundary condition array for the variable (Implicit part) |
[in] | i_massflux | mass flux at interior faces |
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.
[in] | inc | Not an increment flag |
[in] | halo_type | halo type |
[in] | grad | standard gradient |
[out] | grdpa | upwind gradient |
[in] | pvar | values |
[in] | coefa | boundary condition array for the variable (explicit part) |
[in] | coefb | boundary condition array for the variable (implicit part) |
[in] | i_massflux | mass flux at interior faces |
This function assumes the input gradient and pvar values have already been synchronized.
[in] | inc | Not an increment flag |
[in] | halo_type | halo type |
[in] | grad | standard gradient |
[out] | grdpa | upwind gradient |
[in] | pvar | values |
[in] | coefa | boundary condition array for the variable (Explicit part) |
[in] | coefb | boundary condition array for the variable (Implicit part) |
[in] | i_massflux | mass flux at interior faces |
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.
[in] | f_id | field index |
[in] | inc | Not an increment flag |
[in] | halo_type | halo type |
[out] | grdpa | upwind gradient |
[in] | pvar | values |
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | i_massflux | mass flux at interior faces |
Compute the upwind gradient in order to cope with SOLU schemes observed in the litterature.
[in] | f_id | field index |
[in] | inc | Not an increment flag |
[in] | halo_type | halo type |
[in] | coefap | boundary condition array for the variable (explicit part) |
[in] | coefbp | boundary condition array for the variable (implicit part) |
[in] | i_massflux | mass flux at interior faces |
[in] | b_massflux | mass flux at boundary faces |
[in] | pvar | values |
[out] | grdpa | upwind gradient |
void itrgrp | ( | const int *const | f_id, |
const int *const | init, | ||
const int *const | inc, | ||
const int *const | imrgra, | ||
const int *const | iccocg, | ||
const int *const | nswrgp, | ||
const int *const | imligp, | ||
const int *const | iphydp, | ||
const int *const | iwgrp, | ||
const int *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_t | visel[], | ||
cs_real_t | diverg[] | ||
) |
void itrgrv | ( | const int *const | f_id, |
const int *const | init, | ||
const int *const | inc, | ||
const int *const | imrgra, | ||
const int *const | iccocg, | ||
const int *const | nswrgp, | ||
const int *const | imligp, | ||
const int *const | ircflp, | ||
const int *const | iphydp, | ||
const int *const | iwgrp, | ||
const int *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 int *const | f_id, |
const int *const | init, | ||
const int *const | inc, | ||
const int *const | imrgra, | ||
const int *const | iccocg, | ||
const int *const | nswrgp, | ||
const int *const | imligp, | ||
const int *const | iphydp, | ||
const int *const | iwgrp, | ||
const int *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_t | visel[], | ||
cs_real_t | i_massflux[], | ||
cs_real_t | b_massflux[] | ||
) |
void itrmav | ( | const int *const | f_id, |
const int *const | init, | ||
const int *const | inc, | ||
const int *const | imrgra, | ||
const int *const | iccocg, | ||
const int *const | nswrgp, | ||
const int *const | imligp, | ||
const int *const | ircflp, | ||
const int *const | iphydp, | ||
const int *const | iwgrp, | ||
const int *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[] | ||
) |