Subversion Repositories lagranto.ecmwf

Rev

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

Rev Author Line No. Line
21 michaesp 1
      PROGRAM getmima
3 michaesp 2
 
21 michaesp 3
c     ***********************************************************************
4
c     * Get minimum and maximum value of a varaiable on pressure            *
5
c     * or isentropic surface                                               *
6
c     * Michael Sprenger / Spring, summer 2016                              *
7
c     ***********************************************************************
3 michaesp 8
 
21 michaesp 9
      use netcdf
10
 
3 michaesp 11
      implicit none
12
 
21 michaesp 13
c     ----------------------------------------------------------------------
14
c     Declaration of parameters and variables
15
c     ----------------------------------------------------------------------
16
 
17
c     Interpolation method
18
      integer   ipom
19
      parameter (ipom=0)    
20
 
21
c     Flag for timecheck
22
      character*80 timecheck
23
      parameter    (timecheck = 'no' )        
24
 
25
c     netCDF fields      
26
      real,allocatable, dimension (:,:)   :: sp,varf
27
      real,allocatable, dimension (:,:,:) :: field,varl,tt,p3
28
      real,allocatable, dimension (:)     :: ak,bk
29
      integer                                stat
30
      character*80                           cdfnam,varnam
31
      character*80                           vnam(200)
32
      integer                                nvars
33
      integer	                             cdfid,ierr
34
      real	                                 mdv
35
      real                                   stagz
36
 
37
c     Grid description      
38
      real      xmin,xmax       ! Zonal grid extension
39
      real      ymin,ymax       ! Meridional grid extension
40
      integer   nx,ny,nz        ! Grid dimensions
41
      real      dx,dy           ! Horizontal grid resolution
42
      real	    pollon,pollat   ! Pole position
3 michaesp 43
 
21 michaesp 44
c     Output variables
45
      real      min,max                      ! Minimum & maximum value
46
      real      lonmin,lonmax,latmin,latmax  ! Position of min & max
3 michaesp 47
 
21 michaesp 48
c     Auxiliary variables      
49
      integer        iargc
50
      character*(80) arg
51
      real           rd
52
      integer        i,j
53
      integer        minindx,minindy,maxindx,maxindy
54
      character*80   tvar
55
      character*1    mode
56
      character*80   clev
57
      real	         xphys,yphys
58
      real	         level
59
      integer        flag
60
 
61
c     Externals      
3 michaesp 62
      real      lmstolm,phstoph
21 michaesp 63
      external  lmstolm,phstoph
3 michaesp 64
 
21 michaesp 65
c     ----------------------------------------------------------------------
66
c     Preparation - argument handling, grid paraemters, memor allocation
67
c     ----------------------------------------------------------------------
68
 
69
c     Check for sufficient requested arguments
3 michaesp 70
      if (iargc().lt.2) then
71
        print*,'USAGE: getmima filename var ',
72
     >         '[lev in the form Pval or Tval]'
73
        call exit(1)
74
      endif
75
 
21 michaesp 76
c     Read and transform input
3 michaesp 77
      call getarg(1,arg)
78
      cdfnam=trim(arg)
79
      call getarg(2,arg)
80
      varnam=trim(arg)
81
 
82
      if (iargc().eq.3) then
83
        call getarg(3,arg)
84
        mode=arg(1:1)
85
        clev=arg(2:len(arg)) 
86
        call checkchar(clev,".",flag)
87
        if (flag.eq.0) clev=trim(clev)//"."
88
        read(clev,'(f10.2)') level
89
      else
90
        mode='X'
91
        level=0.
92
      endif
93
 
21 michaesp 94
c     Init level type and minimum,maximum
3 michaesp 95
      min= 1.e19
96
      max=-1.e19
97
 
21 michaesp 98
c     Open netCDF file      
99
      call input_open (cdfid,cdfnam)
3 michaesp 100
 
21 michaesp 101
c     Get list of variables on file
102
      call input_getvars (cdfid,vnam,nvars)
103
 
104
C     Get infos about data domain - in particular: dimensions
105
      nx       = 1
106
      ny       = 1
107
      nz       = 1
108
      call input_grid (-cdfid,varnam,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
39 michaesp 109
     >               0.,pollon,pollat,rd,rd,nz,rd,rd,stagz,timecheck)
21 michaesp 110
 
3 michaesp 111
C     Get memory for dynamic arrays
21 michaesp 112
      allocate(sp(nx,ny),stat=stat)
3 michaesp 113
      if (stat.ne.0) print*,'*** error allocating array sp ***'
21 michaesp 114
      allocate(varf(nx,ny),stat=stat)
3 michaesp 115
      if (stat.ne.0) print*,'*** error allocating array varf ***'
21 michaesp 116
      allocate(varl(nx,ny,nz),stat=stat)
3 michaesp 117
      if (stat.ne.0) print*,'*** error allocating array varl ***'
21 michaesp 118
      allocate(tt(nx,ny,nz),stat=stat)
3 michaesp 119
      if (stat.ne.0) print*,'*** error allocating array tt ***'
21 michaesp 120
      allocate(field(nx,ny,nz),stat=stat)
3 michaesp 121
      if (stat.ne.0) print*,'*** error allocating array field ***'
21 michaesp 122
      allocate(ak(nz),stat=stat)
123
      if (stat.ne.0) print*,'*** error allocating array ak ***'
124
      allocate(bk(nz),stat=stat)
125
      if (stat.ne.0) print*,'*** error allocating array bk ***'
126
 
127
c     ----------------------------------------------------------------------
128
c     Load data
129
c     ----------------------------------------------------------------------
130
 
131
c     Get grid info for variable - non-constant part  (P, PS, AK, BK)   
132
      call input_grid                               
133
     >       (cdfid,varnam,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
134
     >        0.,pollon,pollat,varl,sp,nz,ak,bk,stagz,timecheck)
135
 
136
c     Load Variable    
137
      call input_wind (cdfid,varnam,field,0.,stagz,mdv,
138
     >                 xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
23 michaesp 139
 
21 michaesp 140
c     In case of isentropic level, read temperature
141
      if (mode.eq.'T') then
142
            tvar="T"
143
            do i=1,nvars
144
              if (vnam(i).eq."THETA") tvar="THETA"
145
              if (vnam(i).eq."TH")    tvar="TH"
146
            enddo
147
            if (tvar.eq.'T') then
148
                 call input_wind (cdfid,tvar,tt,0.,stagz,mdv,xmin,xmax,
149
     >                            ymin,ymax,dx,dy,nx,ny,nz,timecheck)
150
                 call pottemp(varl,tt,sp,nx,ny,nz,ak,bk)
151
            else
152
                 call input_wind (cdfid,tvar,tt,0.,stagz,mdv,xmin,xmax,
153
     >                            ymin,ymax,dx,dy,nx,ny,nz,timecheck)
154
            endif
3 michaesp 155
      endif
156
 
21 michaesp 157
c     ----------------------------------------------------------------------
158
c     Do interpolation and get minimum, maximum
159
c     ----------------------------------------------------------------------
3 michaesp 160
 
161
C     If required do interpolation on pressure or theta level
162
      if (mode.eq.'P') then
163
        call vipo(field,varl,level,varf,nx,ny,nz,mdv,ipom)
21 michaesp 164
      elseif (mode.eq.'T') then
165
        call thipo(field,varl,level,varf,nx,ny,nz,mdv,ipom)
166
      elseif (mode.eq.'L') then
167
        do i=1,nx
168
           do j=1,ny
169
               varf(i,j) = field(i,j,nint(level))
170
           enddo
171
        enddo
172
      elseif (mode.eq.'X') then
173
         do i=1,nx
174
            do j=1,ny
175
                varf(i,j) = field(i,j,1)
176
            enddo
177
         enddo
3 michaesp 178
      endif
179
 
21 michaesp 180
C     Determine minimum and maximum
181
      do i=1,nx
182
        do j=1,ny
23 michaesp 183
 
21 michaesp 184
            if ((varf(i,j).ne.mdv).and.(varf(i,j).lt.min)) then
185
               min     = varf(i,j)
186
               minindx = i
187
               minindy = j
23 michaesp 188
            endif
189
            if ((varf(i,j).ne.mdv).and.(varf(i,j).gt.max)) then
21 michaesp 190
               max     = varf(i,j)
191
               maxindx = i
192
               maxindy = j
3 michaesp 193
          endif
194
        enddo
21 michaesp 195
      enddo
196
      lonmin = xmin + real(minindx-1) * dx
197
      latmin = ymin + real(minindy-1) * dy
198
      lonmax = xmin + real(maxindx-1) * dx
199
      latmax = ymin + real(maxindy-1) * dy
23 michaesp 200
 
21 michaesp 201
c     Rotate position to true longitude/latitude 
3 michaesp 202
      if ((pollon.ne.0.).or.(pollat.ne.90.)) then
203
        xphys=lmstolm(latmin,lonmin,pollat,pollon)
204
        yphys=phstoph(latmin,lonmin,pollat,pollon)
205
        lonmin=xphys
206
        latmin=yphys
207
        xphys=lmstolm(latmax,lonmax,pollat,pollon)
208
        yphys=phstoph(latmax,lonmax,pollat,pollon)
209
        lonmax=xphys
210
        latmax=yphys
211
      endif
212
 
21 michaesp 213
c     Write output
214
      write(*,101) min,max,lonmin,latmin,nint(level),
215
     >                     lonmax,latmax,nint(level)
216
  101 format(2f10.3,2f8.2,i6,2f8.2,i6)
3 michaesp 217
 
21 michaesp 218
c     Close netCDF file
219
      call input_close(cdfid)
3 michaesp 220
 
221
      end
21 michaesp 222
 
223
c     ----------------------------------------------------------------------
224
C     Interpolates the 3d variable var3d on the pressure surface defined
3 michaesp 225
C     by lev of the variable varl. The interpolated field is
226
C     returned as var.
227
C     ipom determines the way of vertical interpolation:
228
C       ipom=0 is for linear interpolation
229
C       ipom=1 is for logarithmic interpolation
21 michaesp 230
c     ----------------------------------------------------------------------
231
 
232
      subroutine vipo(var3d,varl,lev,var,nx,ny,nz,mdv,ipom)
233
 
3 michaesp 234
      integer   nx,ny,nz,ipom
235
      real      lev,mdv
236
      real      var3d(nx,ny,nz),varl(nx,ny,nz),var(nx,ny)
237
 
238
      integer   i,j,k
239
      real      kind
240
      real      int3dm
241
 
242
      do i=1,nx
243
      do j=1,ny
244
 
245
        do k=1,nz-1
246
          if ((varl(i,j,k).ge.lev).and.(varl(i,j,k+1).le.lev)) then
247
            kind=float(k)+(varl(i,j,k)-lev)/
248
     >                   (varl(i,j,k)-varl(i,j,k+1))
249
            goto 100
250
          endif
251
        enddo
252
 100    continue
253
 
254
        var(i,j)=int3dm(var3d,nx,ny,nz,float(i),float(j),kind,mdv)
255
 
256
      enddo
257
      enddo
258
 
259
      return
260
      end
21 michaesp 261
 
262
c     ----------------------------------------------------------------------
263
C     Interpolates the 3d variable var3d on the isentropic surface defined
264
C     by lev of the variable varl. The interpolated field is
265
C     returned as var.
266
C     ipom determines the way of vertical interpolation:
267
C       ipom=0 is for linear interpolation
268
C       ipom=1 is for logarithmic interpolation
269
c     ----------------------------------------------------------------------
3 michaesp 270
 
271
      subroutine thipo(var3d,th3d,lev,var,nx,ny,nz,mdv,mode)
272
 
273
      integer   nx,ny,nz,mode
274
      real      lev,mdv
275
      real      var3d(nx,ny,nz),th3d(nx,ny,nz),var(nx,ny)
276
 
277
      integer   i,j,k
278
      real      kind
279
      real      int3dm
280
 
281
      do i=1,nx
282
      do j=1,ny
283
 
284
        do k=1,nz-1
285
          if ((th3d(i,j,k).le.lev).and.(th3d(i,j,k+1).ge.lev)) then
286
            kind=float(k)+(th3d(i,j,k)-lev)/
287
     >                   (th3d(i,j,k)-th3d(i,j,k+1))
288
            goto 100
289
          endif
290
        enddo
291
 100    continue
292
 
293
        var(i,j)=int3dm(var3d,nx,ny,nz,float(i),float(j),kind,mdv)
294
 
295
      enddo
296
      enddo
297
 
298
      return
299
      end
21 michaesp 300
 
301
c     ----------------------------------------------------------------------
302
c     Check whether character appears within string
303
c     ----------------------------------------------------------------------
3 michaesp 304
 
305
      subroutine checkchar(string,char,flag)
306
 
307
      character*(*)     string
308
      character*(1)     char
21 michaesp 309
      integer           n,flag
3 michaesp 310
 
311
      flag=0
312
      do n=1,len(string)
313
        if (string(n:n).eq.char) then
314
          flag=n
315
          return
316
        endif
317
      enddo
318
      end
21 michaesp 319
 
320
c     ----------------------------------------------------------------------
321
c     3D interpolation
322
c     ----------------------------------------------------------------------
323
 
324
      real function int3d(ar,n1,n2,n3,rid,rjd,rkd)
3 michaesp 325
 
326
c     Purpose:
327
c        This subroutine interpolates a 3d-array to an arbitrary
328
c        location within the grid.
329
c     Arguments:
330
c        ar        real  input   surface pressure, define as ar(n1,n2,n3)
331
c        n1,n2,n3  int   input   dimensions of ar
332
c        ri,rj,rk  real  input   grid location to be interpolated to
333
c     History:
334
 
335
c     argument declarations
336
      integer   n1,n2,n3
337
      real      ar(n1,n2,n3), rid,rjd,rkd
338
 
339
c     local declarations
340
      integer   i,j,k,ip1,jp1,kp1,ih,jh,kh
341
      real      frac0i,frac0j,frac0k,frac1i,frac1j,frac1k,ri,rj,rk
342
 
343
c     do linear interpolation
344
      ri=amax1(1.,amin1(float(n1),rid))
345
      rj=amax1(1.,amin1(float(n2),rjd))
346
      rk=amax1(1.,amin1(float(n3),rkd))
347
      ih=nint(ri)
348
      jh=nint(rj)
349
      kh=nint(rk)
350
 
351
c     Check for interpolation in i
352
*     if (abs(float(ih)-ri).lt.1.e-3) then
353
*       i  =ih
354
*       ip1=ih
355
*     else
356
        i =min0(int(ri),n1-1)
357
        ip1=i+1
358
*     endif
359
 
360
c     Check for interpolation in j
361
      if (abs(float(jh)-rj).lt.1.e-3) then
362
        j  =jh
363
        jp1=jh
364
      else
365
        j =min0(int(rj),n2-1)
366
        jp1=j+1
367
      endif
368
 
369
c     Check for interpolation in k
370
*     if (abs(float(kh)-rk).lt.1.e-3) then
371
*       k  =kh
372
*       kp1=kh
373
*     else
374
        k =min0(int(rk),n3-1)
375
        kp1=k+1
376
*     endif
377
 
378
      if (k.eq.kp1) then
379
c       no interpolation in k
380
        if ((i.eq.ip1).and.(j.eq.jp1)) then
381
c          no interpolation at all
382
           int3d=ar(i,j,k)
383
c          print *,'int3d 00: ',rid,rjd,rkd,int3d
384
        else
385
c          horizontal interpolation only
386
           frac0i=ri-float(i)
387
           frac0j=rj-float(j)
388
           frac1i=1.-frac0i
389
           frac1j=1.-frac0j
390
           int3d = ar(i  ,j  ,k  ) * frac1i * frac1j
391
     &           + ar(i  ,jp1,k  ) * frac1i * frac0j
392
     &           + ar(ip1,j  ,k  ) * frac0i * frac1j
393
     &           + ar(ip1,jp1,k  ) * frac0i * frac0j
394
c          print *,'int3d 10: ',rid,rjd,rkd,int3d
395
        endif
396
      else 
397
        frac0k=rk-float(k)
398
        frac1k=1.-frac0k
399
        if ((i.eq.ip1).and.(j.eq.jp1)) then
400
c          vertical interpolation only
401
           int3d = ar(i  ,j  ,k  ) * frac1k
402
     &           + ar(i  ,j  ,kp1) * frac0k
403
c          print *,'int3d 01: ',rid,rjd,rkd,int3d
404
        else
405
c          full 3d interpolation
406
           frac0i=ri-float(i)
407
           frac0j=rj-float(j)
408
           frac1i=1.-frac0i
409
           frac1j=1.-frac0j
410
           int3d = ar(i  ,j  ,k  ) * frac1i * frac1j * frac1k
411
     &           + ar(i  ,jp1,k  ) * frac1i * frac0j * frac1k 
412
     &           + ar(ip1,j  ,k  ) * frac0i * frac1j * frac1k
413
     &           + ar(ip1,jp1,k  ) * frac0i * frac0j * frac1k
414
     &           + ar(i  ,j  ,kp1) * frac1i * frac1j * frac0k
415
     &           + ar(i  ,jp1,kp1) * frac1i * frac0j * frac0k 
416
     &           + ar(ip1,j  ,kp1) * frac0i * frac1j * frac0k
417
     &           + ar(ip1,jp1,kp1) * frac0i * frac0j * frac0k
418
c          print *,'int3d 11: ',rid,rjd,rkd,int3d
419
        endif
420
      endif
421
      end
21 michaesp 422
 
423
c     ----------------------------------------------------------------------
424
c     3D interpolation including missing data check
425
c     ----------------------------------------------------------------------
426
 
3 michaesp 427
      real function int3dm(ar,n1,n2,n3,rid,rjd,rkd,misdat)
21 michaesp 428
 
3 michaesp 429
c     Purpose:
430
c        This subroutine interpolates a 3d-array to an arbitrary
431
c        location within the grid. The interpolation includes the 
432
c        testing of the missing data flag 'misdat'. 
433
c     Arguments:
434
c        ar        real  input   surface pressure, define as ar(n1,n2,n3)
435
c        n1,n2,n3  int   input   dimensions of ar
436
c        ri,rj,rk  real  input   grid location to be interpolated to
437
c        misdat    real  input   missing data flag (on if misdat<>0)
438
c     Warning:
439
c        This routine has not yet been seriously tested
440
c     History:
441
 
442
c     argument declarations
443
      integer   n1,n2,n3
444
      real      ar(n1,n2,n3), rid,rjd,rkd, misdat
445
 
446
c     local declarations
447
      integer   i,j,k,ip1,jp1,kp1,ih,jh,kh
448
      real      frac0i,frac0j,frac0k,frac1i,frac1j,frac1k,ri,rj,rk,int3d
449
 
450
c     check if routine without missing data checking can be called instead
451
      if (misdat.eq.0.) then
452
        int3dm=int3d(ar,n1,n2,n3,rid,rjd,rkd)
453
        return
454
      endif
455
 
456
c     do linear interpolation
457
      ri=amax1(1.,amin1(float(n1),rid))
458
      rj=amax1(1.,amin1(float(n2),rjd))
459
      rk=amax1(1.,amin1(float(n3),rkd))
460
      ih=nint(ri)
461
      jh=nint(rj)
462
      kh=nint(rk)
463
 
464
c     Check for interpolation in i
465
*     if (abs(float(ih)-ri).lt.1.e-3) then
466
*       i  =ih
467
*       ip1=ih
468
*     else
469
        i =min0(int(ri),n1-1)
470
        ip1=i+1
471
*     endif
472
 
473
c     Check for interpolation in j
474
*     if (abs(float(jh)-rj).lt.1.e-3) then
475
*       j  =jh
476
*       jp1=jh
477
*     else
478
        j =min0(int(rj),n2-1)
479
        jp1=j+1
480
*     endif
481
 
482
c     Check for interpolation in k
483
*     if (abs(float(kh)-rk).lt.1.e-3) then
484
*       k  =kh
485
*       kp1=kh
486
*     else
487
        k =min0(int(rk),n3-1)
488
        kp1=k+1
489
*     endif
490
 
491
      if (k.eq.kp1) then
492
c       no interpolation in k
493
        if ((i.eq.ip1).and.(j.eq.jp1)) then
494
c          no interpolation at all
495
           if (misdat.eq.ar(i,j,k)) then
496
             int3dm=misdat
497
           else
498
             int3dm=ar(i,j,k)
499
           endif
500
c          print *,'int3dm 00: ',rid,rjd,rkd,int3dm
501
        else
502
c          horizontal interpolation only
503
           if ((misdat.eq.ar(i  ,j  ,k  )).or.
504
     &         (misdat.eq.ar(i  ,jp1,k  )).or.
505
     &         (misdat.eq.ar(ip1,j  ,k  )).or.
506
     &         (misdat.eq.ar(ip1,jp1,k  ))) then
507
             int3dm=misdat
508
           else
509
             frac0i=ri-float(i)
510
             frac0j=rj-float(j)
511
             frac1i=1.-frac0i
512
             frac1j=1.-frac0j
513
             int3dm = ar(i  ,j  ,k  ) * frac1i * frac1j
514
     &              + ar(i  ,jp1,k  ) * frac1i * frac0j
515
     &              + ar(ip1,j  ,k  ) * frac0i * frac1j
516
     &              + ar(ip1,jp1,k  ) * frac0i * frac0j
517
c            print *,'int3dm 10: ',rid,rjd,rkd,int3dm
518
           endif
519
        endif
520
      else 
521
        frac0k=rk-float(k)
522
        frac1k=1.-frac0k
523
        if ((i.eq.ip1).and.(j.eq.jp1)) then
524
c          vertical interpolation only
525
           if ((misdat.eq.ar(i  ,j  ,k  )).or.
526
     &         (misdat.eq.ar(i  ,j  ,kp1))) then
527
             int3dm=misdat
528
           else
529
             int3dm = ar(i  ,j  ,k  ) * frac1k
530
     &              + ar(i  ,j  ,kp1) * frac0k
531
c            print *,'int3dm 01: ',rid,rjd,rkd,int3dm
532
           endif
533
        else
534
c          full 3d interpolation
535
           if ((misdat.eq.ar(i  ,j  ,k  )).or.
536
     &         (misdat.eq.ar(i  ,jp1,k  )).or.
537
     &         (misdat.eq.ar(ip1,j  ,k  )).or.
538
     &         (misdat.eq.ar(ip1,jp1,k  )).or.
539
     &         (misdat.eq.ar(i  ,j  ,kp1)).or.
540
     &         (misdat.eq.ar(i  ,jp1,kp1)).or.
541
     &         (misdat.eq.ar(ip1,j  ,kp1)).or.
542
     &         (misdat.eq.ar(ip1,jp1,kp1))) then
543
             int3dm=misdat
544
           else
545
             frac0i=ri-float(i)
546
             frac0j=rj-float(j)
547
             frac1i=1.-frac0i
548
             frac1j=1.-frac0j
549
             int3dm = ar(i  ,j  ,k  ) * frac1i * frac1j * frac1k
550
     &              + ar(i  ,jp1,k  ) * frac1i * frac0j * frac1k
551
     &              + ar(ip1,j  ,k  ) * frac0i * frac1j * frac1k
552
     &              + ar(ip1,jp1,k  ) * frac0i * frac0j * frac1k
553
     &              + ar(i  ,j  ,kp1) * frac1i * frac1j * frac0k
554
     &              + ar(i  ,jp1,kp1) * frac1i * frac0j * frac0k
555
     &              + ar(ip1,j  ,kp1) * frac0i * frac1j * frac0k
556
     &              + ar(ip1,jp1,kp1) * frac0i * frac0j * frac0k
557
c            print *,'int3dm 11: ',rid,rjd,rkd,int3dm
558
           endif
559
        endif
560
      endif
561
      end
562
 
21 michaesp 563
c     ----------------------------------------------------------------------
564
c     Calculate potential temperature
565
c     ----------------------------------------------------------------------
3 michaesp 566
 
567
      subroutine pottemp(pt,t,sp,ie,je,ke,ak,bk)
568
 
569
c     argument declaration
570
      integer   ie,je,ke
571
      real      pt(ie,je,ke),t(ie,je,ke),sp(ie,je),
572
     >     	ak(ke),bk(ke)
573
 
574
c     variable declaration
575
      integer   i,j,k
576
      real      rdcp,tzero
577
      data      rdcp,tzero /0.286,273.15/
578
 
579
c     statement-functions for the computation of pressure
580
      real      pp,psrf
581
      integer   is
582
      pp(is)=ak(is)+bk(is)*psrf
583
 
584
c     computation of potential temperature
585
      do i=1,ie
586
      do j=1,je
587
        psrf=sp(i,j)
588
        do k=1,ke
589
c     distinction of temperature in K and deg C
590
          if (t(i,j,k).lt.100.) then
591
            pt(i,j,k)=(t(i,j,k)+tzero)*( (1000./pp(k))**rdcp )
592
          else
593
            pt(i,j,k)=t(i,j,k)*( (1000./pp(k))**rdcp )
594
          endif
595
        enddo
596
      enddo
597
      enddo
598
      end
21 michaesp 599
 
600
c     ----------------------------------------------------------------------
601
c     Calculate 3D pressure
602
c     ----------------------------------------------------------------------
3 michaesp 603
 
604
      subroutine pres(pr,sp,ie,je,ke,ak,bk)
21 michaesp 605
 
3 michaesp 606
c     argument declaration
607
      integer  ie,je,ke
608
      real,intent(OUT) :: pr(ie,je,ke)
609
      real,intent(IN)  :: sp(ie,je)
610
      real,intent(IN)  :: ak(ke),bk(ke)
611
 
612
c     variable declaration
613
      integer  i,j,k
614
 
615
c     computation pressure
616
      do i=1,ie
617
      do j=1,je
618
      do k=1,ke
619
        pr(i,j,k)=ak(k)+bk(k)*sp(i,j)
620
      enddo
621
      enddo
622
      enddo
623
      end
21 michaesp 624
 
625
c     ----------------------------------------------------------------------
626
c     Coordinate rotations from COSMO
627
c     ----------------------------------------------------------------------      
628
 
3 michaesp 629
      REAL FUNCTION PHTOPHS (PHI, LAM, POLPHI, POLLAM)
630
C
631
C%Z% Modul %M%, V%I% vom %G%, extrahiert am %H%
632
C
633
C**** PHTOPHS  -   FC:UMRECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE PHI
634
C****                 AUF EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
635
C****                 IM ROTIERTEN SYSTEM. DER NORDPOL DES SYSTEMS HAT
636
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
637
C**   AUFRUF   :   PHI = PHTOPHS (PHI, LAM, POLPHI, POLLAM)
638
C**   ENTRIES  :   KEINE
639
C**   ZWECK    :   UMRECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE PHI AUF
640
C**                EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
641
C**                ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
642
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
643
C**   VERSIONS-
644
C**   DATUM    :   03.05.90
645
C**
646
C**   EXTERNALS:   KEINE
647
C**   EINGABE-
648
C**   PARAMETER:   PHI    REAL BREITE DES PUNKTES IM GEOGR. SYSTEM
649
C**                LAM    REAL LAENGE DES PUNKTES IM GEOGR. SYSTEM
650
C**                POLPHI REAL GEOGR.BREITE DES N-POLS DES ROT. SYSTEMS
651
C**                POLLAM REAL GEOGR.LAENGE DES N-POLS DES ROT. SYSTEMS
652
C**   AUSGABE-
653
C**   PARAMETER:   ROTIERTE BREITE PHIS ALS WERT DER FUNKTION
654
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
655
C**
656
C**   COMMON-
657
C**   BLOECKE  :   KEINE
658
C**
659
C**   FEHLERBE-
660
C**   HANDLUNG :   KEINE
661
C**   VERFASSER:   G. DE MORSIER
662
 
663
      REAL        LAM,PHI,POLPHI,POLLAM
664
 
665
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
666
 
667
      ZSINPOL = SIN(ZPIR18*POLPHI)
668
      ZCOSPOL = COS(ZPIR18*POLPHI)
669
      ZLAMPOL = ZPIR18*POLLAM
670
      ZPHI    = ZPIR18*PHI
671
      ZLAM    = LAM
672
      IF(ZLAM.GT.180.0) ZLAM = ZLAM - 360.0
673
      ZLAM    = ZPIR18*ZLAM
674
      ZARG    = ZCOSPOL*COS(ZPHI)*COS(ZLAM-ZLAMPOL) + ZSINPOL*SIN(ZPHI)
675
 
676
      PHTOPHS = ZRPI18*ASIN(ZARG)
677
 
678
      RETURN
679
      END
21 michaesp 680
 
681
 
3 michaesp 682
      REAL FUNCTION PHSTOPH (PHIS, LAMS, POLPHI, POLLAM)
683
C
684
C%Z% Modul %M%, V%I% vom %G%, extrahiert am %H%
685
C
686
C**** PHSTOPH  -   FC:BERECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE FUER
687
C****                 EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
688
C****                 ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
689
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
690
C**   AUFRUF   :   PHI = PHSTOPH (PHIS, LAMS, POLPHI, POLLAM)
691
C**   ENTRIES  :   KEINE
692
C**   ZWECK    :   BERECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE FUER
693
C**                EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
694
C**                ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
695
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
696
C**   VERSIONS-
697
C**   DATUM    :   03.05.90
698
C**
699
C**   EXTERNALS:   KEINE
700
C**   EINGABE-
701
C**   PARAMETER:   PHIS     REAL   GEOGR. BREITE DES PUNKTES IM ROT.SYS.
702
C**                LAMS     REAL   GEOGR. LAENGE DES PUNKTES IM ROT.SYS.
703
C**                POLPHI   REAL   WAHRE GEOGR. BREITE DES NORDPOLS
704
C**                POLLAM   REAL   WAHRE GEOGR. LAENGE DES NORDPOLS
705
C**   AUSGABE-
706
C**   PARAMETER:   WAHRE GEOGRAPHISCHE BREITE ALS WERT DER FUNKTION
707
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
708
C**
709
C**   COMMON-
710
C**   BLOECKE  :   KEINE
711
C**
712
C**   FEHLERBE-
713
C**   HANDLUNG :   KEINE
714
C**   VERFASSER:   D.MAJEWSKI
715
 
716
      REAL        LAMS,PHIS,POLPHI,POLLAM
717
 
718
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
719
 
720
      SINPOL = SIN(ZPIR18*POLPHI)
721
      COSPOL = COS(ZPIR18*POLPHI)
722
      ZPHIS  = ZPIR18*PHIS
723
      ZLAMS  = LAMS
724
      IF(ZLAMS.GT.180.0) ZLAMS = ZLAMS - 360.0
725
      ZLAMS  = ZPIR18*ZLAMS
726
      ARG     = COSPOL*COS(ZPHIS)*COS(ZLAMS) + SINPOL*SIN(ZPHIS)
727
 
728
      PHSTOPH = ZRPI18*ASIN(ARG)
729
 
730
      RETURN
731
      END
732
      REAL FUNCTION LMTOLMS (PHI, LAM, POLPHI, POLLAM)
733
C
734
C%Z% Modul %M%, V%I% vom %G%, extrahiert am %H%
735
C
736
C**** LMTOLMS  -   FC:UMRECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE LAM
737
C****                 AUF EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
738
C****                 IM ROTIERTEN SYSTEM. DER NORDPOL DES SYSTEMS HAT
739
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
740
C**   AUFRUF   :   LAM = LMTOLMS (PHI, LAM, POLPHI, POLLAM)
741
C**   ENTRIES  :   KEINE
742
C**   ZWECK    :   UMRECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE LAM AUF
743
C**                EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
744
C**                ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
745
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
746
C**   VERSIONS-
747
C**   DATUM    :   03.05.90
748
C**
749
C**   EXTERNALS:   KEINE
750
C**   EINGABE-
751
C**   PARAMETER:   PHI    REAL BREITE DES PUNKTES IM GEOGR. SYSTEM
752
C**                LAM    REAL LAENGE DES PUNKTES IM GEOGR. SYSTEM
753
C**                POLPHI REAL GEOGR.BREITE DES N-POLS DES ROT. SYSTEMS
754
C**                POLLAM REAL GEOGR.LAENGE DES N-POLS DES ROT. SYSTEMS
755
C**   AUSGABE-
756
C**   PARAMETER:   WAHRE GEOGRAPHISCHE LAENGE ALS WERT DER FUNKTION
757
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
758
C**
759
C**   COMMON-
760
C**   BLOECKE  :   KEINE
761
C**
762
C**   FEHLERBE-
763
C**   HANDLUNG :   KEINE
764
C**   VERFASSER:   G. DE MORSIER
765
 
766
      REAL        LAM,PHI,POLPHI,POLLAM
767
 
768
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
769
 
770
      ZSINPOL = SIN(ZPIR18*POLPHI)
771
      ZCOSPOL = COS(ZPIR18*POLPHI)
772
      ZLAMPOL =     ZPIR18*POLLAM
773
      ZPHI    =     ZPIR18*PHI
774
      ZLAM    = LAM
775
      IF(ZLAM.GT.180.0) ZLAM = ZLAM - 360.0
776
      ZLAM    = ZPIR18*ZLAM
777
 
778
      ZARG1   = - SIN(ZLAM-ZLAMPOL)*COS(ZPHI)
779
      ZARG2   = - ZSINPOL*COS(ZPHI)*COS(ZLAM-ZLAMPOL)+ZCOSPOL*SIN(ZPHI)
780
      IF (ABS(ZARG2).LT.1.E-30) THEN
781
        IF (ABS(ZARG1).LT.1.E-30) THEN
782
          LMTOLMS =   0.0
783
        ELSEIF (ZARG1.GT.0.) THEN
784
              LMTOLMS =  90.0
785
            ELSE
786
              LMTOLMS = -90.0
787
            ENDIF
788
      ELSE
789
        LMTOLMS = ZRPI18*ATAN2(ZARG1,ZARG2)
790
      ENDIF
791
 
792
      RETURN
793
      END
21 michaesp 794
 
795
 
3 michaesp 796
      REAL FUNCTION LMSTOLM (PHIS, LAMS, POLPHI, POLLAM)
797
C
798
C%Z% Modul %M%, V%I% vom %G%, extrahiert am %H%
799
C
800
C**** LMSTOLM  -   FC:BERECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE FUER
801
C****                 EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
802
C****                 IM ROTIERTEN SYSTEM. DER NORDPOL DES SYSTEMS HAT
803
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
804
C**   AUFRUF   :   LAM = LMSTOLM (PHIS, LAMS, POLPHI, POLLAM)
805
C**   ENTRIES  :   KEINE
806
C**   ZWECK    :   BERECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE FUER
807
C**                EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
808
C**                IM ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
809
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
810
C**   VERSIONS-
811
C**   DATUM    :   03.05.90
812
C**
813
C**   EXTERNALS:   KEINE
814
C**   EINGABE-
815
C**   PARAMETER:   PHIS     REAL   GEOGR. BREITE DES PUNKTES IM ROT.SYS.
816
C**                LAMS     REAL   GEOGR. LAENGE DES PUNKTES IM ROT.SYS.
817
C**                POLPHI   REAL   WAHRE GEOGR. BREITE DES NORDPOLS
818
C**                POLLAM   REAL   WAHRE GEOGR. LAENGE DES NORDPOLS
819
C**   AUSGABE-
820
C**   PARAMETER:   WAHRE GEOGRAPHISCHE LAENGE ALS WERT DER FUNKTION
821
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
822
C**
823
C**   COMMON-
824
C**   BLOECKE  :   KEINE
825
C**
826
C**   FEHLERBE-
827
C**   HANDLUNG :   KEINE
828
C**   VERFASSER:   D.MAJEWSKI
829
 
830
      REAL        LAMS,PHIS,POLPHI,POLLAM
831
 
832
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
833
 
834
      ZSINPOL = SIN(ZPIR18*POLPHI)
835
      ZCOSPOL = COS(ZPIR18*POLPHI)
836
      ZLAMPOL = ZPIR18*POLLAM
837
      ZPHIS   = ZPIR18*PHIS
838
      ZLAMS   = LAMS
839
      IF(ZLAMS.GT.180.0) ZLAMS = ZLAMS - 360.0
840
      ZLAMS   = ZPIR18*ZLAMS
841
 
842
      ZARG1   = SIN(ZLAMPOL)*(- ZSINPOL*COS(ZLAMS)*COS(ZPHIS)  +
843
     1                          ZCOSPOL*           SIN(ZPHIS)) -
844
     2          COS(ZLAMPOL)*           SIN(ZLAMS)*COS(ZPHIS)
845
      ZARG2   = COS(ZLAMPOL)*(- ZSINPOL*COS(ZLAMS)*COS(ZPHIS)  +
846
     1                          ZCOSPOL*           SIN(ZPHIS)) +
847
     2          SIN(ZLAMPOL)*           SIN(ZLAMS)*COS(ZPHIS)
848
      IF (ABS(ZARG2).LT.1.E-30) THEN
849
        IF (ABS(ZARG1).LT.1.E-30) THEN
850
          LMSTOLM =   0.0
851
        ELSEIF (ZARG1.GT.0.) THEN
852
              LMSTOLAM =  90.0
853
            ELSE
854
              LMSTOLAM = -90.0
855
            ENDIF
856
      ELSE
857
        LMSTOLM = ZRPI18*ATAN2(ZARG1,ZARG2)
858
      ENDIF
859
 
860
      RETURN
861
      END