Subversion Repositories lagranto.wrf

Rev

Rev 2 | 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
110
      integer         ndims
111
      integer         is2d
112
      integer		  levid,lonid,latid
113
      integer      	  nx,ny,nz
114
      integer		  status,i,j,k
115
      character*2	  stag
116
 
117
c     ------------ Get grid info ---------------------------------------
118
      status = NF_INQ_DIMID(fid, 'south_north', LATID)
119
      status = NF_INQ_DIMID(fid, 'west_east',   LONID)
120
      status = NF_INQ_DIMID(fid, 'bottom_top',  LEVID)
121
      if (status .ne. NF_NOERR) print*,"Error in reading grid atts"
122
 
123
      status = NF_INQ_DIMLEN(fid, LONID, nx)
124
      status = NF_INQ_DIMLEN(fid, LATID, ny)
125
      status = NF_INQ_DIMLEN(fid, LEVID, nz)
126
      if (status .ne. NF_NOERR) print*,"Error in reading grid size"
127
 
128
c     ------------ Allocate memory and set staggering mode -------------
129
      if (trim(varname).eq.'U') then
130
          stag = 'X'
131
          is2d = 0
132
 
133
      elseif (trim(varname).eq.'V') then
134
          stag = 'Y'
135
          is2d = 0
136
 
137
      elseif (trim(varname).eq.'W') then
138
          stag = 'Z'
139
          is2d = 0
140
 
141
      elseif ( (trim(varname).eq.'geopot').or.
142
     >         (trim(varname).eq.'Z'     ) )
143
     >then
144
          stag = 'Z'
145
          is2d = 0
146
 
147
      elseif ( (trim(varname).eq.'geopots').or.
148
     >         (trim(varname).eq.'ZB'     ) )
149
     >then
150
          stag = 'Z'
151
          is2d = 1
152
 
153
      elseif ( (trim(varname).eq.'pressure').or.
154
     >         (trim(varname).eq.'P'       ) )
155
     >then
156
	      stag = 'nil'
157
          is2d = 0
158
 
159
	  else if (trim(varname).eq.'TH' .or. trim(varname).eq.'T') then
160
          stag = 'nil'
161
          is2d = 0
162
 
163
      else
164
         status = NF_INQ_VARID (fid, trim(varname), varid)
165
	     status = NF_GET_ATT_TEXT(fid,varid,'stagger',stag)
15 michaesp 166
	     if (status .ne. NF_NOERR) then
167
             print*,'Error in inq:',trim(varname)
168
         endif
2 michaesp 169
         status = NF_INQ_VARNDIMS (fid, varid, ndims)
170
         is2d = 0
171
         if ( ndims.eq.3 ) is2d = 1
172
 
173
      endif
174
 
175
      if ( stag.eq.'X' ) then
176
          allocate( temp (nx+1,ny,nz) )
177
          allocate( temp1(nx+1,ny,nz) )
178
          allocate( temp2(nx+1,ny,nz) )
179
 
180
      elseif ( stag.eq.'Y' ) then
181
          allocate( temp (nx,ny+1,nz) )
182
          allocate( temp1(nx,ny+1,nz) )
183
          allocate( temp2(nx,ny+1,nz) )
184
 
185
      elseif ( stag.eq.'Z' ) then
186
          allocate( temp (nx,ny,nz+1) )
187
          allocate( temp1(nx,ny,nz+1) )
188
          allocate( temp2(nx,ny,nz+1) )
189
 
190
      else
191
          allocate( temp (nx,ny,nz) )
192
          allocate( temp1(nx,ny,nz) )
193
          allocate( temp2(nx,ny,nz) )
194
 
195
      endif
196
 
197
c	  ------------ Read data ------------------------------------------
198
 
199
c	  Zonal wind : temp(nx+1,ny,nk)
200
	  if (trim(varname).eq.'U') then
201
 
202
	     status = NF_INQ_VARID (fid, 'U', varid)
203
	     if (status .ne. NF_NOERR) print*,"Error in inq U"
204
	     status = NF_GET_VAR_REAL (fid, varid, temp)
205
	     if (status .ne. NF_NOERR) print*,"Error in reading U"
206
 
207
c     Meridional wind : temp(nx,ny+1,nk)
208
	  else if (trim(varname).eq.'V') then
209
 
210
	     status = NF_INQ_VARID (fid, 'V', varid)
211
	     if (status .ne. NF_NOERR) print*,"Error in inq V"
212
	     status = NF_GET_VAR_REAL (fid, varid, temp)
213
	     if (status .ne. NF_NOERR) print*,"Error in reading V"
214
 
215
c	  Vertical wind : temp(nx,ny,nz+1)
216
	  else if (trim(varname).eq.'W') then
217
 
218
	     status = NF_INQ_VARID (fid, 'W', varid)
219
	     if (status .ne. NF_NOERR) print*,"Error in inq W"
220
	     status = NF_GET_VAR_REAL (fid, varid, temp)
221
	     if (status .ne. NF_NOERR) print*,"Error in reading U"
222
 
223
c	  Geopotential height (base + perturbation) : temp(nx,ny,nz+1)
224
	  else if ( (trim(varname).eq.'geopot').or.
225
     >          (trim(varname).eq.'Z'     ) )
226
     >then
227
 
228
	     status = NF_INQ_VARID (fid, 'PHB', varid)
229
	     if (status .ne. NF_NOERR) print*,"Error in inq geopot"
230
	     status = NF_GET_VAR_REAL (fid, varid, temp1)
231
	     if (status .ne. NF_NOERR) print*,"Error in reading geopot"
232
 
233
	     status = NF_INQ_VARID (fid, 'PH', varid)
234
	     if (status .ne. NF_NOERR) print*,"Error in inq pgeopot"
235
	     status = NF_GET_VAR_REAL (fid, varid, temp2)
236
	     if (status .ne. NF_NOERR) print*,"Error in reading pgeopot"
237
 
238
	     do k = 1, nz+1
239
	       do j = 1, ny
240
	         do i = 1, nx
241
	           temp(i,j,k) = temp1(i,j,k) + temp2(i,j,k)
242
	         enddo
243
	       enddo
244
	     enddo
245
 
246
c	  surface geopotential: temp(nx,ny,nz+1)
247
	  else if ( (trim(varname).eq.'geopots').or.
248
     >          (trim(varname).eq.'ZB'     ) )
249
     >then
250
 
251
	     status = NF_INQ_VARID (fid, 'PHB', varid)
252
	     if (status .ne. NF_NOERR) print*,"Error in inq sgeopot"
253
	     status = NF_GET_VAR_REAL (fid, varid, temp1)
254
	     if (status .ne. NF_NOERR) print*,"Error in reading sgeopot"
255
 
256
	     status = NF_INQ_VARID (fid, 'PH', varid)
257
	     if (status .ne. NF_NOERR) print*,"Error in inq pgeopot"
258
	     status = NF_GET_VAR_REAL (fid, varid, temp2)
259
	     if (status .ne. NF_NOERR) print*,"Error in reading pgeopot"
260
 
261
	     do k = 1, nz+1
262
	       do j = 1, ny
263
	         do i = 1, nx
264
	           temp(i,j,k) = temp1(i,j,k) + temp2(i,j,k)
265
	         enddo
266
	       enddo
267
	     enddo
268
 
269
c	  Pressure (base + perturbation) : temp(nx,ny,nz)
270
	  elseif ( (trim(varname).eq.'pressure').or.
271
     >         (trim(varname).eq.'P'       ) )
272
     >then
273
 
274
	      status = NF_INQ_VARID (fid, 'PB', varid)
275
	      if (status .ne. NF_NOERR) print*,"Error in inq pb"
276
	      status = NF_GET_VAR_REAL (fid, varid, temp1)
277
	      if (status .ne. NF_NOERR) print*,"Error in reading pb"
278
 
279
       	  status = NF_INQ_VARID (fid, 'P', varid)
280
	      if (status .ne. NF_NOERR) print*,"Error in inq p"
281
	      status = NF_GET_VAR_REAL (fid, varid, temp2)
282
	      if (status .ne. NF_NOERR) print*,"Error in reading p"
283
 
284
	      do k = 1, nz
285
	        do j = 1, ny
286
	          do i = 1, nx
287
	            temp(i,j,k) = temp1(i,j,k) + temp2(i,j,k)
288
	          enddo
289
	        enddo
290
	      enddo
291
 
292
c     Potential temperature: temp(nx,ny,nz)
293
	  else if (trim(varname).eq.'TH' .or. trim(varname).eq.'T') then
294
 
295
	      status = NF_INQ_VARID (fid, 'T', varid)
296
	      if (status .ne. NF_NOERR) print*,"Error in inq T"
297
	      status = NF_GET_VAR_REAL (fid, varid, temp)
