source: LMDZ6/trunk/libf/phylmd/ecrad/radiation/radiation_spectral_definition.F90 @ 4848

Last change on this file since 4848 was 4773, checked in by idelkadi, 11 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: 36.9 KB
Line 
1! radiation_spectral_definition.F90 - Derived type to describe a spectral definition
2!
3! (C) Copyright 2020- 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! License: see the COPYING file for details
15!
16
17module radiation_spectral_definition
18
19  use parkind1,    only : jprb
20
21  implicit none
22
23  public
24
25  real(jprb), parameter :: SolarReferenceTemperature = 5777.0_jprb ! K
26  real(jprb), parameter :: TerrestrialReferenceTemperature = 273.15_jprb ! K
27
28  !---------------------------------------------------------------------
29  ! A derived type describing the contribution of the g points of a
30  ! correlated k-distribution gas-optics model from each part of the
31  ! spectrum. This is used primarily to map the cloud and aerosol
32  ! optical properties on to the gas g points.
33  type spectral_definition_type
34   
35    ! Spectral mapping of g points
36
37    ! Number of wavenumber intervals
38    integer :: nwav = 0
39    ! Number of k terms / g points
40    integer :: ng   = 0
41    ! Start and end wavenumber (cm-1), dimensioned (nwav)
42    real(jprb), allocatable :: wavenumber1(:)
43    real(jprb), allocatable :: wavenumber2(:)
44    ! Fraction of each g point in each wavenumber interval,
45    ! dimensioned (nwav, ng)
46    real(jprb), allocatable :: gpoint_fraction(:,:)
47
48    ! Spectral weighting information for generating mappings to/from
49    ! different spectral grids: this can be in terms of a reference
50    ! temperature (K) to generate a Planck function, or the
51    ! solar_spectral_irradiance (W m-2) if available in the gas-optics
52    ! file.
53    real(jprb) :: reference_temperature = -1.0_jprb
54    real(jprb), allocatable :: solar_spectral_irradiance(:)
55   
56    ! Band information
57
58    ! Number of bands
59    integer :: nband = 0
60    ! Lower and upper bounds of wavenumber bands (cm-1), dimensioned
61    ! (nband)
62    real(jprb), allocatable :: wavenumber1_band(:)
63    real(jprb), allocatable :: wavenumber2_band(:)
64    ! Band (one based) to which each g point belongs
65    integer,    allocatable :: i_band_number(:)
66
67  contains
68    procedure :: read => read_spectral_definition
69    procedure :: allocate_bands_only
70    procedure :: deallocate
71    procedure :: find => find_wavenumber
72    procedure :: calc_mapping
73    procedure :: calc_mapping_from_bands
74    procedure :: calc_mapping_from_wavenumber_bands
75    procedure :: print_mapping_from_bands
76    procedure :: min_wavenumber, max_wavenumber
77
78  end type spectral_definition_type
79
80contains
81
82  !---------------------------------------------------------------------
83  ! Read the description of a spectral definition from a NetCDF
84  ! file of the type used to describe an ecCKD model
85  subroutine read_spectral_definition(this, file)
86
87    use easy_netcdf, only : netcdf_file
88    use yomhook,     only : lhook, dr_hook, jphook
89
90    class(spectral_definition_type), intent(inout) :: this
91    type(netcdf_file),               intent(inout) :: file
92
93    real(jphook) :: hook_handle
94
95    if (lhook) call dr_hook('radiation_spectral_definition:read',0,hook_handle)
96
97    ! Read spectral mapping of g points
98    call file%get('wavenumber1', this%wavenumber1)
99    call file%get('wavenumber2', this%wavenumber2)
100    call file%get('gpoint_fraction', this%gpoint_fraction)
101
102    ! Read band information
103    call file%get('wavenumber1_band', this%wavenumber1_band)
104    call file%get('wavenumber2_band', this%wavenumber2_band)
105    call file%get('band_number', this%i_band_number)
106
107    ! Read spectral weighting information
108    if (file%exists('solar_spectral_irradiance')) then
109      ! This is on the same grid as wavenumber1,2
110      call file%get('solar_spectral_irradiance', &
111           &        this%solar_spectral_irradiance)
112    end if
113    if (file%exists('solar_irradiance')) then
114      ! Shortwave default temperature
115      this%reference_temperature = SolarReferenceTemperature
116    else
117      ! Longwave reference temperature
118      this%reference_temperature = TerrestrialReferenceTemperature
119    end if
120   
121    ! Band number is 0-based: add 1
122    this%i_band_number = this%i_band_number + 1
123
124    this%nwav  = size(this%wavenumber1)
125    this%ng    = size(this%gpoint_fraction, 2);
126    this%nband = size(this%wavenumber1_band)
127
128    if (lhook) call dr_hook('radiation_spectral_definition:read',1,hook_handle)
129
130  end subroutine read_spectral_definition
131
132
133  !---------------------------------------------------------------------
134  ! Store a simple band description by copying over the lower and
135  ! upper wavenumbers of each band
136  subroutine allocate_bands_only(this, wavenumber1, wavenumber2)
137
138    use yomhook,     only : lhook, dr_hook, jphook
139
140    class(spectral_definition_type), intent(inout) :: this
141    real(jprb),        dimension(:), intent(in)    :: wavenumber1, wavenumber2
142
143    real(jphook) :: hook_handle
144
145    if (lhook) call dr_hook('radiation_spectral_definition:allocate_bands_only',0,hook_handle)
146
147    call this%deallocate()
148
149    this%nband = size(wavenumber1)
150    allocate(this%wavenumber1_band(this%nband))
151    allocate(this%wavenumber2_band(this%nband))
152    this%wavenumber1_band = wavenumber1
153    this%wavenumber2_band = wavenumber2
154
155    if (lhook) call dr_hook('radiation_spectral_definition:allocate_bands_only',1,hook_handle)
156
157  end subroutine allocate_bands_only
158
159
160  !---------------------------------------------------------------------
161  ! Deallocate memory inside a spectral definition object
162  subroutine deallocate(this)
163
164    class(spectral_definition_type), intent(inout) :: this
165   
166    this%nwav  = 0
167    this%ng    = 0
168    this%nband = 0
169
170    if (allocated(this%wavenumber1))      deallocate(this%wavenumber1)
171    if (allocated(this%wavenumber2))      deallocate(this%wavenumber2)
172    if (allocated(this%wavenumber1_band)) deallocate(this%wavenumber1_band)
173    if (allocated(this%wavenumber2_band)) deallocate(this%wavenumber2_band)
174    if (allocated(this%gpoint_fraction))  deallocate(this%gpoint_fraction)
175    if (allocated(this%i_band_number))    deallocate(this%i_band_number)
176
177  end subroutine deallocate
178
179
180  !---------------------------------------------------------------------
181  ! Find the index to the highest wavenumber in the spectral
182  ! definition that is lower than or equal to "wavenumber", used for
183  ! implementing look-up tables
184  pure function find_wavenumber(this, wavenumber)
185    class(spectral_definition_type), intent(in) :: this
186    real(jprb),                      intent(in) :: wavenumber ! cm-1
187    integer                                     :: find_wavenumber
188
189    if (wavenumber < this%wavenumber1(1) .or. wavenumber > this%wavenumber2(this%nwav)) then
190      ! Wavenumber not present
191      find_wavenumber = 0
192    else
193      find_wavenumber = 1
194      do while (wavenumber > this%wavenumber2(find_wavenumber) &
195           &    .and. find_wavenumber < this%nwav)
196        find_wavenumber = find_wavenumber + 1
197      end do
198    end if
199  end function find_wavenumber
200
201
202  !---------------------------------------------------------------------
203  ! Compute a mapping matrix "mapping" that can be used in an
204  ! expression y=matmul(mapping,x) where x is a variable containing
205  ! optical properties at each input "wavenumber", and y is this
206  ! variable mapped on to the spectral intervals in the spectral
207  ! definition "this".
208  subroutine calc_mapping(this, wavenumber, mapping, weighting_temperature, use_bands)
209
210    use yomhook,      only : lhook, dr_hook, jphook
211    use radiation_io, only : nulerr, radiation_abort
212
213    class(spectral_definition_type), intent(in)    :: this
214    real(jprb),                      intent(in)    :: wavenumber(:) ! cm-1
215    real(jprb), allocatable,         intent(inout) :: mapping(:,:)
216    real(jprb), optional,            intent(in)    :: weighting_temperature ! K
217    logical,    optional,            intent(in)    :: use_bands
218
219    ! Spectral weights to apply, same length as wavenumber above
220    real(jprb), dimension(:), allocatable :: weight, planck_weight
221
222    ! Wavenumbers (cm-1) marking triangle of influence of a cloud
223    ! spectral point
224    real(jprb) :: wavenum0, wavenum1, wavenum2
225
226    integer    :: nwav ! Number of wavenumbers describing cloud
227
228    ! Indices to wavenumber intervals in spectral definition structure
229    integer    :: isd, isd0, isd1, isd2
230
231    ! Wavenumber index
232    integer    :: iwav
233   
234    ! Loop indices
235    integer    :: jg, jwav, jband
236
237    logical    :: use_bands_local
238
239    real(jphook) :: hook_handle
240
241    if (lhook) call dr_hook('radiation_spectral_definition:calc_mapping',0,hook_handle)
242
243    if (present(use_bands)) then
244      use_bands_local = use_bands
245    else
246      use_bands_local = .false.
247    end if
248
249    nwav = size(wavenumber)
250
251    if (allocated(mapping)) then
252      deallocate(mapping)
253    end if
254   
255    ! Define the mapping matrix
256    if (use_bands_local) then
257      ! Cloud properties per band
258
259      allocate(mapping(this%nband, nwav))
260      allocate(weight(nwav))
261
262      ! Planck weight uses the wavenumbers of the cloud points
263      allocate(planck_weight(nwav))
264      if (present(weighting_temperature)) then
265        if (weighting_temperature > 0.0_jprb) then
266          planck_weight = calc_planck_function_wavenumber(wavenumber, &
267               &                          weighting_temperature)
268        else
269          ! Legacy mode: unweighted average
270          planck_weight = 1.0_jprb
271        end if
272      else
273        planck_weight = calc_planck_function_wavenumber(wavenumber, &
274             &                          this%reference_temperature)
275      end if
276
277      do jband = 1,this%nband
278        weight = 0.0_jprb
279        do jwav = 1,nwav
280          ! Work out wavenumber range for which this cloud wavenumber
281          ! will be applicable
282          if (wavenumber(jwav) >= this%wavenumber1_band(jband) &
283               & .and. wavenumber(jwav) <= this%wavenumber2_band(jband)) then
284            if (jwav > 1) then
285              wavenum1 = max(this%wavenumber1_band(jband), &
286                   &  0.5_jprb*(wavenumber(jwav-1)+wavenumber(jwav)))
287            else
288              wavenum1 = this%wavenumber1_band(jband)
289            end if
290            if (jwav < nwav) then
291              wavenum2 = min(this%wavenumber2_band(jband), &
292                   &  0.5_jprb*(wavenumber(jwav)+wavenumber(jwav+1)))
293            else
294              wavenum2 = this%wavenumber2_band(jband)
295            end if
296            ! This cloud wavenumber is weighted by the wavenumber
297            ! range of its applicability multiplied by the Planck
298            ! function at an appropriate temperature
299            weight(jwav) = (wavenum2-wavenum1) * planck_weight(jwav)
300          end if
301        end do
302        if (sum(weight) <= 0.0_jprb) then
303          ! No cloud wavenumbers lie in the band; interpolate to
304          ! central wavenumber of band instead
305          if (wavenumber(1) >= this%wavenumber2_band(jband)) then
306            ! Band is entirely below first cloudy wavenumber
307            weight(1) = 1.0_jprb
308          else if (wavenumber(nwav) <= this%wavenumber1_band(jband)) then
309            ! Band is entirely above last cloudy wavenumber
310            weight(nwav) = 1.0_jprb
311          else
312            ! Find interpolating points
313            iwav = 2
314            do while (wavenumber(iwav) < this%wavenumber2_band(jband))
315              iwav = iwav+1
316            end do
317            weight(iwav-1) = planck_weight(iwav-1) * (wavenumber(iwav) &
318                 &  - 0.5_jprb*(this%wavenumber2_band(jband)+this%wavenumber1_band(jband)))
319            weight(iwav) = planck_weight(iwav) * (-wavenumber(iwav-1) &
320                 &  + 0.5_jprb*(this%wavenumber2_band(jband)+this%wavenumber1_band(jband)))
321          end if
322        end if
323        mapping(jband,:) = weight / sum(weight)
324      end do
325
326      deallocate(weight)
327      deallocate(planck_weight)
328
329    else
330      ! Cloud properties per g-point
331
332      if (this%ng == 0) then
333        write(nulerr,'(a)') '*** Error: requested cloud/aerosol mapping per g-point but only available per band'
334        call radiation_abort('Radiation configuration error')
335      end if
336
337      allocate(mapping(this%ng, nwav))
338      allocate(weight(this%nwav))
339      allocate(planck_weight(this%nwav))
340
341      if (allocated(this%solar_spectral_irradiance)) then
342        planck_weight = this%solar_spectral_irradiance
343      else
344        planck_weight = calc_planck_function_wavenumber(0.5_jprb &
345             &  * (this%wavenumber1 + this%wavenumber2), &
346             &  this%reference_temperature)
347      end if
348
349      mapping = 0.0_jprb
350      ! Loop over wavenumbers representing cloud
351      do jwav = 1,nwav
352        ! Clear the weights. The weight says for one wavenumber in the
353        ! cloud file, what is its fractional contribution to each of
354        ! the spectral-definition intervals
355        weight = 0.0_jprb
356
357        ! Cloud properties are linearly interpolated between each of
358        ! the nwav cloud points; therefore, the influence of a
359        ! particular cloud point extends as a triangle between
360        ! wavenum0 and wavenum2, peaking at wavenum1
361        wavenum1 = wavenumber(jwav)
362        isd1 = this%find(wavenum1)
363        if (isd1 < 1) then
364          cycle
365        end if
366        if (jwav > 1) then
367          wavenum0 = wavenumber(jwav-1)
368
369          ! Map triangle under (wavenum0,0) to (wavenum1,1) to the
370          ! wavenumbers in this%gpoint_fraction
371          isd0 = this%find(wavenum0)
372          if (isd0 == isd1) then
373            ! Triangle completely within the range
374            ! this%wavenumber1(isd0)-this%wavenumber2(isd0)
375            weight(isd0) = 0.5_jprb*(wavenum1-wavenum0) &
376                 &       / (this%wavenumber2(isd0)-this%wavenumber1(isd0))
377          else
378            if (isd0 >= 1) then
379              ! Left part of triangle
380              weight(isd0) = 0.5_jprb * (this%wavenumber2(isd0)-wavenum0)**2 &
381                   &       / ((this%wavenumber2(isd0)-this%wavenumber1(isd0)) &
382                   &         *(wavenum1-wavenum0))
383            end if
384            ! Right part of triangle (trapezium)
385!            weight(isd1) = 0.5_jprb * (wavenum1-this%wavenumber1(isd1)) &
386!                 &       * (wavenum1 + this%wavenumber1(isd1) - 2.0_jprb*wavenum0) &
387!                 &       / (wavenum1-wavenum0)
388            weight(isd1) = 0.5_jprb * (1.0_jprb &
389                 &  + (this%wavenumber1(isd1)-wavenum1)/(wavenum1-wavenum0)) &
390                 &  * (wavenum1-this%wavenumber1(isd1)) &
391                 &  / (this%wavenumber2(isd1)-this%wavenumber1(isd1))
392            if (isd1-isd0 > 1) then
393              do isd = isd0+1,isd1-1
394                ! Intermediate trapezia
395                weight(isd) = 0.5_jprb * (this%wavenumber1(isd)+this%wavenumber2(isd) &
396                     &                    - 2.0_jprb*wavenum0) &
397                     &      / (wavenum1-wavenum0)
398              end do
399            end if
400          end if
401
402        else
403          ! First cloud wavenumber: all wavenumbers in the spectral
404          ! definition below this will use the first one
405          if (isd1 >= 1) then
406            weight(1:isd1-1) = 1.0_jprb
407            weight(isd1) = (wavenum1-this%wavenumber1(isd1)) &
408                 &       / (this%wavenumber2(isd1)-this%wavenumber1(isd1))
409          end if
410        end if
411
412        if (jwav < nwav) then
413          wavenum2 = wavenumber(jwav+1)
414
415          ! Map triangle under (wavenum1,1) to (wavenum2,0) to the
416          ! wavenumbers in this%gpoint_fraction
417          isd2 = this%find(wavenum2)
418
419          if (isd1 == isd2) then
420            ! Triangle completely within the range
421            ! this%wavenumber1(isd1)-this%wavenumber2(isd1)
422            weight(isd1) = weight(isd1) + 0.5_jprb*(wavenum2-wavenum1) &
423                 &       / (this%wavenumber2(isd1)-this%wavenumber1(isd1))
424          else
425            if (isd2 >= 1 .and. isd2 <= this%nwav) then
426              ! Right part of triangle
427              weight(isd2) = weight(isd2) + 0.5_jprb * (wavenum2-this%wavenumber1(isd2))**2 &
428                   &       / ((this%wavenumber2(isd2)-this%wavenumber1(isd2)) &
429                   &         *(wavenum2-wavenum1))
430            end if
431            ! Left part of triangle (trapezium)
432!            weight(isd1) = weight(isd1) + 0.5_jprb * (this%wavenumber2(isd1)-wavenum1) &
433!                 &       * (wavenum1 + this%wavenumber2(isd1) - 2.0_jprb*wavenum2) &
434!                 &       / (wavenum2-wavenum1)
435            weight(isd1) = weight(isd1) + 0.5_jprb * (1.0_jprb &
436                 &  + (wavenum2-this%wavenumber2(isd1)) / (wavenum2-wavenum1)) &
437                 &  * (this%wavenumber2(isd1)-wavenum1) &
438                 &  / (this%wavenumber2(isd1)-this%wavenumber1(isd1))
439            if (isd2-isd1 > 1) then
440              do isd = isd1+1,isd2-1
441                ! Intermediate trapezia
442                weight(isd) = weight(isd) + 0.5_jprb * (2.0_jprb*wavenum2 &
443                     & - this%wavenumber1(isd) - this%wavenumber2(isd)) &
444                     &      / (wavenum2-wavenum1)
445              end do
446            end if
447          end if
448
449        else
450          ! Last cloud wavenumber: all wavenumbers in the spectral
451          ! definition above this will use the last one
452          if (isd1 <= this%nwav) then
453            weight(isd1+1:this%nwav) = 1.0_jprb
454            weight(isd1) = (this%wavenumber2(isd1)-wavenum1) &
455                 &       / (this%wavenumber2(isd1)-this%wavenumber1(isd1))
456          end if
457        end if
458
459        weight = weight * planck_weight
460
461        do jg = 1,this%ng
462          mapping(jg, jwav) = sum(weight * this%gpoint_fraction(:,jg))
463        end do
464
465      end do
466
467      deallocate(weight)
468      deallocate(planck_weight)
469
470      ! Normalize mapping matrix
471      do jg = 1,this%ng
472        mapping(jg,:) = mapping(jg,:) * (1.0_jprb/sum(mapping(jg,:)))
473      end do
474
475    end if
476
477    if (lhook) call dr_hook('radiation_spectral_definition:calc_mapping',1,hook_handle)
478
479  end subroutine calc_mapping
480
481
482  !---------------------------------------------------------------------
483  ! Under normal operation (if use_fluxes is .false. or not present),
484  ! compute a mapping matrix "mapping" that can be used in an
485  ! expression y=matmul(mapping^T,x) where x is a variable containing
486  ! optical properties in input bands (e.g. albedo in shortwave albedo
487  ! bands), and y is this variable mapped on to the spectral intervals
488  ! in the spectral definition "this". Note that "mapping" is here
489  ! transposed from the convention in the calc_mapping routine.  Under
490  ! the alternative operation (if use_fluxes is present and .true.),
491  ! the mapping works in the reverse sense: if y contains fluxes in
492  ! each ecRad band or g-point, then x=matmul(mapping,y) would return
493  ! fluxes in x averaged to user-supplied "input" bands. In this
494  ! version, the bands are described by their wavelength bounds
495  ! (wavelength_bound, which must be increasing and exclude the end
496  ! points) and the index of the mapping matrix that each band
497  ! corresponds to (i_intervals, which has one more element than
498  ! wavelength_bound and can have duplicated values if an
499  ! albedo/emissivity value is to be associated with more than one
500  ! discontinuous ranges of the spectrum).
501  subroutine calc_mapping_from_bands(this, &
502       &  wavelength_bound, i_intervals, mapping, use_bands, use_fluxes)
503
504    use yomhook,      only : lhook, dr_hook, jphook
505    use radiation_io, only : nulerr, radiation_abort
506
507    class(spectral_definition_type), intent(in)    :: this
508    ! Monotonically increasing wavelength bounds (m) between
509    ! intervals, not including the outer bounds (which are assumed to
510    ! be zero and infinity)
511    real(jprb),                      intent(in)    :: wavelength_bound(:)
512    ! The albedo band indices corresponding to each interval
513    integer,                         intent(in)    :: i_intervals(:)
514    real(jprb), allocatable,         intent(inout) :: mapping(:,:)
515    logical,    optional,            intent(in)    :: use_bands
516    logical,    optional,            intent(in)    :: use_fluxes
517
518    ! Planck function and central wavenumber of each wavenumber
519    ! interval of the spectral definition
520    real(jprb) :: planck(this%nwav)         ! W m-2 (cm-1)-1
521    real(jprb) :: wavenumber_mid(this%nwav) ! cm-1
522
523    real(jprb), allocatable :: mapping_denom(:,:)
524
525    real(jprb) :: wavenumber1_bound, wavenumber2_bound
526
527    ! To work out weights we sample the Planck function at five points
528    ! in the interception between an input interval and a band, and
529    ! use the Trapezium Rule
530    integer, parameter :: nsample = 5
531    integer :: isamp
532    real(jprb), dimension(nsample) :: wavenumber_sample, planck_sample
533    real(jprb), parameter :: weight_sample(nsample) &
534         &        = [0.5_jprb, 1.0_jprb, 1.0_jprb, 1.0_jprb, 0.5_jprb]
535
536    ! Index of input value corresponding to each wavenumber interval
537    integer :: i_input(this%nwav)
538
539    ! Number of albedo/emissivity values that will be provided, some
540    ! of which may span discontinuous intervals in wavelength space
541    integer :: ninput
542
543    ! Number of albedo/emissivity intervals represented, where some
544    ! may be grouped to have the same value of albedo/emissivity (an
545    ! example is in the thermal infrared where classically the IFS has
546    ! ninput=2 and ninterval=3, since only two emissivities are
547    ! provided representing (1) the infrared window, and (2) the
548    ! intervals to each side of the infrared window.
549    integer :: ninterval
550
551    logical    :: use_bands_local, use_fluxes_local
552
553    ! Loop indices
554    integer    :: jg, jband, jin, jint, jwav
555
556    real(jphook) :: hook_handle
557
558    if (lhook) call dr_hook('radiation_spectral_definition:calc_mapping_from_bands',0,hook_handle)
559
560    if (present(use_bands)) then
561      use_bands_local = use_bands
562    else
563      use_bands_local = .false.
564    end if
565
566    if (present(use_fluxes)) then
567      use_fluxes_local = use_fluxes
568    else
569      use_fluxes_local = .false.
570    end if
571
572    ! Count the number of input intervals
573    ninterval = size(i_intervals)
574    ninput    = maxval(i_intervals)
575   
576    if (allocated(mapping)) then
577      deallocate(mapping)
578    end if
579   
580    ! Check wavelength is monotonically increasing
581    if (ninterval > 2) then
582      do jint = 2,ninterval-1
583        if (wavelength_bound(jint) <= wavelength_bound(jint-1)) then
584          write(nulerr, '(a)') '*** Error: wavelength bounds must be monotonically increasing'
585          call radiation_abort()
586        end if
587      end do
588    end if
589
590    ! Define the mapping matrix
591    if (use_bands_local) then
592      ! Require properties per band
593
594      allocate(mapping(ninput, this%nband))
595      mapping = 0.0_jprb
596
597      if (use_fluxes_local) then
598        allocate(mapping_denom(ninput, this%nband))
599        mapping_denom = 0.0_jprb
600      end if
601
602      do jband = 1,this%nband
603        do jint = 1,ninterval
604          if (jint == 1) then
605            ! First input interval in wavelength space: lower
606            ! wavelength bound is 0 m, so infinity cm-1
607            wavenumber2_bound = this%wavenumber2_band(jband)
608          else
609            wavenumber2_bound = min(this%wavenumber2_band(jband), &
610                 &                  0.01_jprb/wavelength_bound(jint-1))
611          end if
612
613          if (jint == ninterval) then
614            ! Final input interval in wavelength space: upper
615            ! wavelength bound is infinity m, so 0 cm-1
616            wavenumber1_bound = this%wavenumber1_band(jband)
617          else
618            wavenumber1_bound = max(this%wavenumber1_band(jband), &
619                 &                  0.01_jprb/wavelength_bound(jint))
620
621          end if
622
623          if (wavenumber2_bound > wavenumber1_bound) then
624            ! Current input interval contributes to current band;
625            ! compute the weight of the contribution in proportion to
626            ! an approximate calculation of the integral of the Planck
627            ! function over the relevant part of the spectrum
628            wavenumber_sample = wavenumber1_bound + [(isamp,isamp=0,nsample-1)] &
629                 &  * (wavenumber2_bound-wavenumber1_bound) / real(nsample-1,jprb)
630            planck_sample = calc_planck_function_wavenumber(wavenumber_sample, &
631                 &                                 this%reference_temperature)
632            mapping(i_intervals(jint),jband) = mapping(i_intervals(jint),jband) &
633                 &  + sum(planck_sample*weight_sample) * (wavenumber2_bound-wavenumber1_bound)
634            if (use_fluxes_local) then
635              ! Compute an equivalent sample containing the entire ecRad band
636              wavenumber_sample = this%wavenumber1_band(jband) + [(isamp,isamp=0,nsample-1)] &
637                   &  * (this%wavenumber2_band(jband)-this%wavenumber1_band(jband)) &
638                   &  / real(nsample-1,jprb)
639              planck_sample = calc_planck_function_wavenumber(wavenumber_sample, &
640                   &                                 this%reference_temperature)
641              mapping_denom(i_intervals(jint),jband) = mapping_denom(i_intervals(jint),jband) &
642                 &  + sum(planck_sample*weight_sample) * (this%wavenumber2_band(jband)-this%wavenumber1_band(jband))
643            end if
644          end if
645
646        end do
647      end do
648
649      if (use_fluxes_local) then
650        mapping = mapping / max(1.0e-12_jprb, mapping_denom)
651        deallocate(mapping_denom)
652      end if
653
654    else
655      ! Require properties per g-point
656
657      if (this%ng == 0) then
658        write(nulerr,'(a)') '*** Error: requested surface mapping per g-point but only available per band'
659        call radiation_abort('Radiation configuration error')
660      end if
661
662      allocate(mapping(ninput,this%ng))
663      mapping = 0.0_jprb
664
665      wavenumber_mid = 0.5_jprb * (this%wavenumber1 + this%wavenumber2)
666      if (allocated(this%solar_spectral_irradiance)) then
667        planck = this%solar_spectral_irradiance
668      else
669        planck = calc_planck_function_wavenumber(wavenumber_mid, &
670             &                       this%reference_temperature)
671      end if
672
673#ifdef USE_COARSE_MAPPING
674      ! In the processing that follows, we assume that the wavenumber
675      ! grid on which the g-points are defined in the spectral
676      ! definition is much finer than the albedo/emissivity intervals
677      ! that the user will provide.  This means that each wavenumber
678      ! is assigned to only one of the albedo/emissivity intervals.
679
680      ! By default set all wavenumbers to use first input
681      ! albedo/emissivity
682      i_input = 1
683     
684      ! All bounded intervals
685      do jint = 2,ninterval-1
686        wavenumber1_bound = 0.01_jprb / wavelength_bound(jint)
687        wavenumber2_bound = 0.01_jprb / wavelength_bound(jint-1)
688        where (wavenumber_mid > wavenumber1_bound &
689             & .and. wavenumber_mid <= wavenumber2_bound)
690          i_input = i_intervals(jint)
691        end where
692      end do
693
694      ! Final interval in wavelength space goes up to wavelength of
695      ! infinity (wavenumber of zero)
696      if (ninterval > 1) then
697        wavenumber2_bound = 0.01_jprb / wavelength_bound(ninterval-1)
698        where (wavenumber_mid <= wavenumber2_bound)
699          i_input = i_intervals(ninterval)
700        end where
701      end if
702
703      do jg = 1,this%ng
704        do jin = 1,ninput
705          mapping(jin,jg) = sum(this%gpoint_fraction(:,jg) * planck, &
706               &                 mask=(i_input==jin))
707          if (use_fluxes_local) then
708            mapping(jin,jg) = mapping(jin,jg) / sum(this%gpoint_fraction(:,jg) * planck)
709          end if
710        end do
711      end do
712
713#else
714
715      ! Loop through all intervals
716      do jint = 1,ninterval
717        ! Loop through the wavenumbers for gpoint_fraction
718        do jwav = 1,this%nwav
719          if (jint == 1) then
720            ! First input interval in wavelength space: lower
721            ! wavelength bound is 0 m, so infinity cm-1
722            wavenumber2_bound = this%wavenumber2(jwav)
723          else
724            wavenumber2_bound = min(this%wavenumber2(jwav), &
725                 &                  0.01_jprb/wavelength_bound(jint-1))
726          end if
727
728          if (jint == ninterval) then
729            ! Final input interval in wavelength space: upper
730            ! wavelength bound is infinity m, so 0 cm-1
731            wavenumber1_bound = this%wavenumber1(jwav)
732          else
733            wavenumber1_bound = max(this%wavenumber1(jwav), &
734                 &                  0.01_jprb/wavelength_bound(jint))
735
736          end if
737
738          if (wavenumber2_bound > wavenumber1_bound) then
739            ! Overlap between input interval and gpoint_fraction
740            ! interval: compute the weight of the contribution in
741            ! proportion to an approximate calculation of the integral
742            ! of the Planck function over the relevant part of the
743            ! spectrum
744            mapping(i_intervals(jint),:) = mapping(i_intervals(jint),:) + this%gpoint_fraction(jwav,:) &
745                 &  * (planck(jwav) * (wavenumber2_bound - wavenumber1_bound) &
746                 &                  / (this%wavenumber2(jwav)-this%wavenumber1(jwav)))
747          end if
748        end do
749      end do
750      if (use_fluxes_local) then
751        do jg = 1,this%ng
752          mapping(:,jg) = mapping(:,jg) / sum(this%gpoint_fraction(:,jg) * planck)
753        end do
754      end if
755
756#endif
757     
758    end if
759
760    if (.not. use_fluxes_local) then
761      ! Normalize mapping matrix
762      do jg = 1,size(mapping,dim=2)
763        mapping(:,jg) = mapping(:,jg) * (1.0_jprb/sum(mapping(:,jg)))
764      end do
765    end if
766
767    if (lhook) call dr_hook('radiation_spectral_definition:calc_mapping_from_bands',1,hook_handle)
768
769  end subroutine calc_mapping_from_bands
770
771
772  !---------------------------------------------------------------------
773  ! As calc_mapping_from_bands but in terms of wavenumber bounds from
774  ! wavenumber1 to wavenumber2
775  subroutine calc_mapping_from_wavenumber_bands(this, &
776       &  wavenumber1, wavenumber2, mapping, use_bands, use_fluxes)
777
778    use yomhook,      only : lhook, dr_hook, jphook
779
780    class(spectral_definition_type), intent(in)    :: this
781    real(jprb), intent(in)    :: wavenumber1(:), wavenumber2(:)
782    real(jprb), allocatable,         intent(inout) :: mapping(:,:)
783    logical,    optional,            intent(in)    :: use_bands
784    logical,    optional,            intent(in)    :: use_fluxes
785
786    ! Monotonically increasing wavelength bounds (m) between
787    ! intervals, not including the outer bounds (which are assumed to
788    ! be zero and infinity)
789    real(jprb) :: wavelength_bound(size(wavenumber1)-1)
790    ! The albedo band indices corresponding to each interval
791    integer    :: i_intervals(size(wavenumber1))
792
793    ! Lower wavelength bound (m) of each band
794    real(jprb) :: wavelength1(size(wavenumber1))
795
796    logical    :: is_band_unassigned(size(wavenumber1))
797
798    ! Number of albedo/emissivity intervals represented, where some
799    ! may be grouped to have the same value of albedo/emissivity (an
800    ! example is in the thermal infrared where classically the IFS has
801    ! ninput=2 and ninterval=3, since only two emissivities are
802    ! provided representing (1) the infrared window, and (2) the
803    ! intervals to each side of the infrared window.
804    integer :: ninterval
805
806    ! Index to next band in order of increasing wavelength
807    integer :: inext
808
809    ! Loop indices
810    integer :: jint
811
812    real(jphook) :: hook_handle
813
814    if (lhook) call dr_hook('radiation_spectral_definition:calc_mapping_from_wavenumber_bands',0,hook_handle)
815
816    wavelength1 = 0.01_jprb / wavenumber2
817    ninterval = size(wavelength1)
818   
819    is_band_unassigned = .true.
820
821    do jint = 1,ninterval
822      inext = minloc(wavelength1, dim=1, mask=is_band_unassigned)
823      if (jint > 1) then
824        wavelength_bound(jint-1) = wavelength1(inext)
825      end if
826      is_band_unassigned(inext) = .false.
827      i_intervals(jint) = inext
828    end do
829
830    call calc_mapping_from_bands(this, wavelength_bound, i_intervals, mapping, use_bands, use_fluxes)
831
832    if (lhook) call dr_hook('radiation_spectral_definition:calc_mapping_from_wavenumber_bands',1,hook_handle)
833
834  end subroutine calc_mapping_from_wavenumber_bands
835
836
837  !---------------------------------------------------------------------
838  ! Print out the mapping computed by calc_mapping_from_bands
839  subroutine print_mapping_from_bands(this, mapping, use_bands)
840
841    use radiation_io, only : nulout
842
843    class(spectral_definition_type), intent(in) :: this
844    real(jprb), allocatable,         intent(in) :: mapping(:,:) ! (ninput,nband/ng)
845    logical,    optional,            intent(in) :: use_bands
846
847    logical :: use_bands_local
848
849    integer :: nin, nout
850    integer :: jin, jout
851
852    if (present(use_bands)) then
853      use_bands_local = use_bands
854    else
855      use_bands_local = .false.
856    end if
857
858    nin = size(mapping,1)
859    nout = size(mapping,2)
860
861    if (nin <= 1) then
862      write(nulout, '(a)') '  All spectral intervals will use the same albedo/emissivity'
863    else if (use_bands_local) then
864      write(nulout, '(a,i0,a,i0,a)') '  Mapping from ', nin, ' values to ', nout, ' bands (wavenumber ranges in cm-1)'
865      if (nout <= 40) then
866        do jout = 1,nout
867          write(nulout,'(i6,a,i6,a)',advance='no') nint(this%wavenumber1_band(jout)), ' to', &
868               &                        nint(this%wavenumber2_band(jout)), ':'
869          do jin = 1,nin
870            write(nulout,'(f5.2)',advance='no') mapping(jin,jout)
871          end do
872          write(nulout, '()')
873        end do
874      else
875        do jout = 1,30
876          write(nulout,'(i6,a,i6,a)',advance='no') nint(this%wavenumber1_band(jout)), ' to', &
877               &                        nint(this%wavenumber2_band(jout)), ':'
878          do jin = 1,nin
879            write(nulout,'(f5.2)',advance='no') mapping(jin,jout)
880          end do
881          write(nulout, '()')
882        end do
883        write(nulout,'(a)') '  ...'
884        write(nulout,'(i6,a,i6,a)',advance='no') nint(this%wavenumber1_band(nout)), ' to', &
885             &                        nint(this%wavenumber2_band(nout)), ':'
886        do jin = 1,nin
887          write(nulout,'(f5.2)',advance='no') mapping(jin,nout)
888        end do
889        write(nulout, '()')
890      end if
891    else
892      write(nulout, '(a,i0,a,i0,a)') '  Mapping from ', nin, ' values to ', nout, ' g-points'
893      if (nout <= 40) then
894        do jout = 1,nout
895          write(nulout,'(i3,a)',advance='no') jout, ':'
896          do jin = 1,nin
897            write(nulout,'(f5.2)',advance='no') mapping(jin,jout)
898          end do
899          write(nulout, '()')
900        end do
901      else
902        do jout = 1,30
903          write(nulout,'(i3,a)',advance='no') jout, ':'
904          do jin = 1,nin
905            write(nulout,'(f5.2)',advance='no') mapping(jin,jout)
906          end do
907          write(nulout, '()')
908        end do
909        write(nulout,'(a)') '  ...'
910        write(nulout,'(i3,a)',advance='no') nout, ':'
911        do jin = 1,nin
912          write(nulout,'(f5.2)',advance='no') mapping(jin,nout)
913        end do
914        write(nulout, '()')
915      end if
916    end if
917
918  end subroutine print_mapping_from_bands
919
920
921  !---------------------------------------------------------------------
922  ! Return the minimum wavenumber of this object in cm-1
923  pure function min_wavenumber(this)
924
925    class(spectral_definition_type), intent(in)    :: this
926    real(jprb) :: min_wavenumber
927
928    if (this%nwav > 0) then
929      min_wavenumber = this%wavenumber1(1)
930    else
931      min_wavenumber = minval(this%wavenumber1_band)
932    end if
933
934  end function min_wavenumber
935
936
937  !---------------------------------------------------------------------
938  ! Return the maximum wavenumber of this object in cm-1
939  pure function max_wavenumber(this)
940
941    class(spectral_definition_type), intent(in)    :: this
942    real(jprb) :: max_wavenumber
943
944    if (this%nwav > 0) then
945      max_wavenumber = this%wavenumber1(this%nwav)
946    else
947      max_wavenumber = maxval(this%wavenumber2_band)
948    end if
949
950  end function max_wavenumber
951
952
953  !---------------------------------------------------------------------
954  ! Return the Planck function (in W m-2 (cm-1)-1) for a given
955  ! wavenumber (cm-1) and temperature (K), ensuring double precision
956  ! for internal calculation.  If temperature is 0 or less then unity
957  ! is returned; since this function is primarily used to weight an
958  ! integral by the Planck function, a temperature of 0 or less means
959  ! no weighting is to be applied.
960  elemental function calc_planck_function_wavenumber(wavenumber, temperature)
961
962    use parkind1,            only : jprb, jprd
963    use radiation_constants, only : SpeedOfLight, BoltzmannConstant, PlanckConstant
964
965    real(jprb), intent(in) :: wavenumber  ! cm-1
966    real(jprb), intent(in) :: temperature ! K
967    real(jprb) :: calc_planck_function_wavenumber
968
969    real(jprd) :: freq ! Hz
970    real(jprd) :: planck_fn_freq ! W m-2 Hz-1
971
972    if (temperature > 0.0_jprd) then
973      freq = 100.0_jprd * real(SpeedOfLight,jprd) * real(wavenumber,jprd)
974      planck_fn_freq = 2.0_jprd * real(PlanckConstant,jprd) * freq**3 &
975           &  / (real(SpeedOfLight,jprd)**2 * (exp(real(PlanckConstant,jprd)*freq &
976           &     /(real(BoltzmannConstant,jprd)*real(temperature,jprd))) - 1.0_jprd))
977      calc_planck_function_wavenumber = real(planck_fn_freq * 100.0_jprd * real(SpeedOfLight,jprd), jprb)
978    else
979      calc_planck_function_wavenumber = 1.0_jprb
980    end if
981
982  end function calc_planck_function_wavenumber
983
984end module radiation_spectral_definition
Note: See TracBrowser for help on using the repository browser.