source: flexpart.git/src/readcommand.f90 @ 78e62dc

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 78e62dc was 78e62dc, checked in by flexpart <>, 9 years ago

New OH parameter in SPECIES files (now 3 instead of 2). New path to OH binariy files.

  • Property mode set to 100644
File size: 23.3 KB
Line 
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
22subroutine 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  !     HSO, 1 July 2014                                                       *
32  !     Added optional namelist input                                          *
33  !                                                                            *
34  !*****************************************************************************
35  !                                                                            *
36  ! Variables:                                                                 *
37  ! bdate                beginning date as Julian date                         *
38  ! ctl                  factor by which time step must be smaller than        *
39  !                      Lagrangian time scale                                 *
40  ! ibdate,ibtime        beginnning date and time (YYYYMMDD, HHMISS)           *
41  ! ideltas [s]          modelling period                                      *
42  ! iedate,ietime        ending date and time (YYYYMMDD, HHMISS)               *
43  ! ifine                reduction factor for vertical wind time step          *
44  ! outputforeachrel     for forward runs it is possible either to create      *
45  !                      one outputfield or several for each releasepoint      *
46  ! iflux                switch to turn on (1)/off (0) flux calculations       *
47  ! iout                 1 for conc. (residence time for backward runs) output,*
48  !                      2 for mixing ratio output, 3 both, 4 for plume        *
49  !                      trajectory output, 5 = options 1 and 4                *
50  ! ipin                 1 continue simulation with dumped particle data, 0 no *
51  ! ipout                0 no particle dump, 1 every output time, 3 only at end*
52  ! itsplit [s]          time constant for particle splitting                  *
53  ! loutaver [s]         concentration output is an average over loutaver      *
54  !                      seconds                                               *
55  ! loutsample [s]       average is computed from samples taken every [s]      *
56  !                      seconds                                               *
57  ! loutstep [s]         time interval of concentration output                 *
58  ! lsynctime [s]        synchronisation time interval for all particles       *
59  ! lagespectra          switch to turn on (1)/off (0) calculation of age      *
60  !                      spectra                                               *
61  ! lconvection          value of either 0 and 1 indicating mixing by          *
62  !                      convection                                            *
63  !                      = 0 .. no convection                                  *
64  !                      + 1 .. parameterisation of mixing by subgrid-scale    *
65  !                              convection = on                               *
66  ! lsubgrid             switch to turn on (1)/off (0) subgrid topography      *
67  !                      parameterization                                      *
68  ! method               method used to compute the particle pseudovelocities  *
69  ! mdomainfill          1 use domain-filling option, 0 not, 2 use strat. O3   *
70  !                                                                            *
71  ! Constants:                                                                 *
72  ! unitcommand          unit connected to file COMMAND                        *
73  !                                                                            *
74  !*****************************************************************************
75
76  use par_mod
77  use com_mod
78
79  implicit none
80
81  real(kind=dp) :: juldate
82  character(len=50) :: line
83  logical :: old
84  integer :: readerror
85
86  namelist /command/ &
87  ldirect, &
88  ibdate,ibtime, &
89  iedate,ietime, &
90  loutstep, &
91  loutaver, &
92  loutsample, &
93  itsplit, &
94  lsynctime, &
95  ctl, &
96  ifine, &
97  iout, &
98  ipout, &
99  lsubgrid, &
100  lconvection, &
101  lagespectra, &
102  ipin, &
103  ioutputforeachrelease, &
104  iflux, &
105  mdomainfill, &
106  ind_source, &
107  ind_receptor, &
108  mquasilag, &
109  nested_output, &
110  linit_cond, &
111  lnetcdfout, &
112  surf_only, &
113  cblflag, &
114  ohfields_path
115
116  ! Presetting namelist command
117  ldirect=0
118  ibdate=20000101
119  ibtime=0
120  iedate=20000102
121  ietime=0
122  loutstep=10800
123  loutaver=10800
124  loutsample=900
125  itsplit=999999999
126  lsynctime=900
127  ctl=-5.0
128  ifine=4
129  iout=3
130  ipout=0
131  lsubgrid=1
132  lconvection=1
133  lagespectra=0
134  ipin=1
135  ioutputforeachrelease=1
136  iflux=1
137  mdomainfill=0
138  ind_source=1
139  ind_receptor=1
140  mquasilag=0
141  nested_output=0
142  linit_cond=0
143  lnetcdfout=0
144  surf_only=0
145  cblflag=0
146  ohfields_path="../../flexin/"
147
148  ! Open the command file and read user options
149  ! Namelist input first: try to read as namelist file
150  !**************************************************************************
151  open(unitcommand,file=path(1)(1:length(1))//'COMMAND',status='old',form='formatted',err=999)
152
153  ! try namelist input (default)
154  read(unitcommand,command,iostat=readerror)
155  close(unitcommand)
156
157  ! distinguish namelist from fixed text input
158  if ((readerror.ne.0).or.(ldirect.eq.0)) then ! parse as text file format
159 
160    open(unitcommand,file=path(1)(1:length(1))//'COMMAND',status='old', err=999)
161
162    ! Check the format of the COMMAND file (either in free format,
163    ! or using formatted mask)
164    ! Use of formatted mask is assumed if line 10 contains the word 'DIRECTION'
165    !**************************************************************************
166
167    call skplin(9,unitcommand)
168    read (unitcommand,901) line
169  901   format (a)
170    if (index(line,'LDIRECT') .eq. 0) then
171      old = .false.
172      write(*,*) 'COMMAND in old short format, please update to namelist format'
173    else
174      old = .true.
175      write(*,*) 'COMMAND in old long format, please update to namelist format'
176    endif
177    rewind(unitcommand)
178
179
180    ! Read parameters
181    !****************
182
183    call skplin(7,unitcommand)
184    if (old) call skplin(1,unitcommand)
185    read(unitcommand,*) ldirect
186    if (old) call skplin(3,unitcommand)
187    read(unitcommand,*) ibdate,ibtime
188    if (old) call skplin(3,unitcommand)
189    read(unitcommand,*) iedate,ietime
190    if (old) call skplin(3,unitcommand)
191    read(unitcommand,*) loutstep
192    if (old) call skplin(3,unitcommand)
193    read(unitcommand,*) loutaver
194    if (old) call skplin(3,unitcommand)
195    read(unitcommand,*) loutsample
196    if (old) call skplin(3,unitcommand)
197    read(unitcommand,*) itsplit
198    if (old) call skplin(3,unitcommand)
199    read(unitcommand,*) lsynctime
200    if (old) call skplin(3,unitcommand)
201    read(unitcommand,*) ctl
202    if (old) call skplin(3,unitcommand)
203    read(unitcommand,*) ifine
204    if (old) call skplin(3,unitcommand)
205    read(unitcommand,*) iout
206    if (old) call skplin(3,unitcommand)
207    read(unitcommand,*) ipout
208    if (old) call skplin(3,unitcommand)
209    read(unitcommand,*) lsubgrid
210    if (old) call skplin(3,unitcommand)
211    read(unitcommand,*) lconvection
212    if (old) call skplin(3,unitcommand)
213    read(unitcommand,*) lagespectra
214    if (old) call skplin(3,unitcommand)
215    read(unitcommand,*) ipin
216    if (old) call skplin(3,unitcommand)
217    read(unitcommand,*) ioutputforeachrelease
218    if (old) call skplin(3,unitcommand)
219    read(unitcommand,*) iflux
220    if (old) call skplin(3,unitcommand)
221    read(unitcommand,*) mdomainfill
222    if (old) call skplin(3,unitcommand)
223    read(unitcommand,*) ind_source
224    if (old) call skplin(3,unitcommand)
225    read(unitcommand,*) ind_receptor
226    if (old) call skplin(3,unitcommand)
227    read(unitcommand,*) mquasilag
228    if (old) call skplin(3,unitcommand)
229    read(unitcommand,*) nested_output
230    if (old) call skplin(3,unitcommand)
231    read(unitcommand,*) linit_cond
232    if (old) call skplin(3,unitcommand)
233    read(unitcommand,*) surf_only
234    if (old) call skplin(3,unitcommand)  !added by mc
235    read(unitcommand,*) cblflag          !added by mc
236
237    close(unitcommand)
238
239  endif ! input format
240
241  ! write command file in namelist format to output directory if requested
242  if (nmlout.and.lroot) then
243    open(unitcommand,file=path(2)(1:length(2))//'COMMAND.namelist',err=1000)
244    write(unitcommand,nml=command)
245    close(unitcommand)
246  endif
247
248  ifine=max(ifine,1)
249
250  ! Determine how Markov chain is formulated (for w or for w/sigw)
251  !***************************************************************
252  if (cblflag.eq.1) then !---- added by mc to properly set parameters for CBL simulations
253    turbswitch=.true.
254    if (lsynctime>maxtl) lsynctime=maxtl  !maxtl defined in com_mod.f90
255    if (ctl.lt.5) then
256      print *,'WARNING: CBL flag active the ratio of TLu/dt has been set to 5'
257      ctl=5.
258    end if
259    if (ifine*ctl.lt.50) then
260      ifine=int(50./ctl)+1
261
262      print *,'WARNING: CBL flag active the ratio of TLW/dt was < 50, ifine has been re-set to',ifine
263!pause
264    endif
265    print *,'WARNING: CBL flag active the ratio of TLW/dt is ',ctl*ifine
266    print *,'WARNING: CBL flag active lsynctime is ',lsynctime
267  else                    !added by mc
268    if (ctl.ge.0.1) then
269      turbswitch=.true.
270    else
271      turbswitch=.false.
272      ifine=1
273    endif
274  endif                   !added by mc
275  fine=1./real(ifine)
276  ctl=1./ctl
277
278  ! Set the switches required for the various options for input/output units
279  !*************************************************************************
280  !AF Set the switches IND_REL and IND_SAMP for the release and sampling
281  !Af switches for the releasefile:
282  !Af IND_REL =  1 : xmass * rho
283  !Af IND_REL =  0 : xmass * 1
284
285  !Af switches for the conccalcfile:
286  !AF IND_SAMP =  0 : xmass * 1
287  !Af IND_SAMP = -1 : xmass / rho
288
289  !AF IND_SOURCE switches between different units for concentrations at the source
290  !Af   NOTE that in backward simulations the release of computational particles
291  !Af   takes place at the "receptor" and the sampling of p[articles at the "source".
292  !Af          1 = mass units
293  !Af          2 = mass mixing ratio units
294  !Af IND_RECEPTOR switches between different units for concentrations at the receptor
295  !Af          1 = mass units
296  !Af          2 = mass mixing ratio units
297
298  if ( ldirect .eq. 1 ) then  ! FWD-Run
299  !Af set release-switch
300     if (ind_source .eq. 1 ) then !mass
301        ind_rel = 0
302     else ! mass mix
303        ind_rel = 1
304     endif
305  !Af set sampling switch
306     if (ind_receptor .eq. 1) then !mass
307        ind_samp = 0
308     else ! mass mix
309        ind_samp = -1
310     endif
311  elseif (ldirect .eq. -1 ) then !BWD-Run
312  !Af set sampling switch
313     if (ind_source .eq. 1 ) then !mass
314        ind_samp = -1
315     else ! mass mix
316        ind_samp = 0
317     endif
318  !Af set release-switch
319     if (ind_receptor .eq. 1) then !mass
320        ind_rel = 1
321     else ! mass mix
322        ind_rel = 0
323     endif
324  endif
325
326  !*************************************************************
327  ! Check whether valid options have been chosen in file COMMAND
328  !*************************************************************
329
330  ! Check options for initial condition output: Switch off for forward runs
331  !************************************************************************
332
333  if (ldirect.eq.1) linit_cond=0
334  if ((linit_cond.lt.0).or.(linit_cond.gt.2)) then
335    write(*,*) ' #### FLEXPART MODEL ERROR! INVALID OPTION    #### '
336    write(*,*) ' #### FOR LINIT_COND IN FILE "COMMAND".       #### '
337    stop
338  endif
339
340  ! Check input dates
341  !******************
342
343  if (iedate.lt.ibdate) then
344    write(*,*) ' #### FLEXPART MODEL ERROR! BEGINNING DATE    #### '
345    write(*,*) ' #### IS LARGER THAN ENDING DATE. CHANGE      #### '
346    write(*,*) ' #### EITHER POINT 2 OR POINT 3 IN FILE       #### '
347    write(*,*) ' #### "COMMAND".                              #### '
348    stop
349  else if (iedate.eq.ibdate) then
350    if (ietime.lt.ibtime) then
351    write(*,*) ' #### FLEXPART MODEL ERROR! BEGINNING TIME    #### '
352    write(*,*) ' #### IS LARGER THAN ENDING TIME. CHANGE      #### '
353    write(*,*) ' #### EITHER POINT 2 OR POINT 3 IN FILE       #### '
354    write(*,*) ' #### "COMMAND".                              #### '
355    stop
356    endif
357  endif
358
359
360  ! Determine kind of dispersion method
361  !************************************
362
363  if (ctl.gt.0.) then
364    method=1
365    mintime=minstep
366  else
367    method=0
368    mintime=lsynctime
369  endif
370
371!  check for netcdf output switch (use for non-namelist input only!)
372  if (iout.ge.8) then
373     lnetcdfout = 1
374     iout = iout - 8
375! #ifndef NETCDF_OUTPUT
376!      print*,'ERROR: netcdf output not activated during compile time but used in COMMAND file!'
377!      print*,'Please recompile with netcdf library or use standard output format.'
378!      stop
379! #endif
380  endif
381
382  ! Check whether a valid option for gridded model output has been chosen
383  !**********************************************************************
384
385  if ((iout.lt.1).or.(iout.gt.5)) then
386    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
387    write(*,*) ' #### IOUT MUST BE 1, 2, 3, 4, OR 5 FOR       #### '
388    write(*,*) ' #### STANDARD FLEXPART OUTPUT OR  9 - 13     #### '
389    write(*,*) ' #### FOR NETCDF OUTPUT                       #### '
390    stop
391  endif
392
393  !AF check consistency between units and volume mixing ratio
394  if ( ((iout.eq.2).or.(iout.eq.3)).and. &
395       (ind_source.gt.1 .or.ind_receptor.gt.1) ) then
396    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
397    write(*,*) ' #### VOLUME MIXING RATIO ONLY SUPPORTED      #### '
398    write(*,*) ' #### FOR MASS UNITS (at the moment)          #### '
399    stop
400  endif
401
402
403  ! For quasilag output for each release is forbidden
404  !*****************************************************************************
405
406  if ((ioutputforeachrelease.eq.1).and.(mquasilag.eq.1)) then
407      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
408      write(*,*) '#### OUTPUTFOREACHRELEASE AND QUASILAGRANGIAN####'
409      write(*,*) '#### MODE IS NOT POSSIBLE   !                ####'
410      stop
411  endif
412
413
414  ! For quasilag backward is forbidden
415  !*****************************************************************************
416
417  if ((ldirect.lt.0).and.(mquasilag.eq.1)) then
418      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
419      write(*,*) '#### FOR BACKWARD RUNS, QUASILAGRANGIAN MODE ####'
420      write(*,*) '#### IS NOT POSSIBLE   !                     ####'
421      stop
422  endif
423
424
425  ! For backward runs one releasefield for all releases makes no sense,
426  ! For quasilag and domainfill ioutputforechrelease is forbidden
427  !*****************************************************************************
428
429  if ((ldirect.lt.0).and.(ioutputforeachrelease.eq.0)) then
430      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
431      write(*,*) '#### FOR BACKWARD RUNS, IOUTPUTFOREACHRLEASE ####'
432      write(*,*) '#### MUST BE SET TO ONE!                     ####'
433      stop
434  endif
435
436
437  ! For backward runs one releasefield for all releases makes no sense,
438  ! and is "forbidden"
439  !*****************************************************************************
440
441  if ((mdomainfill.eq.1).and.(ioutputforeachrelease.eq.1)) then
442      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
443      write(*,*) '#### FOR DOMAIN FILLING RUNS OUTPUT FOR      ####'
444      write(*,*) '#### EACH RELEASE IS FORBIDDEN !             ####'
445      stop
446  endif
447
448
449  ! For domain-filling trajectories, a plume centroid trajectory makes no sense,
450  ! For backward runs, only residence time output (iout=1) or plume trajectories (iout=4),
451  ! or both (iout=5) makes sense; other output options are "forbidden"
452  !*****************************************************************************
453
454  if (ldirect.lt.0) then
455    if ((iout.eq.2).or.(iout.eq.3)) then
456      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
457      write(*,*) '#### FOR BACKWARD RUNS, IOUT MUST BE 1,4,OR 5####'
458      stop
459    endif
460  endif
461
462
463  ! For domain-filling trajectories, a plume centroid trajectory makes no sense,
464  ! and is "forbidden"
465  !*****************************************************************************
466
467  if (mdomainfill.ge.1) then
468    if ((iout.eq.4).or.(iout.eq.5)) then
469      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
470      write(*,*) '#### FOR DOMAIN-FILLING TRAJECTORY OPTION,   ####'
471      write(*,*) '#### IOUT MUST NOT BE SET TO 4 OR 5.         ####'
472      stop
473    endif
474  endif
475
476
477
478  ! Check whether a valid options for particle dump has been chosen
479  !****************************************************************
480
481  if ((ipout.ne.0).and.(ipout.ne.1).and.(ipout.ne.2)) then
482    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
483    write(*,*) ' #### IPOUT MUST BE 1, 2 OR 3!                #### '
484    stop
485  endif
486
487  if(lsubgrid.ne.1.and.verbosity.eq.0) then
488    write(*,*) '             ----------------               '
489    write(*,*) ' INFORMATION: SUBGRIDSCALE TERRAIN EFFECT IS'
490    write(*,*) ' NOT PARAMETERIZED DURING THIS SIMULATION.  '
491    write(*,*) '             ----------------               '
492  endif
493
494
495  ! Check whether convection scheme is either turned on or off
496  !***********************************************************
497
498  if ((lconvection.ne.0).and.(lconvection.ne.1)) then
499    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
500    write(*,*) ' #### LCONVECTION MUST BE SET TO EITHER 1 OR 0#### '
501    stop
502  endif
503
504
505  ! Check whether synchronisation interval is sufficiently short
506  !*************************************************************
507
508  if (lsynctime.gt.(idiffnorm/2)) then
509    write(*,*) ' #### FLEXPART MODEL ERROR! SYNCHRONISATION   #### '
510    write(*,*) ' #### TIME IS TOO LONG. MAKE IT SHORTER.      #### '
511    write(*,*) ' #### MINIMUM HAS TO BE: ', idiffnorm/2
512    stop
513  endif
514
515
516  ! Check consistency of the intervals, sampling periods, etc., for model output
517  !*****************************************************************************
518
519  if (loutaver.eq.0) then
520    write(*,*) ' #### FLEXPART MODEL ERROR! TIME AVERAGE OF   #### '
521    write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE  #### '
522    write(*,*) ' #### ZERO.                                   #### '
523    write(*,*) ' #### CHANGE INPUT IN FILE COMMAND.           #### '
524    stop
525  endif
526
527  if (loutaver.gt.loutstep) then
528    write(*,*) ' #### FLEXPART MODEL ERROR! TIME AVERAGE OF   #### '
529    write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE  #### '
530    write(*,*) ' #### GREATER THAN INTERVAL OF OUTPUT.        #### '
531    write(*,*) ' #### CHANGE INPUT IN FILE COMMAND.           #### '
532    stop
533  endif
534
535  if (loutsample.gt.loutaver) then
536    write(*,*) ' #### FLEXPART MODEL ERROR! SAMPLING TIME OF  #### '
537    write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE  #### '
538    write(*,*) ' #### GREATER THAN TIME AVERAGE OF OUTPUT.    #### '
539    write(*,*) ' #### CHANGE INPUT IN FILE COMMAND.           #### '
540    stop
541  endif
542
543  if (mod(loutaver,lsynctime).ne.0) then
544    write(*,*) ' #### FLEXPART MODEL ERROR! AVERAGING TIME OF #### '
545    write(*,*) ' #### CONCENTRATION FIELD MUST BE A MULTIPLE  #### '
546    write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL         #### '
547    stop
548  endif
549
550  if ((loutaver/lsynctime).lt.2) then
551    write(*,*) ' #### FLEXPART MODEL ERROR! AVERAGING TIME OF #### '
552    write(*,*) ' #### CONCENTRATION FIELD MUST BE AT LEAST    #### '
553    write(*,*) ' #### TWICE THE SYNCHRONISATION INTERVAL      #### '
554    stop
555  endif
556
557  if (mod(loutstep,lsynctime).ne.0) then
558    write(*,*) ' #### FLEXPART MODEL ERROR! INTERVAL BETWEEN  #### '
559    write(*,*) ' #### CONCENTRATION FIELDS MUST BE A MULTIPLE #### '
560    write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL         #### '
561    stop
562  endif
563
564  if ((loutstep/lsynctime).lt.2) then
565    write(*,*) ' #### FLEXPART MODEL ERROR! INTERVAL BETWEEN  #### '
566    write(*,*) ' #### CONCENTRATION FIELDS MUST BE AT LEAST   #### '
567    write(*,*) ' #### TWICE THE SYNCHRONISATION INTERVAL      #### '
568    stop
569  endif
570
571  if (mod(loutsample,lsynctime).ne.0) then
572    write(*,*) ' #### FLEXPART MODEL ERROR! SAMPLING TIME OF  #### '
573    write(*,*) ' #### CONCENTRATION FIELD MUST BE A MULTIPLE  #### '
574    write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL         #### '
575    stop
576  endif
577
578  if (itsplit.lt.loutaver) then
579    write(*,*) ' #### FLEXPART MODEL ERROR! SPLITTING TIME FOR#### '
580    write(*,*) ' #### PARTICLES IS TOO SHORT. PLEASE INCREASE #### '
581    write(*,*) ' #### SPLITTING TIME CONSTANT.                #### '
582    stop
583  endif
584
585  if ((mquasilag.eq.1).and.(iout.ge.4)) then
586    write(*,*) ' #### FLEXPART MODEL ERROR! CONFLICTING       #### '
587    write(*,*) ' #### OPTIONS: IF MQUASILAG=1, PLUME          #### '
588    write(*,*) ' #### TRAJECTORY OUTPUT IS IMPOSSIBLE.        #### '
589    stop
590  endif
591
592  ! Compute modeling time in seconds and beginning date in Julian date
593  !*******************************************************************
594
595  outstep=real(abs(loutstep))
596  if (ldirect.eq.1) then
597    bdate=juldate(ibdate,ibtime)
598    edate=juldate(iedate,ietime)
599    ideltas=nint((edate-bdate)*86400.)
600  else if (ldirect.eq.-1) then
601    loutaver=-1*loutaver
602    loutstep=-1*loutstep
603    loutsample=-1*loutsample
604    lsynctime=-1*lsynctime
605    bdate=juldate(iedate,ietime)
606    edate=juldate(ibdate,ibtime)
607    ideltas=nint((edate-bdate)*86400.)
608  else
609    write(*,*) ' #### FLEXPART MODEL ERROR! DIRECTION IN      #### '
610    write(*,*) ' #### FILE "COMMAND" MUST BE EITHER -1 OR 1.  #### '
611    stop
612  endif
613
614  return
615
616999   write(*,*) ' #### FLEXPART MODEL ERROR! FILE "COMMAND"    #### '
617  write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY       #### '
618  write(*,'(a)') path(1)(1:length(1))
619  stop
620
6211000   write(*,*) ' #### FLEXPART MODEL ERROR! FILE "COMMAND"    #### '
622  write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY       #### '
623  write(*,'(a)') path(2)(1:length(1))
624  stop
625end subroutine readcommand
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG