source: flexpart.git/src/interpol_wind_short_nests.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 4.4 KB
Line 
1subroutine interpol_wind_short_nests(itime,xt,yt,zt)
2  !                                       i   i  i  i
3  !*****************************************************************************
4  !                                                                            *
5  !  This subroutine interpolates the wind data to current trajectory position.*
6  !                                                                            *
7  !    Author: A. Stohl                                                        *
8  !                                                                            *
9  !    16 December 1997                                                        *
10  !                                                                            *
11  !  Revision March 2005 by AST : all output variables in common block         *
12  !                                                                            *
13  !*****************************************************************************
14  !                                                                            *
15  ! Variables:                                                                 *
16  ! u,v,w              wind components                                         *
17  ! itime [s]          current temporal position                               *
18  ! memtime(3) [s]     times of the wind fields in memory                      *
19  ! xt,yt,zt           coordinates position for which wind data shall be       *
20  !                    calculated                                              *
21  !                                                                            *
22  ! Constants:                                                                 *
23  !                                                                            *
24  !*****************************************************************************
25
26  use par_mod
27  use com_mod
28  use interpol_mod
29
30  implicit none
31
32  integer :: itime
33  real :: xt,yt,zt
34
35  ! Auxiliary variables needed for interpolation
36  real :: dz1,dz2,dz
37  real :: u1(2),v1(2),w1(2),uh(2),vh(2),wh(2)
38  integer :: i,m,n,indexh,indzh
39
40
41  !********************************************
42  ! Multilinear interpolation in time and space
43  !********************************************
44
45  ddx=xt-real(ix)
46  ddy=yt-real(jy)
47  rddx=1.-ddx
48  rddy=1.-ddy
49  p1=rddx*rddy
50  p2=ddx*rddy
51  p3=rddx*ddy
52  p4=ddx*ddy
53
54  ! Calculate variables for time interpolation
55  !*******************************************
56
57  dt1=real(itime-memtime(1))
58  dt2=real(memtime(2)-itime)
59  dtt=1./(dt1+dt2)
60
61  ! Determine the level below the current position for u,v
62  !*******************************************************
63
64  do i=2,nz
65    if (height(i).gt.zt) then
66      indz=i-1
67      goto 6
68    endif
69  end do
706   continue
71
72
73  ! Vertical distance to the level below and above current position
74  !****************************************************************
75
76  dz=1./(height(indz+1)-height(indz))
77  dz1=(zt-height(indz))*dz
78  dz2=(height(indz+1)-zt)*dz
79
80
81  !**********************************************************************
82  ! 1.) Bilinear horizontal interpolation
83  ! This has to be done separately for 6 fields (Temporal(2)*Vertical(3))
84  !**********************************************************************
85
86  ! Loop over 2 time steps and 2 levels
87  !************************************
88
89  do m=1,2
90    indexh=memind(m)
91    do n=1,2
92      indzh=indz+n-1
93
94      u1(n)=p1*uun(ix ,jy ,indzh,indexh,ngrid) &
95           +p2*uun(ixp,jy ,indzh,indexh,ngrid) &
96           +p3*uun(ix ,jyp,indzh,indexh,ngrid) &
97           +p4*uun(ixp,jyp,indzh,indexh,ngrid)
98      v1(n)=p1*vvn(ix ,jy ,indzh,indexh,ngrid) &
99           +p2*vvn(ixp,jy ,indzh,indexh,ngrid) &
100           +p3*vvn(ix ,jyp,indzh,indexh,ngrid) &
101           +p4*vvn(ixp,jyp,indzh,indexh,ngrid)
102      w1(n)=p1*wwn(ix ,jy ,indzh,indexh,ngrid) &
103           +p2*wwn(ixp,jy ,indzh,indexh,ngrid) &
104           +p3*wwn(ix ,jyp,indzh,indexh,ngrid) &
105           +p4*wwn(ixp,jyp,indzh,indexh,ngrid)
106
107    end do
108
109
110  !**********************************
111  ! 2.) Linear vertical interpolation
112  !**********************************
113
114    uh(m)=dz2*u1(1)+dz1*u1(2)
115    vh(m)=dz2*v1(1)+dz1*v1(2)
116    wh(m)=dz2*w1(1)+dz1*w1(2)
117  end do
118
119
120  !************************************
121  ! 3.) Temporal interpolation (linear)
122  !************************************
123
124  u=(uh(1)*dt2+uh(2)*dt1)*dtt
125  v=(vh(1)*dt2+vh(2)*dt1)*dtt
126  w=(wh(1)*dt2+wh(2)*dt1)*dtt
127
128end subroutine interpol_wind_short_nests
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG