CCPP SciDoc v7.0.x  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
sfc_land.F90
1
3
8
10 module sfc_land
11
12 use machine, only : kind_phys
13 use funcphys, only : fpvs
14
15 contains
16
21
22!!
23!! \section general General Algorithm
24!! \section detailed Detailed Algorithm
25!! @{
26 subroutine sfc_land_run(im, flag_init, flag_restart, &
27 cpllnd, cpllnd2atm, flag_iter, dry, &
28 t1, q1, prsl1, prslki, ps, tskin, wind, cm, ch, &
29 dlwflx, dswsfc, sfalb, sfcemis, &
30 rd, eps, epsm1, rvrdm1, hvap, cp, con_sbc, &
31 sncovr1_lnd, qsurf_lnd, &
32 evap_lnd, hflx_lnd, ep_lnd, t2mmp_lnd, q2mp_lnd, gflux_lnd, &
33 runoff_lnd, drain_lnd, cmm_lnd, chh_lnd, zvfun_lnd, slc, &
34 sncovr1, qsurf, evap, hflx, ep, t2mmp, q2mp, &
35 gflux, runoff, drain, cmm, chh, zvfun, &
36 errmsg, errflg)
37
38 implicit none
39
40 ! Inputs
41 integer , intent(in) :: im
42 logical , intent(in) :: flag_init
43 logical , intent(in) :: flag_restart
44 logical , intent(in) :: cpllnd
45 logical , intent(in) :: cpllnd2atm
46 logical , intent(in) :: flag_iter(:)
47 logical , intent(in) :: dry(:)
48 real(kind=kind_phys), intent(in) :: t1(:)
49 real(kind=kind_phys), intent(in) :: q1(:)
50 real(kind=kind_phys), intent(in) :: prsl1(:)
51 real(kind=kind_phys), intent(in) :: prslki(:)
52 real(kind=kind_phys), intent(in) :: ps(:)
53 real(kind=kind_phys), intent(in) :: tskin(:)
54 real(kind=kind_phys), intent(in) :: wind(:)
55 real(kind=kind_phys), intent(in) :: cm(:)
56 real(kind=kind_phys), intent(in) :: ch(:)
57 real(kind=kind_phys), intent(in) :: dlwflx(:)
58 real(kind=kind_phys), intent(in) :: dswsfc(:)
59 real(kind=kind_phys), intent(in) :: sfalb(:)
60 real(kind=kind_phys), intent(in) :: sfcemis(:)
61 real(kind=kind_phys), intent(in) :: rd
62 real(kind=kind_phys), intent(in) :: eps
63 real(kind=kind_phys), intent(in) :: epsm1
64 real(kind=kind_phys), intent(in) :: rvrdm1
65 real(kind=kind_phys), intent(in) :: hvap
66 real(kind=kind_phys), intent(in) :: cp
67 real(kind=kind_phys), intent(in) :: con_sbc
68 real(kind=kind_phys), intent(in), optional :: sncovr1_lnd(:)
69 real(kind=kind_phys), intent(in), optional :: qsurf_lnd(:)
70 real(kind=kind_phys), intent(in), optional :: evap_lnd(:)
71 real(kind=kind_phys), intent(in), optional :: hflx_lnd(:)
72 real(kind=kind_phys), intent(in), optional :: ep_lnd(:)
73 real(kind=kind_phys), intent(in), optional :: t2mmp_lnd(:)
74 real(kind=kind_phys), intent(in), optional :: q2mp_lnd(:)
75 real(kind=kind_phys), intent(in), optional :: gflux_lnd(:)
76 real(kind=kind_phys), intent(in), optional :: runoff_lnd(:)
77 real(kind=kind_phys), intent(in), optional :: drain_lnd(:)
78 real(kind=kind_phys), intent(in), optional :: cmm_lnd(:)
79 real(kind=kind_phys), intent(in), optional :: chh_lnd(:)
80 real(kind=kind_phys), intent(in), optional :: zvfun_lnd(:)
81 real(kind=kind_phys), intent(in), optional :: slc(:,:)
82 ! Inputs/Outputs
83 real(kind=kind_phys), intent(inout) :: sncovr1(:)
84 real(kind=kind_phys), intent(inout) :: qsurf(:)
85 real(kind=kind_phys), intent(inout) :: evap(:)
86 real(kind=kind_phys), intent(inout) :: hflx(:)
87 real(kind=kind_phys), intent(inout) :: ep(:)
88 real(kind=kind_phys), intent(inout), optional :: t2mmp(:)
89 real(kind=kind_phys), intent(inout), optional :: q2mp(:)
90 real(kind=kind_phys), intent(inout) :: gflux(:)
91 real(kind=kind_phys), intent(inout) :: runoff(:)
92 real(kind=kind_phys), intent(inout) :: drain(:)
93 real(kind=kind_phys), intent(inout) :: cmm(:)
94 real(kind=kind_phys), intent(inout) :: chh(:)
95 real(kind=kind_phys), intent(inout) :: zvfun(:)
96 ! Outputs
97 character(len=*) , intent(out) :: errmsg
98 integer , intent(out) :: errflg
99
100 ! Constant parameters
101 real(kind=kind_phys), parameter :: &
102 & one = 1.0_kind_phys, &
103 & zero = 0.0_kind_phys, &
104 & qmin = 1.0e-8_kind_phys, &
105 & slc_min = 0.05_kind_phys, & ! estimate dry limit for soil moisture
106 & slc_max = 0.50_kind_phys ! estimate saturated limit for soil moisture
107
108 ! Locals
109 integer :: i
110 real(kind=kind_phys) :: qss, rch, tem, cpinv, hvapi, elocp
111 real(kind=kind_phys) :: available_energy, soil_stress_factor
112 real(kind=kind_phys), dimension(im) :: rho, q0
113
114 ! Initialize CCPP error handling variables
115 errmsg = ''
116 errflg = 0
117
118 cpinv = one/cp
119 hvapi = one/hvap
120 elocp = hvap/cp
121
122 ! Check coupling from component land to atmosphere
123 if (.not. cpllnd2atm) return
124
125 ! Check if it is cold or warm run
126 if (flag_init .and. .not. flag_restart) then
127 ! Calculate fluxes internally
128 do i = 1, im
129 if (dry(i)) then
130 soil_stress_factor = (slc(i,1)-slc_min)/(slc_max-slc_min)
131 soil_stress_factor = min(max(zero,soil_stress_factor),one)
132 available_energy = dswsfc(i)*(one-sfalb(i))+dlwflx(i)*sfcemis(i) - &
133 sfcemis(i)*con_sbc*tskin(i)**4
134 available_energy = min(max(-200.0,available_energy),1000.0) ! set some arbitrary limits
135 q0(i) = max(q1(i), qmin)
136 rho(i) = prsl1(i)/(rd*t1(i)*(one+rvrdm1*q0(i)))
137 qss = fpvs(tskin(i))
138 qss = eps*qss/(ps(i)+epsm1*qss)
139 rch = rho(i)*cp*ch(i)*wind(i)
140 tem = ch(i)*wind(i)
141 sncovr1(i) = zero
142 qsurf(i) = qss
143 hflx(i) = rch*(tskin(i)-t1(i)*prslki(i)) ! first guess hflx [W/m2]
144 evap(i) = elocp*rch*(qss-q0(i)) ! first guess evap [W/m2]
145 evap(i) = evap(i)*soil_stress_factor ! reduce evap for soil moisture stress
146 hflx(i) = min(max(-100.0,hflx(i)),500.0) ! set some arbitrary limits
147 evap(i) = min(max(-100.0,evap(i)),500.0) ! set some arbitrary limits
148 if(evap(i) + hflx(i) /= zero) then
149 hflx(i) = available_energy * hflx(i) / (abs(evap(i)) + abs(hflx(i)))
150 evap(i) = available_energy * evap(i) / (abs(evap(i)) + abs(hflx(i)))
151 else
152 hflx(i) = zero
153 evap(i) = zero
154 end if
155 hflx(i) = min(max(-100.0,hflx(i)),500.0) ! set some arbitrary limits
156 evap(i) = min(max(-100.0,evap(i)),500.0) ! set some arbitrary limits
157 hflx(i) = hflx(i)*(1.0/rho(i))*cpinv ! convert to expected units
158 ep(i) = evap(i)
159 evap(i) = evap(i)*(1.0/rho(i))*hvapi ! convert to expected units
160 t2mmp(i) = tskin(i)
161 q2mp(i) = qsurf(i)
162 gflux(i) = zero
163 drain(i) = zero
164 runoff(i) = zero
165 cmm(i) = cm(i)*wind(i)
166 chh(i) = rho(i)*tem
167 zvfun(i) = one
168 end if
169 enddo
170 else
171 ! Use fluxes from land component model
172 do i = 1, im
173 if (dry(i)) then
174 sncovr1(i) = sncovr1_lnd(i)
175 qsurf(i) = qsurf_lnd(i)
176 hflx(i) = hflx_lnd(i)
177 evap(i) = evap_lnd(i)
178 ep(i) = ep_lnd(i)
179 t2mmp(i) = t2mmp_lnd(i)
180 q2mp(i) = q2mp_lnd(i)
181 gflux(i) = gflux_lnd(i)
182 drain(i) = drain_lnd(i)
183 runoff(i) = runoff_lnd(i)
184 cmm(i) = cmm_lnd(i)
185 chh(i) = chh_lnd(i)
186 zvfun(i) = zvfun_lnd(i)
187 end if
188 enddo
189 endif
190
191 end subroutine sfc_land_run
192
193 end module sfc_land
This module contains the CCPP-compliant GFS land post interstitial codes, which returns updated surfa...
Definition sfc_land.F90:10