source: flexpart.git/src/readcommand.f90 @ 8ee24a5

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

force releaseheights for bkwd depo runs

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