Subversion Repositories lagranto.icon

Rev

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