-
Notifications
You must be signed in to change notification settings - Fork 634
/
Copy pathradi.f90
4844 lines (4190 loc) · 182 KB
/
radi.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
MODULE SPECDATA
USE PRECISION_PARAMETERS
IMPLICIT NONE (TYPE,EXTERNAL)
INTEGER :: I,J
INTEGER, PARAMETER :: NWATERK=183
REAL(EB) :: CPLXREF_WATER(NWATERK,2)
INTEGER, PARAMETER :: NFUELK=94
REAL(EB) :: CPLXREF_FUEL(NFUELK,3)
! HALE, G.M. AND QUERRY, M.R., "OPTICAL CONSTANTS OF WATER IN THE 200-NM TO 200-\MU M WAVELENGTH REGION",
! APPLIED OPTICS, 12(3), 555-63 (1973)
DATA ((CPLXREF_WATER(I,J), J = 1,2), I = 1,10) &
/ 1.0000000E-6_EB, 2.8647890E-6_EB,&
1.0200000E-6_EB, 2.1915636E-6_EB,&
1.0400000E-6_EB, 1.3241691E-6_EB,&
1.0600000E-6_EB, 1.0122254E-6_EB,&
1.0800000E-6_EB, 1.1172677E-6_EB,&
1.1000000E-6_EB, 1.4880987E-6_EB,&
1.1200000E-6_EB, 4.6345919E-6_EB,&
1.1400000E-6_EB, 5.9874090E-6_EB,&
1.1600000E-6_EB, 8.2155782E-6_EB,&
1.1800000E-6_EB, 9.7657473E-6_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =11,20) &
/ 1.2000000E-6_EB, 9.9312684E-6_EB,&
1.2200000E-6_EB, 9.2230290E-6_EB,&
1.2400000E-6_EB, 8.6834937E-6_EB,&
1.2600000E-6_EB, 8.9238177E-6_EB,&
1.2800000E-6_EB, 9.9821980E-6_EB,&
1.3000000E-6_EB, 1.1483029E-5_EB,&
1.3200000E-6_EB, 1.1953809E-5_EB,&
1.3400000E-6_EB, 1.2614780E-5_EB,&
1.3600000E-6_EB, 2.9978533E-5_EB,&
1.3800000E-6_EB, 6.6988316E-5_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =21,30) &
/ 1.4000000E-6_EB, 1.3803508E-4_EB,&
1.4200000E-6_EB, 2.4995602E-4_EB,&
1.4400000E-6_EB, 3.3002369E-4_EB,&
1.4600000E-6_EB, 3.2996003E-4_EB,&
1.4800000E-6_EB, 2.5003560E-4_EB,&
1.5000000E-6_EB, 2.0996516E-4_EB,&
1.5200000E-6_EB, 1.6994565E-4_EB,&
1.5400000E-6_EB, 1.4497583E-4_EB,&
1.5600000E-6_EB, 1.2004433E-4_EB,&
1.5800000E-6_EB, 9.9957388E-5_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =31,40) &
/ 1.6000000E-6_EB, 8.5561825E-5_EB,&
1.6200000E-6_EB, 7.5028952E-5_EB,&
1.6400000E-6_EB, 6.4992643E-5_EB,&
1.6600000E-6_EB, 5.9972898E-5_EB,&
1.6800000E-6_EB, 6.0027012E-5_EB,&
1.7000000E-6_EB, 6.0065211E-5_EB,&
1.7200000E-6_EB, 6.9942231E-5_EB,&
1.7400000E-6_EB, 8.5017388E-5_EB,&
1.7600000E-6_EB, 1.0000023E-4_EB,&
1.7800000E-6_EB, 1.1501809E-4_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =41,50) &
/ 1.8000000E-6_EB, 1.1502128E-4_EB,&
1.8200000E-6_EB, 1.3005838E-4_EB,&
1.8400000E-6_EB, 1.4993669E-4_EB,&
1.8600000E-6_EB, 2.1003200E-4_EB,&
1.8800000E-6_EB, 4.6497435E-4_EB,&
1.9000000E-6_EB, 1.0000183E-3_EB,&
1.9200000E-6_EB, 1.7500423E-3_EB,&
1.9400000E-6_EB, 1.8499391E-3_EB,&
1.9600000E-6_EB, 1.6500261E-3_EB,&
1.9800000E-6_EB, 1.4500559E-3_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =51,60) &
/ 2.0000000E-6_EB, 1.1000790E-3_EB,&
2.0200000E-6_EB, 9.0001961E-4_EB,&
2.0400000E-6_EB, 7.3003417E-4_EB,&
2.0600000E-6_EB, 6.3998112E-4_EB,&
2.0800000E-6_EB, 5.2006742E-4_EB,&
2.1000000E-6_EB, 4.5003447E-4_EB,&
2.1200000E-6_EB, 4.0505888E-4_EB,&
2.1400000E-6_EB, 3.4995785E-4_EB,&
2.1600000E-6_EB, 3.2005422E-4_EB,&
2.1800000E-6_EB, 2.9994500E-4_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =61,70) &
/ 2.2000000E-6_EB, 2.8904129E-4_EB,&
2.2200000E-6_EB, 2.8495578E-4_EB,&
2.2400000E-6_EB, 2.9500960E-4_EB,&
2.2600000E-6_EB, 3.1005293E-4_EB,&
2.2800000E-6_EB, 3.5997028E-4_EB,&
2.3000000E-6_EB, 4.0998313E-4_EB,&
2.3200000E-6_EB, 4.9496551E-4_EB,&
2.3400000E-6_EB, 5.9494505E-4_EB,&
2.3600000E-6_EB, 6.9994116E-4_EB,&
2.3800000E-6_EB, 8.2007768E-4_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =71,80) &
/ 2.4000000E-6_EB, 9.5607557E-4_EB,&
2.4200000E-6_EB, 1.1500727E-3_EB,&
2.4400000E-6_EB, 1.2999617E-3_EB,&
2.4600000E-6_EB, 1.4999176E-3_EB,&
2.4800000E-6_EB, 1.6999912E-3_EB,&
2.5000000E-6_EB, 1.8000424E-3_EB,&
2.5200000E-6_EB, 2.0500716E-3_EB,&
2.5400000E-6_EB, 2.1999478E-3_EB,&
2.5600000E-6_EB, 2.3500946E-3_EB,&
2.5800000E-6_EB, 2.7000302E-3_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =81,90) &
/ 2.6000000E-6_EB, 3.1699367E-3_EB,&
2.6500000E-6_EB, 6.7000889E-3_EB,&
2.7000000E-6_EB, 1.8999997E-2_EB,&
2.7500000E-6_EB, 5.9000926E-2_EB,&
2.8000000E-6_EB, 1.1500027E-1_EB,&
2.8500000E-6_EB, 1.8499960E-1_EB,&
2.9000000E-6_EB, 2.6769861E-1_EB,&
2.9500000E-6_EB, 2.9813700E-1_EB,&
3.0000000E-6_EB, 2.7215495E-1_EB,&
3.0500000E-6_EB, 2.4000020E-1_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =91,100) &
/ 3.1000000E-6_EB, 1.9200142E-1_EB,&
3.1500000E-6_EB, 1.3500032E-1_EB,&
3.2000000E-6_EB, 9.2401540E-2_EB,&
3.2500000E-6_EB, 6.0999713E-2_EB,&
3.3000000E-6_EB, 3.6798931E-2_EB,&
3.3500000E-6_EB, 2.6099958E-2_EB,&
3.4000000E-6_EB, 1.9500046E-2_EB,&
3.4500000E-6_EB, 1.3199993E-2_EB,&
3.5000000E-6_EB, 9.4000888E-3_EB,&
3.6000000E-6_EB, 5.1500311E-3_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =101,110) &
/ 3.7000000E-6_EB, 3.6000769E-3_EB,&
3.8000000E-6_EB, 3.4001225E-3_EB,&
3.9000000E-6_EB, 3.7999516E-3_EB,&
4.0000000E-6_EB, 4.5998962E-3_EB,&
4.1000000E-6_EB, 5.6118033E-3_EB,&
4.2000000E-6_EB, 6.8850428E-3_EB,&
4.3000000E-6_EB, 8.4519233E-3_EB,&
4.4000000E-6_EB, 1.0294142E-2_EB,&
4.5000000E-6_EB, 1.3392888E-2_EB,&
4.6000000E-6_EB, 1.4715466E-2_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =111,120) &
/ 4.7000000E-6_EB, 1.5708593E-2_EB,&
4.8000000E-6_EB, 1.5011494E-2_EB,&
4.9000000E-6_EB, 1.3686529E-2_EB,&
5.0000000E-6_EB, 1.2414086E-2_EB,&
5.1000000E-6_EB, 1.1120156E-2_EB,&
5.2000000E-6_EB, 1.0096790E-2_EB,&
5.3000000E-6_EB, 9.7848459E-3_EB,&
5.4000000E-6_EB, 1.0313240E-2_EB,&
5.5000000E-6_EB, 1.1598416E-2_EB,&
5.6000000E-6_EB, 1.4215720E-2_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =121,130) &
/ 5.7000000E-6_EB, 2.0320903E-2_EB,&
5.8000000E-6_EB, 3.3000777E-2_EB,&
5.9000000E-6_EB, 6.2209688E-2_EB,&
6.0000000E-6_EB, 1.0699987E-1_EB,&
6.1000000E-6_EB, 1.3101555E-1_EB,&
6.2000000E-6_EB, 8.8019050E-2_EB,&
6.3000000E-6_EB, 5.7002139E-2_EB,&
6.4000000E-6_EB, 4.4868962E-2_EB,&
6.5000000E-6_EB, 3.9207820E-2_EB,&
6.6000000E-6_EB, 3.5609327E-2_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =131,140) &
/ 6.7000000E-6_EB, 3.3696285E-2_EB,&
6.8000000E-6_EB, 3.2684059E-2_EB,&
6.9000000E-6_EB, 3.2176355E-2_EB,&
7.0000000E-6_EB, 3.1974228E-2_EB,&
7.1000000E-6_EB, 3.1979003E-2_EB,&
7.2000000E-6_EB, 3.2085637E-2_EB,&
7.3000000E-6_EB, 3.2182721E-2_EB,&
7.4000000E-6_EB, 3.2388031E-2_EB,&
7.5000000E-6_EB, 3.2586975E-2_EB,&
7.6000000E-6_EB, 3.2779552E-2_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =141,150) &
/ 7.7000000E-6_EB, 3.3088313E-2_EB,&
7.8000000E-6_EB, 3.3518031E-2_EB,&
7.9000000E-6_EB, 3.3884883E-2_EB,&
8.0000000E-6_EB, 3.4313806E-2_EB,&
8.2000000E-6_EB, 3.5106397E-2_EB,&
8.4000000E-6_EB, 3.6096341E-2_EB,&
8.6000000E-6_EB, 3.1001791E-2_EB,&
8.8000000E-6_EB, 3.8515496E-2_EB,&
9.0000000E-6_EB, 3.9892186E-2_EB,&
9.2000000E-6_EB, 4.1510792E-2_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =151,160) &
/ 9.4000000E-6_EB, 4.3310835E-2_EB,&
9.6000000E-6_EB, 4.5378257E-2_EB,&
9.8000000E-6_EB, 4.7883356E-2_EB,&
1.0000000E-5_EB, 5.0770427E-2_EB,&
1.0500000E-5_EB, 6.6176625E-2_EB,&
1.1000000E-5_EB, 9.6813952E-2_EB,&
1.1500000E-5_EB, 1.4202987E-1_EB,&
1.2000000E-5_EB, 1.9900734E-1_EB,&
1.2500000E-5_EB, 2.5902467E-1_EB,&
1.3000000E-5_EB, 3.0497270E-1_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =161,170) &
/ 1.3500000E-5_EB, 3.4302267E-1_EB,&
1.4000000E-5_EB, 3.6998750E-1_EB,&
1.4500000E-5_EB, 3.8804760E-1_EB,&
1.5000000E-5_EB, 4.0202539E-1_EB,&
1.5500000E-5_EB, 4.1394609E-1_EB,&
1.6000000E-5_EB, 4.2195159E-1_EB,&
1.6500000E-5_EB, 4.2804722E-1_EB,&
1.7000000E-5_EB, 4.2897828E-1_EB,&
1.7500000E-5_EB, 4.2906183E-1_EB,&
1.8000000E-5_EB, 4.2599412E-1_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =171,180) &
/ 1.8500000E-5_EB, 4.2104440E-1_EB,&
1.9000000E-5_EB, 4.1397792E-1_EB,&
1.9500000E-5_EB, 4.0407849E-1_EB,&
2.0000000E-5_EB, 3.9295355E-1_EB,&
3.0000000E-5_EB, 3.2801834E-1_EB,&
3.8000000E-5_EB, 3.6105890E-1_EB,&
5.0000000E-5_EB, 5.1407047E-1_EB,&
6.0000000E-5_EB, 5.8680428E-1_EB,&
7.0000000E-5_EB, 5.7598174E-1_EB,&
8.0000000E-5_EB, 5.4685638E-1_EB /
DATA ((CPLXREF_WATER(I,J), J = 1,2), I =181,NWATERK) &
/ 9.0000000E-5_EB, 5.3571554E-1_EB,&
1.0000000E-4_EB, 5.3237328E-1_EB,&
2.0000000E-4_EB, 5.0452117E-1_EB /
! HEPTANE PROPERTIES FROM
! L.A. DOMBROVSKY, S.S. SAZHIN, S.V. MIKHALOVSKY, R. WOOD, M.R. HEIKAL
! SPECTRAL PROPERTIES OF DIESEL FUEL PARTICLES
! FUEL, VOL. 82, NO. 1 (2003) PP. 15-22
DATA (CPLXREF_FUEL( 1,J), J=1,3) / 0.7_EB, 1.47_EB, 1.55E-07_EB /
DATA (CPLXREF_FUEL( 2,J), J=1,3) / 0.8_EB, 1.47_EB, 4.33E-07_EB /
DATA (CPLXREF_FUEL( 3,J), J=1,3) / 0.9_EB, 1.47_EB, 9.39E-07_EB /
DATA (CPLXREF_FUEL( 4,J), J=1,3) / 1.0_EB, 1.47_EB, 8.81E-07_EB /
DATA (CPLXREF_FUEL( 5,J), J=1,3) / 1.1_EB, 1.47_EB, 1.19E-06_EB /
DATA (CPLXREF_FUEL( 6,J), J=1,3) / 1.2_EB, 1.47_EB, 2.63E-06_EB /
DATA (CPLXREF_FUEL( 7,J), J=1,3) / 1.3_EB, 1.47_EB, 1.08E-05_EB /
DATA (CPLXREF_FUEL( 8,J), J=1,3) / 1.4_EB, 1.47_EB, 2.24E-05_EB /
DATA (CPLXREF_FUEL( 9,J), J=1,3) / 1.5_EB, 1.46_EB, 3.41E-05_EB /
DATA (CPLXREF_FUEL( 10,J), J=1,3) / 1.6_EB, 1.46_EB, 4.57E-05_EB /
DATA (CPLXREF_FUEL( 11,J), J=1,3) / 1.7_EB, 1.46_EB, 5.73E-05_EB /
DATA (CPLXREF_FUEL( 12,J), J=1,3) / 1.8_EB, 1.46_EB, 6.90E-05_EB /
DATA (CPLXREF_FUEL( 13,J), J=1,3) / 1.9_EB, 1.46_EB, 8.06E-05_EB /
DATA (CPLXREF_FUEL( 14,J), J=1,3) / 2.0_EB, 1.46_EB, 9.22E-05_EB /
DATA (CPLXREF_FUEL( 15,J), J=1,3) / 2.1_EB, 1.46_EB, 1.04E-04_EB /
DATA (CPLXREF_FUEL( 16,J), J=1,3) / 2.2_EB, 1.46_EB, 1.15E-04_EB /
DATA (CPLXREF_FUEL( 17,J), J=1,3) / 2.3_EB, 1.46_EB, 1.27E-04_EB /
DATA (CPLXREF_FUEL( 18,J), J=1,3) / 2.4_EB, 1.45_EB, 1.39E-04_EB /
DATA (CPLXREF_FUEL( 19,J), J=1,3) / 2.5_EB, 1.45_EB, 1.50E-04_EB /
DATA (CPLXREF_FUEL( 20,J), J=1,3) / 2.6_EB, 1.45_EB, 1.62E-04_EB /
DATA (CPLXREF_FUEL( 21,J), J=1,3) / 2.7_EB, 1.45_EB, 1.11E-04_EB /
DATA (CPLXREF_FUEL( 22,J), J=1,3) / 2.8_EB, 1.44_EB, 5.92E-05_EB /
DATA (CPLXREF_FUEL( 23,J), J=1,3) / 2.9_EB, 1.44_EB, 7.45E-05_EB /
DATA (CPLXREF_FUEL( 24,J), J=1,3) / 3.0_EB, 1.43_EB, 9.72E-05_EB /
DATA (CPLXREF_FUEL( 25,J), J=1,3) / 3.1_EB, 1.42_EB, 3.12E-04_EB /
DATA (CPLXREF_FUEL( 26,J), J=1,3) / 3.2_EB, 1.40_EB, 6.09E-04_EB /
DATA (CPLXREF_FUEL( 27,J), J=1,3) / 3.3_EB, 1.17_EB, 5.72E-02_EB /
DATA (CPLXREF_FUEL( 28,J), J=1,3) / 3.4_EB, 1.39_EB, 1.20E-01_EB /
DATA (CPLXREF_FUEL( 29,J), J=1,3) / 3.5_EB, 1.45_EB, 8.24E-02_EB /
DATA (CPLXREF_FUEL( 30,J), J=1,3) / 3.6_EB, 1.51_EB, 1.63E-03_EB /
DATA (CPLXREF_FUEL( 31,J), J=1,3) / 3.7_EB, 1.64_EB, 1.33E-03_EB /
DATA (CPLXREF_FUEL( 32,J), J=1,3) / 3.8_EB, 1.56_EB, 1.02E-03_EB /
DATA (CPLXREF_FUEL( 33,J), J=1,3) / 3.9_EB, 1.52_EB, 6.36E-04_EB /
DATA (CPLXREF_FUEL( 34,J), J=1,3) / 4.0_EB, 1.50_EB, 2.51E-04_EB /
DATA (CPLXREF_FUEL( 35,J), J=1,3) / 4.1_EB, 1.50_EB, 2.59E-04_EB /
DATA (CPLXREF_FUEL( 36,J), J=1,3) / 4.2_EB, 1.49_EB, 3.10E-04_EB /
DATA (CPLXREF_FUEL( 37,J), J=1,3) / 4.3_EB, 1.49_EB, 2.60E-04_EB /
DATA (CPLXREF_FUEL( 38,J), J=1,3) / 4.4_EB, 1.48_EB, 2.11E-04_EB /
DATA (CPLXREF_FUEL( 39,J), J=1,3) / 4.5_EB, 1.48_EB, 1.98E-04_EB /
DATA (CPLXREF_FUEL( 40,J), J=1,3) / 4.6_EB, 1.47_EB, 1.90E-04_EB /
DATA (CPLXREF_FUEL( 41,J), J=1,3) / 4.7_EB, 1.47_EB, 1.58E-04_EB /
DATA (CPLXREF_FUEL( 42,J), J=1,3) / 4.8_EB, 1.47_EB, 1.21E-04_EB /
DATA (CPLXREF_FUEL( 43,J), J=1,3) / 4.9_EB, 1.47_EB, 8.32E-05_EB /
DATA (CPLXREF_FUEL( 44,J), J=1,3) / 5.0_EB, 1.46_EB, 4.58E-05_EB /
DATA (CPLXREF_FUEL( 45,J), J=1,3) / 5.1_EB, 1.46_EB, 6.37E-05_EB /
DATA (CPLXREF_FUEL( 46,J), J=1,3) / 5.2_EB, 1.46_EB, 7.78E-05_EB /
DATA (CPLXREF_FUEL( 47,J), J=1,3) / 5.3_EB, 1.46_EB, 7.66E-05_EB /
DATA (CPLXREF_FUEL( 48,J), J=1,3) / 5.4_EB, 1.45_EB, 7.55E-05_EB /
DATA (CPLXREF_FUEL( 49,J), J=1,3) / 5.5_EB, 1.45_EB, 9.62E-05_EB /
DATA (CPLXREF_FUEL( 50,J), J=1,3) / 5.6_EB, 1.45_EB, 1.17E-04_EB /
DATA (CPLXREF_FUEL( 51,J), J=1,3) / 5.7_EB, 1.45_EB, 1.63E-04_EB /
DATA (CPLXREF_FUEL( 52,J), J=1,3) / 5.8_EB, 1.45_EB, 2.16E-04_EB /
DATA (CPLXREF_FUEL( 53,J), J=1,3) / 5.9_EB, 1.44_EB, 2.68E-04_EB /
DATA (CPLXREF_FUEL( 54,J), J=1,3) / 6.0_EB, 1.44_EB, 3.21E-04_EB /
DATA (CPLXREF_FUEL( 55,J), J=1,3) / 6.1_EB, 1.44_EB, 3.74E-04_EB /
DATA (CPLXREF_FUEL( 56,J), J=1,3) / 6.2_EB, 1.44_EB, 4.36E-04_EB /
DATA (CPLXREF_FUEL( 57,J), J=1,3) / 6.3_EB, 1.43_EB, 5.87E-04_EB /
DATA (CPLXREF_FUEL( 58,J), J=1,3) / 6.4_EB, 1.43_EB, 7.38E-04_EB /
DATA (CPLXREF_FUEL( 59,J), J=1,3) / 6.5_EB, 1.42_EB, 1.35E-03_EB /
DATA (CPLXREF_FUEL( 60,J), J=1,3) / 6.6_EB, 1.41_EB, 6.12E-03_EB /
DATA (CPLXREF_FUEL( 61,J), J=1,3) / 6.7_EB, 1.39_EB, 2.06E-02_EB /
DATA (CPLXREF_FUEL( 62,J), J=1,3) / 6.8_EB, 1.35_EB, 3.51E-02_EB /
DATA (CPLXREF_FUEL( 63,J), J=1,3) / 6.9_EB, 1.37_EB, 2.29E-02_EB /
DATA (CPLXREF_FUEL( 64,J), J=1,3) / 7.0_EB, 1.48_EB, 3.99E-03_EB /
DATA (CPLXREF_FUEL( 65,J), J=1,3) / 7.1_EB, 1.53_EB, 3.24E-03_EB /
DATA (CPLXREF_FUEL( 66,J), J=1,3) / 7.2_EB, 1.48_EB, 2.61E-03_EB /
DATA (CPLXREF_FUEL( 67,J), J=1,3) / 7.3_EB, 1.43_EB, 2.97E-03_EB /
DATA (CPLXREF_FUEL( 68,J), J=1,3) / 7.4_EB, 1.47_EB, 3.33E-03_EB /
DATA (CPLXREF_FUEL( 69,J), J=1,3) / 7.5_EB, 1.50_EB, 2.80E-03_EB /
DATA (CPLXREF_FUEL( 70,J), J=1,3) / 7.6_EB, 1.49_EB, 2.14E-03_EB /
DATA (CPLXREF_FUEL( 71,J), J=1,3) / 7.7_EB, 1.48_EB, 2.28E-03_EB /
DATA (CPLXREF_FUEL( 72,J), J=1,3) / 7.8_EB, 1.47_EB, 2.41E-03_EB /
DATA (CPLXREF_FUEL( 73,J), J=1,3) / 7.9_EB, 1.47_EB, 1.63E-03_EB /
DATA (CPLXREF_FUEL( 74,J), J=1,3) / 8.0_EB, 1.47_EB, 9.84E-04_EB /
DATA (CPLXREF_FUEL( 75,J), J=1,3) / 8.1_EB, 1.46_EB, 9.03E-04_EB /
DATA (CPLXREF_FUEL( 76,J), J=1,3) / 8.2_EB, 1.46_EB, 8.19E-04_EB /
DATA (CPLXREF_FUEL( 77,J), J=1,3) / 8.3_EB, 1.46_EB, 7.13E-04_EB /
DATA (CPLXREF_FUEL( 78,J), J=1,3) / 8.4_EB, 1.46_EB, 6.07E-04_EB /
DATA (CPLXREF_FUEL( 79,J), J=1,3) / 8.5_EB, 1.45_EB, 8.54E-04_EB /
DATA (CPLXREF_FUEL( 80,J), J=1,3) / 8.6_EB, 1.45_EB, 1.10E-03_EB /
DATA (CPLXREF_FUEL( 81,J), J=1,3) / 8.7_EB, 1.45_EB, 1.76E-03_EB /
DATA (CPLXREF_FUEL( 82,J), J=1,3) / 8.8_EB, 1.45_EB, 2.41E-03_EB /
DATA (CPLXREF_FUEL( 83,J), J=1,3) / 8.9_EB, 1.45_EB, 1.58E-03_EB /
DATA (CPLXREF_FUEL( 84,J), J=1,3) / 9.0_EB, 1.45_EB, 5.71E-04_EB /
DATA (CPLXREF_FUEL( 85,J), J=1,3) / 9.1_EB, 1.45_EB, 8.53E-04_EB /
DATA (CPLXREF_FUEL( 86,J), J=1,3) / 9.2_EB, 1.45_EB, 1.28E-03_EB /
DATA (CPLXREF_FUEL( 87,J), J=1,3) / 9.3_EB, 1.45_EB, 1.63E-03_EB /
DATA (CPLXREF_FUEL( 88,J), J=1,3) / 9.4_EB, 1.45_EB, 1.98E-03_EB /
DATA (CPLXREF_FUEL( 89,J), J=1,3) / 9.5_EB, 1.45_EB, 1.69E-03_EB /
DATA (CPLXREF_FUEL( 90,J), J=1,3) / 9.6_EB, 1.45_EB, 1.24E-03_EB /
DATA (CPLXREF_FUEL( 91,J), J=1,3) / 9.7_EB, 1.45_EB, 1.43E-03_EB /
DATA (CPLXREF_FUEL( 92,J), J=1,3) / 9.8_EB, 1.45_EB, 1.69E-03_EB /
DATA (CPLXREF_FUEL( 93,J), J=1,3) / 9.9_EB, 1.45_EB, 1.26E-03_EB /
DATA (CPLXREF_FUEL( 94,J), J=1,3) / 10.0_EB, 1.45_EB, 8.24E-04_EB /
END MODULE SPECDATA
MODULE WSGG_ARRAYS
USE PRECISION_PARAMETERS
REAL(EB), ALLOCATABLE, DIMENSION(:,:) :: WSGG_B1_ARRAY,WSGG_B2_ARRAY,WSGG_D_ARRAY
REAL(EB), ALLOCATABLE, DIMENSION(:,:,:) :: WSGG_C_ARRAY
REAL(EB), ALLOCATABLE, DIMENSION(:) :: WSGG_KAPPAP1_ARRAY,WSGG_KAPPAP2_ARRAY
END MODULE WSGG_ARRAYS
MODULE MIEV
USE PRECISION_PARAMETERS
USE GLOBAL_CONSTANTS, ONLY: NUMBER_SPECTRAL_BANDS,NUMBER_RADIATION_ANGLES,LU_ERR
USE SPECDATA, ONLY: CPLXREF_WATER, NWATERK, CPLXREF_FUEL, NFUELK
USE RADCAL_CALC, ONLY: PLANCK
USE RADCONS
USE MEMORY_FUNCTIONS, ONLY : CHKMEMERR
USE MATH_FUNCTIONS, ONLY : INTERPOLATE1D
IMPLICIT NONE (TYPE,EXTERNAL)
REAL(EB), ALLOCATABLE :: RDMIE(:), LMBDMIE(:),LMBDWGHT(:),REAL_REF_INDX(:),CMPLX_REF_INDX(:)
REAL(EB), ALLOCATABLE :: QSCA(:,:), QABS(:,:), CHI_F(:,:)
PRIVATE
PUBLIC MEAN_CROSS_SECTIONS
CONTAINS
!> \brief Calculates the mean scattering and absorption coefficients for each radiation band and particle size group
!> \param CLASS_NUMBER Index of the Lagrangian particle class
SUBROUTINE MEAN_CROSS_SECTIONS(CLASS_NUMBER)
USE TYPES, ONLY: LAGRANGIAN_PARTICLE_CLASS_TYPE, LAGRANGIAN_PARTICLE_CLASS, TABLES_TYPE, TABLES, SPECIES_MIXTURE, SPECIES, &
SPECIES_MIXTURE_TYPE
USE GLOBAL_CONSTANTS, ONLY: H2O_INDEX
USE MATH_FUNCTIONS, ONLY : INTERPOLATE1D
INTEGER, INTENT(IN) :: CLASS_NUMBER
INTEGER :: NSB,I,J,IBND,IZERO,NX,NLAMBDALOW(1),NLAMBDAHIGH(1),ND,NL_DATA
REAL(EB) :: RMMAX,RMMIN,RDTMP,IB,IBSUM,AVAL,BVAL,ASUM,BSUM,B_WIEN
REAL(EB), ALLOCATABLE, DIMENSION(:,:) :: WQABS,WQSCA
REAL(EB), ALLOCATABLE, DIMENSION(:) :: R50
TYPE (LAGRANGIAN_PARTICLE_CLASS_TYPE), POINTER :: LPC=>NULL()
TYPE (SPECIES_MIXTURE_TYPE), POINTER :: SM2=>NULL()
TYPE (TABLES_TYPE), POINTER :: TA=>NULL()
IF (CLASS_NUMBER > 0) THEN
LPC => LAGRANGIAN_PARTICLE_CLASS(CLASS_NUMBER)
ELSE
SM2 => SPECIES_MIXTURE(-CLASS_NUMBER)
ENDIF
NSB = NUMBER_SPECTRAL_BANDS
! PHYSICAL PARAMETERS
! MINIMUM MEAN RADIUS (M)
RMMIN = 0.5_EB*1.E-6*MIE_MINIMUM_DIAMETER
IF (RMMIN < TWO_EPSILON_EB) RMMIN = 0.5E-6_EB
! MAXIMUM MEAN RADIUS (M)
RMMAX = 0.5_EB*1.E-6*MIE_MAXIMUM_DIAMETER
IF (RMMAX < TWO_EPSILON_EB) THEN
IF (CLASS_NUMBER > 0) RMMAX = 0.5_EB*LPC%DIAMETER
IF (CLASS_NUMBER < 0) RMMAX = 0.5_EB*SM2%MEAN_DIAMETER
! ALLOW INCREASE OF THE MEAN RADIUS
RMMAX = 1.5_EB*RMMAX
ENDIF
! OTHER CONSTANTS
B_WIEN = 2.8977685E-3_EB
! CALCULATE PARAMETERS OF THE PARTICLE GROUP LOOKUP TABLE
DGROUP_A = (LOG(RMMAX)-LOG(RMMIN))/(MIE_NDG-1)
DGROUP_B = LOG(RMMAX)-DGROUP_A*MIE_NDG
! GENERATE THE PARTICLE RADII FOR MIE TABLE (MICRONS)
RDTMP = 0.5_EB*RMMIN
NX = 0
DO WHILE (RDTMP < 2._EB*RMMAX)
NX = NX + 1
RDTMP = RDTMP + MIN(1.0E6_EB,0.2_EB*RDTMP**(1._EB))
ENDDO
NRDMIE = NX
ALLOCATE(RDMIE(1:NRDMIE),STAT=IZERO)
CALL CHKMEMERR('MIEV','RDMIE',IZERO)
RDTMP = 0.5_EB*RMMIN
RDMIE(1) = RDTMP
DO NX = 2, NRDMIE
RDTMP = RDTMP + MIN(3.0E6_EB,0.2_EB*RDTMP**(1.0_EB))
RDMIE(NX) = RDTMP
ENDDO
! RADIATIVE PROPERTIES
NLMBDMIE = NWATERK
IF (CLASS_NUMBER > 0) THEN
IF (LPC%RADIATIVE_PROPERTY_INDEX>0) THEN
TA => TABLES(LPC%RADIATIVE_PROPERTY_INDEX)
NL_DATA = TA%NUMBER_ROWS
ELSEIF (LPC%Y_INDEX==H2O_INDEX) THEN
NL_DATA = NLMBDMIE
ELSEIF (LPC%FUEL) THEN
NL_DATA = NFUELK
ELSE
NL_DATA = 1
ENDIF
ELSE
IF (SM2%SINGLE_SPEC_INDEX > 0) THEN
IF (SM2%SINGLE_SPEC_INDEX==H2O_INDEX) THEN
NL_DATA = NLMBDMIE
ELSEIF (SPECIES(SM2%SINGLE_SPEC_INDEX)%ISFUEL) THEN
NL_DATA = NFUELK
ELSE
NL_DATA = 1
ENDIF
ELSE
NL_DATA = 1
ENDIF
ENDIF
! ALLOCATE ARRAYS
ALLOCATE(QSCA(1:NRDMIE,1:NLMBDMIE),STAT=IZERO)
CALL CHKMEMERR('INIT','QSCA',IZERO)
ALLOCATE(QABS(1:NRDMIE,1:NLMBDMIE),STAT=IZERO)
CALL CHKMEMERR('INIT','QABS',IZERO)
ALLOCATE(CHI_F(1:NRDMIE,1:NLMBDMIE),STAT=IZERO)
CALL CHKMEMERR('INIT','CHI_F',IZERO)
!
ALLOCATE(LMBDMIE(1:NLMBDMIE),STAT=IZERO)
CALL CHKMEMERR('INIT','LMBDMIE',IZERO)
ALLOCATE(LMBDWGHT(1:NLMBDMIE),STAT=IZERO)
CALL CHKMEMERR('INIT','LMBDWGHT',IZERO)
ALLOCATE(REAL_REF_INDX(1:NLMBDMIE),STAT=IZERO)
CALL CHKMEMERR('INIT','REAL_REF_INDX',IZERO)
ALLOCATE(CMPLX_REF_INDX(1:NLMBDMIE),STAT=IZERO)
CALL CHKMEMERR('INIT','CMPLX_REF_INDX',IZERO)
! RADIATIVE PROPERTIES
LMBDMIE(1:NLMBDMIE) = CPLXREF_WATER(1:NLMBDMIE,1)
IF (CLASS_NUMBER > 0) THEN
IF (LPC%RADIATIVE_PROPERTY_INDEX > 0) THEN
DO NX = 1,NLMBDMIE
CALL INTERPOLATE1D(TA%TABLE_DATA(:,1)*1.0E-6_EB,TA%TABLE_DATA(:,2),LMBDMIE(NX),REAL_REF_INDX(NX))
CALL INTERPOLATE1D(TA%TABLE_DATA(:,1)*1.0E-6_EB,TA%TABLE_DATA(:,3),LMBDMIE(NX),CMPLX_REF_INDX(NX))
ENDDO
ELSEIF (LPC%Y_INDEX==H2O_INDEX) THEN
REAL_REF_INDX = 1.33_EB
CMPLX_REF_INDX(1:NLMBDMIE) = CPLXREF_WATER(1:NLMBDMIE,2)
ELSEIF (LPC%FUEL) THEN
DO NX = 1,NLMBDMIE
CALL INTERPOLATE1D(CPLXREF_FUEL(:,1)*1.0E-6_EB,CPLXREF_FUEL(:,2),LMBDMIE(NX),REAL_REF_INDX(NX))
CALL INTERPOLATE1D(CPLXREF_FUEL(:,1)*1.0E-6_EB,CPLXREF_FUEL(:,3),LMBDMIE(NX),CMPLX_REF_INDX(NX))
ENDDO
ELSE
REAL_REF_INDX = LPC%REAL_REFRACTIVE_INDEX
CMPLX_REF_INDX = LPC%COMPLEX_REFRACTIVE_INDEX
ENDIF
ELSE
IF (SM2%SINGLE_SPEC_INDEX > 0) THEN
IF (SM2%SINGLE_SPEC_INDEX==H2O_INDEX) THEN
REAL_REF_INDX = 1.33_EB
CMPLX_REF_INDX(1:NLMBDMIE) = CPLXREF_WATER(1:NLMBDMIE,2)
ELSEIF (SPECIES(SM2%SINGLE_SPEC_INDEX)%ISFUEL) THEN
DO NX = 1,NLMBDMIE
CALL INTERPOLATE1D(CPLXREF_FUEL(:,1)*1.0E-6_EB,CPLXREF_FUEL(:,2),LMBDMIE(NX),REAL_REF_INDX(NX))
CALL INTERPOLATE1D(CPLXREF_FUEL(:,1)*1.0E-6_EB,CPLXREF_FUEL(:,3),LMBDMIE(NX),CMPLX_REF_INDX(NX))
ENDDO
ELSE
REAL_REF_INDX = LPC%REAL_REFRACTIVE_INDEX
CMPLX_REF_INDX = LPC%COMPLEX_REFRACTIVE_INDEX
ENDIF
ELSE
REAL_REF_INDX = LPC%REAL_REFRACTIVE_INDEX
CMPLX_REF_INDX = LPC%COMPLEX_REFRACTIVE_INDEX
ENDIF
ENDIF
CALL MIE_SCATTERING
! GENERATE INTEGRATION WEIGHTS FOR LAMBDA
IF (NLMBDMIE == 1) THEN
LMBDWGHT(1) = 1._EB
ELSE
LMBDWGHT(1) = 0.5_EB*(LMBDMIE(2)-LMBDMIE(1))
DO I = 2, NLMBDMIE-1
LMBDWGHT(I) = 0.5_EB*(LMBDMIE(I+1) - LMBDMIE(I-1))
ENDDO
LMBDWGHT(NLMBDMIE) = 0.5_EB*(LMBDMIE(NLMBDMIE)-LMBDMIE(NLMBDMIE-1))
ENDIF
ALLOCATE(WQABS(0:MIE_NDG,1:NUMBER_SPECTRAL_BANDS))
CALL ChkMemErr('IRAD','WQABS',IZERO)
WQABS = 0._EB
ALLOCATE(WQSCA(0:MIE_NDG,1:NUMBER_SPECTRAL_BANDS))
CALL ChkMemErr('IRAD','WQSCA',IZERO)
WQSCA = 0._EB
ALLOCATE(R50(0:MIE_NDG))
CALL ChkMemErr('IRAD','R50',IZERO)
R50 = 0._EB
! LOOP OVER ALL RADIATION BANDS
BANDLOOP: DO IBND = 1,NSB
!
IF (NSB == 1) THEN
NLAMBDALOW = 1
NLAMBDAHIGH = NLMBDMIE
ELSE
NLAMBDALOW = MINLOC(LMBDMIE,MASK=LMBDMIE>=WL_LOW(IBND)*1.E-6)
NLAMBDAHIGH = MAXLOC(LMBDMIE,MASK=LMBDMIE<=WL_HIGH(IBND)*1.0E-6)
ENDIF
! LOOP OVER ALL PARTICLE SIZE GROUPS
DRGROUPLOOP: DO ND = 1, MIE_NDG
R50(ND) = EXP(DGROUP_A*REAL(ND,EB) + DGROUP_B)
! LOOP OVER WAVELENGTHS
IBSUM = 0._EB
DO J = NLAMBDALOW(1),NLAMBDAHIGH(1)
IB = PLANCK(RADTMP, LMBDMIE(J)*1.0E6_EB)
IBSUM = IBSUM + LMBDWGHT(J)*IB
ASUM = 0._EB
BSUM = 0._EB
! PROPERTIES AT D32
CALL INTERPOLATE1D(RDMIE,QSCA(:,J),R50(ND),AVAL)
CALL INTERPOLATE1D(RDMIE,CHI_F(:,J),R50(ND),BVAL)
BVAL = (1._EB-BVAL)
ASUM = AVAL*BVAL
CALL INTERPOLATE1D(RDMIE,QABS(:,J),R50(ND),BVAL)
BSUM = BVAL
WQSCA(ND,IBND) = WQSCA(ND,IBND) + ASUM*LMBDWGHT(J)*IB
WQABS(ND,IBND) = WQABS(ND,IBND) + BSUM*LMBDWGHT(J)*IB
ENDDO
! NORMALIZE WITH BLACKBODY RADIATION
WQSCA(ND,IBND) = WQSCA(ND,IBND)/IBSUM
WQABS(ND,IBND) = WQABS(ND,IBND)/IBSUM
ENDDO DRGROUPLOOP
ENDDO BANDLOOP
IF (CLASS_NUMBER > 0) THEN
LPC%R50 = R50
LPC%WQSCA = WQSCA
LPC%WQABS = WQABS
ELSE
SM2%R50 = R50
SM2%WQSCA = WQSCA
SM2%WQABS = WQABS
ENDIF
DEALLOCATE(WQABS)
DEALLOCATE(WQSCA)
DEALLOCATE(R50)
DEALLOCATE(RDMIE)
DEALLOCATE(QSCA)
DEALLOCATE(QABS)
DEALLOCATE(CHI_F)
DEALLOCATE(LMBDMIE)
DEALLOCATE(LMBDWGHT)
DEALLOCATE(REAL_REF_INDX)
DEALLOCATE(CMPLX_REF_INDX)
END SUBROUTINE MEAN_CROSS_SECTIONS
SUBROUTINE MIE_SCATTERING
!
! CALCULATES THE SCATTERING AND ABSORPTION CROSS SECTIONS
! AND CALCULATES FORWARD SCATTERING FRACTION BY INTEGRATING
! THE SCATTERING PHASE FUNCTION.
!
! ----------------------------------------------------------------------
! ----------- SPECIFICATIONS FOR SUBROUTINE MIEV0 ---------------
! ----------------------------------------------------------------------
INTEGER MOMDIM
PARAMETER ( MOMDIM = 200)
LOGICAL ANYANG, PERFCT, PRNT( 2 )
INTEGER IPOLZN, NMOM
REAL(EB) GQSC, MIMCUT, PMOM( 0:MOMDIM, 4 ), SPIKE, QE, QS
REAL(EB), ALLOCATABLE :: XXX(:)
COMPLEX(EB) SFORW, SBACK, TFORW( 2 ), TBACK( 2 ), CREFIN
COMPLEX(EB), ALLOCATABLE :: S1(:), S2(:)
! --------------- LOCAL VARIABLES --------------------------------------
INTEGER I, J, K, NX, IZERO
INTEGER NLAMBDA, NRA, NMIEANG2
REAL(EB) STMP, AIJ, FTMP, XX_MAX
REAL(EB) MUMIN1, MUMIN2,THETALIM1, THETALIM2, MUDLOC
REAL(EB) MU1, MU2, NU1, NU2, MUD0LOC, MUDPILOC, MUD1, MUD2, DMUD
REAL(EB),ALLOCATABLE :: XMU1(:), XNU1(:),XMU2(:), &
ANGLE1(:), ANGLE2(:), MUD(:), MUDX(:),PWGHT(:), &
PHSFUN(:), PFOR(:,:), MUD0(:,:), MUDPI(:,:)
! ----------------------------------------------------------------------
NRA = NUMBER_RADIATION_ANGLES
! MIEV-CODE VARIABLES
MIMCUT = 1.E-6_EB
PERFCT = .FALSE.
ANYANG = .TRUE.
! IPOLZN = +1234
IPOLZN = 0
NMOM = 0
PRNT = .FALSE.
! LIMIT FOR XX
XX_MAX = 15000.0_EB
! INTEGRATION LIMITS
THETALIM1 = ACOS(1._EB - 2._EB/REAL(NRA))
MUMIN1 = COS(THETALIM1)
MUMIN2 = MUMIN1**2-(1._EB-MUMIN1**2)
THETALIM2 = ACOS(MUMIN2)
NMIEANG2 = NMIEANG*2
! ALLOCATE LOCAL ARRAYS
ALLOCATE(XXX(1:NRDMIE),STAT=IZERO)
CALL CHKMEMERR('INIT','XXX',IZERO)
ALLOCATE(S1(1:NMIEANG2),STAT=IZERO)
CALL CHKMEMERR('INIT','S1',IZERO)
ALLOCATE(S2(1:NMIEANG2),STAT=IZERO)
CALL CHKMEMERR('INIT','S2',IZERO)
ALLOCATE(XMU1(1:NMIEANG),STAT=IZERO)
CALL CHKMEMERR('INIT','XMU1',IZERO)
ALLOCATE(XNU1(1:NMIEANG),STAT=IZERO)
CALL CHKMEMERR('INIT','XNU1',IZERO)
ALLOCATE(XMU2(1:NMIEANG2),STAT=IZERO)
CALL CHKMEMERR('INIT','XMU2',IZERO)
ALLOCATE(MUD(1:NMIEANG2),STAT=IZERO)
CALL CHKMEMERR('INIT','MUD',IZERO)
ALLOCATE(MUDX(1:NMIEANG2),STAT=IZERO)
CALL CHKMEMERR('INIT','MUDX',IZERO)
ALLOCATE(PWGHT(1:NMIEANG2),STAT=IZERO)
CALL CHKMEMERR('INIT','PWGHT',IZERO)
ALLOCATE(ANGLE1(1:NMIEANG),STAT=IZERO)
CALL CHKMEMERR('INIT','ANGLE1',IZERO)
ALLOCATE(ANGLE2(1:NMIEANG2),STAT=IZERO)
CALL CHKMEMERR('INIT','ANGLE2',IZERO)
ALLOCATE(PHSFUN(1:NMIEANG2),STAT=IZERO)
CALL CHKMEMERR('INIT','PHSFUN',IZERO)
ALLOCATE(PFOR(1:NMIEANG,1:NMIEANG),STAT=IZERO)
CALL CHKMEMERR('INIT','PFOR',IZERO)
ALLOCATE(MUD0(1:NMIEANG,1:NMIEANG),STAT=IZERO)
CALL CHKMEMERR('INIT','MUD0',IZERO)
ALLOCATE(MUDPI(1:NMIEANG,1:NMIEANG),STAT=IZERO)
CALL CHKMEMERR('INIT','MUDPI',IZERO)
! CREATE SOLID ANGLE INTEGRATION ARRAYS
DO I = 1,NMIEANG
ANGLE1(I) = THETALIM1*REAL(I-1)/REAL(NMIEANG-1)
ANGLE1(I) = ANGLE1(I)*(1._EB-0.99*REAL(NMIEANG-I)/REAL(NMIEANG))
ENDDO
DO I = 1,NMIEANG
XMU1(NMIEANG-I+1) = COS(ANGLE1(I))
XNU1(NMIEANG-I+1) = SIN(ANGLE1(I))
ENDDO
! CREATE PHASE FUNCTION INGTEGRATION ARRAYS
DO I = 1,NMIEANG2
ANGLE2(I) = THETALIM2*REAL(I-1)/REAL(NMIEANG2-1)
ANGLE2(I) = ANGLE2(I)*(1._EB-0.99*REAL(NMIEANG2-I)/REAL(NMIEANG2))
ENDDO
DO I = 1,NMIEANG2
XMU2(NMIEANG2-I+1) = COS(ANGLE2(I))
ENDDO
! CALCULATE PHASE FUNCTION INGETRATION LIMITS
DO J = 1,NMIEANG
DO I = 1,NMIEANG
MUD0(I,J) = XMU1(I)*XMU1(J) + XNU1(I)*XNU1(J)
MUDPI(I,J) = XMU1(I)*XMU1(J) - XNU1(I)*XNU1(J)
ENDDO
ENDDO
! CALCULATE PHASE FUNCTION INTEGRATION WEIGHTS
MU1 = 0.7_EB
MU2 = 0.9_EB
NU1 = SQRT(1-MU1**2)
NU2 = SQRT(1-MU2**2)
MUD0LOC = MU1*MU2 + NU1*NU2
MUDPILOC = MU1*MU2 - NU1*NU2
MUD1 = MUDPILOC
MUD2 = MUD0LOC
DMUD = (MUD2-MUD1)/(NMIEANG2-1)
DO I = 1,NMIEANG2
MUD(I) = MUD1+REAL(I-1)*DMUD
ENDDO
MUD(1) = MUD(1) + 0.25_EB*DMUD !EMPIRICAL
MUD(NMIEANG2) = MUD(NMIEANG2) - 0.25_EB*DMUD !EMPIRICAL
DO I = 1,NMIEANG2
MUDX(I) = (MUD(I)-MUD1)/(MUD2-MUD1)
ENDDO
DO I = 1,NMIEANG2
PWGHT(I) = DMUD/SQRT((NU1*NU2)**2-(MUD(I)-MU1*MU2)**2)
ENDDO
PWGHT(2) = 0.5_EB*PWGHT(2)
PWGHT(NMIEANG2-1) = 0.5_EB*PWGHT(NMIEANG2-1)
! LOOP OVER WAVELENGTH
LAMBDALOOP: DO NLAMBDA = 1, NLMBDMIE
!
CREFIN = CMPLX( REAL_REF_INDX(NLAMBDA), CMPLX_REF_INDX(NLAMBDA), EB )
! CHOOSE PERFECTLY REFLECTING SPHERE, IF LARGE REAL INDEX IS GIVEN.
IF (REAL_REF_INDX(NLAMBDA) > 10._EB) PERFCT = .TRUE.
! LOOP OVER PARTICLE RADIUS
RADIUSLOOP: DO NX = 1, NRDMIE
XXX(NX) = MIN(XX_MAX,2._EB*PI*RDMIE(NX)/LMBDMIE(NLAMBDA))
CALL MIEV0( XXX(NX), CREFIN, PERFCT, MIMCUT, ANYANG, &
NMIEANG2, XMU2, NMOM, IPOLZN, MOMDIM, PRNT, &
QE, QS, GQSC, &
PMOM, SFORW, SBACK, S1, S2, TFORW, TBACK, SPIKE )
QSCA(NX,NLAMBDA) = QS
QABS(NX,NLAMBDA) = QE-QS
! CALCULATE SINGLE DROP PHASE FUNCTION
IF (ABS(QS)>TWO_EPSILON_EB) THEN
DO I = 1,NMIEANG2
PHSFUN(I) = 2._EB*(ABS(S1(I))**2 + ABS(S2(I))**2 )
ENDDO
PHSFUN = PHSFUN/(QS*XXX(NX)**2)
ELSE
PHSFUN = 1.0_EB
ENDIF
! CALCULATE THE INNERMOST INTEGRAL OF THE FORWARD SCATTERING FRACTION
PFOR = 0._EB
DO J = 1,NMIEANG
DO I = J,NMIEANG
IF (ABS(MUD0(I,J)-MUDPI(I,J))<=TWO_EPSILON_EB) THEN
CALL INTERPOLATE1D(XMU2,PHSFUN,MUD0(I,J),FTMP)
PFOR(I,J) = PI*FTMP
ELSE
MUD1 = MUDPI(I,J)
MUD2 = MUD0(I,J)
STMP = 0.0_EB
DO K = 1,NMIEANG2
MUDLOC = MUD1+MUDX(K)*(MUD2-MUD1)
CALL INTERPOLATE1D(XMU2,PHSFUN,MUDLOC,FTMP)
STMP = STMP + PWGHT(K)*FTMP
ENDDO
PFOR(I,J) = STMP
ENDIF
ENDDO
ENDDO
! CALCULATE THE TWO OUTER INTEGRALS OF THE FORWARD FRACTION
STMP = 0._EB
DO J = 1,NMIEANG-1
DO I = J+1,NMIEANG-1
AIJ = (XMU1(I+1)-XMU1(I))*(XMU1(J+1)-XMU1(J))/2._EB
STMP = STMP + (2._EB*PFOR(I,J)+PFOR(I+1,J)+PFOR(I,J+1)+2._EB*PFOR(I+1,J+1))*AIJ/3._EB
ENDDO
ENDDO
DO I = 1,NMIEANG-1
AIJ = ((XMU1(I+1)-XMU1(I))**2)/2._EB
STMP = STMP + (PFOR(I,I)+PFOR(I+1,I)+PFOR(I+1,I+1))*AIJ/3._EB
ENDDO
CHI_F(NX,NLAMBDA) = 2._EB*STMP/(4._EB*PI/NRA)
ENDDO RADIUSLOOP
ENDDO LAMBDALOOP
DEALLOCATE(XXX)
DEALLOCATE(S1)
DEALLOCATE(S2)
DEALLOCATE(XMU1)
DEALLOCATE(XNU1)
DEALLOCATE(XMU2)
DEALLOCATE(MUD)
DEALLOCATE(MUDX)
DEALLOCATE(PWGHT)
DEALLOCATE(ANGLE1)
DEALLOCATE(ANGLE2)
DEALLOCATE(PHSFUN)
DEALLOCATE(PFOR)
DEALLOCATE(MUD0)
DEALLOCATE(MUDPI)
END SUBROUTINE MIE_SCATTERING
SUBROUTINE MIEV0( XX, CREFIN, PERFCT, MIMCUT, ANYANG, NUMANG, XMU, &
NMOM, IPOLZN, MOMDIM, PRNT, QEXT, QSCA, GQSC, &
PMOM, SFORW, SBACK, S1, S2, TFORW, TBACK, &
SPIKE )
! MIE SCATTERING FOR A SINGLE PARTICLE AND WAVELENGTH.
! AUTHOR: DR. WARREN J. WISCOMBE (WISCOMBE@CLIMATE.GSFC.NASA.GOV)
! NASA GODDARD SPACE FLIGHT CENTER
! CODE 913
! GREENBELT, MD 20771
! REFERENCES
! ----------
!
! (1) WISCOMBE, W., 1979: MIE SCATTERING CALCULATIONS--ADVANCES
! IN TECHNIQUE AND FAST, VECTOR-SPEED COMPUTER CODES,
! NCAR TECH NOTE TN-140+STR, NATIONAL CENTER FOR
! ATMOSPHERIC RESEARCH, BOULDER, COLORADO (OUT OF PRINT
! BUT AN UPDATED ELECTRONIC VERSION AVAILABLE)
!
! (2) WISCOMBE, W., 1980: IMPROVED MIE SCATTERING ALGORITHMS,
! APPL. OPT. 19, 1505-1509
!
! COMPUTES MIE SCATTERING AND EXTINCTION EFFICIENCIES; ASYMMETRY
! FACTOR; FORWARD- AND BACKSCATTER AMPLITUDE; SCATTERING
! AMPLITUDES VS. SCATTERING ANGLE FOR INCIDENT POLARIZATION PARALLEL
! AND PERPENDICULAR TO THE PLANE OF SCATTERING;
! COEFFICIENTS IN THE LEGENDRE POLYNOMIAL EXPANSIONS OF EITHER THE
! UNPOLARIZED PHASE FUNCTION OR THE POLARIZED PHASE MATRIX;
! SOME QUANTITIES NEEDED IN POLARIZED RADIATIVE TRANSFER; AND
! INFORMATION ABOUT WHETHER OR NOT A RESONANCE HAS BEEN HIT.
!
! INPUT AND OUTPUT VARIABLES ARE DESCRIBED IN FILE MIEV.DOC.
! MANY STATEMENTS ARE ACCOMPANIED BY COMMENTS REFERRING TO
! REFERENCES IN MIEV.DOC, NOTABLY THE NCAR MIE REPORT WHICH IS NOW
! AVAILABLE ELECTRONICALLY AND WHICH IS REFERRED TO USING THE
! SHORTHAND (RN), MEANING EQ. (N) OF THE REPORT.
! CALLING TREE:
!
! MIEV0
! TESTMI
! TSTBAD
! MIPRNT
! ERRMSG
! CKINMI
! WRTBAD
! WRTDIM
! ERRMSG
! SMALL1
! SMALL2
! ERRMSG
! BIGA
! CONFRA
! ERRMSG
! LPCOEF
! LPCO1T
! LPCO2T
! ERRMSG
! MIPRNT
!
! I N P U T V A R I A B L E S
! -----------------------------
!
! ( EVEN IF AN INPUT VARIABLE IS NOT NEEDED FOR A PARTICULAR
! APPLICATION, MAKE SURE IT HAS A LEGITIMATE VALUE THAT CAN
! BE WRITTEN OUT AND READ IN -- NO INDEFINITES, ETC. )
!
! XX MIE SIZE PARAMETER ( 2 * PI * RADIUS / WAVELENGTH )
!
! CREFIN COMPLEX REFRACTIVE INDEX ( IMAG PART CAN BE + OR -,
! BUT INTERNALLY A NEGATIVE IMAGINARY INDEX IS ASSUMED ).
! IF IMAG PART IS - , SCATTERING AMPLITUDES AS IN VAN
! DE HULST ARE RETURNED; IF IMAG PART IS + , COMPLEX
! CONJUGATES OF THOSE SCATTERING AMPLITUDES ARE RETURNED
! (THE LATTER IS THE CONVENTION IN PHYSICS).
! ** NOTE ** IN THE 'PERFECT' CASE, SCATTERING AMPLITUDES
! IN THE VAN DE HULST (REF. 6 ABOVE) CONVENTION WILL
! AUTOMATICALLY BE RETURNED UNLESS IM(CREFIN) IS
! POSITIVE; OTHERWISE, CREFIN PLAYS NO ROLE.
!
! PERFCT TRUE, ASSUME REFRACTIVE INDEX IS INFINITE AND USE
! SPECIAL CASE FORMULAS FOR MIE COEFFICIENTS 'A'
! AND 'B' ( SEE KERKER, M., THE SCATTERING OF
! LIGHT AND OTHER ELECTROMAGNETIC RADIATION, P. 90 ).
! THIS IS SOMETIMES CALLED THE 'TOTALLY REFLECTING',
! SOMETIMES THE 'PERFECTLY CONDUCTING' CASE.
! ( SEE CREFIN FOR ADDITIONAL INFORMATION )
!
! MIMCUT (POSITIVE) VALUE BELOW WHICH IMAGINARY REFRACTIVE
! INDEX IS REGARDED AS ZERO (COMPUTATION PROCEEDS
! FASTER FOR ZERO IMAGINARY INDEX)
!
! ANYANG TRUE, ANY ANGLES WHATSOEVER MAY BE INPUT THROUGH
! XMU. FALSE, THE ANGLES ARE MONOTONE INCREASING
! AND MIRROR SYMMETRIC ABOUT 90 DEGREES (THIS OPTION
! IS ADVANTAGEOUS BECAUSE THE SCATTERING AMPLITUDES
! S1,S2 FOR THE ANGLES BETWEEN 90 AND 180 DEGREES
! ARE EVALUABLE FROM SYMMETRY RELATIONS, AND HENCE
! ARE OBTAINED WITH LITTLE ADDED COMPUTATIONAL COST.)
!
! NUMANG NO. OF ANGLES AT WHICH SCATTERING AMPLITUDES
! S1,S2 ARE TO BE EVALUATED ( SET = 0 TO SKIP
! CALCULATION OF S1,S2 ). MAKE SURE NUMANG DOES
! NOT EXCEED THE PARAMETER MAXANG IN THE PROGRAM.
!
! XMU(N) COSINES OF ANGLES ( N = 1 TO NUMANG ) AT WHICH S1,S2
! ARE TO BE EVALUATED. IF ANYANG = FALSE, THEN
!
! (A) THE ANGLES MUST BE MONOTONE INCREASING AND
! MIRROR SYMMETRIC ABOUT 90 DEGREES (IF 90-A IS
! AN ANGLE, THEN 90+A MUST BE ALSO)
!
! (B) IF NUMANG IS ODD, 90 DEGREES MUST BE AMONG
! THE ANGLES
!
! NMOM HIGHEST LEGENDRE MOMENT PMOM TO CALCULATE,
! NUMBERING FROM ZERO ( NMOM = 0 PREVENTS
! CALCULATION OF PMOM )
!
! IPOLZN POSITIVE, COMPUTE LEGENDRE MOMENTS PMOM FOR THE
! MUELLER MATRIX ELEMENTS DETERMINED BY THE
! DIGITS OF IPOLZN, WITH 1 REFERRING TO M1,
! 2 TO M2, 3 TO S21, AND 4 TO D21 (REF. 3).
! E.G., IF IPOLZN = 14 THEN ONLY MOMENTS FOR
! M1 AND D21 WILL BE RETURNED.
!
! 0, COMPUTE LEGENDRE MOMENTS PMOM FOR THE
! UNPOLARIZED UNNORMALIZED PHASE FUNCTION.
!
! NEGATIVE, COMPUTE LEGENDRE MOMENTS PMOM FOR THE
! SEKERA PHASE QUANTITIES DETERMINED BY THE
! DIGITS OF ABS(IPOLZN), WITH 1 REFERRING TO
! R1, 2 TO R2, 3 TO R3, AND 4 TO R4 (REF. 4).
! E.G., IF IPOLZN = -14 THEN ONLY MOMENTS FOR
! R1 AND R4 WILL BE RETURNED.
!
! ( NOT USED IF NMOM = 0 )
!
! MOMDIM DETERMINES FIRST DIMENSION OF PMOM, WHICH IS DIMENSIONED
! INTERNALLY AS PMOM( 0:MOMDIM, * ) (SECOND DIMENSION MUST
! BE THE LARGER OF UNITY AND THE HIGHEST DIGIT IN
! IPOLZN; IF NOT, SERIOUS ERRORS WILL OCCUR).
! MUST BE GIVEN A VALUE, EVEN IF NMOM = 0. MINIMUM: 1.
!
! PRT(L) PRINT FLAGS (LOGICAL). L = 1 PRINTS S1,S2, THEIR
! SQUARED ABSOLUTE VALUES, AND DEGREE OF POLARIZATION,
! PROVIDED NUMANG IS NON-ZERO. L = 2 PRINTS ALL
! OUTPUT VARIABLES OTHER THAN S1,S2.
!
!
! O U T P U T V A R I A B L E S
! -------------------------------
!
! QEXT (REAL) EXTINCTION EFFICIENCY FACTOR ( REF. 2, EQ. 1A )
!
! QSCA (REAL) SCATTERING EFFICIENCY FACTOR ( REF. 2, EQ. 1B )
!
! GQSC (REAL) ASYMMETRY FACTOR TIMES SCATTERING EFFICIENCY
! ( REF. 2, EQ. 1C ) ( ALLOWS CALCULATION OF RADIATION
! PRESSURE EFFICIENCY FACTOR QPR = QEXT - GQSC )