1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145
|
!**********************************************************************
! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010 *
! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa, *
! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann *
! *
! This file is part of FLEXPART. *
! *
! FLEXPART is free software: you can redistribute it and/or modify *
! it under the terms of the GNU General Public License as published by*
! the Free Software Foundation, either version 3 of the License, or *
! (at your option) any later version. *
! *
! FLEXPART is distributed in the hope that it will be useful, *
! but WITHOUT ANY WARRANTY; without even the implied warranty of *
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
! GNU General Public License for more details. *
! *
! You should have received a copy of the GNU General Public License *
! along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. *
!**********************************************************************
subroutine readdepo
!*****************************************************************************
! *
! Reads dry deposition parameters needed by the procedure of Wesely (1989). *
! Wesely (1989): Parameterization of surface resistances to gaseous *
! dry deposition in regional-scale numerical models. *
! Atmos. Environ. 23, 1293-1304. *
! *
! *
! AUTHOR: Andreas Stohl, 19 May 1995 *
! *
!*****************************************************************************
! *
! Variables: *
! *
! rcl(maxspec,5,9) [s/m] Lower canopy resistance *
! rgs(maxspec,5,9) [s/m] Ground resistance *
! rlu(maxspec,5,9) [s/m] Leaf cuticular resistance *
! rm(maxspec) [s/m] Mesophyll resistance, set in readreleases *
! ri(maxspec) [s/m] Stomatal resistance *
! *
! Constants: *
! *
!*****************************************************************************
use par_mod
use com_mod
implicit none
! FOR THIS SUBROUTINE, numclass=9 IS ASSUMED
!*******************************************
real :: rluh(5,numclass),rgssh(5,numclass),rgsoh(5,numclass)
real :: rclsh(5,numclass),rcloh(5,numclass)
integer :: i,j,ic
! Read deposition constants related with landuse and seasonal category
!*********************************************************************
open(unitwesely,file=path(1)(1:length(1))//'surfdepo.t', &
status='old',err=999)
do i=1,16
read(unitwesely,*)
end do
do i=1,5
read(unitwesely,*)
read(unitwesely,'(8x,13f8.0)') (ri(i,j),j=1,numclass)
read(unitwesely,'(8x,13f8.0)') (rluh(i,j),j=1,numclass)
read(unitwesely,'(8x,13f8.0)') (rac(i,j),j=1,numclass)
read(unitwesely,'(8x,13f8.0)') (rgssh(i,j),j=1,numclass)
read(unitwesely,'(8x,13f8.0)') (rgsoh(i,j),j=1,numclass)
read(unitwesely,'(8x,13f8.0)') (rclsh(i,j),j=1,numclass)
read(unitwesely,'(8x,13f8.0)') (rcloh(i,j),j=1,numclass)
end do
! TEST
! do 31 i=1,5
! ri(i,13)=ri(i,5)
! rluh(i,13)=rluh(i,5)
! rac(i,13)=rac(i,5)
! rgssh(i,13)=rgssh(i,5)
! rgsoh(i,13)=rgsoh(i,5)
! rclsh(i,13)=rclsh(i,5)
! rcloh(i,13)=rcloh(i,5)
!31 continue
! TEST
! Sabine Eckhardt, Dec 06, set resistances of 9999 to 'infinite' (1E25)
do i=1,5
do j=1,numclass
if (ri(i,j).eq.9999.) ri(i,j)=1.E25
if (rluh(i,j).eq.9999.) rluh(i,j)=1.E25
if (rac(i,j).eq.9999.) rac(i,j)=1.E25
if (rgssh(i,j).eq.9999.) rgssh(i,j)=1.E25
if (rgsoh(i,j).eq.9999.) rgsoh(i,j)=1.E25
if (rclsh(i,j).eq.9999.) rclsh(i,j)=1.E25
if (rcloh(i,j).eq.9999.) rcloh(i,j)=1.E25
end do
end do
do i=1,5
do j=1,numclass
ri(i,j)=max(ri(i,j),0.001)
rluh(i,j)=max(rluh(i,j),0.001)
rac(i,j)=max(rac(i,j),0.001)
rgssh(i,j)=max(rgssh(i,j),0.001)
rgsoh(i,j)=max(rgsoh(i,j),0.001)
rclsh(i,j)=max(rclsh(i,j),0.001)
rcloh(i,j)=max(rcloh(i,j),0.001)
end do
end do
close(unitwesely)
! Compute additional parameters
!******************************
do ic=1,nspec
if (reldiff(ic).gt.0.) then ! gas is dry deposited
do i=1,5
do j=1,numclass
rlu(ic,i,j)=rluh(i,j)/(1.e-5*henry(ic)+f0(ic))
rgs(ic,i,j)=1./(henry(ic)/(10.e5*rgssh(i,j))+f0(ic)/ &
rgsoh(i,j))
rcl(ic,i,j)=1./(henry(ic)/(10.e5*rclsh(i,j))+f0(ic)/ &
rcloh(i,j))
end do
end do
endif
end do
return
999 write(*,*) '### FLEXPART ERROR! FILE ###'
write(*,*) '### surfdepo.t DOES NOT EXIST. ###'
stop
end subroutine readdepo
|