Subversion Repositories lagranto.ecmwf

Rev

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

Rev Author Line No. Line
3 michaesp 1
      PROGRAM create_startf
2
 
3
c     **************************************************************
4
c     * Create a <startfile> for <lagrangto>. It can be chosen     *
5
c     * whether to start from an isentropic or an isobaric         *
6
c     * surface. The starting points are equidistantly distributed *
7
c     * Michael Sprenger / Autumn 2004                             *
8
c     **************************************************************
9
 
19 michaesp 10
      use netcdf
11
 
3 michaesp 12
      implicit none
13
 
14
c     --------------------------------------------------------------
15
c     Set parameters
16
c     --------------------------------------------------------------
17
 
18
c     Maximum number of starting positions
19
      integer           nmax
13 michaesp 20
      parameter         (nmax=2000000)
3 michaesp 21
 
22
c     Maximum number of model levels
23
      integer           nlevmax
24
      parameter         (nlevmax=200)
25
 
26
c     Grid constant (distance in km corresponding to 1 deg at the equator)
27
      real              deltat
28
      parameter         (deltat=111.)
29
 
30
c     Mathematical constant (conversion degree -> radian)
31
      real              pi180
32
      parameter         (pi180=3.14159/180.)
33
 
34
c     Numerical epsilon
35
      real              eps
36
      parameter         (eps=0.00001)
37
 
38
c     --------------------------------------------------------------
39
c     Set variables
40
c     --------------------------------------------------------------
41
 
42
c     Filenames and output format
43
      character*80      pfile0,pfile1                  ! P filenames
44
      character*80      sfile0,sfile1                  ! S filenames
45
      character*80      ofile                          ! Output filename
46
      integer           oformat                        ! Output format
47
      real              timeshift                      ! Time shift relative to data files <*0>
48
      real              timeinc                        ! Time increment between input files
49
 
50
c     Horizontal grid
51
      character*80      hmode                          ! Horizontale mode
52
      real              lat1,lat2,lon1,lon2            ! Lat/lon boundaries
53
      real              ds,dlon,dlat                   ! Distance and lat/lon shifts
54
      character*80      hfile                          ! Filename
55
      integer           hn                             ! Number of entries in lat/lon list 
56
      real              latlist(nmax)                  ! List of latitudes
57
      real              lonlist(nmax)                  ! List of longitudes
58
      integer           pn                             ! Number of entries in lat/lon poly
59
      real              latpoly(500)                   ! List of polygon latitudes
60
      real              lonpoly(500)                   ! List of polygon longitudes
61
      real              loninpoly,latinpoly            ! Lon/lat inside polygon
62
      character*80      regionf                        ! Region file
63
      integer           iregion                        ! Region number
64
      real              xcorner(4),ycorner(4)          ! Vertices of region
65
 
66
c     Vertical grid
67
      character*80      vmode                          ! Vertical mode
68
      real              lev1,lev2,levlist(nmax)        ! Single levels, and list of levels
69
      character*80      vfile                          ! Filename
70
      integer           vn                             ! Number of entries
71
 
72
c     Unit of vertical axis
73
      character*80      umode                          ! Unit of vertical axis
74
 
75
c     Flag for 'no time check'
76
      character*80      timecheck                      ! Either 'no' or 'yes'
77
 
78
c     List of all starting positions
79
      integer           start_n                        ! Number of coordinates
80
      real              start_lat(nmax)                ! Latitudes
81
      real              start_lon(nmax)                ! Longitudes
82
      real              start_lev(nmax)                ! Levels (depending on vertical unit)
83
      real              start_pre(nmax)                ! Level in hPa
84
      integer           reftime(6)                     ! Reference time
85
      character*80      vars(10)                       ! Name of output fields (time,lon,lat,p)
86
      real,allocatable, dimension (:,:,:) :: tra       ! Trajectories (ntra,ntim,ncol)
87
      real              latmin,latmax
88
      real              lonmin,lonmax
89
      real              premin,premax
90
 
91
c     Grid description
92
      real              pollon,pollat                  ! Longitude/latitude of pole
93
      real              ak(nlevmax)                    ! Vertical layers and levels
94
      real              bk(nlevmax)
95
      real              xmin,xmax                      ! Zonal grid extension
96
      real              ymin,ymax                      ! Meridional grid extension
97
      integer           nx,ny,nz                       ! Grid dimensions
98
      real              dx,dy                          ! Horizontal grid resolution
99
      real,allocatable, dimension (:,:,:) :: pr        ! 3d pressure
100
      real,allocatable, dimension (:,:)   :: prs       ! surface pressure
101
      real,allocatable, dimension (:,:,:) :: th        ! 3d potential temperature
102
      real,allocatable, dimension (:,:)   :: ths       ! surface poential temperature
103
      real,allocatable, dimension (:,:,:) :: pv        ! 3d potential vorticity
104
      real,allocatable, dimension (:,:)   :: pvs       ! surface potential vorticiy
105
      real,allocatable, dimension (:,:,:) :: in        ! 3d 'dummy' array with vertical indices
106
      character*80      varname                        ! Name of input variable
107
      integer           fid                            ! File identifier
108
      real              stagz                          ! Vertical staggering
109
      real              mdv                            ! Missing data values
110
      real              tstart,tend                    ! Time on P and S file
111
      real              rid,rjd,rkd                    ! Real grid position
112
 
113
c     Auxiliary variable
114
      integer           i,j,k
115
      real              lon,lat
116
      real              rd
117
      integer           stat,flag
118
      real              tmp1,tmp2
119
      real              tfrac,frac
120
      real              radius,dist
121
      character*80      string
122
      character*80      selectstr
123
      character*80      umode_save
19 michaesp 124
      character*80      maskname
125
      real              maskvalue
126
      character*2       maskoper
127
      integer           cdfid,varid,ierr
128
      integer           nmask
129
      integer           indx,indy
130
      character         ch1,ch2
131
      integer           len
3 michaesp 132
      real,allocatable, dimension (:,:,:) :: fld0
133
      real,allocatable, dimension (:,:,:) :: fld1
134
      real,allocatable, dimension (:,:  ) :: sfc0
135
      real,allocatable, dimension (:,:)   :: sfc1
19 michaesp 136
      real,allocatable, dimension (:,:)   :: mask
3 michaesp 137
 
138
c     Externals 
139
      real              int_index3     ! 3d interpolation
140
      external          int_index3
141
      real              sdis           ! Speherical distance
142
      external          sdis
143
      integer           inregion       ! In/out of region
144
      external          inrehion
145
 
146
c     ------------------------------------------------------------------
147
c     Start of program, Read parameters
148
c     ------------------------------------------------------------------
149
 
150
c     Write start message
151
      print*,'========================================================='
152
      print*,'         *** START OF PROGRAM CREATE_STARTF ***'
153
      print*
154
 
155
c     Read parameter file
156
      open(10,file='create_startf.param')
157
 
158
c      Input P and S file
159
       read(10,*) pfile0,pfile1
160
       read(10,*) sfile0,sfile1
161
       read(10,*) ofile
162
 
163
c      Read name of region file
164
       read(10,*) regionf
165
 
166
c      Reference time
167
       do i=1,6
168
          read(10,*) reftime(i)
169
       enddo
170
 
171
c      Time shift relative to data files <pfile0,sfile0> - format (hh.mm)
172
       read(10,*) timeshift
173
 
174
c      Read timeincrement between input files
175
       read(10,*) timeinc
176
 
177
c      Parameters for horizontal grid
178
       read(10,*) hmode
179
       if ( hmode.eq.'file' ) then                  ! from file
180
          read(10,*) hfile
181
       elseif ( hmode.eq.'line' ) then              ! along a line
182
          read(10,*) lon1,lon2,lat1,lat2,hn
183
       elseif ( hmode.eq.'box.eqd' ) then           ! box: 2d equidistant
184
          read(10,*) lon1,lon2,lat1,lat2,ds
185
       elseif ( hmode.eq.'box.grid' ) then          ! box: 2d grid
186
          read(10,*) lon1,lon2,lat1,lat2
187
       elseif ( hmode.eq.'point' ) then             ! single point
188
          read(10,*) lon1,lat1
189
       elseif ( hmode.eq.'shift' ) then             ! centre + shifted
190
          read(10,*) lon1,lat1,dlon,dlat
191
       elseif ( hmode.eq.'polygon.eqd' ) then       ! polygon: 2d equidistant
192
          read(10,*) hfile,ds
193
       elseif ( hmode.eq.'polygon.grid' ) then      ! polygon: 2d grid
194
          read(10,*) hfile
195
       elseif ( hmode.eq.'circle.eqd' ) then        ! circle: 2d equidistant
196
          read(10,*) lon1,lat1,radius,ds 
197
       elseif ( hmode.eq.'circle.grid' ) then       ! circle: 2d grid
198
          read(10,*) lon1,lat1,radius
199
       elseif ( hmode.eq.'region.eqd' ) then        ! region: 2d equidistant
200
          read(10,*) iregion,ds 
201
       elseif ( hmode.eq.'region.grid' ) then       ! iregion: 2d grid
202
          read(10,*) iregion
19 michaesp 203
       elseif ( hmode.eq.'mask.grid' ) then         ! 0/1 masked - grid points
204
          read(10,*) hfile,maskname
205
       elseif ( hmode.eq.'mask.eqd' ) then         ! 0/1 masked - euqidistant points
206
          read(10,*) hfile,maskname,ds
3 michaesp 207
       else
208
          print*,' ERROR: horizontal mode not supported ',trim(hmode)
209
          stop
210
       endif
211
 
212
c      Parameters for vertical grid
213
       read(10,*) vmode
214
       if ( vmode.eq.'file') then                   ! from file
215
          read(10,*) vfile
216
       elseif ( vmode.eq.'level' ) then             ! single level (explicit command)
217
          read(10,*) lev1                         
218
       elseif ( vmode.eq.'list') then               ! a list
219
          read(10,*) vn
220
          read(10,*) (levlist(i),i=1,vn)
221
       elseif ( vmode.eq.'profile') then            ! a profile
222
          read(10,*) lev1,lev2,vn
223
       elseif ( vmode.eq.'grid') then               ! grid points
224
          read(10,*) lev1,lev2
225
       else
226
          print*,' ERROR: vertical mode not supported ',trim(vmode)
227
          stop
228
       endif
229
 
230
c      Read units of vertical axis
231
       read(10,*) umode
232
       if ( ( umode.ne.'hPa'     ).and.
233
     >      ( umode.ne.'hPa,agl' ).and.   
234
     >      ( umode.ne.'K'       ).and.
235
     >      ( umode.ne.'PVU'     ).and.
236
     >      ( umode.ne.'INDEX'   )  ) 
237
     > then  
238
          print*,' ERROR: unit not supported ',trim(umode)
239
          stop
240
       endif
241
 
242
c     Read selection criterion (dummy read)
243
      read(10,*) selectstr
244
 
245
c     Read flag for 'no time check'
246
      read(10,*) timecheck 
247
 
248
c     Close parameter file
249
      close(10)
250
 
251
c     Decide which output format is used (1..4: trajectory format, -1: triple list)
252
      call mode_tra(oformat,ofile)
253
 
254
c     Decide whether all lat/lon/lev coordaintes are read from one file
255
      if ( (hmode.eq.'file').and.(vmode.eq.'nil') ) then
256
         hmode='file3'
257
      elseif ( (hmode.eq.'file').and.(vmode.ne.'nil') ) then
258
         hmode='file2'
259
      endif
260
 
261
c     Convert timeshift (hh.mm) into a fractional time shift
262
      call hhmm2frac(timeshift,tfrac)
263
      if (tfrac.gt.0.) then
264
         tfrac=tfrac/timeinc
265
      else
266
         tfrac=0.
267
      endif
268
 
19 michaesp 269
c     If a mask file is provided, no time interpolation is allowed
270
      if ( (hmode.eq.'mask.grid').or.(hmode.eq.'mask.eqd') ) then
271
         if ( abs(tfrac).gt.eps ) then
272
            print*,' ERROR: no intermediate times allowed for ',
273
     >             trim(hmode)
274
            stop
275
         endif
276
      endif
277
 
3 michaesp 278
c     Read the region coordinates if needed
279
      if ( (hmode.eq.'region.eqd' ).or.
280
     >     (hmode.eq.'region.grid') ) then
281
 
282
         open(10,file=regionf)
283
 
284
 50       read(10,*,end=51) string
285
 
286
          if ( string(1:1).ne.'#' ) then
287
             call regionsplit(string,i,xcorner,ycorner)
288
             if ( i.eq.iregion ) goto 52
289
          endif
290
 
291
          goto 50
292
 
293
 51      close(10)
294
 
295
         print*,' ERROR: region ',iregion,' not found on ',trim(regionf)
296
         stop
297
 
298
 52      continue
299
 
300
      endif
301
 
302
c     Write some status information
303
      print*,'---- INPUT PARAMETERS -----------------------------------'
304
      print*
305
      if ( timeshift.gt.0. ) then
306
         print*,'  P file                     : ',trim(pfile0),
307
     >                                            ' ',
308
     >                                            trim(pfile1)
309
         print*,'  S file                     : ',trim(sfile0),
310
     >                                            ' ',
311
     >                                            trim(sfile1)
312
      else
313
         print*,'  P file                     : ',trim(pfile0)
314
         print*,'  S file                     : ',trim(sfile0)
315
      endif
316
      print*,'  Output file                : ',trim(ofile) 
317
      print*
318
      if (oformat.eq.-1) then
319
         print*,'  Output format              : (lon,lat,lev)-list'
320
      else
321
         print*,'  Output format              : ',oformat
322
      endif
323
      print*
324
      print*,'  Reference time (year)      : ',reftime(1)
325
      print*,'                 (month)     : ',reftime(2)
326
      print*,'                 (day)       : ',reftime(3)
327
      print*,'                 (hour)      : ',reftime(4)
328
      print*,'                 (min)       : ',reftime(5)
329
      print*,'  Time range                 : ',reftime(6)
330
      print*
331
      print*,'  Time shift                 : ',timeshift,' + ',
332
     >                                         trim(pfile0)
333
      print*,'  Region file                : ',trim(regionf)
334
      print*
335
      print*,'  hmode                      : ',trim(hmode)
336
      if ( hmode.eq.'file2' ) then        
337
        print*,'      filename [lat/lon]      : ',trim(hfile)
338
      elseif ( hmode.eq.'file3' ) then        
339
        print*,'      filename [lat/lon/lev]  : ',trim(hfile)
340
      elseif ( hmode.eq.'line' ) then   
341
        write(*,'(a30,4f10.2,i4)') 
342
     >         '      lon1,lon2,lat1,lat2,n   : ',lon1,lon2,lat1,lat2,hn
343
      elseif ( hmode.eq.'box.eqd' ) then     
344
        write(*,'(a30,5f10.2)') 
345
     >         '      lon1,lon2,lat1,lat2,ds  : ',lon1,lon2,lat1,lat2,ds
346
      elseif ( hmode.eq.'box.grid' ) then    
347
        write(*,'(a30,4f10.2)') 
348
     >         '      lon1,lon2,lat1,lat2     : ',lon1,lon2,lat1,lat2 
349
      elseif ( hmode.eq.'point' ) then   
350
        print*,'      lon,lat                 : ',lon1,lat1
351
      elseif ( hmode.eq.'shift' ) then   
352
        write(*,'(a30,4f10.2)') 
353
     >         '      lon,lat,dlon,dlat       : ',lon1,lat1,dlon,dlat
354
      elseif ( hmode.eq.'polygon.eqd' ) then   
355
        write(*,'(a30,a10,f10.2)') 
356
     >         '      hfile, ds               : ',trim(hfile),ds
357
      elseif ( hmode.eq.'polygon.grid' ) then   
358
        write(*,'(a30,a10)') 
359
     >         '      hfile                   : ',trim(hfile)
360
      elseif ( hmode.eq.'circle.eqd' ) then   
361
        write(*,'(a30,4f10.2)') 
362
     >         '      lonc,latc,radius, ds    : ',lon1,lat1,radius,ds
363
      elseif ( hmode.eq.'circle.grid' ) then   
364
        write(*,'(a30,3f10.2)') 
365
     >         '      lonc,latc,radius        : ',lon1,lat1,radius
366
      elseif ( hmode.eq.'region.eqd' ) then   
367
        write(*,'(a30,i4,1f10.2)') 
368
     >         '      iregion, ds             : ',iregion,ds
369
        write(*,'(a30,4f10.2)')
370
     >         '      xcorner                 : ',(xcorner(i),i=1,4)
371
        write(*,'(a30,4f10.2)')
372
     >         '      ycorner                 : ',(ycorner(i),i=1,4)
373
      elseif ( hmode.eq.'region.grid' ) then   
374
        write(*,'(a30,i4)') 
375
     >         '      iregion                 : ',iregion
376
        write(*,'(a30,4f10.2)')
377
     >         '      xcorner                 : ',(xcorner(i),i=1,4)
378
        write(*,'(a30,4f10.2)')
379
     >         '      ycorner                 : ',(ycorner(i),i=1,4)
19 michaesp 380
      elseif ( hmode.eq.'mask.eqd' ) then
381
        write(*,'(a30,a15,a15,f10.2)')
382
     >         '      hfile,variable,ds       : ',trim(hfile),
383
     >                                            trim(maskname),
384
     >                                            ds
385
      elseif ( hmode.eq.'mask.grid' ) then
386
        write(*,'(a30,a15,a15)')
387
     >         '      hfile,variable          : ',trim(hfile),
388
     >                                            trim(maskname)
3 michaesp 389
      endif
19 michaesp 390
 
3 michaesp 391
      print*
392
      print*,'  vmode                      : ',trim(vmode)
393
      if ( vmode.eq.'file') then 
394
         print*,'      filename               : ',trim(vfile)
395
      elseif ( vmode.eq.'level' ) then 
396
         print*,'      level                  : ',lev1                         
397
      elseif ( vmode.eq.'list') then 
398
         print*,'      n                      : ',vn
399
         print*,'      level(i)               : ',(levlist(i),i=1,vn)
400
      elseif ( vmode.eq.'profile') then 
401
         print*,'      lev1,lev2,n            : ',lev1,lev2,vn
402
      elseif ( vmode.eq.'grid') then 
403
         print*,'      lev1,lev2              : ',lev1,lev2
404
      endif
19 michaesp 405
 
3 michaesp 406
      print* 
407
      print*,'  umode                      : ',trim(umode)
408
      print*
409
      print*,'  time check                 : ',trim(timecheck)
410
      print*
411
 
19 michaesp 412
c     <mask.grid> and <mask.eqd>: split variable into name, operator,value
413
      if ( (hmode.eq.'mask.grid').or.(hmode.eq.'mask.eqd') ) then
414
          len = len_trim(maskname)
415
          do i=1,len_trim(maskname)-1
416
            ch1 = maskname(i:i)
417
            ch2 = maskname(i+1:i+1)
418
 
419
            if ( (ch1.eq.'<').and.(ch2.eq.'>') ) then
420
               read(maskname(i+2:len),*) maskvalue
421
               maskname = maskname(1:i-1)
422
               maskoper = 'ne'
423
               goto 90
424
            elseif ( (ch1.eq.'<').and.(ch2.eq.'=') ) then
425
               read(maskname(i+2:len),*) maskvalue
426
               maskname = maskname(1:i-1)
427
               maskoper = 'le'
428
               goto 90
429
            elseif ( (ch1.eq.'>').and.(ch2.eq.'=') ) then
430
               read(maskname(i+2:len),*) maskvalue
431
               maskname = maskname(1:i-1)
432
               maskoper = 'ge'
433
               goto 90
434
            elseif ( (ch1.eq.'=').and.(ch2.eq.'=') ) then
