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_y.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_y(const int m, const int n, const int nt, const double tau, struct cell_var_stru * CV,
31 struct b_f_var * bfv_D, struct b_f_var * bfv_U, const _Bool Transversal)
32{
33 double const eps = config[4]; // the largest value could be seen as zero
34 double const h_y = config[11]; // the length of the initial y spatial grids
35 struct i_f_var ifv_D = {.n_x = 0.0, .n_y = 1.0}, ifv_U = {.n_x = 0.0, .n_y = 1.0};
36 int i, j, data_err;
37
38//===========================
39 for(j = 0; j < m; ++j)
40 for(i = 0; i <= n; ++i)
41 {
42 if(i)
43 {
44 ifv_D.d_rho = CV->t_rho[j][i-1];
45 ifv_D.d_u = CV->t_u[j][i-1];
46 ifv_D.d_v = CV->t_v[j][i-1];
47 ifv_D.d_p = CV->t_p[j][i-1];
48 ifv_D.RHO = CV[nt].RHO[j][i-1] + 0.5*h_y*CV->t_rho[j][i-1];
49 ifv_D.U = CV[nt].U[j][i-1] + 0.5*h_y* CV->t_u[j][i-1];
50 ifv_D.V = CV[nt].V[j][i-1] + 0.5*h_y* CV->t_v[j][i-1];
51 ifv_D.P = CV[nt].P[j][i-1] + 0.5*h_y* CV->t_p[j][i-1];
52 }
53 else
54 {
55 ifv_D.d_rho = bfv_D[j].TRHO;
56 ifv_D.d_u = bfv_D[j].TU;
57 ifv_D.d_v = bfv_D[j].TV;
58 ifv_D.d_p = bfv_D[j].TP;
59 ifv_D.RHO = bfv_D[j].RHO + 0.5*h_y*bfv_D[j].TRHO;
60 ifv_D.U = bfv_D[j].U + 0.5*h_y*bfv_D[j].TU;
61 ifv_D.V = bfv_D[j].V + 0.5*h_y*bfv_D[j].TV;
62 ifv_D.P = bfv_D[j].P + 0.5*h_y*bfv_D[j].TP;
63 }
64 if(i < n)
65 {
66 ifv_U.d_rho = CV->t_rho[j][i];
67 ifv_U.d_u = CV->t_u[j][i];
68 ifv_U.d_v = CV->t_v[j][i];
69 ifv_U.d_p = CV->t_p[j][i];
70 ifv_U.RHO = CV[nt].RHO[j][i] - 0.5*h_y*CV->t_rho[j][i];
71 ifv_U.U = CV[nt].U[j][i] - 0.5*h_y* CV->t_u[j][i];
72 ifv_U.V = CV[nt].V[j][i] - 0.5*h_y* CV->t_v[j][i];
73 ifv_U.P = CV[nt].P[j][i] - 0.5*h_y* CV->t_p[j][i];
74 }
75 else
76 {
77 ifv_U.d_rho = bfv_U[j].TRHO;
78 ifv_U.d_u = bfv_U[j].TU;
79 ifv_U.d_v = bfv_U[j].TV;
80 ifv_U.d_p = bfv_U[j].TP;
81 ifv_U.RHO = bfv_U[j].RHO - 0.5*h_y*bfv_U[j].TRHO;
82 ifv_U.U = bfv_U[j].U - 0.5*h_y*bfv_U[j].TU;
83 ifv_U.V = bfv_U[j].V - 0.5*h_y*bfv_U[j].TV;
84 ifv_U.P = bfv_U[j].P - 0.5*h_y*bfv_U[j].TP;
85 }
86 if(ifv_D.P < eps || ifv_U.P < eps || ifv_D.RHO < eps || ifv_U.RHO < eps)
87 {
88 printf("<0.0 error on [%d, %d, %d] (nt, x, y) - Reconstruction_y\n", nt, j, i);
89 return 1;
90 }
91 if(!isfinite(ifv_D.d_p)|| !isfinite(ifv_U.d_p)|| !isfinite(ifv_D.d_u)|| !isfinite(ifv_U.d_u)|| !isfinite(ifv_D.d_v)|| !isfinite(ifv_U.d_v)|| !isfinite(ifv_D.d_rho)|| !isfinite(ifv_U.d_rho))
92 {
93 printf("NAN or INFinite error on [%d, %d, %d] (nt, x, y) - d_Slope_y\n", nt, j, i);
94 return 1;
95 }
96
97//===========================
98 if (Transversal)
99 {
100 if(i)
101 {
102 ifv_D.t_rho = -CV->s_rho[j][i-1];
103 ifv_D.t_u = - CV->s_u[j][i-1];
104 ifv_D.t_v = - CV->s_v[j][i-1];
105 ifv_D.t_p = - CV->s_p[j][i-1];
106 }
107 else
108 {
109 ifv_D.t_rho = -bfv_D[j].SRHO;
110 ifv_D.t_u = -bfv_D[j].SU;
111 ifv_D.t_v = -bfv_D[j].SV;
112 ifv_D.t_p = -bfv_D[j].SP;
113 }
114 if(i < n)
115 {
116 ifv_U.t_rho = -CV->s_rho[j][i];
117 ifv_U.t_u = - CV->s_u[j][i];
118 ifv_U.t_v = - CV->s_v[j][i];
119 ifv_U.t_p = - CV->s_p[j][i];
120 }
121 else
122 {
123 ifv_U.t_rho = -bfv_U[j].SRHO;
124 ifv_U.t_u = -bfv_U[j].SU;
125 ifv_U.t_v = -bfv_U[j].SV;
126 ifv_U.t_p = -bfv_U[j].SP;
127 }
128 if(!isfinite(ifv_D.t_p)|| !isfinite(ifv_U.t_p)|| !isfinite(ifv_D.t_u)|| !isfinite(ifv_U.t_u)|| !isfinite(ifv_D.t_v)|| !isfinite(ifv_U.t_v)|| !isfinite(ifv_D.t_rho)|| !isfinite(ifv_U.t_rho))
129 {
130 printf("NAN or INFinite error on [%d, %d, %d] (nt, x, y) - t_Slope_y\n", nt, j, i);
131 return 1;
132 }
133 }
134 else
135 {
136 ifv_D.t_rho = -0.0;
137 ifv_D.t_u = -0.0;
138 ifv_D.t_v = -0.0;
139 ifv_D.t_p = -0.0;
140 ifv_U.t_rho = -0.0;
141 ifv_U.t_u = -0.0;
142 ifv_U.t_v = -0.0;
143 ifv_U.t_p = -0.0;
144 }
145//===========================
146
147 data_err = GRP_2D_flux(&ifv_D, &ifv_U, tau);
148 switch (data_err)
149 {
150 case 1:
151 printf("<0.0 error on [%d, %d, %d] (nt, x, y) - STAR_y\n", nt, j, i);
152 return 2;
153 case 2:
154 printf("NAN or INFinite error on [%d, %d, %d] (nt, x, y) - STAR_y\n", nt, j, i);
155 return 2;
156 case 3:
157 printf("NAN or INFinite error on [%d, %d, %d] (nt, x, y) - DIRE_y\n", nt, j, i);
158 return 2;
159 }
160
161 CV->G_rho[j][i] = ifv_D.F_rho;
162 CV->G_u[j][i] = ifv_D.F_u;
163 CV->G_v[j][i] = ifv_D.F_v;
164 CV->G_e[j][i] = ifv_D.F_e;
165
166 CV->rhoIy[j][i] = ifv_D.RHO_int;
167 CV->uIy[j][i] = ifv_D.U_int;
168 CV->vIy[j][i] = ifv_D.V_int;
169 CV->pIy[j][i] = ifv_D.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_y(const int m, const int n, const int nt, const double tau, struct cell_var_stru *CV, struct b_f_var *bfv_D, struct b_f_var *bfv_U, const _Bool Transversal)
This function calculate Eulerian fluxes of 2-D Euler equations in y-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 ** G_rho
Definition: var_struc.h:43
double ** U
Definition: var_struc.h:36
double ** rhoIy
Definition: var_struc.h:41
double ** t_rho
Definition: var_struc.h:39
double ** s_rho
Definition: var_struc.h:38
double ** RHO
Definition: var_struc.h:36
double ** vIy
Definition: var_struc.h:41
double ** s_u
Definition: var_struc.h:38
double ** t_u
Definition: var_struc.h:39
double ** G_v
numerical fluxes at (y_{j-1/2}, t_{n}).
Definition: var_struc.h:43
double ** G_e
Definition: var_struc.h:43
double ** s_v
Definition: var_struc.h:38
double ** t_v
Definition: var_struc.h:39
double ** G_u
Definition: var_struc.h:43
double ** uIy
Definition: var_struc.h:41
double ** V
Definition: var_struc.h:36
double ** s_p
spatial derivatives in coordinate x (slopes).
Definition: var_struc.h:38
double ** pIy
interfacial variable values in coordinate y at t_{n+1}.
Definition: var_struc.h:41
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