Subversion Repositories lagranto.ecmwf

Rev

Rev 29 | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
3 michaesp 1
      PROGRAM trace
2
 
3
C     ********************************************************************
4
C     *                                                                  *
5
C     * Pseudo-lidar plots along trajectories                             *
6
C     *                                                                  *
7
C     * Heini Wernli       first version:       April 1993               *
8
C     * Michael Sprenger   major upgrade:       2008-2009                *
9
C     *                                                                  *
10
C     ********************************************************************
11
 
12
      implicit none
13
 
14
c     --------------------------------------------------------------------
15
c     Declaration of parameters
16
c     --------------------------------------------------------------------
17
 
18
c     Maximum number of levels for input files
19
      integer   nlevmax
20
      parameter (nlevmax=100)
21
 
22
c     Maximum number of input files (dates, length of trajectories)
23
      integer   ndatmax
46 michaesp 24
      parameter (ndatmax=10000)
3 michaesp 25
 
26
c     Numerical epsilon (for float comparison)
27
      real      eps
28
      parameter (eps=0.001)
29
 
30
c     Conversion factors
31
      real      pi180                                   ! deg -> rad
32
      parameter (pi180=3.14159/180.)
33
      real      deg2km                                  ! deg -> km (at equator)
34
      parameter (deg2km=111.)
35
 
36
c     Prefix for primary and secondary fields
37
      character charp
38
      character chars
39
      parameter (charp='P')
40
      parameter (chars='S')
41
 
42
c     --------------------------------------------------------------------
43
c     Declaration of variables
44
c     --------------------------------------------------------------------
45
 
46
c     Input and output format for trajectories (see iotra.f)
47
      integer   inpmode
48
 
49
c     Input parameters
50
      character*80                           inpfile         ! Input trajectory file
51
      character*80                           outfile         ! Output netCDF file
52
      character*80                           outmode         ! Output mode (sum,mean)
53
      integer                                ntra            ! Number of trajectories
54
      integer                                ncol            ! Number of columns (including time, lon, lat, p)
55
      integer                                ntim            ! Number of times per trajectory
56
      integer                                ntrace0         ! Number of trace variables
57
      character*80                           tvar(200)       ! Tracing variable name (only the variable)
58
      character*1                            tfil(200)       ! Filename prefix 
59
      real                                   fac(200)        ! Scaling factor 
60
      integer                                compfl(200)     ! Computation flag (1=compute)
61
      integer                                numdat          ! Number of input files
62
      character*11                           dat(ndatmax)    ! Dates of input files
63
      real                                   timeinc         ! Time increment between input files
64
      real                                   tst             ! Time shift of start relative to first data file
65
      real                                   ten             ! Time shift of end relatiev to first data file  
66
      character*20                           startdate       ! First time/date on trajectory
67
      character*20                           enddate         ! Last time/date on trajectory
68
      character*80                           timecheck       ! Either 'yes' or 'no'
69
      character*80                           intmode         ! Interpolation mode ('normal', 'nearest')
70
	  real                                   pmin,pmax       ! Pressure range for output grid
71
	  integer                                npre            ! Number of pressure levels in output grid
72
	  character*80                           centering       ! Centering around trajectory position ('yes','no')
5 michaesp 73
	  character*80                       direction       ! Direction of lidar (vertical,lat,lon,normal)
74
      character*80                           dumpcoord       ! Dumping coordinates ('yes','no')
3 michaesp 75
 
76
c     Trajectories
5 michaesp 77
      real,allocatable, dimension (:,:,:) :: trainp                ! Input trajectories (ntra,ntim,ncol)
78
      integer                                reftime(6)            ! Reference date
79
      character*80                           varsinp(100)          ! Field names for input trajectory
80
      integer                                fid,fod               ! File identifier for inp and out trajectories
81
      real                                   x0_tra,y0_tra,p0_tra  ! Position of air parcel (physical space)
82
      real                                   reltpos0              ! Relative time of air parcel
83
      real                                   xind,yind,pind        ! Position of air parcel (grid space)
84
      integer                                fbflag                ! Flag for forward (1) or backward (-1) trajectories
3 michaesp 85
 
86
c     Meteorological fields from input file
87
      real,allocatable, dimension (:)     :: spt0,spt1       ! Surface pressure
88
      real,allocatable, dimension (:)     :: p3t0,p3t1       ! 3d-pressure 
89
      real,allocatable, dimension (:)     :: f3t0,f3t1       ! 3d field for tracing 
90
      character*80                           svars(100)      ! List of variables on S file
91
      character*80                           pvars(100)      ! List of variables on P file
92
      integer                                n_svars         ! Number of variables on S file
93
      integer                                n_pvars         ! Number of variables on P file
94
 
95
c     Input grid description
96
      real                                   pollon,pollat   ! Longitude/latitude of pole
97
      real                                   ak(100)         ! Vertical layers and levels
98
      real                                   bk(100) 
99
      real                                   xmin,xmax       ! Zonal grid extension
100
      real                                   ymin,ymax       ! Meridional grid extension
101
      integer                                nx,ny,nz        ! Grid dimensions
102
      real                                   dx,dy           ! Horizontal grid resolution
103
      integer                                hem             ! Flag for hemispheric domain
104
      integer                                per             ! Flag for periodic domain
105
      real                                   stagz           ! Vertical staggering
106
      real                                   mdv             ! Missing data value
107
 
108
c	  Output grid  and fields
109
      real                                   levels(1000)    ! Ouput levels
46 michaesp 110
      real                                   times (10000)    ! Output times
3 michaesp 111
	  real,allocatable, dimension (:,:)   :: out_pos         ! Position of trajectories
112
	  real,allocatable, dimension (:,:)   :: out_val         ! Output lidar field
113
	  real,allocatable, dimension (:,:)   :: out_cnt         ! # output lidar sum ups
114
 
115
 
116
c     Auxiliary variables
13 michaesp 117
      integer                                i,j,k,l,n,m
3 michaesp 118
      real                                   rd
119
      character*80                           filename
120
      real                                   time0,time1,reltpos
121
      integer                                itime0,itime1
122
      integer                                stat
123
      real                                   tstart
124
      integer                                iloaded0,iloaded1
125
      real                                   f0
126
      real                                   frac
127
      real                                   tload,tfrac
128
      integer                                isok
129
      character                              ch
130
      integer                                ind
131
      integer                                ind1,ind2,ind3,ind4,ind5
132
      integer                                ind6,ind7,ind8,ind9,ind0
133
      integer                                noutside
134
      real                                   delta
135
      integer                                itrace0
136
      character*80                           string
137
      character*80                           cdfname
138
      character*80                           varname
139
      real                                   time
140
      character*80                           longname
141
      character*80                           unit
142
      integer                                ind_time
143
      integer                                ind_pre
5 michaesp 144
      real                                   rlat,rlon
145
      real                                   x0,y0,p0
146
      real                                   vx0,vy0,vx1,vy1
147
      real                                   rotation,lon,lat
3 michaesp 148
 
149
c     Externals 
150
      real                                   int_index4
151
      external                               int_index4
152
 
153
c     --------------------------------------------------------------------
154
c     Start of program, Read parameters, get grid parameters
155
c     --------------------------------------------------------------------
156
 
157
c     Write start message
158
      print*,'========================================================='
159
      print*,'              *** START OF PROGRAM LIDAR ***'
160
      print*
161
 
162
c     Read parameters
163
      open(10,file='trace.param')
164
       read(10,*) inpfile
165
       read(10,*) outfile
166
       read(10,*) outmode
167
       read(10,*) startdate
168
       read(10,*) enddate 
169
       read(10,*) fbflag
170
       read(10,*) numdat
171
       if ( fbflag.eq.1) then
172
          do i=1,numdat
173
             read(10,'(a11)') dat(i)
174
          enddo
175
       else
176
          do i=numdat,1,-1
177
             read(10,'(a11)') dat(i)
178
          enddo
179
       endif
180
       read(10,*) timeinc
181
       read(10,*) tst
182
       read(10,*) ten
183
       read(10,*) ntra
184
       read(10,*) ntim
185
       read(10,*) ncol
186
       read(10,*) ntrace0
187
       do i=1,ntrace0
188
          read(10,*) tvar(i), fac(i), compfl(i), tfil(i)
189
       enddo
190
       read(10,*) n_pvars
191
       do i=1,n_pvars
192
          read(10,*) pvars(i)
193
       enddo
194
       read(10,*) n_svars
195
       do i=1,n_svars
196
          read(10,*) svars(i)
197
       enddo
198
       read(10,*) timecheck
199
       read(10,*) intmode
200
       read(10,*) pmin,pmax,npre
201
       read(10,*) centering
5 michaesp 202
       read(10,*) direction
203
       read(10,*) dumpcoord
3 michaesp 204
      close(10)
205
 
5 michaesp 206
c     Check that the direction is ok
207
      if ( ( direction.ne.'vertical' ).and.
208
     >     ( direction.ne.'lat'      ).and.
209
     >     ( direction.ne.'lon'      ).and.
210
     >     ( direction.ne.'normal'   ) ) 
211
     >then 
212
         print*,' ERROR: invalid direction ',trim(direction)
213
         stop
214
      endif
215
 
3 michaesp 216
c     Remove commented tracing fields
217
      itrace0 = 1
218
      do while ( itrace0.le.ntrace0) 
219
         string = tvar(itrace0)
220
         if ( string(1:1).eq.'#' ) then
221
            do i=itrace0,ntrace0-1
222
               tvar(i)   = tvar(i+1)
223
               fac(i)    = fac(i+1)
224
               compfl(i) = compfl(i+1)
225
               tfil(i)   = tfil(i+1)
226
            enddo
227
            ntrace0 = ntrace0 - 1
228
         else
229
            itrace0 = itrace0 + 1
230
         endif
231
      enddo
232
 
233
c     Set the formats of the input  files
234
      call mode_tra(inpmode,inpfile)
235
      if (inpmode.eq.-1) inpmode=1
236
 
237
C     Convert time shifts <tst,ten> from <hh.mm> into fractional time
238
      call hhmm2frac(tst,frac)
239
      tst = frac
240
      call hhmm2frac(ten,frac)
241
      ten = frac
242
 
243
c     Set the time for the first data file (depending on forward/backward mode)
244
      if (fbflag.eq.1) then
245
        tstart = -tst
246
      else
247
        tstart = tst
248
      endif
249
 
250
c     Read the constant grid parameters (nx,ny,nz,xmin,xmax,ymin,ymax,pollon,pollat)
251
c     The negative <-fid> of the file identifier is used as a flag for parameter retrieval  
252
      filename = charp//dat(1)
9 michaesp 253
      varname  = tvar(1)
3 michaesp 254
      call input_open (fid,filename)
255
      call input_grid (-fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
256
     >                 tstart,pollon,pollat,rd,rd,nz,rd,rd,rd,timecheck)
257
      call input_close(fid)
258
 
259
C     Allocate memory for some meteorological arrays
260
      allocate(spt0(nx*ny),stat=stat)
261
      if (stat.ne.0) print*,'*** error allocating array spt0 ***'   ! Surface pressure
262
      allocate(spt1(nx*ny),stat=stat)
263
      if (stat.ne.0) print*,'*** error allocating array spt1 ***'
264
      allocate(p3t0(nx*ny*nz),stat=stat)
265
      if (stat.ne.0) print*,'*** error allocating array p3t0 ***'   ! Pressure
266
      allocate(p3t1(nx*ny*nz),stat=stat)
267
      if (stat.ne.0) print*,'*** error allocating array p3t1 ***'
268
      allocate(f3t0(nx*ny*nz),stat=stat)
269
      if (stat.ne.0) print*,'*** error allocating array p3t0 ***'   ! Lidar field
270
      allocate(f3t1(nx*ny*nz),stat=stat)
271
      if (stat.ne.0) print*,'*** error allocating array p3t1 ***'
272
 
273
c	  Allocate memory for output field
274
	  allocate(out_pos(ntim,npre),stat=stat)
275
      if (stat.ne.0) print*,'*** error allocating array out_pos ***'
276
	  allocate(out_val(ntim,npre),stat=stat)
277
      if (stat.ne.0) print*,'*** error allocating array out_val ***'
278
	  allocate(out_cnt(ntim,npre),stat=stat)
279
      if (stat.ne.0) print*,'*** error allocating array out_cnt ***'
280
 
281
C     Get memory for trajectory arrays
282
      allocate(trainp(ntra,ntim,ncol),stat=stat)
283
      if (stat.ne.0) print*,'*** error allocating array tra      ***' 
284
 
285
c     Set the flags for periodic domains
286
      if ( abs(xmax-xmin-360.).lt.eps ) then
287
         per = 1
288
      elseif ( abs(xmax-xmin-360.+dx).lt.eps ) then
289
         per = 2
290
      else
291
         per = 0
292
      endif
293
 
294
C     Set logical flag for periodic data set (hemispheric or not)
295
      hem = 0
296
      if (per.eq.0.) then
297
         delta=xmax-xmin-360.
298
         if (abs(delta+dx).lt.eps) then               ! Program aborts: arrays must be closed
299
            print*,' ERROR: arrays must be closed... Stop'
300
         else if (abs(delta).lt.eps) then ! Periodic and hemispheric
301
           hem=1
302
           per=360.
303
        endif
304
      else                                            ! Periodic and hemispheric
305
         hem=1
306
      endif
307
 
308
c     Write some status information
309
      print*,'---- INPUT PARAMETERS -----------------------------------'
310
      print*
311
      print*,'  Input trajectory file  : ',trim(inpfile)
312
      print*,'  Format of input file   : ',inpmode
313
      print*,'  Output netCDF    file  : ',trim(outfile)
314
      print*,'  Format of output file  : ',trim(outmode)
315
      print*,'  Forward/backward       : ',fbflag
316
      print*,'  #tra                   : ',ntra
317
      print*,'  #col                   : ',ncol
318
      print*,'  #tim                   : ',ntim
319
      print*,'  No time check          : ',trim(timecheck)
320
      print*,'  Interpolation mode     : ',trim(intmode)
321
      do i=1,ntrace0
322
         if (compfl(i).eq.0) then
323
            print*,'  Tracing field          : ',
324
     >                 trim(tvar(i)), fac(i), ' 0 ', tfil(i)
325
         else
326
            print*,'  Tracing field          : ',
327
     >                trim(tvar(i)),' : online calc not supported'
328
         endif
329
      enddo
330
      print*,'  Output (pmin,pmax,n)   : ',pmin,pmax,npre
331
      print*,'  Centering              : ',trim(centering)
5 michaesp 332
      print*,'  Orientation            : ',trim(direction)
333
      print*,'  Coordinate Dump        : ',trim(dumpcoord)
3 michaesp 334
      print*
335
      print*,'---- INPUT DATA FILES -----------------------------------'
336
      print*
337
      call frac2hhmm(tstart,tload)
338
      print*,'  Time of 1st data file  : ',tload
339
      print*,'  #input files           : ',numdat
340
      print*,'  time increment         : ',timeinc
341
      call frac2hhmm(tst,tload)
342
      print*,'  Shift of start         : ',tload
343
      call frac2hhmm(ten,tload)
344
      print*,'  Shift of end           : ',tload
345
      print*,'  First/last input file  : ',trim(dat(1)),
346
     >                                     ' ... ',
5 michaesp 347
     >                                     trim(dat(numdat)) 
3 michaesp 348
      print*,'  Primary variables      : ',trim(pvars(1))
349
      do i=2,n_pvars
350
         print*,'                         : ',trim(pvars(i))
351
      enddo
352
      if ( n_svars.ge.1 ) then
353
         print*,'  Secondary variables    : ',trim(svars(1))
354
         do i=2,n_svars
355
            print*,'                         : ',trim(svars(i))
356
         enddo
357
      endif
358
      print*
359
      print*,'---- CONSTANT GRID PARAMETERS ---------------------------'
360
      print*
361
      print*,'  xmin,xmax     : ',xmin,xmax
362
      print*,'  ymin,ymax     : ',ymin,ymax
363
      print*,'  dx,dy         : ',dx,dy
364
      print*,'  pollon,pollat : ',pollon,pollat
365
      print*,'  nx,ny,nz      : ',nx,ny,nz
366
      print*,'  per, hem      : ',per,hem
367
      print*
368
 
369
c     --------------------------------------------------------------------
370
c     Load the input trajectories
371
c     --------------------------------------------------------------------
372
 
373
c     Read the input trajectory file
374
      call ropen_tra(fid,inpfile,ntra,ntim,ncol,reftime,varsinp,inpmode)
375
      call read_tra (fid,trainp,ntra,ntim,ncol,inpmode)
376
      call close_tra(fid,inpmode)
377
 
378
c     Check that first four columns correspond to time,lon,lat,p
379
      if ( (varsinp(1).ne.'time' ).or.
380
     >     (varsinp(2).ne.'xpos' ).and.(varsinp(2).ne.'lon' ).or.
381
     >     (varsinp(3).ne.'ypos' ).and.(varsinp(3).ne.'lat' ).or.
382
     >     (varsinp(4).ne.'ppos' ).and.(varsinp(4).ne.'p'   ) )
383
     >then
384
         print*,' ERROR: problem with input trajectories ...'
385
         stop
386
      endif
387
      varsinp(1) = 'time'
388
      varsinp(2) = 'lon'
389
      varsinp(3) = 'lat'
390
      varsinp(4) = 'p'
391
 
392
c     Write some status information of the input trajectories
393
      print*,'---- INPUT TRAJECTORIES ---------------------------------'
394
      print*
395
      print*,' Start date             : ',trim(startdate)
396
      print*,' End date               : ',trim(enddate)
397
      print*,' Reference time (year)  : ',reftime(1)
398
      print*,'                (month) : ',reftime(2)
399
      print*,'                (day)   : ',reftime(3)
400
      print*,'                (hour)  : ',reftime(4)
401
      print*,'                (min)   : ',reftime(5)
402
      print*,' Time range (min)       : ',reftime(6)
403
      do i=1,ncol
404
         print*,' Var                    :',i,trim(varsinp(i))
405
      enddo
406
      print*
407
 
408
c     Check that first time is 0 - otherwise the tracing will produce
409
c     wrong results because later in the code only absolute times are
410
c     considered: <itime0   = int(abs(tfrac-tstart)/timeinc) + 1>. This 
411
c     will be changed in a future version.
412
      if ( abs( trainp(1,1,1) ).gt.eps ) then
413
         print*,' ERROR: First time of trajectory must be 0, i.e. '
414
         print*,'     correspond to the reference date. Otherwise'
415
         print*,'     the tracing will give wrong results... STOP'
416
         stop
417
      endif
5 michaesp 418
 
419
c     If requested, open the coordinate dump file
420
      if ( dumpcoord.eq.'yes' ) then
421
         open(10,file=trim(outfile)//'.coord')
422
      endif
423
 
3 michaesp 424
c     --------------------------------------------------------------------
425
c     Trace the fields (fields available on input files)
426
c     --------------------------------------------------------------------
427
 
428
      print*
429
      print*,'---- LIDAR FROM PRIMARY AND SECONDARY DATA FILES ------'
430
 
431
c     Loop over all tracing fields
432
      do i=1,ntrace0
433
 
434
c	      Skip all fields marked for online calculation
435
          if ( compfl(i).eq.1 ) goto 110
436
 
437
c	      Init the output fields: position and lidar field
438
	      do k=1,ntim
439
	      	do l=1,npre
440
	      		out_pos(k,l) = 0.
441
	      	    out_val(k,l) = 0.
442
	      	    out_cnt(k,l) = 0.
443
	      	 enddo
444
	      enddo
445
 
446
c         Write some status information
447
          print*
448
          print*,' Now lidaring           : ',
449
     >         trim(tvar(i)),compfl(i),' ',trim(tfil(i))
450
 
451
c         Reset flags for load manager
452
          iloaded0 = -1
453
          iloaded1 = -1
454
 
455
c         Reset the counter for fields outside domain
456
          noutside = 0
457
 
458
c         Loop over all times
459
          do j=1,ntim
460
 
461
c            Convert trajectory time from hh.mm to fractional time
462
             call hhmm2frac(trainp(1,j,1),tfrac)
463
 
464
c            Get the times which are needed
465
             itime0   = int(abs(tfrac-tstart)/timeinc) + 1
466
             time0    = tstart + fbflag * real(itime0-1) * timeinc
467
             itime1   = itime0 + 1
468
             time1    = time0 + fbflag * timeinc
469
             if ( itime1.gt.numdat ) then
470
                itime1 = itime0
471
                time1  = time0
472
             endif
473
 
474
c            Load manager: Check whether itime0 can be copied from itime1
475
             if ( itime0.eq.iloaded1 ) then
476
                f3t0     = f3t1
477
                p3t0     = p3t1
478
                spt0     = spt1
479
                iloaded0 = itime0
480
             endif
481
 
482
c            Load manager: Check whether itime1 can be copied from itime0
483
             if ( itime1.eq.iloaded0 ) then
484
                f3t1     = f3t0
485
                p3t1     = p3t0
486
                spt1     = spt0
487
                iloaded1 = itime1
488
             endif
489
 
46 michaesp 490
c            Load manager:  Load first time (tracing variable and grid) 
3 michaesp 491
             if ( itime0.ne.iloaded0 ) then
492
 
493
                filename = tfil(i)//dat(itime0)
494
                call frac2hhmm(time0,tload)
495
                varname  = tvar(i) 
496
                write(*,'(a23,a20,a3,a5,f7.2)') 
497
     >               '    ->  loading          : ',
498
     >               trim(filename),'  ',trim(varname),tload
499
                call input_open (fid,filename)                
500
                call input_wind 
501
     >               (fid,varname,f3t0,tload,stagz,mdv,
502
     >               xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck) 
503
                call input_grid      
504
     >               (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
505
     >                tload,pollon,pollat,p3t0,spt0,nz,ak,bk,stagz,
506
     >                timecheck)
507
                call input_close(fid)
508
 
509
                iloaded0 = itime0
13 michaesp 510
 
3 michaesp 511
             endif
512
 
513
c            Load manager: Load second time (tracing variable and grid)
514
             if ( itime1.ne.iloaded1 ) then
515
 
516
                filename = tfil(i)//dat(itime1)
517
                call frac2hhmm(time1,tload)
518
                varname  = tvar(i) 
519
                write(*,'(a23,a20,a3,a5,f7.2)') 
520
     >               '    ->  loading          : ',
521
     >               trim(filename),'  ',trim(varname),tload
522
                call input_open (fid,filename)
523
                call input_wind 
524
     >               (fid,varname,f3t1,tload,stagz,mdv,
525
     >               xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)      
526
                call input_grid      
527
     >               (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
528
     >               tload,pollon,pollat,p3t1,spt1,nz,ak,bk,stagz,
529
     >               timecheck)
530
                call input_close(fid)
531
 
532
                iloaded1 = itime1
533
 
534
             endif
535
 
536
c            Loop over all trajectories
537
             do k=1,ntra
538
 
5 michaesp 539
c               Set the trajectory position 
540
                x0_tra = trainp(k,j,2)                          ! Longitude
541
                y0_tra = trainp(k,j,3)                          ! Latitude
542
                p0_tra = trainp(k,j,4)                          ! Pressure
3 michaesp 543
 
5 michaesp 544
c               Get rotation angle - orient normal to trajectory
545
                if ( direction.eq.'normal' ) then
546
 
547
                   vx0 = 1.                     
548
                   vy0 = 0.
549
 
550
                   if ( j.lt.ntim ) then
551
                      lat =  0.5 * ( trainp(k,j,3) + trainp(k,j+1,3) )
552
                      vx1 = ( trainp(k,j+1,2) - trainp(k,j,2) ) * 
553
     >                      cos( lat * pi180 )
554
                      vy1 = ( trainp(k,j+1,3) - trainp(k,j,3) )
555
                   else
556
                      lat =  0.5 * ( trainp(k,j,3) + trainp(k,j-1,3) )
557
                      vx1 = ( trainp(k,j,2) - trainp(k,j-1,2) ) * 
558
     >                      cos( lat * pi180 )
559
                      vy1 = ( trainp(k,j,3) - trainp(k,j-1,3) )
560
                   endif
561
 
562
                   if ( vx1.gt.180  ) vx1 = vx1 - 360
563
                   if ( vx1.lt.-180 ) vx1 = vx1 + 360.
564
 
565
                   call getangle (vx0,vy0,vx1,vy1,rotation)
566
                   rotation = -rotation
567
 
568
                else
569
                   rotation = 0.
570
                endif
571
 
3 michaesp 572
c               Set the relative time
573
                call hhmm2frac(trainp(k,j,1),tfrac)
574
                reltpos0 = fbflag * (tfrac-time0)/timeinc    
575
 
5 michaesp 576
c               Loop over pressure profile (or other positions for horizontal mode)
577
	        do l=1,npre
3 michaesp 578
 
5 michaesp 579
c                     Vertical
580
                      if ( direction.eq.'vertical' ) then
581
                        x0 = x0_tra
582
                        y0 = y0_tra
583
                        p0 = pmin + real(l-1)/real(npre-1) * (pmax-pmin)
584
                        if ( centering.eq.'yes' )then
585
                           p0 = p0 + trainp(k,j,4)
586
                        endif
3 michaesp 587
 
5 michaesp 588
c                     Longitude
589
                      elseif ( direction.eq.'lon' ) then
590
                        x0 = pmin + real(l-1)/real(npre-1) * (pmax-pmin)
591
                        y0 = y0_tra
592
                        p0 = p0_tra
593
                        if ( centering.eq.'yes' )then
594
                           x0 = x0 + x0_tra
595
                        endif
596
 
597
c                     Latitude
598
                      elseif ( direction.eq.'lat' ) then
599
                        x0 = x0_tra
600
                        y0 = pmin + real(l-1)/real(npre-1) * (pmax-pmin)
601
                        p0 = p0_tra
602
                        if ( centering.eq.'yes' )then
603
                           y0 = y0 + y0_tra
604
                        endif
3 michaesp 605
 
5 michaesp 606
c                     Normal to trajerctory
607
                      elseif ( direction.eq.'normal' ) then
608
 
609
c                        Set the coordinate in the rotated system
610
                         rlat = pmin + 
611
     >                             real(l-1)/real(npre-1) * (pmax-pmin)
612
                         rlon = 0.
613
 
614
c                        Transform it back to geographical lon/lat
615
                         call getenvir_b (x0_tra,y0_tra,rotation,
616
     >                                              x0,y0,rlon,rlat,1)
617
 
618
c                        Pressure unchanged
619
                         p0 = p0_tra
620
 
621
                      endif
622
 
623
c                     Handle periodic boundaries in zonal direction
624
                      if ( (x0.gt.xmax).and.(per.ne.0) ) x0 = x0 - 360.
625
                      if ( (x0.lt.xmin).and.(per.ne.0) ) x0 = x0 + 360.
626
 
627
c                     Handle pole problems for hemispheric data (taken from caltra.f)
628
                      if ((hem.eq.1).and.(y0.gt.90.)) then
629
                         y0=180.-y0
630
                         x0=x0+per/2.
631
                      endif
632
                      if ((hem.eq.1).and.(y0.lt.-90.)) then
633
                         y0=-180.-y0
634
                         x0=x0+per/2.
635
                      endif
636
                      if (y0.gt.89.99) then
637
                         y0=89.99
638
                      endif     
639
 
640
c                 If requested, dump the lidar coordinates
641
                  if ( (dumpcoord.eq.'yes').and.(i.eq.1) ) then
9 michaesp 642
                     write(10,'(3f10.2)') x0,y0,trainp(k,j,1)
643
                     write(10,'(3f10.2)') x0_tra,y0_tra,5.
5 michaesp 644
                  endif
645
 
3 michaesp 646
C                 Get the index where to interpolate (x0,y0,p0)
647
                  if ( (abs(x0-mdv).gt.eps).and.
648
     >                 (abs(y0-mdv).gt.eps) )
649
     >            then
650
                     call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0,
651
     >                                p3t0,p3t1,spt0,spt1,3,
652
     >                                nx,ny,nz,xmin,ymin,dx,dy,mdv)
653
                  else
654
                     xind = mdv
655
                     yind = mdv
656
                     pind = mdv
657
                  endif
658
 
659
c                 If requested, apply nearest-neighbor interpolation
660
                  if ( intmode.eq.'nearest') then
661
 
662
                     xind = real( nint(xind) )
663
                     yind = real( nint(yind) )
664
                     pind = real( nint(pind) )
665
 
666
                     if ( xind.lt.1.  ) xind = 1.
667
                     if ( xind.gt.nx  ) xind = real(nx)
668
                     if ( yind.lt.1.  ) yind = 1.
669
                     if ( yind.gt.ny  ) yind = real(ny)
670
 
671
                     if ( pind.lt.1.  ) pind = 1.
672
                     if ( pind.gt.nz  ) pind = real(nz)
673
 
674
                  endif
675
 
676
c                 Do the interpolation: everthing is ok
677
                  if ( (xind.ge.1.).and.(xind.le.real(nx)).and.
678
     >                 (yind.ge.1.).and.(yind.le.real(ny)).and.
679
     >                 (pind.ge.1.).and.(pind.le.real(nz)) )
680
     >            then
681
                     f0 = int_index4(f3t0,f3t1,nx,ny,nz,
682
     >                               xind,yind,pind,reltpos0,mdv)
683
 
684
c                 Set to missing data
685
                  else
686
                     f0       = mdv
687
                  endif
688
 
689
c	              Save result to output array
690
                  if (abs(f0-mdv).gt.eps) then
691
                     out_val(j,l) = out_val(j,l) + f0 * fac(i)
692
                     out_cnt(j,l) = out_cnt(j,l) + 1.
693
 
694
	              endif
695
 
696
c              End loop over all pressure levels
697
	           enddo
698
 
5 michaesp 699
c	           Save output - time index
3 michaesp 700
	           ind_time = j
5 michaesp 701
 
702
c                  Save output - space index for 'no centering'
3 michaesp 703
	           if ( centering.eq.'no' ) then
5 michaesp 704
                      if ( direction.eq.'vertical') then 
705
                         ind_pre  = nint( real(npre) *
706
     >          	       ( (p0_tra - pmin)/(pmax-pmin) ) + 1.)
707
                      elseif ( direction.eq.'lon') then 
708
                         ind_pre  = nint( real(npre) *
709
     >          	       ( (x0_tra - pmin)/(pmax-pmin) ) + 1.)
710
                      elseif ( direction.eq.'lat') then 
711
                         ind_pre  = nint( real(npre) *
712
     >          	       ( (y0_tra - pmin)/(pmax-pmin) ) + 1.)
713
                      endif
714
 
715
c                  Save output - space index for 'centering'
3 michaesp 716
	           else
717
	           	  ind_pre  = nint( real(npre) *
718
     >          	       ( (0.            - pmin)/(pmax-pmin) ) + 1.)
719
	           endif
720
 
5 michaesp 721
c                  Update the output array
3 michaesp 722
	           if ( (ind_time.ge.1).and.(ind_time.le.ntim).and.
723
     >              (ind_pre .ge.1).and.(ind_pre .le.npre) )
724
     >         then
725
                    out_pos(ind_time,ind_pre) =
726
     >                          	out_pos(ind_time,ind_pre) + 1.
727
	           endif
728
 
729
c	         End loop over all trajectories
730
             enddo
731
 
732
c	      End loop over all times
733
          enddo
734
 
735
c	      Write the trajectory position to netCDF file - only once
736
	      if ( i.eq.1 ) then
737
	      	  cdfname  = outfile
738
	      	  varname  = 'POSITION'
739
	      	  longname = 'position of trajectory points'
740
	      	  unit     = 'none'
741
	      	  time     = 0.
742
              do k=1,npre
743
              	levels(k) = pmin + real(k-1)/real(npre-1) * (pmax-pmin)
744
              enddo
745
              do k=1,ntim
746
                 times(k) = trainp(1,k,1)
747
              enddo
748
              call writecdf2D_cf
749
     >            (cdfname,varname,longname,unit,out_pos,time,levels,
5 michaesp 750
     >             times,npre,ntim,1,1,direction)
3 michaesp 751
	      endif
752
 
753
c	      If no valid lidar count: set the field to missing data
754
          do k=1,ntim
755
          	do l=1,npre
756
          		if (abs(out_cnt(k,l)).lt.eps) then
757
          			out_val(k,l) = mdv
758
          	    endif
759
          	 enddo
760
          enddo
761
 
762
c	      If requested, calculate the mean of the lidar field
763
	      if ( outmode.eq.'mean' ) then
764
	      	do k=1,ntim
765
	      		do l=1,npre
766
	      			if ( (abs(out_val(k,l)-mdv).gt.eps).and.
767
     >                   (abs(out_cnt(k,l)    ).gt.0. ) )
768
     >              then
769
	      				out_val(k,l) = out_val(k,l) / out_cnt(k,l)
770
	      		    endif
771
	      		 enddo
772
	          enddo
773
	      endif
774
 
775
c	      Write the lidar field and count
776
	      cdfname  = outfile
777
	      if (outmode.eq.'sum' ) then
778
	         varname  = trim(tvar(i))//'_SUM'
779
	      elseif (outmode.eq.'mean' ) then
780
	         varname  = trim(tvar(i))//'_MEAN'
781
	      endif
782
	      longname = 'sum over all '//trim(tvar(i))//' profiles'
783
	      unit     = 'not given'
784
	      time     = 0.
785
          call writecdf2D_cf
786
     >            (cdfname,varname,longname,unit,out_val,time,levels,
5 michaesp 787
     >             times,npre,ntim,0,1,direction)
3 michaesp 788
 
789
	  	  cdfname  = outfile
790
	      varname  = trim(tvar(i))//'_CNT'
791
	      longname = 'counts of all '//trim(tvar(i))//' profiles'
792
	      unit     = 'not given'
793
	      time     = 0.
794
          call writecdf2D_cf
795
     >            (cdfname,varname,longname,unit,out_cnt,time,levels,
5 michaesp 796
     >             times,npre,ntim,0,1,direction)
3 michaesp 797
 
798
c         Exit point for loop over all tracing variables
799
 110      continue
800
 
801
c	   End loop over all lidar variables
802
       enddo
803
 
804
 
805
c     --------------------------------------------------------------------
806
c     Write output to netCDF file
807
c     --------------------------------------------------------------------
808
 
809
c     Write status information
810
      print*
811
      print*,'---- WRITE OUTPUT LIDAR FIELDS --------------------------'
812
      print*
813
 
5 michaesp 814
c     Close coord dump file
815
      print*,' LIDAR written to      : ',trim(outfile)
816
      if ( dumpcoord.eq.'yes' ) then
817
        print*,' Coordinates dumped to : ',trim(outfile)//'.coord'
818
      endif
3 michaesp 819
 
820
c     Write some status information, and end of program message
821
      print*  
822
      print*,'---- STATUS INFORMATION --------------------------------'
823
      print*
824
      print*,' ok'
825
      print*
826
      print*,'              *** END OF PROGRAM LIDAR ***'
827
      print*,'========================================================='
828
 
829
 
830
      end 
831
 
832
 
833
c     ********************************************************************
834
c     * INPUT / OUTPUT SUBROUTINES                                       *
835
c     ********************************************************************
836
 
837
c     --------------------------------------------------------------------
838
c     Subroutines to write 2D CF netcdf output file
839
c     --------------------------------------------------------------------
840
 
841
      subroutine writecdf2D_cf
842
     >          (cdfname,varname,longname,unit,arr,time,levels,times,
5 michaesp 843
     >           npre,ntim,crefile,crevar,direction)
3 michaesp 844
 
845
c     Create and write to the CF netcdf file <cdfname>. The variable
846
c     with name <varname> and with time <time> is written. The data
847
c     are in the two-dimensional array <arr>.  The flags <crefile> and
848
c     <crevar> determine whether the file and/or the variable should
849
c     be created.
850
 
851
      USE netcdf
852
 
853
      IMPLICIT NONE
854
 
855
c     Declaration of input parameters
856
      character*80 cdfname
857
      character*80 varname,longname,unit
858
      integer      npre,ntim
859
      real         arr(ntim,npre)
860
      real         levels(npre)
861
      real         times (ntim)
862
      real         time
863
      integer      crefile,crevar
5 michaesp 864
      character*80 direction
3 michaesp 865
 
866
c     Numerical epsilon
867
      real         eps
868
      parameter    (eps=1.e-5)
869
 
870
c     Local variables
871
      integer      ierr
872
      integer      ncID
873
      integer      LevDimId,    varLevID
874
      integer      TimeDimID,   varTimeID
875
      real         timeindex
876
      integer      i,j
877
      integer      nvars,varids(100)
878
      integer      ndims,dimids(100)
879
      real         timelist(1000)
880
      integer      ntimes
881
      integer      ind
882
      integer      varID
883
 
884
c     Quick an dirty solution for fieldname conflict
885
      if ( varname.eq.'time' ) varname = 'TIME'
886
 
887
c     Initially set error to indicate no errors.
888
      ierr = 0
889
 
890
c     ---- Create the netCDF - skip if <crefile=0> ----------------------
891
      if ( crefile.ne.1 ) goto 100
892
 
893
c     Create the file
894
      ierr = nf90_create(trim(cdfname), NF90_CLOBBER, ncID)
895
 
896
c     Define dimensions
897
      ierr=nf90_def_dim(ncID,'level',npre, LevDimID )
898
      ierr=nf90_def_dim(ncID,'time' ,ntim, TimeDimID)
899
 
5 michaesp 900
c     Define space coordinate 
3 michaesp 901
      ierr = nf90_def_var(ncID,'level',NF90_FLOAT,
902
     >     (/ LevDimID /),varLevID)
5 michaesp 903
      if ( direction.eq.'vertical' ) then
904
        ierr = nf90_put_att(ncID, varLevID, "standard_name","level")
905
        ierr = nf90_put_att(ncID, varLevID, "units"        ,"hPa")
906
      elseif ( direction.eq.'lat' ) then
907
        ierr = nf90_put_att(ncID, varLevID, "standard_name","latitude")
908
        ierr = nf90_put_att(ncID, varLevID, "units"        ,"deg")
909
      elseif ( direction.eq.'lon' ) then
910
        ierr = nf90_put_att(ncID, varLevID, "standard_name","longitude")
911
        ierr = nf90_put_att(ncID, varLevID, "units"        ,"deg")
912
      elseif ( direction.eq.'normal' ) then
913
        ierr = nf90_put_att(ncID, varLevID, "standard_name","normal")
914
        ierr = nf90_put_att(ncID, varLevID, "units"        ,"deg")
915
      endif
3 michaesp 916
 
5 michaesp 917
c     Define time coordinate
3 michaesp 918
      ierr = nf90_def_var(ncID,'time',NF90_FLOAT,
919
     >     (/ TimeDimID /), varTimeID)
920
      ierr = nf90_put_att(ncID, varTimeID, "long_name",    "time")
921
      ierr = nf90_put_att(ncID, varTimeID, "units",       "hours")
922
 
923
c     Write global attributes
924
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'Conventions', 'CF-1.0')
925
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'title',
926
     >     'pseudo-lidar from trajectory file')
927
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'source',
928
     >     'Lagranto Trajectories')
929
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'institution',
930
     >     'ETH Zurich, IACETH')
931
 
932
c     Check whether the definition was successful
933
      ierr = nf90_enddef(ncID)
934
      if (ierr.gt.0) then
935
         print*, 'An error occurred while attempting to ',
936
     >        'finish definition mode.'
937
         stop
938
      endif
939
 
940
c     Write coordinate data
941
      ierr = nf90_put_var(ncID,varLevID  ,levels)
942
      ierr = nf90_put_var(ncID,varTimeID ,times )
943
 
944
c     Close netCDF file
945
      ierr = nf90_close(ncID)
946
 
947
 100  continue
948
 
949
c     ---- Define a new variable - skip if <crevar=0> -----------------------
950
 
951
      if ( crevar.ne.1 ) goto 110
952
 
5 michaesp 953
      print*,'Now defining ',trim(varname)
954
 
3 michaesp 955
c     Open the file for read/write access
956
      ierr = nf90_open  (trim(cdfname), NF90_WRITE  , ncID)
957
 
958
c     Get the IDs for dimensions
959
      ierr = nf90_inq_dimid(ncID,'level', LevDimID )
960
      ierr = nf90_inq_dimid(ncID,'time' , TimeDimID)
961
 
962
c     Enter define mode
963
      ierr = nf90_redef(ncID)
964
 
965
c     Write definition and add attributes
966
      ierr = nf90_def_var(ncID,varname,NF90_FLOAT,
5 michaesp 967
     >                   (/ TimeDimID, LevDimID /),varID)
3 michaesp 968
      ierr = nf90_put_att(ncID, varID, "long_name" , longname )
969
      ierr = nf90_put_att(ncID, varID, "units"     , unit     )
970
      ierr = nf90_put_att(ncID, varID, '_FillValue', -999.99  )
971
 
972
c     Check whether definition was successful
973
      ierr = nf90_enddef(ncID)
974
      if (ierr.gt.0) then
975
         print*, 'An error occurred while attempting to ',
976
     >           'finish definition mode.'
977
         stop
978
      endif
979
 
980
c     Close netCDF file
981
      ierr = nf90_close(ncID)
982
 
983
 110  continue
984
 
985
c     ---- Write data --------------------------------------------------
986
 
987
c     Open the file for read/write access
988
      ierr = nf90_open  (trim(cdfname), NF90_WRITE , ncID)
989
 
990
c     Get the varID
991
      ierr = nf90_inq_varid(ncID,varname, varID )
992
      if (ierr.ne.0) then
993
         print*,'Variable ',trim(varname),' is not defined on ',
994
     >          trim(cdfname)
995
         stop
996
      endif
997
 
998
c     Write data block
999
      ierr = nf90_put_var(ncID,varID,arr,
1000
     >                    start = (/ 1, 1 /),
1001
     >                    count = (/ ntim, npre/) )
1002
 
1003
c     Check whether writing was successful
1004
      ierr = nf90_close(ncID)
1005
      if (ierr.ne.0) then
1006
         write(*,*) trim(nf90_strerror(ierr))
1007
         write(*,*) 'An error occurred while attempting to ',
1008
     >              'close the netcdf file.'
1009
         write(*,*) 'in clscdf_CF'
1010
      endif
1011
 
1012
      end
1013
 
5 michaesp 1014
c     ********************************************************************************
1015
c     * Coordinate rotation - lidar normal to trajectory                             *
1016
c     ********************************************************************************
3 michaesp 1017
 
5 michaesp 1018
c     --------------------------------------------------------------------------------
1019
c     Backward coordinate transformation (Rotated lon/lat -> True lon/lat)
1020
c     --------------------------------------------------------------------------------
3 michaesp 1021
 
5 michaesp 1022
      SUBROUTINE getenvir_b (clon,clat,rotation,
1023
     >                       lon,lat,rlon,rlat,n)
1024
 
1025
      implicit none
1026
 
1027
c     Declaration of input and output parameters
1028
      integer     n
1029
      real        clon,clat,rotation
1030
      real        lon(n), lat(n)
1031
      real        rlon(n),rlat(n)
1032
 
1033
c     Auxiliary variables 
1034
      real         pollon,pollat
1035
      integer      i
1036
      real         rlon1(n),rlat1(n)
1037
 
1038
c     Externals
1039
      real         lmstolm,phstoph
1040
      external     lmstolm,phstoph
1041
 
1042
c     First coordinate transformation (make the local coordinate system parallel to equator)
1043
      pollon=-180.
1044
      pollat=90.+rotation
1045
      do i=1,n
1046
         rlon1(i)=90.+lmstolm(rlat(i),rlon(i)-90.,pollat,pollon)
1047
         rlat1(i)=phstoph(rlat(i),rlon(i)-90.,pollat,pollon)            
1048
      enddo
1049
 
1050
c     Second coordinate transformation (make the local coordinate system parallel to equator)
1051
      pollon=clon-180.
1052
      if (pollon.lt.-180.) pollon=pollon+360.
1053
      pollat=90.-clat
1054
      do i=1,n
1055
         lon(i)=lmstolm(rlat1(i),rlon1(i),pollat,pollon)
1056
         lat(i)=phstoph(rlat1(i),rlon1(i),pollat,pollon)            
1057
      enddo
1058
 
1059
      END
1060
 
1061
c     ---------------------------------------------------------------------
1062
c     Determine the angle between two vectors
1063
c     ---------------------------------------------------------------------
1064
 
1065
      SUBROUTINE getangle (vx1,vy1,vx2,vy2,angle)
1066
 
1067
c     Given two vectors <vx1,vy1> and <vx2,vy2>, determine the angle (in deg) 
1068
c     between the two vectors.
1069
 
1070
      implicit none
1071
 
1072
c     Declaration of subroutine parameters
1073
      real vx1,vy1
1074
      real vx2,vy2
1075
      real angle
1076
 
1077
c     Auxiliary variables and parameters
1078
      real len1,len2,len3
1079
      real val1,val2,val3
1080
      real pi
1081
      parameter (pi=3.14159265359)
1082
 
1083
      len1=sqrt(vx1*vx1+vy1*vy1)
1084
      len2=sqrt(vx2*vx2+vy2*vy2)
1085
 
1086
      if ((len1.gt.0.).and.(len2.gt.0.)) then
1087
         vx1=vx1/len1
1088
         vy1=vy1/len1
1089
         vx2=vx2/len2
1090
         vy2=vy2/len2
3 michaesp 1091
 
5 michaesp 1092
         val1=vx1*vx2+vy1*vy2
1093
         val2=-vy1*vx2+vx1*vy2
3 michaesp 1094
 
5 michaesp 1095
         len3=sqrt(val1*val1+val2*val2)
1096
 
1097
         if ( (val1.ge.0.).and.(val2.ge.0.) ) then
1098
            val3=acos(val1/len3)
1099
         else if ( (val1.lt.0.).and.(val2.ge.0.) ) then
1100
            val3=pi-acos(abs(val1)/len3)
1101
         else if ( (val1.ge.0.).and.(val2.le.0.) ) then
1102
            val3=-acos(val1/len3)
1103
         else if ( (val1.lt.0.).and.(val2.le.0.) ) then
1104
            val3=-pi+acos(abs(val1)/len3)
1105
         endif
1106
      else
1107
         val3=0.
1108
      endif
3 michaesp 1109
 
5 michaesp 1110
      angle=180./pi*val3                                                                                     
3 michaesp 1111
 
5 michaesp 1112
      END
3 michaesp 1113
 
5 michaesp 1114
c     --------------------------------------------------------------------------------
1115
c     Transformation routine: LMSTOLM and PHSTOPH from library gm2em            
1116
c     --------------------------------------------------------------------------------
1117
 
1118
      REAL FUNCTION LMSTOLM (PHIS, LAMS, POLPHI, POLLAM)
1119
C
1120
C**** LMSTOLM  -   FC:BERECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE FUER
1121
C****                 EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
1122
C****                 IM ROTIERTEN SYSTEM. DER NORDPOL DES SYSTEMS HAT
1123
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1124
C**   AUFRUF   :   LAM = LMSTOLM (PHIS, LAMS, POLPHI, POLLAM)
1125
C**   ENTRIES  :   KEINE
1126
C**   ZWECK    :   BERECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE FUER
1127
C**                EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
1128
C**                IM ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
1129
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1130
C**   VERSIONS-
1131
C**   DATUM    :   03.05.90
1132
C**
1133
C**   EXTERNALS:   KEINE
1134
C**   EINGABE-
1135
C**   PARAMETER:   PHIS     REAL   GEOGR. BREITE DES PUNKTES IM ROT.SYS.
1136
C**                LAMS     REAL   GEOGR. LAENGE DES PUNKTES IM ROT.SYS.
1137
C**                POLPHI   REAL   WAHRE GEOGR. BREITE DES NORDPOLS
1138
C**                POLLAM   REAL   WAHRE GEOGR. LAENGE DES NORDPOLS
1139
C**   AUSGABE-
1140
C**   PARAMETER:   WAHRE GEOGRAPHISCHE LAENGE ALS WERT DER FUNKTION
1141
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
1142
C**
1143
C**   COMMON-
1144
C**   BLOECKE  :   KEINE
1145
C**
1146
C**   FEHLERBE-
1147
C**   HANDLUNG :   KEINE
1148
C**   VERFASSER:   D.MAJEWSKI
1149
 
1150
      REAL        LAMS,PHIS,POLPHI,POLLAM
1151
 
1152
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
1153
 
1154
      ZSINPOL = SIN(ZPIR18*POLPHI)
1155
      ZCOSPOL = COS(ZPIR18*POLPHI)
1156
      ZLAMPOL = ZPIR18*POLLAM
1157
      ZPHIS   = ZPIR18*PHIS
1158
      ZLAMS   = LAMS
1159
      IF(ZLAMS.GT.180.0) ZLAMS = ZLAMS - 360.0
1160
      ZLAMS   = ZPIR18*ZLAMS
1161
 
1162
      ZARG1   = SIN(ZLAMPOL)*(- ZSINPOL*COS(ZLAMS)*COS(ZPHIS)  +
1163
     1                          ZCOSPOL*           SIN(ZPHIS)) -
1164
     2          COS(ZLAMPOL)*           SIN(ZLAMS)*COS(ZPHIS)
1165
      ZARG2   = COS(ZLAMPOL)*(- ZSINPOL*COS(ZLAMS)*COS(ZPHIS)  +
1166
     1                          ZCOSPOL*           SIN(ZPHIS)) +
1167
     2          SIN(ZLAMPOL)*           SIN(ZLAMS)*COS(ZPHIS)
1168
      IF (ABS(ZARG2).LT.1.E-30) THEN
1169
        IF (ABS(ZARG1).LT.1.E-30) THEN
1170
          LMSTOLM =   0.0
1171
        ELSEIF (ZARG1.GT.0.) THEN
1172
              LMSTOLAM =  90.0
1173
            ELSE
1174
              LMSTOLAM = -90.0
1175
            ENDIF
1176
      ELSE
1177
        LMSTOLM = ZRPI18*ATAN2(ZARG1,ZARG2)
1178
      ENDIF
1179
 
1180
      RETURN
1181
      END
1182
 
1183
 
1184
      REAL FUNCTION PHSTOPH (PHIS, LAMS, POLPHI, POLLAM)
1185
C
1186
C**** PHSTOPH  -   FC:BERECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE FUER
1187
C****                 EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
1188
C****                 ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
1189
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1190
C**   AUFRUF   :   PHI = PHSTOPH (PHIS, LAMS, POLPHI, POLLAM)
1191
C**   ENTRIES  :   KEINE
1192
C**   ZWECK    :   BERECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE FUER
1193
C**                EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
1194
C**                ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
1195
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1196
C**   VERSIONS-
1197
C**   DATUM    :   03.05.90
1198
C**
1199
C**   EXTERNALS:   KEINE
1200
C**   EINGABE-
1201
C**   PARAMETER:   PHIS     REAL   GEOGR. BREITE DES PUNKTES IM ROT.SYS.
1202
C**                LAMS     REAL   GEOGR. LAENGE DES PUNKTES IM ROT.SYS.
1203
C**                POLPHI   REAL   WAHRE GEOGR. BREITE DES NORDPOLS
1204
C**                POLLAM   REAL   WAHRE GEOGR. LAENGE DES NORDPOLS
1205
C**   AUSGABE-
1206
C**   PARAMETER:   WAHRE GEOGRAPHISCHE BREITE ALS WERT DER FUNKTION
1207
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
1208
C**
1209
C**   COMMON-
1210
C**   BLOECKE  :   KEINE
1211
C**
1212
C**   FEHLERBE-
1213
C**   HANDLUNG :   KEINE
1214
C**   VERFASSER:   D.MAJEWSKI
1215
 
