ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_Radiation.H
Go to the documentation of this file.
1 #ifndef ERF_RADIATION_H
2 #define ERF_RADIATION_H
3 
4 /*
5  * RTE-RRTMGP radiation model interface to ERF
6  * The original code is developed by RobertPincus, and the code is open source available at:
7  * https://github.com/earth-system-radiation/rte-rrtmgp
8  * Please reference to the following paper,
9  * https://agupubs.onlinelibrary.wiley.com/doi/10.1029/2019MS001621
10  * NOTE: we use the C++ version of RTE-RRTMGP, which is the implementation of the original Fortran
11  * code using C++ KOKKOS for CUDA, HiP and SYCL application by E3SM ECP team, the C++ version
12  * of the rte-rrtmgp code is located at:
13  * https://github.com/E3SM-Project/rte-rrtmgp
14  * The RTE-RRTMGP uses BSD-3-Clause Open Source License, if you want to make changes,
15  * and modifications to the code, please refer to BSD-3-Clause Open Source License.
16  */
17 
18 #include <ctime>
19 #include <string>
20 #include <vector>
21 #include <memory>
22 
23 #include <mo_gas_concentrations.h>
24 
25 #include <Kokkos_Core.hpp>
26 
27 #include <AMReX_ParmParse.H>
28 #include <AMReX_FArrayBox.H>
29 #include <AMReX_Geometry.H>
30 #include <AMReX_TableData.H>
31 #include <AMReX_MultiFabUtil.H>
32 #include <AMReX_PlotFileUtil.H>
33 
34 #include <ERF_RadiationInterface.H>
35 #include <ERF_RRTMGP_Interface.H>
36 #include <ERF_RRTMGP_Utils.H>
37 #include <ERF_OrbCosZenith.H>
38 #include <ERF_Constants.H>
39 #include <ERF_IndexDefines.H>
40 #include <ERF_DataStruct.H>
41 #include <ERF_EOS.H>
42 
43 #include <ERF_LandSurface.H>
44 
45 class Radiation : public IRadiation {
46 public:
47 
48  #include <ERF_Kokkos.H>
49 
50  // Constructor
51  Radiation (const int& lev,
52  SolverChoice& sc);
53 
54  // Destructor
56  {
57  // Note that Kokkos is now finalized in main.cpp
58  }
59 
60  virtual
61  void
62  Init (const amrex::Geometry& geom,
63  const amrex::BoxArray& ba,
64  amrex::MultiFab* cons_in) override
65  {
66  // Ensure the boxes span klo -> khi
67  int klo = geom.Domain().smallEnd(2);
68  int khi = geom.Domain().bigEnd(2);
69 
70  // Reset vector of offsets for columnar data
71  m_nlay = geom.Domain().length(2);
72 
73  m_ncol = 0;
74  m_col_offsets.clear();
75  m_col_offsets.resize(int(ba.size()));
76  for (amrex::MFIter mfi(*cons_in); mfi.isValid(); ++mfi) {
77  const amrex::Box& vbx = mfi.validbox();
78  AMREX_ALWAYS_ASSERT_WITH_MESSAGE((klo == vbx.smallEnd(2)) &&
79  (khi == vbx.bigEnd(2)),
80  "Vertical decomposition with radiation is not allowed.");
81  int nx = vbx.length(0);
82  int ny = vbx.length(1);
83  m_col_offsets[mfi.index()] = m_ncol;
84  m_ncol += nx * ny;
85  }
86  };
87 
88  virtual
89  void
90  Run (int& level,
91  int& step,
92  amrex::Real& time,
93  const amrex::Real& dt,
94  const amrex::BoxArray& ba,
95  amrex::Geometry& geom,
96  amrex::MultiFab* cons_in,
97  amrex::iMultiFab* lmask,
98  amrex::MultiFab* t_surf,
99  amrex::MultiFab* lsm_fluxes,
100  amrex::MultiFab* lsm_zenith,
101  amrex::Vector<amrex::MultiFab*>& lsm_input_ptrs,
102  amrex::Vector<amrex::MultiFab*>& lsm_output_ptrs,
103  amrex::MultiFab* qheating_rates,
104  amrex::MultiFab* rad_fluxes,
105  amrex::MultiFab* z_phys,
106  amrex::MultiFab* lat_ptr,
107  amrex::MultiFab* lon_ptr) override
108  {
109  set_grids(level, step, time, dt, ba, geom,
110  cons_in, lmask, t_surf, lsm_fluxes,
111  lsm_zenith, lsm_input_ptrs, qheating_rates,
112  rad_fluxes, z_phys, lat_ptr, lon_ptr);
113  rad_run_impl(lsm_output_ptrs);
114  }
115 
116  // Set the grid info for columnar data in KOKKOS
117  void
118  set_grids (int& level,
119  int& step,
120  amrex::Real& time,
121  const amrex::Real& dt,
122  const amrex::BoxArray& ba,
123  amrex::Geometry& geom,
124  amrex::MultiFab* cons_in,
125  amrex::iMultiFab* lmask,
126  amrex::MultiFab* t_surf,
127  amrex::MultiFab* lsm_fluxes,
128  amrex::MultiFab* lsm_zenith,
129  amrex::Vector<amrex::MultiFab*>& lsm_input_ptrs,
130  amrex::MultiFab* qheating_rates,
131  amrex::MultiFab* rad_fluxes,
132  amrex::MultiFab* z_phys,
133  amrex::MultiFab* lat,
134  amrex::MultiFab* lon);
135 
136  // Initialize the temporary variables
137  void
138  alloc_buffers ();
139 
140  // Clear the temporary variables
141  void
142  dealloc_buffers ();
143 
144  // Fill KOKKOS Views from AMReX MultiFabs
145  void
146  mf_to_kokkos_buffers(amrex::iMultiFab* lmask,
147  amrex::MultiFab* t_surf,
148  amrex::Vector<amrex::MultiFab*>& lsm_input_ptrs);
149 
150  // Fill AMReX MultiFabs from KOKKOS Views
151  void
152  kokkos_buffers_to_mf (amrex::Vector<amrex::MultiFab*>& lsm_output_ptrs);
153 
154  // Write the rrtmgp fluxes
155  void
157 
158  // Initialize the implementation
159  void
160  initialize_impl ();
161 
162  // Run the implementation
163  void
164  run_impl ();
165 
166  // Finalize the implementation
167  void
168  finalize_impl (amrex::Vector<amrex::MultiFab*>& lsm_output_ptrs);
169 
170  // Wrapper for implementation steps
171  void
172  rad_run_impl (amrex::Vector<amrex::MultiFab*>& lsm_output_ptrs)
173  {
174  if (m_update_rad) {
175  amrex::Print() << "Radiation advancing level " << m_lev << " at (YY-MM-DD SS) " << m_orbital_year << '-'
176  << m_orbital_mon << '-' << m_orbital_day << ' ' << m_orbital_sec << " ...";
177  this->initialize_impl();
178  this->run_impl();
179  this->finalize_impl(lsm_output_ptrs);
180  amrex::Print() << "DONE\n";
181  }
182  }
183 
184  // Get names input varnames for lsm
185  virtual
186  amrex::Vector<std::string>
188  {
189  return m_lsm_input_names;
190  }
191 
192  // Get names output varnames for lsm
193  virtual
194  amrex::Vector<std::string>
196  {
197  return m_lsm_output_names;
198  }
199 
200  // Populate datalog structures
201  void
203 
204  // Write datalog
205  virtual
206  void
207  WriteDataLog (const amrex::Real &time) override;
208 
209 private:
210 
211  // Process interface vars from ERF/AMReX
212  //===================================================================================
213 
214  // Grid level
215  int m_lev;
216 
217  // Step number
218  int m_step;
219 
220  // Current time
222 
223  // Timestep at given level
225 
226  // Geometry at given level
227  amrex::Geometry m_geom;
228 
229  // Boxarray at given level
230  amrex::BoxArray m_ba;
231 
232  // Are we updating radiation?
233  bool m_update_rad = false;
234 
235  // Are we writing out the fluxes?
236  bool m_rad_write_fluxes = false;
237 
238  bool m_first_step = true;
239 
240  // Do we have moisture and cold comps?
241  bool m_moist = false;
242  bool m_ice = false;
243 
244  // Do we have a land surface model?
245  bool m_lsm = false;
246 
247  // List of input parameter names
248  amrex::Vector<std::string> m_lsm_input_names = {"t_sfc" , "sfc_emis" ,
249  "sfc_alb_dir_vis", "sfc_alb_dir_nir",
250  "sfc_alb_dif_vis", "sfc_alb_dif_nir"};
251 
252  // List of output parameter names
253  amrex::Vector<std::string> m_lsm_output_names = {"cos_zenith_angle", "sw_flux_dn", "lw_flux_dn"};
254 
255  // T surf if no LSM is available
257 
258  // Pointer to the CC conserved vars
259  amrex::MultiFab* m_cons_in = nullptr;
260 
261  // Pointer to the radiation source terms
262  amrex::MultiFab* m_qheating_rates = nullptr;
263 
264  // Pointer to the radiation fluxes
265  amrex::MultiFab* m_rad_fluxes = nullptr;
266 
267  // Pointer to the terrain heights
268  amrex::MultiFab* m_z_phys = nullptr;
269 
270  // Pointer to latitude and longitude
271  amrex::MultiFab* m_lat = nullptr;
272  amrex::MultiFab* m_lon = nullptr;
273 
274  // Constant lat/lon if the above MFs are not valid
275  amrex::Real m_lat_cons = 39.809860;
276  amrex::Real m_lon_cons = -98.555183;
277 
278  // Pointer to output data for LSM
279  amrex::MultiFab* m_lsm_fluxes = nullptr;
280  amrex::MultiFab* m_lsm_zenith = nullptr;
281 
282  // Holds output from KOKKOS views used in the datalog
283  amrex::MultiFab datalog_mf;
284 
285  // Path, data file, and coefficient file for K-distribution
286  std::string rrtmgp_file_path = ".";
287  std::string rrtmgp_coeffs_sw = "rrtmgp-data-sw-g224-2018-12-04.nc";
288  std::string rrtmgp_coeffs_lw = "rrtmgp-data-lw-g256-2018-12-04.nc";
289  std::string rrtmgp_cloud_optics_sw = "rrtmgp-cloud-optics-coeffs-sw.nc";
290  std::string rrtmgp_cloud_optics_lw = "rrtmgp-cloud-optics-coeffs-lw.nc";
295 
296  // Active gases
297  int m_ngas = 8;
298  const std::vector<std::string> m_gas_names = {"H2O", "CO2", "O3", "N2O",
299  "CO" , "CH4", "O2", "N2" };
300  const std::vector<amrex::Real> m_mol_weight_gas = {18.01528, 44.00950, 47.9982, 44.0128,
301  28.01010, 16.04246, 31.9980, 28.0134}; // g/mol
302 
303  // Prescribed greenhouse gas surface concentrations in moles / moles air
304  amrex::Real m_co2vmr = 388.717e-6;
305  amrex::Vector<amrex::Real> m_o3vmr;
306  amrex::Real m_n2ovmr = 323.141e-9;
308  amrex::Real m_ch4vmr = 1807.851e-9;
309  amrex::Real m_o2vmr = 0.209448;
311  //amrex::Real m_f11vmr = 768.7644e-12;
312  //amrex::Real m_f12vmr = 531.2820e-12;
313 
316  std::vector<std::string> gas_names_offset;
317 
318  GasConcsK<amrex::Real, layout_t, KokkosDefaultDevice> m_gas_concs;
319 
320  // Process interface vars modeled after EAMXX
321  //===================================================================================
322 
323  // Keep track of number of columns and levels
324  int m_ncol;
325  int m_nlay;
326 
327  // Offsets for MultiFab <-> KOKKOS transfer
328  amrex::Vector<int> m_col_offsets;
329 
330  // Whether we use aerosol forcing in radiation
331  bool m_do_aerosol_rad = true;
332 
333  // Whether we do extra aerosol forcing calls
334  bool m_extra_clnsky_diag = false;
336 
337  // The orbital year, used for zenith angle calculations:
338  // If > 0, use constant orbital year for duration of simulation
339  // If < 0, use year from timestamp for orbital parameters
340  int m_orbital_year = -9999;
341  int m_orbital_mon = -9999;
342  int m_orbital_day = -9999;
343  int m_orbital_sec = -9999;
344 
345  // Orbital parameters, used for zenith angle calculations.
346  // If >= 0, bypass computation based on orbital year and use fixed parameters
347  // If < 0, compute based on orbital year, specified above
348  bool m_fixed_orbital_year = false;
349  amrex::Real m_orbital_eccen = -9999.; // Eccentricity
350  amrex::Real m_orbital_obliq = -9999.; // Obliquity
351  amrex::Real m_orbital_mvelp = -9999.; // Vernal Equinox Mean Longitude of Perihelion
352 
353  // Value for prescribing an invariant solar constant (i.e. total solar irradiance
354  // at TOA). Used for idealized experiments such as RCE. This is only used when a
355  // positive value is supplied.
357 
358  // Fixed solar zenith angle to use for shortwave calculations
359  // This is only used if a positive value is supplied
361 
362  // Dimensions to be read from lookup data
367 
368  // Rad frequency in number of steps
370 
371  // Whether or not to do subcolumn sampling of cloud state for MCICA
372  bool m_do_subcol_sampling = true;
373 
374  // 1d size (1 or nlay)
376 
377  // 1d size (ncol)
392 
393  // 2d size (ncol, nlay)
410 
411  // 2d size (ncol, nlay+1)
435 
436  // 3d size (ncol, nlay+1, nswbands)
441 
442  // 3d size (ncol, nlay+1, nlwbands)
445 
446  // 2d size (ncol, nswbands)
449 
450  // 2d size (ncol, nlwbands)
452 
453  // 3d size (ncol, nlay, n[sw,lw]bands)
458 
459  // 3d size (ncol, nlay, n[sw,lw]bnds)
462 
463  // 3d size (ncol, nlay, n[sw,lw]gpts)
466 };
467 #endif
Kokkos::View< RealT *, KokkosDefaultDevice > real1d_k
Definition: ERF_Kokkos.H:18
Kokkos::View< RealT ***, layout_t, KokkosDefaultDevice > real3d_k
Definition: ERF_Kokkos.H:20
Kokkos::View< RealT **, layout_t, KokkosDefaultDevice > real2d_k
Definition: ERF_Kokkos.H:19
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Definition: ERF_RadiationInterface.H:14
Definition: ERF_Radiation.H:45
real3d_k sw_bnd_flux_dn
Definition: ERF_Radiation.H:438
amrex::MultiFab * m_lon
Definition: ERF_Radiation.H:272
real2d_k lw_flux_up
Definition: ERF_Radiation.H:418
virtual amrex::Vector< std::string > get_lsm_output_varnames() override
Definition: ERF_Radiation.H:195
real3d_k aero_tau_sw
Definition: ERF_Radiation.H:454
int m_o3_size
Definition: ERF_Radiation.H:314
real3d_k sw_bnd_flux_dir
Definition: ERF_Radiation.H:439
virtual void Run(int &level, int &step, amrex::Real &time, const amrex::Real &dt, const amrex::BoxArray &ba, amrex::Geometry &geom, amrex::MultiFab *cons_in, amrex::iMultiFab *lmask, amrex::MultiFab *t_surf, amrex::MultiFab *lsm_fluxes, amrex::MultiFab *lsm_zenith, amrex::Vector< amrex::MultiFab * > &lsm_input_ptrs, amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs, amrex::MultiFab *qheating_rates, amrex::MultiFab *rad_fluxes, amrex::MultiFab *z_phys, amrex::MultiFab *lat_ptr, amrex::MultiFab *lon_ptr) override
Definition: ERF_Radiation.H:90
real2d_k sw_clnsky_flux_dn
Definition: ERF_Radiation.H:427
std::string rrtmgp_coeffs_file_sw
Definition: ERF_Radiation.H:291
std::string rrtmgp_coeffs_sw
Definition: ERF_Radiation.H:287
real2d_k d_tint
Definition: ERF_Radiation.H:412
real2d_k lw_clnclrsky_flux_dn
Definition: ERF_Radiation.H:430
int m_rad_freq_in_steps
Definition: ERF_Radiation.H:369
real2d_k lwp
Definition: ERF_Radiation.H:404
real3d_k cld_tau_lw_gpt
Definition: ERF_Radiation.H:465
real1d_k lw_src
Definition: ERF_Radiation.H:391
real2d_k eff_radius_qi
Definition: ERF_Radiation.H:403
void dealloc_buffers()
Definition: ERF_Radiation.cpp:334
real1d_k m_gas_mol_weights
Definition: ERF_Radiation.H:315
amrex::Real m_lon_cons
Definition: ERF_Radiation.H:276
void initialize_impl()
Definition: ERF_Radiation.cpp:1031
int m_nswbands
Definition: ERF_Radiation.H:363
bool m_do_aerosol_rad
Definition: ERF_Radiation.H:331
std::string rrtmgp_coeffs_lw
Definition: ERF_Radiation.H:288
bool m_moist
Definition: ERF_Radiation.H:241
std::string rrtmgp_cloud_optics_file_lw
Definition: ERF_Radiation.H:294
amrex::MultiFab * m_cons_in
Definition: ERF_Radiation.H:259
amrex::MultiFab datalog_mf
Definition: ERF_Radiation.H:283
Radiation(const int &lev, SolverChoice &sc)
Definition: ERF_Radiation.cpp:20
real2d_k sw_heating
Definition: ERF_Radiation.H:406
bool m_lsm
Definition: ERF_Radiation.H:245
bool m_do_subcol_sampling
Definition: ERF_Radiation.H:372
real3d_k lw_bnd_flux_dn
Definition: ERF_Radiation.H:444
bool m_rad_write_fluxes
Definition: ERF_Radiation.H:236
real3d_k sw_bnd_flux_up
Definition: ERF_Radiation.H:437
real2d_k qv_lay
Definition: ERF_Radiation.H:398
std::string rrtmgp_cloud_optics_file_sw
Definition: ERF_Radiation.H:293
amrex::Vector< std::string > m_lsm_output_names
Definition: ERF_Radiation.H:253
void set_grids(int &level, int &step, amrex::Real &time, const amrex::Real &dt, const amrex::BoxArray &ba, amrex::Geometry &geom, amrex::MultiFab *cons_in, amrex::iMultiFab *lmask, amrex::MultiFab *t_surf, amrex::MultiFab *lsm_fluxes, amrex::MultiFab *lsm_zenith, amrex::Vector< amrex::MultiFab * > &lsm_input_ptrs, amrex::MultiFab *qheating_rates, amrex::MultiFab *rad_fluxes, amrex::MultiFab *z_phys, amrex::MultiFab *lat, amrex::MultiFab *lon)
Definition: ERF_Radiation.cpp:134
void run_impl()
Definition: ERF_Radiation.cpp:1042
real2d_k sw_clnclrsky_flux_dn_dir
Definition: ERF_Radiation.H:422
real1d_k sfc_flux_dif_vis
Definition: ERF_Radiation.H:385
int m_step
Definition: ERF_Radiation.H:218
real2d_k lw_clnclrsky_flux_up
Definition: ERF_Radiation.H:429
real3d_k cld_tau_sw_gpt
Definition: ERF_Radiation.H:464
real1d_k lat
Definition: ERF_Radiation.H:387
void populateDatalogMF()
Definition: ERF_Radiation.cpp:790
real2d_k qi_lay
Definition: ERF_Radiation.H:400
real2d_k sw_clrsky_flux_up
Definition: ERF_Radiation.H:423
amrex::Real m_orbital_mvelp
Definition: ERF_Radiation.H:351
real2d_k t_lev
Definition: ERF_Radiation.H:414
amrex::Vector< amrex::Real > m_o3vmr
Definition: ERF_Radiation.H:305
real1d_k o3_lay
Definition: ERF_Radiation.H:375
std::string rrtmgp_cloud_optics_sw
Definition: ERF_Radiation.H:289
real1d_k sfc_alb_dif_vis
Definition: ERF_Radiation.H:381
int m_nlwgpts
Definition: ERF_Radiation.H:366
real1d_k mu0
Definition: ERF_Radiation.H:378
real2d_k cldfrac_tot
Definition: ERF_Radiation.H:401
real3d_k cld_tau_lw_bnd
Definition: ERF_Radiation.H:461
real1d_k sfc_flux_dir_nir
Definition: ERF_Radiation.H:384
amrex::MultiFab * m_z_phys
Definition: ERF_Radiation.H:268
real2d_k sw_clnclrsky_flux_up
Definition: ERF_Radiation.H:420
real3d_k aero_g_sw
Definition: ERF_Radiation.H:456
real2d_k sw_flux_up
Definition: ERF_Radiation.H:415
std::string rrtmgp_cloud_optics_lw
Definition: ERF_Radiation.H:290
amrex::MultiFab * m_rad_fluxes
Definition: ERF_Radiation.H:265
real3d_k sw_bnd_flux_dif
Definition: ERF_Radiation.H:440
real2d_k sfc_alb_dif
Definition: ERF_Radiation.H:448
amrex::Real m_co2vmr
Definition: ERF_Radiation.H:304
bool m_extra_clnsky_diag
Definition: ERF_Radiation.H:334
real2d_k emis_sfc
Definition: ERF_Radiation.H:451
real1d_k sfc_alb_dif_nir
Definition: ERF_Radiation.H:382
amrex::MultiFab * m_lsm_fluxes
Definition: ERF_Radiation.H:279
virtual void Init(const amrex::Geometry &geom, const amrex::BoxArray &ba, amrex::MultiFab *cons_in) override
Definition: ERF_Radiation.H:62
real2d_k lw_clnsky_flux_dn
Definition: ERF_Radiation.H:434
amrex::BoxArray m_ba
Definition: ERF_Radiation.H:230
real2d_k p_lay
Definition: ERF_Radiation.H:395
real2d_k r_lay
Definition: ERF_Radiation.H:394
int m_ncol
Definition: ERF_Radiation.H:324
real2d_k sw_flux_dn_dir
Definition: ERF_Radiation.H:417
const std::vector< amrex::Real > m_mol_weight_gas
Definition: ERF_Radiation.H:300
amrex::Real m_lat_cons
Definition: ERF_Radiation.H:275
amrex::Real m_fixed_total_solar_irradiance
Definition: ERF_Radiation.H:356
void kokkos_buffers_to_mf(amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
Definition: ERF_Radiation.cpp:638
std::vector< std::string > gas_names_offset
Definition: ERF_Radiation.H:316
real2d_k sw_clrsky_flux_dn
Definition: ERF_Radiation.H:424
amrex::Real m_dt
Definition: ERF_Radiation.H:224
virtual amrex::Vector< std::string > get_lsm_input_varnames() override
Definition: ERF_Radiation.H:187
int m_orbital_mon
Definition: ERF_Radiation.H:341
real3d_k aero_ssa_sw
Definition: ERF_Radiation.H:455
real2d_k sw_clnsky_flux_dn_dir
Definition: ERF_Radiation.H:428
amrex::Real m_o2vmr
Definition: ERF_Radiation.H:309
void rad_run_impl(amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
Definition: ERF_Radiation.H:172
GasConcsK< amrex::Real, layout_t, KokkosDefaultDevice > m_gas_concs
Definition: ERF_Radiation.H:318
amrex::MultiFab * m_qheating_rates
Definition: ERF_Radiation.H:262
real2d_k qc_lay
Definition: ERF_Radiation.H:399
virtual void WriteDataLog(const amrex::Real &time) override
Definition: ERF_Radiation.cpp:896
real2d_k lw_clrsky_flux_up
Definition: ERF_Radiation.H:431
amrex::Real m_n2ovmr
Definition: ERF_Radiation.H:306
amrex::Real m_rad_t_sfc
Definition: ERF_Radiation.H:256
std::string rrtmgp_file_path
Definition: ERF_Radiation.H:286
real2d_k sw_flux_dn
Definition: ERF_Radiation.H:416
amrex::MultiFab * m_lsm_zenith
Definition: ERF_Radiation.H:280
void mf_to_kokkos_buffers(amrex::iMultiFab *lmask, amrex::MultiFab *t_surf, amrex::Vector< amrex::MultiFab * > &lsm_input_ptrs)
Definition: ERF_Radiation.cpp:441
int m_nlwbands
Definition: ERF_Radiation.H:364
real2d_k lw_clrsky_flux_dn
Definition: ERF_Radiation.H:432
amrex::Real m_orbital_eccen
Definition: ERF_Radiation.H:349
void finalize_impl(amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
Definition: ERF_Radiation.cpp:1301
amrex::Real m_covmr
Definition: ERF_Radiation.H:307
int m_orbital_sec
Definition: ERF_Radiation.H:343
real2d_k z_del
Definition: ERF_Radiation.H:397
real3d_k aero_tau_lw
Definition: ERF_Radiation.H:457
real2d_k lw_flux_dn
Definition: ERF_Radiation.H:419
~Radiation()
Definition: ERF_Radiation.H:55
real2d_k sw_clrsky_flux_dn_dir
Definition: ERF_Radiation.H:425
void write_rrtmgp_fluxes()
Definition: ERF_Radiation.cpp:750
std::string rrtmgp_coeffs_file_lw
Definition: ERF_Radiation.H:292
int m_lev
Definition: ERF_Radiation.H:215
real1d_k sfc_flux_dif_nir
Definition: ERF_Radiation.H:386
int m_ngas
Definition: ERF_Radiation.H:297
real2d_k lw_heating
Definition: ERF_Radiation.H:407
real1d_k sfc_alb_dir_nir
Definition: ERF_Radiation.H:380
real2d_k lw_clnsky_flux_up
Definition: ERF_Radiation.H:433
amrex::Real m_ch4vmr
Definition: ERF_Radiation.H:308
bool m_update_rad
Definition: ERF_Radiation.H:233
real2d_k sfc_alb_dir
Definition: ERF_Radiation.H:447
bool m_ice
Definition: ERF_Radiation.H:242
amrex::Real m_fixed_solar_zenith_angle
Definition: ERF_Radiation.H:360
real1d_k lon
Definition: ERF_Radiation.H:388
int m_orbital_day
Definition: ERF_Radiation.H:342
real1d_k sfc_emis
Definition: ERF_Radiation.H:389
amrex::Real m_orbital_obliq
Definition: ERF_Radiation.H:350
amrex::Real m_n2vmr
Definition: ERF_Radiation.H:310
real3d_k cld_tau_sw_bnd
Definition: ERF_Radiation.H:460
bool m_fixed_orbital_year
Definition: ERF_Radiation.H:348
real2d_k sw_clnsky_flux_up
Definition: ERF_Radiation.H:426
int m_nlay
Definition: ERF_Radiation.H:325
real3d_k lw_bnd_flux_up
Definition: ERF_Radiation.H:443
int m_nswgpts
Definition: ERF_Radiation.H:365
bool m_first_step
Definition: ERF_Radiation.H:238
real1d_k sfc_alb_dir_vis
Definition: ERF_Radiation.H:379
const std::vector< std::string > m_gas_names
Definition: ERF_Radiation.H:298
real2d_k sw_clrsky_heating
Definition: ERF_Radiation.H:408
bool m_extra_clnclrsky_diag
Definition: ERF_Radiation.H:335
real2d_k t_lay
Definition: ERF_Radiation.H:396
real1d_k sfc_flux_dir_vis
Definition: ERF_Radiation.H:383
real2d_k iwp
Definition: ERF_Radiation.H:405
amrex::Geometry m_geom
Definition: ERF_Radiation.H:227
amrex::MultiFab * m_lat
Definition: ERF_Radiation.H:271
real2d_k p_lev
Definition: ERF_Radiation.H:413
real2d_k lw_clrsky_heating
Definition: ERF_Radiation.H:409
real2d_k eff_radius_qc
Definition: ERF_Radiation.H:402
amrex::Vector< int > m_col_offsets
Definition: ERF_Radiation.H:328
real1d_k t_sfc
Definition: ERF_Radiation.H:390
void alloc_buffers()
Definition: ERF_Radiation.cpp:205
int m_orbital_year
Definition: ERF_Radiation.H:340
real2d_k sw_clnclrsky_flux_dn
Definition: ERF_Radiation.H:421
amrex::Vector< std::string > m_lsm_input_names
Definition: ERF_Radiation.H:248
amrex::Real m_time
Definition: ERF_Radiation.H:221
Definition: ERF_DataStruct.H:128