source: branches/FLEXPART_9.1.3/src/readcommand.f90 @ 13

Last change on this file since 13 was 13, checked in by saeck, 11 years ago

update to wetdepo.f90

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