Subversion Repositories lagranto.wrf

Rev

Rev 15 | Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
2 michaesp 1
c     ************************************************************
2
c     * This package provides input routines to read if grid     *
3
c     * informations or field from WRF output			 *
4
c     *		                                                 *
5
c     * Author: Sebastian Schemm, ETH, 2013 			 *
6
c     * V.1.0 only support for ideal channel implemented         *
7
c     *      for input_grid_wrf				         *
8
c     *								 *
9
c     ************************************************************
10
 
11
c     ------------------------------------------------------------
12
c     Open input file
13
c     ------------------------------------------------------------
14
 
15
      subroutine input_open (fid,filename)
16
 
17
c     Open the input file with filename <filename> and return the
18
c     file identifier <fid> for further reference. 
19
 
20
      implicit none
21
 
22
      include 'netcdf.inc'
23
 
24
c     Declaration of subroutine parameters
25
      integer      fid              ! File identifier
26
      character*80 filename         ! Filename
27
 
28
c     Declaration of auxiliary variables
29
      integer      ierr
30
 
31
c     Open IVE netcdf file
32
      fid = ncopn(filename,NCWRITE,ierr)
33
      if (ierr.ne.0) goto 900
34
 
35
c     Exception handling
36
      return
37
 
38
 900  print*,'Cannot open input file ',trim(filename)
39
      stop
40
 
41
      end
42
 
43
c     ------------------------------------------------------------
44
c     Read information about the grid
45
c     ------------------------------------------------------------
46
 
47
	  subroutine input_grid_wrf(fid,
48
     >				  xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz)
49
 
50
  	  implicit none
51
 
52
      include 'netcdf.inc'
53
 
54
c     Declaration of subroutine parameters
55
      integer      fid                 ! File identifier
56
	  integer	   latid,lonid,levid
57
      real         xmin,xmax,ymin,ymax ! Domain size
58
      real         dx,dy               ! Horizontal resolution
59
      integer      nx,ny,nz            ! Grid dimensions
60
	  integer	   status
61
 
62
 
63
C	  get grid info
64
	  status = NF_INQ_DIMID(fid, 'south_north', LATID)
65
	  status = NF_INQ_DIMID(fid, 'west_east',   LONID)
66
	  status = NF_INQ_DIMID(fid, 'bottom_top',  LEVID)
67
	  if (status .ne. NF_NOERR) print*,"Error in reading grid atts"
68
 
69
	  status = NF_INQ_DIMLEN(fid, LONID, nx)
70
	  status = NF_INQ_DIMLEN(fid, LATID, ny)
71
	  status = NF_INQ_DIMLEN(fid, LEVID, nz)
72
	  if (status .ne. NF_NOERR) print*,"Error in reading grid size"
73
 
74
	  status = NF_GET_ATT_REAL (fid, NF_GLOBAL, 'DX', dx)
75
	  status = NF_GET_ATT_REAL (fid, NF_GLOBAL, 'DY', dy)
76
	  if (status .ne. NF_NOERR) print*,"Error in reading dx, dy"
77
 
78
c	  setup a pseudo grid
79
      xmin = 0.0
80
      ymin = 0.0
81
	  xmax = xmin+real(nx-1)*dx
82
      ymax = ymin+real(ny-1)*dy
83
 
84
	  end
85
 
86
c     ------------------------------------------------------------
87
c     Read variables
88
c     ------------------------------------------------------------
89
 
90
      subroutine input_var_wrf(fid,varname,field,n1,n2,n3)
91
 
92
      implicit none
93
 
94
      include 'netcdf.inc'
95
 
96
c     Declaration of subroutine parameters
97
      integer		  fid, varid
98
      integer      	  n1,n2,n3
99
      character*80	  varname
100
      real            field(n1,n2,n3)
101
 
102
c     Parmeters
103
      real            gearth         ! Earth's acceleration
104
      parameter       (gearth=9.81)
105
      real            addtheta       ! Offset for potential temperature
106
      parameter       (addtheta=300.)
107
 
108
c     Auxiliary variables
109
      real,allocatable,dimension (:,:,:)  :: temp,temp1,temp2
21 michaesp 110
      real,allocatable,dimension (:,:  )  :: temp3
2 michaesp 111
      integer         ndims
112
      integer         is2d
113
      integer		  levid,lonid,latid
114
      integer      	  nx,ny,nz
115
      integer		  status,i,j,k
116
      character*2	  stag
117
 
118
c     ------------ Get grid info ---------------------------------------
119
      status = NF_INQ_DIMID(fid, 'south_north', LATID)
