-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmiesubs.F
executable file
·1586 lines (1586 loc) · 59.3 KB
/
miesubs.F
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 miev0 ( xx, crefin, perfct, mimcut, anyang,
$ numang, xmu, nmom, ipolzn, momdim, prnt,
$ qext, qsca, gqsc, pmom, sforw, sback, s1,
$ s2, tforw, tback )
c
c computes mie scattering and extinction efficiencies; asymmetry
c factor; forward- and backscatter amplitude; scattering
c amplitudes for incident polarization parallel and perpendicular
c to the plane of scattering, as functions of scattering angle;
c coefficients in the legendre polynomial expansions of either the
c unpolarized phase function or the polarized phase matrix;
c and some quantities needed in polarized radiative transfer.
c
c calls : biga, ckinmi, small1, small2, testmi, miprnt,
c lpcoef, errmsg
c
c i n t e r n a l v a r i a b l e s
c -----------------------------------
c
c an,bn mie coefficients little-a-sub-n, little-b-sub-n
c ( ref. 1, eq. 16 )
c anm1,bnm1 mie coefficients little-a-sub-(n-1),
c little-b-sub-(n-1); used in -gqsc- sum
c anp coeffs. in s+ expansion ( ref. 2, p. 1507 )
c bnp coeffs. in s- expansion ( ref. 2, p. 1507 )
c anpm coeffs. in s+ expansion ( ref. 2, p. 1507 )
c when mu is replaced by - mu
c bnpm coeffs. in s- expansion ( ref. 2, p. 1507 )
c when mu is replaced by - mu
c calcmo(k) true, calculate moments for k-th phase quantity
c (derived from -ipolzn-; used only in 'lpcoef')
c cbiga(n) bessel function ratio capital-a-sub-n (ref. 2, eq. 2)
c ( complex version )
c cior complex index of refraction with negative
c imaginary part (van de hulst convention)
c cioriv 1 / cior
c coeff ( 2n + 1 ) / ( n ( n + 1 ) )
c fn floating point version of index in loop performing
c mie series summation
c lita,litb(n) mie coefficients -an-, -bn-, saved in arrays for
c use in calculating legendre moments *pmom*
c maxtrm max. possible no. of terms in mie series
c mm + 1 and - 1, alternately.
c mim magnitude of imaginary refractive index
c mre real part of refractive index
c maxang max. possible value of input variable -numang-
c nangd2 (numang+1)/2 ( no. of angles in 0-90 deg; anyang=f )
c noabs true, sphere non-absorbing (determined by -mimcut-)
c np1dn ( n + 1 ) / n
c npquan highest-numbered phase quantity for which moments are
c to be calculated (the largest digit in -ipolzn-
c if ipolzn .ne. 0)
c ntrm no. of terms in mie series
c pass1 true on first entry, false thereafter; for self-test
c pin(j) angular function little-pi-sub-n ( ref. 2, eq. 3 )
c at j-th angle
c pinm1(j) little-pi-sub-(n-1) ( see -pin- ) at j-th angle
c psinm1 ricatti-bessel function psi-sub-(n-1), argument -xx-
c psin ricatti-bessel function psi-sub-n of argument -xx-
c ( ref. 1, p. 11 ff. )
c rbiga(n) bessel function ratio capital-a-sub-n (ref. 2, eq. 2)
c ( real version, for when imag refrac index = 0 )
c rioriv 1 / mre
c rn 1 / n
c rtmp (real) temporary variable
c sp(j) s+ for j-th angle ( ref. 2, p. 1507 )
c sm(j) s- for j-th angle ( ref. 2, p. 1507 )
c sps(j) s+ for (numang+1-j)-th angle ( anyang=false )
c sms(j) s- for (numang+1-j)-th angle ( anyang=false )
c taun angular function little-tau-sub-n ( ref. 2, eq. 4 )
c at j-th angle
c tcoef n ( n+1 ) ( 2n+1 ) (for summing tforw,tback series)
c twonp1 2n + 1
c yesang true if scattering amplitudes are to be calculated
c zetnm1 ricatti-bessel function zeta-sub-(n-1) of argument
c -xx- ( ref. 2, eq. 17 )
c zetn ricatti-bessel function zeta-sub-n of argument -xx-
c
c ----------------------------------------------------------------------
c -------- i / o specifications for subroutine miev0 -----------------
c ----------------------------------------------------------------------
implicit none
logical anyang, perfct, prnt(*)
integer ipolzn, momdim, numang, nmom
real*8 gqsc, mimcut, pmom( 0:momdim, * ), qext, qsca,
$ xmu(*), xx
complex*16 crefin, sforw, sback, s1(*), s2(*), tforw(*),
$ tback(*)
integer maxang,mxang2,maxtrm
real*8 onethr
c ----------------------------------------------------------------------
c
parameter ( maxang = 501, mxang2 = maxang/2 + 1 )
c
c ** note -- maxtrm = 10100 is neces-
c ** sary to do some of the test probs,
c ** but 1100 is sufficient for most
c ** conceivable applications
parameter ( maxtrm = 1100 )
parameter ( onethr = 1./3. )
c
logical anysav, calcmo(4), noabs, ok, pass1, persav, yesang
integer npquan
integer i,j,n,nmosav,iposav,numsav,ntrm,nangd2
real*8 mim, mimsav, mre, mm, np1dn
real*8 rioriv,xmusav,xxsav,sq,fn,rn,twonp1,tcoef, coeff
real*8 xinv,psinm1,chinm1,psin,chin,rtmp,taun
real*8 rbiga( maxtrm ), pin( maxang ), pinm1( maxang )
complex*16 an, bn, anm1, bnm1, anp, bnp, anpm, bnpm, cresav,
$ cior, cioriv, ctmp, zet, zetnm1, zetn
complex*16 cbiga( maxtrm ), lita( maxtrm ), litb( maxtrm ),
$ sp( maxang ), sm( maxang ), sps( mxang2 ), sms( mxang2 )
equivalence ( cbiga, rbiga )
save pass1
sq( ctmp ) = dble( ctmp )**2 + dimag( ctmp )**2
data pass1 / .true. /
c
c
if ( pass1 ) then
c ** save certain user input values
xxsav = xx
cresav = crefin
mimsav = mimcut
persav = perfct
anysav = anyang
nmosav = nmom
iposav = ipolzn
numsav = numang
xmusav = xmu( 1 )
c ** reset input values for test case
xx = 10.0
crefin = ( 1.5, - 0.1 )
perfct = .false.
mimcut = 0.0
anyang = .true.
numang = 1
xmu( 1 )= - 0.7660444
nmom = 1
ipolzn = - 1
c
end if
c ** check input and calculate
c ** certain variables from input
c
10 call ckinmi( numang, maxang, xx, perfct, crefin, momdim,
$ nmom, ipolzn, anyang, xmu, calcmo, npquan )
c
if ( perfct .and. xx .le. 0.1 ) then
c ** use totally-reflecting
c ** small-particle limit
c
call small1 ( xx, numang, xmu, qext, qsca, gqsc, sforw,
$ sback, s1, s2, tforw, tback, lita, litb )
ntrm = 2
go to 200
end if
c
if ( .not.perfct ) then
c
cior = crefin
if ( dimag( cior ) .gt. 0.0 ) cior = dconjg( cior )
mre = dble( cior )
mim = - dimag( cior )
noabs = mim .le. mimcut
cioriv = 1.0 / cior
rioriv = 1.0 / mre
c
if ( xx * dmax1( 1.d0, cdabs(cior) ) .le. 0.d1 ) then
c
c ** use general-refractive-index
c ** small-particle limit
c ** ( ref. 2, p. 1508 )
c
call small2 ( xx, cior, .not.noabs, numang, xmu, qext,
$ qsca, gqsc, sforw, sback, s1, s2, tforw,
$ tback, lita, litb )
ntrm = 2
go to 200
end if
c
end if
c
nangd2 = ( numang + 1 ) / 2
yesang = numang .gt. 0
c ** estimate number of terms in mie series
c ** ( ref. 2, p. 1508 )
if ( xx.le.8.0 ) then
ntrm = xx + 4. * xx**onethr + 1.
else if ( xx.lt.4200. ) then
ntrm = xx + 4.05 * xx**onethr + 2.
else
ntrm = xx + 4. * xx**onethr + 2.
end if
if ( ntrm+1 .gt. maxtrm )
$ call errmsg( 'miev0--parameter maxtrm too small', .true. )
c
c ** calculate logarithmic derivatives of
c ** j-bessel-fcn., big-a-sub-(1 to ntrm)
if ( .not.perfct )
$ call biga( cior, xx, ntrm, noabs, yesang, rbiga, cbiga )
c
c ** initialize ricatti-bessel functions
c ** (psi,chi,zeta)-sub-(0,1) for upward
c ** recurrence ( ref. 1, eq. 19 )
xinv = 1.0 / xx
psinm1 = dsin( xx )
chinm1 = dcos( xx )
psin = psinm1 * xinv - chinm1
chin = chinm1 * xinv + psinm1
zetnm1 = dcmplx( psinm1, chinm1 )
zetn = dcmplx( psin, chin )
c ** initialize previous coeffi-
c ** cients for -gqsc- series
anm1 = ( 0.0, 0.0 )
bnm1 = ( 0.0, 0.0 )
c ** initialize angular function little-pi
c ** and sums for s+, s- ( ref. 2, p. 1507 )
if ( anyang ) then
do 60 j = 1, numang
pinm1( j ) = 0.0
pin( j ) = 1.0
sp ( j ) = ( 0.0, 0.0 )
sm ( j ) = ( 0.0, 0.0 )
60 continue
else
do 70 j = 1, nangd2
pinm1( j ) = 0.0
pin( j ) = 1.0
sp ( j ) = ( 0.0, 0.0 )
sm ( j ) = ( 0.0, 0.0 )
sps( j ) = ( 0.0, 0.0 )
sms( j ) = ( 0.0, 0.0 )
70 continue
end if
c ** initialize mie sums for efficiencies, etc.
qsca = 0.0
gqsc = 0.0
sforw = ( 0., 0. )
sback = ( 0., 0. )
tforw( 1 ) = ( 0., 0. )
tback( 1 ) = ( 0., 0. )
c
c
c --------- loop to sum mie series -----------------------------------
c
mm = + 1.0
do 100 n = 1, ntrm
c ** compute various numerical coefficients
fn = n
rn = 1.0 / fn
np1dn = 1.0 + rn
twonp1 = 2 * n + 1
coeff = twonp1 / ( fn * ( n + 1 ) )
tcoef = twonp1 * ( fn * ( n + 1 ) )
c
c ** calculate mie series coefficients
if ( perfct ) then
c ** totally-reflecting case
c
an = ( ( fn*xinv ) * psin - psinm1 ) /
$ ( ( fn*xinv ) * zetn - zetnm1 )
bn = psin / zetn
c
else if ( noabs ) then
c ** no-absorption case
c
an = ( ( rioriv*rbiga(n) + ( fn*xinv ) ) * psin - psinm1 )
$ / ( ( rioriv*rbiga(n) + ( fn*xinv ) ) * zetn - zetnm1 )
bn = ( ( mre * rbiga(n) + ( fn*xinv ) ) * psin - psinm1 )
$ / ( ( mre * rbiga(n) + ( fn*xinv ) ) * zetn - zetnm1 )
else
c ** absorptive case
c
an = ( ( cioriv * cbiga(n) + ( fn*xinv ) ) * psin - psinm1 )
$ /( ( cioriv * cbiga(n) + ( fn*xinv ) ) * zetn - zetnm1 )
bn = ( ( cior * cbiga(n) + ( fn*xinv ) ) * psin - psinm1 )
$ /( ( cior * cbiga(n) + ( fn*xinv ) ) * zetn - zetnm1 )
qsca = qsca + twonp1 * ( sq( an ) + sq( bn ) )
c
end if
c ** save mie coefficients for *pmom* calculation
lita( n ) = an
litb( n ) = bn
c ** increment mie sums for non-angle-
c ** dependent quantities
c
sforw = sforw + twonp1 * ( an + bn )
tforw( 1 ) = tforw( 1 ) + tcoef * ( an - bn )
sback = sback + ( mm * twonp1 ) * ( an - bn )
tback( 1 ) = tback( 1 ) + ( mm * tcoef ) * ( an + bn )
gqsc = gqsc + ( fn - rn ) * dble( anm1 * dconjg( an )
$ + bnm1 * dconjg( bn ) )
$ + coeff * dble( an * dconjg( bn ) )
c
if ( yesang ) then
c ** put mie coefficients in form
c ** needed for computing s+, s-
c ** ( ref. 2, p. 1507 )
anp = coeff * ( an + bn )
bnp = coeff * ( an - bn )
c ** increment mie sums for s+, s-
c ** while upward recursing
c ** angular functions little pi
c ** and little tau
if ( anyang ) then
c ** arbitrary angles
c
c ** vectorizable loop
do 80 j = 1, numang
rtmp = ( xmu( j ) * pin( j ) ) - pinm1( j )
taun = fn * rtmp - pinm1( j )
sp( j ) = sp( j ) + anp * ( pin( j ) + taun )
sm( j ) = sm( j ) + bnp * ( pin( j ) - taun )
pinm1( j ) = pin( j )
pin( j ) = ( xmu( j ) * pin( j ) ) + np1dn * rtmp
80 continue
c
else
c ** angles symmetric about 90 degrees
anpm = mm * anp
bnpm = mm * bnp
c ** vectorizable loop
do 90 j = 1, nangd2
rtmp = ( xmu( j ) * pin( j ) ) - pinm1( j )
taun = fn * rtmp - pinm1( j )
sp ( j ) = sp ( j ) + anp * ( pin( j ) + taun )
sms( j ) = sms( j ) + bnpm * ( pin( j ) + taun )
sm ( j ) = sm ( j ) + bnp * ( pin( j ) - taun )
sps( j ) = sps( j ) + anpm * ( pin( j ) - taun )
pinm1( j ) = pin( j )
pin( j ) = ( xmu( j ) * pin( j ) ) + np1dn * rtmp
90 continue
c
end if
end if
c ** update relevant quantities for next
c ** pass through loop
mm = - mm
anm1 = an
bnm1 = bn
c ** upward recurrence for ricatti-bessel
c ** functions ( ref. 1, eq. 17 )
c
zet = ( twonp1 * xinv ) * zetn - zetnm1
zetnm1 = zetn
zetn = zet
psinm1 = psin
psin = dble( zetn )
100 continue
c
c ---------- end loop to sum mie series --------------------------------
c
c
qext = 2. / xx**2 * dble( sforw )
if ( perfct .or. noabs ) then
qsca = qext
else
qsca = 2. / xx**2 * qsca
end if
c
gqsc = 4. / xx**2 * gqsc
sforw = 0.5 * sforw
sback = 0.5 * sback
tforw( 2 ) = 0.5 * ( sforw + 0.25 * tforw( 1 ) )
tforw( 1 ) = 0.5 * ( sforw - 0.25 * tforw( 1 ) )
tback( 2 ) = 0.5 * ( sback + 0.25 * tback( 1 ) )
tback( 1 ) = 0.5 * ( - sback + 0.25 * tback( 1 ) )
c
if ( yesang ) then
c ** recover scattering amplitudes
c ** from s+, s- ( ref. 1, eq. 11 )
if ( anyang ) then
c ** vectorizable loop
do 110 j = 1, numang
s1( j ) = 0.5 * ( sp( j ) + sm( j ) )
s2( j ) = 0.5 * ( sp( j ) - sm( j ) )
110 continue
c
else
c ** vectorizable loop
do 120 j = 1, nangd2
s1( j ) = 0.5 * ( sp( j ) + sm( j ) )
s2( j ) = 0.5 * ( sp( j ) - sm( j ) )
120 continue
c ** vectorizable loop
do 130 j = 1, nangd2
s1( numang+1 - j ) = 0.5 * ( sps( j ) + sms( j ) )
s2( numang+1 - j ) = 0.5 * ( sps( j ) - sms( j ) )
130 continue
end if
c
end if
c ** calculate legendre moments
200 if ( nmom.gt.0 )
$ call lpcoef ( ntrm, nmom, ipolzn, momdim, calcmo, npquan,
$ lita, litb, pmom )
c
if ( dimag(crefin) .gt. 0.0 ) then
c ** take complex conjugates
c ** of scattering amplitudes
sforw = dconjg( sforw )
sback = dconjg( sback )
do 210 i = 1, 2
tforw( i ) = dconjg( tforw(i) )
tback( i ) = dconjg( tback(i) )
210 continue
c
do 220 j = 1, numang
s1( j ) = dconjg( s1(j) )
s2( j ) = dconjg( s2(j) )
220 continue
c
end if
c
if ( pass1 ) then
c ** compare test case results with
c ** correct answers and abort if bad
c
call testmi ( qext, qsca, gqsc, sforw, sback, s1, s2,
$ tforw, tback, pmom, momdim, ok )
if ( .not. ok ) then
prnt(1) = .false.
prnt(2) = .false.
call miprnt( prnt, xx, perfct, crefin, numang, xmu, qext,
$ qsca, gqsc, nmom, ipolzn, momdim, calcmo,
$ pmom, sforw, sback, tforw, tback, s1, s2 )
call errmsg( 'miev0 -- self-test failed', .true. )
end if
c ** restore user input values
xx = xxsav
crefin = cresav
mimcut = mimsav
perfct = persav
anyang = anysav
nmom = nmosav
ipolzn = iposav
numang = numsav
xmu(1) = xmusav
pass1 = .false.
go to 10
c
end if
c
if ( prnt(1) .or. prnt(2) )
$ call miprnt( prnt, xx, perfct, crefin, numang, xmu, qext,
$ qsca, gqsc, nmom, ipolzn, momdim, calcmo,
$ pmom, sforw, sback, tforw, tback, s1, s2 )
c
return
c
end
subroutine ckinmi( numang, maxang, xx, perfct, crefin, momdim,
$ nmom, ipolzn, anyang, xmu, calcmo, npquan )
c
c check for bad input to 'miev0' and calculate -calcmo,npquan-
c
implicit none
logical perfct, anyang, calcmo(*)
integer numang, maxang, momdim, nmom, ipolzn, npquan
real*8 xx, xmu(*)
integer i,l,j,ip
complex*16 crefin
c
character*4 string
logical inperr
c
inperr = .false.
c
if ( numang.gt.maxang ) then
call errmsg( 'miev0--parameter maxang too small', .true. )
inperr = .true.
end if
if ( numang.lt.0 ) call wrtbad( 'numang', inperr )
if ( xx.lt.0. ) call wrtbad( 'xx', inperr )
if ( .not.perfct .and. dble(crefin).le.0. )then
print *,'crefin=',crefin
call wrtbad( 'crefin', inperr )
endif
if ( momdim.lt.1 ) call wrtbad( 'momdim', inperr )
c
if ( nmom.ne.0 ) then
if ( nmom.lt.0 .or. nmom.gt.momdim ) call wrtbad('nmom',inperr)
if ( iabs(ipolzn).gt.4444 ) call wrtbad( 'ipolzn', inperr )
npquan = 0
do 5 l = 1, 4
calcmo( l ) = .false.
5 continue
if ( ipolzn.ne.0 ) then
c ** parse out -ipolzn- into its digits
c ** to find which phase quantities are
c ** to have their moments calculated
c
write( string, '(i4)' ) iabs(ipolzn)
do 10 j = 1, 4
ip = ichar( string(j:j) ) - ichar( '0' )
if ( ip.ge.1 .and. ip.le.4 ) calcmo( ip ) = .true.
if ( ip.eq.0 .or. (ip.ge.5 .and. ip.le.9) )
$ call wrtbad( 'ipolzn', inperr )
npquan = max0( npquan, ip )
10 continue
end if
end if
c
if ( anyang ) then
c ** allow for slight imperfections in
c ** computation of cosine
do 20 i = 1, numang
if ( xmu(i).lt.-1.00001 .or. xmu(i).gt.1.00001 )
$ call wrtbad( 'xmu', inperr )
20 continue
else
do 22 i = 1, ( numang + 1 ) / 2
if ( xmu(i).lt.-0.00001 .or. xmu(i).gt.1.00001 )
$ call wrtbad( 'xmu', inperr )
22 continue
end if
c
if ( inperr )
$ call errmsg( 'miev0--input error(s). aborting...', .true. )
c
if ( xx.gt.20000.0 .or. dble(crefin).gt.10.0 .or.
$ dabs( dimag(crefin) ).gt.10.0 ) call errmsg(
$ 'miev0--xx or crefin outside tested range', .false. )
c
return
end
subroutine lpcoef ( ntrm, nmom, ipolzn, momdim, calcmo, npquan,
$ a, b, pmom )
c
c calculate legendre polynomial expansion coefficients (also
c called moments) for phase quantities ( ref. 5 formulation )
c
c input: ntrm number terms in mie series
c nmom, ipolzn, momdim 'miev0' arguments
c calcmo flags calculated from -ipolzn-
c npquan defined in 'miev0'
c a, b mie series coefficients
c
c output: pmom legendre moments ('miev0' argument)
c
c *** notes ***
c
c (1) eqs. 2-5 are in error in dave, appl. opt. 9,
c 1888 (1970). eq. 2 refers to m1, not m2; eq. 3 refers to
c m2, not m1. in eqs. 4 and 5, the subscripts on the second
c term in square brackets should be interchanged.
c
c (2) the general-case logic in this subroutine works correctly
c in the two-term mie series case, but subroutine 'lpco2t'
c is called instead, for speed.
c
c (3) subroutine 'lpco1t', to do the one-term case, is never
c called within the context of 'miev0', but is included for
c complete generality.
c
c (4) some improvement in speed is obtainable by combining the
c 310- and 410-loops, if moments for both the third and fourth
c phase quantities are desired, because the third phase quantity
c is the real part of a complex series, while the fourth phase
c quantity is the imaginary part of that very same series. but
c most users are not interested in the fourth phase quantity,
c which is related to circular polarization, so the present
c scheme is usually more efficient.
c
implicit none
logical calcmo(*)
integer ipolzn, momdim, nmom, ntrm, npquan
real*8 pmom( 0:momdim, * )
complex*16 a(*), b(*)
c
c ** specification of local variables
c
c am(m) numerical coefficients a-sub-m-super-l
c in dave, eqs. 1-15, as simplified in ref. 5.
c
c bi(i) numerical coefficients b-sub-i-super-l
c in dave, eqs. 1-15, as simplified in ref. 5.
c
c bidel(i) 1/2 bi(i) times factor capital-del in dave
c
c cm,dm() arrays c and d in dave, eqs. 16-17 (mueller form),
c calculated using recurrence derived in ref. 5
c
c cs,ds() arrays c and d in ref. 4, eqs. a5-a6 (sekera form),
c calculated using recurrence derived in ref. 5
c
c c,d() either -cm,dm- or -cs,ds-, depending on -ipolzn-
c
c evenl true for even-numbered moments; false otherwise
c
c idel 1 + little-del in dave
c
c maxtrm max. no. of terms in mie series
c
c maxmom max. no. of non-zero moments
c
c nummom number of non-zero moments
c
c recip(k) 1 / k
c
integer maxtrm,maxmom,mxmom2,maxrcp
parameter ( maxtrm = 1102, maxmom = 2*maxtrm, mxmom2 = maxmom/2,
$ maxrcp = 4*maxtrm + 2 )
real*8 am( 0:maxtrm ), bi( 0:mxmom2 ), bidel( 0:mxmom2 ),
$ recip( maxrcp )
complex*16 cm( maxtrm ), dm( maxtrm ), cs( maxtrm ), ds( maxtrm ),
$ c( maxtrm ), d( maxtrm )
integer k,j,l,nummom,ld2,idel,m,i,mmax,imax
real*8 sum
equivalence ( c, cm ), ( d, dm )
logical pass1, evenl
save pass1, recip
data pass1 / .true. /
c
c
if ( pass1 ) then
c
do 1 k = 1, maxrcp
recip( k ) = 1.0 / k
1 continue
pass1 = .false.
c
end if
c
do 5 j = 1, max0( 1, npquan )
do 5 l = 0, nmom
pmom( l, j ) = 0.0
5 continue
c
if ( ntrm.eq.1 ) then
call lpco1t ( nmom, ipolzn, momdim, calcmo, a, b, pmom )
return
else if ( ntrm.eq.2 ) then
call lpco2t ( nmom, ipolzn, momdim, calcmo, a, b, pmom )
return
end if
c
if ( ntrm+2 .gt. maxtrm )
$ call errmsg( 'lpcoef--parameter maxtrm too small', .true. )
c
c ** calculate mueller c, d arrays
cm( ntrm+2 ) = ( 0., 0. )
dm( ntrm+2 ) = ( 0., 0. )
cm( ntrm+1 ) = ( 1. - recip( ntrm+1 ) ) * b( ntrm )
dm( ntrm+1 ) = ( 1. - recip( ntrm+1 ) ) * a( ntrm )
cm( ntrm ) = ( recip(ntrm) + recip(ntrm+1) ) * a( ntrm )
$ + ( 1. - recip(ntrm) ) * b( ntrm-1 )
dm( ntrm ) = ( recip(ntrm) + recip(ntrm+1) ) * b( ntrm )
$ + ( 1. - recip(ntrm) ) * a( ntrm-1 )
c
do 10 k = ntrm-1, 2, -1
cm( k ) = cm( k+2 ) - ( 1. + recip(k+1) ) * b( k+1 )
$ + ( recip(k) + recip(k+1) ) * a( k )
$ + ( 1. - recip(k) ) * b( k-1 )
dm( k ) = dm( k+2 ) - ( 1. + recip(k+1) ) * a( k+1 )
$ + ( recip(k) + recip(k+1) ) * b( k )
$ + ( 1. - recip(k) ) * a( k-1 )
10 continue
cm( 1 ) = cm( 3 ) + 1.5 * ( a( 1 ) - b( 2 ) )
dm( 1 ) = dm( 3 ) + 1.5 * ( b( 1 ) - a( 2 ) )
c
if ( ipolzn.ge.0 ) then
c
do 20 k = 1, ntrm + 2
c( k ) = ( 2*k - 1 ) * cm( k )
d( k ) = ( 2*k - 1 ) * dm( k )
20 continue
c
else
c ** compute sekera c and d arrays
cs( ntrm+2 ) = ( 0., 0. )
ds( ntrm+2 ) = ( 0., 0. )
cs( ntrm+1 ) = ( 0., 0. )
ds( ntrm+1 ) = ( 0., 0. )
c
do 30 k = ntrm, 1, -1
cs( k ) = cs( k+2 ) + ( 2*k + 1 ) * ( cm( k+1 ) - b( k ) )
ds( k ) = ds( k+2 ) + ( 2*k + 1 ) * ( dm( k+1 ) - a( k ) )
30 continue
c
do 40 k = 1, ntrm + 2
c( k ) = ( 2*k - 1 ) * cs( k )
d( k ) = ( 2*k - 1 ) * ds( k )
40 continue
c
end if
c
c
if( ipolzn.lt.0 ) nummom = min0( nmom, 2*ntrm - 2 )
if( ipolzn.ge.0 ) nummom = min0( nmom, 2*ntrm )
if ( nummom .gt. maxmom )
$ call errmsg( 'lpcoef--parameter maxtrm too small', .true. )
c
c ** loop over moments
do 500 l = 0, nummom
ld2 = l / 2
evenl = mod( l,2 ) .eq. 0
c ** calculate numerical coefficients
c ** a-sub-m and b-sub-i in dave
c ** double-sums for moments
if( l.eq.0 ) then
c
idel = 1
do 60 m = 0, ntrm
am( m ) = 2.0 * recip( 2*m + 1 )
60 continue
bi( 0 ) = 1.0
c
else if( evenl ) then
c
idel = 1
do 70 m = ld2, ntrm
am( m ) = ( 1. + recip( 2*m-l+1 ) ) * am( m )
70 continue
do 75 i = 0, ld2-1
bi( i ) = ( 1. - recip( l-2*i ) ) * bi( i )
75 continue
bi( ld2 ) = ( 2. - recip( l ) ) * bi( ld2-1 )
c
else
c
idel = 2
do 80 m = ld2, ntrm
am( m ) = ( 1. - recip( 2*m+l+2 ) ) * am( m )
80 continue
do 85 i = 0, ld2
bi( i ) = ( 1. - recip( l+2*i+1 ) ) * bi( i )
85 continue
c
end if
c ** establish upper limits for sums
c ** and incorporate factor capital-
c ** del into b-sub-i
mmax = ntrm - idel
if( ipolzn.ge.0 ) mmax = mmax + 1
imax = min0( ld2, mmax - ld2 )
if( imax.lt.0 ) go to 600
do 90 i = 0, imax
bidel( i ) = bi( i )
90 continue
if( evenl ) bidel( 0 ) = 0.5 * bidel( 0 )
c
c ** perform double sums just for
c ** phase quantities desired by user
if( ipolzn.eq.0 ) then
c
do 110 i = 0, imax
c ** vectorizable loop (cray)
sum = 0.0
do 100 m = ld2, mmax - i
sum = sum + am( m ) *
$ ( dble( c(m-i+1) * dconjg( c(m+i+idel) ) )
$ + dble( d(m-i+1) * dconjg( d(m+i+idel) ) ) )
100 continue
pmom( l,1 ) = pmom( l,1 ) + bidel( i ) * sum
110 continue
pmom( l,1 ) = 0.5 * pmom( l,1 )
go to 500
c
end if
c
if ( calcmo(1) ) then
do 160 i = 0, imax
c ** vectorizable loop (cray)
sum = 0.0
do 150 m = ld2, mmax - i
sum = sum + am( m ) *
$ dble( c(m-i+1) * dconjg( c(m+i+idel) ) )
150 continue
pmom( l,1 ) = pmom( l,1 ) + bidel( i ) * sum
160 continue
end if
c
c
if ( calcmo(2) ) then
do 210 i = 0, imax
c ** vectorizable loop (cray)
sum = 0.0
do 200 m = ld2, mmax - i
sum = sum + am( m ) *
$ dble( d(m-i+1) * dconjg( d(m+i+idel) ) )
200 continue
pmom( l,2 ) = pmom( l,2 ) + bidel( i ) * sum
210 continue
end if
c
c
if ( calcmo(3) ) then
do 310 i = 0, imax
c ** vectorizable loop (cray)
sum = 0.0
do 300 m = ld2, mmax - i
sum = sum + am( m ) *
$ ( dble( c(m-i+1) * dconjg( d(m+i+idel) ) )
$ + dble( c(m+i+idel) * dconjg( d(m-i+1) ) ) )
300 continue
pmom( l,3 ) = pmom( l,3 ) + bidel( i ) * sum
310 continue
pmom( l,3 ) = 0.5 * pmom( l,3 )
end if
c
c
if ( calcmo(4) ) then
do 410 i = 0, imax
c ** vectorizable loop (cray)
sum = 0.0
do 400 m = ld2, mmax - i
sum = sum + am( m ) *
$ ( dimag( c(m-i+1) * dconjg( d(m+i+idel) ) )
$ + dimag( c(m+i+idel) * dconjg( d(m-i+1) ) ))
400 continue
pmom( l,4 ) = pmom( l,4 ) + bidel( i ) * sum
410 continue
pmom( l,4 ) = - 0.5 * pmom( l,4 )
end if
c
500 continue
c
c
600 return
end
subroutine lpco1t ( nmom, ipolzn, momdim, calcmo, a, b, pmom )
c
c calculate legendre polynomial expansion coefficients (also
c called moments) for phase quantities in special case where
c no. terms in mie series = 1
c
c input: nmom, ipolzn, momdim 'miev0' arguments
c calcmo flags calculated from -ipolzn-
c a(1), b(1) mie series coefficients
c
c output: pmom legendre moments
c
implicit none
logical calcmo(*)
integer ipolzn, momdim, nmom,nummom,l
real*8 pmom( 0:momdim, * ),sq,a1sq,b1sq
complex*16 a(*), b(*), ctmp, a1b1c
sq( ctmp ) = dble( ctmp )**2 + dimag( ctmp )**2
c
c
a1sq = sq( a(1) )
b1sq = sq( b(1) )
a1b1c = a(1) * dconjg( b(1) )
c
if( ipolzn.lt.0 ) then
c
if( calcmo(1) ) pmom( 0,1 ) = 2.25 * b1sq
if( calcmo(2) ) pmom( 0,2 ) = 2.25 * a1sq
if( calcmo(3) ) pmom( 0,3 ) = 2.25 * dble( a1b1c )
if( calcmo(4) ) pmom( 0,4 ) = 2.25 *dimag( a1b1c )
c
else
c
nummom = min0( nmom, 2 )
c ** loop over moments
do 100 l = 0, nummom
c
if( ipolzn.eq.0 ) then
if( l.eq.0 ) pmom( l,1 ) = 1.5 * ( a1sq + b1sq )
if( l.eq.1 ) pmom( l,1 ) = 1.5 * dble( a1b1c )
if( l.eq.2 ) pmom( l,1 ) = 0.15 * ( a1sq + b1sq )
go to 100
end if
c
if( calcmo(1) ) then
if( l.eq.0 ) pmom( l,1 ) = 2.25 * ( a1sq + b1sq / 3. )
if( l.eq.1 ) pmom( l,1 ) = 1.5 * dble( a1b1c )
if( l.eq.2 ) pmom( l,1 ) = 0.3 * b1sq
end if
c
if( calcmo(2) ) then
if( l.eq.0 ) pmom( l,2 ) = 2.25 * ( b1sq + a1sq / 3. )
if( l.eq.1 ) pmom( l,2 ) = 1.5 * dble( a1b1c )
if( l.eq.2 ) pmom( l,2 ) = 0.3 * a1sq
end if
c
if( calcmo(3) ) then
if( l.eq.0 ) pmom( l,3 ) = 3.0 * dble( a1b1c )
if( l.eq.1 ) pmom( l,3 ) = 0.75 * ( a1sq + b1sq )
if( l.eq.2 ) pmom( l,3 ) = 0.3 * dble( a1b1c )
end if
c
if( calcmo(4) ) then
if( l.eq.0 ) pmom( l,4 ) = - 1.5 * dimag( a1b1c )
if( l.eq.1 ) pmom( l,4 ) = 0.0
if( l.eq.2 ) pmom( l,4 ) = 0.3 * dimag( a1b1c )
end if
c
100 continue
c
end if
c
return
end
subroutine lpco2t ( nmom, ipolzn, momdim, calcmo, a, b, pmom )
c
c calculate legendre polynomial expansion coefficients (also
c called moments) for phase quantities in special case where
c no. terms in mie series = 2
c
c input: nmom, ipolzn, momdim 'miev0' arguments
c calcmo flags calculated from -ipolzn-
c a(1-2), b(1-2) mie series coefficients
c
c output: pmom legendre moments
c
implicit none
logical calcmo(*)
integer ipolzn, momdim, nmom,l,nummom
real*8 pmom( 0:momdim, * ),sq,a1sq,b1sq,pm1,pm2,a2sq,b2sq
complex*16 a(*), b(*)
complex*16 a2c, b2c, ctmp, ca, cac, cat, cb, cbc, cbt, cg, ch
sq( ctmp ) = dble( ctmp )**2 + dimag( ctmp )**2
c
c
ca = 3. * a(1) - 5. * b(2)
cat= 3. * b(1) - 5. * a(2)
cac = dconjg( ca )
a2sq = sq( a(2) )
b2sq = sq( b(2) )
a2c = dconjg( a(2) )
b2c = dconjg( b(2) )
c
if( ipolzn.lt.0 ) then
c ** loop over sekera moments
nummom = min0( nmom, 2 )
do 50 l = 0, nummom
c
if( calcmo(1) ) then
if( l.eq.0 ) pmom( l,1 ) = 0.25 * ( sq(cat) +
$ (100./3.) * b2sq )
if( l.eq.1 ) pmom( l,1 ) = (5./3.) * dble( cat * b2c )
if( l.eq.2 ) pmom( l,1 ) = (10./3.) * b2sq
end if
c
if( calcmo(2) ) then
if( l.eq.0 ) pmom( l,2 ) = 0.25 * ( sq(ca) +
$ (100./3.) * a2sq )
if( l.eq.1 ) pmom( l,2 ) = (5./3.) * dble( ca * a2c )
if( l.eq.2 ) pmom( l,2 ) = (10./3.) * a2sq
end if
c
if( calcmo(3) ) then
if( l.eq.0 ) pmom( l,3 ) = 0.25 * dble( cat*cac +
$ (100./3.)*b(2)*a2c )
if( l.eq.1 ) pmom( l,3 ) = 5./6. * dble( b(2)*cac +
$ cat*a2c )
if( l.eq.2 ) pmom( l,3 ) = 10./3. * dble( b(2) * a2c )
end if
c
if( calcmo(4) ) then
if( l.eq.0 ) pmom( l,4 ) = -0.25 * dimag( cat*cac +
$ (100./3.)*b(2)*a2c )
if( l.eq.1 ) pmom( l,4 ) = -5./6. * dimag( b(2)*cac +
$ cat*a2c )
if( l.eq.2 ) pmom( l,4 ) = -10./3. * dimag( b(2) * a2c )
end if
c
50 continue
c
else
c
cb = 3. * b(1) + 5. * a(2)
cbt= 3. * a(1) + 5. * b(2)
cbc = dconjg( cb )
cg = ( cbc*cbt + 10.*( cac*a(2) + b2c*cat) ) / 3.
ch = 2.*( cbc*a(2) + b2c*cbt )
c
c ** loop over mueller moments
nummom = min0( nmom, 4 )
do 100 l = 0, nummom
c
if( ipolzn.eq.0 .or. calcmo(1) ) then
if( l.eq.0 ) pm1 = 0.25 * sq(ca) + sq(cb) / 12.
$ + (5./3.) * dble(ca*b2c) + 5.*b2sq
if( l.eq.1 ) pm1 = dble( cb * ( cac/6. + b2c ) )
if( l.eq.2 ) pm1 = sq(cb)/30. + (20./7.) * b2sq
$ + (2./3.) * dble( ca * b2c )
if( l.eq.3 ) pm1 = (2./7.) * dble( cb * b2c )
if( l.eq.4 ) pm1 = (40./63.) * b2sq
if ( calcmo(1) ) pmom( l,1 ) = pm1
end if
c
if( ipolzn.eq.0 .or. calcmo(2) ) then
if( l.eq.0 ) pm2 = 0.25*sq(cat) + sq(cbt) / 12.
$ + (5./3.) * dble(cat*a2c) + 5.*a2sq
if( l.eq.1 ) pm2 = dble( cbt * ( dconjg(cat)/6. + a2c) )
if( l.eq.2 ) pm2 = sq(cbt)/30. + (20./7.) * a2sq
$ + (2./3.) * dble( cat * a2c )
if( l.eq.3 ) pm2 = (2./7.) * dble( cbt * a2c )
if( l.eq.4 ) pm2 = (40./63.) * a2sq
if ( calcmo(2) ) pmom( l,2 ) = pm2
end if
c
if( ipolzn.eq.0 ) then
pmom( l,1 ) = 0.5 * ( pm1 + pm2 )
go to 100
end if
c