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