source: flexpart.git/src/interpol_all_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: 7.7 KB
Line 
1subroutine interpol_all_nests(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  !  Version for interpolating nested grids.                                   *
8  !                                                                            *
9  !    Author: A. Stohl                                                        *
10  !                                                                            *
11  !    9 February 1999                                                         *
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  ! itime [s]          current temporal position                               *
23  ! memtime(3) [s]     times of the wind fields in memory                      *
24  ! xt,yt,zt           coordinates position for which wind data shall be       *
25  !                    calculated                                              *
26  !                                                                            *
27  ! Constants:                                                                 *
28  !                                                                            *
29  !*****************************************************************************
30
31  use par_mod
32  use com_mod
33  use interpol_mod
34  use hanna_mod
35
36  implicit none
37
38  integer :: itime
39  real :: xt,yt,zt
40
41  ! Auxiliary variables needed for interpolation
42  real :: ust1(2),wst1(2),oli1(2),oliaux
43  real :: y1(2),y2(2),y3(2),rho1(2),rhograd1(2)
44  real :: usl,vsl,wsl,usq,vsq,wsq,xaux
45  integer :: i,m,n,indexh
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
73  !*****************************************
74  ! 1. Interpolate u*, w* and Obukhov length
75  !*****************************************
76
77  ! a) Bilinear horizontal interpolation
78
79  do m=1,2
80    indexh=memind(m)
81
82    ust1(m)=p1*ustarn(ix ,jy ,1,indexh,ngrid) &
83         + p2*ustarn(ixp,jy ,1,indexh,ngrid) &
84         + p3*ustarn(ix ,jyp,1,indexh,ngrid) &
85         + p4*ustarn(ixp,jyp,1,indexh,ngrid)
86    wst1(m)=p1*wstarn(ix ,jy ,1,indexh,ngrid) &
87         + p2*wstarn(ixp,jy ,1,indexh,ngrid) &
88         + p3*wstarn(ix ,jyp,1,indexh,ngrid) &
89         + p4*wstarn(ixp,jyp,1,indexh,ngrid)
90    oli1(m)=p1*olin(ix ,jy ,1,indexh,ngrid) &
91         + p2*olin(ixp,jy ,1,indexh,ngrid) &
92         + p3*olin(ix ,jyp,1,indexh,ngrid) &
93         + p4*olin(ixp,jyp,1,indexh,ngrid)
94  end do
95
96  ! b) Temporal interpolation
97
98  ust=(ust1(1)*dt2+ust1(2)*dt1)*dtt
99  wst=(wst1(1)*dt2+wst1(2)*dt1)*dtt
100  oliaux=(oli1(1)*dt2+oli1(2)*dt1)*dtt
101
102  if (oliaux.ne.0.) then
103    ol=1./oliaux
104  else
105    ol=99999.
106  endif
107
108
109  !*****************************************************
110  ! 2. Interpolate vertical profiles of u,v,w,rho,drhodz
111  !*****************************************************
112
113
114  ! Determine the level below the current position
115  !***********************************************
116
117  do i=2,nz
118    if (height(i).gt.zt) then
119      indz=i-1
120      indzp=i
121      goto 6
122    endif
123  end do
1246   continue
125
126  !**************************************
127  ! 1.) Bilinear horizontal interpolation
128  ! 2.) Temporal interpolation (linear)
129  !**************************************
130
131  ! Loop over 2 time steps and indz levels
132  !***************************************
133
134  do n=indz,indz+1
135    usl=0.
136    vsl=0.
137    wsl=0.
138    usq=0.
139    vsq=0.
140    wsq=0.
141    do m=1,2
142      indexh=memind(m)
143      y1(m)=p1*uun(ix ,jy ,n,indexh,ngrid) &
144           +p2*uun(ixp,jy ,n,indexh,ngrid) &
145           +p3*uun(ix ,jyp,n,indexh,ngrid) &
146           +p4*uun(ixp,jyp,n,indexh,ngrid)
147      y2(m)=p1*vvn(ix ,jy ,n,indexh,ngrid) &
148           +p2*vvn(ixp,jy ,n,indexh,ngrid) &
149           +p3*vvn(ix ,jyp,n,indexh,ngrid) &
150           +p4*vvn(ixp,jyp,n,indexh,ngrid)
151      y3(m)=p1*wwn(ix ,jy ,n,indexh,ngrid) &
152           +p2*wwn(ixp,jy ,n,indexh,ngrid) &
153           +p3*wwn(ix ,jyp,n,indexh,ngrid) &
154           +p4*wwn(ixp,jyp,n,indexh,ngrid)
155      rhograd1(m)=p1*drhodzn(ix ,jy ,n,indexh,ngrid) &
156           +p2*drhodzn(ixp,jy ,n,indexh,ngrid) &
157           +p3*drhodzn(ix ,jyp,n,indexh,ngrid) &
158           +p4*drhodzn(ixp,jyp,n,indexh,ngrid)
159      rho1(m)=p1*rhon(ix ,jy ,n,indexh,ngrid) &
160           +p2*rhon(ixp,jy ,n,indexh,ngrid) &
161           +p3*rhon(ix ,jyp,n,indexh,ngrid) &
162           +p4*rhon(ixp,jyp,n,indexh,ngrid)
163
164     usl=usl+uun(ix ,jy ,n,indexh,ngrid)+uun(ixp,jy ,n,indexh,ngrid) &
165          +uun(ix ,jyp,n,indexh,ngrid)+uun(ixp,jyp,n,indexh,ngrid)
166     vsl=vsl+vvn(ix ,jy ,n,indexh,ngrid)+vvn(ixp,jy ,n,indexh,ngrid) &
167          +vvn(ix ,jyp,n,indexh,ngrid)+vvn(ixp,jyp,n,indexh,ngrid)
168     wsl=wsl+wwn(ix ,jy ,n,indexh,ngrid)+wwn(ixp,jy ,n,indexh,ngrid) &
169          +wwn(ix ,jyp,n,indexh,ngrid)+wwn(ixp,jyp,n,indexh,ngrid)
170
171    usq=usq+uun(ix ,jy ,n,indexh,ngrid)*uun(ix ,jy ,n,indexh,ngrid)+ &
172         uun(ixp,jy ,n,indexh,ngrid)*uun(ixp,jy ,n,indexh,ngrid)+ &
173         uun(ix ,jyp,n,indexh,ngrid)*uun(ix ,jyp,n,indexh,ngrid)+ &
174         uun(ixp,jyp,n,indexh,ngrid)*uun(ixp,jyp,n,indexh,ngrid)
175    vsq=vsq+vvn(ix ,jy ,n,indexh,ngrid)*vvn(ix ,jy ,n,indexh,ngrid)+ &
176         vvn(ixp,jy ,n,indexh,ngrid)*vvn(ixp,jy ,n,indexh,ngrid)+ &
177         vvn(ix ,jyp,n,indexh,ngrid)*vvn(ix ,jyp,n,indexh,ngrid)+ &
178         vvn(ixp,jyp,n,indexh,ngrid)*vvn(ixp,jyp,n,indexh,ngrid)
179    wsq=wsq+wwn(ix ,jy ,n,indexh,ngrid)*wwn(ix ,jy ,n,indexh,ngrid)+ &
180         wwn(ixp,jy ,n,indexh,ngrid)*wwn(ixp,jy ,n,indexh,ngrid)+ &
181         wwn(ix ,jyp,n,indexh,ngrid)*wwn(ix ,jyp,n,indexh,ngrid)+ &
182         wwn(ixp,jyp,n,indexh,ngrid)*wwn(ixp,jyp,n,indexh,ngrid)
183    end do
184    uprof(n)=(y1(1)*dt2+y1(2)*dt1)*dtt
185    vprof(n)=(y2(1)*dt2+y2(2)*dt1)*dtt
186    wprof(n)=(y3(1)*dt2+y3(2)*dt1)*dtt
187    rhoprof(n)=(rho1(1)*dt2+rho1(2)*dt1)*dtt
188    rhogradprof(n)=(rhograd1(1)*dt2+rhograd1(2)*dt1)*dtt
189    indzindicator(n)=.false.
190
191  ! Compute standard deviations
192  !****************************
193
194    xaux=usq-usl*usl/8.
195    if (xaux.lt.eps) then
196      usigprof(n)=0.
197    else
198      usigprof(n)=sqrt(xaux/7.)
199    endif
200
201    xaux=vsq-vsl*vsl/8.
202    if (xaux.lt.eps) then
203      vsigprof(n)=0.
204    else
205      vsigprof(n)=sqrt(xaux/7.)
206    endif
207
208
209    xaux=wsq-wsl*wsl/8.
210    if (xaux.lt.eps) then
211      wsigprof(n)=0.
212    else
213      wsigprof(n)=sqrt(xaux/7.)
214    endif
215
216  end do
217
218end subroutine interpol_all_nests
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG