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