source: flexpart.git/src/readcommand.f90 @ 027e844

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

Fixed an inconsistency (serial vs parallel) with domain-filling option

  • Property mode set to 100644
File size: 24.4 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         write(*,*) ' #### FLEXPART WET DEPOSITION BACKWARD MODE    #### '
335         write(*,*) ' #### Releaseheight is forced to 0 - 20km      #### '
336         write(*,*) ' #### Release is performed above ground lev    #### '
337         WETBKDEP=.true.
338         allocate(xscav_frac1(maxpart,maxspec))
339     case (4)  ! 4 .. dry deposition in outputfield
340         ind_rel = 4
341         write(*,*) ' #### FLEXPART DRY DEPOSITION BACKWARD MODE    #### '
342         write(*,*) ' #### Releaseheight is forced to 0 - 2*href    #### '
343         write(*,*) ' #### Release is performed above ground lev    #### '
344         DRYBKDEP=.true.
345         allocate(xscav_frac1(maxpart,maxspec))
346     end select
347  endif
348
349  !*************************************************************
350  ! Check whether valid options have been chosen in file COMMAND
351  !*************************************************************
352
353  ! Check options for initial condition output: Switch off for forward runs
354  !************************************************************************
355
356  if (ldirect.eq.1) linit_cond=0
357  if ((linit_cond.lt.0).or.(linit_cond.gt.2)) then
358    write(*,*) ' #### FLEXPART MODEL ERROR! INVALID OPTION    #### '
359    write(*,*) ' #### FOR LINIT_COND IN FILE "COMMAND".       #### '
360    stop
361  endif
362
363  ! Check input dates
364  !******************
365
366  if (iedate.lt.ibdate) then
367    write(*,*) ' #### FLEXPART MODEL ERROR! BEGINNING DATE    #### '
368    write(*,*) ' #### IS LARGER THAN ENDING DATE. CHANGE      #### '
369    write(*,*) ' #### EITHER POINT 2 OR POINT 3 IN FILE       #### '
370    write(*,*) ' #### "COMMAND".                              #### '
371    stop
372  else if (iedate.eq.ibdate) then
373    if (ietime.lt.ibtime) then
374    write(*,*) ' #### FLEXPART MODEL ERROR! BEGINNING TIME    #### '
375    write(*,*) ' #### IS LARGER THAN ENDING TIME. CHANGE      #### '
376    write(*,*) ' #### EITHER POINT 2 OR POINT 3 IN FILE       #### '
377    write(*,*) ' #### "COMMAND".                              #### '
378    stop
379    endif
380  endif
381
382
383  ! Determine kind of dispersion method
384  !************************************
385
386  if (ctl.gt.0.) then
387    method=1
388    mintime=minstep
389  else
390    method=0
391    mintime=lsynctime
392  endif
393
394!  check for netcdf output switch (use for non-namelist input only!)
395  if (iout.ge.8) then
396     lnetcdfout = 1
397     iout = iout - 8
398! #ifndef NETCDF_OUTPUT
399!      print*,'ERROR: netcdf output not activated during compile time but used in COMMAND file!'
400!      print*,'Please recompile with netcdf library or use standard output format.'
401!      stop
402! #endif
403  endif
404
405  ! Check whether a valid option for gridded model output has been chosen
406  !**********************************************************************
407
408  if ((iout.lt.1).or.(iout.gt.6)) then
409    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
410    write(*,*) ' #### IOUT MUST BE 1, 2, 3, 4 OR 5 FOR        #### '
411    write(*,*) ' #### STANDARD FLEXPART OUTPUT OR  9 - 13     #### '
412    write(*,*) ' #### FOR NETCDF OUTPUT                       #### '
413    stop
414  endif
415
416  !AF check consistency between units and volume mixing ratio
417  if ( ((iout.eq.2).or.(iout.eq.3)).and. &
418       (ind_source.gt.1 .or.ind_receptor.gt.1) ) then
419    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
420    write(*,*) ' #### VOLUME MIXING RATIO ONLY SUPPORTED      #### '
421    write(*,*) ' #### FOR MASS UNITS (at the moment)          #### '
422    stop
423  endif
424
425
426  ! For quasilag output for each release is forbidden
427  !*****************************************************************************
428
429  if ((ioutputforeachrelease.eq.1).and.(mquasilag.eq.1)) then
430      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
431      write(*,*) '#### OUTPUTFOREACHRELEASE AND QUASILAGRANGIAN####'
432      write(*,*) '#### MODE IS NOT POSSIBLE   !                ####'
433      stop
434  endif
435
436
437  ! For quasilag backward is forbidden
438  !*****************************************************************************
439
440  if ((ldirect.lt.0).and.(mquasilag.eq.1)) then
441      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
442      write(*,*) '#### FOR BACKWARD RUNS, QUASILAGRANGIAN MODE ####'
443      write(*,*) '#### IS NOT POSSIBLE   !                     ####'
444      stop
445  endif
446
447
448  ! For backward runs one releasefield for all releases makes no sense,
449  ! For quasilag and domainfill ioutputforechrelease is forbidden
450  !*****************************************************************************
451
452  if ((ldirect.lt.0).and.(ioutputforeachrelease.eq.0)) then
453      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
454      write(*,*) '#### FOR BACKWARD RUNS, IOUTPUTFOREACHRLEASE ####'
455      write(*,*) '#### MUST BE SET TO ONE!                     ####'
456      stop
457  endif
458
459
460  ! For backward runs one releasefield for all releases makes no sense,
461  ! and is "forbidden"
462  !*****************************************************************************
463
464  if ((mdomainfill.eq.1).and.(ioutputforeachrelease.eq.1)) then
465      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
466      write(*,*) '#### FOR DOMAIN FILLING RUNS OUTPUT FOR      ####'
467      write(*,*) '#### EACH RELEASE IS FORBIDDEN !             ####'
468      stop
469  endif
470
471
472  ! For domain-filling trajectories, a plume centroid trajectory makes no sense,
473  ! For backward runs, only residence time output (iout=1) or plume trajectories (iout=4),
474  ! or both (iout=5) makes sense; other output options are "forbidden"
475  !*****************************************************************************
476
477  if (ldirect.lt.0) then
478    if ((iout.eq.2).or.(iout.eq.3)) then
479      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
480      write(*,*) '#### FOR BACKWARD RUNS, IOUT MUST BE 1,4,OR 5####'
481      stop
482    endif
483  endif
484
485
486  ! For domain-filling trajectories, a plume centroid trajectory makes no sense,
487  ! and is "forbidden"
488  !*****************************************************************************
489
490  if (mdomainfill.ge.1) then
491    if ((iout.eq.4).or.(iout.eq.5)) then
492      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
493      write(*,*) '#### FOR DOMAIN-FILLING TRAJECTORY OPTION,   ####'
494      write(*,*) '#### IOUT MUST NOT BE SET TO 4 OR 5.         ####'
495      stop
496    endif
497  endif
498
499
500
501  ! Check whether a valid options for particle dump has been chosen
502  !****************************************************************
503
504  if ((ipout.ne.0).and.(ipout.ne.1).and.(ipout.ne.2)) then
505    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
506    write(*,*) ' #### IPOUT MUST BE 1, 2 OR 3!                #### '
507    stop
508  endif
509
510  if(lsubgrid.ne.1.and.verbosity.eq.0) then
511    write(*,*) '             ----------------               '
512    write(*,*) ' INFORMATION: SUBGRIDSCALE TERRAIN EFFECT IS'
513    write(*,*) ' NOT PARAMETERIZED DURING THIS SIMULATION.  '
514    write(*,*) '             ----------------               '
515  endif
516
517
518  ! Check whether convection scheme is either turned on or off
519  !***********************************************************
520
521  if ((lconvection.ne.0).and.(lconvection.ne.1)) then
522    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
523    write(*,*) ' #### LCONVECTION MUST BE SET TO EITHER 1 OR 0#### '
524    stop
525  endif
526
527
528  ! Check whether synchronisation interval is sufficiently short
529  !*************************************************************
530
531  if (lsynctime.gt.(idiffnorm/2)) then
532    write(*,*) ' #### FLEXPART MODEL ERROR! SYNCHRONISATION   #### '
533    write(*,*) ' #### TIME IS TOO LONG. MAKE IT SHORTER.      #### '
534    write(*,*) ' #### MINIMUM HAS TO BE: ', idiffnorm/2
535    stop
536  endif
537
538
539  ! Check consistency of the intervals, sampling periods, etc., for model output
540  !*****************************************************************************
541
542  if (loutaver.eq.0) then
543    write(*,*) ' #### FLEXPART MODEL ERROR! TIME AVERAGE OF   #### '
544    write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE  #### '
545    write(*,*) ' #### ZERO.                                   #### '
546    write(*,*) ' #### CHANGE INPUT IN FILE COMMAND.           #### '
547    stop
548  endif
549
550  if (loutaver.gt.loutstep) then
551    write(*,*) ' #### FLEXPART MODEL ERROR! TIME AVERAGE OF   #### '
552    write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE  #### '
553    write(*,*) ' #### GREATER THAN INTERVAL OF OUTPUT.        #### '
554    write(*,*) ' #### CHANGE INPUT IN FILE COMMAND.           #### '
555    stop
556  endif
557
558  if (loutsample.gt.loutaver) then
559    write(*,*) ' #### FLEXPART MODEL ERROR! SAMPLING TIME OF  #### '
560    write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE  #### '
561    write(*,*) ' #### GREATER THAN TIME AVERAGE OF OUTPUT.    #### '
562    write(*,*) ' #### CHANGE INPUT IN FILE COMMAND.           #### '
563    stop
564  endif
565
566  if (mod(loutaver,lsynctime).ne.0) then
567    write(*,*) ' #### FLEXPART MODEL ERROR! AVERAGING TIME OF #### '
568    write(*,*) ' #### CONCENTRATION FIELD MUST BE A MULTIPLE  #### '
569    write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL         #### '
570    stop
571  endif
572
573  if ((loutaver/lsynctime).lt.2) then
574    write(*,*) ' #### FLEXPART MODEL ERROR! AVERAGING TIME OF #### '
575    write(*,*) ' #### CONCENTRATION FIELD MUST BE AT LEAST    #### '
576    write(*,*) ' #### TWICE THE SYNCHRONISATION INTERVAL      #### '
577    stop
578  endif
579
580  if (mod(loutstep,lsynctime).ne.0) then
581    write(*,*) ' #### FLEXPART MODEL ERROR! INTERVAL BETWEEN  #### '
582    write(*,*) ' #### CONCENTRATION FIELDS MUST BE A MULTIPLE #### '
583    write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL         #### '
584    stop
585  endif
586
587  if ((loutstep/lsynctime).lt.2) then
588    write(*,*) ' #### FLEXPART MODEL ERROR! INTERVAL BETWEEN  #### '
589    write(*,*) ' #### CONCENTRATION FIELDS MUST BE AT LEAST   #### '
590    write(*,*) ' #### TWICE THE SYNCHRONISATION INTERVAL      #### '
591    stop
592  endif
593
594  if (mod(loutsample,lsynctime).ne.0) then
595    write(*,*) ' #### FLEXPART MODEL ERROR! SAMPLING TIME OF  #### '
596    write(*,*) ' #### CONCENTRATION FIELD MUST BE A MULTIPLE  #### '
597    write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL         #### '
598    stop
599  endif
600
601  if (itsplit.lt.loutaver) then
602    write(*,*) ' #### FLEXPART MODEL ERROR! SPLITTING TIME FOR#### '
603    write(*,*) ' #### PARTICLES IS TOO SHORT. PLEASE INCREASE #### '
604    write(*,*) ' #### SPLITTING TIME CONSTANT.                #### '
605    stop
606  endif
607
608  if ((mquasilag.eq.1).and.(iout.ge.4)) then
609    write(*,*) ' #### FLEXPART MODEL ERROR! CONFLICTING       #### '
610    write(*,*) ' #### OPTIONS: IF MQUASILAG=1, PLUME          #### '
611    write(*,*) ' #### TRAJECTORY OUTPUT IS IMPOSSIBLE.        #### '
612    stop
613  endif
614
615  ! Compute modeling time in seconds and beginning date in Julian date
616  !*******************************************************************
617
618  outstep=real(abs(loutstep))
619  if (ldirect.eq.1) then
620    bdate=juldate(ibdate,ibtime)
621    edate=juldate(iedate,ietime)
622    ideltas=nint((edate-bdate)*86400.)
623  else if (ldirect.eq.-1) then
624    loutaver=-1*loutaver
625    loutstep=-1*loutstep
626    loutsample=-1*loutsample
627    lsynctime=-1*lsynctime
628    bdate=juldate(iedate,ietime)
629    edate=juldate(ibdate,ibtime)
630    ideltas=nint((edate-bdate)*86400.)
631  else
632    write(*,*) ' #### FLEXPART MODEL ERROR! DIRECTION IN      #### '
633    write(*,*) ' #### FILE "COMMAND" MUST BE EITHER -1 OR 1.  #### '
634    stop
635  endif
636
637  return
638
639999   write(*,*) ' #### FLEXPART MODEL ERROR! FILE "COMMAND"    #### '
640  write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY       #### '
641  write(*,'(a)') path(1)(1:length(1))
642  stop
643
6441000   write(*,*) ' #### FLEXPART MODEL ERROR! FILE "COMMAND"    #### '
645  write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY       #### '
646  write(*,'(a)') path(2)(1:length(1))
647  stop
648end subroutine readcommand
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG