-
Notifications
You must be signed in to change notification settings - Fork 9
/
evap_mono.F90
126 lines (106 loc) · 4.09 KB
/
evap_mono.F90
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
! Include shortname defintions, so that the F77 code does not have to be modified to
! reference the CARMA structure.
#include "carma_globaer.h"
!! This routine calculates particle source terms <evappe> due to total
!! evaporation from bin <ibin> group <ig> into a monodisperse
!! distribution.
!!
!! Distinct evaporation of cores has not been treated.
!!
!! @author Andy Ackerman
!! @version Aug-2001
subroutine evap_mono(carma,cstate,iz,ibin,ig,iavg,ieto,igto,rc)
! types
use carma_precision_mod
use carma_enums_mod
use carma_constants_mod
use carma_types_mod
use carmastate_mod
use carma_mod
implicit none
type(carma_type), intent(in) :: carma !! the carma object
type(carmastate_type), intent(inout) :: cstate !! the carma state object
integer, intent(in) :: iz !! z index
integer, intent(in) :: ibin !! bin index
integer, intent(in) :: ig !! group index
integer, intent(in) :: iavg
integer, intent(in) :: ieto
integer, intent(in) :: igto
integer, intent(inout) :: rc !! return code, negative indicates failure
! Local declarations
integer :: ic
integer :: iecore
integer :: ie2cn
integer :: jbin
logical :: conserve_mass
real(kind=f) :: factor
real(kind=f) :: fracmass
! Define option to conserve mass or number when a choice must be made
! during monodisperse total evaporation beyond CN grid -- should be done in setupaer()
conserve_mass = .true.
! Set automatic flag for total evaporation used in gasexchange()
totevap(ibin,ig) = .true.
! Possibly put all of core mass into largest, smallest, or
! smallest nucelated CN bin
if( too_big .or. too_small .or. nuc_small )then
if( too_big )then
jbin = NBIN
elseif( too_small )then
jbin = 1
else
jbin = 1
endif
if( conserve_mass )then
factor = coreavg/rmass(jbin,igto)
else
factor = ONE
endif
! First the CN number concentration element
evappe(jbin,ieto) = evappe(jbin,ieto) + factor*evdrop
! Now the CN cores
do ic = 2, ncore(ig)
iecore = icorelem(ic,ig)
ie2cn = ievp2elem(iecore)
! It is possible to have coremasses in the particle
! that don't participate in nucleation. If there is no
! evp2elem defined, then skip this part.
!
! NOTE: This cam up in a sensitivity test where there were
! two nucleation cores (dust and sulfates) and there was
! a desire to turn off one of the nucleation mechanisms for
! sensitivity testing.
if (ie2cn .gt. 0) then
evappe(jbin,ie2cn) = evappe(jbin,ie2cn) + &
factor*evcore(ic)*rmass(jbin,igto)
end if
enddo
else
! Partition core mass between two CN bins, conserving total core mass
! and number. The number will be subdivided into bins <iavg> and <iavg>-1.
if( iavg .le. 1 .or. iavg .gt. NBIN )then
if (do_print) write(LUNOPRT, *) "evap_mono: bad iavg = , ", iavg
rc = RC_ERROR
return
endif
fracmass = ( rmass(iavg,igto) - coreavg ) / diffmass(iavg,igto,iavg-1,igto)
! fracmass = max( 0._f, min( ONE, fracmass ) )
! First the CN number concentration element
evappe(iavg-1,ieto) = evappe(iavg-1,ieto) + evdrop*fracmass
evappe(iavg,ieto) = evappe(iavg,ieto) + evdrop*( ONE - fracmass )
! Now the cores
do ic = 2, ncore(ig)
iecore = icorelem(ic,ig)
ie2cn = ievp2elem(iecore)
! It is possible to have coremasses in the particle
! that don't participate in nucleation. If there is no
! evp2elem defined, then skip this part.
if (ie2cn .gt. 0) then
evappe(iavg-1,ie2cn) = evappe(iavg-1,ie2cn) + &
rmass(iavg-1,igto)*evcore(ic)*fracmass
evappe(iavg,ie2cn) = evappe(iavg,ie2cn) + &
rmass(iavg,igto)*evcore(ic)*( ONE - fracmass )
end if
enddo
endif
return
end