1 | ! radiation_aerosol.F90 - Derived type describing aerosol |
---|
2 | ! |
---|
3 | ! (C) Copyright 2014- ECMWF. |
---|
4 | ! |
---|
5 | ! This software is licensed under the terms of the Apache Licence Version 2.0 |
---|
6 | ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0. |
---|
7 | ! |
---|
8 | ! In applying this licence, ECMWF does not waive the privileges and immunities |
---|
9 | ! granted to it by virtue of its status as an intergovernmental organisation |
---|
10 | ! nor does it submit to any jurisdiction. |
---|
11 | ! |
---|
12 | ! Author: Robin Hogan |
---|
13 | ! Email: r.j.hogan@ecmwf.int |
---|
14 | ! |
---|
15 | ! Modifications |
---|
16 | ! 2018-04-15 R. Hogan Add "direct" option |
---|
17 | ! 2019-01-14 R. Hogan Added out_of_physical_bounds routine |
---|
18 | |
---|
19 | module radiation_aerosol |
---|
20 | |
---|
21 | use parkind1, only : jprb |
---|
22 | |
---|
23 | implicit none |
---|
24 | public |
---|
25 | |
---|
26 | !--------------------------------------------------------------------- |
---|
27 | ! Type describing the aerosol content in the atmosphere |
---|
28 | type aerosol_type |
---|
29 | ! The mass mixing ratio of config%n_aerosol_types different |
---|
30 | ! aerosol types dimensioned |
---|
31 | ! (ncol,istartlev:iendlev,config%n_aerosol_types), where ncol is |
---|
32 | ! the number of columns, istartlev:iendlev is the range of model |
---|
33 | ! levels where aerosols are present |
---|
34 | real(jprb), allocatable, dimension(:,:,:) :: & |
---|
35 | & mixing_ratio ! mass mixing ratio (kg/kg) |
---|
36 | |
---|
37 | ! Alternatively, if is_direct=true, the optical properties are |
---|
38 | ! provided directly and are dimensioned |
---|
39 | ! (nband,istartlev:iendlev,ncol) |
---|
40 | real(jprb), allocatable, dimension(:,:,:) :: & |
---|
41 | & od_sw, ssa_sw, g_sw, & ! Shortwave optical properties |
---|
42 | & od_lw, ssa_lw, g_lw ! Longwave optical properties |
---|
43 | |
---|
44 | ! Range of levels in which the aerosol properties are provided |
---|
45 | integer :: istartlev, iendlev |
---|
46 | |
---|
47 | ! Are the optical properties going to be provided directly by the |
---|
48 | ! user? |
---|
49 | logical :: is_direct = .false. |
---|
50 | |
---|
51 | contains |
---|
52 | procedure :: allocate => allocate_aerosol_arrays |
---|
53 | procedure :: allocate_direct => allocate_aerosol_arrays_direct |
---|
54 | procedure :: deallocate => deallocate_aerosol_arrays |
---|
55 | procedure :: out_of_physical_bounds |
---|
56 | end type aerosol_type |
---|
57 | |
---|
58 | contains |
---|
59 | |
---|
60 | !--------------------------------------------------------------------- |
---|
61 | ! Allocate array for describing aerosols, although in the offline |
---|
62 | ! code these are allocated when they are read from the NetCDF file |
---|
63 | subroutine allocate_aerosol_arrays(this, ncol, istartlev, iendlev, ntype) |
---|
64 | |
---|
65 | use yomhook, only : lhook, dr_hook, jphook |
---|
66 | |
---|
67 | class(aerosol_type), intent(inout) :: this |
---|
68 | integer, intent(in) :: ncol ! Number of columns |
---|
69 | integer, intent(in) :: istartlev, iendlev ! Level range |
---|
70 | integer, intent(in) :: ntype ! Number of aerosol types |
---|
71 | real(jphook) :: hook_handle |
---|
72 | |
---|
73 | if (lhook) call dr_hook('radiation_aerosol:allocate',0,hook_handle) |
---|
74 | |
---|
75 | allocate(this%mixing_ratio(ncol,istartlev:iendlev,ntype)) |
---|
76 | this%is_direct = .false. |
---|
77 | this%istartlev = istartlev |
---|
78 | this%iendlev = iendlev |
---|
79 | |
---|
80 | if (lhook) call dr_hook('radiation_aerosol:allocate',1,hook_handle) |
---|
81 | |
---|
82 | end subroutine allocate_aerosol_arrays |
---|
83 | |
---|
84 | |
---|
85 | !--------------------------------------------------------------------- |
---|
86 | ! Allocate arrays for describing aerosol optical properties |
---|
87 | subroutine allocate_aerosol_arrays_direct(this, config, & |
---|
88 | & ncol, istartlev, iendlev) |
---|
89 | |
---|
90 | use yomhook, only : lhook, dr_hook, jphook |
---|
91 | use radiation_config, only : config_type |
---|
92 | |
---|
93 | class(aerosol_type), intent(inout) :: this |
---|
94 | type(config_type), intent(in) :: config |
---|
95 | integer, intent(in) :: ncol ! Number of columns |
---|
96 | integer, intent(in) :: istartlev, iendlev ! Level range |
---|
97 | |
---|
98 | real(jphook) :: hook_handle |
---|
99 | |
---|
100 | if (lhook) call dr_hook('radiation_aerosol:allocate_direct',0,hook_handle) |
---|
101 | |
---|
102 | this%is_direct = .true. |
---|
103 | this%istartlev = istartlev |
---|
104 | this%iendlev = iendlev |
---|
105 | |
---|
106 | if (config%do_sw) then |
---|
107 | allocate(this%od_sw (config%n_bands_sw,istartlev:iendlev,ncol)) |
---|
108 | allocate(this%ssa_sw(config%n_bands_sw,istartlev:iendlev,ncol)) |
---|
109 | allocate(this%g_sw (config%n_bands_sw,istartlev:iendlev,ncol)) |
---|
110 | end if |
---|
111 | |
---|
112 | if (config%do_lw) then |
---|
113 | allocate(this%od_lw (config%n_bands_lw,istartlev:iendlev,ncol)) |
---|
114 | allocate(this%ssa_lw(config%n_bands_lw,istartlev:iendlev,ncol)) |
---|
115 | allocate(this%g_lw (config%n_bands_lw,istartlev:iendlev,ncol)) |
---|
116 | ! If longwave scattering by aerosol is not to be represented, |
---|
117 | ! then the user may wish to just provide absorption optical |
---|
118 | ! depth in od_lw, in which case we must set the following two |
---|
119 | ! variables to zero |
---|
120 | this%ssa_lw = 0.0_jprb |
---|
121 | this%g_lw = 0.0_jprb |
---|
122 | end if |
---|
123 | |
---|
124 | if (lhook) call dr_hook('radiation_aerosol:allocate_direct',1,hook_handle) |
---|
125 | |
---|
126 | end subroutine allocate_aerosol_arrays_direct |
---|
127 | |
---|
128 | |
---|
129 | !--------------------------------------------------------------------- |
---|
130 | ! Deallocate arrays |
---|
131 | subroutine deallocate_aerosol_arrays(this) |
---|
132 | |
---|
133 | use yomhook, only : lhook, dr_hook, jphook |
---|
134 | |
---|
135 | class(aerosol_type), intent(inout) :: this |
---|
136 | |
---|
137 | real(jphook) :: hook_handle |
---|
138 | |
---|
139 | if (lhook) call dr_hook('radiation_aerosol:deallocate',0,hook_handle) |
---|
140 | |
---|
141 | if (allocated(this%mixing_ratio)) deallocate(this%mixing_ratio) |
---|
142 | if (allocated(this%od_sw)) deallocate(this%od_sw) |
---|
143 | if (allocated(this%ssa_sw)) deallocate(this%ssa_sw) |
---|
144 | if (allocated(this%g_sw)) deallocate(this%g_sw) |
---|
145 | if (allocated(this%od_lw)) deallocate(this%od_lw) |
---|
146 | if (allocated(this%ssa_lw)) deallocate(this%ssa_lw) |
---|
147 | if (allocated(this%g_lw)) deallocate(this%g_lw) |
---|
148 | |
---|
149 | if (lhook) call dr_hook('radiation_aerosol:deallocate',1,hook_handle) |
---|
150 | |
---|
151 | end subroutine deallocate_aerosol_arrays |
---|
152 | |
---|
153 | |
---|
154 | !--------------------------------------------------------------------- |
---|
155 | ! Return .true. if variables are out of a physically sensible range, |
---|
156 | ! optionally only considering columns between istartcol and iendcol |
---|
157 | function out_of_physical_bounds(this, istartcol, iendcol, do_fix) result(is_bad) |
---|
158 | |
---|
159 | use yomhook, only : lhook, dr_hook, jphook |
---|
160 | use radiation_check, only : out_of_bounds_3d |
---|
161 | |
---|
162 | class(aerosol_type), intent(inout) :: this |
---|
163 | integer, optional,intent(in) :: istartcol, iendcol |
---|
164 | logical, optional,intent(in) :: do_fix |
---|
165 | logical :: is_bad |
---|
166 | |
---|
167 | logical :: do_fix_local |
---|
168 | |
---|
169 | real(jphook) :: hook_handle |
---|
170 | |
---|
171 | if (lhook) call dr_hook('radiation_aerosol:out_of_physical_bounds',0,hook_handle) |
---|
172 | |
---|
173 | if (present(do_fix)) then |
---|
174 | do_fix_local = do_fix |
---|
175 | else |
---|
176 | do_fix_local = .false. |
---|
177 | end if |
---|
178 | |
---|
179 | is_bad = out_of_bounds_3d(this%mixing_ratio, 'aerosol%mixing_ratio', & |
---|
180 | & 0.0_jprb, 1.0_jprb, do_fix_local, i1=istartcol, i2=iendcol) & |
---|
181 | & .or. out_of_bounds_3d(this%od_sw, 'aerosol%od_sw', & |
---|
182 | & 0.0_jprb, 100.0_jprb, do_fix_local, k1=istartcol, k2=iendcol) & |
---|
183 | & .or. out_of_bounds_3d(this%od_lw, 'aerosol%od_lw', & |
---|
184 | & 0.0_jprb, 100.0_jprb, do_fix_local, k1=istartcol, k2=iendcol) & |
---|
185 | & .or. out_of_bounds_3d(this%ssa_sw, 'aerosol%ssa_sw', & |
---|
186 | & 0.0_jprb, 1.0_jprb, do_fix_local, k1=istartcol, k2=iendcol) & |
---|
187 | & .or. out_of_bounds_3d(this%ssa_lw, 'aerosol%ssa_lw', & |
---|
188 | & 0.0_jprb, 1.0_jprb, do_fix_local, k1=istartcol, k2=iendcol) & |
---|
189 | & .or. out_of_bounds_3d(this%g_sw, 'aerosol%g_sw', & |
---|
190 | & 0.0_jprb, 1.0_jprb, do_fix_local, k1=istartcol, k2=iendcol) & |
---|
191 | & .or. out_of_bounds_3d(this%g_lw, 'aerosol%g_lw', & |
---|
192 | & 0.0_jprb, 1.0_jprb, do_fix_local, k1=istartcol, k2=iendcol) |
---|
193 | |
---|
194 | if (lhook) call dr_hook('radiation_aerosol:out_of_physical_bounds',1,hook_handle) |
---|
195 | |
---|
196 | end function out_of_physical_bounds |
---|
197 | |
---|
198 | end module radiation_aerosol |
---|