435
               read(maskname(i+2:len),*) maskvalue
436
               maskname = maskname(1:i-1)
437
               maskoper = 'eq'
438
               goto 90
439
            elseif ( ch1.eq.'=' ) then
440
               read(maskname(i+1:len),*) maskvalue
441
               maskname = maskname(1:i-1)
442
               maskoper = 'eq'
443
               goto 90
444
            elseif ( ch1.eq.'>' ) then
445
               read(maskname(i+1:len),*) maskvalue
446
               maskname = maskname(1:i-1)
447
               maskoper = 'gt'
448
               goto 90
449
            elseif ( ch1.eq.'<' ) then
450
               read(maskname(i+1:len),*) maskvalue
451
               maskname = maskname(1:i-1)
452
               maskoper = 'lt'
453
               goto 90
454
             endif
455
 
456
          enddo
457
      endif
458
 
459
  90  continue
460
 
3 michaesp 461
c     ------------------------------------------------------------------
462
c     Read grid parameters from inital files
463
c     ------------------------------------------------------------------
464
 
465
c     Get the time of the first and second data file
466
      tstart = -timeshift                              ! Format hh.mm
467
      call hhmm2frac(tstart,frac)
468
      frac   = frac + timeinc          
469
      call frac2hhmm(frac,tend)                        ! Format hh.mm                  
470
 
471
c     Convert timeshift (hh.mm) into a fractional time shift
472
      tfrac=real(int(timeshift))+
473
     >          100.*(timeshift-real(int(timeshift)))/60.
474
      if (tfrac.gt.0.) then
475
         tfrac=tfrac/timeinc
476
      else
477
         tfrac=0.
478
      endif
479
 
480
c     Read the constant grid parameters (nx,ny,nz,xmin,xmax,ymin,ymax,
481
c     pollon,pollat) The negative <-fid> of the file identifier is used 
482
c     as a flag for parameter retrieval
483
      varname  = 'U'
484
      nx       = 1
485
      ny       = 1
486
      nz       = 1
487
      call input_open (fid,pfile0)
488
      call input_grid (-fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
489
     >                 tstart,pollon,pollat,rd,rd,nz,rd,rd,rd,timecheck)
490
      call input_close(fid)
491
 
492
c     Check whether region coordinates are within the domain
493
      if ( (hmode.eq.'region.eqd' ).or.
494
     >     (hmode.eq.'region.grid') ) then
495
 
496
         do i=1,4
497
            if ( (xcorner(i).lt.xmin).or.
498
     >           (ycorner(i).lt.ymin).or.
499
     >           (xcorner(i).gt.xmax).or.
500
     >           (ycorner(i).gt.ymax) )
501
     >      then
502
               print*,' ERROR: region not included in data domain...'
503
               print*,'     ',trim(string)
504
               print*,'     ',(xcorner(j),j=1,4)
505
               print*,'     ',(ycorner(j),j=1,4)
506
               stop
507
            endif
508
 
509
         enddo
510
 
511
      endif
512
 
513
C     Check if the number of levels is too large
514
      if (nz.gt.nlevmax) goto 993
515
 
516
c     Allocate memory for 3d arrays: pressure, theta, pv
517
      allocate(pr(nx,ny,nz),stat=stat)
518
      if (stat.ne.0) print*,'*** error allocating array pr ***'
519
      allocate(prs(nx,ny),stat=stat)
520
      if (stat.ne.0) print*,'*** error allocating array prs **'
521
      allocate(th(nx,ny,nz),stat=stat)
522
      if (stat.ne.0) print*,'*** error allocating array th ***'
523
      allocate(ths(nx,ny),stat=stat)
524
      if (stat.ne.0) print*,'*** error allocating array ths **'
525
      allocate(pv(nx,ny,nz),stat=stat)
526
      if (stat.ne.0) print*,'*** error allocating array pv ***'
527
      allocate(pvs(nx,ny),stat=stat)
528
      if (stat.ne.0) print*,'*** error allocating array pvs **'
529
      allocate(in(nx,ny,nz),stat=stat)
530
      if (stat.ne.0) print*,'*** error allocating array in ***'
19 michaesp 531
      allocate(mask(nx,ny),stat=stat)
532
      if (stat.ne.0) print*,'*** error allocating array mask ***'
3 michaesp 533
 
534
c     Allocate memory for temporary arrays for time interpolation
535
      allocate(fld0(nx,ny,nz),stat=stat)
536
      if (stat.ne.0) print*,'*** error allocating array tmp0 ***'
537
      allocate(fld1(nx,ny,nz),stat=stat)
538
      if (stat.ne.0) print*,'*** error allocating array tmp1 ***'
539
      allocate(sfc0(nx,ny),stat=stat)
540
      if (stat.ne.0) print*,'*** error allocating array sfc0 ***'
541
      allocate(sfc1(nx,ny),stat=stat)
542
      if (stat.ne.0) print*,'*** error allocating array sfc1 ***'
543
 
544
c     ------ Index -----------------------------------------------------
545
 
546
c     Init the dummy array with vertical index
547
      do i=1,nx
548
         do j=1,ny
549
            do k=1,nz
550
               in(i,j,k) = real(k)
551
            enddo
552
         enddo
553
      enddo
554
 
555
c     ------ Pressure --------------------------------------------------
556
 
557
c     Read pressure from first data file (pfile0) on U-grid; we have to set
558
c     mdv explicitely, because it's not read from netCDF
559
      call input_open (fid,pfile0)
560
      varname='U'
15 michaesp 561
      stagz  = -0.5
3 michaesp 562
      call input_grid                                      
563
     >     (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
564
     >      tstart,pollon,pollat,fld0,sfc0,nz,ak,bk,stagz,timecheck)
565
      mdv = -999.99 
566
      call input_close(fid)
567
 
568
c     Read or set pressure for second data file (pfile1)
569
      if ( timeshift.ne.0.) then
570
         call input_open (fid,pfile1)
571
         varname='U'
572
         call input_grid                                      
573
     >        (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
574
     >        tend,pollon,pollat,fld1,sfc1,nz,ak,bk,stagz,timecheck)
575
         call input_close(fid)
576
      else
577
         do i=1,nx
578
            do j=1,ny
579
               do k=1,nz
580
                  fld1(i,j,k) = fld0(i,j,k)
581
               enddo
582
               sfc1(i,j) = sfc0(i,j)
583
            enddo
584
         enddo
585
      endif
586
 
587
c     Time interpolation to get the final pressure field
588
      do i=1,nx
589
         do j=1,ny
590
            do k=1,nz
591
               pr(i,j,k) = (1.-tfrac) * fld0(i,j,k) +
592
     >                     tfrac      * fld1(i,j,k)
593
            enddo
594
            prs(i,j) = (1.-tfrac) * sfc0(i,j) +
595
     >                 tfrac      * sfc1(i,j)
596
         enddo
597
      enddo
598
 
599
c     ------ Potential temperature -------------------------------------
600
 
601
       if ( (umode.eq.'K').or.(umode.eq.'PVU') ) then
602
 
603
c         Read potential temperature from first data file <sfile0>
604
          call input_open (fid,sfile0)
605
          varname='TH'                                      ! Theta                                      
606
          call input_wind
