CTBTO FLEXPART WO4 (2015-10-15)
 All Classes Files Functions Variables
readcommand.F90
Go to the documentation of this file.
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 subroutine readcommand
23 
24  !*****************************************************************************
25  ! *
26  ! This routine reads the user specifications for the current model run. *
27  ! *
28  ! Author: A. Stohl *
29  ! *
30  ! 18 May 1996 *
31  ! *
32  !*****************************************************************************
33  ! *
34  ! Variables: *
35  ! bdate beginning date as Julian date *
36  ! ctl factor by which time step must be smaller than *
37  ! Lagrangian time scale *
38  ! ibdate,ibtime beginnning date and time (YYYYMMDD, HHMISS) *
39  ! ideltas [s] modelling period *
40  ! iedate,ietime ending date and time (YYYYMMDD, HHMISS) *
41  ! ifine reduction factor for vertical wind time step *
42  ! outputforeachrel for forward runs it is possible either to create *
43  ! one outputfield or several for each releasepoint *
44  ! iflux switch to turn on (1)/off (0) flux calculations *
45  ! iout 1 for conc. (residence time for backward runs) output,*
46  ! 2 for mixing ratio output, 3 both, 4 for plume *
47  ! trajectory output, 5 = options 1 and 4 *
48  ! ipin 1 continue simulation with dumped particle data, 0 no *
49  ! ipout 0 no particle dump, 1 every output time, 3 only at end*
50  ! itsplit [s] time constant for particle splitting *
51  ! loutaver [s] concentration output is an average over loutaver *
52  ! seconds *
53  ! loutsample [s] average is computed from samples taken every [s] *
54  ! seconds *
55  ! loutstep [s] time interval of concentration output *
56  ! lsynctime [s] synchronisation time interval for all particles *
57  ! lagespectra switch to turn on (1)/off (0) calculation of age *
58  ! spectra *
59  ! lconvection value of either 0 and 1 indicating mixing by *
60  ! convection *
61  ! = 0 .. no convection *
62  ! + 1 .. parameterisation of mixing by subgrid-scale *
63  ! convection = on *
64  ! lsubgrid switch to turn on (1)/off (0) subgrid topography *
65  ! parameterization *
66  ! method method used to compute the particle pseudovelocities *
67  ! mdomainfill 1 use domain-filling option, 0 not, 2 use strat. O3 *
68  ! *
69  ! Constants: *
70  ! unitcommand unit connected to file COMMAND *
71  ! *
72  !*****************************************************************************
73 
74  use par_mod
75  use com_mod
76 
77  implicit none
78 
79  real(kind=dp) :: juldate
80  character(len=50) :: line
81  logical :: old
82 
83 
84  ! Open the command file and read user options
85  !********************************************
86 
87 
88  open(unitcommand,file=path(1)(1:length(1))//'COMMAND',status='old', &
89  err=999)
90 
91  ! Check the format of the COMMAND file (either in free format,
92  ! or using formatted mask)
93  ! Use of formatted mask is assumed if line 10 contains the word 'DIRECTION'
94  !**************************************************************************
95 
96  call skplin(9,unitcommand)
97  read (unitcommand,901) line
98 901 format (a)
99  if (index(line,'LDIRECT') .eq. 0) then
100  old = .false.
101  else
102  old = .true.
103  endif
104  rewind(unitcommand)
105 
106  ! Read parameters
107  !****************
108 
109  call skplin(7,unitcommand)
110  if (old) call skplin(1,unitcommand)
111 
112  read(unitcommand,*) ldirect
113  if (old) call skplin(3,unitcommand)
114  read(unitcommand,*) ibdate,ibtime
115  if (old) call skplin(3,unitcommand)
116  read(unitcommand,*) iedate,ietime
117  if (old) call skplin(3,unitcommand)
118  read(unitcommand,*) loutstep
119  if (old) call skplin(3,unitcommand)
120  read(unitcommand,*) loutaver
121  if (old) call skplin(3,unitcommand)
122  read(unitcommand,*) loutsample
123  if (old) call skplin(3,unitcommand)
124  read(unitcommand,*) itsplit
125  if (old) call skplin(3,unitcommand)
126  read(unitcommand,*) lsynctime
127  if (old) call skplin(3,unitcommand)
128  read(unitcommand,*) ctl
129  if (old) call skplin(3,unitcommand)
130  read(unitcommand,*) ifine
131  if (old) call skplin(3,unitcommand)
132  read(unitcommand,*) iout
133  if (old) call skplin(3,unitcommand)
134  read(unitcommand,*) ipout
135  if (old) call skplin(3,unitcommand)
136  read(unitcommand,*) lsubgrid
137  if (old) call skplin(3,unitcommand)
138  read(unitcommand,*) lconvection
139  if (old) call skplin(3,unitcommand)
140  read(unitcommand,*) lagespectra
141  if (old) call skplin(3,unitcommand)
142  read(unitcommand,*) ipin
143  if (old) call skplin(3,unitcommand)
144  read(unitcommand,*) ioutputforeachrelease
145  if (old) call skplin(3,unitcommand)
146  read(unitcommand,*) iflux
147  if (old) call skplin(3,unitcommand)
148  read(unitcommand,*) mdomainfill
149  if (old) call skplin(3,unitcommand)
150  read(unitcommand,*) ind_source
151  if (old) call skplin(3,unitcommand)
152  read(unitcommand,*) ind_receptor
153  if (old) call skplin(3,unitcommand)
154  read(unitcommand,*) mquasilag
155  if (old) call skplin(3,unitcommand)
156  read(unitcommand,*) nested_output
157  if (old) call skplin(3,unitcommand)
158  read(unitcommand,*) linit_cond
159  if (old) call skplin(3,unitcommand)
160  read(unitcommand,*) preprocessed_metdata
161  close(unitcommand)
162 
163  ifine=max(ifine,1)
164 
165 
166  ! Determine how Markov chain is formulated (for w or for w/sigw)
167  !***************************************************************
168 
169  if (ctl.ge.0.1) then
170  turbswitch=.true.
171  else
172  turbswitch=.false.
173  ifine=1
174  endif
175  fine=1./real(ifine)
176  ctl=1./ctl
177 
178  ! Set the switches required for the various options for input/output units
179  !*************************************************************************
180  !AF Set the switches IND_REL and IND_SAMP for the release and sampling
181  !Af switches for the releasefile:
182  !Af IND_REL = 1 : xmass * rho
183  !Af IND_REL = 0 : xmass * 1
184 
185  !Af switches for the conccalcfile:
186  !AF IND_SAMP = 0 : xmass * 1
187  !Af IND_SAMP = -1 : xmass / rho
188 
189  !AF IND_SOURCE switches between different units for concentrations at the source
190  !Af NOTE that in backward simulations the release of computational particles
191  !Af takes place at the "receptor" and the sampling of p[articles at the "source".
192  !Af 1 = mass units
193  !Af 2 = mass mixing ratio units
194  !Af IND_RECEPTOR switches between different units for concentrations at the receptor
195  !Af 1 = mass units
196  !Af 2 = mass mixing ratio units
197 
198  if ( ldirect .eq. 1 ) then ! FWD-Run
199  !Af set release-switch
200  if (ind_source .eq. 1 ) then !mass
201  ind_rel = 0
202  else ! mass mix
203  ind_rel = 1
204  endif
205  !Af set sampling switch
206  if (ind_receptor .eq. 1) then !mass
207  ind_samp = 0
208  else ! mass mix
209  ind_samp = -1
210  endif
211  elseif (ldirect .eq. -1 ) then !BWD-Run
212  !Af set sampling switch
213  if (ind_source .eq. 1 ) then !mass
214  ind_samp = -1
215  else ! mass mix
216  ind_samp = 0
217  endif
218  !Af set release-switch
219  if (ind_receptor .eq. 1) then !mass
220  ind_rel = 1
221  else ! mass mix
222  ind_rel = 0
223  endif
224  endif
225 
226  !*************************************************************
227  ! Check whether valid options have been chosen in file COMMAND
228  !*************************************************************
229 
230  ! Check options for initial condition output: Switch off for forward runs
231  !************************************************************************
232 
233  if (ldirect.eq.1) linit_cond=0
234  if ((linit_cond.lt.0).or.(linit_cond.gt.2)) then
235  write(*,*) ' #### FLEXPART MODEL ERROR! INVALID OPTION #### '
236  write(*,*) ' #### FOR LINIT_COND IN FILE "COMMAND". #### '
237  stop
238  endif
239 
240  ! Check input dates
241  !******************
242 
243  if (iedate.lt.ibdate) then
244  write(*,*) ' #### FLEXPART MODEL ERROR! BEGINNING DATE #### '
245  write(*,*) ' #### IS LARGER THAN ENDING DATE. CHANGE #### '
246  write(*,*) ' #### EITHER POINT 2 OR POINT 3 IN FILE #### '
247  write(*,*) ' #### "COMMAND". #### '
248  stop
249  else if (iedate.eq.ibdate) then
250  if (ietime.lt.ibtime) then
251  write(*,*) ' #### FLEXPART MODEL ERROR! BEGINNING TIME #### '
252  write(*,*) ' #### IS LARGER THAN ENDING TIME. CHANGE #### '
253  write(*,*) ' #### EITHER POINT 2 OR POINT 3 IN FILE #### '
254  write(*,*) ' #### "COMMAND". #### '
255  stop
256  endif
257  endif
258 
259 
260  ! Determine kind of dispersion method
261  !************************************
262 
263  if (ctl.gt.0.) then
264  method=1
265  mintime=minstep
266  else
267  method=0
268  mintime=lsynctime
269  endif
270 
271  ! Check whether a valid option for gridded model output has been chosen
272  !**********************************************************************
273 
274  if ((iout.lt.1).or.(iout.gt.5)) then
275  write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND: #### '
276  write(*,*) ' #### IOUT MUST BE 1, 2, 3, 4, OR 5! #### '
277  stop
278  endif
279 
280  !AF check consistency between units and volume mixing ratio
281  if ( ((iout.eq.2).or.(iout.eq.3)).and. &
282  (ind_source.gt.1 .or.ind_receptor.gt.1) ) then
283  write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND: #### '
284  write(*,*) ' #### VOLUME MIXING RATIO ONLY SUPPORTED #### '
285  write(*,*) ' #### FOR MASS UNITS (at the moment) #### '
286  stop
287  endif
288 
289 
290 
291  ! For quasilag output for each release is forbidden
292  !*****************************************************************************
293 
294  if ((ioutputforeachrelease.eq.1).and.(mquasilag.eq.1)) then
295  write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND: ####'
296  write(*,*) '#### OUTPUTFOREACHRELEASE AND QUASILAGRANGIAN####'
297  write(*,*) '#### MODE IS NOT POSSIBLE ! ####'
298  stop
299  endif
300 
301 
302  ! For quasilag backward is forbidden
303  !*****************************************************************************
304 
305  if ((ldirect.lt.0).and.(mquasilag.eq.1)) then
306  write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND: ####'
307  write(*,*) '#### FOR BACKWARD RUNS, QUASILAGRANGIAN MODE ####'
308  write(*,*) '#### IS NOT POSSIBLE ! ####'
309  stop
310  endif
311 
312 
313  ! For backward runs one releasefield for all releases makes no sense,
314  ! For quasilag and domainfill ioutputforechrelease is forbidden
315  !*****************************************************************************
316 #if defined WITH_CTBTO_PATCHES
317 ! if ((ldirect.lt.0).and.(ioutputforeachrelease.eq.0)) then
318 ! write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND: ####'
319 ! write(*,*) '#### FOR BACKWARD RUNS, IOUTPUTFOREACHRLEASE ####'
320 ! write(*,*) '#### MUST BE SET TO ONE! ####'
321 ! stop
322 ! endif
323 #else
324  if ((ldirect.lt.0).and.(ioutputforeachrelease.eq.0)) then
325  write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND: ####'
326  write(*,*) '#### FOR BACKWARD RUNS, IOUTPUTFOREACHRLEASE ####'
327  write(*,*) '#### MUST BE SET TO ONE! ####'
328  stop
329  endif
330 #endif
331 
332 
333  ! For backward runs one releasefield for all releases makes no sense,
334  ! and is "forbidden"
335  !*****************************************************************************
336 
337  if ((mdomainfill.eq.1).and.(ioutputforeachrelease.eq.1)) then
338  write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND: ####'
339  write(*,*) '#### FOR DOMAIN FILLING RUNS OUTPUT FOR ####'
340  write(*,*) '#### EACH RELEASE IS FORBIDDEN ! ####'
341  stop
342  endif
343 
344 
345  ! For domain-filling trajectories, a plume centroid trajectory makes no sense,
346  ! For backward runs, only residence time output (iout=1) or plume trajectories (iout=4),
347  ! or both (iout=5) makes sense; other output options are "forbidden"
348  !*****************************************************************************
349 
350  if (ldirect.lt.0) then
351  if ((iout.eq.2).or.(iout.eq.3)) then
352  write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND: ####'
353  write(*,*) '#### FOR BACKWARD RUNS, IOUT MUST BE 1,4,OR 5####'
354  stop
355  endif
356  endif
357 
358 
359  ! For domain-filling trajectories, a plume centroid trajectory makes no sense,
360  ! and is "forbidden"
361  !*****************************************************************************
362 
363  if (mdomainfill.ge.1) then
364  if ((iout.eq.4).or.(iout.eq.5)) then
365  write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND: ####'
366  write(*,*) '#### FOR DOMAIN-FILLING TRAJECTORY OPTION, ####'
367  write(*,*) '#### IOUT MUST NOT BE SET TO 4 OR 5. ####'
368  stop
369  endif
370  endif
371 
372 
373 
374  ! Check whether a valid options for particle dump has been chosen
375  !****************************************************************
376 
377  if ((ipout.ne.0).and.(ipout.ne.1).and.(ipout.ne.2)) then
378  write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND: #### '
379  write(*,*) ' #### IPOUT MUST BE 1, 2 OR 3! #### '
380  stop
381  endif
382 
383  if(lsubgrid.ne.1) then
384  write(*,*) ' ---------------- '
385  write(*,*) ' INFORMATION: SUBGRIDSCALE TERRAIN EFFECT IS'
386  write(*,*) ' NOT PARAMETERIZED DURING THIS SIMULATION. '
387  write(*,*) ' ---------------- '
388  endif
389 
390 
391  ! Check whether convection scheme is either turned on or off
392  !***********************************************************
393 
394  if ((lconvection.ne.0).and.(lconvection.ne.1)) then
395  write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND: #### '
396  write(*,*) ' #### LCONVECTION MUST BE SET TO EITHER 1 OR 0#### '
397  stop
398  endif
399 
400 
401  ! Check whether synchronisation interval is sufficiently short
402  !*************************************************************
403 
404  if (lsynctime.gt.(idiffnorm/2)) then
405  write(*,*) ' #### FLEXPART MODEL ERROR! SYNCHRONISATION #### '
406  write(*,*) ' #### TIME IS TOO LONG. MAKE IT SHORTER. #### '
407  write(*,*) ' #### MINIMUM HAS TO BE: ', idiffnorm/2
408  stop
409  endif
410 
411 
412  ! Check consistency of the intervals, sampling periods, etc., for model output
413  !*****************************************************************************
414 
415  if (loutaver.eq.0) then
416  write(*,*) ' #### FLEXPART MODEL ERROR! TIME AVERAGE OF #### '
417  write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE #### '
418  write(*,*) ' #### ZERO. #### '
419  write(*,*) ' #### CHANGE INPUT IN FILE COMMAND. #### '
420  stop
421  endif
422 
423  if (loutaver.gt.loutstep) then
424  write(*,*) ' #### FLEXPART MODEL ERROR! TIME AVERAGE OF #### '
425  write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE #### '
426  write(*,*) ' #### GREATER THAN INTERVAL OF OUTPUT. #### '
427  write(*,*) ' #### CHANGE INPUT IN FILE COMMAND. #### '
428  stop
429  endif
430 
431  if (loutsample.gt.loutaver) then
432  write(*,*) ' #### FLEXPART MODEL ERROR! SAMPLING TIME OF #### '
433  write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE #### '
434  write(*,*) ' #### GREATER THAN TIME AVERAGE OF OUTPUT. #### '
435  write(*,*) ' #### CHANGE INPUT IN FILE COMMAND. #### '
436  stop
437  endif
438 
439  if (mod(loutaver,lsynctime).ne.0) then
440  write(*,*) ' #### FLEXPART MODEL ERROR! AVERAGING TIME OF #### '
441  write(*,*) ' #### CONCENTRATION FIELD MUST BE A MULTIPLE #### '
442  write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL #### '
443  stop
444  endif
445 
446  if ((loutaver/lsynctime).lt.2) then
447  write(*,*) ' #### FLEXPART MODEL ERROR! AVERAGING TIME OF #### '
448  write(*,*) ' #### CONCENTRATION FIELD MUST BE AT LEAST #### '
449  write(*,*) ' #### TWICE THE SYNCHRONISATION INTERVAL #### '
450  stop
451  endif
452 
453  if (mod(loutstep,lsynctime).ne.0) then
454  write(*,*) ' #### FLEXPART MODEL ERROR! INTERVAL BETWEEN #### '
455  write(*,*) ' #### CONCENTRATION FIELDS MUST BE A MULTIPLE #### '
456  write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL #### '
457  stop
458  endif
459 
460  if ((loutstep/lsynctime).lt.2) then
461  write(*,*) ' #### FLEXPART MODEL ERROR! INTERVAL BETWEEN #### '
462  write(*,*) ' #### CONCENTRATION FIELDS MUST BE AT LEAST #### '
463  write(*,*) ' #### TWICE THE SYNCHRONISATION INTERVAL #### '
464  stop
465  endif
466 
467  if (mod(loutsample,lsynctime).ne.0) then
468  write(*,*) ' #### FLEXPART MODEL ERROR! SAMPLING TIME OF #### '
469  write(*,*) ' #### CONCENTRATION FIELD MUST BE A MULTIPLE #### '
470  write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL #### '
471  stop
472  endif
473 
474  if (itsplit.lt.loutaver) then
475  write(*,*) ' #### FLEXPART MODEL ERROR! SPLITTING TIME FOR#### '
476  write(*,*) ' #### PARTICLES IS TOO SHORT. PLEASE INCREASE #### '
477  write(*,*) ' #### SPLITTING TIME CONSTANT. #### '
478  stop
479  endif
480 
481  if ((mquasilag.eq.1).and.(iout.ge.4)) then
482  write(*,*) ' #### FLEXPART MODEL ERROR! CONFLICTING #### '
483  write(*,*) ' #### OPTIONS: IF MQUASILAG=1, PLUME #### '
484  write(*,*) ' #### TRAJECTORY OUTPUT IS IMPOSSIBLE. #### '
485  stop
486  endif
487 
488  ! Compute modeling time in seconds and beginning date in Julian date
489  !*******************************************************************
490 
491  outstep=real(abs(loutstep))
492  if (ldirect.eq.1) then
493  bdate=juldate(ibdate,ibtime)
494  edate=juldate(iedate,ietime)
495  ideltas=nint((edate-bdate)*86400.)
496  else if (ldirect.eq.-1) then
497  loutaver=-1*loutaver
498  loutstep=-1*loutstep
499  loutsample=-1*loutsample
500  lsynctime=-1*lsynctime
501  bdate=juldate(iedate,ietime)
502  edate=juldate(ibdate,ibtime)
503  ideltas=nint((edate-bdate)*86400.)
504  else
505  write(*,*) ' #### FLEXPART MODEL ERROR! DIRECTION IN #### '
506  write(*,*) ' #### FILE "COMMAND" MUST BE EITHER -1 OR 1. #### '
507  stop
508  endif
509 
510  return
511 
512 999 write(*,*) ' #### FLEXPART MODEL ERROR! FILE "COMMAND" #### '
513  write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY #### '
514  write(*,'(a)') path(1)(1:length(1))
515  stop
516 
517 end subroutine readcommand
subroutine skplin(nlines, iunit)
Definition: skplin.f90:22
subroutine readcommand
Definition: readcommand.F90:22
real(kind=dp) function juldate(yyyymmdd, hhmiss)
Definition: juldate.f90:22