1 #include "../src/cpp_macros.h"
5 program compressible_flow_wall
9 integer,
parameter :: n_vars = 2+ndim
10 integer,
parameter :: box_size = 8
13 integer,
parameter :: i_rho = 1
14 integer,
parameter :: i_mom(NDIM) = [2, 3]
15 integer,
parameter :: i_e = 4
18 real(dp),
parameter :: k_boltzmann = 1.380649e-23_dp
19 real(dp),
parameter :: dalton = 1.66053906660e-27_dp
21 real(dp) :: euler_gamma = 1.4_dp
22 real(dp) :: cfl_number = 0.5_dp
23 integer :: time_integrator = af_heuns_method
24 real(dp) :: rho_0, u_0, p_0, T_0, sound_speed, mach_number
25 real(dp) :: molecular_mass, domain_size(NDIM)
26 real(dp) :: wedge_angle, wedge_normal(2)
27 integer :: variables(n_vars)
28 integer :: cc_rho, cc_mom(2), cc_e
29 integer :: fluxes(n_vars)
30 integer :: grid_size(NDIM)
32 real(dp) :: u_init(1, n_vars)
34 real(dp) :: dt, dt_lim, time, end_time, dt_output
35 character(len=100) :: fname
36 integer :: n, output_cnt
39 character(len=10) :: var_names(n_vars)
44 domain_size = [0.45_dp, 0.3_dp]
45 grid_size = [12, 8] * box_size
48 molecular_mass = 28.9_dp
51 wedge_angle = 15.0_dp * acos(-1.0_dp) / 180.0_dp
52 wedge_normal = [-sin(wedge_angle), cos(wedge_angle)]
55 rho_0 = p_0 * molecular_mass * dalton / (k_boltzmann * t_0)
56 sound_speed = sqrt(euler_gamma * p_0 / rho_0)
58 u_0 = mach_number * sound_speed
60 u_init(1, :) = [rho_0, u_0, 0.0_dp, p_0]
61 call to_conservative(1, n_vars, u_init)
63 var_names(i_rho) =
"rho"
64 var_names(i_mom) = [
"mom_x",
"mom_y"]
68 call af_add_cc_variable(tree, var_names(n), ix=variables(n), n_copies=2)
69 call af_add_fc_variable(tree,
"flux", ix=fluxes(n))
72 call af_set_cc_methods(tree, variables(i_rho), bc_gas_rho)
73 call af_set_cc_methods(tree, variables(i_mom(1)), bc_gas_mom_x)
74 call af_set_cc_methods(tree, variables(i_mom(2)), bc_gas_mom_y)
75 call af_set_cc_methods(tree, variables(i_e), bc_gas_e)
77 call af_add_cc_variable(tree,
"lsf", ix=i_lsf)
78 call af_set_cc_methods(tree, i_lsf, funcval=set_lsf_box)
80 cc_rho = variables(i_rho)
81 cc_mom = variables(i_mom)
86 call af_init(tree, box_size, domain_size, grid_size)
87 call af_loop_box(tree, set_init_conds)
88 call af_gc_tree(tree, variables)
94 dt_output = end_time / 100
97 iteration = iteration + 1
98 if (output_cnt * dt_output <= time)
then
99 write(fname,
"(A,I0)")
"output/compressible_flow_wall_" &
100 // dimname //
"_", output_cnt
101 call af_write_silo(tree, trim(fname), output_cnt, time, &
102 add_vars = write_primitives, add_names=[
"xVel",
"yVel",
"pres"])
103 output_cnt = output_cnt + 1
106 call af_advance(tree, dt, dt_lim, time, variables, &
107 time_integrator, forward_euler)
108 dt = cfl_number * dt_lim
109 if (time > end_time)
exit
112 print *,
"Number of time steps: ", iteration
113 print *,
"Number of grid cells: ", af_num_leaves_used(tree) * box_size**ndim
114 call af_destroy(tree)
118 subroutine forward_euler(tree, dt, dt_stiff, dt_lim, time, s_deriv, n_prev, &
119 s_prev, w_prev, s_out, i_step, n_steps)
120 type(af_t),
intent(inout) :: tree
121 real(dp),
intent(in) :: dt
122 real(dp),
intent(in) :: dt_stiff
123 real(dp),
intent(inout) :: dt_lim
124 real(dp),
intent(in) :: time
125 integer,
intent(in) :: s_deriv
126 integer,
intent(in) :: n_prev
127 integer,
intent(in) :: s_prev(n_prev)
128 real(dp),
intent(in) :: w_prev(n_prev)
129 integer,
intent(in) :: s_out
130 integer,
intent(in) :: i_step, n_steps
131 real(dp) :: dummy_dt(0)
133 call flux_generic_tree(tree, n_vars, variables, s_deriv, fluxes, dt_lim, &
134 max_wavespeed, get_fluxes, flux_dummy_modify, line_modify, &
135 to_primitive, to_conservative, af_limiter_vanleer_t)
136 call flux_update_densities(tree, dt, n_vars, variables, n_vars, variables, fluxes, &
137 s_deriv, n_prev, s_prev, w_prev, s_out, flux_dummy_source, 0, dummy_dt, &
139 end subroutine forward_euler
142 subroutine set_box_mask(box, mask)
143 type(box_t),
intent(in) :: box
144 logical,
intent(out) :: mask(DTIMES(box%n_cell))
148 mask = (box%cc(dtimes(1:nc), i_lsf) > 0.0_dp)
149 end subroutine set_box_mask
151 subroutine set_init_conds(box)
152 type(box_t),
intent(inout) :: box
157 box%cc(ijk, variables) = u_init(1, :)
158 if (box%cc(ijk, i_lsf) <= 0)
then
159 box%cc(ijk, i_mom) = 0.0_dp
162 end subroutine set_init_conds
165 subroutine to_primitive(n_values, n_vars, u)
166 integer,
intent(in) :: n_values, n_vars
167 real(dp),
intent(inout) :: u(n_values, n_vars)
169 u(:, i_mom(1)) = u(:, i_mom(1))/u(:, i_rho)
170 u(:, i_mom(2)) = u(:, i_mom(2))/u(:, i_rho)
171 u(:, i_e) = (euler_gamma-1.0_dp) * (u(:, i_e) - &
172 0.5_dp*u(:, i_rho)* sum(u(:, i_mom(:))**2, dim=2))
173 end subroutine to_primitive
176 subroutine to_conservative(n_values, n_vars, u)
177 integer,
intent(in) :: n_values, n_vars
178 real(dp),
intent(inout) :: u(n_values, n_vars)
179 real(dp) :: kin_en(n_values)
184 kin_en = 0.5_dp * u(:, i_rho) * sum(u(:, i_mom(:))**2, dim=2)
187 inv_fac = 1/(euler_gamma - 1.0_dp)
188 u(:, i_e) = u(:, i_e) * inv_fac + kin_en
192 u(:, i_mom(i)) = u(:, i_rho) * u(:, i_mom(i))
194 end subroutine to_conservative
197 subroutine max_wavespeed(n_values, n_var, flux_dim, u, w)
198 integer,
intent(in) :: n_values
199 integer,
intent(in) :: n_var
200 integer,
intent(in) :: flux_dim
201 real(dp),
intent(in) :: u(n_values, n_var)
202 real(dp),
intent(out) :: w(n_values)
203 real(dp) :: sound_speeds(n_values)
205 sound_speeds = sqrt(euler_gamma * u(:, i_e) / u(:, i_rho))
206 w = sound_speeds + abs(u(:, i_mom(flux_dim)))
207 end subroutine max_wavespeed
210 subroutine get_fluxes(n_values, n_var, flux_dim, u, flux, box, line_ix, s_deriv)
211 integer,
intent(in) :: n_values
212 integer,
intent(in) :: n_var
213 integer,
intent(in) :: flux_dim
214 real(dp),
intent(in) :: u(n_values, n_var)
215 real(dp),
intent(out) :: flux(n_values, n_var)
216 type(box_t),
intent(in) :: box
217 integer,
intent(in) :: line_ix(NDIM-1)
218 integer,
intent(in) :: s_deriv
219 real(dp) :: E(n_values), inv_fac
226 flux(:, i_rho) = u(:, i_rho) * u(:, i_mom(flux_dim))
230 flux(:, i_mom(i)) = u(:, i_rho) * &
231 u(:, i_mom(i)) * u(:, i_mom(flux_dim))
235 flux(:, i_mom(flux_dim)) = flux(:, i_mom(flux_dim)) + u(:, i_e)
238 inv_fac = 1/(euler_gamma-1.0_dp)
239 e = u(:, i_e) * inv_fac + 0.5_dp * u(:, i_rho) * &
240 sum(u(:, i_mom(:))**2, dim=2)
243 flux(:, i_e) = u(:, i_mom(flux_dim)) * (e + u(:, i_e))
245 end subroutine get_fluxes
248 subroutine line_modify(n_cc, n_var, cc_line, flux_dim, box, line_ix, s_deriv)
249 integer,
intent(in) :: n_cc
250 integer,
intent(in) :: n_var
251 real(dp),
intent(inout) :: cc_line(-1:n_cc-2, n_var)
252 integer,
intent(in) :: flux_dim
253 type(box_t),
intent(in) :: box
254 integer,
intent(in) :: line_ix(NDIM-1)
255 integer,
intent(in) :: s_deriv
257 real(dp) :: lsf(0:box%n_cell+1)
258 real(dp) :: mom(2), mom_par(2), mom_perp(2)
262 call flux_get_line_1cc(box, i_lsf, flux_dim, line_ix, lsf)
264 if (all(lsf > 0))
return
267 if (lsf(i) * lsf(i+1) <= 0)
then
270 cc_line(i+1, :) = cc_line(i, :)
271 cc_line(i+2, :) = cc_line(i-1, :)
278 mom = cc_line(i, i_mom)
279 mom_perp = dot_product(mom, wedge_normal) * wedge_normal
280 mom_par = mom - mom_perp
281 cc_line(i+1, i_mom) = mom_par - mom_perp
283 mom = cc_line(i-1, i_mom)
284 mom_perp = dot_product(mom, wedge_normal) * wedge_normal
285 mom_par = mom - mom_perp
286 cc_line(i+2, i_mom) = mom_par - mom_perp
288 cc_line(i, :) = cc_line(i+1, :)
289 cc_line(i-1, :) = cc_line(i+2, :)
291 mom = cc_line(i+1, i_mom)
292 mom_perp = dot_product(mom, wedge_normal) * wedge_normal
293 mom_par = mom - mom_perp
294 cc_line(i, i_mom) = mom_par - mom_perp
296 mom = cc_line(i+2, i_mom)
297 mom_perp = dot_product(mom, wedge_normal) * wedge_normal
298 mom_par = mom - mom_perp
299 cc_line(i-1, i_mom) = mom_par - mom_perp
304 end subroutine line_modify
307 subroutine write_primitives(box, new_vars, n_var)
308 type(box_t),
intent(in) :: box
309 integer,
intent(in) :: n_var
310 real(dp) :: new_vars(DTIMES(0:box%n_cell+1),n_var)
313 new_vars(dtimes(:), 1) = box%cc(dtimes(:), cc_mom(1))/box%cc(dtimes(:), cc_rho)
315 new_vars(dtimes(:), 2) = box%cc(dtimes(:), cc_mom(2))/box%cc(dtimes(:), cc_rho)
317 new_vars(dtimes(:), 3) = (euler_gamma-1.0_dp)*(box%cc(dtimes(:), cc_e) - &
318 sum(box%cc(dtimes(:), cc_mom(:))**2, dim=ndim+1) &
319 /(2.0_dp*box%cc(dtimes(:), cc_rho)))
320 end subroutine write_primitives
323 subroutine bc_gas_rho(box, nb, iv, coords, bc_val, bc_type)
324 type(box_t),
intent(in) :: box
325 integer,
intent(in) :: nb
326 integer,
intent(in) :: iv
327 real(dp),
intent(in) :: coords(NDIM, box%n_cell**(NDIM-1))
328 real(dp),
intent(out) :: bc_val(box%n_cell**(NDIM-1))
329 integer,
intent(out) :: bc_type
331 if (nb == af_neighb_lowx)
then
332 bc_type = af_bc_dirichlet
333 bc_val = u_init(1, i_rho)
335 call af_bc_neumann_zero(box, nb, iv, coords, bc_val, bc_type)
337 end subroutine bc_gas_rho
340 subroutine bc_gas_mom_x(box, nb, iv, coords, bc_val, bc_type)
341 type(box_t),
intent(in) :: box
342 integer,
intent(in) :: nb
343 integer,
intent(in) :: iv
344 real(dp),
intent(in) :: coords(NDIM, box%n_cell**(NDIM-1))
345 real(dp),
intent(out) :: bc_val(box%n_cell**(NDIM-1))
346 integer,
intent(out) :: bc_type
348 if (nb == af_neighb_lowx)
then
349 bc_type = af_bc_dirichlet
350 bc_val = u_init(1, i_mom(1))
352 call af_bc_neumann_zero(box, nb, iv, coords, bc_val, bc_type)
354 end subroutine bc_gas_mom_x
357 subroutine bc_gas_mom_y(box, nb, iv, coords, bc_val, bc_type)
358 type(box_t),
intent(in) :: box
359 integer,
intent(in) :: nb
360 integer,
intent(in) :: iv
361 real(dp),
intent(in) :: coords(NDIM, box%n_cell**(NDIM-1))
362 real(dp),
intent(out) :: bc_val(box%n_cell**(NDIM-1))
363 integer,
intent(out) :: bc_type
365 if (nb == af_neighb_lowx)
then
366 bc_type = af_bc_dirichlet
367 bc_val = u_init(1, i_mom(2))
368 else if (nb == af_neighb_lowy)
then
369 bc_type = af_bc_dirichlet
371 else if (nb == af_neighb_highy)
then
372 bc_type = af_bc_dirichlet
375 call af_bc_neumann_zero(box, nb, iv, coords, bc_val, bc_type)
377 end subroutine bc_gas_mom_y
380 subroutine bc_gas_e(box, nb, iv, coords, bc_val, bc_type)
381 type(box_t),
intent(in) :: box
382 integer,
intent(in) :: nb
383 integer,
intent(in) :: iv
384 real(dp),
intent(in) :: coords(NDIM, box%n_cell**(NDIM-1))
385 real(dp),
intent(out) :: bc_val(box%n_cell**(NDIM-1))
386 integer,
intent(out) :: bc_type
388 if (nb == af_neighb_lowx)
then
389 bc_type = af_bc_dirichlet
390 bc_val = u_init(1, i_e)
392 call af_bc_neumann_zero(box, nb, iv, coords, bc_val, bc_type)
394 end subroutine bc_gas_e
396 subroutine set_lsf_box(box, iv)
397 type(box_t),
intent(inout) :: box
398 integer,
intent(in) :: iv
400 real(dp) :: rr(NDIM), norm_dr
403 norm_dr = norm2(box%dr)
406 rr = af_r_cc(box, [ijk])
407 box%cc(ijk, iv) = get_lsf(rr)
409 end subroutine set_lsf_box
412 real(dp) function get_lsf(rr)
413 real(dp),
intent(in) :: rr(NDIM)
417 get_lsf = cos(wedge_angle) * (rr(2) - 0.0_dp) - &
418 sin(wedge_angle) * (rr(1) - 0.15_dp)
Module which contains all Afivo modules, so that a user does not have to include them separately.