ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_RadConstants.H
Go to the documentation of this file.
1 //
2 // This module contains constexprants that are specific to the radiative transfer
3 // code used in the RRTMG model.
4 //
5 #ifndef ERF_RAD_CONSTANTS_H
6 #define ERF_RAD_CONSTANTS_H
7 
8 #include "rrtmgp_const.h"
9 using yakl::fortran::parallel_for;
10 using yakl::fortran::SimpleBounds;
11 
12 class RadConstants {
13  public:
14  // number of shorwave spectral intervals
15  static constexpr int nswbands = 14;
16  static constexpr int nbndsw = 14;
17 
18  // Wavenumbers of band boundaries
19  //
20  // Note: Currently rad_solar_var extends the lowest band down to
21  // 100 cm^-1 if it is too high to cover the far-IR. Any changes meant
22  // to affect IR solar variability should take note of this.
23 
24  static constexpr real wavenum_low[] = // in cm^-1
25  {2600., 3250., 4000., 4650., 5150., 6150., 7700.,
26  8050., 12850., 16000., 22650., 29000., 38000., 820.};
27 
28  static constexpr real wavenum_high[] = // in cm^-1
29  {3250., 4000., 4650., 5150., 6150., 7700., 8050.,
30  12850.,16000., 22650., 29000., 38000., 50000., 2600.};
31 
32  // Solar irradiance at 1 A.U. in W/m^2 assumed by radiation code
33  // Rescaled so that sum is precisely 1368.22 and fractional amounts sum to 1.0
34  static constexpr real solar_ref_band_irradiance[] =
35  {12.11, 20.3600000000001, 23.73,
36  22.43, 55.63, 102.93, 24.29,
37  345.74, 218.19, 347.20,
38  129.49, 50.15, 3.08, 12.89};
39 
40  // parameterizations related to radiation (when microphysics model doesn't provide)
41  static constexpr real icesize_table_min_temp = 180.;
42  static constexpr real retab[] = { 5.92779, 6.26422, 6.61973, 6.99539, 7.39234,
43  7.81177, 8.25496, 8.72323, 9.21800, 9.74075, 10.2930,
44  10.8765, 11.4929, 12.1440, 12.8317, 13.5581, 14.2319,
45  15.0351, 15.8799, 16.7674, 17.6986, 18.6744, 19.6955,
46  20.7623, 21.8757, 23.0364, 24.2452, 25.5034, 26.8125,
47  27.7895, 28.6450, 29.4167, 30.1088, 30.7306, 31.2943,
48  31.8151, 32.3077, 32.7870, 33.2657, 33.7540, 34.2601,
49  34.7892, 35.3442, 35.9255, 36.5316, 37.1602, 37.8078,
50  38.4720, 39.1508, 39.8442, 40.5552, 41.2912, 42.0635,
51  42.8876, 43.7863, 44.7853, 45.9170, 47.2165, 48.7221,
52  50.4710, 52.4980, 54.8315, 57.4898, 60.4785, 63.7898,
53  65.5604, 71.2885, 75.4113, 79.7368, 84.2351, 88.8833,
54  93.6658, 98.5739, 103.603, 108.752, 114.025, 119.424,
55  124.954, 130.630, 136.457, 142.446, 148.608, 154.956,
56  161.503, 168.262, 175.248, 182.473, 189.952, 197.699,
57  205.728, 214.055, 222.694, 231.661, 240.971, 250.639};
58 
59  // These are indices to the band for diagnostic output
60  static constexpr int idx_sw_diag = 10; // index to sw visible band
61  static constexpr int idx_nir_diag = 8; // index to sw near infrared (778-1240 nm) band
62  static constexpr int idx_uv_diag = 11; // index to sw uv (345-441 nm) band
63 
64  static constexpr int rrtmg_sw_cloudsim_band = 9; // rrtmg band for .67 micron
65 
66  // Number of evenly spaced intervals in rh
67  // The globality of this mesh may not be necessary
68  // Perhaps it could be specific to the aerosol
69  // But it is difficult to see how refined it must be
70  // for lookup. This value was found to be sufficient
71  // for Sulfate and probably necessary to resolve the
72  // high variation near rh = 1. Alternative methods
73  // were found to be too slow.
74  // Optimal approach would be for cam to specify size of aerosol
75  // based on each aerosol's characteristics. Radiation
76  // should know nothing about hygroscopic growth!
77  static constexpr int nrh = 1000;
78 
79  // LONGWAVE DATA
80 
81  // These are indices to the band for diagnostic output
82  static constexpr int idx_lw_diag = 7; // index to (H20 window) LW band
83  static constexpr int rrtmg_lw_cloudsim_band = 6; // rrtmg band for 10.5 micron
84 
85  // number of lw bands
86  static constexpr int nlwbands = 16;
87  static constexpr int nbndlw = 16;
88 
89  static constexpr real wavenumber1_longwave[] = // Longwave spectral band limits (cm-1)
90  {10., 350., 500., 630., 700., 820., 980., 1080.,
91  1180., 1390., 1480., 1800., 2080., 2250., 2390., 2600.};
92 
93  static constexpr real wavenumber2_longwave[] = // Longwave spectral band limits (cm-1)
94  {350., 500., 630., 700., 820., 980., 1080., 1180.,
95  1390., 1480., 1800., 2080., 2250., 2390., 2600., 3250.};
96 
97  //GASES TREATED BY RADIATION (line spectrae)
98 
99  // gasses required by radiation
100  static constexpr int nradgas = 8;
101  static constexpr const char* gaslist[] =
102  {"H2O","O3", "O2", "CO2", "N2O", "CH4", "CFC11", "CFC12"};
103 
104  // use enum unit
105  enum Units {
110  centimeter
111  };
112 
113  // provide Solar Irradiance for each band in RRTMG
114  inline
115  static void get_solar_band_fraction_irrad (real1d& fractional_irradiance)
116  {
117  real tsi = 0; // total solar irradiance
118  for(auto i = 0; i < nswbands; ++i)
119  tsi = tsi + solar_ref_band_irradiance[i];
120 
121  for(auto i = 0; i < nswbands; ++i)
122  fractional_irradiance(i) = solar_ref_band_irradiance[i] / tsi;
123  }
124 
125  // provide Total Solar Irradiance assumed by RRTMG
126  inline
127  static void get_ref_total_solar_irrad (real& tsi)
128  {
129  tsi = 0;
130  for(auto i = 0; i < nswbands; ++i)
131  tsi = tsi + solar_ref_band_irradiance[i];
132  }
133 
134  // solar irradiance in each band (W/m^2)
135  inline
136  static void get_ref_solar_band_irrad (real1d& band_irrad )
137  {
138  for(auto i = 0; i < nswbands; ++i)
139  band_irrad(i) = solar_ref_band_irradiance[i];
140  }
141 
142  // number of solar (shortwave) bands in the rrtmg code
143  inline
144  static void get_number_sw_bands (int& number_of_bands)
145  {
146  number_of_bands = nswbands;
147  }
148 
149  // provide spectral boundaries of each longwave band
150  inline
151  static void get_lw_spectral_boundaries (real1d& low_boundaries,
152  real1d& high_boundaries,
153  Units units)
154  {
155  switch (units) {
156  case inv_cm:
157  for(auto i = 0; i < nlwbands; ++i) {
158  low_boundaries(i) = wavenumber1_longwave[i];
159  high_boundaries(i) = wavenumber2_longwave[i];
160  }
161  case meter:
162  for(auto i = 0; i < nlwbands; ++i) {
163  low_boundaries(i) = 1.e-2/wavenumber2_longwave[i];
164  high_boundaries(i) = 1.e-2/wavenumber1_longwave[i];
165  }
166  case nanometer:
167  for(auto i = 0; i < nlwbands; ++i) {
168  low_boundaries(i) = 1.e7/wavenumber2_longwave[i];
169  high_boundaries(i) = 1.e7/wavenumber1_longwave[i];
170  }
171  case micrometer:
172  for(auto i = 0; i < nlwbands; ++i) {
173  low_boundaries(i) = 1.e4/wavenumber2_longwave[i];
174  high_boundaries(i) = 1.e4/wavenumber1_longwave[i];
175  }
176  case centimeter:
177  for(auto i = 0; i < nlwbands; ++i) {
178  low_boundaries(i) = 1./wavenumber2_longwave[i];
179  high_boundaries(i) = 1./wavenumber1_longwave[i];
180  }
181  default:
182  amrex::Print() << "get_lw_spectral_boundaries: spectral units not acceptable\n";
183  }
184  }
185 
186  inline
187  static void get_lw_spectral_midpoints (real1d& band_midpoints, Units units)
188  {
189  real1d lower_boundaries("lower_boundaries", nlwbands);
190  real1d upper_boundaries("upper_boundaries", nlwbands);
191  //int iband;
192 
193  // Get band boundaries
194  get_lw_spectral_boundaries(lower_boundaries, upper_boundaries, units);
195 
196  // Get band midpoints
197  for(auto iband = 0; iband < nlwbands; ++iband) {
198  band_midpoints(iband) = 0.5*(lower_boundaries(iband) + upper_boundaries(iband));
199  }
200 
201  }
202 
203  // provide spectral boundaries of each shortwave band
204  inline
205  static void get_sw_spectral_boundaries (real1d& low_boundaries, real1d& high_boundaries, Units units)
206  {
207  switch (units) {
208  case inv_cm:
209  for(auto i = 0; i < nswbands; ++i) {
210  low_boundaries(i) = wavenum_low[i];
211  high_boundaries(i) = wavenum_high[i];
212  }
213  case meter:
214  for(auto i = 0; i < nswbands; ++i) {
215  low_boundaries(i) = 1.e-2/wavenum_high[i];
216  high_boundaries(i) = 1.e-2/wavenum_low[i];
217  }
218  case nanometer:
219  for(auto i = 0; i < nswbands; ++i) {
220  low_boundaries(i) = 1.e7/wavenum_high[i];
221  high_boundaries(i) = 1.e7/wavenum_low[i];
222  }
223  case micrometer:
224  for(auto i = 0; i < nswbands; ++i) {
225  low_boundaries(i) = 1.e4/wavenum_high[i];
226  high_boundaries(i) = 1.e4/wavenum_low[i];
227  }
228  case centimeter:
229  for(auto i = 0; i < nswbands; ++i) {
230  low_boundaries(i) = 1./wavenum_high[i];
231  high_boundaries(i) = 1./wavenum_low[i];
232  }
233  default:
234  amrex::Print() << "rad_constants.F90: spectral units not acceptable\n";
235  }
236  }
237 
238  inline
239  static void get_sw_spectral_midpoints (real1d& band_midpoints, Units units)
240  {
241  real1d lower_boundaries("lower_boundaries", nswbands);
242  real1d upper_boundaries("upper_boundaries", nswbands);
243 
244  // Get band boundaries
245  get_sw_spectral_boundaries(lower_boundaries, upper_boundaries, units);
246 
247  //Get band midpoints
248  for(auto iband = 0; iband < nswbands; ++iband) {
249  band_midpoints(iband) = 0.5 * (lower_boundaries(iband) + upper_boundaries(iband));
250  }
251  }
252 
253  // return the index in the gaslist array of the specified gasname
254  inline
255  static int rad_gas_index (std::string gasname)
256  {
257  for(auto igas = 0; igas < nradgas; ++igas) {
258  if (gaslist[igas] == gasname) return igas;
259  }
260  amrex::Abort("rad_gas_index: can not find gas with name " + gasname);
261  return -1;
262  }
263 };
264 #endif
Definition: ERF_RadConstants.H:12
static constexpr real wavenum_high[]
Definition: ERF_RadConstants.H:28
static constexpr real solar_ref_band_irradiance[]
Definition: ERF_RadConstants.H:34
static void get_lw_spectral_midpoints(real1d &band_midpoints, Units units)
Definition: ERF_RadConstants.H:187
static void get_lw_spectral_boundaries(real1d &low_boundaries, real1d &high_boundaries, Units units)
Definition: ERF_RadConstants.H:151
static constexpr int nbndsw
Definition: ERF_RadConstants.H:16
Units
Definition: ERF_RadConstants.H:105
@ centimeter
Definition: ERF_RadConstants.H:110
@ meter
Definition: ERF_RadConstants.H:107
@ nanometer
Definition: ERF_RadConstants.H:108
@ inv_cm
Definition: ERF_RadConstants.H:106
@ micrometer
Definition: ERF_RadConstants.H:109
static int rad_gas_index(std::string gasname)
Definition: ERF_RadConstants.H:255
static void get_ref_total_solar_irrad(real &tsi)
Definition: ERF_RadConstants.H:127
static constexpr int idx_lw_diag
Definition: ERF_RadConstants.H:82
static constexpr int nlwbands
Definition: ERF_RadConstants.H:86
static constexpr int idx_nir_diag
Definition: ERF_RadConstants.H:61
static void get_sw_spectral_boundaries(real1d &low_boundaries, real1d &high_boundaries, Units units)
Definition: ERF_RadConstants.H:205
static constexpr real wavenumber2_longwave[]
Definition: ERF_RadConstants.H:93
static constexpr int nbndlw
Definition: ERF_RadConstants.H:87
static constexpr int idx_sw_diag
Definition: ERF_RadConstants.H:60
static constexpr int nrh
Definition: ERF_RadConstants.H:77
static void get_solar_band_fraction_irrad(real1d &fractional_irradiance)
Definition: ERF_RadConstants.H:115
static constexpr real icesize_table_min_temp
Definition: ERF_RadConstants.H:41
static constexpr real wavenum_low[]
Definition: ERF_RadConstants.H:24
static constexpr int idx_uv_diag
Definition: ERF_RadConstants.H:62
static constexpr const char * gaslist[]
Definition: ERF_RadConstants.H:101
static constexpr int nswbands
Definition: ERF_RadConstants.H:15
static constexpr int nradgas
Definition: ERF_RadConstants.H:100
static constexpr int rrtmg_sw_cloudsim_band
Definition: ERF_RadConstants.H:64
static void get_sw_spectral_midpoints(real1d &band_midpoints, Units units)
Definition: ERF_RadConstants.H:239
static constexpr real wavenumber1_longwave[]
Definition: ERF_RadConstants.H:89
static void get_ref_solar_band_irrad(real1d &band_irrad)
Definition: ERF_RadConstants.H:136
static constexpr int rrtmg_lw_cloudsim_band
Definition: ERF_RadConstants.H:83
static constexpr real retab[]
Definition: ERF_RadConstants.H:42
static void get_number_sw_bands(int &number_of_bands)
Definition: ERF_RadConstants.H:144