source: flexpart.git/src/readdepo.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

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