source: branches/flexpart91_hasod/src_parallel/readcommand.f90 @ 8

Last change on this file since 8 was 8, checked in by hasod, 11 years ago

Added parallel version of Flexpart91

File size: 21.9 KB
RevLine 
[8]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=0
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  ! Failsafe detection with ldirect initial value change
157  if ((readerror.ne.0).or.(ldirect.ne.0)) then
158
159    open(unitcommand,file=path(1)(1:length(1))//'COMMAND',status='old', &
160         err=999)
161
162    ! Check the format of the COMMAND file (either in free format,
163    ! or using formatted mask)
164    ! Use of formatted mask is assumed if line 10 contains the word 'DIRECTION'
165    !**************************************************************************
166
167    call skplin(9,unitcommand)
168    read (unitcommand,901) line
169  901   format (a)
170    if (index(line,'LDIRECT') .eq. 0) then
171      old = .false.
172    else
173      old = .true.
174    endif
175    rewind(unitcommand)
176
177    ! Read parameters
178    !****************
179
180    call skplin(7,unitcommand)
181    if (old) call skplin(1,unitcommand)
182
183    read(unitcommand,*) ldirect
184    if (old) call skplin(3,unitcommand)
185    read(unitcommand,*) ibdate,ibtime
186    if (old) call skplin(3,unitcommand)
187    read(unitcommand,*) iedate,ietime
188    if (old) call skplin(3,unitcommand)
189    read(unitcommand,*) loutstep
190    if (old) call skplin(3,unitcommand)
191    read(unitcommand,*) loutaver
192    if (old) call skplin(3,unitcommand)
193    read(unitcommand,*) loutsample
194    if (old) call skplin(3,unitcommand)
195    read(unitcommand,*) itsplit
196    if (old) call skplin(3,unitcommand)
197    read(unitcommand,*) lsynctime
198    if (old) call skplin(3,unitcommand)
199    read(unitcommand,*) ctl
200    if (old) call skplin(3,unitcommand)
201    read(unitcommand,*) ifine
202    if (old) call skplin(3,unitcommand)
203    read(unitcommand,*) iout
204    if (old) call skplin(3,unitcommand)
205    read(unitcommand,*) ipout
206    if (old) call skplin(3,unitcommand)
207    read(unitcommand,*) lsubgrid
208    if (old) call skplin(3,unitcommand)
209    read(unitcommand,*) lconvection
210    if (old) call skplin(3,unitcommand)
211    read(unitcommand,*) lagespectra
212    if (old) call skplin(3,unitcommand)
213    read(unitcommand,*) ipin
214    if (old) call skplin(3,unitcommand)
215    read(unitcommand,*) ioutputforeachrelease
216    if (old) call skplin(3,unitcommand)
217    read(unitcommand,*) iflux
218    if (old) call skplin(3,unitcommand)
219    read(unitcommand,*) mdomainfill
220    if (old) call skplin(3,unitcommand)
221    read(unitcommand,*) ind_source
222    if (old) call skplin(3,unitcommand)
223    read(unitcommand,*) ind_receptor
224    if (old) call skplin(3,unitcommand)
225    read(unitcommand,*) mquasilag
226    if (old) call skplin(3,unitcommand)
227    read(unitcommand,*) nested_output
228    if (old) call skplin(3,unitcommand)
229    read(unitcommand,*) linit_cond
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',status='new',err=999)
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
349  lncdfout = .false.
350  if (iout.ge.8) then
351     lncdfout = .true.
352     iout = iout - 8
353#ifndef NETCDF_OUTPUT
354     print*,'ERROR: netcdf output not activated during compile time'
355#endif
356  endif
357
358  ! Check whether a valid option for gridded model output has been chosen
359  !**********************************************************************
360
361  if ((iout.lt.1).or.(iout.gt.5)) then
362    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
363    write(*,*) ' #### IOUT MUST BE 1, 2, 3, 4, OR 5 FOR       #### '
364    write(*,*) ' #### STANDARD FLEXPART OUTPUT OR  9 - 13     #### '
365    write(*,*) ' #### FOR NETCDF OUTPUT                       #### '
366    stop
367  endif
368
369  !AF check consistency between units and volume mixing ratio
370  if ( ((iout.eq.2).or.(iout.eq.3)).and. &
371       (ind_source.gt.1 .or.ind_receptor.gt.1) ) then
372    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
373    write(*,*) ' #### VOLUME MIXING RATIO ONLY SUPPORTED      #### '
374    write(*,*) ' #### FOR MASS UNITS (at the moment)          #### '
375    stop
376  endif
377
378
379
380  ! For quasilag output for each release is forbidden
381  !*****************************************************************************
382
383  if ((ioutputforeachrelease.eq.1).and.(mquasilag.eq.1)) then
384      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
385      write(*,*) '#### OUTPUTFOREACHRELEASE AND QUASILAGRANGIAN####'
386      write(*,*) '#### MODE IS NOT POSSIBLE   !                ####'
387      stop
388  endif
389
390
391  ! For quasilag backward is forbidden
392  !*****************************************************************************
393
394  if ((ldirect.lt.0).and.(mquasilag.eq.1)) then
395      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
396      write(*,*) '#### FOR BACKWARD RUNS, QUASILAGRANGIAN MODE ####'
397      write(*,*) '#### IS NOT POSSIBLE   !                     ####'
398      stop
399  endif
400
401
402  ! For backward runs one releasefield for all releases makes no sense,
403  ! For quasilag and domainfill ioutputforechrelease is forbidden
404  !*****************************************************************************
405
406  if ((ldirect.lt.0).and.(ioutputforeachrelease.eq.0)) then
407      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
408      write(*,*) '#### FOR BACKWARD RUNS, IOUTPUTFOREACHRLEASE ####'
409      write(*,*) '#### MUST BE SET TO ONE!                     ####'
410      stop
411  endif
412
413
414  ! For backward runs one releasefield for all releases makes no sense,
415  ! and is "forbidden"
416  !*****************************************************************************
417
418  if ((mdomainfill.eq.1).and.(ioutputforeachrelease.eq.1)) then
419      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
420      write(*,*) '#### FOR DOMAIN FILLING RUNS OUTPUT FOR      ####'
421      write(*,*) '#### EACH RELEASE IS FORBIDDEN !             ####'
422      stop
423  endif
424
425
426  ! For domain-filling trajectories, a plume centroid trajectory makes no sense,
427  ! For backward runs, only residence time output (iout=1) or plume trajectories (iout=4),
428  ! or both (iout=5) makes sense; other output options are "forbidden"
429  !*****************************************************************************
430
431  if (ldirect.lt.0) then
432    if ((iout.eq.2).or.(iout.eq.3)) then
433      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
434      write(*,*) '#### FOR BACKWARD RUNS, IOUT MUST BE 1,4,OR 5####'
435      stop
436    endif
437  endif
438
439
440  ! For domain-filling trajectories, a plume centroid trajectory makes no sense,
441  ! and is "forbidden"
442  !*****************************************************************************
443
444  if (mdomainfill.ge.1) then
445    if ((iout.eq.4).or.(iout.eq.5)) then
446      write(*,*) '#### FLEXPART MODEL ERROR! FILE COMMAND:     ####'
447      write(*,*) '#### FOR DOMAIN-FILLING TRAJECTORY OPTION,   ####'
448      write(*,*) '#### IOUT MUST NOT BE SET TO 4 OR 5.         ####'
449      stop
450    endif
451  endif
452
453
454
455  ! Check whether a valid options for particle dump has been chosen
456  !****************************************************************
457
458  if ((ipout.ne.0).and.(ipout.ne.1).and.(ipout.ne.2)) then
459    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
460    write(*,*) ' #### IPOUT MUST BE 1, 2 OR 3!                #### '
461    stop
462  endif
463
464  if(lsubgrid.ne.1) then
465    write(*,*) '             ----------------               '
466    write(*,*) ' INFORMATION: SUBGRIDSCALE TERRAIN EFFECT IS'
467    write(*,*) ' NOT PARAMETERIZED DURING THIS SIMULATION.  '
468    write(*,*) '             ----------------               '
469  endif
470
471
472  ! Check whether convection scheme is either turned on or off
473  !***********************************************************
474
475  if ((lconvection.ne.0).and.(lconvection.ne.1)) then
476    write(*,*) ' #### FLEXPART MODEL ERROR! FILE COMMAND:     #### '
477    write(*,*) ' #### LCONVECTION MUST BE SET TO EITHER 1 OR 0#### '
478    stop
479  endif
480
481
482  ! Check whether synchronisation interval is sufficiently short
483  !*************************************************************
484
485  if (lsynctime.gt.(idiffnorm/2)) then
486    write(*,*) ' #### FLEXPART MODEL ERROR! SYNCHRONISATION   #### '
487    write(*,*) ' #### TIME IS TOO LONG. MAKE IT SHORTER.      #### '
488    write(*,*) ' #### MINIMUM HAS TO BE: ', idiffnorm/2
489    stop
490  endif
491
492
493  ! Check consistency of the intervals, sampling periods, etc., for model output
494  !*****************************************************************************
495
496  if (loutaver.eq.0) then
497    write(*,*) ' #### FLEXPART MODEL ERROR! TIME AVERAGE OF   #### '
498    write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE  #### '
499    write(*,*) ' #### ZERO.                                   #### '
500    write(*,*) ' #### CHANGE INPUT IN FILE COMMAND.           #### '
501    stop
502  endif
503
504  if (loutaver.gt.loutstep) then
505    write(*,*) ' #### FLEXPART MODEL ERROR! TIME AVERAGE OF   #### '
506    write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE  #### '
507    write(*,*) ' #### GREATER THAN INTERVAL OF OUTPUT.        #### '
508    write(*,*) ' #### CHANGE INPUT IN FILE COMMAND.           #### '
509    stop
510  endif
511
512  if (loutsample.gt.loutaver) then
513    write(*,*) ' #### FLEXPART MODEL ERROR! SAMPLING TIME OF  #### '
514    write(*,*) ' #### CONCENTRATION FIELD OUTPUT MUST NOT BE  #### '
515    write(*,*) ' #### GREATER THAN TIME AVERAGE OF OUTPUT.    #### '
516    write(*,*) ' #### CHANGE INPUT IN FILE COMMAND.           #### '
517    stop
518  endif
519
520  if (mod(loutaver,lsynctime).ne.0) then
521    write(*,*) ' #### FLEXPART MODEL ERROR! AVERAGING TIME OF #### '
522    write(*,*) ' #### CONCENTRATION FIELD MUST BE A MULTIPLE  #### '
523    write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL         #### '
524    stop
525  endif
526
527  if ((loutaver/lsynctime).lt.2) then
528    write(*,*) ' #### FLEXPART MODEL ERROR! AVERAGING TIME OF #### '
529    write(*,*) ' #### CONCENTRATION FIELD MUST BE AT LEAST    #### '
530    write(*,*) ' #### TWICE THE SYNCHRONISATION INTERVAL      #### '
531    stop
532  endif
533
534  if (mod(loutstep,lsynctime).ne.0) then
535    write(*,*) ' #### FLEXPART MODEL ERROR! INTERVAL BETWEEN  #### '
536    write(*,*) ' #### CONCENTRATION FIELDS MUST BE A MULTIPLE #### '
537    write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL         #### '
538    stop
539  endif
540
541  if ((loutstep/lsynctime).lt.2) then
542    write(*,*) ' #### FLEXPART MODEL ERROR! INTERVAL BETWEEN  #### '
543    write(*,*) ' #### CONCENTRATION FIELDS MUST BE AT LEAST   #### '
544    write(*,*) ' #### TWICE THE SYNCHRONISATION INTERVAL      #### '
545    stop
546  endif
547
548  if (mod(loutsample,lsynctime).ne.0) then
549    write(*,*) ' #### FLEXPART MODEL ERROR! SAMPLING TIME OF  #### '
550    write(*,*) ' #### CONCENTRATION FIELD MUST BE A MULTIPLE  #### '
551    write(*,*) ' #### OF THE SYNCHRONISATION INTERVAL         #### '
552    stop
553  endif
554
555  if (itsplit.lt.loutaver) then
556    write(*,*) ' #### FLEXPART MODEL ERROR! SPLITTING TIME FOR#### '
557    write(*,*) ' #### PARTICLES IS TOO SHORT. PLEASE INCREASE #### '
558    write(*,*) ' #### SPLITTING TIME CONSTANT.                #### '
559    stop
560  endif
561
562  if ((mquasilag.eq.1).and.(iout.ge.4)) then
563    write(*,*) ' #### FLEXPART MODEL ERROR! CONFLICTING       #### '
564    write(*,*) ' #### OPTIONS: IF MQUASILAG=1, PLUME          #### '
565    write(*,*) ' #### TRAJECTORY OUTPUT IS IMPOSSIBLE.        #### '
566    stop
567  endif
568
569  ! Compute modeling time in seconds and beginning date in Julian date
570  !*******************************************************************
571
572  outstep=real(abs(loutstep))
573  if (ldirect.eq.1) then
574    bdate=juldate(ibdate,ibtime)
575    edate=juldate(iedate,ietime)
576    ideltas=nint((edate-bdate)*86400.)
577  else if (ldirect.eq.-1) then
578    loutaver=-1*loutaver
579    loutstep=-1*loutstep
580    loutsample=-1*loutsample
581    lsynctime=-1*lsynctime
582    bdate=juldate(iedate,ietime)
583    edate=juldate(ibdate,ibtime)
584    ideltas=nint((edate-bdate)*86400.)
585  else
586    write(*,*) ' #### FLEXPART MODEL ERROR! DIRECTION IN      #### '
587    write(*,*) ' #### FILE "COMMAND" MUST BE EITHER -1 OR 1.  #### '
588    stop
589  endif
590
591  return
592
593999   write(*,*) ' #### FLEXPART MODEL ERROR! FILE "COMMAND"    #### '
594  write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY       #### '
595  write(*,'(a)') path(1)(1:length(1))
596  stop
597
598end subroutine readcommand
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG