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

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

Added code, makefile for dev branch

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