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 | |
---|
22 | subroutine 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 | !***************************************************************************** |
---|
34 | ! * |
---|
35 | ! Variables: * |
---|
36 | ! decaytime(maxtable) half time for radiological decay * |
---|
37 | ! specname(maxtable) names of chemical species, radionuclides * |
---|
38 | ! wetscava, wetscavb Parameters for determining scavenging coefficient * |
---|
39 | ! ohreact OH reaction rate * |
---|
40 | ! id_spec SPECIES number as referenced in RELEASE file * |
---|
41 | ! id_pos position where SPECIES data shall be stored * |
---|
42 | ! * |
---|
43 | ! Constants: * |
---|
44 | ! * |
---|
45 | !***************************************************************************** |
---|
46 | |
---|
47 | use par_mod |
---|
48 | use com_mod |
---|
49 | |
---|
50 | implicit none |
---|
51 | |
---|
52 | integer :: i, pos_spec,j |
---|
53 | integer :: idow,ihour,id_spec |
---|
54 | character(len=3) :: aspecnumb |
---|
55 | logical :: spec_found |
---|
56 | |
---|
57 | ! Open the SPECIES file and read species names and properties |
---|
58 | !************************************************************ |
---|
59 | specnum(pos_spec)=id_spec |
---|
60 | write(aspecnumb,'(i3.3)') specnum(pos_spec) |
---|
61 | open(unitspecies,file= & |
---|
62 | path(1)(1:length(1))//'SPECIES/SPECIES_'//aspecnumb,status='old', & |
---|
63 | err=998) |
---|
64 | !write(*,*) 'reading SPECIES',specnum(pos_spec) |
---|
65 | |
---|
66 | ASSSPEC=.FALSE. |
---|
67 | |
---|
68 | do i=1,6 |
---|
69 | read(unitspecies,*) |
---|
70 | end do |
---|
71 | |
---|
72 | read(unitspecies,'(a10)',end=22) species(pos_spec) |
---|
73 | ! write(*,*) species(pos_spec) |
---|
74 | read(unitspecies,'(f18.1)',end=22) decay(pos_spec) |
---|
75 | ! write(*,*) decay(pos_spec) |
---|
76 | read(unitspecies,'(e18.1)',end=22) weta(pos_spec) |
---|
77 | ! write(*,*) weta(pos_spec) |
---|
78 | read(unitspecies,'(f18.2)',end=22) wetb(pos_spec) |
---|
79 | ! write(*,*) wetb(pos_spec) |
---|
80 | read(unitspecies,'(f18.1)',end=22) reldiff(pos_spec) |
---|
81 | ! write(*,*) reldiff(pos_spec) |
---|
82 | read(unitspecies,'(e18.1)',end=22) henry(pos_spec) |
---|
83 | ! write(*,*) henry(pos_spec) |
---|
84 | read(unitspecies,'(f18.1)',end=22) f0(pos_spec) |
---|
85 | ! write(*,*) f0(pos_spec) |
---|
86 | read(unitspecies,'(e18.1)',end=22) density(pos_spec) |
---|
87 | ! write(*,*) density(pos_spec) |
---|
88 | read(unitspecies,'(e18.1)',end=22) dquer(pos_spec) |
---|
89 | ! write(*,*) dquer(pos_spec) |
---|
90 | read(unitspecies,'(e18.1)',end=22) dsigma(pos_spec) |
---|
91 | ! write(*,*) dsigma(pos_spec) |
---|
92 | read(unitspecies,'(f18.2)',end=22) dryvel(pos_spec) |
---|
93 | ! write(*,*) dryvel(pos_spec) |
---|
94 | read(unitspecies,'(f18.2)',end=22) weightmolar(pos_spec) |
---|
95 | ! write(*,*) weightmolar(pos_spec) |
---|
96 | read(unitspecies,'(e18.1)',end=22) ohreact(pos_spec) |
---|
97 | ! write(*,*) ohreact(pos_spec) |
---|
98 | read(unitspecies,'(i18)',end=22) spec_ass(pos_spec) |
---|
99 | ! write(*,*) spec_ass(pos_spec) |
---|
100 | read(unitspecies,'(f18.2)',end=22) kao(pos_spec) |
---|
101 | ! write(*,*) kao(pos_spec) |
---|
102 | i=pos_spec |
---|
103 | |
---|
104 | if ((weta(pos_spec).gt.0).and.(henry(pos_spec).le.0)) then |
---|
105 | if (dquer(pos_spec).le.0) goto 996 ! no particle, no henry set |
---|
106 | endif |
---|
107 | |
---|
108 | if (spec_ass(pos_spec).gt.0) then |
---|
109 | spec_found=.FALSE. |
---|
110 | do j=1,pos_spec-1 |
---|
111 | if (spec_ass(pos_spec).eq.specnum(j)) then |
---|
112 | spec_ass(pos_spec)=j |
---|
113 | spec_found=.TRUE. |
---|
114 | ASSSPEC=.TRUE. |
---|
115 | endif |
---|
116 | end do |
---|
117 | if (spec_found.eqv..FALSE.) then |
---|
118 | goto 997 |
---|
119 | endif |
---|
120 | endif |
---|
121 | |
---|
122 | if (dsigma(i).eq.1.) dsigma(i)=1.0001 ! avoid floating exception |
---|
123 | if (dsigma(i).eq.0.) dsigma(i)=1.0001 ! avoid floating exception |
---|
124 | |
---|
125 | if ((reldiff(i).gt.0.).and.(density(i).gt.0.)) then |
---|
126 | write(*,*) '#### FLEXPART MODEL ERROR! FILE "SPECIES" ####' |
---|
127 | write(*,*) '#### IS CORRUPT. SPECIES CANNOT BE BOTH ####' |
---|
128 | write(*,*) '#### PARTICLE AND GAS. ####' |
---|
129 | write(*,*) '#### SPECIES NUMBER',aspecnumb |
---|
130 | stop |
---|
131 | endif |
---|
132 | 20 continue |
---|
133 | |
---|
134 | |
---|
135 | ! Read in daily and day-of-week variation of emissions, if available |
---|
136 | !******************************************************************* |
---|
137 | |
---|
138 | do j=1,24 ! initialize everything to no variation |
---|
139 | area_hour(i,j)=1. |
---|
140 | point_hour(i,j)=1. |
---|
141 | end do |
---|
142 | do j=1,7 |
---|
143 | area_dow(i,j)=1. |
---|
144 | point_dow(i,j)=1. |
---|
145 | end do |
---|
146 | |
---|
147 | read(unitspecies,*,end=22) |
---|
148 | do j=1,24 ! 24 hours, starting with 0-1 local time |
---|
149 | read(unitspecies,*) ihour,area_hour(i,j),point_hour(i,j) |
---|
150 | end do |
---|
151 | read(unitspecies,*) |
---|
152 | do j=1,7 ! 7 days of the week, starting with Monday |
---|
153 | read(unitspecies,*) idow,area_dow(i,j),point_dow(i,j) |
---|
154 | end do |
---|
155 | |
---|
156 | 22 close(unitspecies) |
---|
157 | |
---|
158 | return |
---|
159 | |
---|
160 | 996 write(*,*) '#####################################################' |
---|
161 | write(*,*) '#### FLEXPART MODEL ERROR! #### ' |
---|
162 | write(*,*) '#### WET DEPOSITION SWITCHED ON, BUT NO HENRYS #### ' |
---|
163 | write(*,*) '#### CONSTANT IS SET ####' |
---|
164 | write(*,*) '#### PLEASE MODIFY SPECIES DESCR. FILE! #### ' |
---|
165 | write(*,*) '#####################################################' |
---|
166 | stop |
---|
167 | |
---|
168 | |
---|
169 | 997 write(*,*) '#####################################################' |
---|
170 | write(*,*) '#### FLEXPART MODEL ERROR! #### ' |
---|
171 | write(*,*) '#### THE ASSSOCIATED SPECIES HAS TO BE DEFINED #### ' |
---|
172 | write(*,*) '#### BEFORE THE ONE WHICH POINTS AT IT #### ' |
---|
173 | write(*,*) '#### PLEASE CHANGE ORDER IN RELEASES OR ADD #### ' |
---|
174 | write(*,*) '#### THE ASSOCIATED SPECIES IN RELEASES #### ' |
---|
175 | write(*,*) '#####################################################' |
---|
176 | stop |
---|
177 | |
---|
178 | |
---|
179 | 998 write(*,*) '#####################################################' |
---|
180 | write(*,*) '#### FLEXPART MODEL ERROR! #### ' |
---|
181 | write(*,*) '#### THE SPECIES FILE FOR SPECIES ', id_spec |
---|
182 | write(*,*) '#### CANNOT BE FOUND: CREATE FILE' |
---|
183 | write(*,*) '#### ',path(1)(1:length(1)),'SPECIES/SPECIES_',aspecnumb |
---|
184 | write(*,*) '#####################################################' |
---|
185 | stop |
---|
186 | |
---|
187 | |
---|
188 | end subroutine readspecies |
---|