File: etatrafo.f

package info (click to toggle)
flextra 5.0-18
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 928 kB
  • sloc: fortran: 7,018; makefile: 61; sh: 17
file content (90 lines) | stat: -rw-r--r-- 4,539 bytes parent folder | download | duplicates (7)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
      subroutine etatrafo(xt,yt,zt,itime1,itime2,itime,indexf,ngrid,
     +psint)
***********************************************************************
*                                                                     *
*             TRAJECTORY MODEL SUBROUTINE ETATRAFO                    *
*                                                                     *
***********************************************************************
*                                                                     *
*             AUTHOR:      G. WOTAWA                                  *
*             DATE:        1994-04-06                                 *
*             LAST UPDATE: ----------                                 *
*                                                                     *
***********************************************************************
*                                                                     *
* DESCRIPTION: This subroutine transforms the vertical coordinate     *
*              from z-coordinate system [m] to eta-coordinate         *
*              system [ECMWF]                                         *
*              Remark: Just call after initialization of a new        *
*              trajectory (first call of getwind)                     *
*                                                                     *
***********************************************************************
*                                                                     *
* INPUT:                                                              *
*                                                                     *
* xt                         x-coordinate of point [grid units]       *
* yt                         y-coordinate of point [grid units]       *
* zt                         z-coordinate of point [m]                *
* itime1                     time [s] of first windfield              *
* itime2                     time [s] of second windfield             *
* itime                      time [s] of calculation                  *
* indexf                     time index of field xx                   *
* psint                      surface pressure at point (xt,yt) [Pa]   *
*                                                                     *
***********************************************************************
*                                                                     *
* OUTPUT:                                                             *
*                                                                     *
* zt                         z-coordinate of point [eta ECMWF]        *
*                                                                     *
***********************************************************************
*
      include 'includepar'
      include 'includecom'

      integer itime1,itime2,itime,indexf,k,ngrid
      real xt,yt,zt,psint,fract,pp1,pp2,tv
      real tlev(nuvzmax),zzlev(nwzmax)
*
* calculate interpolated vertical temperature profile on model levels
*
      do 10 k=1,nuvz
       if (ngrid.gt.0) then
         call levlininterpoln(ttn,maxnests,nxmaxn,nymaxn,nuvzmax,ngrid,
     +   nxn,nyn,memind,xt,yt,k,itime1,itime2,itime,indexf,tlev(k))
       else
         call levlininterpol(tt,nxmax,nymax,nuvzmax,nx,ny,memind,
     +   xt,yt,k,itime1,itime2,itime,indexf,tlev(k))
       endif
10     continue
*
* calculate geometric height on model levels
*
      zzlev(1)=0.
      if (zt.lt.zzlev(1)) then
        write(*,*) ' TRAJECTORY MODEL: NOTICE - STARTING POINT OUT'//
     &             ' OF MODEL DOMAIN'
        write(*,*) ' (VERTICAL) HAS BEEN DETECTED --> IT IS SET TO'//
     &             ' THE BOTTOM OF THE MODEL ...'
        zt=uvheight(1)
        return
      endif
      do 15 k=2,nwz
        pp1=akm(k-1)+bkm(k-1)*psint
        pp2=akm(k)+bkm(k)*psint
        tv=tlev(k-1)      !! NO HUMIDITY INFORMATION AVAILABLE IN FLEXTRA
        zzlev(k)=zzlev(k-1)+r_air/ga*log(pp1/pp2)*tv
        if(zt.lt.zzlev(k)) goto 20
15    continue
      write(*,*) ' TRAJECTORY MODEL: NOTICE - STARTING POINT OUT'//
     &           ' OF MODEL DOMAIN'
      write(*,*) ' (VERTICAL) HAS BEEN DETECTED --> IT IS SET TO'//
     &           ' THE TOP OF THE MODEL ...'
      zt=uvheight(nuvz)
      return
20    fract=(zt-zzlev(k-1))/(zzlev(k)-zzlev(k-1))
      zt=wheight(k-1)*(1.-fract)+wheight(k)*fract
      if(zt.lt.uvheight(1)) zt=uvheight(1)
      if(zt.gt.uvheight(nuvz)) zt=uvheight(nuvz)
      return
      end