1 | !********************************************************************** |
---|
2 | ! Copyright 2016 * |
---|
3 | ! Andreas Stohl, Massimo Cassiani, Petra Seibert, A. Frank, * |
---|
4 | ! Gerhard Wotawa, Caroline Forster, Sabine Eckhardt, John Burkhart, * |
---|
5 | ! Harald Sodemann, Ignacio Pisso * |
---|
6 | ! * |
---|
7 | ! This file is part of FLEXPART-NorESM * |
---|
8 | ! * |
---|
9 | ! FLEXPART-NorESM is free software: you can redistribute it * |
---|
10 | ! and/or modify * |
---|
11 | ! it under the terms of the GNU General Public License as published by* |
---|
12 | ! the Free Software Foundation, either version 3 of the License, or * |
---|
13 | ! (at your option) any later version. * |
---|
14 | ! * |
---|
15 | ! FLEXPART-NorESM is distributed in the hope that it will be useful, * |
---|
16 | ! but WITHOUT ANY WARRANTY; without even the implied warranty of * |
---|
17 | ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * |
---|
18 | ! GNU General Public License for more details. * |
---|
19 | ! * |
---|
20 | ! You should have received a copy of the GNU General Public License * |
---|
21 | ! along with FLEXPART-NorESM. * |
---|
22 | ! If not, see <http://www.gnu.org/licenses/>. * |
---|
23 | !********************************************************************** |
---|
24 | |
---|
25 | |
---|
26 | subroutine read_delta_ps_intime(indj,index) |
---|
27 | |
---|
28 | !* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * |
---|
29 | !* * |
---|
30 | !* Read pressure fileds at time t-dt and t+dt * |
---|
31 | !* note: part of netcdf reading strucutre adapted from the routines by * |
---|
32 | !* Fast and Easter in FLEXPART-WRF * |
---|
33 | !* * |
---|
34 | !* * |
---|
35 | !* Author: * |
---|
36 | !* M. Cassiani 2016 * |
---|
37 | !* * |
---|
38 | !* * |
---|
39 | !* * |
---|
40 | !* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * |
---|
41 | |
---|
42 | |
---|
43 | |
---|
44 | use noresm_variables |
---|
45 | use par_mod |
---|
46 | use com_mod |
---|
47 | use conv_mod |
---|
48 | use cmapf_mod, only: stlmbr,stcm2p |
---|
49 | |
---|
50 | implicit none |
---|
51 | |
---|
52 | include 'netcdf.inc' |
---|
53 | |
---|
54 | !integer maxdim,maxvar,maxtime !moved to par_mod.f90 |
---|
55 | real(kind=4) :: eps |
---|
56 | |
---|
57 | !parameter(maxdim=9,maxvar=60,maxtime=12, !moved in in common block par_mod |
---|
58 | parameter(eps=0.0001) |
---|
59 | |
---|
60 | integer :: gotGrid |
---|
61 | real(kind=4) :: sizesouth,sizenorth,xauxa,pint |
---|
62 | |
---|
63 | !c--------------------------------------------- |
---|
64 | real(kind=4) :: xaux1,xaux2,yaux1,yaux2 |
---|
65 | real(kind=dp) :: xaux1in,xaux2in,yaux1in,yaux2in |
---|
66 | |
---|
67 | integer :: nyfieldin,nxfieldin |
---|
68 | |
---|
69 | !C-------------- extra function to calculate the dew point ------ |
---|
70 | real(kind=4) :: dewpoint ! this is a function |
---|
71 | |
---|
72 | !c---------------varibles used to calculate the time interval --------- |
---|
73 | integer :: date_aid(maxtime),idumb,itime ,dimtimenum,indextime, & |
---|
74 | time_interval |
---|
75 | real(kind=dp) :: jul,juldate,dumb !variable jul used for date in days, juldate is a function |
---|
76 | |
---|
77 | !C-------------- some declaration for netcdf reading ------ |
---|
78 | real(kind=4) :: duma |
---|
79 | real(kind=4), allocatable, dimension(:) :: duma_alloc |
---|
80 | integer :: ierr !error code message |
---|
81 | integer :: idiagaa !flag |
---|
82 | integer :: id_var,id_dim(maxdim) |
---|
83 | integer :: nvar_exp_in_grid_atm_nc, & |
---|
84 | nvar_exp_in_meteo_field |
---|
85 | integer :: i,ii,j,k, iatt, idimid_unlim, idum, iret, ivtype |
---|
86 | integer :: ix,jy,kz |
---|
87 | integer :: lenatt, lendim(maxdim) |
---|
88 | integer :: natts_tot, ncid, ndims_tot, nvars_tot |
---|
89 | integer :: sizetype |
---|
90 | character*110 :: fnamenc |
---|
91 | character*80 :: dimname(maxdim) |
---|
92 | character*80 :: attname |
---|
93 | character*160 :: varname,vartype |
---|
94 | character*160 :: varnamev(maxvar) |
---|
95 | character*160 :: units(maxvar) |
---|
96 | character*160 :: vartypev(maxvar) |
---|
97 | integer :: ndimsv(maxvar) |
---|
98 | integer :: dimidsv(maxvar,maxdim) |
---|
99 | character*1000 :: dumch1000 |
---|
100 | |
---|
101 | !integer istart(maxdim),icount(maxdim) |
---|
102 | integer :: xtype,xtypev(maxvar) |
---|
103 | integer :: ndims |
---|
104 | integer :: dimids(maxdim),dimidsm(maxvar,maxdim) |
---|
105 | integer :: LENDIM_EXP(maxdim),LENDIM_MAX(maxdim) |
---|
106 | integer :: natts |
---|
107 | integer :: varid |
---|
108 | !C-------------------------------------------------------------------------------- |
---|
109 | integer :: index !index of the time locaction with respect to wftime(indj) in readwind: 1 for indj-1 & 2 for indj+1 |
---|
110 | |
---|
111 | real(kind=4) :: plev1,pmean,tv,fu,hlev1,ff10m,fflev1 |
---|
112 | |
---|
113 | logical :: hflswitch,strswitch |
---|
114 | integer :: indj,n |
---|
115 | |
---|
116 | !*************************************************************************************************** |
---|
117 | |
---|
118 | |
---|
119 | 9100 format( / '** read_deltap_intime -- ', a /& |
---|
120 | 'file = ', a ) |
---|
121 | 9110 format( / '** read_deltap_intime -- ', a, 1x, i8 /& |
---|
122 | 'file = ', a ) |
---|
123 | 9120 format( / '** read_deltap_intime -- ', a, 2(1x,i8) /& |
---|
124 | 'file = ', a ) |
---|
125 | 9030 format( a, 2i6, 2(2x,a) ) |
---|
126 | 9031 format( 1i4,1a20,2i4,1a30) |
---|
127 | |
---|
128 | idiagaa=0 !set it to zero supress some diagnostic in files 73,74 see opening below! |
---|
129 | !if (idiagaa.eq.1) then |
---|
130 | !open(71,file='..\options\data_type.txt') ! for testing: mc |
---|
131 | !open(72,file='..\options\list_windfield.txt') ! for testing: mc |
---|
132 | !open(73,file='..\options\list_global_att.txt') ! for testing: mc |
---|
133 | !open(74,file='..\options\list_variable_att.txt') ! for testing: mc |
---|
134 | !open(75,file='..\options\seq_diagnostict.txt') ! for testing: mc |
---|
135 | !continue |
---|
136 | !end if |
---|
137 | |
---|
138 | gotgrid=0 |
---|
139 | ierr=0 |
---|
140 | |
---|
141 | !c-------------- open 4D (3D plus dtime dimension) wind meteo file |
---|
142 | if (indj.le.1) goto 100 |
---|
143 | fnamenc=path(3)(1:length(3))& |
---|
144 | //trim(wfname(indj)) |
---|
145 | |
---|
146 | iret = nf_open(fnamenc,NF_NOWRITE,ncid) |
---|
147 | if (iret .ne. nf_noerr) then |
---|
148 | write(*,9100) 'error doing open', fnamenc |
---|
149 | ierr = -1 |
---|
150 | !stop |
---|
151 | end if |
---|
152 | !c-------------- get information on dimensions --------- |
---|
153 | iret = nf_inq( ncid, & |
---|
154 | ndims_tot, nvars_tot, natts_tot, idimid_unlim ) |
---|
155 | if (iret .ne. nf_noerr) then |
---|
156 | write(*,9100) 'error inquiring dimensions', fnamenc |
---|
157 | ierr = -2 |
---|
158 | stop |
---|
159 | end if |
---|
160 | |
---|
161 | |
---|
162 | !c-------------inquiring abouot dimensions name ------- |
---|
163 | dimtimenum=0 |
---|
164 | do i = 1, min(ndims_tot,maxdim) |
---|
165 | iret = nf_inq_dim( ncid, i, dimname(i), lendim(i) ) |
---|
166 | if (iret .ne. nf_noerr) then |
---|
167 | write(*,9110) 'error inquiring dimensions for dim#', i, fnamenc |
---|
168 | ierr = -2 |
---|
169 | stop |
---|
170 | end if |
---|
171 | if (dimname(i).eq.'time') dimtimenum=i |
---|
172 | end do |
---|
173 | !c----------------------------------------------------------- |
---|
174 | |
---|
175 | !C------------inquiring about global attributes ------------- |
---|
176 | |
---|
177 | if (idiagaa .gt. 0) then |
---|
178 | write(73,*) |
---|
179 | write(73,*) 'attribute #, name, type, value' |
---|
180 | end if |
---|
181 | |
---|
182 | do 3401 iatt = 1, natts_tot |
---|
183 | iret = nf_inq_attname( ncid, nf_global, iatt, attname) |
---|
184 | if (iret .ne. nf_noerr) goto 3601 |
---|
185 | iret = nf_inq_att( ncid, nf_global, attname, ivtype, lenatt ) |
---|
186 | if (iret .ne. nf_noerr) goto 3601 |
---|
187 | if (ivtype .eq. 2) then |
---|
188 | iret = nf_get_att_text( ncid, nf_global, attname, dumch1000 ) |
---|
189 | if (iret .ne. nf_noerr) goto 3601 |
---|
190 | i = max(1,min(1000,lenatt)) |
---|
191 | if (idiagaa .gt. 0) write(73,91010) & |
---|
192 | iatt, attname(1:40), ivtype, lenatt, dumch1000(1:i) |
---|
193 | else if (ivtype .eq. 4) then |
---|
194 | iret = nf_get_att_int( ncid, nf_global, attname, idum ) |
---|
195 | if (iret .ne. nf_noerr) goto 3601 |
---|
196 | if (idiagaa .gt. 0) write(73,91020) & |
---|
197 | iatt, attname(1:40), ivtype, lenatt, idum |
---|
198 | else if ((ivtype .eq. 5) .and. (lenatt .eq. 1)) then |
---|
199 | iret = nf_get_att_real( ncid, nf_global, attname, duma ) |
---|
200 | if (iret .ne. nf_noerr) goto 3601 |
---|
201 | if (idiagaa .gt. 0) write(73,91030) & |
---|
202 | iatt, attname(1:40), ivtype, lenatt, duma |
---|
203 | else if ((ivtype .eq. 5) .and. (lenatt .gt. 1)) then |
---|
204 | allocate(duma_alloc(lenatt)) |
---|
205 | iret = nf_get_att_real( ncid, nf_global, attname, duma_alloc ) |
---|
206 | if (iret .ne. nf_noerr) goto 3601 |
---|
207 | if (idiagaa .gt. 0) then |
---|
208 | write(73,91010) iatt, attname(1:40), ivtype, lenatt |
---|
209 | write(73,91040) (duma_alloc(i), i=1,lenatt) |
---|
210 | end if |
---|
211 | deallocate(duma_alloc) |
---|
212 | else |
---|
213 | if (idiagaa .gt. 0) write(73,'(i4,1x,a,2(1x,i6))') & |
---|
214 | iatt, attname(1:40), ivtype, lenatt |
---|
215 | goto 3401 |
---|
216 | endif |
---|
217 | |
---|
218 | 3401 continue |
---|
219 | 91010 format( i4, 1x, a, 2(1x,i6), 1x, a ) |
---|
220 | 91020 format( i4, 1x, a, 2(1x,i6), 1x, i10 ) |
---|
221 | 91030 format( i4, 1x, a, 2(1x,i6), 1x, 1pe12.4 ) |
---|
222 | 91040 format(( 12x, 5(1pe12.4) )) |
---|
223 | |
---|
224 | 91050 format(a,1x,2(1x,i5), 12(1x,e15.8)) |
---|
225 | goto 3901 |
---|
226 | |
---|
227 | 3601 write(*,9110) 'error inquiring global attribute', iatt, fnamenc |
---|
228 | write(73,9110) 'error inquiring global attribute', iatt, fnamenc |
---|
229 | stop |
---|
230 | |
---|
231 | 3901 continue |
---|
232 | |
---|
233 | |
---|
234 | !c----------------some checks on the varibles-------------------- |
---|
235 | do ii=1, nvars_tot |
---|
236 | varid=ii |
---|
237 | iret = nf_inq_var ( ncid, varid, varname, xtype, & |
---|
238 | ndims, dimids, natts) |
---|
239 | if (iret .ne. nf_noerr) then |
---|
240 | write (*,*) 'error in nf_inq_var' |
---|
241 | stop |
---|
242 | end if |
---|
243 | |
---|
244 | do iatt = 1, natts |
---|
245 | iret = nf_inq_attname( ncid, varid, iatt, attname) |
---|
246 | if (iret .ne. nf_noerr) goto 3602 |
---|
247 | iret = nf_inq_att( ncid, varid, attname, ivtype, lenatt ) |
---|
248 | if (iret .ne. nf_noerr) goto 3602 |
---|
249 | if (ivtype .eq. 2) then |
---|
250 | dumch1000='' |
---|
251 | iret = nf_get_att_text( ncid, varid, attname, dumch1000 ) |
---|
252 | if (iret .ne. nf_noerr) goto 3602 |
---|
253 | i = max(1,min(1000,lenatt)) |
---|
254 | if (attname.eq.'units') then |
---|
255 | units(ii)=dumch1000 |
---|
256 | end if |
---|
257 | if (idiagaa .gt. 0) write(74,91010) & |
---|
258 | iatt, attname(1:40), ivtype, lenatt, dumch1000(1:i) |
---|
259 | else if (ivtype .eq. 4) then |
---|
260 | iret = nf_get_att_int( ncid, varid, attname, idum ) |
---|
261 | if (iret .ne. nf_noerr) goto 3602 |
---|
262 | if (idiagaa .gt. 0) write(74,91020) & |
---|
263 | iatt, attname(1:40), ivtype, lenatt, idum |
---|
264 | else if ((ivtype .eq. 5) .and. (lenatt .eq. 1)) then |
---|
265 | iret = nf_get_att_real( ncid, varid, attname, duma ) |
---|
266 | if (iret .ne. nf_noerr) goto 3602 |
---|
267 | if (idiagaa .gt. 0) write(74,91030) & |
---|
268 | iatt, attname(1:40), ivtype, lenatt, duma |
---|
269 | else if ((ivtype .eq. 5) .and. (lenatt .gt. 1)) then |
---|
270 | allocate( duma_alloc(lenatt) ) |
---|
271 | iret = nf_get_att_real( ncid, varid, attname, duma_alloc ) |
---|
272 | if (iret .ne. nf_noerr) goto 3602 |
---|
273 | if (idiagaa .gt. 0) then |
---|
274 | write(74,91010) iatt, attname(1:40), ivtype, lenatt |
---|
275 | write(74,91040) (duma_alloc(i), i=1,lenatt) |
---|
276 | end if |
---|
277 | deallocate( duma_alloc ) |
---|
278 | else |
---|
279 | if (idiagaa .gt. 0) write(74,'(i4,1x,a,2(1x,i6))') & |
---|
280 | iatt, attname(1:40), ivtype, lenatt |
---|
281 | exit |
---|
282 | endif |
---|
283 | |
---|
284 | end do |
---|
285 | goto 3902 |
---|
286 | 3602 write(*,9110) 'error inquiring variable attribute', varname, & |
---|
287 | iatt, attname, fnamenc |
---|
288 | write(74,9110) 'error inquiring variable attribute', varname, & |
---|
289 | iatt, attname, fnamenc |
---|
290 | stop |
---|
291 | 3902 continue |
---|
292 | |
---|
293 | !if (idiagaa.eq.1) then |
---|
294 | !write(75,*)'start',ncid, varid, varname, xtype, & !more diagnostic for testing mc |
---|
295 | !ndims, dimids, natts,'fine' |
---|
296 | !end if |
---|
297 | |
---|
298 | xtypev(ii)=xtype |
---|
299 | varnamev(ii)=varname |
---|
300 | ndimsv(ii)=ndims |
---|
301 | do j=1,maxdim |
---|
302 | dimidsm(ii,j)=dimids(j) |
---|
303 | end do |
---|
304 | |
---|
305 | !if (idiagaa.eq.1) then |
---|
306 | !write(71,9031)i,varnamev(ii),xtypev(ii),ndims,units(ii) !WRITE diagnostic to file BY M C |
---|
307 | !end if |
---|
308 | |
---|
309 | |
---|
310 | end do |
---|
311 | |
---|
312 | !C-------------- find the time -------------- |
---|
313 | |
---|
314 | varname='date' !read date |
---|
315 | call check_variable(varname,fnamenc,maxdim,nf_int, & |
---|
316 | id_var,ndims,id_dim,ierr,ncid) |
---|
317 | if (ierr.lt.0) goto 100 |
---|
318 | lendim_exp=0 |
---|
319 | if (ndims.eq.0) then !this shoudl happen if there is only one time per file |
---|
320 | ndims=1 |
---|
321 | lendim_exp(1)=1 |
---|
322 | else |
---|
323 | do j = 1, ndims |
---|
324 | lendim_exp(j) = lendim(id_dim(j)) |
---|
325 | end do |
---|
326 | end if |
---|
327 | call allocatedumarray(ndims,lendim_exp,maxdim,nf_int) |
---|
328 | do j=1, ndims |
---|
329 | istart(j) = 1 !ndims |
---|
330 | icount(j) = lendim_exp(j) |
---|
331 | end do |
---|
332 | !here we assume taht ndims=1 for date variable |
---|
333 | iret = & |
---|
334 | nf_get_vara_int( ncid, id_var, istart, icount, dumarray1D_int) |
---|
335 | |
---|
336 | if (iret .ne. nf_noerr) then |
---|
337 | write(*,9100) 'error inquiring var value ' // varname, fnamenc |
---|
338 | ierr = -5 |
---|
339 | goto 100 |
---|
340 | end if |
---|
341 | |
---|
342 | do i=1,lendim_exp(1) !one dimension only expected! |
---|
343 | date_aid(i)=dumarray1D_int(i) |
---|
344 | end do |
---|
345 | |
---|
346 | |
---|
347 | varname='datesec' !read date in seconds |
---|
348 | call check_variable(varname,fnamenc,maxdim,nf_int, & |
---|
349 | id_var,ndims,id_dim,ierr,ncid) |
---|
350 | if (ierr.lt.0) goto 100 |
---|
351 | lendim_exp=0 |
---|
352 | if (ndims.eq.0) then |
---|
353 | ndims=1 |
---|
354 | lendim_exp(1)=1 |
---|
355 | else |
---|
356 | do j = 1, ndims |
---|
357 | lendim_exp(j) = lendim(id_dim(j)) |
---|
358 | end do |
---|
359 | end if |
---|
360 | call allocatedumarray(ndims,lendim_exp,maxdim,nf_int) |
---|
361 | do j=1, ndims |
---|
362 | istart(j) = 1 !ndims |
---|
363 | icount(j) = lendim_exp(j) |
---|
364 | end do |
---|
365 | iret = & |
---|
366 | nf_get_vara_int( ncid, id_var, istart, icount, dumarray1D_int) |
---|
367 | if (iret .ne. nf_noerr) then |
---|
368 | write(*,9100) 'error inquiring var value ' // varname, fnamenc |
---|
369 | ierr = -5 |
---|
370 | goto 100 |
---|
371 | end if |
---|
372 | |
---|
373 | do i=1,lendim_exp(1) !one dimension only expected |
---|
374 | jul=juldate(date_aid(i),0) |
---|
375 | jul=jul+real(dumarray1D_int(i),kind=dp)/86400._dp |
---|
376 | dumb= (jul-bdate)*86400._dp |
---|
377 | idumb=nint(dumb) |
---|
378 | if (idumb.eq.wftime(indj)) exit |
---|
379 | end do |
---|
380 | itime=i !selected right file and record.... |
---|
381 | time_interval=wftime(indj+1)-wftime(indj-1) !for calculating time derivative of pressure see below |
---|
382 | !c-------- finish operation of finding the time -------------- |
---|
383 | print *,wfname(indj),wftime(indj),itime,n !write diagnostic to file for test reason by, by mc |
---|
384 | |
---|
385 | !c---------------- Inquiring about varnames and real fields --------------------- |
---|
386 | nvar_exp_in_meteo_field=1 |
---|
387 | varnamev(1)='PS' !horizontal wind m/s |
---|
388 | |
---|
389 | !---------------- check and load file contents |
---|
390 | do i=1,nvar_exp_in_meteo_field !variable_loop |
---|
391 | |
---|
392 | varname=varnamev(i) |
---|
393 | call check_variable(varname,fnamenc,maxdim,nf_float, & |
---|
394 | id_var,ndims,id_dim,ierr,ncid) |
---|
395 | if (ierr.ne.0) goto 100 |
---|
396 | lendim_exp=0 |
---|
397 | indextime=0 |
---|
398 | do j = 1, ndims |
---|
399 | lendim_exp(j) = lendim(id_dim(j)) |
---|
400 | if (id_dim(j).eq.dimtimenum) indextime=j !this indextime select whcih is the dimension associated with time |
---|
401 | end do |
---|
402 | if (dimtimenum.eq.0.and.indextime.eq.0) then !this means that there is no time dimension in the file and therefore we create an extra time dimension of leght 1. |
---|
403 | !notimedim=1 |
---|
404 | ndims=ndims+1 |
---|
405 | lendim_exp(ndims)=1 |
---|
406 | indextime=ndims |
---|
407 | end if |
---|
408 | |
---|
409 | if (indextime.ne.ndims) stop 'ERROR readwind the time is expected & |
---|
410 | to be always the last dimension in the stored variables' |
---|
411 | |
---|
412 | call allocatedumarray(ndims,lendim_exp,maxdim,nf_float) |
---|
413 | |
---|
414 | do j=1, ndims |
---|
415 | istart(j) = 1 !ndims |
---|
416 | icount(j) = lendim_exp(j) |
---|
417 | end do |
---|
418 | |
---|
419 | if (ndims.eq.1) then |
---|
420 | iret = & |
---|
421 | nf_get_vara_real( ncid, id_var, istart, icount, dumarray1D_real) |
---|
422 | if (iret .ne. nf_noerr) then |
---|
423 | write(*,9100) 'error inquiring var value ' // varname, fnamenc |
---|
424 | ierr = -5 |
---|
425 | goto 100 |
---|
426 | end if |
---|
427 | else if (ndims.eq.2) then |
---|
428 | iret = & |
---|
429 | nf_get_vara_real( ncid, id_var, istart, icount, dumarray2D_real) |
---|
430 | if (iret .ne. nf_noerr) then |
---|
431 | write(*,9100) 'error inquiring var value ' // varname, fnamenc |
---|
432 | ierr = -5 |
---|
433 | goto 100 |
---|
434 | end if |
---|
435 | else if (ndims.eq.3) then |
---|
436 | iret = & |
---|
437 | nf_get_vara_real( ncid, id_var, istart, icount, dumarray3D_real) |
---|
438 | if (iret .ne. nf_noerr) then |
---|
439 | write(*,9100) 'error inquiring var value ' // varname, fnamenc |
---|
440 | ierr = -5 |
---|
441 | if ((varname.eq.'TAUx').or.(varname.eq.'TAUY').or.(varname.eq.& |
---|
442 | 'SHFLX')) then |
---|
443 | print *,'************************************' |
---|
444 | print *,'* no stress available **************' |
---|
445 | print *,'************************************' |
---|
446 | cycle ! skyp the rest and go to the next variable |
---|
447 | end if |
---|
448 | goto 100 |
---|
449 | end if |
---|
450 | else if (ndims.eq.4) then |
---|
451 | iret = & |
---|
452 | nf_get_vara_real( ncid, id_var, istart, icount, dumarray4D_real) |
---|
453 | if (iret .ne. nf_noerr) then |
---|
454 | write(*,9100) 'error inquiring var value ' // varname, fnamenc |
---|
455 | ierr = -5 |
---|
456 | goto 100 |
---|
457 | end if |
---|
458 | end if |
---|
459 | !!if (idiagaa.eq.1) then |
---|
460 | !!write(75,*) varname |
---|
461 | !!end if |
---|
462 | |
---|
463 | do jy=0,nymin1 |
---|
464 | do ix=0, nxfield-1 |
---|
465 | if (varname.eq.'PS') then |
---|
466 | ps_tplus1_and_min1(ix,jy,index)=dumarray3D_real(ix+1,jy+1,itime)/time_interval |
---|
467 | end if |
---|
468 | end do |
---|
469 | end do |
---|
470 | |
---|
471 | end do |
---|
472 | |
---|
473 | iret = nf_close( ncid ) |
---|
474 | |
---|
475 | return |
---|
476 | |
---|
477 | |
---|
478 | |
---|
479 | 100 continue |
---|
480 | print *,'error reading for read_delta_ps_intime.f90',varname,'erros code',ierr,'or indj < 1:',indj |
---|
481 | print *,'delta_ps_in_time is set to be zero' |
---|
482 | do jy=0,nymin1 |
---|
483 | do ix=0,nxfield-1 |
---|
484 | ps_tplus1_and_min1(ix,jy,index)=0. |
---|
485 | end do |
---|
486 | end do |
---|
487 | |
---|
488 | |
---|
489 | |
---|
490 | |
---|
491 | return |
---|
492 | end subroutine read_delta_ps_intime |
---|
493 | |
---|
494 | |
---|