607
     >         (fid,varname,fld0,tstart,stagz,mdv,
608
     >         xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
609
          call input_close(fid)
610
 
611
c         Read or set potential temperature for second data file (sfile1)
612
          if ( timeshift.ne.0.) then
613
             call input_open (fid,sfile1)
614
             varname='TH' 
615
             call input_wind
616
     >            (fid,varname,fld1,tend,stagz,mdv,
617
     >            xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
618
             call input_close(fid)
619
          else
620
             do i=1,nx
621
                do j=1,ny
622
                   do k=1,nz
623
                      fld1(i,j,k) = fld0(i,j,k)
624
                   enddo
625
                enddo
626
             enddo
627
          endif
628
 
629
c         Time interpolation to get the final potential temperature field
630
          do i=1,nx
631
             do j=1,ny
632
                do k=1,nz
633
                   th(i,j,k) = (1.-tfrac) * fld0(i,j,k) +
634
     >                         tfrac      * fld1(i,j,k)
635
                enddo
636
             enddo
637
          enddo
638
 
639
c         Set the surface potential temperature
640
          do i=1,nx                            
641
             do j=1,ny
642
                ths(i,j)=th(i,j,1)
643
             enddo
644
          enddo
645
 
646
       endif
647
 
648
c     ------ Potential vorticity -----------------------------------------
649
 
650
       if ( (umode.eq.'PVU') ) then
651
 
652
c         Read potential vorticity from first data file <sfile0>
653
          call input_open (fid,sfile0)
654
          varname='PV'            
655
          call input_wind
656
     >         (fid,varname,fld0,tstart,stagz,mdv,
657
     >         xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
658
          call input_close(fid)
659
 
660
c         Read or set potential vorticity for second data file (sfile1)
661
          if ( timeshift.ne.0.) then
662
             call input_open (fid,sfile1)
663
             varname='PV' 
664
             call input_wind
665
     >            (fid,varname,fld1,tend,stagz,mdv,
666
     >            xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
667
             call input_close(fid)
668
          else
669
             do i=1,nx
670
                do j=1,ny
671
                   do k=1,nz
672
                      fld1(i,j,k) = fld0(i,j,k)
673
                   enddo
674
                enddo
675
             enddo
676
          endif
677
 
678
c         Time interpolation to get the final potential vorticity field
679
          do i=1,nx
680
             do j=1,ny
681
                do k=1,nz
682
                   pv(i,j,k) = (1.-tfrac) * fld0(i,j,k) +
683
     >                         tfrac      * fld1(i,j,k)
684
                enddo
685
             enddo
686
          enddo
687
 
688
c         Set the surface potential vorticity
689
          do i=1,nx                          
690
             do j=1,ny
691
                pvs(i,j)=pv(i,j,1)
692
             enddo
693
          enddo
694
       endif
695
 
19 michaesp 696
c     ------ Load 0/1 label (mask) ------------------------------------
697
      if ( (hmode.eq.'mask.grid').or.(hmode.eq.'mask.eqd') ) then
698
 
699
c         Explicit netCDF calls - no consistency check for grids!
700
          ierr = NF90_OPEN(hfile,nf90_nowrite, cdfid)
701
          IF ( ierr /= nf90_NoErr ) PRINT *,NF90_STRERROR(ierr)
702
          ierr = NF90_INQ_VARID(cdfid,maskname,varid)
703
          IF(ierr /= nf90_NoErr) PRINT *,NF90_STRERROR(ierr)
704
          ierr = nf90_get_var(cdfid,varid,mask)
705
          IF(ierr /= nf90_NoErr) PRINT *,NF90_STRERROR(ierr)
706
          ierr = NF90_CLOSE(cdfid)
707
          IF( ierr /= nf90_NoErr) PRINT *,NF90_STRERROR(ierr)
708
 
709
          nmask = 0
710
          do i=1,nx
711
            do j=1,ny
712
 
713
                if ( maskoper.eq.'eq' ) then
714
                   if ( abs(mask(i,j)-maskvalue).lt.eps ) then
715
                      nmask = nmask + 1
716
                      mask(i,j) = 1.
717
                   else
718
                      mask(i,j) = 0.
719
                   endif
720
                endif
721
 
722
                if ( maskoper.eq.'ne' ) then
723
                   if ( abs(mask(i,j)-maskvalue).gt.eps ) then
724
                      nmask = nmask + 1
725
                      mask(i,j) = 1.
726
                   else
727
                      mask(i,j) = 0.
728
                   endif
729
                endif
730
 
731
                if ( maskoper.eq.'gt' ) then
732
                   if ( mask(i,j).gt.maskvalue ) then
733
                      nmask = nmask + 1
734
                      mask(i,j) = 1.
735
                   else
736
                      mask(i,j) = 0.
737
                   endif
738
                endif
739
 
740
                if ( maskoper.eq.'lt' ) then
741
                   if ( mask(i,j).lt.maskvalue ) then
742
                      nmask = nmask + 1
743
                      mask(i,j) = 1.
744
                   else
745
                      mask(i,j) = 0.
746
                   endif
747
                endif
748
 
749
                if ( maskoper.eq.'ge' ) then
750
                   if ( mask(i,j).ge.maskvalue ) then
751
                      nmask = nmask + 1
752
                      mask(i,j) = 1.
753
                   else
754
                      mask(i,j) = 0.
755
                   endif
756
                endif
757
 
758
                if ( maskoper.eq.'le' ) then
759
                   if ( mask(i,j).le.maskvalue ) then
760
                      nmask = nmask + 1
761
                      mask(i,j) = 1.
762
                   else
763
                      mask(i,j) = 0.
764
                   endif
765
                endif
766
 
767
             enddo
768
          enddo
769
 
770
          print*
771
          print*,'  Mask read from file -> ',nmask,' 1-points'
772
          print*
773
 
774
      endif
775
 
3 michaesp 776
c     Write some status information
777
      print*,'---- CONSTANT GRID PARAMETERS ---------------------------'
778
      print*
779
      print*,'  xmin,xmax        : ',xmin,xmax
780
      print*,'  ymin,ymax        : ',ymin,ymax
781
      print*,'  dx,dy            : ',dx,dy
782
      print*,'  pollon,pollat    : ',pollon,pollat
783
      print*,'  nx,ny,nz         : ',nx,ny,nz
784
      print*
785
      print*,'  Pressure loaded  : ',trim(pfile0),' ',trim(pfile1)
786
      if ( (umode.eq.'K').or.(umode.eq.'PVU') ) then
787
         print*,'  Theta loaded     : ',trim(sfile0),' ',trim(sfile1)
788
      endif
789
      if ( (umode.eq.'PVU') ) then
790
         print*,'  PV loaded        : ',trim(sfile0),' ',trim(sfile1)
791
      endif
792
      print*
793
 
794
c     ------------------------------------------------------------------
795
c     Determine the expanded list of starting coordinates
796
c     ------------------------------------------------------------------
797
 
798
c     Write some status information
799
      print*,'---- EXPAND LIST OF STARTING POSITIONS  -----------------'
800
      print*
801
 
802
c     ------ Read lat/lon/lev from <hfile> -----------------------------
803
      if ( hmode.eq.'file3' ) then
804
         start_n = 0
805
         open(10,file=hfile)
806
 100       continue
807
              start_n = start_n + 1
808
              read(10,*,end=101) start_lon(start_n),
809
     >                           start_lat(start_n),
810
     >                           start_lev(start_n)
811
              goto 100
812
 101       continue
813
           start_n = start_n - 1
814
         close(10)
815
         goto 400
816
      endif
817
 
818
c     ------ Get lat/lon (horizontal) coordinates ---------------------
819
 
820
c     Read lat/lon from <hfile> 
821
      if ( hmode.eq.'file2' ) then
822
         hn = 0
823
         open(10,file=hfile)
824
 200       continue
825
              hn = hn + 1
826
              read(10,*,end=201) lonlist(hn),
827
     >                           latlist(hn)
828
              goto 200
829
 201       continue
830
           hn = hn - 1
831
         close(10)
832
      endif
833
 
834
c     Get lat/lon along a line (linear in lat/lon space)
835
      if ( hmode.eq.'line' ) then
836
         do i=1,hn
837
            lonlist(i) = lon1 + real(i-1)/real(hn-1)*(lon2-lon1)
838
            latlist(i) = lat1 + real(i-1)/real(hn-1)*(lat2-lat1)
839
         enddo
840
      endif
841
 
842
c     Lat/lon box: equidistant
843
      if ( hmode.eq.'box.eqd' ) then
844
         hn  = 0
845
         lat = lat1
846
         do while ( lat.le.lat2 )      
847
           lon = lon1
848
           do while ( lon.le.lon2 ) 
849
             hn          = hn+1
850
             lonlist(hn) = lon
851
             latlist(hn) = lat
852
             lon         = lon+ds/(deltat*cos(pi180*lat))
853
           enddo
854
           lat = lat+ds/deltat
855
        enddo
856
      endif
857
 
858
c     Lat/lon box: grid
859
      if ( hmode.eq.'box.grid' ) then
860
         hn = 0
861
         do j=1,ny
862
            do i=1,nx
863
               lon = xmin + real(i-1) * dx
864
               lat = ymin + real(j-1) * dy
865
               if ( (lon.ge.lon1).and.(lon.le.lon2).and.
866
     >              (lat.ge.lat1).and.(lat.le.lat2) )
867
     >         then
868
                  hn          = hn+1
869
                  lonlist(hn) = lon
870
                  latlist(hn) = lat
871
               endif
872
            enddo
873
         enddo
874
      endif
875
 
876
c     Get single starting point
877
      if ( hmode.eq.'point' ) then
878
         hn          = 1
879
         lonlist(hn) = lon1
880
         latlist(hn) = lat1
881
      endif
882
 
883
c     Get shifted and central starting point
884
      if ( hmode.eq.'shift' ) then
885
         hn         = 5
886
         lonlist(1) = lon1
887
         latlist(1) = lat1
888
         lonlist(2) = lon1+dlon
889
         latlist(2) = lat1
890
         lonlist(3) = lon1-dlon
891
         latlist(3) = lat1
892
         lonlist(4) = lon1
893
         latlist(4) = lat1+dlat
894
         lonlist(5) = lon1
895
         latlist(5) = lat1-dlat
896
      endif
897
 
898
c     Lat/lon polygon: grid
899
      if ( hmode.eq.'polygon.grid' ) then
900
 
901
c        Read list of polygon coordinates
902
         pn = 0
903
         open(10,file=hfile)
904
           read(10,*) loninpoly,latinpoly
905
 210       continue
906
              pn = pn + 1
907
              read(10,*,end=211) lonpoly(pn),
908
     >                           latpoly(pn)
909
 
910
              print*,pn,lonpoly(pn),latpoly(pn)
911
 
912
              goto 210
913
 211       continue
914
           pn = pn - 1
915
         close(10)
916
 
917
c        Define the polygon boundaries
918
         call DefSPolyBndry(latpoly,lonpoly,pn,latinpoly,loninpoly)
919
 
920
c        Get the grid points inside the polygon
921
         hn = 0
922
         do j=1,ny
923
            do i=1,nx
924
               lon = xmin + real(i-1) * dx
925
               lat = ymin + real(j-1) * dy
926
 
927
               call LctPtRelBndry(lat,lon,flag)
928
 
929
               if ( (flag.eq.1).or.(flag.eq.2) ) then
930
                  hn          = hn+1
931
                  lonlist(hn) = lon
932
                  latlist(hn) = lat
933
               endif
934
 
935
            enddo
936
         enddo
937
 
938
      endif
939
 
940
c     Lat/lon polygon: equidistant
941
      if ( hmode.eq.'polygon.eqd' ) then
942
 
943
c        Read list of polygon coordinates
944
         pn = 0
945
 
946
         open(10,file=hfile)
947
           read(10,*) loninpoly,latinpoly
948
 220       continue
949
              pn = pn + 1
950
              read(10,*,end=221) lonpoly(pn),
951
     >                           latpoly(pn)
952
              goto 220
953
 221       continue
954
           pn = pn - 1
955
         close(10)
956
 
957
 
958
c        Define the polygon boundaries
959
         call DefSPolyBndry(latpoly,lonpoly,pn,latinpoly,loninpoly)
960
 
961
c        Get the grid points inside the polygon
962
         hn  = 0
963
         lat = -90.
964
         do while ( lat.le.90. )      
965
           lon = -180.
966
           do while ( lon.lt.180. )
967
 
968
              call LctPtRelBndry(lat,lon,flag)
969
 
970
               if ( (flag.eq.1).or.(flag.eq.2) ) then
971
                  hn          = hn+1
972
                  lonlist(hn) = lon
973
                  latlist(hn) = lat
974
 
975
               endif
976
 
977
               lon = lon+ds/(deltat*cos(pi180*lat))
978
           enddo
979
           lat = lat+ds/deltat
980
 
981
        enddo
982
 
983
      endif
984
 
985
c     Circle: equidistant
986
      if ( hmode.eq.'circle.eqd' ) then
987
         hn  = 0
988
         lat = ymin
989
         do while ( lat.le.ymax )      
990
           lon = xmin
991
           do while ( lon.le.xmax ) 
992
              dist = sdis(lon1,lat1,lon,lat)
993
              if ( dist.le.radius ) then
994
                 hn          = hn+1
995
                 lonlist(hn) = lon
996
                 latlist(hn) = lat
997
              endif
998
              lon = lon+ds/(deltat*cos(pi180*lat))
999
           enddo
1000
           lat = lat+ds/deltat
1001
        enddo
1002
      endif
1003
 
1004
c     Circle: grid
1005
      if ( hmode.eq.'circle.grid' ) then
1006
         hn = 0
1007
         do j=1,ny
1008
            do i=1,nx
1009
               lon = xmin + real(i-1) * dx
1010
               lat = ymin + real(j-1) * dy
1011
               dist = sdis(lon1,lat1,lon,lat)
1012
               if ( dist.le.radius ) then
1013
                  hn          = hn+1
1014
                  lonlist(hn) = lon
1015
                  latlist(hn) = lat
1016
               endif
1017
            enddo
1018
         enddo
1019
 
1020
      endif
1021
 
1022
c     Region: equidistant
1023
      if ( hmode.eq.'region.eqd' ) then
1024
         hn  = 0
1025
         lat = ymin
1026
         do while ( lat.le.ymax )      
1027
           lon = xmin
1028
           do while ( lon.le.xmax ) 
1029
              flag = inregion(lon,lat,xcorner,ycorner)
1030
              if ( flag.eq.1 ) then
1031
                 hn          = hn+1
1032
                 lonlist(hn) = lon
1033
                 latlist(hn) = lat
1034
              endif
1035
              lon = lon+ds/(deltat*cos(pi180*lat))
1036
           enddo
1037
           lat = lat+ds/deltat
1038
        enddo
1039
      endif
1040
 
1041
c     Region: grid
1042
      if ( hmode.eq.'region.grid' ) then
1043
         hn = 0
1044
         do j=1,ny
1045
            do i=1,nx
1046
               lon  = xmin + real(i-1) * dx
1047
               lat  = ymin + real(j-1) * dy
1048
               flag = inregion(lon,lat,xcorner,ycorner)
1049
               if ( flag.eq.1 ) then
1050
                  hn          = hn+1
1051
                  lonlist(hn) = lon
1052
                  latlist(hn) = lat
1053
               endif
1054
            enddo
1055
         enddo
1056
      endif
1057
 
19 michaesp 1058
c     mask: grid
1059
      if ( hmode.eq.'mask.grid' ) then
1060
         hn  = 0
1061
         do i=1,nx
1062
          do j=1,ny
1063
            if ( mask(i,j).gt.0.5 ) then
1064
               hn = hn +1
1065
               lonlist(hn) = xmin + real(i-1) * dx
1066
               latlist(hn) = ymin + real(j-1) * dy
1067
            endif
1068
          enddo
1069
         enddo
1070
      endif
1071
 
1072
c     mask: equidistant
1073
      if ( hmode.eq.'mask.eqd' ) then
1074
         hn  = 0
1075
         lat = -90.+dy
1076
         do while ( lat.le.ymax )
1077
           lon = -180.
1078
           do while ( lon.le.xmax )
1079
              indx = nint( (lon-xmin)/dx + 1. )
1080
              indy = nint( (lat-ymin)/dy + 1. )
1081
              if ( indx.lt. 1 ) indx =  1
1082
              if ( indy.lt. 1 ) indy =  1
1083
              if ( indx.gt.nx ) indx = nx
1084
              if ( indy.gt.ny ) indy = ny
1085
              if ( mask(indx,indy).gt.0.5 ) then
1086
                 hn          = hn+1
1087
                 lonlist(hn) = lon
1088
                 latlist(hn) = lat
1089
              endif
1090
              lon = lon+ds/(deltat*cos(pi180*lat))
1091
           enddo
1092
           lat = lat+ds/deltat
1093
        enddo
1094
      endif
1095
 
1096
 
3 michaesp 1097
c     ------ Get lev (vertical) coordinates -------------------------
1098
 
1099
c     Read level list from file
1100
      if ( vmode.eq.'file' ) then
1101
         vn = 0
1102
         open(10,file=vfile)
1103
 300       continue
1104
              vn = vn + 1
1105
              read(10,*,end=301) levlist(vn)
1106
              goto 300
1107
 301       continue
1108
           vn = vn - 1
1109
         close(10)
1110
      endif
1111
 
1112
c     Get single starting level
1113
      if ( vmode.eq.'level' ) then
1114
         vn          = 1
1115
         levlist(vn) = lev1
1116
      endif
1117
 
1118
c     Get level profile
1119
      if ( vmode.eq.'profile' ) then
1120
         do i=1,vn
1121
            levlist(i) = lev1 + real(i-1)/real(vn-1)*(lev2-lev1)
1122
         enddo
1123
      endif
1124
 
1125
c     Get all grid points in a layer: at the moment set the list of levels to 
1126
c     all indices from 1 to nz; later the correct subset of indices will be chosen
1127
      if ( vmode.eq.'grid' ) then
1128
         vn = nz
1129
         do i=1,vn
1130
            levlist(i) = real(i)
1131
         enddo
1132
         umode_save = umode
1133
         umode      = 'INDEX'
1134
 
1135
      endif
1136
 
1137
c     ------ Compile the complete list of starting positions  ------
1138
 
1139
c     Get all starting points in specified vertical coordinate system
1140
      start_n = 0
1141
      do i=1,vn
1142
         do j=1,hn
1143
 
1144
            start_n = start_n + 1
1145
            start_lon(start_n) = lonlist(j)
1146
            start_lat(start_n) = latlist(j)
1147
            start_lev(start_n) = levlist(i)
1148
 
1149
         enddo
1150
      enddo
1151
 
1152
c     ------ Exit point of this section
1153
 400  continue
1154
 
1155
c     Write status information
1156
      print*,'  # expanded points : ', start_n
1157
      print*
1158
 
1159
c     ------------------------------------------------------------------
1160
c     Transform starting levels into pressure
1161
c     ------------------------------------------------------------------
1162
 
1163
c     Write some status information
1164
      print*,'---- STARTING POSITIONS ---------------------------------'
1165
      print*
1166
 
1167
c     Vertical mode <hPa,asl> or simply <hPa>
1168
      if ( (umode.eq.'hPa,asl').or.(umode.eq.'hPa') ) then
1169
 
1170
         do i=1,start_n
1171
            start_pre(i) = start_lev(i)
1172
         enddo
1173
 
1174
c     Vertical mode <hPa,agl>
1175
      elseif ( umode.eq.'hPa,agl' ) then
1176
 
1177
         do i=1,start_n
1178
            call get_index3(rid,rjd,rkd,start_lon(i),start_lat(i),1050.,
1179
     >                      3,pr,prs,nx,ny,nz,xmin,ymin,dx,dy)
1180
            tmp1 = int_index3 (prs,nx,ny,1,rid,rjd,1,mdv)
1181
            start_pre(i) = tmp1 - start_lev(i) 
1182
         enddo
1183
 
1184
c     Vertical mode <K>
1185
      elseif ( umode.eq.'K' ) then
1186
 
1187
         do i=1,start_n
1188
            call get_index3(rid,rjd,rkd,start_lon(i),start_lat(i),
1189
     >                start_lev(i),1,th,ths,nx,ny,nz,xmin,ymin,dx,dy)
1190
            tmp1 = int_index3 (pr,nx,ny,nz,rid,rjd,rkd,mdv)
1191
            start_pre(i) = tmp1
1192
         enddo
1193
 
1194
c     Vertical mode <PVU>
1195
      elseif ( umode.eq.'PVU' ) then
1196
 
1197
         do i=1,start_n
1198
            call get_index3(rid,rjd,rkd,start_lon(i),start_lat(i),
1199
     >               start_lev(i),2,pv,pvs,nx,ny,nz,xmin,ymin,dx,dy)
1200
            tmp1 = int_index3 (pr,nx,ny,nz,rid,rjd,rkd,mdv)
1201
            start_pre(i) = tmp1
1202
         enddo
1203
 
1204
c     Vertical mode <INDEX>
1205
      elseif ( umode.eq.'INDEX' ) then
1206
 
1207
         do i=1,start_n
1208
            call get_index3(rid,rjd,rkd,start_lon(i),start_lat(i),
1209
     >               1050.,2,pv,pvs,nx,ny,nz,xmin,ymin,dx,dy)
1210
            rkd = start_lev(i)
1211
            tmp1 = int_index3 (pr,nx,ny,nz,rid,rjd,rkd,mdv)
1212
            start_pre(i) = tmp1
1213
         enddo
1214
 
1215
      endif
1216
 
1217
c     ------------------------------------------------------------------
1218
c     Remove invalid points from the list
1219
c     ------------------------------------------------------------------
1220
 
1221
c     Select the correct subset if <vmode=grid>: starting points outside the layer
1222
c     will receive a <mdv> vertical pressure and will be removed
1223
      if ( vmode.eq.'grid' ) then
1224
 
1225
         do i=1,start_n
1226
 
1227
c           Get the pressure at the grid point
1228
            if ( ( umode_save.eq.'hPa'      ).or.
1229
     >            (umode_save.eq.'hPa,asl') ) 
1230
     >      then
1231
 
1232
               call get_index3(rid,rjd,rkd,start_lon(i),start_lat(i),
1233
     >                         start_pre(i),3,pr,prs,nx,ny,nz,xmin,
1234
     >                         ymin,dx,dy)
1235
               tmp1 = int_index3 (pr,nx,ny,nz,rid,rjd,rkd,mdv)
1236
 
1237
c           Get pressure AGL at grid point
1238
            elseif ( umode_save.eq.'hPa,agl' ) then
1239
 
1240
               call get_index3(rid,rjd,rkd,start_lon(i),
1241
     >                         start_lat(i),start_pre(i),3,pr,prs,
1242
     >                         nx,ny,nz,xmin,ymin,dx,dy)
1243
               tmp1 = int_index3 (pr,nx,ny,nz,rid,rjd,rkd,mdv)
1244
               call get_index3(rid,rjd,rkd,start_lon(i),
1245
     >                         start_lat(i),1050.,3,pr,prs,nx,ny,
1246
     >                         nz,xmin,ymin,dx,dy)
1247
               tmp2 = int_index3 (prs,nx,ny,1,rid,rjd,1,mdv)
1248
               tmp1 = tmp2 - tmp1
1249
 
1250
c           Get potential temperature at grid point
1251
            elseif ( umode_save.eq.'K' ) then
1252
 
1253
               call get_index3(rid,rjd,rkd,start_lon(i),start_lat(i),
1254
     >                         start_pre(i),3,pr,prs,nx,ny,nz,
1255
     >                         xmin,ymin,dx,dy)
1256
               tmp1 = int_index3 (th,nx,ny,nz,rid,rjd,rkd,mdv)
1257
 
1258
c           Get potential vorticity at the grid point
1259
            elseif ( umode_save.eq.'PVU' ) then
1260
 
1261
               call get_index3(rid,rjd,rkd,start_lon(i),start_lat(i),
1262
     >                         start_pre(i),3,pr,prs,nx,ny,nz,xmin,
1263
     >                         ymin,dx,dy)
1264
               tmp1 = int_index3 (pv,nx,ny,nz,rid,rjd,rkd,mdv)
1265
 
1266
c           Get vertical index at the grid point
1267
            elseif ( umode_save.eq.'INDEX' ) then
1268
 
1269
              call get_index3(rid,rjd,rkd,start_lon(i),start_lat(i),
1270
     >                        start_pre(i),3,pr,prs,nx,ny,nz,
1271
     >                        xmin,ymin,dx,dy)
1272
              tmp1 = int_index3 (in,nx,ny,nz,rid,rjd,rkd,mdv)
1273
 
1274
           endif
1275
 
1276
c          Remove points outside layer
1277
           if ( ( tmp1.lt.lev1).or.(tmp1.gt.lev2) ) then
1278
              start_pre(i) = mdv
1279
           endif
1280
 
1281
         enddo         
1282
 
1283
      endif
1284
 
1285
c     Check whether the starting levels are valid (in data domain)
1286
      do i=1,start_n
1287
 
1288
         call get_index3(rid,rjd,rkd,start_lon(i),start_lat(i),1050.,
1289
     >                   3,pr,prs,nx,ny,nz,xmin,ymin,dx,dy)
1290
         tmp1 = int_index3 (prs,nx,ny, 1,rid,rjd,real( 1),mdv)           ! Surface
1291
         tmp2 = int_index3 (pr ,nx,ny,nz,rid,rjd,real(nz),mdv)           ! Top of domain
1292
 
1293
         if ( (start_pre(i).gt.tmp1).or.
1294
     >        (start_pre(i).lt.tmp2).or.
1295
     >        (start_lon(i).lt.xmin).or. 
1296
     >        (start_lon(i).gt.xmax).or. 
1297
     >        (start_lat(i).lt.ymin).or.
1298
     >        (start_lat(i).gt.ymax) )
1299
     >   then
1300
            start_pre(i) = mdv
1301
         endif
1302
 
1303
      enddo
1304
 
1305
c     Remove all starting points outside the domain
1306
      i        = 1
1307
      do while ( i.le.start_n )
1308
 
1309
         if ( abs(start_pre(i)-mdv).lt.eps ) then
1310
 
1311
            if ( vmode.ne.'grid') then
1312
             print*,'  Outside ', start_lon(i),start_lat(i),start_lev(i)
1313
            endif
1314
 
1315
            do j=i,start_n
1316
               start_lon(j) = start_lon(j+1)
1317
               start_lat(j) = start_lat(j+1)
1318
               start_pre(j) = start_pre(j+1)
1319
               start_lev(j) = start_lev(j+1)
1320
            enddo
1321
            start_n = start_n - 1
1322
 
1323
         else
1324
 
1325
            i = i + 1
1326
 
1327
         endif
1328
 
1329
      enddo
1330
 
1331
c     Write some status information
1332
      latmin = start_lat(1)
1333
      latmax = start_lat(1)
1334
      lonmin = start_lon(1)
1335
      lonmax = start_lon(1)
1336
      premin = start_pre(1)
1337
      premax = start_pre(1)
1338
      do i=1,start_n
1339
         if (start_lat(i).lt.latmin) latmin = start_lat(i)
1340
         if (start_lat(i).gt.latmax) latmax = start_lat(i)
1341
         if (start_lon(i).lt.lonmin) lonmin = start_lon(i)
1342
         if (start_lon(i).gt.lonmax) lonmax = start_lon(i)
1343
         if (start_pre(i).lt.premin) premin = start_pre(i)
1344
         if (start_pre(i).gt.premax) premax = start_pre(i)
1345
      enddo
1346
      print*,'  min(lat),max(lat) : ', latmin,latmax
1347
      print*,'  min(lon),max(lon) : ', lonmin,lonmax
1348
      print*,'  min(pre),max(pre) : ', premin,premax
1349
      print*
1350
      print*,'  # starting points : ', start_n
1351
      print*
1352
 
1353
 
1354
c     ------------------------------------------------------------------
1355
c     Write starting positions to output file
1356
c     ------------------------------------------------------------------
1357
 
1358
c     Output as a trajectory file (with only one time == 0)
1359
      if (oformat.ne.-1) then
1360
 
1361
         allocate(tra(start_n,1,5),stat=stat)
1362
 
1363
         vars(1)  ='time'
1364
         vars(2)  ='lon'
1365
         vars(3)  ='lat'
1366
         vars(4)  ='p'
1367
         vars(5)  ='level'
1368
         call wopen_tra(fid,ofile,start_n,1,5,reftime,vars,oformat)
1369
 
1370
         do i=1,start_n
1371
            tra(i,1,1) = 0.
1372
            tra(i,1,2) = start_lon(i)
1373
            tra(i,1,3) = start_lat(i)
1374
            tra(i,1,4) = start_pre(i)
1375
            tra(i,1,5) = start_lev(i)
1376
         enddo
1377
         call write_tra(fid,tra,start_n,1,5,oformat)
1378
 
1379
         call close_tra(fid,oformat)
1380
 
1381
c     Output as a triple list (corresponding to <startf> file)
1382
      else
1383
 
1384
         fid = 10
1385
         open(fid,file=ofile)
1386
          do i=1,start_n
1387
             write(fid,'(3f10.3)') start_lon(i),start_lat(i),
1388
     >                             start_pre(i) 
1389
          enddo
1390
         close(fid)
1391
 
1392
      endif
1393
 
1394
c     Write some status information, and end of program message
1395
      print*  
1396
      print*,'---- STATUS INFORMATION --------------------------------'
1397
      print*
1398
      print*,'ok'
1399
      print*
1400
      print*,'       *** END OF PROGRAM CREATE_STARTF ***'
1401
      print*,'========================================================='
1402
 
1403
c     ------------------------------------------------------------------
1404
c     Exception handling
1405
c     ------------------------------------------------------------------
1406
 
1407
      stop
1408
 
1409
 993  write(*,*) '*** ERROR: problems with array size'
1410
      call exit(1)
1411
 
1412
      end
1413
 
1414
c     --------------------------------------------------------------------------
1415
c     Split a region string and get corners of the domain
1416
c     --------------------------------------------------------------------------
1417
 
1418
      subroutine regionsplit(string,iregion,xcorner,ycorner)
1419
 
1420
c     The region string comes either as <lonw,lone,lats,latn> or as <lon1,lat1,
1421
c     lon2,lat2,lon3,lat3,lon4,lat4>: split it into ints components and get the
1422
c     four coordinates for the region
1423
 
1424
      implicit none
1425
 
1426
c     Declaration of subroutine parameters
1427
      character*80    string
1428
      real            xcorner(4),ycorner(4)
1429
      integer         iregion
1430
 
1431
c     Local variables
1432
      integer         i,n
1433
      integer         il,ir
1434
      real            subfloat (80)
1435
      integer         stat
1436
      integer         len
1437
 
1438
c     ------- Split the string
1439
      i    = 1
1440
      n    = 0
1441
      stat = 0
1442
      il   = 1
1443
      len  = len_trim(string)
1444
 
1445
 100  continue
1446
 
1447
c     Find start of a substring
1448
      do while ( stat.eq.0 )
1449
         if ( string(i:i).ne.' ' ) then
1450
            stat = 1
1451
            il   = i
1452
         else
1453
            i = i + 1
1454
         endif
1455
      enddo
1456
 
1457
c     Find end of substring
1458
      do while ( stat.eq.1 )         
1459
         if ( ( string(i:i).eq.' ' ) .or. ( i.eq.len ) ) then
1460
            stat = 2
1461
            ir   = i
1462
         else
1463
            i    = i + 1
1464
         endif
1465
      enddo
1466
 
1467
c     Convert the substring into a number
1468
      if ( stat.eq.2 ) then
1469
         n = n + 1
1470
         read(string(il:ir),*) subfloat(n)
1471
         stat = 0
1472
      endif
1473
 
1474
      if ( i.lt.len ) goto 100
1475
 
1476
 
1477
c     -------- Get the region number
1478
 
1479
      iregion = nint(subfloat(1))
1480
 
1481
c     -------- Get the corners of the region
1482
 
1483
      if ( n.eq.5 ) then     ! lonw(2),lone(3),lats(4),latn(5)
1484
 
1485
         xcorner(1) = subfloat(2)
1486
         ycorner(1) = subfloat(4)
1487
 
1488
         xcorner(2) = subfloat(3)
1489
         ycorner(2) = subfloat(4)
1490
 
1491
         xcorner(3) = subfloat(3)
1492
         ycorner(3) = subfloat(5)
1493
 
1494
         xcorner(4) = subfloat(2)
1495
         ycorner(4) = subfloat(5)
1496
 
1497
      elseif ( n.eq.9 ) then     ! lon1,lat1,lon2,lat2,lon3,lon4,lat4
1498
 
1499
         xcorner(1) = subfloat(2)
1500
         ycorner(1) = subfloat(3)
1501
 
1502
         xcorner(2) = subfloat(4)
1503
         ycorner(2) = subfloat(5)
1504
 
1505
         xcorner(3) = subfloat(6)
1506
         ycorner(3) = subfloat(7)
1507
 
1508
         xcorner(4) = subfloat(8)
1509
         ycorner(4) = subfloat(9)
1510
 
1511
      else
1512
 
1513
         print*,' ERROR: invalid region specification '
1514
         print*,'     ',trim(string)
1515
         stop
1516
 
1517
      endif
1518
 
1519
 
1520
      end
1521
 
1522
c     --------------------------------------------------------------------------
1523
c     Decide whether lat/lon point is in or out of region
1524
c     --------------------------------------------------------------------------
1525
 
1526
      integer function inregion (lon,lat,xcorner,ycorner)
1527
 
1528
c     Decide whether point (lon/lat) is in the region specified by <xcorner(1..4),
1529
c     ycorner(1..4).
1530
 
1531
      implicit none
1532
 
1533
c     Declaration of subroutine parameters
1534
      real    lon,lat
1535
      real    xcorner(4),ycorner(4)
1536
 
1537
c     Local variables
1538
      integer flag
1539
      real    xmin,xmax,ymin,ymax
1540
      integer i
1541
 
1542
c     Reset the flag
1543
      flag = 0
1544
 
1545
c     Set some boundaries
1546
      xmax = xcorner(1)
1547
      xmin = xcorner(1)
1548
      ymax = ycorner(1)
1549
      ymin = ycorner(1)
1550
      do i=2,4
1551
        if (xcorner(i).lt.xmin) xmin = xcorner(i)
1552
        if (xcorner(i).gt.xmax) xmax = xcorner(i)
1553
        if (ycorner(i).lt.ymin) ymin = ycorner(i)
1554
        if (ycorner(i).gt.ymax) ymax = ycorner(i)
1555
      enddo
1556
 
1557
c     Do the tests - set flag=1 if all tests pased
1558
      if (lon.lt.xmin) goto 970
1559
      if (lon.gt.xmax) goto 970
1560
      if (lat.lt.ymin) goto 970
1561
      if (lat.gt.ymax) goto 970
1562
 
1563
      if ((lon-xcorner(1))*(ycorner(2)-ycorner(1))-
1564
     >    (lat-ycorner(1))*(xcorner(2)-xcorner(1)).gt.0.) goto 970
1565
      if ((lon-xcorner(2))*(ycorner(3)-ycorner(2))-
1566
     >    (lat-ycorner(2))*(xcorner(3)-xcorner(2)).gt.0.) goto 970
1567
      if ((lon-xcorner(3))*(ycorner(4)-ycorner(3))-
1568
     >    (lat-ycorner(3))*(xcorner(4)-xcorner(3)).gt.0.) goto 970
1569
      if ((lon-xcorner(4))*(ycorner(1)-ycorner(4))-
1570
     >    (lat-ycorner(4))*(xcorner(1)-xcorner(4)).gt.0.) goto 970
1571
 
1572
      flag = 1
1573
 
1574
c     Return the value
1575
 970  continue
1576
 
1577
      inregion = flag
1578
 
1579
      return
1580
 
1581
      end
1582
 
1583
c     --------------------------------------------------------------------------
1584
c     Spherical distance between lat/lon points                                                       
1585
c     --------------------------------------------------------------------------
1586
 
1587
      real function sdis(xp,yp,xq,yq)
1588
c
1589
c     calculates spherical distance (in km) between two points given
1590
c     by their spherical coordinates (xp,yp) and (xq,yq), respectively.
1591
c
1592
      real      re
1593
      parameter (re=6370.)
1594
      real      pi180
1595
      parameter (pi180=3.14159/180.)
1596
      real      xp,yp,xq,yq,arg
1597
 
1598
      arg=sin(pi180*yp)*sin(pi180*yq)+
1599
     >    cos(pi180*yp)*cos(pi180*yq)*cos(pi180*(xp-xq))
1600
      if (arg.lt.-1.) arg=-1.
1601
      if (arg.gt.1.) arg=1.
1602
 
1603
      sdis=re*acos(arg)
1604
 
1605
      end
1606
 
1607
 
1608
c     ****************************************************************
1609
c     * Given some spherical polygon S and some point X known to be  *
1610
c     * located inside S, these routines will determine if an arbit- *
1611
c     * -rary point P lies inside S, outside S, or on its boundary.  *
1612
c     * The calling program must first call DefSPolyBndry to define  *
1613
c     * the boundary of S and the point X. Any subsequent call to    *
1614
c     * subroutine LctPtRelBndry will determine if some point P lies *
1615
c     * inside or outside S, or on its boundary. (Usually            *
1616
c     * DefSPolyBndry is called once, then LctPrRelBndry is called   *
1617
c     * many times).                                                 *
1618
c     *                                                              * 
1619
c     * REFERENCE:            Bevis, M. and Chatelain, J.-L. (1989)  * 
1620
c     *                       Maflaematical Geology, vol 21.         *
1621
c     * VERSION 1.0                                                  *
1622
c     ****************************************************************
1623
 
1624
      Subroutine DefSPolyBndry(vlat,vlon,nv,xlat, xlon)
1625
 
1626
c     ****************************************************************
1627
c     * This mmn entry point is used m define ~e spheric~ polygon S  *
1628
c     * and the point X.                                             *
1629
c     * ARGUMENTS:                                                   *
1630
c     * vlat,vlon (sent) ... vectors containing the latitude and     * 
1631
c     *                      longitude of each vertex of the         *
1632
c     *                      spherical polygon S. The ith.vertex is  *
1633
c     *                      located at [vlat(i),vlon(i)].           *
1634
c     * nv        (sent) ... the number of vertices and sides in the *
1635
c     *                      spherical polygon S                     *
1636
c     * xlat,xlon (sent) ... latitude and longitude of some point X  *
1637
c     *                      located inside S. X must not be located *
1638
c     *                      on any great circle that includes two   *
1639
c     *                      vertices of S.                          *
1640
c     *                                                              *
1641
c     * UNITS AND SIGN CONVENTION:                                   *
1642
c     *  Latitudes and longitudes are specified in degrees.          *
1643
c     *  Latitudes are positive to the north and negative to the     *
1644
c     *  south.                                                      *
1645
c     *  Longitudes are positive to the east and negative to the     *
1646
c     *  west.                                                       *
1647
c     *                                                              * 
1648
c     * VERTEX ENUMERATION:                                          * 
1649
c     * The vertices of S should be numbered sequentially around the *
1650
c     * border of the spherical polygon. Vertex 1 lies between vertex*
1651
c     * nv and vertex 2. Neighbouring vertices must be seperated by  *
1652
c     * less than 180 degrees. (In order to generate a polygon side  *
1653
c     * whose arc length equals or exceeds 180 degrees simply        *
1654
c     * introduce an additional (pseudo)vertex). Having chosen       *
1655
c     * vertex 1, the user may number the remaining vertices in      *
1656
c     * either direction. However if the user wishes to use the      *
1657
c     * subroutine SPA to determine the area of the polygon S (Bevis *
1658
c     * & Cambareri, 1987, Math. Geol., v.19, p. 335-346) then he or *
1659
c     * she must follow the convention whereby in moving around the  *
1660
c     * polygon border in the direction of increasing vertex number  *
1661
c     * clockwise bends occur at salient vertices. A vertex is       *
1662
c     * salient if the interior angle is less than 180 degrees.      *
1663
c     * (In the case of a convex polygon this convention implies     *
1664
c     * that vertices are numbered in clockwise sequence).           *
1665
c     ****************************************************************
1666
 
1667
      implicit none
1668
 
1669
      integer mxnv,nv
1670
 
1671
c     ----------------------------------------------------------------
1672
c     Edit next statement to increase maximum number of vertices that 
1673
c     may be used to define the spherical polygon S               
1674
c     The value of parameter mxnv in subroutine LctPtRelBndry must match
1675
c     that of parameter mxnv in this subroutine, as assigned above.
1676
c     ----------------------------------------------------------------
1677
      parameter (mxnv=500)
1678
 
1679
      real  vlat(nv),vlon(nv),xlat,xlon,dellon
1680
      real  tlonv(mxnv),vlat_c(mxnv),vlon_c(mxnv),xlat_c,xlon_c
1681
      integer i,ibndry,nv_c,ip
1682
 
1683
      data ibndry/0/
1684
 
1685
      common/spolybndry/vlat_c,vlon_c,nv_c,xlat_c,xlon_c,tlonv,ibndry
1686
 
1687
      if (nv.gt.mxnv) then
1688
         print *,'nv exceeds maximum allowed value'
1689
         print *,'adjust parameter mxnv in subroutine DefSPolyBndry'
1690
         stop
1691
      endif
1692
 
1693
      ibndry=1                  ! boundary defined at least once (flag)
1694
      nv_c=nv                   ! copy for named common
1695
      xlat_c=xlat               ! . . . .
1696
      xlon_c=xlon               !
1697
 
1698
      do i=1,nv
1699
         vlat_c(i)=vlat(i)      ! "
1700
         vlon_c(i)=vlon(i)      !
1701
 
1702
         call TrnsfmLon(xlat,xlon,vlat(i),vlon(i),tlonv(i))
1703
 
1704
         if (i.gt.1) then
1705
            ip=i-1
1706
         else
1707
            ip=nv
1708
         endif
1709
 
1710
         if ((vlat(i).eq.vlat(ip)).and.(vlon(i).eq.vlon(ip))) then
1711
            print *,'DefSPolyBndry detects user error:'
1712
            print *,'vertices ',i,' and ',ip,' are not distinct'
1713
            print*,'lat ',i,ip,vlat(i),vlat(ip)
1714
            print*,'lon ',i,ip,vlon(i),vlon(ip)            
1715
            stop
1716
         endif
1717
 
1718
         if (tlonv(i).eq.tlonv(ip)) then
1719
            print *,'DefSPolyBndry detects user error:'
1720
            print *,'vertices ',i,' & ',ip,' on same gt. circle as X'
1721
            stop
1722
         endif
1723
 
1724
         if (vlat(i).eq.(-vlat(ip))) then
1725
            dellon=vlon(i)-vlon(ip)
1726
            if (dellon.gt.+180.) dellon=dellon-360.
1727
            if (dellon.lt.-180.) dellon=dellon-360.
1728
            if ((dellon.eq.+180.0).or.(dellon.eq.-180.0)) then
1729
               print *,'DefSPolyBndry detects user error:'
1730
               print *,'vertices ',i,' and ',ip,' are antipodal'
1731
               stop
1732
            endif
1733
         endif
1734
      enddo
1735
 
1736
      return
1737
 
1738
      end
1739
 
1740
 
1741
c     ****************************************************************
1742
 
1743
      Subroutine LctPtRelBndry(plat,plon,location)
1744
 
1745
c     ****************************************************************
1746
 
1747
c     ****************************************************************
1748
c     * This routine is used to see if some point P is located       *
1749
c     * inside, outside or on the boundary of the spherical polygon  *
1750
c     * S previously defined by a call to subroutine DefSPolyBndry.  *
1751
c     * There is a single restriction on point P: it must not be     *
1752
c     * antipodal to the point X defined in the call to DefSPolyBndry*
1753
c     * (ie.P and X cannot be seperated by exactly 180 degrees).     *
1754
c     * ARGUMENTS:                                                   *  
1755
c     * plat,plon (sent)... the latitude and longitude of point P    *
1756
c     * location (returned)... specifies the location of P:          *
1757
c     *                        location=0 implies P is outside of S  *
1758
c     *                        location=1 implies P is inside of S   *
1759
c     *                        location=2 implies P on boundary of S *
1760
c     *                        location=3 implies user error (P is   *
1761
c     *                                     antipodal to X)          *
1762
c     * UNFfS AND SIGN CONVENTION:                                   * 
1763
c     *  Latitudes and longitudes are specified in degrees.          *
1764
c     *  Latitudes are positive to the north and negative to the     *
1765
c     *  south.                                                      *    
1766
c     *  Longitudes are positive to the east and negative to the     *
1767
c     *  west.                                                       *
1768
c     ****************************************************************
1769
 
1770
      implicit none
1771
 
1772
      integer mxnv
1773
 
1774
c     ----------------------------------------------------------------
1775
c     The statement below must match that in subroutine DefSPolyBndry
1776
c     ----------------------------------------------------------------
1777
 
1778
      parameter (mxnv=500)
1779
 
1780
      real tlonv(mxnv),vlat_c(mxnv),vlon_c(mxnv),xlat_c,xlon_c
1781
      real plat,plon,vAlat,vAlon,vBlat,vBlon,tlonA,tlonB,tlonP
1782
      real tlon_X,tlon_P,tlon_B,dellon
1783
      integer i,ibndry,nv_c,location,icross,ibrngAB,ibrngAP,ibrngPB
1784
      integer ibrng_BX,ibrng_BP,istrike
1785
 
1786
      common/spolybndry/vlat_c,vlon_c,nv_c,xlat_c,xlon_c,tlonv,ibndry
1787
 
1788
      if (ibndry.eq.0) then     ! user has never defined the bndry
1789
         print*,'Subroutine LctPtRelBndry detects user error:'
1790
         print*,'Subroutine DefSPolyBndry must be called before'
1791
         print*,'subroutine LctPtRelBndry can be called'
1792
         stop
1793
      endif
1794
 
1795
      if (plat.eq.(-xlat_c)) then
1796
         dellon=plon-xlon_c
1797
         if (dellon.lt.(-180.)) dellon=dellon+360.
1798
         if (dellon.gt.+180.) dellon=dellon-360.
1799
         if ((dellon.eq.+180.0).or.(dellon.eq.-180.)) then
1800
            print*,'Warning: LctPtRelBndry detects case P antipodal
1801
     >           to X'
1802
            print*,'location of P relative to S is undetermined'
1803
            location=3
1804
            return
1805
         endif
1806
      endif 
1807
 
1808
      location=0                ! default ( P is outside S)
1809
      icross=0                  ! initialize counter
1810
 
1811
      if ((plat.eq.xlat_c).and.(plon.eq.xlon_c)) then
1812
         location=1
1813
         return
1814
      endif
1815
 
1816
 
1817
      call TrnsfmLon (xlat_c,xlon_c,plat,plon,tlonP)
1818
 
1819
      do i=1,nv_c              ! start of loop over sides of S 
1820
 
1821
         vAlat=vlat_c(i)
1822
         vAlon=vlon_c(i)
1823
         tlonA=tlonv(i)
1824
 
1825
         if (i.lt.nv_c) then
1826
            vBlat=vlat_c(i+1)
1827
            vBlon=vlon_c(i+1)
1828
            tlonB=tlonv(i+1)
1829
         else
1830
            vBlat=vlat_c(1)
1831
            vBlon=vlon_c(1)
1832
            tlonB=tlonv(1)
1833
         endif
1834
 
1835
         istrike=0
1836
 
1837
         if (tlonP.eq.tlonA) then
1838
            istrike=1
1839
         else
1840
            call EastOrWest(tlonA,tlonB,ibrngAB)
1841
            call EastOrWest(tlonA,tlonP,ibrngAP)
1842
            call EastOrWest(tlonP,tlonB,ibrngPB)
1843
 
1844
 
1845
            if((ibrngAP.eq.ibrngAB).and.(ibrngPB.eq.ibrngAB)) istrike=1
1846
         endif
1847
 
1848
 
1849
         if (istrike.eq.1) then
1850
 
1851
            if ((plat.eq.vAlat).and.(plon.eq.vAlon)) then
1852
               location=2       ! P lies on a vertex of S
1853
               return
1854
            endif
1855
            call TrnsfmLon(vAlat,vAlon,xlat_c,xlon_c,tlon_X)
1856
            call TrnsfmLon(vAlat,vAlon,vBlat,vBlon,tlon_B)
1857
            call TrnsfmLon(vAlat,vAlon,plat,plon,tlon_P)
1858
 
1859
            if (tlon_P.eq.tlon_B) then
1860
               location=2       ! P lies on side of S
1861
               return 
1862
            else
1863
               call EastOrWest(tlon_B,tlon_X,ibrng_BX)
1864
               call EastOrWest(tlon_B,tlon_P,ibrng_BP)
1865
               if(ibrng_BX.eq.(-ibrng_BP)) icross=icross+1
1866
            endif
1867
 
1868
         endif
1869
      enddo                     ! end of loop over the sides of S
1870
 
1871
 
1872
c     if the arc XP crosses the boundary S an even number of times then P
1873
c     is in S
1874
 
1875
      if (mod(icross,2).eq.0) location=1
1876
 
1877
      return
1878
 
1879
      end
1880
 
1881
 
1882
c     ****************************************************************
1883
 
1884
      subroutine TrnsfmLon(plat,plon,qlat,qlon,tranlon)
1885
 
1886
c     ****************************************************************
1887
c     * This subroutine is required by subroutines DefSPolyBndry &   *
1888
c     * LctPtRelBndry. It finds the 'longitude' of point Q in a      *
1889
c     * geographic coordinate system for which point P acts as a     *
1890
c     * 'north pole'. SENT: plat,plon,qlat,qlon, in degrees.         *
1891
c     * RETURNED: tranlon, in degrees.                               *
1892
c     ****************************************************************
1893
 
1894
      implicit none
1895
 
1896
      real pi,dtr,plat,plon,qlat,qlon,tranlon,t,b
1897
      parameter (pi=3.141592654,dtr=pi/180.0)
1898
 
1899
      if (plat.eq.90.) then
1900
         tranlon=qlon
1901
      else
1902
         t=sin((qlon-plon)*dtr)*cos(qlat*dtr)
1903
         b=sin(dtr*qlat)*cos(plat*dtr)-cos(qlat*dtr)*sin(plat*dtr)
1904
     >    *cos((qlon-plon)*dtr)
1905
         tranlon=atan2(t,b)/dtr
1906
      endif
1907
 
1908
      return
1909
      end
1910
 
1911
c     ****************************************************************
1912
 
1913
      subroutine EastOrWest(clon,dlon,ibrng)
1914
 
1915
c     ****************************************************************
1916
c     * This subroutine is required by subroutine LctPtRelBndry.     *
1917
c     * This routine determines if in travelling the shortest path   *
1918
c     * from point C (at longitude clon) to point D (at longitude    *
1919
c     * dlon) one is heading east, west or neither.                  *
1920
c     * SENT: clon,dlon; in degrees. RETURNED: ibrng                 *
1921
c     * (1=east,-1=west, 0=neither).                                 *
1922
c     ****************************************************************
1923
 
1924
      implicit none
1925
      real clon,dlon,del
1926
      integer ibrng
1927
      del=dlon-clon
1928
      if (del.gt.180.) del=del-360.
1929
      if (del.lt.-180.) del=del+360.
1930
      if ((del.gt.0.0).and.(del.ne.180.)) then
1931
         ibrng=-1               ! (D is west of C)
1932
      elseif ((del.lt.0.0).and.(del.ne.-180.)) then
1933
         ibrng=+1               ! (D is east of C)
1934
      else
1935
         ibrng=0                ! (D north or south of C)
1936
      endif
1937
      return
1938
      end