source: flexpart.git/src/interpol_rain.f90

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file was 92fab65, checked in by Ignacio Pisso <ip@…>, 4 years ago

add SPDX-License-Identifier to all .f90 files

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