source: trunk/src/readcommand.f90 @ 20

Last change on this file since 20 was 20, checked in by igpis, 10 years ago

move version 9.1.8 form branches to trunk. Contributions from HSO, saeck, pesei, NIK, RT, XKF, IP and others

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