source: flexpart.git/src/readcommand.f90 @ 2eefa58

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

Added Ronas changes for inversion output

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