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

Last change on this file since 5418 was 4853, checked in by idelkadi, 9 months ago

Ecrad update in LMDZ / Implementation of Ecrad double call in LMDZ

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