1 | ! radiation_mcica_sw.F90 - Monte-Carlo Independent Column Approximation shortwave solver |
---|
2 | ! |
---|
3 | ! (C) Copyright 2015- 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 | ! 2017-04-11 R. Hogan Receive albedos at g-points |
---|
17 | ! 2017-04-22 R. Hogan Store surface fluxes at all g-points |
---|
18 | ! 2017-10-23 R. Hogan Renamed single-character variables |
---|
19 | |
---|
20 | #include "ecrad_config.h" |
---|
21 | |
---|
22 | module radiation_mcica_sw |
---|
23 | |
---|
24 | public |
---|
25 | |
---|
26 | contains |
---|
27 | |
---|
28 | ! Provides elemental function "delta_eddington" |
---|
29 | #include "radiation_delta_eddington.h" |
---|
30 | |
---|
31 | !--------------------------------------------------------------------- |
---|
32 | ! Shortwave Monte Carlo Independent Column Approximation |
---|
33 | ! (McICA). This implementation performs a clear-sky and a cloudy-sky |
---|
34 | ! calculation, and then weights the two to get the all-sky fluxes |
---|
35 | ! according to the total cloud cover. This method reduces noise for |
---|
36 | ! low cloud cover situations, and exploits the clear-sky |
---|
37 | ! calculations that are usually performed for diagnostic purposes |
---|
38 | ! simultaneously. The cloud generator has been carefully written |
---|
39 | ! such that the stochastic cloud field satisfies the prescribed |
---|
40 | ! overlap parameter accounting for this weighting. |
---|
41 | subroutine solver_mcica_sw(nlev,istartcol,iendcol, & |
---|
42 | & config, single_level, cloud, & |
---|
43 | & od, ssa, g, od_cloud, ssa_cloud, g_cloud, & |
---|
44 | & albedo_direct, albedo_diffuse, incoming_sw, & |
---|
45 | & flux) |
---|
46 | |
---|
47 | use parkind1, only : jprb |
---|
48 | use yomhook, only : lhook, dr_hook, jphook |
---|
49 | |
---|
50 | use radiation_io, only : nulerr, radiation_abort |
---|
51 | use radiation_config, only : config_type |
---|
52 | use radiation_single_level, only : single_level_type |
---|
53 | use radiation_cloud, only : cloud_type |
---|
54 | use radiation_flux, only : flux_type |
---|
55 | use radiation_two_stream, only : calc_ref_trans_sw |
---|
56 | use radiation_adding_ica_sw, only : adding_ica_sw |
---|
57 | use radiation_cloud_generator, only: cloud_generator |
---|
58 | |
---|
59 | implicit none |
---|
60 | |
---|
61 | ! Inputs |
---|
62 | integer, intent(in) :: nlev ! number of model levels |
---|
63 | integer, intent(in) :: istartcol, iendcol ! range of columns to process |
---|
64 | type(config_type), intent(in) :: config |
---|
65 | type(single_level_type), intent(in) :: single_level |
---|
66 | type(cloud_type), intent(in) :: cloud |
---|
67 | |
---|
68 | ! Gas and aerosol optical depth, single-scattering albedo and |
---|
69 | ! asymmetry factor at each shortwave g-point |
---|
70 | real(jprb), intent(in), dimension(config%n_g_sw, nlev, istartcol:iendcol) :: & |
---|
71 | & od, ssa, g |
---|
72 | |
---|
73 | ! Cloud and precipitation optical depth, single-scattering albedo and |
---|
74 | ! asymmetry factor in each shortwave band |
---|
75 | real(jprb), intent(in), dimension(config%n_bands_sw,nlev,istartcol:iendcol) :: & |
---|
76 | & od_cloud, ssa_cloud, g_cloud |
---|
77 | |
---|
78 | ! Direct and diffuse surface albedos, and the incoming shortwave |
---|
79 | ! flux into a plane perpendicular to the incoming radiation at |
---|
80 | ! top-of-atmosphere in each of the shortwave g points |
---|
81 | real(jprb), intent(in), dimension(config%n_g_sw,istartcol:iendcol) :: & |
---|
82 | & albedo_direct, albedo_diffuse, incoming_sw |
---|
83 | |
---|
84 | ! Output |
---|
85 | type(flux_type), intent(inout):: flux |
---|
86 | |
---|
87 | ! Local variables |
---|
88 | |
---|
89 | ! Cosine of solar zenith angle |
---|
90 | real(jprb) :: cos_sza |
---|
91 | |
---|
92 | ! Diffuse reflectance and transmittance for each layer in clear |
---|
93 | ! and all skies |
---|
94 | real(jprb), dimension(config%n_g_sw, nlev) :: ref_clear, trans_clear, reflectance, transmittance |
---|
95 | |
---|
96 | ! Fraction of direct beam scattered by a layer into the upwelling |
---|
97 | ! or downwelling diffuse streams, in clear and all skies |
---|
98 | real(jprb), dimension(config%n_g_sw, nlev) :: ref_dir_clear, trans_dir_diff_clear, ref_dir, trans_dir_diff |
---|
99 | |
---|
100 | ! Transmittance for the direct beam in clear and all skies |
---|
101 | real(jprb), dimension(config%n_g_sw, nlev) :: trans_dir_dir_clear, trans_dir_dir |
---|
102 | |
---|
103 | ! Fluxes per g point |
---|
104 | real(jprb), dimension(config%n_g_sw, nlev+1) :: flux_up, flux_dn_diffuse, flux_dn_direct |
---|
105 | |
---|
106 | ! Combined gas+aerosol+cloud optical depth, single scattering |
---|
107 | ! albedo and asymmetry factor |
---|
108 | real(jprb), dimension(config%n_g_sw) :: od_total, ssa_total, g_total |
---|
109 | |
---|
110 | ! Combined scattering optical depth |
---|
111 | real(jprb) :: scat_od |
---|
112 | |
---|
113 | ! Optical depth scaling from the cloud generator, zero indicating |
---|
114 | ! clear skies |
---|
115 | real(jprb), dimension(config%n_g_sw,nlev) :: od_scaling |
---|
116 | |
---|
117 | ! Modified optical depth after McICA scaling to represent cloud |
---|
118 | ! inhomogeneity |
---|
119 | real(jprb), dimension(config%n_g_sw) :: od_cloud_new |
---|
120 | |
---|
121 | ! Total cloud cover output from the cloud generator |
---|
122 | real(jprb) :: total_cloud_cover |
---|
123 | |
---|
124 | ! Temporary storage for more efficient summation |
---|
125 | #ifdef DWD_REDUCTION_OPTIMIZATIONS |
---|
126 | real(jprb), dimension(nlev+1,3) :: sum_aux |
---|
127 | #else |
---|
128 | real(jprb) :: sum_up, sum_dn_diff, sum_dn_dir |
---|
129 | #endif |
---|
130 | |
---|
131 | ! Number of g points |
---|
132 | integer :: ng |
---|
133 | |
---|
134 | ! Loop indices for level, column and g point |
---|
135 | integer :: jlev, jcol, jg |
---|
136 | |
---|
137 | real(jphook) :: hook_handle |
---|
138 | |
---|
139 | if (lhook) call dr_hook('radiation_mcica_sw:solver_mcica_sw',0,hook_handle) |
---|
140 | |
---|
141 | if (.not. config%do_clear) then |
---|
142 | write(nulerr,'(a)') '*** Error: shortwave McICA requires clear-sky calculation to be performed' |
---|
143 | call radiation_abort() |
---|
144 | end if |
---|
145 | |
---|
146 | ng = config%n_g_sw |
---|
147 | |
---|
148 | ! Loop through columns |
---|
149 | do jcol = istartcol,iendcol |
---|
150 | ! Only perform calculation if sun above the horizon |
---|
151 | if (single_level%cos_sza(jcol) > 0.0_jprb) then |
---|
152 | cos_sza = single_level%cos_sza(jcol) |
---|
153 | |
---|
154 | ! Clear-sky calculation - first compute clear-sky reflectance, |
---|
155 | ! transmittance etc at each model level |
---|
156 | if (.not. config%do_sw_delta_scaling_with_gases) then |
---|
157 | ! Delta-Eddington scaling has already been performed to the |
---|
158 | ! aerosol part of od, ssa and g |
---|
159 | call calc_ref_trans_sw(ng*nlev, & |
---|
160 | & cos_sza, od(:,:,jcol), ssa(:,:,jcol), g(:,:,jcol), & |
---|
161 | & ref_clear, trans_clear, & |
---|
162 | & ref_dir_clear, trans_dir_diff_clear, & |
---|
163 | & trans_dir_dir_clear) |
---|
164 | else |
---|
165 | ! Apply delta-Eddington scaling to the aerosol-gas mixture |
---|
166 | do jlev = 1,nlev |
---|
167 | od_total = od(:,jlev,jcol) |
---|
168 | ssa_total = ssa(:,jlev,jcol) |
---|
169 | g_total = g(:,jlev,jcol) |
---|
170 | call delta_eddington(od_total, ssa_total, g_total) |
---|
171 | call calc_ref_trans_sw(ng, & |
---|
172 | & cos_sza, od_total, ssa_total, g_total, & |
---|
173 | & ref_clear(:,jlev), trans_clear(:,jlev), & |
---|
174 | & ref_dir_clear(:,jlev), trans_dir_diff_clear(:,jlev), & |
---|
175 | & trans_dir_dir_clear(:,jlev) ) |
---|
176 | end do |
---|
177 | end if |
---|
178 | |
---|
179 | ! Use adding method to compute fluxes |
---|
180 | call adding_ica_sw(ng, nlev, incoming_sw(:,jcol), & |
---|
181 | & albedo_diffuse(:,jcol), albedo_direct(:,jcol), spread(cos_sza,1,ng), & |
---|
182 | & ref_clear, trans_clear, ref_dir_clear, trans_dir_diff_clear, & |
---|
183 | & trans_dir_dir_clear, flux_up, flux_dn_diffuse, flux_dn_direct) |
---|
184 | |
---|
185 | ! Sum over g-points to compute and save clear-sky broadband |
---|
186 | ! fluxes. Note that the built-in "sum" function is very slow, |
---|
187 | ! and before being replaced by the alternatives below |
---|
188 | ! accounted for around 40% of the total cost of this routine. |
---|
189 | #ifdef DWD_REDUCTION_OPTIMIZATIONS |
---|
190 | ! Optimized summation for the NEC architecture |
---|
191 | sum_aux(:,:) = 0.0_jprb |
---|
192 | do jg = 1,ng |
---|
193 | do jlev = 1,nlev+1 |
---|
194 | sum_aux(jlev,1) = sum_aux(jlev,1) + flux_up(jg,jlev) |
---|
195 | sum_aux(jlev,2) = sum_aux(jlev,2) + flux_dn_direct(jg,jlev) |
---|
196 | sum_aux(jlev,3) = sum_aux(jlev,3) + flux_dn_diffuse(jg,jlev) |
---|
197 | end do |
---|
198 | end do |
---|
199 | flux%sw_up_clear(jcol,:) = sum_aux(:,1) |
---|
200 | flux%sw_dn_clear(jcol,:) = sum_aux(:,2) + sum_aux(:,3) |
---|
201 | if (allocated(flux%sw_dn_direct_clear)) then |
---|
202 | flux%sw_dn_direct_clear(jcol,:) = sum_aux(:,2) |
---|
203 | end if |
---|
204 | #else |
---|
205 | ! Optimized summation for the x86-64 architecture |
---|
206 | do jlev = 1,nlev+1 |
---|
207 | sum_up = 0.0_jprb |
---|
208 | sum_dn_diff = 0.0_jprb |
---|
209 | sum_dn_dir = 0.0_jprb |
---|
210 | !$omp simd reduction(+:sum_up, sum_dn_diff, sum_dn_dir) |
---|
211 | do jg = 1,ng |
---|
212 | sum_up = sum_up + flux_up(jg,jlev) |
---|
213 | sum_dn_diff = sum_dn_diff + flux_dn_diffuse(jg,jlev) |
---|
214 | sum_dn_dir = sum_dn_dir + flux_dn_direct(jg,jlev) |
---|
215 | end do |
---|
216 | flux%sw_up_clear(jcol,jlev) = sum_up |
---|
217 | flux%sw_dn_clear(jcol,jlev) = sum_dn_diff + sum_dn_dir |
---|
218 | if (allocated(flux%sw_dn_direct_clear)) then |
---|
219 | flux%sw_dn_direct_clear(jcol,jlev) = sum_dn_dir |
---|
220 | end if |
---|
221 | end do |
---|
222 | #endif |
---|
223 | |
---|
224 | ! Store spectral downwelling fluxes at surface |
---|
225 | do jg = 1,ng |
---|
226 | flux%sw_dn_diffuse_surf_clear_g(jg,jcol) = flux_dn_diffuse(jg,nlev+1) |
---|
227 | flux%sw_dn_direct_surf_clear_g(jg,jcol) = flux_dn_direct(jg,nlev+1) |
---|
228 | end do |
---|
229 | |
---|
230 | ! Do cloudy-sky calculation |
---|
231 | call cloud_generator(ng, nlev, config%i_overlap_scheme, & |
---|
232 | & single_level%iseed(jcol), & |
---|
233 | & config%cloud_fraction_threshold, & |
---|
234 | & cloud%fraction(jcol,:), cloud%overlap_param(jcol,:), & |
---|
235 | & config%cloud_inhom_decorr_scaling, cloud%fractional_std(jcol,:), & |
---|
236 | & config%pdf_sampler, od_scaling, total_cloud_cover, & |
---|
237 | & use_beta_overlap=config%use_beta_overlap, & |
---|
238 | & use_vectorizable_generator=config%use_vectorizable_generator) |
---|
239 | |
---|
240 | ! Store total cloud cover |
---|
241 | flux%cloud_cover_sw(jcol) = total_cloud_cover |
---|
242 | |
---|
243 | if (total_cloud_cover >= config%cloud_fraction_threshold) then |
---|
244 | ! Total-sky calculation |
---|
245 | do jlev = 1,nlev |
---|
246 | ! Compute combined gas+aerosol+cloud optical properties |
---|
247 | if (cloud%fraction(jcol,jlev) >= config%cloud_fraction_threshold) then |
---|
248 | do jg = 1,ng |
---|
249 | od_cloud_new(jg) = od_scaling(jg,jlev) & |
---|
250 | & * od_cloud(config%i_band_from_reordered_g_sw(jg),jlev,jcol) |
---|
251 | od_total(jg) = od(jg,jlev,jcol) + od_cloud_new(jg) |
---|
252 | ssa_total(jg) = 0.0_jprb |
---|
253 | g_total(jg) = 0.0_jprb |
---|
254 | |
---|
255 | ! In single precision we need to protect against the |
---|
256 | ! case that od_total > 0.0 and ssa_total > 0.0 but |
---|
257 | ! od_total*ssa_total == 0 due to underflow |
---|
258 | if (od_total(jg) > 0.0_jprb) then |
---|
259 | scat_od = ssa(jg,jlev,jcol)*od(jg,jlev,jcol) & |
---|
260 | & + ssa_cloud(config%i_band_from_reordered_g_sw(jg),jlev,jcol) & |
---|
261 | & * od_cloud_new(jg) |
---|
262 | ssa_total(jg) = scat_od / od_total(jg) |
---|
263 | if (scat_od > 0.0_jprb) then |
---|
264 | g_total(jg) = (g(jg,jlev,jcol)*ssa(jg,jlev,jcol)*od(jg,jlev,jcol) & |
---|
265 | & + g_cloud(config%i_band_from_reordered_g_sw(jg),jlev,jcol) & |
---|
266 | & * ssa_cloud(config%i_band_from_reordered_g_sw(jg),jlev,jcol) & |
---|
267 | & * od_cloud_new(jg)) & |
---|
268 | & / scat_od |
---|
269 | end if |
---|
270 | end if |
---|
271 | end do |
---|
272 | |
---|
273 | ! Apply delta-Eddington scaling to the cloud-aerosol-gas |
---|
274 | ! mixture |
---|
275 | if (config%do_sw_delta_scaling_with_gases) then |
---|
276 | call delta_eddington(od_total, ssa_total, g_total) |
---|
277 | end if |
---|
278 | |
---|
279 | ! Compute cloudy-sky reflectance, transmittance etc at |
---|
280 | ! each model level |
---|
281 | call calc_ref_trans_sw(ng, & |
---|
282 | & cos_sza, od_total, ssa_total, g_total, & |
---|
283 | & reflectance(:,jlev), transmittance(:,jlev), & |
---|
284 | & ref_dir(:,jlev), trans_dir_diff(:,jlev), & |
---|
285 | & trans_dir_dir(:,jlev)) |
---|
286 | |
---|
287 | else |
---|
288 | ! Clear-sky layer: copy over clear-sky values |
---|
289 | do jg = 1,ng |
---|
290 | reflectance(jg,jlev) = ref_clear(jg,jlev) |
---|
291 | transmittance(jg,jlev) = trans_clear(jg,jlev) |
---|
292 | ref_dir(jg,jlev) = ref_dir_clear(jg,jlev) |
---|
293 | trans_dir_diff(jg,jlev) = trans_dir_diff_clear(jg,jlev) |
---|
294 | trans_dir_dir(jg,jlev) = trans_dir_dir_clear(jg,jlev) |
---|
295 | end do |
---|
296 | end if |
---|
297 | end do |
---|
298 | |
---|
299 | ! Use adding method to compute fluxes for an overcast sky |
---|
300 | call adding_ica_sw(ng, nlev, incoming_sw(:,jcol), & |
---|
301 | & albedo_diffuse(:,jcol), albedo_direct(:,jcol), spread(cos_sza,1,ng), & |
---|
302 | & reflectance, transmittance, ref_dir, trans_dir_diff, & |
---|
303 | & trans_dir_dir, flux_up, flux_dn_diffuse, flux_dn_direct) |
---|
304 | |
---|
305 | ! Store overcast broadband fluxes |
---|
306 | #ifdef DWD_REDUCTION_OPTIMIZATIONS |
---|
307 | sum_aux(:,:) = 0.0_jprb |
---|
308 | do jg = 1,ng |
---|
309 | do jlev = 1,nlev+1 |
---|
310 | sum_aux(jlev,1) = sum_aux(jlev,1) + flux_up(jg,jlev) |
---|
311 | sum_aux(jlev,2) = sum_aux(jlev,2) + flux_dn_direct(jg,jlev) |
---|
312 | sum_aux(jlev,3) = sum_aux(jlev,3) + flux_dn_diffuse(jg,jlev) |
---|
313 | end do |
---|
314 | end do |
---|
315 | flux%sw_up(jcol,:) = sum_aux(:,1) |
---|
316 | flux%sw_dn(jcol,:) = sum_aux(:,2) + sum_aux(:,3) |
---|
317 | if (allocated(flux%sw_dn_direct)) then |
---|
318 | flux%sw_dn_direct(jcol,:) = sum_aux(:,2) |
---|
319 | end if |
---|
320 | #else |
---|
321 | do jlev = 1,nlev+1 |
---|
322 | sum_up = 0.0_jprb |
---|
323 | sum_dn_diff = 0.0_jprb |
---|
324 | sum_dn_dir = 0.0_jprb |
---|
325 | !$omp simd reduction(+:sum_up, sum_dn_diff, sum_dn_dir) |
---|
326 | do jg = 1,ng |
---|
327 | sum_up = sum_up + flux_up(jg,jlev) |
---|
328 | sum_dn_diff = sum_dn_diff + flux_dn_diffuse(jg,jlev) |
---|
329 | sum_dn_dir = sum_dn_dir + flux_dn_direct(jg,jlev) |
---|
330 | end do |
---|
331 | flux%sw_up(jcol,jlev) = sum_up |
---|
332 | flux%sw_dn(jcol,jlev) = sum_dn_diff + sum_dn_dir |
---|
333 | if (allocated(flux%sw_dn_direct)) then |
---|
334 | flux%sw_dn_direct(jcol,jlev) = sum_dn_dir |
---|
335 | end if |
---|
336 | end do |
---|
337 | #endif |
---|
338 | |
---|
339 | ! Cloudy flux profiles currently assume completely overcast |
---|
340 | ! skies; perform weighted average with clear-sky profile |
---|
341 | do jlev = 1, nlev+1 |
---|
342 | flux%sw_up(jcol,jlev) = total_cloud_cover *flux%sw_up(jcol,jlev) & |
---|
343 | & + (1.0_jprb - total_cloud_cover)*flux%sw_up_clear(jcol,jlev) |
---|
344 | flux%sw_dn(jcol,jlev) = total_cloud_cover *flux%sw_dn(jcol,jlev) & |
---|
345 | & + (1.0_jprb - total_cloud_cover)*flux%sw_dn_clear(jcol,jlev) |
---|
346 | if (allocated(flux%sw_dn_direct)) then |
---|
347 | flux%sw_dn_direct(jcol,jlev) = total_cloud_cover *flux%sw_dn_direct(jcol,jlev) & |
---|
348 | & + (1.0_jprb - total_cloud_cover)*flux%sw_dn_direct_clear(jcol,jlev) |
---|
349 | end if |
---|
350 | end do |
---|
351 | ! Likewise for surface spectral fluxes |
---|
352 | do jg = 1,ng |
---|
353 | flux%sw_dn_diffuse_surf_g(jg,jcol) = flux_dn_diffuse(jg,nlev+1) |
---|
354 | flux%sw_dn_direct_surf_g(jg,jcol) = flux_dn_direct(jg,nlev+1) |
---|
355 | flux%sw_dn_diffuse_surf_g(jg,jcol) = total_cloud_cover *flux%sw_dn_diffuse_surf_g(jg,jcol) & |
---|
356 | & + (1.0_jprb - total_cloud_cover)*flux%sw_dn_diffuse_surf_clear_g(jg,jcol) |
---|
357 | flux%sw_dn_direct_surf_g(jg,jcol) = total_cloud_cover *flux%sw_dn_direct_surf_g(jg,jcol) & |
---|
358 | & + (1.0_jprb - total_cloud_cover)*flux%sw_dn_direct_surf_clear_g(jg,jcol) |
---|
359 | end do |
---|
360 | |
---|
361 | else |
---|
362 | ! No cloud in profile and clear-sky fluxes already |
---|
363 | ! calculated: copy them over |
---|
364 | do jlev = 1, nlev+1 |
---|
365 | flux%sw_up(jcol,jlev) = flux%sw_up_clear(jcol,jlev) |
---|
366 | flux%sw_dn(jcol,jlev) = flux%sw_dn_clear(jcol,jlev) |
---|
367 | if (allocated(flux%sw_dn_direct)) then |
---|
368 | flux%sw_dn_direct(jcol,jlev) = flux%sw_dn_direct_clear(jcol,jlev) |
---|
369 | end if |
---|
370 | end do |
---|
371 | do jg = 1,ng |
---|
372 | flux%sw_dn_diffuse_surf_g(jg,jcol) = flux%sw_dn_diffuse_surf_clear_g(jg,jcol) |
---|
373 | flux%sw_dn_direct_surf_g(jg,jcol) = flux%sw_dn_direct_surf_clear_g(jg,jcol) |
---|
374 | end do |
---|
375 | |
---|
376 | end if ! Cloud is present in profile |
---|
377 | |
---|
378 | else |
---|
379 | ! Set fluxes to zero if sun is below the horizon |
---|
380 | do jlev = 1, nlev+1 |
---|
381 | flux%sw_up(jcol,jlev) = 0.0_jprb |
---|
382 | flux%sw_dn(jcol,jlev) = 0.0_jprb |
---|
383 | if (allocated(flux%sw_dn_direct)) then |
---|
384 | flux%sw_dn_direct(jcol,jlev) = 0.0_jprb |
---|
385 | end if |
---|
386 | flux%sw_up_clear(jcol,jlev) = 0.0_jprb |
---|
387 | flux%sw_dn_clear(jcol,jlev) = 0.0_jprb |
---|
388 | if (allocated(flux%sw_dn_direct_clear)) then |
---|
389 | flux%sw_dn_direct_clear(jcol,jlev) = 0.0_jprb |
---|
390 | end if |
---|
391 | end do |
---|
392 | do jg = 1,ng |
---|
393 | flux%sw_dn_diffuse_surf_g(jg,jcol) = 0.0_jprb |
---|
394 | flux%sw_dn_direct_surf_g(jg,jcol) = 0.0_jprb |
---|
395 | flux%sw_dn_diffuse_surf_clear_g(jg,jcol) = 0.0_jprb |
---|
396 | flux%sw_dn_direct_surf_clear_g(jg,jcol) = 0.0_jprb |
---|
397 | end do |
---|
398 | end if ! Sun above horizon |
---|
399 | |
---|
400 | end do ! Loop over columns |
---|
401 | |
---|
402 | if (lhook) call dr_hook('radiation_mcica_sw:solver_mcica_sw',1,hook_handle) |
---|
403 | |
---|
404 | end subroutine solver_mcica_sw |
---|
405 | |
---|
406 | end module radiation_mcica_sw |
---|