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