This file is the header file of several Riemann solvers and GRP solvers. 更多...
#include "../include/var_struc.h"
宏定义 | |
#define | Riemann_solver_exact_single Riemann_solver_exact_Ben |
Which solver is chosen as the exact Riemann solver for single-component flow. 更多... | |
函数 | |
double | Riemann_solver_exact (double *U_star, double *P_star, const double gammaL, const double gammaR, const double u_L, const double u_R, const double p_L, const double p_R, const double c_L, const double c_R, _Bool *CRW, const double eps, const double tol, int N) |
EXACT RIEMANN SOLVER FOR Two-Component γ-Law Gas 更多... | |
double | Riemann_solver_starPU (double *U_star, double *P_star, const double GammaL, const double GammaR, const double UL, const double UR, const double PL, const double PR, const double CL, const double CR, _Bool *CRW, const double eps, const double TOLPRE, const int NRITER) |
double | Riemann_solver_exact_Ben (double *U_star, double *P_star, const double gamma, const double u_L, const double u_R, const double p_L, const double p_R, const double c_L, const double c_R, _Bool *CRW, const double eps, const double tol, const int N) |
EXACT RIEMANN SOLVER FOR A γ-Law Gas 更多... | |
double | Riemann_solver_exact_Toro (double *U_star, double *P_star, const double gamma, const double U_l, const double U_r, const double P_l, const double P_r, const double c_l, const double c_r, _Bool *CRW, const double eps, const double tol, const int N) |
EXACT RIEMANN SOLVER FOR THE EULER EQUATIONS 更多... | |
void | linear_GRP_solver_LAG (double *D, double *U, const struct i_f_var *ifv_L, const struct i_f_var *ifv_R, const double eps, const double atc) |
void | linear_GRP_solver_Edir (double *D, double *U, const struct i_f_var *ifv_L, const struct i_f_var *ifv_R, const double eps, const double atc) |
void | linear_GRP_solver_Edir_Q1D (double *wave_speed, double *D, double *U, double *U_star, const struct i_f_var *ifv_L, const struct i_f_var *ifv_R, const double eps, const double atc) |
A Quasi-1D direct Eulerian GRP solver for unsteady compressible inviscid two-component flow in two space dimension. 更多... | |
void | linear_GRP_solver_Edir_G2D (double *wave_speed, double *D, double *U, double *U_star, const struct i_f_var *ifv_L, const struct i_f_var *ifv_R, const double eps, const double atc) |
A Genuinely-2D direct Eulerian GRP solver for unsteady compressible inviscid two-component flow in two space dimension. 更多... | |
void | AcousticRLagTangent (double *dire, double *U_star, const struct i_f_var *ifv_L, const struct i_f_var *ifv_R, const double r, const double M, const double eps) |
void | GRPsolverRLag (double *wave_speed, double *dire, double *U_star, const struct i_f_var *ifv_L, const struct i_f_var *ifv_R, const double r, const double M, const double eps, const double atc) |
void | HLL_2D_solver (double *F, double *lambda_max, const struct i_f_var *ifv_L, const struct i_f_var *ifv_R) |
A HLL approxiamate Riemann solver for unsteady compressible inviscid single-component flow in two space dimension. 更多... | |
void | Roe_solver (double *F, double *lambda_max, const struct i_f_var *ifv_L, const struct i_f_var *ifv_R, const double delta) |
An approxiamate Riemann solver of Roe for unsteady compressible inviscid single-component flow in one space dimension. 更多... | |
void | Roe_2D_solver (double *F, double *lambda_max, const struct i_f_var *ifv_L, const struct i_f_var *ifv_R, const double delta) |
An approxiamate Riemann solver of Roe for unsteady compressible inviscid single-component flow in two space dimension. 更多... | |
void | Roe_HLL_solver (double *V_mk, double *F, double *lambda_max, const struct i_f_var *ifv_L, const struct i_f_var *ifv_R, const double delta) |
An approxiamate Riemann solver hybridizing the Roe flux and the HLL flux for unsteady compressible inviscid single-component flow in two space dimension. 更多... | |
This file is the header file of several Riemann solvers and GRP solvers.
This header file declares functions in the folder 'Riemann_solver'.
在文件 riemann_solver.h 中定义.
#define Riemann_solver_exact_single Riemann_solver_exact_Ben |
Which solver is chosen as the exact Riemann solver for single-component flow.
在文件 riemann_solver.h 第 17 行定义.
void AcousticRLagTangent | ( | double * | dire, |
double * | U_star, | ||
const struct i_f_var * | ifv_L, | ||
const struct i_f_var * | ifv_R, | ||
const double | r, | ||
const double | M, | ||
const double | eps | ||
) |
void GRPsolverRLag | ( | double * | wave_speed, |
double * | dire, | ||
double * | U_star, | ||
const struct i_f_var * | ifv_L, | ||
const struct i_f_var * | ifv_R, | ||
const double | r, | ||
const double | M, | ||
const double | eps, | ||
const double | atc | ||
) |
void HLL_2D_solver | ( | double * | F, |
double * | lambda_max, | ||
const struct i_f_var * | ifv_L, | ||
const struct i_f_var * | ifv_R | ||
) |
A HLL approxiamate Riemann solver for unsteady compressible inviscid single-component flow in two space dimension.
[out] | F | All four fluxes. |
[out] | lambda_max | Maximum characteristic velocity. |
[in] | ifv_L | Left States (rho_L, u_L, v_L, p_L, gamma, n_x, n_y). |
[in] | ifv_R | Right States (rho_R, u_R, v_R, p_R).
|
在文件 hll_2D_solver.c 第 25 行定义.
void linear_GRP_solver_Edir | ( | double * | D, |
double * | U, | ||
const struct i_f_var * | ifv_L, | ||
const struct i_f_var * | ifv_R, | ||
const double | eps, | ||
const double | atc | ||
) |
void linear_GRP_solver_Edir_G2D | ( | double * | wave_speed, |
double * | D, | ||
double * | U, | ||
double * | U_star, | ||
const struct i_f_var * | ifv_L, | ||
const struct i_f_var * | ifv_R, | ||
const double | eps, | ||
const double | atc | ||
) |
A Genuinely-2D direct Eulerian GRP solver for unsteady compressible inviscid two-component flow in two space dimension.
[out] | wave_speed | the velocity of left and right waves. |
[out] | D | the temporal derivative of fluid variables. [rho, u, v, p, phi, z_a]_t |
[out] | U | the intermediate Riemann solutions at t-axis. [rho_mid, u_mid, v_mid, p_mid, phi_mid, z_a_mid] |
[out] | U_star | the Riemann solutions in star region. [rho_star_L, u_star, rho_star_R, p_star, c_star_L, c_star_R] |
[in] | ifv_L | Left States (rho/u/v/p/phi/z, d_, t_, gammaL). |
[in] | ifv_R | Right States (rho/u/v/p/phi/z, d_, t_, gammaR).
|
[in] | eps | the largest value could be seen as zero. |
[in] | atc | Parameter that determines the solver type.
|
在文件 linear_grp_solver_Edir_G2D.c 第 48 行定义.
void linear_GRP_solver_Edir_Q1D | ( | double * | wave_speed, |
double * | D, | ||
double * | U, | ||
double * | U_star, | ||
const struct i_f_var * | ifv_L, | ||
const struct i_f_var * | ifv_R, | ||
const double | eps, | ||
const double | atc | ||
) |
A Quasi-1D direct Eulerian GRP solver for unsteady compressible inviscid two-component flow in two space dimension.
[out] | wave_speed | the velocity of left and right waves. |
[out] | D | the temporal derivative of fluid variables. [rho, u, v, p, phi, z_a]_t |
[out] | U | the intermediate Riemann solutions at t-axis. [rho_mid, u_mid, v_mid, p_mid, phi_mid, z_a_mid] |
[out] | U_star | the Riemann solutions in star region. [rho_star_L, u_star, rho_star_R, p_star, c_star_L, c_star_R] |
[in] | ifv_L | Left States (rho/u/v/p/phi/z, d_, t_, gammaL). |
[in] | ifv_R | Right States (rho/u/v/p/phi/z, d_, t_, gammaR).
|
[in] | eps | the largest value could be seen as zero. |
[in] | atc | Parameter that determines the solver type.
|
在文件 linear_grp_solver_Edir_Q1D.c 第 38 行定义.
void linear_GRP_solver_LAG | ( | double * | D, |
double * | U, | ||
const struct i_f_var * | ifv_L, | ||
const struct i_f_var * | ifv_R, | ||
const double | eps, | ||
const double | atc | ||
) |
double Riemann_solver_exact | ( | double * | U_star, |
double * | P_star, | ||
const double | gammaL, | ||
const double | gammaR, | ||
const double | u_L, | ||
const double | u_R, | ||
const double | p_L, | ||
const double | p_R, | ||
const double | c_L, | ||
const double | c_R, | ||
_Bool * | CRW, | ||
const double | eps, | ||
const double | tol, | ||
const int | N | ||
) |
EXACT RIEMANN SOLVER FOR Two-Component γ-Law Gas
The purpose of this function is to solve the Riemann problem exactly, for the time dependent one dimensional Euler equations for two-component γ-law gas.
[out] | U_star,P_star | Velocity/Pressure in star region. |
[in] | u_L,p_L,c_L | Initial Velocity/Pressure/sound_speed on left state. |
[in] | u_R,p_R,c_R | Initial Velocity/Pressure/sound_speed on right state. |
[in] | gammaL,gammaR | Ratio of specific heats. |
[out] | CRW | Centred Rarefaction Wave (CRW) Indicator of left and right waves.
|
[in] | eps | The largest value can be seen as zero. |
[in] | tol | Condition value of 'gap' at the end of the iteration. |
[in] | N | Maximum iteration step. |
在文件 riemann_solver_exact_Ben.c 第 30 行定义.
double Riemann_solver_exact_Ben | ( | double * | U_star, |
double * | P_star, | ||
const double | gamma, | ||
const double | u_L, | ||
const double | u_R, | ||
const double | p_L, | ||
const double | p_R, | ||
const double | c_L, | ||
const double | c_R, | ||
_Bool * | CRW, | ||
const double | eps, | ||
const double | tol, | ||
const int | N | ||
) |
EXACT RIEMANN SOLVER FOR A γ-Law Gas
The purpose of this function is to solve the Riemann problem exactly, for the time dependent one dimensional Euler equations for a γ-law gas.
[out] | U_star,P_star | Velocity/Pressure in star region. |
[in] | u_L,p_L,c_L | Initial Velocity/Pressure/sound_speed on left state. |
[in] | u_R,p_R,c_R | Initial Velocity/Pressure/sound_speed on right state. |
[in] | gamma | Ratio of specific heats. |
[out] | CRW | Centred Rarefaction Wave (CRW) Indicator of left and right waves.
|
[in] | eps | The largest value can be seen as zero. |
[in] | tol | Condition value of 'gap' at the end of the iteration. |
[in] | N | Maximum iteration step. |
在文件 riemann_solver_exact_Ben.c 第 230 行定义.
double Riemann_solver_exact_Toro | ( | double * | U_star, |
double * | P_star, | ||
const double | gamma, | ||
const double | U_l, | ||
const double | U_r, | ||
const double | P_l, | ||
const double | P_r, | ||
const double | c_l, | ||
const double | c_r, | ||
_Bool * | CRW, | ||
const double | eps, | ||
const double | tol, | ||
const int | N | ||
) |
EXACT RIEMANN SOLVER FOR THE EULER EQUATIONS
The purpose of this function is to solve the Riemann problem exactly, for the time dependent one dimensional Euler equations for an ideal gas.
[out] | U_star,P_star | Velocity/Pressure in star region. |
[in] | U_l,P_l,c_l | Initial Velocity/Pressure/sound_speed on left state. |
[in] | U_r,P_r,c_r | Initial Velocity/Pressure/sound_speed on right state. |
[in] | gamma | Ratio of specific heats. |
[out] | CRW | Centred Rarefaction Wave (CRW) Indicator of left and right waves.
|
[in] | eps | The largest value can be seen as zero. |
[in] | tol | Condition value of 'gap' at the end of the iteration. |
[in] | N | Maximum iteration step. |
在文件 riemann_solver_exact_Toro.c 第 36 行定义.
double Riemann_solver_starPU | ( | double * | U_star, |
double * | P_star, | ||
const double | GammaL, | ||
const double | GammaR, | ||
const double | UL, | ||
const double | UR, | ||
const double | PL, | ||
const double | PR, | ||
const double | CL, | ||
const double | CR, | ||
_Bool * | CRW, | ||
const double | eps, | ||
const double | TOLPRE, | ||
const int | NRITER | ||
) |
void Roe_2D_solver | ( | double * | F, |
double * | lambda_max, | ||
const struct i_f_var * | ifv_L, | ||
const struct i_f_var * | ifv_R, | ||
const double | delta | ||
) |
An approxiamate Riemann solver of Roe for unsteady compressible inviscid single-component flow in two space dimension.
[out] | F | All four fluxes. |
[out] | lambda_max | Maximum characteristic velocity. |
[in] | ifv_L | Left States (rho_L, u_L, v_L, p_L, gamma, n_x, n_y). |
[in] | ifv_R | Right States (rho_R, u_R, v_R, p_R).
|
[in] | delta | Parameter to modify the modulus of the eigenvalues. |
在文件 roe_2D_solver.c 第 26 行定义.
void Roe_HLL_solver | ( | double * | V_mk, |
double * | F, | ||
double * | lambda_max, | ||
const struct i_f_var * | ifv_L, | ||
const struct i_f_var * | ifv_R, | ||
const double | delta | ||
) |
An approxiamate Riemann solver hybridizing the Roe flux and the HLL flux for unsteady compressible inviscid single-component flow in two space dimension.
[out] | V_mk | ? |
[out] | F | All four fluxes. |
[out] | lambda_max | Maximum characteristic velocity. |
[in] | ifv_L | Left States (rho_L, u_L, v_L, p_L, gamma). |
[in] | ifv_R | Right States (rho_R, u_R, v_R, p_R).
|
[in] | delta | Parameter to modify the modulus of the eigenvalues. |
在文件 roe_hll_solver.c 第 28 行定义.
void Roe_solver | ( | double * | F, |
double * | lambda_max, | ||
const struct i_f_var * | ifv_L, | ||
const struct i_f_var * | ifv_R, | ||
const double | delta | ||
) |
An approxiamate Riemann solver of Roe for unsteady compressible inviscid single-component flow in one space dimension.
[out] | F | All three fluxes. |
[out] | lambda_max | Maximum characteristic velocity. |
[in] | ifv_L | Left States (rho_L, u_L, p_L, gamma). |
[in] | ifv_R | Right States (rho_R, u_R, p_R).
|
[in] | delta | Parameter to modify the modulus of the eigenvalues. (NOT USED!) |
在文件 roe_solver.c 第 25 行定义.