120
      status = NF_INQ_DIMID(fid, 'west_east',   LONID)
121
      status = NF_INQ_DIMID(fid, 'bottom_top',  LEVID)
122
      if (status .ne. NF_NOERR) print*,"Error in reading grid atts"
123
 
124
      status = NF_INQ_DIMLEN(fid, LONID, nx)
125
      status = NF_INQ_DIMLEN(fid, LATID, ny)
126
      status = NF_INQ_DIMLEN(fid, LEVID, nz)
127
      if (status .ne. NF_NOERR) print*,"Error in reading grid size"
128
 
129
c     ------------ Allocate memory and set staggering mode -------------
130
      if (trim(varname).eq.'U') then
131
          stag = 'X'
132
          is2d = 0
133
 
134
      elseif (trim(varname).eq.'V') then
135
          stag = 'Y'
136
          is2d = 0
137
 
138
      elseif (trim(varname).eq.'W') then
139
          stag = 'Z'
140
          is2d = 0
141
 
142
      elseif ( (trim(varname).eq.'geopot').or.
143
     >         (trim(varname).eq.'Z'     ) )
144
     >then
145
          stag = 'Z'
146
          is2d = 0
147
 
148
      elseif ( (trim(varname).eq.'geopots').or.
149
     >         (trim(varname).eq.'ZB'     ) )
150
     >then
21 michaesp 151
          stag = 'nil'
2 michaesp 152
          is2d = 1
153
 
154
      elseif ( (trim(varname).eq.'pressure').or.
155
     >         (trim(varname).eq.'P'       ) )
156
     >then
157
	      stag = 'nil'
158
          is2d = 0
159
 
160
	  else if (trim(varname).eq.'TH' .or. trim(varname).eq.'T') then
161
          stag = 'nil'
162
          is2d = 0
163
 
164
      else
165
         status = NF_INQ_VARID (fid, trim(varname), varid)
166
	     status = NF_GET_ATT_TEXT(fid,varid,'stagger',stag)
15 michaesp 167
	     if (status .ne. NF_NOERR) then
168
             print*,'Error in inq:',trim(varname)
169
         endif
2 michaesp 170
         status = NF_INQ_VARNDIMS (fid, varid, ndims)
171
         is2d = 0
172
         if ( ndims.eq.3 ) is2d = 1
173
 
174
      endif
175
 
176
      if ( stag.eq.'X' ) then
177
          allocate( temp (nx+1,ny,nz) )
178
          allocate( temp1(nx+1,ny,nz) )
179
          allocate( temp2(nx+1,ny,nz) )
180
 
181
      elseif ( stag.eq.'Y' ) then
182
          allocate( temp (nx,ny+1,nz) )
183
          allocate( temp1(nx,ny+1,nz) )
184
          allocate( temp2(nx,ny+1,nz) )
185
 
186
      elseif ( stag.eq.'Z' ) then
187
          allocate( temp (nx,ny,nz+1) )
188
          allocate( temp1(nx,ny,nz+1) )
189
          allocate( temp2(nx,ny,nz+1) )
190
      else
191
          allocate( temp (nx,ny,nz) )
192
          allocate( temp1(nx,ny,nz) )
193
          allocate( temp2(nx,ny,nz) )
21 michaesp 194
          allocate( temp3(nx,ny   ) )
2 michaesp 195
 
196
      endif
197
 
198
c	  ------------ Read data ------------------------------------------
199
 
200
c	  Zonal wind : temp(nx+1,ny,nk)
201
	  if (trim(varname).eq.'U') then
202
 
203
	     status = NF_INQ_VARID (fid, 'U', varid)
204
	     if (status .ne. NF_NOERR) print*,"Error in inq U"
205
	     status = NF_GET_VAR_REAL (fid, varid, temp)
206
	     if (status .ne. NF_NOERR) print*,"Error in reading U"
207
 
208
c     Meridional wind : temp(nx,ny+1,nk)
209
	  else if (trim(varname).eq.'V') then
210
 
211
	     status = NF_INQ_VARID (fid, 'V', varid)
212
	     if (status .ne. NF_NOERR) print*,"Error in inq V"
213
	     status = NF_GET_VAR_REAL (fid, varid, temp)
214
	     if (status .ne. NF_NOERR) print*,"Error in reading V"
215
 
216
c	  Vertical wind : temp(nx,ny,nz+1)
217
	  else if (trim(varname).eq.'W') then
218
 
219
	     status = NF_INQ_VARID (fid, 'W', varid)
220
	     if (status .ne. NF_NOERR) print*,"Error in inq W"
