CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
rrtmg_sw_post.F90
1
2
5 contains
6
14 subroutine rrtmg_sw_post_run (im, levr, levs, ltp, nday, lm, kd, lsswr, &
15 swhtr, sfcalb1, sfcalb2, sfcalb3, sfcalb4, htswc, htsw0, &
16 nirbmdi, nirdfdi, visbmdi, visdfdi, nirbmui, nirdfui, visbmui,&
17 visdfui, sfcdsw, sfcnsw, htrsw, swhc, scmpsw, sfcfsw, topfsw, &
18 errmsg, errflg)
19
20 use machine, only: kind_phys
23
24 implicit none
25
26 integer, intent(in) :: im, levr, levs, &
27 ltp, nday, lm, kd
28 logical, intent(in) :: lsswr, swhtr
29 real(kind=kind_phys), dimension(:), intent(in) :: sfcalb1, sfcalb2, &
30 sfcalb3, sfcalb4
31 real(kind=kind_phys), dimension(:,:), intent(in) :: htswc, htsw0
32
33 real(kind=kind_phys), dimension(:), intent(inout) :: nirbmdi, nirdfdi, &
34 visbmdi, visdfdi, &
35 nirbmui, nirdfui, &
36 visbmui, visdfui, &
37 sfcdsw, sfcnsw
38 real(kind=kind_phys), dimension(:,:), intent(inout) :: htrsw, swhc
39
40 type(cmpfsw_type), dimension(:), intent(inout) :: scmpsw
41 type(sfcfsw_type), dimension(:), intent(inout) :: sfcfsw
42 type(topfsw_type), dimension(:), intent(inout) :: topfsw
43
44 character(len=*), intent(out) :: errmsg
45 integer, intent(out) :: errflg
46 ! Local variables
47 integer :: i, k1, k
48
49 ! Initialize CCPP error handling variables
50 errmsg = ''
51 errflg = 0
52
53 if (lsswr) then
54 if (nday > 0) then
55 do k = 1, lm
56 k1 = k + kd
57 htrsw(1:im,k) = htswc(1:im,k1)
58 enddo
59 ! We are assuming that radiative tendencies are from bottom to top
60 ! --- repopulate the points above levr i.e. LM
61 if (lm < levs) then
62 do k = lm+1, levs
63 htrsw(1:im,k) = htrsw(1:im,lm)
64 enddo
65 endif
66
67 if (swhtr) then
68 do k = 1, lm
69 k1 = k + kd
70 swhc(1:im,k) = htsw0(1:im,k1)
71 enddo
72 ! --- repopulate the points above levr i.e. LM
73 if (lm < levs) then
74 do k = lm+1, levs
75 swhc(1:im,k) = swhc(1:im,lm)
76 enddo
77 endif
78 endif
79
80! --- surface down and up spectral component fluxes
81! Save two spectral bands' surface downward and upward fluxes for
82! output.
83
84 do i=1,im
85 nirbmdi(i) = scmpsw(i)%nirbm
86 nirdfdi(i) = scmpsw(i)%nirdf
87 visbmdi(i) = scmpsw(i)%visbm
88 visdfdi(i) = scmpsw(i)%visdf
89
90 nirbmui(i) = scmpsw(i)%nirbm * sfcalb1(i)
91 nirdfui(i) = scmpsw(i)%nirdf * sfcalb2(i)
92 visbmui(i) = scmpsw(i)%visbm * sfcalb3(i)
93 visdfui(i) = scmpsw(i)%visdf * sfcalb4(i)
94 enddo
95
96 else ! if_nday_block
97
98 htrsw(:,:) = 0.0
99
100 sfcfsw = sfcfsw_type( 0.0, 0.0, 0.0, 0.0 )
101 topfsw = topfsw_type( 0.0, 0.0, 0.0 )
102 scmpsw = cmpfsw_type( 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 )
103
104 do i=1,im
105 nirbmdi(i) = 0.0
106 nirdfdi(i) = 0.0
107 visbmdi(i) = 0.0
108 visdfdi(i) = 0.0
109
110 nirbmui(i) = 0.0
111 nirdfui(i) = 0.0
112 visbmui(i) = 0.0
113 visdfui(i) = 0.0
114 enddo
115
116 if (swhtr) then
117 swhc(:,:) = 0
118 endif
119
120 endif ! end_if_nday
121
122! --- radiation fluxes for other physics processes
123 do i=1,im
124 sfcnsw(i) = sfcfsw(i)%dnfxc - sfcfsw(i)%upfxc
125 sfcdsw(i) = sfcfsw(i)%dnfxc
126 enddo
127
128 endif ! end_if_lsswr
129
130 end subroutine rrtmg_sw_post_run
132 end module rrtmg_sw_post
subroutine rrtmg_sw_post_run(im, levr, levs, ltp, nday, lm, kd, lsswr, swhtr, sfcalb1, sfcalb2, sfcalb3, sfcalb4, htswc, htsw0, nirbmdi, nirdfdi, visbmdi, visdfdi, nirbmui, nirdfui, visbmui, visdfui, sfcdsw, sfcnsw, htrsw, swhc, scmpsw, sfcfsw, topfsw, errmsg, errflg)
This module is for specifying the band structures and program parameters used by the RRTMG-SW scheme.
Definition radsw_param.f:62
This module contains RRTMG-SW scheme post.
derived type for special components of surface SW fluxes
derived type for SW fluxes at surface
Definition radsw_param.f:89
derived type for SW fluxes at TOA
Definition radsw_param.f:79