source: flexpart.git/src/readspecies.f90 @ 43c8684

10.4.1_peseiFPv9.3.1FPv9.3.1b_testingFPv9.3.2GFS_025bugfixes+enhancementsdevfp9.3.1-20161214-nc4grib2nc4_repairrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 43c8684 was 43c8684, checked in by Ignacio Pisso <Ignacio.Pisso@…>, 9 years ago

add verbosity option to readspecies.f90

  • Property mode set to 100644
File size: 11.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 readspecies(id_spec,pos_spec)
23
24  !*****************************************************************************
25  !                                                                            *
26  !     This routine reads names and physical constants of chemical species/   *
27  !     radionuclides given in the parameter pos_spec                          *
28  !                                                                            *
29  !   Author: A. Stohl                                                         *
30  !                                                                            *
31  !   11 July 1996                                                             *
32  !
33  !   Changes:                                                                 *
34  !   N. Kristiansen, 31.01.2013: Including parameters for in-cloud scavenging *
35  !                                                                            *
36  !   HSO, 13 August 2013
37  !   added optional namelist input
38  !                                                                            *
39  !*****************************************************************************
40  !                                                                            *
41  ! Variables:                                                                 *
42  ! decaytime(maxtable)  half time for radiological decay                      *
43  ! specname(maxtable)   names of chemical species, radionuclides              *
44  ! weta, wetb           Parameters for determining below-cloud scavenging     *
45  ! weta_in              Parameter for determining in-cloud scavenging         *
46  ! wetb_in              Parameter for determining in-cloud scavenging         *
47  ! wetc_in              Parameter for determining in-cloud scavenging         *
48  ! wetd_in              Parameter for determining in-cloud scavenging         *
49  ! ohreact              OH reaction rate                                      *
50  ! id_spec              SPECIES number as referenced in RELEASE file          *
51  ! id_pos               position where SPECIES data shall be stored           *
52  !                                                                            *
53  ! Constants:                                                                 *
54  !                                                                            *
55  !*****************************************************************************
56
57  use par_mod
58  use com_mod
59
60  implicit none
61
62  integer :: i, pos_spec,j
63  integer :: idow,ihour,id_spec
64  character(len=3) :: aspecnumb
65  logical :: spec_found
66
67  character(len=16) :: pspecies
68  real :: pdecay, pweta, pwetb, preldiff, phenry, pf0, pdensity, pdquer
69  real :: pdsigma, pdryvel, pweightmolar, pohreact, pspec_ass, pkao
70  real :: pweta_in, pwetb_in, pwetc_in, pwetd_in
71  integer :: readerror
72
73  ! declare namelist
74  namelist /species_params/ &
75   pspecies, pdecay, pweta, pwetb, &
76   pweta_in, pwetb_in, pwetc_in, pwetd_in, &
77   preldiff, phenry, pf0, pdensity, pdquer, &
78   pdsigma, pdryvel, pweightmolar, pohreact, pspec_ass, pkao
79
80  pspecies=" "
81  pdecay=-999.9
82  pweta=-9.9E-09
83  pwetb=0.0
84  pweta_in=-9.9E-09
85  pwetb_in=-9.9E-09
86  pwetc_in=-9.9E-09
87  pwetd_in=-9.9E-09
88  preldiff=-9.9
89  phenry=0.0
90  pf0=0.0
91  pdensity=-9.9E09
92  pdquer=0.0
93  pdsigma=0.0
94  pdryvel=-9.99
95  pohreact=-9.9E-09
96  pspec_ass=-9
97  pkao=-99.99
98  pweightmolar=-789.0 ! read failure indicator value
99
100  ! Open the SPECIES file and read species names and properties
101  !************************************************************
102  specnum(pos_spec)=id_spec
103  write(aspecnumb,'(i3.3)') specnum(pos_spec)
104  open(unitspecies,file=path(1)(1:length(1))//'SPECIES/SPECIES_'//aspecnumb,status='old',form='formatted',err=998)
105 
106  if (verbosity.gt.0) then
107    write(*,*) 'reading SPECIES',specnum(pos_spec)
108  endif
109
110  ASSSPEC=.FALSE.
111
112  ! try namelist input
113  read(unitspecies,species_params,iostat=readerror)
114  close(unitspecies)
115   
116  if ((pweightmolar.eq.-789.0).or.(readerror.ne.0)) then ! no namelist found
117
118    readerror=1
119
120    open(unitspecies,file=path(1)(1:length(1))//'SPECIES/SPECIES_'//aspecnumb,status='old',err=998)
121
122    do i=1,6
123      read(unitspecies,*)
124    end do
125
126    read(unitspecies,'(a10)',end=22) species(pos_spec)
127  !  write(*,*) species(pos_spec)
128    read(unitspecies,'(f18.1)',end=22) decay(pos_spec)
129  !  write(*,*) decay(pos_spec)
130    read(unitspecies,'(e18.1)',end=22) weta(pos_spec)
131  !  write(*,*) weta(pos_spec)
132    read(unitspecies,'(f18.2)',end=22) wetb(pos_spec)
133  !  write(*,*) wetb(pos_spec)
134
135  !*** NIK 31.01.2013: including in-cloud scavening parameters
136   read(unitspecies,'(e18.1)',end=22) weta_in(pos_spec)
137  !  write(*,*) weta_in(pos_spec)
138   read(unitspecies,'(f18.2)',end=22) wetb_in(pos_spec)
139  !  write(*,*) wetb_in(pos_spec)
140   read(unitspecies,'(f18.2)',end=22) wetc_in(pos_spec)
141  !  write(*,*) wetc_in(pos_spec)
142   read(unitspecies,'(f18.2)',end=22) wetd_in(pos_spec)
143  !  write(*,*) wetd_in(pos_spec)
144
145    read(unitspecies,'(f18.1)',end=22) reldiff(pos_spec)
146  !  write(*,*) reldiff(pos_spec)
147    read(unitspecies,'(e18.1)',end=22) henry(pos_spec)
148  !  write(*,*) henry(pos_spec)
149    read(unitspecies,'(f18.1)',end=22) f0(pos_spec)
150  !  write(*,*) f0(pos_spec)
151    read(unitspecies,'(e18.1)',end=22) density(pos_spec)
152  !  write(*,*) density(pos_spec)
153    read(unitspecies,'(e18.1)',end=22) dquer(pos_spec)
154  !  write(*,*) dquer(pos_spec)
155    read(unitspecies,'(e18.1)',end=22) dsigma(pos_spec)
156  !  write(*,*) dsigma(pos_spec)
157    read(unitspecies,'(f18.2)',end=22) dryvel(pos_spec)
158  !  write(*,*) dryvel(pos_spec)
159    read(unitspecies,'(f18.2)',end=22) weightmolar(pos_spec)
160  !  write(*,*) weightmolar(pos_spec)
161    read(unitspecies,'(e18.1)',end=22) ohreact(pos_spec)
162  !  write(*,*) ohreact(pos_spec)
163    read(unitspecies,'(i18)',end=22) spec_ass(pos_spec)
164  !  write(*,*) spec_ass(pos_spec)
165    read(unitspecies,'(f18.2)',end=22) kao(pos_spec)
166  !       write(*,*) kao(pos_spec)
167
168    pspecies=species(pos_spec)
169    pdecay=decay(pos_spec)
170    pweta=weta(pos_spec)
171    pwetb=wetb(pos_spec)
172    pweta_in=weta_in(pos_spec)
173    pwetb_in=wetb_in(pos_spec)
174    pwetc_in=wetc_in(pos_spec)
175    pwetd_in=wetd_in(pos_spec)
176    preldiff=reldiff(pos_spec)
177    phenry=henry(pos_spec)
178    pf0=f0(pos_spec)
179    pdensity=density(pos_spec)
180    pdquer=dquer(pos_spec)
181    pdsigma=dsigma(pos_spec)
182    pdryvel=dryvel(pos_spec)
183    pweightmolar=weightmolar(pos_spec)
184    pohreact=ohreact(pos_spec)
185    pspec_ass=spec_ass(pos_spec)
186    pkao=kao(pos_spec)
187
188  else
189
190    species(pos_spec)=pspecies
191    decay(pos_spec)=pdecay
192    weta(pos_spec)=pweta
193    wetb(pos_spec)=pwetb
194    weta_in(pos_spec)=pweta_in
195    wetb_in(pos_spec)=pwetb_in
196    wetc_in(pos_spec)=pwetc_in
197    wetd_in(pos_spec)=pwetd_in
198    reldiff(pos_spec)=preldiff
199    henry(pos_spec)=phenry
200    f0(pos_spec)=pf0
201    density(pos_spec)=pdensity
202    dquer(pos_spec)=pdquer
203    dsigma(pos_spec)=pdsigma
204    dryvel(pos_spec)=pdryvel
205    weightmolar(pos_spec)=pweightmolar
206    ohreact(pos_spec)=pohreact
207    spec_ass(pos_spec)=pspec_ass
208    kao(pos_spec)=pkao
209
210  endif
211
212  i=pos_spec
213
214  if ((weta(pos_spec).gt.0).and.(henry(pos_spec).le.0)) then
215   if (dquer(pos_spec).le.0) goto 996 ! no particle, no henry set
216  endif
217
218  if (spec_ass(pos_spec).gt.0) then
219    spec_found=.FALSE.
220    do j=1,pos_spec-1
221      if (spec_ass(pos_spec).eq.specnum(j)) then
222        spec_ass(pos_spec)=j
223        spec_found=.TRUE.
224        ASSSPEC=.TRUE.
225      endif
226    end do
227    if (spec_found.eqv..FALSE.) then
228      goto 997
229    endif
230  endif
231
232  if (dsigma(i).eq.1.) dsigma(i)=1.0001   ! avoid floating exception
233  if (dsigma(i).eq.0.) dsigma(i)=1.0001   ! avoid floating exception
234
235  if ((reldiff(i).gt.0.).and.(density(i).gt.0.)) then
236    write(*,*) '#### FLEXPART MODEL ERROR! FILE "SPECIES"    ####'
237    write(*,*) '#### IS CORRUPT. SPECIES CANNOT BE BOTH      ####'
238    write(*,*) '#### PARTICLE AND GAS.                       ####'
239    write(*,*) '#### SPECIES NUMBER',aspecnumb
240    stop
241  endif
24220   continue
243
244  if (readerror.ne.0) then ! text format input
245
246    ! Read in daily and day-of-week variation of emissions, if available
247    !*******************************************************************
248    ! HSO: This is not yet implemented as namelist parameters
249
250    do j=1,24           ! initialize everything to no variation
251      area_hour(i,j)=1.
252      point_hour(i,j)=1.
253    end do
254    do j=1,7
255      area_dow(i,j)=1.
256      point_dow(i,j)=1.
257    end do
258
259    read(unitspecies,*,end=22)
260    do j=1,24     ! 24 hours, starting with 0-1 local time
261      read(unitspecies,*) ihour,area_hour(i,j),point_hour(i,j)
262    end do
263    read(unitspecies,*)
264    do j=1,7      ! 7 days of the week, starting with Monday
265      read(unitspecies,*) idow,area_dow(i,j),point_dow(i,j)
266    end do
267
268  endif
269
27022 close(unitspecies)
271
272  ! namelist output if requested
273  if (nmlout.eqv..true.) then
274    open(unitspecies,file=path(2)(1:length(2))//'SPECIES_'//aspecnumb//'.namelist',access='append',status='new',err=1000)
275    write(unitspecies,nml=species_params)
276    close(unitspecies)
277  endif
278
279  return
280
281996   write(*,*) '#####################################################'
282  write(*,*) '#### FLEXPART MODEL ERROR!                      #### '
283  write(*,*) '#### WET DEPOSITION SWITCHED ON, BUT NO HENRYS  #### '
284  write(*,*) '#### CONSTANT IS SET                            ####'
285  write(*,*) '#### PLEASE MODIFY SPECIES DESCR. FILE!        #### '
286  write(*,*) '#####################################################'
287  stop
288
289
290997   write(*,*) '#####################################################'
291  write(*,*) '#### FLEXPART MODEL ERROR!                      #### '
292  write(*,*) '#### THE ASSSOCIATED SPECIES HAS TO BE DEFINED  #### '
293  write(*,*) '#### BEFORE THE ONE WHICH POINTS AT IT          #### '
294  write(*,*) '#### PLEASE CHANGE ORDER IN RELEASES OR ADD     #### '
295  write(*,*) '#### THE ASSOCIATED SPECIES IN RELEASES         #### '
296  write(*,*) '#####################################################'
297  stop
298
299
300998   write(*,*) '#####################################################'
301  write(*,*) '#### FLEXPART MODEL ERROR!                      #### '
302  write(*,*) '#### THE SPECIES FILE FOR SPECIES ', id_spec
303  write(*,*) '#### CANNOT BE FOUND: CREATE FILE'
304  write(*,*) '#### ',path(1)(1:length(1)),'SPECIES/SPECIES_',aspecnumb
305  write(*,*) '#####################################################'
306  stop
307
3081000 write(*,*) ' #### FLEXPART MODEL ERROR! FILE "SPECIES_',aspecnumb,'.namelist'
309  write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY       #### '
310  write(*,'(a)') path(2)(1:length(2))
311  stop
312
313end subroutine readspecies
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG