CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
rrtmg_lw_post.F90
1
2
5 contains
6
10
14 subroutine rrtmg_lw_post_run (im, levs, ltp, lm, kd, lslwr, lwhtr, &
15 tsfa, htlwc, htlw0, sfcflw, tsflw, sfcdlw, htrlw, lwhc, &
16 errmsg, errflg)
17
18 use machine, only: kind_phys
20
21 implicit none
22
23 integer, intent(in) :: im, levs, ltp, lm, kd
24 logical, intent(in) :: lslwr, lwhtr
25 real(kind=kind_phys), dimension(im), intent(in) :: tsfa
26 real(kind=kind_phys), dimension(im, LM+LTP), intent(in) :: htlwc
27 real(kind=kind_phys), dimension(im, LM+LTP), intent(in) :: htlw0
28
29 type(sfcflw_type), dimension(im), intent(in) :: sfcflw
30
31 real(kind=kind_phys), dimension(im), intent(inout) :: tsflw, sfcdlw
32 real(kind=kind_phys), dimension(im, levs), intent(inout) :: htrlw, lwhc
33 character(len=*), intent(out) :: errmsg
34 integer, intent(out) :: errflg
35
36 ! local variables
37 integer :: k1, k
38
39 ! Initialize CCPP error handling variables
40 errmsg = ''
41 errflg = 0
42
43 if (lslwr) then
44! Save calculation results
45! Save surface air temp for diurnal adjustment at model t-steps
46
47 tsflw(:) = tsfa(:)
48
49 do k = 1, lm
50 k1 = k + kd
51 htrlw(1:im,k) = htlwc(1:im,k1)
52 enddo
53 ! --- repopulate the points above levr
54 if (lm < levs) then
55 do k = lm+1, levs
56 htrlw(1:im,k) = htrlw(1:im,lm)
57 enddo
58 endif
59
60 if (lwhtr) then
61 do k = 1, lm
62 k1 = k + kd
63 lwhc(1:im,k) = htlw0(1:im,k1)
64 enddo
65 ! --- repopulate the points above levr
66 if (lm < levs) then
67 do k = lm+1, levs
68 lwhc(1:im,k) = lwhc(1:im,lm)
69 enddo
70 endif
71 endif
72
73! --- radiation fluxes for other physics processes
74 sfcdlw(:) = sfcflw(:)%dnfxc
75
76 endif ! end_if_lslwr
77
78 end subroutine rrtmg_lw_post_run
79
81 end module rrtmg_lw_post
subroutine rrtmg_lw_post_run(im, levs, ltp, lm, kd, lslwr, lwhtr, tsfa, htlwc, htlw0, sfcflw, tsflw, sfcdlw, htrlw, lwhc, errmsg, errflg)
This module contains LW band parameters set up.
Definition radlw_param.f:61
This module contains code executed after RRTMG-LW scheme.
derived type for LW fluxes at surface
Definition radlw_param.f:87