source: branches/petra/src/readreleases.f90 @ 36

Last change on this file since 36 was 36, checked in by pesei, 9 years ago

Implement switch for incremental deposition, see ticket:113 and many small changes, see changelog.txt

  • Property svn:executable set to *
File size: 24.2 KB
Line 
1!**********************************************************************
2! Copyright 1998-2015                                                 *
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: Added optional namelist input
36  !     PS, 2/2015: access= -> position=
37  !
38  !                                                                            *
39  !*****************************************************************************
40  !                                                                            *
41  ! Variables:                                                                 *
42  ! decay               decay constant of species                              *
43  ! dquer [um]          mean particle diameters                                *
44  ! dsigma              e.g. dsigma=10 or dsigma=0.1 means that 68% of the mass*
45  !                     are between 0.1*dquer and 10*dquer                     *
46  ! ireleasestart, ireleaseend [s] starting time and ending time of each       *
47  !                     release                                                *
48  ! kindz               1: zpoint is in m agl, 2: zpoint is in m asl, 3: zpoint*
49  !                     is in hPa                                              *
50  ! npart               number of particles to be released                     *
51  ! nspec               number of species to be released                       *
52  ! density [kg/m3]     density of the particles                               *
53  ! rm [s/m]            Mesophyll resistance                                   *
54  ! species             name of species                                        *
55  ! xmass               total mass of each species                             *
56  ! xpoint1,ypoint1     geograf. coordinates of lower left corner of release   *
57  !                     area                                                   *
58  ! xpoint2,ypoint2     geograf. coordinates of upper right corner of release  *
59  !                     area                                                   *
60  ! weta, wetb          parameters to determine the below-cloud scavenging     *
61  ! weta_in, wetb_in    parameters to determine the in-cloud scavenging        *
62  ! wetc_in, wetd_in    parameters to determine the in-cloud scavenging        *
63  ! zpoint1,zpoint2     height range, over which release takes place           *
64  ! num_min_discrete    if less, release cannot be randomized and happens at   *
65  !                     time mid-point of release interval                     *
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.eqv..true.) then
129    open(unitreleasesout,file=path(2)(1:length(2))//'RELEASES.namelist',position='append',status='new',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  allocate(specnum_rel(nspec),stat=stat)
229  if (stat.ne.0) write(*,*)'ERROR: could not allocate specnum_rel'
230  specnum_rel=specnum_rel2
231  deallocate(specnum_rel2)
232
233  !allocate memory for numpoint releaspoints
234  allocate(ireleasestart(numpoint),stat=stat)
235  if (stat.ne.0) write(*,*)'ERROR: could not allocate ireleasestart'
236  allocate(ireleaseend(numpoint),stat=stat)
237  if (stat.ne.0) write(*,*)'ERROR: could not allocate ireleaseend'
238  allocate(xpoint1(numpoint),stat=stat)
239  if (stat.ne.0) write(*,*)'ERROR: could not allocate xpoint1'
240  allocate(xpoint2(numpoint),stat=stat)
241  if (stat.ne.0) write(*,*)'ERROR: could not allocate xpoint2'
242  allocate(ypoint1(numpoint),stat=stat)
243  if (stat.ne.0) write(*,*)'ERROR: could not allocate ypoint1'
244  allocate(ypoint2(numpoint),stat=stat)
245  if (stat.ne.0) write(*,*)'ERROR: could not allocate ypoint2'
246  allocate(zpoint1(numpoint),stat=stat)
247  if (stat.ne.0) write(*,*)'ERROR: could not allocate zpoint1'
248  allocate(zpoint2(numpoint),stat=stat)
249  if (stat.ne.0) write(*,*)'ERROR: could not allocate zpoint2'
250  allocate(kindz(numpoint),stat=stat)
251  if (stat.ne.0) write(*,*)'ERROR: could not allocate kindz'
252  allocate(xmass(numpoint,maxspec),stat=stat)
253  if (stat.ne.0) write(*,*)'ERROR: could not allocate xmass'
254  allocate(rho_rel(numpoint),stat=stat)
255  if (stat.ne.0) write(*,*)'ERROR: could not allocate rho_rel'
256  allocate(npart(numpoint),stat=stat)
257  if (stat.ne.0) write(*,*)'ERROR: could not allocate npart'
258  allocate(xmasssave(numpoint),stat=stat)
259  if (stat.ne.0) write(*,*)'ERROR: could not allocate xmasssave'
260
261  write (*,*) 'Releasepoints : ', numpoint
262
263  do i=1,numpoint
264    xmasssave(i)=0.
265  end do
266
267  !now save the information
268  DEP=.false.
269  DRYDEP=.false.
270  WETDEP=.false.
271  OHREA=.false.
272  do i=1,maxspec
273    DRYDEPSPEC(i)=.false.
274  end do
275
276  if (readerror.ne.0) then
277    ! Skip first 11 lines (file header)
278    !**********************************
279
280    call skplin(11,unitreleases)
281
282    ! Assign species-specific parameters needed for physical processes
283    !*************************************************************************
284
285    read(unitreleases,*,err=998) nspec
286    if (nspec.gt.maxspec) goto 994
287    if (old) call skplin(2,unitreleases)
288  endif
289
290  ! namelist output
291  if (nmlout.eqv..true.) then
292    write(unitreleasesout,nml=releases_ctrl)
293  endif
294
295  do i=1,nspec
296    if (readerror.ne.0) then
297      read(unitreleases,*,err=998) specnum_rel(i)
298      if (old) call skplin(2,unitreleases)
299      call readspecies(specnum_rel(i),i)
300    else
301      call readspecies(specnum_rel(i),i)
302    endif
303
304    ! For backward runs, only 1 species is allowed
305    !*********************************************
306
307    !if ((ldirect.lt.0).and.(nspec.gt.1)) then
308    !write(*,*) '#####################################################'
309    !write(*,*) '#### FLEXPART MODEL SUBROUTINE READRELEASES:     ####'
310    !write(*,*) '#### FOR BACKWARD RUNS, ONLY 1 SPECIES IS ALLOWED####'
311    !write(*,*) '#####################################################'
312    !  stop
313    !endif
314
315    ! Molecular weight
316    !*****************
317
318    if (((iout.eq.2).or.(iout.eq.3)).and.(weightmolar(i).lt.0.)) then
319      write(*,*) 'For mixing ratio output, valid molar weight'
320      write(*,*) 'must be specified for all simulated species.'
321      write(*,*) 'Check table SPECIES or choose concentration'
322      write(*,*) 'output instead if molar weight is not known.'
323      stop
324    endif
325
326    ! Radioactive decay
327    !******************
328
329    decay(i)=0.693147/decay(i) !conversion half life to decay constant
330
331
332    ! Dry deposition of gases
333    !************************
334
335    if (reldiff(i).gt.0.) rm(i)=1./(henry(i)/3000.+100.*f0(i))    ! mesophyll resistance
336
337    ! Dry deposition of particles
338    !****************************
339
340    vsetaver(i)=0.
341    cunningham(i)=0.
342    dquer(i)=dquer(i)*1000000.         ! Conversion m to um
343    if (density(i).gt.0.) then         ! Additional parameters
344      call part0(dquer(i),dsigma(i),density(i),fracth,schmih,cun,vsh)
345      do j=1,ni
346        fract(i,j)=fracth(j)
347        schmi(i,j)=schmih(j)
348        vset(i,j)=vsh(j)
349        cunningham(i)=cunningham(i)+cun*fract(i,j)
350        vsetaver(i)=vsetaver(i)-vset(i,j)*fract(i,j)
351      end do
352      write(*,*) 'Average settling velocity: ',i,vsetaver(i)
353    endif
354
355    ! Dry deposition for constant deposition velocity
356    !************************************************
357
358    dryvel(i)=dryvel(i)*0.01         ! conversion to m/s
359
360    ! Check if wet deposition or OH reaction shall be calculated
361    !***********************************************************
362    if (weta(i).gt.0.)  then
363      WETDEP=.true.
364      write (*,*) 'Below-cloud scavenging: ON'
365      write (*,*) 'Below-cloud scavenging coefficients: ',weta(i),i
366    else
367      write (*,*) 'Below-cloud scavenging: OFF'
368    endif
369   
370    ! NIK 31.01.2013 + 10.12.2013
371    if (weta_in(i).gt.0.)  then
372      WETDEP=.true.
373      write (*,*) 'In-cloud scavenging: ON'
374      write (*,*) 'In-cloud scavenging coefficients: ',weta_in(i),wetb_in(i), wetc_in(i), wetd_in(i),i
375    else
376      write (*,*) 'In-cloud scavenging: OFF'
377    endif
378
379    if (ohreact(i).gt.0) then
380      OHREA=.true.
381      write (*,*) 'OHreaction: ON (',ohreact(i),i,')'
382    endif
383
384    if ((reldiff(i).gt.0.).or.(density(i).gt.0.).or.(dryvel(i).gt.0.)) then
385      DRYDEP=.true.
386      DRYDEPSPEC(i)=.true.
387    endif
388
389  end do
390
391  if (WETDEP.or.DRYDEP) DEP=.true.
392
393  ! Read specifications for each release point
394  !*******************************************
395  numpoints=numpoint
396  numpoint=0
397  numpartmax=0
398  releaserate=0.
399101   numpoint=numpoint+1
400
401  if (readerror.lt.1) then ! reading namelist format
402
403    if (numpoint.gt.numpoints) goto 250
404    zkind = 1
405    mass = 0
406    parts = 0
407    comment = ' '
408    read(unitreleases,release,iostat=readerror)
409    id1=idate1
410    it1=itime1
411    id2=idate2
412    it2=itime2
413    xpoint1(numpoint)=lon1
414    xpoint2(numpoint)=lon2
415    ypoint1(numpoint)=lat1
416    ypoint2(numpoint)=lat2
417    zpoint1(numpoint)=z1
418    zpoint2(numpoint)=z2
419    kindz(numpoint)=zkind
420    do i=1,nspec
421      xmass(numpoint,i)=mass(i)
422    end do
423    npart(numpoint)=parts
424    compoint(min(1001,numpoint))=comment
425
426    ! namelist output
427    if (nmlout.eqv..true.) then
428      write(unitreleasesout,nml=release)
429    endif
430
431  else
432
433    read(unitreleases,*,end=250)
434    read(unitreleases,*,err=998,end=250) id1,it1
435    if (old) call skplin(2,unitreleases)
436    read(unitreleases,*,err=998) id2,it2
437    if (old) call skplin(2,unitreleases)
438    read(unitreleases,*,err=998) xpoint1(numpoint)
439    if (old) call skplin(2,unitreleases)
440    read(unitreleases,*,err=998) ypoint1(numpoint)
441    if (old) call skplin(2,unitreleases)
442    read(unitreleases,*,err=998) xpoint2(numpoint)
443    if (old) call skplin(2,unitreleases)
444    read(unitreleases,*,err=998) ypoint2(numpoint)
445    if (old) call skplin(2,unitreleases)
446    read(unitreleases,*,err=998) kindz(numpoint)
447    if (old) call skplin(2,unitreleases)
448    read(unitreleases,*,err=998) zpoint1(numpoint)
449    if (old) call skplin(2,unitreleases)
450    read(unitreleases,*,err=998) zpoint2(numpoint)
451    if (old) call skplin(2,unitreleases)
452    read(unitreleases,*,err=998) npart(numpoint)
453    if (old) call skplin(2,unitreleases)
454    do i=1,nspec
455      read(unitreleases,*,err=998) xmass(numpoint,i)
456      if (old) call skplin(2,unitreleases)
457      mass(i)=xmass(numpoint,i)
458    end do
459    !save compoint only for the first 1000 release points
460    if (numpoint.le.1000) then
461      read(unitreleases,'(a40)',err=998) compoint(numpoint)(1:40)
462      comment=compoint(numpoint)(1:40)
463    else
464      read(unitreleases,'(a40)',err=998) compoint(1001)(1:40)
465      comment=compoint(1001)(1:40)
466    endif
467    if (old) call skplin(1,unitreleases)
468
469    ! namelist output
470    if (nmlout.eqv..true.) then
471      idate1=id1
472      itime1=it1
473      idate2=id2
474      itime2=it2
475      lon1=xpoint1(numpoint)
476      lon2=xpoint2(numpoint)
477      lat1=ypoint1(numpoint)
478      lat2=ypoint2(numpoint)
479      z1=zpoint1(numpoint)
480      z2=zpoint2(numpoint)
481      zkind=kindz(numpoint)
482      parts=npart(numpoint)
483      write(unitreleasesout,nml=release)
484    endif
485
486    if (numpoint.le.1000) then
487      if((xpoint1(numpoint).eq.0.).and.(ypoint1(numpoint).eq.0.).and. &
488           (xpoint2(numpoint).eq.0.).and.(ypoint2(numpoint).eq.0.).and. &
489           (compoint(numpoint)(1:8).eq.'        ')) goto 250
490    else
491      if((xpoint1(numpoint).eq.0.).and.(ypoint1(numpoint).eq.0.).and. &
492           (xpoint2(numpoint).eq.0.).and.(ypoint2(numpoint).eq.0.)) goto 250
493    endif
494
495  endif ! if namelist format
496
497  ! If a release point contains no particles, stop and issue error message
498  !***********************************************************************
499
500  if (npart(numpoint).eq.0) then
501    write(*,*) 'FLEXPART MODEL ERROR'
502    write(*,*) 'RELEASES file is corrupt.'
503    write(*,*) 'At least for one release point, there are zero'
504    write(*,*) 'particles released. Make changes to RELEASES.'
505    stop
506  endif
507
508  ! Check whether x coordinates of release point are within model domain
509  !*********************************************************************
510
511   if (xpoint1(numpoint).lt.xlon0) &
512        xpoint1(numpoint)=xpoint1(numpoint)+360.
513   if (xpoint1(numpoint).gt.xlon0+(nxmin1)*dx) &
514        xpoint1(numpoint)=xpoint1(numpoint)-360.
515   if (xpoint2(numpoint).lt.xlon0) &
516        xpoint2(numpoint)=xpoint2(numpoint)+360.
517   if (xpoint2(numpoint).gt.xlon0+(nxmin1)*dx) &
518        xpoint2(numpoint)=xpoint2(numpoint)-360.
519
520  ! Determine relative beginning and ending times of particle release
521  !******************************************************************
522
523  jul1=juldate(id1,it1)
524  jul2=juldate(id2,it2)
525  julm=(jul1+jul2)/2.
526  if (jul1.gt.jul2) then
527    write(*,*) 'FLEXPART MODEL ERROR'
528    write(*,*) 'Release stops before it begins.'
529    write(*,*) 'Make changes to file RELEASES.'
530    stop
531  endif
532  if (mdomainfill.eq.0) then   ! no domain filling
533    if (ldirect.eq.1) then
534      if ((jul1.lt.bdate).or.(jul2.gt.edate)) then
535        write(*,*) 'FLEXPART MODEL ERROR'
536        write(*,*) 'Release starts before simulation begins or ends'
537        write(*,*) 'after simulation stops.'
538        write(*,*) 'Make files COMMAND and RELEASES consistent.'
539        stop
540      endif
541      if (npart(numpoint).gt.num_min_discrete) then
542        ireleasestart(numpoint)=int((jul1-bdate)*86400.)
543        ireleaseend(numpoint)=int((jul2-bdate)*86400.)
544      else
545        ireleasestart(numpoint)=int((julm-bdate)*86400.)
546        ireleaseend(numpoint)=int((julm-bdate)*86400.)
547      endif
548    else if (ldirect.eq.-1) then
549      if ((jul1.lt.edate).or.(jul2.gt.bdate)) then
550        write(*,*) 'FLEXPART MODEL ERROR'
551        write(*,*) 'Release starts before simulation begins or ends'
552        write(*,*) 'after simulation stops.'
553        write(*,*) 'Make files COMMAND and RELEASES consistent.'
554        stop
555      endif
556      if (npart(numpoint).gt.num_min_discrete) then
557        ireleasestart(numpoint)=int((jul1-bdate)*86400.)
558        ireleaseend(numpoint)=int((jul2-bdate)*86400.)
559      else
560        ireleasestart(numpoint)=int((julm-bdate)*86400.)
561        ireleaseend(numpoint)=int((julm-bdate)*86400.)
562      endif
563    endif
564  endif
565
566  ! Determine the release rate (particles per second) and total number
567  ! of particles released during the simulation
568  !*******************************************************************
569
570  if (ireleasestart(numpoint).ne.ireleaseend(numpoint)) then
571    releaserate=releaserate+real(npart(numpoint))/ &
572         real(ireleaseend(numpoint)-ireleasestart(numpoint))
573  else
574    releaserate=99999999
575  endif
576  numpartmax=numpartmax+npart(numpoint)
577  goto 101
578
579250   close(unitreleases)
580
581  if (nmlout.eqv..true.) then
582    close(unitreleasesout)
583  endif
584
585  write (*,*) 'Particles allocated (maxpart)  : ',maxpart
586  write (*,*) 'Particles released (numpartmax): ',numpartmax
587  numpoint=numpoint-1
588
589  if (ioutputforeachrelease.eq.1) then
590    maxpointspec_act=numpoint
591  else
592    maxpointspec_act=1
593  endif
594
595  ! Check, whether the total number of particles may exceed totally allowed
596  ! number of particles at some time during the simulation
597  !************************************************************************
598
599  if (releaserate.gt. &
600       0.99*real(maxpart)/real(lage(nageclass))) then
601    if (numpartmax.gt.maxpart) then
602  write(*,*) '#####################################################'
603  write(*,*) '#### FLEXPART MODEL SUBROUTINE READRELEASES:     ####'
604  write(*,*) '####                                             ####'
605  write(*,*) '####WARNING - TOTAL NUMBER OF PARTICLES SPECIFIED####'
606  write(*,*) '#### IN FILE "RELEASES" MAY AT SOME POINT DURING ####'
607  write(*,*) '#### THE SIMULATION EXCEED THE MAXIMUM ALLOWED   ####'
608  write(*,*) '#### NUMBER (MAXPART).IF RELEASES DO NOT OVERLAP,####'
609  write(*,*) '#### FLEXPART CAN POSSIBLY COMPLETE SUCCESSFULLY.####'
610  write(*,*) '#### HOWEVER, FLEXPART MAY HAVE TO STOP          ####'
611  write(*,*) '#### AT SOME TIME DURING THE SIMULATION. PLEASE  ####'
612  write(*,*) '#### MAKE SURE THAT YOUR SETTINGS ARE CORRECT.   ####'
613  write(*,*) '#####################################################'
614      write(*,*) 'Maximum release rate may be: ',releaserate, &
615           ' particles per second'
616      write(*,*) 'Maximum allowed release rate is: ', &
617           real(maxpart)/real(lage(nageclass)),' particles per second'
618      write(*,*) &
619           'Total number of particles released during the simulation is: ', &
620           numpartmax
621      write(*,*) 'Maximum allowed number of particles is: ',maxpart
622    endif
623  endif
624
625  return
626
627994   write(*,*) '#####################################################'
628  write(*,*) '#### FLEXPART MODEL SUBROUTINE READRELEASES:     ####'
629  write(*,*) '####                                             ####'
630  write(*,*) '#### ERROR - MAXIMUM NUMBER OF EMITTED SPECIES IS####'
631  write(*,*) '#### TOO LARGE. PLEASE REDUCE NUMBER OF SPECIES. ####'
632  write(*,*) '#####################################################'
633  stop
634
635998   write(*,*) '#####################################################'
636  write(*,*) '#### FLEXPART MODEL SUBROUTINE READRELEASES:     ####'
637  write(*,*) '####                                             ####'
638  write(*,*) '#### FATAL ERROR - FILE "RELEASES" IS            ####'
639  write(*,*) '#### CORRUPT. PLEASE CHECK YOUR INPUTS FOR       ####'
640  write(*,*) '#### MISTAKES OR GET A NEW "RELEASES"-           ####'
641  write(*,*) '#### FILE ...                                    ####'
642  write(*,*) '#####################################################'
643  stop
644
645
646999   write(*,*) '#####################################################'
647  write(*,*) '   FLEXPART MODEL SUBROUTINE READRELEASES: '
648  write(*,*)
649  write(*,*) 'FATAL ERROR - FILE CONTAINING PARTICLE RELEASE POINTS'
650  write(*,*) 'POINTS IS NOT AVAILABLE OR YOU ARE NOT'
651  write(*,*) 'PERMITTED FOR ANY ACCESS'
652  write(*,*) '#####################################################'
653  stop
654
6551000 write(*,*) ' #### FLEXPART MODEL ERROR! FILE "RELEASES"    #### '
656  write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY       #### '
657  write(*,'(a)') path(2)(1:length(2))
658  stop
659
660end subroutine readreleases
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG