forked from AstroSnow/PIP
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathrelax_prom.f90
279 lines (219 loc) · 8.85 KB
/
relax_prom.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
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
subroutine relax_prom
use parameters,only:pi
use globalvar,only:ix,jx,kx,U_h,U_m,flag_bnd,beta,flag_b_stg,flag_divb &
,dtout,&
flag_mhd,flag_mpi,my_rank,flag_pip,gm,beta,tend,mpi_pos,mpi_siz,&
x,y,z,dx,dy,dz,n_fraction,gra,flag_grav,scl_height,margin
use scheme_rot,only:pv2cq_mhd,pv2cq_hd
use model_rot, only:set_coordinate,setcq
use matrix_rot,only:inverse_tridiagonal
use boundary_rot,only:PIPbnd ! here all the boundary routines are called
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 f_n,f_p,f_p_n,f_p_p,start(3),end(3)
double precision pr_val
real*4 :: lenx,lenz,znull,zaxis
integer :: num,nx,nz,i,j,opt,div_x,div_y,start_x(2),start_y(2),o_height(1)
character*40 :: filename
double precision, allocatable:: b_mag(:,:)
real*4 :: bx(1000,1000),by(1000,1000),bz(1000,1000)
double precision, allocatable:: p_tot(:),ro_tot(:),temp(:),y_1(:)
double precision t0,t1,trw,trp,b_norm
!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
else
f_n=n_fraction
f_p=1.0d0-n_fraction
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
!----------------------------------------
! load the file
num=1 ! 1 is the x-point, 3 is the no-x-point
if(num.eq.0) stop 'Quit'
write (filename,fmt='("shear_",i1,".dat")') num
if (my_rank.eq.1) print *,'Reading data from '//filename
open(unit=1,file=filename,form='unformatted',status='old')
read (1) opt
if (my_rank.eq.3) print *,'File type: opt=',opt
if(.not.(opt.eq.1)) stop 'Unknown file type'
read (1) nx, nz ! read array dimensions
if (my_rank.eq.1) print 10,nx ,nz
if (my_rank.eq.0) print 20,ix,jx
10 format(" nx =",i4,", nz =",i4)
20 format(" ix=",i4,", jx=",i4)
! read in the grid
if (my_rank.eq.2) print *, 'Start read'
read (1) lenx,lenz,znull,zaxis
if (my_rank.eq.1) print *, lenx,lenz,znull,zaxis
if (my_rank.eq.2) print *, 'loaded units'
read (1) ((bx(i,j),i=1,nx),j=1,nz)
if (my_rank.eq.3) print *, 'loaded b_x'
read (1) ((by(i,j),i=1,nx),j=1,nz)
if (my_rank.eq.0) print *, 'loaded B_y'
read (1) ((bz(i,j),i=1,nx),j=1,nz)
if (my_rank.eq.1) print *, 'loaded B_z'
close(1)
if (my_rank.eq.2) print *,'OK'
! NEED TO SET UP B FOR DIFFERENT CORES!!
div_x=ix-2*margin(1)
div_y=jx-2*margin(2)
if (my_rank.eq.3) print *, div_x,div_y, my_rank
start_x(1)=1+mpi_pos(1)*div_x
start_x(2)=(mpi_pos(1)+1)*div_x+2*margin(1)
start_y(1)=1+mpi_pos(2)*div_y
start_y(2)=(mpi_pos(2)+1)*div_y+2*margin(2)
if (my_rank.eq.0) print *, ix,jx,kx
if (my_rank.eq.1) print *, start_x,start_y,my_rank
! deal with Margin!!
!-----------------------------------------------------
!Set coordinate (uniform grid)--------------------------
!!set lower and upper coordinate
start(1)=-2.5d0 ;end(1)=2.5d0
start(2)=0.0d0 ;end(2)=5.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)=2
if (flag_bnd(2) .eq.-1) flag_bnd(2)=2
if (flag_bnd(3) .eq.-1) flag_bnd(3)=4
if (flag_bnd(4) .eq.-1) flag_bnd(4)=2
if (flag_bnd(5) .eq.-1) flag_bnd(5)=1
if (flag_bnd(6) .eq.-1) flag_bnd(6)=1
!-------------------------------------------------
!!!========================================================
!density of lower fluid is unity
nz=1000
nx=1000
flag_grav=1
gra(:,:,:,:)=0.0d0
scl_height=1.d0
gra(:,:,:,2)=-1.0d0/(gm*scl_height)
if (mpi_pos(2).eq.0) then
gra(:,1:margin(2),:,2)=1.0d0/(gm*scl_height)
endif
if (mpi_pos(2).eq.(mpi_siz(2)-1)) then
gra(:,jx-margin(2)+1:jx,:,2)=1.0d0/(gm*scl_height)
endif
allocate(p_tot(nz+2*margin(2)),ro_tot(nz+2*margin(2)) &
,temp(nz+2*margin(2)),y_1(nz+2*margin(2)))
t0 = 0.02d0
t1 = 1.0d0
trp = 0.3d0
trw = 0.09d0
y_1=((/ (DBLE(i),i=1,nz+2*margin(2))/)-1.d0)*dy(1)
temp(:)=t0 + (t1 - t0) * 0.5d0*(tanh((y_1(:)-trp)/trw)+1.d0)
! pr_val=1.d0/gm
! do i=margin(2)+1,nz+2*margin(2)
! if (i .eq. margin(2)+1) then
! p_tot(i)=pr_val
! else
! p_tot(i)=p_tot(i-1) - dy(1) * p_tot(i-1)/temp(i-1)
! endif
! enddo
! p_tot(1:margin(2))=p_tot(margin(2)*2:margin(2)+1:-1)
! ro_tot(:)=p_tot(:) *gm/temp(:)
!--------------------
! read in relaxed background
write (filename,fmt='("pr_n_ro.dat")')
if (my_rank.eq.1) print *,'Reading data from '//filename
open(unit=1,file=filename,form='unformatted',status='old')
read (1) (p_tot(i),i=1,nz+2*margin(2))
read (1) (ro_tot(i),i=1,nz+2*margin(2))
close(1)
! if (my_rank.eq.3) print *, p_tot
!---------------------
! if (my_rank.eq.0) print *, ro_tot
! print *, 2*margin(2)+start_y
! if (my_rank.eq.0) print *, p_tot
! MAKW MPI READY
! P_m= pr_val*f_p_p
P_m=spread(spread(p_tot(start_y(1):start_y(2)) &
*f_p_p,1,ix),3,kx)
! ro_m= 1.d0*f_p
ro_m= spread(spread(ro_tot(start_y(1):start_y(2)) &
*f_p ,1,ix),3,kx)
! P_h= pr_val*f_p_n
P_h= spread(spread(p_tot(start_y(1):start_y(2)) &
*f_p_n ,1,ix),3,kx)
! ro_h= 1.d0*f_n
ro_h= spread(spread(ro_tot(start_y(1):start_y(2)) &
*f_n ,1,ix),3,kx)
o_height(1)=maxloc(abs(bz(499,:)),1)
print *, o_height
b_norm=sqrt(p_tot(328)*2.d0/beta)/abs(bz(499,o_height(1)))
!sqrt(p_tot(164)/2.d0/beta)/bz(499,164)
! if (my_rank.eq.0) print *, b_norm,bz(499,328),maxloc(abs(bz(499,:)))
! STOP
vx_h=0.0d0;vy_h=0.0d0;vz_h=0.0d0
vx_m=0.0d0;vy_m=0.0d0;vz_m=0.0d0
! B_x=0.0d0;B_y=0.0d0;B_z=0.0d0
allocate(b_mag(nx+2*margin(1),nz+2*margin(2)))
b_mag(margin(1)+1:nx+margin(1),margin(2)+1:nz+margin(2))= &
bx
b_mag(1:margin(1),:)=-b_mag(2*margin(1):margin(1)+1:-1,:)
b_mag(nx+margin(1)+1:nx+2*margin(1),:)= &
-b_mag(nx+1:nx+margin(1):-1,:)
b_mag(:,1:margin(2))=-b_mag(:,2*margin(2):margin(2)+1:-1)
b_mag(:,nz+margin(2)+1:nz+2*margin(2))= &
b_mag(:,nz+margin(2):nz+1:-1)
B_x=DBLE(spread(b_mag(start_x(1):start_x(2),start_y(1) &
:start_y(2)),3,kx))*b_norm
if (my_rank.eq.3) print *, 'done bx', my_rank
b_mag(margin(1)+1:nx+margin(1),margin(2)+1:nz+margin(2))= &
by
b_mag(1:margin(1),:)=b_mag(2*margin(1):margin(1)+1:-1,:)
b_mag(nx+margin(1)+1:nx+2*margin(1),:)= &
b_mag(nx+1:nx+margin(1):-1,:)
b_mag(:,1:margin(2))=b_mag(:,2*margin(2):margin(2)+1:-1)
b_mag(:,nz+margin(2)+1:nz+2*margin(2))= &
-b_mag(:,nz+margin(2):nz+1:-1)
B_y=DBLE(spread(b_mag(start_x(1):start_x(2),start_y(1) &
:start_y(2)),3,kx))*b_norm
if (my_rank.eq.2) print *, 'done by',my_rank
b_mag(margin(1)+1:nx+margin(1),margin(2)+1:nz+margin(2))= &
bz
b_mag(1:margin(1),:)=-b_mag(2*margin(1):margin(1)+1:-1,:)
b_mag(nx+margin(1)+1:nx+2*margin(1),:)= &
-b_mag(nx+1:nx+margin(1):-1,:)
b_mag(:,1:margin(2))=-b_mag(:,2*margin(2):margin(2)+1:-1)
b_mag(:,nz+margin(2)+1:nz+2*margin(2))= &
b_mag(:,nz+margin(2):nz+1:-1)
B_z=DBLE(spread(b_mag(start_x(1):start_x(2),start_y(1) &
:start_y(2)),3,kx))*b_norm
if (my_rank.eq.1) print *, 'done bz',my_rank
! print *, B_z(200,1:12,1)
! STOP
! B_x=0.0d0;B_y=0.0d0;B_z=0.0d0
!!!=================
! Add mass
ro_m=ro_m+25.d0*f_p*spread(0.5d0*(1.d0-dtanh((abs(spread(x,2,jx))-0.1d0 )/0.1d0)) &
*0.5d0*(1.d0-dtanh((abs(spread(y-1.1d0,1,ix))-0.3d0 )/0.05d0)),3,kx)*ro_tot(328)
ro_h=ro_h+25.d0*f_n*spread(0.5d0*(1.d0-dtanh((abs(spread(x,2,jx))-0.1d0 )/0.1d0)) &
*0.5d0*(1.d0-dtanh((abs(spread(y-1.1d0,1,ix))-0.3d0 )/0.05d0)),3,kx)*ro_tot(328)
deallocate(ro_tot,temp,y_1)
!!!========================================================
!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 default output period and time duration--------------------------
if(tend.lt.0.0) then
tend=40.0d0
dtout=tend/40.0
if(flag_mpi.eq.0 .or.my_rank.eq.0) print *,"TEND",dtout,tend
endif
!---------------------------------------------------------------------
end subroutine relax_prom