298
	      if (status .ne. NF_NOERR) print*,"Error in reading T"
299
 
300
c	  Any other field (2D + 3D): temp(:,:,:), depending on staggering
301
	  else
302
 
303
	      status = NF_INQ_VARID (fid, trim(varname), varid)
304
	      if (status .ne. NF_NOERR) then
305
	         print*,"Error in inq:",trim(varname)
306
	 	  endif
307
	 	  status = NF_GET_VAR_REAL (fid, varid, temp)
308
	 	  if (status .ne. NF_NOERR) then
309
	 	     print*,"Error in reading:",trim(varname)
310
          endif
311
 
312
      endif
313
 
314
c     ------------ Destaggering in X, Y and Z direction ---------------
315
      if (trim(stag).eq.'X') then
316
         do k = 1, n3
317
	  	 do j = 1, n2
318
	   	 do i = 1, n1
319
	    	 field(i,j,k) = 0.5*(temp(i,j,k) + temp(i+1,j,k))
320
	   	 enddo
321
	  	 enddo
322
	 	 enddo
323
 
324
	  elseif (trim(stag).eq.'Y') then
325
	 	 do k = 1, n3
326
	  	 do j = 1, n2
327
	   	 do i = 1, n1
328
	    	 field(i,j,k) = 0.5*(temp(i,j,k) + temp(i,j+1,k))
329
	   	 enddo
330
	  	 enddo
331
	 	 enddo
332
 
333
	  elseif (trim(stag).eq.'Z') then
334
	     do k = 1, n3
335
	  	 do j = 1, n2
336
	   	 do i = 1, n1
337
	    	 field(i,j,k) = 0.5*(temp(i,j,k) + temp(i,j,k+1))
338
	   	 enddo
339
	  	 enddo
340
	 	 enddo
341
 
342
	  else
343
	  	 do k = 1, n3
344
	  	 do j = 1, n2
345
	   	 do i = 1, n1
346
	    	 field(i,j,k) = temp(i,j,k)
347
	   	 enddo
348
	  	 enddo
349
	 	 enddo
350
 
351
	  endif
352
 
353
c     ---------- Change units -----------------------------------------
354
 
355
c	  Pressure Pa -> hPa
356
	  if (trim(varname).eq.'pressure' .or. trim(varname).eq.'P') then
357
	    do k = 1, n3
358
	      do j = 1, n2
359
	        do i = 1, n1
360
	          field(i,j,k)=0.01 * field(i,j,k)
361
	        enddo
362
	      enddo
363
	    enddo
364
 
365
c     Potential temperature + 300
366
      else if (trim(varname).eq.'TH' .or. trim(varname).eq.'T') then
367
	    do k = 1, n3
368
	      do j = 1, n2
369
	        do i = 1, n1
370
	          field(i,j,k)=field(i,j,k) + addtheta
371
	        enddo
372
	      enddo
373
	    enddo
374
 
375
c     Geopotential -> Geopotential Height (3D + Surface)
376
	  else if ( (trim(varname).eq.'geopot' ).or.
377
     >          (trim(varname).eq.'Z'      ).or.
378
     >          (trim(varname).eq.'geopots').or.
379
     >          (trim(varname).eq.'ZB'     ) )
380
     >then
381
     	do k = 1, n3
382
	      do j = 1, n2
383
	        do i = 1, n1
384
	          field(i,j,k)=field(i,j,k)/gearth
385
	        enddo
386
	      enddo
387
	    enddo
388
 
389
	  endif
390
 
391
c     ---------- Copy lowest level for 2d tracing ---------------------
392
      if ( is2d.eq.1 ) then
393
         do i=1,n1
394
           do j=1,n2
395
             do k=1,n3
396
                field(i,j,k) = field(i,j,1)
397
             enddo
398
           enddo
399
         enddo
400
      endif
401
 
402
 
403
	  end
404
 
405
c     ------------------------------------------------------------
406
c     Close input file
407
c     ------------------------------------------------------------
408
 
409
      subroutine input_close(fid)
410
 
411
c     Close the input file with file identifier <fid>.
412
 
413
      implicit none
414
 
415
      include 'netcdf.inc'
416
 
417
c     Declaration of subroutine parameters
418
      integer fid
419
 
420
c     Auxiliary variables
421
      integer ierr
422
 
423
c     Close file
424
      call ncclos(fid,ierr)
425
 
426
      end
427
 
428
c     ------------------------------------------------------------
429
c     ------------------------------------------------------------