CCPP SciDoc v7.0.x  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 sfcdswc, 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 sfcdswc
39 real(kind=kind_phys), dimension(:,:), intent(inout) :: htrsw, swhc
40
41 type(cmpfsw_type), dimension(:), intent(inout) :: scmpsw
42 type(sfcfsw_type), dimension(:), intent(inout) :: sfcfsw
43 type(topfsw_type), dimension(:), intent(inout) :: topfsw
44
45 character(len=*), intent(out) :: errmsg
46 integer, intent(out) :: errflg
47 ! Local variables
48 integer :: i, k1, k
49
50 ! Initialize CCPP error handling variables
51 errmsg = ''
52 errflg = 0
53
54 if (lsswr) then
55 if (nday > 0) then
56 do k = 1, lm
57 k1 = k + kd
58 htrsw(1:im,k) = htswc(1:im,k1)
59 enddo
60 ! We are assuming that radiative tendencies are from bottom to top
61 ! --- repopulate the points above levr i.e. LM
62 if (lm < levs) then
63 do k = lm+1, levs
64 htrsw(1:im,k) = htrsw(1:im,lm)
65 enddo
66 endif
67
68 if (swhtr) then
69 do k = 1, lm
70 k1 = k + kd
71 swhc(1:im,k) = htsw0(1:im,k1)
72 enddo
73 ! --- repopulate the points above levr i.e. LM
74 if (lm < levs) then
75 do k = lm+1, levs
76 swhc(1:im,k) = swhc(1:im,lm)
77 enddo
78 endif
79 endif
80
81! --- surface down and up spectral component fluxes
82! Save two spectral bands' surface downward and upward fluxes for
83! output.
84
85 do i=1,im
86 nirbmdi(i) = scmpsw(i)%nirbm
87 nirdfdi(i) = scmpsw(i)%nirdf
88 visbmdi(i) = scmpsw(i)%visbm
89 visdfdi(i) = scmpsw(i)%visdf
90
91 nirbmui(i) = scmpsw(i)%nirbm * sfcalb1(i)
92 nirdfui(i) = scmpsw(i)%nirdf * sfcalb2(i)
93 visbmui(i) = scmpsw(i)%visbm * sfcalb3(i)
94 visdfui(i) = scmpsw(i)%visdf * sfcalb4(i)
95 enddo
96
97 else ! if_nday_block
98
99 htrsw(:,:) = 0.0
100
101 sfcfsw = sfcfsw_type( 0.0, 0.0, 0.0, 0.0 )
102 topfsw = topfsw_type( 0.0, 0.0, 0.0 )
103 scmpsw = cmpfsw_type( 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 )
104
105 do i=1,im
106 nirbmdi(i) = 0.0
107 nirdfdi(i) = 0.0
108 visbmdi(i) = 0.0
109 visdfdi(i) = 0.0
110
111 nirbmui(i) = 0.0
112 nirdfui(i) = 0.0
113 visbmui(i) = 0.0
114 visdfui(i) = 0.0
115 enddo
116
117 if (swhtr) then
118 swhc(:,:) = 0
119 endif
120
121 endif ! end_if_nday
122
123! --- radiation fluxes for other physics processes
124 do i=1,im
125 sfcnsw(i) = sfcfsw(i)%dnfxc - sfcfsw(i)%upfxc
126 sfcdsw(i) = sfcfsw(i)%dnfxc
127 sfcdswc(i)= sfcfsw(i)%dnfx0
128 enddo
129
130 endif ! end_if_lsswr
131
132 end subroutine rrtmg_sw_post_run
134 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, sfcdswc, 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