source: flexpart.git/src/interpol_all.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: 8.4 KB
Line 
1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
3
4subroutine interpol_all(itime,xt,yt,zt)
5  !                          i   i  i  i
6  !*****************************************************************************
7  !                                                                            *
8  !  This subroutine interpolates everything that is needed for calculating the*
9  !  dispersion.                                                               *
10  !                                                                            *
11  !    Author: A. Stohl                                                        *
12  !                                                                            *
13  !    16 December 1997                                                        *
14  !                                                                            *
15  !  Revision March 2005 by AST : all output variables in common block cal-    *
16  !                               culation of standard deviation done in this  *
17  !                               routine rather than subroutine call in order *
18  !                               to save computation time                     *
19  !                                                                            *
20  !*****************************************************************************
21  !                                                                            *
22  ! Variables:                                                                 *
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  !                    culated                                                 *
27  !                                                                            *
28  ! Constants:                                                                 *
29  !                                                                            *
30  !*****************************************************************************
31
32  use par_mod
33  use com_mod
34  use interpol_mod
35  use hanna_mod
36
37  implicit none
38
39  integer :: itime
40  real :: xt,yt,zt
41
42  ! Auxiliary variables needed for interpolation
43  real :: ust1(2),wst1(2),oli1(2),oliaux
44  real :: y1(2),y2(2),y3(2),rho1(2),rhograd1(2)
45  real :: usl,vsl,wsl,usq,vsq,wsq,xaux
46  integer :: i,m,n,indexh
47  real,parameter :: eps=1.0e-30
48
49
50  !********************************************
51  ! Multilinear interpolation in time and space
52  !********************************************
53
54  ! Determine the lower left corner and its distance to the current position
55  !*************************************************************************
56
57  ddx=xt-real(ix)
58  ddy=yt-real(jy)
59  rddx=1.-ddx
60  rddy=1.-ddy
61  p1=rddx*rddy
62  p2=ddx*rddy
63  p3=rddx*ddy
64  p4=ddx*ddy
65
66  ! Calculate variables for time interpolation
67  !*******************************************
68
69  dt1=real(itime-memtime(1))
70  dt2=real(memtime(2)-itime)
71  dtt=1./(dt1+dt2)
72
73
74  !*****************************************
75  ! 1. Interpolate u*, w* and Obukhov length
76  !*****************************************
77
78  ! a) Bilinear horizontal interpolation
79
80  do m=1,2
81    indexh=memind(m)
82
83    ust1(m)=p1*ustar(ix ,jy ,1,indexh) &
84         + p2*ustar(ixp,jy ,1,indexh) &
85         + p3*ustar(ix ,jyp,1,indexh) &
86         + p4*ustar(ixp,jyp,1,indexh)
87    wst1(m)=p1*wstar(ix ,jy ,1,indexh) &
88         + p2*wstar(ixp,jy ,1,indexh) &
89         + p3*wstar(ix ,jyp,1,indexh) &
90         + p4*wstar(ixp,jyp,1,indexh)
91    oli1(m)=p1*oli(ix ,jy ,1,indexh) &
92         + p2*oli(ixp,jy ,1,indexh) &
93         + p3*oli(ix ,jyp,1,indexh) &
94         + p4*oli(ixp,jyp,1,indexh)
95  end do
96
97  ! b) Temporal interpolation
98
99  ust=(ust1(1)*dt2+ust1(2)*dt1)*dtt
100  wst=(wst1(1)*dt2+wst1(2)*dt1)*dtt
101  oliaux=(oli1(1)*dt2+oli1(2)*dt1)*dtt
102
103  if (oliaux.ne.0.) then
104    ol=1./oliaux
105  else
106    ol=99999.
107  endif
108
109
110  !*****************************************************
111  ! 2. Interpolate vertical profiles of u,v,w,rho,drhodz
112  !*****************************************************
113
114
115  ! Determine the level below the current position
116  !***********************************************
117
118  do i=2,nz
119    if (height(i).gt.zt) then
120      indz=i-1
121      indzp=i
122      goto 6
123    endif
124  end do
1256   continue
126
127  !**************************************
128  ! 1.) Bilinear horizontal interpolation
129  ! 2.) Temporal interpolation (linear)
130  !**************************************
131
132  ! Loop over 2 time steps and indz levels
133  !***************************************
134
135  do n=indz,indzp
136    usl=0.
137    vsl=0.
138    wsl=0.
139    usq=0.
140    vsq=0.
141    wsq=0.
142    do m=1,2
143      indexh=memind(m)
144      if (ngrid.lt.0) then
145        y1(m)=p1*uupol(ix ,jy ,n,indexh) &
146             +p2*uupol(ixp,jy ,n,indexh) &
147             +p3*uupol(ix ,jyp,n,indexh) &
148             +p4*uupol(ixp,jyp,n,indexh)
149        y2(m)=p1*vvpol(ix ,jy ,n,indexh) &
150             +p2*vvpol(ixp,jy ,n,indexh) &
151             +p3*vvpol(ix ,jyp,n,indexh) &
152             +p4*vvpol(ixp,jyp,n,indexh)
153        usl=usl+uupol(ix ,jy ,n,indexh)+uupol(ixp,jy ,n,indexh) &
154             +uupol(ix ,jyp,n,indexh)+uupol(ixp,jyp,n,indexh)
155        vsl=vsl+vvpol(ix ,jy ,n,indexh)+vvpol(ixp,jy ,n,indexh) &
156             +vvpol(ix ,jyp,n,indexh)+vvpol(ixp,jyp,n,indexh)
157
158        usq=usq+uupol(ix ,jy ,n,indexh)*uupol(ix ,jy ,n,indexh)+ &
159             uupol(ixp,jy ,n,indexh)*uupol(ixp,jy ,n,indexh)+ &
160             uupol(ix ,jyp,n,indexh)*uupol(ix ,jyp,n,indexh)+ &
161             uupol(ixp,jyp,n,indexh)*uupol(ixp,jyp,n,indexh)
162        vsq=vsq+vvpol(ix ,jy ,n,indexh)*vvpol(ix ,jy ,n,indexh)+ &
163             vvpol(ixp,jy ,n,indexh)*vvpol(ixp,jy ,n,indexh)+ &
164             vvpol(ix ,jyp,n,indexh)*vvpol(ix ,jyp,n,indexh)+ &
165             vvpol(ixp,jyp,n,indexh)*vvpol(ixp,jyp,n,indexh)
166      else
167        y1(m)=p1*uu(ix ,jy ,n,indexh) &
168             +p2*uu(ixp,jy ,n,indexh) &
169             +p3*uu(ix ,jyp,n,indexh) &
170             +p4*uu(ixp,jyp,n,indexh)
171        y2(m)=p1*vv(ix ,jy ,n,indexh) &
172             +p2*vv(ixp,jy ,n,indexh) &
173             +p3*vv(ix ,jyp,n,indexh) &
174             +p4*vv(ixp,jyp,n,indexh)
175        usl=usl+uu(ix ,jy ,n,indexh)+uu(ixp,jy ,n,indexh) &
176             +uu(ix ,jyp,n,indexh)+uu(ixp,jyp,n,indexh)
177        vsl=vsl+vv(ix ,jy ,n,indexh)+vv(ixp,jy ,n,indexh) &
178             +vv(ix ,jyp,n,indexh)+vv(ixp,jyp,n,indexh)
179
180        usq=usq+uu(ix ,jy ,n,indexh)*uu(ix ,jy ,n,indexh)+ &
181             uu(ixp,jy ,n,indexh)*uu(ixp,jy ,n,indexh)+ &
182             uu(ix ,jyp,n,indexh)*uu(ix ,jyp,n,indexh)+ &
183             uu(ixp,jyp,n,indexh)*uu(ixp,jyp,n,indexh)
184        vsq=vsq+vv(ix ,jy ,n,indexh)*vv(ix ,jy ,n,indexh)+ &
185             vv(ixp,jy ,n,indexh)*vv(ixp,jy ,n,indexh)+ &
186             vv(ix ,jyp,n,indexh)*vv(ix ,jyp,n,indexh)+ &
187             vv(ixp,jyp,n,indexh)*vv(ixp,jyp,n,indexh)
188      endif
189      y3(m)=p1*ww(ix ,jy ,n,indexh) &
190           +p2*ww(ixp,jy ,n,indexh) &
191           +p3*ww(ix ,jyp,n,indexh) &
192           +p4*ww(ixp,jyp,n,indexh)
193      rhograd1(m)=p1*drhodz(ix ,jy ,n,indexh) &
194           +p2*drhodz(ixp,jy ,n,indexh) &
195           +p3*drhodz(ix ,jyp,n,indexh) &
196           +p4*drhodz(ixp,jyp,n,indexh)
197      rho1(m)=p1*rho(ix ,jy ,n,indexh) &
198           +p2*rho(ixp,jy ,n,indexh) &
199           +p3*rho(ix ,jyp,n,indexh) &
200           +p4*rho(ixp,jyp,n,indexh)
201      wsl=wsl+ww(ix ,jy ,n,indexh)+ww(ixp,jy ,n,indexh) &
202           +ww(ix ,jyp,n,indexh)+ww(ixp,jyp,n,indexh)
203      wsq=wsq+ww(ix ,jy ,n,indexh)*ww(ix ,jy ,n,indexh)+ &
204           ww(ixp,jy ,n,indexh)*ww(ixp,jy ,n,indexh)+ &
205           ww(ix ,jyp,n,indexh)*ww(ix ,jyp,n,indexh)+ &
206           ww(ixp,jyp,n,indexh)*ww(ixp,jyp,n,indexh)
207    end do
208    uprof(n)=(y1(1)*dt2+y1(2)*dt1)*dtt
209    vprof(n)=(y2(1)*dt2+y2(2)*dt1)*dtt
210    wprof(n)=(y3(1)*dt2+y3(2)*dt1)*dtt
211    rhoprof(n)=(rho1(1)*dt2+rho1(2)*dt1)*dtt
212    rhogradprof(n)=(rhograd1(1)*dt2+rhograd1(2)*dt1)*dtt
213    indzindicator(n)=.false.
214
215  ! Compute standard deviations
216  !****************************
217
218    xaux=usq-usl*usl/8.
219    if (xaux.lt.eps) then
220      usigprof(n)=0.
221    else
222      usigprof(n)=sqrt(xaux/7.)
223    endif
224
225    xaux=vsq-vsl*vsl/8.
226    if (xaux.lt.eps) then
227      vsigprof(n)=0.
228    else
229      vsigprof(n)=sqrt(xaux/7.)
230    endif
231
232
233    xaux=wsq-wsl*wsl/8.
234    if (xaux.lt.eps) then
235      wsigprof(n)=0.
236    else
237      wsigprof(n)=sqrt(xaux/7.)
238    endif
239
240  end do
241
242
243end subroutine interpol_all
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG