source: branches/flexpart91_hasod/src_parallel/readreleases.f90

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

add verbose mode to version 9.1.7.1

  • Property svn:executable set to *
File size: 20.8 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 readreleases
23
24  !*****************************************************************************
25  !                                                                            *
26  !     This routine reads the release point specifications for the current    *
27  !     model run. Several release points can be used at the same time.        *
28  !                                                                            *
29  !     Author: A. Stohl                                                       *
30  !                                                                            *
31  !     18 May 1996                                                            *
32  !                                                                            *
33  !     Update: 29 January 2001                                                *
34  !     Release altitude can be either in magl or masl                         *
35  !     HSO, 12 August 2013
36  !     Added optional namelist input
37  !                                                                            *
38  !*****************************************************************************
39  !                                                                            *
40  ! Variables:                                                                 *
41  ! decay               decay constant of species                              *
42  ! dquer [um]          mean particle diameters                                *
43  ! dsigma              e.g. dsigma=10 or dsigma=0.1 means that 68% of the mass*
44  !                     are between 0.1*dquer and 10*dquer                     *
45  ! ireleasestart, ireleaseend [s] starting time and ending time of each       *
46  !                     release                                                *
47  ! kindz               1: zpoint is in m agl, 2: zpoint is in m asl, 3: zpoint*
48  !                     is in hPa                                              *
49  ! npart               number of particles to be released                     *
50  ! nspec               number of species to be released                       *
51  ! density [kg/m3]     density of the particles                               *
52  ! rm [s/m]            Mesophyll resistance                                   *
53  ! species             name of species                                        *
54  ! xmass               total mass of each species                             *
55  ! xpoint1,ypoint1     geograf. coordinates of lower left corner of release   *
56  !                     area                                                   *
57  ! xpoint2,ypoint2     geograf. coordinates of upper right corner of release  *
58  !                     area                                                   *
59  ! weta, wetb          parameters to determine the wet scavenging coefficient *
60  ! zpoint1,zpoint2     height range, over which release takes place           *
61  !                                                                            *
62  !*****************************************************************************
63
64  use point_mod
65  use xmass_mod
66  use par_mod
67  use com_mod
68
69  implicit none
70
71  integer :: numpartmax,i,j,id1,it1,id2,it2,idum,stat
72  integer :: specnum_rel(maxspec)
73  real :: vsh(ni),fracth(ni),schmih(ni),releaserate,xdum,cun
74  real(kind=dp) :: jul1,jul2,juldate
75  character(len=50) :: line
76  logical :: old
77
78  ! help variables for namelist reading
79  integer :: numpoints, parts, readerror
80  integer*2 :: zkind
81  integer :: idate1, itime1, idate2, itime2
82  real :: lon1,lon2,lat1,lat2,z1,z2 !,mass(nspec)
83  real :: mass(45)
84  character*40 :: comment
85
86  ! declare namelists
87  namelist /releases_ctrl/ &
88    nspec, &
89    specnum_rel
90
91  namelist /release/ &
92    idate1, itime1, &
93    idate2, itime2, &
94    lon1, lon2, &
95    lat1, lat2, &
96    z1, z2, &
97    zkind, &
98    mass, &
99    parts, &
100    comment
101
102  numpoint=0
103
104  ! presetting namelist releases_ctrl
105  nspec = -1  ! use negative value to determine failed namelist input
106  specnum_rel(1) = 0
107
108  !sec, read release to find how many releasepoints should be allocated
109  open(unitreleases,file=path(1)(1:length(1))//'RELEASES',status='old', &
110       err=999)
111
112  ! check if namelist input provided
113  read(unitreleases,releases_ctrl,iostat=readerror)
114
115  if ((readerror.ne.0).or.(nspec.lt.0)) then
116
117  ! no namelist format, close file and allow reopening in old format
118  rewind(unitreleases)
119
120  readerror=1 ! indicates old format
121
122  ! Check the format of the RELEASES file (either in free format,
123  ! or using a formatted mask)
124  ! Use of formatted mask is assumed if line 10 contains the word 'DIRECTION'
125  !**************************************************************************
126
127  call skplin(12,unitreleases)
128  read (unitreleases,901) line
129901   format (a)
130  if (index(line,'Total') .eq. 0) then
131    old = .false.
132  else
133    old = .true.
134  endif
135  rewind(unitreleases)
136
137  ! Skip first 11 lines (file header)
138  !**********************************
139
140  call skplin(11,unitreleases)
141  read(unitreleases,*,err=998) nspec
142  if (old) call skplin(2,unitreleases)
143  do i=1,nspec
144    read(unitreleases,*,err=998) specnum_rel(1)
145    if (old) call skplin(2,unitreleases)
146  end do
147
148100   numpoint=numpoint+1
149  read(unitreleases,*,end=25)
150  read(unitreleases,*,err=998,end=25) idum,idum
151  if (old) call skplin(2,unitreleases)
152  read(unitreleases,*,err=998) idum,idum
153  if (old) call skplin(2,unitreleases)
154  read(unitreleases,*,err=998) xdum
155  if (old) call skplin(2,unitreleases)
156  read(unitreleases,*,err=998) xdum
157  if (old) call skplin(2,unitreleases)
158  read(unitreleases,*,err=998) xdum
159  if (old) call skplin(2,unitreleases)
160  read(unitreleases,*,err=998) xdum
161  if (old) call skplin(2,unitreleases)
162  read(unitreleases,*,err=998) idum
163  if (old) call skplin(2,unitreleases)
164  read(unitreleases,*,err=998) xdum
165  if (old) call skplin(2,unitreleases)
166  read(unitreleases,*,err=998) xdum
167  if (old) call skplin(2,unitreleases)
168  read(unitreleases,*,err=998) idum
169  if (old) call skplin(2,unitreleases)
170  do i=1,nspec
171    read(unitreleases,*,err=998) xdum
172    if (old) call skplin(2,unitreleases)
173  end do
174  !save compoint only for the first 1000 release points
175  read(unitreleases,'(a40)',err=998) compoint(1)(1:40)
176  if (old) call skplin(1,unitreleases)
177
178  goto 100
179
18025   numpoint=numpoint-1
181
182  else
183    readerror=0
184    do while (readerror.eq.0)
185      idate1=-1
186      read(unitreleases,release,iostat=readerror)
187      if ((idate1.lt.0).or.(readerror.ne.0)) then
188        readerror=1
189      else
190        numpoint=numpoint+1
191      endif
192    end do
193    readerror=0
194  endif ! if namelist input
195
196  rewind(unitreleases)
197
198  !allocate memory for numpoint releaspoint
199  allocate(ireleasestart(numpoint),stat=stat)
200       
201  if (stat.ne.0) write(*,*)'ERROR: could not allocate releasepoint'
202  allocate(ireleaseend(numpoint),stat=stat)
203  if (stat.ne.0) write(*,*)'ERROR: could not allocate releasepoint'
204  allocate(xpoint1(numpoint),stat=stat)
205  if (stat.ne.0) write(*,*)'ERROR: could not allocate releasepoint'
206  allocate(xpoint2(numpoint),stat=stat)
207  if (stat.ne.0) write(*,*)'ERROR: could not allocate releasepoint'
208  allocate(ypoint1(numpoint),stat=stat)
209  if (stat.ne.0) write(*,*)'ERROR: could not allocate releasepoint'
210  allocate(ypoint2(numpoint),stat=stat)
211  if (stat.ne.0) write(*,*)'ERROR: could not allocate releasepoint'
212  allocate(zpoint1(numpoint),stat=stat)
213  if (stat.ne.0) write(*,*)'ERROR: could not allocate releasepoint'
214  allocate(zpoint2(numpoint),stat=stat)
215  if (stat.ne.0) write(*,*)'ERROR: could not allocate releasepoint'
216  allocate(kindz(numpoint),stat=stat)
217  if (stat.ne.0) write(*,*)'ERROR: could not allocate releasepoint'
218  allocate(xmass(numpoint,maxspec),stat=stat)
219  if (stat.ne.0) write(*,*)'ERROR: could not allocate releasepoint'
220  allocate(rho_rel(numpoint),stat=stat)
221  if (stat.ne.0) write(*,*)'ERROR: could not allocate releasepoint'
222  allocate(npart(numpoint),stat=stat)
223  if (stat.ne.0) write(*,*)'ERROR: could not allocate releasepoint'
224  allocate(xmasssave(numpoint),stat=stat)
225  if (stat.ne.0) write(*,*)'ERROR: could not allocate releasepoint'
226
227  write (*,*) numpoint,' releasepoints allocated.'
228
229  do i=1,numpoint
230    xmasssave(i)=0.
231  end do
232
233  !now save the information
234  DEP=.false.
235  DRYDEP=.false.
236  WETDEP=.false.
237  OHREA=.false.
238  do i=1,maxspec
239    DRYDEPSPEC(i)=.false.
240  end do
241
242  if (readerror.ne.0) then
243  ! Skip first 11 lines (file header)
244  !**********************************
245
246    call skplin(11,unitreleases)
247
248  ! Assign species-specific parameters needed for physical processes
249  !*************************************************************************
250
251    read(unitreleases,*,err=998) nspec
252    if (nspec.gt.maxspec) goto 994
253  if (old) call skplin(2,unitreleases)
254  endif
255  do i=1,nspec
256    if (readerror.ne.0) then
257      read(unitreleases,*,err=998) specnum_rel(1)
258      if (old) call skplin(2,unitreleases)
259      call readspecies(specnum_rel(1),i)
260    else
261      call readspecies(specnum_rel(i),i)
262    endif
263
264  ! For backward runs, only 1 species is allowed
265  !*********************************************
266
267  !if ((ldirect.lt.0).and.(nspec.gt.1)) then
268  !write(*,*) '#####################################################'
269  !write(*,*) '#### FLEXPART MODEL SUBROUTINE READRELEASES:     ####'
270  !write(*,*) '#### FOR BACKWARD RUNS, ONLY 1 SPECIES IS ALLOWED####'
271  !write(*,*) '#####################################################'
272  !  stop
273  !endif
274
275  ! Molecular weight
276  !*****************
277
278    if (((iout.eq.2).or.(iout.eq.3)).and. &
279         (weightmolar(i).lt.0.)) then
280      write(*,*) 'For mixing ratio output, valid molar weight'
281      write(*,*) 'must be specified for all simulated species.'
282      write(*,*) 'Check table SPECIES or choose concentration'
283      write(*,*) 'output instead if molar weight is not known.'
284      stop
285    endif
286
287
288  ! Radioactive decay
289  !******************
290
291    decay(i)=0.693147/decay(i) !conversion half life to decay constant
292
293
294  ! Dry deposition of gases
295  !************************
296
297    if (reldiff(i).gt.0.) &
298         rm(i)=1./(henry(i)/3000.+100.*f0(i))    ! mesophyll resistance
299
300  ! Dry deposition of particles
301  !****************************
302
303    vsetaver(i)=0.
304    cunningham(i)=0.
305    dquer(i)=dquer(i)*1000000.         ! Conversion m to um
306    if (density(i).gt.0.) then                  ! Additional parameters
307     call part0(dquer(i),dsigma(i),density(i),fracth,schmih,cun,vsh)
308      do j=1,ni
309        fract(i,j)=fracth(j)
310        schmi(i,j)=schmih(j)
311        vset(i,j)=vsh(j)
312        cunningham(i)=cunningham(i)+cun*fract(i,j)
313        vsetaver(i)=vsetaver(i)-vset(i,j)*fract(i,j)
314      end do
315      write(*,*) 'Average setting velocity: ',i,vsetaver(i)
316    endif
317
318  ! Dry deposition for constant deposition velocity
319  !************************************************
320
321    dryvel(i)=dryvel(i)*0.01         ! conversion to m/s
322
323  ! Check if wet deposition or OH reaction shall be calculated
324  !***********************************************************
325    if (weta(i).gt.0.)  then
326      WETDEP=.true.
327      write (*,*) 'Wetdeposition switched on: ',weta(i),i
328    endif
329    if (ohreact(i).gt.0) then
330      OHREA=.true.
331      write (*,*) 'OHreaction switched on: ',ohreact(i),i
332    endif
333
334
335    if ((reldiff(i).gt.0.).or.(density(i).gt.0.).or. &
336         (dryvel(i).gt.0.)) then
337      DRYDEP=.true.
338      DRYDEPSPEC(i)=.true.
339    endif
340
341  end do
342
343    if (WETDEP.or.DRYDEP) DEP=.true.
344
345  ! Read specifications for each release point
346  !*******************************************
347  numpoints=numpoint
348  numpoint=0
349  numpartmax=0
350  releaserate=0.
3511000   numpoint=numpoint+1
352
353  if (readerror.eq.0) then ! reading namelist format
354
355  if (numpoint.gt.numpoints) goto 250
356  zkind = 1
357  mass = 0
358  parts = 0
359  comment = ' '
360  read(unitreleases,release,iostat=readerror)
361  id1=idate1
362  it1=itime1
363  id2=idate2
364  it2=itime2
365  xpoint1(numpoint)=lon1
366  xpoint2(numpoint)=lon2
367  ypoint1(numpoint)=lat1
368  ypoint2(numpoint)=lat2
369  zpoint1(numpoint)=z1
370  zpoint2(numpoint)=z2
371  kindz(numpoint)=zkind
372  do i=1,nspec
373    xmass(numpoint,i)=mass(i)
374  end do
375  npart(numpoint)=parts
376  compoint(min(1001,numpoint))=comment
377
378  else
379
380  read(unitreleases,*,end=250)
381  read(unitreleases,*,err=998,end=250) id1,it1
382  if (old) call skplin(2,unitreleases)
383  read(unitreleases,*,err=998) id2,it2
384  if (old) call skplin(2,unitreleases)
385  read(unitreleases,*,err=998) xpoint1(numpoint)
386  if (old) call skplin(2,unitreleases)
387  read(unitreleases,*,err=998) ypoint1(numpoint)
388  if (old) call skplin(2,unitreleases)
389  read(unitreleases,*,err=998) xpoint2(numpoint)
390  if (old) call skplin(2,unitreleases)
391  read(unitreleases,*,err=998) ypoint2(numpoint)
392  if (old) call skplin(2,unitreleases)
393  read(unitreleases,*,err=998) kindz(numpoint)
394  if (old) call skplin(2,unitreleases)
395  read(unitreleases,*,err=998) zpoint1(numpoint)
396  if (old) call skplin(2,unitreleases)
397  read(unitreleases,*,err=998) zpoint2(numpoint)
398  if (old) call skplin(2,unitreleases)
399  read(unitreleases,*,err=998) npart(numpoint)
400  if (old) call skplin(2,unitreleases)
401  do i=1,nspec
402    read(unitreleases,*,err=998) xmass(numpoint,i)
403    if (old) call skplin(2,unitreleases)
404  end do
405  !save compoint only for the first 1000 release points
406  if (numpoint.le.1000) then
407    read(unitreleases,'(a40)',err=998) compoint(numpoint)(1:40)
408  else
409    read(unitreleases,'(a40)',err=998) compoint(1001)(1:40)
410  endif
411  if (old) call skplin(1,unitreleases)
412  if (numpoint.le.1000) then
413    if((xpoint1(numpoint).eq.0.).and.(ypoint1(numpoint).eq.0.).and. &
414         (xpoint2(numpoint).eq.0.).and.(ypoint2(numpoint).eq.0.).and. &
415         (compoint(numpoint)(1:8).eq.'        ')) goto 250
416  else
417    if((xpoint1(numpoint).eq.0.).and.(ypoint1(numpoint).eq.0.).and. &
418         (xpoint2(numpoint).eq.0.).and.(ypoint2(numpoint).eq.0.)) goto 250
419  endif
420
421  endif ! if namelist format
422
423  ! If a release point contains no particles, stop and issue error message
424  !***********************************************************************
425
426  if (npart(numpoint).eq.0) then
427    write(*,*) 'FLEXPART MODEL ERROR'
428    write(*,*) 'RELEASES file is corrupt.'
429    write(*,*) 'At least for one release point, there are zero'
430    write(*,*) 'particles released. Make changes to RELEASES.'
431    stop
432  endif
433
434  ! Check whether x coordinates of release point are within model domain
435  !*********************************************************************
436
437   if (xpoint1(numpoint).lt.xlon0) &
438        xpoint1(numpoint)=xpoint1(numpoint)+360.
439   if (xpoint1(numpoint).gt.xlon0+(nxmin1)*dx) &
440        xpoint1(numpoint)=xpoint1(numpoint)-360.
441   if (xpoint2(numpoint).lt.xlon0) &
442        xpoint2(numpoint)=xpoint2(numpoint)+360.
443   if (xpoint2(numpoint).gt.xlon0+(nxmin1)*dx) &
444        xpoint2(numpoint)=xpoint2(numpoint)-360.
445
446  ! Determine relative beginning and ending times of particle release
447  !******************************************************************
448
449  jul1=juldate(id1,it1)
450  jul2=juldate(id2,it2)
451  if (jul1.gt.jul2) then
452    write(*,*) 'FLEXPART MODEL ERROR'
453    write(*,*) 'Release stops before it begins.'
454    write(*,*) 'Make changes to file RELEASES.'
455    stop
456  endif
457  if (mdomainfill.eq.0) then   ! no domain filling
458    if (ldirect.eq.1) then
459      if ((jul1.lt.bdate).or.(jul2.gt.edate)) then
460        write(*,*) 'FLEXPART MODEL ERROR'
461        write(*,*) 'Release starts before simulation begins or ends'
462        write(*,*) 'after simulation stops.'
463        write(*,*) 'Make files COMMAND and RELEASES consistent.'
464        stop
465      endif
466      ireleasestart(numpoint)=int((jul1-bdate)*86400.)
467      ireleaseend(numpoint)=int((jul2-bdate)*86400.)
468    else if (ldirect.eq.-1) then
469      if ((jul1.lt.edate).or.(jul2.gt.bdate)) then
470        write(*,*) 'FLEXPART MODEL ERROR'
471        write(*,*) 'Release starts before simulation begins or ends'
472        write(*,*) 'after simulation stops.'
473        write(*,*) 'Make files COMMAND and RELEASES consistent.'
474        stop
475      endif
476      ireleasestart(numpoint)=int((jul1-bdate)*86400.)
477      ireleaseend(numpoint)=int((jul2-bdate)*86400.)
478    endif
479  endif
480
481
482  ! Check, whether the total number of particles may exceed totally allowed
483  ! number of particles at some time during the simulation
484  !************************************************************************
485
486  ! Determine the release rate (particles per second) and total number
487  ! of particles released during the simulation
488  !*******************************************************************
489
490  if (ireleasestart(numpoint).ne.ireleaseend(numpoint)) then
491    releaserate=releaserate+real(npart(numpoint))/ &
492         real(ireleaseend(numpoint)-ireleasestart(numpoint))
493  else
494    releaserate=99999999
495  endif
496  numpartmax=numpartmax+npart(numpoint)
497  goto 1000
498
499
500250   close(unitreleases)
501
502  write (*,*) ' Particles allocated for this run: ',maxpart, ', released in simulation: ',  numpartmax
503  numpoint=numpoint-1
504
505  if (ioutputforeachrelease.eq.1) then
506    maxpointspec_act=numpoint
507  else
508    maxpointspec_act=1
509  endif
510
511  if (releaserate.gt. &
512       0.99*real(maxpart)/real(lage(nageclass))) then
513    if (numpartmax.gt.maxpart) then
514  write(*,*) '#####################################################'
515  write(*,*) '#### FLEXPART MODEL SUBROUTINE READRELEASES:     ####'
516  write(*,*) '####                                             ####'
517  write(*,*) '####WARNING - TOTAL NUMBER OF PARTICLES SPECIFIED####'
518  write(*,*) '#### IN FILE "RELEASES" MAY AT SOME POINT DURING ####'
519  write(*,*) '#### THE SIMULATION EXCEED THE MAXIMUM ALLOWED   ####'
520  write(*,*) '#### NUMBER (MAXPART).IF RELEASES DO NOT OVERLAP,####'
521  write(*,*) '#### FLEXPART CAN POSSIBLY COMPLETE SUCCESSFULLY.####'
522  write(*,*) '#### HOWEVER, FLEXPART MAY HAVE TO STOP          ####'
523  write(*,*) '#### AT SOME TIME DURING THE SIMULATION. PLEASE  ####'
524  write(*,*) '#### MAKE SURE THAT YOUR SETTINGS ARE CORRECT.   ####'
525  write(*,*) '#####################################################'
526      write(*,*) 'Maximum release rate may be: ',releaserate, &
527           ' particles per second'
528      write(*,*) 'Maximum allowed release rate is: ', &
529           real(maxpart)/real(lage(nageclass)),' particles per second'
530      write(*,*) &
531           'Total number of particles released during the simulation is: ', &
532           numpartmax
533      write(*,*) 'Maximum allowed number of particles is: ',maxpart
534    endif
535  endif
536
537  return
538
539994   write(*,*) '#####################################################'
540  write(*,*) '#### FLEXPART MODEL SUBROUTINE READRELEASES:     ####'
541  write(*,*) '####                                             ####'
542  write(*,*) '#### ERROR - MAXIMUM NUMBER OF EMITTED SPECIES IS####'
543  write(*,*) '#### TOO LARGE. PLEASE REDUCE NUMBER OF SPECIES. ####'
544  write(*,*) '#####################################################'
545  stop
546
547998   write(*,*) '#####################################################'
548  write(*,*) '#### FLEXPART MODEL SUBROUTINE READRELEASES:     ####'
549  write(*,*) '####                                             ####'
550  write(*,*) '#### FATAL ERROR - FILE "RELEASES" IS            ####'
551  write(*,*) '#### CORRUPT. PLEASE CHECK YOUR INPUTS FOR       ####'
552  write(*,*) '#### MISTAKES OR GET A NEW "RELEASES"-           ####'
553  write(*,*) '#### FILE ...                                    ####'
554  write(*,*) '#####################################################'
555  stop
556
557
558999   write(*,*) '#####################################################'
559  write(*,*) '   FLEXPART MODEL SUBROUTINE READRELEASES: '
560  write(*,*)
561  write(*,*) 'FATAL ERROR - FILE CONTAINING PARTICLE RELEASE POINTS'
562  write(*,*) 'POINTS IS NOT AVAILABLE OR YOU ARE NOT'
563  write(*,*) 'PERMITTED FOR ANY ACCESS'
564  write(*,*) '#####################################################'
565  stop
566
567end subroutine readreleases
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG