source: flexpart.git/src/initialize.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: 7.6 KB
Line 
1subroutine initialize(itime,ldt,up,vp,wp, &
2       usigold,vsigold,wsigold,xt,yt,zt,icbt)
3  !                        i    i   o  o  o
4  !        o       o       o    i  i  i   o
5  !*****************************************************************************
6  !                                                                            *
7  !  Calculation of trajectories utilizing a zero-acceleration scheme. The time*
8  !  step is determined by the Courant-Friedrichs-Lewy (CFL) criterion. This   *
9  !  means that the time step must be so small that the displacement within    *
10  !  this time step is smaller than 1 grid distance. Additionally, a temporal  *
11  !  CFL criterion is introduced: the time step must be smaller than the time  *
12  !  interval of the wind fields used for interpolation.                       *
13  !  For random walk simulations, these are the only time step criteria.       *
14  !  For the other options, the time step is also limited by the Lagrangian    *
15  !  time scale.                                                               *
16  !                                                                            *
17  !     Author: A. Stohl                                                       *
18  !                                                                            *
19  !     16 December 1997                                                       *
20  !                                                                            *
21  !  Literature:                                                               *
22  !                                                                            *
23  !*****************************************************************************
24  !                                                                            *
25  ! Variables:                                                                 *
26  ! h [m]              Mixing height                                           *
27  ! lwindinterv [s]    time interval between two wind fields                   *
28  ! itime [s]          current temporal position                               *
29  ! ldt [s]            Suggested time step for next integration                *
30  ! ladvance [s]       Total integration time period                           *
31  ! rannumb(maxrand)   normally distributed random variables                   *
32  ! up,vp,wp           random velocities due to turbulence                     *
33  ! usig,vsig,wsig     uncertainties of wind velocities due to interpolation   *
34  ! usigold,vsigold,wsigold  like usig, etc., but for the last time step       *
35  ! xt,yt,zt           Next time step's spatial position of trajectory         *
36  !                                                                            *
37  !                                                                            *
38  ! Constants:                                                                 *
39  ! cfl                factor, by which the time step has to be smaller than   *
40  !                    the spatial CFL-criterion                               *
41  ! cflt               factor, by which the time step has to be smaller than   *
42  !                    the temporal CFL-criterion                              *
43  !                                                                            *
44  !*****************************************************************************
45
46  use par_mod
47  use com_mod
48  use interpol_mod
49  use hanna_mod
50  use random_mod, only: ran3
51
52  implicit none
53
54  integer :: itime
55  integer :: ldt,nrand
56  integer(kind=2) :: icbt
57  real :: zt,dz,dz1,dz2,up,vp,wp,usigold,vsigold,wsigold
58  real(kind=dp) :: xt,yt
59  save idummy
60
61  integer :: idummy = -7
62
63  icbt=1           ! initialize particle to "no reflection"
64
65  nrand=int(ran3(idummy)*real(maxrand-1))+1
66
67
68  !******************************
69  ! 2. Interpolate necessary data
70  !******************************
71
72  ! Compute maximum mixing height around particle position
73  !*******************************************************
74
75  ix=int(xt)
76  jy=int(yt)
77  ixp=ix+1
78  jyp=jy+1
79
80  h=max(hmix(ix ,jy,1,memind(1)), &
81       hmix(ixp,jy ,1,memind(1)), &
82       hmix(ix ,jyp,1,memind(1)), &
83       hmix(ixp,jyp,1,memind(1)), &
84       hmix(ix ,jy ,1,memind(2)), &
85       hmix(ixp,jy ,1,memind(2)), &
86       hmix(ix ,jyp,1,memind(2)), &
87       hmix(ixp,jyp,1,memind(2)))
88
89  zeta=zt/h
90
91
92  !*************************************************************
93  ! If particle is in the PBL, interpolate once and then make a
94  ! time loop until end of interval is reached
95  !*************************************************************
96
97  if (zeta.le.1.) then
98
99    call interpol_all(itime,real(xt),real(yt),zt)
100
101
102  ! Vertical interpolation of u,v,w,rho and drhodz
103  !***********************************************
104
105  ! Vertical distance to the level below and above current position
106  ! both in terms of (u,v) and (w) fields
107  !****************************************************************
108
109    dz1=zt-height(indz)
110    dz2=height(indzp)-zt
111    dz=1./(dz1+dz2)
112
113    u=(dz1*uprof(indzp)+dz2*uprof(indz))*dz
114    v=(dz1*vprof(indzp)+dz2*vprof(indz))*dz
115    w=(dz1*wprof(indzp)+dz2*wprof(indz))*dz
116
117
118  ! Compute the turbulent disturbances
119
120  ! Determine the sigmas and the timescales
121  !****************************************
122
123    if (turbswitch) then
124      call hanna(zt)
125    else
126      call hanna1(zt)
127    endif
128
129
130  ! Determine the new diffusivity velocities
131  !*****************************************
132
133    if (nrand+2.gt.maxrand) nrand=1
134    up=rannumb(nrand)*sigu
135    vp=rannumb(nrand+1)*sigv
136    wp=rannumb(nrand+2)
137    if (.not.turbswitch) then     ! modified by mc
138      wp=wp*sigw
139    else if (cblflag.eq.1) then   ! modified by mc
140      if(-h/ol.gt.5) then
141!if (ol.lt.0.) then
142!if (ol.gt.0.) then !by mc : only for test correct is lt.0
143        call initialize_cbl_vel(idummy,zt,ust,wst,h,sigw,wp,ol)
144      else
145        wp=wp*sigw
146      end if
147    end if
148
149
150  ! Determine time step for next integration
151  !*****************************************
152
153    if (turbswitch) then
154      ldt=int(min(tlw,h/max(2.*abs(wp*sigw),1.e-5), &
155           0.5/abs(dsigwdz),600.)*ctl)
156    else
157      ldt=int(min(tlw,h/max(2.*abs(wp),1.e-5),600.)*ctl)
158    endif
159    ldt=max(ldt,mintime)
160
161
162    usig=(usigprof(indzp)+usigprof(indz))/2.
163    vsig=(vsigprof(indzp)+vsigprof(indz))/2.
164    wsig=(wsigprof(indzp)+wsigprof(indz))/2.
165
166  else
167
168
169
170  !**********************************************************
171  ! For all particles that are outside the PBL, make a single
172  ! time step. Only horizontal turbulent disturbances are
173  ! calculated. Vertical disturbances are reset.
174  !**********************************************************
175
176
177  ! Interpolate the wind
178  !*********************
179
180    call interpol_wind(itime,real(xt),real(yt),zt)
181
182
183  ! Compute everything for above the PBL
184
185  ! Assume constant turbulent perturbations
186  !****************************************
187
188    ldt=abs(lsynctime)
189
190    if (nrand+1.gt.maxrand) nrand=1
191    up=rannumb(nrand)*0.3
192    vp=rannumb(nrand+1)*0.3
193    nrand=nrand+2
194    wp=0.
195    sigw=0.
196
197  endif
198
199  !****************************************************************
200  ! Add mesoscale random disturbances
201  ! This is done only once for the whole lsynctime interval to save
202  ! computation time
203  !****************************************************************
204
205
206  ! It is assumed that the average interpolation error is 1/2 sigma
207  ! of the surrounding points, autocorrelation time constant is
208  ! 1/2 of time interval between wind fields
209  !****************************************************************
210
211  if (nrand+2.gt.maxrand) nrand=1
212  usigold=rannumb(nrand)*usig
213  vsigold=rannumb(nrand+1)*vsig
214  wsigold=rannumb(nrand+2)*wsig
215
216end subroutine initialize
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG