-
Notifications
You must be signed in to change notification settings - Fork 0
/
GasDensity.f90
688 lines (613 loc) · 19.8 KB
/
GasDensity.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
module projections
DOUBLE PRECISION,PARAMETER ::pi_n = atan(1.d0)*4.d0
DOUBLE PRECISION,SAVE ::aline,argaline
DOUBLE PRECISION,SAVE ::apitc
DOUBLE PRECISION,SAVE ::aolin
CONTAINS
FUNCTION R(th)
IMPLICIT NONE
DOUBLE PRECISION ::R(2,2)
DOUBLE PRECISION ::th
R(1,1) = cos(th)
R(1,2) = -sin(th)
R(2,1) = sin(th)
R(2,2) = cos(th)
ENDFUNCTION
FUNCTION C(th)
IMPLICIT NONE
DOUBLE PRECISION ::C(2,2)
DOUBLE PRECISION ::th
C = 0.d0
C(1,1) = cos(th)
C(2,2) = 1.d0
ENDFUNCTION
FUNCTION INVC(th)
IMPLICIT NONE
DOUBLE PRECISION ::INVC(2,2)
DOUBLE PRECISION ::th
INVC = 0.d0
INVC(1,1) = 1.d0/cos(th)
INVC(2,2) = 1.d0
ENDFUNCTION
endmodule
MODULE argument
IMPLICIT NONE
DOUBLE PRECISION ::zmax = 30.d0
CHARACTER(len=225) ::filename
INTEGER ::frames(2)
INTEGER ::drawmode = 1 !1 for draw static
!density, 2 for movie,
LOGICAL ::toproject = .false.
LOGICAL ::rauto = .true.
LOGICAL ::drawcir = .false.
LOGICAL ::drawstellar=.false.
LOGICAL ::givefnm = .false.
LOGICAL ::drawq = .false.
LOGICAL ::save = .false.
LOGICAL ::cut = .false.
LOGICAL ::velocity = .false.
CONTAINS
SUBROUTINE readarg
USE projections,only:argaline
IMPLICIT NONE
CHARACTER(len=32) ::arg
INTEGER ::i
INTEGER ::ioerr
!read in rotating parameters and gas filename by stdin
if(iargc().ne.0)then
DO i = 1, iargc()
CALL getarg(i,arg)
SELECT CASE(arg)
CASE('--help','-h')
write(6,'(a)')'Draw gas density distribution. '
write(6,'(a)')'usage: GasDensity.exe [option] '
write(6,'(a)')'options: '
write(6,'(a)')' -c, --circle To draw circles. '
write(6,'(a)')' -i [file name.h5] Read file.'
write(6,'(a)')' -m [start] [end] Read files.'
write(6,'(a)')' -p, --project [degree] To project the density.'
write(6,'(a)')' -r, --zauto Automatic find scale of z'
write(6,'(a)')' -z, --zmax [scale of z] Specified the scale of z'
write(6,'(a)')' -h, --help Show this help page '
write(6,'(a)')' -s, --stellar Draw stellar contour. '
write(6,'(a)')' -q, Draw instability.'
write(6,'(a)')' -v, Version information.'
write(6,'(a)')' --save Save as png file. '
write(6,'(a)')' --cut Make a cut at CO. '
write(6,'(a)')' --velocity Draw velocity contour.'
STOP
CASE('--circle','-c')
drawcir = .true.
CASE('-i')
CALL getarg(i+1,arg)
READ(arg,*)filename
givefnm = .true.
CASE('-m')
CALL getarg(i+1,arg)
READ(arg,*)frames(1)
CALL getarg(i+2,arg)
READ(arg,*)frames(2)
drawmode = 2
CASE('--project','-p')
toproject = .true.
CALL getarg(i+1,arg)
READ(arg,*,iostat=ioerr)argaline
if(ioerr.ne.0)argaline = -3.d0
CASE('--stellar','-s')
drawstellar = .true.
CASE('--save')
save = .true.
CASE('--zmax','-z')
CALL getarg(i+1,arg)
READ(arg,*)zmax
rauto = .false.
CASE('-q')
drawq = .true.
CASE('--cut')
cut = .true.
CASE('-v')
write(6,'(a)')'Compiled at: '//__DATE__//' '//__TIME__
STOP
CASE('--velocity')
velocity = .true.
ENDSELECT
ENDDO
ENDIF
ENDSUBROUTINE
ENDMODULE
PROGRAM density1
USE PLOTTING
USE STELLARDISK_MODEL
USE STELLARDISK,only:FindSpiral,pi_n=>pi,sigma1
USE projections,only:argaline
USE io
USE argument
USE math
IMPLICIT NONE
TYPE gastyp
TYPE(typintplt2) ::Intpl_Den
TYPE(typintplt2) ::Intpl_Mom_x
TYPE(typintplt2) ::Intpl_Mom_y
CHARACTER(len=32) ::filename
DOUBLE PRECISION,ALLOCATABLE::density(:,:)
DOUBLE PRECISION,ALLOCATABLE::momentum(:,:,:)
DOUBLE PRECISION,ALLOCATABLE::x(:)
DOUBLE PRECISION,ALLOCATABLE::y(:)
DOUBLE PRECISION ::dx,dy
ENDTYPE
type(gastyp) ::gas
INTEGER ::i,j,k,l
CHARACTER(len=32) ::arg
DOUBLE PRECISION ::domain= 12.d0,dx,dy,r,th,pf(2),pi(2)
DOUBLE PRECISION,ALLOCATABLE ::density(:,:,:),xcoord(:),ycoord(:)
DOUBLE PRECISION,ALLOCATABLE ::GasDensity(:,:)
DOUBLE PRECISION ::limit = 100.d0
DOUBLE PRECISION ::d
DOUBLE PRECISION ::intb(4)
INTEGER,PARAMETER ::n=600
INTEGER ::PGBEG
type(typspiral) ::spiral
!!Reading from arguments
CALL readarg
gas.filename = filename
!!check if gas reading files are set:
if(.not.(givefnm.or.drawmode.eq.2))then
write(0,'(a)')'No input file or files set'
STOP
ENDIF
!preparing spiral
CALL stdpara.readstd
CALL spiral.init(500,12.d0,stdpara,2)
CALL spiral.readw(2)
CALL FindSpiral(spiral)
!set up grid
ALLOCATE(density(n,n,4))!1 for stellar, 2 for gas, 3 for q of gas,
!4 for approaching velocity
ALLOCATE(xcoord(n))
ALLOCATE(ycoord(n))
dx = domain/dble(n)*2.d0
dy = domain/dble(n)*2.d0
DO i = 1, n
xcoord(i) = dx*0.5d0 - domain + dble(i-1)*dx
ycoord(i) = dy*0.5d0 - domain + dble(i-1)*dy
ENDDO
!!Filling stellar density
!$OMP PARALLEL SHARED(density,spiral,gas) PRIVATE(j,r,th,pi,pf,d,k,l,intb)
!$OMP DO
DO i = 1, n
DO j = 1, n
!read coordinate
pf = (/xcoord(i),ycoord(j)/)
pi = pf
!project to original coordinate if nessesary
if(toproject)CALL projection(pi,pf)
!Fill in Stellar Density
r = sqrt(pi(1)**2+pi(2)**2)
th = atan2(pi(2),pi(1))
IF(r.lt.spiral.fortoone)then
d = sigma1(r,th,spiral)
ELSE
d = 0.d0
ENDIF
!remove wired pixel
if(isnan(d))d = 0.d0
density(i,j,1) = d
ENDDO
ENDDO
!$OMP END DO
!$OMP END PARALLEL
FORALL(i = 1:n,j = 1:n,density(i,j,1)<0)
density(i,j,1) = 0.d0
ENDFORALL
points(1,:) = (/0.0,10.636/)
points(2,:) = (/0.0,-10.636/)
points(3,:) = (/0.0,4.727/)
points(4,:) = (/0.0,-4.727/)
CALL dprojection(points)
!!looping to draw movie
!read in gas density
SELECT CASE(drawmode)
CASE(1) ! plot static density map
CALL ReadGasDensity
CALL FillGasDensity
IF(.not.drawq)THEN !plot density only
IF (PGBEG(0,'/xs',1,1) .NE. 1) STOP
!!set page
CALL PGENV(-real(domain),real(domain),-real(domain),real(domain),1,0)
CALL PlotGasDensity(gas.filename)
ELSE !plot density and q
IF (PGBEG(0,'/xs',1,1) .NE. 1) STOP
!!set coordinate
CALL PGSUBP(2,1)
CALL PGENV(-real(domain),real(domain),-real(domain),real(domain),1,0)
CALL PlotGasDensity(gas.filename)
CALL PGENV(-real(domain),real(domain),-real(domain),real(domain),1,0)
CALL Plotq(gas.filename)
ENDIF
CASE(2) ! plot movie
IF (PGBEG(0,'/xs',1,1) .NE. 1) STOP
!!set plot layout
IF(.not.drawq)THEN
CALL PGENV(-real(domain),real(domain),-real(domain),real(domain),1,0)
ELSE
CALL PGSUBP(2,1)
CALL PGENV(-real(domain),real(domain),-real(domain),real(domain),1,0)
ENDIF
DO i = frames(1),frames(2)
write(gas.filename,'("M",I0.4,".h5")')i
print *,gas.filename
CALL ReadGasDensity
CALL FillGasDensity
IF(rauto)THEN
CALL meshplot(density(:,:,2),n,domain,maxval(density(:,:,2)),zmin_in=0.d0)
ELSE
CALL meshplot(density(:,:,2),n,domain,zmax,zmin_in=0.d0)
ENDIF
CALL PGETXT
CALL PGPTXT(5.,8.,0.,0.,gas.filename)
CALL PlotCircle
CALL PGBBUF
CALL PlotStellarDensity(0)
CALL PGEBUF
IF(drawq)THEN
CALL PGPANL(2,1)
CALL Plotq(gas.filename)
CALL PGPANL(1,1)
ELSE
ENDIF
ENDDO
ENDSELECT
CALL PGCLOS
CALL PrintGasDensity
!CALL StellarGasPlot(density,n,domain,4)
!!close h5 files writing
!CALL h5write(xcoord,N,'density.h5','xcoord')
!CALL h5write(ycoord,N,'density.h5','ycoord')
!CALL h5write(density(:,:,2),N,N,'density.h5','density')
DEALLOCATE(xcoord)
DEALLOCATE(ycoord)
DEALLOCATE(density)
!CALL PhaseIntegrate
!CALL ENDSTELLARDISK
DEALLOCATE(gas.density)
DEALLOCATE(gas.x)
DEALLOCATE(gas.y)
CALL gas.Intpl_Den.free
CALL gas.Intpl_Mom_x.free
CALL gas.Intpl_Mom_y.free
STOP
CONTAINS
SUBROUTINE ReadGasDensity
IMPLICIT NONE
INTEGER ::fsize(2)
fsize = h5size2d(gas.filename,'density')
IF(.not.ALLOCATED(gas.density))ALLOCATE(gas.density(fsize(1),fsize(2)))
IF(.not.ALLOCATED(gas.momentum))ALLOCATE(gas.momentum(fsize(1),fsize(2),2))
IF(.not.ALLOCATED(gas.x))ALLOCATE(gas.x(fsize(1)))
IF(.not.ALLOCATED(gas.y))ALLOCATE(gas.y(fsize(2)))
CALL h5read(gas.density,fsize(1),fsize(2),gas.filename,'density')
CALL h5read(gas.momentum(:,:,1),fsize(1),fsize(2),gas.filename,'momx')
CALL h5read(gas.momentum(:,:,2),fsize(1),fsize(2),gas.filename,'momy')
CALL h5read(gas.x,fsize(1),gas.filename,'x')
CALL h5read(gas.y,fsize(2),gas.filename,'y')
gas.dx = gas.x(2)-gas.x(1)
gas.dy = gas.y(2)-gas.y(1)
CALL gas.Intpl_Den.init(gas.density(:,:),gas.x,gas.y)
CALL gas.Intpl_Mom_x.init(gas.momentum(:,:,1),gas.x,gas.y)
CALL gas.Intpl_Mom_y.init(gas.momentum(:,:,2),gas.x,gas.y)
ENDSUBROUTINE
SUBROUTINE PlotGasDensity(title)
USE argument
USE plotting
IMPLICIT NONE
CHARACTER(len=32) ::title
INTEGER ::i
IF(rauto)THEN
zmax = maxval(density(:,:,2))*1.1d0
write(6,*)achar(27)//'[33m Plotting z scale :',zmax,achar(27)//'[0m'
ELSE
write(6,*)achar(27)//'[33m Plotting z scale :',zmax,achar(27)//'[0m'
ENDIF
!!plot gas
CALL meshplot(density(:,:,2),n,domain,zmax,zmin_in=0.d0)
!!plot velocity
IF(velocity.and.toproject)THEN
CALL PGSAVE
CALL PGSCI(0)
CALL contourplot(density(:,:,4),n,domain,20)
CALL PGUNSA
ENDIF
!!plot circle
CALL PlotCircle
!!plot stellar contour
CALL PlotStellarDensity(0)
CALL PGLAB('kpc','kpc',title)
ENDSUBROUTINE
SUBROUTINE PrintGasDensity
USE argument
USE plotting
IMPLICIT NONE
INTEGER ::i
INTEGER ::PGOPEN
CHARACTER(len=225) ::CH
CALL GETCWD(CH)
CH = TRIM(CH(INDEX(CH,'/',.true.)+1:LEN(CH)))
IF (PGOPEN(TRIM(CH)//'_GasDensity.png/png') .NE. 1) STOP
!!plot gas
CALL PGENV(-real(domain),real(domain),-real(domain),real(domain),1,0)
CALL meshplot(density(:,:,2),n,domain,zmax,zmin_in=0.d0)
!!plot circle
CALL PlotCircle
!!plot stellar contour
CALL PlotStellarDensity(0)
CALL PGLAB('kpc','kpc','Gas Density')
CALL PGCLOS
IF(drawq)THEN
IF (PGOPEN(TRIM(CH)//'_GasDensity_Q.png/png') .NE. 1) STOP
CALL PGENV(-real(domain),real(domain),-real(domain),real(domain),1,0)
CALL plotq(CH)
CALL PGCLOS
ENDIF
ENDSUBROUTINE
SUBROUTINE plotq(title)
USE argument
USE plotting
IMPLICIT NONE
CHARACTER(*) ::title
INTEGER ::i
!IF(rauto)THEN
! zmax = maxval(density(:,:,3))*1.1d0
! write(6,*)achar(27)//'[33m Plotting z scale :',zmax,achar(27)//'[0m'
!ELSE
! write(6,*)achar(27)//'[33m Plotting z scale :',zmax,achar(27)//'[0m'
!ENDIF
!!plot q of gas
CALL meshplot(density(:,:,3),n,domain, 3.d0,zmin_in=0.d0)
!!plot circle
CALL PlotCircle
!!plot stellar contour
CALL PlotStellarDensity(1)
CALL PGLAB('kpc','kpc','-log(Q)')
ENDSUBROUTINE
SUBROUTINE FillGasDensity
USE STELLARDISK, only:kappa,GravConst,pi_n=>pi,Omega
USE galaxy, only:gasdensity
IMPLICIT NONE
DOUBLE PRECISION, EXTERNAL ::vprojection
DOUBLE PRECISION ::v(2)
!!Filling stellar density
!$OMP PARALLEL SHARED(density,spiral,gas) PRIVATE(j,r,th,pi,pf,d,k,l,intb,v)
!$OMP DO
DO i = 1, n
DO j = 1, n
!read coordinate
pf = (/xcoord(i),ycoord(j)/)
pi = pf
!project to original coordinate if nessesary
if(toproject)CALL projection(pi,pf)
r = sqrt(pi(1)**2+pi(2)**2)
th = atan2(pi(2),pi(1))
!Fill in Gas Density
!find interpolating grid
!!if pixel is not in simulation grid
IF(pi(1)<minval(gas.x).or.pi(1)>maxval(gas.x).or.pi(2)<minval(gas.y).or.pi(2)>maxval(gas.y))then
density(i,j,:) = 0.d0
ELSE
!! gas density
density(i,j,2) = gas.Intpl_Den.find(pi(1),pi(2))
! Toomre Q
density(i,j,3) = &
kappa(r,spiral)*8.d0/pi_n/GravConst/(density(i,j,2))/1d6
density(i,j,3) = -log(density(i,j,3))
! approaching velocity
v(1) = gas.Intpl_Mom_x.find(pi(1),pi(2))/density(i,j,2)-r*Omega(r,spiral)*cos(th)
v(2) = gas.Intpl_Mom_y.find(pi(1),pi(2))/density(i,j,2)+r*Omega(r,spiral)*sin(th)
density(i,j,4) = vprojection(v)
ENDIF
IF(cut)THEN
IF(r>spiral.co)THEN
density(i,j,2) = 0.d0
density(i,j,3) = 0.d0
ENDIF
ENDIF
ENDDO
ENDDO
!$OMP END DO
!$OMP END PARALLEL
ENDSUBROUTINE
SUBROUTINE PlotStellarDensity(color)
USE argument,only:drawstellar
USE plotting
IMPLICIT NONE
INTEGER ::color
IF(drawstellar)then
write(6,*)achar(27)//'[33m Drawing stellar contour',achar(27)//'[0m'
!CALL PGSCI(0)
CALL PGSAVE
CALL PGSCI(color)
CALL contourplot(density(:,:,1),n,domain,4)
CALL PGUNSA
ENDIF
ENDSUBROUTINE
END PROGRAM
SUBROUTINE projection(pi,pf)
USE projections
IMPLICIT NONE
DOUBLE PRECISION ::pi(2),pf(2)
!!BLAS
CHARACTER(1) ::TRANS
DOUBLE PRECISION ::ALPHA = 1.d0
DOUBLE PRECISION ::A(2,2)
DOUBLE PRECISION ::X(2),Y(2)
DOUBLE PRECISION ::BETA = 0.d0
INTEGER ::M = 2
INTEGER ::N = 2
INTEGER ::LDA = 2
INTEGER ::INCX = 1
INTEGER ::INCY = 1
CALL set_angles
X = pf
A = R(-aolin)
TRANS = 'n'
CALL DGEMV (TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
X = Y
A = INVC(apitc)
TRANS = 'n'
CALL DGEMV (TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
X = Y
A = R(-aline)
TRANS = 'n'
CALL DGEMV (TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
pi = y
ENDSUBROUTINE
SUBROUTINE dprojection(p)
USE projections
USE plotting,only:points
IMPLICIT NONE
DOUBLE PRECISION ::pi(2),pf(2)
REAL ::p(4,2)
INTEGER ::i
!!BLAS
CHARACTER(1) ::TRANS
DOUBLE PRECISION ::ALPHA = 1.d0
DOUBLE PRECISION ::A(2,2)
DOUBLE PRECISION ::X(2),Y(2)
DOUBLE PRECISION ::BETA = 0.d0
INTEGER ::M = 2
INTEGER ::N = 2
INTEGER ::LDA = 2
INTEGER ::INCX = 1
INTEGER ::INCY = 1
call set_angles
DO i = 1,4
X = points(i,:)
A = R(aolin)
TRANS = 'n'
CALL DGEMV (TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
p(i,:) = y
ENDDO
ENDSUBROUTINE
FUNCTION vprojection(p)
USE projections
IMPLICIT NONE
DOUBLE PRECISION ::p(2)
DOUBLE PRECISION ::vprojection
!!BLAS
CHARACTER(1) ::TRANS
DOUBLE PRECISION ::ALPHA = 1.d0
DOUBLE PRECISION ::A(2,2)
DOUBLE PRECISION ::X(2),Y(2)
DOUBLE PRECISION ::BETA = 0.d0
INTEGER ::M = 2
INTEGER ::N = 2
INTEGER ::LDA = 2
INTEGER ::INCX = 1
INTEGER ::INCY = 1
CALL set_angles
X = p
A = R(-aolin)
TRANS = 'n'
CALL DGEMV (TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
p = y
!! now velocity have changed to the coordinate of vertical line of node.
vprojection = p(1)*sin(apitc)
ENDFUNCTION
SUBROUTINE set_angles
USE projections
IMPLICIT NONE
!observed angle of line of node 28 degree
aolin = -20.0d0/180.d0*pi_n
!aolin = -28.3d0/180.d0*pi_n
!aolin = -32.d0/180.d0*pi_n
!inclination angle is 55 degree ?
apitc = 55.d0/180.d0*pi_n
!tuned angle of line of node
!aline = 30.d0/180.d0*pi_n
aline = argaline/180.d0*pi_n
ENDSUBROUTINE
SUBROUTINE PhaseIntegrate(spiral)
USE STELLARDISK_MODEL
USE STELLARDISK,only:k3sqrt,pi
IMPLICIT NONE
type(typspiral),TARGET ::spiral
DOUBLE PRECISION ::ri,rf
DOUBLE PRECISION ::RE,AE,RR,ERR,ANS
INTEGER ::IFLAG
ri = 0.d0
rf = 4.d0
RR = 2.d0
RE = 1d-8
AE = 1d-8
CALL DFZERO(F,ri,rf,RR,RE,AE,IFLAG)
ri = 0.d0
rf = min(5.1d0,rf)
ERR = 1.d-8
CALL DGAUS8(g,ri,rf,ERR,ANS,IFLAG)
print *,'phase before zero',ans/pi,'to',rf
CONTAINS
FUNCTION F(r)
IMPLICIT NONE
DOUBLE PRECISION ::r,F
F = k3sqrt(r,spiral)
ENDFUNCTION
FUNCTION g(r)
IMPLICIT NONE
DOUBLE PRECISION ::r,G
g = REAL(sqrt(k3sqrt(r,spiral)))
ENDFUNCTION
ENDSUBROUTINE
SUBROUTINE SOLVE(rho,phi,N,L)
IMPLICIT NONE
DOUBLE PRECISION ::rho(2*n,2*n),phi(2*n,2*n),L
REAL,ALLOCATABLE ::rrho(:,:)
real ::A,B,C,D,ELMBDA,PERTRB
real,ALLOCATABLE::BDA(:),BDB(:),BDC(:),BDD(:),w(:)
INTEGER ::IDIMF,MBDCND,NBDCND
INTEGER ::IERROR,wdim,n,i
IDIMF = 2*n
A = -real(l)
B = real(l)
MBDCND = 1
C = -real(l)
D = real(l)
NBDCND = 1
ELMBDA = 0.0
wdim = 13 + INT(LOG(dble(N))/log(2.))*N + 4*N
wdim = wdim*4
ALLOCATE(rrho(2*n,2*n))
ALLOCATE(BDA(2*n))
ALLOCATE(BDB(2*n))
ALLOCATE(BDC(2*n))
ALLOCATE(BDD(2*n))
ALLOCATE(W(wdim))
!print *,n,size(rho),size(rrho)
rrho = real(rho)
BDA(:) = 0.
BDB = BDA
BDC = BDA
BDD = BDA
CALL HSTCRT (A, B, 2*N, MBDCND, BDA, BDB, C, D, 2*N, NBDCND, BDC, BDD, ELMBDA, rrho,IDIMF, PERTRB, IERROR, W)
phi = dble(rrho)
DEALLOCATE(rrho)
DEALLOCATE(BDA)
DEALLOCATE(BDB)
DEALLOCATE(BDC)
DEALLOCATE(BDD)
DEALLOCATE(W)
ENDSUBROUTINE
SUBROUTINE PlotCircle
USE argument
IMPLICIT NONE
IF(.not.toproject.and.drawcir)THEN
write(6,*)achar(27)//'[33m Drawing circles',achar(27)//'[0m'
CALL PGSAVE
CALL PGSFS(2)
CALL PGSCI(0)
CALL PGCIRC(0.,0.,1.26)
CALL PGCIRC(0.,0.,2.36)
CALL PGCIRC(0.,0.,4.72)
CALL PGCIRC(0.,0.,10.636)
CALL PGCIRC(0.,0.,8.83)
CALL PGUNSA
ENDIF
ENDSUBROUTINE