source: trunk/src/readcommand.f90 @ 28

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