Subversion Repositories lagranto.icon

Rev

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

Rev Author Line No. Line
3 michaesp 1
      PROGRAM density
2
 
3
      use netcdf
4
 
5
      implicit none
6
 
7
c     ---------------------------------------------------------------------
8
c     Declaration of variables
9
c     ---------------------------------------------------------------------
10
 
11
c     Parameter and working arrays
12
      real                                    radius
13
      character*80                            runit
14
      integer                                 nx,ny
15
      integer                                 nlonlat
16
      real                                    dlonlat
17
      real                                    xmin,ymin,dx,dy
18
      real                                    clon,clat
19
      integer                                 ntime,nfield,ntra
20
      character*80                            inpfile
21
      character*80                            outfile
22
      character*80                            mode
23
      real                                    param
24
      integer                                 opts,npts
25
      integer                                 step
26
      character*80                            gridtype
27
      character*80                            field
28
      integer                                 crefile,crevar
29
      real                                    pollon,pollat,polgam
30
      real,allocatable,    dimension (:,:) :: cnt,res,fld,area
31
      real,allocatable,    dimension (:)   :: traj
32
      real,allocatable,    dimension (:)   :: olon,olat,otim,ofld
33
      real,allocatable,    dimension (:)   :: nlon,nlat,ntim,nfld
34
 
35
c     Output format
36
      character*80                            outformat
37
 
38
c     Physical and mathematical constants
39
      real                                    pi180
40
      parameter                               (pi180=3.14159/180.)
41
      real                                    deltay
42
      parameter                               (deltay=111.)
43
      real                                    eps
44
      parameter                               (eps=0.001)
45
 
46
c     Input trajectories (see iotra.f)
47
      integer                                 inpmode
48
      real,allocatable, dimension (:,:,:) ::  trainp     
49
      integer                                 reftime(6)      
50
      character*80                            varsinp(100)   
51
      integer,allocatable, dimension (:) ::   sel_flag
52
      character*80                            sel_file
53
      character*80                            sel_format
54
 
55
c     Auxiliary variables
56
      character*80                            cdfname,varname
57
      integer                                 i,j,k
58
      integer                                 stat
59
      integer,allocatable, dimension (:,:) :: connect0
60
      integer                                 connectval0
61
      integer,allocatable, dimension (:,:) :: connect1
62
      integer                                 connectval1
63
      integer,allocatable, dimension (:,:) :: connect2
64
      integer                                 connectval2
65
      real                                    slat
66
      integer                                 ipre
67
      real                                    addvalue
68
      real                                    xmax,ymax
69
      real ,allocatable, dimension (:)  ::    odist,ndist
70
      real                                    dt
71
      integer                                 fid
72
      integer                                 dynamic_grid
73
      real                                    ycen,xcen
74
      integer                                 indx,indy
75
      character*80                            unit
76
      real                                    rlon0,rlat0,rlon,rlat
77
      real                                    lon,lat
78
      real                                    crot
79
      integer                                 count
80
      character*80                            longname, varunit
81
      real                                    time
82
      integer                                 ind
83
      integer                                 ifield
84
      real                                    hhmm,frac
85
      integer                                 ierr,ncID
86
      integer                                 rotflag
87
 
88
c     External functions
89
      real         lmstolm,lmtolms
90
      real         phstoph,phtophs
91
      external     lmstolm,lmtolms,phstoph,phtophs
92
 
93
      real         sdis
94
      external     sdis
95
 
96
c     ---------------------------------------------------------------------
97
c     Preparations
98
c     ---------------------------------------------------------------------
99
 
100
c     Write start message
101
      print*,'========================================================='
102
      print*,'              *** START OF PROGRAM DENSITY ***'
103
      print*
104
 
105
c     Read input parameters
106
      open(10,file='density.param')
107
       read(10,*) inpfile
108
       read(10,*) outfile
109
       read(10,*) field
110
       read(10,*) ntime,nfield,ntra
111
       read(10,*) gridtype
112
       if ( gridtype.eq.'latlon' ) then 
113
          read(10,*) nx,ny,xmin,ymin,dx,dy
114
       elseif ( gridtype.eq.'centered') then
115
          read(10,*) clon,clat,nlonlat,dlonlat
116
       else
117
          print*,' ERROR: unsupported grid type ',trim(gridtype)
118
          stop
119
       endif
120
       read(10,*) radius,runit
121
       read(10,*) mode
122
       read(10,*) param
123
       read(10,*) step
124
       read(10,*) sel_file
125
       read(10,*) sel_format
126
       read(10,*) crefile
127
       read(10,*) crevar
128
       read(10,*) pollon,pollat
129
      close(10)
130
 
131
c	  Set the flag for trajectory rotation
132
      if ( abs(pollat-90.).gt.eps ) then
133
          rotflag = 1
134
      else
135
          rotflag = 0
136
      endif
137
 
138
c     Get the grid parameters if <crefile=0>
139
      if ( crefile.eq.0 ) then
140
 
141
           ierr = nf90_open  (trim(outfile), NF90_NOWRITE  , ncID)
142
 
143
           ierr = nf90_get_att(ncID, NF90_GLOBAL, 'grid'   ,gridtype ) 
144
           ierr = nf90_get_att(ncID, NF90_GLOBAL, 'clon'   ,clon     )
145
           ierr = nf90_get_att(ncID, NF90_GLOBAL, 'clat'   ,clat     )
146
           ierr = nf90_get_att(ncID, NF90_GLOBAL, 'nlonlat',nlonlat  )
147
           ierr = nf90_get_att(ncID, NF90_GLOBAL, 'dlonlat',dlonlat  )
148
           ierr = nf90_get_att(ncID, NF90_GLOBAL, 'nx'     ,nx       )
149
           ierr = nf90_get_att(ncID, NF90_GLOBAL, 'ny'     ,ny       )
150
           ierr = nf90_get_att(ncID, NF90_GLOBAL, 'dx'     ,dx       )
151
           ierr = nf90_get_att(ncID, NF90_GLOBAL, 'dy'     ,dy       )
152
           ierr = nf90_get_att(ncID, NF90_GLOBAL, 'xmin'   ,xmin     )
153
           ierr = nf90_get_att(ncID, NF90_GLOBAL, 'ymin'   ,ymin     )
154
 
155
           ierr = nf90_close(ncID)
156
 
157
           print*,'**** GRID PARAMETERS IMPORTED ',
158
     >            'FROM NETCDF FILE!!!! ****'
159
           print*
160
 
161
      endif
162
 
163
c     Check for consistency
164
      if ( (step.ne.0).and.(mode.ne.'keep') ) then
165
         print*," ERROR: interpolation is only possible for all",
166
     >                   ' time steps... Stop'
167
         stop
168
      endif
169
 
170
c     Set the number of times (just code aesthetics)
171
      opts=ntime
172
 
173
c     Set grid parameters for centered grid
174
      if ( gridtype.eq.'centered' ) then
175
         nx   = nlonlat
176
         ny   = nlonlat
177
         dx   = dlonlat
178
         dy   = dlonlat
179
         xmin = - real(nlonlat-1)/2. * dx
180
         xmax = + real(nlonlat-1)/2. * dx
181
         ymin = - real(nlonlat-1)/2. * dy
182
         ymax = + real(nlonlat-1)/2. * dy
183
      endif
184
 
185
c     Set the flag for dynamic grid adjustment
186
      if ( (nx.eq.0).or.(ny.eq.0) ) then
187
         dynamic_grid = 1
188
      else
189
         dynamic_grid = 0
190
      endif
191
 
192
c     Print status information
193
      print*,'---- INPUT PARAMETERS -----------------------------------'
194
      print* 
195
      print*,'Input                : ',trim(inpfile)
196
      print*,'Output               : ',trim(outfile)
197
      print*,'Field                : ',trim(field)
198
      print*,'Trajectory           : ',ntime,nfield,ntra
199
      print*,'Grid type            : ',trim(gridtype)
200
      if ( dynamic_grid.eq.1 ) then
201
         print*,'Grid                 : dynamic (see below)'
202
      elseif ( gridtype.eq.'latlon' ) then
203
         print*,'Grid   nlon,nlat     : ',nx,ny
204
         print*,'       lonmin,latmin : ',xmin,ymin
205
         print*,'       dlon,dlat     : ',dx,dy
206
      elseif ( gridtype.eq.'centered' ) then
207
         print*,'Grid   clon,clat     : ',clon,clat
208
         print*,'       nlonlat       : ',nlonlat
209
         print*,'       dlonlat       : ',dlonlat
210
      endif
211
      if ( rotflag.eq.0 ) then
212
          print*,'Pole lon/lat         : ',pollon,pollat
213
      else
214
          print*,'Pole lon/lat         : ',pollon,pollat,' (ROTATION)'
215
      endif
216
      print*,'Filter radius        : ',radius,' ',trim(runit)
217
      print*,'Mode                 : ',trim(mode)
218
      if ( ( mode.eq.'time'  ).or.
219
     >     ( mode.eq.'space' ).or.
220
     >     ( mode.eq.'grid' ) )
221
     >then
222
         print*,'Parameter            : ',param
223
      endif
224
      if ( step.eq.0 ) then
225
         print*,'Time step            : all'
226
      elseif (step.gt.0) then
227
         print*,'Time step            : ',step
228
      endif
229
      print*,'Selection file       : ',trim(sel_file)
230
      print*,'Selection format     : ',trim(sel_file)
231
      print*,'Flag <crefile>       : ',crefile
232
      print*,'Flag <crevar>        : ',crevar
233
 
234
c     Check whether mode is valid
235
      if ((mode.ne.'keep'  ).and.
236
     >    (mode.ne.'time'  ).and.
237
     >    (mode.ne.'space' ).and.
238
     >    (mode.ne.'grid'  ))  
239
     >then
240
         print*,' ERROR: Invalid mode ',trim(mode)
241
         stop
242
      endif
243
 
244
c     Allocate memory for old and new (reparameterised) trajectory
245
      allocate(olon(ntime),stat=stat)
246
      if (stat.ne.0) print*,'*** error allocating array olon ***'
247
      allocate(olat(ntime),stat=stat)
248
      if (stat.ne.0) print*,'*** error allocating array olat ***'
249
      allocate(otim(ntime),stat=stat)
250
      if (stat.ne.0) print*,'*** error allocating array otim ***'
251
      allocate(nlon(1000*ntime),stat=stat)
252
      if (stat.ne.0) print*,'*** error allocating array nlon ***'
253
      allocate(nlat(1000*ntime),stat=stat)
254
      if (stat.ne.0) print*,'*** error allocating array nlat ***'
255
      allocate(ntim(1000*ntime),stat=stat)
256
      if (stat.ne.0) print*,'*** error allocating array ntim ***'
257
      allocate(odist(ntime),stat=stat)
258
      if (stat.ne.0) print*,'*** error allocating array odist ***'
259
      allocate(ndist(1000*ntime),stat=stat)
260
      if (stat.ne.0) print*,'*** error allocating array ndist ***'
261
      allocate(ofld(ntime),stat=stat)
262
      if (stat.ne.0) print*,'*** error allocating array ofld ***'
263
      allocate(nfld(1000*ntime),stat=stat)
264
      if (stat.ne.0) print*,'*** error allocating array nfld ***'
265
 
266
c     Allocate memory for complete trajectory set
267
      allocate(trainp(ntra,ntime,nfield),stat=stat)
268
      if (stat.ne.0) print*,'*** error allocating array trainp ***'
269
      allocate(sel_flag(ntra),stat=stat)
270
      if (stat.ne.0) print*,'*** error allocating array sel_flag ***'
271
 
272
c     Allocate memory for auxiliary fields
273
      allocate(traj(nfield),stat=stat)
274
      if (stat.ne.0) print*,'*** error allocating array traj ***'
275
 
276
c     Set the format of the input file
277
      call mode_tra(inpmode,inpfile)
278
      if (inpmode.eq.-1) inpmode=1
279
 
280
c     Read the input trajectory file
281
      call ropen_tra(fid,inpfile,ntra,ntime,nfield,
282
     >                   reftime,varsinp,inpmode)
283
      call read_tra (fid,trainp,ntra,ntime,nfield,inpmode)
284
      call close_tra(fid,inpmode)
285
 
286
c     Check that first four columns correspond to time,lon,lat,p
287
      if ( (varsinp(1).ne.'time' ).or.
288
     >     (varsinp(2).ne.'xpos' ).and.(varsinp(2).ne.'lon' ).or.
289
     >     (varsinp(3).ne.'ypos' ).and.(varsinp(3).ne.'lat' ).or.
290
     >     (varsinp(4).ne.'zpos' ).and.(varsinp(4).ne.'z'   ) )
291
     >then
292
         print*,' ERROR: problem with input trajectories ...'
293
         stop
294
      endif
295
      varsinp(1) = 'TIME'
296
      varsinp(2) = 'lon'
297
      varsinp(3) = 'lat'
298
      varsinp(4) = 'z'
299
 
300
c     Get the index of the field (if needed)
301
      if ( field.ne.'nil' ) then
302
         ifield = 0
303
         do i=1,nfield
304
            if ( varsinp(i).eq.field ) ifield = i
305
         enddo
306
         if ( ifield.eq.0 ) then
307
            print*,' ERROR: field ',trim(field),' not found... Stop'
308
            stop
309
         endif
310
      endif
311
 
312
c     Rotate coordinates - if pollat<>90
313
      do i=1,ntra
314
        do j=1,ntime
315
 
316
            lon = trainp(i,j,2)
317
            lat = trainp(i,j,3)
318
 
319
            if ( rotflag.eq.1 ) then
320
               rlon = lmtolms(lat,lon,pollat,pollon)
321
               rlat = phtophs(lat,lon,pollat,pollon)
322
            else
323
               rlon = lon
324
               rlat = lat
325
            endif
326
 
327
            trainp(i,j,2) = rlon
328
            trainp(i,j,3) = rlat
329
 
330
        enddo
331
      enddo
332
 
333
c     Write some status information of the input trajectories
334
      print*
335
      print*,'---- INPUT TRAJECTORIES ---------------------------------'
336
      print*
337
      print*,' Reference time (year)  : ',reftime(1)
338
      print*,'                (month) : ',reftime(2)
339
      print*,'                (day)   : ',reftime(3)
340
      print*,'                (hour)  : ',reftime(4)
341
      print*,'                (min)   : ',reftime(5)
342
      print*,' Time range (min)       : ',reftime(6)
343
      do i=1,nfield
344
         if ( i.ne.ifield ) then
345
            print*,' Var                    :',i,trim(varsinp(i))
346
         else
347
            print*,' Var                    :',i,trim(varsinp(i)),
348
     >                                        '       [ gridding ]'
349
         endif
350
      enddo
351
      print*,' List of selected times'
352
      do i=1,ntime
353
         if ( (step.eq.0).or.(step.eq.i) ) then
354
            print*,'     ',i,'  -> ',trainp(1,i,1)
355
         endif
356
      enddo
357
      print*
358
 
359
c     Select flag: all trajectories are selected
360
      if ( sel_file.eq.'nil' ) then
361
 
362
         do i=1,ntra
363
            sel_flag(i) = 1
364
         enddo
365
 
366
c     Select flag: index file
367
      elseif ( sel_format.eq.'index' ) then
368
 
369
         do i=1,ntra
370
            sel_flag(i) = 0
371
         enddo
372
 
373
         open(10,file=sel_file)
374
 142      read(10,*,end=141) ind
375
          sel_flag(ind) = 1
376
          goto 142 
377
 141     continue
378
         close(10)
379
 
380
c     Select flag: boolean file
381
      elseif ( sel_format.eq.'boolean' ) then
382
 
383
         open(10,file=sel_file)
384
          do i=1,ntra
385
            read(10,*) ind
386
            if ( ind.eq.1 ) sel_flag(i) = ind
387
          enddo
388
         close(10)
389
 
390
      endif
391
 
392
c     Write status information
393
      if ( sel_file.eq.'nil' ) then
394
          print*,' Selected trajectories  : all ',ntra          
395
       else
396
          count = 0
397
          do i=1,ntra
398
             if ( sel_flag(i).eq.1 ) count = count + 1
399
          enddo
400
          print*,' #selected trajectories : ',count,
401
     >            ' [ ',real(count)/real(ntra) * 100.,' % ] '
402
       endif
403
       print*
404
 
405
c     ---------------------------------------------------------------------
406
c     Coordinate transformations and grid adjustment
407
c     ---------------------------------------------------------------------
408
 
409
c     Transform from lat/lon to rotated lat/lon, if requested
410
      if ( gridtype.eq.'centered') then
411
 
412
         crot = 0.
413
 
414
         pollon=clon-180.
415
         if (pollon.lt.-180.) pollon=pollon+360.
416
         pollat=90.-clat
417
         do i=1,ntra
418
            do j=1,ntime
419
 
420
               if ( sel_flag(i).eq.1 ) then
421
 
422
c                Get lat/lon coordinates for trajectory point
423
                 lon = trainp(i,j,2)
424
                 lat = trainp(i,j,3)
425
 
426
c                First Rotation
427
                 pollon=clon-180.
428
                 if (pollon.lt.-180.) pollon=pollon+360.
429
                 pollat=90.-clat
430
                 rlon0=lmtolms(lat,lon,pollat,pollon)
431
                 rlat0=phtophs(lat,lon,pollat,pollon)            
432
 
433
c                Second rotation
434
                 pollon=-180.
435
                 pollat=90.+crot
436
                 rlon=90.+lmtolms(rlat0,rlon0-90.,pollat,pollon)
437
                 rlat=phtophs(rlat0,rlon0-90.,pollat,pollon)   
438
 
439
c                Get rotated latitude and longitude
440
 100             if (rlon.lt.xmin) then
441
                  rlon=rlon+360.
442
                  goto 100
443
                 endif
444
 102             if (rlon.gt.(xmin+real(nx-1)*dx)) then
445
                  rlon=rlon-360.
446
                  goto 102
447
                 endif
448
 
449
c                Set the new trajectory coordinates
450
                 trainp(i,j,2) = rlon
451
                 trainp(i,j,3) = rlat
452
 
453
              endif
454
 
455
            enddo
456
         enddo
457
 
458
      endif
459
 
460
c     Dynamic grid adjustment
461
      if ( dynamic_grid.eq.1 ) then
462
 
463
c        Get the grid parameters
464
         xmin =  180.
465
         ymin =   90.
466
         xmax = -180.
467
         ymax =  -90.
468
 
469
         do i=1,ntra
470
 
471
            if ( sel_flag(i).eq.1 ) then
472
 
473
              if ( step.eq.0 ) then
474
               do j=1,ntime
475
                  if ( trainp(i,j,2).lt.xmin) xmin =  trainp(i,j,2)
476
                  if ( trainp(i,j,2).gt.xmax) xmax =  trainp(i,j,2)
477
                  if ( trainp(i,j,3).lt.ymin) ymin =  trainp(i,j,3)
478
                  if ( trainp(i,j,3).gt.ymax) ymax =  trainp(i,j,3)
479
               enddo
480
              else
481
                if ( trainp(i,step,2).lt.xmin) xmin =  trainp(i,step,2)
482
                if ( trainp(i,step,2).gt.xmax) xmax =  trainp(i,step,2)
483
                if ( trainp(i,step,3).lt.ymin) ymin =  trainp(i,step,3)
484
                if ( trainp(i,step,3).gt.ymax) ymax =  trainp(i,step,3)
485
              endif
486
 
487
            endif
488
 
489
         enddo
490
 
491
c        Get first guess for "optimal" grid
492
         nx = 400
493
         ny = 400
494
         dx = (xmax - xmin)/real(nx-1)
495
         dy = (ymax - ymin)/real(ny-1)
496
 
497
c        Make the grid spacing equal in zonal and meridional direction
498
         if ( dx.gt.dy ) then
499
 
500
            dy = dx
501
            ny = (ymax - ymin)/dy + 1
502
            if (ny.lt.nx/2)              ny = nx / 2
503
            if ( real(ny)*dy .ge. 180. ) ny = 180./dy + 1
504
            ycen = 0.5* (ymin+ymax)
505
            ymin = ycen - 0.5 * real(ny/2) * dy
506
            if (ymin.le.-90.) ymin = -90.
507
 
508
         else
509
 
510
            dx = dy
511
            nx = (xmax - xmin)/dx + 1
512
            if (nx.lt.ny/2)              nx = ny / 2
513
            if ( real(nx)*dx .ge. 360. ) nx = 360./dx + 1
514
            xcen = 0.5* (xmin+xmax)
515
            xmin = xcen - 0.5 * real(nx/2) * dx
516
            if (xmin.le.-180.) xmin = -180.
517
 
518
         endif
519
 
520
c        Write information
521
         print*
522
         print*,'---- DYNAMIC GRID ADJUSTMENT',
523
     >          ' ----------------------------'  
524
         print*
525
         print*,'Grid   nlon,nlat     : ',nx,ny
526
         print*,'       lonmin,latmin : ',xmin,ymin
527
         print*,'       dlon,dlat     : ',dx,dy
528
         print*
529
 
530
c     Write grid information for rotated grid (if not already done
531
      elseif ( gridtype.eq.'centered') then
532
 
533
         print*
534
         print*,'---- GRID PARAMETERS -------',
535
     >          ' ----------------------------'  
536
         print*
537
         print*,'Grid   nlon,nlat     : ',nx,ny
538
         print*,'       lonmin,latmin : ',xmin,ymin
539
         print*,'       dlon,dlat     : ',dx,dy
540
         print*
541
 
542
 
543
      endif
544
 
545
c     Set the grid boundaries
546
      xmax=xmin+real(nx-1)*dx
547
      ymax=ymin+real(ny-1)*dy
548
 
549
c     Allocate memory for output array and auxiliary gridding array 
550
      allocate(cnt(nx,ny),stat=stat)
551
      if (stat.ne.0) print*,'*** error allocating array cnt  ***'
552
      allocate(res(nx,ny),stat=stat)
553
      if (stat.ne.0) print*,'*** error allocating array res  ***'
554
      allocate(fld(nx,ny),stat=stat)
555
      if (stat.ne.0) print*,'*** error allocating array fld  ***'
556
      allocate(area(nx,ny),stat=stat)
557
      if (stat.ne.0) print*,'*** error allocating array area ***'
558
 
559
      allocate(connect0(nx,ny),stat=stat)
560
      if (stat.ne.0) print*,'*** error allocating array connect0 ***'
561
      allocate(connect1(nx,ny),stat=stat)
562
      if (stat.ne.0) print*,'*** error allocating array connect1 ***'
563
      allocate(connect2(nx,ny),stat=stat)
564
      if (stat.ne.0) print*,'*** error allocating array connect2 ***'
565
 
566
 
567
c     Init the output array
568
      do i=1,nx
569
         do j=1,ny
570
            connect0(i,j) = 0
571
            connect1(i,j) = 0
572
            connect2(i,j) = 0
573
            cnt(i,j)      = 0.
574
            res(i,j)      = 0.
575
            fld(i,j)      = 0.
576
         enddo
577
      enddo  
578
 
579
c     ---------------------------------------------------------------------
580
c     Gridding
581
c     ---------------------------------------------------------------------
582
 
583
c     Write some status information 
584
      print*,'---- GRIDDING -------------------------------------------'
585
      print*
586
 
587
c     Loop over all entries of sampling table
588
      connectval0 = 0
589
      connectval1 = 0
590
      connectval2 = 0
591
      count       = 0
592
 
593
      do i=1,ntra
594
 
595
         if (mod(i,100).eq.0) print*,i,' of ',ntra
596
 
597
c        Skip all trajectories which are not selected
598
         if ( sel_flag(i).eq.0 ) goto 300
599
 
600
c        ------- Read a complete trajectory ---------------------------
601
         do j=1,ntime
602
            otim(j) = trainp(i,j,1)
603
            olon(j) = trainp(i,j,2)
604
            olat(j) = trainp(i,j,3)
605
            if ( field.ne.'nil' ) then
606
               ofld(j) =trainp(i,j,ifield)
607
            endif
608
         enddo
609
 
610
c        -------- Convert hh.m time into fractional time --------------
611
         do j=1,ntime
612
            hhmm    = otim(j)
613
            call hhmm2frac (hhmm,frac)
614
            otim(j) = frac
615
         enddo
616
 
617
c        -------- Interpolation ---------------------------------------
618
 
619
c        Keep the trajectory points as they are
620
         if ( ( mode.eq.'keep').and.(step.eq.0) ) then
621
            npts=opts
622
            do j=1,opts
623
               ntim(j)=otim(j)
624
               nlon(j)=olon(j)
625
               nlat(j)=olat(j)
626
               if ( field.ne.'nil' ) then
627
                  nfld(j)=ofld(j)
628
               endif
629
            enddo
630
 
631
c        Select a single time step
632
         elseif ( ( mode.eq.'keep').and.(step.gt.0) ) then
633
            npts    = 1
634
            ntim(1) = otim(step)
635
            nlon(1) = olon(step)
636
            nlat(1) = olat(step)
637
            if ( field.ne.'nil' ) then
638
               nfld(1) = ofld(step)
639
            endif
640
 
641
c        Perform a reparameterisation in time
642
         else if ( (mode.eq.'time').and.(step.eq.0) ) then
643
 
644
c           Get the new number of trajectory points
645
            npts=nint(abs(otim(opts)-otim(1))/param)+1
646
 
647
c           Handle date line problem
648
            do j=2,opts
649
               if ( (olon(j-1)-olon(j)).gt.180. ) then
650
                  olon(j) = olon(j) + 360.
651
               else if ( (olon(j-1)-olon(j)).lt.-180. ) then
652
                  olon(j) = olon(j) - 360.
653
               endif
654
            enddo
655
 
656
c           Cubic spline fitting
657
            call curvefit(otim,olon,opts,ntim,nlon,npts)
658
            call curvefit(otim,olat,opts,ntim,nlat,npts)
659
            if ( field.ne.'nil' ) then
660
               call curvefit(otim,ofld,opts,ntim,nfld,npts)
661
            endif
662
 
663
c           Reverse date line handling
664
            do j=1,npts
665
               if ( nlon(j).gt.xmax ) then
666
                  nlon(j) = nlon(j) -360.
667
               else if ( nlon(j).lt.xmin ) then
668
                  nlon(j) = nlon(j) +360.
669
               endif
670
            enddo
671
 
672
c        Perform a reparameterisation with equally spaced gridpoint
673
         elseif ( (mode.eq.'space').and.(step.eq.0) ) then
674
 
675
c           Calculate the distance and spacing
676
            odist(1) = 0.
677
            unit     = 'km'
678
            do j=2,ntime
679
               odist(j)=odist(j-1) + 
680
     >                  sdis(olon(j-1),olat(j-1),olon(j),olat(j),unit)
681
            enddo
682
 
683
c           Determine the new number of trajectory points
684
            npts=nint(odist(ntime)/param)+1
685
            if (npts.eq.0) then
686
               npts=1.
687
            endif
688
 
689
c           Handle date line problem
690
            do j=2,opts
691
               if ( (olon(j-1)-olon(j)).gt.180. ) then
692
                  olon(j) = olon(j) + 360.
693
               else if ( (olon(j-1)-olon(j)).lt.-180. ) then
694
                  olon(j) = olon(j) - 360.
695
               endif
696
            enddo
697
 
698
c           Cubic spline fitting
699
            call curvefit(odist,olon,opts,ndist,nlon,npts)
700
            call curvefit(odist,olat,opts,ndist,nlat,npts)
701
            call curvefit(odist,otim,opts,ndist,ntim,npts)
702
            if ( field.ne.'nil' ) then
703
               call curvefit(odist,ofld,opts,ndist,nfld,npts)
704
            endif
705
 
706
c           Reverse date line handling
707
            do j=1,npts
708
               if ( nlon(j).gt.xmax ) then
709
                  nlon(j) = nlon(j) -360.
710
               else if ( nlon(j).lt.xmin ) then
711
                  nlon(j) = nlon(j) +360.
712
               endif
713
            enddo
714
 
715
c        Perform a reparameterisation with equally spaced gridpoint
716
         elseif ( (mode.eq.'grid').and.(step.eq.0) ) then
717
 
718
c           Calculate the distance and spacing
719
            odist(1) = 0.
720
            unit     = 'deg'
721
            do j=2,ntime
722
               odist(j)=odist(j-1) + 
723
     >                  sdis(olon(j-1),olat(j-1),olon(j),olat(j),unit)
724
            enddo
725
 
726
c           Determine the new number of trajectory points
727
            npts=nint(odist(ntime)/param)+1
728
            if (npts.eq.0) then
729
               npts=1.
730
            endif
731
 
732
c           Handle date line problem
733
            do j=2,opts
734
               if ( (olon(j-1)-olon(j)).gt.180. ) then
735
                  olon(j) = olon(j) + 360.
736
               else if ( (olon(j-1)-olon(j)).lt.-180. ) then
737
                  olon(j) = olon(j) - 360.
738
               endif
739
            enddo
740
 
741
c           Cubic spline fitting
742
            call curvefit(odist,olon,opts,ndist,nlon,npts)
743
            call curvefit(odist,olat,opts,ndist,nlat,npts)
744
            call curvefit(odist,otim,opts,ndist,ntim,npts)
745
            if ( field.ne.'nil' ) then
746
               call curvefit(odist,ofld,opts,ndist,nfld,npts)
747
            endif
748
 
749
c           Reverse date line handling
750
            do j=1,npts
751
               if ( nlon(j).gt.xmax ) then
752
                  nlon(j) = nlon(j) -360.
753
               else if ( nlon(j).lt.xmin ) then
754
                  nlon(j) = nlon(j) +360.
755
               endif
756
            enddo
757
 
758
         endif
759
 
760
c        -------- Do the gridding -------------------------------------
761
 
762
c        Gridding of trajectory
763
         do j=1,npts
764
 
765
c           Check whether point is in data domain
766
	    if ( (nlon(j).gt.xmin).and.(nlon(j).lt.xmax).and.
767
     >           (nlat(j).gt.ymin).and.(nlat(j).lt.ymax))
768
     >      then
769
 
770
c              Increase counter for gridded points
771
               count = count + 1
772
 
773
c              ----------------- Gridding: simple count -----------------
774
               connectval0 = connectval0+1
775
 
776
               addvalue    = 1.
777
 
778
               call  gridding1
779
     >              (nlat(j),nlon(j),addvalue,
780
     >               radius,runit,connect0,connectval0,
781
     >               cnt,nx,ny,xmin,ymin,dx,dy)
782
 
783
c              ----------------- Gridding: residence time ---------------
784
               connectval1 = connectval1+1
785
 
786
               if ( ntime.eq.1 ) then
787
                  addvalue = 0.
788
               elseif ( j.eq.1 )  then
789
                  addvalue=abs(ntim(2)-ntim(1))
790
               else
791
                  addvalue=abs(ntim(j)-ntim(j-1))
792
               endif
793
 
794
               call  gridding1
795
     >              (nlat(j),nlon(j),addvalue,
796
     >               radius,runit,connect1,connectval1,
797
     >               res,nx,ny,xmin,ymin,dx,dy)
798
 
799
 
800
c              --------------- Gridding: field -------------------------
801
               if ( field.ne.'nil' ) then
802
 
803
                   connectval2 = connectval2+1
804
 
805
                   addvalue    = nfld(j)
806
 
807
                   call  gridding1
808
     >                  (nlat(j),nlon(j),addvalue,
809
     >                  radius,runit,connect2,connectval2,
810
     >                  fld,nx,ny,xmin,ymin,dx,dy)
811
 
812
               endif
813
 
814
	    endif
815
 
816
         enddo
817
 
818
c        Exit point for loop over all trajectories
819
 300     continue
820
 
821
      enddo
822
 
823
c     Write status information
824
      print*
825
      print*,' # gridded points       : ',count
826
 
827
c     ---------------------------------------------------------------------
828
c     Unit conversions and output to netCDF file
829
c     ---------------------------------------------------------------------
830
 
831
c     Write some status information 
832
      print*
833
      print*,'---- WRITE OUTPUT ---------------------------------------'
834
      print*
835
 
836
c     Area (in km^2)
837
      do i=1,nx	         
838
         do j=1,ny	
839
            slat=ymin+real(j-1)*dy
840
            if (abs(abs(slat)-90.).gt.eps) then
841
               area(i,j) = dy*dx*cos(pi180*slat)*deltay**2
842
            else
843
               area(i,j) = 0.
844
            endif
845
         enddo
846
      enddo
847
 
848
c     Normalise gridded field
849
      if ( field.ne.'nil' ) then
850
         do i=1,nx
851
            do j=1,ny
852
               if ( cnt(i,j).gt.0. ) then
853
                  fld(i,j) = fld(i,j) / cnt(i,j)
854
               endif
855
            enddo
856
         enddo
857
      endif
858
 
859
c     Set the time for the output netCDF files - if a composite is
860
c     calculatd, then the time is set to 
861
      if ( step.eq.0 ) then
862
         time = -999.
863
         print*,'   ... COMPOSITE OVER ALL TRAJECTORY TIMES (-999)'
864
         print*
865
      else
866
         time = trainp(1,step,1)
867
      endif
868
 
869
c     Write output to CF netCDF
870
      cdfname  = outfile
871
 
872
      varname  = 'COUNT'
873
      longname = 'trajectory counts'
874
      varunit  = 'counts per grid point'
875
      call  writecdf2D_cf (cdfname,varname,longname,varunit,gridtype,
876
     >       clon,clat,nlonlat,dlonlat,cnt,time,dx,dy,xmin,ymin,nx,
877
     >       ny,crefile,crefile,1,pollon,pollat)
878
      write(*,'(a8,a10,a5,a10,a10,f7.2,a2)') 
879
     >     '    ... ',trim(varname),' -> ',trim(cdfname),
880
     >     ' [ time = ',time,' ]'    
881
 
882
      varname  = 'RESIDENCE'
883
      longname = 'residence time'
884
      varunit  = 'hours per grid point'
885
      call  writecdf2D_cf (cdfname,varname,longname,varunit,gridtype,
886
     >       clon,clat,nlonlat,dlonlat,res,time,dx,dy,xmin,ymin,nx,
887
     >       ny,0,crefile,1,pollon,pollat)
888
      write(*,'(a8,a10,a5,a10,a10,f7.2,a2)') 
889
     >     '    ... ',trim(varname),' -> ',trim(cdfname),
890
     >     ' [ time = ',time,' ]'    
891
 
892
      varname  = 'AREA'
893
      longname = 'area corresponding to grid points'
894
      varunit  = 'square kilometers'
895
      call  writecdf2D_cf (cdfname,varname,longname,varunit,gridtype,
896
     >       clon,clat,nlonlat,dlonlat,area,time,dx,dy,xmin,ymin,nx,
897
     >       ny,0,crefile,1,pollon,pollat)
898
      write(*,'(a8,a10,a5,a10,a10,f7.2,a2)') 
899
     >     '    ... ',trim(varname),' -> ',trim(cdfname),
900
     >     ' [ time = ',time,' ]'    
901
 
902
      if ( field.ne.'nil' ) then
903
         varname  = field
904
         longname = field
905
         varunit  = 'as on trajectory file'
906
         call  writecdf2D_cf (cdfname,varname,longname,varunit,gridtype,
907
     >       clon,clat,nlonlat,dlonlat,fld,time,dx,dy,xmin,ymin,nx,
908
     >       ny,0,crevar,1,pollon,pollat)
909
 
910
         write(*,'(a8,a10,a5,a10,a10,f7.2,a2)') 
911
     >        '    ... ',trim(varname),' -> ',trim(cdfname),
912
     >        ' [ time = ',time,' ]'    
913
      endif
914
 
915
c     Write status information
916
      print*
917
      print*,'              *** END OF PROGRAM DENSITY **'
918
      print*,'========================================================='
919
 
920
      end
921
 
922
c     ********************************************************************
923
c     * GRIDDING SUBROUTINES                                             *
924
c     ********************************************************************
925
 
926
c     ---------------------------------------------------------------------
927
c     Gridding of one single data point (smoothing in km, deg, gridp)
928
c     ---------------------------------------------------------------------
929
 
930
      subroutine gridding1 (lat,lon,addval,radius,unit,
931
     >                      connect,connectval,
932
     >                      out,nx,ny,xmin,ymin,dx,dy)
933
 
934
      implicit none
935
 
936
c     Declaration of subroutine parameters
937
      real         lat,lon
938
      integer      nx,ny
939
      real         xmin,ymin,dx,dy
940
      real         out(nx,ny)
941
      real         radius
942
      character*80 unit
943
      integer      connectval
944
      integer      connect(nx,ny)
945
      real         addval
946
 
947
c     Auxiliary variables
948
      integer   i,j,k
949
      integer   mu,md,nr,nl,n,m
950
      integer   stackx(nx*ny),stacky(nx*ny)
951
      integer   tab_x(nx*ny),tab_y(nx*ny)
952
      real      tab_r(nx*ny)
953
      integer   sp
954
      real      lat2,lon2
955
      real      dist,sum
956
      real      xmax
957
      integer   periodic
958
      integer   test
959
 
960
c     Numerical epsilon
961
      real      eps
962
      parameter (eps=0.01)
963
 
964
c     Externals
965
      real      sdis,weight
966
      external  sdis,weight
967
 
968
c     Check whether lat/lon point is valid
969
      xmax=xmin+real(nx-1)*dx
970
      if (lon.lt.xmin-eps) lon=lon+360.
971
      if (lon.gt.xmax+eps) lon=lon-360.
972
      if (abs(lat-90).lt.eps) lat=90.
973
      if (abs(lat+90).lt.eps) lat=-90.
974
      if ((abs(lat).gt.(90.+eps)).or.
975
     >    (lon.lt.xmin-eps).or.(lon.gt.xmax+eps)) then
976
         print*,'Invalid lat/lon point ',lat,lon
977
         return
978
      endif
979
 
980
c     Set flag for periodic domain
981
      if (abs(xmax-xmin-360.).lt.eps) then
982
         periodic=1
983
      else if (abs(xmax-xmin-360+dx).lt.eps) then
984
         periodic=2
985
      else
986
         periodic=0
987
      endif
988
 
989
c     Get indices of one coarse grid point within search radius
990
      i=nint((lon-xmin)/dx)+1
991
      if ((i.eq.nx).and.(periodic.eq.1)) i=1
992
      j=nint((lat-ymin)/dy)+1
993
      lat2=ymin+real(j-1)*dy
994
      lon2=xmin+real(i-1)*dx
995
      dist=sdis(lon,lat,lon2,lat2,unit)
996
      if (dist.gt.radius) then
997
         print*,'1: Search radius is too small...'
998
         stop
999
      endif
1000
 
1001
c     Get connected points
1002
      k=0
1003
      stackx(1)=i
1004
      stacky(1)=j
1005
      sp=1
1006
      do while (sp.ne.0) 
1007
 
1008
c        Get an element from stack
1009
         n=stackx(sp)
1010
         m=stacky(sp)
1011
         sp=sp-1
1012
 
1013
c        Get distance from reference point
1014
         lat2=ymin+real(m-1)*dy
1015
         lon2=xmin+real(n-1)*dx
1016
         dist=sdis(lon,lat,lon2,lat2,unit)
1017
 
1018
c        Check whether distance is smaller than search radius: connected
1019
         if (dist.lt.radius) then
1020
 
1021
c           Make entry in filter mask
1022
            k=k+1
1023
            tab_x(k)=n
1024
            tab_y(k)=m
1025
            tab_r(k)=weight(dist,radius)
1026
 
1027
c           Mark this point as visited
1028
            connect(n,m)=connectval
1029
 
1030
c           Get coordinates of neighbouring points
1031
            nr=n+1
1032
            if ((nr.gt.nx)  .and.(periodic.eq.0)) nr=nx
1033
            if ((nr.gt.nx-1).and.(periodic.eq.1)) nr=1
1034
            if ((nr.gt.nx)  .and.(periodic.eq.2)) nr=1
1035
            nl=n-1
1036
            if ((nl.lt.1).and.(periodic.eq.0)) nl=1
1037
            if ((nl.lt.1).and.(periodic.eq.1)) nl=nx-1
1038
            if ((nl.lt.1).and.(periodic.eq.2)) nl=nx
1039
            mu=m+1
1040
            if (mu.gt.ny) mu=ny
1041
            md=m-1
1042
            if (md.lt.1) md=1
1043
 
1044
c           Update stack
1045
            if (connect(nr,m).ne.connectval) then
1046
               connect(nr,m)=connectval
1047
               sp=sp+1
1048
               stackx(sp)=nr
1049
               stacky(sp)=m
1050
            endif
1051
            if (connect(nl,m).ne.connectval) then
1052
               connect(nl,m)=connectval
1053
               sp=sp+1
1054
               stackx(sp)=nl
1055
               stacky(sp)=m
1056
            endif
1057
            if (connect(n,mu).ne.connectval) then
1058
               connect(n,mu)=connectval
1059
               sp=sp+1
1060
               stackx(sp)=n
1061
               stacky(sp)=mu
1062
            endif
1063
            if (connect(n,md).ne.connectval) then
1064
               connect(n,md)=connectval
1065
               sp=sp+1
1066
               stackx(sp)=n
1067
               stacky(sp)=md
1068
            endif
1069
         endif
1070
 
1071
      end do
1072
 
1073
      if (k.ge.1) then
1074
         sum=0.
1075
         do i=1,k
1076
            sum=sum+tab_r(i)
1077
         enddo
1078
         do i=1,k
1079
            out(tab_x(i),tab_y(i))=out(tab_x(i),tab_y(i))+
1080
     >                             addval*tab_r(i)/sum
1081
 
1082
            if ((tab_x(i).eq.1).and.(periodic.eq.1)) then
1083
               out(nx,tab_y(i))=out(nx,tab_y(i))+
1084
     >                             addval*tab_r(i)/sum
1085
            endif
1086
         enddo
1087
      else
1088
         print*,'2: Search radius is too small...'
1089
         stop
1090
      endif
1091
 
1092
      end
1093
 
1094
 
1095
c     ----------------------------------------------------------------------
1096
c     Get spherical distance between lat/lon points
1097
c     ----------------------------------------------------------------------
1098
 
1099
      real function sdis(xp,yp,xq,yq,unit)
1100
 
1101
c     Calculates spherical distance (in km) between two points given
1102
c     by their spherical coordinates (xp,yp) and (xq,yq), respectively.
1103
 
1104
      real         re
1105
      parameter    (re=6370.)
1106
      real         xp,yp,xq,yq,arg
1107
      character*80 unit
1108
      real         dlon
1109
 
1110
      if ( unit.eq.'km' ) then
1111
 
5 michaesp 1112
         arg=sin(yp*pi180)*sin(yq*pi180)+cos(yp*pi180)*cos(yq*pi180)*
1113
     >       cos((xp-xq)*pi180)
1114
c         arg=sind(yp)*sind(yq)+cosd(yp)*cosd(yq)*cosd(xp-xq)
3 michaesp 1115
         if (arg.lt.-1.) arg=-1.
1116
         if (arg.gt.1.) arg=1.
1117
         sdis=re*acos(arg)
1118
 
1119
      elseif ( unit.eq.'deg' ) then
1120
 
1121
         dlon = xp-xq
1122
         if ( dlon.gt. 180. ) dlon = dlon - 360.
1123
         if ( dlon.lt.-180. ) dlon = dlon + 360.
1124
         sdis = sqrt( dlon**2 + (yp-yq)**2 )
1125
 
1126
      endif
1127
 
1128
 
1129
c     Quick and dirty trick to avoid zero distances
1130
      if (sdis.eq.0.) sdis=0.1
1131
 
1132
      end
1133
 
1134
c     ----------------------------------------------------------------------
1135
c     Weight function for the filter mask
1136
c     ----------------------------------------------------------------------
1137
 
1138
      real function weight (r,radius)
1139
 
1140
c     Attribute to each distanc r its corresponding weight in the filter mask
1141
 
1142
      implicit none
1143
 
1144
c     Declaration of subroutine parameters
1145
      real r
1146
      real radius
1147
 
1148
c     Simple 0/1 mask
1149
      if (r.lt.radius) then
1150
         weight=exp(-r/radius)
1151
      else
1152
         weight=0.
1153
      endif
1154
 
1155
      end
1156
 
1157
 
1158
c     ********************************************************************
1159
c     * REPARAMETERIZATION SUBROUTINES                                   *
1160
c     ********************************************************************
1161
 
1162
c     -------------------------------------------------------------
1163
c     Interpolation of the trajectory with a natural cubic spline
1164
c     -------------------------------------------------------------
1165
 
1166
      SUBROUTINE curvefit (time,lon,n,
1167
     >                     sptime,splon,spn)
1168
 
1169
c     Given the curve <time,lon> with <n> data points, fit a
1170
c     cubic spline to this curve. The new curve is returned in 
1171
c     <sptime,splon,spn> with <spn> data points. The parameter
1172
c     <spn> specifies on entry the number of spline interpolated points
1173
c     along the curve.
1174
 
1175
      implicit none
1176
 
1177
c     Declaration of subroutine parameters
1178
      integer n
1179
      real time(n),lon(n)
1180
      integer spn
1181
      real sptime(spn),splon(spn)
1182
 
1183
c     Auxiliary variables
1184
      real y2ax(n)
1185
      real dt
1186
      real s
1187
      integer i
1188
      real order
1189
 
1190
c     Determine whether the input array is ascending or descending
1191
      if (time(1).gt.time(n)) then
1192
         order=-1.
1193
      else
1194
         order= 1.
1195
      endif
1196
 
1197
c     Bring the time array into ascending order
1198
      do i=1,n
1199
         time(i)=order*time(i)
1200
      enddo
1201
 
1202
c     Prepare the (natural) cubic spline interpolation
1203
      call spline (time,lon,n,1.e30,1.e30,y2ax)
1204
      dt=(time(n)-time(1))/real(spn-1)
1205
      do i=1,spn
1206
         sptime(i)=time(1)+real(i-1)*dt
1207
      enddo
1208
 
1209
c     Do the spline interpolation
1210
      do i=1,spn
1211
         call splint(time,lon,y2ax,n,sptime(i),s)
1212
         splon(i)=s
1213
      enddo
1214
 
1215
c     Change the time arrays back
1216
      do i=1,spn
1217
         sptime(i)=order*sptime(i)
1218
      enddo
1219
      do i=1,n
1220
         time(i)=order*time(i)
1221
      enddo
1222
 
1223
      return
1224
      end
1225
 
1226
c     -------------------------------------------------------------
1227
c     Basic routines for spline interpolation (Numerical Recipes)
1228
c     -------------------------------------------------------------
1229
 
1230
      SUBROUTINE spline(x,y,n,yp1,ypn,y2)
1231
      INTEGER n,NMAX
1232
      REAL yp1,ypn,x(n),y(n),y2(n)
1233
      PARAMETER (NMAX=500)
1234
      INTEGER i,k
1235
      REAL p,qn,sig,un,u(NMAX)
1236
      if (yp1.gt..99e30) then
1237
        y2(1)=0.
1238
        u(1)=0.
1239
      else
1240
        y2(1)=-0.5
1241
        u(1)=(3./(x(2)-x(1)))*((y(2)-y(1))/(x(2)-x(1))-yp1)
1242
      endif
1243
      do 11 i=2,n-1
1244
        sig=(x(i)-x(i-1))/(x(i+1)-x(i-1))
1245
        p=sig*y2(i-1)+2.
1246
        y2(i)=(sig-1.)/p
1247
        u(i)=(6.*((y(i+1)-y(i))/(x(i+
1248
     *1)-x(i))-(y(i)-y(i-1))/(x(i)-x(i-1)))/(x(i+1)-x(i-1))-sig*
1249
     *u(i-1))/p
1250
11    continue
1251
      if (ypn.gt..99e30) then
1252
        qn=0.
1253
        un=0.
1254
      else
1255
        qn=0.5
1256
        un=(3./(x(n)-x(n-1)))*(ypn-(y(n)-y(n-1))/(x(n)-x(n-1)))
1257
      endif
1258
      y2(n)=(un-qn*u(n-1))/(qn*y2(n-1)+1.)
1259
      do 12 k=n-1,1,-1
1260
        y2(k)=y2(k)*y2(k+1)+u(k)
1261
12    continue
1262
      return
1263
      END
1264
 
1265
      SUBROUTINE splint(xa,ya,y2a,n,x,y)
1266
      INTEGER n
1267
      REAL x,y,xa(n),y2a(n),ya(n)
1268
      INTEGER k,khi,klo
1269
      REAL a,b,h
1270
      klo=1
1271
      khi=n
1272
1     if (khi-klo.gt.1) then
1273
        k=(khi+klo)/2
1274
        if(xa(k).gt.x)then
1275
          khi=k
1276
        else
1277
          klo=k
1278
        endif
1279
      goto 1
1280
      endif
1281
      h=xa(khi)-xa(klo)
5 michaesp 1282
c      if (h.eq.0.) pause 'bad xa input in splint'
1283
      if (h.eq.0.) then
1284
        write(*,*) 'bad xa input in splint'
1285
        read(*,*)
1286
      end if
3 michaesp 1287
      a=(xa(khi)-x)/h
1288
      b=(x-xa(klo))/h
1289
      y=a*ya(klo)+b*ya(khi)+((a**3-a)*y2a(klo)+(b**3-b)*y2a(khi))*(h**
1290
     *2)/6.
1291
      return
1292
      END
1293
 
1294
c     ********************************************************************
1295
c     * INPUT / OUTPUT SUBROUTINES                                       *
1296
c     ********************************************************************
1297
 
1298
 
1299
c     --------------------------------------------------------------------
1300
c     Subroutines to write the CF netcdf output file
1301
c     --------------------------------------------------------------------
1302
 
1303
      subroutine writecdf2D_cf 
1304
     >          (cdfname,varname,longname,unit,gridtype,clon,clat,
1305
     >           nlonlat,dlonlat,arr,time,dx,dy,xmin,ymin,nx,ny,
1306
     >           crefile,crevar,cretime,pollon,pollat)
1307
 
1308
c     Create and write to the CF netcdf file <cdfname>. The variable
1309
c     with name <varname> and with time <time> is written. The data
1310
c     are in the two-dimensional array <arr>. The list <dx,dy,xmin,
1311
c     ymin,nx,ny> specifies the output grid. The flags <crefile> and
1312
c     <crevar> determine whether the file and/or the variable should
1313
c     be created; correspondingly for the unlimited dimension <time>
1314
c     with the flag <cretime>.
1315
 
1316
      USE netcdf
1317
 
1318
      IMPLICIT NONE
1319
 
1320
c     Declaration of input parameters
1321
      character*80 cdfname
1322
      character*80 varname,longname,unit
1323
      integer      nx,ny
1324
      real         arr(nx,ny)
1325
      real         dx,dy,xmin,ymin
1326
      real         time
1327
      integer      crefile,crevar,cretime
1328
      character*80 gridtype
1329
      real         clon,clat
1330
      integer      nlonlat
1331
      real         dlonlat
1332
      real         pollon,pollat
1333
 
1334
c     Local variables
1335
      integer      ierr
1336
      integer      ncID
1337
      integer      LonDimId,    varLonID
1338
      integer      LatDimID,    varLatID
1339
      integer      TimeDimID,   varTimeID
1340
      real         longitude(nx)
1341
      real         latitude (ny)
1342
      real         timeindex
1343
      integer      i
1344
      integer      nvars,varids(100)
1345
      integer      ndims,dimids(100)
1346
      real         timelist(1000)
1347
      integer      ntimes
1348
      integer      ind
1349
      integer      varID
1350
 
1351
c     Quick an dirty solution for fieldname conflict
1352
      if ( varname.eq.'time' ) varname = 'TIME'
1353
 
1354
c     Initially set error to indicate no errors.
1355
      ierr = 0
1356
 
1357
c     ---- Create the netCDF - skip if <crefile=0> ----------------------
1358
      if ( crefile.ne.1 ) goto 100
1359
 
1360
c     Create the file 
1361
      ierr = nf90_create(trim(cdfname), NF90_CLOBBER, ncID)
1362
 
1363
c     Define dimensions 
1364
      ierr=nf90_def_dim(ncID,'longitude',nx            , LonDimID )
1365
      ierr=nf90_def_dim(ncID,'latitude' ,ny            , LatDimID )
1366
      ierr=nf90_def_dim(ncID,'time'     ,nf90_unlimited, TimeDimID)
1367
 
1368
c     Define coordinate Variables 
1369
      ierr = nf90_def_var(ncID,'longitude',NF90_FLOAT,
1370
     >     (/ LonDimID /),varLonID)
1371
      ierr = nf90_put_att(ncID, varLonID, "standard_name","longitude")
1372
      ierr = nf90_put_att(ncID, varLonID, "units"      ,"degree_east")
1373
 
1374
      ierr = nf90_def_var(ncID,'latitude',NF90_FLOAT,
1375
     >     (/ LatDimID /),varLatID)
1376
      ierr = nf90_put_att(ncID, varLatID, "standard_name", "latitude")
1377
      ierr = nf90_put_att(ncID, varLatID, "units"    ,"degree_north")
1378
 
1379
      ierr = nf90_def_var(ncID,'time',NF90_FLOAT, 
1380
     >     (/ TimeDimID /), varTimeID)
1381
      ierr = nf90_put_att(ncID, varTimeID, "axis",            "T")
1382
      ierr = nf90_put_att(ncID, varTimeID, "calendar", "standard")
1383
      ierr = nf90_put_att(ncID, varTimeID, "long_name",    "time")
1384
      ierr = nf90_put_att(ncID, varTimeID, "units",       "hours")
1385
 
1386
c     Write global attributes 
1387
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'Conventions', 'CF-1.0')
1388
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'title',  
1389
     >     'Trajectory Densities')
1390
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'source', 
1391
     >     'Lagranto Trajectories')
1392
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'institution', 
1393
     >     'ETH Zurich, IACETH')
1394
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'grid',trim(gridtype) ) 
1395
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'clon',clon )
1396
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'clat',clat )
1397
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'nlonlat',nlonlat )
1398
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'dlonlat',dlonlat )
1399
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'nx',nx )
1400
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'ny',ny )
1401
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'dx',dx )
1402
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'dy',dy )
1403
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'xmin',xmin )
1404
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'ymin',ymin )
1405
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'pollon',pollon )
1406
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'pollat',pollat )
1407
 
1408
c     Write coordinate data 
1409
      do i = 1,nx+1
1410
         longitude(i) = xmin + real(i-1) * dx
1411
      enddo
1412
      do i = 1,ny+1
1413
         latitude(i)  = ymin + real(i-1) * dy
1414
      enddo
1415
 
1416
c     Check whether the definition was successful
1417
      ierr = nf90_enddef(ncID)
1418
      if (ierr.gt.0) then
1419
         print*, 'An error occurred while attempting to ', 
1420
     >        'finish definition mode.'
1421
         stop
1422
      endif
1423
 
1424
c     Write coordinate data  
1425
      ierr = nf90_put_var(ncID,varLonID ,longitude)
1426
      ierr = nf90_put_var(ncID,varLatID ,latitude )
1427
 
1428
c     Close netCDF file 
1429
      ierr = nf90_close(ncID)
1430
 
1431
 100  continue
1432
 
1433
c     ---- Define a new variable - skip if <crevar=0> -----------------------
1434
 
1435
      if ( crevar.ne.1 ) goto 110
1436
 
1437
c     Open the file for read(write access
1438
      ierr = nf90_open  (trim(cdfname), NF90_WRITE  , ncID)
1439
 
1440
c     Get the IDs for dimensions
1441
      ierr = nf90_inq_dimid(ncID,'longitude', LonDimID )
1442
      ierr = nf90_inq_dimid(ncID,'latitude' , LatDimID )
1443
      ierr = nf90_inq_dimid(ncID,'time'     , TimeDimID)
1444
 
1445
c     Enter define mode
1446
      ierr = nf90_redef(ncID)
1447
 
1448
c     Write definition and add attributes
1449
      ierr = nf90_def_var(ncID,varname,NF90_FLOAT,
1450
     >                   (/ LonDimID, LatDimID, TimeDimID /),varID)
1451
      ierr = nf90_put_att(ncID, varID, "long_name" , longname )
1452
      ierr = nf90_put_att(ncID, varID, "units"     , unit     ) 
1453
      ierr = nf90_put_att(ncID, varID, '_FillValue', -999.99  ) 
1454
 
1455
c     Check whether definition was successful
1456
      ierr = nf90_enddef(ncID)
1457
      if (ierr.gt.0) then
1458
         print*, 'An error occurred while attempting to ', 
1459
     >           'finish definition mode.'
1460
         stop
1461
      endif
1462
 
1463
c     Close netCDF file 
1464
      ierr = nf90_close(ncID)
1465
 
1466
 110  continue
1467
 
1468
c     ---- Create a new time (unlimited dimension) - skip if <cretime=0> ------
1469
 
1470
      if ( cretime.ne.1 ) goto 120
1471
 
1472
c     Open the file for read/write access
1473
      ierr = nf90_open  (trim(cdfname), NF90_WRITE, ncID)
1474
 
1475
c     Get the list of times on the netCDF file
1476
      ierr = nf90_inq_dimid(ncID,'time', TimeDimID)
1477
      if ( ierr.ne.0 ) then
1478
         print*,'Time dimension is not defined on ',trim(cdfname),
1479
     >          ' .... Stop'
1480
         stop
1481
      endif
1482
      ierr = nf90_inquire_dimension(ncID, TimeDimID, len = ntimes)
1483
      ierr = nf90_inq_varid(ncID,'time', varTimeID)
1484
      if ( ierr.ne.0 ) then
1485
         print*,'Variable time is not defined on ',trim(cdfname),
1486
     >          ' ... Stop'
1487
         stop
1488
      endif
1489
      ierr = nf90_get_var(ncID,varTimeID,timelist(1:ntimes))
1490
 
1491
c     Decide whether a new time must be written
1492
      ind = 0
1493
      do i=1,ntimes
1494
         if ( time.eq.timelist(i) ) ind = i
1495
      enddo
1496
 
1497
c     Extend the time list if required 
1498
      if ( ind.eq.0 ) then
1499
         ntimes           = ntimes + 1
1500
         timelist(ntimes) = time
1501
         ierr = nf90_put_var(ncID,varTimeID,timelist(1:ntimes))
1502
      endif
1503
 
1504
c     Close netCDF file 
1505
      ierr = nf90_close(ncID)
1506
 
1507
 120  continue
1508
 
1509
c     ---- Write data --------------------------------------------------
1510
 
1511
c     Open the file for read/write access
1512
      ierr = nf90_open  (trim(cdfname), NF90_WRITE , ncID)
1513
 
1514
c     Get the varID
1515
      ierr = nf90_inq_varid(ncID,varname, varID )
1516
      if (ierr.ne.0) then
1517
         print*,'Variable ',trim(varname),' is not defined on ',
1518
     >          trim(cdfname)
1519
         stop
1520
      endif
1521
 
1522
c     Get the time index
1523
      ierr = nf90_inq_dimid(ncID,'time', TimeDimID)
1524
      if ( ierr.ne.0 ) then
1525
         print*,'Time dimension is not defined on ',trim(cdfname),
1526
     >          ' .... Stop'
1527
         stop
1528
      endif
1529
      ierr = nf90_inquire_dimension(ncID, TimeDimID, len = ntimes)
1530
      ierr = nf90_inq_varid(ncID,'time', varTimeID)
1531
      if ( ierr.ne.0 ) then
1532
         print*,'Variable time is not defined on ',trim(cdfname),
1533
     >          ' ... Stop'
1534
         stop
1535
      endif
1536
      ierr = nf90_get_var(ncID,varTimeID,timelist(1:ntimes))
1537
      ind = 0
1538
      do i=1,ntimes
1539
         if ( time.eq.timelist(i) ) ind = i
1540
      enddo
1541
      if (ind.eq.0) then
1542
         print*,'Time',time,' is not defined on the netCDF file',
1543
     >          trim(cdfname),' ... Stop'
1544
         stop
1545
      endif
1546
 
1547
c     Write data block      
1548
      ierr = nf90_put_var(ncID,varID,arr,
1549
     >                    start = (/ 1, 1, ind /), 
1550
     >                    count = (/ nx, ny, 1 /) )
1551
 
1552
c     Check whether writing was successful 
1553
      ierr = nf90_close(ncID)
1554
      if (ierr.ne.0) then
1555
         write(*,*) trim(nf90_strerror(ierr))
1556
         write(*,*) 'An error occurred while attempting to ', 
1557
     >              'close the netcdf file.'
1558
         write(*,*) 'in clscdf_CF'
1559
      endif
1560
 
1561
      end
1562
 
1563
 
1564
c     ********************************************************************************
1565
c     * Transformation routine: LMSTOLM and PHSTOPH from library gm2em               *
1566
c     ********************************************************************************
1567
 
1568
      REAL FUNCTION LMSTOLM (PHIS, LAMS, POLPHI, POLLAM)
1569
C
1570
C**** LMSTOLM  -   FC:BERECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE FUER
1571
C****                 EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
1572
C****                 IM ROTIERTEN SYSTEM. DER NORDPOL DES SYSTEMS HAT
1573
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1574
C**   AUFRUF   :   LAM = LMSTOLM (PHIS, LAMS, POLPHI, POLLAM)
1575
C**   ENTRIES  :   KEINE
1576
C**   ZWECK    :   BERECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE FUER
1577
C**                EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
1578
C**                IM ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
1579
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1580
C**   VERSIONS-
1581
C**   DATUM    :   03.05.90
1582
C**
1583
C**   EXTERNALS:   KEINE
1584
C**   EINGABE-
1585
C**   PARAMETER:   PHIS     REAL   GEOGR. BREITE DES PUNKTES IM ROT.SYS.
1586
C**                LAMS     REAL   GEOGR. LAENGE DES PUNKTES IM ROT.SYS.
1587
C**                POLPHI   REAL   WAHRE GEOGR. BREITE DES NORDPOLS
1588
C**                POLLAM   REAL   WAHRE GEOGR. LAENGE DES NORDPOLS
1589
C**   AUSGABE-
1590
C**   PARAMETER:   WAHRE GEOGRAPHISCHE LAENGE ALS WERT DER FUNKTION
1591
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
1592
C**
1593
C**   COMMON-
1594
C**   BLOECKE  :   KEINE
1595
C**
1596
C**   FEHLERBE-
1597
C**   HANDLUNG :   KEINE
1598
C**   VERFASSER:   D.MAJEWSKI
1599
 
1600
      REAL        LAMS,PHIS,POLPHI,POLLAM
1601
 
1602
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
1603
 
1604
      ZSINPOL = SIN(ZPIR18*POLPHI)
1605
      ZCOSPOL = COS(ZPIR18*POLPHI)
1606
      ZLAMPOL = ZPIR18*POLLAM
1607
      ZPHIS   = ZPIR18*PHIS
1608
      ZLAMS   = LAMS
1609
      IF(ZLAMS.GT.180.0) ZLAMS = ZLAMS - 360.0
1610
      ZLAMS   = ZPIR18*ZLAMS
1611
 
1612
      ZARG1   = SIN(ZLAMPOL)*(- ZSINPOL*COS(ZLAMS)*COS(ZPHIS)  +
1613
     1                          ZCOSPOL*           SIN(ZPHIS)) -
1614
     2          COS(ZLAMPOL)*           SIN(ZLAMS)*COS(ZPHIS)
1615
      ZARG2   = COS(ZLAMPOL)*(- ZSINPOL*COS(ZLAMS)*COS(ZPHIS)  +
1616
     1                          ZCOSPOL*           SIN(ZPHIS)) +
1617
     2          SIN(ZLAMPOL)*           SIN(ZLAMS)*COS(ZPHIS)
1618
      IF (ABS(ZARG2).LT.1.E-30) THEN
1619
        IF (ABS(ZARG1).LT.1.E-30) THEN
1620
          LMSTOLM =   0.0
1621
        ELSEIF (ZARG1.GT.0.) THEN
1622
              LMSTOLAM =  90.0
1623
            ELSE
1624
              LMSTOLAM = -90.0
1625
            ENDIF
1626
      ELSE
1627
        LMSTOLM = ZRPI18*ATAN2(ZARG1,ZARG2)
1628
      ENDIF
1629
 
1630
      RETURN
1631
      END
1632
 
1633
 
1634
      REAL FUNCTION PHSTOPH (PHIS, LAMS, POLPHI, POLLAM)
1635
C
1636
C**** PHSTOPH  -   FC:BERECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE FUER
1637
C****                 EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
1638
C****                 ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
1639
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1640
C**   AUFRUF   :   PHI = PHSTOPH (PHIS, LAMS, POLPHI, POLLAM)
1641
C**   ENTRIES  :   KEINE
1642
C**   ZWECK    :   BERECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE FUER
1643
C**                EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
1644
C**                ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
1645
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1646
C**   VERSIONS-
1647
C**   DATUM    :   03.05.90
1648
C**
1649
C**   EXTERNALS:   KEINE
1650
C**   EINGABE-
1651
C**   PARAMETER:   PHIS     REAL   GEOGR. BREITE DES PUNKTES IM ROT.SYS.
1652
C**                LAMS     REAL   GEOGR. LAENGE DES PUNKTES IM ROT.SYS.
1653
C**                POLPHI   REAL   WAHRE GEOGR. BREITE DES NORDPOLS
1654
C**                POLLAM   REAL   WAHRE GEOGR. LAENGE DES NORDPOLS
1655
C**   AUSGABE-
1656
C**   PARAMETER:   WAHRE GEOGRAPHISCHE BREITE ALS WERT DER FUNKTION
1657
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
1658
C**
1659
C**   COMMON-
1660
C**   BLOECKE  :   KEINE
1661
C**
1662
C**   FEHLERBE-
1663
C**   HANDLUNG :   KEINE
1664
C**   VERFASSER:   D.MAJEWSKI
1665
 
1666
      REAL        LAMS,PHIS,POLPHI,POLLAM
1667
 
1668
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
1669
 
1670
      SINPOL = SIN(ZPIR18*POLPHI)
1671
      COSPOL = COS(ZPIR18*POLPHI)
1672
      ZPHIS  = ZPIR18*PHIS
1673
      ZLAMS  = LAMS
1674
      IF(ZLAMS.GT.180.0) ZLAMS = ZLAMS - 360.0
1675
      ZLAMS  = ZPIR18*ZLAMS
1676
      ARG     = COSPOL*COS(ZPHIS)*COS(ZLAMS) + SINPOL*SIN(ZPHIS)
1677
 
1678
      PHSTOPH = ZRPI18*ASIN(ARG)
1679
 
1680
      RETURN
1681
      END
1682
 
1683
 
1684
      REAL FUNCTION LMTOLMS (PHI, LAM, POLPHI, POLLAM)
1685
C
1686
C%Z% Modul %M%, V%I% vom %G%, extrahiert am %H%
1687
C
1688
C**** LMTOLMS  -   FC:UMRECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE LAM
1689
C****                 AUF EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
1690
C****                 IM ROTIERTEN SYSTEM. DER NORDPOL DES SYSTEMS HAT
1691
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1692
C**   AUFRUF   :   LAM = LMTOLMS (PHI, LAM, POLPHI, POLLAM)
1693
C**   ENTRIES  :   KEINE
1694
C**   ZWECK    :   UMRECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE LAM AUF
1695
C**                EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
1696
C**                ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
1697
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1698
C**   VERSIONS-
1699
C**   DATUM    :   03.05.90
1700
C**
1701
C**   EXTERNALS:   KEINE
1702
C**   EINGABE-
1703
C**   PARAMETER:   PHI    REAL BREITE DES PUNKTES IM GEOGR. SYSTEM
1704
C**                LAM    REAL LAENGE DES PUNKTES IM GEOGR. SYSTEM
1705
C**                POLPHI REAL GEOGR.BREITE DES N-POLS DES ROT. SYSTEMS
1706
C**                POLLAM REAL GEOGR.LAENGE DES N-POLS DES ROT. SYSTEMS
1707
C**   AUSGABE-
1708
C**   PARAMETER:   WAHRE GEOGRAPHISCHE LAENGE ALS WERT DER FUNKTION
1709
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
1710
C**
1711
C**   COMMON-
1712
C**   BLOECKE  :   KEINE
1713
C**
1714
C**   FEHLERBE-
1715
C**   HANDLUNG :   KEINE
1716
C**   VERFASSER:   G. DE MORSIER
1717
 
1718
      REAL        LAM,PHI,POLPHI,POLLAM
1719
 
1720
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
1721
 
1722
      ZSINPOL = SIN(ZPIR18*POLPHI)
1723
      ZCOSPOL = COS(ZPIR18*POLPHI)
1724
      ZLAMPOL =     ZPIR18*POLLAM
1725
      ZPHI    =     ZPIR18*PHI
1726
      ZLAM    = LAM
1727
      IF(ZLAM.GT.180.0) ZLAM = ZLAM - 360.0
1728
      ZLAM    = ZPIR18*ZLAM
1729
 
1730
      ZARG1   = - SIN(ZLAM-ZLAMPOL)*COS(ZPHI)
1731
      ZARG2   = - ZSINPOL*COS(ZPHI)*COS(ZLAM-ZLAMPOL)+ZCOSPOL*SIN(ZPHI)
1732
      IF (ABS(ZARG2).LT.1.E-30) THEN
1733
        IF (ABS(ZARG1).LT.1.E-30) THEN
1734
          LMTOLMS =   0.0
1735
        ELSEIF (ZARG1.GT.0.) THEN
1736
              LMTOLMS =  90.0
1737
            ELSE
1738
              LMTOLMS = -90.0
1739
            ENDIF
1740
      ELSE
1741
        LMTOLMS = ZRPI18*ATAN2(ZARG1,ZARG2)
1742
      ENDIF
1743
 
1744
      RETURN
1745
      END
1746
 
1747
 
1748
      REAL FUNCTION PHTOPHS (PHI, LAM, POLPHI, POLLAM)
1749
C
1750
C%Z% Modul %M%, V%I% vom %G%, extrahiert am %H%
1751
C
1752
C**** PHTOPHS  -   FC:UMRECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE PHI
1753
C****                 AUF EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
1754
C****                 IM ROTIERTEN SYSTEM. DER NORDPOL DES SYSTEMS HAT
1755
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1756
C**   AUFRUF   :   PHI = PHTOPHS (PHI, LAM, POLPHI, POLLAM)
1757
C**   ENTRIES  :   KEINE
1758
C**   ZWECK    :   UMRECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE PHI AUF
1759
C**                EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
1760
C**                ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
1761
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1762
C**   VERSIONS-
1763
C**   DATUM    :   03.05.90
1764
C**
1765
C**   EXTERNALS:   KEINE
1766
C**   EINGABE-
1767
C**   PARAMETER:   PHI    REAL BREITE DES PUNKTES IM GEOGR. SYSTEM
1768
C**                LAM    REAL LAENGE DES PUNKTES IM GEOGR. SYSTEM
1769
C**                POLPHI REAL GEOGR.BREITE DES N-POLS DES ROT. SYSTEMS
1770
C**                POLLAM REAL GEOGR.LAENGE DES N-POLS DES ROT. SYSTEMS
1771
C**   AUSGABE-
1772
C**   PARAMETER:   ROTIERTE BREITE PHIS ALS WERT DER FUNKTION
1773
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
1774
C**
1775
C**   COMMON-
1776
C**   BLOECKE  :   KEINE
1777
C**
1778
C**   FEHLERBE-
1779
C**   HANDLUNG :   KEINE
1780
C**   VERFASSER:   G. DE MORSIER
1781
 
1782
      REAL        LAM,PHI,POLPHI,POLLAM
1783
 
1784
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
1785
 
1786
      ZSINPOL = SIN(ZPIR18*POLPHI)
1787
      ZCOSPOL = COS(ZPIR18*POLPHI)
1788
      ZLAMPOL = ZPIR18*POLLAM
1789
      ZPHI    = ZPIR18*PHI
1790
      ZLAM    = LAM
1791
      IF(ZLAM.GT.180.0) ZLAM = ZLAM - 360.0
1792
      ZLAM    = ZPIR18*ZLAM
1793
      ZARG    = ZCOSPOL*COS(ZPHI)*COS(ZLAM-ZLAMPOL) + ZSINPOL*SIN(ZPHI)
1794
 
1795
      PHTOPHS = ZRPI18*ASIN(ZARG)
1796
 
1797
      RETURN
1798
      END