CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
radiation_tools.F90
1
3
6 use machine, only: &
7 kind_phys ! Working type
8 implicit none
9
10 real(kind_phys) :: &
11 rrtmgp_minp, & ! Minimum pressure allowed in RRTMGP
12 rrtmgp_mint ! Minimum temperature allowed in RRTMGP
13contains
14
16 subroutine cmp_tlev(nCol,nLev,minP,p_lay,t_lay,p_lev,tsfc,t_lev)
17 ! Inputs
18 integer, intent(in) :: &
19 nCol,nLev
20 real(kind_phys),intent(in) :: &
21 minP
22 real(kind_phys),dimension(nCol),intent(in) :: &
23 tsfc
24 real(kind_phys),dimension(nCol,nLev),intent(in) :: &
25 p_lay,t_lay
26 real(kind_phys),dimension(nCol,nLev+1),intent(in) :: &
27 p_lev
28
29 ! Outputs
30 real(kind_phys),dimension(nCol,nLev+1),intent(out) :: &
31 t_lev
32
33 ! Local
34 integer :: iCol,iLay, iSFC, iTOA
35 logical :: top_at_1
36 real(kind_phys), dimension(nCol,nLev) :: tem2da, tem2db
37
38 top_at_1 = (p_lev(1,1) .lt. p_lev(1, nlev))
39 if (top_at_1) then
40 isfc = nlev
41 itoa = 1
42 else
43 isfc = 1
44 itoa = nlev
45 endif
46
47 if (itoa .eq. 1) then
48 tem2da(1:ncol,2:isfc) = log(p_lay(1:ncol,2:isfc))
49 tem2db(1:ncol,2:isfc) = log(p_lev(1:ncol,2:isfc))
50 do icol = 1, ncol
51 tem2da(icol,1) = log(p_lay(icol,1) )
52 tem2db(icol,1) = log(max(minp, p_lev(icol,1)) )
53 tem2db(icol,isfc) = log(p_lev(icol,isfc) )
54 enddo
55 t_lev(1:ncol,1) = t_lay(1:ncol,itoa)
56 do ilay = 2, isfc
57 do icol = 1, ncol
58 t_lev(icol,ilay) = t_lay(icol,ilay) + (t_lay(icol,ilay-1) - t_lay(icol,ilay))&
59 * (tem2db(icol,ilay) - tem2da(icol,ilay)) &
60 / (tem2da(icol,ilay-1) - tem2da(icol,ilay))
61 enddo
62 enddo
63 t_lev(1:ncol,isfc+1) = tsfc(1:ncol)
64 else
65 tem2da(1:ncol,2:itoa) = log(p_lay(1:ncol,2:itoa))
66 tem2db(1:ncol,2:itoa) = log(p_lev(1:ncol,2:itoa))
67 do icol = 1, ncol
68 tem2da(icol,1) = log(p_lay(icol,1))
69 tem2db(icol,1) = log(p_lev(icol,1))
70 tem2db(icol,itoa) = log(max(minp, p_lev(icol,itoa)) )
71 enddo
72
73 t_lev(1:ncol,1) = tsfc(1:ncol)
74 do ilay = 1, itoa-1
75 do icol = 1, ncol
76 t_lev(icol,ilay+1) = t_lay(icol,ilay) + (t_lay(icol,ilay+1) - t_lay(icol,ilay))&
77 * (tem2db(icol,ilay+1) - tem2da(icol,ilay)) &
78 / (tem2da(icol,ilay+1) - tem2da(icol,ilay))
79 enddo
80 enddo
81 t_lev(1:ncol,itoa+1) = t_lay(1:ncol,itoa)
82 endif
83
84 end subroutine cmp_tlev
85
87 subroutine check_error_msg(routine_name, error_msg)
88 character(len=*), intent(in) :: &
89 error_msg, routine_name
90
91 if(error_msg /= "") then
92 print*,"ERROR("//trim(routine_name)//"): "
93 print*,trim(error_msg)
94 return
95 end if
96 end subroutine check_error_msg
97
98end module radiation_tools
This module contains tools for radiation.