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