Line data Source code
1 :
2 : module radheat
3 : !-----------------------------------------------------------------------
4 : !
5 : ! Purpose: Provide an interface to convert shortwave and longwave
6 : ! radiative heating terms into net heating.
7 : !
8 : ! This module provides a hook to allow incorporating additional
9 : ! radiative terms (eUV heating and nonLTE longwave cooling).
10 : !
11 : ! Original version: B.A. Boville
12 : !-----------------------------------------------------------------------
13 :
14 : use shr_kind_mod, only: r8 => shr_kind_r8
15 : use ppgrid, only: pcols, pver
16 : use physics_types, only: physics_state, physics_ptend, physics_ptend_init
17 :
18 : use physics_buffer, only : physics_buffer_desc
19 :
20 : implicit none
21 : private
22 : save
23 :
24 : ! Public interfaces
25 : public &
26 : radheat_readnl, &!
27 : radheat_register, &!
28 : radheat_init, &!
29 : radheat_timestep_init, &!
30 : radheat_tend ! return net radiative heating
31 :
32 : public :: radheat_disable_waccm ! disable waccm heating in the upper atm
33 :
34 : !===============================================================================
35 : contains
36 : !===============================================================================
37 :
38 1536 : subroutine radheat_readnl(nlfile)
39 :
40 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
41 :
42 : ! No options for this version of radheat; this is just a stub.
43 :
44 1536 : end subroutine radheat_readnl
45 :
46 : !================================================================================================
47 :
48 1536 : subroutine radheat_register
49 :
50 : ! No options for this version of radheat; this is just a stub.
51 :
52 1536 : end subroutine radheat_register
53 :
54 : !================================================================================================
55 :
56 1536 : subroutine radheat_init(pref_mid)
57 :
58 : use pmgrid, only: plev
59 : use physics_buffer, only : physics_buffer_desc
60 :
61 : real(r8), intent(in) :: pref_mid(plev)
62 :
63 1536 : end subroutine radheat_init
64 :
65 : !================================================================================================
66 :
67 370944 : subroutine radheat_timestep_init (state, pbuf2d)
68 1536 : use physics_types,only : physics_state
69 : use ppgrid, only : begchunk, endchunk
70 : use physics_buffer, only : physics_buffer_desc
71 :
72 : type(physics_state), intent(in):: state(begchunk:endchunk)
73 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
74 :
75 :
76 370944 : end subroutine radheat_timestep_init
77 :
78 : !================================================================================================
79 :
80 62675424 : subroutine radheat_tend(state, pbuf, ptend, qrl, qrs, fsns, &
81 : fsnt, flns, flnt, asdir, net_flx)
82 : #if ( defined OFFLINE_DYN )
83 : use metdata, only: met_rlx, met_srf_feedback
84 : #endif
85 : !-----------------------------------------------------------------------
86 : ! Compute net radiative heating from qrs and qrl, and the associated net
87 : ! boundary flux.
88 : !-----------------------------------------------------------------------
89 :
90 : ! Arguments
91 : type(physics_state), intent(in) :: state ! Physics state variables
92 :
93 : type(physics_buffer_desc), pointer :: pbuf(:)
94 : type(physics_ptend), intent(out) :: ptend ! indivdual parameterization tendencie
95 : real(r8), intent(in) :: qrl(pcols,pver) ! longwave heating
96 : real(r8), intent(in) :: qrs(pcols,pver) ! shortwave heating
97 : real(r8), intent(in) :: fsns(pcols) ! Surface solar absorbed flux
98 : real(r8), intent(in) :: fsnt(pcols) ! Net column abs solar flux at model top
99 : real(r8), intent(in) :: flns(pcols) ! Srf longwave cooling (up-down) flux
100 : real(r8), intent(in) :: flnt(pcols) ! Net outgoing lw flux at model top
101 : real(r8), intent(in) :: asdir(pcols) ! shortwave, direct albedo
102 : real(r8), intent(out) :: net_flx(pcols)
103 :
104 :
105 : ! Local variables
106 : integer :: i, k
107 : integer :: ncol
108 : !-----------------------------------------------------------------------
109 :
110 1492272 : ncol = state%ncol
111 :
112 1492272 : call physics_ptend_init(ptend,state%psetcols, 'radheat', ls=.true.)
113 :
114 : #if ( defined OFFLINE_DYN )
115 : ptend%s(:ncol,:) = 0._r8
116 : do k = 1,pver
117 : if (met_rlx(k) < 1._r8 .or. met_srf_feedback) then
118 : ptend%s(:ncol,k) = (qrs(:ncol,k) + qrl(:ncol,k))
119 : endif
120 : enddo
121 : #else
122 2318817168 : ptend%s(:ncol,:) = (qrs(:ncol,:) + qrl(:ncol,:))
123 : #endif
124 :
125 24917472 : do i = 1, ncol
126 24917472 : net_flx(i) = fsnt(i) - fsns(i) - flnt(i) + flns(i)
127 : end do
128 :
129 370944 : end subroutine radheat_tend
130 :
131 : !================================================================================================
132 0 : subroutine radheat_disable_waccm()
133 0 : end subroutine radheat_disable_waccm
134 : end module radheat
|