1 | ! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt |
---|
2 | ! SPDX-License-Identifier: GPL-3.0-or-later |
---|
3 | |
---|
4 | subroutine readwind(indj,n,uuh,vvh,wwh) |
---|
5 | |
---|
6 | !********************************************************************** |
---|
7 | ! * |
---|
8 | ! TRAJECTORY MODEL SUBROUTINE READWIND * |
---|
9 | ! * |
---|
10 | !********************************************************************** |
---|
11 | ! * |
---|
12 | ! AUTHOR: G. WOTAWA * |
---|
13 | ! DATE: 1997-08-05 * |
---|
14 | ! LAST UPDATE: 2000-10-17, Andreas Stohl * |
---|
15 | ! * |
---|
16 | !********************************************************************** |
---|
17 | ! Changes, Bernd C. Krueger, Feb. 2001: |
---|
18 | ! Variables tth and qvh (on eta coordinates) in common block |
---|
19 | !********************************************************************** |
---|
20 | ! * |
---|
21 | ! DESCRIPTION: * |
---|
22 | ! * |
---|
23 | ! READING OF ECMWF METEOROLOGICAL FIELDS FROM INPUT DATA FILES. THE * |
---|
24 | ! INPUT DATA FILES ARE EXPECTED TO BE AVAILABLE IN GRIB CODE * |
---|
25 | ! * |
---|
26 | ! INPUT: * |
---|
27 | ! indj indicates number of the wind field to be read in * |
---|
28 | ! n temporal index for meteorological fields (1 to 3)* |
---|
29 | ! * |
---|
30 | ! IMPORTANT VARIABLES FROM COMMON BLOCK: * |
---|
31 | ! * |
---|
32 | ! wfname File name of data to be read in * |
---|
33 | ! nx,ny,nuvz,nwz expected field dimensions * |
---|
34 | ! nlev_ec number of vertical levels ecmwf model * |
---|
35 | ! uu,vv,ww wind fields * |
---|
36 | ! tt,qv temperature and specific humidity * |
---|
37 | ! ps surface pressure * |
---|
38 | ! * |
---|
39 | !********************************************************************** |
---|
40 | |
---|
41 | use par_mod |
---|
42 | use com_mod |
---|
43 | |
---|
44 | implicit none |
---|
45 | |
---|
46 | real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
47 | real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax) |
---|
48 | real :: wwh(0:nxmax-1,0:nymax-1,nwzmax) |
---|
49 | integer :: indj,i,j,k,n,levdiff2,ifield,iumax,iwmax,lunit |
---|
50 | |
---|
51 | ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING |
---|
52 | |
---|
53 | ! dimension of isec2 at least (22+n), where n is the number of parallels or |
---|
54 | ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid |
---|
55 | |
---|
56 | ! dimension of zsec2 at least (10+nn), where nn is the number of vertical |
---|
57 | ! coordinate parameters |
---|
58 | |
---|
59 | integer :: isec0(2),isec1(56),isec2(22+nxmax+nymax),isec3(2) |
---|
60 | integer :: isec4(64),inbuff(jpack),ilen,ierr,iword |
---|
61 | !integer iswap |
---|
62 | real :: zsec2(60+2*nuvzmax),zsec3(2),zsec4(jpunp) |
---|
63 | real :: xaux,yaux,xaux0,yaux0 |
---|
64 | real,parameter :: eps=1.e-4 |
---|
65 | real :: ewss(0:nxmax-1,0:nymax-1),nsss(0:nxmax-1,0:nymax-1) |
---|
66 | real :: plev1,pmean,tv,fu,hlev1,ff10m,fflev1 |
---|
67 | |
---|
68 | character(len=1) :: yoper = 'D' |
---|
69 | logical :: hflswitch,strswitch |
---|
70 | |
---|
71 | hflswitch=.false. |
---|
72 | strswitch=.false. |
---|
73 | levdiff2=nlev_ec-nwz+1 |
---|
74 | iumax=0 |
---|
75 | iwmax=0 |
---|
76 | |
---|
77 | ! |
---|
78 | ! OPENING OF DATA FILE (GRIB CODE) |
---|
79 | ! |
---|
80 | 5 call pbopen(lunit,path(3)(1:length(3))//wfname(indj),'r',ierr) |
---|
81 | if(ierr.lt.0) goto 999 |
---|
82 | |
---|
83 | ifield=0 |
---|
84 | 10 ifield=ifield+1 |
---|
85 | ! |
---|
86 | ! GET NEXT FIELDS |
---|
87 | ! |
---|
88 | call pbgrib(lunit,inbuff,jpack,ilen,ierr) |
---|
89 | if(ierr.eq.-1) goto 50 ! EOF DETECTED |
---|
90 | if(ierr.lt.-1) goto 888 ! ERROR DETECTED |
---|
91 | |
---|
92 | ierr=1 |
---|
93 | |
---|
94 | ! Check whether we are on a little endian or on a big endian computer |
---|
95 | !******************************************************************** |
---|
96 | |
---|
97 | !if (inbuff(1).eq.1112101447) then ! little endian, swap bytes |
---|
98 | ! iswap=1+ilen/4 |
---|
99 | ! call swap32(inbuff,iswap) |
---|
100 | !else if (inbuff(1).ne.1196575042) then ! big endian |
---|
101 | ! stop 'subroutine gridcheck: corrupt GRIB data' |
---|
102 | !endif |
---|
103 | |
---|
104 | |
---|
105 | call gribex(isec0,isec1,isec2,zsec2,isec3,zsec3,isec4, & |
---|
106 | zsec4,jpunp,inbuff,jpack,iword,yoper,ierr) |
---|
107 | if (ierr.ne.0) goto 888 ! ERROR DETECTED |
---|
108 | |
---|
109 | if(ifield.eq.1) then |
---|
110 | |
---|
111 | ! CHECK GRID SPECIFICATIONS |
---|
112 | |
---|
113 | if(isec2(2).ne.nxfield) stop 'READWIND: NX NOT CONSISTENT' |
---|
114 | if(isec2(3).ne.ny) stop 'READWIND: NY NOT CONSISTENT' |
---|
115 | if(isec2(12)/2-1.ne.nlev_ec) & |
---|
116 | stop 'READWIND: VERTICAL DISCRETIZATION NOT CONSISTENT' |
---|
117 | xaux=real(isec2(5))/1000.+real(nxshift)*dx |
---|
118 | yaux=real(isec2(7))/1000. |
---|
119 | xaux0=xlon0 |
---|
120 | yaux0=ylat0 |
---|
121 | if(xaux.lt.0.) xaux=xaux+360. |
---|
122 | if(yaux.lt.0.) yaux=yaux+360. |
---|
123 | if(xaux0.lt.0.) xaux0=xaux0+360. |
---|
124 | if(yaux0.lt.0.) yaux0=yaux0+360. |
---|
125 | if(abs(xaux-xaux0).gt.eps) & |
---|
126 | stop 'READWIND: LOWER LEFT LONGITUDE NOT CONSISTENT' |
---|
127 | if(abs(yaux-yaux0).gt.eps) & |
---|
128 | stop 'READWIND: LOWER LEFT LATITUDE NOT CONSISTENT' |
---|
129 | endif |
---|
130 | |
---|
131 | do j=0,nymin1 |
---|
132 | do i=0,nxfield-1 |
---|
133 | k=isec1(8) |
---|
134 | if(isec1(6).eq.130) tth(i,j,nlev_ec-k+2,n)= &!! TEMPERATURE |
---|
135 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
136 | if(isec1(6).eq.131) uuh(i,j,nlev_ec-k+2)= &!! U VELOCITY |
---|
137 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
138 | if(isec1(6).eq.132) vvh(i,j,nlev_ec-k+2)= &!! V VELOCITY |
---|
139 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
140 | if(isec1(6).eq.133) then !! SPEC. HUMIDITY |
---|
141 | qvh(i,j,nlev_ec-k+2,n)=zsec4(nxfield*(ny-j-1)+i+1) |
---|
142 | if (qvh(i,j,nlev_ec-k+2,n) .lt. 0.) & |
---|
143 | qvh(i,j,nlev_ec-k+2,n) = 0. |
---|
144 | ! this is necessary because the gridded data may contain |
---|
145 | ! spurious negative values |
---|
146 | endif |
---|
147 | if(isec1(6).eq.134) ps(i,j,1,n)= &!! SURF. PRESS. |
---|
148 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
149 | |
---|
150 | if(isec1(6).eq.135) wwh(i,j,nlev_ec-k+1)= &!! W VELOCITY |
---|
151 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
152 | if(isec1(6).eq.141) sd(i,j,1,n)= &!! SNOW DEPTH |
---|
153 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
154 | if(isec1(6).eq.151) msl(i,j,1,n)= &!! SEA LEVEL PRESS. |
---|
155 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
156 | if(isec1(6).eq.164) tcc(i,j,1,n)= &!! CLOUD COVER |
---|
157 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
158 | if(isec1(6).eq.165) u10(i,j,1,n)= &!! 10 M U VELOCITY |
---|
159 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
160 | if(isec1(6).eq.166) v10(i,j,1,n)= &!! 10 M V VELOCITY |
---|
161 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
162 | if(isec1(6).eq.167) tt2(i,j,1,n)= &!! 2 M TEMPERATURE |
---|
163 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
164 | if(isec1(6).eq.168) td2(i,j,1,n)= &!! 2 M DEW POINT |
---|
165 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
166 | if(isec1(6).eq.142) then !! LARGE SCALE PREC. |
---|
167 | lsprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1) |
---|
168 | if (lsprec(i,j,1,n).lt.0.) lsprec(i,j,1,n)=0. |
---|
169 | endif |
---|
170 | if(isec1(6).eq.143) then !! CONVECTIVE PREC. |
---|
171 | convprec(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1) |
---|
172 | if (convprec(i,j,1,n).lt.0.) convprec(i,j,1,n)=0. |
---|
173 | endif |
---|
174 | if(isec1(6).eq.146) sshf(i,j,1,n)= &!! SENS. HEAT FLUX |
---|
175 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
176 | if((isec1(6).eq.146).and.(zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) & |
---|
177 | hflswitch=.true. ! Heat flux available |
---|
178 | if(isec1(6).eq.176) then !! SOLAR RADIATION |
---|
179 | ssr(i,j,1,n)=zsec4(nxfield*(ny-j-1)+i+1) |
---|
180 | if (ssr(i,j,1,n).lt.0.) ssr(i,j,1,n)=0. |
---|
181 | endif |
---|
182 | if(isec1(6).eq.180) ewss(i,j)= &!! EW SURFACE STRESS |
---|
183 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
184 | if(isec1(6).eq.181) nsss(i,j)= &!! NS SURFACE STRESS |
---|
185 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
186 | if(((isec1(6).eq.180).or.(isec1(6).eq.181)).and. & |
---|
187 | (zsec4(nxfield*(ny-j-1)+i+1).ne.0.)) strswitch=.true. ! stress available |
---|
188 | if(isec1(6).eq.129) oro(i,j)= &!! ECMWF OROGRAPHY |
---|
189 | zsec4(nxfield*(ny-j-1)+i+1)/ga |
---|
190 | if(isec1(6).eq.160) excessoro(i,j)= &!! STANDARD DEVIATION OF OROGRAPHY |
---|
191 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
192 | if(isec1(6).eq.172) lsm(i,j)= &!! ECMWF LAND SEA MASK |
---|
193 | zsec4(nxfield*(ny-j-1)+i+1) |
---|
194 | if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1) |
---|
195 | if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1) |
---|
196 | end do |
---|
197 | end do |
---|
198 | |
---|
199 | goto 10 !! READ NEXT LEVEL OR PARAMETER |
---|
200 | ! |
---|
201 | ! CLOSING OF INPUT DATA FILE |
---|
202 | ! |
---|
203 | 50 call pbclose(lunit,ierr) !! FINNISHED READING / CLOSING GRIB FILE |
---|
204 | |
---|
205 | if(levdiff2.eq.0) then |
---|
206 | iwmax=nlev_ec+1 |
---|
207 | do i=0,nxmin1 |
---|
208 | do j=0,nymin1 |
---|
209 | wwh(i,j,nlev_ec+1)=0. |
---|
210 | end do |
---|
211 | end do |
---|
212 | endif |
---|
213 | |
---|
214 | |
---|
215 | ! For global fields, assign the leftmost data column also to the rightmost |
---|
216 | ! data column; if required, shift whole grid by nxshift grid points |
---|
217 | !************************************************************************* |
---|
218 | |
---|
219 | if (xglobal) then |
---|
220 | call shift_field_0(ewss,nxfield,ny) |
---|
221 | call shift_field_0(nsss,nxfield,ny) |
---|
222 | call shift_field_0(oro,nxfield,ny) |
---|
223 | call shift_field_0(excessoro,nxfield,ny) |
---|
224 | call shift_field_0(lsm,nxfield,ny) |
---|
225 | call shift_field(ps,nxfield,ny,1,1,2,n) |
---|
226 | call shift_field(sd,nxfield,ny,1,1,2,n) |
---|
227 | call shift_field(msl,nxfield,ny,1,1,2,n) |
---|
228 | call shift_field(tcc,nxfield,ny,1,1,2,n) |
---|
229 | call shift_field(u10,nxfield,ny,1,1,2,n) |
---|
230 | call shift_field(v10,nxfield,ny,1,1,2,n) |
---|
231 | call shift_field(tt2,nxfield,ny,1,1,2,n) |
---|
232 | call shift_field(td2,nxfield,ny,1,1,2,n) |
---|
233 | call shift_field(lsprec,nxfield,ny,1,1,2,n) |
---|
234 | call shift_field(convprec,nxfield,ny,1,1,2,n) |
---|
235 | call shift_field(sshf,nxfield,ny,1,1,2,n) |
---|
236 | call shift_field(ssr,nxfield,ny,1,1,2,n) |
---|
237 | call shift_field(tth,nxfield,ny,nuvzmax,nuvz,2,n) |
---|
238 | call shift_field(qvh,nxfield,ny,nuvzmax,nuvz,2,n) |
---|
239 | call shift_field(uuh,nxfield,ny,nuvzmax,nuvz,1,1) |
---|
240 | call shift_field(vvh,nxfield,ny,nuvzmax,nuvz,1,1) |
---|
241 | call shift_field(wwh,nxfield,ny,nwzmax,nwz,1,1) |
---|
242 | endif |
---|
243 | |
---|
244 | do i=0,nxmin1 |
---|
245 | do j=0,nymin1 |
---|
246 | surfstr(i,j,1,n)=sqrt(ewss(i,j)**2+nsss(i,j)**2) |
---|
247 | end do |
---|
248 | end do |
---|
249 | |
---|
250 | if ((.not.hflswitch).or.(.not.strswitch)) then |
---|
251 | write(*,*) 'WARNING: No flux data contained in GRIB file ', & |
---|
252 | wfname(indj) |
---|
253 | |
---|
254 | ! CALCULATE USTAR AND SSHF USING THE PROFILE METHOD |
---|
255 | ! As ECMWF has increased the model resolution, such that now the first model |
---|
256 | ! level is at about 10 m (where 10-m wind is given), use the 2nd ECMWF level |
---|
257 | ! (3rd model level in FLEXPART) for the profile method |
---|
258 | !*************************************************************************** |
---|
259 | |
---|
260 | do i=0,nxmin1 |
---|
261 | do j=0,nymin1 |
---|
262 | plev1=akz(3)+bkz(3)*ps(i,j,1,n) |
---|
263 | pmean=0.5*(ps(i,j,1,n)+plev1) |
---|
264 | tv=tth(i,j,3,n)*(1.+0.61*qvh(i,j,3,n)) |
---|
265 | fu=-r_air*tv/ga/pmean |
---|
266 | hlev1=fu*(plev1-ps(i,j,1,n)) ! HEIGTH OF FIRST MODEL LAYER |
---|
267 | ff10m= sqrt(u10(i,j,1,n)**2+v10(i,j,1,n)**2) |
---|
268 | fflev1=sqrt(uuh(i,j,3)**2+vvh(i,j,3)**2) |
---|
269 | call pbl_profile(ps(i,j,1,n),td2(i,j,1,n),hlev1, & |
---|
270 | tt2(i,j,1,n),tth(i,j,3,n),ff10m,fflev1, & |
---|
271 | surfstr(i,j,1,n),sshf(i,j,1,n)) |
---|
272 | if(sshf(i,j,1,n).gt.200.) sshf(i,j,1,n)=200. |
---|
273 | if(sshf(i,j,1,n).lt.-400.) sshf(i,j,1,n)=-400. |
---|
274 | end do |
---|
275 | end do |
---|
276 | endif |
---|
277 | |
---|
278 | |
---|
279 | ! Assign 10 m wind to model level at eta=1.0 to have one additional model |
---|
280 | ! level at the ground |
---|
281 | ! Specific humidity is taken the same as at one level above |
---|
282 | ! Temperature is taken as 2 m temperature |
---|
283 | !************************************************************************** |
---|
284 | |
---|
285 | do i=0,nxmin1 |
---|
286 | do j=0,nymin1 |
---|
287 | uuh(i,j,1)=u10(i,j,1,n) |
---|
288 | vvh(i,j,1)=v10(i,j,1,n) |
---|
289 | qvh(i,j,1,n)=qvh(i,j,2,n) |
---|
290 | tth(i,j,1,n)=tt2(i,j,1,n) |
---|
291 | end do |
---|
292 | end do |
---|
293 | |
---|
294 | if(iumax.ne.nuvz-1) stop 'READWIND: NUVZ NOT CONSISTENT' |
---|
295 | if(iwmax.ne.nwz) stop 'READWIND: NWZ NOT CONSISTENT' |
---|
296 | |
---|
297 | return |
---|
298 | 888 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD #### ' |
---|
299 | write(*,*) ' #### ',wfname(indj),' #### ' |
---|
300 | write(*,*) ' #### IS NOT GRIB FORMAT !!! #### ' |
---|
301 | stop 'Execution terminated' |
---|
302 | 999 write(*,*) ' #### FLEXPART MODEL ERROR! WINDFIELD #### ' |
---|
303 | write(*,*) ' #### ',wfname(indj),' #### ' |
---|
304 | write(*,*) ' #### CANNOT BE OPENED !!! #### ' |
---|
305 | stop 'Execution terminated' |
---|
306 | |
---|
307 | end subroutine readwind |
---|