-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathstepND_mhd.f90
284 lines (266 loc) · 9.43 KB
/
stepND_mhd.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
280
281
282
283
284
!------------------------------------------------------------------------------!
! NDSPMHD: A Smoothed Particle (Magneto)Hydrodynamics code for (astrophysical) !
! fluid dynamics simulations in 1, 2 and 3 spatial dimensions. !
! !
! (c) 2002-2015 Daniel Price !
! !
! http://users.monash.edu.au/~dprice/ndspmhd !
! daniel.price@monash.edu -or- dprice@cantab.net (forwards to current address) !
! !
! NDSPMHD comes with ABSOLUTELY NO WARRANTY. !
! This is free software; and you are welcome to redistribute !
! it under the terms of the GNU General Public License !
! (see LICENSE file for details) and the provision that !
! this notice remains intact. If you modify this file, please !
! note section 2a) of the GPLv2 states that: !
! !
! a) You must cause the modified files to carry prominent notices !
! stating that you changed the files and the date of any change. !
! !
! ChangeLog: !
!------------------------------------------------------------------------------!
!!------------------------------------------------------------------------
!! Computes one timestep
!! Change this subroutine to change the timestepping algorithm
!!
!! This subroutine is a slightly modified version of the
!! predictor-corrector scheme. Algorithm goes as follows:
!!
!! Predictor:
!!
!! v^{1/2} = v^0 + dt/2*f^(-1/2)
!! x^{1/2} = x^0 + dt/2*v^{1/2} (note used updated v)
!! rho^{1/2} = rho^0 + dt/2*drhodt^{-1/2}
!!
!! --> calculate f^{1/2}, drhodt^{1/2} using x^{1/2} and v^{1/2}
!!
!! Corrector:
!!
!! v^* = v^0 + dt/2*f^{1/2}
!! x^* = x^0 + dt/2*v^* (note uses updated v)
!! rho^* = rho^0 + dt/2*drhodt^*
!!
!! v^1 = 2*v^* - v^0 = v^0 + dt*f^{1/2}
!! x^1 = 2*x^* - x^0 = x^0 + dt*v^*
!! rho^1 = 2*rho^* - rho^0 = rho^0 + dt*drhodt^{1/2}
!!
!! Energy, smoothing length, alpha and magnetic field follow density
!!
!! Gives good results (and good stability) on wave and shock-type problems
!! with a reasonable courant number (C_cour = 0.4 seems to be enough).
!! On these problems seems to do much better than leapfrog/symplectic,
!! mainly related to the derivatives which depend on velocity (for which the
!! other two methods are not reversible/symplectic).
!!
!! However, results for disks and orbits are slightly better with leapfrog
!!
!!---------------------------------------------------------------------------
SUBROUTINE step
USE dimen_mhd
USE debug
USE loguns
USE bound
USE derivB
USE eos
USE hterms
USE options
USE part
USE part_in
USE rates
USE timestep
USE setup_params
USE xsph
USE kernels, only:eps
USE utils, only:cross_product3D
!
!--define local variables
!
IMPLICIT NONE
INTEGER :: i
REAL :: hdt, maxdivB
REAL :: dtrhoi
REAL, DIMENSION(ndimV,SIZE(rho)) :: gradpsiprev
REAL, DIMENSION(SIZE(divB)) :: divBprev
real, dimension(ndimV) :: vcrossB
!
!--allow for tracing flow
!
IF (trace) WRITE(iprint,*) ' Entering subroutine step'
!
!--Mid-point Predictor step
!
hdt = 0.5*dt
DO i=1,npart
xin(:,i) = x(:,i)
velin(:,i) = vel(:,i)
Bevolin(:,i) = Bevol(:,i)
rhoin(i) = rho(i)
hhin(i) = hh(i)
enin(i) = en(i)
alphain(:,i) = alpha(:,i)
psiin(i) = psi(i)
ENDDO
!
!--if doing divergence correction then do correction to magnetic field
!
nsubsteps_divB = -1
IF (idivBzero.GE.10 .and. MOD(nsteps,10).EQ.0) THEN
CALL divBcorrect(npart,ntotal)
Bevolin(:,1:ntotal) = Bevol(:,1:ntotal)
IF (any(ibound.ne.0)) WRITE(iprint,*) 'Warning: boundaries not correct'
ELSEIF (idivBzero.GE.1) THEN
maxdivB = MAXVAL(ABS(divB(1:npart)))
print*,'max div B = ',maxdivB
nsubsteps_divB = -1
ENDIF
DO i=1,npart
IF (itype(i).EQ.itypebnd .OR. itype(i).EQ.itypebnd2) THEN ! fixed particles
vel(:,i) = velin(:,i)
IF (icty.GE.1) rho(i) = rhoin(i)
if (imhd.lt.0) then
call cross_product3D(vel(:,i),Bconst(:),vcrossB)
Bevol(:,i) = Bevolin(:,i) + hdt*vcrossB(:)
else
Bevol(:,i) = Bevolin(:,i)
endif
IF (iener.NE.0) en(i) = enin(i)
hh(i) = hhin(i)
x(:,i) = xin(:,i) + hdt*(vel(1:ndim,i) + xsphfac*xsphterm(1:ndim,i))
alpha(:,i) = alphain(:,i)
psi(i) = psiin(i)
ELSE
vel(:,i) = (velin(:,i) + hdt*force(:,i))/(1.+damp)
!--vertical damping
if (ndim.ge.3) vel(3,i) = vel(3,i)/(1.+dampz)
!--radial damping
if (dampr.gt.tiny(dampr) .and.ndim.ge.2) call radialdamping(vel(1:2,i),x(1:2,i),dampr)
x(:,i) = xin(:,i) + hdt*(vel(1:ndim,i) + xsphfac*xsphterm(1:ndim,i))
IF (imhd.NE.0) THEN
IF (nsubsteps_divB.LT.0) THEN
Bevol(:,i) = Bevolin(:,i) + hdt*gradpsi(:,i)
psi(i) = (psiin(i) + hdt*dpsidt(i))/(1.+eps)
ENDIF
Bevol(:,i) = Bevol(:,i) + hdt*dBevoldt(:,i)
ENDIF
rho(i) = rhoin(i) + hdt*drhodt(i)
IF (ihvar.EQ.1) THEN
! hh(i) = hfact(pmass(i)/rho(i))**dndim ! my version
hh(i) = hhin(i)*(rhoin(i)/rho(i))**dndim ! Joe's
ELSEIF (ihvar.EQ.2 .OR. ihvar.EQ.3) THEN
hh(i) = hhin(i) + hdt*dhdt(i)
ENDIF
IF (iener.NE.0) en(i) = enin(i) + hdt*dendt(i)
IF (ANY(iavlim.NE.0)) alpha(:,i) = alphain(:,i) + hdt*daldt(:,i)
ENDIF
ENDDO
! WHERE (alpha(:,:).GT.0.5)
! alpha(:,:) = 1.0
! END WHERE
!
!--calculate all derivatives
!
call derivs
IF (nsubsteps_divB.GE.0) THEN
DO i=1,npart
Bevol(:,i) = Bevolin(:,i)
psi(i) = psiin(i)
ENDDO
ENDIF
dtrho = huge(dtrho)
DO i=1,npart
dtrhoi = abs(rho(i)/(drhodt(i) + 1.e-8))
dtrho = min(dtrho,dtrhoi)
ENDDO
if (C_rho*dtrho/dtcourant .lt. C_cour) then
write(iprint,*) 'dtrho equiv courant number = ',C_rho*dtrho/dtcourant
endif
!
!--do substepping on corrector for div B correction
!
! IF (imhd.GT.0 .and. idivBzero.GE.2 .and. nsubsteps_divB.GE.0) THEN
! CALL substep_divB(2,dt,nsubsteps_divB,Bevol,psi,divB,gradpsi, &
! x,hh,pmass,itype,npart,ntotal)
! ENDIF
!
!--Mid-point Corrector step
!
DO i=1,npart
IF (itype(i).EQ.itypebnd .or. itype(i).EQ.itypebnd2) THEN
vel(:,i) = velin(:,i)
IF (icty.GE.1) rho(i) = rhoin(i)
if (imhd.lt.0) then
call cross_product3D(vel(:,i),Bconst(:),vcrossB)
Bevol(:,i) = Bevolin(:,i) + dt*vcrossB(:)
else
Bevol(:,i) = Bevolin(:,i)
endif
IF (iener.NE.0) en(i) = enin(i)
x(:,i) = xin(:,i) + dt*(vel(1:ndim,i) + xsphfac*xsphterm(1:ndim,i))
alpha(:,i) = alphain(:,i)
hh(i) = hhin(i)
psi(i) = psiin(i)
ELSE
vel(:,i) = (velin(:,i) + dt*force(:,i))/(1.+damp)
!--vertical damping
if (ndim.ge.3) vel(3,i) = vel(3,i)/(1.+dampz)
!--radial damping
if (dampr.gt.tiny(dampr) .and.ndim.ge.2) call radialdamping(vel(1:2,i),x(1:2,i),dampr)
x(:,i) = xin(:,i) + dt*(vel(1:ndim,i) + xsphfac*xsphterm(1:ndim,i))
IF (ihvar.EQ.2 .OR. (ihvar.EQ.3 .and. itsdensity.eq.0)) THEN
hh(i) = hhin(i) + dt*dhdt(i)
IF (hh(i).LE.0.) THEN
WRITE(iprint,*) 'step: hh -ve ',i,hh(i)
CALL quit
ENDIF
ENDIF
IF (icty.GE.1) THEN
rho(i) = rhoin(i) + dt*drhodt(i)
IF (rho(i).LE.0.) THEN
WRITE(iprint,*) 'step: rho -ve ',i,rho(i)
CALL quit
ENDIF
ENDIF
IF (iener.NE.0) en(i) = enin(i) + dt*dendt(i)
IF (ANY(iavlim.NE.0)) alpha(:,i) = alphain(:,i) + dt*daldt(:,i)
IF (imhd.NE.0) THEN
IF (nsubsteps_divB.LT.0) THEN
Bevol(:,i) = Bevolin(:,i) + dt*gradpsi(:,i)
psi(i) = (psiin(i) + dt*dpsidt(i))/(1.+eps)
ENDIF
Bevol(:,i) = Bevol(:,i) + dt*dBevoldt(:,i)
ENDIF
ENDIF
ENDDO
!
!--if doing divergence correction then do correction to magnetic field
!
! IF (idivBzero.EQ.1 .and. MOD(nsteps,10.EQ.0) .AND. ALL(ibound.EQ.0)) CALL divBcorrect
!
IF (ANY(ibound.NE.0)) CALL boundary ! inflow/outflow/periodic boundary conditions
!
!--set new timestep from courant/forces condition
!
dt = min(C_force*dtforce,C_cour*dtcourant,C_force*dtdrag,C_force*dtvisc)
!
!--this is just so that alpha looks right in the output
! (overwritten in predictor step, but shows where alpha is 1 in rates)
!
! WHERE (alpha(:,:).GT.0.5)
! alpha(:,:) = 1.0
! END WHERE
IF (trace) WRITE (iprint,*) ' Exiting subroutine step'
RETURN
END SUBROUTINE step
subroutine radialdamping(vxy,xy,dampr)
implicit none
real, intent(in), dimension(2) :: xy
real, intent(inout), dimension(2) :: vxy
real, intent(in) :: dampr
real, dimension(2) :: rhat
real :: vrad,vradnew
rhat(:) = xy(:)/sqrt(dot_product(xy,xy))
vrad = dot_product(vxy,rhat)
vradnew = vrad/(1. + dampr)
vxy(:) = vxy(:) + (vradnew - vrad)*rhat(:)
return
end subroutine radialdamping