-
Notifications
You must be signed in to change notification settings - Fork 11
/
optool_guts.f90
2081 lines (1857 loc) · 74.2 KB
/
optool_guts.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
subroutine MeerhoffMie(rad,lam,e1,e2,Csca,Cext,F11,F12,F33,F34,na)
implicit none
integer na,nangle
real*8 rad,lam,e1,e2,csca,cext,F11(na),F12(na),F33(na),F34(na)
integer nparts,develop,nsubr(1),ngaur(1),idis(1),ndis,i
real*8 delta,cutoff,thmin,thmax,step,wavel(1),Rem(1),fImm(1)
real*8 par1(1),par2(1),par3(1),rdis(1,300),nwrdis(1,300)
real*8 F(4,6000),theta(6000)
nangle = na
nparts = 1
develop = 0
delta = 1d-8
cutoff = 1d-8
thmin = 180d0*(real(1)-0.5d0)/real(na)
thmax = 180d0*(real(na)-0.5d0)/real(na)
step = (thmax-thmin)/real(na-1)
wavel = lam
Rem = e1
fImm = e2
nsubr = 1
ngaur = 1
idis = 0
par1 = rad
par2 = 0d0
par3 = 0d0
ndis = 1
rdis(1,1) =rad
nwrdis(1,1)=1d0
call mie(nparts, develop, delta, cutoff, thmin, thmax, step, &
wavel, Rem, fImm, nsubr , ngaur, idis, &
par1, par2, par3, ndis, rdis, nwrdis, &
nangle, theta, F ,Cext, Csca)
do i=1,nangle
F11(i) = F(1,nangle-i+1)
F12(i) = F(2,nangle-i+1)
F33(i) = F(3,nangle-i+1)
F34(i) = F(4,nangle-i+1)
enddo
return
end subroutine MeerhoffMie
subroutine mie(nparts, develop, delta, cutoff, thmin, thmax, step, &
wavel, Rem, fImm, nsubr , ngaur, idis, &
par1, par2, par3, ndis, rdis, nwrdis, &
nangle, theta, F ,outCext, outCsca)
implicit double precision (a-h,o-z)
implicit integer (i-n)
parameter ( NDang=6000, NDcoef=6000, NDpart = 4, nrunit=88 )
parameter ( NDdis=300 )
double precision lambda, nr, ni, miec, nwrdis,outCext,outCsca
integer develop
double complex m
character*60 s1, s2, s3, s4, s5, s6
character*10 scfile
dimension F(4,NDang), &
miec(13),u(NDang),w(NDang),coefs(4,4,0:NDcoef), &
wavel(*), Rem(*), fImm(*), nsubr(*), &
ngaur(*), idis(*), &
par1(*), par2(*), par3(*), &
rdis(1,*), nwrdis(1,*), &
scfile(NDpart), theta(NDang)
pi = dacos(-1.D0)
radtod = 180.D0/pi
scfile(1) = 'mie1_sc '
scfile(2) = 'mie2_sc '
scfile(3) = 'mie3_sc '
scfile(4) = 'mie4_sc '
!************************************************************************
!* Start a loop over different 'particles' referring to the different *
!* cases specified in the input *
!************************************************************************
do iparts=1,nparts
!************************************************************************
!* Determine the integration bounds for the radius integration *
!************************************************************************
call rminmax( idis(iparts), nsubr(iparts), ngaur(iparts), &
par1(iparts), par2(iparts), par3(iparts), cutoff, &
rmin, rmax, &
iparts, ndis, rdis, nwrdis )
lambda = wavel(iparts)
nr = Rem(iparts)
ni = fImm(iparts)
m = dcmplx(nr,-ni)
!************************************************************************
!* Calculate the scattering matrix : this is most of the work ! *
!************************************************************************
call scatmat( u, w, m, wavel(iparts), idis(iparts), &
thmin, thmax, step, develop, &
nsubr(iparts), ngaur(iparts), rmin, rmax, &
par1(iparts), par2(iparts), par3(iparts), &
delta, F, miec, nangle, &
ndis, rdis, nwrdis, iparts )
!************************************************************************
!* Test if the number of coefficients needed fits in the array coefs *
!* The number of coefficients is equal to the number of integration *
!* points in scattering angle !! (think about this !) *
!************************************************************************
ncoef = nangle
if ((ncoef .gt. NDcoef) .and. (develop .eq. 1)) then
write(*,*) ' main: too many coefficients needed :',ncoef
write(*,*) ' maximum array size NDcoef = ',NDcoef
write(*,*) ' Therefore I cannot do the expansion.'
develop=0
call clossc(nrunit)
endif
!************************************************************************
!* If required, calculate the expansion coefficients *
!************************************************************************
if (develop.eq.1) then
call devel(ncoef,nangle,u,w,F,coefs)
npunt = idnint((thmax-thmin)/step)+1
if (npunt .gt. NDang) then
write(*,*) ' main: too many angles for table npunt=',npunt
write(*,*) ' maximum array size NDang = ',NDang
npunt=61
write(*,*) ' I have set the number of angles to ',npunt
endif
endif
outCsca=miec(1)
outCext=miec(2)
if (nangle .ge. 1) then
cosbar = coefs(1,1,1)/3.D0
else
cosbar = 0.D0
endif
do j=1,nangle
i = nangle-j+1
xp = -100.d0*F(2,i)/F(1,i)
theta1 = radtod*dacos(u(i))
enddo
!************************************************************************
!* Print the coefficients on output. *
!* Evaluate the expansion in GSF at an equidistant set of scattering *
!* angles and print the resulting scattering matrix. *
!************************************************************************
if (develop.eq.1) then
call expand( ncoef, npunt, coefs, u, F, thmin, thmax, step )
do i=1,npunt
xp = -100.d0*F(2,i)/F(1,i)
theta1 = radtod*dacos(u(i))
enddo
endif
enddo
!************************************************************************
!* End of loop over 'particles'. *
!************************************************************************
return
end subroutine mie
subroutine anbn( m, x, nmax, psi, chi, d, an, bn )
!************************************************************************
!* Calculate the Mie coefficients an and bn. *
!************************************************************************
implicit double precision (a-h,o-z)
parameter( NDn=10000000 )
double complex m, zn, znm1, save, perm
double complex an(nmax), bn(nmax), D(nmax)
dimension psi(0:nmax), chi(0:nmax)
perm = 1.D0/m
perx = 1.D0/x
xn = 0.D0
do n=1, nmax
zn = dcmplx( psi(n), chi(n))
znm1 = dcmplx( psi(n-1), chi(n-1))
xn = dble(n)*perx
save = D(n)*perm+xn
an(n)= (save*psi(n)-psi(n-1)) / (save*zn-znm1)
save = m*D(n)+xn
bn(n)= (save*psi(n)-psi(n-1)) / (save*zn-znm1)
enddo
return
end subroutine anbn
subroutine clossc(iunit)
!************************************************************************
!* On entry : *
!* iunit number of the unit to be closed *
!* On exit : *
!* The file is closed. *
!* VLD January 9, 1989 *
!************************************************************************
!c close(unit=iunit,err=999)
return
!c 999 write(7,*) ' clossc: error in closing file with unit number',iunit
!c stop 'in clossc error in closing file'
end subroutine clossc
subroutine devel(ncoef,nangle,u,w,F,coefs)
!************************************************************************
!* Calculate the expansion coefficients of the scattering matrix in *
!* generalized spherical functions by numerical integration over the *
!* scattering angle. *
!************************************************************************
implicit double precision (a-h,o-z)
parameter( NDcoef=3000, NDang = 6000 )
dimension u(NDang),w(NDang),F(4,NDang)
dimension coefs(4,4,0:NDcoef), P00(NDang,2), P02(NDang,2) &
, P22(NDang,2), P2m2(NDang,2)
!************************************************************************
!* Initialization *
!************************************************************************
qroot6 = -0.25D0*dsqrt(6.D0)
do j=0,NDcoef
do i=1,4
do ii=1,4
coefs(ii,i,j)=0.D0
enddo
enddo
enddo
!************************************************************************
!* Multiply the scattering matrix F with the weights w for all angles *
!* We do this here because otherwise it should be done for each l *
!************************************************************************
do k=1,4
do i=1, nangle
F(k,i) = w(i)*F(k,i)
enddo
enddo
!************************************************************************
!* Start loop over the coefficient index l *
!* first update generalized spherical functions, then calculate coefs. *
!* lold and lnew are pointer-like indices used in recurrence *
!************************************************************************
lnew = 1
lold = 2
do l=0, ncoef
if (l .eq. 0) then
!************************************************************************
!* Adding paper Eq. (77) with m=0 *
!************************************************************************
do i=1, nangle
P00(i,lold) = 1.D0
P00(i,lnew) = 0.D0
P02(i,lold) = 0.D0
P22(i,lold) = 0.D0
P2m2(i,lold) = 0.D0
P02(i,lnew) = 0.D0
P22(i,lnew) = 0.D0
P2m2(i,lnew) = 0.D0
enddo
else
fac1 = (2.D0*l-1.d0)/dble(l)
fac2 = dble(l-1.d0)/dble(l)
!************************************************************************
!* Adding paper Eq. (81) with m=0 *
!************************************************************************
do i=1, nangle
P00(i,lold) = fac1*u(i)*P00(i,lnew) - fac2*P00(i,lold)
enddo
endif
if (l .eq. 2) then
!************************************************************************
!* Adding paper Eqs. (78) and (80) with m=2 *
!* sql4 contains the factor dsqrt(l*l-4) needed in *
!* the recurrence Eqs. (81) and (82) *
!************************************************************************
do i=1, nangle
P02(i,lold) = qroot6*(1.D0-u(i)*u(i))
P22(i,lold) = 0.25D0*(1.D0+u(i))*(1.D0+u(i))
P2m2(i,lold) = 0.25D0*(1.D0-u(i))*(1.D0-u(i))
P02(i,lnew) = 0.D0
P22(i,lnew) = 0.D0
P2m2(i,lnew) = 0.D0
enddo
sql41 = 0.D0
else if (l .gt. 2) then
!************************************************************************
!* Adding paper Eq. (82) with m=0 and m=2 *
!************************************************************************
sql4 = sql41
sql41 = dsqrt(dble(l*l)-4.d0)
twol1 = 2.D0*dble(l)-1.d0
tmp1 = twol1/sql41
tmp2 = sql4/sql41
denom = (dble(l)-1.d0)*(dble(l*l)-4.d0)
fac1 = twol1*(dble(l)-1.d0)*dble(l)/denom
fac2 = 4.D0*twol1/denom
fac3 = dble(l)*((dble(l)-1.d0)*(dble(l)-1.d0)-4.d0)/denom
do i=1, nangle
P02(i,lold) = tmp1*u(i)*P02(i,lnew) - tmp2*P02(i,lold)
P22(i,lold) = (fac1*u(i)-fac2)*P22(i,lnew) &
- fac3*P22(i,lold)
P2m2(i,lold)= (fac1*u(i)+fac2)*P2m2(i,lnew) &
- fac3*P2m2(i,lold)
enddo
endif
!************************************************************************
!* Switch indices so that lnew indicates the function with *
!* the present index value l, this mechanism prevents swapping *
!* of entire arrays. *
!************************************************************************
itmp = lnew
lnew = lold
lold = itmp
!************************************************************************
!* Now calculate the coefficients by integration over angle *
!* See de Haan et al. (1987) Eqs. (68)-(73). *
!* Remember for Mie scattering : F11 = F22 and F33 = F44 *
!************************************************************************
alfap = 0.D0
alfam = 0.D0
do i=1, nangle
coefs(1,1,l) = coefs(1,1,l) + P00(i,lnew)*F(1,i)
alfap = alfap + P22(i,lnew)*(F(1,i)+F(3,i))
alfam = alfam + P2m2(i,lnew)*(F(1,i)-F(3,i))
coefs(4,4,l) = coefs(4,4,l) + P00(i,lnew)*F(3,i)
coefs(1,2,l) = coefs(1,2,l) + P02(i,lnew)*F(2,i)
coefs(3,4,l) = coefs(3,4,l) + P02(i,lnew)*F(4,i)
enddo
!************************************************************************
!* Multiply with trivial factors like 0.5D0*(2*l+1) *
!************************************************************************
fl = dble(l)+0.5D0
coefs(1,1,l) = fl*coefs(1,1,l)
coefs(2,2,l) = fl*0.5D0*(alfap+alfam)
coefs(3,3,l) = fl*0.5D0*(alfap-alfam)
coefs(4,4,l) = fl*coefs(4,4,l)
coefs(1,2,l) = fl*coefs(1,2,l)
coefs(3,4,l) = fl*coefs(3,4,l)
coefs(2,1,l) = coefs(1,2,l)
coefs(4,3,l) = -coefs(3,4,l)
enddo
!************************************************************************
!* End of loop over index l *
!************************************************************************
!************************************************************************
!* Remove the weight factor from the scattering matrix *
!************************************************************************
do k=1, 4
do i=1, nangle
F(k,i) = F(k,i)/w(i)
enddo
enddo
return
end subroutine devel
subroutine expand(ncoef,nangle,coefs,u,F,thmin,thmax,step)
! ***********************************************************
implicit double precision (a-h,o-z)
parameter( pi=3.141592653589793238462643D0, radfac= pi/180.D0 )
parameter( NDcoef=3000, NDang=6000 )
dimension u(NDang),F(4,NDang)
dimension coefs(4,4,0:NDcoef), P00(NDang,2), P02(NDang,2)
!
! **********************************************************
!
! initialize
do i=1, nangle
u(i) = dcos(radfac*(thmin+dble(i-1)*step))
enddo
qroot6 = -0.25D0*dsqrt(6.D0)
!************************************************************************
!* Set scattering matrix F to zero *
!************************************************************************
do k=1,4
do i=1, nangle
F(k,i) = 0.D0
enddo
enddo
!************************************************************************
!* Start loop over the coefficient index l *
!* first update generalized spherical functions, then calculate coefs. *
!* lold and lnew are pointer-like indices used in recurrence *
!************************************************************************
lnew = 1
lold = 2
do l=0, ncoef
if (l .eq. 0) then
!************************************************************************
!* Adding paper Eqs. (76) and (77) with m=0 *
!************************************************************************
do i=1, nangle
P00(i,lold) = 1.D0
P00(i,lnew) = 0.D0
P02(i,lold) = 0.D0
P02(i,lnew) = 0.D0
enddo
else
fac1 = (2.D0*dble(l)-1.d0)/dble(l)
fac2 = (dble(l)-1.d0)/dble(l)
!************************************************************************
!* Adding paper Eq. (81) with m=0 *
!************************************************************************
do i=1, nangle
P00(i,lold) = fac1*u(i)*P00(i,lnew) - fac2*P00(i,lold)
enddo
endif
if (l .eq. 2) then
!************************************************************************
!* Adding paper Eq. (78) *
!* sql4 contains the factor dsqrt((l+1)*(l+1)-4) needed in *
!* the recurrence Eqs. (81) and (82) *
!************************************************************************
do i=1, nangle
P02(i,lold) = qroot6*(1.D0-u(i)*u(i))
P02(i,lnew) = 0.D0
enddo
sql41 = 0.D0
else if (l .gt. 2) then
!************************************************************************
!* Adding paper Eq. (82) with m=0 *
!************************************************************************
sql4 = sql41
sql41 = dsqrt(dble(l*l)-4.d0)
tmp1 = (2.D0*dble(l)-1.d0)/sql41
tmp2 = sql4/sql41
do i=1, nangle
P02(i,lold) = tmp1*u(i)*P02(i,lnew) - tmp2*P02(i,lold)
enddo
endif
!************************************************************************
!* Switch indices so that lnew indicates the function with *
!* the present index value l, this mechanism prevents swapping *
!* of entire arrays. *
!************************************************************************
itmp = lnew
lnew = lold
lold = itmp
!************************************************************************
!* Now add the l-th term to the scattering matrix. *
!* See de Haan et al. (1987) Eqs. (68)-(73). *
!* Remember for Mie scattering : F11 = F22 and F33 = F44 *
!************************************************************************
do i=1, nangle
F(1,i) = F(1,i) + coefs(1,1,l)*P00(i,lnew)
F(2,i) = F(2,i) + coefs(1,2,l)*P02(i,lnew)
F(3,i) = F(3,i) + coefs(4,4,l)*P00(i,lnew)
F(4,i) = F(4,i) + coefs(3,4,l)*P02(i,lnew)
enddo
enddo
!************************************************************************
!* End of loop over index l *
!************************************************************************
return
end subroutine expand
subroutine fichid( m, x, nchi, nmax, nD, psi, chi, D )
!************************************************************************
!* Calculate functions psi(x) chi(x) and D(z) where z = mx. *
!* On entry, the following should be supplied : *
!* m : complex index of refraction *
!* x : sizeparameter *
!* nchi : starting index for backwards recurrence of chi *
!* nmax : number of chi, psi and D that must be available *
!* nD : starting index for backwards recurrence of D *
!* On exit, the desired functions are returned through psi, chi and D *
!************************************************************************
implicit double precision (a-h,o-z)
double complex D,m,z, perz, zn1
dimension psi(0:nchi), chi(0:nmax+1), D(nd)
z = m*x
perz = 1.D0/z
perx = 1.D0/x
sinx = dsin(x)
cosx = dcos(x)
psi(nchi)=0d0
!************************************************************************
!* (mis-) use the array psi to calculate the functions rn(x)
!* De Rooij and van der Stap Eq. (A6)
!************************************************************************
do n=nchi-1, 0, -1
psi(n) = 1.D0 / (dble(2*n+1)/x - psi(n+1))
enddo
!************************************************************************
!* Calculate functions D(z) by backward recurrence
!* De Rooij and van der Stap Eq. (A11)
!************************************************************************
D(nD) = dcmplx(0.D0,0.D0)
do n=nD - 1, 1, -1
zn1 = dble(n+1)*perz
D(n) = zn1 - 1.D0/(D(n+1)+zn1)
enddo
!************************************************************************
!* De Rooij and van der Stap Eqs. (A3) and (A1)
!* De Rooij and van der Stap Eq. (A5) where psi(n) was equal to r(n)
!* and Eq. (A2)
!************************************************************************
psi(0) = sinx
psi1 = psi(0)*perx - cosx
if (dabs(psi1) .gt. 1.d-4) then
psi(1) = psi1
do n=2,nmax
psi(n) = psi(n)*psi(n-1)
enddo
else
do n=1,nmax
psi(n) = psi(n)*psi(n-1)
enddo
endif
!************************************************************************
!* De Rooij and van der Stap Eqs. (A4) and (A2)
!************************************************************************
chi(0) = cosx
chi(1) = chi(0)*perx + sinx
do n=1, nmax-1
chi(n+1) = dble(2*n+1)*chi(n)*perx - chi(n-1)
enddo
!* DEBUG
! write(7,*) ' fichid: x = ',x
! write(7,12)
! do n=0, nchi
! write(7,11) n,psi(n),chi(n),D(n)
! enddo
!11 format(i4,4e24.14)
!12 format(' n',t20,'psi(n)',t44,'chi(n)',t68 &
! ,'Re(D(n))',t92,'Im(D(n))',/ &
! ,' ----------------------------------------------------' &
! ,'-----------------------------------------------------')
!* END DEBUG
return
end subroutine fichid
subroutine gaulegCP(ndim,ngauss,a,b,x,w)
!************************************************************************
!* Given the lower and upper limits of integration a and b, and given *
!* the number of Gauss-Legendre points ngauss, this routine returns *
!* through array x the abscissas and through array w the weights of *
!* the Gauss-Legendre quadrature formula. Eps is the desired accuracy *
!* of the abscissas. This routine is documented further in : *
!* W.H. Press et al. 'Numerical Recipes' Cambridge Univ. Pr. (1987) *
!* page 125 ISBN 0-521-30811-9 *
!************************************************************************
implicit double precision (a-h,o-z)
parameter( eps = 1.d-14 )
double precision x(ndim), w(ndim)
double precision a,b,xm,xl,z,p1,p2,p3,pp,z1,pi
pi = 4.D0*datan(1.D0)
m = (ngauss+1)/2
xm = 0.5D0*(a+b)
xl = 0.5D0*(b-a)
do i=1,m
!* THIS IS A REALLY CLEVER ESTIMATE :
z = dcos(pi*(dble(i)-0.25D0)/(dble(ngauss)+0.5D0))
1 continue
p1 = 1.D0
p2 = 0.D0
do j=1,ngauss
p3 = p2
p2 = p1
p1 = ((dble(2*j)-1.d0)*z*p2-(dble(j)-1.d0)*p3)/dble(j)
enddo
pp = ngauss*(z*p1-p2)/(z*z-1.d0)
z1 = z
z = z1-p1/pp
if (dabs(z-z1) .gt. eps) goto 1
x(i) = xm-xl*z
x(ngauss+1-i)= xm+xl*z
w(i) = 2.D0*xl/((1.D0-z*z)*pp*pp)
!****** write(7,*) ' gaulegCP: ',i,' x=',x(i),' w=',w(i)
w(ngauss+1-i)= w(i)
enddo
return
end subroutine gaulegCP
subroutine gausspt(ndim,ngauss,a,b,x,w)
!
! ***********************************************************
!
! put the gauss-legendre points of order ndim in array x,
! the weights in array w. the points are transformed to the
! interval [a,b]
!
! ***********************************************************
implicit double precision (a-h,o-z)
dimension x(ndim),w(ndim)
! ***********************************************************
!
! find starting values
!
gn = 0.5D0/dble(ngauss)
extra = 1.0D0/(.4D0*dble(ngauss)*dble(ngauss)+5.0D0)
xz = -gn
nt = 0
nteken = 0
5 pnm2 = 1.0D0
pnm1 = xz
do i=2,ngauss
pnm1xz = pnm1*xz
pn = 2.0D0*pnm1xz-pnm2-(pnm1xz-pnm2)/dble(i)
pnm2 = pnm1
pnm1 = pn
enddo
mteken=1
if (pn .le. 0.0D0) mteken=-1
if ((mteken+nteken) .eq. 0) then
nt = nt+1
x(nt) = xz
endif
nteken = mteken
if ((1.0D0-xz) .le. extra) goto 30
xz = xz+(1.D0-xz*xz)*gn+extra
goto 5
30 continue
! ***********************************************************
!
! determine zero's and weights
do i=1,nt
xz = x(i)
delta2 = 1.D0
35 pnm2 = 1.0D0
pnm1 = xz
pnm1af = 1.0D0
z = .5D0 + 1.5D0*xz*xz
do k=2,ngauss
pnm1xz = pnm1*xz
pn = 2.0D0*pnm1xz-pnm2-(pnm1xz-pnm2)/dble(k)
pnaf = xz*pnm1af+k*pnm1
z = z+(dble(k)+0.5D0)*pn*pn
pnm2 = pnm1
pnm1 = pn
pnm1af = pnaf
enddo
delta1 = pn/pnaf
xz = xz-delta1
if(delta1.lt.0.0D0) delta1=-delta1
if( (delta1.ge.delta2) .and. (delta2.lt.1.d-14) ) goto 50
delta2 = delta1
goto 35
50 x(i) = xz
w(i) = 1.0D0/z
enddo
! ***********************************************************
!
! transform to the interval [a,b]
nghalf = ngauss/2
ngp1 = ngauss+1
ntp1 = nt+1
apb = a+b
bmag2 = (b-a)/2.0D0
do i=1,nghalf
x(ngp1-i) = b-bmag2*(1.0D0-x(ntp1-i))
w(ngp1-i) = bmag2*w(ntp1-i)
enddo
if (nghalf .ne. nt) then
x(nt) = apb/2.0D0
w(nt) = w(1)*bmag2
endif
do i=1,nghalf
x(i) = apb-x(ngp1-i)
w(i) = w(ngp1-i)
enddo
return
end subroutine gausspt
subroutine pitau(u,nmax,pi,tau)
! ***********************************************************
! calculates pi,n(u) and tau,n(u) with upward recursion
!
! ***********************************************************
implicit double precision (a-h,o-z)
dimension pi(nmax),tau(nmax)
!
! starting values:
pi(1) = 1.D0
pi(2) = 3.D0*u
delta = 3.D0*u*u-1.d0
tau(1)= u
tau(2)= 2.D0*delta-1.d0
! upward recursion:
do n=2, nmax-1
pi(n+1) = dble(n+1)/dble(n)*delta + u*pi(n)
delta = u*pi(n+1) - pi(n)
tau(n+1)= dble(n+1)*delta - pi(n)
enddo
return
end subroutine pitau
subroutine rminmax( idis, nsubr, ngaur, par1, par2, par3, cutoff &
, rmin, rmax, iparts, ndis, rdis, nwrdis )
!************************************************************************
!* Find the integration bounds rmin and rmax for the integration over *
!* a size distribution. These bounds are chosen such that the size *
!* distribution falls below the user specified cutoff. It is essential *
!* that the size distribution is normalized such that the integral *
!* over all r is equal to one ! *
!* This is programmed rather clumsy and will in the future be changed *
!************************************************************************
implicit double precision (a-h,o-z)
parameter( eps = 1.d-10, NDpart = 4, NDdis=300 )
double precision nwithr, rdis(NDpart,NDdis), nwrdis(NDpart,NDdis)
dimension r(1), nwithr(1)
!*
if (idis.eq.0) then
rmin= par1
rmax= par1
else
goto (10,20,30,40,50,60,70,80,90) idis
write(*,*) ' rminmax: illegal size distribution index :',idis
stop 'in rminmax illegal size distribution index'
!************************************************************************
10 sef = 1.D0/dsqrt(par2+3.D0)
ref = 1.D0/(sef*sef*par2)
rref= ref
goto 100
20 ref = par1
sef = dsqrt(par2)
rref= ref
goto 100
30 sef = dsqrt(par3)
ref = dmax1(par1,par2)+sef
rref= dmin1(par1,par2)
goto 100
40 sef = dsqrt(dexp(dlog(par2)**2)-1.d0)
ref = par1*(1.D0+sef*sef)**0.4D0
rref= ref
goto 100
50 ref = par1
sef = dsqrt(ref)
rref= ref
goto 100
60 rmin= par2
rmax= par3
goto 999
70 ref = par2
sef = 2.D0*ref
rref=0.5D0*ref
goto 100
80 ref = (par1/(par2*par3))**par3
sef = 2.D0*ref
rref= 0.5D0*ref
goto 100
90 rmin = par1
rmax = par2
goto 999
!************************************************************************
!* search for a value of r such that the size distribution
!* is less than the cutoff. Start the search at ref+sef which *
!* guarantees that such a value will be found on the TAIL of the *
!* distribution. *
!************************************************************************
100 r(1) = ref+sef
r0 = ref
200 call sizedis( idis, par1, par2, par3, r, 1, nwithr, &
iparts, ndis, rdis, nwrdis )
if (nwithr(1) .gt. cutoff) then
r0 = r(1)
r(1) = 2.D0*r(1)
goto 200
endif
r1 = r(1)
!************************************************************************
!* Now the size distribution assumes the cutoff value somewhere *
!* between r0 and r1 Use bisection to find the corresponding r *
!************************************************************************
300 r(1) = 0.5D0*(r0+r1)
call sizedis( idis, par1, par2, par3, r, 1, nwithr, &
iparts, ndis, rdis, nwrdis )
if (nwithr(1) .gt. cutoff) then
r0 = r(1)
else
r1 = r(1)
endif
if ((r1-r0) .gt. eps) goto 300
rmax = 0.5D0*(r0+r1)
!************************************************************************
!* Search for a value of r on the low end of the size distribution *
!* such that the distribution falls below the cutoff. There is no *
!* guarantee that such a value exists, so use an extra test to see if *
!* the search comes very near to r = 0 *
!************************************************************************
r1 = rref
r0 = 0.D0
400 r(1) = 0.5D0*r1
call sizedis( idis, par1, par2, par3, r, 1, nwithr, &
iparts, ndis, rdis, nwrdis )
if (nwithr(1) .gt. cutoff) then
r1 = r(1)
if (r1 .gt. eps) goto 400
else
r0 = r(1)
endif
!************************************************************************
!* Possibly the size distribution goes through cutoff between r0 *
!* and r1 try to find the exact value of r where this happens by *
!* bisection. *
!* In case there is no solution, the algorithm will terminate soon. *
!************************************************************************
500 r(1) = 0.5D0*(r0+r1)
call sizedis( idis, par1, par2, par3, r, 1, nwithr, &
iparts, ndis, rdis, nwrdis )
if (nwithr(1) .gt. cutoff) then
r1 = r(1)
else
r0 = r(1)
endif
if ((r1-r0) .gt. eps) goto 500
if(r1 .le. eps) then
rmin = 0.D0
else
rmin = 0.5D0*(r0+r1)
endif
endif
999 continue!write(7,*) ' rminmax: found rmin = ',rmin,' rmax = ',rmax
return
end subroutine rminmax
subroutine scatmat( u, wth, m, lambda, idis, &
thmin, thmax, step, develop, &
nsub, ngauss, rmin, rmax, &
par1, par2, par3, &
delta, F, miec, nangle, &
ndis, rdis, nwrdis, iparts )
!************************************************************************
!* Calculate the scattering matrix of an ensemble of homogenous *
!* spheres. On entry, the following must be supplied : *
!* m : complex index of refraction *
!* lambda : wavelength *
!* idis : index of the size distribution *
!* nsub : number of subintervals for integration over r *
!* ngauss : number of Gauss points used per subinterval *
!* rmin : lower bound for integration over r *
!* rmax : upper bound for integration over r *
!* par1,2,3 : parameters of the size distribution *
!* delta : cutoff used in truncation of the Mie sum *
!* thmin : minimum scattering angle in degrees *
!* thmax : maximum scattering angle in degrees *
!* step : step in scattering angle in degrees *
!* develop : expansion in GSF (1) or not (0) *
!* On exit, the following results are returned : *
!* u : cosines of scattering angles *
!* wth : Gaussian weights associated with u *
!* F : scattering matrix for all cosines in u *
!* miec : array containing cross sections etc. *
!* nangle : the number of scattering angles *
!* When develop=1 u contains Gauss points, when develop=0 u contains *
!* cosines of equidistant angles between thmin and thmax with step. *
!************************************************************************
implicit double precision (a-h,o-z)
parameter( NDn=10000000, NDr=1000, NDang=6000, NDdis=300,NDpart=4)
double complex m, ci, Splusf, Sminf, cSplusf
double complex cSminf, Splusb, Sminb, cSplusb, cSminb
double complex, allocatable :: an(:), bn(:),D(:)
double precision,allocatable :: pi(:), tau(:), fi(:), chi(:)
double precision,allocatable :: facf(:), facb(:)
double precision lambda, nwithr, miec, numpar, thmin, thmax, step, nwrdis
integer develop
dimension u(NDang), wth(NDang), F(4,NDang),rdis(NDpart, NDdis)
dimension miec(13), nwrdis(NDpart, NDdis),r(NDr), w(NDr), nwithr(NDr)
logical symth
!************************************************************************
!* Initialization *
!************************************************************************
do j=1,NDang
do k=1,4
F(k,j)=0.D0
enddo
enddo
Csca = 0.D0
Cext = 0.D0
numpar= 0.D0
G = 0.D0
reff = 0.D0
nfou = 0
fac90 = 1.D0
ci = dcmplx(0.D0,1.D0)
call tstsym( develop, thmin, thmax, step, symth )
!************************************************************************
!* Constants *
!************************************************************************
pie = dacos(-1.d0)
radfac= pie/180.D0
rtox = 2.D0*pie/lambda
fakt = lambda*lambda/(2.D0*pie)
!* nfac is the number of precomputed factors (2n+1)/(n*(n+1))
nfac = 0
!************************************************************************
!* distinguish between distribution or not *
!************************************************************************
if (idis.eq.0) then
w(1) = 1.D0
r(1) = rmin
nwithr(1)= 1.D0
nsub = 1
ngauss = 1
dr = 0.D0
else
dr = (rmax-rmin)/dble(nsub)
!* call gausspt( ngauss, ngauss, (rmax-dr) , rmax ,r ,w )
call gaulegCP( ngauss, ngauss, (rmax-dr) , rmax ,r ,w )
call sizedis( idis, par1, par2, par3, r, ngauss, nwithr, &
iparts, ndis, rdis, nwrdis )
endif
!************************************************************************
!* Start integration over radius r with largest radius ! *
!************************************************************************
do l=nsub,1,-1
do i=ngauss,1,-1
sw = nwithr(i)*w(i)
x = rtox*r(i)
nmax = x + 4.05D0*x**(1.D0/3.D0) + 20
nfi = nmax+60
zabs = x*cdabs(m)
nD = zabs + 4.05D0*zabs**(1.D0/3.D0) + 70
allocate(an(max(nD,nfi,nmax)))
allocate(bn(max(nD,nfi,nmax)))
allocate(pi(max(nD,nfi,nmax)))
allocate(tau(max(nD,nfi,nmax)))
allocate(fi(0:max(nD,nfi,nmax)))
allocate(chi(0:max(nD,nfi,nmax)))
allocate(D(max(nD,nfi,nmax)))
allocate(facf(max(nD,nfi,nmax)))
allocate(facb(max(nD,nfi,nmax)))
if ((nD.gt.NDn) .or. (nfi.gt.NDn)) then
write(*,*) ' scatmat: estimated number of Mie-terms:',nD
write(*,*) ' for particle sizeparameter :',x
write(*,*) ' maximum NDn is only : ',NDn
stop 'in scatmat no room for all Mie terms'
endif
call fichid( m, x, nfi, nmax, nD, fi, chi, D )
call anbn( m, x, nmax, fi, chi, D, an, bn )
!************************************************************************
!* Precompute the factor (2n+1)/(n*(n+1)) needed in Mie sum over n *
!************************************************************************
if (nmax .gt. nfac) then
do n=nfac+1, nmax
facf(n) = dble(2*n+1)/dble(n*(n+1))
facb(n) = facf(n)
if (mod(n,2) .eq. 1) facb(n) = -facb(n)
enddo
nfac = nmax
endif
!************************************************************************
!* Calculate extinction and scattering cross section *
!* Use the convergence criterion to determine the number of terms that *
!* will later be used in the mie sum for the scattering matrix itself *
!************************************************************************
Cextsum = 0.D0
Cscasum = 0.D0
nstop = nmax
do n=1, nmax
aux = (2.D0*dble(n)+1.D0) * &
dabs(dble(an(n)*conjg(an(n)) + bn(n)*conjg(bn(n))))
Cscasum = Cscasum + aux
Cextsum = Cextsum + (2.D0*n+1.D0)*dble(an(n)+bn(n))
if (aux .lt. delta) then
nstop = n
goto 53
endif
enddo
53 nfou = nstop
if (nfou .ge. nmax) then
write(*,*) ' WARNING from scatmat : Mie sum not converged for', &
' scattering cross section'
write(*,*) ' radius r = ',r(i),' sizeparameter x = ',x, &
' sizedistribution nr. ',idis
write(*,*) ' Re(m) = ',dble(m),' Im(m) = ',dimag(m)
write(*,*) ' a priori estimate of number of Mie terms:',nmax
write(*,*) ' term ',nmax,' for Csca was ',aux
write(*,*) ' should have been less than ',delta
write(*,*) ' the a priori estimate will be used'
endif
!************************************************************************
!* Only for the first run through the loop set points in u= dcos(th) *
!************************************************************************
if ((l.eq.nsub) .and. (i.eq.ngauss)) then
!************************************************************************
!* In case of expansion in GSF : set Gauss points for dcos(th) *
!************************************************************************