Subversion Repositories lagranto.icon

Rev

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