source: flexpart.git/src/readcommand.f90 @ 3481cc1

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

move license from headers to a different file

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