Subversion Repositories lagranto.um

Rev

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