Subversion Repositories lagranto.ecmwf

Rev

Rev 13 | Go to most recent revision | 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,
109
     >               0.,pollon,pollat,rd,rd,nz,rd,rd,rd,stagz,timecheck)
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)
139
 
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
183
            if ((varf(i,j).ne.mdv).and.(varf(i,j).lt.min)) then
184
               min     = varf(i,j)
185
               minindx = i
186
               minindy = j
187
            elseif ((varf(i,j).ne.mdv).and.(varf(i,j).gt.max)) then
188
               max     = varf(i,j)
189
               maxindx = i
190
               maxindy = j
3 michaesp 191
          endif
192
        enddo
21 michaesp 193
      enddo
194
      lonmin = xmin + real(minindx-1) * dx
195
      latmin = ymin + real(minindy-1) * dy
196
      lonmax = xmin + real(maxindx-1) * dx
197
      latmax = ymin + real(maxindy-1) * dy
3 michaesp 198
 
21 michaesp 199
c     Rotate position to true longitude/latitude 
3 michaesp 200
      if ((pollon.ne.0.).or.(pollat.ne.90.)) then
201
        xphys=lmstolm(latmin,lonmin,pollat,pollon)
202
        yphys=phstoph(latmin,lonmin,pollat,pollon)
203
        lonmin=xphys
204
        latmin=yphys
205
        xphys=lmstolm(latmax,lonmax,pollat,pollon)
206
        yphys=phstoph(latmax,lonmax,pollat,pollon)
207
        lonmax=xphys
208
        latmax=yphys
209
      endif
210
 
21 michaesp 211
c     Write output
212
      write(*,101) min,max,lonmin,latmin,nint(level),
213
     >                     lonmax,latmax,nint(level)
214
  101 format(2f10.3,2f8.2,i6,2f8.2,i6)
3 michaesp 215
 
21 michaesp 216
c     Close netCDF file
217
      call input_close(cdfid)
3 michaesp 218
 
219
      end
21 michaesp 220
 
221
c     ----------------------------------------------------------------------
222
C     Interpolates the 3d variable var3d on the pressure surface defined
3 michaesp 223
C     by lev of the variable varl. The interpolated field is
224
C     returned as var.
225
C     ipom determines the way of vertical interpolation:
226
C       ipom=0 is for linear interpolation
227
C       ipom=1 is for logarithmic interpolation
21 michaesp 228
c     ----------------------------------------------------------------------
229
 
230
      subroutine vipo(var3d,varl,lev,var,nx,ny,nz,mdv,ipom)
231
 
3 michaesp 232
      integer   nx,ny,nz,ipom
233
      real      lev,mdv
234
      real      var3d(nx,ny,nz),varl(nx,ny,nz),var(nx,ny)
235
 
236
      integer   i,j,k
237
      real      kind
238
      real      int3dm
239
 
240
      do i=1,nx
241
      do j=1,ny
242
 
243
        do k=1,nz-1
244
          if ((varl(i,j,k).ge.lev).and.(varl(i,j,k+1).le.lev)) then
245
            kind=float(k)+(varl(i,j,k)-lev)/
246
     >                   (varl(i,j,k)-varl(i,j,k+1))
247
            goto 100
248
          endif
249
        enddo
250
 100    continue
251
 
252
        var(i,j)=int3dm(var3d,nx,ny,nz,float(i),float(j),kind,mdv)
253
 
254
      enddo
255
      enddo
256
 
257
      return
258
      end
21 michaesp 259
 
260
c     ----------------------------------------------------------------------
261
C     Interpolates the 3d variable var3d on the isentropic surface defined
262
C     by lev of the variable varl. The interpolated field is
263
C     returned as var.
264
C     ipom determines the way of vertical interpolation:
265
C       ipom=0 is for linear interpolation
266
C       ipom=1 is for logarithmic interpolation
267
c     ----------------------------------------------------------------------
3 michaesp 268
 
