source: flexpart.git/src/interpol_rain.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: 5.7 KB
Line 
1subroutine interpol_rain(yy1,yy2,yy3,nxmax,nymax,nzmax,nx, &
2       ny,iwftouse,xt,yt,level,itime1,itime2,itime,yint1,yint2,yint3)
3  !                          i   i   i    i    i     i   i
4  !i    i    i  i    i     i      i      i     o     o     o
5  !****************************************************************************
6  !                                                                           *
7  !  Interpolation of meteorological fields on 2-d model layers.              *
8  !  In horizontal direction bilinear interpolation interpolation is used.    *
9  !  Temporally a linear interpolation is used.                               *
10  !  Three fields are interpolated at the same time.                          *
11  !                                                                           *
12  !  This is a special version of levlininterpol to save CPU time.            *
13  !                                                                           *
14  !  1 first time                                                             *
15  !  2 second time                                                            *
16  !                                                                           *
17  !                                                                           *
18  !     Author: A. Stohl                                                      *
19  !                                                                           *
20  !     30 August 1996                                                        *
21  !                                                                           *
22  !****************************************************************************
23  !                                                                           *
24  ! Variables:                                                                *
25  !                                                                           *
26  ! dt1,dt2              time differences between fields and current position *
27  ! dz1,dz2              z distance between levels and current position       *
28  ! height(nzmax)        heights of the model levels                          *
29  ! indexh               help variable                                        *
30  ! indz                 the level closest to the current trajectory position *
31  ! indzh                help variable                                        *
32  ! itime                current time                                         *
33  ! itime1               time of the first wind field                         *
34  ! itime2               time of the second wind field                        *
35  ! ix,jy                x,y coordinates of lower left subgrid point          *
36  ! level                level at which interpolation shall be done           *
37  ! iwftouse             points to the place of the wind field                *
38  ! nx,ny                actual field dimensions in x,y and z direction       *
39  ! nxmax,nymax,nzmax    maximum field dimensions in x,y and z direction      *
40  ! xt                   current x coordinate                                 *
41  ! yint                 the final interpolated value                         *
42  ! yt                   current y coordinate                                 *
43  ! yy(0:nxmax,0:nymax,nzmax,3) meteorological field used for interpolation   *
44  ! zt                   current z coordinate                                 *
45  !                                                                           *
46  !****************************************************************************
47  use par_mod, only: numwfmem
48
49  implicit none
50
51  integer :: nx,ny,nxmax,nymax,nzmax,memind(numwfmem),m,ix,jy,ixp,jyp
52  integer :: itime,itime1,itime2,level,indexh
53  real :: yy1(0:nxmax-1,0:nymax-1,nzmax,numwfmem)
54  real :: yy2(0:nxmax-1,0:nymax-1,nzmax,numwfmem)
55  real :: yy3(0:nxmax-1,0:nymax-1,nzmax,numwfmem)
56  real :: ddx,ddy,rddx,rddy,dt1,dt2,dt,y1(2),y2(2),y3(2)
57  real :: xt,yt,yint1,yint2,yint3,p1,p2,p3,p4
58  integer :: iwftouse
59
60
61
62  ! If point at border of grid -> small displacement into grid
63  !***********************************************************
64
65  if (xt.ge.real(nx-1)) xt=real(nx-1)-0.00001
66  if (yt.ge.real(ny-1)) yt=real(ny-1)-0.00001
67
68
69
70  !**********************************************************************
71  ! 1.) Bilinear horizontal interpolation
72  ! This has to be done separately for 2 fields (Temporal)
73  !*******************************************************
74
75  ! Determine the lower left corner and its distance to the current position
76  !*************************************************************************
77
78  ix=int(xt)
79  jy=int(yt)
80  ixp=ix+1
81  jyp=jy+1
82  ddx=xt-real(ix)
83  ddy=yt-real(jy)
84  rddx=1.-ddx
85  rddy=1.-ddy
86  p1=rddx*rddy
87  p2=ddx*rddy
88  p3=rddx*ddy
89  p4=ddx*ddy
90
91
92  ! Loop over 2 time steps
93  !***********************
94
95!  do m=1,2
96    indexh=iwftouse
97
98    y1(1)=p1*yy1(ix ,jy ,level,indexh) &
99         + p2*yy1(ixp,jy ,level,indexh) &
100         + p3*yy1(ix ,jyp,level,indexh) &
101         + p4*yy1(ixp,jyp,level,indexh)
102    y2(1)=p1*yy2(ix ,jy ,level,indexh) &
103         + p2*yy2(ixp,jy ,level,indexh) &
104         + p3*yy2(ix ,jyp,level,indexh) &
105         + p4*yy2(ixp,jyp,level,indexh)
106    y3(1)=p1*yy3(ix ,jy ,level,indexh) &
107         + p2*yy3(ixp,jy ,level,indexh) &
108         + p3*yy3(ix ,jyp,level,indexh) &
109         + p4*yy3(ixp,jyp,level,indexh)
110!  end do
111
112
113  !************************************
114  ! 2.) Temporal interpolation (linear) - skip to be consistent with clouds
115  !************************************
116
117!  dt1=real(itime-itime1)
118!  dt2=real(itime2-itime)
119!  dt=dt1+dt2
120
121!  yint1=(y1(1)*dt2+y1(2)*dt1)/dt
122!  yint2=(y2(1)*dt2+y2(2)*dt1)/dt
123!  yint3=(y3(1)*dt2+y3(2)*dt1)/dt
124
125   yint1=y1(1)
126   yint2=y2(1)
127   yint3=y3(1)
128
129end subroutine interpol_rain
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG