CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
sfc_ocean.F
1
4
7 module sfc_ocean
8 implicit none
9 private
10 public :: sfc_ocean_run
11
12 contains
13
14
23 subroutine sfc_ocean_run &
24!...................................
25! --- inputs:
26 & ( im, hvap, cp, rd, eps, epsm1, rvrdm1, ps, u1, v1, t1, q1, &
27 & tskin, cm, ch, lseaspray, fm, fm10, &
28 & usfco, vsfco, icplocn2atm, &
29 & prsl1, prslki, wet, use_lake_model, wind, &, ! --- inputs
30 & flag_iter, use_med_flux, dqsfc_med, dtsfc_med, &
31 & qsurf, cmm, chh, gflux, evap, hflx, ep, & ! --- outputs
32 & errmsg, errflg &
33 & )
34
35! ===================================================================== !
36! description: !
37! !
38! usage: !
39! !
40! call sfc_ocean !
41! inputs: !
42! ( im, ps, u1, v1, t1, q1, tskin, cm, ch, lseaspray, fm, fm10, !
43! usfco, vsfco, icplocn2atm, !
44! prsl1, prslki, wet, use_lake_model, wind, flag_iter, !
45! use_med_flux, !
46! outputs: !
47! qsurf, cmm, chh, gflux, evap, hflx, ep ) !
48! !
49! !
50! subprograms/functions called: fpvs !
51! !
52! !
53! program history log: !
54! 2005 -- created from the original progtm to account for !
55! ocean only !
56! oct 2006 -- h. wei added cmm and chh to the output !
57! apr 2009 -- y.-t. hou modified to match the modified gbphys.f !
58! reformatted the code and added program documentation !
59! sep 2009 -- s. moorthi removed rcl and made pa as pressure unit !
60! and furthur reformatted the code !
61! dec 2021 -- u. turuncoglu added support for receiving fluxes !
62! from mediator !
63! !
64! !
65! ==================== defination of variables ==================== !
66! !
67! inputs: size !
68! im - integer, horizontal dimension 1 !
69! ps - real, surface pressure im !
70! u1, v1 - real, u/v component of surface layer wind (m/s) im !
71! usfco - real, u component of surface ocean current (m/s) im !
72! vsfco - real, v component of surface ocean current (m/s) im !
73! icplocn2atm - integer, =1 if usfco and vsfco are used in the 1 !
74! computation of air-sea fluxes !
75! t1 - real, surface layer mean temperature ( k ) im !
76! q1 - real, surface layer mean specific humidity im !
77! tskin - real, ground surface skin temperature ( k ) im !
78! cm - real, surface exchange coeff for momentum (m/s) im !
79! ch - real, surface exchange coeff heat & moisture(m/s) im !
80! lseaspray- logical, .t. for parameterization for sea spray 1 !
81! fm - real, a stability profile function for momentum im !
82! fm10 - real, a stability profile function for momentum im !
83! at 10m !
84! prsl1 - real, surface layer mean pressure im !
85! prslki - real, im !
86! wet - logical, =T if any ocean/lak, =F otherwise im !
87! wind - real, wind speed (m/s) im !
88! flag_iter- logical, im !
89! use_med_flux - logical, =T to use fluxes coming from med 1 !
90! dqsfc_med- real, latent heat flux im !
91! dtsfc_med- real, sensible heat flux im !
92! !
93! outputs: !
94! qsurf - real, specific humidity at sfc im !
95! cmm - real, im !
96! chh - real, im !
97! gflux - real, ground heat flux (zero for ocean) im !
98! evap - real, evaporation from latent heat flux im !
99! hflx - real, sensible heat flux im !
100! ep - real, potential evaporation im !
101! !
102! ===================================================================== !
103!
104 use machine , only : kind_phys
105 use funcphys, only : fpvs
106!
107 implicit none
108
109! --- constant parameters:
110 real (kind=kind_phys), parameter :: one = 1.0_kind_phys, &
111 & zero = 0.0_kind_phys, qmin = 1.0e-8_kind_phys
112! --- inputs:
113 integer, intent(in) :: im
114 real (kind=kind_phys), intent(in) :: hvap, cp, rd, &
115 & eps, epsm1, rvrdm1
116
117 real (kind=kind_phys), dimension(:), intent(in) :: ps, u1, v1, &
118 & t1, q1, tskin, cm, ch, fm, fm10, prsl1, prslki, wind, &
119 & usfco, vsfco
120
121! For sea spray effect
122 logical, intent(in) :: lseaspray
123!
124 logical, dimension(:), intent(in) :: flag_iter, wet
125 integer, dimension(:), intent(in) :: use_lake_model
126 integer, intent(in) :: icplocn2atm
127!
128 logical, intent(in) :: use_med_flux
129
130! To receive fluxes from mediator
131 real (kind=kind_phys), dimension(:), intent(in), optional :: &
132 & dqsfc_med, dtsfc_med
133
134! --- outputs:
135 real (kind=kind_phys), dimension(:), intent(inout) :: qsurf, &
136 & cmm, chh, gflux, evap, hflx, ep
137
138 character(len=*), intent(out) :: errmsg
139 integer, intent(out) :: errflg
140
141! --- locals:
142
143 real (kind=kind_phys) :: qss, rch, tem,
144 & elocp, cpinv, hvapi
145 real (kind=kind_phys), dimension(im) :: rho, q0
146 real (kind=kind_phys), dimension(im) :: windrel
147
148 integer :: i
149
150 logical :: flag(im)
151!
152! parameters for sea spray effect
153!
154 real (kind=kind_phys) :: f10m, u10m, v10m, ws10, ru10, qss1,
155 & bb1, hflxs, evaps, ptem
156!
157! real (kind=kind_phys), parameter :: alps=0.5, bets=0.5, gams=0.1,
158! real (kind=kind_phys), parameter :: alps=0.5, bets=0.5, gams=0.0,
159! real (kind=kind_phys), parameter :: alps=1.0, bets=1.0, gams=0.2,
160 real (kind=kind_phys), parameter :: alps=0.75,bets=0.75,gams=0.15,
161 & ws10cr=30., conlf=7.2e-9, consf=6.4e-8
162!
163!======================================================================================================
164!===> ... begin here
165!
166! -- ... initialize CCPP error handling variables
167 errmsg = ''
168 errflg = 0
169
170 cpinv = one/cp
171 hvapi = one/hvap
172 elocp = hvap/cp
173!
175 do i = 1, im
176 flag(i) = (wet(i) .and. flag_iter(i) .and. use_lake_model(i)/=1)
178! ps is in pascals, wind is wind speed,
179! rho is density, qss is sat. hum. at surface
180
181 if ( flag(i) ) then
182 if (use_med_flux) then
183 q0(i) = max( q1(i), qmin )
184 rho(i) = prsl1(i) / (rd*t1(i)*(one + rvrdm1*q0(i)))
185
186 if (icplocn2atm == 0) then
187 tem = ch(i) * wind(i)
188 cmm(i) = cm(i) * wind(i)
189 else if (icplocn2atm ==1) then
190 windrel(i)=sqrt( (u1(i)-usfco(i))**2+(v1(i)-vsfco(i))**2 )
191 tem = ch(i) * windrel(i)
192 cmm(i) = cm(i) * windrel(i)
193 endif
194 chh(i) = rho(i) * tem
195
196 hflx(i) = dtsfc_med(i)
197 evap(i) = dqsfc_med(i)
198
199 qsurf(i) = q1(i) + dqsfc_med(i) / (hvap*chh(i))
200 gflux(i) = zero
201 else
202 q0(i) = max( q1(i), qmin )
203 rho(i) = prsl1(i) / (rd*t1(i)*(one + rvrdm1*q0(i)))
204
205 qss = fpvs( tskin(i) )
206 qss = eps*qss / (ps(i) + epsm1*qss)
207
208! --- ... rcp = rho cp ch v
209
210 if (icplocn2atm == 0) then
211 rch = rho(i) * cp * ch(i) * wind(i)
212 tem = ch(i) * wind(i)
213 cmm(i) = cm(i) * wind(i)
214 else if (icplocn2atm ==1) then
215 windrel(i)=sqrt( (u1(i)-usfco(i))**2+(v1(i)-vsfco(i))**2 )
216 rch = rho(i) * cp * ch(i) * windrel(i)
217 tem = ch(i) * windrel(i)
218 cmm(i) = cm(i) * windrel(i)
219 endif
220 chh(i) = rho(i) * tem
221
223
224 hflx(i) = rch * (tskin(i) - t1(i) * prslki(i))
225
226 evap(i) = elocp * rch * (qss - q0(i))
227
228 qsurf(i) = qss
229 gflux(i) = zero
230 endif
231 endif
232 enddo
233!
235!
236 do i=1,im
237 if(lseaspray .and. flag(i)) then
238 f10m = fm10(i) / fm(i)
239 u10m = f10m * u1(i)
240 v10m = f10m * v1(i)
241 ws10 = sqrt(u10m*u10m + v10m*v10m)
242 ws10 = max(ws10,1.)
243 ws10 = min(ws10,ws10cr)
244 tem = .015 * ws10 * ws10
245 ru10 = 1. - .087 * log(10./tem)
246 qss1 = fpvs(t1(i))
247 qss1 = eps * qss1 / (prsl1(i) + epsm1 * qss1)
248 tem = rd * cp * t1(i) * t1(i)
249 tem = 1. + eps * hvap * hvap * qss1 / tem
250 bb1 = 1. / tem
251 evaps = conlf * (ws10**5.4) * ru10 * bb1
252 evaps = evaps * rho(i) * hvap * (qss1 - q0(i))
253 evap(i) = evap(i) + alps * evaps
254 hflxs = consf * (ws10**3.4) * ru10
255 hflxs = hflxs * rho(i) * cp * (tskin(i) - t1(i))
256 ptem = alps - gams
257 hflx(i) = hflx(i) + bets * hflxs - ptem * evaps
258 endif
259 enddo
260!
261 do i = 1, im
262 if ( flag(i) ) then
263 tem = 1.0 / rho(i)
264 hflx(i) = hflx(i) * tem * cpinv
265 evap(i) = evap(i) * tem * hvapi
266 ep(i) = evap(i)
267 endif
268 enddo
269!
270 return
271!...................................
272 end subroutine sfc_ocean_run
273!-----------------------------------
275 end module sfc_ocean
subroutine, public sfc_ocean_run(im, hvap, cp, rd, eps, epsm1, rvrdm1, ps, u1, v1, t1, q1, tskin, cm, ch, lseaspray, fm, fm10, usfco, vsfco, icplocn2atm, prsl1, prslki, wet, use_lake_model, wind,,
Definition sfc_ocean.F:30
This module contains the CCPP-compliant GFS near-surface sea temperature scheme when the model is ini...
Definition sfc_ocean.F:7