initialize.f90 9.02 KB
Newer Older
Matthias Langer's avatar
 
Matthias Langer committed
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
!**********************************************************************
! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010         *
! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa,             *
! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann   *
!                                                                     *
! This file is part of FLEXPART.                                      *
!                                                                     *
! FLEXPART is free software: you can redistribute it and/or modify    *
! it under the terms of the GNU General Public License as published by*
! the Free Software Foundation, either version 3 of the License, or   *
! (at your option) any later version.                                 *
!                                                                     *
! FLEXPART is distributed in the hope that it will be useful,         *
! but WITHOUT ANY WARRANTY; without even the implied warranty of      *
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
! GNU General Public License for more details.                        *
!                                                                     *
! You should have received a copy of the GNU General Public License   *
! along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
!**********************************************************************

subroutine initialize(itime,ldt,up,vp,wp, &
       usigold,vsigold,wsigold,xt,yt,zt,icbt)
  !                        i    i   o  o  o
  !        o       o       o    i  i  i   o
  !*****************************************************************************
  !                                                                            *
  !  Calculation of trajectories utilizing a zero-acceleration scheme. The time*
  !  step is determined by the Courant-Friedrichs-Lewy (CFL) criterion. This   *
  !  means that the time step must be so small that the displacement within    *
  !  this time step is smaller than 1 grid distance. Additionally, a temporal  *
  !  CFL criterion is introduced: the time step must be smaller than the time  *
  !  interval of the wind fields used for interpolation.                       *
  !  For random walk simulations, these are the only time step criteria.       *
  !  For the other options, the time step is also limited by the Lagrangian    *
  !  time scale.                                                               *
  !                                                                            *
  !     Author: A. Stohl                                                       *
  !                                                                            *
  !     16 December 1997                                                       *
  !                                                                            *
  !  Literature:                                                               *
  !                                                                            *
  !*****************************************************************************
  !                                                                            *
  ! Variables:                                                                 *
  ! h [m]              Mixing height                                           *
  ! lwindinterv [s]    time interval between two wind fields                   *
  ! itime [s]          current temporal position                               *
  ! ldt [s]            Suggested time step for next integration                *
  ! ladvance [s]       Total integration time period                           *
  ! rannumb(maxrand)   normally distributed random variables                   *
  ! up,vp,wp           random velocities due to turbulence                     *
  ! usig,vsig,wsig     uncertainties of wind velocities due to interpolation   *
  ! usigold,vsigold,wsigold  like usig, etc., but for the last time step       *
  ! xt,yt,zt           Next time step's spatial position of trajectory         *
  !                                                                            *
  !                                                                            *
  ! Constants:                                                                 *
  ! cfl                factor, by which the time step has to be smaller than   *
  !                    the spatial CFL-criterion                               *
  ! cflt               factor, by which the time step has to be smaller than   *
  !                    the temporal CFL-criterion                              *
  !                                                                            *
  !*****************************************************************************

  use par_mod
  use com_mod
  use interpol_mod
  use hanna_mod
71
  use random_mod, only: ran3
Matthias Langer's avatar
 
Matthias Langer committed
72
73
74
75
76
77

  implicit none

  integer :: itime
  integer :: ldt,nrand
  integer(kind=2) :: icbt
78
  real :: zt,dz,dz1,dz2,up,vp,wp,usigold,vsigold,wsigold
Matthias Langer's avatar
 
Matthias Langer committed
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
  real(kind=dp) :: xt,yt
  save idummy

  integer :: idummy = -7

  icbt=1           ! initialize particle to "no reflection"

  nrand=int(ran3(idummy)*real(maxrand-1))+1


  !******************************
  ! 2. Interpolate necessary data
  !******************************

  ! Compute maximum mixing height around particle position
  !*******************************************************

  ix=int(xt)
  jy=int(yt)
  ixp=ix+1
  jyp=jy+1

101
  h=max(hmix(ix ,jy,1,memind(1)), &
Matthias Langer's avatar
 
Matthias Langer committed
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
       hmix(ixp,jy ,1,memind(1)), &
       hmix(ix ,jyp,1,memind(1)), &
       hmix(ixp,jyp,1,memind(1)), &
       hmix(ix ,jy ,1,memind(2)), &
       hmix(ixp,jy ,1,memind(2)), &
       hmix(ix ,jyp,1,memind(2)), &
       hmix(ixp,jyp,1,memind(2)))

  zeta=zt/h


  !*************************************************************
  ! If particle is in the PBL, interpolate once and then make a
  ! time loop until end of interval is reached
  !*************************************************************

  if (zeta.le.1.) then

    call interpol_all(itime,real(xt),real(yt),zt)


  ! Vertical interpolation of u,v,w,rho and drhodz
  !***********************************************

  ! Vertical distance to the level below and above current position
  ! both in terms of (u,v) and (w) fields
  !****************************************************************

    dz1=zt-height(indz)
    dz2=height(indzp)-zt
    dz=1./(dz1+dz2)

    u=(dz1*uprof(indzp)+dz2*uprof(indz))*dz
    v=(dz1*vprof(indzp)+dz2*vprof(indz))*dz
    w=(dz1*wprof(indzp)+dz2*wprof(indz))*dz


  ! Compute the turbulent disturbances

  ! Determine the sigmas and the timescales
  !****************************************

    if (turbswitch) then
      call hanna(zt)
    else
      call hanna1(zt)
    endif


  ! Determine the new diffusivity velocities
  !*****************************************

    if (nrand+2.gt.maxrand) nrand=1
    up=rannumb(nrand)*sigu
    vp=rannumb(nrand+1)*sigv
    wp=rannumb(nrand+2)
158
    if (.not.turbswitch) then     ! modified by mc
159
      wp=wp*sigw
pesei's avatar
pesei committed
160
    else if (iflagcbl.eq.1) then   ! modified by mc
161
162
163
164
165
166
167
      if(-h/ol.gt.5) then
!if (ol.lt.0.) then
!if (ol.gt.0.) then !by mc : only for test correct is lt.0
        call initialize_cbl_vel(idummy,zt,ust,wst,h,sigw,wp,ol)
      else
        wp=wp*sigw
      end if
168
    end if
Matthias Langer's avatar
 
Matthias Langer committed
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237


  ! Determine time step for next integration
  !*****************************************

    if (turbswitch) then
      ldt=int(min(tlw,h/max(2.*abs(wp*sigw),1.e-5), &
           0.5/abs(dsigwdz),600.)*ctl)
    else
      ldt=int(min(tlw,h/max(2.*abs(wp),1.e-5),600.)*ctl)
    endif
    ldt=max(ldt,mintime)


    usig=(usigprof(indzp)+usigprof(indz))/2.
    vsig=(vsigprof(indzp)+vsigprof(indz))/2.
    wsig=(wsigprof(indzp)+wsigprof(indz))/2.

  else



  !**********************************************************
  ! For all particles that are outside the PBL, make a single
  ! time step. Only horizontal turbulent disturbances are
  ! calculated. Vertical disturbances are reset.
  !**********************************************************


  ! Interpolate the wind
  !*********************

    call interpol_wind(itime,real(xt),real(yt),zt)


  ! Compute everything for above the PBL

  ! Assume constant turbulent perturbations
  !****************************************

    ldt=abs(lsynctime)

    if (nrand+1.gt.maxrand) nrand=1
    up=rannumb(nrand)*0.3
    vp=rannumb(nrand+1)*0.3
    nrand=nrand+2
    wp=0.
    sigw=0.

  endif

  !****************************************************************
  ! Add mesoscale random disturbances
  ! This is done only once for the whole lsynctime interval to save
  ! computation time
  !****************************************************************


  ! It is assumed that the average interpolation error is 1/2 sigma
  ! of the surrounding points, autocorrelation time constant is
  ! 1/2 of time interval between wind fields
  !****************************************************************

  if (nrand+2.gt.maxrand) nrand=1
  usigold=rannumb(nrand)*usig
  vsigold=rannumb(nrand+1)*vsig
  wsigold=rannumb(nrand+2)*wsig

end subroutine initialize