1216
      REAL        LAMS,PHIS,POLPHI,POLLAM
1217
 
1218
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
1219
 
1220
      SINPOL = SIN(ZPIR18*POLPHI)
1221
      COSPOL = COS(ZPIR18*POLPHI)
1222
      ZPHIS  = ZPIR18*PHIS
1223
      ZLAMS  = LAMS
1224
      IF(ZLAMS.GT.180.0) ZLAMS = ZLAMS - 360.0
1225
      ZLAMS  = ZPIR18*ZLAMS
1226
      ARG     = COSPOL*COS(ZPHIS)*COS(ZLAMS) + SINPOL*SIN(ZPHIS)
1227
 
1228
      PHSTOPH = ZRPI18*ASIN(ARG)
1229
 
1230
      RETURN
1231
      END
1232
 
1233
 
1234
      REAL FUNCTION LMTOLMS (PHI, LAM, POLPHI, POLLAM)
1235
C
1236
C%Z% Modul %M%, V%I% vom %G%, extrahiert am %H%
1237
C
1238
C**** LMTOLMS  -   FC:UMRECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE LAM
1239
C****                 AUF EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
1240
C****                 IM ROTIERTEN SYSTEM. DER NORDPOL DES SYSTEMS HAT
1241
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1242
C**   AUFRUF   :   LAM = LMTOLMS (PHI, LAM, POLPHI, POLLAM)
1243
C**   ENTRIES  :   KEINE
1244
C**   ZWECK    :   UMRECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE LAM AUF
1245
C**                EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
1246
C**                ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
1247
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1248
C**   VERSIONS-
1249
C**   DATUM    :   03.05.90
1250
C**
1251
C**   EXTERNALS:   KEINE
1252
C**   EINGABE-
1253
C**   PARAMETER:   PHI    REAL BREITE DES PUNKTES IM GEOGR. SYSTEM
1254
C**                LAM    REAL LAENGE DES PUNKTES IM GEOGR. SYSTEM
1255
C**                POLPHI REAL GEOGR.BREITE DES N-POLS DES ROT. SYSTEMS
1256
C**                POLLAM REAL GEOGR.LAENGE DES N-POLS DES ROT. SYSTEMS
1257
C**   AUSGABE-
1258
C**   PARAMETER:   WAHRE GEOGRAPHISCHE LAENGE ALS WERT DER FUNKTION
1259
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
1260
C**
1261
C**   COMMON-
1262
C**   BLOECKE  :   KEINE
1263
C**
1264
C**   FEHLERBE-
1265
C**   HANDLUNG :   KEINE
1266
C**   VERFASSER:   G. DE MORSIER
1267
 
1268
      REAL        LAM,PHI,POLPHI,POLLAM
1269
 
1270
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
1271
 
1272
      ZSINPOL = SIN(ZPIR18*POLPHI)
1273
      ZCOSPOL = COS(ZPIR18*POLPHI)
1274
      ZLAMPOL =     ZPIR18*POLLAM
1275
      ZPHI    =     ZPIR18*PHI
1276
      ZLAM    = LAM
1277
      IF(ZLAM.GT.180.0) ZLAM = ZLAM - 360.0
1278
      ZLAM    = ZPIR18*ZLAM
1279
 
1280
      ZARG1   = - SIN(ZLAM-ZLAMPOL)*COS(ZPHI)
1281
      ZARG2   = - ZSINPOL*COS(ZPHI)*COS(ZLAM-ZLAMPOL)+ZCOSPOL*SIN(ZPHI)
1282
      IF (ABS(ZARG2).LT.1.E-30) THEN
1283
        IF (ABS(ZARG1).LT.1.E-30) THEN
1284
          LMTOLMS =   0.0
1285
        ELSEIF (ZARG1.GT.0.) THEN
1286
              LMTOLMS =  90.0
1287
            ELSE
1288
              LMTOLMS = -90.0
1289
            ENDIF
1290
      ELSE
1291
        LMTOLMS = ZRPI18*ATAN2(ZARG1,ZARG2)
1292
      ENDIF
1293
 
1294
      RETURN
1295
      END
1296
 
1297
 
1298
      REAL FUNCTION PHTOPHS (PHI, LAM, POLPHI, POLLAM)
1299
C
1300
C%Z% Modul %M%, V%I% vom %G%, extrahiert am %H%
1301
C
1302
C**** PHTOPHS  -   FC:UMRECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE PHI
1303
C****                 AUF EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
1304
C****                 IM ROTIERTEN SYSTEM. DER NORDPOL DES SYSTEMS HAT
1305
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1306
C**   AUFRUF   :   PHI = PHTOPHS (PHI, LAM, POLPHI, POLLAM)
1307
C**   ENTRIES  :   KEINE
1308
C**   ZWECK    :   UMRECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE PHI AUF
1309
C**                EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
1310
C**                ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
1311
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
1312
C**   VERSIONS-
1313
C**   DATUM    :   03.05.90
1314
C**
1315
C**   EXTERNALS:   KEINE
1316
C**   EINGABE-
1317
C**   PARAMETER:   PHI    REAL BREITE DES PUNKTES IM GEOGR. SYSTEM
1318
C**                LAM    REAL LAENGE DES PUNKTES IM GEOGR. SYSTEM
1319
C**                POLPHI REAL GEOGR.BREITE DES N-POLS DES ROT. SYSTEMS
1320
C**                POLLAM REAL GEOGR.LAENGE DES N-POLS DES ROT. SYSTEMS
1321
C**   AUSGABE-
1322
C**   PARAMETER:   ROTIERTE BREITE PHIS ALS WERT DER FUNKTION
1323
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
1324
C**
1325
C**   COMMON-
1326
C**   BLOECKE  :   KEINE
1327
C**
1328
C**   FEHLERBE-
1329
C**   HANDLUNG :   KEINE
1330
C**   VERFASSER:   G. DE MORSIER
1331
 
1332
      REAL        LAM,PHI,POLPHI,POLLAM
1333
 
1334
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
1335
 
1336
      ZSINPOL = SIN(ZPIR18*POLPHI)
1337
      ZCOSPOL = COS(ZPIR18*POLPHI)
1338
      ZLAMPOL = ZPIR18*POLLAM
1339
      ZPHI    = ZPIR18*PHI
1340
      ZLAM    = LAM
1341
      IF(ZLAM.GT.180.0) ZLAM = ZLAM - 360.0
1342
      ZLAM    = ZPIR18*ZLAM
1343
      ZARG    = ZCOSPOL*COS(ZPHI)*COS(ZLAM-ZLAMPOL) + ZSINPOL*SIN(ZPHI)
1344
 
1345
      PHTOPHS = ZRPI18*ASIN(ZARG)
1346
 
1347
      RETURN
1348
      END
3 michaesp 1349