source: LMDZ6/trunk/libf/phylmd/ecrad/radiation/radiation_ecckd_gas.F90 @ 5327

Last change on this file since 5327 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: 4.2 KB
Line 
1! radiation_ecckd_gas.F90 - type representing a single ecCKD gas
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_ecckd_gas
20
21  use parkind1, only : jprb
22  use radiation_gas_constants
23
24  implicit none
25
26  public
27
28  ! Concentration dependence of individual gases
29  enum, bind(c)
30    enumerator :: IConcDependenceNone = 0, &
31         &        IConcDependenceLinear, &
32         &        IConcDependenceLUT, &
33         &        IConcDependenceRelativeLinear
34  end enum
35
36  !---------------------------------------------------------------------
37  ! This derived type describes a correlated k-distribution
38  ! representation of an individual gas (including composite gases)
39  type ckd_gas_type
40
41    ! Code identifying the gas, from the codes in the
42    ! radiation_gas_constants module
43    integer :: i_gas_code = -1
44
45    ! One of the IConcDependence* enumerators
46    integer :: i_conc_dependence
47
48    ! Molar absorption coefficient in m2 mol-1. If
49    ! i_conc_dependence==IConcDependenceNone then it is the absorption
50    ! cross section per mole of dry air.  If
51    ! conc_dependence==IConcDependenceLinear|IConcDependenceRelativeLinear,
52    ! it is the absorption cross section per mole of the gas in
53    ! question. It is dimensioned (g_point,pressure,temperature).
54    real(jprb), allocatable :: molar_abs(:,:,:)
55   
56    ! If i_conc_dependence==IConcDependenceLUT then we have an
57    ! additional dimension for concentration. It is dimensioned
58    ! (g_point,pressure,temperature,conc)
59    real(jprb), allocatable :: molar_abs_conc(:,:,:,:)
60
61    ! If i_conc_dependence==IConcDependenceRelativeLinear then the
62    ! following reference concentration is subtracted from the actual
63    ! concentration before the result is multiplied by the mass
64    ! absorption coefficient
65    real(jprb) :: reference_mole_frac = 0.0_jprb
66
67    ! Mole fraction coordinate variable if
68    ! i_conc_dependence==IConcDependenceLUT
69    real(jprb) :: log_mole_frac1 = 0.0_jprb, d_log_mole_frac = 1.0_jprb
70    integer    :: n_mole_frac = 0
71
72  contains
73
74    procedure :: read => read_ckd_gas
75!    procedure :: deallocate => deallocate_ckd_gas
76
77  end type ckd_gas_type
78
79contains
80
81  !---------------------------------------------------------------------
82  ! Read information about the representation of a single gas from a
83  ! NetCDF file, identifying it with code i_gas_code
84  subroutine read_ckd_gas(this, file, gas_name, i_gas_code)
85
86#ifdef EASY_NETCDF_READ_MPI
87    use easy_netcdf_read_mpi, only : netcdf_file
88#else
89    use easy_netcdf,          only : netcdf_file
90#endif
91
92    class(ckd_gas_type), intent(inout) :: this
93    type(netcdf_file),   intent(inout) :: file
94    character(len=*),    intent(in)    :: gas_name
95    integer,             intent(in)    :: i_gas_code
96   
97    ! Local storage for mole fraction coordinate variable
98    real(jprb), allocatable :: mole_fraction(:)
99
100    this%i_gas_code = i_gas_code
101
102    call file%get(gas_name // "_conc_dependence_code", this%i_conc_dependence)
103    if (this%i_conc_dependence == IConcDependenceLut) then
104      call file%get(gas_name // "_molar_absorption_coeff", &
105           &        this%molar_abs_conc)
106      call file%get(gas_name // "_mole_fraction", mole_fraction)
107      this%log_mole_frac1  = log(mole_fraction(1))
108      this%n_mole_frac     = size(mole_fraction)
109      this%d_log_mole_frac = (log(mole_fraction(size(mole_fraction))) &
110           &                  - this%log_mole_frac1) / (this%n_mole_frac-1)
111      deallocate(mole_fraction)
112    else
113      call file%get(gas_name // "_molar_absorption_coeff", &
114           &        this%molar_abs)
115    end if
116
117    if (this%i_conc_dependence == IConcDependenceRelativeLinear) then
118      call file%get(gas_name // "_reference_mole_fraction", &
119           &        this%reference_mole_frac)
120    end if
121
122  end subroutine read_ckd_gas
123
124end module radiation_ecckd_gas
Note: See TracBrowser for help on using the repository browser.