source: LMDZ6/branches/contrails/libf/phylmd/ecrad/radiation/radiation_single_level.F90 @ 5456

Last change on this file since 5456 was 4773, checked in by idelkadi, 13 months ago
  • Update of Ecrad in LMDZ The same organization of the Ecrad offline version is retained in order to facilitate the updating of Ecrad in LMDZ and the comparison between online and offline results. version 1.6.1 of Ecrad (https://github.com/lguez/ecrad.git)
  • Implementation of the double call of Ecrad in LMDZ


File size: 16.5 KB
Line 
1! radiation_single_level.F90 - Derived type for single-level fields
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!   2019-01-14  R. Hogan  Added out_of_physical_bounds routine
17!   2019-01-18  R. Hogan  Added weighted albedo mapping
18
19module radiation_single_level
20
21  use parkind1, only : jprb
22
23  implicit none
24  public
25
26  !---------------------------------------------------------------------
27  ! Derived type to contain variables that don't vary with height;
28  ! mostly they are surface properties
29  type single_level_type
30    ! Note that skin temperature is only defined if
31    ! is_simple_surface==.true.
32    real(jprb), allocatable, dimension(:) :: &
33         &   cos_sza, &       ! (ncol) Cosine of solar zenith angle
34         &   skin_temperature ! (ncol) Skin temperature (K)
35
36    ! Shortwave albedo: if sw_albedo_direct is not allocated then
37    ! sw_albedo will be used for both direct and diffuse solar
38    ! radiation; otherwise, sw_albedo will be used for diffuse
39    ! radiation and sw_albedo_direct for direct radiation.
40    real(jprb), allocatable, dimension(:,:) :: &
41         &   sw_albedo, &        ! (ncol,nalbedobands)
42         &   sw_albedo_direct    ! (ncol,nalbedobands)
43   
44    ! If is_simple_surface==.true., we provide longwave emissivity
45    ! coarser spectral intervals, while if
46    ! is_simple_surface==.false. it will be generated by the surface
47    ! radiation scheme and may optionally be on the full spectral grid
48    ! of the atmospheric model.
49    real(jprb), allocatable, dimension(:,:) :: &
50         &   lw_emissivity       ! (ncol,nemissbands) If
51    ! is_simple_surface==.false. then we specify the emission instead
52    ! of the skin temperature, and again this may be on the emissivity
53    ! bands or the full spectral grid
54    real(jprb), allocatable, dimension(:,:) :: &
55         &   lw_emission         ! (ncol,nemissbands)
56
57    ! Incoming solar irradiance at the Earth is specified as the total
58    ! across the shortwave spectrum. Note that this needs to be
59    ! adjusted for Earth-Sun distance, solar cycle etc. To get
60    ! normalized fluxes out, simply set this to 1.0.
61    real(jprb) :: solar_irradiance = 1366.0_jprb ! W m-2
62
63    ! In addition to the effect of the solar cycle on total solar
64    ! irradiance (which the user is expected to input via
65    ! solar_irradiance), it can change the balance between different
66    ! parts of the solar spectrum, e.g. more UV at solar maximum. This
67    ! variable provides the scaling to be applied to the spectral
68    ! amplitude of the solar cycle (assuming it has been provided),
69    ! e.g. 1.0 for solar maximum, -1.0 for solar maximum and 0.0 for
70    ! a mean solar spectrum.
71    real(jprb) :: spectral_solar_cycle_multiplier = 0.0_jprb
72   
73    ! If config%use_spectral_solar_irradiance==true then this will be
74    ! scaled by spectral_solar_scaling
75    real(jprb), allocatable, dimension(:) :: &
76         &   spectral_solar_scaling ! (n_bands_sw)
77 
78    ! Seed for random number generator in McICA; it is expected that
79    ! the host model generates this from the model time, longitude
80    ! and latitude, in order that the model is deterministic
81    integer, allocatable, dimension(:) :: iseed ! (ncol)
82
83    ! Is the underlying surface a simple single "flat" tile? If so we
84    ! describe it with skin_temperature and lw_emissivity. Otherwise
85    ! we describe it with lw_emission and lw_emissivity coming out of
86    ! the surface radiation scheme.
87    logical :: is_simple_surface = .true.
88
89    ! If is_simple_surface==.false., do we use the full number of g
90    ! points for dimensioning sw_albedo, sw_albedo_direct and
91    ! lw_emission?
92    !logical :: use_full_spectrum_sw = .false.
93    !logical :: use_full_spectrum_lw = .true.
94
95  contains
96    procedure :: allocate   => allocate_single_level
97    procedure :: deallocate => deallocate_single_level
98    procedure :: init_seed_simple
99    procedure :: get_albedos
100    procedure :: out_of_physical_bounds
101
102  end type single_level_type
103
104contains
105 
106  !---------------------------------------------------------------------
107  ! Allocate the arrays of a single-level type
108  subroutine allocate_single_level(this, ncol, nalbedobands, nemisbands, &
109       &                           use_sw_albedo_direct, is_simple_surface)
110
111    use yomhook, only : lhook, dr_hook, jphook
112
113    class(single_level_type), intent(inout) :: this
114    integer,                  intent(in)    :: ncol, nalbedobands, nemisbands
115    logical,        optional, intent(in)    :: use_sw_albedo_direct
116    logical,        optional, intent(in)    :: is_simple_surface
117
118    real(jphook) :: hook_handle
119
120    if (lhook) call dr_hook('radiation_single_level:allocate',0,hook_handle)
121
122    call this%deallocate()
123
124    if (present(is_simple_surface)) then
125      this%is_simple_surface = is_simple_surface
126    else
127      this%is_simple_surface = .true.
128    end if
129
130    allocate(this%cos_sza(ncol))
131
132    if (this%is_simple_surface) then
133      allocate(this%skin_temperature(ncol))
134    else
135      allocate(this%lw_emission(ncol, nemisbands))
136    end if
137    allocate(this%lw_emissivity(ncol, nemisbands))
138
139    allocate(this%sw_albedo(ncol, nalbedobands))
140
141    if (present(use_sw_albedo_direct)) then
142      if (use_sw_albedo_direct) then
143        allocate(this%sw_albedo_direct(ncol, nalbedobands))
144      end if
145    end if
146
147    allocate(this%iseed(ncol))
148
149    if (lhook) call dr_hook('radiation_single_level:allocate',1,hook_handle)
150
151  end subroutine allocate_single_level
152
153
154  !---------------------------------------------------------------------
155  ! Deallocate the arrays of a single-level type
156  subroutine deallocate_single_level(this)
157
158    use yomhook, only : lhook, dr_hook, jphook
159
160    class(single_level_type), intent(inout) :: this
161
162    real(jphook) :: hook_handle
163
164    if (lhook) call dr_hook('radiation_single_level:deallocate',0,hook_handle)
165
166    if (allocated(this%cos_sza)) then
167      deallocate(this%cos_sza)
168    end if
169    if (allocated(this%skin_temperature)) then
170      deallocate(this%skin_temperature)
171    end if
172    if (allocated(this%sw_albedo)) then
173      deallocate(this%sw_albedo)
174    end if
175    if (allocated(this%sw_albedo_direct)) then
176      deallocate(this%sw_albedo_direct)
177    end if
178    if (allocated(this%lw_emissivity)) then
179      deallocate(this%lw_emissivity)
180    end if
181    if (allocated(this%lw_emission)) then
182      deallocate(this%lw_emission)
183    end if
184    if (allocated(this%spectral_solar_scaling)) then
185      deallocate(this%spectral_solar_scaling)
186    end if
187    if (allocated(this%iseed)) then
188      deallocate(this%iseed)
189    end if
190
191    if (lhook) call dr_hook('radiation_single_level:deallocate',1,hook_handle)
192
193  end subroutine deallocate_single_level
194
195
196  !---------------------------------------------------------------------
197  ! Unimaginative initialization of random-number seeds
198  subroutine init_seed_simple(this, istartcol, iendcol)
199    class(single_level_type), intent(inout) :: this
200    integer, intent(in)                     :: istartcol, iendcol
201
202    integer :: jcol
203
204    if (.not. allocated(this%iseed)) then
205      allocate(this%iseed(istartcol:iendcol))
206    end if
207
208    do jcol = istartcol,iendcol
209      this%iseed(jcol) = jcol
210    end do
211  end subroutine init_seed_simple
212
213
214  !---------------------------------------------------------------------
215  ! Extract the shortwave and longwave surface albedos in each g-point
216  subroutine get_albedos(this, istartcol, iendcol, config, &
217       &                 sw_albedo_direct, sw_albedo_diffuse, lw_albedo)
218
219    use radiation_config, only : config_type
220    use radiation_io,     only : nulerr, radiation_abort
221    use yomhook,          only : lhook, dr_hook, jphook
222
223    class(single_level_type), intent(in) :: this
224    type(config_type),        intent(in) :: config
225    integer,                  intent(in) :: istartcol, iendcol
226
227    ! The longwave albedo of the surface in each longwave g-point;
228    ! note that these are weighted averages of the values from
229    ! individual tiles
230    real(jprb), intent(out), optional &
231         &  :: lw_albedo(config%n_g_lw, istartcol:iendcol)
232
233    ! Direct and diffuse shortwave surface albedo in each shortwave
234    ! g-point; note that these are weighted averages of the values
235    ! from individual tiles
236    real(jprb), intent(out), dimension(config%n_g_sw, istartcol:iendcol) &
237         &  :: sw_albedo_direct, sw_albedo_diffuse
238
239    ! Temporary storage of albedo in ecRad bands
240    real(jprb) :: sw_albedo_band(istartcol:iendcol, config%n_bands_sw)
241    real(jprb) :: lw_albedo_band(istartcol:iendcol, config%n_bands_lw)
242
243    ! Number of albedo bands
244    integer :: nalbedoband
245
246    ! Loop indices for ecRad bands and albedo bands
247    integer :: jband, jalbedoband, jcol
248
249    real(jphook) :: hook_handle
250
251    if (lhook) call dr_hook('radiation_single_level:get_albedos',0,hook_handle)
252
253    if (config%do_sw) then
254      ! Albedos/emissivities are stored in single_level in their own
255      ! spectral intervals and with column as the first dimension
256      if (config%use_canopy_full_spectrum_sw) then
257        ! Albedos provided in each g point
258        if (size(this%sw_albedo,2) /= config%n_g_sw) then
259          write(nulerr,'(a,i0,a)') '*** Error: single_level%sw_albedo does not have the expected ', &
260               &  config%n_g_sw, ' spectral intervals'
261          call radiation_abort()
262        end if
263        sw_albedo_diffuse = transpose(this%sw_albedo(istartcol:iendcol,:))
264        if (allocated(this%sw_albedo_direct)) then
265          sw_albedo_direct = transpose(this%sw_albedo_direct(istartcol:iendcol,:))
266        end if
267      else if (.not. config%do_nearest_spectral_sw_albedo) then
268        ! Albedos averaged accurately to ecRad spectral bands
269        nalbedoband = size(config%sw_albedo_weights,1)
270        if (size(this%sw_albedo,2) /= nalbedoband) then
271          write(nulerr,'(a,i0,a)') '*** Error: single_level%sw_albedo does not have the expected ', &
272               &  nalbedoband, ' bands'
273          call radiation_abort()
274        end if
275
276        sw_albedo_band = 0.0_jprb
277        do jband = 1,config%n_bands_sw
278          do jalbedoband = 1,nalbedoband
279            if (config%sw_albedo_weights(jalbedoband,jband) /= 0.0_jprb) then
280              do jcol = istartcol,iendcol
281                sw_albedo_band(jcol,jband) &
282                    &  = sw_albedo_band(jcol,jband) &
283                    &  + config%sw_albedo_weights(jalbedoband,jband) &
284                    &    * this%sw_albedo(jcol, jalbedoband)
285              end do
286            end if
287          end do
288        end do
289
290        sw_albedo_diffuse = transpose(sw_albedo_band(istartcol:iendcol, &
291             &                              config%i_band_from_reordered_g_sw))
292        if (allocated(this%sw_albedo_direct)) then
293          sw_albedo_band = 0.0_jprb
294          do jband = 1,config%n_bands_sw
295            do jalbedoband = 1,nalbedoband
296              if (config%sw_albedo_weights(jalbedoband,jband) /= 0.0_jprb) then
297                sw_albedo_band(istartcol:iendcol,jband) &
298                     &  = sw_albedo_band(istartcol:iendcol,jband) &
299                     &  + config%sw_albedo_weights(jalbedoband,jband) &
300                     &    * this%sw_albedo_direct(istartcol:iendcol, jalbedoband)
301              end if
302            end do
303          end do
304          sw_albedo_direct = transpose(sw_albedo_band(istartcol:iendcol, &
305               &                             config%i_band_from_reordered_g_sw))
306        else
307          sw_albedo_direct = sw_albedo_diffuse
308        end if
309      else
310        ! Albedos mapped less accurately to ecRad spectral bands
311        if (maxval(config%i_albedo_from_band_sw) > size(this%sw_albedo,2)) then
312          write(nulerr,'(a,i0,a)') '*** Error: single_level%sw_albedo has fewer than required ', &
313               &  maxval(config%i_albedo_from_band_sw), ' bands'
314          call radiation_abort()
315        end if     
316        sw_albedo_diffuse = transpose(this%sw_albedo(istartcol:iendcol, &
317             &  config%i_albedo_from_band_sw(config%i_band_from_reordered_g_sw)))
318        if (allocated(this%sw_albedo_direct)) then
319          sw_albedo_direct = transpose(this%sw_albedo_direct(istartcol:iendcol, &
320               &  config%i_albedo_from_band_sw(config%i_band_from_reordered_g_sw)))
321        else
322          sw_albedo_direct = sw_albedo_diffuse
323        end if
324      end if
325    end if
326
327    if (config%do_lw .and. present(lw_albedo)) then
328      if (config%use_canopy_full_spectrum_lw) then
329        if (config%n_g_lw /= size(this%lw_emissivity,2)) then
330          write(nulerr,'(a,i0,a)') '*** Error: single_level%lw_emissivity does not have the expected ', &
331               &  config%n_g_lw, ' spectral intervals'
332          call radiation_abort()
333        end if
334        lw_albedo = 1.0_jprb - transpose(this%lw_emissivity(istartcol:iendcol,:))
335      else if (.not. config%do_nearest_spectral_lw_emiss) then
336        ! Albedos averaged accurately to ecRad spectral bands
337        nalbedoband = size(config%lw_emiss_weights,1)
338        if (nalbedoband /= size(this%lw_emissivity,2)) then
339          write(nulerr,'(a,i0,a)') '*** Error: single_level%lw_emissivity does not have the expected ', &
340               &  nalbedoband, ' bands'
341          call radiation_abort()
342        end if
343        lw_albedo_band = 0.0_jprb
344        do jband = 1,config%n_bands_lw
345          do jalbedoband = 1,nalbedoband
346            if (config%lw_emiss_weights(jalbedoband,jband) /= 0.0_jprb) then
347              do jcol = istartcol,iendcol
348                lw_albedo_band(jcol,jband) &
349                    &  = lw_albedo_band(jcol,jband) &
350                    &  + config%lw_emiss_weights(jalbedoband,jband) &
351                    &    * (1.0_jprb-this%lw_emissivity(jcol, jalbedoband))
352              end do
353            end if
354          end do
355        end do
356
357        lw_albedo = transpose(lw_albedo_band(istartcol:iendcol, &
358             &                config%i_band_from_reordered_g_lw))
359      else
360        if (maxval(config%i_emiss_from_band_lw) > size(this%lw_emissivity,2)) then
361          write(nulerr,'(a,i0,a)') '*** Error: single_level%lw_emissivity has fewer than required ', &
362               &  maxval(config%i_emiss_from_band_lw), ' bands'
363          call radiation_abort()
364        end if
365        lw_albedo = 1.0_jprb - transpose(this%lw_emissivity(istartcol:iendcol, &
366             &  config%i_emiss_from_band_lw(config%i_band_from_reordered_g_lw)))
367      end if
368    end if
369
370    if (lhook) call dr_hook('radiation_single_level:get_albedos',1,hook_handle)
371
372  end subroutine get_albedos
373
374
375  !---------------------------------------------------------------------
376  ! Return .true. if the contents of a single_level structure are out
377  ! of a physically sensible range, optionally only considering only
378  ! columns between istartcol and iendcol
379  function out_of_physical_bounds(this, istartcol, iendcol, do_fix) result(is_bad)
380
381    use yomhook,          only : lhook, dr_hook, jphook
382    use radiation_check,  only : out_of_bounds_1d, out_of_bounds_2d
383
384    class(single_level_type), intent(inout) :: this
385    integer,         optional,intent(in) :: istartcol, iendcol
386    logical,         optional,intent(in) :: do_fix
387    logical                              :: is_bad
388
389    logical    :: do_fix_local
390
391    real(jphook) :: hook_handle
392
393    if (lhook) call dr_hook('radiation_single_level:out_of_physical_bounds',0,hook_handle)
394
395    if (present(do_fix)) then
396      do_fix_local = do_fix
397    else
398      do_fix_local = .false.
399    end if
400
401    is_bad =    out_of_bounds_1d(this%cos_sza, "cos_sza", -1.0_jprb, 1.0_jprb, &
402         &                       do_fix_local, i1=istartcol, i2=iendcol) &
403         & .or. out_of_bounds_1d(this%skin_temperature, "skin_temperature", 173.0_jprb, 373.0_jprb, &
404         &                       do_fix_local, i1=istartcol, i2=iendcol) &
405         & .or. out_of_bounds_2d(this%sw_albedo, "sw_albedo", 0.0_jprb, 1.0_jprb, &
406         &                       do_fix_local, i1=istartcol, i2=iendcol) &
407         & .or. out_of_bounds_2d(this%sw_albedo_direct, "sw_albedo", 0.0_jprb, 1.0_jprb, &
408         &                       do_fix_local, i1=istartcol, i2=iendcol) &
409         & .or. out_of_bounds_2d(this%lw_emissivity, "lw_emissivity", 0.0_jprb, 1.0_jprb, &
410         &                       do_fix_local, i1=istartcol, i2=iendcol)
411
412    if (lhook) call dr_hook('radiation_single_level:out_of_physical_bounds',1,hook_handle)
413
414  end function out_of_physical_bounds
415
416end module radiation_single_level
Note: See TracBrowser for help on using the repository browser.