source: flexpart.git/src/readcommand.f90

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file was 92fab65, checked in by Ignacio Pisso <ip@…>, 4 years ago

add SPDX-License-Identifier to all .f90 files

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