source: flexpart.git/src/readcommand.f90 @ 77778f8

univie
Last change on this file since 77778f8 was 77778f8, checked in by pesei <petra seibert at univie ac at>, 6 years ago

Introduce changelog.txt, update version string, make makefile more useful for use outside NILU, improve error msg in readcomma$

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