269
      subroutine thipo(var3d,th3d,lev,var,nx,ny,nz,mdv,mode)
270
 
271
      integer   nx,ny,nz,mode
272
      real      lev,mdv
273
      real      var3d(nx,ny,nz),th3d(nx,ny,nz),var(nx,ny)
274
 
275
      integer   i,j,k
276
      real      kind
277
      real      int3dm
278
 
279
      do i=1,nx
280
      do j=1,ny
281
 
282
        do k=1,nz-1
283
          if ((th3d(i,j,k).le.lev).and.(th3d(i,j,k+1).ge.lev)) then
284
            kind=float(k)+(th3d(i,j,k)-lev)/
285
     >                   (th3d(i,j,k)-th3d(i,j,k+1))
286
            goto 100
287
          endif
288
        enddo
289
 100    continue
290
 
291
        var(i,j)=int3dm(var3d,nx,ny,nz,float(i),float(j),kind,mdv)
292
 
293
      enddo
294
      enddo
295
 
296
      return
297
      end
21 michaesp 298
 
299
c     ----------------------------------------------------------------------
300
c     Check whether character appears within string
301
c     ----------------------------------------------------------------------
3 michaesp 302
 
303
      subroutine checkchar(string,char,flag)
304
 
305
      character*(*)     string
306
      character*(1)     char
21 michaesp 307
      integer           n,flag
3 michaesp 308
 
309
      flag=0
310
      do n=1,len(string)
311
        if (string(n:n).eq.char) then
312
          flag=n
313
          return
314
        endif
315
      enddo
316
      end
21 michaesp 317
 
318
c     ----------------------------------------------------------------------
319
c     3D interpolation
320
c     ----------------------------------------------------------------------
321
 
322
      real function int3d(ar,n1,n2,n3,rid,rjd,rkd)
3 michaesp 323
 
324
c     Purpose:
325
c        This subroutine interpolates a 3d-array to an arbitrary
326
c        location within the grid.
327
c     Arguments:
328
c        ar        real  input   surface pressure, define as ar(n1,n2,n3)
329
c        n1,n2,n3  int   input   dimensions of ar
330
c        ri,rj,rk  real  input   grid location to be interpolated to
331
c     History:
332
 
333
c     argument declarations
334
      integer   n1,n2,n3
335
      real      ar(n1,n2,n3), rid,rjd,rkd
336
 
337
c     local declarations
338
      integer   i,j,k,ip1,jp1,kp1,ih,jh,kh
339
      real      frac0i,frac0j,frac0k,frac1i,frac1j,frac1k,ri,rj,rk
340
 
341
c     do linear interpolation
342
      ri=amax1(1.,amin1(float(n1),rid))
343
      rj=amax1(1.,amin1(float(n2),rjd))
344
      rk=amax1(1.,amin1(float(n3),rkd))
345
      ih=nint(ri)
346
      jh=nint(rj)
347
      kh=nint(rk)
348
 
349
c     Check for interpolation in i
350
*     if (abs(float(ih)-ri).lt.1.e-3) then
351
*       i  =ih
352
*       ip1=ih
353
*     else
354
        i =min0(int(ri),n1-1)
355
        ip1=i+1
356
*     endif
357
 
358
c     Check for interpolation in j
359
      if (abs(float(jh)-rj).lt.1.e-3) then
360
        j  =jh
361
        jp1=jh
362
      else
363
        j =min0(int(rj),n2-1)
364
        jp1=j+1
365
      endif
366
 
367
c     Check for interpolation in k
368
*     if (abs(float(kh)-rk).lt.1.e-3) then
369
*       k  =kh
370
*       kp1=kh
371
*     else
372
        k =min0(int(rk),n3-1)
373
        kp1=k+1
374
*     endif
375
 
376
      if (k.eq.kp1) then
377
c       no interpolation in k
378
        if ((i.eq.ip1).and.(j.eq.jp1)) then
379
c          no interpolation at all
380
           int3d=ar(i,j,k)
381
c          print *,'int3d 00: ',rid,rjd,rkd,int3d
382
        else
383
c          horizontal interpolation only
384
           frac0i=ri-float(i)
385
           frac0j=rj-float(j)
386
           frac1i=1.-frac0i
387
           frac1j=1.-frac0j
388
           int3d = ar(i  ,j  ,k  ) * frac1i * frac1j
389
     &           + ar(i  ,jp1,k  ) * frac1i * frac0j
390
     &           + ar(ip1,j  ,k  ) * frac0i * frac1j
391
     &           + ar(ip1,jp1,k  ) * frac0i * frac0j
392
c          print *,'int3d 10: ',rid,rjd,rkd,int3d
393
        endif
394
      else 
395
        frac0k=rk-float(k)
396
        frac1k=1.-frac0k
397
        if ((i.eq.ip1).and.(j.eq.jp1)) then
398
c          vertical interpolation only
399
           int3d = ar(i  ,j  ,k  ) * frac1k
400
     &           + ar(i  ,j  ,kp1) * frac0k
401
c          print *,'int3d 01: ',rid,rjd,rkd,int3d
402
        else
403
c          full 3d interpolation
404
           frac0i=ri-float(i)
405
           frac0j=rj-float(j)
406
           frac1i=1.-frac0i
407
           frac1j=1.-frac0j
408
           int3d = ar(i  ,j  ,k  ) * frac1i * frac1j * frac1k
409
     &           + ar(i  ,jp1,k  ) * frac1i * frac0j * frac1k 
410
     &           + ar(ip1,j  ,k  ) * frac0i * frac1j * frac1k
411
     &           + ar(ip1,jp1,k  ) * frac0i * frac0j * frac1k
412
     &           + ar(i  ,j  ,kp1) * frac1i * frac1j * frac0k
413
     &           + ar(i  ,jp1,kp1) * frac1i * frac0j * frac0k 
414
     &           + ar(ip1,j  ,kp1) * frac0i * frac1j * frac0k
415
     &           + ar(ip1,jp1,kp1) * frac0i * frac0j * frac0k
416
c          print *,'int3d 11: ',rid,rjd,rkd,int3d
417
        endif
418
      endif
419
      end
21 michaesp 420
 
421
c     ----------------------------------------------------------------------
422
c     3D interpolation including missing data check
423
c     ----------------------------------------------------------------------
424
 
3 michaesp 425
      real function int3dm(ar,n1,n2,n3,rid,rjd,rkd,misdat)
21 michaesp 426
 
3 michaesp 427
c     Purpose:
428
c        This subroutine interpolates a 3d-array to an arbitrary
429
c        location within the grid. The interpolation includes the 
430
c        testing of the missing data flag 'misdat'. 
431
c     Arguments:
432
c        ar        real  input   surface pressure, define as ar(n1,n2,n3)
433
c        n1,n2,n3  int   input   dimensions of ar
434
c        ri,rj,rk  real  input   grid location to be interpolated to
435
c        misdat    real  input   missing data flag (on if misdat<>0)
436
c     Warning:
437
c        This routine has not yet been seriously tested
438
c     History:
439
 
440
c     argument declarations
441
      integer   n1,n2,n3
442
      real      ar(n1,n2,n3), rid,rjd,rkd, misdat
443
 
444
c     local declarations
445
      integer   i,j,k,ip1,jp1,kp1,ih,jh,kh
446
      real      frac0i,frac0j,frac0k,frac1i,frac1j,frac1k,ri,rj,rk,int3d
447
 
448
c     check if routine without missing data checking can be called instead
449
      if (misdat.eq.0.) then
450
        int3dm=int3d(ar,n1,n2,n3,rid,rjd,rkd)
451
        return
452
      endif
453
 
454
c     do linear interpolation
455
      ri=amax1(1.,amin1(float(n1),rid))
456
      rj=amax1(1.,amin1(float(n2),rjd))
457
      rk=amax1(1.,amin1(float(n3),rkd))
458
      ih=nint(ri)
459
      jh=nint(rj)
460
      kh=nint(rk)
461
 
462
c     Check for interpolation in i
463
*     if (abs(float(ih)-ri).lt.1.e-3) then
464
*       i  =ih
465
*       ip1=ih
466
*     else
467
        i =min0(int(ri),n1-1)
468
        ip1=i+1
469
*     endif
470
 
471
c     Check for interpolation in j
472
*     if (abs(float(jh)-rj).lt.1.e-3) then
473
*       j  =jh
474
*       jp1=jh
475
*     else
476
        j =min0(int(rj),n2-1)
477
        jp1=j+1
478
*     endif
479
 
480
c     Check for interpolation in k
481
*     if (abs(float(kh)-rk).lt.1.e-3) then
482
*       k  =kh
483
*       kp1=kh
484
*     else
485
        k =min0(int(rk),n3-1)
486
        kp1=k+1
487
*     endif
488
 
489
      if (k.eq.kp1) then
490
c       no interpolation in k
491
        if ((i.eq.ip1).and.(j.eq.jp1)) then
492
c          no interpolation at all
493
           if (misdat.eq.ar(i,j,k)) then
494
             int3dm=misdat
495
           else
496
             int3dm=ar(i,j,k)
497
           endif
498
c          print *,'int3dm 00: ',rid,rjd,rkd,int3dm
499
        else
500
c          horizontal interpolation only
501
           if ((misdat.eq.ar(i  ,j  ,k  )).or.
502
     &         (misdat.eq.ar(i  ,jp1,k  )).or.
503
     &         (misdat.eq.ar(ip1,j  ,k  )).or.
504
     &         (misdat.eq.ar(ip1,jp1,k  ))) then
505
             int3dm=misdat
506
           else
507
             frac0i=ri-float(i)
508
             frac0j=rj-float(j)
509
             frac1i=1.-frac0i
510
             frac1j=1.-frac0j
511
             int3dm = ar(i  ,j  ,k  ) * frac1i * frac1j
512
     &              + ar(i  ,jp1,k  ) * frac1i * frac0j
513
     &              + ar(ip1,j  ,k  ) * frac0i * frac1j
514
     &              + ar(ip1,jp1,k  ) * frac0i * frac0j
515
c            print *,'int3dm 10: ',rid,rjd,rkd,int3dm
516
           endif
517
        endif
518
      else 
519
        frac0k=rk-float(k)
520
        frac1k=1.-frac0k
521
        if ((i.eq.ip1).and.(j.eq.jp1)) then
522
c          vertical interpolation only
523
           if ((misdat.eq.ar(i  ,j  ,k  )).or.
524
     &         (misdat.eq.ar(i  ,j  ,kp1))) then
525
             int3dm=misdat
526
           else
527
             int3dm = ar(i  ,j  ,k  ) * frac1k
528
     &              + ar(i  ,j  ,kp1) * frac0k
529
c            print *,'int3dm 01: ',rid,rjd,rkd,int3dm
530
           endif
531
        else
532
c          full 3d interpolation
533
           if ((misdat.eq.ar(i  ,j  ,k  )).or.
534
     &         (misdat.eq.ar(i  ,jp1,k  )).or.
535
     &         (misdat.eq.ar(ip1,j  ,k  )).or.
536
     &         (misdat.eq.ar(ip1,jp1,k  )).or.
537
     &         (misdat.eq.ar(i  ,j  ,kp1)).or.
538
     &         (misdat.eq.ar(i  ,jp1,kp1)).or.
539
     &         (misdat.eq.ar(ip1,j  ,kp1)).or.
540
     &         (misdat.eq.ar(ip1,jp1,kp1))) then
541
             int3dm=misdat
542
           else
543
             frac0i=ri-float(i)
544
             frac0j=rj-float(j)
545
             frac1i=1.-frac0i
546
             frac1j=1.-frac0j
