source: branches/petra/src/readcommand.f90 @ 36

Last change on this file since 36 was 36, checked in by pesei, 9 years ago

Implement switch for incremental deposition, see ticket:113 and many small changes, see changelog.txt

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