494 & ( plyr,plvl,tlyr,tlvl,qlyr,olyr, &
495 & gasvmr_co2,gasvmr_n2o,gasvmr_ch4,gasvmr_o2,gasvmr_co, &
496 & gasvmr_cfc11,gasvmr_cfc12,gasvmr_cfc22,gasvmr_ccl4, &
497 & icseed, aeraod, aerssa, aerasy, &
498 & sfcalb_nir_dir, sfcalb_nir_dif, &
499 & sfcalb_uvis_dir, sfcalb_uvis_dif, &
500 & dzlyr,delpin,de_lgth,alpha, &
501 & cosz,solcon,nday,idxday, &
502 & npts, nlay, nlp1, lprnt, inc_minor_gas, iswcliq, iswcice, &
503 & isubcsw, iovr, top_at_1, iswmode, cld_cf, lsswr, iovr_rand,&
504 & iovr_maxrand, iovr_max, iovr_dcorr, iovr_exp, iovr_exprand,&
505 & hswc,topflx,sfcflx,cldtau, &
506 & hsw0,hswb,flxprf,fdncmp, &
507 & cld_lwp, cld_ref_liq, cld_iwp, cld_ref_ice, &
508 & cld_rwp,cld_ref_rain, cld_swp, cld_ref_snow, &
509 & cld_od, cld_ssa, cld_asy, errmsg, errflg &
688 integer,
intent(in) :: npts, nlay, nlp1, nday, iswcliq, iswcice, &
689 isubcsw, iovr, iswmode, iovr_dcorr, iovr_exp, iovr_exprand, &
690 iovr_rand, iovr_maxrand, iovr_max
692 integer,
dimension(:),
intent(in) :: idxday
693 integer,
dimension(:),
intent(in),
optional :: icseed
695 logical,
intent(in) :: lprnt, lsswr, inc_minor_gas, top_at_1
697 real (kind=kind_phys),
dimension(:,:),
intent(in) :: &
699 real (kind=kind_phys),
dimension(:,:),
intent(in) :: &
700 & plyr, tlyr, qlyr, olyr, dzlyr, delpin
702 real (kind=kind_phys),
dimension(:),
intent(in):: sfcalb_nir_dir
703 real (kind=kind_phys),
dimension(:),
intent(in):: sfcalb_nir_dif
704 real (kind=kind_phys),
dimension(:),
intent(in):: sfcalb_uvis_dir
705 real (kind=kind_phys),
dimension(:),
intent(in):: sfcalb_uvis_dif
707 real(kind=kind_phys),
dimension(:,:),
intent(in)::gasvmr_co2
708 real(kind=kind_phys),
dimension(:,:),
intent(in)::gasvmr_n2o
709 real(kind=kind_phys),
dimension(:,:),
intent(in)::gasvmr_ch4
710 real(kind=kind_phys),
dimension(:,:),
intent(in)::gasvmr_o2
711 real(kind=kind_phys),
dimension(:,:),
intent(in)::gasvmr_co
712 real(kind=kind_phys),
dimension(:,:),
intent(in)::gasvmr_cfc11
713 real(kind=kind_phys),
dimension(:,:),
intent(in)::gasvmr_cfc12
714 real(kind=kind_phys),
dimension(:,:),
intent(in)::gasvmr_cfc22
715 real(kind=kind_phys),
dimension(:,:),
intent(in)::gasvmr_ccl4
717 real (kind=kind_phys),
dimension(:,:),
intent(in):: cld_cf
718 real (kind=kind_phys),
dimension(:,:),
intent(in),
optional:: &
719 & cld_lwp, cld_ref_liq, cld_iwp, cld_ref_ice, &
720 & cld_rwp, cld_ref_rain, cld_swp, cld_ref_snow, &
721 & cld_od, cld_ssa, cld_asy
723 real(kind=kind_phys),
dimension(:,:,:),
intent(in)::aeraod
724 real(kind=kind_phys),
dimension(:,:,:),
intent(in)::aerssa
725 real(kind=kind_phys),
dimension(:,:,:),
intent(in)::aerasy
727 real (kind=kind_phys),
intent(in) :: cosz(npts), solcon, &
729 real (kind=kind_phys),
dimension(npts,nlay),
intent(in) :: alpha
732 real (kind=kind_phys),
dimension(:,:),
intent(inout) :: hswc
733 real (kind=kind_phys),
dimension(:,:),
intent(inout) :: &
736 type (
topfsw_type),
dimension(:),
intent(inout) :: topflx
737 type (
sfcfsw_type),
dimension(:),
intent(inout) :: sfcflx
739 character(len=*),
intent(out) :: errmsg
740 integer,
intent(out) :: errflg
743 real (kind=kind_phys),
dimension(:,:,:),
optional, &
744 &
intent(inout) :: hswb
746 real (kind=kind_phys),
dimension(:,:),
optional, &
747 &
intent(inout) :: hsw0
749 & intent(inout) :: flxprf
751 & intent(inout) :: fdncmp
754 real (kind=kind_phys),
dimension(nlay,ngptsw) :: cldfmc, &
757 real (kind=kind_phys),
dimension(nlp1,nbdsw):: fxupc, fxdnc, &
760 real (kind=kind_phys),
dimension(nlay,nbdsw) :: &
761 & tauae, ssaae, asyae, taucw, ssacw, asycw
763 real (kind=kind_phys),
dimension(ngptsw) :: sfluxzen
765 real (kind=kind_phys),
dimension(nlay) :: cldfrc, delp, &
766 & pavel, tavel, coldry, colmol, h2ovmr, o3vmr, temcol, &
767 & cliqp, reliq, cicep, reice, cdat1, cdat2, cdat3, cdat4, &
768 & cfrac, fac00, fac01, fac10, fac11, forfac, forfrac, &
769 & selffac, selffrac, rfdelp, dz
771 real (kind=kind_phys),
dimension(nlp1) :: fnet, flxdc, flxuc, &
774 real (kind=kind_phys),
dimension(2) :: albbm, albdf, sfbmc, &
775 & sfbm0, sfdfc, sfdf0
777 real (kind=kind_phys) :: cosz1, sntz1, tem0, tem1, tem2, s0fac, &
778 & ssolar, zcf0, zcf1, ftoau0, ftoauc, ftoadc, &
779 & fsfcu0, fsfcuc, fsfcd0, fsfcdc, suvbfc, suvbf0, delgth
780 real (kind=kind_phys),
dimension(nlay) :: alph
784 real (kind=kind_phys) :: colamt(nlay,maxgas)
786 integer,
dimension(npts) :: ipseed
787 integer,
dimension(nlay) :: indfor, indself, jp, jt, jt1
789 integer :: i, ib, ipt, j1, k, kk, laytrop, mb, ig
790 integer :: inflgsw, iceflgsw, liqflgsw
791 integer :: irng, permuteseed
809 if (.not. lsswr)
return
810 if (nday <= 0)
return
812 lhswb =
present ( hswb )
813 lhsw0 =
present ( hsw0 )
814 lflxprf=
present ( flxprf )
815 lfdncmp=
present ( fdncmp )
828 sfcflx =
sfcfsw_type( f_zero, f_zero, f_zero, f_zero )
832 flxprf =
profsw_type( f_zero, f_zero, f_zero, f_zero )
836 fdncmp =
cmpfsw_type(f_zero,f_zero,f_zero,f_zero,f_zero,f_zero)
848 if (iswcliq > 0)
then
849 if ( .not.
present(cld_lwp) .or. .not.
present(cld_ref_liq) .or. &
850 & .not.
present(cld_iwp) .or. .not.
present(cld_ref_ice) .or. &
851 & .not.
present(cld_rwp) .or. .not.
present(cld_ref_rain) .or. &
852 & .not.
present(cld_swp) .or. .not.
present(cld_ref_snow) )
then
853 write(errmsg,
'(*(a))') &
854 &
'Logic error: iswcliq>0 requires the following', &
855 &
' optional arguments to be present:', &
856 &
' cld_lwp, cld_ref_liq, cld_iwp, cld_ref_ice,', &
857 &
' cld_rwp, cld_ref_rain, cld_swp, cld_ref_snow'
862 if ( .not.
present(cld_od) .or. .not.
present(cld_ssa) .or. &
863 & .not.
present(cld_asy))
then
864 write(errmsg,
'(*(a))') &
865 &
'Logic error: iswcliq<=0 requires the following', &
866 &
' optional arguments to be present:', &
867 &
' cld_od, cld_ssa, cld_asy'
876 if ( isubcsw == 1 )
then
878 ipseed(i) = ipsdsw0 + i
880 elseif ( isubcsw == 2 )
then
882 ipseed(i) = icseed(i)
887 write(0,*)
' In radsw, isubcsw, ipsdsw0,ipseed =', &
888 & isubcsw, ipsdsw0, ipseed
893 lab_do_ipt :
do ipt = 1, nday
898 sntz1 = f_one / cosz(j1)
899 ssolar = s0fac * cosz(j1)
900 if (iovr == iovr_dcorr) delgth = de_lgth(j1)
903 albbm(1) = sfcalb_nir_dir(j1)
904 albdf(1) = sfcalb_nir_dif(j1)
905 albbm(2) = sfcalb_uvis_dir(j1)
906 albdf(2) = sfcalb_uvis_dif(j1)
914 tem2 = 1.0e-20 * 1.0e3 * con_avgd
918 pavel(k) = plyr(j1,kk)
919 tavel(k) = tlyr(j1,kk)
920 delp(k) = delpin(j1,kk)
922 if (iovr == iovr_exp .or. iovr == iovr_exprand) alph(k) = alpha(j1,k)
935 h2ovmr(k)= max(f_zero,qlyr(j1,kk)*amdw/(f_one-qlyr(j1,kk)))
936 o3vmr(k)= max(f_zero,olyr(j1,kk)*amdo3)
938 tem0 = (f_one - h2ovmr(k))*con_amd + h2ovmr(k)*con_amw
939 coldry(k) = tem2 * delp(k) / (tem1*tem0*(f_one + h2ovmr(k)))
940 temcol(k) = 1.0e-12 * coldry(k)
942 colamt(k,1) = max(f_zero, coldry(k)*h2ovmr(k))
943 colamt(k,2) = max(temcol(k), coldry(k)*gasvmr_co2(j1,kk))
944 colamt(k,3) = max(f_zero, coldry(k)*o3vmr(k))
945 colmol(k) = coldry(k) + colamt(k,1)
951 if (inc_minor_gas)
then
954 colamt(k,4) = max(temcol(k), coldry(k)*gasvmr_n2o(j1,kk))
955 colamt(k,5) = max(temcol(k), coldry(k)*gasvmr_ch4(j1,kk))
956 colamt(k,6) = max(temcol(k), coldry(k)*gasvmr_o2(j1,kk))
961 colamt(k,4) = temcol(k)
962 colamt(k,5) = temcol(k)
963 colamt(k,6) = temcol(k)
973 tauae(k,ib) = aeraod(j1,kk,ib)
974 ssaae(k,ib) = aerssa(j1,kk,ib)
975 asyae(k,ib) = aerasy(j1,kk,ib)
980 if (iswcliq > 0)
then
983 cfrac(k) = cld_cf(j1,kk)
984 cliqp(k) = cld_lwp(j1,kk)
985 reliq(k) = cld_ref_liq(j1,kk)
986 cicep(k) = cld_iwp(j1,kk)
987 reice(k) = cld_ref_ice(j1,kk)
988 cdat1(k) = cld_rwp(j1,kk)
989 cdat2(k) = cld_ref_rain(j1,kk)
990 cdat3(k) = cld_swp(j1,kk)
991 cdat4(k) = cld_ref_snow(j1,kk)
996 cfrac(k) = cld_cf(j1,kk)
997 cdat1(k) = cld_od(j1,kk)
998 cdat2(k) = cld_ssa(j1,kk)
999 cdat3(k) = cld_asy(j1,kk)
1005 tem1 = 100.0 * con_g
1006 tem2 = 1.0e-20 * 1.0e3 * con_avgd
1009 pavel(k) = plyr(j1,k)
1010 tavel(k) = tlyr(j1,k)
1011 delp(k) = delpin(j1,k)
1013 if (iovr == iovr_exp .or. iovr == iovr_exprand) alph(k) = alpha(j1,k)
1021 h2ovmr(k)= max(f_zero,qlyr(j1,k)*amdw/(f_one-qlyr(j1,k)))
1022 o3vmr(k)= max(f_zero,olyr(j1,k)*amdo3)
1024 tem0 = (f_one - h2ovmr(k))*con_amd + h2ovmr(k)*con_amw
1025 coldry(k) = tem2 * delp(k) / (tem1*tem0*(f_one + h2ovmr(k)))
1026 temcol(k) = 1.0e-12 * coldry(k)
1028 colamt(k,1) = max(f_zero, coldry(k)*h2ovmr(k))
1029 colamt(k,2) = max(temcol(k), coldry(k)*gasvmr_co2(j1,k))
1030 colamt(k,3) = max(f_zero, coldry(k)*o3vmr(k))
1031 colmol(k) = coldry(k) + colamt(k,1)
1037 write(0,*)
' pavel=',pavel
1038 write(0,*)
' tavel=',tavel
1039 write(0,*)
' delp=',delp
1040 write(0,*)
' h2ovmr=',h2ovmr*1000
1041 write(0,*)
' o3vmr=',o3vmr*1000000
1048 if (inc_minor_gas)
then
1050 colamt(k,4) = max(temcol(k), coldry(k)*gasvmr_n2o(j1,k))
1051 colamt(k,5) = max(temcol(k), coldry(k)*gasvmr_ch4(j1,k))
1052 colamt(k,6) = max(temcol(k), coldry(k)*gasvmr_o2(j1,k))
1057 colamt(k,4) = temcol(k)
1058 colamt(k,5) = temcol(k)
1059 colamt(k,6) = temcol(k)
1068 tauae(k,ib) = aeraod(j1,k,ib)
1069 ssaae(k,ib) = aerssa(j1,k,ib)
1070 asyae(k,ib) = aerasy(j1,k,ib)
1074 if (iswcliq > 0)
then
1076 cfrac(k) = cld_cf(j1,k)
1077 cliqp(k) = cld_lwp(j1,k)
1078 reliq(k) = cld_ref_liq(j1,k)
1079 cicep(k) = cld_iwp(j1,k)
1080 reice(k) = cld_ref_ice(j1,k)
1081 cdat1(k) = cld_rwp(j1,k)
1082 cdat2(k) = cld_ref_rain(j1,k)
1083 cdat3(k) = cld_swp(j1,k)
1084 cdat4(k) = cld_ref_snow(j1,k)
1088 cfrac(k) = cld_cf(j1,k)
1089 cdat1(k) = cld_od(j1,k)
1090 cdat2(k) = cld_ssa(j1,k)
1091 cdat3(k) = cld_asy(j1,k)
1104 if (iovr == iovr_rand)
then
1106 zcf0 = zcf0 * (f_one - cfrac(k))
1108 else if (iovr == iovr_maxrand)
then
1110 if (cfrac(k) > ftiny)
then
1111 zcf1 = min( zcf1, f_one-cfrac(k) )
1112 elseif (zcf1 < f_one)
then
1118 else if (iovr >= 2)
then
1120 zcf0 = min( zcf0, f_one-cfrac(k) )
1124 if (zcf0 <= ftiny) zcf0 = f_zero
1125 if (zcf0 > oneminus) zcf0 = f_one
1131 if (zcf1 > f_zero)
then
1135 & ( cfrac,cliqp,reliq,cicep,reice,cdat1,cdat2,cdat3,cdat4, &
1136 & zcf1, nlay, ipseed(j1), dz, delgth, alph, iswcliq, iswcice,&
1139 & taucw, ssacw, asycw, cldfrc, cldfmc &
1148 cldtau(j1,kk) = taucw(k,10)
1152 cldtau(j1,k) = taucw(k,10)
1172 & ( pavel,tavel,h2ovmr, nlay,nlp1, &
1174 & laytrop,jp,jt,jt1,fac00,fac01,fac10,fac11, &
1175 & selffac,selffrac,indself,forfac,forfrac,indfor &
1182 & ( colamt,colmol,fac00,fac01,fac10,fac11,jp,jt,jt1,laytrop, &
1183 & forfac,forfrac,indfor,selffac,selffrac,indself, nlay, &
1185 & sfluxzen, taug, taur &
1194 if ( isubcsw <= 0 )
then
1198 & ( ssolar,cosz1,sntz1,albbm,albdf,sfluxzen,cldfrc, &
1199 & zcf1,zcf0,taug,taur,tauae,ssaae,asyae,taucw,ssacw,asycw, &
1200 & nlay, nlp1, iswmode, &
1202 & fxupc,fxdnc,fxup0,fxdn0, &
1203 & ftoauc,ftoau0,ftoadc,fsfcuc,fsfcu0,fsfcdc,fsfcd0, &
1204 & sfbmc,sfdfc,sfbm0,sfdf0,suvbfc,suvbf0 &
1211 & ( ssolar,cosz1,sntz1,albbm,albdf,sfluxzen,cldfmc, &
1212 & zcf1,zcf0,taug,taur,tauae,ssaae,asyae,taucw,ssacw,asycw, &
1213 & nlay, nlp1, iswmode, &
1215 & fxupc,fxdnc,fxup0,fxdn0, &
1216 & ftoauc,ftoau0,ftoadc,fsfcuc,fsfcu0,fsfcdc,fsfcd0, &
1217 & sfbmc,sfdfc,sfbm0,sfdf0,suvbfc,suvbf0 &
1230 flxuc(k) = flxuc(k) + fxupc(k,ib)
1231 flxdc(k) = flxdc(k) + fxdnc(k,ib)
1237 if ( lhsw0 .or. lflxprf )
then
1243 flxu0(k) = flxu0(k) + fxup0(k,ib)
1244 flxd0(k) = flxd0(k) + fxdn0(k,ib)
1252 rfdelp(k) = heatfac / delp(k)
1257 fdncmp(j1)%uvbf0 = suvbf0
1258 fdncmp(j1)%uvbfc = suvbfc
1261 fdncmp(j1)%nirbm = sfbmc(1)
1262 fdncmp(j1)%nirdf = sfdfc(1)
1263 fdncmp(j1)%visbm = sfbmc(2)
1264 fdncmp(j1)%visdf = sfdfc(2)
1269 topflx(j1)%upfxc = ftoauc
1270 topflx(j1)%dnfxc = ftoadc
1271 topflx(j1)%upfx0 = ftoau0
1273 sfcflx(j1)%upfxc = fsfcuc
1274 sfcflx(j1)%dnfxc = fsfcdc
1275 sfcflx(j1)%upfx0 = fsfcu0
1276 sfcflx(j1)%dnfx0 = fsfcd0
1282 fnet(1) = flxdc(1) - flxuc(1)
1286 fnet(k) = flxdc(k) - flxuc(k)
1287 hswc(j1,kk) = (fnet(k)-fnet(k-1)) * rfdelp(k-1)
1295 flxprf(j1,kk)%upfxc = flxuc(k)
1296 flxprf(j1,kk)%dnfxc = flxdc(k)
1297 flxprf(j1,kk)%upfx0 = flxu0(k)
1298 flxprf(j1,kk)%dnfx0 = flxd0(k)
1305 fnet(1) = flxd0(1) - flxu0(1)
1309 fnet(k) = flxd0(k) - flxu0(k)
1310 hsw0(j1,kk) = (fnet(k)-fnet(k-1)) * rfdelp(k-1)
1318 fnet(1) = fxdnc(1,mb) - fxupc(1,mb)
1322 fnet(k) = fxdnc(k,mb) - fxupc(k,mb)
1323 hswb(j1,kk,mb) = (fnet(k) - fnet(k-1)) * rfdelp(k-1)
1332 fnet(1) = flxdc(1) - flxuc(1)
1335 fnet(k) = flxdc(k) - flxuc(k)
1336 hswc(j1,k-1) = (fnet(k)-fnet(k-1)) * rfdelp(k-1)
1343 flxprf(j1,k)%upfxc = flxuc(k)
1344 flxprf(j1,k)%dnfxc = flxdc(k)
1345 flxprf(j1,k)%upfx0 = flxu0(k)
1346 flxprf(j1,k)%dnfx0 = flxd0(k)
1353 fnet(1) = flxd0(1) - flxu0(1)
1356 fnet(k) = flxd0(k) - flxu0(k)
1357 hsw0(j1,k-1) = (fnet(k)-fnet(k-1)) * rfdelp(k-1)
1365 fnet(1) = fxdnc(1,mb) - fxupc(1,mb)
1368 fnet(k+1) = fxdnc(k+1,mb) - fxupc(k+1,mb)
1369 hswb(j1,k,mb) = (fnet(k+1) - fnet(k)) * rfdelp(k)
1568 & ( cfrac,cliqp,reliq,cicep,reice,cdat1,cdat2,cdat3,cdat4, &
1569 & cf1, nlay, ipseed, dz, delgth, alpha, iswcliq, iswcice, &
1570 & isubcsw, iovr, taucw, ssacw, asycw, cldfrc, cldfmc &
1656 integer,
intent(in) :: nlay, ipseed, iswcliq, iswcice, isubcsw, &
1658 real (kind=kind_phys),
intent(in) :: cf1, delgth
1660 real (kind=kind_phys),
dimension(nlay),
intent(in) :: cliqp, &
1661 & reliq, cicep, reice, cdat1, cdat2, cdat3, cdat4, cfrac, dz
1662 real (kind=kind_phys),
dimension(nlay),
intent(in) :: alpha
1665 real (kind=kind_phys),
dimension(nlay,ngptsw),
intent(out) :: &
1667 real (kind=kind_phys),
dimension(nlay,nbdsw),
intent(out) :: &
1668 & taucw, ssacw, asycw
1669 real (kind=kind_phys),
dimension(nlay),
intent(out) :: cldfrc
1672 real (kind=kind_phys),
dimension(nblow:nbhgh) :: tauliq, tauice, &
1673 & ssaliq, ssaice, ssaran, ssasnw, asyliq, asyice, &
1675 real (kind=kind_phys),
dimension(nlay) :: cldf
1677 real (kind=kind_phys) :: dgeice, factor, fint, tauran, tausnw, &
1678 & cldliq, refliq, cldice, refice, cldran, cldsnw, refsnw, &
1679 & extcoliq, ssacoliq, asycoliq, extcoice, ssacoice, asycoice,&
1682 logical :: lcloudy(nlay,ngptsw)
1683 integer :: ia, ib, ig, jb, k, index
1690 taucw(k,ib) = f_zero
1692 asycw(k,ib) = f_zero
1698 lab_if_iswcliq :
if (iswcliq > 0)
then
1700 lab_do_k :
do k = 1, nlay
1701 lab_if_cld :
if (cfrac(k) > ftiny)
then
1719 dgesnw = 1.0315 * refsnw
1721 tauran = cldran * a0r
1729 if (cldsnw>f_zero .and. refsnw>10.0_kind_phys)
then
1731 tausnw = cldsnw*1.09087*(a0s + a1s/dgesnw)
1736 do ib = nblow, nbhgh
1737 ssaran(ib) = tauran * (f_one - b0r(ib))
1738 ssasnw(ib) = tausnw * (f_one - (b0s(ib)+b1s(ib)*dgesnw))
1739 asyran(ib) = ssaran(ib) * c0r(ib)
1740 asysnw(ib) = ssasnw(ib) * c0s(ib)
1750 if ( cldliq <= f_zero )
then
1751 do ib = nblow, nbhgh
1757 factor = refliq - 1.5
1758 index = max( 1, min( 57, int( factor ) ))
1759 fint = factor - float(index)
1761 if ( iswcliq == 1 )
then
1762 do ib = nblow, nbhgh
1763 extcoliq = max(f_zero, extliq1(index,ib) &
1764 & + fint*(extliq1(index+1,ib)-extliq1(index,ib)) )
1765 ssacoliq = max(f_zero, min(f_one, ssaliq1(index,ib) &
1766 & + fint*(ssaliq1(index+1,ib)-ssaliq1(index,ib)) ))
1768 asycoliq = max(f_zero, min(f_one, asyliq1(index,ib) &
1769 & + fint*(asyliq1(index+1,ib)-asyliq1(index,ib)) ))
1772 tauliq(ib) = cldliq * extcoliq
1773 ssaliq(ib) = tauliq(ib) * ssacoliq
1774 asyliq(ib) = ssaliq(ib) * asycoliq
1776 elseif ( iswcliq == 2 )
then
1777 do ib = nblow, nbhgh
1778 extcoliq = max(f_zero, extliq2(index,ib) &
1779 & + fint*(extliq2(index+1,ib)-extliq2(index,ib)) )
1780 ssacoliq = max(f_zero, min(f_one, ssaliq2(index,ib) &
1781 & + fint*(ssaliq2(index+1,ib)-ssaliq2(index,ib)) ))
1783 asycoliq = max(f_zero, min(f_one, asyliq2(index,ib) &
1784 & + fint*(asyliq2(index+1,ib)-asyliq2(index,ib)) ))
1787 tauliq(ib) = cldliq * extcoliq
1788 ssaliq(ib) = tauliq(ib) * ssacoliq
1789 asyliq(ib) = ssaliq(ib) * asycoliq
1796 if ( cldice <= f_zero )
then
1797 do ib = nblow, nbhgh
1807 if ( iswcice == 1 )
then
1808 refice = min(130.0_kind_phys,max(13.0_kind_phys,refice))
1810 do ib = nblow, nbhgh
1813 extcoice = max(f_zero, abari(ia)+bbari(ia)/refice )
1814 ssacoice = max(f_zero, min(f_one, &
1815 & f_one-cbari(ia)-dbari(ia)*refice ))
1816 asycoice = max(f_zero, min(f_one, &
1817 & ebari(ia)+fbari(ia)*refice ))
1820 tauice(ib) = cldice * extcoice
1821 ssaice(ib) = tauice(ib) * ssacoice
1822 asyice(ib) = ssaice(ib) * asycoice
1827 elseif ( iswcice == 2 )
then
1828 refice = min(131.0_kind_phys,max(5.0_kind_phys,refice))
1830 factor = (refice - 2.0) / 3.0
1831 index = max( 1, min( 42, int( factor ) ))
1832 fint = factor - float(index)
1834 do ib = nblow, nbhgh
1835 extcoice = max(f_zero, extice2(index,ib) &
1836 & + fint*(extice2(index+1,ib)-extice2(index,ib)) )
1837 ssacoice = max(f_zero, min(f_one, ssaice2(index,ib) &
1838 & + fint*(ssaice2(index+1,ib)-ssaice2(index,ib)) ))
1839 asycoice = max(f_zero, min(f_one, asyice2(index,ib) &
1840 & + fint*(asyice2(index+1,ib)-asyice2(index,ib)) ))
1843 tauice(ib) = cldice * extcoice
1844 ssaice(ib) = tauice(ib) * ssacoice
1845 asyice(ib) = ssaice(ib) * asycoice
1851 elseif ( iswcice == 3 )
then
1852 dgeice = max( 5.0, min( 140.0, 1.0315*refice ))
1854 factor = (dgeice - 2.0) / 3.0
1855 index = max( 1, min( 45, int( factor ) ))
1856 fint = factor - float(index)
1858 do ib = nblow, nbhgh
1859 extcoice = max(f_zero, extice3(index,ib) &
1860 & + fint*(extice3(index+1,ib)-extice3(index,ib)) )
1861 ssacoice = max(f_zero, min(f_one, ssaice3(index,ib) &
1862 & + fint*(ssaice3(index+1,ib)-ssaice3(index,ib)) ))
1863 asycoice = max(f_zero, min(f_one, asyice3(index,ib) &
1864 & + fint*(asyice3(index+1,ib)-asyice3(index,ib)) ))
1869 tauice(ib) = cldice * extcoice
1870 ssaice(ib) = tauice(ib) * ssacoice
1871 asyice(ib) = ssaice(ib) * asycoice
1879 taucw(k,ib) = tauliq(jb)+tauice(jb)+tauran+tausnw
1880 ssacw(k,ib) = ssaliq(jb)+ssaice(jb)+ssaran(jb)+ssasnw(jb)
1881 asycw(k,ib) = asyliq(jb)+asyice(jb)+asyran(jb)+asysnw(jb)
1890 if (cfrac(k) > ftiny)
then
1892 taucw(k,ib) = cdat1(k)
1893 ssacw(k,ib) = cdat1(k) * cdat2(k)
1894 asycw(k,ib) = ssacw(k,ib) * cdat3(k)
1899 endif lab_if_iswcliq
1904 if ( isubcsw > 0 )
then
1907 where (cldf(:) < ftiny)
1915 & ( cldf, nlay, ipseed, dz, delgth, alpha, iovr, &
1922 if ( lcloudy(k,ig) )
then
1923 cldfmc(k,ig) = f_one
1925 cldfmc(k,ig) = f_zero
1933 cldfrc(k) = cfrac(k) / cf1
2434 & ( ssolar,cosz,sntz,albbm,albdf,sfluxzen,cldfrc, &
2435 & cf1,cf0,taug,taur,tauae,ssaae,asyae,taucw,ssacw,asycw, &
2436 & nlay, nlp1, iswmode, &
2437 & fxupc,fxdnc,fxup0,fxdn0, &
2438 & ftoauc,ftoau0,ftoadc,fsfcuc,fsfcu0,fsfcdc,fsfcd0, &
2439 & sfbmc,sfdfc,sfbm0,sfdf0,suvbfc,suvbf0 &
2533 real (kind=kind_phys),
parameter :: zcrit = 0.9999995
2534 real (kind=kind_phys),
parameter :: zsr3 = sqrt(3.0)
2535 real (kind=kind_phys),
parameter :: od_lo = 0.06
2536 real (kind=kind_phys),
parameter :: eps1 = 1.0e-8
2539 integer,
intent(in) :: nlay, nlp1, iswmode
2541 real (kind=kind_phys),
dimension(nlay,ngptsw),
intent(in) :: &
2543 real (kind=kind_phys),
dimension(nlay,nbdsw),
intent(in) :: &
2544 & taucw, ssacw, asycw, tauae, ssaae, asyae
2546 real (kind=kind_phys),
dimension(ngptsw),
intent(in) :: sfluxzen
2547 real (kind=kind_phys),
dimension(nlay),
intent(in) :: cldfrc
2549 real (kind=kind_phys),
dimension(2),
intent(in) :: albbm, albdf
2551 real (kind=kind_phys),
intent(in) :: cosz, sntz, cf1, cf0, ssolar
2554 real (kind=kind_phys),
dimension(nlp1,nbdsw),
intent(out) :: &
2555 & fxupc, fxdnc, fxup0, fxdn0
2557 real (kind=kind_phys),
dimension(2),
intent(out) :: sfbmc, sfdfc, &
2560 real (kind=kind_phys),
intent(out) :: suvbfc, suvbf0, ftoadc, &
2561 & ftoauc, ftoau0, fsfcuc, fsfcu0, fsfcdc, fsfcd0
2564 real (kind=kind_phys),
dimension(nlay) :: ztaus, zssas, zasys, &
2567 real (kind=kind_phys),
dimension(nlp1) :: zrefb, zrefd, ztrab, &
2568 & ztrad, ztdbt, zldbt, zfu, zfd
2570 real (kind=kind_phys) :: ztau1, zssa1, zasy1, ztau0, zssa0, &
2571 & zasy0, zasy3, zssaw, zasyw, zgam1, zgam2, zgam3, zgam4, &
2572 & zc0, zc1, za1, za2, zb1, zb2, zrk, zrk2, zrp, zrp1, zrm1, &
2573 & zrpp, zrkg1, zrkg3, zrkg4, zexp1, zexm1, zexp2, zexm2, &
2574 & zexp3, zexp4, zden1, ze1r45, ftind, zsolar, zrefb1, &
2575 & zrefd1, ztrab1, ztrad1, ztdbt0, zr1, zr2, zr3, zr4, zr5, &
2576 & zt1, zt2, zt3, zf1, zf2, zrpp1
2578 integer :: ib, ibd, jb, jg, k, kp, itind
2585 fxdnc(k,ib) = f_zero
2586 fxupc(k,ib) = f_zero
2587 fxdn0(k,ib) = f_zero
2588 fxup0(k,ib) = f_zero
2616 lab_do_jg :
do jg = 1, ngptsw
2622 zsolar = ssolar * sfluxzen(jg)
2631 zrefb(1) = albbm(ibd)
2632 zrefd(1) = albdf(ibd)
2634 zrefb(1) = 0.5 * (albbm(1) + albbm(2))
2635 zrefd(1) = 0.5 * (albdf(1) + albdf(2))
2654 ztau0 = max( ftiny, taur(k,jg)+taug(k,jg)+tauae(k,ib) )
2655 zssa0 = taur(k,jg) + tauae(k,ib)*ssaae(k,ib)
2656 zasy0 = asyae(k,ib)*ssaae(k,ib)*tauae(k,ib)
2657 zssaw = min( oneminus, zssa0 / ztau0 )
2658 zasyw = zasy0 / max( ftiny, zssa0 )
2669 ztau1 = (f_one - za2) * ztau0
2670 zssa1 = (zssaw - za2) / (f_one - za2)
2672 zasy1 = zasyw / (f_one + zasyw)
2673 zasy3 = 0.75 * zasy1
2681 if ( iswmode == 1 )
then
2682 zgam1 = 1.75 - zssa1 * (f_one + zasy3)
2683 zgam2 =-0.25 + zssa1 * (f_one - zasy3)
2684 zgam3 = 0.5 - zasy3 * cosz
2685 elseif ( iswmode == 2 )
then
2686 zgam1 = 2.0 - zssa1 * (1.25 + zasy3)
2687 zgam2 = 0.75* zssa1 * (f_one- zasy1)
2688 zgam3 = 0.5 - zasy3 * cosz
2689 elseif ( iswmode == 3 )
then
2690 zgam1 = zsr3 * (2.0 - zssa1 * (1.0 + zasy1)) * 0.5
2691 zgam2 = zsr3 * zssa1 * (1.0 - zasy1) * 0.5
2692 zgam3 = (1.0 - zsr3 * zasy1 * cosz) * 0.5
2694 zgam4 = f_one - zgam3
2699 if ( zssaw >= zcrit )
then
2700 za1 = zgam1 * cosz - zgam3
2706 zb1 = min( ztau1*sntz , 500.0 )
2707 if ( zb1 <= od_lo )
then
2708 zb2 = f_one - zb1 + 0.5*zb1*zb1
2710 ftind = zb1 / (bpade + zb1)
2711 itind = ftind*ntbmx + 0.5
2712 zb2 = exp_tbl(itind)
2716 zrefb(kp) = max(f_zero, min(f_one, &
2717 & (za2 - za1*(f_one - zb2))/(f_one + za2) ))
2718 ztrab(kp) = max(f_zero, min(f_one, f_one-zrefb(kp) ))
2721 zrefd(kp) = max(f_zero, min(f_one, za2/(f_one + za2) ))
2722 ztrad(kp) = max(f_zero, min(f_one, f_one-zrefd(kp) ))
2725 za1 = zgam1*zgam4 + zgam2*zgam3
2726 za2 = zgam1*zgam3 + zgam2*zgam4
2727 zrk = (zgam1 - zgam2) * (zgam1 + zgam2)
2728 if (zrk > eps1)
then
2738 zrpp1= f_one - zrp*zrp
2739 zrpp = sign( max(flimit, abs(zrpp1)), zrpp1 )
2744 zr1 = zrm1 * (za2 + zrkg3)
2745 zr2 = zrp1 * (za2 - zrkg3)
2746 zr3 = zrk2 * (zgam3 - za2*cosz)
2748 zr5 = zrpp * (zrk - zgam1)
2750 zt1 = zrp1 * (za1 + zrkg4)
2751 zt2 = zrm1 * (za1 - zrkg4)
2752 zt3 = zrk2 * (zgam4 + za1*cosz)
2757 zb1 = min( zrk*ztau1, 500.0 )
2758 if ( zb1 <= od_lo )
then
2759 zexm1 = f_one - zb1 + 0.5*zb1*zb1
2761 ftind = zb1 / (bpade + zb1)
2762 itind = ftind*ntbmx + 0.5
2763 zexm1 = exp_tbl(itind)
2765 zexp1 = f_one / zexm1
2767 zb2 = min( sntz*ztau1, 500.0 )
2768 if ( zb2 <= od_lo )
then
2769 zexm2 = f_one - zb2 + 0.5*zb2*zb2
2771 ftind = zb2 / (bpade + zb2)
2772 itind = ftind*ntbmx + 0.5
2773 zexm2 = exp_tbl(itind)
2775 zexp2 = f_one / zexm2
2776 ze1r45 = zr4*zexp1 + zr5*zexm1
2780 if (abs(ze1r45) <= eps1)
then
2784 zden1 = zssa1 / ze1r45
2785 zrefb(kp) = max(f_zero, min(f_one, &
2786 & (zr1*zexp1 - zr2*zexm1 - zr3*zexm2)*zden1 ))
2787 ztrab(kp) = max(f_zero, min(f_one, zexm2*(f_one &
2788 & - (zt1*zexp1 - zt2*zexm1 - zt3*zexp2)*zden1) ))
2792 if (ze1r45 >= f_zero)
then
2793 zden1 = zr4 / max(eps1, ze1r45*zrkg1)
2795 zden1 = zr4 / min(-eps1, ze1r45*zrkg1)
2797 zrefd(kp) = max(f_zero, min(f_one, &
2798 & zgam2*(zexp1 - zexm1)*zden1 ))
2799 ztrad(kp) = max(f_zero, min(f_one, zrk2*zden1 ))
2806 if ( zr1 <= od_lo )
then
2807 zexp3 = f_one - zr1 + 0.5*zr1*zr1
2809 ftind = zr1 / (bpade + zr1)
2810 itind = max(0, min(ntbmx, int(0.5+ntbmx*ftind) ))
2811 zexp3 = exp_tbl(itind)
2814 ztdbt(k) = zexp3 * ztdbt(kp)
2821 if ( zr1 <= od_lo )
then
2822 zexp4 = f_one - zr1 + 0.5*zr1*zr1
2824 ftind = zr1 / (bpade + zr1)
2825 itind = max(0, min(ntbmx, int(0.5+ntbmx*ftind) ))
2826 zexp4 = exp_tbl(itind)
2830 ztdbt0 = zexp4 * ztdbt0
2836 & ( zrefb,zrefd,ztrab,ztrad,zldbt,ztdbt, &
2844 fxup0(k,ib) = fxup0(k,ib) + zsolar*zfu(k)
2845 fxdn0(k,ib) = fxdn0(k,ib) + zsolar*zfd(k)
2850 zb2 = zsolar*(zfd(1) - ztdbt0)
2853 sfbm0(ibd) = sfbm0(ibd) + zb1
2854 sfdf0(ibd) = sfdf0(ibd) + zb2
2858 sfbm0(1) = sfbm0(1) + zf1
2859 sfdf0(1) = sfdf0(1) + zf2
2860 sfbm0(2) = sfbm0(2) + zf1
2861 sfdf0(2) = sfdf0(2) + zf2
2876 if ( cf1 > eps )
then
2884 zc0 = f_one - cldfrc(k)
2886 if ( zc1 > ftiny )
then
2888 ztau0 = ztaus(k) + taucw(k,ib)
2889 zssa0 = zssas(k) + ssacw(k,ib)
2890 zasy0 = zasys(k) + asycw(k,ib)
2891 zssaw = min(oneminus, zssa0 / ztau0)
2892 zasyw = zasy0 / max(ftiny, zssa0)
2898 ztau1 = (f_one - za2) * ztau0
2899 zssa1 = (zssaw - za2) / (f_one - za2)
2901 zasy1 = zasyw / (f_one + zasyw)
2902 zasy3 = 0.75 * zasy1
2911 if ( iswmode == 1 )
then
2912 zgam1 = 1.75 - zssa1 * (f_one + zasy3)
2913 zgam2 =-0.25 + zssa1 * (f_one - zasy3)
2914 zgam3 = 0.5 - zasy3 * cosz
2915 elseif ( iswmode == 2 )
then
2916 zgam1 = 2.0 - zssa1 * (1.25 + zasy3)
2917 zgam2 = 0.75* zssa1 * (f_one- zasy1)
2918 zgam3 = 0.5 - zasy3 * cosz
2919 elseif ( iswmode == 3 )
then
2920 zgam1 = zsr3 * (2.0 - zssa1 * (1.0 + zasy1)) * 0.5
2921 zgam2 = zsr3 * zssa1 * (1.0 - zasy1) * 0.5
2922 zgam3 = (1.0 - zsr3 * zasy1 * cosz) * 0.5
2924 zgam4 = f_one - zgam3
2934 if ( zssaw >= zcrit )
then
2935 za1 = zgam1 * cosz - zgam3
2941 zb1 = min( ztau1*sntz , 500.0 )
2942 if ( zb1 <= od_lo )
then
2943 zb2 = f_one - zb1 + 0.5*zb1*zb1
2945 ftind = zb1 / (bpade + zb1)
2946 itind = ftind*ntbmx + 0.5
2947 zb2 = exp_tbl(itind)
2951 zrefb(kp) = max(f_zero, min(f_one, &
2952 & (za2 - za1*(f_one - zb2))/(f_one + za2) ))
2953 ztrab(kp) = max(f_zero, min(f_one, f_one-zrefb(kp)))
2956 zrefd(kp) = max(f_zero, min(f_one, za2 / (f_one+za2) ))
2957 ztrad(kp) = max(f_zero, min(f_one, f_one - zrefd(kp) ))
2960 za1 = zgam1*zgam4 + zgam2*zgam3
2961 za2 = zgam1*zgam3 + zgam2*zgam4
2962 zrk = (zgam1 - zgam2) * (zgam1 + zgam2)
2963 if (zrk > eps1)
then
2973 zrpp1= f_one - zrp*zrp
2974 zrpp = sign( max(flimit, abs(zrpp1)), zrpp1 )
2979 zr1 = zrm1 * (za2 + zrkg3)
2980 zr2 = zrp1 * (za2 - zrkg3)
2981 zr3 = zrk2 * (zgam3 - za2*cosz)
2983 zr5 = zrpp * (zrk - zgam1)
2985 zt1 = zrp1 * (za1 + zrkg4)
2986 zt2 = zrm1 * (za1 - zrkg4)
2987 zt3 = zrk2 * (zgam4 + za1*cosz)
2992 zb1 = min( zrk*ztau1, 500.0 )
2993 if ( zb1 <= od_lo )
then
2994 zexm1 = f_one - zb1 + 0.5*zb1*zb1
2996 ftind = zb1 / (bpade + zb1)
2997 itind = ftind*ntbmx + 0.5
2998 zexm1 = exp_tbl(itind)
3000 zexp1 = f_one / zexm1
3002 zb2 = min( ztau1*sntz, 500.0 )
3003 if ( zb2 <= od_lo )
then
3004 zexm2 = f_one - zb2 + 0.5*zb2*zb2
3006 ftind = zb2 / (bpade + zb2)
3007 itind = ftind*ntbmx + 0.5
3008 zexm2 = exp_tbl(itind)
3010 zexp2 = f_one / zexm2
3011 ze1r45 = zr4*zexp1 + zr5*zexm1
3015 if ( abs(ze1r45) <= eps1 )
then
3019 zden1 = zssa1 / ze1r45
3020 zrefb(kp) = max(f_zero, min(f_one, &
3021 & (zr1*zexp1-zr2*zexm1-zr3*zexm2)*zden1 ))
3022 ztrab(kp) = max(f_zero, min(f_one, zexm2*(f_one - &
3023 & (zt1*zexp1-zt2*zexm1-zt3*zexp2)*zden1) ))
3027 if (ze1r45 >= f_zero)
then
3028 zden1 = zr4 / max(eps1, ze1r45*zrkg1)
3030 zden1 = zr4 / min(-eps1, ze1r45*zrkg1)
3032 zrefd(kp) = max(f_zero, min(f_one, &
3033 & zgam2*(zexp1 - zexm1)*zden1 ))
3034 ztrad(kp) = max(f_zero, min(f_one, zrk2*zden1 ))
3040 zrefb(kp) = zc0*zrefb1 + zc1*zrefb(kp)
3041 zrefd(kp) = zc0*zrefd1 + zc1*zrefd(kp)
3042 ztrab(kp) = zc0*ztrab1 + zc1*ztrab(kp)
3043 ztrad(kp) = zc0*ztrad1 + zc1*ztrad(kp)
3050 if ( zr1 <= od_lo )
then
3051 zexp3 = f_one - zr1 + 0.5*zr1*zr1
3053 ftind = zr1 / (bpade + zr1)
3054 itind = max(0, min(ntbmx, int(0.5+ntbmx*ftind) ))
3055 zexp3 = exp_tbl(itind)
3058 zldbt(kp) = zc0*zldbt(kp) + zc1*zexp3
3059 ztdbt(k) = zldbt(kp) * ztdbt(kp)
3065 if ( zr1 <= od_lo )
then
3066 zexp4 = f_one - zr1 + 0.5*zr1*zr1
3068 ftind = zr1 / (bpade + zr1)
3069 itind = max(0, min(ntbmx, int(0.5+ntbmx*ftind) ))
3070 zexp4 = exp_tbl(itind)
3073 ztdbt0 = (zc0*zldbt0(k) + zc1*zexp4) * ztdbt0
3078 ztdbt(k) = zldbt(kp) * ztdbt(kp)
3081 ztdbt0 = zldbt0(k) * ztdbt0
3090 & ( zrefb,zrefd,ztrab,ztrad,zldbt,ztdbt, &
3098 fxupc(k,ib) = fxupc(k,ib) + zsolar*zfu(k)
3099 fxdnc(k,ib) = fxdnc(k,ib) + zsolar*zfd(k)
3105 zb2 = zsolar*(zfd(1) - ztdbt0)
3108 sfbmc(ibd) = sfbmc(ibd) + zb1
3109 sfdfc(ibd) = sfdfc(ibd) + zb2
3113 sfbmc(1) = sfbmc(1) + zf1
3114 sfdfc(1) = sfdfc(1) + zf2
3115 sfbmc(2) = sfbmc(2) + zf1
3116 sfdfc(2) = sfdfc(2) + zf2
3128 ftoadc = ftoadc + fxdn0(nlp1,ib)
3129 ftoau0 = ftoau0 + fxup0(nlp1,ib)
3130 fsfcu0 = fsfcu0 + fxup0(1,ib)
3131 fsfcd0 = fsfcd0 + fxdn0(1,ib)
3135 ibd = nuvb - nblow + 1
3136 suvbf0 = fxdn0(1,ibd)
3138 if ( cf1 <= eps )
then
3141 fxupc(k,ib) = fxup0(k,ib)
3142 fxdnc(k,ib) = fxdn0(k,ib)
3161 fxupc(k,ib) = cf1*fxupc(k,ib) + cf0*fxup0(k,ib)
3162 fxdnc(k,ib) = cf1*fxdnc(k,ib) + cf0*fxdn0(k,ib)
3167 ftoauc = ftoauc + fxupc(nlp1,ib)
3168 fsfcuc = fsfcuc + fxupc(1,ib)
3169 fsfcdc = fsfcdc + fxdnc(1,ib)
3173 suvbfc = fxdnc(1,ibd)
3176 sfbmc(1) = cf1*sfbmc(1) + cf0*sfbm0(1)
3177 sfbmc(2) = cf1*sfbmc(2) + cf0*sfbm0(2)
3178 sfdfc(1) = cf1*sfdfc(1) + cf0*sfdf0(1)
3179 sfdfc(2) = cf1*sfdfc(2) + cf0*sfdf0(2)
3230 & ( ssolar,cosz,sntz,albbm,albdf,sfluxzen,cldfmc, &
3231 & cf1,cf0,taug,taur,tauae,ssaae,asyae,taucw,ssacw,asycw, &
3232 & nlay, nlp1, iswmode, &
3233 & fxupc,fxdnc,fxup0,fxdn0, &
3234 & ftoauc,ftoau0,ftoadc,fsfcuc,fsfcu0,fsfcdc,fsfcd0, &
3235 & sfbmc,sfdfc,sfbm0,sfdf0,suvbfc,suvbf0 &
3329 real (kind=kind_phys),
parameter :: zcrit = 0.9999995
3330 real (kind=kind_phys),
parameter :: zsr3 = sqrt(3.0)
3331 real (kind=kind_phys),
parameter :: od_lo = 0.06
3332 real (kind=kind_phys),
parameter :: eps1 = 1.0e-8
3335 integer,
intent(in) :: nlay, nlp1, iswmode
3337 real (kind=kind_phys),
dimension(nlay,ngptsw),
intent(in) :: &
3338 & taug, taur, cldfmc
3339 real (kind=kind_phys),
dimension(nlay,nbdsw),
intent(in) :: &
3340 & taucw, ssacw, asycw, tauae, ssaae, asyae
3342 real (kind=kind_phys),
dimension(ngptsw),
intent(in) :: sfluxzen
3344 real (kind=kind_phys),
dimension(2),
intent(in) :: albbm, albdf
3346 real (kind=kind_phys),
intent(in) :: cosz, sntz, cf1, cf0, ssolar
3349 real (kind=kind_phys),
dimension(nlp1,nbdsw),
intent(out) :: &
3350 & fxupc, fxdnc, fxup0, fxdn0
3352 real (kind=kind_phys),
dimension(2),
intent(out) :: sfbmc, sfdfc, &
3355 real (kind=kind_phys),
intent(out) :: suvbfc, suvbf0, ftoadc, &
3356 & ftoauc, ftoau0, fsfcuc, fsfcu0, fsfcdc, fsfcd0
3359 real (kind=kind_phys),
dimension(nlay) :: ztaus, zssas, zasys, &
3362 real (kind=kind_phys),
dimension(nlp1) :: zrefb, zrefd, ztrab, &
3363 & ztrad, ztdbt, zldbt, zfu, zfd
3365 real (kind=kind_phys) :: ztau1, zssa1, zasy1, ztau0, zssa0, &
3366 & zasy0, zasy3, zssaw, zasyw, zgam1, zgam2, zgam3, zgam4, &
3367 & za1, za2, zb1, zb2, zrk, zrk2, zrp, zrp1, zrm1, zrpp, &
3368 & zrkg1, zrkg3, zrkg4, zexp1, zexm1, zexp2, zexm2, zden1, &
3369 & zexp3, zexp4, ze1r45, ftind, zsolar, ztdbt0, zr1, zr2, &
3370 & zr3, zr4, zr5, zt1, zt2, zt3, zf1, zf2, zrpp1
3372 integer :: ib, ibd, jb, jg, k, kp, itind
3380 fxdnc(k,ib) = f_zero
3381 fxupc(k,ib) = f_zero
3382 fxdn0(k,ib) = f_zero
3383 fxup0(k,ib) = f_zero
3411 lab_do_jg :
do jg = 1, ngptsw
3417 zsolar = ssolar * sfluxzen(jg)
3426 zrefb(1) = albbm(ibd)
3427 zrefd(1) = albdf(ibd)
3429 zrefb(1) = 0.5 * (albbm(1) + albbm(2))
3430 zrefd(1) = 0.5 * (albdf(1) + albdf(2))
3448 ztau0 = max( ftiny, taur(k,jg)+taug(k,jg)+tauae(k,ib) )
3449 zssa0 = taur(k,jg) + tauae(k,ib)*ssaae(k,ib)
3450 zasy0 = asyae(k,ib)*ssaae(k,ib)*tauae(k,ib)
3451 zssaw = min( oneminus, zssa0 / ztau0 )
3452 zasyw = zasy0 / max( ftiny, zssa0 )
3463 ztau1 = (f_one - za2) * ztau0
3464 zssa1 = (zssaw - za2) / (f_one - za2)
3466 zasy1 = zasyw / (f_one + zasyw)
3467 zasy3 = 0.75 * zasy1
3475 if ( iswmode == 1 )
then
3476 zgam1 = 1.75 - zssa1 * (f_one + zasy3)
3477 zgam2 =-0.25 + zssa1 * (f_one - zasy3)
3478 zgam3 = 0.5 - zasy3 * cosz
3479 elseif ( iswmode == 2 )
then
3480 zgam1 = 2.0 - zssa1 * (1.25 + zasy3)
3481 zgam2 = 0.75* zssa1 * (f_one- zasy1)
3482 zgam3 = 0.5 - zasy3 * cosz
3483 elseif ( iswmode == 3 )
then
3484 zgam1 = zsr3 * (2.0 - zssa1 * (1.0 + zasy1)) * 0.5
3485 zgam2 = zsr3 * zssa1 * (1.0 - zasy1) * 0.5
3486 zgam3 = (1.0 - zsr3 * zasy1 * cosz) * 0.5
3488 zgam4 = f_one - zgam3
3492 if ( zssaw >= zcrit )
then
3493 za1 = zgam1 * cosz - zgam3
3499 zb1 = min( ztau1*sntz , 500.0 )
3500 if ( zb1 <= od_lo )
then
3501 zb2 = f_one - zb1 + 0.5*zb1*zb1
3503 ftind = zb1 / (bpade + zb1)
3504 itind = ftind*ntbmx + 0.5
3505 zb2 = exp_tbl(itind)
3509 zrefb(kp) = max(f_zero, min(f_one, &
3510 & (za2 - za1*(f_one - zb2))/(f_one + za2) ))
3511 ztrab(kp) = max(f_zero, min(f_one, f_one-zrefb(kp) ))
3514 zrefd(kp) = max(f_zero, min(f_one, za2/(f_one + za2) ))
3515 ztrad(kp) = max(f_zero, min(f_one, f_one-zrefd(kp) ))
3518 za1 = zgam1*zgam4 + zgam2*zgam3
3519 za2 = zgam1*zgam3 + zgam2*zgam4
3520 zrk = (zgam1 - zgam2) * (zgam1 + zgam2)
3521 if (zrk > eps1)
then
3531 zrpp1= f_one - zrp*zrp
3532 zrpp = sign( max(flimit, abs(zrpp1)), zrpp1 )
3537 zr1 = zrm1 * (za2 + zrkg3)
3538 zr2 = zrp1 * (za2 - zrkg3)
3539 zr3 = zrk2 * (zgam3 - za2*cosz)
3541 zr5 = zrpp * (zrk - zgam1)
3543 zt1 = zrp1 * (za1 + zrkg4)
3544 zt2 = zrm1 * (za1 - zrkg4)
3545 zt3 = zrk2 * (zgam4 + za1*cosz)
3550 zb1 = min( zrk*ztau1, 500.0 )
3551 if ( zb1 <= od_lo )
then
3552 zexm1 = f_one - zb1 + 0.5*zb1*zb1
3554 ftind = zb1 / (bpade + zb1)
3555 itind = ftind*ntbmx + 0.5
3556 zexm1 = exp_tbl(itind)
3558 zexp1 = f_one / zexm1
3560 zb2 = min( sntz*ztau1, 500.0 )
3561 if ( zb2 <= od_lo )
then
3562 zexm2 = f_one - zb2 + 0.5*zb2*zb2
3564 ftind = zb2 / (bpade + zb2)
3565 itind = ftind*ntbmx + 0.5
3566 zexm2 = exp_tbl(itind)
3568 zexp2 = f_one / zexm2
3569 ze1r45 = zr4*zexp1 + zr5*zexm1
3573 if (abs(ze1r45) <= eps1)
then
3577 zden1 = zssa1 / ze1r45
3578 zrefb(kp) = max(f_zero, min(f_one, &
3579 & (zr1*zexp1 - zr2*zexm1 - zr3*zexm2)*zden1 ))
3580 ztrab(kp) = max(f_zero, min(f_one, zexm2*(f_one &
3581 & - (zt1*zexp1 - zt2*zexm1 - zt3*zexp2)*zden1) ))
3585 if (ze1r45 >= f_zero)
then
3586 zden1 = zr4 / max(eps1, ze1r45*zrkg1)
3588 zden1 = zr4 / min(-eps1, ze1r45*zrkg1)
3590 zrefd(kp) = max(f_zero, min(f_one, &
3591 & zgam2*(zexp1 - zexm1)*zden1 ))
3592 ztrad(kp) = max(f_zero, min(f_one, zrk2*zden1 ))
3599 if ( zr1 <= od_lo )
then
3600 zexp3 = f_one - zr1 + 0.5*zr1*zr1
3602 ftind = zr1 / (bpade + zr1)
3603 itind = max(0, min(ntbmx, int(0.5+ntbmx*ftind) ))
3604 zexp3 = exp_tbl(itind)
3607 ztdbt(k) = zexp3 * ztdbt(kp)
3614 if ( zr1 <= od_lo )
then
3615 zexp4 = f_one - zr1 + 0.5*zr1*zr1
3617 ftind = zr1 / (bpade + zr1)
3618 itind = max(0, min(ntbmx, int(0.5+ntbmx*ftind) ))
3619 zexp4 = exp_tbl(itind)
3623 ztdbt0 = zexp4 * ztdbt0
3629 & ( zrefb,zrefd,ztrab,ztrad,zldbt,ztdbt, &
3637 fxup0(k,ib) = fxup0(k,ib) + zsolar*zfu(k)
3638 fxdn0(k,ib) = fxdn0(k,ib) + zsolar*zfd(k)
3643 zb2 = zsolar*(zfd(1) - ztdbt0)
3646 sfbm0(ibd) = sfbm0(ibd) + zb1
3647 sfdf0(ibd) = sfdf0(ibd) + zb2
3651 sfbm0(1) = sfbm0(1) + zf1
3652 sfdf0(1) = sfdf0(1) + zf2
3653 sfbm0(2) = sfbm0(2) + zf1
3654 sfdf0(2) = sfdf0(2) + zf2
3669 if ( cf1 > eps )
then
3677 if ( cldfmc(k,jg) > ftiny )
then
3679 ztau0 = ztaus(k) + taucw(k,ib)
3680 zssa0 = zssas(k) + ssacw(k,ib)
3681 zasy0 = zasys(k) + asycw(k,ib)
3682 zssaw = min(oneminus, zssa0 / ztau0)
3683 zasyw = zasy0 / max(ftiny, zssa0)
3689 ztau1 = (f_one - za2) * ztau0
3690 zssa1 = (zssaw - za2) / (f_one - za2)
3692 zasy1 = zasyw / (f_one + zasyw)
3693 zasy3 = 0.75 * zasy1
3696 if ( iswmode == 1 )
then
3697 zgam1 = 1.75 - zssa1 * (f_one + zasy3)
3698 zgam2 =-0.25 + zssa1 * (f_one - zasy3)
3699 zgam3 = 0.5 - zasy3 * cosz
3700 elseif ( iswmode == 2 )
then
3701 zgam1 = 2.0 - zssa1 * (1.25 + zasy3)
3702 zgam2 = 0.75* zssa1 * (f_one- zasy1)
3703 zgam3 = 0.5 - zasy3 * cosz
3704 elseif ( iswmode == 3 )
then
3705 zgam1 = zsr3 * (2.0 - zssa1 * (1.0 + zasy1)) * 0.5
3706 zgam2 = zsr3 * zssa1 * (1.0 - zasy1) * 0.5
3707 zgam3 = (1.0 - zsr3 * zasy1 * cosz) * 0.5
3709 zgam4 = f_one - zgam3
3714 if ( zssaw >= zcrit )
then
3715 za1 = zgam1 * cosz - zgam3
3721 zb1 = min( ztau1*sntz , 500.0 )
3722 if ( zb1 <= od_lo )
then
3723 zb2 = f_one - zb1 + 0.5*zb1*zb1
3725 ftind = zb1 / (bpade + zb1)
3726 itind = ftind*ntbmx + 0.5
3727 zb2 = exp_tbl(itind)
3731 zrefb(kp) = max(f_zero, min(f_one, &
3732 & (za2 - za1*(f_one - zb2))/(f_one + za2) ))
3733 ztrab(kp) = max(f_zero, min(f_one, f_one-zrefb(kp)))
3736 zrefd(kp) = max(f_zero, min(f_one, za2 / (f_one+za2) ))
3737 ztrad(kp) = max(f_zero, min(f_one, f_one - zrefd(kp) ))
3740 za1 = zgam1*zgam4 + zgam2*zgam3
3741 za2 = zgam1*zgam3 + zgam2*zgam4
3742 zrk = (zgam1 - zgam2) * (zgam1 + zgam2)
3743 if (zrk > eps1)
then
3753 zrpp1= f_one - zrp*zrp
3754 zrpp = sign( max(flimit, abs(zrpp1)), zrpp1 )
3759 zr1 = zrm1 * (za2 + zrkg3)
3760 zr2 = zrp1 * (za2 - zrkg3)
3761 zr3 = zrk2 * (zgam3 - za2*cosz)
3763 zr5 = zrpp * (zrk - zgam1)
3765 zt1 = zrp1 * (za1 + zrkg4)
3766 zt2 = zrm1 * (za1 - zrkg4)
3767 zt3 = zrk2 * (zgam4 + za1*cosz)
3772 zb1 = min( zrk*ztau1, 500.0 )
3773 if ( zb1 <= od_lo )
then
3774 zexm1 = f_one - zb1 + 0.5*zb1*zb1
3776 ftind = zb1 / (bpade + zb1)
3777 itind = ftind*ntbmx + 0.5
3778 zexm1 = exp_tbl(itind)
3780 zexp1 = f_one / zexm1
3782 zb2 = min( ztau1*sntz, 500.0 )
3783 if ( zb2 <= od_lo )
then
3784 zexm2 = f_one - zb2 + 0.5*zb2*zb2
3786 ftind = zb2 / (bpade + zb2)
3787 itind = ftind*ntbmx + 0.5
3788 zexm2 = exp_tbl(itind)
3790 zexp2 = f_one / zexm2
3791 ze1r45 = zr4*zexp1 + zr5*zexm1
3795 if ( abs(ze1r45) <= eps1 )
then
3799 zden1 = zssa1 / ze1r45
3800 zrefb(kp) = max(f_zero, min(f_one, &
3801 & (zr1*zexp1-zr2*zexm1-zr3*zexm2)*zden1 ))
3802 ztrab(kp) = max(f_zero, min(f_one, zexm2*(f_one - &
3803 & (zt1*zexp1-zt2*zexm1-zt3*zexp2)*zden1) ))
3807 if (ze1r45 >= f_zero)
then
3808 zden1 = zr4 / max(eps1, ze1r45*zrkg1)
3810 zden1 = zr4 / min(-eps1, ze1r45*zrkg1)
3812 zrefd(kp) = max(f_zero, min(f_one, &
3813 & zgam2*(zexp1 - zexm1)*zden1 ))
3814 ztrad(kp) = max(f_zero, min(f_one, zrk2*zden1 ))
3822 if ( zr1 <= od_lo )
then
3823 zexp3 = f_one - zr1 + 0.5*zr1*zr1
3825 ftind = zr1 / (bpade + zr1)
3826 itind = max(0, min(ntbmx, int(0.5+ntbmx*ftind) ))
3827 zexp3 = exp_tbl(itind)
3831 ztdbt(k) = zexp3 * ztdbt(kp)
3837 if ( zr1 <= od_lo )
then
3838 zexp4 = f_one - zr1 + 0.5*zr1*zr1
3840 ftind = zr1 / (bpade + zr1)
3841 itind = max(0, min(ntbmx, int(0.5+ntbmx*ftind) ))
3842 zexp4 = exp_tbl(itind)
3845 ztdbt0 = zexp4 * ztdbt0
3850 ztdbt(k) = zldbt(kp) * ztdbt(kp)
3853 ztdbt0 = zldbt0(k) * ztdbt0
3862 & ( zrefb,zrefd,ztrab,ztrad,zldbt,ztdbt, &
3870 fxupc(k,ib) = fxupc(k,ib) + zsolar*zfu(k)
3871 fxdnc(k,ib) = fxdnc(k,ib) + zsolar*zfd(k)
3877 zb2 = zsolar*(zfd(1) - ztdbt0)
3880 sfbmc(ibd) = sfbmc(ibd) + zb1
3881 sfdfc(ibd) = sfdfc(ibd) + zb2
3885 sfbmc(1) = sfbmc(1) + zf1
3886 sfdfc(1) = sfdfc(1) + zf2
3887 sfbmc(2) = sfbmc(2) + zf1
3888 sfdfc(2) = sfdfc(2) + zf2
3900 ftoadc = ftoadc + fxdn0(nlp1,ib)
3901 ftoau0 = ftoau0 + fxup0(nlp1,ib)
3902 fsfcu0 = fsfcu0 + fxup0(1,ib)
3903 fsfcd0 = fsfcd0 + fxdn0(1,ib)
3907 ibd = nuvb - nblow + 1
3908 suvbf0 = fxdn0(1,ibd)
3910 if ( cf1 <= eps )
then
3913 fxupc(k,ib) = fxup0(k,ib)
3914 fxdnc(k,ib) = fxdn0(k,ib)
3932 ftoauc = ftoauc + fxupc(nlp1,ib)
3933 fsfcuc = fsfcuc + fxupc(1,ib)
3934 fsfcdc = fsfcdc + fxdnc(1,ib)
3938 suvbfc = fxdnc(1,ibd)
subroutine, public rrtmg_sw_run(plyr, plvl, tlyr, tlvl, qlyr, olyr, gasvmr_co2, gasvmr_n2o, gasvmr_ch4, gasvmr_o2, gasvmr_co, gasvmr_cfc11, gasvmr_cfc12, gasvmr_cfc22, gasvmr_ccl4, icseed, aeraod, aerssa, aerasy, sfcalb_nir_dir, sfcalb_nir_dif, sfcalb_uvis_dir, sfcalb_uvis_dif, dzlyr, delpin, de_lgth, alpha, cosz, solcon, nday, idxday, npts, nlay, nlp1, lprnt, inc_minor_gas, iswcliq, iswcice, isubcsw, iovr, top_at_1, iswmode, cld_cf, lsswr, iovr_rand, iovr_maxrand, iovr_max, iovr_dcorr, iovr_exp, iovr_exprand, hswc, topflx, sfcflx, cldtau, hsw0, hswb, flxprf, fdncmp, cld_lwp, cld_ref_liq, cld_iwp, cld_ref_ice, cld_rwp, cld_ref_rain, cld_swp, cld_ref_snow, cld_od, cld_ssa, cld_asy, errmsg, errflg)