source: branches/flexpart91_hasod/src_parallel/interpol_wind_short_nests.f90 @ 8

Last change on this file since 8 was 8, checked in by hasod, 11 years ago

Added parallel version of Flexpart91

File size: 5.8 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  !                                       i   i  i  i
24  !*****************************************************************************
25  !                                                                            *
26  !  This subroutine interpolates the wind data to current trajectory position.*
27  !                                                                            *
28  !    Author: A. Stohl                                                        *
29  !                                                                            *
30  !    16 December 1997                                                        *
31  !                                                                            *
32  !  Revision March 2005 by AST : all output variables in common block         *
33  !                                                                            *
34  !*****************************************************************************
35  !                                                                            *
36  ! Variables:                                                                 *
37  ! u,v,w              wind components                                         *
38  ! itime [s]          current temporal position                               *
39  ! memtime(3) [s]     times of the wind fields in memory                      *
40  ! xt,yt,zt           coordinates position for which wind data shall be       *
41  !                    calculated                                              *
42  !                                                                            *
43  ! Constants:                                                                 *
44  !                                                                            *
45  !*****************************************************************************
46
47  use par_mod
48  use com_mod
49  use interpol_mod
50
51  implicit none
52
53  integer :: itime
54  real :: xt,yt,zt
55
56  ! Auxiliary variables needed for interpolation
57  real :: dz1,dz2,dz
58  real :: u1(2),v1(2),w1(2),uh(2),vh(2),wh(2)
59  integer :: i,m,n,indexh,indzh
60
61
62  !********************************************
63  ! Multilinear interpolation in time and space
64  !********************************************
65
66  ddx=xt-real(ix)
67  ddy=yt-real(jy)
68  rddx=1.-ddx
69  rddy=1.-ddy
70  p1=rddx*rddy
71  p2=ddx*rddy
72  p3=rddx*ddy
73  p4=ddx*ddy
74
75  ! Calculate variables for time interpolation
76  !*******************************************
77
78  dt1=real(itime-memtime(1))
79  dt2=real(memtime(2)-itime)
80  dtt=1./(dt1+dt2)
81
82  ! Determine the level below the current position for u,v
83  !*******************************************************
84
85  do i=2,nz
86    if (height(i).gt.zt) then
87      indz=i-1
88      goto 6
89    endif
90  end do
916   continue
92
93
94  ! Vertical distance to the level below and above current position
95  !****************************************************************
96
97  dz=1./(height(indz+1)-height(indz))
98  dz1=(zt-height(indz))*dz
99  dz2=(height(indz+1)-zt)*dz
100
101
102  !**********************************************************************
103  ! 1.) Bilinear horizontal interpolation
104  ! This has to be done separately for 6 fields (Temporal(2)*Vertical(3))
105  !**********************************************************************
106
107  ! Loop over 2 time steps and 2 levels
108  !************************************
109
110  do m=1,2
111    indexh=memind(m)
112    do n=1,2
113      indzh=indz+n-1
114
115      u1(n)=p1*uun(ix ,jy ,indzh,indexh,ngrid) &
116           +p2*uun(ixp,jy ,indzh,indexh,ngrid) &
117           +p3*uun(ix ,jyp,indzh,indexh,ngrid) &
118           +p4*uun(ixp,jyp,indzh,indexh,ngrid)
119      v1(n)=p1*vvn(ix ,jy ,indzh,indexh,ngrid) &
120           +p2*vvn(ixp,jy ,indzh,indexh,ngrid) &
121           +p3*vvn(ix ,jyp,indzh,indexh,ngrid) &
122           +p4*vvn(ixp,jyp,indzh,indexh,ngrid)
123      w1(n)=p1*wwn(ix ,jy ,indzh,indexh,ngrid) &
124           +p2*wwn(ixp,jy ,indzh,indexh,ngrid) &
125           +p3*wwn(ix ,jyp,indzh,indexh,ngrid) &
126           +p4*wwn(ixp,jyp,indzh,indexh,ngrid)
127
128    end do
129
130
131  !**********************************
132  ! 2.) Linear vertical interpolation
133  !**********************************
134
135    uh(m)=dz2*u1(1)+dz1*u1(2)
136    vh(m)=dz2*v1(1)+dz1*v1(2)
137    wh(m)=dz2*w1(1)+dz1*w1(2)
138  end do
139
140
141  !************************************
142  ! 3.) Temporal interpolation (linear)
143  !************************************
144
145  u=(uh(1)*dt2+uh(2)*dt1)*dtt
146  v=(vh(1)*dt2+vh(2)*dt1)*dtt
147  w=(wh(1)*dt2+wh(2)*dt1)*dtt
148
149end subroutine interpol_wind_short_nests
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG