source: LMDZ6/trunk/libf/phylmd/ecrad/radiation/radiation_aerosol_optics_description.F90 @ 5325

Last change on this file since 5325 was 4853, checked in by idelkadi, 8 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: 13.8 KB
Line 
1! radiation_aerosol_optics_description.F90 - Type to store aerosol optics metadata
2!
3! (C) Copyright 2022- 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
16#include "ecrad_config.h"
17
18module radiation_aerosol_optics_description
19
20  use parkind1,      only : jprb
21
22  implicit none
23  public
24
25  !---------------------------------------------------------------------
26  ! This type holds the metadata from an aerosol optical property
27  ! file, enabling the user to request the index to the aerosol type
28  ! with particular properties.  Note that string information is held
29  ! in the form of 2D arrays of single characters, so comparison to
30  ! character strings requires the to_string helper function at the
31  ! end of this file.
32  type aerosol_optics_description_type
33
34    ! Two-character code describing the aerosol family, dimensioned
35    ! (2,naer), e.g.
36    !   SS: Sea salt
37    !   OM: Organic matter
38    !   SU: Sulfate
39    !   OB: Secondary organic biogenic
40    !   OA: Secondary organic anthropogenic
41    !   AM: Fine-mode ammonium sulfate
42    !   NI: Nitrate
43    !   DD: Desert dust
44    !   BC: Black carbon
45    character(len=1), allocatable :: code_phobic(:,:)
46    character(len=1), allocatable :: code_philic(:,:)
47
48    ! Size bin, typically 1-2 or 1-3 in order fine to coarse, or zero
49    ! if no division by size is used, dimensioned (naer)
50    integer, allocatable :: bin_phobic(:)
51    integer, allocatable :: bin_philic(:)
52
53    ! Character string characterizing the optical model, e.g. OPAC,
54    ! GACP, GLOMAP, Dubovik2002 etc.
55    character(len=1), allocatable :: optical_model_phobic(:,:)
56    character(len=1), allocatable :: optical_model_philic(:,:)
57
58    ! The user can call preferred_optical_model to specify that a
59    ! certain optical model for a certain aerosol family is to be
60    ! preferred when get_index is called
61    logical, allocatable :: is_preferred_phobic(:)
62    logical, allocatable :: is_preferred_philic(:)
63
64    ! Verbosity level
65    integer :: iverbose
66   
67  contains
68    procedure :: read
69    procedure :: preferred_optical_model
70    procedure :: get_index
71
72  end type aerosol_optics_description_type
73
74contains
75
76  !---------------------------------------------------------------------
77  ! Read optical property file file_name into an
78  ! aerosol_optics_description_type object
79  subroutine read(this, file_name, iverbose)
80
81    use yomhook,              only : lhook, dr_hook, jphook
82#ifdef EASY_NETCDF_READ_MPI
83    use easy_netcdf_read_mpi, only : netcdf_file
84#else
85    use easy_netcdf,          only : netcdf_file
86#endif
87
88    class(aerosol_optics_description_type), intent(inout) :: this
89    character(len=*), intent(in)              :: file_name
90    integer, intent(in), optional             :: iverbose
91   
92    ! The NetCDF file containing the aerosol optics data
93    type(netcdf_file)  :: file
94
95    real(jphook)       :: hook_handle
96
97    if (lhook) call dr_hook('radiation_aerosol_optics_description:load',0,hook_handle)
98
99    ! Open the aerosol scattering file and configure the way it is
100    ! read
101    call file%open(trim(file_name), iverbose=iverbose)
102
103    ! Read metadata variables
104    call file%get('code_hydrophilic', this%code_philic)
105    call file%get('code_hydrophobic', this%code_phobic)
106    call file%get('bin_hydrophilic',  this%bin_philic)
107    call file%get('bin_hydrophobic',  this%bin_phobic)
108    call file%get('optical_model_hydrophilic', this%optical_model_philic)
109    call file%get('optical_model_hydrophobic', this%optical_model_phobic)
110
111    ! Allocate logical arrays of the appropriate size and set to FALSE
112    allocate(this%is_preferred_philic(size(this%bin_philic)))
113    allocate(this%is_preferred_phobic(size(this%bin_phobic)))
114    this%is_preferred_philic = .false.
115    this%is_preferred_phobic = .false.
116
117    call file%close()
118
119    if (present(iverbose)) then
120      this%iverbose = iverbose
121    else
122      this%iverbose = 3
123    end if
124   
125    if (lhook) call dr_hook('radiation_aerosol_optics_description:load',1,hook_handle)
126
127  end subroutine read
128
129  !---------------------------------------------------------------------
130  ! Specify the preferred optical model for a particular aerosol
131  ! family, e.g. "call
132  ! aer_desc%preferred_optical_model('DD','Woodward2001')" would mean
133  ! that subsequent calls to get_index in which the optical model is
134  ! not specified would return the Woodward model rather than the
135  ! first matching model in the file.  The check is only done on the
136  ! first len(optical_model_str) characters, so "Woodward" and
137  ! "Woodward2001" would both match the Woodward2001 model.
138  subroutine preferred_optical_model(this, code_str, optical_model_str)
139
140    use yomhook,              only : lhook, dr_hook, jphook
141    use radiation_io,         only : nulout, nulerr, radiation_abort
142   
143    class(aerosol_optics_description_type), intent(inout) :: this
144    character(len=2), intent(in) :: code_str
145    character(len=*), intent(in) :: optical_model_str
146
147    ! Aerosol loop counter
148    integer :: ja
149
150    logical :: is_found, is_philic, is_phobic
151
152    real(jphook)         :: hook_handle
153
154    if (lhook) call dr_hook('radiation_aerosol_optics_description:preferred_optical_model',0,hook_handle)
155
156    ! Check for empty string
157    if (optical_model_str == ' ') then
158      if (lhook) call dr_hook('radiation_aerosol_optics_description:preferred_optical_model',1,hook_handle)
159      return
160    end if
161
162    is_found  = .false.
163    is_philic = .false.
164    is_phobic = .false.
165   
166    ! Loop over hydrophilic types
167    do ja = 1,size(this%bin_philic)
168      ! Check if we have a match
169      if (to_string(this%code_philic(:,ja)) == code_str &
170           &  .and. trim(to_string(this%optical_model_philic(:,ja))) &
171           &          == optical_model_str) then
172        this%is_preferred_philic(ja) = .true.
173        is_found  = .true.
174        is_philic = .true.
175      end if
176    end do
177    ! Repeat for the hydrophobic types
178    do ja = 1,size(this%bin_phobic)
179      if (to_string(this%code_phobic(:,ja)) == code_str &
180           &  .and. trim(to_string(this%optical_model_phobic(:,ja))) &
181           &          == optical_model_str) then
182        this%is_preferred_phobic(ja) = .true.
183        is_found  = .true.
184        is_phobic = .true.
185      end if
186    end do
187
188    if (.not. is_found) then
189      write(nulerr,'(a,a2,a,a,a)') '*** Error: Preferred "', code_str ,'" aerosol optical model "', &
190           &  trim(optical_model_str), '" not found in file'
191      call radiation_abort()
192    else if (this%iverbose > 2) then
193      write(nulout,'(a,a2,a,a,a)',advance='no') 'Preferred "', code_str, '" aerosol optical model set to "', &
194           &  trim(optical_model_str), '" ('
195      if (is_phobic) then
196        write(nulout,'(a)',advance='no') ' hydrophobic'
197      end if
198      if (is_philic) then
199        write(nulout,'(a)',advance='no') ' hydrophilic'
200      end if
201      write(nulout,'(a)') ' )'
202    end if
203   
204    if (lhook) call dr_hook('radiation_aerosol_optics_description:preferred_optical_model',1,hook_handle)
205
206  end subroutine preferred_optical_model
207
208 
209  !---------------------------------------------------------------------
210  ! Return the index to the aerosol optical properties corresponding
211  ! to the aerosol family in code_str (e.g. SS, DD etc), whether or
212  ! not the requested aerosol is hydrophilic in the logical
213  ! lhydrophilic, and optionally the size bin ibin and optical model
214  ! in optical_model_str. The return value may be used to populate the
215  ! radiation_config%i_aerosol_map vector, where a positive number is
216  ! a hydrophobic index, a negative number is a hydrophilic index and
217  ! zero indicates that the aerosol type was not found in the file.
218  ! This is a valid entry in i_aerosol_map meaning the aerosol is
219  ! ignored, but the calling routine to get_index might wish to throw
220  ! a warning or error. This routine works by assigning a score based
221  ! on the closeness of the match.
222  function get_index(this, code_str, lhydrophilic, ibin, optical_model_str)
223   
224    use yomhook,              only : lhook, dr_hook, jphook
225    use radiation_io,         only : nulout
226
227    class(aerosol_optics_description_type), intent(in) :: this
228    character(len=2), intent(in) :: code_str
229    logical, intent(in) :: lhydrophilic
230    integer, intent(in), optional :: ibin
231    character(len=*), intent(in), optional :: optical_model_str
232
233    ! Score of the currently selected aerosol index, and the score of
234    ! the current one under consideration
235    integer :: score, current_score
236
237    ! Loop index for aerosol type
238    integer :: ja
239
240    ! Return value
241    integer :: get_index
242
243    ! Issue a warning if there is more than one equal match
244    logical :: is_ambiguous
245
246    real(jphook) :: hook_handle
247
248    if (lhook) call dr_hook('radiation_aerosol_optics_description:get_index',0,hook_handle)
249
250    ! Initial values
251    get_index = 0
252    score = 0
253    is_ambiguous = .false.
254
255    if (lhydrophilic) then
256      ! Loop over hydrophilic aerosol types
257      do ja = 1,size(this%bin_philic)
258        current_score = 0
259        if (to_string(this%code_philic(:,ja)) == code_str) then
260          ! Aerosol code matches
261          if (present(ibin) .and. this%bin_philic(ja) > 0) then
262            if (ibin > 0) then
263              if (ibin == this%bin_philic(ja)) then
264                ! Requested bin number matches
265                current_score = 4
266              else
267                ! Requested bin number does not match
268                current_score = -1
269              end if
270            else
271              ! Bin number is zero: no request
272              current_score = 2
273            end if
274          else
275            ! No bin number present
276            current_score = 2
277          end if
278          if (present(optical_model_str)) then
279            if (trim(to_string(this%optical_model_philic(:,ja))) &
280                 &  == optical_model_str) then
281              ! Requested optical model matches
282              if (current_score >= 0) then
283                current_score = current_score + 4
284              end if
285            else
286              ! Requested optical model does not match
287              current_score = -1
288            end if
289          else if (current_score >= 0) then
290            ! No requested optical model
291            current_score = current_score + 2
292          end if
293          if (current_score > 0 .and. this%is_preferred_philic(ja)) then
294            current_score = current_score + 1
295          end if
296          if (current_score > score) then
297            ! Better score than any existing aerosol type
298            get_index = -ja
299            score = current_score
300            is_ambiguous = .false.
301          else if (current_score > 0 .and. current_score == score) then
302            is_ambiguous = .true.
303          end if
304        end if
305      end do
306    else
307      ! Loop over hydrophobic aerosol types
308      do ja = 1,size(this%bin_phobic)
309        current_score = 0
310        if (to_string(this%code_phobic(:,ja)) == code_str) then
311          ! Aerosol code matches
312          if (present(ibin) .and. this%bin_phobic(ja) > 0) then
313            if (ibin > 0) then
314              if (ibin == this%bin_phobic(ja)) then
315                ! Requested bin number matches
316                current_score = 4
317              else
318                ! Requested bin number does not match
319                current_score = -1
320              end if
321            else
322              ! Bin number is zero: no request
323              current_score = 2
324            end if
325          else
326            ! No bin number requested or present
327            current_score = 2
328          end if
329          if (present(optical_model_str)) then
330            if (trim(to_string(this%optical_model_phobic(:,ja))) &
331                 &  == optical_model_str) then
332              ! Requested optical model matches
333              if (current_score >= 0) then
334                current_score = current_score + 4
335              end if
336            else
337              ! Requested optical model does not match
338              current_score = -1
339            end if
340          else if (current_score >= 0) then
341            ! No requested optical model
342            current_score = current_score + 2
343          end if
344          if (current_score > 0 .and. this%is_preferred_phobic(ja)) then
345            current_score = current_score + 1
346          end if
347          if (current_score > score) then
348            ! Better score than any existing aerosol type
349            get_index = ja
350            score = current_score
351            is_ambiguous = .false.
352          else if (current_score > 0 .and. current_score == score) then
353            is_ambiguous = .true.
354          end if         
355        end if
356      end do
357    end if
358
359    if (is_ambiguous) then
360      write(nulout,'(a,a2,a,l,a)') 'Warning: radiation_aerosol_optics_description:get_index("', &
361           &  code_str, '",', lhydrophilic, &
362           &  ',...) does not unambiguously identify an aerosol optical property index'
363    end if
364   
365    if (lhook) call dr_hook('radiation_aerosol_optics_description:get_index',1,hook_handle)
366
367  end function get_index
368
369  !---------------------------------------------------------------------
370  ! Utility function to convert an array of single characters to a
371  ! character string (yes Fortran's string handling is a bit
372  ! rubbish). We set NULL characters (ASCII code 0) returned from the
373  ! NetCDF library to spaces, so that TRIM can remove them.
374  pure function to_string(arr) result(str)
375    character, intent(in)  :: arr(:)
376    character(len=size(arr)) :: str
377    integer :: jc
378    do jc = 1,size(arr)
379      if (ichar(arr(jc)) == 0) then
380        ! Replace NULL character with a space
381        str(jc:jc) = ' '
382      else
383        str(jc:jc) = arr(jc)
384      end if
385    end do
386  end function to_string
387
388end module radiation_aerosol_optics_description
Note: See TracBrowser for help on using the repository browser.