forked from HiroyukiTsujino/JRA55-do
-
Notifications
You must be signed in to change notification settings - Fork 0
/
bulk_shift_zrough.F90
291 lines (251 loc) · 11.8 KB
/
bulk_shift_zrough.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
! -*-F90-*-
!------------------------ bulk-ncar-shift.F90 ----------------------------
! Information:
! Shift 2m temperature and specific humidity to those at 10m.
!
subroutine bulk_shift(tmptrg,sphtrg, &
& sat,qar,wdv,tau,slp,sst,zrough,sphsrf,&
& sens, evap, &
& aexl,ice,num_data,altu,altt,altq,alt_target)
implicit none
integer(4), intent(in) :: num_data
real(8), intent(out) :: tmptrg(num_data), sphtrg(num_data) ! temperature and specific humidity at target altitude alt_target
real(8), intent(in) :: sat(num_data) ! Observed Air Temperature at height altt [degC]
real(8), intent(in) :: qar(num_data) ! specific humidity at height altq
real(8), intent(in) :: wdv(num_data) ! Observed Wind Speed [cm/s] at height altu - not used, and call from shift_2m_to_10m_hybrid sets it to zero
real(8), intent(in) :: slp(num_data) ! Observed sea level pressure [hPa]
real(8), intent(in) :: tau(num_data) ! wind stress [Pa]
real(8), intent(in) :: sst(num_data) ! [degC]
real(8), intent(in) :: zrough(num_data) ! roughness length
! sensible heat and evaporation : negative for the "gain" on the atmospheric side
real(8), intent(in) :: sens(num_data), evap(num_data)
real(8), intent(in) :: ice(num_data) ! sea ice concentration (0 or 1, no fractional values)
real(8), intent(in) :: aexl(num_data) ! land mask (1.0=ocean)
real(8), intent(in) :: altu, altt, altq ! original altitudes of wind, temperature, specific humidity [m], defined in `namelist.shift_height_rough*`
real(8), intent(in) :: alt_target ! altitude to shift temperature & specific humidity to [m], defined in `namelist.shift_height_rough*`
real(8), intent(out) :: sphsrf(num_data) ! surface saturated specific humidity
real(8) :: aexl_tmp(num_data)
! [cgs]
real(8), parameter :: grav = 981.0d0 ! [cm/s2]
! ...
real(8), parameter :: tab = 273.15d0 ! 0 degC expressed in Kelvin
! [MKS]
real(8), parameter :: grav_mks = grav * 1.0d-2
real(8), parameter :: cpa_mks = 1004.67d0 ! specific heat of dry air at constant pressure [J/Kg/K]
real(8), parameter :: gasr = 287.04d0 ! gas constant of dry air [J/kg/K]
real(8), parameter :: sphmin = 5.0d-5 ! minimum of specific humidity
! parameters to calculate saturation water vapor presure
real(8), parameter :: ali1=0.7859d0, ali2=0.03477d0, ali3=0.00412d0, ali4=1.0d0
real(8), parameter :: bli1=1.d-6, bli2=4.5d0, bli3=0.0006d0
real(8), parameter :: cli1=2.839d6, cli2=3.6d0, cli3=35.0d0
real(8), parameter :: ewa = 0.62197d0 !Molecular Weight Ratio Water/DryAir
! work variables
real(8) :: es, qs
real(8) :: ewp1, ewp2
!
real(8) :: dqr, dtemp
!
real(8) :: qatmos, satmos, slpres, tssurf
real(8) :: rhoair
!
! defined at 10 [m]
!
real(8) :: cdn10, cen10, ctn10
real(8) :: cdn10_rt
real(8) :: wv10
!
! defined at velocity level
!
real(8) :: qatmosu, satmosu
real(8) :: cdn, cen, ctn
real(8) :: cd_rt
real(8) :: wv
!
real(8), parameter :: karman = 0.4 ! von Karman constant
real(8) :: stab
!
real(8) :: ustar, tstar, qstar, bstar
real(8) :: zetau, zetat, zetaq, zetatrg
real(8) :: psi_mu, psi_hu
real(8) :: psi_mt, psi_ht
real(8) :: psi_mq, psi_hq
real(8) :: psi_mtrg, psi_htrg
!
real(8) :: x, xx, x2, tv
integer(4) :: m
!
! 入力(単位に注意! 用意するデータの単位を変更するか、Input (Be careful of the unit! Change the unit of the prepared data
! このサブルーチンで変更すること)or change it with this subroutine.)
! TBL : Model SST [degC]
! SAT,SATA: Observed Air Temperature [degC]
! DWT,DWTA: Observed Specific Humidity [g/g]
! WDV,WDVA: Observed Wind Speed [cm/s]
! SLP,SLPA: Observed sea level pressure [hPa]
! 出力 output
! QSN : Sensible Heat Flux [g/s3] (W/m2での1000倍の値)(1000 times the value in W / m2)
! QLA : Latent Heat Flux [g/s3] (W/m2での1000倍の値)(1000 times the value in W / m2)
!
! ES : 飽和水蒸気圧 (hPa) (!! MKS !!) Saturated water vapor pressure
! QS : 海面における飽和比湿 (g/kg) Saturated specific humidity at sea level
! WV : 海上風速 (cm/s) offshore wind speed
! CDN,CEN: バルク係数(無次元)Bulk coefficient (dimensionless)
! DQR,DTEMP: 大気海洋間の比湿差、温度差 Specific humidity difference and temperature difference between atmosphere and ocean
!
!----------------------------------------------------------------------
aexl_tmp(1:num_data) = aexl(1:num_data)
! aexl_tmp(1:num_data) = 0.0d0
!$omp parallel
!$omp do private(slpres,qatmos,satmos, &
!$omp & tssurf, dtemp, es, qs, dqr, wv, tv, wv10, &
!$omp & cdn10, cdn10_rt, cen10, stab, ctn10, cdn, ctn, cen, &
!$omp & cd_rt, ustar, tstar, qstar, bstar, zetau, x2, x, xx, &
!$omp & psi_mu, psi_hu, zetat, psi_mt, psi_ht, &
!$omp & zetatrg, psi_mtrg, psi_htrg, &
!$omp & zetaq, psi_mq, psi_hq, satmosu, qatmosu, rhoair)
do m = 1, num_data
slpres = slp(m)
qatmos = qar(m)
satmos = sat(m) + tab
tssurf = sst(m)
! calculate surface saturated specific humidity
if (ice(m) == 0.0d0) then
es = 9.8d-1 * 6.1078d0 * 1.0d1**(7.5d0 * tssurf / (2.373d2 + tssurf))
qs = ewa * es / (slpres - 3.78d-1 * es)
else
ewp1 = 10.D0**((ali1 + ali2 * tssurf) / (1.0d0 + ali3 * tssurf)) ! [hPa]
ewp2 = ewp1 * (1.0d0 + bli1 * slpres * (bli2 + bli3 * tssurf**2)) ! [hPa]
qs = ewa * ewp2 / (slpres - (1.0d0 - ewa) * ewp2)
end if
sphsrf(m) = qs
if (aexl_tmp(m) == 1.0d0) then ! ocean
dtemp = tssurf + tab - satmos ! SST - SAT
dqr = qs - qatmos
tmptrg(m) = tssurf - dtemp * log(alt_target / zrough(m)) / log(altt / zrough(m))
sphtrg(m) = qs - dqr * log(alt_target / zrough(m)) / log(altq / zrough(m))
else ! land
tv = satmos * (1.0d0 + 0.6078d0 * qatmos) ! virtual temperature
rhoair = (slpres * 1.0d2) / gasr / tv
! L-Y is:
! Large, W. G. and Yeager, S. (2004).
! Diurnal to decadal global forcing for ocean and sea-ice models: The data sets and flux climatologies.
! Technical Note NCAR/TN-460+STR, NCAR.
! http://dx.doi.org/10.5065/D6KK98Q6
ustar = sqrt(tau(m)/rhoair) ! L-Y eqn. 7a
tstar = sens(m) / rhoair / cpa_mks / ustar ! L-Y eqn. 7b
qstar = evap(m) / rhoair / ustar ! L-Y eqn. 7c
! VVVVVVVVV the code below is unused VVVVVVVVV
bstar = grav_mks * &
& ( tstar / tv + qstar / (qatmos + 1.0d0 / 0.6078d0))
!
! velocity level
!
zetau = karman * bstar * altu / (ustar * ustar) ! L-Y eqn. 8a
zetau = sign(min(abs(zetau),10.d0),zetau) ! undocumented NCAR
x2 = sqrt(abs(1.0d0 - 16.0d0 * zetau)) ! L-Y eqn. 8b (squared)
x2 = max(x2,1.0d0) ! undocumented NCAR
x = sqrt(x2)
if (zetau > 0.0d0) then
psi_mu = - 5.0d0 * zetau ! L-Y eqn. 8c
psi_hu = - 5.0d0 * zetau ! L-Y eqn. 8c
else
psi_mu = log((1.0d0 + 2.0d0 * x + x2) &
& * (1.0d0 + x2) / 8.0d0) &
& - 2.0d0 * (atan(x) - atan(1.0d0)) ! L-Y eqn. 8d
psi_hu = 2.0d0 * log((1.0d0 + x2) / 2.0d0) ! L-Y eqn. 8e
end if
!
! temperature level
!
zetat = karman * bstar * altt / (ustar * ustar) ! L-Y eqn. 8a
zetat = sign(min(abs(zetat),10.d0),zetat) ! undocumented NCAR
x2 = sqrt(abs(1.0d0 - 16.0d0 * zetat)) ! L-Y eqn. 8b (squared)
x2 = max(x2,1.0d0) ! undocumented NCAR
x = sqrt(x2)
if (zetat > 0.0d0) then
psi_mt = - 5.0d0 * zetat ! L-Y eqn. 8c
psi_ht = - 5.0d0 * zetat ! L-Y eqn. 8c
else
psi_mt = log((1.0d0 + 2.0d0 * x + x2) &
& * (1.0d0 + x2) / 8.0d0) &
& - 2.0d0 * (atan(x) - atan(1.0d0)) ! L-Y eqn. 8d
psi_ht = 2.0d0 * log((1.0d0 + x2) / 2.0d0) ! L-Y eqn. 8e
end if
!
! water vapor level
!
zetaq = karman * bstar * altq / (ustar * ustar) ! L-Y eqn. 8a
zetaq = sign(min(abs(zetaq),10.d0),zetaq) ! undocumented NCAR
x2 = sqrt(abs(1.0d0 - 16.0d0 * zetaq)) ! L-Y eqn. 8b (squared)
x2 = max(x2,1.0d0) ! undocumented NCAR
x = sqrt(x2)
if (zetaq > 0.0d0) then
psi_mq = - 5.0d0 * zetaq ! L-Y eqn. 8c
psi_hq = - 5.0d0 * zetaq ! L-Y eqn. 8c
else
psi_mq = log((1.0d0 + 2.0d0 * x + x2) &
& * (1.0d0 + x2) / 8.0d0) &
& - 2.0d0 * (atan(x) - atan(1.0d0)) ! L-Y eqn. 8d
psi_hq = 2.0d0 * log((1.0d0 + x2) / 2.0d0) ! L-Y eqn. 8e
end if
!
! target level
!
zetatrg = karman * bstar * alt_target / (ustar * ustar) ! L-Y eqn. 8a
zetatrg = sign(min(abs(zetatrg),10.d0),zetatrg) ! undocumented NCAR
x2 = sqrt(abs(1.0d0 - 16.0d0 * zetatrg)) ! L-Y eqn. 8b (squared)
x2 = max(x2,1.0d0) ! undocumented NCAR
x = sqrt(x2)
if (zetatrg > 0.0d0) then
psi_mtrg = - 5.0d0 * zetatrg ! L-Y eqn. 8c
psi_htrg = - 5.0d0 * zetatrg ! L-Y eqn. 8c
else
psi_mtrg = log((1.0d0 + 2.0d0 * x + x2) &
& * (1.0d0 + x2) / 8.0d0) &
& - 2.0d0 * (atan(x) - atan(1.0d0)) ! L-Y eqn. 8d
psi_htrg = 2.0d0 * log((1.0d0 + x2) / 2.0d0) ! L-Y eqn. 8e
end if
! ^^^^^^^^^ the code above is unused ^^^^^^^^^
! neglect stability: set psi_* to zero in L-Y eqn. 9bc
tmptrg(m) = satmos - tab - tstar * (log(altt / alt_target)) / karman ! L-Y eqn. 9b, neglecting psi_h terms
sphtrg(m) = qatmos - qstar * (log(altq / alt_target)) / karman ! L-Y eqn. 9c, neglecting psi_h terms
!tmptrg(m) = satmos - tab - tstar * &
! & (log(altt / alt_target) + psi_htrg - psi_ht) / karman ! L-Y eqn. 9b
!sphtrg(m) = qatmos - qstar * &
! & (log(altq / alt_target) + psi_htrg - psi_hq) / karman ! L-Y eqn. 9c
end if
sphtrg(m) = max(sphtrg(m),sphmin)
! if (sphtrg(m) < 0.0d0) then
! if (sphtrg(m) < sphmin) then
!
! !write(6,*) ' !!!!! WARNING !!!!! '
! !write(6,'(1a,i8,1a,i8,1a,e10.3,1a,f8.3,f5.1)') ' humidity at ', m, '/', num_data, &
! ! & ' is negative ', real(sphtrg(m),4), ' temp = ', real(tmptrg(m),4), real(aexl(m),4)
! write(6,'(1a,i8,1a,i8,1a,e8.3,1a,f8.3,f5.1)') ' humidity at ', m, '/', num_data, &
! & ' is replaced with ', real(sphmin,4), ' temp = ', real(tmptrg(m),4), real(aexl(m),4)
! sphtrg(m) = sphmin
!
! write(15,*)
!
! if (aexl(m) == 1.0d0) then
! write(15,*) 'Ocean', m, ice(m)
! write(15,*) 'qar ', real(sphtrg(m),4), real(qatmos,4), real(qs,4)
! write(15,*) 'sat ', real(tmptrg(m),4), real(satmos-tab,4),real(tssurf,4)
! else
! write(15,*) 'Land', m, ice(m)
! write(15,*) 'qar ', real(sphtrg(m),4), real(qatmos,4), real(qs,4)
! write(15,*) 'sat ', real(tmptrg(m),4), real(satmos-tab,4),real(tssurf,4)
! write(15,*) 'qstar ', real(qstar,4)
! write(15,*) 'tstar ', real(tstar,4)
! write(15,*) 'bstar ', real(bstar,4)
! write(15,*) 'ustar ', real(ustar,4)
! write(15,*) 'psi(10), psi(2) ', psi_htrg, psi_hq
! write(15,*) 'factor ', (log(altq / alt_target) + psi_htrg - psi_hq) / karman
! write(15,*) 'zetaq ', real(zetaq,4)
! write(15,*) 'zetatrg ', real(zetatrg,4)
! end if
! write(15,*) 'sens, evap ', real(sens(m),4), real(evap(m),4)
!
! end if
end do
!$omp end parallel
end subroutine bulk_shift