547
             int3dm = ar(i  ,j  ,k  ) * frac1i * frac1j * frac1k
548
     &              + ar(i  ,jp1,k  ) * frac1i * frac0j * frac1k
549
     &              + ar(ip1,j  ,k  ) * frac0i * frac1j * frac1k
550
     &              + ar(ip1,jp1,k  ) * frac0i * frac0j * frac1k
551
     &              + ar(i  ,j  ,kp1) * frac1i * frac1j * frac0k
552
     &              + ar(i  ,jp1,kp1) * frac1i * frac0j * frac0k
553
     &              + ar(ip1,j  ,kp1) * frac0i * frac1j * frac0k
554
     &              + ar(ip1,jp1,kp1) * frac0i * frac0j * frac0k
555
c            print *,'int3dm 11: ',rid,rjd,rkd,int3dm
556
           endif
557
        endif
558
      endif
559
      end
560
 
21 michaesp 561
c     ----------------------------------------------------------------------
562
c     Calculate potential temperature
563
c     ----------------------------------------------------------------------
3 michaesp 564
 
565
      subroutine pottemp(pt,t,sp,ie,je,ke,ak,bk)
566
 
567
c     argument declaration
568
      integer   ie,je,ke
569
      real      pt(ie,je,ke),t(ie,je,ke),sp(ie,je),
570
     >     	ak(ke),bk(ke)
571
 
572
c     variable declaration
573
      integer   i,j,k
574
      real      rdcp,tzero
575
      data      rdcp,tzero /0.286,273.15/
576
 
577
c     statement-functions for the computation of pressure
578
      real      pp,psrf
579
      integer   is
580
      pp(is)=ak(is)+bk(is)*psrf
581
 
582
c     computation of potential temperature
583
      do i=1,ie
584
      do j=1,je
585
        psrf=sp(i,j)
586
        do k=1,ke
587
c     distinction of temperature in K and deg C
588
          if (t(i,j,k).lt.100.) then
589
            pt(i,j,k)=(t(i,j,k)+tzero)*( (1000./pp(k))**rdcp )
590
          else
591
            pt(i,j,k)=t(i,j,k)*( (1000./pp(k))**rdcp )
592
          endif
593
        enddo
594
      enddo
595
      enddo
596
      end
21 michaesp 597
 
598
c     ----------------------------------------------------------------------
599
c     Calculate 3D pressure
600
c     ----------------------------------------------------------------------
3 michaesp 601
 
602
      subroutine pres(pr,sp,ie,je,ke,ak,bk)
21 michaesp 603
 
3 michaesp 604
c     argument declaration
605
      integer  ie,je,ke
606
      real,intent(OUT) :: pr(ie,je,ke)
607
      real,intent(IN)  :: sp(ie,je)
608
      real,intent(IN)  :: ak(ke),bk(ke)
609
 
610
c     variable declaration
611
      integer  i,j,k
612
 
613
c     computation pressure
614
      do i=1,ie
615
      do j=1,je
616
      do k=1,ke
617
        pr(i,j,k)=ak(k)+bk(k)*sp(i,j)
618
      enddo
619
      enddo
620
      enddo
621
      end
21 michaesp 622
 
623
c     ----------------------------------------------------------------------
624
c     Coordinate rotations from COSMO
625
c     ----------------------------------------------------------------------      
626
 
3 michaesp 627
      REAL FUNCTION PHTOPHS (PHI, LAM, POLPHI, POLLAM)
628
C
629
C%Z% Modul %M%, V%I% vom %G%, extrahiert am %H%
630
C
631
C**** PHTOPHS  -   FC:UMRECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE PHI
632
C****                 AUF EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
633
C****                 IM ROTIERTEN SYSTEM. DER NORDPOL DES SYSTEMS HAT
634
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
635
C**   AUFRUF   :   PHI = PHTOPHS (PHI, LAM, POLPHI, POLLAM)
636
C**   ENTRIES  :   KEINE
637
C**   ZWECK    :   UMRECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE PHI AUF
638
C**                EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
639
C**                ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
640
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
641
C**   VERSIONS-
642
C**   DATUM    :   03.05.90
643
C**
644
C**   EXTERNALS:   KEINE
645
C**   EINGABE-
646
C**   PARAMETER:   PHI    REAL BREITE DES PUNKTES IM GEOGR. SYSTEM
647
C**                LAM    REAL LAENGE DES PUNKTES IM GEOGR. SYSTEM
648
C**                POLPHI REAL GEOGR.BREITE DES N-POLS DES ROT. SYSTEMS
649
C**                POLLAM REAL GEOGR.LAENGE DES N-POLS DES ROT. SYSTEMS
650
C**   AUSGABE-
651
C**   PARAMETER:   ROTIERTE BREITE PHIS ALS WERT DER FUNKTION
652
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
653
C**
654
C**   COMMON-
655
C**   BLOECKE  :   KEINE
656
C**
657
C**   FEHLERBE-
658
C**   HANDLUNG :   KEINE
659
C**   VERFASSER:   G. DE MORSIER
660
 
661
      REAL        LAM,PHI,POLPHI,POLLAM
662
 
663
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
664
 
665
      ZSINPOL = SIN(ZPIR18*POLPHI)
666
      ZCOSPOL = COS(ZPIR18*POLPHI)
667
      ZLAMPOL = ZPIR18*POLLAM
668
      ZPHI    = ZPIR18*PHI
669
      ZLAM    = LAM
670
      IF(ZLAM.GT.180.0) ZLAM = ZLAM - 360.0
671
      ZLAM    = ZPIR18*ZLAM
672
      ZARG    = ZCOSPOL*COS(ZPHI)*COS(ZLAM-ZLAMPOL) + ZSINPOL*SIN(ZPHI)
673
 
674
      PHTOPHS = ZRPI18*ASIN(ZARG)
675
 
676
      RETURN
677
      END
21 michaesp 678
 
679
 
3 michaesp 680
      REAL FUNCTION PHSTOPH (PHIS, LAMS, POLPHI, POLLAM)
681
C
682
C%Z% Modul %M%, V%I% vom %G%, extrahiert am %H%
683
C
684
C**** PHSTOPH  -   FC:BERECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE FUER
685
C****                 EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
686
C****                 ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
687
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
688
C**   AUFRUF   :   PHI = PHSTOPH (PHIS, LAMS, POLPHI, POLLAM)
689
C**   ENTRIES  :   KEINE
690
C**   ZWECK    :   BERECHNUNG DER WAHREN GEOGRAPHISCHEN BREITE FUER
691
C**                EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
692
C**                ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
693
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
694
C**   VERSIONS-
695
C**   DATUM    :   03.05.90
696
C**
697
C**   EXTERNALS:   KEINE
698
C**   EINGABE-
699
C**   PARAMETER:   PHIS     REAL   GEOGR. BREITE DES PUNKTES IM ROT.SYS.
700
C**                LAMS     REAL   GEOGR. LAENGE DES PUNKTES IM ROT.SYS.
701
C**                POLPHI   REAL   WAHRE GEOGR. BREITE DES NORDPOLS
702
C**                POLLAM   REAL   WAHRE GEOGR. LAENGE DES NORDPOLS
703
C**   AUSGABE-
704
C**   PARAMETER:   WAHRE GEOGRAPHISCHE BREITE ALS WERT DER FUNKTION
705
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
706
C**
707
C**   COMMON-
708
C**   BLOECKE  :   KEINE
709
C**
710
C**   FEHLERBE-
711
C**   HANDLUNG :   KEINE
712
C**   VERFASSER:   D.MAJEWSKI
713
 
714
      REAL        LAMS,PHIS,POLPHI,POLLAM
715
 
716
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
717
 
718
      SINPOL = SIN(ZPIR18*POLPHI)
719
      COSPOL = COS(ZPIR18*POLPHI)
720
      ZPHIS  = ZPIR18*PHIS
721
      ZLAMS  = LAMS
722
      IF(ZLAMS.GT.180.0) ZLAMS = ZLAMS - 360.0
723
      ZLAMS  = ZPIR18*ZLAMS
724
      ARG     = COSPOL*COS(ZPHIS)*COS(ZLAMS) + SINPOL*SIN(ZPHIS)
725
 
726
      PHSTOPH = ZRPI18*ASIN(ARG)
727
 
728
      RETURN
729
      END
730
      REAL FUNCTION LMTOLMS (PHI, LAM, POLPHI, POLLAM)
731
C
732
C%Z% Modul %M%, V%I% vom %G%, extrahiert am %H%
733
C
734
C**** LMTOLMS  -   FC:UMRECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE LAM
735
C****                 AUF EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
736
C****                 IM ROTIERTEN SYSTEM. DER NORDPOL DES SYSTEMS HAT
737
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
738
C**   AUFRUF   :   LAM = LMTOLMS (PHI, LAM, POLPHI, POLLAM)
739
C**   ENTRIES  :   KEINE
740
C**   ZWECK    :   UMRECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE LAM AUF
741
C**                EINEM PUNKT MIT DEN KOORDINATEN (PHIS, LAMS) IM
742
C**                ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
743
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
744
C**   VERSIONS-
745
C**   DATUM    :   03.05.90
746
C**
747
C**   EXTERNALS:   KEINE
748
C**   EINGABE-
749
C**   PARAMETER:   PHI    REAL BREITE DES PUNKTES IM GEOGR. SYSTEM
750
C**                LAM    REAL LAENGE DES PUNKTES IM GEOGR. SYSTEM
751
C**                POLPHI REAL GEOGR.BREITE DES N-POLS DES ROT. SYSTEMS
752
C**                POLLAM REAL GEOGR.LAENGE DES N-POLS DES ROT. SYSTEMS
753
C**   AUSGABE-
754
C**   PARAMETER:   WAHRE GEOGRAPHISCHE LAENGE ALS WERT DER FUNKTION
755
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
756
C**
757
C**   COMMON-
758
C**   BLOECKE  :   KEINE
759
C**
760
C**   FEHLERBE-
761
C**   HANDLUNG :   KEINE
762
C**   VERFASSER:   G. DE MORSIER
763
 
764
      REAL        LAM,PHI,POLPHI,POLLAM
765
 
766
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
767
 
768
      ZSINPOL = SIN(ZPIR18*POLPHI)
769
      ZCOSPOL = COS(ZPIR18*POLPHI)
770
      ZLAMPOL =     ZPIR18*POLLAM
771
      ZPHI    =     ZPIR18*PHI
772
      ZLAM    = LAM
773
      IF(ZLAM.GT.180.0) ZLAM = ZLAM - 360.0
774
      ZLAM    = ZPIR18*ZLAM
775
 
776
      ZARG1   = - SIN(ZLAM-ZLAMPOL)*COS(ZPHI)
777
      ZARG2   = - ZSINPOL*COS(ZPHI)*COS(ZLAM-ZLAMPOL)+ZCOSPOL*SIN(ZPHI)
778
      IF (ABS(ZARG2).LT.1.E-30) THEN
779
        IF (ABS(ZARG1).LT.1.E-30) THEN
780
          LMTOLMS =   0.0
781
        ELSEIF (ZARG1.GT.0.) THEN
782
              LMTOLMS =  90.0
783
            ELSE
784
              LMTOLMS = -90.0
785
            ENDIF
786
      ELSE
787
        LMTOLMS = ZRPI18*ATAN2(ZARG1,ZARG2)
788
      ENDIF
789
 
790
      RETURN
791
      END
21 michaesp 792
 
793
 
3 michaesp 794
      REAL FUNCTION LMSTOLM (PHIS, LAMS, POLPHI, POLLAM)
795
C
796
C%Z% Modul %M%, V%I% vom %G%, extrahiert am %H%
797
C
798
C**** LMSTOLM  -   FC:BERECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE FUER
799
C****                 EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
800
C****                 IM ROTIERTEN SYSTEM. DER NORDPOL DES SYSTEMS HAT
801
C****                 DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
802
C**   AUFRUF   :   LAM = LMSTOLM (PHIS, LAMS, POLPHI, POLLAM)
803
C**   ENTRIES  :   KEINE
804
C**   ZWECK    :   BERECHNUNG DER WAHREN GEOGRAPHISCHEN LAENGE FUER
805
C**                EINEN PUNKT MIT DEN KOORDINATEN (PHIS, LAMS)
806
C**                IM ROTIERTEN SYSTEM. DER NORDPOL DIESES SYSTEMS HAT
807
C**                DIE WAHREN KOORDINATEN (POLPHI, POLLAM)
808
C**   VERSIONS-
809
C**   DATUM    :   03.05.90
810
C**
811
C**   EXTERNALS:   KEINE
812
C**   EINGABE-
813
C**   PARAMETER:   PHIS     REAL   GEOGR. BREITE DES PUNKTES IM ROT.SYS.
814
C**                LAMS     REAL   GEOGR. LAENGE DES PUNKTES IM ROT.SYS.
815
C**                POLPHI   REAL   WAHRE GEOGR. BREITE DES NORDPOLS
816
C**                POLLAM   REAL   WAHRE GEOGR. LAENGE DES NORDPOLS
817
C**   AUSGABE-
818
C**   PARAMETER:   WAHRE GEOGRAPHISCHE LAENGE ALS WERT DER FUNKTION
819
C**                ALLE WINKEL IN GRAD (NORDEN>0, OSTEN>0)
820
C**
821
C**   COMMON-
822
C**   BLOECKE  :   KEINE
823
C**
824
C**   FEHLERBE-
825
C**   HANDLUNG :   KEINE
826
C**   VERFASSER:   D.MAJEWSKI
827
 
828
      REAL        LAMS,PHIS,POLPHI,POLLAM
829
 
830
      DATA        ZRPI18 , ZPIR18  / 57.2957795 , 0.0174532925 /
831
 
832
      ZSINPOL = SIN(ZPIR18*POLPHI)
833
      ZCOSPOL = COS(ZPIR18*POLPHI)
834
      ZLAMPOL = ZPIR18*POLLAM
835
      ZPHIS   = ZPIR18*PHIS
836
      ZLAMS   = LAMS
837
      IF(ZLAMS.GT.180.0) ZLAMS = ZLAMS - 360.0
838
      ZLAMS   = ZPIR18*ZLAMS
839
 
840
      ZARG1   = SIN(ZLAMPOL)*(- ZSINPOL*COS(ZLAMS)*COS(ZPHIS)  +
841
     1                          ZCOSPOL*           SIN(ZPHIS)) -
842
     2          COS(ZLAMPOL)*           SIN(ZLAMS)*COS(ZPHIS)
843
      ZARG2   = COS(ZLAMPOL)*(- ZSINPOL*COS(ZLAMS)*COS(ZPHIS)  +
844
     1                          ZCOSPOL*           SIN(ZPHIS)) +
845
     2          SIN(ZLAMPOL)*           SIN(ZLAMS)*COS(ZPHIS)
846
      IF (ABS(ZARG2).LT.1.E-30) THEN
847
        IF (ABS(ZARG1).LT.1.E-30) THEN
848
          LMSTOLM =   0.0
849
        ELSEIF (ZARG1.GT.0.) THEN
850
              LMSTOLAM =  90.0
851
            ELSE
852
              LMSTOLAM = -90.0
853
            ENDIF
854
      ELSE
855
        LMSTOLM = ZRPI18*ATAN2(ZARG1,ZARG2)
856
      ENDIF
857
 
858
      RETURN
859
      END