-
Notifications
You must be signed in to change notification settings - Fork 11
/
MODULE_INITIAL_CONDITION.f90
113 lines (94 loc) · 4.26 KB
/
MODULE_INITIAL_CONDITION.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
!=================================================================================================
!> CARTESIAN ADAPTIVE FIVE EQUATION MODEL
!> AUTHOR: VAN-DAT THANG
!> E-MAIL: [email protected]
!> E-MAIL: [email protected]
!> SOURCE CODE LINK: https://github.com/dattv/2D_CARFIVE
!=================================================================================================
MODULE MODULE_INITIALCONDITION
use MODULE_PRECISION
use MODULE_NODECOORD
use MODULE_CONSTANTS
use MODULE_QUADTREE
use MODULE_GENERICMETHOD
use MODULE_CFDMAINDATA
use MODULE_STIFFNESSGAS
contains
!==================================================================================================
subroutine initial_2D_simle_check_grid(nelm, tree)
implicit none
integer(ip), intent(in) :: nelm
type(quadtree), dimension(:), intent(inout), target :: tree
integer(ip) :: i
!call loop_on_quadtree_array(1, nelm, tree, initial_2d_dambreak_single_level)
call loop_on_quadtree_array(1, nelm, tree, initial_2d_01_simple)
return
end subroutine initial_2D_simle_check_grid
!==================================================================================================
subroutine initial_condition(nelm, tree)
implicit none
integer(ip), intent(in) :: nelm
type(quadtree), dimension(:), intent(inout), target :: tree
call loop_on_quadtree_array(1, nelm, tree, initial_2d_dambreak_single_level)
return
end subroutine initial_condition
!==================================================================================================
subroutine initial_2d_dambreak_single_level(tree)
implicit none
type(quadtree), pointer, intent(inout) :: tree
real(rp) :: rho1_a1, rho2_a2, rho_u, rho_v, p, a1, a2, rho_e
real(rp) :: gamma_mix, pi_mix, u_vel, v_vel, rho_mix, rho
! body
rho = 1.e3_rp
if (tree%pts(5)%coord(1) <= half) then
a1 = one - tolerance; a2 = tolerance
rho1_a1 = rho*a1; rho2_a2 = rho*a2; rho_u = 1500._rp; rho_v = zero; p = 1.e9_rp
else
a2 = one - tolerance; a1 = tolerance
rho1_a1 = rho*a1; rho2_a2 = rho*a2; rho_u = 1500._rp; rho_v = zero; p = 1.e9_rp
endif
! ===> NONE DIMENSIONAL <======================================================================
rho1_a1 = rho1_a1/Rho_inf
rho2_a2 = rho2_a2/Rho_inf
rho_u = rho_u/ u_inf
rho_v = rho_v/ u_inf
p = p/p_inf
gamma_mix = compute_gamma_mixture(a1, a2, matInfo(1)%gamma, matInfo(2)%gamma)
pi_mix = compute_pi_mixture(a1, a2, matInfo(1)%pi, matInfo(2)%pi, matInfo(1)%gamma, matInfo(2)%gamma)
rho_mix = rho1_a1 + rho2_a2
u_vel = rho_u/rho_mix
v_vel = rho_v/rho_mix
rho_E = compute_rhoE(P, gamma_mix, pi_mix) + half*(rho_mix*(u_vel**2 + v_vel**2))
tree%data%u(:) = (/rho1_a1, rho2_a2, rho_u, rho_v, rho_E, a1/)
tree%data%w(:) = (/rho1_a1/a1, rho2_a2/(one - a1), u_vel, v_vel, p, a1/)
return
end subroutine initial_2d_dambreak_single_level
!==================================================================================================
subroutine initial_2d_01_simple(tree)
implicit none
type(quadtree), pointer, intent(inout) :: tree
if (tree%pts(5)%coord(1) >= half ) then
tree%data%u(:) = one
tree%data%w(:) = one
else
tree%data%u(:) = zero
tree%data%w(:) = zero
endif
return
end subroutine initial_2d_01_simple
!==================================================================================================
subroutine initial_2d_02_simple(tree)
implicit none
type(quadtree), pointer, intent(inout) :: tree
if (tree%pts(5)%coord(1) >= 0.7_rp) then
tree%data%u(:) = one
tree%data%w(:) = one
else
tree%data%u(:) = zero
tree%data%w(:) = zero
endif
return
end subroutine initial_2d_02_simple
!==================================================================================================
!==================================================================================================
END MODULE MODULE_INITIALCONDITION