source: branches/petra/src/interpol_rain.f90 @ 37

Last change on this file since 37 was 37, checked in by pesei, 9 years ago

Wet dep quick fix and other small changes. Wet depo quick fix not final yet.

File size: 8.9 KB
Line 
1!**********************************************************************
2! Copyright 1998-2015                                                 *
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_rain(yy1,yy2,yy3,nxmax,nymax,nzmax,nx,ny, &
23       memind,xt,yt,level,itime1,itime2,itime,yint1,yint2,yint3, &
24       icbot,ictop)
25  !                          i   i   i    i    i     i   i
26  !i    i    i  i    i     i      i      i     o     o     o
27  !****************************************************************************
28  !                                                                           *
29  !  Interpolation of meteorological fields on 2-d model layers.              *
30  !  In horizontal direction bilinear interpolation interpolation is used.    *
31  !  Temporally a linear interpolation is used.                               *
32  !  Three fields are interpolated at the same time.                          *
33  !                                                                           *
34  !  This is a special version of levlininterpol to save CPU time.            *
35  !                                                                           *
36  !  1 first time                                                             *
37  !  2 second time                                                            *
38  !                                                                           *
39  !                                                                           *
40  !     Author: A. Stohl                                                      *
41  !                                                                           *
42  !     30 August 1996                                                        *
43  !****************************************************************************
44  ! Petra Seibert, 2011/2012-2015:
45  !  Add interpolation of cloud bottom and cloud thickness
46  !  for fix to SE's new wet scavenging scheme
47  !                                                                           *
48  !****************************************************************************
49  !                                                                           *
50  ! Variables:                                                                *
51  !                                                                           *
52  ! dt1,dt2              time differences between fields and current position *
53  ! dz1,dz2              z distance between levels and current position       *
54  ! height(nzmax)        heights of the model levels                          *
55  ! indexh               help variable                                        *
56  ! indz                 the level closest to the current trajectory position *
57  ! indzh                help variable                                        *
58  ! itime                current time                                         *
59  ! itime1               time of the first wind field                         *
60  ! itime2               time of the second wind field                        *
61  ! ix,jy                x,y coordinates of lower left subgrid point          *
62  ! level                level at which interpolation shall be done           *
63  ! memind(3)            points to the places of the wind fields              *
64  ! nx,ny                actual field dimensions in x,y and z direction       *
65  ! nxmax,nymax,nzmax    maximum field dimensions in x,y and z direction      *
66  ! xt                   current x coordinate                                 *
67  ! yint                 the final interpolated value                         *
68  ! yt                   current y coordinate                                 *
69  ! yy(0:nxmax,0:nymax,nzmax,3) meteorological field used for interpolation   *
70  ! zt                   current z coordinate                                 *
71  !                                                                           *
72  !****************************************************************************
73 
74  use com_mod, only: icloudbot, icloudthck
75  use par_mod, only: icmv
76
77  implicit none
78
79  integer :: nx,ny,nxmax,nymax,nzmax,memind(2),m,ix,jy,ixp,jyp
80  integer :: itime,itime1,itime2,level,indexh,ip1,ip2,ip3,ip4
81  integer :: icbot,icthck,ictop,ipsum
82  real :: yy1(0:nxmax-1,0:nymax-1,nzmax,2)
83  real :: yy2(0:nxmax-1,0:nymax-1,nzmax,2)
84  real :: yy3(0:nxmax-1,0:nymax-1,nzmax,2)
85  real :: ddx,ddy,rddx,rddy,dt1,dt2,dt,y1(2),y2(2),y3(2),cbot(2),cthck(2)
86  real :: xt,yt,yint1,yint2,yint3,p1,p2,p3,p4
87
88
89  ! If point at border of grid -> small displacement into grid
90  !***********************************************************
91
92  if (xt.ge.real(nx-1)) xt=real(nx-1)-0.00001
93  if (yt.ge.real(ny-1)) yt=real(ny-1)-0.00001
94
95
96
97  !**********************************************************************
98  ! 1.) Bilinear horizontal interpolation
99  ! This has to be done separately for 2 fields (Temporal)
100  !*******************************************************
101
102  ! Determine the lower left corner and its distance to the current position
103  !*************************************************************************
104
105  ix=int(xt)
106  jy=int(yt)
107  ixp=ix+1
108  jyp=jy+1
109  ddx=xt-real(ix)
110  ddy=yt-real(jy)
111  rddx=1.-ddx
112  rddy=1.-ddy
113  p1=rddx*rddy
114  p2=ddx*rddy
115  p3=rddx*ddy
116  p4=ddx*ddy
117
118
119  ! Loop over 2 time steps
120  !***********************
121
122  memind_loop: &
123  do m=1,2
124
125    indexh=memind(m)
126
127    y1(m)=p1*yy1(ix ,jy ,level,indexh) &
128         + p2*yy1(ixp,jy ,level,indexh) &
129         + p3*yy1(ix ,jyp,level,indexh) &
130         + p4*yy1(ixp,jyp,level,indexh)
131    y2(m)=p1*yy2(ix ,jy ,level,indexh) &
132         + p2*yy2(ixp,jy ,level,indexh) &
133         + p3*yy2(ix ,jyp,level,indexh) &
134         + p4*yy2(ixp,jyp,level,indexh)
135    y3(m)=p1*yy3(ix ,jy ,level,indexh) &
136         + p2*yy3(ixp,jy ,level,indexh) &
137         + p3*yy3(ix ,jyp,level,indexh) &
138         + p4*yy3(ixp,jyp,level,indexh)
139
140
141! PS clouds:
142    ip1=1
143    ip2=1
144    ip3=1
145    ip4=1
146    if (icloudbot(ix ,jy ,indexh) .eq. icmv) ip1=0
147    if (icloudbot(ixp,jy ,indexh) .eq. icmv) ip2=0
148    if (icloudbot(ix ,jyp,indexh) .eq. icmv) ip3=0
149    if (icloudbot(ixp,jyp,indexh) .eq. icmv) ip4=0
150    ipsum = ip1+ip2+ip3+ip4
151    if (ipsum .eq. 0) then
152      cbot(m)=icmv
153    else
154      cbot(m)=(ip1*p1*icloudbot(ix ,jy ,indexh) &
155             + ip2*p2*icloudbot(ixp,jy ,indexh) &
156             + ip3*p3*icloudbot(ix ,jyp,indexh) &
157             + ip4*p4*icloudbot(ixp,jyp,indexh)) / ipsum
158    endif
159
160    ip1=1
161    ip2=1
162    ip3=1
163    ip4=1
164    if (icloudthck(ix ,jy ,indexh) .eq. icmv) ip1=0
165    if (icloudthck(ixp,jy ,indexh) .eq. icmv) ip2=0
166    if (icloudthck(ix ,jyp,indexh) .eq. icmv) ip3=0
167    if (icloudthck(ixp,jyp,indexh) .eq. icmv) ip4=0
168    ipsum = ip1+ip2+ip3+ip4
169    if (ipsum .eq. 0) then
170      cthck(m)=icmv
171    else
172      cthck(m)=(ip1*p1*icloudthck(ix ,jy ,indexh) &
173              + ip2*p2*icloudthck(ixp,jy ,indexh) &
174              + ip3*p3*icloudthck(ix ,jyp,indexh) &
175              + ip4*p4*icloudthck(ixp,jyp,indexh)) / ipsum
176    endif
177! PS end clouds
178
179  enddo memind_loop
180
181
182  !************************************
183  ! 2.) Temporal interpolation (linear)
184  !************************************
185
186  dt1=real(itime-itime1)
187  dt2=real(itime2-itime)
188  dt=dt1+dt2
189
190  yint1=(y1(1)*dt2+y1(2)*dt1)/dt
191  yint2=(y2(1)*dt2+y2(2)*dt1)/dt
192  yint3=(y3(1)*dt2+y3(2)*dt1)/dt
193
194! PS clouds:
195  icbot = nint( (cbot(1)*dt2 + cbot(2)*dt1)/dt )
196  if (nint(cbot(1)) .eq. icmv) icbot=cbot(2)
197  if (nint(cbot(2)) .eq. icmv) icbot=cbot(1)
198
199  icthck = nint( (cthck(1)*dt2 + cthck(2)*dt1)/dt )
200  if (nint(cthck(1)) .eq. icmv) icthck=cthck(2)
201  if (nint(cthck(2)) .eq. icmv) icthck=cthck(1)
202
203  if (icbot .ne. icmv .and. icthck .ne. icmv) then
204    ictop = icbot + icthck ! convert cloud thickness to cloud top
205  else
206    icbot=icmv
207    ictop=icmv
208  endif
209! PS end clouds
210
211
212end subroutine interpol_rain
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG