HydroCODE_2D 0.1
This is a implementation of fully explict forward Euler scheme for 2-D Euler equations of motion on Eulerian coordinate
flux_generator_x.c
浏览该文件的文档.
1
6#include <stdio.h>
7#include <math.h>
8
9#include "../include/var_struc.h"
10#include "../include/flux_calc.h"
11
12
30int flux_generator_x(const int m, const int n, const int nt, const double tau, struct cell_var_stru * CV,
31 struct b_f_var * bfv_L, struct b_f_var * bfv_R, const _Bool Transversal)
32{
33 double const eps = config[4]; // the largest value could be seen as zero
34 double const h_x = config[10]; // the length of the initial x spatial grids
35 struct i_f_var ifv_L = {.n_x = 1.0, .n_y = 0.0}, ifv_R = {.n_x = 1.0, .n_y = 0.0};
36 int i, j, data_err;
37
38//===========================
39 for(i = 0; i < n; ++i)
40 for(j = 0; j <= m; ++j)
41 {
42 if(j)
43 {
44 ifv_L.d_rho = CV->s_rho[j-1][i];
45 ifv_L.d_u = CV->s_u[j-1][i];
46 ifv_L.d_v = CV->s_v[j-1][i];
47 ifv_L.d_p = CV->s_p[j-1][i];
48 ifv_L.RHO = CV[nt].RHO[j-1][i] + 0.5*h_x*CV->s_rho[j-1][i];
49 ifv_L.U = CV[nt].U[j-1][i] + 0.5*h_x* CV->s_u[j-1][i];
50 ifv_L.V = CV[nt].V[j-1][i] + 0.5*h_x* CV->s_v[j-1][i];
51 ifv_L.P = CV[nt].P[j-1][i] + 0.5*h_x* CV->s_p[j-1][i];
52 }
53 else
54 {
55 ifv_L.d_rho = bfv_L[i].SRHO;
56 ifv_L.d_u = bfv_L[i].SU;
57 ifv_L.d_v = bfv_L[i].SV;
58 ifv_L.d_p = bfv_L[i].SP;
59 ifv_L.RHO = bfv_L[i].RHO + 0.5*h_x*bfv_L[i].SRHO;
60 ifv_L.U = bfv_L[i].U + 0.5*h_x*bfv_L[i].SU;
61 ifv_L.V = bfv_L[i].V + 0.5*h_x*bfv_L[i].SV;
62 ifv_L.P = bfv_L[i].P + 0.5*h_x*bfv_L[i].SP;
63 }
64 if(j < m)
65 {
66 ifv_R.d_rho = CV->s_rho[j][i];
67 ifv_R.d_u = CV->s_u[j][i];
68 ifv_R.d_v = CV->s_v[j][i];
69 ifv_R.d_p = CV->s_p[j][i];
70 ifv_R.RHO = CV[nt].RHO[j][i] - 0.5*h_x*CV->s_rho[j][i];
71 ifv_R.U = CV[nt].U[j][i] - 0.5*h_x* CV->s_u[j][i];
72 ifv_R.V = CV[nt].V[j][i] - 0.5*h_x* CV->s_v[j][i];
73 ifv_R.P = CV[nt].P[j][i] - 0.5*h_x* CV->s_p[j][i];
74 }
75 else
76 {
77 ifv_R.d_rho = bfv_R[i].SRHO;
78 ifv_R.d_u = bfv_R[i].SU;
79 ifv_R.d_v = bfv_R[i].SV;
80 ifv_R.d_p = bfv_R[i].SP;
81 ifv_R.RHO = bfv_R[i].RHO - 0.5*h_x*bfv_R[i].SRHO;
82 ifv_R.U = bfv_R[i].U - 0.5*h_x*bfv_R[i].SU;
83 ifv_R.V = bfv_R[i].V - 0.5*h_x*bfv_R[i].SV;
84 ifv_R.P = bfv_R[i].P - 0.5*h_x*bfv_R[i].SP;
85 }
86 if(ifv_L.P < eps || ifv_R.P < eps || ifv_L.RHO < eps || ifv_R.RHO < eps)
87 {
88 printf("<0.0 error on [%d, %d, %d] (nt, x, y) - Reconstruction_x\n", nt, j, i);
89 return 1;
90 }
91 if(!isfinite(ifv_L.d_p)|| !isfinite(ifv_R.d_p)|| !isfinite(ifv_L.d_u)|| !isfinite(ifv_R.d_u)|| !isfinite(ifv_L.d_v)|| !isfinite(ifv_R.d_v)|| !isfinite(ifv_L.d_rho)|| !isfinite(ifv_R.d_rho))
92 {
93 printf("NAN or INFinite error on [%d, %d, %d] (nt, x, y) - d_Slope_x\n", nt, j, i);
94 return 1;
95 }
96
97//===========================
98 if (Transversal)
99 {
100 if(j)
101 {
102 ifv_L.t_rho = CV->t_rho[j-1][i];
103 ifv_L.t_u = CV->t_u[j-1][i];
104 ifv_L.t_v = CV->t_v[j-1][i];
105 ifv_L.t_p = CV->t_p[j-1][i];
106 }
107 else
108 {
109 ifv_L.t_rho = bfv_L[i].TRHO;
110 ifv_L.t_u = bfv_L[i].TU;
111 ifv_L.t_v = bfv_L[i].TV;
112 ifv_L.t_p = bfv_L[i].TP;
113 }
114 if(j < m)
115 {
116 ifv_R.t_rho = CV->t_rho[j][i];
117 ifv_R.t_u = CV->t_u[j][i];
118 ifv_R.t_v = CV->t_v[j][i];
119 ifv_R.t_p = CV->t_p[j][i];
120 }
121 else
122 {
123 ifv_R.t_rho = bfv_R[i].TRHO;
124 ifv_R.t_u = bfv_R[i].TU;
125 ifv_R.t_v = bfv_R[i].TV;
126 ifv_R.t_p = bfv_R[i].TP;
127 }
128 if(!isfinite(ifv_L.t_p)|| !isfinite(ifv_R.t_p)|| !isfinite(ifv_L.t_u)|| !isfinite(ifv_R.t_u)|| !isfinite(ifv_L.t_v)|| !isfinite(ifv_R.t_v)|| !isfinite(ifv_L.t_rho)|| !isfinite(ifv_R.t_rho))
129 {
130 printf("NAN or INFinite error on [%d, %d, %d] (nt, x, y) - t_Slope_x\n", nt, j, i);
131 return 1;
132 }
133 }
134 else
135 {
136 ifv_L.t_rho = 0.0;
137 ifv_L.t_u = 0.0;
138 ifv_L.t_v = 0.0;
139 ifv_L.t_p = 0.0;
140 ifv_R.t_rho = 0.0;
141 ifv_R.t_u = 0.0;
142 ifv_R.t_v = 0.0;
143 ifv_R.t_p = 0.0;
144 }
145//===========================
146
147 data_err = GRP_2D_flux(&ifv_L, &ifv_R, tau);
148 switch (data_err)
149 {
150 case 1:
151 printf("<0.0 error on [%d, %d, %d] (nt, x, y) - STAR_x\n", nt, j, i);
152 return 2;
153 case 2:
154 printf("NAN or INFinite error on [%d, %d, %d] (nt, x, y) - STA_x\n", nt, j, i);
155 return 2;
156 case 3:
157 printf("NAN or INFinite error on [%d, %d, %d] (nt, x, y) - DIRE_x\n", nt, j, i);
158 return 2;
159 }
160
161 CV->F_rho[j][i] = ifv_L.F_rho;
162 CV->F_u[j][i] = ifv_L.F_u;
163 CV->F_v[j][i] = ifv_L.F_v;
164 CV->F_e[j][i] = ifv_L.F_e;
165
166 CV->rhoIx[j][i] = ifv_L.RHO_int;
167 CV->uIx[j][i] = ifv_L.U_int;
168 CV->vIx[j][i] = ifv_L.V_int;
169 CV->pIx[j][i] = ifv_L.P_int;
170 }
171 return 0;
172}
int GRP_2D_flux(struct i_f_var *ifv, struct i_f_var *ifv_R, const double tau)
This function calculate Eulerian fluxes of 2-D Euler equations by 2-D GRP solver.
Definition: flux_solver.c:24
int flux_generator_x(const int m, const int n, const int nt, const double tau, struct cell_var_stru *CV, struct b_f_var *bfv_L, struct b_f_var *bfv_R, const _Bool Transversal)
This function calculate Eulerian fluxes of 2-D Euler equations in x-direction by 2-D GRP solver.
Fluid VARiables at Boundary.
Definition: var_struc.h:63
double SRHO
Definition: var_struc.h:65
double V
Definition: var_struc.h:64
double SV
spatial derivatives in coordinate x (slopes).
Definition: var_struc.h:65
double TRHO
Definition: var_struc.h:66
double RHO
Definition: var_struc.h:64
double U
Definition: var_struc.h:64
double P
Definition: var_struc.h:64
double SU
Definition: var_struc.h:65
double SP
Definition: var_struc.h:65
double TV
spatial derivatives in coordinate y (slopes).
Definition: var_struc.h:66
double TU
Definition: var_struc.h:66
double TP
Definition: var_struc.h:66
pointer structure of VARiables on STRUctural computational grid CELLs.
Definition: var_struc.h:35
double ** U
Definition: var_struc.h:36
double ** t_rho
Definition: var_struc.h:39
double ** F_v
numerical fluxes at (x_{j-1/2}, t_{n}).
Definition: var_struc.h:42
double ** pIx
interfacial variable values in coordinate x at t_{n+1}.
Definition: var_struc.h:40
double ** s_rho
Definition: var_struc.h:38
double ** RHO
Definition: var_struc.h:36
double ** s_u
Definition: var_struc.h:38
double ** F_e
Definition: var_struc.h:42
double ** t_u
Definition: var_struc.h:39
double ** F_u
Definition: var_struc.h:42
double ** F_rho
Definition: var_struc.h:42
double ** s_v
Definition: var_struc.h:38
double ** t_v
Definition: var_struc.h:39
double ** rhoIx
Definition: var_struc.h:40
double ** V
Definition: var_struc.h:36
double ** s_p
spatial derivatives in coordinate x (slopes).
Definition: var_struc.h:38
double ** uIx
Definition: var_struc.h:40
double ** vIx
Definition: var_struc.h:40
double ** P
Definition: var_struc.h:36
double ** t_p
spatial derivatives in coordinate y (slopes).
Definition: var_struc.h:39
Interfacial Fluid VARiables.
Definition: var_struc.h:47
double t_rho
Definition: var_struc.h:53
double V
variable values at t_{n}.
Definition: var_struc.h:49
double V_int
interfacial variables at t_{n+1}.
Definition: var_struc.h:50
double RHO
Definition: var_struc.h:49
double t_v
tangential spatial derivatives OR spatial derivatives in Lagrangian coordinate ξ
Definition: var_struc.h:53
double d_p
Definition: var_struc.h:52
double U_int
Definition: var_struc.h:50
double U
Definition: var_struc.h:49
double t_p
Definition: var_struc.h:53
double P
Definition: var_struc.h:49
double F_rho
Definition: var_struc.h:51
double RHO_int
Definition: var_struc.h:50
double n_x
Definition: var_struc.h:48
double d_rho
Definition: var_struc.h:52
double F_e
Definition: var_struc.h:51
double d_u
Definition: var_struc.h:52
double t_u
Definition: var_struc.h:53
double F_u
Definition: var_struc.h:51
double F_v
interfacial fluxes at t_{n+1/2}.
Definition: var_struc.h:51
double d_v
normal spatial derivatives.
Definition: var_struc.h:52
double P_int
Definition: var_struc.h:50
double config[]
Initial configuration data array.
Definition: hydrocode.c:115