source: flexpart.git/src/readcommand.f90 @ 660bcb7

NetCDF
Last change on this file since 660bcb7 was 660bcb7, checked in by Anne Fouilloux <annefou@…>, 9 years ago

NETCDF outputs from svn trunk (hasod: ADD: netcdf module file)
I did not take changes in advance.f90; it will be added later.
I also changed OPENs where status was set to new and access is set to APPEND (had problems on abel.uio.no with intel compilers 2013.sp1.3)
It should contains technical changes only and results should be identical.

  • Property mode set to 100644
File size: 22.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
114  ! Presetting namelist command
115  ldirect=0
116  ibdate=20000101
117  ibtime=0
118  iedate=20000102
119  ietime=0
120  loutstep=10800
121  loutaver=10800
122  loutsample=900
123  itsplit=999999999
124  lsynctime=900
125  ctl=-5.0
126  ifine=4
127  iout=3
128  ipout=0
129  lsubgrid=1
130  lconvection=1
131  lagespectra=0
132  ipin=1
133  ioutputforeachrelease=1
134  iflux=1
135  mdomainfill=0
136  ind_source=1
137  ind_receptor=1
138  mquasilag=0
139  nested_output=0
140  linit_cond=0
141  lnetcdfout=0
142  surf_only=0
143
144  ! Open the command file and read user options
145  ! Namelist input first: try to read as namelist file
146  !**************************************************************************
147  open(unitcommand,file=path(1)(1:length(1))//'COMMAND',status='old',form='formatted',err=999)
148
149  ! try namelist input (default)
150  read(unitcommand,command,iostat=readerror)
151  close(unitcommand)
152
153  ! distinguish namelist from fixed text input
154  if ((readerror.ne.0).or.(ldirect.eq.0)) then ! parse as text file format
155 
156    open(unitcommand,file=path(1)(1:length(1))//'COMMAND',status='old', err=999)
157
158    ! Check the format of the COMMAND file (either in free format,
159    ! or using formatted mask)
160    ! Use of formatted mask is assumed if line 10 contains the word 'DIRECTION'
161    !**************************************************************************
162
163    call skplin(9,unitcommand)
164    read (unitcommand,901) line
165  901   format (a)
166    if (index(line,'LDIRECT') .eq. 0) then
167      old = .false.
168      write(*,*) 'COMMAND in old short format, please update to namelist format'
169    else
170      old = .true.
171      write(*,*) 'COMMAND in old long format, please update to namelist format'
172    endif
173    rewind(unitcommand)
174
175
176    ! Read parameters
177    !****************
178
179    call skplin(7,unitcommand)
180    if (old) call skplin(1,unitcommand)
181    read(unitcommand,*) ldirect
182    if (old) call skplin(3,unitcommand)
183    read(unitcommand,*) ibdate,ibtime
184    if (old) call skplin(3,unitcommand)
185    read(unitcommand,*) iedate,ietime
186    if (old) call skplin(3,unitcommand)
187    read(unitcommand,*) loutstep
188    if (old) call skplin(3,unitcommand)
189    read(unitcommand,*) loutaver
190    if (old) call skplin(3,unitcommand)
191    read(unitcommand,*) loutsample
192    if (old) call skplin(3,unitcommand)
193    read(unitcommand,*) itsplit
194    if (old) call skplin(3,unitcommand)
195    read(unitcommand,*) lsynctime
196    if (old) call skplin(3,unitcommand)
197    read(unitcommand,*) ctl
198    if (old) call skplin(3,unitcommand)
199    read(unitcommand,*) ifine
200    if (old) call skplin(3,unitcommand)
201    read(unitcommand,*) iout
202    if (old) call skplin(3,unitcommand)
203    read(unitcommand,*) ipout
204    if (old) call skplin(3,unitcommand)
205    read(unitcommand,*) lsubgrid
206    if (old) call skplin(3,unitcommand)
207    read(unitcommand,*) lconvection
208    if (old) call skplin(3,unitcommand)
209    read(unitcommand,*) lagespectra
210    if (old) call skplin(3,unitcommand)
211    read(unitcommand,*) ipin
212    if (old) call skplin(3,unitcommand)
213    read(unitcommand,*) ioutputforeachrelease
214    if (old) call skplin(3,unitcommand)
215    read(unitcommand,*) iflux
216    if (old) call skplin(3,unitcommand)
217    read(unitcommand,*) mdomainfill
218    if (old) call skplin(3,unitcommand)
219    read(unitcommand,*) ind_source
220    if (old) call skplin(3,unitcommand)
221    read(unitcommand,*) ind_receptor
222    if (old) call skplin(3,unitcommand)
223    read(unitcommand,*) mquasilag
224    if (old) call skplin(3,unitcommand)
225    read(unitcommand,*) nested_output
226    if (old) call skplin(3,unitcommand)
227    read(unitcommand,*) linit_cond
228    if (old) call skplin(3,unitcommand)
229    read(unitcommand,*) surf_only
230    close(unitcommand)
231
232  endif ! input format
233
234  ! write command file in namelist format to output directory if requested
235  if (nmlout.eqv..true.) then
236    open(unitcommand,file=path(2)(1:length(2))//'COMMAND.namelist',err=1000)
237    write(unitcommand,nml=command)
238    close(unitcommand)
239  endif
240
241  ifine=max(ifine,1)
242
243  ! Determine how Markov chain is formulated (for w or for w/sigw)
244  !***************************************************************
245
246  if (ctl.ge.0.1) then
247    turbswitch=.true.
248  else
249    turbswitch=.false.
250    ifine=1
251  endif
252  fine=1./real(ifine)
253  ctl=1./ctl
254
255  ! Set the switches required for the various options for input/output units
256  !*************************************************************************
257  !AF Set the switches IND_REL and IND_SAMP for the release and sampling
258  !Af switches for the releasefile:
259  !Af IND_REL =  1 : xmass * rho
260  !Af IND_REL =  0 : xmass * 1
261
262  !Af switches for the conccalcfile:
263  !AF IND_SAMP =  0 : xmass * 1
264  !Af IND_SAMP = -1 : xmass / rho
265
266  !AF IND_SOURCE switches between different units for concentrations at the source
267  !Af   NOTE that in backward simulations the release of computational particles
268  !Af   takes place at the "receptor" and the sampling of p[articles at the "source".
269  !Af          1 = mass units
270  !Af          2 = mass mixing ratio units
271  !Af IND_RECEPTOR switches between different units for concentrations at the receptor
272  !Af          1 = mass units
273  !Af          2 = mass mixing ratio units
274
275  if ( ldirect .eq. 1 ) then  ! FWD-Run
276  !Af set release-switch
277     if (ind_source .eq. 1 ) then !mass
278        ind_rel = 0
279     else ! mass mix
280        ind_rel = 1
281     endif
282  !Af set sampling switch
283     if (ind_receptor .eq. 1) then !mass
284        ind_samp = 0
285     else ! mass mix
286        ind_samp = -1
287     endif
288  elseif (ldirect .eq. -1 ) then !BWD-Run
289  !Af set sampling switch
290     if (ind_source .eq. 1 ) then !mass
291        ind_samp = -1
292     else ! mass mix
293        ind_samp = 0
294     endif
295  !Af set release-switch
296     if (ind_receptor .eq. 1) then !mass
297        ind_rel = 1
298     else ! mass mix
299        ind_rel = 0
300     endif
301  endif
302
303  !*************************************************************
304  ! Check whether valid options have been chosen in file COMMAND
305  !*************************************************************
306
307  ! Check options for initial condition output: Switch off for forward runs
308  !************************************************************************
309
310  if (ldirect.eq.1) linit_cond=0
311  if ((linit_cond.lt.0).or.(linit_cond.gt.2)) then
312    write(*,*) ' #### FLEXPART MODEL ERROR! INVALID OPTION    #### '
313    write(*,*) ' #### FOR LINIT_COND IN FILE "COMMAND".       #### '
314    stop
315  endif
316
317  ! Check input dates
318  !******************
319
320  if (iedate.lt.ibdate) then
321    write(*,*) ' #### FLEXPART MODEL ERROR! BEGINNING DATE    #### '
322    write(*,*) ' #### IS LARGER THAN ENDING DATE. CHANGE      #### '
323    write(*,*) ' #### EITHER POINT 2 OR POINT 3 IN FILE       #### '
324    write(*,*) ' #### "COMMAND".                              #### '
325    stop
326  else if (iedate.eq.ibdate) then
327    if (ietime.lt.ibtime) then
328    write(*,*) ' #### FLEXPART MODEL ERROR! BEGINNING TIME    #### '
329    write(*,*) ' #### IS LARGER THAN ENDING TIME. CHANGE      #### '
330    write(*,*) ' #### EITHER POINT 2 OR POINT 3 IN FILE       #### '
331    write(*,*) ' #### "COMMAND".                              #### '
332    stop
333    endif
334  endif
335
336
337  ! Determine kind of dispersion method
338  !************************************
339
340  if (ctl.gt.0.) then
341    method=1
342    mintime=minstep
343  else
344    method=0
345    mintime=lsynctime
346  endif
347
348!  check for netcdf output switch (use for non-namelist input only!)
349  if (iout.ge.8) then
350     lnetcdfout = 1
351     iout = iout - 8
352  endif
353
354  ! Check whether a valid option for gridded model output has been chosen
355  !**********************************************************************
356
357  if ((iout.lt.1).or.(iout.gt.5)) then
358    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
359    write(*,*) ' #### IOUT MUST BE 1, 2, 3, 4, OR 5 FOR       #### '
360    write(*,*) ' #### STANDARD FLEXPART OUTPUT OR  9 - 13     #### '
361    write(*,*) ' #### FOR NETCDF OUTPUT                       #### '
362    stop
363  endif
364
365  !AF check consistency between units and volume mixing ratio
366  if ( ((iout.eq.2).or.(iout.eq.3)).and. &
367       (ind_source.gt.1 .or.ind_receptor.gt.1) ) then
368    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
369    write(*,*) ' #### VOLUME MIXING RATIO ONLY SUPPORTED      #### '
370    write(*,*) ' #### FOR MASS UNITS (at the moment)          #### '
371    stop
372  endif
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