HydroCODE_1D 0.1
This is a implementation of fully explict forward Euler scheme for 1-D Euler equations of motion on Lagrangian/Eulerian coordinate
riemann_solver.h 文件参考

This file is the header file of several Riemann solvers and GRP solvers. 更多...

riemann_solver.h 的引用(Include)关系图:
此图展示该文件直接或间接的被哪些文件引用了:

浏览源代码.

宏定义

#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)
 A Lagrangian GRP solver for unsteady compressible inviscid two-component flow in one space dimension. 更多...
 
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)
 A direct Eulerian GRP solver for unsteady compressible inviscid flow in one space dimension. 更多...
 
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)
 
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)
 
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)
 
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)
 
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)
 
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)
 

详细描述

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 中定义.

宏定义说明

◆ Riemann_solver_exact_single

#define Riemann_solver_exact_single   Riemann_solver_exact_Ben

Which solver is chosen as the exact Riemann solver for single-component flow.

在文件 riemann_solver.h17 行定义.

函数说明

◆ AcousticRLagTangent()

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 
)

◆ GRPsolverRLag()

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 
)

◆ HLL_2D_solver()

void HLL_2D_solver ( double *  F,
double *  lambda_max,
const struct i_f_var ifv_L,
const struct i_f_var ifv_R 
)

◆ linear_GRP_solver_Edir()

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 
)

A direct Eulerian GRP solver for unsteady compressible inviscid flow in one space dimension.

参数
[out]Dthe temporal derivative of fluid variables.
[rho, u, p]_t
[out]Uthe intermediate Riemann solutions at t-axis.
[rho_mid, u_mid, p_mid]
[in]ifv_LLeft States (rho_L, u_L, p_L, s_rho_L, s_u_L, s_p_L, gamma).
[in]ifv_RRight States (rho_R, u_R, p_R, s_rho_R, s_u_R, s_p_R).
  • s_rho, s_u, s_p: x-spatial derivatives.
  • gamma: the constant of the perfect gas.
[in]epsthe largest value could be seen as zero.
[in]atcParameter that determines the solver type.
  • INFINITY: acoustic approximation
    • ifv_.s_ = -0.0: exact Riemann solver
  • eps: 1D GRP solver(nonlinear + acoustic case)
  • -0.0: 1D GRP solver(only nonlinear case)
参见
Theory is found in Reference [1].
[1] M. Ben-Artzi, J. Li & G. Warnecke, A direct Eulerian GRP scheme for compressible fluid flows. Journal of Computational Physics, 218.1: 19-43, 2006.

在文件 linear_grp_solver_Edir.c33 行定义.

这是这个函数的调用关系图:

◆ linear_GRP_solver_Edir_G2D()

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 
)

◆ linear_GRP_solver_Edir_Q1D()

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 
)

◆ linear_GRP_solver_LAG()

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 
)

A Lagrangian GRP solver for unsteady compressible inviscid two-component flow in one space dimension.

参数
[out]Dthe temporal derivative of fluid variables in the Star Region.
[rho_L, u, p, rho_R]_t
[out]Uthe Riemann solutions in the Star Region.
[rho_star_L, u_star, p_star, rho_star_R]
[in]ifv_LLeft States (rho_L, u_L, p_L, s_rho_L, s_u_L, s_p_L, gammaL).
[in]ifv_RRight States (rho_R, u_R, p_R, s_rho_R, s_u_R, s_p_R, gammaR).
  • s_rho, s_u, s_p: ξ-Lagrangian spatial derivatives.
  • gamma: the constant of the perfect gas.
[in]epsthe largest value could be seen as zero.
[in]atcParameter that determines the solver type.
  • INFINITY: acoustic approximation
  • eps: GRP solver(nonlinear + acoustic case)
  • -0.0: GRP solver(only nonlinear case)
参见
Theory is found in Reference [1].
[1] M. Ben-Artzi & J. Falcovitz, A second-order Godunov-type scheme for compressible fluid dynamics, Journal of Computational Physics, 55.1: 1-32, 1984

在文件 linear_grp_solver_LAG.c32 行定义.

函数调用图:
这是这个函数的调用关系图:

◆ Riemann_solver_exact()

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_starVelocity/Pressure in star region.
[in]u_L,p_L,c_LInitial Velocity/Pressure/sound_speed on left state.
[in]u_R,p_R,c_RInitial Velocity/Pressure/sound_speed on right state.
[in]gammaL,gammaRRatio of specific heats.
[out]CRWCentred Rarefaction Wave (CRW) Indicator of left and right waves.
  • true: CRW
  • false: Shock wave
[in]epsThe largest value can be seen as zero.
[in]tolCondition value of 'gap' at the end of the iteration.
[in]NMaximum iteration step.
返回
gap: Relative pressure change after the last iteration.

在文件 riemann_solver_exact_Ben.c30 行定义.

这是这个函数的调用关系图:

◆ Riemann_solver_exact_Ben()

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_starVelocity/Pressure in star region.
[in]u_L,p_L,c_LInitial Velocity/Pressure/sound_speed on left state.
[in]u_R,p_R,c_RInitial Velocity/Pressure/sound_speed on right state.
[in]gammaRatio of specific heats.
[out]CRWCentred Rarefaction Wave (CRW) Indicator of left and right waves.
  • true: CRW
  • false: Shock wave
[in]epsThe largest value can be seen as zero.
[in]tolCondition value of 'gap' at the end of the iteration.
[in]NMaximum iteration step.
返回
gap: Relative pressure change after the last iteration.

在文件 riemann_solver_exact_Ben.c230 行定义.

◆ Riemann_solver_exact_Toro()

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_starVelocity/Pressure in star region.
[in]U_l,P_l,c_lInitial Velocity/Pressure/sound_speed on left state.
[in]U_r,P_r,c_rInitial Velocity/Pressure/sound_speed on right state.
[in]gammaRatio of specific heats.
[out]CRWCentred Rarefaction Wave (CRW) Indicator of left and right waves.
  • true: CRW
  • false: Shock wave
[in]epsThe largest value can be seen as zero.
[in]tolCondition value of 'gap' at the end of the iteration.
[in]NMaximum iteration step.
返回
gap: Relative pressure change after the last iteration.
Bug:
Some program errors in this exact Riemann solver!
作者
E. F. Toro
日期
February 1st 1999
参见
Theory is found in Chapter 4 of Reference [1].
[1] E. F. Toro, "Riemann Solvers and Numerical Methods for Fluid Dynamics". Springer-Verlag, Second Edition, 1999

在文件 riemann_solver_exact_Toro.c36 行定义.

◆ Riemann_solver_starPU()

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 
)

◆ Roe_2D_solver()

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 
)

◆ Roe_HLL_solver()

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 
)

◆ Roe_solver()

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 
)