[e200b7a] | 1 | !********************************************************************** |
---|
| 2 | ! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010 * |
---|
| 3 | ! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa, * |
---|
| 4 | ! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann * |
---|
| 5 | ! * |
---|
| 6 | ! This file is part of FLEXPART. * |
---|
| 7 | ! * |
---|
| 8 | ! FLEXPART is free software: you can redistribute it and/or modify * |
---|
| 9 | ! it under the terms of the GNU General Public License as published by* |
---|
| 10 | ! the Free Software Foundation, either version 3 of the License, or * |
---|
| 11 | ! (at your option) any later version. * |
---|
| 12 | ! * |
---|
| 13 | ! FLEXPART is distributed in the hope that it will be useful, * |
---|
| 14 | ! but WITHOUT ANY WARRANTY; without even the implied warranty of * |
---|
| 15 | ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * |
---|
| 16 | ! GNU General Public License for more details. * |
---|
| 17 | ! * |
---|
| 18 | ! You should have received a copy of the GNU General Public License * |
---|
| 19 | ! along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. * |
---|
| 20 | !********************************************************************** |
---|
| 21 | |
---|
| 22 | program flexpart |
---|
| 23 | |
---|
| 24 | !***************************************************************************** |
---|
| 25 | ! * |
---|
| 26 | ! This is the Lagrangian Particle Dispersion Model FLEXPART. * |
---|
| 27 | ! The main program manages the reading of model run specifications, etc. * |
---|
| 28 | ! All actual computing is done within subroutine timemanager. * |
---|
| 29 | ! * |
---|
| 30 | ! Author: A. Stohl * |
---|
| 31 | ! * |
---|
| 32 | ! 18 May 1996 * |
---|
| 33 | ! * |
---|
| 34 | !***************************************************************************** |
---|
| 35 | ! * |
---|
| 36 | ! Variables: * |
---|
| 37 | ! * |
---|
| 38 | ! Constants: * |
---|
| 39 | ! * |
---|
| 40 | !***************************************************************************** |
---|
| 41 | |
---|
| 42 | use point_mod |
---|
| 43 | use par_mod |
---|
| 44 | use com_mod |
---|
| 45 | use conv_mod |
---|
[8a65cb0] | 46 | use netcdf_output_mod, only: writeheader_netcdf |
---|
| 47 | use random_mod, only: gasdev1 |
---|
[e200b7a] | 48 | |
---|
| 49 | implicit none |
---|
| 50 | |
---|
| 51 | integer :: i,j,ix,jy,inest |
---|
| 52 | integer :: idummy = -320 |
---|
[f13406c] | 53 | character(len=256) :: inline_options !pathfile, flexversion, arg2 |
---|
[8a65cb0] | 54 | |
---|
| 55 | |
---|
| 56 | ! Initialize arrays in com_mod |
---|
| 57 | !***************************** |
---|
[b0434e1] | 58 | call com_mod_allocate_part(maxpart) |
---|
| 59 | |
---|
[8a65cb0] | 60 | |
---|
[e200b7a] | 61 | ! Generate a large number of random numbers |
---|
| 62 | !****************************************** |
---|
| 63 | |
---|
| 64 | do i=1,maxrand-1,2 |
---|
| 65 | call gasdev1(idummy,rannumb(i),rannumb(i+1)) |
---|
| 66 | end do |
---|
| 67 | call gasdev1(idummy,rannumb(maxrand),rannumb(maxrand-1)) |
---|
| 68 | |
---|
[b4d29ce] | 69 | ! FLEXPART version string |
---|
[fddc6ec] | 70 | flexversion_major = '10' ! Major version number, also used for species file names |
---|
| 71 | flexversion='Version '//trim(flexversion_major)//'.0beta (2015-05-01)' |
---|
[8a65cb0] | 72 | verbosity=0 |
---|
| 73 | |
---|
[f13406c] | 74 | ! Read the pathnames where input/output files are stored |
---|
| 75 | !******************************************************* |
---|
| 76 | |
---|
[8a65cb0] | 77 | inline_options='none' |
---|
[f13406c] | 78 | select case (iargc()) |
---|
[8a65cb0] | 79 | case (2) |
---|
[f13406c] | 80 | call getarg(1,arg1) |
---|
| 81 | pathfile=arg1 |
---|
| 82 | call getarg(2,arg2) |
---|
| 83 | inline_options=arg2 |
---|
[8a65cb0] | 84 | case (1) |
---|
[f13406c] | 85 | call getarg(1,arg1) |
---|
| 86 | pathfile=arg1 |
---|
| 87 | if (arg1(1:1).eq.'-') then |
---|
[b4d29ce] | 88 | write(pathfile,'(a11)') './pathnames' |
---|
| 89 | inline_options=arg1 |
---|
[f13406c] | 90 | endif |
---|
[8a65cb0] | 91 | case (0) |
---|
[f13406c] | 92 | write(pathfile,'(a11)') './pathnames' |
---|
| 93 | end select |
---|
| 94 | |
---|
[e200b7a] | 95 | ! Print the GPL License statement |
---|
| 96 | !******************************************************* |
---|
[4fbe7a5] | 97 | print*,'Welcome to FLEXPART ', trim(flexversion) |
---|
[b4d29ce] | 98 | print*,'FLEXPART is free software released under the GNU General Public License.' |
---|
[414a5e5] | 99 | |
---|
[8a65cb0] | 100 | if (inline_options(1:1).eq.'-') then |
---|
| 101 | if (trim(inline_options).eq.'-v'.or.trim(inline_options).eq.'-v1') then |
---|
| 102 | print*, 'Verbose mode 1: display detailed information during run' |
---|
[b4d29ce] | 103 | verbosity=1 |
---|
[8a65cb0] | 104 | endif |
---|
| 105 | if (trim(inline_options).eq.'-v2') then |
---|
| 106 | print*, 'Verbose mode 2: display more detailed information during run' |
---|
[b4d29ce] | 107 | verbosity=2 |
---|
| 108 | endif |
---|
[8a65cb0] | 109 | if (trim(inline_options).eq.'-i') then |
---|
| 110 | print*, 'Info mode: provide detailed run specific information and stop' |
---|
[b4d29ce] | 111 | verbosity=1 |
---|
| 112 | info_flag=1 |
---|
| 113 | endif |
---|
[8a65cb0] | 114 | if (trim(inline_options).eq.'-i2') then |
---|
| 115 | print*, 'Info mode: provide more detailed run specific information and stop' |
---|
| 116 | verbosity=2 |
---|
| 117 | info_flag=1 |
---|
[b4d29ce] | 118 | endif |
---|
| 119 | endif |
---|
| 120 | |
---|
[f13406c] | 121 | if (verbosity.gt.0) then |
---|
[b4d29ce] | 122 | write(*,*) 'call readpaths' |
---|
[f13406c] | 123 | endif |
---|
| 124 | call readpaths(pathfile) |
---|
| 125 | |
---|
[8a65cb0] | 126 | if (verbosity.gt.1) then !show clock info |
---|
| 127 | !print*,'length(4)',length(4) |
---|
[f13406c] | 128 | !count=0,count_rate=1000 |
---|
[8a65cb0] | 129 | CALL SYSTEM_CLOCK(count_clock0, count_rate, count_max) |
---|
[f13406c] | 130 | !WRITE(*,*) 'SYSTEM_CLOCK',count, count_rate, count_max |
---|
| 131 | !WRITE(*,*) 'SYSTEM_CLOCK, count_clock0', count_clock0 |
---|
| 132 | !WRITE(*,*) 'SYSTEM_CLOCK, count_rate', count_rate |
---|
| 133 | !WRITE(*,*) 'SYSTEM_CLOCK, count_max', count_max |
---|
[8a65cb0] | 134 | endif |
---|
[e200b7a] | 135 | |
---|
| 136 | ! Read the user specifications for the current model run |
---|
| 137 | !******************************************************* |
---|
| 138 | |
---|
[f13406c] | 139 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 140 | write(*,*) 'call readcommand' |
---|
[f13406c] | 141 | endif |
---|
[e200b7a] | 142 | call readcommand |
---|
[f13406c] | 143 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 144 | write(*,*) ' ldirect=', ldirect |
---|
| 145 | write(*,*) ' ibdate,ibtime=',ibdate,ibtime |
---|
[b4d29ce] | 146 | write(*,*) ' iedate,ietime=', iedate,ietime |
---|
[8a65cb0] | 147 | if (verbosity.gt.1) then |
---|
[b4d29ce] | 148 | CALL SYSTEM_CLOCK(count_clock, count_rate, count_max) |
---|
| 149 | write(*,*) 'SYSTEM_CLOCK',(count_clock - count_clock0)/real(count_rate) !, count_rate, count_max |
---|
| 150 | endif |
---|
[8a65cb0] | 151 | endif |
---|
[e200b7a] | 152 | |
---|
| 153 | ! Read the age classes to be used |
---|
| 154 | !******************************** |
---|
[f13406c] | 155 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 156 | write(*,*) 'call readageclasses' |
---|
[f13406c] | 157 | endif |
---|
[e200b7a] | 158 | call readageclasses |
---|
| 159 | |
---|
[8a65cb0] | 160 | if (verbosity.gt.1) then |
---|
[b4d29ce] | 161 | CALL SYSTEM_CLOCK(count_clock, count_rate, count_max) |
---|
| 162 | write(*,*) 'SYSTEM_CLOCK',(count_clock - count_clock0)/real(count_rate) !, count_rate, count_max |
---|
[f13406c] | 163 | endif |
---|
[e200b7a] | 164 | |
---|
| 165 | ! Read, which wind fields are available within the modelling period |
---|
| 166 | !****************************************************************** |
---|
| 167 | |
---|
[f13406c] | 168 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 169 | write(*,*) 'call readavailable' |
---|
[f13406c] | 170 | endif |
---|
[e200b7a] | 171 | call readavailable |
---|
| 172 | |
---|
[b0434e1] | 173 | ! If nested wind fields are used, allocate arrays |
---|
| 174 | !************************************************ |
---|
| 175 | |
---|
| 176 | if (verbosity.gt.0) then |
---|
| 177 | write(*,*) 'call com_mod_allocate_nests' |
---|
| 178 | endif |
---|
| 179 | call com_mod_allocate_nests |
---|
| 180 | |
---|
[e200b7a] | 181 | ! Read the model grid specifications, |
---|
| 182 | ! both for the mother domain and eventual nests |
---|
| 183 | !********************************************** |
---|
[f13406c] | 184 | |
---|
| 185 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 186 | write(*,*) 'call gridcheck' |
---|
[f13406c] | 187 | endif |
---|
[8a65cb0] | 188 | |
---|
[e200b7a] | 189 | call gridcheck |
---|
[f13406c] | 190 | |
---|
[8a65cb0] | 191 | if (verbosity.gt.1) then |
---|
[b4d29ce] | 192 | CALL SYSTEM_CLOCK(count_clock, count_rate, count_max) |
---|
| 193 | write(*,*) 'SYSTEM_CLOCK',(count_clock - count_clock0)/real(count_rate) !, count_rate, count_max |
---|
[f13406c] | 194 | endif |
---|
| 195 | |
---|
| 196 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 197 | write(*,*) 'call gridcheck_nests' |
---|
[f13406c] | 198 | endif |
---|
[e200b7a] | 199 | call gridcheck_nests |
---|
| 200 | |
---|
| 201 | ! Read the output grid specifications |
---|
| 202 | !************************************ |
---|
| 203 | |
---|
[f13406c] | 204 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 205 | write(*,*) 'call readoutgrid' |
---|
[f13406c] | 206 | endif |
---|
| 207 | |
---|
[e200b7a] | 208 | call readoutgrid |
---|
| 209 | |
---|
[f13406c] | 210 | if (nested_output.eq.1) then |
---|
[b4d29ce] | 211 | call readoutgrid_nest |
---|
[f13406c] | 212 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 213 | write(*,*) '# readoutgrid_nest' |
---|
[f13406c] | 214 | endif |
---|
| 215 | endif |
---|
[e200b7a] | 216 | |
---|
| 217 | ! Read the receptor points for which extra concentrations are to be calculated |
---|
| 218 | !***************************************************************************** |
---|
| 219 | |
---|
[f13406c] | 220 | if (verbosity.eq.1) then |
---|
[8a65cb0] | 221 | print*,'call readreceptors' |
---|
[f13406c] | 222 | endif |
---|
[e200b7a] | 223 | call readreceptors |
---|
| 224 | |
---|
| 225 | ! Read the physico-chemical species property table |
---|
| 226 | !************************************************* |
---|
| 227 | !SEC: now only needed SPECIES are read in readreleases.f |
---|
| 228 | !call readspecies |
---|
| 229 | |
---|
| 230 | |
---|
| 231 | ! Read the landuse inventory |
---|
| 232 | !*************************** |
---|
| 233 | |
---|
[f13406c] | 234 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 235 | print*,'call readlanduse' |
---|
[f13406c] | 236 | endif |
---|
[e200b7a] | 237 | call readlanduse |
---|
| 238 | |
---|
| 239 | ! Assign fractional cover of landuse classes to each ECMWF grid point |
---|
| 240 | !******************************************************************** |
---|
| 241 | |
---|
[f13406c] | 242 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 243 | print*,'call assignland' |
---|
[f13406c] | 244 | endif |
---|
[e200b7a] | 245 | call assignland |
---|
| 246 | |
---|
| 247 | ! Read the coordinates of the release locations |
---|
| 248 | !********************************************** |
---|
| 249 | |
---|
[f13406c] | 250 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 251 | print*,'call readreleases' |
---|
[f13406c] | 252 | endif |
---|
[e200b7a] | 253 | call readreleases |
---|
| 254 | |
---|
| 255 | ! Read and compute surface resistances to dry deposition of gases |
---|
| 256 | !**************************************************************** |
---|
| 257 | |
---|
[f13406c] | 258 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 259 | print*,'call readdepo' |
---|
[f13406c] | 260 | endif |
---|
[e200b7a] | 261 | call readdepo |
---|
| 262 | |
---|
| 263 | ! Convert the release point coordinates from geografical to grid coordinates |
---|
| 264 | !*************************************************************************** |
---|
| 265 | |
---|
[f13406c] | 266 | call coordtrafo |
---|
| 267 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 268 | print*,'call coordtrafo' |
---|
[f13406c] | 269 | endif |
---|
[e200b7a] | 270 | |
---|
| 271 | ! Initialize all particles to non-existent |
---|
| 272 | !***************************************** |
---|
| 273 | |
---|
[f13406c] | 274 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 275 | print*,'Initialize all particles to non-existent' |
---|
[f13406c] | 276 | endif |
---|
[e200b7a] | 277 | do j=1,maxpart |
---|
| 278 | itra1(j)=-999999999 |
---|
| 279 | end do |
---|
| 280 | |
---|
| 281 | ! For continuation of previous run, read in particle positions |
---|
| 282 | !************************************************************* |
---|
| 283 | |
---|
| 284 | if (ipin.eq.1) then |
---|
[f13406c] | 285 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 286 | print*,'call readpartpositions' |
---|
[f13406c] | 287 | endif |
---|
[e200b7a] | 288 | call readpartpositions |
---|
| 289 | else |
---|
[f13406c] | 290 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 291 | print*,'numpart=0, numparticlecount=0' |
---|
[f13406c] | 292 | endif |
---|
[e200b7a] | 293 | numpart=0 |
---|
| 294 | numparticlecount=0 |
---|
| 295 | endif |
---|
| 296 | |
---|
| 297 | ! Calculate volume, surface area, etc., of all output grid cells |
---|
| 298 | ! Allocate fluxes and OHfield if necessary |
---|
| 299 | !*************************************************************** |
---|
| 300 | |
---|
[f13406c] | 301 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 302 | print*,'call outgrid_init' |
---|
[f13406c] | 303 | endif |
---|
[e200b7a] | 304 | call outgrid_init |
---|
| 305 | if (nested_output.eq.1) call outgrid_init_nest |
---|
| 306 | |
---|
| 307 | ! Read the OH field |
---|
| 308 | !****************** |
---|
| 309 | |
---|
[f13406c] | 310 | if (OHREA.eqv..TRUE.) then |
---|
| 311 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 312 | print*,'call readOHfield' |
---|
[f13406c] | 313 | endif |
---|
[b4d29ce] | 314 | call readOHfield |
---|
[f13406c] | 315 | endif |
---|
[e200b7a] | 316 | |
---|
| 317 | ! Write basic information on the simulation to a file "header" |
---|
| 318 | ! and open files that are to be kept open throughout the simulation |
---|
| 319 | !****************************************************************** |
---|
| 320 | |
---|
[8a65cb0] | 321 | if (lnetcdfout.eq.1) then |
---|
| 322 | call writeheader_netcdf(lnest=.false.) |
---|
| 323 | else |
---|
| 324 | call writeheader |
---|
| 325 | end if |
---|
| 326 | |
---|
| 327 | if (nested_output.eq.1) then |
---|
| 328 | if (lnetcdfout.eq.1) then |
---|
| 329 | call writeheader_netcdf(lnest=.true.) |
---|
| 330 | else |
---|
| 331 | call writeheader_nest |
---|
| 332 | endif |
---|
| 333 | endif |
---|
| 334 | |
---|
[f13406c] | 335 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 336 | print*,'call writeheader' |
---|
[f13406c] | 337 | endif |
---|
| 338 | |
---|
[e200b7a] | 339 | call writeheader |
---|
[8a65cb0] | 340 | ! FLEXPART 9.2 ticket ?? write header in ASCII format |
---|
[f13406c] | 341 | call writeheader_txt |
---|
| 342 | !if (nested_output.eq.1) call writeheader_nest |
---|
| 343 | if (nested_output.eq.1.and.surf_only.ne.1) call writeheader_nest |
---|
| 344 | if (nested_output.eq.1.and.surf_only.eq.1) call writeheader_nest_surf |
---|
| 345 | if (nested_output.ne.1.and.surf_only.eq.1) call writeheader_surf |
---|
| 346 | |
---|
| 347 | !open(unitdates,file=path(2)(1:length(2))//'dates') |
---|
| 348 | |
---|
| 349 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 350 | print*,'call openreceptors' |
---|
[f13406c] | 351 | endif |
---|
[e200b7a] | 352 | call openreceptors |
---|
| 353 | if ((iout.eq.4).or.(iout.eq.5)) call openouttraj |
---|
| 354 | |
---|
| 355 | ! Releases can only start and end at discrete times (multiples of lsynctime) |
---|
| 356 | !*************************************************************************** |
---|
| 357 | |
---|
[f13406c] | 358 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 359 | print*,'discretize release times' |
---|
[f13406c] | 360 | endif |
---|
[e200b7a] | 361 | do i=1,numpoint |
---|
[b4d29ce] | 362 | ireleasestart(i)=nint(real(ireleasestart(i))/real(lsynctime))*lsynctime |
---|
| 363 | ireleaseend(i)=nint(real(ireleaseend(i))/real(lsynctime))*lsynctime |
---|
[e200b7a] | 364 | end do |
---|
| 365 | |
---|
| 366 | ! Initialize cloud-base mass fluxes for the convection scheme |
---|
| 367 | !************************************************************ |
---|
| 368 | |
---|
[f13406c] | 369 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 370 | print*,'Initialize cloud-base mass fluxes for the convection scheme' |
---|
[f13406c] | 371 | endif |
---|
| 372 | |
---|
[e200b7a] | 373 | do jy=0,nymin1 |
---|
| 374 | do ix=0,nxmin1 |
---|
| 375 | cbaseflux(ix,jy)=0. |
---|
| 376 | end do |
---|
| 377 | end do |
---|
| 378 | do inest=1,numbnests |
---|
| 379 | do jy=0,nyn(inest)-1 |
---|
| 380 | do ix=0,nxn(inest)-1 |
---|
| 381 | cbasefluxn(ix,jy,inest)=0. |
---|
| 382 | end do |
---|
| 383 | end do |
---|
| 384 | end do |
---|
| 385 | |
---|
| 386 | ! Calculate particle trajectories |
---|
| 387 | !******************************** |
---|
| 388 | |
---|
[f13406c] | 389 | if (verbosity.gt.0) then |
---|
[8a65cb0] | 390 | if (verbosity.gt.1) then |
---|
| 391 | CALL SYSTEM_CLOCK(count_clock, count_rate, count_max) |
---|
| 392 | write(*,*) 'SYSTEM_CLOCK',(count_clock - count_clock0)/real(count_rate) !, count_rate, count_max |
---|
| 393 | endif |
---|
| 394 | if (info_flag.eq.1) then |
---|
| 395 | print*, 'info only mode (stop)' |
---|
| 396 | stop |
---|
| 397 | endif |
---|
| 398 | print*,'call timemanager' |
---|
[f13406c] | 399 | endif |
---|
| 400 | |
---|
[e200b7a] | 401 | call timemanager |
---|
| 402 | |
---|
[8a65cb0] | 403 | ! NIK 16.02.2005 |
---|
| 404 | write(*,*) '**********************************************' |
---|
| 405 | write(*,*) 'Total number of occurences of below-cloud scavenging', tot_blc_count |
---|
| 406 | write(*,*) 'Total number of occurences of in-cloud scavenging', tot_inc_count |
---|
| 407 | write(*,*) '**********************************************' |
---|
| 408 | |
---|
| 409 | write(*,*) 'CONGRATULATIONS: YOU HAVE SUCCESSFULLY COMPLETED A FLE& |
---|
| 410 | &XPART MODEL RUN!' |
---|
[e200b7a] | 411 | |
---|
| 412 | end program flexpart |
---|