source: flexpart.git/src/readcommand.f90 @ adead08

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since adead08 was 25b4532, checked in by Ignacio Pisso <ip@…>, 5 years ago

turbulence factors can change for different runs with the same executable

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