-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathpostshock.f90
2056 lines (1578 loc) · 58.5 KB
/
postshock.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
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!$Id: postshock.f90,v 1.11 2007-08-14 12:05:58 bgiacoma Exp $
!!$ Copyright (C) 2005 B. Giacomazzo, L. Rezzolla
subroutine postshock(Vs,Byb,Bzb,vxb,ahead,behind)
use type
use global
use wavecheck !contains wave_error
use eos_param
implicit none
real(DP),intent(IN)::Vs,Byb,Bzb,vxb
real(DP),intent(IN),dimension(7)::ahead !value of primitives ahead the shock
real(DP),intent(OUT),dimension(7)::behind
!Given the values of the post-shock By, Bz and v^x and the value of the shock velocity Vs,
!returns the values of the remaining postshock quantities using the Bt-method.
real(DP)::vyb,vzb,Db,pb,Wb,v2b,Vb,Ea,normBa
real(DP)::rhoa,pa,vxa,vya,vza,Bya,Bza,etaa,b2a,ha,v2a,Wa,Ws,j,Da,taua,Va
real(DP)::rhob,etab,b2b,F,B2
wave_error=.false.
!values ahead the shock
rhoa = ahead(1) !mass density
pa = ahead(2) !Totale pressure (gas pressure+ magnetic pressure)
vxa = ahead(3) !velocity
vya = ahead(4)
vza = ahead(5)
Bya = ahead(6) !magnetic field
Bza = ahead(7)
normBa=sqrt(Bya**2+Bza**2) !Bt
etaa=Bx*vxa+Bya*vya+Bza*vza !B^j v_j
v2a=vxa**2+vya**2+vza**2 !v^2
if (v2a>=1.0e0_dp) stop 'error in postshock.f90: v^2a>c'
Wa=1.0e0_dp/sqrt(1.0e0_dp-v2a) !Lorentz factor
b2a=(Bx**2+Bya**2+Bza**2)/Wa**2+etaa**2 !(b_\mu b^\mu)
!EoS
call eos_enthalpy(pa-0.5e0_dp*b2a,rhoa,gamma,ha)
!ha=1.0e0_dp+gamma*(pa-0.5e0_dp*b2a)/(rhoa*(gamma-1.0e0_dp)) !specifc gas enthalpy
taua=(rhoa*ha+b2a)*Wa**2-pa-Da
Da=rhoa*Wa
Ea=taua-Wa**2*etaa**2
Va=1/Da
Ws=1.0e0_dp/sqrt(1-Vs**2) !Lorents factor of the shock wave
j=Ws*rhoa*Wa*(Vs-vxa) !J
F=j/Ws
!-----------------------------------------------------!
!--------------COMPUTE THE POST-SHOCK VALUES----------!
!-----------------------------------------------------!
B2=Bx**2+Byb**2+Bzb**2 !B^2
Vb = Va + (vxa - vxb)/F
if (Bx==0.0e0_dp) then
stop 'riemann.f90: postshock: no solution for Bx=0 within the Bt-method'
else
vyb=(Bya*F*Va - Byb*(F*Va + vxa - vxb) + Bx*vya)/Bx
vzb=(Bza*F*Va - Bzb*(F*Va + vxa - vxb) + Bx*vza)/Bx
end if
if (((vxb**2+vyb**2+vzb**2)>1.0e0_dp).OR.(abs(Vs)>1.0e0_dp)) then
wave_error=.true.
return
end if
Wb=1.0e0_dp/sqrt(1.0e0_dp-vxb**2-vyb**2-vzb**2)
Db=1.0e0_dp/Vb
rhob=Db*sqrt(1.0e0_dp-vxb**2-vyb**2-vzb**2)
etab=vxb*Bx+vyb*Byb+vzb*Bzb
b2b=B2/Wb**2+etab**2
!This equation is obtained using the invariance of (iii) (see Anile, page 260)
if (eos_ideal) then
pb = (b2b*gamma*rhoa*Wa*(etab*j*Wb - Bx*rhob*Ws) + &
2*(-1 + gamma)*rhob*(-(etab*j*rhoa*Wa*Wb) + &
etaa*ha*j*rhob*Wa*Wb + Bx*rhoa*rhob*(Wa - ha*Wb)*Ws))/ &
(2.*gamma*rhoa*Wa*(etab*j*Wb - Bx*rhob*Ws))
elseif (eos_meliani) then
pb = (b2b*rhoa*Wa*(etab*j*Wb - Bx*rhob*Ws)* &
(j*(2*etaa*ha + etab*(-1 + gamma)*rhoa)*Wa*Wb - &
Bx*rhoa*((-1 + gamma)*rhob*Wa + 2*ha*Wb)*Ws) - &
2*(-1 + gamma)*(j*(etab*rhoa + etaa*ha*rhob)*Wa*Wb - &
Bx*rhoa*rhob*(Wa + ha*Wb)*Ws)* &
(etab*j*rhoa*Wa*Wb - &
rhob*(etaa*ha*j*Wa*Wb + Bx*rhoa*(Wa - ha*Wb)*Ws)))/ &
(2.*rhoa*Wa*(etab*j*Wb - Bx*rhob*Ws)* &
(j*(2*etaa*ha + etab*(-1 + gamma)*rhoa)*Wa*Wb - &
Bx*rhoa*((-1 + gamma)*rhob*Wa + 2*ha*Wb)*Ws))
else
write(*,*) 'postshock.f90:postshock: Error in the EOS'
STOP
end if
behind(1)=rhob
behind(2)=pb
behind(3)=vxb
behind(4)=vyb
behind(5)=vzb
behind(6)=Byb
behind(7)=Bzb
end subroutine postshock
!---------------------------------------------------------------!
!-------------------------VELOCITY------------------------------!
!---------------------------------------------------------------!
subroutine velocity(unk1,unk2,ahead,switchLR,vx,Vs,switchPB,behind)
use type
use global
use accmod !contains accuracy
use wavecheck
use interfaces,only:velocity_eqn,velocity_df,postshock
implicit none
real(DP),intent(IN)::unk1,unk2
real(DP),dimension(7),intent(IN)::ahead
character(len=2),intent(IN)::switchLR
real(DP),intent(OUT)::vx,Vs
character(len=1),intent(IN)::switchPB
real(DP),dimension(7),intent(OUT),OPTIONAL::behind
!given the post-shock value of By and Bz it computes
! the postshock v^x and the shock velocity Vs
!(It's used only in the Bt-method)
real(DP)::vx1,vx2,f1,f2,vxl,vxh,dvx,dvxold,Vs1,Vs2,d,df,f,temp,ftol,step,vxinit
real(DP)::xmin,xmax,step_min,cf,vx_max,vx_min
integer(I4B)::j
logical::upperbound,lowerbound
real(DP),dimension(7)::plotstate
!This routine is needed only if you use the norm of Bt
if (switchPB=='P') stop 'postshock.f90: velocity: it has no sense to use this routine with the p-method'
cf=1.0e0_dp !relaxation factor
ftol=1.0e-10_dp !accuracy used in this subroutine
if (accuracy<ftol) ftol=accuracy !to be consistent with the solution of the equations at the CD
!It will search for the value of vx in the interval [xmin,xmax]
vx=ahead(3)
xmin=-0.9e0_dp
xmax=0.9e0_dp
select case(initial_data)
case (13)
!For the generic Alfven test (from HLLE)
if (switchLR=='LS') then
vx=0.03885e0_dp
xmin=0.0388e0_dp
xmax=0.0395e0_dp
else if (switchLR=='RS') then
vx=0.03885e0_dp
xmin=0.0388e0_dp
xmax=0.0395e0_dp
end if
case (7)
!For KO: Collision (from HLLE)
if (switchLR=='LS') then
vx=0.0e0_dp
xmin=-1.0e-3_dp
xmax=1.0e-3_dp
else if (switchLR=='RS') then
vx=0.0e0_dp
xmin=-1.0e-3_dp
xmax=1.0e-3_dp
end if
case (8)
!For Balsara 1 (from HLLE)
if (switchLR=='LS') then
vx=0.2555e0_dp
xmin=0.25e0_dp
xmax=0.26e0_dp
else if (switchLR=='RS') then
vx=0.25545e0_dp
xmin=0.25540e0_dp
xmax=0.25546e0_dp
end if
case (9)
!For Balsara 2
if (switchLR=='RS') then
vx=0.676980015e0_dp
xmin=0.675
xmax=0.677
end if
case (10)
!For Balsara 3 (from HLLE)
if (switchLR=='RS') then
vx=0.953e0_dp
xmin=0.95e0_dp
xmax=0.954e0_dp
end if
case(11)
!For Balsara 4 (from HLLE)
if (switchLR=='LS') then
vx=0.0e0_dp
xmin=-1.0e-4_dp
xmax=5.0e-4_dp
else if (switchLR=='RS') then
vx=0.0e0_dp
xmin=-2.0e-4_dp
xmax=5.0e-4_dp
end if
case (12)
!For Balsara 5 (from HLLE)
if (switchLR=='RS') then
vx=-0.454303238E-01
xmin=-0.0455
xmax=-0.0454
end if
case default
!change these values if you are using your own initial condition
! (and if Bx is different from zero)
if (switchLR=='LS') then
vx=0.0e0_dp
xmin=-0.99e0_dp
xmax=0.99e0_dp
else if (switchLR=='RS') then
vx=0.0e0_dp
xmin=-0.99e0_dp
xmax=0.99e0_dp
end if
end select
vxinit=vx
!step used to search the right interval in which the solution is contained
! (i.e. it must be f(xmin)*f(xmax)<0, and the function f must be defined in that interval)
step=1.0e-1_dp*(xmax-xmin)
step_min=1.0e-5*step
!-------------------------------------------------------------------------------------!
!-----------------------ONLY FOR DEBUG------------------------------------------------!
!---------------------------PLOT------------------------------------------------------!
!-------------------------------------------------------------------------------------!
!This plot the function whose roots give the solution for the velocity
!Try this to get an initial guess about the interval in which the root is present
!switchLR can be set equal to RS or LS
if ((switchLR=="noplot").AND.(niter>-1)) then
open(UNIT=1234,FILE='vel_func.dat',STATUS='REPLACE',ACTION='WRITE')
print *
print *,'vmax = ',xmax
print *,'vmin = ',xmin
print *,'switchLR = ',switchLR
wave_error=.false.
vx=xmin
write(1234,*)
do while (vx<=xmax)
call velocity_eqn(unk1,unk2,ahead,vx,switchLR,f1,Vs1,plotstate)
print *,vx,wave_error
!if wave_error==true then f is not defined for that value
if (.NOT.(wave_error)) write(1234,FMT='(13E30.18)') vx,f1,df,Vs1,unk1,unk2,plotstate
wave_error=.false.
vx=vx+1.0D-4*(xmax-xmin)
end do
close (1234)
stop 'function plotted in vel_func.dat'
end if
!-------------------------------------------------------------------------------------!
!-------------------------------------------------------------------------------------!
!-------------------------------------------------------------------------------------!
!-----------------------!
!Bracketing the solution!
!-----------------------!
2 if(upperbound .and. lowerbound) then
step=step/2.0e0_dp
if (step<step_min) then !STEP TOO SMALL
wave_error=.true.
stop 'error in velocity'
!return
end if
if (verbose) then
print *,'STEP=',step,'switchLR=',switchLR,'vx=',vx
print *
end if
lowerbound=.false.
upperbound=.false.
vx=vxinit
end if
call velocity_eqn(unk1,unk2,ahead,vx,switchLR,f2,Vs)
if ((abs(f2)<ftol).AND.(.NOT.(wave_error))) goto 1000
if (wave_error) then
vx=vx-step !Go to the left
if (vx<xmin) then
vx=xmax-(xmin-vx)
if (lowerbound) then
!reduce the step
!It cannot find a valid interval
upperbound=.true.
goto 2
end if
lowerbound=.true.
end if
goto 2
end if
vx2=vx !In this point the function is defined
lowerbound=.false.
upperbound=.false.
if (verbose) then
print *
print *,'vx2=',vx2,'f2=',f2
end if
!Now start from vx2 and go to the left, if the function changes its sign
!then we have found the right interval
!else we have to reduce the step
do
21 vx=vx-step !Go to the left
if (vx<xmin) then
wave_error=.true.
goto 22
end if
call velocity_eqn(unk1,unk2,ahead,vx,switchLR,f1,Vs)
if ((abs(f1)<ftol).AND.(.NOT.(wave_error))) goto 1000
22 if (wave_error) then !f is not defined go to the right
vx1=vx
vx=vx2
4 vx=vx+step
if (vx>xmax) then
wave_error=.true.
goto 5
end if
call velocity_eqn(unk1,unk2,ahead,vx,switchLR,f1,Vs)
if ((abs(f1)<ftol).AND.(.NOT.(wave_error))) goto 1000
5 if (wave_error) then
step=step/2.0e0_dp
if (step<step_min) then !STEP TOO SMALL
wave_error=.true.
stop 'error in velocity'
return
end if
if (verbose) then
print *,'reducing the step and tring again...'
print *,'STEP=',step
print *
end if
vx=vx2
goto 21
end if
if ((f1*f2)<0.0e0_dp) then
vx1=vx
goto 10
else
goto 4
end if
end if
if ((f1*f2)<0.0e0_dp) then
vx1=vx
goto 10
end if
end do
10 if (verbose) then
print *
print *,'Byb=',unk1,'Bzb=',unk2
print *,'I will try to find the solution in the following interval:'
print *,'vx_min=',vx1,' f1=',f1
print *,'vx_max=',vx2,' f2=',f2
print *
end if
!Find the post-shock velocity
!Newton-Raphson and bisection method
if (abs(f1)<=ftol) then
if (veryverbose) print *,'postshock.f90: pressure: f1 is zero'
vx=vx1
Vs=Vs1
if (present(behind)) call postshock(Vs,unk1,unk2,vx,ahead,behind)
return
end if
if (abs(f2)<=ftol) then
if (veryverbose) print *,'postshock.f90: pressure: f2 is zero'
vx=vx2
Vs=Vs2
if (present(behind)) call postshock(Vs,unk1,unk2,vx,ahead,behind)
return
end if
if (f1*f2>0.0e0_dp) stop 'f1 and f2 have the same sign; there is no solution!'
if (f1 < 0.0) then !Orient the search so that f(Vsl)<0.
vxl=vx1
vxh=vx2
else
vxh=vx1
vxl=vx2
end if
vx=0.5e0_dp*(vx1+vx2) !Initialize the guess for root
dvxold=abs(vx2-vx1) !the "stepsize before last",
dvx=dvxold !and the last step.
call velocity_df(unk1,unk2,ahead,vx,switchLR,f,df,Vs)
do j=1,200 !200 is the maximum number of iterations
if (((vx-vxh)*df-f)*((vx-vxl)*df-f) > 0.0 .or. &
abs(2.0_DP*f) > abs(dvxold*df) ) then
!Bisect if Newton out of range, or not decreasing fast enough.
dvxold=dvx
dvx=0.5_DP*(vxh-vxl)
vx=vxl+dvx
if (vxl == vx) goto 1000 !Change in root is negligible.
else !Newton step acceptable. Take it.
dvxold=dvx
dvx=f/df
temp=vx
vx=vx-dvx
if (temp == vx) goto 1000
end if
if (abs(dvx) < ftol) goto 1000 !Convergence criterion.
call velocity_df(unk1,unk2,ahead,vx,switchLR,f,df,Vs)
!One new function evaluation per iteration.
if (f < 0.0) then !Maintain the bracket on the root.
vxl=vx
else
vxh=vx
end if
end do
stop 'I cannot find post-shock velocity'
999 call velocity_eqn(unk1,unk2,ahead,vx,switchLR,f1,Vs)
if (wave_error) stop 'postshock.f90: velocity: no solution for the value vx used'
1000 if (present(behind)) call postshock(Vs,unk1,unk2,vx,ahead,behind)
end subroutine velocity
!---------------------------------------------------------------!
!-------------------------VELOCITY_DF---------------------------!
!---------------------------------------------------------------!
subroutine velocity_df(unk1,unk2,ahead,vx,switchLR,f,df,Vs)
use type
use global
use wavecheck
use interfaces,only:velocity_eqn
implicit none
real(DP),intent(IN)::unk1,unk2,vx
real(DP),dimension(7),intent(IN)::ahead
character(len=2),intent(IN)::switchLR
real(DP),intent(OUT)::f,df
real(DP),intent(OUT),OPTIONAL::Vs
!Compute the first derivative of f(vx)
real(DP)::delta,f2,dump
if (veryverbose) then
print *
print *,'postshock.f90: velocity_df: I''m computing the derivative of f(vx)'
print *
end if
call velocity_eqn(unk1,unk2,ahead,vx,switchLR,f,Vs)
if (wave_error) return
!First derivative
delta=1.0D-8*vx
if (delta<=1.0D-10) delta=1.0D-10
call velocity_eqn(unk1,unk2,ahead,vx-delta,switchLR,f2,dump)
if (wave_error) then
print *,'vx=',vx,'delta=',delta,'switchLR=',switchLR
stop 'postshock.f90: velocity_df: error computing the derivative'
end if
df=(f-f2)/delta
end subroutine velocity_df
!---------------------------------------------------------------!
!-------------------------VELOCITY_EQN--------------------------!
!---------------------------------------------------------------!
subroutine velocity_eqn(unk1,unk2,ahead,vxb,switchLR,f,Vs,stateb)
use type
use global
use wavecheck!This module contains wave_error
use interfaces,only:ContactVelocity,shockvelocity,postshock
implicit none
real(DP),intent(IN)::unk1,unk2
real(DP),dimension(7),intent(IN)::ahead
real(DP),intent(IN)::vxb
character(len=2),intent(IN)::switchLR
real(DP),intent(OUT)::f
real(DP),intent(OUT)::Vs
real(DP),dimension(7),OPTIONAL,intent(OUT)::stateb
!The zero of the equation (4.16) gives the right value for vxb
!Used only within the Bt-method
real(DP)::rhoa,pa,vxa,vya,vza,Bya,Bza,etaa,v2a,Wa,b2a,ha,taua,Da,Ws,j
real(DP)::rhob,pb,vyb,vzb,Byb,Bzb,etab,v2b,Wb,b2b,hb,taub,Db
real(DP)::Ea,Eb,dump,Va,Vb,normBa,B2,bigF
real(DP),dimension(7)::behind
Byb=unk1
Bzb=unk2
rhoa = ahead(1) !mass density
pa = ahead(2) !Totale pressure (gas pressure+ magnetic pressure)
vxa = ahead(3) !velocity
vya = ahead(4)
vza = ahead(5)
Bya = ahead(6) !magnetic field
Bza = ahead(7)
if (abs(vxa-vxb)<1.0e-15_dp) then
!It's not a shock
wave_error=.true.
return
end if
normBa=sqrt(Bya**2+Bza**2) !Bt ahead the shock
etaa=Bx*vxa+Bya*vya+Bza*vza !B^j v_j
v2a=vxa**2+vya**2+vza**2 !v^2
if (v2a>=1.0e0_dp) stop 'error in postshock.f90: pressure_eqn: v^2a>c'
Wa=1.0e0_dp/sqrt(1.0e0_dp-v2a) !Lorentz factor
b2a=(Bx**2+Bya**2+Bza**2)/Wa**2+etaa**2 !(b_\mu b^\mu)
!EoS
call eos_enthalpy(pa-0.5e0_dp*b2a,rhoa,gamma,ha)
!ha=1.0e0_dp+gamma*(pa-0.5e0_dp*b2a)/(rhoa*(gamma-1.0e0_dp)) !specific gas enthalpy
Da=rhoa*Wa
Va=1.0e0_dp/Da
taua=(rhoa*ha+b2a)*Wa**2-pa-Da
Ea=taua-Wa**2*etaa**2
!Compute the shock velocity from normB, vxb and the postshock values
wave_error=.false.
call ContactVelocity(Byb,Bzb,vxb,Vs,ahead,dump,switchLR,'B',behind)
if (wave_error) then
if (veryverbose) print *,'postshock.f90: velocity_eqn: wave_error'
return
end if
Ws=1.0e0_dp/sqrt(1-Vs**2) !Lorentz factor for the shock
j=Ws*rhoa*Wa*(Vs-vxa) !J
bigF=j/Ws
rhob = behind(1) !mass density
pb = behind(2) !total pressure
vyb = behind(4) !tangential components of velocity
vzb = behind(5)
Byb = behind(6) !magnetic field
Bzb = behind(7)
if (present(stateb)) stateb=behind
B2=Bx**2+Byb**2+Bzb**2
etab=Bx*vxb+Byb*vyb+Bzb*vzb
v2b=vxb**2+vyb**2+vzb**2
if (v2b>=1.0e0_dp) then
if (verbose) print *,'velocity_eqn: v2b=',v2b
wave_error=.true.
return
end if
Wb=1.0e0_dp/sqrt(1.0e0_dp-v2b)
!Compare j before and after the shock (it's only a check)
if (veryverbose) print *,'ja=',j,'jb=',Ws*rhob*Wb*(Vs-vxb)
b2b=(Bx**2+Byb**2+Bzb**2)/Wb**2+etab**2 !this depends on normB
!EoS
call eos_enthalpy(pb-0.5e0_dp*b2b,rhob,gamma,hb)
!hb=1.0e0_dp+gamma*(pb-0.5e0_dp*b2b)/(rhob*(gamma-1.0e0_dp))
Db=rhob*Wb
Vb=1.0e0_dp/Db
taub=(rhob*hb+b2b)*Wb**2-pb-Db
Eb=taub-Wb**2*etab**2
!equation (4.16)
f = bigF*(((Ea+pa+Da)*vya-etaa*Bya)/Da-((Eb+pb+Db)*vyb-etab*Byb)/Db) + Bx*(etaa*vya-etab*vyb)+ &
Bx*(Bya/Wa**2-Byb/Wb**2)
if (.NOT.(abs(f)>=0.0e0_dp)) then
if (verbose) print *,'velocity_eqn: f=',f
wave_error=.true.
return
end if
end subroutine velocity_eqn
!---------------------------------------------------------------!
!---------------------CONTACT VELOCITY--------------------------!
!---------------------------------------------------------------!
subroutine ContactVelocity(unk1,unk2,vxb,Vs,ahead,vx,switchLR,switchPB,solution)
use type
use accmod !This module contains accuracy
use global !This module contains Bx and gamma
use wavecheck!This module contains wave_error
use odeswitch!This module contains switch
use interfaces,only: zeroeqn,postshock,eos_cs2,eos_enthalpy,postshock_pressure
implicit none
real(DP),intent(IN)::unk1,unk2,vxb
real(DP),intent(OUT)::Vs
real(DP),dimension(7),intent(IN)::ahead
real(DP),intent(OUT)::vx
character(len=2),intent(IN)::switchLR
character(len=1),intent(IN)::switchPB
real(DP),dimension(7),intent(OUT),OPTIONAL::solution
!Given vxb this subroutine returns Vs and
!the value of the postshock vx
!It can be used with both the Bt-method (switchPB=B) and the p-method (switchPB=P)
integer(I4B)::j
real(DP)::dVs,f,df,Vsinit,temp,Vs1,Vs2,f1,f2
real(DP),dimension(7)::behind
logical::errcheck,secondtry,lowerbound,upperbound
real(DP)::Vsl,Vsh,dVsold,step
real(DP)::b2,pgas,wtot,leftAlfven,rightAlfven,cdvel,xmax,xmin,B2big,W2inv,eta,enthalpy
real(DP)::step_min,ftol,cf !Minimum step allowed, accuracy and relaxation factor
ftol=1.0e-10_dp !accuracy used in this subroutine
if (accuracy<ftol) ftol=accuracy !to be consistent with the solution of the equations at the CD
cf=0.5e0_dp
wave_error=.false.
switch=switchLR
!Compute Alfven velocities (left and right-going)
B2big=Bx**2+ahead(6)**2+ahead(7)**2 !B^2
W2inv=1.0e0_dp-ahead(3)**2-ahead(4)**2-ahead(5)**2
eta=Bx*ahead(3)+ahead(6)*ahead(4)+ahead(7)*ahead(5) !B_j v^j
b2=B2big*W2inv+eta**2 !b_mu b^\mu
pgas=ahead(2)-0.5e0_dp*b2 !gas pressure
call eos_enthalpy(pgas,ahead(1),gamma,enthalpy)
wtot=ahead(1)*enthalpy + b2 !total relativistic enthalpy
leftAlfven =ahead(3)+Bx*W2inv/(eta-sqrt(wtot)) !left-going Alfven velocity
rightAlfven=ahead(3)+Bx*W2inv/(eta+sqrt(wtot)) !right-going Alfven velocity
!vx
if ((switchLR=='LF').OR.(switchLR=='LS'))then
cdvel=min(ahead(3),vxb) !because J should be always < 0
if (switchPB=='P') cdvel=ahead(3)
else
cdvel=max(ahead(3),vxb) !because J should be always > 0
if (switchPB=='P') cdvel=ahead(3)
end if
if (switchPB=='B') then
if(abs(vxb-ahead(3))<epsilon(vxb)) then !No change in v^x (there is no discontinuity)
if (present(solution)) solution=ahead
vx=ahead(3)
if ((switchLR=='LS').OR.(switchLR=='RS')) then
Vs=vx
else if (switchLR=='LF') then
Vs=leftAlfven
else if (switchLR=='RF') then
Vs=rightAlfven
else
stop 'riemann.f90: ContactVelocity: error in switchLR'
end if
return
end if
end if
if (switchLR=='LF') then
xmax=leftAlfven
xmin=-0.9999999e0_dp
!if (initial_data==0) xmin=0.99 !for user defined Riemann problem (if needed)
else if (switchLR=='LS') then
xmax=cdvel
xmin=leftAlfven
!This is due to the presence of more the one root in the interval (xmin,xmax)
if (initial_data==13) xmin=-0.2e0_dp !for test alfven 1
if (initial_data==8) xmin=0.02e0_dp !for Balsara 1
if (initial_data==11) xmin=-0.150e0_dp !for Balsara 4
!if (initial_data==0) xmin=??? !for user defined Riemann problem (if needed)
else if (switchLR=='RS') then
xmax=rightAlfven
xmin=cdvel
!This is due to the presence of more the one root in the interval (xmin,xmax)
if (initial_data==9) xmin=0.8e0_dp !for Balsara 2
if (initial_data==8) xmin=0.38e0_dp !for Balsara 1
if (initial_data==11) xmax=0.16e0_dp !for Balsara 4
if (initial_data==12) xmin=0.36e0_dp !for Balsara 5
!if (initial_data==10) xmax=0.995e0_dp !for Balsara 4 !Only for Meliani EOS
!if (initial_data==0) xmin=??? !for user defined Riemann problem (if needed)
else if (switchLR=='RF') then
xmax=0.9999999e0_dp
xmin=rightAlfven
!if (initial_data==0) xmin=0.999999 !for user defined Riemann problem (if needed)
!This is due to the presence of more than one root in the interval (xmin,xmax)
if (initial_data==10) xmin=0.994e0_dp !for Balsara 3
else
stop 'riemann.f90: ContactVelocity: error in switchLR'
end if
call zeroeqn(vxb,unk1,unk2,xmax,ahead,f,df,switchLR,switchPB,errcheck)
if (veryverbose) print *,'ContactVelocity: f(xmax)=',f,'xmax=',xmax,'errcheck=',errcheck
if ((abs(f)<=ftol).AND.(.NOT.(errcheck))) then
Vs=xmax
write(*,*) vxb,unk1,unk2,xmax,ahead,f,df,switchLR,switchPB,errcheck
stop 'postshock.f90: contact_velocity: Vs=vmax; this is strange!'
if (veryverbose) print *,'Vs=',Vs
return
end if
call zeroeqn(vxb,unk1,unk2,xmin,ahead,f,df,switchLR,switchPB,errcheck)
if (veryverbose) print *,'ContactVelocity: f(xmin)=',f ,'xmin=',xmin,'errcheck=',errcheck
if ((abs(f)<=ftol).AND.(.NOT.(errcheck))) then
Vs=xmin
write(*,*) vxb,unk1,unk2,xmin,ahead,f,df,switchLR,switchPB,errcheck
stop 'postshock.f90: contact_velocity: Vs=vmin; this is strange!'
if (veryverbose) print *,'Vs=',Vs
return
end if
secondtry=.false.
Vs=0.5e0_dp*(xmax+xmin)!First guess
Vsinit=Vs
step=1.0D-3*abs(xmax-xmin) !Initial step
step_min=1.0e-4_dp*step !Minimum step allowed
upperbound=.false.
lowerbound=.false.
!-----------------------ONLY FOR DEBUG-----------------------------!
!--------------------------PLOT------------------------------------!
!This plots the function whose roots give the value of the shock velocity
!This can be used to have an idea of what is going wrong if the routine doesn't work and
! also to check for the presence of more than one root in the interval
!substitute noplot with one of the following: RF, RS, LF, LS
! (respectively right-fast, right-slow, left-fast, left-slow shock)
if (switchLR=='noplot') then
open(UNIT=123,FILE='funcd.dat',STATUS='REPLACE',ACTION='WRITE')
print *
print *,'cdvel =',cdvel
print *,'xmax =',xmax
print *,'xmin =',xmin
!!$ xmin=-0.99999e0_dp
!!$ xmax=0.99999e0_dp
errcheck=.false.
Vs=xmin
write(123,*)
do while (Vs<=xmax)
call zeroeqn(vxb,unk1,unk2,Vs,ahead,f,df,switchLR,switchPB,errcheck)
!if errcheck==true then f is not defined for that value
if (.NOT.(errcheck)) write(123,FMT='(3E30.18)') Vs,f,vxb
!!$ write(*,*) Vs,f
errcheck=.false.
Vs=Vs+1.0D-4*(xmax-xmin)
end do
close(123)
stop 'function plotted in funcd.dat'
end if
!------------------------------------------------------------------!
!------------------------------------------------------------------!
Vs=Vsinit
!-----------------------!
!Bracketing the solution!
!-----------------------!
2 if(upperbound .and. lowerbound) then
step=step/2.0e0_dp
if (step<step_min) then !STEP TOO SMALL
wave_error=.true.
!rarefaction_pmin=ahead(2)
shock=.true.
return
end if
if (veryverbose) then
print *,'STEP=',step,'switchLR=',switchLR,'vxb=',vxb
print *
end if
lowerbound=.false.
upperbound=.false.
Vs=Vsinit
end if
wave_error=.false.
call zeroeqn(vxb,unk1,unk2,Vs,ahead,f2,df,switchLR,switchPB,errcheck)
if ((abs(f2)<ftol).AND.(.NOT.(errcheck))) goto 1000
if (errcheck) then
Vs=Vs-step !Go to the left
if (Vs<xmin) then
Vs=xmax-(xmin-Vs)
if (lowerbound) then
!reduce the step
!It cannot find a valid interval
upperbound=.true.
goto 2
end if
lowerbound=.true.
end if
goto 2
end if
Vs2=Vs !In this point the function is defined
lowerbound=.false.
upperbound=.false.
if (veryverbose) then
print *
print *,'Vs2=',Vs2,'f2=',f2,'vxb=',vxb
end if
!Now start from Vs2 and go to the left, if the function changes its sign
!then we have found the right interval
!else we have to change interval
do
Vs=Vs-step !Go to the left
if (Vs<xmin) then
errcheck=.true.
goto 22
end if
call zeroeqn(vxb,unk1,unk2,Vs,ahead,f1,df,switchLR,switchPB,errcheck)
if ((abs(f1)<ftol).AND.(.NOT.(errcheck))) goto 1000
22 if (errcheck) then !f is not defined go to the right
Vs1=Vs
Vs=Vs2
4 Vs=Vs+step
if (Vs>xmax) then
errcheck=.true.
goto 5
end if
call zeroeqn(vxb,unk1,unk2,Vs,ahead,f1,df,switchLR,switchPB,errcheck)
if ((abs(f1)<ftol).AND.(.NOT.(errcheck))) goto 1000
5 if (errcheck) then !wrong interval. chose a different starting point
step=step/2.0e0_dp
if (step<step_min) then !STEP TOO SMALL
wave_error=.true.
!rarefaction_pmin=ahead(2)
shock=.true.
return
end if
if (veryverbose) then
print *,'wrong interval chose a different starting point'
print *,'STEP=',step
print *
end if
Vs=Vs1-step
if (Vs<xmin) Vs=xmax-(xmin-Vs)
goto 2
end if
if ((f1*f2)<0.0e0_dp) then
Vs1=Vs
goto 10
else
goto 4
end if
end if
if ((f1*f2)<0.0e0_dp) then
Vs1=Vs
goto 10
end if
end do
!Find the shock velocity
!Newton-Raphson and bisection method
10 if (veryverbose) then
print *,'Finding Root for'
print *,'vxb=',vxb
print *,'in the following interval'
print *,'Vs1=',Vs1,'f1=',f1
print *,'Vs2=',Vs2,'f2=',f2
end if
if (f1==0.0e0_dp) then