!********************************************************************** ! 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 . * !********************************************************************** program flexpart !***************************************************************************** ! * ! This is the Lagrangian Particle Dispersion Model FLEXPART. * ! The main program manages the reading of model run specifications, etc. * ! All actual computing is done within subroutine timemanager. * ! * ! Author: A. Stohl * ! * ! 18 May 1996 * ! * !***************************************************************************** ! * ! Variables: * ! * ! Constants: * ! * !***************************************************************************** use point_mod use par_mod use com_mod use conv_mod implicit none integer :: i,j,ix,jy,inest integer :: idummy = -320 ! character(len=256) :: pathfile, flexversion, arg2 ! Generate a large number of random numbers !****************************************** do i=1,maxrand-1,2 call gasdev1(idummy,rannumb(i),rannumb(i+1)) end do call gasdev1(idummy,rannumb(maxrand),rannumb(maxrand-1)) ! flexversion='Version 9.1.6.2 (Build 2013-09-04)' !verbosity=0 ! Read the pathnames where input/output files are stored !******************************************************* select case (iargc()) case (2) call getarg(1,arg1) pathfile=arg1 call getarg(2,arg2) if (trim(arg2).eq.'-v') then print*, 'verbose mode: additional information will be displayed' verbosity=1 endif case (1) call getarg(1,arg1) pathfile=arg1 if (trim(arg1).eq.'-i') then write(*,*) trim(flexversion) !add a function with more info stop endif if (trim(arg1).eq.'-v') then print*, 'verbose mode: additional information will be displayed' verbosity=1 write(pathfile,'(a11)') './pathnames' endif case (0) write(pathfile,'(a11)') './pathnames' end select ! Print the GPL License statement !******************************************************* ! print*,'Welcome to FLEXPART Version 9.1 (Build 20121029)' print*,'Welcome to FLEXPART', trim(flexversion) print*,'FLEXPART is free software released under the GNU Genera'// & 'l Public License.' call readpaths(pathfile) if (verbosity.eq.1) then print*,'length(4)',length(4) !count=0,count_rate=1000 CALL SYSTEM_CLOCK(count_clock0, count_rate, count_max) !WRITE(*,*) 'SYSTEM_CLOCK',count, count_rate, count_max WRITE(*,*) 'SYSTEM_CLOCK, count_clock0', count_clock0 WRITE(*,*) 'SYSTEM_CLOCK, count_rate', count_rate WRITE(*,*) 'SYSTEM_CLOCK, count_max', count_max endif ! Read the user specifications for the current model run !******************************************************* call readcommand ! Read the age classes to be used !******************************** call readageclasses ! Read, which wind fields are available within the modelling period !****************************************************************** call readavailable ! Read the model grid specifications, ! both for the mother domain and eventual nests !********************************************** call gridcheck call gridcheck_nests ! Read the output grid specifications !************************************ call readoutgrid if (nested_output.eq.1) call readoutgrid_nest ! Read the receptor points for which extra concentrations are to be calculated !***************************************************************************** call readreceptors if (verbosity.eq.1) then print*,'readreceptors' endif ! Read the physico-chemical species property table !************************************************* !SEC: now only needed SPECIES are read in readreleases.f !call readspecies ! Read the landuse inventory !*************************** if (verbosity.eq.1) then print*,'readlanduse' endif call readlanduse ! Assign fractional cover of landuse classes to each ECMWF grid point !******************************************************************** if (verbosity.eq.1) then print*,'assignland' endif call assignland ! Read the coordinates of the release locations !********************************************** if (verbosity.eq.1) then print*,'readreleases' endif call readreleases ! Read and compute surface resistances to dry deposition of gases !**************************************************************** if (verbosity.eq.1) then print*,'readdepo' endif call readdepo ! Convert the release point coordinates from geografical to grid coordinates !*************************************************************************** if (verbosity.eq.1) then print*,'coordtrafo' endif call coordtrafo ! Initialize all particles to non-existent !***************************************** if (verbosity.eq.1) then print*,'Initialize all particles to non-existent' endif do j=1,maxpart itra1(j)=-999999999 end do ! For continuation of previous run, read in particle positions !************************************************************* if (ipin.eq.1) then if (verbosity.eq.1) then print*,'readpartpositions' endif call readpartpositions else numpart=0 numparticlecount=0 endif ! Calculate volume, surface area, etc., of all output grid cells ! Allocate fluxes and OHfield if necessary !*************************************************************** call outgrid_init if (nested_output.eq.1) call outgrid_init_nest ! Read the OH field !****************** if (OHREA.eqv..TRUE.) & call readOHfield ! Write basic information on the simulation to a file "header" ! and open files that are to be kept open throughout the simulation !****************************************************************** call writeheader call writeheader_txt !if (nested_output.eq.1) call writeheader_nest if (nested_output.eq.1.and.surf_only.ne.1) call writeheader_nest if (nested_output.eq.1.and.surf_only.eq.1) call writeheader_nest_surf if (nested_output.ne.1.and.surf_only.eq.1) call writeheader_surf !open(unitdates,file=path(2)(1:length(2))//'dates') call openreceptors if ((iout.eq.4).or.(iout.eq.5)) call openouttraj ! Releases can only start and end at discrete times (multiples of lsynctime) !*************************************************************************** do i=1,numpoint ireleasestart(i)=nint(real(ireleasestart(i))/ & real(lsynctime))*lsynctime ireleaseend(i)=nint(real(ireleaseend(i))/ & real(lsynctime))*lsynctime end do ! Initialize cloud-base mass fluxes for the convection scheme !************************************************************ do jy=0,nymin1 do ix=0,nxmin1 cbaseflux(ix,jy)=0. end do end do do inest=1,numbnests do jy=0,nyn(inest)-1 do ix=0,nxn(inest)-1 cbasefluxn(ix,jy,inest)=0. end do end do end do ! Calculate particle trajectories !******************************** if (verbosity.eq.1) then print*,'call timemanager' CALL SYSTEM_CLOCK(count_clock, count_rate, count_max) WRITE(*,*) 'SYSTEM_CLOCK',(count_clock - count_clock0)/count_rate !, count_rate, count_max endif call timemanager write(*,*) 'CONGRATULATIONS: YOU HAVE SUCCESSFULLY COMPLETED A FLE& &XPART MODEL RUN!' end program flexpart