Skip to content
Snippets Groups Projects
Commit 1b3c86a7 authored by Harald Sodemann's avatar Harald Sodemann
Browse files

Initial import

- namelist input for COMMAND
- pathnames optionally as command line argument
- conversion utility from COMMAND to COMMAND namelist


git-svn-id: http://flexpart.flexpart.eu:8088/svn/FlexPart90/branches/command2nml@7 ef8cc7e1-21b7-489e-abab-c1baa636049d
parents
No related branches found
No related tags found
No related merge requests found
COMMAND 0 → 100644
**********************************************
COMMAND Input file for FLEXPART RDF
RDF releases
*********************************************
1
20121001 000000
20121004 120000
10800 tout
10800 tavg
600
999999999
600 SYNC
5.0 CTL
4 IFINE
5 IOUT
2 IPOUT
1 LSUBGRID
1 LCONVECTION
1 LAGESPECTRA
0 IPIN
0 OUTPUTFOREACHRELEASE
0 IFLUX
0 MDOMAINFILL
1 IND_SOURCE
1 IND_RECEPTOR
0 MQUASILAG
0 NESTED_OUTPUT
0 LINIT_COND
&COMMAND
LDIRECT= 1,
IBDATE= 20121001,
IBTIME= 0,
IEDATE= 20121004,
IETIME= 120000,
LOUTSTEP= 10800,
LOUTAVER= 10800,
LOUTSAMPLE= 600,
ITSPLIT= 999999999,
LSYNCTIME= 600,
CTL= 5.0000000 ,
IFINE= 4,
IOUT= 5,
IPOUT= 2,
LSUBGRID= 1,
LCONVECTION= 1,
LAGESPECTRA= 1,
IPIN= 0,
IOUTPUTFOREACHRELEASE= 0,
IFLUX= 0,
MDOMAINFILL= 0,
IND_SOURCE= 1,
IND_RECEPTOR= 1,
MQUASILAG= 0,
NESTED_OUTPUT= 0,
LINIT_COND= 0,
/
Makefile 0 → 100644
SHELL = /bin/bash
MAIN = command2nml
#
FC = gfortran
LIBPATH2 = /usr/local/lib
FFLAGS = -O3 -m64
LDFLAGS = $(FFLAGS) -L$(LIBPATH2) -lm
#
MODOBJS = \
par_mod.o com_mod.o
OBJECTS = \
skplin.o command2nml.o
$(MAIN): $(MODOBJS) $(OBJECTS)
$(FC) *.o -o $(MAIN) $(LDFLAGS)
$(OBJECTS): $(MODOBJS)
%.o: %.f90
$(FC) -c $(FFLAGS) $<
clean:
rm *.o *.mod
README 0 → 100644
command2nml
---------------------
Converts FLEXPART COMMAND files in pre-9.1 version to namelist format
Call sequence:
command2nml COMMAND_old COMMAND_new
HS, 14.05.2013
This diff is collapsed.
program command2nml
!*****************************************************************************
! *
! This program reads the command file in any known format and writes *
! it in namelist format to the output argument
! *
! Author: Harald Sodemann
! 29 Oct 2012 *
!
! Input argument: COMMAND file
! Output argument: COMMAND file in namelist format
! *
!*****************************************************************************
! *
! Variables: *
! bdate beginning date as Julian date *
! ctl factor by which time step must be smaller than *
! Lagrangian time scale *
! ibdate,ibtime beginnning date and time (YYYYMMDD, HHMISS) *
! ideltas [s] modelling period *
! iedate,ietime ending date and time (YYYYMMDD, HHMISS) *
! ifine reduction factor for vertical wind time step *
! outputforeachrel for forward runs it is possible either to create *
! one outputfield or several for each releasepoint *
! iflux switch to turn on (1)/off (0) flux calculations *
! iout 1 for conc. (residence time for backward runs) output,*
! 2 for mixing ratio output, 3 both, 4 for plume *
! trajectory output, 5 = options 1 and 4 *
! ipin 1 continue simulation with dumped particle data, 0 no *
! ipout 0 no particle dump, 1 every output time, 3 only at end*
! itsplit [s] time constant for particle splitting *
! loutaver [s] concentration output is an average over loutaver *
! seconds *
! loutsample [s] average is computed from samples taken every [s] *
! seconds *
! loutstep [s] time interval of concentration output *
! lsynctime [s] synchronisation time interval for all particles *
! lagespectra switch to turn on (1)/off (0) calculation of age *
! spectra *
! lconvection value of either 0 and 1 indicating mixing by *
! convection *
! = 0 .. no convection *
! + 1 .. parameterisation of mixing by subgrid-scale *
! convection = on *
! lsubgrid switch to turn on (1)/off (0) subgrid topography *
! parameterization *
! method method used to compute the particle pseudovelocities *
! mdomainfill 1 use domain-filling option, 0 not, 2 use strat. O3 *
! *
! Constants: *
! 10 unit connected to file COMMAND *
! *
!*****************************************************************************
use com_mod
implicit none
real(kind=dp) :: juldate
character(len=50) :: line
logical :: old
integer :: readerror
character(256) :: infile
character(256) :: outfile
namelist /command/ &
ldirect, &
ibdate,ibtime, &
iedate,ietime, &
loutstep, &
loutaver, &
loutsample, &
itsplit, &
lsynctime, &
ctl, &
ifine, &
iout, &
ipout, &
lsubgrid, &
lconvection, &
lagespectra, &
ipin, &
ioutputforeachrelease, &
iflux, &
mdomainfill, &
ind_source, &
ind_receptor, &
mquasilag, &
nested_output, &
linit_cond
! Presetting namelist command
ldirect=1
ibdate=20000101
ibtime=0
iedate=20000102
ietime=0
loutstep=10800
loutaver=10800
loutsample=900
itsplit=999999999
lsynctime=900
ctl=-5.0
ifine=4
iout=3
ipout=0
lsubgrid=1
lconvection=1
lagespectra=0
ipin=1
ioutputforeachrelease=0
iflux=1
mdomainfill=0
ind_source=1
ind_receptor=1
mquasilag=0
nested_output=0
linit_cond=0
print*,'command2nml V1.0 converts FLEXPART COMMAND files to namelist format'
select case (iargc())
case (2)
call getarg(1,infile)
call getarg(2,outfile)
case default
print*,'USAGE: command2nml COMMAND.input COMMAND.namelist.output'
stop
end select
! Open the command file and read user options
! Namelist input first: try to read as namelist file
!**************************************************************************
open(10,file=trim(infile),status='old', form='formatted',iostat=readerror)
! If fail, check if file does not exist
if (readerror.ne.0) then
print*,'***ERROR: file COMMAND not found at ',trim(infile)
stop
endif
read(10,command,iostat=readerror)
close(10)
! If error in namelist format, try to open with old input code
if (readerror.ne.0) then
open(10,file=trim(infile),status='old', err=999)
! Check the format of the COMMAND file (either in free format,
! or using formatted mask)
! Use of formatted mask is assumed if line 10 contains the word 'DIRECTION'
!**************************************************************************
call skplin(9,10)
read (10,901) line
901 format (a)
if (index(line,'LDIRECT') .eq. 0) then
old = .false.
else
old = .true.
endif
rewind(10)
! Read parameters
!****************
call skplin(7,10)
if (old) call skplin(1,10)
read(10,*) ldirect
if (old) call skplin(3,10)
read(10,*) ibdate,ibtime
if (old) call skplin(3,10)
read(10,*) iedate,ietime
if (old) call skplin(3,10)
read(10,*) loutstep
if (old) call skplin(3,10)
read(10,*) loutaver
if (old) call skplin(3,10)
read(10,*) loutsample
if (old) call skplin(3,10)
read(10,*) itsplit
if (old) call skplin(3,10)
read(10,*) lsynctime
if (old) call skplin(3,10)
read(10,*) ctl
if (old) call skplin(3,10)
read(10,*) ifine
if (old) call skplin(3,10)
read(10,*) iout
if (old) call skplin(3,10)
read(10,*) ipout
if (old) call skplin(3,10)
read(10,*) lsubgrid
if (old) call skplin(3,10)
read(10,*) lconvection
if (old) call skplin(3,10)
read(10,*) lagespectra
if (old) call skplin(3,10)
read(10,*) ipin
if (old) call skplin(3,10)
read(10,*) ioutputforeachrelease
if (old) call skplin(3,10)
read(10,*) iflux
if (old) call skplin(3,10)
read(10,*) mdomainfill
if (old) call skplin(3,10)
read(10,*) ind_source
if (old) call skplin(3,10)
read(10,*) ind_receptor
if (old) call skplin(3,10)
read(10,*) mquasilag
if (old) call skplin(3,10)
read(10,*) nested_output
if (old) call skplin(3,10)
read(10,*) linit_cond
close(10)
endif ! input format
print*,'Input file read from ',trim(infile)
! write command file in namelist format to output directory if requested
open(11,file=trim(outfile),status='replace',err=998)
write(11,nml=command)
close(11)
print*,'Output file successfully created at ',trim(outfile)
stop
998 print*,' ERROR: Output file not found at ',trim(outfile)
stop
999 print*,' ERROR: Input file "COMMAND" not found at ',trim(infile)
stop
end program command2nml
!**********************************************************************
! 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/>. *
!**********************************************************************
!*******************************************************************************
! Include file for calculation of particle trajectories (Program FLEXPART) *
! This file contains the parameter statements used in FLEXPART *
! *
! Author: A. Stohl *
! *
! 1997 *
! *
! Last update 10 August 2000 *
! *
!*******************************************************************************
module par_mod
implicit none
!****************************************************************
! Parameter defining KIND parameter for "double precision"
!****************************************************************
integer,parameter :: dp=selected_real_kind(P=15)
!***********************************************************
! Number of directories/files used for FLEXPART input/output
!***********************************************************
integer,parameter :: numpath=4
! numpath Number of different pathnames for input/output files
!*****************************
! Physical and other constants
!*****************************
real,parameter :: pi=3.14159265, r_earth=6.371e6, r_air=287.05, ga=9.81
real,parameter :: cpa=1004.6, kappa=0.286, pi180=pi/180., vonkarman=0.4
! pi number "pi"
! pi180 pi/180.
! r_earth radius of earth [m]
! r_air individual gas constant for dry air [J/kg/K]
! ga gravity acceleration of earth [m/s**2]
! cpa specific heat for dry air
! kappa exponent of formula for potential temperature
! vonkarman von Karman constant
real,parameter :: karman=0.40, href=15., convke=2.0
real,parameter :: hmixmin=100., hmixmax=4500., turbmesoscale=0.16
real,parameter :: d_trop=50., d_strat=0.1
! karman Karman's constant
! href [m] Reference height for dry deposition
! konvke Relative share of kinetic energy used for parcel lifting
! hmixmin,hmixmax Minimum and maximum allowed PBL height
! turbmesoscale the factor by which standard deviations of winds at grid
! points surrounding the particle positions are scaled to
! yield the scales for the mesoscale wind velocity fluctuations
! d_trop [m2/s] Turbulent diffusivity for horizontal components in the troposphere
! d_strat [m2/s] Turbulent diffusivity for vertical component in the stratosphere
real,parameter :: xmwml=18.016/28.960
! xmwml ratio of molar weights of water vapor and dry air
!****************************************************
! Constants related to the stratospheric ozone tracer
!****************************************************
real,parameter :: ozonescale=60., pvcrit=2.0
! ozonescale ppbv O3 per PV unit
! pvcrit PV level of the tropopause
!********************
! Some time constants
!********************
integer,parameter :: idiffnorm=10800, idiffmax=2*idiffnorm, minstep=1
! idiffnorm [s] normal time interval between two wind fields
! idiffmax [s] maximum time interval between two wind fields
! minstep [s] minimum time step to be used within FLEXPART
!*****************************************************************
! Parameters for polar stereographic projection close to the poles
!*****************************************************************
real,parameter :: switchnorth=75., switchsouth=-75.
! switchnorth use polar stereographic grid north of switchnorth
! switchsouth use polar stereographic grid south of switchsouth
!*********************************************
! Maximum dimensions of the input mother grids
!*********************************************
integer,parameter :: nxmax=361,nymax=181,nuvzmax=92,nwzmax=92,nzmax=92
!integer,parameter :: nxmax=361,nymax=181,nuvzmax=26,nwzmax=26,nzmax=26
!integer,parameter :: nxmax=721,nymax=361,nuvzmax=64,nwzmax=64,nzmax=64
!integer,parameter :: nxmax=1201,nymax=235,nuvzmax=58,nwzmax=58,nzmax=58
integer,parameter :: nxshift=359 ! for ECMWF
!integer,parameter :: nxshift=0 ! for GFS
integer,parameter :: nconvlevmax = nuvzmax-1
integer,parameter :: na = nconvlevmax+1
! nxmax,nymax maximum dimension of wind fields in x and y
! direction, respectively
! nuvzmax,nwzmax maximum dimension of (u,v) and (w) wind fields in z
! direction (for fields on eta levels)
! nzmax maximum dimension of wind fields in z direction
! for the transformed Cartesian coordinates
! nxshift for global grids (in x), the grid can be shifted by
! nxshift grid points, in order to accomodate nested
! grids, and output grids overlapping the domain "boundary"
! nxshift must not be negative; "normal" setting would be 0
! ntracermax maximum number of tracer species in convection
! nconvlevmax maximum number of levels for convection
! na parameter used in Emanuel's convect subroutine
!*********************************************
! Maximum dimensions of the nested input grids
!*********************************************
integer,parameter :: maxnests=0, nxmaxn=0, nymaxn=0
!integer,parameter :: maxnests=1,nxmaxn=251,nymaxn=151
! maxnests maximum number of nested grids
! nxmaxn,nymaxn maximum dimension of nested wind fields in
! x and y direction, respectively
!*********************************
! Parmaters for GRIB file decoding
!*********************************
integer,parameter :: jpack=4*nxmax*nymax, jpunp=4*jpack
! jpack,jpunp maximum dimensions needed for GRIB file decoding
!**************************************
! Maximum dimensions of the output grid
!**************************************
!integer,parameter :: maxageclass=1,maxzgrid=10,nclassunc=1
integer,parameter :: maxageclass=1,nclassunc=1
! nclassunc number of classes used to calculate the uncertainty
! of the output
! maxageclass maximum number of age classes used for output
! Sabine Eckhardt, June, 2008
! the dimensions of the OUTGRID are now set dynamically during runtime
! maxxgrid,maxygrid,maxzgrid maximum dimensions in x,y,z direction
! maxxgridn,maxygridn maximum dimension of the nested grid
!integer maxxgrid,maxygrid,maxzgrid,maxxgridn,maxygridn
!integer,parameter :: maxxgrid=361,maxygrid=181,maxxgridn=0,maxygridn=0)
integer,parameter :: maxreceptor=200
! maxreceptor maximum number of receptor points
!**************************************************
! Maximum number of particles, species, and similar
!**************************************************
integer,parameter :: maxpart=50000
integer,parameter :: maxspec=1
! maxpart Maximum number of particles
! maxspec Maximum number of chemical species per release
! maxpoint is also set dynamically during runtime
! maxpoint Maximum number of release locations
! ---------
! Sabine Eckhardt: change of landuse inventary numclass=13
! ---------
integer,parameter :: maxwf=50000, maxtable=1000, numclass=13, ni=11
! maxwf maximum number of wind fields to be used for simulation
! maxtable Maximum number of chemical species that can be
! tabulated for FLEXPART
! numclass Number of landuse classes available to FLEXPART
! ni Number of diameter classes of particles
!**************************************************************************
! dimension of the OH field
!**************************************************************************
integer,parameter :: maxxOH=72, maxyOH=46, maxzOH=7
!**************************************************************************
! Maximum number of particles to be released in a single atmospheric column
! for the domain-filling trajectories option
!**************************************************************************
integer,parameter :: maxcolumn=3000
!*********************************
! Dimension of random number field
!*********************************
integer,parameter :: maxrand=2000000
! maxrand number of random numbers used
!*****************************************************
! Number of clusters to be used for plume trajectories
!*****************************************************
integer,parameter :: ncluster=5
!************************************
! Unit numbers for input/output files
!************************************
integer,parameter :: unitpath=1, unitcommand=1, unitageclasses=1, unitgrid=1
integer,parameter :: unitavailab=1, unitreleases=88, unitpartout=93
integer,parameter :: unitpartin=93, unitflux=98, unitouttraj=96
integer,parameter :: unitvert=1, unitoro=1, unitpoin=1, unitreceptor=1
integer,parameter :: unitoutgrid=97, unitoutgridppt=99, unitoutinfo=1
integer,parameter :: unitspecies=1, unitoutrecept=91, unitoutreceptppt=92
integer,parameter :: unitlsm=1, unitsurfdata=1, unitland=1, unitwesely=1
integer,parameter :: unitOH=1
integer,parameter :: unitdates=94, unitheader=90, unitshortpart=95
integer,parameter :: unitboundcond=89
end module par_mod
!**********************************************************************
! 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 skplin(nlines,iunit)
! i i
!*****************************************************************************
! *
! This routine reads nlines from unit iunit and discards them *
! *
! Authors: Petra Seibert *
! *
! 31 Dec 1998 *
! *
!*****************************************************************************
! *
! Variables: *
! *
! iunit unit number from which lines are to be skipped *
! nlines number of lines to be skipped *
! *
!*****************************************************************************
implicit none
integer :: i,iunit, nlines
do i=1,nlines
read(iunit,*)
end do
end subroutine skplin
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment