source: branches/flexpart91_hasod/src_parallel/interpol_misslev.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: 7.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_misslev(n)
23  !                            i
24  !*****************************************************************************
25  !                                                                            *
26  !  This subroutine interpolates u,v,w, density and density gradients.        *
27  !                                                                            *
28  !    Author: A. Stohl                                                        *
29  !                                                                            *
30  !    16 December 1997                                                        *
31  !    Update: 2 March 1999                                                    *
32  !                                                                            *
33  !  Revision March 2005 by AST : all output variables in common block cal-    *
34  !                               culation of standard deviation done in this  *
35  !                               routine rather than subroutine call in order *
36  !                               to save computation time                     *
37  !                                                                            *
38  !*****************************************************************************
39  !                                                                            *
40  ! Variables:                                                                 *
41  ! n                  level                                                   *
42  !                                                                            *
43  ! Constants:                                                                 *
44  !                                                                            *
45  !*****************************************************************************
46
47  use par_mod
48  use com_mod
49  use interpol_mod
50  use hanna_mod
51
52  implicit none
53
54  ! Auxiliary variables needed for interpolation
55  real :: y1(2),y2(2),y3(2),rho1(2),rhograd1(2)
56  real :: usl,vsl,wsl,usq,vsq,wsq,xaux
57  integer :: m,n,indexh
58  real,parameter :: eps=1.0e-30
59
60
61  !********************************************
62  ! Multilinear interpolation in time and space
63  !********************************************
64
65
66  !**************************************
67  ! 1.) Bilinear horizontal interpolation
68  ! 2.) Temporal interpolation (linear)
69  !**************************************
70
71  ! Loop over 2 time steps
72  !***********************
73
74  usl=0.
75  vsl=0.
76  wsl=0.
77  usq=0.
78  vsq=0.
79  wsq=0.
80  do m=1,2
81    indexh=memind(m)
82    if (ngrid.lt.0) then
83      y1(m)=p1*uupol(ix ,jy ,n,indexh) &
84           +p2*uupol(ixp,jy ,n,indexh) &
85           +p3*uupol(ix ,jyp,n,indexh) &
86           +p4*uupol(ixp,jyp,n,indexh)
87      y2(m)=p1*vvpol(ix ,jy ,n,indexh) &
88           +p2*vvpol(ixp,jy ,n,indexh) &
89           +p3*vvpol(ix ,jyp,n,indexh) &
90           +p4*vvpol(ixp,jyp,n,indexh)
91        usl=usl+uupol(ix ,jy ,n,indexh)+uupol(ixp,jy ,n,indexh) &
92             +uupol(ix ,jyp,n,indexh)+uupol(ixp,jyp,n,indexh)
93        vsl=vsl+vvpol(ix ,jy ,n,indexh)+vvpol(ixp,jy ,n,indexh) &
94             +vvpol(ix ,jyp,n,indexh)+vvpol(ixp,jyp,n,indexh)
95
96        usq=usq+uupol(ix ,jy ,n,indexh)*uupol(ix ,jy ,n,indexh)+ &
97             uupol(ixp,jy ,n,indexh)*uupol(ixp,jy ,n,indexh)+ &
98             uupol(ix ,jyp,n,indexh)*uupol(ix ,jyp,n,indexh)+ &
99             uupol(ixp,jyp,n,indexh)*uupol(ixp,jyp,n,indexh)
100        vsq=vsq+vvpol(ix ,jy ,n,indexh)*vvpol(ix ,jy ,n,indexh)+ &
101             vvpol(ixp,jy ,n,indexh)*vvpol(ixp,jy ,n,indexh)+ &
102             vvpol(ix ,jyp,n,indexh)*vvpol(ix ,jyp,n,indexh)+ &
103             vvpol(ixp,jyp,n,indexh)*vvpol(ixp,jyp,n,indexh)
104    else
105      y1(m)=p1*uu(ix ,jy ,n,indexh) &
106           +p2*uu(ixp,jy ,n,indexh) &
107           +p3*uu(ix ,jyp,n,indexh) &
108           +p4*uu(ixp,jyp,n,indexh)
109      y2(m)=p1*vv(ix ,jy ,n,indexh) &
110           +p2*vv(ixp,jy ,n,indexh) &
111           +p3*vv(ix ,jyp,n,indexh) &
112           +p4*vv(ixp,jyp,n,indexh)
113      usl=usl+uu(ix ,jy ,n,indexh)+uu(ixp,jy ,n,indexh) &
114           +uu(ix ,jyp,n,indexh)+uu(ixp,jyp,n,indexh)
115      vsl=vsl+vv(ix ,jy ,n,indexh)+vv(ixp,jy ,n,indexh) &
116           +vv(ix ,jyp,n,indexh)+vv(ixp,jyp,n,indexh)
117
118      usq=usq+uu(ix ,jy ,n,indexh)*uu(ix ,jy ,n,indexh)+ &
119           uu(ixp,jy ,n,indexh)*uu(ixp,jy ,n,indexh)+ &
120           uu(ix ,jyp,n,indexh)*uu(ix ,jyp,n,indexh)+ &
121           uu(ixp,jyp,n,indexh)*uu(ixp,jyp,n,indexh)
122      vsq=vsq+vv(ix ,jy ,n,indexh)*vv(ix ,jy ,n,indexh)+ &
123           vv(ixp,jy ,n,indexh)*vv(ixp,jy ,n,indexh)+ &
124           vv(ix ,jyp,n,indexh)*vv(ix ,jyp,n,indexh)+ &
125           vv(ixp,jyp,n,indexh)*vv(ixp,jyp,n,indexh)
126    endif
127    y3(m)=p1*ww(ix ,jy ,n,indexh) &
128         +p2*ww(ixp,jy ,n,indexh) &
129         +p3*ww(ix ,jyp,n,indexh) &
130         +p4*ww(ixp,jyp,n,indexh)
131    rhograd1(m)=p1*drhodz(ix ,jy ,n,indexh) &
132         +p2*drhodz(ixp,jy ,n,indexh) &
133         +p3*drhodz(ix ,jyp,n,indexh) &
134         +p4*drhodz(ixp,jyp,n,indexh)
135    rho1(m)=p1*rho(ix ,jy ,n,indexh) &
136         +p2*rho(ixp,jy ,n,indexh) &
137         +p3*rho(ix ,jyp,n,indexh) &
138         +p4*rho(ixp,jyp,n,indexh)
139    wsl=wsl+ww(ix ,jy ,n,indexh)+ww(ixp,jy ,n,indexh) &
140         +ww(ix ,jyp,n,indexh)+ww(ixp,jyp,n,indexh)
141    wsq=wsq+ww(ix ,jy ,n,indexh)*ww(ix ,jy ,n,indexh)+ &
142         ww(ixp,jy ,n,indexh)*ww(ixp,jy ,n,indexh)+ &
143         ww(ix ,jyp,n,indexh)*ww(ix ,jyp,n,indexh)+ &
144         ww(ixp,jyp,n,indexh)*ww(ixp,jyp,n,indexh)
145  end do
146  uprof(n)=(y1(1)*dt2+y1(2)*dt1)*dtt
147  vprof(n)=(y2(1)*dt2+y2(2)*dt1)*dtt
148  wprof(n)=(y3(1)*dt2+y3(2)*dt1)*dtt
149  rhoprof(n)=(rho1(1)*dt2+rho1(2)*dt1)*dtt
150  rhogradprof(n)=(rhograd1(1)*dt2+rhograd1(2)*dt1)*dtt
151  indzindicator(n)=.false.
152
153
154  ! Compute standard deviations
155  !****************************
156
157  xaux=usq-usl*usl/8.
158  if (xaux.lt.eps) then
159    usigprof(n)=0.
160  else
161    usigprof(n)=sqrt(xaux/7.)
162  endif
163
164  xaux=vsq-vsl*vsl/8.
165  if (xaux.lt.eps) then
166    vsigprof(n)=0.
167  else
168    vsigprof(n)=sqrt(xaux/7.)
169  endif
170
171
172  xaux=wsq-wsl*wsl/8.
173  if (xaux.lt.eps) then
174    wsigprof(n)=0.
175  else
176    wsigprof(n)=sqrt(xaux/7.)
177  endif
178
179
180end subroutine interpol_misslev
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG