forked from AstroSnow/PIP
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathshock_tube_ion.f90
258 lines (243 loc) · 8.93 KB
/
shock_tube_ion.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
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
subroutine shock_tube_ion
use parameters,only:pi
use globalvar,only:ix,jx,kx,U_h,U_m,flag_bnd,col,beta,flag_b_stg,dtout,&
flag_mhd,flag_mpi,my_rank,flag_pip,gm,beta,tend,flag_ir,&
x,y,z,dx,dy,dz,n_fraction,debug_direction,debug_option,T0,n0,Nexcite,&
f_p_ini,f_p_p_ini,n0fac,Gm_rec_ref
use scheme_rot,only:pv2cq_mhd,pv2cq_hd,get_Te_MHD
use model_rot, only:set_coordinate,setcq
use PIP_rot, only:get_col_ion_coeff
implicit none
double precision :: ro_h(1:ix,1:jx,1:kx),ro_m(1:ix,1:jx,1:kx)
double precision :: vx_h(1:ix,1:jx,1:kx),vx_m(1:ix,1:jx,1:kx)
double precision :: vy_h(1:ix,1:jx,1:kx),vy_m(1:ix,1:jx,1:kx)
double precision :: vz_h(1:ix,1:jx,1:kx),vz_m(1:ix,1:jx,1:kx)
double precision :: P_h (1:ix,1:jx,1:kx),P_m (1:ix,1:jx,1:kx)
double precision :: B_x (1:ix,1:jx,1:kx)
double precision :: B_y (1:ix,1:jx,1:kx)
double precision :: B_z (1:ix,1:jx,1:kx)
double precision ::mask(ix,jx,kx)
double precision ::v_para(ix,jx,kx),v_perp(ix,jx,kx)
double precision ::b_para(ix,jx,kx),b_perp(ix,jx,kx)
double precision ::v_b_para(ix,jx,kx),b_b_para(ix,jx,kx)
double precision ::v_theta(ix,jx,kx),v_phi(ix,jx,kx)
double precision ::b_theta(ix,jx,kx),b_phi(ix,jx,kx)
double precision ::Te_p(ix,jx,kx),Gm_ion0(ix,jx,kx),Gm_rec0(ix,jx,kx)
double precision::Eion(6) !Energy to ionise
double precision f_n,f_p,f_p_n,f_p_p,start(3),end(3),B0
double precision theta_p,phi_p,tmp,v_L(8),v_R(8),wtr,ioneq,Te_0,tfac
double precision,parameter::kbhat=1.38064852,mehat=9.10938356,hhat=6.62607004
double precision,parameter::kboltz=1.38064852e-23 !Boltzmann Constant [m^2 kg s^-2 K^-1]
integer i,j,k
!Find the equilibrium neutral fraction
if(flag_IR .eq. 4) then
print*,'Calculating LTE excitation state'
allocate(Nexcite(ix,jx,kx,6)) !Allocate the fractional array
Eion=[13.6,3.4,1.51,0.85,0.54,0.0] !in eV
Eion=Eion/13.6*2.18e-18 !Convert to joules (to be dimensionally correct)
Nexcite(:,:,:,6)=n0
Nexcite(:,:,:,1)=(2.0/n0/2.d0*(2.0*pi*mehat*kbhat*T0/hhat/hhat*1.0e14)**(3.0/2.0)*exp(-Eion(1)/kboltz/T0))
Nexcite(:,:,:,2)=(2.0/n0/8.d0*(2.0*pi*mehat*kbhat*T0/hhat/hhat*1.0e14)**(3.0/2.0)*exp(-Eion(2)/kboltz/T0))
Nexcite(:,:,:,3)=(2.0/n0/18.d0*(2.0*pi*mehat*kbhat*T0/hhat/hhat*1.0e14)**(3.0/2.0)*exp(-Eion(3)/kboltz/T0))
Nexcite(:,:,:,4)=(2.0/n0/32.d0*(2.0*pi*mehat*kbhat*T0/hhat/hhat*1.0e14)**(3.0/2.0)*exp(-Eion(4)/kboltz/T0))
Nexcite(:,:,:,5)=(2.0/n0/50.d0*(2.0*pi*mehat*kbhat*T0/hhat/hhat*1.0e14)**(3.0/2.0)*exp(-Eion(5)/kboltz/T0))
Nexcite(:,:,:,1:5)=n0/Nexcite(:,:,:,1:5)
Nexcite=Nexcite/n0
Nexcite=Nexcite/sum(Nexcite(1,1,1,:))
f_n=sum(Nexcite(1,1,1,1:5))
f_p=Nexcite(1,1,1,6)
f_p_n=f_n/(f_n+2.0d0*f_p)
f_p_p=2.0d0*f_p/(f_n+2.0d0*f_p)
f_p_ini=f_p
f_p_p_ini=f_p_p
n0fac=Nexcite(1,1,1,6)
! print*,Nexcite(1,1,1,:)
! print*,f_n,f_p
! stop
else
Te_0=T0/1.1604e4
ioneq=(2.6e-19/dsqrt(Te_0))/(2.91e-14/(0.232+13.6/Te_0)*(13.6/Te_0)**0.39*dexp(-13.6/Te_0))
f_n=ioneq/(ioneq+1.0d0)
f_p=1.0d0-f_n
f_p_n=f_n/(f_n+2.0d0*f_p)
f_p_p=2.0d0*f_p/(f_n+2.0d0*f_p)
endif
if (my_rank.eq.0) print*,'Neutral fraction = ',f_n
! !set ionization fraction-----------------
if(flag_pip.eq.0) then
f_n=1.0d0
f_p=1.0d0
f_p_n=1.0d0
f_p_p=1.0d0
endif
!----------------------------------------
!Set coordinate (uniform grid)--------------------------
!!set lower and upper coordinate
start(1)=0.0d0/f_p ;end(1)=600.d0!30000.0d0 !400.0d0/f_p
start(2)=-1.0d0 ;end(2)=1.0d0
start(3)=-1.0d0 ;end(3)=1.0d0
call set_coordinate(start,end)
!---------------------------------------
!!default boundary condition----------------------
if (flag_bnd(1) .eq.-1) flag_bnd(1)=3
if (flag_bnd(2) .eq.-1) flag_bnd(2)=10
if (flag_bnd(3) .eq.-1) flag_bnd(3)=10
if (flag_bnd(4) .eq.-1) flag_bnd(4)=10
if (flag_bnd(5) .eq.-1) flag_bnd(5)=10
if (flag_bnd(6) .eq.-1) flag_bnd(6)=10
!-------------------------------------------------
!!!========================================================
!parameters---------------------
! B0=sqrt(2.0*gm/beta)
B0=1.0d0
! B0=0.0d0
!/sqrt(2.0)
!-------------------------
if(debug_direction.eq.1) then
mask=spread(spread(pi*x,2,jx),3,kx)
phi_p=0.0d0
theta_p=pi/2.0d0
tmp=0.0
else if(debug_direction.eq.2) then
mask=spread(spread(pi*y,1,ix),3,kx)
phi_p=pi/2.0d0
theta_p=pi/2.0d0
tmp=0.0
else if(debug_direction.eq.3) then
mask=spread(spread(pi*z,1,jx),1,ix)
phi_p=0.0d0
theta_p=0.0d0
tmp=0.0
else if(debug_direction.eq.4) then
do k=1,kx;do j=1,jx; do i=1,ix
if(x(i)+y(j).gt.0)then
mask(i,j,k)=1
else
mask(i,j,k)=-1
endif
enddo;enddo;enddo
phi_p=pi/4.0d0
theta_p=pi/2.0d0
tmp=pi/2.0d0
else if(debug_direction.eq.5) then
do k=1,kx;do j=1,jx; do i=1,ix
if(x(i)+z(k).gt.0)then
mask(i,j,k)=1
else
mask(i,j,k)=-1
endif
enddo;enddo;enddo
phi_p=0.0d0
theta_p=pi/4.0d0
tmp=0.0
else if(debug_direction.eq.6) then
do k=1,kx;do j=1,jx; do i=1,ix
if(y(j)+z(k).gt.0)then
mask(i,j,k)=1
else
mask(i,j,k)=-1
endif
enddo;enddo;enddo
phi_p=pi/2.0d0
theta_p=pi/4.0d0
tmp=0.0
else if(debug_direction.eq.7) then
do k=1,kx;do j=1,jx; do i=1,ix
if(x(i)+y(j)+z(k).gt.0)then
mask(i,j,k)=1
else
mask(i,j,k)=-1
endif
enddo;enddo;enddo
phi_p=pi/4.0d0
theta_p=pi/4.0d0
tmp=0.0
endif
select case(debug_option)
case (2) !Switch-off shock (Falle et al.1998)
v_l=(/1.368d0,1.769d0,0.269d0,1.0d0,0.0d0,1.0d0,0.0d0,0.0d0/)
v_r=(/1.0d0,1.0d0,0.0d0,0.0d0,0.0d0,1.0d0,1.0d0,0.0d0/)
case default
! v_l=(/1.0d0,1.0d0,0.0d0,0.0d0,0.0d0,B0*0.75d0,B0,0.0d0/)
! v_r=(/0.125d0,0.1d0,0.0d0,0.0d0,0.0d0,B0*0.75d0,-B0,0.0d0/)
v_l=(/1.0d0,beta*B0**2/2.d0,0.0d0,0.0d0,0.0d0,B0*0.3d0,B0,0.0d0/)
v_r=(/1.0d0,beta*B0**2/2.d0,0.0d0,0.0d0,0.0d0,B0*0.3d0,-B0,0.0d0/)
! Density, pressure, velocity * 3, magnetic field *3
end select
where(mask<0)
ro_h=f_n*v_l(1)
p_h=f_p_n*v_l(2)
ro_m=f_p*v_l(1)
p_m=f_p_p*v_l(2)
vx_h=v_l(3)
vy_h=v_l(4)
vz_h=v_l(5)
vx_m=v_l(3)
vy_m=v_l(4)
vz_m=v_l(5)
b_para=v_l(6)
b_perp=v_l(7)
elsewhere
ro_h=f_n*v_r(1)
p_h=f_p_n*v_r(2)
ro_m=f_p*v_r(1)
p_m=f_p_p*v_r(2)
vx_h=v_r(3)
vy_h=v_r(4)
vz_h=v_r(5)
vx_m=v_r(3)
vy_m=v_r(4)
vz_m=v_r(5)
b_para=v_r(6)
b_perp=v_r(7)
end where
! wtr=0.001d0
wtr=dx(1)*10.d0
!0.001d0
do k=1,kx;do j=1,jx;do i=1,ix
ro_h(i,j,k)=f_n*(v_l(1)+(v_r(1)-v_l(1))*(1.0d0+tanh(mask(i,j,k)/wtr))*0.5d0)
ro_m(i,j,k)=f_p*(v_l(1)+(v_r(1)-v_l(1))*(1.0d0+tanh(mask(i,j,k)/wtr))*0.5d0)
p_h(i,j,k)=f_p_n*(v_l(2)+(v_r(2)-v_l(2))*(1.0d0+tanh(mask(i,j,k)/wtr))*0.5d0)
p_m(i,j,k)=f_p_p*(v_l(2)+(v_r(2)-v_l(2))*(1.0d0+tanh(mask(i,j,k)/wtr))*0.5d0)
vx_h(i,j,k)=(v_l(3)+(v_r(3)-v_l(3))*(1.0d0+tanh(mask(i,j,k)/wtr))*0.5d0)
vx_m(i,j,k)=(v_l(3)+(v_r(3)-v_l(3))*(1.0d0+tanh(mask(i,j,k)/wtr))*0.5d0)
vy_h(i,j,k)=(v_l(4)+(v_r(4)-v_l(4))*(1.0d0+tanh(mask(i,j,k)/wtr))*0.5d0)
vy_m(i,j,k)=(v_l(4)+(v_r(4)-v_l(4))*(1.0d0+tanh(mask(i,j,k)/wtr))*0.5d0)
vz_h(i,j,k)=(v_l(5)+(v_r(5)-v_l(5))*(1.0d0+tanh(mask(i,j,k)/wtr))*0.5d0)
vz_m(i,j,k)=(v_l(5)+(v_r(5)-v_l(5))*(1.0d0+tanh(mask(i,j,k)/wtr))*0.5d0)
b_para(i,j,k)=(v_l(6)+(v_r(6)-v_l(6))*(1.0d0+tanh(mask(i,j,k)/wtr))*0.5d0)
b_perp(i,j,k)=(v_l(7)+(v_r(7)-v_l(7))*(1.0d0+tanh(mask(i,j,k)/wtr))*0.5d0)
enddo;enddo;enddo
b_theta=b_perp*sin(tmp)
b_phi =b_perp*cos(tmp)
b_x=b_para*sin(theta_p)*cos(phi_p)+ &
b_theta*cos(theta_p)*cos(phi_p)-b_phi*sin(phi_p)
b_y=b_para*sin(theta_p)*sin(phi_p)+ &
b_theta*cos(theta_p)*sin(phi_p)+b_phi*cos(phi_p)
b_z=b_para*cos(theta_p)-b_theta*sin(theta_p)
!!!========================================================
!convert PV2cq and set that value to global variable 'U_h' and/or 'U_m'
call setcq(ro_m,vx_m,vy_m,vz_m,p_m,B_x,B_y,B_z, &
ro_h,vx_h,vy_h,vz_h,p_h)
!---------------------------------------------------------------------
!Set the reference recombination rate
! if (flag_IR .eq. 4) then
! call get_Te_MHD(U_m,Te_p)
! tfac=beta/2.0d0*f_p_p_ini*5.0d0/6.0d0/f_p_ini
!print*,'Te_p',Te_p
!print*,'T0',T0
!print*,'tfac',tfac
!print*,U_m(:,:,:,1)
!print*,n0
!print*,n0fac
! call get_col_ion_coeff(Te_p*T0/tfac,U_m(:,:,:,1)*n0/n0fac,Gm_ion0,Gm_rec0)
! Gm_rec_ref=Gm_rec0(1,1,1)
! endif
!set default output period and time duration--------------------------
if(tend.lt.0.0) then
tend=0.4d0
dtout=tend/10.0
if(flag_mpi.eq.0 .or.my_rank.eq.0) print *,"TEND",dtout,tend
endif
!---------------------------------------------------------------------
end subroutine shock_tube_ion