source: flexpart.git/src/readcommand.f90 @ 6985a98

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 6985a98 was 6985a98, checked in by Sabine <sabine.eckhardt@…>, 7 years ago

compiles after merge scavenging into test dev

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