-
Notifications
You must be signed in to change notification settings - Fork 23
/
disort.f
7432 lines (5589 loc) · 216 KB
/
disort.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 DISORT( NLYR, DTAUC, SSALB, corint, NMOM, PMOM, TEMPER,
& WVNMLO, WVNMHI, USRTAU, NTAU, UTAU, NSTR, USRANG, NUMU, UMU,
& NPHI, PHI, IBCND, FBEAM, UMU0, PHI0, FISOT, LAMBER, ALBEDO,
& BTEMP, TTEMP, TEMIS, PLANK, ONLYFL, ACCUR, PRNT, HEADER,
& MAXCLY, MAXULV, MAXUMU, MAXPHI, MAXMOM, RFLDIR, RFLDN, FLUP,
& DFDT, UAVG, UU, ALBMED, TRNMED )
c *******************************************************************
c Plane-parallel discrete ordinates radiative transfer program
c (all SINGLE precision)
c ( see DISORT.doc for complete documentation )
c *******************************************************************
c
c +------------------------------------------------------------------+
c Calling Tree (omitting calls to ERRMSG):
c (routines in parentheses are not in this file)
c
c DISORT-+-(R1MACH)
c +-SLFTST-+-(TSTBAD)
c +-ZEROIT
c +-CHEKIN-+-(WRTBAD)
c | +-(WRTDIM)
c | +-DREF
c +-ZEROAL
c +-SETDIS-+-QGAUSN-+-(R1MACH)
c +-PRTINP
c +-ALBTRN-+-LEPOLY
c | +-ZEROIT
c | +-SOLEIG-+-ASYMTX-+-(R1MACH)
c | +-TERPEV
c | +-SETMTX-+-ZEROIT
c | +-(SGBCO)
c | +-SOLVE1-+-ZEROIT
c | | +-(SGBSL)
c | +-ALTRIN
c | +-SPALTR
c | +-PRALTR
c +-PLKAVG-+-(R1MACH)
c +-LEPOLY
c +-SURFAC-+-QGAUSN-+-(R1MACH)
c | +-BDREF
c | +-ZEROIT
c +-SOLEIG-+-ASYMTX-+-(R1MACH)
c +-UPBEAM-+-(SGECO)
c | +-(SGESL)
c +-UPISOT-+-(SGECO)
c | +-(SGESL)
c +-TERPEV
c +-TERPSO
c +-SETMTX-+-ZEROIT
c +-SOLVE0-+-ZEROIT
c | +-(SGBCO)
c | +-(SGBSL)
c +-FLUXES--ZEROIT
c +-ZEROIT
c +-USRINT
c +-CMPINT
c +-PRAVIN
c +-ZEROIT
c +-RATIO--(R1MACH)
c +-INTCOR-+-SINSCA
c | +-SECSCA-+-XIFUNC
c +-PRTINT
c
c *** Intrinsic Functions used in DISORT package which take
c non-negligible amount of time:
c
c EXP : Called by- ALBTRN, ALTRIN, CMPINT, FLUXES, SETDIS,
c SETMTX, SPALTR, USRINT, PLKAVG
c
c SQRT : Called by- ASYMTX, SOLEIG
c
c +-------------------------------------------------------------------+
c
c Index conventions (for all DO-loops and all variable descriptions):
c
c IU : for user polar angles
c
c IQ,JQ,KQ : for computational polar angles ('quadrature angles')
c
c IQ/2 : for half the computational polar angles (just the ones
c in either 0-90 degrees, or 90-180 degrees)
c
c J : for user azimuthal angles
c
c K,L : for Legendre expansion coefficients or, alternatively,
c subscripts of associated Legendre polynomials
c
c LU : for user levels
c
c LC : for computational layers (each having a different
c single-scatter albedo and/or phase function)
c
c LEV : for computational levels
c
c MAZIM : for azimuthal components in Fourier cosine expansion
c of intensity and phase function
c
c +------------------------------------------------------------------+
c
c I N T E R N A L V A R I A B L E S
c
c AMB(IQ/2,IQ/2) First matrix factor in reduced eigenvalue problem
c of Eqs. SS(12), STWJ(8E), STWL(23f)
c (used only in SOLEIG)
c
c APB(IQ/2,IQ/2) Second matrix factor in reduced eigenvalue problem
c of Eqs. SS(12), STWJ(8E), STWL(23f)
c (used only in SOLEIG)
c
c ARRAY(IQ,IQ) Scratch matrix for SOLEIG, UPBEAM and UPISOT
c (see each subroutine for definition)
c
c B() Right-hand side vector of Eq. SC(5) going into
c SOLVE0,1; returns as solution vector
c vector L, the constants of integration
c
c BDR(IQ/2,0:IQ/2) Bottom-boundary bidirectional reflectivity for a
c given azimuthal component. First index always
c refers to a computational angle. Second index:
c if zero, refers to incident beam angle UMU0;
c if non-zero, refers to a computational angle.
c
c BEM(IQ/2) Bottom-boundary directional emissivity at compu-
c tational angles.
c
c BPLANK Intensity emitted from bottom boundary
c
c CBAND() Matrix of left-hand side of the linear system
c Eq. SC(5), scaled by Eq. SC(12); in banded
c form required by LINPACK solution routines
c
c CC(IQ,IQ) C-sub-IJ in Eq. SS(5)
c
c CMU(IQ) Computational polar angles (Gaussian)
c
c CWT(IQ) Quadrature weights corresponding to CMU
c
c CORINT When set TRUE, correct intensities for
c delta-scaling effects (see Nakajima and Tanaka,
c 1988). When FALSE, intensities are not corrected.
c In general, CORINT should be set true when beam
c source is present (FBEAM is not zero) and DELTAM
c is TRUE in a problem including scattering.
c However, execution is faster when CORINT is FALSE,
c and intensities outside the aureole may still be
c accurate enough. When CORINT is TRUE, it is
c important to have a sufficiently high order of
c Legendre approximation of the phase function. This
c is because the intensities are corrected by
c calculating the single-scattered radiation, for
c which an adequate representation of the phase
c function is crucial. In case of a low order
c Legendre approximation of an otherwise highly
c anisotropic phase function, the intensities might
c actually be more accurate when CORINT is FALSE.
c When only fluxes are calculated (ONLYFL is TRUE),
c or there is no beam source (FBEAM=0.0), or there
c is no scattering (SSALB=0.0 for all layers) CORINT
c is set FALSE by the code.
c
c DELM0 Kronecker delta, delta-sub-M0, where M = MAZIM
c is the number of the Fourier component in the
c azimuth cosine expansion
c
c DELTAM TRUE, use delta-M method ( see Wiscombe, 1977 );
c FALSE, do not use delta-M method. In general, for
c a given number of streams, intensities and
c fluxes will be more accurate for phase functions
c with a large forward peak if DELTAM is set true.
c Intensities close to the forward scattering
c direction are often less accurate, however, when
c the delta-M method is applied. The intensity
c correction of Nakajima and Tanaka is used to
c improve the accuracy of the intensities.
c
c DITHER Small quantity subtracted from single-scattering
c albedos of unity, in order to avoid using special
c case formulas; prevents an eigenvalue of exactly
c zero from occurring, which would cause an
c immediate overflow
c
c DTAUCP(LC) Computational-layer optical depths (delta-M-scaled
c if DELTAM = TRUE, otherwise equal to DTAUC)
c
c EMU(IU) Bottom-boundary directional emissivity at user
c angles.
c
c EVAL(IQ) Temporary storage for eigenvalues of Eq. SS(12)
c
c EVECC(IQ,IQ) Complete eigenvectors of SS(7) on return from
c SOLEIG; stored permanently in GC
c
c EXPBEA(LC) Transmission of direct beam in delta-M optical
c depth coordinates
c
c FLYR(LC) Separated fraction in delta-M method
c
c GL(K,LC) Phase function Legendre polynomial expansion
c coefficients, calculated from PMOM by
c including single-scattering albedo, factor
c 2K+1, and (if DELTAM=TRUE) the delta-M
c scaling
c
c GC(IQ,IQ,LC) Eigenvectors at polar quadrature angles,
c g in Eq. SC(1)
c
c GU(IU,IQ,LC) Eigenvectors interpolated to user polar angles
c ( g in Eqs. SC(3) and S1(8-9), i.e.
c G without the L factor )
c
c IPVT(LC*IQ) Integer vector of pivot indices for LINPACK
c routines
c
c KK(IQ,LC) Eigenvalues of coeff. matrix in Eq. SS(7)
c
c KCONV Counter in azimuth convergence test
c
c LAYRU(LU) Computational layer in which user output level
c UTAU(LU) is located
c
c LL(IQ,LC) Constants of integration L in Eq. SC(1),
c obtained by solving scaled version of Eq. SC(5)
c
c LYRCUT TRUE, radiation is assumed zero below layer
c NCUT because of almost complete absorption
c
c NAZ Number of azimuthal components considered
c
c NCUT Computational layer number in which absorption
c optical depth first exceeds ABSCUT
c
c OPRIM(LC) Single scattering albedo after delta-M scaling
c
c PASS1 TRUE on first entry, FALSE thereafter
c
c PKAG(0:LC) Integrated Planck function for internal emission
c
c PRNTU0(L) logical flag to trigger printing of azimuthally-
c averaged intensities:
c L quantities printed
c -- ------------------
c 1 azimuthally-averaged intensities at user
c levels and computational polar angles
c 2 azimuthally-averaged intensities at user
c levels and user polar angles
c
c PSI0(IQ) Sum just after square bracket in Eq. SD(9)
c
c PSI1(IQ) Sum in Eq. STWL(31d)
c
c RMU(IU,0:IQ) Bottom-boundary bidirectional reflectivity for a
c given azimuthal component. First index always
c refers to a user angle. Second index:
c if zero, refers to incident beam angle UMU0;
c if non-zero, refers to a computational angle.
c
c SQT(k) Square root of k (used only in LEPOLY for
c computing associated Legendre polynomials)
c
c TAUC(0:LC) Cumulative optical depth (un-delta-M-scaled)
c
c TAUCPR(0:LC) Cumulative optical depth (delta-M-scaled if
c DELTAM = TRUE, otherwise equal to TAUC)
c
c TPLANK Intensity emitted from top boundary
c
c UUM(IU,LU) Expansion coefficients when the intensity
c (u-super-M) is expanded in Fourier cosine series
c in azimuth angle
c
c U0C(IQ,LU) Azimuthally-averaged intensity at quadrature
c angle
c
c U0U(IU,LU) If ONLYFL = FALSE, azimuthally-averaged intensity
c at user angles and user levels
c
c If ONLYFL = TRUE and MAXUMU.GE.NSTR,
c azimuthally-averaged intensity at computational
c (Gaussian quadrature) angles and user levels;
c the corresponding quadrature angle cosines are
c returned in UMU. If MAXUMU.LT.NSTR, U0U will be
c zeroed, and UMU, NUMU will not be set.
c
c UTAUPR(LU) Optical depths of user output levels in delta-M
c coordinates; equal to UTAU(LU) if no delta-M
c
c WK() scratch array
c
c XR0(LC) X-sub-zero in expansion of thermal source func-
c tion preceding Eq. SS(14)(has no mu-dependence);
c b-sub-zero in Eq. STWL(24d)
c
c XR1(LC) X-sub-one in expansion of thermal source func-
c tion; see Eqs. SS(14-16); b-sub-one in STWL(24d)
c
c YLM0(L) Normalized associated Legendre polynomial
c of subscript L at the beam angle (not saved
c as function of superscipt M)
c
c YLMC(L,IQ) Normalized associated Legendre polynomial
c of subscript L at the computational angles
c (not saved as function of superscipt M)
c
c YLMU(L,IU) Normalized associated Legendre polynomial
c of subscript L at the user angles
c (not saved as function of superscipt M)
c
c Z() scratch array used in SOLVE0, ALBTRN to solve
c a linear system for the constants of integration
c
c Z0(IQ) Solution vectors Z-sub-zero of Eq. SS(16)
c
c Z0U(IU,LC) Z-sub-zero in Eq. SS(16) interpolated to user
c angles from an equation derived from SS(16)
c
c Z1(IQ) Solution vectors Z-sub-one of Eq. SS(16)
c
c Z1U(IU,LC) Z-sub-one in Eq. SS(16) interpolated to user
c angles from an equation derived from SS(16)
c
c ZBEAM(IU,LC) Particular solution for beam source
c
c ZJ(IQ) Right-hand side vector X-sub-zero in
c Eq. SS(19), also the solution vector
c Z-sub-zero after solving that system
c
c ZZ(IQ,LC) Permanent storage for the beam source vectors ZJ
c
c ZPLK0(IQ,LC) Permanent storage for the thermal source
c vectors Z0 obtained by solving Eq. SS(16)
c
c ZPLK1(IQ,LC) Permanent storage for the thermal source
c vectors Z1 obtained by solving Eq. SS(16)
c
c +-------------------------------------------------------------------+
c
c LOCAL SYMBOLIC DIMENSIONS (have big effect on storage requirements):
c
c MXCLY = Max no. of computational layers
c MXULV = Max no. of output levels
c MXCMU = Max no. of computation polar angles
c MXUMU = Max no. of output polar angles
c MXPHI = Max no. of output azimuthal angles
c MXSQT = Max no. of square roots of integers (for LEPOLY)
c +-------------------------------------------------------------------+
use params, only: mxly, nstrms, kr
c .. Parameters ..
INTEGER MXCLY, MXULV, MXCMU, MXUMU, MXPHI, MI, MI9M2, NNLYRI,
& MXSQT
PARAMETER ( MXCLY=mxly, MXULV = mxly+1,
& MXCMU = nstrms, MXUMU = nstrms, MXPHI = nstrms,
& MI = MXCMU / 2, MI9M2 = 9*MI - 2,
& NNLYRI = MXCMU*MXCLY, MXSQT = 1000 )
c ..
c .. Scalar Arguments ..
CHARACTER HEADER*127
LOGICAL LAMBER, ONLYFL, PLANK, USRANG, USRTAU
INTEGER IBCND, MAXCLY, MAXMOM, MAXPHI, MAXULV, MAXUMU, NLYR,
& NMOM, NPHI, NSTR, NTAU, NUMU
REAL(KR) ACCUR, ALBEDO, BTEMP, FBEAM, FISOT, PHI0, TEMIS, TTEMP,
& UMU0, WVNMHI, WVNMLO
c ..
c .. Array Arguments ..
LOGICAL PRNT( 5 )
REAL(KR) ALBMED( MAXUMU ), DFDT( MAXULV ), DTAUC( MAXCLY ),
& FLUP( MAXULV ), PHI( MAXPHI ), PMOM( 0:MAXMOM, MAXCLY ),
& RFLDIR( MAXULV ), RFLDN( MAXULV ), SSALB( MAXCLY ),
& TEMPER( 0:MAXCLY ), TRNMED( MAXUMU ), UAVG( MAXULV ),
& UMU( MAXUMU ), UTAU( MAXULV ),
& UU( MAXUMU, MAXULV, MAXPHI )
c ..
c .. Local Scalars ..
LOGICAL COMPAR, CORINT, DELTAM, LYRCUT, PASS1
INTEGER IQ, IU, J, KCONV, L, LC, LEV, LU, MAZIM, NAZ, NCOL,
& NCOS, NCUT, NN, NS
REAL(KR) ANGCOS(1), AZERR, AZTERM, BPLANK, COSPHI, DELM0, DITHER,
& DUM, PI, RPD, SGN, TPLANK
c ..
c .. Local Arrays ..
LOGICAL PRNTU0( 2 )
INTEGER IPVT( NNLYRI ), LAYRU( MXULV )
REAL(KR) AMB( MI, MI ), APB( MI, MI ), ARRAY( MXCMU, MXCMU ),
& B( NNLYRI ), BDR( MI, 0:MI ), BEM( MI ),
& CBAND( MI9M2, NNLYRI ), CC( MXCMU, MXCMU ),
& CMU( MXCMU ), CWT( MXCMU ), DTAUCP( MXCLY ),
& EMU( MXUMU ), EVAL( MI ), EVECC( MXCMU, MXCMU ),
& EXPBEA( 0:MXCLY ), FLDIR( MXULV ), FLDN( MXULV ),
& FLYR( MXCLY ), GC( MXCMU, MXCMU, MXCLY ),
& GL( 0:MXCMU, MXCLY ), GU( MXUMU, MXCMU, MXCLY ),
& KK( MXCMU, MXCLY ), LL( MXCMU, MXCLY ), OPRIM( MXCLY ),
& PHASA( MXCLY ), PHAST( MXCLY ), PHASM( MXCLY ),
& PHIRAD( MXPHI ), PKAG( 0:MXCLY ), PSI0( MXCMU ),
& PSI1( MXCMU ), RMU( MXUMU, 0:MI ), SQT( MXSQT ),
& TAUC( 0:MXCLY ), TAUCPR( 0:MXCLY ), U0C( MXCMU, MXULV ),
& U0U( MXUMU, MXULV ), UTAUPR( MXULV ),
& UUM( MXUMU, MXULV ), WK( MXCMU ), XR0( MXCLY ),
& XR1( MXCLY ), YLM0( 0:MXCMU ), YLMC( 0:MXCMU, MXCMU ),
& YLMU( 0:MXCMU, MXUMU ), Z( NNLYRI ), Z0( MXCMU ),
& Z0U( MXUMU, MXCLY ), Z1( MXCMU ), Z1U( MXUMU, MXCLY ),
& ZBEAM( MXUMU, MXCLY )
REAL(KR) ZJ( MXCMU ), ZPLK0( MXCMU, MXCLY ),
& ZPLK1( MXCMU, MXCLY ), ZZ( MXCMU, MXCLY )
REAL(KR) WKD( MXCMU )
c ..
c .. External Functions ..
REAL(KR) PLKAVG, R1MACH, RATIO
EXTERNAL PLKAVG, R1MACH, RATIO
c ..
c .. External Subroutines ..
EXTERNAL ALBTRN, CHEKIN, CMPINT, FLUXES, INTCOR, LEPOLY, PRAVIN,
& PRTINP, PRTINT, SETDIS, SETMTX, SLFTST, SOLEIG, SOLVE0,
& SURFAC, TERPEV, TERPSO, UPBEAM, UPISOT, USRINT, ZEROAL,
& ZEROIT
c ..
c .. Intrinsic Functions ..
INTRINSIC ABS, ASIN, COS, FLOAT, LEN, MAX, SQRT
c ..
SAVE DITHER, PASS1, PI, RPD, SQT
DATA PASS1 / .TRUE. /, PRNTU0 / 2*.FALSE. /
DELTAM = .TRUE.
c CORINT = .TRUE. ! set by calling routine (pjr)
IF( PASS1 ) THEN
PI = 2.*ASIN( 1.0 )
DITHER = 10.*R1MACH( 4 )
c ** Must dither more on high (14-digit)
c ** precision machine
c mfl
c print *, 'DITHER = ', DITHER, ', ', R1MACH( 4 )
IF( DITHER.LT.1.E-10 ) DITHER = 10.*DITHER
c print *, 'DITHER2 = ', DITHER
RPD = PI / 180.0
DO 10 NS = 1, MXSQT
SQT( NS ) = SQRT( FLOAT( NS ) )
10 CONTINUE
c ** Set input values for self-test.
c ** Be sure SLFTST sets all print flags off.
COMPAR = .FALSE.
CALL SLFTST( CORINT, ACCUR, ALBEDO, BTEMP, DELTAM, DTAUC( 1 ),
& FBEAM, FISOT, IBCND, LAMBER, NLYR, PLANK, NPHI,
& NUMU, NSTR, NTAU, ONLYFL, PHI( 1 ), PHI0, NMOM,
& PMOM( 0,1 ), PRNT, PRNTU0, SSALB( 1 ), TEMIS,
& TEMPER( 0 ), TTEMP, UMU( 1 ), USRANG, USRTAU,
& UTAU( 1 ), UMU0, WVNMHI, WVNMLO, COMPAR, DUM,
& DUM, DUM, DUM )
END IF
20 CONTINUE
c IF( .NOT.PASS1 .AND. len(HEADER).ne.0 )
c & WRITE( *,'(//,1X,100(''*''),/,A,/,1X,100(''*''))' )
c & ' DISORT: '//HEADER
c ** Calculate cumulative optical depth
c ** and dither single-scatter albedo
c ** to improve numerical behavior of
c ** eigenvalue/vector computation
CALL ZEROIT( TAUC, MXCLY + 1 )
DO 30 LC = 1, NLYR
IF( SSALB( LC ).EQ.1.0 ) SSALB( LC ) = 1.0 - DITHER
TAUC( LC ) = TAUC( LC - 1 ) + DTAUC( LC )
30 CONTINUE
c ** Check input dimensions and variables
CALL CHEKIN( NLYR, DTAUC, SSALB, NMOM, PMOM, TEMPER, WVNMLO,
& WVNMHI, USRTAU, NTAU, UTAU, NSTR, USRANG,
& NUMU, UMU, NPHI, PHI, IBCND, FBEAM, UMU0,
& PHI0, FISOT, LAMBER, ALBEDO, BTEMP, TTEMP,
& TEMIS, PLANK, ONLYFL, DELTAM, CORINT, ACCUR,
& TAUC, MAXCLY, MAXULV, MAXUMU, MAXPHI, MAXMOM,
& MXCLY, MXULV, MXUMU, MXCMU, MXPHI, MXSQT )
c ** Zero internal and output arrays
CALL ZEROAL( MXCLY, EXPBEA(1), FLYR, OPRIM, PHASA, PHAST, PHASM,
& TAUCPR(1), XR0, XR1,
& MXCMU, CMU, CWT, PSI0, PSI1, WK, Z0, Z1, ZJ,
& MXCMU+1, YLM0,
& MXCMU**2, ARRAY, CC, EVECC,
& (MXCMU+1)*MXCLY, GL,
& (MXCMU+1)*MXCMU, YLMC,
& (MXCMU+1)*MXUMU, YLMU,
& MXCMU*MXCLY, KK, LL, ZZ, ZPLK0, ZPLK1,
& MXCMU**2*MXCLY, GC,
& MXULV, LAYRU, UTAUPR,
& MXUMU*MXCMU*MXCLY, GU,
& MXUMU*MXCLY, Z0U, Z1U, ZBEAM,
& MI, EVAL,
& MI**2, AMB, APB,
& NNLYRI, IPVT, Z,
& MAXULV, RFLDIR, RFLDN, FLUP, UAVG, DFDT,
& MAXUMU, ALBMED, TRNMED,
& MXUMU*MXULV, U0U,
& MAXUMU*MAXULV*MAXPHI, UU )
c ** Perform various setup operations
CALL SETDIS( CMU, CWT, DELTAM, DTAUC, DTAUCP, EXPBEA, FBEAM, FLYR,
& GL, IBCND, LAYRU, LYRCUT, MAXMOM, MAXUMU, MXCMU,
& NCUT, NLYR, NTAU, NN, NSTR, PLANK, NUMU, ONLYFL,
& CORINT, OPRIM, PMOM, SSALB, TAUC, TAUCPR, UTAU,
& UTAUPR, UMU, UMU0, USRTAU, USRANG )
if(nstr.lt.0) return
c ** Print input information
IF( PRNT( 1 ) )
& CALL PRTINP( NLYR, DTAUC, DTAUCP, SSALB, NMOM, PMOM, TEMPER,
& WVNMLO, WVNMHI, NTAU, UTAU, NSTR, NUMU, UMU,
& NPHI, PHI, IBCND, FBEAM, UMU0, PHI0, FISOT,
& LAMBER, ALBEDO, BTEMP, TTEMP, TEMIS, DELTAM,
& PLANK, ONLYFL, CORINT, ACCUR, FLYR, LYRCUT,
& OPRIM, TAUC, TAUCPR, MAXMOM, PRNT( 5 ) )
c ** Handle special case for getting albedo
c ** and transmissivity of medium for many
c ** beam angles at once
IF( IBCND.EQ.1 ) THEN
CALL ALBTRN( ALBEDO, AMB, APB, ARRAY, B, BDR, CBAND, CC, CMU,
& CWT, DTAUCP, EVAL, EVECC, GL, GC, GU, IPVT, KK,
& LL, NLYR, NN, NSTR, NUMU, PRNT, TAUCPR, UMU, U0U,
& WK, YLMC, YLMU, Z, WKD, MI, MI9M2, MAXUMU, MXCMU,
& MXUMU, NNLYRI, SQT, ALBMED, TRNMED )
RETURN
END IF
c ** Calculate Planck functions
IF( .NOT.PLANK ) THEN
BPLANK = 0.0
TPLANK = 0.0
CALL ZEROIT( PKAG, MXCLY + 1 )
ELSE
TPLANK = TEMIS*PLKAVG( WVNMLO, WVNMHI, TTEMP )
BPLANK = PLKAVG( WVNMLO, WVNMHI, BTEMP )
DO 40 LEV = 0, NLYR
PKAG( LEV ) = PLKAVG( WVNMLO, WVNMHI, TEMPER( LEV ) )
40 CONTINUE
END IF
c ======== BEGIN LOOP TO SUM AZIMUTHAL COMPONENTS OF INTENSITY =======
c (EQ STWJ 5, STWL 6)
KCONV = 0
NAZ = NSTR - 1
c ** Azimuth-independent case
IF( FBEAM.EQ.0.0 .OR. ABS(1.-UMU0).LT.1.E-5 .OR. ONLYFL .OR.
& ( NUMU.EQ.1 .AND. ABS(1.-UMU(1)).LT.1.E-5 ) .OR.
& ( NUMU.EQ.1 .AND. ABS(1.+UMU(1)).LT.1.E-5 ) .OR.
& ( NUMU.EQ.2 .AND. ABS(1.+UMU(1)).LT.1.E-5 .AND.
& ABS(1.-UMU(2)).LT.1.E-5 ) )
& NAZ = 0
DO 180 MAZIM = 0, NAZ
IF( MAZIM.EQ.0 ) DELM0 = 1.0
IF( MAZIM.GT.0 ) DELM0 = 0.0
c ** Get normalized associated Legendre
c ** polynomials for
c ** (a) incident beam angle cosine
c ** (b) computational and user polar angle
c ** cosines
IF( FBEAM.GT.0.0 ) THEN
NCOS = 1
ANGCOS(1) = -UMU0
CALL LEPOLY( NCOS, MAZIM, MXCMU, NSTR-1, ANGCOS, SQT, YLM0 )
END IF
IF( .NOT.ONLYFL .AND. USRANG )
& CALL LEPOLY( NUMU, MAZIM, MXCMU, NSTR-1, UMU, SQT, YLMU )
CALL LEPOLY( NN, MAZIM, MXCMU, NSTR-1, CMU, SQT, YLMC )
c ** Get normalized associated Legendre polys.
c ** with negative arguments from those with
c ** positive arguments; Dave/Armstrong Eq. (15),
c ** STWL(59)
SGN = -1.0
DO 70 L = MAZIM, NSTR - 1
SGN = -SGN
DO 60 IQ = NN + 1, NSTR
YLMC( L, IQ ) = SGN*YLMC( L, IQ - NN )
60 CONTINUE
70 CONTINUE
c ** Specify users bottom reflectivity
c ** and emissivity properties
IF( .NOT.LYRCUT )
& CALL SURFAC( ALBEDO, DELM0, CMU, FBEAM, LAMBER, MI, MAZIM,
& MXUMU, NN, NUMU, ONLYFL, PI, UMU, UMU0,
& USRANG, WVNMLO, WVNMHI, BDR, EMU, BEM, RMU )
c =================== BEGIN LOOP ON COMPUTATIONAL LAYERS =============
DO 80 LC = 1, NCUT
c ** Solve eigenfunction problem in Eq. STWJ(8B),
c ** STWL(23f); return eigenvalues and eigenvectors
CALL SOLEIG( AMB, APB, ARRAY, CMU, CWT, GL( 0,LC ), MI,
& MAZIM, MXCMU, NN, NSTR, YLMC, CC, EVECC, EVAL,
& KK( 1,LC ), GC( 1,1,LC ), WKD )
c ** Calculate particular solutions of
c ** Eq. SS(18), STWL(24a) for incident
c ** beam source
IF( FBEAM.GT.0.0 )
& CALL UPBEAM( ARRAY, CC, CMU, DELM0, FBEAM, GL( 0,LC ),
& IPVT, MAZIM, MXCMU, NN, NSTR, PI, UMU0, WK,
& YLM0, YLMC, ZJ, ZZ( 1,LC ) )
c ** Calculate particular solutions of Eq.
c ** SS(15), STWL(25) for thermal emission
c ** source
c
IF( PLANK .AND. MAZIM.EQ.0 ) THEN
XR1( LC ) = 0.0
IF( DTAUCP( LC ).GT.0.0 ) XR1( LC ) =
& ( PKAG( LC ) - PKAG( LC-1 ) ) / DTAUCP( LC )
XR0( LC ) = PKAG( LC-1 ) - XR1( LC )*TAUCPR( LC-1 )
CALL UPISOT( ARRAY, CC, CMU, IPVT, MXCMU, NN, NSTR,
& OPRIM( LC ), WK, XR0( LC ), XR1( LC ),
& Z0, Z1, ZPLK0( 1,LC ), ZPLK1( 1,LC ) )
END IF
IF( .NOT.ONLYFL .AND. USRANG ) THEN
c ** Interpolate eigenvectors
c ** to user angles
CALL TERPEV( CWT, EVECC, GL( 0,LC ), GU( 1,1,LC ), MAZIM,
& MXCMU, MXUMU, NN, NSTR, NUMU, WK, YLMC,
& YLMU )
c ** Interpolate source terms
c ** to user angles
CALL TERPSO( CWT, DELM0, FBEAM, GL( 0,LC ), MAZIM, MXCMU,
& PLANK, NUMU, NSTR, OPRIM( LC ), PI, YLM0,
& YLMC, YLMU, PSI0, PSI1, XR0( LC ),
& XR1( LC ), Z0, Z1, ZJ, ZBEAM( 1,LC ),
& Z0U( 1,LC ), Z1U( 1,LC ) )
END IF
80 CONTINUE
c =================== END LOOP ON COMPUTATIONAL LAYERS ===============
c ** Set coefficient matrix of equations combining
c ** boundary and layer interface conditions
CALL SETMTX( BDR, CBAND, CMU, CWT, DELM0, DTAUCP, GC, KK,
& LAMBER, LYRCUT, MI, MI9M2, MXCMU, NCOL, NCUT,
& NNLYRI, NN, NSTR, TAUCPR, WK )
c ** Solve for constants of integration in homo-
c ** geneous solution (general boundary conditions)
CALL SOLVE0( B, BDR, BEM, BPLANK, CBAND, CMU, CWT, EXPBEA,
& FBEAM, FISOT, IPVT, LAMBER, LL, LYRCUT, MAZIM, MI,
& MI9M2, MXCMU, NCOL, NCUT, NN, NSTR, NNLYRI, PI,
& TPLANK, TAUCPR, UMU0, Z, ZZ, ZPLK0, ZPLK1 )
c ** Compute upward and downward fluxes
IF( MAZIM.EQ.0 )
& CALL FLUXES( CMU, CWT, FBEAM, GC, KK, LAYRU, LL, LYRCUT,
& MAXULV, MXCMU, MXULV, NCUT, NN, NSTR, NTAU,
& PI, PRNT, PRNTU0( 1 ), SSALB, TAUCPR, UMU0,
& UTAU, UTAUPR, XR0, XR1, ZZ, ZPLK0, ZPLK1,
& DFDT, FLUP, FLDN, FLDIR, RFLDIR, RFLDN, UAVG,
& U0C )
IF( ONLYFL ) THEN
IF( MAXUMU.GE.NSTR ) THEN
c ** Save azimuthal-avg intensities
c ** at quadrature angles
DO 100 LU = 1, NTAU
DO 90 IQ = 1, NSTR
U0U( IQ, LU ) = U0C( IQ, LU )
90 CONTINUE
100 CONTINUE
END IF
GO TO 190
END IF
CALL ZEROIT( UUM, MXUMU*MXULV )
IF( USRANG ) THEN
c ** Compute azimuthal intensity
c ** components at user angles
CALL USRINT( BPLANK, CMU, CWT, DELM0, DTAUCP, EMU, EXPBEA,
& FBEAM, FISOT, GC, GU, KK, LAMBER, LAYRU, LL,
& LYRCUT, MAZIM, MXCMU, MXULV, MXUMU, NCUT, NLYR,
& NN, NSTR, PLANK, NUMU, NTAU, PI, RMU, TAUCPR,
& TPLANK, UMU, UMU0, UTAUPR, WK, ZBEAM, Z0U, Z1U,
& ZZ, ZPLK0, ZPLK1, UUM )
ELSE
c ** Compute azimuthal intensity
c ** components at quadrature angles
CALL CMPINT( FBEAM, GC, KK, LAYRU, LL, LYRCUT, MAZIM, MXCMU,
& MXULV, MXUMU, NCUT, NN, NSTR, PLANK, NTAU,
& TAUCPR, UMU0, UTAUPR, ZZ, ZPLK0, ZPLK1, UUM )
END IF
IF( MAZIM.EQ.0 ) THEN
c ** Save azimuthally averaged intensities
DO 130 LU = 1, NTAU
DO 120 IU = 1, NUMU
U0U( IU, LU ) = UUM( IU, LU )
DO 110 J = 1, NPHI
UU( IU, LU, J ) = UUM( IU, LU )
110 CONTINUE
120 CONTINUE
130 CONTINUE
c ** Print azimuthally averaged intensities
c ** at user angles
IF( PRNTU0( 2 ) )
& CALL PRAVIN( UMU, NUMU, MXUMU, UTAU, NTAU, U0U )
IF( NAZ.GT.0 ) THEN
CALL ZEROIT( PHIRAD, MXPHI )
DO 140 J = 1, NPHI
PHIRAD( J ) = RPD*( PHI( J ) - PHI0 )
140 CONTINUE
END IF
ELSE
c ** Increment intensity by current
c ** azimuthal component (Fourier
c ** cosine series); Eq SD(2), STWL(6)
AZERR = 0.0
DO 170 J = 1, NPHI
COSPHI = COS( MAZIM*PHIRAD( J ) )
DO 160 LU = 1, NTAU
DO 150 IU = 1, NUMU
AZTERM = UUM( IU, LU )*COSPHI
UU( IU, LU, J ) = UU( IU, LU, J ) + AZTERM
AZERR = MAX( AZERR,
& RATIO( ABS(AZTERM), ABS(UU(IU,LU,J)) ) )
150 CONTINUE
160 CONTINUE
170 CONTINUE
IF( AZERR.LE.ACCUR ) KCONV = KCONV + 1
IF( KCONV.GE.2 ) GO TO 190
END IF
180 CONTINUE
c =================== END LOOP ON AZIMUTHAL COMPONENTS ===============
190 CONTINUE
c ** Apply Nakajima/Tanaka intensity
c ** corrections
IF( CORINT )
& CALL INTCOR( DITHER, FBEAM, FLYR, LAYRU, LYRCUT, MAXMOM,
& MAXULV, MAXUMU, NMOM, NCUT, NPHI, NSTR, NTAU,
& NUMU, OPRIM, PHASA, PHAST, PHASM, PHIRAD, PI,
& RPD, PMOM, SSALB, DTAUC, TAUC, TAUCPR, UMU,
& UMU0, UTAU, UTAUPR, UU )
c ** Print intensities
IF( PRNT( 3 ) .AND. .NOT.ONLYFL )
& CALL PRTINT( UU, UTAU, NTAU, UMU, NUMU, PHI, NPHI, MAXULV,
& MAXUMU )
IF( PASS1 ) THEN
c ** Compare test case results with
c ** correct answers and abort if bad
COMPAR = .TRUE.
CALL SLFTST( CORINT, ACCUR, ALBEDO, BTEMP, DELTAM, DTAUC( 1 ),
& FBEAM, FISOT, IBCND, LAMBER, NLYR, PLANK, NPHI,
& NUMU, NSTR, NTAU, ONLYFL, PHI( 1 ), PHI0, NMOM,
& PMOM( 0,1 ), PRNT, PRNTU0, SSALB( 1 ), TEMIS,
& TEMPER( 0 ), TTEMP, UMU( 1 ), USRANG, USRTAU,
& UTAU( 1 ), UMU0, WVNMHI, WVNMLO, COMPAR,
& FLUP( 1 ), RFLDIR( 1 ), RFLDN( 1 ), UU( 1,1,1 ) )
PASS1 = .FALSE.
GO TO 20
END IF
RETURN
END
SUBROUTINE ASYMTX( AA, EVEC, EVAL, M, IA, IEVEC, IER, WK )
c ====== S I N G L E P R E C I S I O N V E R S I O N ======
c Solves eigenfunction problem for real asymmetric matrix
c for which it is known a priori that the eigenvalues are real.
c This is an adaptation of a subroutine EIGRF in the IMSL
c library to use real instead of complex arithmetic, accounting
c for the known fact that the eigenvalues and eigenvectors in
c the discrete ordinate solution are real. Other changes include
c putting all the called subroutines in-line, deleting the
c performance index calculation, updating many DO-loops
c to Fortran77, and in calculating the machine precision
c TOL instead of specifying it in a data statement.
c EIGRF is based primarily on EISPACK routines. The matrix is
c first balanced using the Parlett-Reinsch algorithm. Then
c the Martin-Wilkinson algorithm is applied.
c There is a statement 'J = WK( I )' that converts a REAL
c variable to an integer variable, that seems dangerous
c to us in principle, but seems to work fine in practice.
c References:
c Dongarra, J. and C. Moler, EISPACK -- A Package for Solving
c Matrix Eigenvalue Problems, in Cowell, ed., 1984:
c Sources and Development of Mathematical Software,
c Prentice-Hall, Englewood Cliffs, NJ
c Parlett and Reinsch, 1969: Balancing a Matrix for Calculation
c of Eigenvalues and Eigenvectors, Num. Math. 13, 293-304
c Wilkinson, J., 1965: The Algebraic Eigenvalue Problem,
c Clarendon Press, Oxford
c I N P U T V A R I A B L E S:
c
c AA : input asymmetric matrix, destroyed after solved
c
c M : order of AA
c
c IA : first dimension of AA
c
c IEVEC : first dimension of EVEC
c
c
c O U T P U T V A R I A B L E S:
c
c EVEC : (unnormalized) eigenvectors of AA
c ( column J corresponds to EVAL(J) )
c
c EVAL : (unordered) eigenvalues of AA ( dimension at least M )
c
c IER : if .NE. 0, signals that EVAL(IER) failed to converge;
c in that case eigenvalues IER+1,IER+2,...,M are
c correct but eigenvalues 1,...,IER are set to zero.
c
c
c S C R A T C H V A R I A B L E S:
c
c WK : work area ( dimension at least 2*M )
c
c Called by- SOLEIG
c Calls- R1MACH, ERRMSG
c +-------------------------------------------------------------------+
c .. Scalar Arguments ..
use params, only: kr
INTEGER IA, IER, IEVEC, M
c ..
c .. Array Arguments ..
REAL(KR) AA( IA, M ), EVAL( M ), EVEC( IEVEC, M ), WK( * )
c ..
c .. Local Scalars ..
LOGICAL NOCONV, NOTLAS
INTEGER I, II, IN, J, K, KA, KKK, L, LB, LLL, N, N1, N2
REAL(KR) C1, C2, C3, C4, C5, C6, COL, DISCRI, F, G, H, ONE,
& P, Q, R, REPL, RNORM, ROW, S, SCALE, SGN, T, TOL, UU, VV, W,
& X, Y, Z, ZERO
c ..
c .. External Functions ..
REAL(KR) R1MACH
EXTERNAL R1MACH
c ..
c .. External Subroutines ..
EXTERNAL ERRMSG
c ..
c .. Intrinsic Functions ..
INTRINSIC ABS, MIN, SIGN, SQRT
c ..
DATA C1 / 0.4375 / , C2 / 0.5 / , C3 / 0.75 / , C4 / 0.95 / ,
& C5 / 16.0 / , C6 / 256.0 / , ZERO / 0.0 / , ONE / 1.0 /
c 2017 new inititializations
DATA P /0.0_kr/, R /0.0_kr/, Q /0.0_kr/
IER = 0
TOL = R1MACH( 4 )
IF( M.LT.1 .OR. IA.LT.M .OR. IEVEC.LT.M )
& call errmsg(0, 'ASYMTX--bad input variable(s)')
c ** Handle 1x1 and 2x2 special cases
IF( M.EQ.1 ) THEN
EVAL( 1 ) = AA( 1,1 )
EVEC( 1,1 ) = ONE
RETURN
ELSE IF( M.EQ.2 ) THEN
DISCRI = ( AA( 1,1 ) - AA( 2,2 ) )**2 + 4.*AA( 1,2 )*AA( 2,1 )
IF( DISCRI .LT. ZERO )
& call errmsg(0, 'ASYMTX--complex evals in 2x2 case')
SGN = ONE
IF( AA( 1,1 ) .LT. AA( 2,2 ) ) SGN = - ONE
EVAL( 1 ) = 0.5*( AA( 1,1 ) + AA( 2,2 ) + SGN*SQRT( DISCRI ) )