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