221
	     status = NF_GET_VAR_REAL (fid, varid, temp)
21 michaesp 222
	     if (status .ne. NF_NOERR) print*,"Error in reading W"
2 michaesp 223
 
224
c	  Geopotential height (base + perturbation) : temp(nx,ny,nz+1)
225
	  else if ( (trim(varname).eq.'geopot').or.
226
     >          (trim(varname).eq.'Z'     ) )
227
     >then
228
 
229
	     status = NF_INQ_VARID (fid, 'PHB', varid)
230
	     if (status .ne. NF_NOERR) print*,"Error in inq geopot"
231
	     status = NF_GET_VAR_REAL (fid, varid, temp1)
232
	     if (status .ne. NF_NOERR) print*,"Error in reading geopot"
233
 
234
	     status = NF_INQ_VARID (fid, 'PH', varid)
235
	     if (status .ne. NF_NOERR) print*,"Error in inq pgeopot"
236
	     status = NF_GET_VAR_REAL (fid, varid, temp2)
237
	     if (status .ne. NF_NOERR) print*,"Error in reading pgeopot"
238
 
239
	     do k = 1, nz+1
240
	       do j = 1, ny
241
	         do i = 1, nx
242
	           temp(i,j,k) = temp1(i,j,k) + temp2(i,j,k)
243
	         enddo
244
	       enddo
245
	     enddo
246
 
247
c	  surface geopotential: temp(nx,ny,nz+1)
248
	  else if ( (trim(varname).eq.'geopots').or.
249
     >          (trim(varname).eq.'ZB'     ) )
250
     >then
251
 
21 michaesp 252
c	     status = NF_INQ_VARID (fid, 'PHB', varid)
253
c	     if (status .ne. NF_NOERR) print*,"Error in inq sgeopot"
254
c	     status = NF_GET_VAR_REAL (fid, varid, temp1)
255
c	     if (status .ne. NF_NOERR) print*,"Error in reading sgeopot"
256
c
257
c	     status = NF_INQ_VARID (fid, 'PH', varid)
258
c	     if (status .ne. NF_NOERR) print*,"Error in inq pgeopot"
259
c	     status = NF_GET_VAR_REAL (fid, varid, temp2)
260
c	     if (status .ne. NF_NOERR) print*,"Error in reading pgeopot"
261
c
262
c	     do k = 1, nz+1
263
c	       do j = 1, ny
264
c	         do i = 1, nx
265
c	           temp(i,j,k) = temp1(i,j,k) + temp2(i,j,k)
266
c	         enddo
267
c	       enddo
268
c	     enddo
2 michaesp 269
 
21 michaesp 270
	     status = NF_INQ_VARID (fid, 'HGT', varid)
271
	     if (status .ne. NF_NOERR) print*,"Error in inq HGT"
272
	     status = NF_GET_VAR_REAL (fid, varid, temp3)
273
	     if (status .ne. NF_NOERR) print*,"Error in reading HGT"
274
	     do k = 1, nz
2 michaesp 275
	       do j = 1, ny
276
	         do i = 1, nx
21 michaesp 277
	           temp(i,j,k) = temp3(i,j) * gearth
2 michaesp 278
	         enddo
279
	       enddo
280
	     enddo
281
 
282
c	  Pressure (base + perturbation) : temp(nx,ny,nz)
283
	  elseif ( (trim(varname).eq.'pressure').or.
284
     >         (trim(varname).eq.'P'       ) )
285
     >then
286
 
287
	      status = NF_INQ_VARID (fid, 'PB', varid)
288
	      if (status .ne. NF_NOERR) print*,"Error in inq pb"
289
	      status = NF_GET_VAR_REAL (fid, varid, temp1)
290
	      if (status .ne. NF_NOERR) print*,"Error in reading pb"
291
 
292
       	  status = NF_INQ_VARID (fid, 'P', varid)
293
	      if (status .ne. NF_NOERR) print*,"Error in inq p"
294
	      status = NF_GET_VAR_REAL (fid, varid, temp2)
295
	      if (status .ne. NF_NOERR) print*,"Error in reading p"
296
 
297
	      do k = 1, nz
298
	        do j = 1, ny
299
	          do i = 1, nx
300
	            temp(i,j,k) = temp1(i,j,k) + temp2(i,j,k)
301
	          enddo
302
	        enddo
303
	      enddo
304
 
305
c     Potential temperature: temp(nx,ny,nz)
306
	  else if (trim(varname).eq.'TH' .or. trim(varname).eq.'T') then
307
 
308
	      status = NF_INQ_VARID (fid, 'T', varid)
309
	      if (status .ne. NF_NOERR) print*,"Error in inq T"
310
	      status = NF_GET_VAR_REAL (fid, varid, temp)
311
	      if (status .ne. NF_NOERR) print*,"Error in reading T"
312
 
313
c	  Any other field (2D + 3D): temp(:,:,:), depending on staggering
314
	  else
315
 
316
	      status = NF_INQ_VARID (fid, trim(varname), varid)
317
	      if (status .ne. NF_NOERR) then
318
	         print*,"Error in inq:",trim(varname)
319
	 	  endif
320
	 	  status = NF_GET_VAR_REAL (fid, varid, temp)
321
	 	  if (status .ne. NF_NOERR) then
322
	 	     print*,"Error in reading:",trim(varname)
323
          endif
324
 
325
      endif
326
 
327
c     ------------ Destaggering in X, Y and Z direction ---------------
328
      if (trim(stag).eq.'X') then
329
         do k = 1, n3
330
	  	 do j = 1, n2
331
	   	 do i = 1, n1
332
	    	 field(i,j,k) = 0.5*(temp(i,j,k) + temp(i+1,j,k))
333
	   	 enddo
334
	  	 enddo
335
	 	 enddo
336
 
337
	  elseif (trim(stag).eq.'Y') then
338
	 	 do k = 1, n3
339
	  	 do j = 1, n2
340
	   	 do i = 1, n1
341
	    	 field(i,j,k) = 0.5*(temp(i,j,k) + temp(i,j+1,k))
342
	   	 enddo
343
	  	 enddo
344
	 	 enddo
345
 
346
	  elseif (trim(stag).eq.'Z') then
347
	     do k = 1, n3
348
	  	 do j = 1, n2
349
	   	 do i = 1, n1
350
	    	 field(i,j,k) = 0.5*(temp(i,j,k) + temp(i,j,k+1))
351
	   	 enddo
352
	  	 enddo
353
	 	 enddo
354
 
355
	  else
356
	  	 do k = 1, n3
357
	  	 do j = 1, n2
358
	   	 do i = 1, n1
359
	    	 field(i,j,k) = temp(i,j,k)
360
	   	 enddo
361
	  	 enddo
362
	 	 enddo
363
 
364
	  endif
365
 
366
c     ---------- Change units -----------------------------------------
367
 
368
c	  Pressure Pa -> hPa
369
	  if (trim(varname).eq.'pressure' .or. trim(varname).eq.'P') then
370
	    do k = 1, n3
371
	      do j = 1, n2
372
	        do i = 1, n1
373
	          field(i,j,k)=0.01 * field(i,j,k)
374
	        enddo
375
	      enddo
376
	    enddo
377
 
378
c     Potential temperature + 300
379
      else if (trim(varname).eq.'TH' .or. trim(varname).eq.'T') then
380
	    do k = 1, n3
381
	      do j = 1, n2
382
	        do i = 1, n1
383
	          field(i,j,k)=field(i,j,k) + addtheta
384
	        enddo
385
	      enddo
386
	    enddo
387
 
388
c     Geopotential -> Geopotential Height (3D + Surface)
389
	  else if ( (trim(varname).eq.'geopot' ).or.
390
     >          (trim(varname).eq.'Z'      ).or.
391
     >          (trim(varname).eq.'geopots').or.
392
     >          (trim(varname).eq.'ZB'     ) )
393
     >then
394
     	do k = 1, n3
395
	      do j = 1, n2
396
	        do i = 1, n1
397
	          field(i,j,k)=field(i,j,k)/gearth
398
	        enddo
399
	      enddo
400
	    enddo
401
 
402
	  endif
403
 
404
c     ---------- Copy lowest level for 2d tracing ---------------------
405
      if ( is2d.eq.1 ) then
406
         do i=1,n1
407
           do j=1,n2
408
             do k=1,n3
409
                field(i,j,k) = field(i,j,1)
410
             enddo
411
           enddo
412
         enddo
413
      endif
414
 
415
 
416
	  end
417
 
418
c     ------------------------------------------------------------
419
c     Close input file
420
c     ------------------------------------------------------------
421
 
422
      subroutine input_close(fid)
423
 
424
c     Close the input file with file identifier <fid>.
425
 
426
      implicit none
427
 
428
      include 'netcdf.inc'
429
 
430
c     Declaration of subroutine parameters
431
      integer fid
432
 
433
c     Auxiliary variables
434
      integer ierr
435
 
436
c     Close file
437
      call ncclos(fid,ierr)
438
 
439
      end
440
 
441
c     ------------------------------------------------------------