1 | ! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt |
---|
2 | ! SPDX-License-Identifier: GPL-3.0-or-later |
---|
3 | |
---|
4 | subroutine readdepo |
---|
5 | |
---|
6 | !***************************************************************************** |
---|
7 | ! * |
---|
8 | ! Reads dry deposition parameters needed by the procedure of Wesely (1989). * |
---|
9 | ! Wesely (1989): Parameterization of surface resistances to gaseous * |
---|
10 | ! dry deposition in regional-scale numerical models. * |
---|
11 | ! Atmos. Environ. 23, 1293-1304. * |
---|
12 | ! * |
---|
13 | ! * |
---|
14 | ! AUTHOR: Andreas Stohl, 19 May 1995 * |
---|
15 | ! * |
---|
16 | !***************************************************************************** |
---|
17 | ! * |
---|
18 | ! Variables: * |
---|
19 | ! * |
---|
20 | ! rcl(maxspec,5,9) [s/m] Lower canopy resistance * |
---|
21 | ! rgs(maxspec,5,9) [s/m] Ground resistance * |
---|
22 | ! rlu(maxspec,5,9) [s/m] Leaf cuticular resistance * |
---|
23 | ! rm(maxspec) [s/m] Mesophyll resistance, set in readreleases * |
---|
24 | ! ri(maxspec) [s/m] Stomatal resistance * |
---|
25 | ! * |
---|
26 | ! Constants: * |
---|
27 | ! * |
---|
28 | !***************************************************************************** |
---|
29 | |
---|
30 | use par_mod |
---|
31 | use com_mod |
---|
32 | |
---|
33 | implicit none |
---|
34 | |
---|
35 | ! FOR THIS SUBROUTINE, numclass=9 IS ASSUMED |
---|
36 | !******************************************* |
---|
37 | |
---|
38 | real :: rluh(5,numclass),rgssh(5,numclass),rgsoh(5,numclass) |
---|
39 | real :: rclsh(5,numclass),rcloh(5,numclass) |
---|
40 | integer :: i,j,ic |
---|
41 | |
---|
42 | |
---|
43 | ! Read deposition constants related with landuse and seasonal category |
---|
44 | !********************************************************************* |
---|
45 | open(unitwesely,file=path(1)(1:length(1))//'surfdepo.t', & |
---|
46 | status='old',err=999) |
---|
47 | |
---|
48 | do i=1,16 |
---|
49 | read(unitwesely,*) |
---|
50 | end do |
---|
51 | do i=1,5 |
---|
52 | read(unitwesely,*) |
---|
53 | read(unitwesely,'(8x,13f8.0)') (ri(i,j),j=1,numclass) |
---|
54 | read(unitwesely,'(8x,13f8.0)') (rluh(i,j),j=1,numclass) |
---|
55 | read(unitwesely,'(8x,13f8.0)') (rac(i,j),j=1,numclass) |
---|
56 | read(unitwesely,'(8x,13f8.0)') (rgssh(i,j),j=1,numclass) |
---|
57 | read(unitwesely,'(8x,13f8.0)') (rgsoh(i,j),j=1,numclass) |
---|
58 | read(unitwesely,'(8x,13f8.0)') (rclsh(i,j),j=1,numclass) |
---|
59 | read(unitwesely,'(8x,13f8.0)') (rcloh(i,j),j=1,numclass) |
---|
60 | end do |
---|
61 | |
---|
62 | ! TEST |
---|
63 | ! do 31 i=1,5 |
---|
64 | ! ri(i,13)=ri(i,5) |
---|
65 | ! rluh(i,13)=rluh(i,5) |
---|
66 | ! rac(i,13)=rac(i,5) |
---|
67 | ! rgssh(i,13)=rgssh(i,5) |
---|
68 | ! rgsoh(i,13)=rgsoh(i,5) |
---|
69 | ! rclsh(i,13)=rclsh(i,5) |
---|
70 | ! rcloh(i,13)=rcloh(i,5) |
---|
71 | !31 continue |
---|
72 | ! TEST |
---|
73 | ! Sabine Eckhardt, Dec 06, set resistances of 9999 to 'infinite' (1E25) |
---|
74 | do i=1,5 |
---|
75 | do j=1,numclass |
---|
76 | if (ri(i,j).eq.9999.) ri(i,j)=1.E25 |
---|
77 | if (rluh(i,j).eq.9999.) rluh(i,j)=1.E25 |
---|
78 | if (rac(i,j).eq.9999.) rac(i,j)=1.E25 |
---|
79 | if (rgssh(i,j).eq.9999.) rgssh(i,j)=1.E25 |
---|
80 | if (rgsoh(i,j).eq.9999.) rgsoh(i,j)=1.E25 |
---|
81 | if (rclsh(i,j).eq.9999.) rclsh(i,j)=1.E25 |
---|
82 | if (rcloh(i,j).eq.9999.) rcloh(i,j)=1.E25 |
---|
83 | end do |
---|
84 | end do |
---|
85 | |
---|
86 | |
---|
87 | |
---|
88 | do i=1,5 |
---|
89 | do j=1,numclass |
---|
90 | ri(i,j)=max(ri(i,j),0.001) |
---|
91 | rluh(i,j)=max(rluh(i,j),0.001) |
---|
92 | rac(i,j)=max(rac(i,j),0.001) |
---|
93 | rgssh(i,j)=max(rgssh(i,j),0.001) |
---|
94 | rgsoh(i,j)=max(rgsoh(i,j),0.001) |
---|
95 | rclsh(i,j)=max(rclsh(i,j),0.001) |
---|
96 | rcloh(i,j)=max(rcloh(i,j),0.001) |
---|
97 | end do |
---|
98 | end do |
---|
99 | close(unitwesely) |
---|
100 | |
---|
101 | |
---|
102 | ! Compute additional parameters |
---|
103 | !****************************** |
---|
104 | |
---|
105 | do ic=1,nspec |
---|
106 | if (reldiff(ic).gt.0.) then ! gas is dry deposited |
---|
107 | do i=1,5 |
---|
108 | do j=1,numclass |
---|
109 | rlu(ic,i,j)=rluh(i,j)/(1.e-5*henry(ic)+f0(ic)) |
---|
110 | rgs(ic,i,j)=1./(henry(ic)/(10.e5*rgssh(i,j))+f0(ic)/ & |
---|
111 | rgsoh(i,j)) |
---|
112 | rcl(ic,i,j)=1./(henry(ic)/(10.e5*rclsh(i,j))+f0(ic)/ & |
---|
113 | rcloh(i,j)) |
---|
114 | end do |
---|
115 | end do |
---|
116 | endif |
---|
117 | end do |
---|
118 | |
---|
119 | |
---|
120 | return |
---|
121 | |
---|
122 | |
---|
123 | 999 write(*,*) '### FLEXPART ERROR! FILE ###' |
---|
124 | write(*,*) '### surfdepo.t DOES NOT EXIST. ###' |
---|
125 | stop |
---|
126 | |
---|
127 | end subroutine readdepo |
---|