source: flexpart.git/src/interpol_wind.f90

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

add SPDX-License-Identifier to all .f90 files

  • Property mode set to 100644
File size: 7.8 KB
Line 
1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
3
4subroutine interpol_wind(itime,xt,yt,zt)
5  !                           i   i  i  i
6  !*****************************************************************************
7  !                                                                            *
8  !  This subroutine interpolates the wind data to current trajectory position.*
9  !                                                                            *
10  !    Author: A. Stohl                                                        *
11  !                                                                            *
12  !    16 December 1997                                                        *
13  !                                                                            *
14  !  Revision March 2005 by AST : all output variables in common block cal-    *
15  !                               culation of standard deviation done in this  *
16  !                               routine rather than subroutine call in order *
17  !                               to save computation time                     *
18  !                                                                            *
19  !*****************************************************************************
20  !                                                                            *
21  ! Variables:                                                                 *
22  ! u,v,w              wind components                                         *
23  ! itime [s]          current temporal position                               *
24  ! memtime(3) [s]     times of the wind fields in memory                      *
25  ! xt,yt,zt           coordinates position for which wind data shall be       *
26  !                    calculated                                              *
27  !                                                                            *
28  ! Constants:                                                                 *
29  !                                                                            *
30  !*****************************************************************************
31
32  use par_mod
33  use com_mod
34  use interpol_mod
35
36  implicit none
37
38  integer :: itime
39  real :: xt,yt,zt
40
41  ! Auxiliary variables needed for interpolation
42  real :: dz1,dz2,dz
43  real :: u1(2),v1(2),w1(2),uh(2),vh(2),wh(2)
44  real :: usl,vsl,wsl,usq,vsq,wsq,xaux
45  integer :: i,m,n,indexh,indzh
46  real,parameter :: eps=1.0e-30
47
48
49  !********************************************
50  ! Multilinear interpolation in time and space
51  !********************************************
52
53  ! Determine the lower left corner and its distance to the current position
54  !*************************************************************************
55
56  ddx=xt-real(ix)
57  ddy=yt-real(jy)
58  rddx=1.-ddx
59  rddy=1.-ddy
60  p1=rddx*rddy
61  p2=ddx*rddy
62  p3=rddx*ddy
63  p4=ddx*ddy
64
65  ! Calculate variables for time interpolation
66  !*******************************************
67
68  dt1=real(itime-memtime(1))
69  dt2=real(memtime(2)-itime)
70  dtt=1./(dt1+dt2)
71
72  ! Determine the level below the current position for u,v
73  !*******************************************************
74
75  do i=2,nz
76    if (height(i).gt.zt) then
77      indz=i-1
78      goto 6
79    endif
80  end do
816   continue
82
83
84  ! Vertical distance to the level below and above current position
85  !****************************************************************
86
87  dz=1./(height(indz+1)-height(indz))
88  dz1=(zt-height(indz))*dz
89  dz2=(height(indz+1)-zt)*dz
90
91  !**********************************************************************
92  ! 1.) Bilinear horizontal interpolation
93  ! This has to be done separately for 6 fields (Temporal(2)*Vertical(3))
94  !**********************************************************************
95
96  ! Loop over 2 time steps and 2 levels
97  !************************************
98
99  usl=0.
100  vsl=0.
101  wsl=0.
102  usq=0.
103  vsq=0.
104  wsq=0.
105  do m=1,2
106    indexh=memind(m)
107    do n=1,2
108      indzh=indz+n-1
109
110      if (ngrid.lt.0) then
111        u1(n)=p1*uupol(ix ,jy ,indzh,indexh) &
112             +p2*uupol(ixp,jy ,indzh,indexh) &
113             +p3*uupol(ix ,jyp,indzh,indexh) &
114             +p4*uupol(ixp,jyp,indzh,indexh)
115        v1(n)=p1*vvpol(ix ,jy ,indzh,indexh) &
116             +p2*vvpol(ixp,jy ,indzh,indexh) &
117             +p3*vvpol(ix ,jyp,indzh,indexh) &
118             +p4*vvpol(ixp,jyp,indzh,indexh)
119        usl=usl+uupol(ix ,jy ,indzh,indexh)+ &
120             uupol(ixp,jy ,indzh,indexh) &
121             +uupol(ix ,jyp,indzh,indexh)+uupol(ixp,jyp,indzh,indexh)
122        vsl=vsl+vvpol(ix ,jy ,indzh,indexh)+ &
123             vvpol(ixp,jy ,indzh,indexh) &
124             +vvpol(ix ,jyp,indzh,indexh)+vvpol(ixp,jyp,indzh,indexh)
125
126        usq=usq+uupol(ix ,jy ,indzh,indexh)* &
127             uupol(ix ,jy ,indzh,indexh)+ &
128             uupol(ixp,jy ,indzh,indexh)*uupol(ixp,jy ,indzh,indexh)+ &
129             uupol(ix ,jyp,indzh,indexh)*uupol(ix ,jyp,indzh,indexh)+ &
130             uupol(ixp,jyp,indzh,indexh)*uupol(ixp,jyp,indzh,indexh)
131        vsq=vsq+vvpol(ix ,jy ,indzh,indexh)* &
132             vvpol(ix ,jy ,indzh,indexh)+ &
133             vvpol(ixp,jy ,indzh,indexh)*vvpol(ixp,jy ,indzh,indexh)+ &
134             vvpol(ix ,jyp,indzh,indexh)*vvpol(ix ,jyp,indzh,indexh)+ &
135             vvpol(ixp,jyp,indzh,indexh)*vvpol(ixp,jyp,indzh,indexh)
136      else
137        u1(n)=p1*uu(ix ,jy ,indzh,indexh) &
138             +p2*uu(ixp,jy ,indzh,indexh) &
139             +p3*uu(ix ,jyp,indzh,indexh) &
140             +p4*uu(ixp,jyp,indzh,indexh)
141        v1(n)=p1*vv(ix ,jy ,indzh,indexh) &
142             +p2*vv(ixp,jy ,indzh,indexh) &
143             +p3*vv(ix ,jyp,indzh,indexh) &
144             +p4*vv(ixp,jyp,indzh,indexh)
145        usl=usl+uu(ix ,jy ,indzh,indexh)+uu(ixp,jy ,indzh,indexh) &
146             +uu(ix ,jyp,indzh,indexh)+uu(ixp,jyp,indzh,indexh)
147        vsl=vsl+vv(ix ,jy ,indzh,indexh)+vv(ixp,jy ,indzh,indexh) &
148             +vv(ix ,jyp,indzh,indexh)+vv(ixp,jyp,indzh,indexh)
149
150        usq=usq+uu(ix ,jy ,indzh,indexh)*uu(ix ,jy ,indzh,indexh)+ &
151             uu(ixp,jy ,indzh,indexh)*uu(ixp,jy ,indzh,indexh)+ &
152             uu(ix ,jyp,indzh,indexh)*uu(ix ,jyp,indzh,indexh)+ &
153             uu(ixp,jyp,indzh,indexh)*uu(ixp,jyp,indzh,indexh)
154        vsq=vsq+vv(ix ,jy ,indzh,indexh)*vv(ix ,jy ,indzh,indexh)+ &
155             vv(ixp,jy ,indzh,indexh)*vv(ixp,jy ,indzh,indexh)+ &
156             vv(ix ,jyp,indzh,indexh)*vv(ix ,jyp,indzh,indexh)+ &
157             vv(ixp,jyp,indzh,indexh)*vv(ixp,jyp,indzh,indexh)
158      endif
159      w1(n)=p1*ww(ix ,jy ,indzh,indexh) &
160           +p2*ww(ixp,jy ,indzh,indexh) &
161           +p3*ww(ix ,jyp,indzh,indexh) &
162           +p4*ww(ixp,jyp,indzh,indexh)
163      wsl=wsl+ww(ix ,jy ,indzh,indexh)+ww(ixp,jy ,indzh,indexh) &
164           +ww(ix ,jyp,indzh,indexh)+ww(ixp,jyp,indzh,indexh)
165      wsq=wsq+ww(ix ,jy ,indzh,indexh)*ww(ix ,jy ,indzh,indexh)+ &
166           ww(ixp,jy ,indzh,indexh)*ww(ixp,jy ,indzh,indexh)+ &
167           ww(ix ,jyp,indzh,indexh)*ww(ix ,jyp,indzh,indexh)+ &
168           ww(ixp,jyp,indzh,indexh)*ww(ixp,jyp,indzh,indexh)
169    end do
170
171
172  !**********************************
173  ! 2.) Linear vertical interpolation
174  !**********************************
175
176    uh(m)=dz2*u1(1)+dz1*u1(2)
177    vh(m)=dz2*v1(1)+dz1*v1(2)
178    wh(m)=dz2*w1(1)+dz1*w1(2)
179  end do
180
181
182  !************************************
183  ! 3.) Temporal interpolation (linear)
184  !************************************
185
186  u=(uh(1)*dt2+uh(2)*dt1)*dtt
187  v=(vh(1)*dt2+vh(2)*dt1)*dtt
188  w=(wh(1)*dt2+wh(2)*dt1)*dtt
189
190
191  ! Compute standard deviations
192  !****************************
193
194  xaux=usq-usl*usl/16.
195  if (xaux.lt.eps) then
196    usig=0.
197  else
198    usig=sqrt(xaux/15.)
199  endif
200
201  xaux=vsq-vsl*vsl/16.
202  if (xaux.lt.eps) then
203    vsig=0.
204  else
205    vsig=sqrt(xaux/15.)
206  endif
207
208
209  xaux=wsq-wsl*wsl/16.
210  if (xaux.lt.eps) then
211    wsig=0.
212  else
213    wsig=sqrt(xaux/15.)
214  endif
215
216end subroutine interpol_wind
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG