source: branches/jerome/src_flexwrf_v3.1/interpol_wind_short_nests.f90 @ 16

Last change on this file since 16 was 16, checked in by jebri, 10 years ago

sources for flexwrf v3.1

File size: 6.4 KB
Line 
1!**********************************************************************
2! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010         *
3! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa,             *
4! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann   *
5!                                                                     *
6! This file is part of FLEXPART.                                      *
7!                                                                     *
8! FLEXPART is free software: you can redistribute it and/or modify    *
9! it under the terms of the GNU General Public License as published by*
10! the Free Software Foundation, either version 3 of the License, or   *
11! (at your option) any later version.                                 *
12!                                                                     *
13! FLEXPART is distributed in the hope that it will be useful,         *
14! but WITHOUT ANY WARRANTY; without even the implied warranty of      *
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
16! GNU General Public License for more details.                        *
17!                                                                     *
18! You should have received a copy of the GNU General Public License   *
19! along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
20!**********************************************************************
21
22subroutine interpol_wind_short_nests(itime,xt,yt,zt, &
23    uprof,vprof,wprof, usigprof,vsigprof,wsigprof, &
24    rhoprof,rhogradprof, tkeprof,pttprof, &
25    u,v,w,usig,vsig,wsig,pvi, &
26    p1,p2,p3,p4,ddx,ddy,rddx,rddy,dtt,dt1,dt2, &
27    ix,jy,ixp,jyp,ngrid,indz,indzp, depoindicator, &
28    indzindicator)
29  !                                      i   i  i  i
30  !*****************************************************************************
31  !                                                                            *
32  !  This subroutine interpolates the wind data to current trajectory position.*
33  !                                                                            *
34  !    Author: A. Stohl                                                        *
35  !                                                                            *
36  !    16 December 1997                                                        *
37  !                                                                            *
38  !  Revision March 2005 by AST : all output variables in common block         *
39  !                                                                            *
40  !*****************************************************************************
41  !                                                                            *
42  ! Variables:                                                                 *
43  ! u,v,w              wind components                                         *
44  ! itime [s]          current temporal position                               *
45  ! memtime(3) [s]     times of the wind fields in memory                      *
46  ! xt,yt,zt           coordinates position for which wind data shall be       *
47  !                    calculated                                              *
48  !                                                                            *
49  ! Constants:                                                                 *
50  !                                                                            *
51  !*****************************************************************************
52
53  use par_mod
54  use com_mod
55!  use interpol_mod
56
57  implicit none
58
59  integer :: itime
60  real :: xt,yt,zt
61
62  ! Auxiliary variables needed for interpolation
63  real :: dz1,dz2,dz
64  real :: u1(2),v1(2),w1(2),uh(2),vh(2),wh(2)
65  integer :: i,m,n,indexh,indzh
66
67  real :: uprof(nzmax),vprof(nzmax),wprof(nzmax)
68  real :: usigprof(nzmax),vsigprof(nzmax),wsigprof(nzmax)
69  real :: rhoprof(nzmax),rhogradprof(nzmax)
70  real :: tkeprof(nzmax),pttprof(nzmax)
71  real :: u,v,w,usig,vsig,wsig,pvi
72
73  real :: p1,p2,p3,p4,ddx,ddy,rddx,rddy,dtt,dt1,dt2
74  integer :: ix,jy,ixp,jyp,ngrid,indz,indzp
75  logical :: depoindicator(maxspec)
76  logical :: indzindicator(nzmax)
77
78
79  !********************************************
80  ! Multilinear interpolation in time and space
81  !********************************************
82
83  ddx=xt-real(ix)
84  ddy=yt-real(jy)
85  rddx=1.-ddx
86  rddy=1.-ddy
87  p1=rddx*rddy
88  p2=ddx*rddy
89  p3=rddx*ddy
90  p4=ddx*ddy
91
92  ! Calculate variables for time interpolation
93  !*******************************************
94
95  dt1=real(itime-memtime(1))
96  dt2=real(memtime(2)-itime)
97  dtt=1./(dt1+dt2)
98
99  ! Determine the level below the current position for u,v
100  !*******************************************************
101
102  do i=2,nz
103    if (height(i).gt.zt) then
104      indz=i-1
105      goto 6
106    endif
107  end do
1086   continue
109
110
111  ! Vertical distance to the level below and above current position
112  !****************************************************************
113
114  dz=1./(height(indz+1)-height(indz))
115  dz1=(zt-height(indz))*dz
116  dz2=(height(indz+1)-zt)*dz
117
118
119  !**********************************************************************
120  ! 1.) Bilinear horizontal interpolation
121  ! This has to be done separately for 6 fields (Temporal(2)*Vertical(3))
122  !**********************************************************************
123
124  ! Loop over 2 time steps and 2 levels
125  !************************************
126
127  do m=1,2
128    indexh=memind(m)
129    do n=1,2
130      indzh=indz+n-1
131
132      u1(n)=p1*uun(ix ,jy ,indzh,indexh,ngrid) &
133           +p2*uun(ixp,jy ,indzh,indexh,ngrid) &
134           +p3*uun(ix ,jyp,indzh,indexh,ngrid) &
135           +p4*uun(ixp,jyp,indzh,indexh,ngrid)
136      v1(n)=p1*vvn(ix ,jy ,indzh,indexh,ngrid) &
137           +p2*vvn(ixp,jy ,indzh,indexh,ngrid) &
138           +p3*vvn(ix ,jyp,indzh,indexh,ngrid) &
139           +p4*vvn(ixp,jyp,indzh,indexh,ngrid)
140      w1(n)=p1*wwn(ix ,jy ,indzh,indexh,ngrid) &
141           +p2*wwn(ixp,jy ,indzh,indexh,ngrid) &
142           +p3*wwn(ix ,jyp,indzh,indexh,ngrid) &
143           +p4*wwn(ixp,jyp,indzh,indexh,ngrid)
144
145    end do
146
147
148  !**********************************
149  ! 2.) Linear vertical interpolation
150  !**********************************
151
152    uh(m)=dz2*u1(1)+dz1*u1(2)
153    vh(m)=dz2*v1(1)+dz1*v1(2)
154    wh(m)=dz2*w1(1)+dz1*w1(2)
155  end do
156
157
158  !************************************
159  ! 3.) Temporal interpolation (linear)
160  !************************************
161
162  u=(uh(1)*dt2+uh(2)*dt1)*dtt
163  v=(vh(1)*dt2+vh(2)*dt1)*dtt
164  w=(wh(1)*dt2+wh(2)*dt1)*dtt
165
166end subroutine interpol_wind_short_nests
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG