-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathrrtmg_lw_init.f90
2632 lines (2406 loc) · 105 KB
/
rrtmg_lw_init.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
! path: $Source: /storm/rc1/cvsroot/rc/rrtmg_lw/src/rrtmg_lw_init.f90,v $
! author: $Author: mike $
! revision: $Revision: 1.2 $
! created: $Date: 2007/08/22 19:20:03 $
!
module rrtmg_lw_init
! --------------------------------------------------------------------------
! | |
! | Copyright 2002-2007, Atmospheric & Environmental Research, Inc. (AER). |
! | This software may be used, copied, or redistributed as long as it is |
! | not sold and this copyright notice is reproduced on each copy made. |
! | This model is provided as is without any express or implied warranties. |
! | (http://www.rtweb.aer.com/) |
! | |
! --------------------------------------------------------------------------
use shr_kind_mod, only: r8 => shr_kind_r8
! use parkind, only : jpim, jprb
use rrlw_wvn
use rrtmg_lw_setcoef, only: lwatmref, lwavplank
implicit none
contains
! **************************************************************************
subroutine rrtmg_lw_ini
! **************************************************************************
!
! Original version: Michael J. Iacono; July, 1998
! First revision for NCAR CCM: September, 1998
! Second revision for RRTM_V3.0: September, 2002
!
! This subroutine performs calculations necessary for the initialization
! of the longwave model. Lookup tables are computed for use in the LW
! radiative transfer, and input absorption coefficient data for each
! spectral band are reduced from 256 g-point intervals to 140.
! **************************************************************************
use parrrtm, only : mg, nbndlw, ngptlw
use rrlw_tbl, only: ntbl, tblint, pade, bpade, tau_tbl, exp_tbl, tfn_tbl
use rrlw_vsn, only: hvrini, hnamini
! ------- Local -------
integer :: itr, ibnd, igc, ig, ind, ipr
integer :: igcsm, iprsm
real(kind=r8) :: wtsum, wtsm(mg) !
real(kind=r8) :: tfn !
! ------- Definitions -------
! Arrays for 10000-point look-up tables:
! TAU_TBL Clear-sky optical depth (used in cloudy radiative transfer)
! EXP_TBL Exponential lookup table for ransmittance
! TFN_TBL Tau transition function; i.e. the transition of the Planck
! function from that for the mean layer temperature to that for
! the layer boundary temperature as a function of optical depth.
! The "linear in tau" method is used to make the table.
! PADE Pade approximation constant (= 0.278)
! BPADE Inverse of the Pade approximation constant
!
hvrini = '$Revision: 1.2 $'
! Initialize model data
call lwdatinit
call lwcmbdat ! g-point interval reduction data
call lwcldpr ! cloud optical properties
call lwatmref ! reference MLS profile
call lwavplank ! Planck function
call lw_kgb01 ! molecular absorption coefficients
call lw_kgb02
call lw_kgb03
call lw_kgb04
call lw_kgb05
call lw_kgb06
call lw_kgb07
call lw_kgb08
call lw_kgb09
call lw_kgb10
call lw_kgb11
call lw_kgb12
call lw_kgb13
call lw_kgb14
call lw_kgb15
call lw_kgb16
! Compute lookup tables for transmittance, tau transition function,
! and clear sky tau (for the cloudy sky radiative transfer). Tau is
! computed as a function of the tau transition function, transmittance
! is calculated as a function of tau, and the tau transition function
! is calculated using the linear in tau formulation at values of tau
! above 0.01. TF is approximated as tau/6 for tau < 0.01. All tables
! are computed at intervals of 0.001. The inverse of the constant used
! in the Pade approximation to the tau transition function is set to b.
tau_tbl(0) = 0.0_r8
tau_tbl(ntbl) = 1.e10_r8
exp_tbl(0) = 1.0_r8
exp_tbl(ntbl) = 0.0_r8
tfn_tbl(0) = 0.0_r8
tfn_tbl(ntbl) = 1.0_r8
bpade = 1.0_r8 / pade
do itr = 1, ntbl-1
tfn = float(itr) / float(ntbl)
tau_tbl(itr) = bpade * tfn / (1._r8 - tfn)
exp_tbl(itr) = exp(-tau_tbl(itr))
if (tau_tbl(itr) .lt. 0.06_r8) then
tfn_tbl(itr) = tau_tbl(itr)/6._r8
else
tfn_tbl(itr) = 1._r8-2._r8*((1._r8/tau_tbl(itr))-(exp_tbl(itr)/(1.-exp_tbl(itr))))
endif
enddo
! Perform g-point reduction from 16 per band (256 total points) to
! a band dependant number (140 total points) for all absorption
! coefficient input data and Planck fraction input data.
! Compute relative weighting for new g-point combinations.
igcsm = 0
do ibnd = 1,nbndlw
iprsm = 0
if (ngc(ibnd).lt.mg) then
do igc = 1,ngc(ibnd)
igcsm = igcsm + 1
wtsum = 0._r8
do ipr = 1, ngn(igcsm)
iprsm = iprsm + 1
wtsum = wtsum + wt(iprsm)
enddo
wtsm(igc) = wtsum
enddo
do ig = 1, ng(ibnd)
ind = (ibnd-1)*mg + ig
rwgt(ind) = wt(ig)/wtsm(ngm(ind))
enddo
else
do ig = 1, ng(ibnd)
igcsm = igcsm + 1
ind = (ibnd-1)*mg + ig
rwgt(ind) = 1.0_r8
enddo
endif
enddo
! Reduce g-points for absorption coefficient data in each LW spectral band.
call cmbgb1
call cmbgb2
call cmbgb3
call cmbgb4
call cmbgb5
call cmbgb6
call cmbgb7
call cmbgb8
call cmbgb9
call cmbgb10
call cmbgb11
call cmbgb12
call cmbgb13
call cmbgb14
call cmbgb15
call cmbgb16
end subroutine rrtmg_lw_ini
!***************************************************************************
subroutine lwdatinit
!***************************************************************************
! --------- Modules ----------
use parrrtm, only : maxxsec, maxinpx
use rrlw_con, only: heatfac, grav, planck, boltz, &
clight, avogad, alosmt, gascon, radcn1, radcn2
use shr_const_mod, only: shr_const_avogad
use physconst, only: cday, gravit, cpair
use rrlw_vsn
save
! Longwave spectral band limits (wavenumbers)
wavenum1(:) = (/ 10._r8, 350._r8, 500._r8, 630._r8, 700._r8, 820._r8, &
980._r8,1080._r8,1180._r8,1390._r8,1480._r8,1800._r8, &
2080._r8,2250._r8,2390._r8,2600._r8/)
wavenum2(:) = (/350._r8, 500._r8, 630._r8, 700._r8, 820._r8, 980._r8, &
1080._r8,1180._r8,1390._r8,1480._r8,1800._r8,2080._r8, &
2250._r8,2390._r8,2600._r8,3250._r8/)
delwave(:) = (/340._r8, 150._r8, 130._r8, 70._r8, 120._r8, 160._r8, &
100._r8, 100._r8, 210._r8, 90._r8, 320._r8, 280._r8, &
170._r8, 130._r8, 220._r8, 650._r8/)
! Spectral band information
ng(:) = (/16,16,16,16,16,16,16,16,16,16,16,16,16,16,16,16/)
nspa(:) = (/1,1,9,9,9,1,9,1,9,1,1,9,9,1,9,9/)
nspb(:) = (/1,1,5,5,5,0,1,1,1,1,1,0,0,1,0,0/)
! Use constants set in CAM for consistency
grav = gravit
avogad = shr_const_avogad * 1.e-3_r8
! Heatfac is the factor by which one must multiply delta-flux/
! delta-pressure, with flux in w/m-2 and pressure in mbar, to get
! the heating rate in units of degrees/day. It is equal to
! (g)x(#sec/day)x(1e-5)/(specific heat of air at const. p)
! = (9.8066)(86400)(1e-5)/(1.004)
! heatfac = 8.4391_r8
! Modified values for consistency with CAM:
! = (9.80616)(86400)(1e-5)/(1.00464)
! heatfac = 8.43339130434_r8
! Calculate heatfac directly from CAM constants:
heatfac = grav * cday * 1.e-5_r8 / (cpair * 1.e-3_r8)
! nxmol - number of cross-sections input by user
! ixindx(i) - index of cross-section molecule corresponding to Ith
! cross-section specified by user
! = 0 -- not allowed in rrtm
! = 1 -- ccl4
! = 2 -- cfc11
! = 3 -- cfc12
! = 4 -- cfc22
nxmol = 4
ixindx(1) = 1
ixindx(2) = 2
ixindx(3) = 3
ixindx(4) = 4
ixindx(5:maxinpx) = 0
! Constants from NIST 01/11/2002
! grav = 9.8066_r8
planck = 6.62606876e-27_r8
boltz = 1.3806503e-16_r8
clight = 2.99792458e+10_r8
! avogad = 6.02214199e+23_r8
alosmt = 2.6867775e+19_r8
gascon = 8.31447200e+07_r8
radcn1 = 1.191042722e-12_r8
radcn2 = 1.4387752_r8
!
! units are generally cgs
!
! The first and second radiation constants are taken from NIST.
! They were previously obtained from the relations:
! radcn1 = 2.*planck*clight*clight*1.e-07
! radcn2 = planck*clight/boltz
end subroutine lwdatinit
!***************************************************************************
subroutine lwcmbdat
!***************************************************************************
save
! ------- Definitions -------
! Arrays for the g-point reduction from 256 to 140 for the 16 LW bands:
! This mapping from 256 to 140 points has been carefully selected to
! minimize the effect on the resulting fluxes and cooling rates, and
! caution should be used if the mapping is modified. The full 256
! g-point set can be restored with ngptlw=256, ngc=16*16, ngn=256*1., etc.
! ngptlw The total number of new g-points
! ngc The number of new g-points in each band
! ngs The cumulative sum of new g-points for each band
! ngm The index of each new g-point relative to the original
! 16 g-points for each band.
! ngn The number of original g-points that are combined to make
! each new g-point in each band.
! ngb The band index for each new g-point.
! wt RRTM weights for 16 g-points.
! ------- Data statements -------
ngc(:) = (/10,12,16,14,16,8,12,8,12,6,8,8,4,2,2,2/)
ngs(:) = (/10,22,38,52,68,76,88,96,108,114,122,130,134,136,138,140/)
ngm(:) = (/1,2,3,3,4,4,5,5,6,6,7,7,8,8,9,10, & ! band 1
1,2,3,4,5,6,7,8,9,9,10,10,11,11,12,12, & ! band 2
1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 3
1,2,3,4,5,6,7,8,9,10,11,12,13,14,14,14, & ! band 4
1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16, & ! band 5
1,1,2,2,3,3,4,4,5,5,6,6,7,7,8,8, & ! band 6
1,1,2,2,3,4,5,6,7,8,9,10,11,11,12,12, & ! band 7
1,1,2,2,3,3,4,4,5,5,6,6,7,7,8,8, & ! band 8
1,2,3,4,5,6,7,8,9,9,10,10,11,11,12,12, & ! band 9
1,1,2,2,3,3,4,4,5,5,5,5,6,6,6,6, & ! band 10
1,2,3,3,4,4,5,5,6,6,7,7,7,8,8,8, & ! band 11
1,2,3,4,5,5,6,6,7,7,7,7,8,8,8,8, & ! band 12
1,1,1,2,2,2,3,3,3,3,4,4,4,4,4,4, & ! band 13
1,1,1,1,1,1,1,1,2,2,2,2,2,2,2,2, & ! band 14
1,1,1,1,1,1,1,1,2,2,2,2,2,2,2,2, & ! band 15
1,1,1,1,2,2,2,2,2,2,2,2,2,2,2,2/) ! band 16
ngn(:) = (/1,1,2,2,2,2,2,2,1,1, & ! band 1
1,1,1,1,1,1,1,1,2,2,2,2, & ! band 2
1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 3
1,1,1,1,1,1,1,1,1,1,1,1,1,3, & ! band 4
1,1,1,1,1,1,1,1,1,1,1,1,1,1,1,1, & ! band 5
2,2,2,2,2,2,2,2, & ! band 6
2,2,1,1,1,1,1,1,1,1,2,2, & ! band 7
2,2,2,2,2,2,2,2, & ! band 8
1,1,1,1,1,1,1,1,2,2,2,2, & ! band 9
2,2,2,2,4,4, & ! band 10
1,1,2,2,2,2,3,3, & ! band 11
1,1,1,1,2,2,4,4, & ! band 12
3,3,4,6, & ! band 13
8,8, & ! band 14
8,8, & ! band 15
4,12/) ! band 16
ngb(:) = (/1,1,1,1,1,1,1,1,1,1, & ! band 1
2,2,2,2,2,2,2,2,2,2,2,2, & ! band 2
3,3,3,3,3,3,3,3,3,3,3,3,3,3,3,3, & ! band 3
4,4,4,4,4,4,4,4,4,4,4,4,4,4, & ! band 4
5,5,5,5,5,5,5,5,5,5,5,5,5,5,5,5, & ! band 5
6,6,6,6,6,6,6,6, & ! band 6
7,7,7,7,7,7,7,7,7,7,7,7, & ! band 7
8,8,8,8,8,8,8,8, & ! band 8
9,9,9,9,9,9,9,9,9,9,9,9, & ! band 9
10,10,10,10,10,10, & ! band 10
11,11,11,11,11,11,11,11, & ! band 11
12,12,12,12,12,12,12,12, & ! band 12
13,13,13,13, & ! band 13
14,14, & ! band 14
15,15, & ! band 15
16,16/) ! band 16
wt(:) = (/ 0.1527534276_r8, 0.1491729617_r8, 0.1420961469_r8, &
0.1316886544_r8, 0.1181945205_r8, 0.1019300893_r8, &
0.0832767040_r8, 0.0626720116_r8, 0.0424925000_r8, &
0.0046269894_r8, 0.0038279891_r8, 0.0030260086_r8, &
0.0022199750_r8, 0.0014140010_r8, 0.0005330000_r8, &
0.0000750000_r8/)
end subroutine lwcmbdat
!***************************************************************************
subroutine cmbgb1
!***************************************************************************
!
! Original version: MJIacono; July 1998
! Revision for GCMs: MJIacono; September 1998
! Revision for RRTMG: MJIacono, September 2002
! Revision for F90 reformatting: MJIacono, June 2006
!
! The subroutines CMBGB1->CMBGB16 input the absorption coefficient
! data for each band, which are defined for 16 g-points and 16 spectral
! bands. The data are combined with appropriate weighting following the
! g-point mapping arrays specified in RRTMINIT. Plank fraction data
! in arrays FRACREFA and FRACREFB are combined without weighting. All
! g-point reduced data are put into new arrays for use in RRTM.
!
! band 1: 10-350 cm-1 (low key - h2o; low minor - n2)
! (high key - h2o; high minor - n2)
! note: previous versions of rrtm band 1:
! 10-250 cm-1 (low - h2o; high - h2o)
!***************************************************************************
use parrrtm, only : mg, nbndlw, ngptlw, ng1
use rrlw_kg01, only: fracrefao, fracrefbo, kao, kbo, kao_mn2, kbo_mn2, &
selfrefo, forrefo, &
fracrefa, fracrefb, ka, kb, ka_mn2, kb_mn2, &
selfref, forref
! ------- Local -------
integer :: jt, jp, igc, ipr, iprsm
real(kind=r8) :: sumk, sumk1, sumk2, sumf1, sumf2
do jt = 1,5
do jp = 1,13
iprsm = 0
do igc = 1,ngc(1)
sumk = 0.
do ipr = 1, ngn(igc)
iprsm = iprsm + 1
sumk = sumk + kao(jt,jp,iprsm)*rwgt(iprsm)
enddo
ka(jt,jp,igc) = sumk
enddo
enddo
do jp = 13,59
iprsm = 0
do igc = 1,ngc(1)
sumk = 0.
do ipr = 1, ngn(igc)
iprsm = iprsm + 1
sumk = sumk + kbo(jt,jp,iprsm)*rwgt(iprsm)
enddo
kb(jt,jp,igc) = sumk
enddo
enddo
enddo
do jt = 1,10
iprsm = 0
do igc = 1,ngc(1)
sumk = 0.
do ipr = 1, ngn(igc)
iprsm = iprsm + 1
sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm)
enddo
selfref(jt,igc) = sumk
enddo
enddo
do jt = 1,4
iprsm = 0
do igc = 1,ngc(1)
sumk = 0.
do ipr = 1, ngn(igc)
iprsm = iprsm + 1
sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm)
enddo
forref(jt,igc) = sumk
enddo
enddo
do jt = 1,19
iprsm = 0
do igc = 1,ngc(1)
sumk1 = 0.
sumk2 = 0.
do ipr = 1, ngn(igc)
iprsm = iprsm + 1
sumk1 = sumk1 + kao_mn2(jt,iprsm)*rwgt(iprsm)
sumk2 = sumk2 + kbo_mn2(jt,iprsm)*rwgt(iprsm)
enddo
ka_mn2(jt,igc) = sumk1
kb_mn2(jt,igc) = sumk2
enddo
enddo
iprsm = 0
do igc = 1,ngc(1)
sumf1 = 0.
sumf2 = 0.
do ipr = 1, ngn(igc)
iprsm = iprsm + 1
sumf1= sumf1+ fracrefao(iprsm)
sumf2= sumf2+ fracrefbo(iprsm)
enddo
fracrefa(igc) = sumf1
fracrefb(igc) = sumf2
enddo
end subroutine cmbgb1
!***************************************************************************
subroutine cmbgb2
!***************************************************************************
!
! band 2: 350-500 cm-1 (low key - h2o; high key - h2o)
!
! note: previous version of rrtm band 2:
! 250 - 500 cm-1 (low - h2o; high - h2o)
!***************************************************************************
use parrrtm, only : mg, nbndlw, ngptlw, ng2
use rrlw_kg02, only: fracrefao, fracrefbo, kao, kbo, selfrefo, forrefo, &
fracrefa, fracrefb, ka, kb, selfref, forref
! ------- Local -------
integer :: jt, jp, igc, ipr, iprsm
real(kind=r8) :: sumk, sumf1, sumf2
do jt = 1,5
do jp = 1,13
iprsm = 0
do igc = 1,ngc(2)
sumk = 0.
do ipr = 1, ngn(ngs(1)+igc)
iprsm = iprsm + 1
sumk = sumk + kao(jt,jp,iprsm)*rwgt(iprsm+16)
enddo
ka(jt,jp,igc) = sumk
enddo
enddo
do jp = 13,59
iprsm = 0
do igc = 1,ngc(2)
sumk = 0.
do ipr = 1, ngn(ngs(1)+igc)
iprsm = iprsm + 1
sumk = sumk + kbo(jt,jp,iprsm)*rwgt(iprsm+16)
enddo
kb(jt,jp,igc) = sumk
enddo
enddo
enddo
do jt = 1,10
iprsm = 0
do igc = 1,ngc(2)
sumk = 0.
do ipr = 1, ngn(ngs(1)+igc)
iprsm = iprsm + 1
sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm+16)
enddo
selfref(jt,igc) = sumk
enddo
enddo
do jt = 1,4
iprsm = 0
do igc = 1,ngc(2)
sumk = 0.
do ipr = 1, ngn(ngs(1)+igc)
iprsm = iprsm + 1
sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm+16)
enddo
forref(jt,igc) = sumk
enddo
enddo
iprsm = 0
do igc = 1,ngc(2)
sumf1 = 0.
sumf2 = 0.
do ipr = 1, ngn(ngs(1)+igc)
iprsm = iprsm + 1
sumf1= sumf1+ fracrefao(iprsm)
sumf2= sumf2+ fracrefbo(iprsm)
enddo
fracrefa(igc) = sumf1
fracrefb(igc) = sumf2
enddo
end subroutine cmbgb2
!***************************************************************************
subroutine cmbgb3
!***************************************************************************
!
! band 3: 500-630 cm-1 (low key - h2o,co2; low minor - n2o)
! (high key - h2o,co2; high minor - n2o)
!
! old band 3: 500-630 cm-1 (low - h2o,co2; high - h2o,co2)
!***************************************************************************
use parrrtm, only : mg, nbndlw, ngptlw, ng3
use rrlw_kg03, only: fracrefao, fracrefbo, kao, kbo, kao_mn2o, kbo_mn2o, &
selfrefo, forrefo, &
fracrefa, fracrefb, ka, kb, ka_mn2o, kb_mn2o, &
selfref, forref
! ------- Local -------
integer :: jn, jt, jp, igc, ipr, iprsm
real(kind=r8) :: sumk, sumf
do jn = 1,9
do jt = 1,5
do jp = 1,13
iprsm = 0
do igc = 1,ngc(3)
sumk = 0.
do ipr = 1, ngn(ngs(2)+igc)
iprsm = iprsm + 1
sumk = sumk + kao(jn,jt,jp,iprsm)*rwgt(iprsm+32)
enddo
ka(jn,jt,jp,igc) = sumk
enddo
enddo
enddo
enddo
do jn = 1,5
do jt = 1,5
do jp = 13,59
iprsm = 0
do igc = 1,ngc(3)
sumk = 0.
do ipr = 1, ngn(ngs(2)+igc)
iprsm = iprsm + 1
sumk = sumk + kbo(jn,jt,jp,iprsm)*rwgt(iprsm+32)
enddo
kb(jn,jt,jp,igc) = sumk
enddo
enddo
enddo
enddo
do jn = 1,9
do jt = 1,19
iprsm = 0
do igc = 1,ngc(3)
sumk = 0.
do ipr = 1, ngn(ngs(2)+igc)
iprsm = iprsm + 1
sumk = sumk + kao_mn2o(jn,jt,iprsm)*rwgt(iprsm+32)
enddo
ka_mn2o(jn,jt,igc) = sumk
enddo
enddo
enddo
do jn = 1,5
do jt = 1,19
iprsm = 0
do igc = 1,ngc(3)
sumk = 0.
do ipr = 1, ngn(ngs(2)+igc)
iprsm = iprsm + 1
sumk = sumk + kbo_mn2o(jn,jt,iprsm)*rwgt(iprsm+32)
enddo
kb_mn2o(jn,jt,igc) = sumk
enddo
enddo
enddo
do jt = 1,10
iprsm = 0
do igc = 1,ngc(3)
sumk = 0.
do ipr = 1, ngn(ngs(2)+igc)
iprsm = iprsm + 1
sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm+32)
enddo
selfref(jt,igc) = sumk
enddo
enddo
do jt = 1,4
iprsm = 0
do igc = 1,ngc(3)
sumk = 0.
do ipr = 1, ngn(ngs(2)+igc)
iprsm = iprsm + 1
sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm+32)
enddo
forref(jt,igc) = sumk
enddo
enddo
do jp = 1,9
iprsm = 0
do igc = 1,ngc(3)
sumf = 0.
do ipr = 1, ngn(ngs(2)+igc)
iprsm = iprsm + 1
sumf = sumf + fracrefao(iprsm,jp)
enddo
fracrefa(igc,jp) = sumf
enddo
enddo
do jp = 1,5
iprsm = 0
do igc = 1,ngc(3)
sumf = 0.
do ipr = 1, ngn(ngs(2)+igc)
iprsm = iprsm + 1
sumf = sumf + fracrefbo(iprsm,jp)
enddo
fracrefb(igc,jp) = sumf
enddo
enddo
end subroutine cmbgb3
!***************************************************************************
subroutine cmbgb4
!***************************************************************************
!
! band 4: 630-700 cm-1 (low key - h2o,co2; high key - o3,co2)
!
! old band 4: 630-700 cm-1 (low - h2o,co2; high - o3,co2)
!***************************************************************************
use parrrtm, only : mg, nbndlw, ngptlw, ng4
use rrlw_kg04, only: fracrefao, fracrefbo, kao, kbo, selfrefo, forrefo, &
fracrefa, fracrefb, ka, kb, selfref, forref
! ------- Local -------
integer :: jn, jt, jp, igc, ipr, iprsm
real(kind=r8) :: sumk, sumf
do jn = 1,9
do jt = 1,5
do jp = 1,13
iprsm = 0
do igc = 1,ngc(4)
sumk = 0.
do ipr = 1, ngn(ngs(3)+igc)
iprsm = iprsm + 1
sumk = sumk + kao(jn,jt,jp,iprsm)*rwgt(iprsm+48)
enddo
ka(jn,jt,jp,igc) = sumk
enddo
enddo
enddo
enddo
do jn = 1,5
do jt = 1,5
do jp = 13,59
iprsm = 0
do igc = 1,ngc(4)
sumk = 0.
do ipr = 1, ngn(ngs(3)+igc)
iprsm = iprsm + 1
sumk = sumk + kbo(jn,jt,jp,iprsm)*rwgt(iprsm+48)
enddo
kb(jn,jt,jp,igc) = sumk
enddo
enddo
enddo
enddo
do jt = 1,10
iprsm = 0
do igc = 1,ngc(4)
sumk = 0.
do ipr = 1, ngn(ngs(3)+igc)
iprsm = iprsm + 1
sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm+48)
enddo
selfref(jt,igc) = sumk
enddo
enddo
do jt = 1,4
iprsm = 0
do igc = 1,ngc(4)
sumk = 0.
do ipr = 1, ngn(ngs(3)+igc)
iprsm = iprsm + 1
sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm+48)
enddo
forref(jt,igc) = sumk
enddo
enddo
do jp = 1,9
iprsm = 0
do igc = 1,ngc(4)
sumf = 0.
do ipr = 1, ngn(ngs(3)+igc)
iprsm = iprsm + 1
sumf = sumf + fracrefao(iprsm,jp)
enddo
fracrefa(igc,jp) = sumf
enddo
enddo
do jp = 1,5
iprsm = 0
do igc = 1,ngc(4)
sumf = 0.
do ipr = 1, ngn(ngs(3)+igc)
iprsm = iprsm + 1
sumf = sumf + fracrefbo(iprsm,jp)
enddo
fracrefb(igc,jp) = sumf
enddo
enddo
end subroutine cmbgb4
!***************************************************************************
subroutine cmbgb5
!***************************************************************************
!
! band 5: 700-820 cm-1 (low key - h2o,co2; low minor - o3, ccl4)
! (high key - o3,co2)
!
! old band 5: 700-820 cm-1 (low - h2o,co2; high - o3,co2)
!***************************************************************************
use parrrtm, only : mg, nbndlw, ngptlw, ng5
use rrlw_kg05, only: fracrefao, fracrefbo, kao, kbo, kao_mo3, ccl4o, &
selfrefo, forrefo, &
fracrefa, fracrefb, ka, kb, ka_mo3, ccl4, &
selfref, forref
! ------- Local -------
integer :: jn, jt, jp, igc, ipr, iprsm
real(kind=r8) :: sumk, sumf
do jn = 1,9
do jt = 1,5
do jp = 1,13
iprsm = 0
do igc = 1,ngc(5)
sumk = 0.
do ipr = 1, ngn(ngs(4)+igc)
iprsm = iprsm + 1
sumk = sumk + kao(jn,jt,jp,iprsm)*rwgt(iprsm+64)
enddo
ka(jn,jt,jp,igc) = sumk
enddo
enddo
enddo
enddo
do jn = 1,5
do jt = 1,5
do jp = 13,59
iprsm = 0
do igc = 1,ngc(5)
sumk = 0.
do ipr = 1, ngn(ngs(4)+igc)
iprsm = iprsm + 1
sumk = sumk + kbo(jn,jt,jp,iprsm)*rwgt(iprsm+64)
enddo
kb(jn,jt,jp,igc) = sumk
enddo
enddo
enddo
enddo
do jn = 1,9
do jt = 1,19
iprsm = 0
do igc = 1,ngc(5)
sumk = 0.
do ipr = 1, ngn(ngs(4)+igc)
iprsm = iprsm + 1
sumk = sumk + kao_mo3(jn,jt,iprsm)*rwgt(iprsm+64)
enddo
ka_mo3(jn,jt,igc) = sumk
enddo
enddo
enddo
do jt = 1,10
iprsm = 0
do igc = 1,ngc(5)
sumk = 0.
do ipr = 1, ngn(ngs(4)+igc)
iprsm = iprsm + 1
sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm+64)
enddo
selfref(jt,igc) = sumk
enddo
enddo
do jt = 1,4
iprsm = 0
do igc = 1,ngc(5)
sumk = 0.
do ipr = 1, ngn(ngs(4)+igc)
iprsm = iprsm + 1
sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm+64)
enddo
forref(jt,igc) = sumk
enddo
enddo
do jp = 1,9
iprsm = 0
do igc = 1,ngc(5)
sumf = 0.
do ipr = 1, ngn(ngs(4)+igc)
iprsm = iprsm + 1
sumf = sumf + fracrefao(iprsm,jp)
enddo
fracrefa(igc,jp) = sumf
enddo
enddo
do jp = 1,5
iprsm = 0
do igc = 1,ngc(5)
sumf = 0.
do ipr = 1, ngn(ngs(4)+igc)
iprsm = iprsm + 1
sumf = sumf + fracrefbo(iprsm,jp)
enddo
fracrefb(igc,jp) = sumf
enddo
enddo
iprsm = 0
do igc = 1,ngc(5)
sumk = 0.
do ipr = 1, ngn(ngs(4)+igc)
iprsm = iprsm + 1
sumk = sumk + ccl4o(iprsm)*rwgt(iprsm+64)
enddo
ccl4(igc) = sumk
enddo
end subroutine cmbgb5
!***************************************************************************
subroutine cmbgb6
!***************************************************************************
!
! band 6: 820-980 cm-1 (low key - h2o; low minor - co2)
! (high key - nothing; high minor - cfc11, cfc12)
!
! old band 6: 820-980 cm-1 (low - h2o; high - nothing)
!***************************************************************************
use parrrtm, only : mg, nbndlw, ngptlw, ng6
use rrlw_kg06, only: fracrefao, kao, kao_mco2, cfc11adjo, cfc12o, &
selfrefo, forrefo, &
fracrefa, ka, ka_mco2, cfc11adj, cfc12, &
selfref, forref
! ------- Local -------
integer :: jt, jp, igc, ipr, iprsm
real(kind=r8) :: sumk, sumf, sumk1, sumk2
do jt = 1,5
do jp = 1,13
iprsm = 0
do igc = 1,ngc(6)
sumk = 0.
do ipr = 1, ngn(ngs(5)+igc)
iprsm = iprsm + 1
sumk = sumk + kao(jt,jp,iprsm)*rwgt(iprsm+80)
enddo
ka(jt,jp,igc) = sumk
enddo
enddo
enddo
do jt = 1,19
iprsm = 0
do igc = 1,ngc(6)
sumk = 0.
do ipr = 1, ngn(ngs(5)+igc)
iprsm = iprsm + 1
sumk = sumk + kao_mco2(jt,iprsm)*rwgt(iprsm+80)
enddo
ka_mco2(jt,igc) = sumk
enddo
enddo
do jt = 1,10
iprsm = 0
do igc = 1,ngc(6)
sumk = 0.
do ipr = 1, ngn(ngs(5)+igc)
iprsm = iprsm + 1
sumk = sumk + selfrefo(jt,iprsm)*rwgt(iprsm+80)
enddo
selfref(jt,igc) = sumk
enddo
enddo
do jt = 1,4
iprsm = 0
do igc = 1,ngc(6)
sumk = 0.
do ipr = 1, ngn(ngs(5)+igc)
iprsm = iprsm + 1
sumk = sumk + forrefo(jt,iprsm)*rwgt(iprsm+80)
enddo
forref(jt,igc) = sumk
enddo
enddo
iprsm = 0
do igc = 1,ngc(6)
sumf = 0.
sumk1= 0.
sumk2= 0.
do ipr = 1, ngn(ngs(5)+igc)
iprsm = iprsm + 1
sumf = sumf + fracrefao(iprsm)
sumk1= sumk1+ cfc11adjo(iprsm)*rwgt(iprsm+80)
sumk2= sumk2+ cfc12o(iprsm)*rwgt(iprsm+80)
enddo
fracrefa(igc) = sumf
cfc11adj(igc) = sumk1
cfc12(igc) = sumk2
enddo
end subroutine cmbgb6
!***************************************************************************
subroutine cmbgb7
!***************************************************************************
!
! band 7: 980-1080 cm-1 (low key - h2o,o3; low minor - co2)
! (high key - o3; high minor - co2)
!
! old band 7: 980-1080 cm-1 (low - h2o,o3; high - o3)
!***************************************************************************
use parrrtm, only : mg, nbndlw, ngptlw, ng7
use rrlw_kg07, only: fracrefao, fracrefbo, kao, kbo, kao_mco2, kbo_mco2, &
selfrefo, forrefo, &
fracrefa, fracrefb, ka, kb, ka_mco2, kb_mco2, &
selfref, forref
! ------- Local -------
integer :: jn, jt, jp, igc, ipr, iprsm
real(kind=r8) :: sumk, sumf
do jn = 1,9
do jt = 1,5
do jp = 1,13