|
programmer's documentation
|

Go to the source code of this file.
Functions | |
| static void | cs_slope_test (const cs_real_t pi, const cs_real_t pj, const cs_real_t distf, const cs_real_t srfan, const cs_real_3_t i_face_normal, const cs_real_3_t gradi, const cs_real_3_t gradj, const cs_real_3_t grdpai, const cs_real_3_t grdpaj, const cs_real_t i_massflux, double *testij, double *tesqck) |
| Compute slope test criteria at internal face between cell i and j. More... | |
| static void | cs_i_compute_quantities (const int ircflp, const double pnd, const cs_real_3_t cell_ceni, const cs_real_3_t cell_cenj, const cs_real_3_t i_face_cog, const cs_real_3_t dijpf, const cs_real_3_t gradi, const cs_real_3_t gradj, const cs_real_t pi, const cs_real_t pj, cs_real_t *recoi, cs_real_t *recoj, cs_real_t *pip, cs_real_t *pjp) |
| Reconstruct values in I' and J'. More... | |
| static void | cs_i_relax_c_val (const double relaxp, const cs_real_t pia, const cs_real_t pja, const cs_real_t recoi, const cs_real_t recoj, const cs_real_t pi, const cs_real_t pj, cs_real_t *pir, cs_real_t *pjr, cs_real_t *pipr, cs_real_t *pjpr) |
| Compute relaxed values at cell i and j. More... | |
| static void | cs_i_no_relax_c_val (const cs_real_t pi, const cs_real_t pj, const cs_real_t pip, const cs_real_t pjp, cs_real_t *pir, cs_real_t *pjr, cs_real_t *pipr, cs_real_t *pjpr) |
| Copy reconstructed or not cell values for consistency with relaxation case. More... | |
| static void | cs_upwind_f_val (const cs_real_t pi, const cs_real_t pj, const cs_real_t pir, const cs_real_t pjr, cs_real_t *pifri, cs_real_t *pifrj, cs_real_t *pjfri, cs_real_t *pjfrj) |
| Prepare value at face ij by using an upwind scheme. More... | |
| static void | cs_centered_f_val (const double pnd, const cs_real_t pip, const cs_real_t pjp, const cs_real_t pipr, const cs_real_t pjpr, cs_real_t *pifri, cs_real_t *pifrj, cs_real_t *pjfri, cs_real_t *pjfrj) |
| Prepare value at face ij by using a centered scheme. More... | |
| static void | cs_solu_f_val (const cs_real_3_t cell_ceni, const cs_real_3_t cell_cenj, const cs_real_3_t i_face_cog, const cs_real_3_t gradi, const cs_real_3_t gradj, const cs_real_t pi, const cs_real_t pj, const cs_real_t pir, const cs_real_t pjr, cs_real_t *pifri, cs_real_t *pifrj, cs_real_t *pjfri, cs_real_t *pjfrj) |
| Prepare value at face ij by using a Second Order Linear Upwind scheme. More... | |
| static void | cs_blend_f_val (const double blencp, const cs_real_t pi, const cs_real_t pj, const cs_real_t pir, const cs_real_t pjr, cs_real_t *pifri, cs_real_t *pifrj, cs_real_t *pjfri, cs_real_t *pjfrj) |
| Blend face values for a centered or SOLU scheme with face values for an upwind scheme. More... | |
| static void | cs_i_conv_flux (const int iconvp, const cs_real_t pi, const cs_real_t pj, const cs_real_t pifri, const cs_real_t pifrj, const cs_real_t pjfri, const cs_real_t pjfrj, const cs_real_t i_massflux, cs_real_2_t fluxij) |
| Add convective fluxes (substracting the mass accumulation from them) to fluxes at face ij. More... | |
| static void | cs_i_conv_flux_cons (const int iconvp, const cs_real_t pifri, const cs_real_t pifrj, const cs_real_t pjfri, const cs_real_t pjfrj, const cs_real_t i_massflux, const cs_real_t xcppi, const cs_real_t xcppj, cs_real_2_t fluxij) |
| Add convective fluxes (cons. formulation) to fluxes at face ij. More... | |
| static void | cs_i_diff_flux (const int idiffp, const cs_real_t pip, const cs_real_t pjp, const cs_real_t pipr, const cs_real_t pjpr, const cs_real_t i_visc, cs_real_2_t fluxij) |
| Add diffusive fluxes to fluxes at face ij. More... | |
| static void | cs_i_cd_steady_upwind (const int ircflp, const cs_real_t relaxp, const cs_real_t weight, const cs_real_3_t cell_ceni, const cs_real_3_t cell_cenj, const cs_real_3_t i_face_cog, const cs_real_3_t dijpf, const cs_real_3_t gradi, const cs_real_3_t gradj, const cs_real_t pi, const cs_real_t pj, const cs_real_t pia, const cs_real_t pja, cs_real_t *pifri, cs_real_t *pifrj, cs_real_t *pjfri, cs_real_t *pjfrj, cs_real_t *pip, cs_real_t *pjp, cs_real_t *pipr, cs_real_t *pjpr) |
| Handle preparation of internal face values for the fluxes computation in case of a steady algorithm and a pure upwind flux. More... | |
| static void | cs_i_cd_unsteady_upwind (const int ircflp, const cs_real_t weight, const cs_real_3_t cell_ceni, const cs_real_3_t cell_cenj, const cs_real_3_t i_face_cog, const cs_real_3_t dijpf, const cs_real_3_t gradi, const cs_real_3_t gradj, const cs_real_t pi, const cs_real_t pj, cs_real_t *pifri, cs_real_t *pifrj, cs_real_t *pjfri, cs_real_t *pjfrj, cs_real_t *pip, cs_real_t *pjp, cs_real_t *pipr, cs_real_t *pjpr) |
| Handle preparation of internal face values for the fluxes computation in case of an unsteady algorithm and a pure upwind flux. More... | |
| static void | cs_i_cd_steady (const int ircflp, const int ischcp, const double relaxp, const double blencp, const cs_real_t weight, const cs_real_3_t cell_ceni, const cs_real_3_t cell_cenj, const cs_real_3_t i_face_cog, const cs_real_3_t dijpf, const cs_real_3_t gradi, const cs_real_3_t gradj, const cs_real_t pi, const cs_real_t pj, const cs_real_t pia, const cs_real_t pja, cs_real_t *pifri, cs_real_t *pifrj, cs_real_t *pjfri, cs_real_t *pjfrj, cs_real_t *pip, cs_real_t *pjp, cs_real_t *pipr, cs_real_t *pjpr) |
| Handle preparation of internal face values for the fluxes computation in case of a steady algorithm and without enabling slope tests. More... | |
| static void | cs_i_cd_unsteady (const int ircflp, const int ischcp, const double blencp, const cs_real_t weight, const cs_real_3_t cell_ceni, const cs_real_3_t cell_cenj, const cs_real_3_t i_face_cog, const cs_real_3_t dijpf, const cs_real_3_t gradi, const cs_real_3_t gradj, const cs_real_t pi, const cs_real_t pj, cs_real_t *pifri, cs_real_t *pifrj, cs_real_t *pjfri, cs_real_t *pjfrj, cs_real_t *pip, cs_real_t *pjp, cs_real_t *pipr, cs_real_t *pjpr) |
| Handle preparation of internal face values for the fluxes computation in case of a unsteady algorithm and without enabling slope tests. More... | |
| static void | cs_i_cd_steady_slope_test (bool *upwind_switch, const int ircflp, const int ischcp, const double relaxp, const double blencp, const cs_real_t weight, const cs_real_t i_dist, const cs_real_t i_face_surf, const cs_real_3_t cell_ceni, const cs_real_3_t cell_cenj, const cs_real_3_t i_face_normal, const cs_real_3_t i_face_cog, const cs_real_3_t dijpf, const cs_real_t i_massflux, const cs_real_3_t gradi, const cs_real_3_t gradj, const cs_real_3_t grdpai, const cs_real_3_t grdpaj, const cs_real_t pi, const cs_real_t pj, const cs_real_t pia, const cs_real_t pja, cs_real_t *pifri, cs_real_t *pifrj, cs_real_t *pjfri, cs_real_t *pjfrj, cs_real_t *pip, cs_real_t *pjp, cs_real_t *pipr, cs_real_t *pjpr) |
| Handle preparation of internal face values for the fluxes computation in case of a steady algorithm and using slope tests. More... | |
| static void | cs_i_cd_unsteady_slope_test (bool *upwind_switch, const int ircflp, const int ischcp, const double blencp, const cs_real_t weight, const cs_real_t i_dist, const cs_real_t i_face_surf, const cs_real_3_t cell_ceni, const cs_real_3_t cell_cenj, const cs_real_3_t i_face_normal, const cs_real_3_t i_face_cog, const cs_real_3_t dijpf, const cs_real_t i_massflux, const cs_real_3_t gradi, const cs_real_3_t gradj, const cs_real_3_t grdpai, const cs_real_3_t grdpaj, const cs_real_t pi, const cs_real_t pj, cs_real_t *pifri, cs_real_t *pifrj, cs_real_t *pjfri, cs_real_t *pjfrj, cs_real_t *pip, cs_real_t *pjp, cs_real_t *pipr, cs_real_t *pjpr) |
| Handle preparation of internal face values for the fluxes computation in case of a unsteady algorithm and using slope tests. More... | |
| static void | cs_b_compute_quantities (const cs_real_3_t diipb, const cs_real_3_t gradi, const int ircflp, cs_real_t *recoi) |
| Reconstruct values in I' at boundary cell i. More... | |
| static void | cs_b_relax_c_val (const double relaxp, const cs_real_t pi, const cs_real_t pia, const cs_real_t recoi, cs_real_t *pir, cs_real_t *pipr) |
| Compute relaxed values at boundary cell i. More... | |
| static void | cs_b_no_relax_c_val (const cs_real_t pi, const cs_real_t recoi, cs_real_t *pir, cs_real_t *pipr) |
| Copy values at bounadry cell i for consistency with relaxation case. More... | |
| static void | cs_b_imposed_conv_flux (const int iconvp, const int inc, const int ifaccp, const cs_int_t bc_type, const int icvfli, const cs_real_t pi, const cs_real_t pir, const cs_real_t pipr, const cs_real_t coefap, const cs_real_t coefbp, const cs_real_t coface, const cs_real_t cofbce, const cs_real_t b_massflux, cs_real_t *flux) |
| Add convective flux (substracting the mass accumulation from it) to flux at boundary face. The convective flux can be either an upwind flux or an imposed value. More... | |
| static void | cs_b_imposed_conv_flux_cons (const int iconvp, const int inc, const int ifaccp, const cs_int_t bc_type, const int icvfli, const cs_real_t pir, const cs_real_t pipr, const cs_real_t coefap, const cs_real_t coefbp, const cs_real_t coface, const cs_real_t cofbce, const cs_real_t xcpp, const cs_real_t b_massflux, cs_real_t *flux) |
| Add convective flux (conservative formulation) to flux at boundary face. The convective flux can be either an upwind flux or an imposed value. More... | |
| static void | cs_b_upwind_flux (const int iconvp, const int inc, const int ifaccp, const int bc_type, const cs_real_t pi, const cs_real_t pir, const cs_real_t pipr, const cs_real_t coefap, const cs_real_t coefbp, const cs_real_t b_massflux, cs_real_t *flux) |
| Add convective flux (substracting the mass accumulation from it) to flux at boundary face. The convective flux is a pure upwind flux. More... | |
| static void | cs_b_upwind_flux_cons (const int iconvp, const int inc, const int ifaccp, const cs_int_t bc_type, const cs_real_t pir, const cs_real_t pipr, const cs_real_t coefap, const cs_real_t coefbp, const cs_real_t b_massflux, const cs_real_t xcpp, cs_real_t *flux) |
| Add convective flux (conservative formulation) to flux at boundary face. The convective flux is a pure upwind flux. More... | |
| static void | cs_b_diff_flux (const int idiffp, const int inc, const cs_real_t pipr, const cs_real_t cofafp, const cs_real_t cofbfp, const cs_real_t b_visc, cs_real_t *flux) |
| Add diffusive flux to flux at boundary face. More... | |
| static void | cs_b_cd_steady (const int ircflp, const double relaxp, const cs_real_3_t diipb, const cs_real_3_t gradi, const cs_real_t pi, const cs_real_t pia, cs_real_t *pir, cs_real_t *pipr) |
| Handle preparation of boundary face values for the flux computation in case of a steady algorithm. More... | |
| static void | cs_b_cd_unsteady (const int ircflp, const cs_real_3_t diipb, const cs_real_3_t gradi, const cs_real_t pi, cs_real_t *pir, cs_real_t *pipr) |
| Handle preparation of boundary face values for the flux computation in case of an unsteady algorithm. More... | |
| static void | cs_slope_test_gradient (const int f_id, const int inc, const cs_halo_type_t halo_type, cs_real_3_t *grad, cs_real_3_t *grdpa, 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 | bilsc2 (const cs_int_t *const idtvar, const cs_int_t *const f_id, const cs_var_cal_opt_t *const var_cal_opt, const cs_int_t *const icvflb, const cs_int_t *const inc, const cs_int_t *const iccocg, const cs_int_t *const ifaccp, cs_real_t pvar[], const cs_real_t pvara[], const cs_int_t bc_type[], const cs_int_t icvfli[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_t rhs[]) |
| void | bilsc4 (const cs_int_t *const idtvar, const cs_int_t *const f_id, const cs_var_cal_opt_t *const var_cal_opt, const cs_int_t *const icvflb, const cs_int_t *const inc, const cs_int_t *const ifaccp, const cs_int_t *const ivisep, cs_real_3_t pvar[], const cs_real_3_t pvara[], const cs_int_t bc_type[], const cs_int_t icvfli[], const cs_real_3_t coefav[], const cs_real_33_t coefbv[], const cs_real_3_t cofafv[], const cs_real_33_t cofbfv[], const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t secvif[], cs_real_3_t rhs[]) |
| void | bilsct (const cs_int_t *const idtvar, const cs_int_t *const f_id, const cs_var_cal_opt_t *const var_cal_opt, const cs_int_t *const inc, const cs_int_t *const iccocg, const cs_int_t *const ifaccp, cs_real_t pvar[], const cs_real_t pvara[], const cs_int_t bc_type[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t xcpp[], cs_real_t rhs[]) |
| void | diften (const cs_int_t *const idtvar, const cs_int_t *const f_id, const cs_var_cal_opt_t *const var_cal_opt, const cs_int_t *const inc, const cs_int_t *const iccocg, cs_real_t pvar[], const cs_real_t pvara[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t viscel[], const cs_real_2_t weighf[], const cs_real_t weighb[], cs_real_t rhs[]) |
| void | diftnv (const cs_int_t *const idtvar, const cs_int_t *const f_id, const cs_var_cal_opt_t *const var_cal_opt, const cs_int_t *const inc, const cs_int_t *const ifaccp, const cs_int_t *const ivisep, cs_real_3_t pvar[], const cs_real_3_t pvara[], const cs_int_t bc_type[], const cs_real_3_t coefav[], const cs_real_33_t coefbv[], const cs_real_3_t cofafv[], const cs_real_33_t cofbfv[], const cs_real_33_t i_visc[], const cs_real_t b_visc[], const cs_real_t secvif[], cs_real_3_t rhs[]) |
| void | itrmas (const cs_int_t *const f_id, const cs_int_t *const init, const cs_int_t *const inc, const cs_int_t *const imrgra, const cs_int_t *const iccocg, const cs_int_t *const nswrgp, const cs_int_t *const imligp, const cs_int_t *const iphydp, const cs_int_t *const iwarnp, const cs_real_t *const epsrgp, const cs_real_t *const climgp, const cs_real_t *const extrap, cs_real_3_t frcxt[], cs_real_t pvar[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t viselx[], const cs_real_t visely[], const cs_real_t viselz[], cs_real_t i_massflux[], cs_real_t b_massflux[]) |
| void | itrmav (const cs_int_t *const init, const cs_int_t *const inc, const cs_int_t *const imrgra, const cs_int_t *const iccocg, const cs_int_t *const nswrgp, const cs_int_t *const imligp, const cs_int_t *const ircflp, const cs_int_t *const iphydp, const cs_int_t *const iwarnp, const cs_real_t *const epsrgp, const cs_real_t *const climgp, const cs_real_t *const extrap, cs_real_3_t frcxt[], cs_real_t pvar[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t viscel[], const cs_real_2_t weighf[], const cs_real_t weighb[], cs_real_t i_massflux[], cs_real_t b_massflux[]) |
| void | itrgrp (const cs_int_t *const f_id, const cs_int_t *const init, const cs_int_t *const inc, const cs_int_t *const imrgra, const cs_int_t *const iccocg, const cs_int_t *const nswrgp, const cs_int_t *const imligp, const cs_int_t *const iphydp, const cs_int_t *const iwarnp, const cs_real_t *const epsrgp, const cs_real_t *const climgp, const cs_real_t *const extrap, cs_real_3_t frcxt[], cs_real_t pvar[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t viselx[], const cs_real_t visely[], const cs_real_t viselz[], cs_real_t diverg[]) |
| void | itrgrv (const cs_int_t *const init, const cs_int_t *const inc, const cs_int_t *const imrgra, const cs_int_t *const iccocg, const cs_int_t *const nswrgp, const cs_int_t *const imligp, const cs_int_t *const ircflp, const cs_int_t *const iphydp, const cs_int_t *const iwarnp, const cs_real_t *const epsrgp, const cs_real_t *const climgp, const cs_real_t *const extrap, cs_real_3_t frcxt[], cs_real_t pvar[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t viscel[], const cs_real_2_t weighf[], const cs_real_t weighb[], cs_real_t diverg[]) |
| void | cs_convection_diffusion_scalar (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int icvflb, int inc, int iccocg, int ifaccp, cs_real_t *restrict pvar, const cs_real_t *restrict pvara, const cs_int_t bc_type[], const cs_int_t icvfli[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_t *restrict rhs) |
Add the explicit part of the convection/diffusion terms of a standard transport equation of a scalar field . More... | |
| void | cs_convection_diffusion_vector (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int icvflb, int inc, int ifaccp, int ivisep, cs_real_3_t *restrict pvar, const cs_real_3_t *restrict pvara, const cs_int_t bc_type[], const cs_int_t icvfli[], const cs_real_3_t coefav[], const cs_real_33_t coefbv[], const cs_real_3_t cofafv[], const cs_real_33_t cofbfv[], const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t secvif[], cs_real_3_t *restrict rhs) |
Add the explicit part of the convection/diffusion terms of a transport equation of a vector field . More... | |
| void | cs_convection_diffusion_thermal (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int inc, int iccocg, int ifaccp, cs_real_t *restrict pvar, const cs_real_t *restrict pvara, const cs_int_t bc_type[], const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t xcpp[], cs_real_t *restrict rhs) |
Add the explicit part of the convection/diffusion terms of a transport equation of a scalar field such as the temperature. More... | |
| void | cs_anisotropic_diffusion_scalar (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int inc, int iccocg, cs_real_t *restrict pvar, const cs_real_t *restrict pvara, const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t *restrict viscel, const cs_real_2_t weighf[], const cs_real_t weighb[], cs_real_t *restrict rhs) |
Add the explicit part of the diffusion terms with a symmetric tensor diffusivity for a transport equation of a scalar field . More... | |
| void | cs_anisotropic_diffusion_vector (int idtvar, int f_id, const cs_var_cal_opt_t var_cal_opt, int inc, int ifaccp, int ivisep, cs_real_3_t *restrict pvar, const cs_real_3_t *restrict pvara, const cs_int_t bc_type[], const cs_real_3_t coefav[], const cs_real_33_t coefbv[], const cs_real_3_t cofafv[], const cs_real_33_t cofbfv[], const cs_real_33_t i_visc[], const cs_real_t b_visc[], const cs_real_t secvif[], cs_real_3_t *restrict rhs) |
Add the explicit part of the diffusion terms with a symmetric tensorial diffusivity for a transport equation of a vector field . More... | |
| void | cs_face_diffusion_potential (const int f_id, const cs_mesh_t *m, cs_mesh_quantities_t *fvq, int init, int inc, int imrgra, int iccocg, int nswrgp, int imligp, int iphydp, int iwarnp, double epsrgp, double climgp, double extrap, cs_real_3_t *restrict frcxt, cs_real_t *restrict pvar, const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t viselx[], const cs_real_t visely[], const cs_real_t viselz[], cs_real_t *restrict i_massflux, cs_real_t *restrict b_massflux) |
| Update the face mass flux with the face pressure (or pressure increment, or pressure double increment) gradient. More... | |
| void | cs_face_anisotropic_diffusion_potential (const cs_mesh_t *m, cs_mesh_quantities_t *fvq, int init, int inc, int imrgra, int iccocg, int nswrgp, int imligp, int ircflp, int iphydp, int iwarnp, double epsrgp, double climgp, double extrap, cs_real_3_t *restrict frcxt, cs_real_t *restrict pvar, const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t *restrict viscel, const cs_real_2_t weighf[], const cs_real_t weighb[], cs_real_t *restrict i_massflux, cs_real_t *restrict b_massflux) |
Add the explicit part of the pressure gradient term to the mass flux in case of anisotropic diffusion of the pressure field . More... | |
| void | cs_diffusion_potential (const int f_id, const cs_mesh_t *m, cs_mesh_quantities_t *fvq, int init, int inc, int imrgra, int iccocg, int nswrgp, int imligp, int iphydp, int iwarnp, double epsrgp, double climgp, double extrap, cs_real_3_t *restrict frcxt, cs_real_t *restrict pvar, const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t viselx[], const cs_real_t visely[], const cs_real_t viselz[], cs_real_t *restrict diverg) |
| Update the cell mass flux divergence with the face pressure (or pressure increment, or pressure double increment) gradient. More... | |
| void | cs_anisotropic_diffusion_potential (const cs_mesh_t *m, cs_mesh_quantities_t *fvq, int init, int inc, int imrgra, int iccocg, int nswrgp, int imligp, int ircflp, int iphydp, int iwarnp, double epsrgp, double climgp, double extrap, cs_real_3_t *restrict frcxt, cs_real_t *restrict pvar, const cs_real_t coefap[], const cs_real_t coefbp[], const cs_real_t cofafp[], const cs_real_t cofbfp[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_6_t *restrict viscel, const cs_real_2_t weighf[], const cs_real_t weighb[], cs_real_t *restrict diverg) |
| Add the explicit part of the divergence of the mass flux due to the pressure gradient (routine analog to cs_anisotropic_diffusion_scalar). More... | |
| void bilsc2 | ( | const cs_int_t *const | idtvar, |
| const cs_int_t *const | f_id, | ||
| const cs_var_cal_opt_t *const | var_cal_opt, | ||
| const cs_int_t *const | icvflb, | ||
| const cs_int_t *const | inc, | ||
| const cs_int_t *const | iccocg, | ||
| const cs_int_t *const | ifaccp, | ||
| cs_real_t | pvar[], | ||
| const cs_real_t | pvara[], | ||
| const cs_int_t | bc_type[], | ||
| const cs_int_t | icvfli[], | ||
| const cs_real_t | coefap[], | ||
| const cs_real_t | coefbp[], | ||
| const cs_real_t | cofafp[], | ||
| const cs_real_t | cofbfp[], | ||
| const cs_real_t | i_massflux[], | ||
| const cs_real_t | b_massflux[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| cs_real_t | rhs[] | ||
| ) |
| void bilsc4 | ( | const cs_int_t *const | idtvar, |
| const cs_int_t *const | f_id, | ||
| const cs_var_cal_opt_t *const | var_cal_opt, | ||
| const cs_int_t *const | icvflb, | ||
| const cs_int_t *const | inc, | ||
| const cs_int_t *const | ifaccp, | ||
| const cs_int_t *const | ivisep, | ||
| cs_real_3_t | pvar[], | ||
| const cs_real_3_t | pvara[], | ||
| const cs_int_t | bc_type[], | ||
| const cs_int_t | icvfli[], | ||
| const cs_real_3_t | coefav[], | ||
| const cs_real_33_t | coefbv[], | ||
| const cs_real_3_t | cofafv[], | ||
| const cs_real_33_t | cofbfv[], | ||
| const cs_real_t | i_massflux[], | ||
| const cs_real_t | b_massflux[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| const cs_real_t | secvif[], | ||
| cs_real_3_t | rhs[] | ||
| ) |
| void bilsct | ( | const cs_int_t *const | idtvar, |
| const cs_int_t *const | f_id, | ||
| const cs_var_cal_opt_t *const | var_cal_opt, | ||
| const cs_int_t *const | inc, | ||
| const cs_int_t *const | iccocg, | ||
| const cs_int_t *const | ifaccp, | ||
| cs_real_t | pvar[], | ||
| const cs_real_t | pvara[], | ||
| const cs_int_t | bc_type[], | ||
| const cs_real_t | coefap[], | ||
| const cs_real_t | coefbp[], | ||
| const cs_real_t | cofafp[], | ||
| const cs_real_t | cofbfp[], | ||
| const cs_real_t | i_massflux[], | ||
| const cs_real_t | b_massflux[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| const cs_real_t | xcpp[], | ||
| cs_real_t | rhs[] | ||
| ) |
| void cs_anisotropic_diffusion_potential | ( | const cs_mesh_t * | m, |
| cs_mesh_quantities_t * | fvq, | ||
| int | init, | ||
| int | inc, | ||
| int | imrgra, | ||
| int | iccocg, | ||
| int | nswrgp, | ||
| int | imligp, | ||
| int | ircflp, | ||
| int | iphydp, | ||
| int | iwarnp, | ||
| double | epsrgp, | ||
| double | climgp, | ||
| double | extrap, | ||
| cs_real_3_t *restrict | frcxt, | ||
| cs_real_t *restrict | pvar, | ||
| const cs_real_t | coefap[], | ||
| const cs_real_t | coefbp[], | ||
| const cs_real_t | cofafp[], | ||
| const cs_real_t | cofbfp[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| cs_real_6_t *restrict | viscel, | ||
| const cs_real_2_t | weighf[], | ||
| const cs_real_t | weighb[], | ||
| cs_real_t *restrict | diverg | ||
| ) |
Add the explicit part of the divergence of the mass flux due to the pressure gradient (routine analog to cs_anisotropic_diffusion_scalar).
More precisely, the divergence of the mass flux side
is updated as follows:
| [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] | iwarnp | verbosity |
| [in] | epsrgp | relative precision for the gradient reconstruction |
| [in] | climgp | clipping coeffecient for the computation of the gradient |
| [in] | extrap | coefficient for extrapolation 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 | at interior faces for the r.h.s. |
| [in] | b_visc | at border faces for the r.h.s. |
| [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] | 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] | iwarnp | verbosity |
| [in] | epsrgp | relative precision for the gradient reconstruction |
| [in] | climgp | clipping coeffecient for the computation of the gradient |
| [in] | extrap | coefficient for extrapolation 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 | at interior faces for the r.h.s. |
| [in] | b_visc | at border faces for the r.h.s. |
| [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:
has already been initialized before calling cs_anisotropic_diffusion_scalar!| [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 | at interior faces for the r.h.s. |
| [in] | b_visc | at border faces for the r.h.s. |
| [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_vector | ( | int | idtvar, |
| int | f_id, | ||
| const cs_var_cal_opt_t | var_cal_opt, | ||
| int | inc, | ||
| int | ifaccp, | ||
| int | ivisep, | ||
| cs_real_3_t *restrict | pvar, | ||
| const cs_real_3_t *restrict | pvara, | ||
| const cs_int_t | bc_type[], | ||
| const cs_real_3_t | coefav[], | ||
| const cs_real_33_t | coefbv[], | ||
| const cs_real_3_t | cofafv[], | ||
| const cs_real_33_t | cofbfv[], | ||
| const cs_real_33_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| const cs_real_t | secvif[], | ||
| cs_real_3_t *restrict | rhs | ||
| ) |
Add the explicit part of the diffusion terms with a symmetric tensorial diffusivity for a transport equation of a vector field
.
More precisely, the right hand side
is updated as follows:
Warning:
has already been initialized before calling diftnv!| [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] | ifaccp | indicator
|
| [in] | ivisep | indicator to take ![]()
|
| [in] | pvar | solved variable (current time step) |
| [in] | pvara | solved variable (previous time step) |
| [in] | bc_type | boundary condition type |
| [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 | at interior faces for the r.h.s. |
| [in] | b_visc | at border faces for the r.h.s. |
| [in] | secvif | secondary viscosity at interior faces |
| [in,out] | rhs | right hand side |
|
inlinestatic |
Handle preparation of boundary face values for the flux computation in case of a steady algorithm.
| [in] | ircflp | recontruction flag |
| [in] | relaxp | relaxation coefficient |
| [in] | diipb | distance I'I' |
| [in] | gradi | gradient at cell i |
| [in] | pi | value at cell i |
| [in] | pia | old value at cell i |
| [out] | pir | relaxed value at cell i |
| [out] | pipr | relaxed reconstructed value at cell i |
|
inlinestatic |
Handle preparation of boundary face values for the flux computation in case of an unsteady algorithm.
| [in] | ircflp | recontruction flag |
| [in] | diipb | distance I'I' |
| [in] | gradi | gradient at cell i |
| [in] | pi | value at cell i |
| [in] | pir | value at cell i (for consistancy) |
| [out] | pipr | reconstructed value at cell i |
|
inlinestatic |
Reconstruct values in I' at boundary cell i.
| [in] | diipb | distance I'I' |
| [in] | gradi | gradient at cell i |
| [in] | ircflp | recontruction flag |
| [out] | recoi | reconstruction at cell i |
|
inlinestatic |
Add diffusive flux to flux at boundary face.
| [in] | idiffp | diffusion flag |
| [in] | inc | Not an increment flag |
| [in] | pipr | relaxed reconstructed value at cell i |
| [in] | cofafp | explicit boundary coefficient for diffusion operator |
| [in] | cofbfp | implicit boundary coefficient for diffusion operator |
| [in] | b_visc | mass flux at boundary face |
| [in,out] | flux | flux at boundary face |
|
inlinestatic |
Add convective flux (substracting the mass accumulation from it) to flux at boundary face. The convective flux can be either an upwind flux or an imposed value.
| [in] | iconvp | convection flag |
| [in] | inc | Not an increment flag |
| [in] | ifaccp | indicator
|
| [in] | bc_type | type of boundary face |
| [in] | icvfli | imposed convective flux flag |
| [in] | pi | value at cell i |
| [in] | pir | relaxed value at cell i |
| [in] | pipr | relaxed reconstructed value at cell i |
| [in] | coefap | explicit boundary coefficient for convection operator |
| [in] | coefbp | implicit boundary coefficient for convection operator |
| [in] | coface | explicit imposed convective flux value (0 otherwise). |
| [in] | cofbce | implicit part of imp. conv. flux value |
| [in] | b_massflux | mass flux at boundary face |
| [in,out] | flux | flux at boundary face |
|
inlinestatic |
Add convective flux (conservative formulation) to flux at boundary face. The convective flux can be either an upwind flux or an imposed value.
| [in] | iconvp | convection flag |
| [in] | inc | Not an increment flag |
| [in] | ifaccp | indicator
|
| [in] | bc_type | type of boundary face |
| [in] | icvfli | imposed convective flux flag |
| [in] | pir | relaxed value at cell i |
| [in] | pipr | relaxed reconstructed value at cell i |
| [in] | coefap | explicit boundary coefficient for convection operator |
| [in] | coefbp | implicit boundary coefficient for convection operator |
| [in] | coface | explicit imposed convective flux value (0 otherwise). |
| [in] | cofbce | implicit part of imp. conv. flux value |
| [in] | xcpp | specific heat value if the scalar is the temperature, 1 otherwise |
| [in] | b_massflux | mass flux at boundary face |
| [in,out] | flux | flux at boundary face |
|
inlinestatic |
Copy values at bounadry cell i for consistency with relaxation case.
| [in] | pi | value at cell i |
| [in] | recoi | reconstruction at cell i |
| [out] | pir | value at cell i |
| [out] | pipr | reconstructed value at cell i |
|
inlinestatic |
Compute relaxed values at boundary cell i.
| [in] | relaxp | relaxation coefficient |
| [in] | pi | value at cell i |
| [in] | pia | old value at cell i |
| [in] | recoi | reconstruction at cell i |
| [out] | pir | relaxed value at cell i |
| [out] | pipr | relaxed reconstructed value at cell i |
|
inlinestatic |
Add convective flux (substracting the mass accumulation from it) to flux at boundary face. The convective flux is a pure upwind flux.
| [in] | iconvp | convection flag |
| [in] | inc | Not an increment flag |
| [in] | ifaccp | indicator
|
| [in] | bc_type | type of boundary face |
| [in] | pi | value at cell i |
| [in] | pir | relaxed value at cell i |
| [in] | pipr | relaxed reconstructed value at cell i |
| [in] | coefap | explicit boundary coefficient for convection operator |
| [in] | coefbp | implicit boundary coefficient for convection operator |
| [in] | b_massflux | mass flux at boundary face |
| [in,out] | flux | flux at boundary face |
|
inlinestatic |
Add convective flux (conservative formulation) to flux at boundary face. The convective flux is a pure upwind flux.
| [in] | iconvp | convection flag |
| [in] | inc | Not an increment flag |
| [in] | ifaccp | indicator
|
| [in] | bc_type | type of boundary face |
| [in] | pir | relaxed value at cell i |
| [in] | pipr | relaxed reconstructed value at cell i |
| [in] | coefap | explicit boundary coefficient for convection operator |
| [in] | coefbp | implicit boundary coefficient for convection operator |
| [in] | b_massflux | mass flux at boundary face |
| [in] | xcpp | specific heat value if the scalar is the temperature, 1 otherwise |
| [in,out] | flux | flux at boundary face |
|
inlinestatic |
Blend face values for a centered or SOLU scheme with face values for an upwind scheme.
| [in] | blencp | proportion of centered or SOLU scheme, (1-blencp) is the proportion of upwind. |
| [in] | pi | value at cell i |
| [in] | pj | value at cell j |
| [in] | pir | relaxed value at cell i |
| [in] | pjr | relaxed value at cell j |
| [out] | pifri | contribution of i to flux from i to j |
| [out] | pifrj | contribution of i to flux from j to i |
| [out] | pjfri | contribution of j to flux from i to j |
| [out] | pjfrj | contribution of j to flux from j to i |
|
inlinestatic |
Prepare value at face ij by using a centered scheme.
| [in] | pnd | weight |
| [in] | pip | reconstructed value at cell i |
| [in] | pjp | reconstructed value at cell j |
| [in] | pipr | relaxed reconstructed value at cell i |
| [in] | pjpr | relaxed reconstructed value at cell j |
| [out] | pifri | contribution of i to flux from i to j |
| [out] | pifrj | contribution of i to flux from j to i |
| [out] | pjfri | contribution of j to flux from i to j |
| [out] | pjfrj | contribution of j to flux from j to i |
| void cs_convection_diffusion_scalar | ( | int | idtvar, |
| int | f_id, | ||
| const cs_var_cal_opt_t | var_cal_opt, | ||
| int | icvflb, | ||
| int | inc, | ||
| int | iccocg, | ||
| int | ifaccp, | ||
| cs_real_t *restrict | pvar, | ||
| const cs_real_t *restrict | pvara, | ||
| const cs_int_t | bc_type[], | ||
| const cs_int_t | icvfli[], | ||
| const cs_real_t | coefap[], | ||
| const cs_real_t | coefbp[], | ||
| const cs_real_t | cofafp[], | ||
| const cs_real_t | cofbfp[], | ||
| const cs_real_t | i_massflux[], | ||
| const cs_real_t | b_massflux[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| cs_real_t *restrict | rhs | ||
| ) |
Add the explicit part of the convection/diffusion terms of a standard transport equation of a scalar field
.
More precisely, the right hand side
is updated as follows:
Warning:
has already been initialized before calling bilsc2!| [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] | ifaccp | indicator
|
| [in] | pvar | solved variable (current time step) |
| [in] | pvara | solved variable (previous time step) |
| [in] | bc_type | boundary condition type |
| [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 | at interior faces for the r.h.s. |
| [in] | b_visc | at border faces for the r.h.s. |
| [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 | ifaccp, | ||
| cs_real_t *restrict | pvar, | ||
| const cs_real_t *restrict | pvara, | ||
| const cs_int_t | bc_type[], | ||
| const cs_real_t | coefap[], | ||
| const cs_real_t | coefbp[], | ||
| const cs_real_t | cofafp[], | ||
| const cs_real_t | cofbfp[], | ||
| const cs_real_t | i_massflux[], | ||
| const cs_real_t | b_massflux[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| const cs_real_t | xcpp[], | ||
| cs_real_t *restrict | rhs | ||
| ) |
Add the explicit part of the convection/diffusion terms of a transport equation of a scalar field
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] | ifaccp | indicator
|
| [in] | pvar | solved variable (current time step) |
| [in] | pvara | solved variable (previous time step) |
| [in] | bc_type | boundary condition type |
| [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 | at interior faces for the r.h.s. |
| [in] | b_visc | at border faces for the r.h.s. |
| [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 | ifaccp, | ||
| int | ivisep, | ||
| cs_real_3_t *restrict | pvar, | ||
| const cs_real_3_t *restrict | pvara, | ||
| const cs_int_t | bc_type[], | ||
| const cs_int_t | icvfli[], | ||
| const cs_real_3_t | coefav[], | ||
| const cs_real_33_t | coefbv[], | ||
| const cs_real_3_t | cofafv[], | ||
| const cs_real_33_t | cofbfv[], | ||
| const cs_real_t | i_massflux[], | ||
| const cs_real_t | b_massflux[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| const cs_real_t | secvif[], | ||
| cs_real_3_t *restrict | rhs | ||
| ) |
Add the explicit part of the convection/diffusion terms of a transport equation of a vector field
.
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:
has already been initialized before calling bilsc!| [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] | ifaccp | indicator
|
| [in] | ivisep | indicator to take ![]()
|
| [in] | pvar | solved velocity (current time step) |
| [in] | pvara | solved velocity (previous time step) |
| [in] | bc_type | boundary condition type |
| [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 | at interior faces for the r.h.s. |
| [in] | b_visc | at border faces for the r.h.s. |
| [in] | secvif | secondary viscosity at interior 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 | iwarnp, | ||
| double | epsrgp, | ||
| double | climgp, | ||
| double | extrap, | ||
| cs_real_3_t *restrict | frcxt, | ||
| cs_real_t *restrict | pvar, | ||
| const cs_real_t | coefap[], | ||
| const cs_real_t | coefbp[], | ||
| const cs_real_t | cofafp[], | ||
| const cs_real_t | cofbfp[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| const cs_real_t | viselx[], | ||
| const cs_real_t | visely[], | ||
| const cs_real_t | viselz[], | ||
| cs_real_t *restrict | diverg | ||
| ) |
Update the cell mass flux divergence with the face pressure (or pressure increment, or pressure double increment) gradient.
| [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] | epsrgp | relative precision for the gradient reconstruction |
| [in] | climgp | clipping coeffecient for the computation of the gradient |
| [in] | extrap | coefficient for extrapolation 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 | at interior faces for the r.h.s. |
| [in] | b_visc | at border faces for the r.h.s. |
| [in] | viselx | viscosity by cell, dir x |
| [in] | visely | viscosity by cell, dir y |
| [in] | viselz | viscosity by cell, dir z |
| [in,out] | diverg | mass flux divergence |
| void cs_face_anisotropic_diffusion_potential | ( | const cs_mesh_t * | m, |
| cs_mesh_quantities_t * | fvq, | ||
| int | init, | ||
| int | inc, | ||
| int | imrgra, | ||
| int | iccocg, | ||
| int | nswrgp, | ||
| int | imligp, | ||
| int | ircflp, | ||
| int | iphydp, | ||
| int | iwarnp, | ||
| double | epsrgp, | ||
| double | climgp, | ||
| double | extrap, | ||
| cs_real_3_t *restrict | frcxt, | ||
| cs_real_t *restrict | pvar, | ||
| const cs_real_t | coefap[], | ||
| const cs_real_t | coefbp[], | ||
| const cs_real_t | cofafp[], | ||
| const cs_real_t | cofbfp[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| cs_real_6_t *restrict | viscel, | ||
| const cs_real_2_t | weighf[], | ||
| const cs_real_t | weighb[], | ||
| cs_real_t *restrict | i_massflux, | ||
| cs_real_t *restrict | b_massflux | ||
| ) |
Add the explicit part of the pressure gradient term to the mass flux in case of anisotropic diffusion of the pressure field
.
More precisely, the mass flux side
is updated as follows:
| [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] | iwarnp | verbosity |
| [in] | epsrgp | relative precision for the gradient reconstruction |
| [in] | climgp | clipping coeffecient for the computation of the gradient |
| [in] | extrap | coefficient for extrapolation 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 | at interior faces for the r.h.s. |
| [in] | b_visc | at border faces for the r.h.s. |
| [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_diffusion_potential | ( | const int | f_id, |
| const cs_mesh_t * | m, | ||
| cs_mesh_quantities_t * | fvq, | ||
| int | init, | ||
| int | inc, | ||
| int | imrgra, | ||
| int | iccocg, | ||
| int | nswrgp, | ||
| int | imligp, | ||
| int | iphydp, | ||
| int | iwarnp, | ||
| double | epsrgp, | ||
| double | climgp, | ||
| double | extrap, | ||
| cs_real_3_t *restrict | frcxt, | ||
| cs_real_t *restrict | pvar, | ||
| const cs_real_t | coefap[], | ||
| const cs_real_t | coefbp[], | ||
| const cs_real_t | cofafp[], | ||
| const cs_real_t | cofbfp[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| const cs_real_t | viselx[], | ||
| const cs_real_t | visely[], | ||
| const cs_real_t | viselz[], | ||
| cs_real_t *restrict | i_massflux, | ||
| cs_real_t *restrict | b_massflux | ||
| ) |
Update the face mass flux with the face pressure (or pressure increment, or pressure double increment) gradient.
| [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] | epsrgp | relative precision for the gradient reconstruction |
| [in] | climgp | clipping coeffecient for the computation of the gradient |
| [in] | extrap | coefficient for extrapolation 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 | at interior faces for the r.h.s. |
| [in] | b_visc | at border faces for the r.h.s. |
| [in] | viselx | viscosity by cell, dir x |
| [in] | visely | viscosity by cell, dir y |
| [in] | viselz | viscosity by cell, dir z |
| [in,out] | i_massflux | mass flux at interior faces |
| [in,out] | b_massflux | mass flux at boundary faces |
|
inlinestatic |
Handle preparation of internal face values for the fluxes computation in case of a steady algorithm and without enabling slope tests.
| [in] | ircflp | recontruction flag |
| [in] | ischcp | second order convection scheme flag |
| [in] | relaxp | relaxation coefficient |
| [in] | blencp | proportion of centered or SOLU scheme, (1-blencp) is the proportion of upwind. |
| [in] | weight | geometrical weight |
| [in] | cell_ceni | center of gravity coordinates of cell i |
| [in] | cell_cenj | center of gravity coordinates of cell i |
| [in] | i_face_cog | center of gravity coordinates of face ij |
| [in] | dijpf | distance I'J' |
| [in] | gradi | gradient at cell i |
| [in] | gradj | gradient at cell j |
| [in] | pi | value at cell i |
| [in] | pj | value at cell j |
| [in] | pia | old value at cell i |
| [in] | pja | old value at cell j |
| [out] | pifri | contribution of i to flux from i to j |
| [out] | pifrj | contribution of i to flux from j to i |
| [out] | pjfri | contribution of j to flux from i to j |
| [out] | pjfrj | contribution of j to flux from j to i |
| [out] | pip | reconstructed value at cell i |
| [out] | pjp | reconstructed value at cell j |
| [out] | pipr | relaxed reconstructed value at cell i |
| [out] | pjpr | relaxed reconstructed value at cell j |
|
inlinestatic |
Handle preparation of internal face values for the fluxes computation in case of a steady algorithm and using slope tests.
| [out] | upwind_switch | slope test result |
| [in] | ircflp | recontruction flag |
| [in] | ischcp | second order convection scheme flag |
| [in] | relaxp | relaxation coefficient |
| [in] | blencp | proportion of centered or SOLU scheme, (1-blencp) is the proportion of upwind. |
| [in] | weight | geometrical weight |
| [in] | i_dist | distance IJ.Nij |
| [in] | i_face_surf | face surface |
| [in] | cell_ceni | center of gravity coordinates of cell i |
| [in] | cell_cenj | center of gravity coordinates of cell i |
| [in] | i_face_normal | face normal |
| [in] | i_face_cog | center of gravity coordinates of face ij |
| [in] | dijpf | distance I'J' |
| [in] | i_massflux | mass flux at face ij |
| [in] | gradi | gradient at cell i |
| [in] | gradj | gradient at cell j |
| [in] | grdpai | upwind gradient at cell i |
| [in] | grdpaj | upwind gradient at cell j |
| [in] | pi | value at cell i |
| [in] | pj | value at cell j |
| [in] | pia | old value at cell i |
| [in] | pja | old value at cell j |
| [out] | pifri | contribution of i to flux from i to j |
| [out] | pifrj | contribution of i to flux from j to i |
| [out] | pjfri | contribution of j to flux from i to j |
| [out] | pjfrj | contribution of j to flux from j to i |
| [out] | pip | reconstructed value at cell i |
| [out] | pjp | reconstructed value at cell j |
| [out] | pipr | relaxed reconstructed value at cell i |
| [out] | pjpr | relaxed reconstructed value at cell j |
|
inlinestatic |
Handle preparation of internal face values for the fluxes computation in case of a steady algorithm and a pure upwind flux.
| [in] | ircflp | recontruction flag |
| [in] | relaxp | relaxation coefficient |
| [in] | weight | geometrical weight |
| [in] | cell_ceni | center of gravity coordinates of cell i |
| [in] | cell_cenj | center of gravity coordinates of cell i |
| [in] | i_face_cog | center of gravity coordinates of face ij |
| [in] | dijpf | distance I'J' |
| [in] | gradi | gradient at cell i |
| [in] | gradj | gradient at cell j |
| [in] | pi | value at cell i |
| [in] | pj | value at cell j |
| [in] | pia | old value at cell i |
| [in] | pja | old value at cell j |
| [out] | pifri | contribution of i to flux from i to j |
| [out] | pifrj | contribution of i to flux from j to i |
| [out] | pjfri | contribution of j to flux from i to j |
| [out] | pjfrj | contribution of j to flux from j to i |
| [out] | pip | reconstructed value at cell i |
| [out] | pjp | reconstructed value at cell j |
| [out] | pipr | relaxed reconstructed value at cell i |
| [out] | pjpr | relaxed reconstructed value at cell j |
|
inlinestatic |
Handle preparation of internal face values for the fluxes computation in case of a unsteady algorithm and without enabling slope tests.
| [in] | ircflp | recontruction flag |
| [in] | ischcp | second order convection scheme flag |
| [in] | blencp | proportion of centered or SOLU scheme, (1-blencp) is the proportion of upwind. |
| [in] | weight | geometrical weight |
| [in] | cell_ceni | center of gravity coordinates of cell i |
| [in] | cell_cenj | center of gravity coordinates of cell i |
| [in] | i_face_cog | center of gravity coordinates of face ij |
| [in] | dijpf | distance I'J' |
| [in] | gradi | gradient at cell i |
| [in] | gradj | gradient at cell j |
| [in] | pi | value at cell i |
| [in] | pj | value at cell j |
| [out] | pifri | contribution of i to flux from i to j |
| [out] | pifrj | contribution of i to flux from j to i |
| [out] | pjfri | contribution of j to flux from i to j |
| [out] | pjfrj | contribution of j to flux from j to i |
| [out] | pip | reconstructed value at cell i |
| [out] | pjp | reconstructed value at cell j |
| [out] | pipr | relaxed reconstructed value at cell i |
| [out] | pjpr | relaxed reconstructed value at cell j |
|
inlinestatic |
Handle preparation of internal face values for the fluxes computation in case of a unsteady algorithm and using slope tests.
| [out] | upwind_switch | slope test result |
| [in] | ircflp | recontruction flag |
| [in] | ischcp | second order convection scheme flag |
| [in] | blencp | proportion of centered or SOLU scheme, (1-blencp) is the proportion of upwind. |
| [in] | weight | geometrical weight |
| [in] | i_dist | distance IJ.Nij |
| [in] | i_face_surf | face surface |
| [in] | cell_ceni | center of gravity coordinates of cell i |
| [in] | cell_cenj | center of gravity coordinates of cell i |
| [in] | i_face_normal | face normal |
| [in] | i_face_cog | center of gravity coordinates of face ij |
| [in] | dijpf | distance I'J' |
| [in] | i_massflux | mass flux at face ij |
| [in] | gradi | gradient at cell i |
| [in] | gradj | gradient at cell j |
| [in] | grdpai | upwind gradient at cell i |
| [in] | grdpaj | upwind gradient at cell j |
| [in] | pi | value at cell i |
| [in] | pj | value at cell j |
| [out] | pifri | contribution of i to flux from i to j |
| [out] | pifrj | contribution of i to flux from j to i |
| [out] | pjfri | contribution of j to flux from i to j |
| [out] | pjfrj | contribution of j to flux from j to i |
| [out] | pip | reconstructed value at cell i |
| [out] | pjp | reconstructed value at cell j |
| [out] | pipr | relaxed reconstructed value at cell i |
| [out] | pjpr | relaxed reconstructed value at cell j |
|
inlinestatic |
Handle preparation of internal face values for the fluxes computation in case of an unsteady algorithm and a pure upwind flux.
| [in] | ircflp | recontruction flag |
| [in] | weight | geometrical weight |
| [in] | cell_ceni | center of gravity coordinates of cell i |
| [in] | cell_cenj | center of gravity coordinates of cell i |
| [in] | i_face_cog | center of gravity coordinates of face ij |
| [in] | dijpf | distance I'J' |
| [in] | gradi | gradient at cell i |
| [in] | gradj | gradient at cell j |
| [in] | pi | value at cell i |
| [in] | pj | value at cell j |
| [out] | pifri | contribution of i to flux from i to j |
| [out] | pifrj | contribution of i to flux from j to i |
| [out] | pjfri | contribution of j to flux from i to j |
| [out] | pjfrj | contribution of j to flux from j to i |
| [out] | pip | reconstructed value at cell i |
| [out] | pjp | reconstructed value at cell j |
| [out] | pipr | relaxed reconstructed value at cell i |
| [out] | pjpr | relaxed reconstructed value at cell j |
|
inlinestatic |
Reconstruct values in I' and J'.
| [in] | ircflp | recontruction flag |
| [in] | pnd | geometrical weight |
| [in] | cell_ceni | center of gravity coordinates of cell i |
| [in] | cell_cenj | center of gravity coordinates of cell i |
| [in] | i_face_cog | center of gravity coordinates of face ij |
| [in] | dijpf | distance I'J' |
| [in] | gradi | gradient at cell i |
| [in] | gradj | gradient at cell j |
| [in] | pi | value at cell i |
| [in] | pj | value at cell j |
| [out] | recoi | reconstruction at cell i |
| [out] | recoj | reconstruction at cell j |
| [out] | pip | reconstructed value at cell i |
| [out] | pjp | reconstructed value at cell j |
|
inlinestatic |
Add convective fluxes (substracting the mass accumulation from them) to fluxes at face ij.
| [in] | iconvp | convection flag |
| [in] | pi | value at cell i |
| [in] | pj | value at cell j |
| [in] | pifri | contribution of i to flux from i to j |
| [in] | pifrj | contribution of i to flux from j to i |
| [in] | pjfri | contribution of j to flux from i to j |
| [in] | pjfrj | contribution of j to flux from j to i |
| [in] | i_massflux | mass flux at face ij |
| [in,out] | fluxij | fluxes at face ij |
|
inlinestatic |
Add convective fluxes (cons. formulation) to fluxes at face ij.
| [in] | iconvp | convection flag |
| [in] | pifri | contribution of i to flux from i to j |
| [in] | pifrj | contribution of i to flux from j to i |
| [in] | pjfri | contribution of j to flux from i to j |
| [in] | pjfrj | contribution of j to flux from j to i |
| [in] | i_massflux | mass flux at face ij |
| [in] | xcppi | specific heat value if the scalar is the temperature, 1 otherwise at cell i |
| [in] | xcppj | specific heat value if the scalar is the temperature, 1 otherwise at cell j |
| [in,out] | fluxij | fluxes at face ij |
|
inlinestatic |
Add diffusive fluxes to fluxes at face ij.
| [in] | idiffp | diffusion flag |
| [in] | pip | reconstructed value at cell i |
| [in] | pjp | reconstructed value at cell j |
| [in] | pipr | relaxed reconstructed value at cell i |
| [in] | pjpr | relaxed reconstructed value at cell j |
| [in] | i_visc | diffusion coefficient (divided by IJ) at face ij |
| [in,out] | fluxij | fluxes at face ij |
|
inlinestatic |
Copy reconstructed or not cell values for consistency with relaxation case.
| [in] | pi | value at cell i |
| [in] | pj | value at cell j |
| [in] | pip | reconstructed value at cell i |
| [in] | pjp | reconstructed value at cell j |
| [out] | pir | value at cell i |
| [out] | pjr | value at cell j |
| [out] | pipr | reconstructed value at cell i |
| [out] | pjpr | reconstructed value at cell j |
|
inlinestatic |
Compute relaxed values at cell i and j.
| [in] | relaxp | relaxation coefficient |
| [in] | pia | old value at cell i |
| [in] | pja | old value at cell j |
| [in] | recoi | reconstruction at cell i |
| [in] | recoj | reconstruction at cell j |
| [in] | pi | value at cell i |
| [in] | pj | value at cell j |
| [out] | pir | relaxed value at cell i |
| [out] | pjr | relaxed value at cell j |
| [out] | pipr | relaxed reconstructed value at cell i |
| [out] | pjpr | relaxed reconstructed value at cell j |
|
inlinestatic |
Compute slope test criteria at internal face between cell i and j.
| [in] | pi | value at cell i |
| [in] | pj | value at cell j |
| [in] | distf | distance IJ.Nij |
| [in] | srfan | face surface |
| [in] | i_face_normal | face normal |
| [in] | gradi | gradient at cell i |
| [in] | gradj | gradient at cell j |
| [in] | grdpai | upwind gradient at cell i |
| [in] | grdpaj | upwind gradient at cell j |
| [in] | i_massflux | mass flux at face (from i to j) |
| [out] | testij | value of slope test first criterium |
| [out] | tesqck | value of slope test second criterium |
|
inlinestatic |
Compute the upwind gradient used in the slope tests.
| [in] | f_id | field index |
| [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 |
|
inlinestatic |
Prepare value at face ij by using a Second Order Linear Upwind scheme.
| [in] | cell_ceni | center of gravity coordinates of cell i |
| [in] | cell_cenj | center of gravity coordinates of cell i |
| [in] | i_face_cog | center of gravity coordinates of face ij |
| [in] | gradi | gradient at cell i |
| [in] | gradj | gradient at cell j |
| [in] | pi | value at cell i |
| [in] | pj | value at cell j |
| [in] | pir | relaxed value at cell i |
| [in] | pjr | relaxed value at cell j |
| [out] | pifri | contribution of i to flux from i to j |
| [out] | pifrj | contribution of i to flux from j to i |
| [out] | pjfri | contribution of j to flux from i to j |
| [out] | pjfrj | contribution of j to flux from j to i |
|
inlinestatic |
Prepare value at face ij by using an upwind scheme.
| [in] | pi | value at cell i |
| [in] | pj | value at cell j |
| [in] | pir | relaxed value at cell i |
| [in] | pjr | relaxed value at cell j |
| [out] | pifri | contribution of i to flux from i to j |
| [out] | pifrj | contribution of i to flux from j to i |
| [out] | pjfri | contribution of j to flux from i to j |
| [out] | pjfrj | contribution of j to flux from j to i |
| void diften | ( | const cs_int_t *const | idtvar, |
| const cs_int_t *const | f_id, | ||
| const cs_var_cal_opt_t *const | var_cal_opt, | ||
| const cs_int_t *const | inc, | ||
| const cs_int_t *const | iccocg, | ||
| cs_real_t | pvar[], | ||
| const cs_real_t | pvara[], | ||
| const cs_real_t | coefap[], | ||
| const cs_real_t | coefbp[], | ||
| const cs_real_t | cofafp[], | ||
| const cs_real_t | cofbfp[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| cs_real_6_t | viscel[], | ||
| const cs_real_2_t | weighf[], | ||
| const cs_real_t | weighb[], | ||
| cs_real_t | rhs[] | ||
| ) |
| void diftnv | ( | const cs_int_t *const | idtvar, |
| const cs_int_t *const | f_id, | ||
| const cs_var_cal_opt_t *const | var_cal_opt, | ||
| const cs_int_t *const | inc, | ||
| const cs_int_t *const | ifaccp, | ||
| const cs_int_t *const | ivisep, | ||
| cs_real_3_t | pvar[], | ||
| const cs_real_3_t | pvara[], | ||
| const cs_int_t | bc_type[], | ||
| const cs_real_3_t | coefav[], | ||
| const cs_real_33_t | coefbv[], | ||
| const cs_real_3_t | cofafv[], | ||
| const cs_real_33_t | cofbfv[], | ||
| const cs_real_33_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| const cs_real_t | secvif[], | ||
| cs_real_3_t | rhs[] | ||
| ) |
| void itrgrp | ( | const cs_int_t *const | f_id, |
| const cs_int_t *const | init, | ||
| const cs_int_t *const | inc, | ||
| const cs_int_t *const | imrgra, | ||
| const cs_int_t *const | iccocg, | ||
| const cs_int_t *const | nswrgp, | ||
| const cs_int_t *const | imligp, | ||
| const cs_int_t *const | iphydp, | ||
| const cs_int_t *const | iwarnp, | ||
| const cs_real_t *const | epsrgp, | ||
| const cs_real_t *const | climgp, | ||
| const cs_real_t *const | extrap, | ||
| cs_real_3_t | frcxt[], | ||
| cs_real_t | pvar[], | ||
| const cs_real_t | coefap[], | ||
| const cs_real_t | coefbp[], | ||
| const cs_real_t | cofafp[], | ||
| const cs_real_t | cofbfp[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| const cs_real_t | viselx[], | ||
| const cs_real_t | visely[], | ||
| const cs_real_t | viselz[], | ||
| cs_real_t | diverg[] | ||
| ) |
| void itrgrv | ( | const cs_int_t *const | init, |
| const cs_int_t *const | inc, | ||
| const cs_int_t *const | imrgra, | ||
| const cs_int_t *const | iccocg, | ||
| const cs_int_t *const | nswrgp, | ||
| const cs_int_t *const | imligp, | ||
| const cs_int_t *const | ircflp, | ||
| const cs_int_t *const | iphydp, | ||
| const cs_int_t *const | iwarnp, | ||
| const cs_real_t *const | epsrgp, | ||
| const cs_real_t *const | climgp, | ||
| const cs_real_t *const | extrap, | ||
| cs_real_3_t | frcxt[], | ||
| cs_real_t | pvar[], | ||
| const cs_real_t | coefap[], | ||
| const cs_real_t | coefbp[], | ||
| const cs_real_t | cofafp[], | ||
| const cs_real_t | cofbfp[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| cs_real_6_t | viscel[], | ||
| const cs_real_2_t | weighf[], | ||
| const cs_real_t | weighb[], | ||
| cs_real_t | diverg[] | ||
| ) |
| void itrmas | ( | const cs_int_t *const | f_id, |
| const cs_int_t *const | init, | ||
| const cs_int_t *const | inc, | ||
| const cs_int_t *const | imrgra, | ||
| const cs_int_t *const | iccocg, | ||
| const cs_int_t *const | nswrgp, | ||
| const cs_int_t *const | imligp, | ||
| const cs_int_t *const | iphydp, | ||
| const cs_int_t *const | iwarnp, | ||
| const cs_real_t *const | epsrgp, | ||
| const cs_real_t *const | climgp, | ||
| const cs_real_t *const | extrap, | ||
| cs_real_3_t | frcxt[], | ||
| cs_real_t | pvar[], | ||
| const cs_real_t | coefap[], | ||
| const cs_real_t | coefbp[], | ||
| const cs_real_t | cofafp[], | ||
| const cs_real_t | cofbfp[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| const cs_real_t | viselx[], | ||
| const cs_real_t | visely[], | ||
| const cs_real_t | viselz[], | ||
| cs_real_t | i_massflux[], | ||
| cs_real_t | b_massflux[] | ||
| ) |
| void itrmav | ( | const cs_int_t *const | init, |
| const cs_int_t *const | inc, | ||
| const cs_int_t *const | imrgra, | ||
| const cs_int_t *const | iccocg, | ||
| const cs_int_t *const | nswrgp, | ||
| const cs_int_t *const | imligp, | ||
| const cs_int_t *const | ircflp, | ||
| const cs_int_t *const | iphydp, | ||
| const cs_int_t *const | iwarnp, | ||
| const cs_real_t *const | epsrgp, | ||
| const cs_real_t *const | climgp, | ||
| const cs_real_t *const | extrap, | ||
| cs_real_3_t | frcxt[], | ||
| cs_real_t | pvar[], | ||
| const cs_real_t | coefap[], | ||
| const cs_real_t | coefbp[], | ||
| const cs_real_t | cofafp[], | ||
| const cs_real_t | cofbfp[], | ||
| const cs_real_t | i_visc[], | ||
| const cs_real_t | b_visc[], | ||
| cs_real_6_t | viscel[], | ||
| const cs_real_2_t | weighf[], | ||
| const cs_real_t | weighb[], | ||
| cs_real_t | i_massflux[], | ||
| cs_real_t | b_massflux[] | ||
| ) |
1.8.7