source: flexpart.git/src/readcommand.f90 @ 20963b1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 20963b1 was 20963b1, checked in by Espen Sollum ATMOS <eso@…>, 6 years ago

Backwards deposition for the MPI version

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