source: flexpart.git/src/readreleases.f90 @ 54cbd6c

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 54cbd6c was 05cf28d, checked in by Espen Sollum ATMOS <eso@…>, 8 years ago

Updated checks and warning messages for wet deposition parameters given in SPECIES files

  • Property mode set to 100644
File size: 24.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_gas, wetb_gas  parameters for below-cloud scavenging (gas)            *
60  ! crain_aero, csnow_aero  parameters for below-cloud scavenging (aerosol)    *
61  ! ccn_aero, in_aero   parameters for in-cloud scavenging (aerosol)           *
62  ! zpoint1,zpoint2     height range, over which release takes place           *
63  ! num_min_discrete    if less, release cannot be randomized and happens at   *
64  !                     time mid-point of release interval                     *
65  ! lroot               true if serial version, or if MPI and root process     *
66  !                                                                            *
67  !*****************************************************************************
68
69  use point_mod
70  use xmass_mod
71  use par_mod
72  use com_mod
73
74  implicit none
75
76  integer :: numpartmax,i,j,id1,it1,id2,it2,idum,stat
77  integer,parameter :: num_min_discrete=100
78  real :: vsh(ni),fracth(ni),schmih(ni),releaserate,xdum,cun
79  real(kind=dp) :: jul1,jul2,julm,juldate
80  character(len=50) :: line
81  logical :: old
82
83  ! help variables for namelist reading
84  integer :: numpoints, parts, readerror
85  integer*2 :: zkind
86  integer :: idate1, itime1, idate2, itime2
87  real :: lon1,lon2,lat1,lat2,z1,z2
88  character*40 :: comment
89  integer,parameter :: unitreleasesout=2
90  real,allocatable, dimension (:) :: mass
91  integer,allocatable, dimension (:) :: specnum_rel,specnum_rel2
92
93  ! declare namelists
94  namelist /releases_ctrl/ &
95    nspec, &
96    specnum_rel
97
98  namelist /release/ &
99    idate1, itime1, &
100    idate2, itime2, &
101    lon1, lon2, &
102    lat1, lat2, &
103    z1, z2, &
104    zkind, &
105    mass, &
106    parts, &
107    comment
108
109  numpoint=0
110
111  ! allocate with maxspec for first input loop
112  allocate(mass(maxspec),stat=stat)
113  if (stat.ne.0) write(*,*)'ERROR: could not allocate mass'
114  allocate(specnum_rel(maxspec),stat=stat)
115  if (stat.ne.0) write(*,*)'ERROR: could not allocate specnum_rel'
116
117  ! presetting namelist releases_ctrl
118  nspec = -1  ! use negative value to determine failed namelist input
119  specnum_rel = 0
120
121  !sec, read release to find how many releasepoints should be allocated
122  open(unitreleases,file=path(1)(1:length(1))//'RELEASES',status='old',form='formatted',err=999)
123
124  ! check if namelist input provided
125  read(unitreleases,releases_ctrl,iostat=readerror)
126
127  ! prepare namelist output if requested
128  if (nmlout.and.lroot) then
129    open(unitreleasesout,file=path(2)(1:length(2))//'RELEASES.namelist',access='append',status='replace',err=1000)
130  endif
131
132  if ((readerror.ne.0).or.(nspec.lt.0)) then
133
134    ! no namelist format, close file and allow reopening in old format
135    close(unitreleases)
136    open(unitreleases,file=path(1)(1:length(1))//'RELEASES',status='old',err=999)
137
138    readerror=1 ! indicates old format
139
140    ! Check the format of the RELEASES file (either in free format,
141    ! or using a formatted mask)
142    ! Use of formatted mask is assumed if line 10 contains the word 'DIRECTION'
143    !**************************************************************************
144
145    call skplin(12,unitreleases)
146    read (unitreleases,901) line
147901 format (a)
148    if (index(line,'Total') .eq. 0) then
149      old = .false.
150    else
151      old = .true.
152    endif
153    rewind(unitreleases)
154
155
156    ! Skip first 11 lines (file header)
157    !**********************************
158
159    call skplin(11,unitreleases)
160
161    read(unitreleases,*,err=998) nspec
162    if (old) call skplin(2,unitreleases)
163    do i=1,nspec
164      read(unitreleases,*,err=998) specnum_rel(i)
165      if (old) call skplin(2,unitreleases)
166    end do
167
168    numpoint=0
169100 numpoint=numpoint+1
170    read(unitreleases,*,end=25)
171    read(unitreleases,*,err=998,end=25) idum,idum
172    if (old) call skplin(2,unitreleases)
173    read(unitreleases,*,err=998) idum,idum
174    if (old) call skplin(2,unitreleases)
175    read(unitreleases,*,err=998) xdum
176    if (old) call skplin(2,unitreleases)
177    read(unitreleases,*,err=998) xdum
178    if (old) call skplin(2,unitreleases)
179    read(unitreleases,*,err=998) xdum
180    if (old) call skplin(2,unitreleases)
181    read(unitreleases,*,err=998) xdum
182    if (old) call skplin(2,unitreleases)
183    read(unitreleases,*,err=998) idum
184    if (old) call skplin(2,unitreleases)
185    read(unitreleases,*,err=998) xdum
186    if (old) call skplin(2,unitreleases)
187    read(unitreleases,*,err=998) xdum
188    if (old) call skplin(2,unitreleases)
189    read(unitreleases,*,err=998) idum
190    if (old) call skplin(2,unitreleases)
191    do i=1,nspec
192      read(unitreleases,*,err=998) xdum
193      if (old) call skplin(2,unitreleases)
194    end do
195    !save compoint only for the first 1000 release points
196    read(unitreleases,'(a40)',err=998) compoint(1)(1:40)
197    if (old) call skplin(1,unitreleases)
198
199    goto 100
200
20125  numpoint=numpoint-1
202
203  else
204
205    readerror=0
206    do while (readerror.eq.0)
207      idate1=-1
208      read(unitreleases,release,iostat=readerror)
209      if ((idate1.lt.0).or.(readerror.ne.0)) then
210        readerror=1
211      else
212        numpoint=numpoint+1
213      endif
214    end do
215    readerror=0
216  endif ! if namelist input
217
218  rewind(unitreleases)
219
220  ! allocate arrays of matching size for number of species (namelist output)
221  deallocate(mass)
222  allocate(mass(nspec),stat=stat)
223  if (stat.ne.0) write(*,*)'ERROR: could not allocate mass'
224  allocate(specnum_rel2(nspec),stat=stat)
225  if (stat.ne.0) write(*,*)'ERROR: could not allocate specnum_rel2'
226  specnum_rel2=specnum_rel(1:nspec)
227  deallocate(specnum_rel)
228! eso: BUG, crashes here for nspec=12 and maxspec=6,
229! TODO: catch error and exit
230  allocate(specnum_rel(nspec),stat=stat)
231  if (stat.ne.0) write(*,*)'ERROR: could not allocate specnum_rel'
232  specnum_rel=specnum_rel2
233  deallocate(specnum_rel2)
234
235  !allocate memory for numpoint releaspoints
236  allocate(ireleasestart(numpoint),stat=stat)
237  if (stat.ne.0) write(*,*)'ERROR: could not allocate ireleasestart'
238  allocate(ireleaseend(numpoint),stat=stat)
239  if (stat.ne.0) write(*,*)'ERROR: could not allocate ireleaseend'
240  allocate(xpoint1(numpoint),stat=stat)
241  if (stat.ne.0) write(*,*)'ERROR: could not allocate xpoint1'
242  allocate(xpoint2(numpoint),stat=stat)
243  if (stat.ne.0) write(*,*)'ERROR: could not allocate xpoint2'
244  allocate(ypoint1(numpoint),stat=stat)
245  if (stat.ne.0) write(*,*)'ERROR: could not allocate ypoint1'
246  allocate(ypoint2(numpoint),stat=stat)
247  if (stat.ne.0) write(*,*)'ERROR: could not allocate ypoint2'
248  allocate(zpoint1(numpoint),stat=stat)
249  if (stat.ne.0) write(*,*)'ERROR: could not allocate zpoint1'
250  allocate(zpoint2(numpoint),stat=stat)
251  if (stat.ne.0) write(*,*)'ERROR: could not allocate zpoint2'
252  allocate(kindz(numpoint),stat=stat)
253  if (stat.ne.0) write(*,*)'ERROR: could not allocate kindz'
254  allocate(xmass(numpoint,maxspec),stat=stat)
255  if (stat.ne.0) write(*,*)'ERROR: could not allocate xmass'
256  allocate(rho_rel(numpoint),stat=stat)
257  if (stat.ne.0) write(*,*)'ERROR: could not allocate rho_rel'
258  allocate(npart(numpoint),stat=stat)
259  if (stat.ne.0) write(*,*)'ERROR: could not allocate npart'
260  allocate(xmasssave(numpoint),stat=stat)
261  if (stat.ne.0) write(*,*)'ERROR: could not allocate xmasssave'
262
263  if (lroot) write (*,*) 'Releasepoints : ', numpoint
264
265  do i=1,numpoint
266    xmasssave(i)=0.
267  end do
268
269  !now save the information
270  DEP=.false.
271  DRYDEP=.false.
272  WETDEP=.false.
273  OHREA=.false.
274  do i=1,maxspec
275    DRYDEPSPEC(i)=.false.
276  end do
277
278  if (readerror.ne.0) then
279    ! Skip first 11 lines (file header)
280    !**********************************
281
282    call skplin(11,unitreleases)
283
284    ! Assign species-specific parameters needed for physical processes
285    !*************************************************************************
286
287    read(unitreleases,*,err=998) nspec
288    if (nspec.gt.maxspec) goto 994
289    if (old) call skplin(2,unitreleases)
290  endif
291
292  ! namelist output
293  if (nmlout.and.lroot) then
294    write(unitreleasesout,nml=releases_ctrl)
295  endif
296
297  do i=1,nspec
298    if (readerror.ne.0) then
299      read(unitreleases,*,err=998) specnum_rel(i)
300      if (old) call skplin(2,unitreleases)
301      call readspecies(specnum_rel(i),i)
302    else
303      call readspecies(specnum_rel(i),i)
304    endif
305
306    ! For backward runs, only 1 species is allowed
307    !*********************************************
308
309    !if ((ldirect.lt.0).and.(nspec.gt.1)) then
310    !write(*,*) '#####################################################'
311    !write(*,*) '#### FLEXPART MODEL SUBROUTINE READRELEASES:     ####'
312    !write(*,*) '#### FOR BACKWARD RUNS, ONLY 1 SPECIES IS ALLOWED####'
313    !write(*,*) '#####################################################'
314    !  stop
315    !endif
316
317    ! Molecular weight
318    !*****************
319
320    if (((iout.eq.2).or.(iout.eq.3)).and.(weightmolar(i).lt.0.)) then
321      write(*,*) 'For mixing ratio output, valid molar weight'
322      write(*,*) 'must be specified for all simulated species.'
323      write(*,*) 'Check table SPECIES or choose concentration'
324      write(*,*) 'output instead if molar weight is not known.'
325      stop
326    endif
327
328    ! Radioactive decay
329    !******************
330
331    decay(i)=0.693147/decay(i) !conversion half life to decay constant
332
333
334    ! Dry deposition of gases
335    !************************
336
337    if (reldiff(i).gt.0.) rm(i)=1./(henry(i)/3000.+100.*f0(i))    ! mesophyll resistance
338
339    ! Dry deposition of particles
340    !****************************
341
342    vsetaver(i)=0.
343    cunningham(i)=0.
344    dquer(i)=dquer(i)*1000000.         ! Conversion m to um
345    if (density(i).gt.0.) then         ! Additional parameters
346      call part0(dquer(i),dsigma(i),density(i),fracth,schmih,cun,vsh)
347      do j=1,ni
348        fract(i,j)=fracth(j)
349        schmi(i,j)=schmih(j)
350        vset(i,j)=vsh(j)
351        cunningham(i)=cunningham(i)+cun*fract(i,j)
352        vsetaver(i)=vsetaver(i)-vset(i,j)*fract(i,j)
353      end do
354      if (lroot) write(*,*) 'Average settling velocity: ',i,vsetaver(i)
355    endif
356
357    ! Dry deposition for constant deposition velocity
358    !************************************************
359
360    dryvel(i)=dryvel(i)*0.01         ! conversion to m/s
361
362  ! Check if wet deposition or OH reaction shall be calculated
363  !***********************************************************
364
365  ! ESO 04.2016 check for below-cloud scavenging (gas or aerosol)
366    if ((dquer(i).le.0..and.(weta_gas(i).gt.0. .or. wetb_gas(i).gt.0.)) .or. &
367         &(dquer(i).gt.0. .and. (crain_aero(i) .gt. 0. .or. csnow_aero(i).gt.0.)))  then
368      WETDEP=.true.
369      if (lroot) then
370        write (*,*) '  Below-cloud scavenging: ON'
371!  write (*,*) 'Below-cloud scavenging coefficients: ',weta(i),i
372      end if
373    else
374      if (lroot) write (*,*) '  Below-cloud scavenging: OFF'
375    endif
376   
377! NIK 31.01.2013 + 10.12.2013 + 15.02.2015
378    if (dquer(i).gt.0..and.(ccn_aero(i).gt.0. .or. in_aero(i).gt.0.))  then
379      WETDEP=.true.
380      if (lroot) then
381        write (*,*) '  In-cloud scavenging: ON'
382!        write (*,*) 'In-cloud scavenging coefficients: ',&
383!           &ccn_aero(i),in_aero(i),i !,wetc_in(i), wetd_in(i),i
384      end if
385    else
386      if (lroot) write (*,*) '  In-cloud scavenging: OFF'
387    endif
388
389    if (ohcconst(i).gt.0.) then
390      OHREA=.true.
391      if (lroot) write (*,*) '  OHreaction switched on: ',ohcconst(i),i
392    endif
393
394    if ((reldiff(i).gt.0.).or.(density(i).gt.0.).or.(dryvel(i).gt.0.)) then
395      DRYDEP=.true.
396      DRYDEPSPEC(i)=.true.
397    endif
398
399  end do
400
401  if (WETDEP.or.DRYDEP) DEP=.true.
402
403  ! Read specifications for each release point
404  !*******************************************
405  numpoints=numpoint
406  numpoint=0
407  numpartmax=0
408  releaserate=0.
409101   numpoint=numpoint+1
410
411  if (readerror.lt.1) then ! reading namelist format
412
413    if (numpoint.gt.numpoints) goto 250
414    zkind = 1
415    mass = 0
416    parts = 0
417    comment = ' '
418    read(unitreleases,release,iostat=readerror)
419    id1=idate1
420    it1=itime1
421    id2=idate2
422    it2=itime2
423    xpoint1(numpoint)=lon1
424    xpoint2(numpoint)=lon2
425    ypoint1(numpoint)=lat1
426    ypoint2(numpoint)=lat2
427    zpoint1(numpoint)=z1
428    zpoint2(numpoint)=z2
429    kindz(numpoint)=zkind
430    do i=1,nspec
431      xmass(numpoint,i)=mass(i)
432    end do
433    npart(numpoint)=parts
434    compoint(min(1001,numpoint))=comment
435
436    ! namelist output
437    if (nmlout.and.lroot) then
438      write(unitreleasesout,nml=release)
439    endif
440
441  else
442
443    read(unitreleases,*,end=250)
444    read(unitreleases,*,err=998,end=250) id1,it1
445    if (old) call skplin(2,unitreleases)
446    read(unitreleases,*,err=998) id2,it2
447    if (old) call skplin(2,unitreleases)
448    read(unitreleases,*,err=998) xpoint1(numpoint)
449    if (old) call skplin(2,unitreleases)
450    read(unitreleases,*,err=998) ypoint1(numpoint)
451    if (old) call skplin(2,unitreleases)
452    read(unitreleases,*,err=998) xpoint2(numpoint)
453    if (old) call skplin(2,unitreleases)
454    read(unitreleases,*,err=998) ypoint2(numpoint)
455    if (old) call skplin(2,unitreleases)
456    read(unitreleases,*,err=998) kindz(numpoint)
457    if (old) call skplin(2,unitreleases)
458    read(unitreleases,*,err=998) zpoint1(numpoint)
459    if (old) call skplin(2,unitreleases)
460    read(unitreleases,*,err=998) zpoint2(numpoint)
461    if (old) call skplin(2,unitreleases)
462    read(unitreleases,*,err=998) npart(numpoint)
463    if (old) call skplin(2,unitreleases)
464    do i=1,nspec
465      read(unitreleases,*,err=998) xmass(numpoint,i)
466      if (old) call skplin(2,unitreleases)
467      mass(i)=xmass(numpoint,i)
468    end do
469    !save compoint only for the first 1000 release points
470    if (numpoint.le.1000) then
471      read(unitreleases,'(a40)',err=998) compoint(numpoint)(1:40)
472      comment=compoint(numpoint)(1:40)
473    else
474      read(unitreleases,'(a40)',err=998) compoint(1001)(1:40)
475      comment=compoint(1001)(1:40)
476    endif
477    if (old) call skplin(1,unitreleases)
478
479    ! namelist output
480    if (nmlout.and.lroot) then
481      idate1=id1
482      itime1=it1
483      idate2=id2
484      itime2=it2
485      lon1=xpoint1(numpoint)
486      lon2=xpoint2(numpoint)
487      lat1=ypoint1(numpoint)
488      lat2=ypoint2(numpoint)
489      z1=zpoint1(numpoint)
490      z2=zpoint2(numpoint)
491      zkind=kindz(numpoint)
492      parts=npart(numpoint)
493      write(unitreleasesout,nml=release)
494    endif
495
496    if (numpoint.le.1000) then
497      if((xpoint1(numpoint).eq.0.).and.(ypoint1(numpoint).eq.0.).and. &
498           (xpoint2(numpoint).eq.0.).and.(ypoint2(numpoint).eq.0.).and. &
499           (compoint(numpoint)(1:8).eq.'        ')) goto 250
500    else
501      if((xpoint1(numpoint).eq.0.).and.(ypoint1(numpoint).eq.0.).and. &
502           (xpoint2(numpoint).eq.0.).and.(ypoint2(numpoint).eq.0.)) goto 250
503    endif
504
505  endif ! if namelist format
506
507  ! If a release point contains no particles, stop and issue error message
508  !***********************************************************************
509
510  if (npart(numpoint).eq.0) then
511    write(*,*) 'FLEXPART MODEL ERROR'
512    write(*,*) 'RELEASES file is corrupt.'
513    write(*,*) 'At least for one release point, there are zero'
514    write(*,*) 'particles released. Make changes to RELEASES.'
515    stop
516  endif
517
518  ! Check whether x coordinates of release point are within model domain
519  !*********************************************************************
520
521   if (xpoint1(numpoint).lt.xlon0) &
522        xpoint1(numpoint)=xpoint1(numpoint)+360.
523   if (xpoint1(numpoint).gt.xlon0+(nxmin1)*dx) &
524        xpoint1(numpoint)=xpoint1(numpoint)-360.
525   if (xpoint2(numpoint).lt.xlon0) &
526        xpoint2(numpoint)=xpoint2(numpoint)+360.
527   if (xpoint2(numpoint).gt.xlon0+(nxmin1)*dx) &
528        xpoint2(numpoint)=xpoint2(numpoint)-360.
529
530  ! Determine relative beginning and ending times of particle release
531  !******************************************************************
532
533  jul1=juldate(id1,it1)
534  jul2=juldate(id2,it2)
535  julm=(jul1+jul2)/2.
536  if (jul1.gt.jul2) then
537    write(*,*) 'FLEXPART MODEL ERROR'
538    write(*,*) 'Release stops before it begins.'
539    write(*,*) 'Make changes to file RELEASES.'
540    stop
541  endif
542  if (mdomainfill.eq.0) then   ! no domain filling
543    if (ldirect.eq.1) then
544      if ((jul1.lt.bdate).or.(jul2.gt.edate)) then
545        write(*,*) 'FLEXPART MODEL ERROR'
546        write(*,*) 'Release starts before simulation begins or ends'
547        write(*,*) 'after simulation stops.'
548        write(*,*) 'Make files COMMAND and RELEASES consistent.'
549        stop
550      endif
551      if (npart(numpoint).gt.num_min_discrete) then
552        ireleasestart(numpoint)=int((jul1-bdate)*86400.)
553        ireleaseend(numpoint)=int((jul2-bdate)*86400.)
554      else
555        ireleasestart(numpoint)=int((julm-bdate)*86400.)
556        ireleaseend(numpoint)=int((julm-bdate)*86400.)
557      endif
558    else if (ldirect.eq.-1) then
559      if ((jul1.lt.edate).or.(jul2.gt.bdate)) then
560        write(*,*) 'FLEXPART MODEL ERROR'
561        write(*,*) 'Release starts before simulation begins or ends'
562        write(*,*) 'after simulation stops.'
563        write(*,*) 'Make files COMMAND and RELEASES consistent.'
564        stop
565      endif
566      if (npart(numpoint).gt.num_min_discrete) then
567        ireleasestart(numpoint)=int((jul1-bdate)*86400.)
568        ireleaseend(numpoint)=int((jul2-bdate)*86400.)
569      else
570        ireleasestart(numpoint)=int((julm-bdate)*86400.)
571        ireleaseend(numpoint)=int((julm-bdate)*86400.)
572      endif
573    endif
574  endif
575
576  ! Determine the release rate (particles per second) and total number
577  ! of particles released during the simulation
578  !*******************************************************************
579
580  if (ireleasestart(numpoint).ne.ireleaseend(numpoint)) then
581    releaserate=releaserate+real(npart(numpoint))/ &
582         real(ireleaseend(numpoint)-ireleasestart(numpoint))
583  else
584    releaserate=99999999
585  endif
586  numpartmax=numpartmax+npart(numpoint)
587  goto 101
588
589250   close(unitreleases)
590
591  if (nmlout.and.lroot) then
592    close(unitreleasesout)
593  endif
594
595  if (lroot) write (*,*) 'Particles allocated (maxpart)  : ',maxpart
596  if (lroot) write (*,*) 'Particles released (numpartmax): ',numpartmax
597  numpoint=numpoint-1
598
599  if (ioutputforeachrelease.eq.1) then
600    maxpointspec_act=numpoint
601  else
602    maxpointspec_act=1
603  endif
604
605  ! Check, whether the total number of particles may exceed totally allowed
606  ! number of particles at some time during the simulation
607  !************************************************************************
608
609  if (releaserate.gt. &
610       0.99*real(maxpart)/real(lage(nageclass))) then
611    if (numpartmax.gt.maxpart.and.lroot) then
612  write(*,*) '#####################################################'
613  write(*,*) '#### FLEXPART MODEL SUBROUTINE READRELEASES:     ####'
614  write(*,*) '####                                             ####'
615  write(*,*) '####WARNING - TOTAL NUMBER OF PARTICLES SPECIFIED####'
616  write(*,*) '#### IN FILE "RELEASES" MAY AT SOME POINT DURING ####'
617  write(*,*) '#### THE SIMULATION EXCEED THE MAXIMUM ALLOWED   ####'
618  write(*,*) '#### NUMBER (MAXPART).IF RELEASES DO NOT OVERLAP,####'
619  write(*,*) '#### FLEXPART CAN POSSIBLY COMPLETE SUCCESSFULLY.####'
620  write(*,*) '#### HOWEVER, FLEXPART MAY HAVE TO STOP          ####'
621  write(*,*) '#### AT SOME TIME DURING THE SIMULATION. PLEASE  ####'
622  write(*,*) '#### MAKE SURE THAT YOUR SETTINGS ARE CORRECT.   ####'
623  write(*,*) '#####################################################'
624      write(*,*) 'Maximum release rate may be: ',releaserate, &
625           ' particles per second'
626      write(*,*) 'Maximum allowed release rate is: ', &
627           real(maxpart)/real(lage(nageclass)),' particles per second'
628      write(*,*) &
629           'Total number of particles released during the simulation is: ', &
630           numpartmax
631      write(*,*) 'Maximum allowed number of particles is: ',maxpart
632    endif
633  endif
634
635  return
636
637994   write(*,*) '#####################################################'
638  write(*,*) '#### FLEXPART MODEL SUBROUTINE READRELEASES:     ####'
639  write(*,*) '####                                             ####'
640  write(*,*) '#### ERROR - MAXIMUM NUMBER OF EMITTED SPECIES IS####'
641  write(*,*) '#### TOO LARGE. PLEASE REDUCE NUMBER OF SPECIES. ####'
642  write(*,*) '#####################################################'
643  stop
644
645998   write(*,*) '#####################################################'
646  write(*,*) '#### FLEXPART MODEL SUBROUTINE READRELEASES:     ####'
647  write(*,*) '####                                             ####'
648  write(*,*) '#### FATAL ERROR - FILE "RELEASES" IS            ####'
649  write(*,*) '#### CORRUPT. PLEASE CHECK YOUR INPUTS FOR       ####'
650  write(*,*) '#### MISTAKES OR GET A NEW "RELEASES"-           ####'
651  write(*,*) '#### FILE ...                                    ####'
652  write(*,*) '#####################################################'
653  stop
654
655
656999   write(*,*) '#####################################################'
657  write(*,*) '   FLEXPART MODEL SUBROUTINE READRELEASES: '
658  write(*,*)
659  write(*,*) 'FATAL ERROR - FILE CONTAINING PARTICLE RELEASE POINTS'
660  write(*,*) 'POINTS IS NOT AVAILABLE OR YOU ARE NOT'
661  write(*,*) 'PERMITTED FOR ANY ACCESS'
662  write(*,*) '#####################################################'
663  stop
664
6651000 write(*,*) ' #### FLEXPART MODEL ERROR! FILE "RELEASES"    #### '
666  write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY       #### '
667  write(*,'(a)') path(2)(1:length(2))
668  stop
669
670end subroutine readreleases
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG