source: flexpart.git/src/readcommand.f90 @ a9cf4b1

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

Made enabling netCDF output during compilation optional

  • 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
395!*******************************
396  if (iout.ge.8) then
397     lnetcdfout = 1
398     iout = iout - 8
399#ifndef USE_NCF
400     write(*,*) 'ERROR: netcdf output not activated during compile time but used in COMMAND file!'
401     write(*,*) 'Please recompile with netcdf library (`make [...] ncf=yes`) or use standard output format.'
402     stop
403#endif
404  endif
405
406  ! Check whether a valid option for gridded model output has been chosen
407  !**********************************************************************
408
409  if ((iout.lt.1).or.(iout.gt.6)) then
410    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
411    write(*,*) ' #### IOUT MUST BE 1, 2, 3, 4 OR 5 FOR        #### '
412    write(*,*) ' #### STANDARD FLEXPART OUTPUT OR  9 - 13     #### '
413    write(*,*) ' #### FOR NETCDF OUTPUT                       #### '
414    stop
415  endif
416
417  !AF check consistency between units and volume mixing ratio
418  if ( ((iout.eq.2).or.(iout.eq.3)).and. &
419       (ind_source.gt.1 .or.ind_receptor.gt.1) ) then
420    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
421    write(*,*) ' #### VOLUME MIXING RATIO ONLY SUPPORTED      #### '
422    write(*,*) ' #### FOR MASS UNITS (at the moment)          #### '
423    stop
424  endif
425
426
427  ! For quasilag output for each release is forbidden
428  !*****************************************************************************
429
430  if ((ioutputforeachrelease.eq.1).and.(mquasilag.eq.1)) then
431      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
432      write(*,*) '#### OUTPUTFOREACHRELEASE AND QUASILAGRANGIAN####'
433      write(*,*) '#### MODE IS NOT POSSIBLE   !                ####'
434      stop
435  endif
436
437
438  ! For quasilag backward is forbidden
439  !*****************************************************************************
440
441  if ((ldirect.lt.0).and.(mquasilag.eq.1)) then
442      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
443      write(*,*) '#### FOR BACKWARD RUNS, QUASILAGRANGIAN MODE ####'
444      write(*,*) '#### IS NOT POSSIBLE   !                     ####'
445      stop
446  endif
447
448
449  ! For backward runs one releasefield for all releases makes no sense,
450  ! For quasilag and domainfill ioutputforechrelease is forbidden
451  !*****************************************************************************
452
453  if ((ldirect.lt.0).and.(ioutputforeachrelease.eq.0)) then
454      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
455      write(*,*) '#### FOR BACKWARD RUNS, IOUTPUTFOREACHRLEASE ####'
456      write(*,*) '#### MUST BE SET TO ONE!                     ####'
457      stop
458  endif
459
460
461  ! For backward runs one releasefield for all releases makes no sense,
462  ! and is "forbidden"
463  !*****************************************************************************
464
465  if ((mdomainfill.eq.1).and.(ioutputforeachrelease.eq.1)) then
466      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
467      write(*,*) '#### FOR DOMAIN FILLING RUNS OUTPUT FOR      ####'
468      write(*,*) '#### EACH RELEASE IS FORBIDDEN !             ####'
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)) then
506    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
507    write(*,*) ' #### IPOUT MUST BE 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