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