ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
Radiation Class Reference

#include <ERF_Radiation.H>

Collaboration diagram for Radiation:

Public Member Functions

 Radiation (const int &lev, SolverChoice &sc)
 
 ~Radiation ()=default
 
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::MultiFab *lsm_fluxes, amrex::MultiFab *lsm_zenith, amrex::MultiFab *qheating_rates, amrex::MultiFab *z_phys, amrex::MultiFab *lat, amrex::MultiFab *lon)
 
void alloc_buffers ()
 
void dealloc_buffers ()
 
void mf_to_yakl_buffers ()
 
void yakl_buffers_to_mf ()
 
void write_rrtmgp_fluxes ()
 
void initialize_impl ()
 
void run_impl ()
 
void finalize_impl ()
 
void rad_run_impl ()
 

Private Attributes

int m_lev
 
int m_step
 
amrex::Real m_time
 
amrex::Real m_dt
 
amrex::Geometry m_geom
 
amrex::BoxArray m_ba
 
bool m_update_rad = false
 
bool m_rad_write_fluxes = false
 
bool m_moist = false
 
bool m_ice = false
 
amrex::MultiFab * m_cons_in = nullptr
 
amrex::MultiFab * m_qheating_rates = nullptr
 
amrex::MultiFab * m_z_phys = nullptr
 
amrex::MultiFab * m_lat = nullptr
 
amrex::MultiFab * m_lon = nullptr
 
amrex::Real m_lat_cons = 39.809860
 
amrex::Real m_lon_cons = -98.555183
 
amrex::MultiFab * m_lsm_fluxes = nullptr
 
amrex::MultiFab * m_lsm_zenith = nullptr
 
std::string rrtmgp_file_path = "."
 
std::string rrtmgp_coeffs_sw = "rrtmgp-data-sw-g224-2018-12-04.nc"
 
std::string rrtmgp_coeffs_lw = "rrtmgp-data-lw-g224-2018-12-04.nc"
 
std::string rrtmgp_cloud_optics_sw = "rrtmgp-cloud-optics-coeffs-sw.nc"
 
std::string rrtmgp_cloud_optics_lw = "rrtmgp-cloud-optics-coeffs-lw.nc"
 
std::string rrtmgp_coeffs_file_sw
 
std::string rrtmgp_coeffs_file_lw
 
std::string rrtmgp_cloud_optics_file_sw
 
std::string rrtmgp_cloud_optics_file_lw
 
int m_ngas = 8
 
const std::vector< std::string > m_gas_names
 
const std::vector< amrex::Real > m_mol_weight_gas
 
amrex::Real m_co2vmr = 388.717e-6
 
amrex::Vector< amrex::Real > m_o3vmr
 
amrex::Real m_n2ovmr = 323.141e-9
 
amrex::Real m_covmr = 1.0e-7
 
amrex::Real m_ch4vmr = 1807.851e-9
 
amrex::Real m_o2vmr = 0.209448
 
amrex::Real m_n2vmr = 0.7906
 
int m_o3_size
 
real1d m_gas_mol_weights
 
string1dv gas_names_yakl_offset
 
GasConcs m_gas_concs
 
int m_ncol
 
int m_nlay
 
amrex::Vector< int > m_col_offsets
 
bool m_do_aerosol_rad = true
 
bool m_extra_clnsky_diag = false
 
bool m_extra_clnclrsky_diag = false
 
int m_orbital_year = -9999
 
int m_orbital_mon = -9999
 
int m_orbital_day = -9999
 
int m_orbital_sec = -9999
 
bool m_fixed_orbital_year = false
 
amrex::Real m_orbital_eccen = -9999.
 
amrex::Real m_orbital_obliq = -9999.
 
amrex::Real m_orbital_mvelp = -9999.
 
amrex::Real m_fixed_total_solar_irradiance = -9999.
 
amrex::Real m_fixed_solar_zenith_angle = -9999.
 
int m_nswbands = 14
 
int m_nlwbands = 16
 
int m_nswgpts = 112
 
int m_nlwgpts = 128
 
int m_rad_freq_in_steps = 1
 
bool m_do_subcol_sampling = true
 
real1d o3_lay
 
real1d cosine_zenith
 
real1d mu0
 
real1d sfc_alb_dir_vis
 
real1d sfc_alb_dir_nir
 
real1d sfc_alb_dif_vis
 
real1d sfc_alb_dif_nir
 
real1d sfc_flux_dir_vis
 
real1d sfc_flux_dir_nir
 
real1d sfc_flux_dif_vis
 
real1d sfc_flux_dif_nir
 
real1d lat
 
real1d lon
 
real2d r_lay
 
real2d p_lay
 
real2d t_lay
 
real2d z_del
 
real2d p_del
 
real2d qv_lay
 
real2d qc_lay
 
real2d qi_lay
 
real2d cldfrac_tot
 
real2d eff_radius_qc
 
real2d eff_radius_qi
 
real2d tmp2d
 
real2d lwp
 
real2d iwp
 
real2d sw_heating
 
real2d lw_heating
 
real2d d_tint
 
real2d p_lev
 
real2d t_lev
 
real2d sw_flux_up
 
real2d sw_flux_dn
 
real2d sw_flux_dn_dir
 
real2d lw_flux_up
 
real2d lw_flux_dn
 
real2d sw_clnclrsky_flux_up
 
real2d sw_clnclrsky_flux_dn
 
real2d sw_clnclrsky_flux_dn_dir
 
real2d sw_clrsky_flux_up
 
real2d sw_clrsky_flux_dn
 
real2d sw_clrsky_flux_dn_dir
 
real2d sw_clnsky_flux_up
 
real2d sw_clnsky_flux_dn
 
real2d sw_clnsky_flux_dn_dir
 
real2d lw_clnclrsky_flux_up
 
real2d lw_clnclrsky_flux_dn
 
real2d lw_clrsky_flux_up
 
real2d lw_clrsky_flux_dn
 
real2d lw_clnsky_flux_up
 
real2d lw_clnsky_flux_dn
 
real3d sw_bnd_flux_up
 
real3d sw_bnd_flux_dn
 
real3d sw_bnd_flux_dir
 
real3d sw_bnd_flux_dif
 
real3d lw_bnd_flux_up
 
real3d lw_bnd_flux_dn
 
real2d sfc_alb_dir
 
real2d sfc_alb_dif
 
real3d aero_tau_sw
 
real3d aero_ssa_sw
 
real3d aero_g_sw
 
real3d aero_tau_lw
 
real3d cld_tau_sw_bnd
 
real3d cld_tau_lw_bnd
 
real3d cld_tau_sw_gpt
 
real3d cld_tau_lw_gpt
 

Constructor & Destructor Documentation

◆ Radiation()

Radiation::Radiation ( const int &  lev,
SolverChoice sc 
)
25 {
26  // Initialize YAKL
27  if (!yakl::isInitialized()) { yakl::init(); }
28 
29  // Check if we have a valid moisture model
30  if (sc.moisture_type != MoistureType::None) { m_moist = true; }
31 
32  // Check if we have a moisture model with ice
33  if (sc.moisture_type == MoistureType::SAM) { m_ice = true; }
34 
35  // Construct parser object for following reads
36  ParmParse pp("erf");
37 
38  // Radiation timestep, as a number of atm steps
39  pp.query("rad_freq_in_steps", m_rad_freq_in_steps);
40 
41  // Flag to write fluxes to plt file
42  pp.query("rad_write_fluxes", m_rad_write_fluxes);
43 
44  // Do MCICA subcolumn sampling
45  pp.query("rad_do_subcol_sampling", m_do_subcol_sampling);
46 
47  // Determine orbital year. If orbital_year is negative, use current year
48  // from timestamp for orbital year; if positive, use provided orbital year
49  // for duration of simulation.
50  m_fixed_orbital_year = pp.query("rad_orbital_year", m_orbital_year);
51 
52  // Get orbital parameters from inputs file
53  pp.query("rad_orbital_eccentricity", m_orbital_eccen);
54  pp.query("rad_orbital_obliquity" , m_orbital_obliq);
55  pp.query("rad_orbital_mvelp" , m_orbital_mvelp);
56 
57  // Get a constant lat/lon for idealized simulations
58  pp.query("rad_cons_lat", m_lat_cons);
59  pp.query("rad_cons_lon", m_lon_cons);
60 
61  // Value for prescribing an invariant solar constant (i.e. total solar irradiance at
62  // TOA). Used for idealized experiments such as RCE. Disabled when value is less than 0.
63  pp.query("fixed_total_solar_irradiance", m_fixed_total_solar_irradiance);
64 
65  // Determine whether or not we are using a fixed solar zenith angle (positive value)
66  pp.query("fixed_solar_zenith_angle", m_fixed_solar_zenith_angle);
67 
68  // Get prescribed surface values of greenhouse gases
69  pp.query("co2vmr", m_co2vmr);
70  pp.queryarr("o3vmr" , m_o3vmr );
71  pp.query("n2ovmr", m_n2ovmr);
72  pp.query("covmr" , m_covmr );
73  pp.query("ch4vmr", m_ch4vmr);
74  pp.query("o2vmr" , m_o2vmr );
75  pp.query("n2vmr" , m_n2vmr );
76 
77  // Required aerosol optical properties from SPA
78  pp.query("rad_do_aerosol", m_do_aerosol_rad);
79 
80  // Whether we do extra clean/clear sky calculations
81  pp.query("rad_extra_clnclrsky_diag", m_extra_clnclrsky_diag);
82  pp.query("rad_extra_clnsky_diag" , m_extra_clnsky_diag);
83 
84  // Parse the band and gauss pt sizes
85  pp.query("nswbands", m_nswbands);
86  pp.query("nlwbands", m_nlwbands);
87  pp.query("nswgpts" , m_nswgpts );
88  pp.query("nlwgpts" , m_nlwgpts );
89 
90  // Parse path and file names
91  pp.query("rrtmgp_file_path" , rrtmgp_file_path);
92  pp.query("rrtmgp_coeffs_sw" , rrtmgp_coeffs_sw );
93  pp.query("rrtmgp_coeffs_lw" , rrtmgp_coeffs_lw );
94  pp.query("rrtmgp_cloud_optics_sw", rrtmgp_cloud_optics_sw);
95  pp.query("rrtmgp_cloud_optics_lw", rrtmgp_cloud_optics_lw);
96 
97  // Append file names to path
102 
103  // Output for user
104  if (lev == 0) {
105  Print() << "Radiation interface constructed:\n";
106  Print() << "========================================================\n";
107  Print() << "Coeff SW file: " << rrtmgp_coeffs_file_sw << "\n";
108  Print() << "Coeff LW file: " << rrtmgp_coeffs_file_lw << "\n";
109  Print() << "Cloud SW file: " << rrtmgp_cloud_optics_file_sw << "\n";
110  Print() << "Cloud LW file: " << rrtmgp_cloud_optics_file_lw << "\n";
111  Print() << "========================================================\n";
112  }
113 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:219
std::string rrtmgp_coeffs_file_sw
Definition: ERF_Radiation.H:173
std::string rrtmgp_coeffs_sw
Definition: ERF_Radiation.H:169
int m_rad_freq_in_steps
Definition: ERF_Radiation.H:251
amrex::Real m_lon_cons
Definition: ERF_Radiation.H:161
int m_nswbands
Definition: ERF_Radiation.H:245
bool m_do_aerosol_rad
Definition: ERF_Radiation.H:212
std::string rrtmgp_coeffs_lw
Definition: ERF_Radiation.H:170
bool m_moist
Definition: ERF_Radiation.H:143
std::string rrtmgp_cloud_optics_file_lw
Definition: ERF_Radiation.H:176
bool m_do_subcol_sampling
Definition: ERF_Radiation.H:254
bool m_rad_write_fluxes
Definition: ERF_Radiation.H:140
std::string rrtmgp_cloud_optics_file_sw
Definition: ERF_Radiation.H:175
amrex::Real m_orbital_mvelp
Definition: ERF_Radiation.H:232
amrex::Vector< amrex::Real > m_o3vmr
Definition: ERF_Radiation.H:187
std::string rrtmgp_cloud_optics_sw
Definition: ERF_Radiation.H:171
int m_nlwgpts
Definition: ERF_Radiation.H:248
std::string rrtmgp_cloud_optics_lw
Definition: ERF_Radiation.H:172
amrex::Real m_co2vmr
Definition: ERF_Radiation.H:186
bool m_extra_clnsky_diag
Definition: ERF_Radiation.H:215
amrex::Real m_lat_cons
Definition: ERF_Radiation.H:160
amrex::Real m_fixed_total_solar_irradiance
Definition: ERF_Radiation.H:237
amrex::Real m_o2vmr
Definition: ERF_Radiation.H:191
amrex::Real m_n2ovmr
Definition: ERF_Radiation.H:188
std::string rrtmgp_file_path
Definition: ERF_Radiation.H:168
int m_nlwbands
Definition: ERF_Radiation.H:246
amrex::Real m_orbital_eccen
Definition: ERF_Radiation.H:230
amrex::Real m_covmr
Definition: ERF_Radiation.H:189
std::string rrtmgp_coeffs_file_lw
Definition: ERF_Radiation.H:174
amrex::Real m_ch4vmr
Definition: ERF_Radiation.H:190
bool m_ice
Definition: ERF_Radiation.H:144
amrex::Real m_fixed_solar_zenith_angle
Definition: ERF_Radiation.H:241
amrex::Real m_orbital_obliq
Definition: ERF_Radiation.H:231
amrex::Real m_n2vmr
Definition: ERF_Radiation.H:192
bool m_fixed_orbital_year
Definition: ERF_Radiation.H:229
int m_nswgpts
Definition: ERF_Radiation.H:247
bool m_extra_clnclrsky_diag
Definition: ERF_Radiation.H:216
int m_orbital_year
Definition: ERF_Radiation.H:221
MoistureType moisture_type
Definition: ERF_DataStruct.H:677
Here is the call graph for this function:

◆ ~Radiation()

Radiation::~Radiation ( )
default

Member Function Documentation

◆ alloc_buffers()

void Radiation::alloc_buffers ( )
187 {
188  // 1d size (m_ngas)
189  m_gas_mol_weights = real1d("m_gas_mol_weights", m_ngas);
190  realHost1d m_gas_mol_weights_h("m_gas_mol_weights_h", m_ngas);
191  gas_names_yakl_offset.clear();
192  parallel_for(m_ngas, YAKL_LAMBDA (int igas)
193  {
194  m_gas_mol_weights_h(igas) = m_mol_weight_gas[igas-1];
195  gas_names_yakl_offset.push_back(m_gas_names[igas-1]);
196  });
197  m_gas_mol_weights_h.deep_copy_to(m_gas_mol_weights);
198 
199  // 1d size (1 or nlay)
200  m_o3_size = m_o3vmr.size();
201  AMREX_ASSERT_WITH_MESSAGE(((m_o3_size==1) || (m_o3_size==m_nlay)), "O3 VMR array must be length 1 or nlay");
202  o3_lay = real1d("o3_lay", m_o3_size);
203  realHost1d o3_lay_h("o3_lay_h", m_o3_size);
204  parallel_for(m_o3_size, YAKL_LAMBDA (int io3)
205  {
206  o3_lay_h(io3) = m_o3vmr[io3-1];
207  });
208  o3_lay_h.deep_copy_to(o3_lay);
209 
210  // 1d size (ncol)
211  cosine_zenith = real1d("cosine_zenith" , m_ncol);
212  mu0 = real1d("mu0" , m_ncol);
213  sfc_alb_dir_vis = real1d("sfc_alb_dir_vis" , m_ncol);
214  sfc_alb_dir_nir = real1d("sfc_alb_dir_nir" , m_ncol);
215  sfc_alb_dif_vis = real1d("sfc_alb_dif_vis" , m_ncol);
216  sfc_alb_dif_nir = real1d("sfc_alb_dif_nir" , m_ncol);
217  sfc_flux_dir_vis = real1d("sfc_flux_dir_vis", m_ncol);
218  sfc_flux_dir_nir = real1d("sfc_flux_dir_nir", m_ncol);
219  sfc_flux_dif_vis = real1d("sfc_flux_dif_vis", m_ncol);
220  sfc_flux_dif_nir = real1d("sfc_flux_dif_nir", m_ncol);
221  lat = real1d("lat" , m_ncol);
222  lon = real1d("lon" , m_ncol);;
223 
224  // 2d size (ncol, nlay)
225  r_lay = real2d("r_lay" , m_ncol, m_nlay);
226  p_lay = real2d("p_lay" , m_ncol, m_nlay);
227  t_lay = real2d("t_lay" , m_ncol, m_nlay);
228  z_del = real2d("z_del" , m_ncol, m_nlay);
229  p_del = real2d("p_del" , m_ncol, m_nlay);
230  qv_lay = real2d("qv" , m_ncol, m_nlay);
231  qc_lay = real2d("qc" , m_ncol, m_nlay);
232  qi_lay = real2d("qi" , m_ncol, m_nlay);
233  cldfrac_tot = real2d("cldfrac_tot" , m_ncol, m_nlay);
234  eff_radius_qc = real2d("eff_radius_qc", m_ncol, m_nlay);
235  eff_radius_qi = real2d("eff_radius_qi", m_ncol, m_nlay);
236  tmp2d = real2d("tmp2d" , m_ncol, m_nlay);
237  lwp = real2d("lwp" , m_ncol, m_nlay);
238  iwp = real2d("iwp" , m_ncol, m_nlay);
239  sw_heating = real2d("sw_heating" , m_ncol, m_nlay);
240  lw_heating = real2d("lw_heating" , m_ncol, m_nlay);
241 
242  // 2d size (ncol, nlay+1)
243  d_tint = real2d("d_tint" , m_ncol, m_nlay+1);
244  p_lev = real2d("p_lev" , m_ncol, m_nlay+1);
245  t_lev = real2d("t_lev" , m_ncol, m_nlay+1);
246  sw_flux_up = real2d("sw_flux_up" , m_ncol, m_nlay+1);
247  sw_flux_dn = real2d("sw_flux_dn" , m_ncol, m_nlay+1);
248  sw_flux_dn_dir = real2d("sw_flux_dn_dir" , m_ncol, m_nlay+1);
249  lw_flux_up = real2d("sw_flux_up" , m_ncol, m_nlay+1);
250  lw_flux_dn = real2d("sw_flux_dn" , m_ncol, m_nlay+1);
251  sw_clnclrsky_flux_up = real2d("sw_clnclrsky_flux_up" , m_ncol, m_nlay+1);
252  sw_clnclrsky_flux_dn = real2d("sw_clnclrsky_flux_dn" , m_ncol, m_nlay+1);
253  sw_clnclrsky_flux_dn_dir = real2d("sw_clnclrsky_flux_dn_dir", m_ncol, m_nlay+1);
254  sw_clrsky_flux_up = real2d("sw_clrsky_flux_up" , m_ncol, m_nlay+1);
255  sw_clrsky_flux_dn = real2d("sw_clrsky_flux_dn" , m_ncol, m_nlay+1);
256  sw_clrsky_flux_dn_dir = real2d("sw_clrsky_flux_dn_dir" , m_ncol, m_nlay+1);
257  sw_clnsky_flux_up = real2d("sw_clnsky_flux_up" , m_ncol, m_nlay+1);
258  sw_clnsky_flux_dn = real2d("sw_clnsky_flux_dn" , m_ncol, m_nlay+1);
259  sw_clnsky_flux_dn_dir = real2d("sw_clnsky_flux_dn_dir" , m_ncol, m_nlay+1);
260  lw_clnclrsky_flux_up = real2d("lw_clnclrsky_flux_up" , m_ncol, m_nlay+1);
261  lw_clnclrsky_flux_dn = real2d("lw_clnclrsky_flux_dn" , m_ncol, m_nlay+1);
262  lw_clrsky_flux_up = real2d("lw_clrsky_flux_up" , m_ncol, m_nlay+1);
263  lw_clrsky_flux_dn = real2d("lw_clrsky_flux_dn" , m_ncol, m_nlay+1);
264  lw_clnsky_flux_up = real2d("lw_clnsky_flux_up" , m_ncol, m_nlay+1);
265  lw_clnsky_flux_dn = real2d("lw_clnsky_flux_dn" , m_ncol, m_nlay+1);
266 
267  // 3d size (ncol, nlay+1, nswbands)
268  sw_bnd_flux_up = real3d("sw_bnd_flux_up" , m_ncol, m_nlay+1, m_nswbands);
269  sw_bnd_flux_dn = real3d("sw_bnd_flux_dn" , m_ncol, m_nlay+1, m_nswbands);
270  sw_bnd_flux_dir = real3d("sw_bnd_flux_dir", m_ncol, m_nlay+1, m_nswbands);
271  sw_bnd_flux_dif = real3d("sw_bnd_flux_dif", m_ncol, m_nlay+1, m_nswbands);
272 
273  // 3d size (ncol, nlay+1, nlwbands)
274  lw_bnd_flux_up = real3d("lw_bnd_flux_up" , m_ncol, m_nlay+1, m_nlwbands);
275  lw_bnd_flux_dn = real3d("lw_bnd_flux_up" , m_ncol, m_nlay+1, m_nlwbands);
276 
277  // 2d size (ncol, nswbands)
278  sfc_alb_dir = real2d("sfc_alb_dir", m_ncol, m_nswbands);
279  sfc_alb_dif = real2d("sfc_alb_dif", m_ncol, m_nswbands);
280 
281  // 3d size (ncol, nlay, n[sw,lw]bands)
282  aero_tau_sw = real3d("aero_tau_sw", m_ncol, m_nlay, m_nswbands);
283  aero_ssa_sw = real3d("aero_ssa_sw", m_ncol, m_nlay, m_nswbands);
284  aero_g_sw = real3d("aero_g_sw" , m_ncol, m_nlay, m_nswbands);
285  aero_tau_lw = real3d("aero_tau_lw", m_ncol, m_nlay, m_nlwbands);
286 
287  // 3d size (ncol, nlay, n[sw,lw]bnds)
288  cld_tau_sw_bnd = real3d("cld_tau_sw_bnd", m_ncol, m_nlay, m_nswbands);
289  cld_tau_lw_bnd = real3d("cld_tau_lw_bnd", m_ncol, m_nlay, m_nlwbands);
290 
291  // 3d size (ncol, nlay, n[sw,lw]gpts)
292  cld_tau_sw_gpt = real3d("cld_tau_sw_gpt", m_ncol, m_nlay, m_nswgpts);
293  cld_tau_lw_gpt = real3d("cld_tau_lw_gpt", m_ncol, m_nlay, m_nlwgpts);
294 }
real2d z_del
Definition: ERF_Radiation.H:277
real2d qv_lay
Definition: ERF_Radiation.H:279
real1d sfc_flux_dif_nir
Definition: ERF_Radiation.H:269
real2d sw_clrsky_flux_up
Definition: ERF_Radiation.H:303
real1d sfc_alb_dif_vis
Definition: ERF_Radiation.H:264
int m_o3_size
Definition: ERF_Radiation.H:196
real2d sw_clrsky_flux_dn_dir
Definition: ERF_Radiation.H:305
real2d lw_flux_dn
Definition: ERF_Radiation.H:299
real1d sfc_flux_dir_vis
Definition: ERF_Radiation.H:266
real2d lw_clrsky_flux_dn
Definition: ERF_Radiation.H:312
real2d sw_clnclrsky_flux_dn
Definition: ERF_Radiation.H:301
real2d sw_clnsky_flux_dn
Definition: ERF_Radiation.H:307
real2d sw_clnsky_flux_up
Definition: ERF_Radiation.H:306
real3d cld_tau_lw_bnd
Definition: ERF_Radiation.H:338
real3d aero_ssa_sw
Definition: ERF_Radiation.H:332
real2d qi_lay
Definition: ERF_Radiation.H:281
real2d sw_flux_dn_dir
Definition: ERF_Radiation.H:297
real1d sfc_flux_dif_vis
Definition: ERF_Radiation.H:268
real2d lw_clrsky_flux_up
Definition: ERF_Radiation.H:311
real1d lat
Definition: ERF_Radiation.H:270
real2d iwp
Definition: ERF_Radiation.H:287
real2d t_lev
Definition: ERF_Radiation.H:294
real2d sw_clnclrsky_flux_dn_dir
Definition: ERF_Radiation.H:302
real2d lw_clnsky_flux_dn
Definition: ERF_Radiation.H:314
real2d eff_radius_qc
Definition: ERF_Radiation.H:283
real3d sw_bnd_flux_dif
Definition: ERF_Radiation.H:320
real2d sw_heating
Definition: ERF_Radiation.H:288
real2d sfc_alb_dir
Definition: ERF_Radiation.H:327
real3d aero_tau_sw
Definition: ERF_Radiation.H:331
real2d r_lay
Definition: ERF_Radiation.H:274
real1d o3_lay
Definition: ERF_Radiation.H:257
real2d sfc_alb_dif
Definition: ERF_Radiation.H:328
int m_ncol
Definition: ERF_Radiation.H:205
real2d sw_clnclrsky_flux_up
Definition: ERF_Radiation.H:300
const std::vector< amrex::Real > m_mol_weight_gas
Definition: ERF_Radiation.H:182
real3d cld_tau_lw_gpt
Definition: ERF_Radiation.H:342
real1d sfc_flux_dir_nir
Definition: ERF_Radiation.H:267
real1d lon
Definition: ERF_Radiation.H:271
real3d lw_bnd_flux_up
Definition: ERF_Radiation.H:323
real3d cld_tau_sw_gpt
Definition: ERF_Radiation.H:341
real3d lw_bnd_flux_dn
Definition: ERF_Radiation.H:324
real2d lw_flux_up
Definition: ERF_Radiation.H:298
real1d mu0
Definition: ERF_Radiation.H:261
real3d sw_bnd_flux_up
Definition: ERF_Radiation.H:317
real2d cldfrac_tot
Definition: ERF_Radiation.H:282
real3d aero_tau_lw
Definition: ERF_Radiation.H:334
real3d sw_bnd_flux_dn
Definition: ERF_Radiation.H:318
real1d cosine_zenith
Definition: ERF_Radiation.H:260
real2d p_lay
Definition: ERF_Radiation.H:275
real2d lw_clnclrsky_flux_up
Definition: ERF_Radiation.H:309
real2d lw_heating
Definition: ERF_Radiation.H:289
real2d lw_clnclrsky_flux_dn
Definition: ERF_Radiation.H:310
real3d aero_g_sw
Definition: ERF_Radiation.H:333
real1d sfc_alb_dif_nir
Definition: ERF_Radiation.H:265
real1d sfc_alb_dir_vis
Definition: ERF_Radiation.H:262
int m_ngas
Definition: ERF_Radiation.H:179
real2d p_del
Definition: ERF_Radiation.H:278
real2d t_lay
Definition: ERF_Radiation.H:276
real2d sw_flux_up
Definition: ERF_Radiation.H:295
real2d p_lev
Definition: ERF_Radiation.H:293
real2d lw_clnsky_flux_up
Definition: ERF_Radiation.H:313
real1d sfc_alb_dir_nir
Definition: ERF_Radiation.H:263
int m_nlay
Definition: ERF_Radiation.H:206
real2d lwp
Definition: ERF_Radiation.H:286
const std::vector< std::string > m_gas_names
Definition: ERF_Radiation.H:180
real2d sw_clnsky_flux_dn_dir
Definition: ERF_Radiation.H:308
string1dv gas_names_yakl_offset
Definition: ERF_Radiation.H:198
real3d cld_tau_sw_bnd
Definition: ERF_Radiation.H:337
real2d eff_radius_qi
Definition: ERF_Radiation.H:284
real2d tmp2d
Definition: ERF_Radiation.H:285
real2d qc_lay
Definition: ERF_Radiation.H:280
real2d sw_flux_dn
Definition: ERF_Radiation.H:296
real2d sw_clrsky_flux_dn
Definition: ERF_Radiation.H:304
real2d d_tint
Definition: ERF_Radiation.H:292
real1d m_gas_mol_weights
Definition: ERF_Radiation.H:197
real3d sw_bnd_flux_dir
Definition: ERF_Radiation.H:319

◆ dealloc_buffers()

void Radiation::dealloc_buffers ( )
298 {
299  // 1d size (m_ngas)
300  m_gas_mol_weights.deallocate();
301 
302  // 1d size (1 or nlay)
303  o3_lay.deallocate();
304 
305  // 1d size (ncol)
306  cosine_zenith.deallocate();
307  mu0.deallocate();
308  sfc_alb_dir_vis.deallocate();
309  sfc_alb_dir_nir.deallocate();
310  sfc_alb_dif_vis.deallocate();
311  sfc_alb_dif_nir.deallocate();
312  sfc_flux_dir_vis.deallocate();
313  sfc_flux_dir_nir.deallocate();
314  sfc_flux_dif_vis.deallocate();
315  sfc_flux_dif_nir.deallocate();
316  lat.deallocate();
317  lon.deallocate();
318 
319  // 2d size (ncol, nlay)
320  r_lay.deallocate();
321  p_lay.deallocate();
322  t_lay.deallocate();
323  z_del.deallocate();
324  p_del.deallocate();
325  qv_lay.deallocate();
326  qc_lay.deallocate();
327  qi_lay.deallocate();
328  cldfrac_tot.deallocate();
329  eff_radius_qc.deallocate();
330  eff_radius_qi.deallocate();
331  tmp2d.deallocate();
332  lwp.deallocate();
333  iwp.deallocate();
334 
335  sw_heating.deallocate();
336  lw_heating.deallocate();
337 
338  // 2d size (ncol, nlay+1)
339  d_tint.deallocate();
340  p_lev.deallocate();
341  t_lev.deallocate();
342 
343  sw_flux_up.deallocate();
344  sw_flux_dn.deallocate();
345  sw_flux_dn_dir.deallocate();
346  lw_flux_up.deallocate();
347  lw_flux_dn.deallocate();
348  sw_clnclrsky_flux_up.deallocate();
349  sw_clnclrsky_flux_dn.deallocate();
350  sw_clnclrsky_flux_dn_dir.deallocate();
351  sw_clrsky_flux_up.deallocate();
352  sw_clrsky_flux_dn.deallocate();
353  sw_clrsky_flux_dn_dir.deallocate();
354  sw_clnsky_flux_up.deallocate();
355  sw_clnsky_flux_dn.deallocate();
356  sw_clnsky_flux_dn_dir.deallocate();
357  lw_clnclrsky_flux_up.deallocate();
358  lw_clnclrsky_flux_dn.deallocate();
359  lw_clrsky_flux_up.deallocate();
360  lw_clrsky_flux_dn.deallocate();
361  lw_clnsky_flux_up.deallocate();
362  lw_clnsky_flux_dn.deallocate();
363 
364  // 3d size (ncol, nlay+1, nswbands)
365  sw_bnd_flux_up.deallocate();
366  sw_bnd_flux_dn.deallocate();
367  sw_bnd_flux_dir.deallocate();
368  sw_bnd_flux_dif.deallocate();
369 
370  // 3d size (ncol, nlay+1, nlwbands)
371  lw_bnd_flux_up.deallocate();
372  lw_bnd_flux_dn.deallocate();
373 
374  // 2d size (ncol, nswbands)
375  sfc_alb_dir.deallocate();
376  sfc_alb_dif.deallocate();
377 
378  // 3d size (ncol, nlay, n[sw,lw]bands)
379  aero_tau_sw.deallocate();
380  aero_ssa_sw.deallocate();
381  aero_g_sw.deallocate();
382  aero_tau_lw.deallocate();
383 
384  // 3d size (ncol, nlay, n[sw,lw]bnds)
385  cld_tau_sw_bnd.deallocate();
386  cld_tau_lw_bnd.deallocate();
387 
388  // 3d size (ncol, nlay, n[sw,lw]gpts)
389  cld_tau_sw_gpt.deallocate();
390  cld_tau_lw_gpt.deallocate();
391 }

◆ finalize_impl()

void Radiation::finalize_impl ( )
838 {
839  // Finish rrtmgp
840  m_gas_concs.reset();
842 
843  // Fill the AMReX MFs from YAKL Arrays
845 
846  // Write fluxes if requested
848 
849  // Deallocate the buffer arrays
850  dealloc_buffers();
851 }
void dealloc_buffers()
Definition: ERF_Radiation.cpp:297
void write_rrtmgp_fluxes()
Definition: ERF_Radiation.cpp:568
void yakl_buffers_to_mf()
Definition: ERF_Radiation.cpp:505
GasConcs m_gas_concs
Definition: ERF_Radiation.H:199
void rrtmgp_finalize()
Definition: ERF_RRTMGP_Interface.cpp:257

Referenced by rad_run_impl().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ initialize_impl()

void Radiation::initialize_impl ( )
604 {
605  // Call API to initialize
610 }
void rrtmgp_initialize(GasConcs &gas_concs, const std::string &coefficients_file_sw, const std::string &coefficients_file_lw, const std::string &cloud_optics_file_sw, const std::string &cloud_optics_file_lw)
Definition: ERF_RRTMGP_Interface.cpp:234

Referenced by rad_run_impl().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ mf_to_yakl_buffers()

void Radiation::mf_to_yakl_buffers ( )
396 {
397  bool moist = m_moist;
398  bool ice = m_ice;
399  int ncol = m_ncol;
400  int nlay = m_nlay;
401  Real dz = m_geom.CellSize(2);
402  Real cons_lat = m_lat_cons;
403  Real cons_lon = m_lon_cons;
404  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
405  const auto& vbx = mfi.validbox();
406  const int nx = vbx.length(0);
407  const int imin = vbx.smallEnd(0);
408  const int jmin = vbx.smallEnd(1);
409  const int offset = m_col_offsets[mfi.index()];
410  const Array4<const Real>& cons_arr = m_cons_in->const_array(mfi);
411  const Array4<const Real>& z_arr = (m_z_phys) ? m_z_phys->const_array(mfi) :
412  Array4<const Real>{};
413  const Array4<const Real>& lat_arr = (m_lat) ? m_lat->const_array(mfi) :
414  Array4<const Real>{};
415  const Array4<const Real>& lon_arr = (m_lon) ? m_lon->const_array(mfi) :
416  Array4<const Real>{};
417  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
418  {
419  // map [i,j,k] 0-based to [icol, ilay] 1-based
420  const int icol = (j-jmin)*nx + (i-imin) + 1 + offset;
421  const int ilay = k+1;
422 
423  // EOS input (at CC)
424  Real r = cons_arr(i,j,k,Rho_comp);
425  Real rt = cons_arr(i,j,k,RhoTheta_comp);
426  Real qv = (moist) ? cons_arr(i,j,k,RhoQ1_comp)/r : 0.0;
427  Real qc = (moist) ? cons_arr(i,j,k,RhoQ2_comp)/r : 0.0;
428  Real qi = (ice) ? cons_arr(i,j,k,RhoQ3_comp)/r : 0.0;
429 
430  // EOS avg to z-face
431  Real r_lo = cons_arr(i,j,k-1,Rho_comp);
432  Real rt_lo = cons_arr(i,j,k-1,RhoTheta_comp);
433  Real qv_lo = (moist) ? cons_arr(i,j,k-1,RhoQ1_comp)/r_lo : 0.0;
434  Real r_avg = 0.5 * (r + r_lo);
435  Real rt_avg = 0.5 * (rt + rt_lo);
436  Real qv_avg = 0.5 * (qv + qv_lo);
437 
438  // Buffers at CC
439  r_lay(icol,ilay) = r;
440  p_lay(icol,ilay) = getPgivenRTh(rt, qv);
441  t_lay(icol,ilay) = getTgivenRandRTh(r, rt, qv);
442  z_del(icol,ilay) = (z_arr) ? 0.25 * ( (z_arr(i ,j ,k+1) - z_arr(i ,j ,k))
443  + (z_arr(i+1,j ,k+1) - z_arr(i+1,j ,k))
444  + (z_arr(i ,j+1,k+1) - z_arr(i ,j+1,k))
445  + (z_arr(i+1,j ,k+1) - z_arr(i+1,j ,k)) ) : dz;
446  qv_lay(icol,ilay) = qv;
447  qc_lay(icol,ilay) = qc;
448  qi_lay(icol,ilay) = qi;
449  cldfrac_tot(icol,ilay) = ((qc+qi)>1.0e-5) ? 1. : 0.;
450 
451  // NOTE: These are populated in 'mixing_ratio_to_cloud_mass'
452  lwp(icol,ilay) = 0.0;
453  iwp(icol,ilay) = 0.0;
454 
455  // NOTE: These would be populated from P3 (we use the constants in p3_main_impl.hpp)
456  eff_radius_qc(icol,ilay) = (qc>1.0e-12) ? 10.0e-6 : 0.0;
457  eff_radius_qi(icol,ilay) = (qi>1.0e-12) ? 25.0e-6 : 0.0;
458 
459  // Buffers on z-faces (nlay+1)
460  p_lev(icol,ilay) = getPgivenRTh(rt_avg, qv_avg);
461  t_lev(icol,ilay) = getTgivenRandRTh(r_avg, rt_avg, qv_avg);
462  if (ilay==nlay) {
463  Real r_hi = cons_arr(i,j,k+1,Rho_comp);
464  Real rt_hi = cons_arr(i,j,k+1,RhoTheta_comp);
465  Real qv_hi = (moist) ? cons_arr(i,j,k+1,RhoQ1_comp)/r_hi : 0.0;
466  r_avg = 0.5 * (r + r_hi);
467  rt_avg = 0.5 * (rt + rt_hi);
468  qv_avg = 0.5 * (qv + qv_hi);
469  p_lev(icol,ilay+1) = getPgivenRTh(rt_avg, qv_avg);
470  t_lev(icol,ilay+1) = getTgivenRandRTh(r_avg, rt_avg, qv_avg);
471  }
472 
473  // 1D data structures
474  if (k==0) {
475  lat(icol) = (m_lat) ? lat_arr(i,j,0) : cons_lat;
476  lon(icol) = (m_lon) ? lon_arr(i,j,0) : cons_lon;
477  }
478 
479  });
480  }
481 
482  // Separate YAKL kernel for derived quantities
483  parallel_for(SimpleBounds<2>(ncol, nlay), YAKL_LAMBDA (int icol, int ilay)
484  {
485  p_del(icol,ilay) = p_lev(icol,ilay+1) - p_lev(icol,ilay);
486  });
487 
488  // TODO: Fill properly
489  // No LSM, so follow EAMXX dummy atmos and set constants
490  yakl::memset(mu0, 0.86);
491  yakl::memset(sfc_alb_dir_vis, 0.06);
492  yakl::memset(sfc_alb_dir_nir, 0.06);
493  yakl::memset(sfc_alb_dif_vis, 0.06);
494  yakl::memset(sfc_alb_dif_nir, 0.06);
495 
496  // TODO: Fill properly
497  yakl::memset(aero_tau_sw, 0.0);
498  yakl::memset(aero_ssa_sw, 0.0);
499  yakl::memset(aero_g_sw , 0.0);
500  yakl::memset(aero_tau_lw, 0.0);
501 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=0.)
Definition: ERF_EOS.H:84
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:46
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define RhoQ3_comp
Definition: ERF_IndexDefines.H:44
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
amrex::MultiFab * m_lon
Definition: ERF_Radiation.H:157
amrex::MultiFab * m_cons_in
Definition: ERF_Radiation.H:147
amrex::MultiFab * m_z_phys
Definition: ERF_Radiation.H:153
amrex::Geometry m_geom
Definition: ERF_Radiation.H:131
amrex::MultiFab * m_lat
Definition: ERF_Radiation.H:156
amrex::Vector< int > m_col_offsets
Definition: ERF_Radiation.H:209
@ qv
Definition: ERF_Kessler.H:36
@ qc
Definition: ERF_SatAdj.H:36
Here is the call graph for this function:

◆ rad_run_impl()

void Radiation::rad_run_impl ( )
inline
103  {
104  if (m_update_rad) {
105  amrex::Print() << "Advancing radiation at level: " << m_lev << " ...";
106  this->initialize_impl();
107  this->run_impl();
108  this->finalize_impl();
109  amrex::Print() << "DONE\n";
110  }
111  }
void initialize_impl()
Definition: ERF_Radiation.cpp:603
void run_impl()
Definition: ERF_Radiation.cpp:614
void finalize_impl()
Definition: ERF_Radiation.cpp:837
int m_lev
Definition: ERF_Radiation.H:119
bool m_update_rad
Definition: ERF_Radiation.H:137
Here is the call graph for this function:

◆ run_impl()

void Radiation::run_impl ( )
615 {
616  // Local copies
617  const auto ncol = m_ncol;
618  const auto nlay = m_nlay;
619  const auto nswbands = m_nswbands;
620 
621  // Compute orbital parameters; these are used both for computing
622  // the solar zenith angle and also for computing total solar
623  // irradiance scaling (tsi_scaling).
624  real obliqr, lambm0, mvelpp;
625  int orbital_year = m_orbital_year;
626  real eccen = m_orbital_eccen;
627  real obliq = m_orbital_obliq;
628  real mvelp = m_orbital_mvelp;
629  if (eccen >= 0 && obliq >= 0 && mvelp >= 0) {
630  // fixed orbital parameters forced with orbital_year == ORB_UNDEF_INT
631  orbital_year = ORB_UNDEF_INT;
632  }
633  orbital_params(orbital_year, eccen, obliq,
634  mvelp, obliqr, lambm0, mvelpp);
635 
636  // Use the orbital parameters to calculate the solar declination and eccentricity factor
637  real delta, eccf;
638  // TODO: Generalize the days per month.
639  // Want day + fraction; calday 1 == Jan 1 0Z
640  static constexpr real dpm = (365.0/12.0);
641  real calday = (m_orbital_mon-1.0)*dpm + std::min(m_orbital_day-1.0,dpm-1.0) + m_orbital_sec/86400.0;
642  orbital_decl(calday, eccen, mvelpp, lambm0, obliqr, delta, eccf);
643 
644  // Overwrite eccf if using a fixed solar constant.
645  auto fixed_total_solar_irradiance = m_fixed_total_solar_irradiance;
646  if (fixed_total_solar_irradiance >= 0){
647  eccf = fixed_total_solar_irradiance/1360.9;
648  }
649 
650  // Precompute volume mixing ratio (VMR) for all gases
651  //
652  // H2O is obtained from qv.
653  // All other comps are set to constants for now
654  for (int igas(0); igas < m_ngas; ++igas) {
655  auto name = m_gas_names[igas];
656  auto gas_mol_weight = m_mol_weight_gas[igas];
657  if (name == "H2O") {
658  parallel_for(SimpleBounds<2>(ncol, nlay), YAKL_LAMBDA (int icol, int ilay)
659  {
660  tmp2d(icol,ilay) = qv_lay(icol,ilay) * mwdair/ gas_mol_weight;
661  });
662  } else if (name == "CO2") {
663  yakl::memset(tmp2d, m_co2vmr);
664  } else if (name == "O3") {
665  if (m_o3_size==1) {
666  yakl::memset(tmp2d, m_o3vmr[0] );
667  } else {
668  parallel_for(SimpleBounds<2>(ncol, nlay), YAKL_LAMBDA (int icol, int ilay)
669  {
670  tmp2d(icol,ilay) = o3_lay(ilay);
671  });
672  }
673  } else if (name == "N2O") {
674  yakl::memset(tmp2d, m_n2ovmr);
675  } else if (name == "CO") {
676  yakl::memset(tmp2d, m_covmr );
677  } else if (name == "CH4") {
678  yakl::memset(tmp2d, m_ch4vmr);
679  } else if (name == "O2") {
680  yakl::memset(tmp2d, m_o2vmr );
681  } else if (name == "N2") {
682  yakl::memset(tmp2d, m_n2vmr );
683  } else {
684  Abort("Radiation: Unknown gas component.");
685  }
686 
687  // Populate GasConcs object
688  m_gas_concs.set_vmr(name, tmp2d);
689  }
690 
691 
692  // TODO: No LSM so leaving comment for code
693  // Calculate T_int from longwave flux up from the surface, assuming
694  // blackbody emission with emissivity of 1.
695 
696 
697  // Determine the cosine zenith angle.
698  // This must be done on HOST and copied to device.
699  realHost1d h_mu0("h_mu0", ncol);
700  if (m_fixed_solar_zenith_angle > 0) {
701  yakl::memset(h_mu0, m_fixed_solar_zenith_angle);
702  } else {
703  realHost1d h_lat("h_lat", ncol);
704  realHost1d h_lon("h_lon", ncol);
705  lat.deep_copy_to(h_lat);
706  lon.deep_copy_to(h_lon);
707  parallel_for(ncol, YAKL_LAMBDA (int icol)
708  {
709  // Convert lat/lon to radians
710  real dt = real(m_dt);
711  real lat_col = h_lat(icol)*PI/180.0;
712  real lon_col = h_lon(icol)*PI/180.0;
713  real lcalday = calday;
714  real ldelta = delta;
715  h_mu0(icol) = orbital_cos_zenith(lcalday, lat_col, lon_col, ldelta, m_rad_freq_in_steps * dt);
716  });
717  }
718  h_mu0.deep_copy_to(mu0);
719 
720  // Compute layer cloud mass (per unit area), populates lwp/iwp
723 
724  // Convert to g/m2 (needed by RRTMGP)
725  parallel_for(SimpleBounds<2>(ncol, nlay), YAKL_LAMBDA (int icol, int ilay)
726  {
727  lwp(icol,ilay) *= 1.e3;
728  iwp(icol,ilay) *= 1.e3;
729  });
730 
731  // Compute band-by-band surface_albedos. This is needed since
732  // the AD passes broadband albedos, but rrtmgp require band-by-band.
737 
738  // Run RRTMGP driver
740  p_lay, t_lay, p_lev, t_lev,
741  m_gas_concs,
756 
757 #if 0
758  // UNIT TEST
759  //================================================================================
760  yakl::memset(mu0, 0.86);
761  yakl::memset(sfc_alb_dir_vis, 0.06);
762  yakl::memset(sfc_alb_dir_nir, 0.06);
763  yakl::memset(sfc_alb_dif_vis, 0.06);
764  yakl::memset(sfc_alb_dif_nir, 0.06);
765 
766  // Generate some fake liquid and ice water data. We pick values to be midway between
767  // the min and max of the valid lookup table values for effective radii
768  real rel_val = 0.5 * (rrtmgp::cloud_optics_sw.get_min_radius_liq()
769  + rrtmgp::cloud_optics_sw.get_max_radius_liq());
770  real rei_val = 0.5 * (rrtmgp::cloud_optics_sw.get_min_radius_ice()
771  + rrtmgp::cloud_optics_sw.get_max_radius_ice());
772 
773  // Restrict clouds to troposphere (> 100 hPa = 100*100 Pa) and not very close to the ground (< 900 hPa), and
774  // put them in 2/3 of the columns since that's roughly the total cloudiness of earth.
775  // Set sane values for liquid and ice water path.
776  // NOTE: these "sane" values are in g/m2!
777  parallel_for( SimpleBounds<2>(nlay,ncol) , YAKL_LAMBDA (int ilay, int icol)
778  {
779  cldfrac_tot(icol,ilay) = (p_lay(icol,ilay) > 100._wp * 100._wp) &&
780  (p_lay(icol,ilay) < 900._wp * 100._wp) &&
781  (icol%3 != 0);
782  // Ice and liquid will overlap in a few layers
783  lwp(icol,ilay) = (cldfrac_tot(icol,ilay) && t_lay(icol,ilay) > 263._wp) ? 10._wp : 0._wp;
784  iwp(icol,ilay) = (cldfrac_tot(icol,ilay) && t_lay(icol,ilay) < 273._wp) ? 10._wp : 0._wp;
785  eff_radius_qc(icol,ilay) = (lwp(icol,ilay) > 0._wp) ? rel_val : 0._wp;
786  eff_radius_qi(icol,ilay) = (iwp(icol,ilay) > 0._wp) ? rei_val : 0._wp;
787  });
788 
793 
794  yakl::memset(aero_tau_sw, 0.0);
795  yakl::memset(aero_ssa_sw, 0.0);
796  yakl::memset(aero_g_sw , 0.0);
797  yakl::memset(aero_tau_lw, 0.0);
798 
800  p_lay, t_lay, p_lev, t_lev,
801  m_gas_concs,
815  1.0);
816  //================================================================================
817 #endif
818 
819  // Update heating tendency
822 
823  // Compute surface fluxes
824  const int kbot = nlay + 1; // Should this be 1 for our layout?
825  parallel_for(SimpleBounds<3>(ncol, nlay+1, nswbands), YAKL_LAMBDA (int icol, int ilay, int ibnd)
826  {
827  sw_bnd_flux_dif(icol,ilay,ibnd) = sw_bnd_flux_dn(icol,ilay,ibnd) - sw_bnd_flux_dir(icol,ilay,ibnd);
828  });
829  rrtmgp::compute_broadband_surface_fluxes(ncol, kbot, nswbands,
833 }
static constexpr int ORB_UNDEF_INT
Definition: ERF_Constants.H:96
constexpr amrex::Real mwdair
Definition: ERF_Constants.H:64
constexpr amrex::Real PI
Definition: ERF_Constants.H:6
double real
Definition: ERF_OrbCosZenith.H:9
real orbital_cos_zenith(real &jday, real &lat, real &lon, real &declin, real dt_avg, real uniform_angle, real constant_zenith_angle_deg)
Definition: ERF_OrbCosZenith.cpp:4
void orbital_params(int &iyear_AD, real &eccen, real &obliq, real &mvelp, real &obliqr, real &lambm0, real &mvelpp)
Definition: ERF_OrbCosZenith.cpp:139
void orbital_decl(real &calday, real &eccen, real &mvelpp, real &lambm0, real &obliqr, real &delta, real &eccf)
Definition: ERF_OrbCosZenith.cpp:512
amrex::Real m_dt
Definition: ERF_Radiation.H:128
int m_orbital_mon
Definition: ERF_Radiation.H:222
int m_orbital_sec
Definition: ERF_Radiation.H:224
int m_orbital_day
Definition: ERF_Radiation.H:223
CloudOptics cloud_optics_sw
Definition: ERF_RRTMGP_Interface.cpp:34
void compute_band_by_band_surface_albedos(const int ncol, const int nswbands, real1d &sfc_alb_dir_vis, real1d &sfc_alb_dir_nir, real1d &sfc_alb_dif_vis, real1d &sfc_alb_dif_nir, real2d &sfc_alb_dir, real2d &sfc_alb_dif)
Definition: ERF_RRTMGP_Interface.cpp:268
void rrtmgp_main(const int ncol, const int nlay, real2d &p_lay, real2d &t_lay, real2d &p_lev, real2d &t_lev, GasConcs &gas_concs, real2d &sfc_alb_dir, real2d &sfc_alb_dif, real1d &mu0, real2d &lwp, real2d &iwp, real2d &rel, real2d &rei, real2d &cldfrac, real3d &aer_tau_sw, real3d &aer_ssa_sw, real3d &aer_asm_sw, real3d &aer_tau_lw, real3d &cld_tau_sw_bnd, real3d &cld_tau_lw_bnd, real3d &cld_tau_sw_gpt, real3d &cld_tau_lw_gpt, real2d &sw_flux_up, real2d &sw_flux_dn, real2d &sw_flux_dn_dir, real2d &lw_flux_up, real2d &lw_flux_dn, real2d &sw_clnclrsky_flux_up, real2d &sw_clnclrsky_flux_dn, real2d &sw_clnclrsky_flux_dn_dir, real2d &sw_clrsky_flux_up, real2d &sw_clrsky_flux_dn, real2d &sw_clrsky_flux_dn_dir, real2d &sw_clnsky_flux_up, real2d &sw_clnsky_flux_dn, real2d &sw_clnsky_flux_dn_dir, real2d &lw_clnclrsky_flux_up, real2d &lw_clnclrsky_flux_dn, real2d &lw_clrsky_flux_up, real2d &lw_clrsky_flux_dn, real2d &lw_clnsky_flux_up, real2d &lw_clnsky_flux_dn, real3d &sw_bnd_flux_up, real3d &sw_bnd_flux_dn, real3d &sw_bnd_flux_dn_dir, real3d &lw_bnd_flux_up, real3d &lw_bnd_flux_dn, const real tsi_scaling, const bool extra_clnclrsky_diag, const bool extra_clnsky_diag)
Definition: ERF_RRTMGP_Interface.cpp:391
void compute_heating_rate(yakl::Array< T, 2, myMem, myStyle > const &flux_up, yakl::Array< T, 2, myMem, myStyle > const &flux_dn, yakl::Array< T, 2, myMem, myStyle > const &rho, yakl::Array< T, 2, myMem, myStyle > const &dz, yakl::Array< T, 2, myMem, myStyle > &heating_rate)
Definition: ERF_RRTMGP_Utils.H:73
void compute_broadband_surface_fluxes(const int ncol, const int ktop, const int nswbands, real3d &sw_bnd_flux_dir, real3d &sw_bnd_flux_dif, real1d &sfc_flux_dir_vis, real1d &sfc_flux_dir_nir, real1d &sfc_flux_dif_vis, real1d &sfc_flux_dif_nir)
Definition: ERF_RRTMGP_Interface.cpp:328
void mixing_ratio_to_cloud_mass(yakl::Array< T, 2, myMem, myStyle > const &mixing_ratio, yakl::Array< T, 2, myMem, myStyle > const &cloud_fraction, yakl::Array< T, 2, myMem, myStyle > const &dp, yakl::Array< T, 2, myMem, myStyle > &cloud_mass)
Definition: ERF_RRTMGP_Utils.H:23

Referenced by rad_run_impl().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_grids()

void Radiation::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::MultiFab *  lsm_fluxes,
amrex::MultiFab *  lsm_zenith,
amrex::MultiFab *  qheating_rates,
amrex::MultiFab *  z_phys,
amrex::MultiFab *  lat,
amrex::MultiFab *  lon 
)
130 {
131  // Set data members that may change
132  m_lev = level;
133  m_step = step;
134  m_time = time;
135  m_dt = dt;
136  m_geom = geom;
137  m_cons_in = cons_in;
138  m_lsm_fluxes = lsm_fluxes;
139  m_qheating_rates = qheating_rates;
140  m_z_phys = z_phys;
141  m_lat = lat;
142  m_lon = lon;
143 
144  // Update the day and month
145  time_t timestamp = time_t(int(time));
146  struct tm *timeinfo = localtime(&timestamp);
147  if (m_fixed_orbital_year) {
148  m_orbital_mon = timeinfo->tm_mon + 1;
149  m_orbital_day = timeinfo->tm_mday;
150  m_orbital_sec = timeinfo->tm_hour*3600 + timeinfo->tm_min*60 + timeinfo->tm_sec;
151  } else {
152  m_orbital_year = timeinfo->tm_year + 1900;
153  m_orbital_mon = timeinfo->tm_mon + 1;
154  m_orbital_day = timeinfo->tm_mday;
155  m_orbital_sec = timeinfo->tm_hour*3600 + timeinfo->tm_min*60 + timeinfo->tm_sec;
156  }
157 
158  // Only allocate and proceed if we are going to update radiation
159  m_update_rad = false;
160  if (m_rad_freq_in_steps > 0) { m_update_rad = ( (m_step == 0) || (m_step % m_rad_freq_in_steps == 0) ); }
161 
162  if (m_update_rad) {
163  // Reset vector of offsets for columnar data
164  m_nlay = geom.Domain().length(2);
165 
166  m_ncol = 0;
167  m_col_offsets.clear();
168  m_col_offsets.resize(int(ba.size()));
169  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
170  const auto& vbx = mfi.validbox();
171  int nx = vbx.length(0);
172  int ny = vbx.length(1);
173  m_col_offsets[mfi.index()] = m_ncol;
174  m_ncol += nx * ny;
175  }
176 
177  // Allocate the buffer arrays
178  alloc_buffers();
179 
180  // Fill the YAKL Arrays from AMReX MFs
182  }
183 }
int m_step
Definition: ERF_Radiation.H:122
void mf_to_yakl_buffers()
Definition: ERF_Radiation.cpp:395
amrex::MultiFab * m_lsm_fluxes
Definition: ERF_Radiation.H:164
amrex::MultiFab * m_qheating_rates
Definition: ERF_Radiation.H:150
void alloc_buffers()
Definition: ERF_Radiation.cpp:186
amrex::Real m_time
Definition: ERF_Radiation.H:125

◆ write_rrtmgp_fluxes()

void Radiation::write_rrtmgp_fluxes ( )
569 {
570  int n_fluxes = 5;
571  MultiFab mf_flux(m_cons_in->boxArray(), m_cons_in->DistributionMap(), n_fluxes, 0);
572 
573  for (MFIter mfi(mf_flux); mfi.isValid(); ++mfi) {
574  const auto& vbx = mfi.validbox();
575  const int nx = vbx.length(0);
576  const int imin = vbx.smallEnd(0);
577  const int jmin = vbx.smallEnd(1);
578  const int offset = m_col_offsets[mfi.index()];
579  const Array4<Real>& dst_arr = mf_flux.array(mfi);
580  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
581  {
582  // map [i,j,k] 0-based to [icol, ilay] 1-based
583  const int icol = (j-jmin)*nx + (i-imin) + 1 + offset;
584  const int ilay = k+1;
585 
586  // SW and LW fluxes
587  dst_arr(i,j,k,0) = sw_flux_up(icol,ilay);
588  dst_arr(i,j,k,1) = sw_flux_dn(icol,ilay);
589  dst_arr(i,j,k,2) = sw_flux_dn_dir(icol,ilay);
590  dst_arr(i,j,k,3) = lw_flux_up(icol,ilay);
591  dst_arr(i,j,k,4) = lw_flux_dn(icol,ilay);
592  });
593  }
594 
595 
596  std::string plotfilename = amrex::Concatenate("plt_rad", m_step, 5);
597  Vector<std::string> flux_names = {"sw_flux_up", "sw_flux_dn", "sw_flux_dir",
598  "lw_flux_up", "lw_flux_dn"};
599  WriteSingleLevelPlotfile(plotfilename, mf_flux, flux_names, m_geom, m_time, m_step);
600 }
Here is the call graph for this function:

◆ yakl_buffers_to_mf()

void Radiation::yakl_buffers_to_mf ( )
506 {
507  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
508  const auto& vbx = mfi.validbox();
509  const auto& sbx = makeSlab(vbx,2,vbx.smallEnd(2));
510  const int nx = vbx.length(0);
511  const int imin = vbx.smallEnd(0);
512  const int jmin = vbx.smallEnd(1);
513  const int offset = m_col_offsets[mfi.index()];
514  const Array4<Real>& q_arr = m_qheating_rates->array(mfi);
515  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
516  {
517  // map [i,j,k] 0-based to [icol, ilay] 1-based
518  const int icol = (j-jmin)*nx + (i-imin) + 1 + offset;
519  const int ilay = k+1;
520 
521  // Temperature heating rate for SW and LW
522  q_arr(i,j,k,0) = sw_heating(icol,ilay);
523  q_arr(i,j,k,1) = lw_heating(icol,ilay);
524 
525  // Convert the rates for theta_d
526  Real exner = getExnergivenP(Real(p_lay(icol,ilay)), R_d/Cp_d);
527  q_arr(i,j,k,0) *= exner;
528  q_arr(i,j,k,1) *= exner;
529 
530  /*
531  if (i==0 && j==0) {
532  AllPrint() << "Qsrcs: " << IntVect(i,j,k) << ' '
533  << IntVect(icol,0,ilay) << ' '
534  << q_arr(i,j,k,0) << ' '
535  << q_arr(i,j,k,1) << ' '
536  << sw_heating(icol,ilay) << ' '
537  << lw_heating(icol,ilay) << ' '
538  << exner << "\n";
539  }
540  */
541 
542  });
543  if (m_lsm_fluxes) {
544  const Array4<Real>& lsm_arr = m_lsm_fluxes->array(mfi);
545  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
546  {
547  // map [i,j,k] 0-based to [icol, ilay] 1-based
548  const int icol = (j-jmin)*nx + (i-imin) + 1 + offset;
549 
550  // SW fluxes for LSM
551  lsm_arr(i,j,k,0) = sfc_flux_dir_vis(icol);
552  lsm_arr(i,j,k,1) = sfc_flux_dir_nir(icol);
553  lsm_arr(i,j,k,2) = sfc_flux_dif_vis(icol);
554  lsm_arr(i,j,k,3) = sfc_flux_dif_nir(icol);
555 
556  // Net SW flux for LSM
557  lsm_arr(i,j,k,4) = sfc_flux_dir_vis(icol) + sfc_flux_dir_nir(icol)
558  + sfc_flux_dif_vis(icol) + sfc_flux_dif_nir(icol);
559 
560  // LW flux for LSM (at bottom surface)
561  lsm_arr(i,j,k,5) = lw_flux_dn(icol,1);
562  });
563  }
564  }
565 }
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:12
constexpr amrex::Real R_d
Definition: ERF_Constants.H:10
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getExnergivenP(const amrex::Real P, const amrex::Real rdOcp)
Definition: ERF_EOS.H:144
Here is the call graph for this function:

Member Data Documentation

◆ aero_g_sw

real3d Radiation::aero_g_sw
private

◆ aero_ssa_sw

real3d Radiation::aero_ssa_sw
private

◆ aero_tau_lw

real3d Radiation::aero_tau_lw
private

◆ aero_tau_sw

real3d Radiation::aero_tau_sw
private

◆ cld_tau_lw_bnd

real3d Radiation::cld_tau_lw_bnd
private

◆ cld_tau_lw_gpt

real3d Radiation::cld_tau_lw_gpt
private

◆ cld_tau_sw_bnd

real3d Radiation::cld_tau_sw_bnd
private

◆ cld_tau_sw_gpt

real3d Radiation::cld_tau_sw_gpt
private

◆ cldfrac_tot

real2d Radiation::cldfrac_tot
private

◆ cosine_zenith

real1d Radiation::cosine_zenith
private

◆ d_tint

real2d Radiation::d_tint
private

◆ eff_radius_qc

real2d Radiation::eff_radius_qc
private

◆ eff_radius_qi

real2d Radiation::eff_radius_qi
private

◆ gas_names_yakl_offset

string1dv Radiation::gas_names_yakl_offset
private

◆ iwp

real2d Radiation::iwp
private

◆ lat

real1d Radiation::lat
private

◆ lon

real1d Radiation::lon
private

◆ lw_bnd_flux_dn

real3d Radiation::lw_bnd_flux_dn
private

◆ lw_bnd_flux_up

real3d Radiation::lw_bnd_flux_up
private

◆ lw_clnclrsky_flux_dn

real2d Radiation::lw_clnclrsky_flux_dn
private

◆ lw_clnclrsky_flux_up

real2d Radiation::lw_clnclrsky_flux_up
private

◆ lw_clnsky_flux_dn

real2d Radiation::lw_clnsky_flux_dn
private

◆ lw_clnsky_flux_up

real2d Radiation::lw_clnsky_flux_up
private

◆ lw_clrsky_flux_dn

real2d Radiation::lw_clrsky_flux_dn
private

◆ lw_clrsky_flux_up

real2d Radiation::lw_clrsky_flux_up
private

◆ lw_flux_dn

real2d Radiation::lw_flux_dn
private

◆ lw_flux_up

real2d Radiation::lw_flux_up
private

◆ lw_heating

real2d Radiation::lw_heating
private

◆ lwp

real2d Radiation::lwp
private

◆ m_ba

amrex::BoxArray Radiation::m_ba
private

◆ m_ch4vmr

amrex::Real Radiation::m_ch4vmr = 1807.851e-9
private

◆ m_co2vmr

amrex::Real Radiation::m_co2vmr = 388.717e-6
private

◆ m_col_offsets

amrex::Vector<int> Radiation::m_col_offsets
private

◆ m_cons_in

amrex::MultiFab* Radiation::m_cons_in = nullptr
private

◆ m_covmr

amrex::Real Radiation::m_covmr = 1.0e-7
private

◆ m_do_aerosol_rad

bool Radiation::m_do_aerosol_rad = true
private

◆ m_do_subcol_sampling

bool Radiation::m_do_subcol_sampling = true
private

◆ m_dt

amrex::Real Radiation::m_dt
private

◆ m_extra_clnclrsky_diag

bool Radiation::m_extra_clnclrsky_diag = false
private

◆ m_extra_clnsky_diag

bool Radiation::m_extra_clnsky_diag = false
private

◆ m_fixed_orbital_year

bool Radiation::m_fixed_orbital_year = false
private

◆ m_fixed_solar_zenith_angle

amrex::Real Radiation::m_fixed_solar_zenith_angle = -9999.
private

◆ m_fixed_total_solar_irradiance

amrex::Real Radiation::m_fixed_total_solar_irradiance = -9999.
private

◆ m_gas_concs

GasConcs Radiation::m_gas_concs
private

◆ m_gas_mol_weights

real1d Radiation::m_gas_mol_weights
private

◆ m_gas_names

const std::vector<std::string> Radiation::m_gas_names
private
Initial value:
= {"H2O", "CO2", "O3", "N2O",
"CO" , "CH4", "O2", "N2" }

◆ m_geom

amrex::Geometry Radiation::m_geom
private

◆ m_ice

bool Radiation::m_ice = false
private

◆ m_lat

amrex::MultiFab* Radiation::m_lat = nullptr
private

◆ m_lat_cons

amrex::Real Radiation::m_lat_cons = 39.809860
private

◆ m_lev

int Radiation::m_lev
private

Referenced by rad_run_impl().

◆ m_lon

amrex::MultiFab* Radiation::m_lon = nullptr
private

◆ m_lon_cons

amrex::Real Radiation::m_lon_cons = -98.555183
private

◆ m_lsm_fluxes

amrex::MultiFab* Radiation::m_lsm_fluxes = nullptr
private

◆ m_lsm_zenith

amrex::MultiFab* Radiation::m_lsm_zenith = nullptr
private

◆ m_moist

bool Radiation::m_moist = false
private

◆ m_mol_weight_gas

const std::vector<amrex::Real> Radiation::m_mol_weight_gas
private
Initial value:
= {18.01528, 44.00950, 47.9982, 44.0128,
28.01010, 16.04246, 31.9980, 28.0134}

◆ m_n2ovmr

amrex::Real Radiation::m_n2ovmr = 323.141e-9
private

◆ m_n2vmr

amrex::Real Radiation::m_n2vmr = 0.7906
private

◆ m_ncol

int Radiation::m_ncol
private

◆ m_ngas

int Radiation::m_ngas = 8
private

◆ m_nlay

int Radiation::m_nlay
private

◆ m_nlwbands

int Radiation::m_nlwbands = 16
private

◆ m_nlwgpts

int Radiation::m_nlwgpts = 128
private

◆ m_nswbands

int Radiation::m_nswbands = 14
private

◆ m_nswgpts

int Radiation::m_nswgpts = 112
private

◆ m_o2vmr

amrex::Real Radiation::m_o2vmr = 0.209448
private

◆ m_o3_size

int Radiation::m_o3_size
private

◆ m_o3vmr

amrex::Vector<amrex::Real> Radiation::m_o3vmr
private

◆ m_orbital_day

int Radiation::m_orbital_day = -9999
private

◆ m_orbital_eccen

amrex::Real Radiation::m_orbital_eccen = -9999.
private

◆ m_orbital_mon

int Radiation::m_orbital_mon = -9999
private

◆ m_orbital_mvelp

amrex::Real Radiation::m_orbital_mvelp = -9999.
private

◆ m_orbital_obliq

amrex::Real Radiation::m_orbital_obliq = -9999.
private

◆ m_orbital_sec

int Radiation::m_orbital_sec = -9999
private

◆ m_orbital_year

int Radiation::m_orbital_year = -9999
private

◆ m_qheating_rates

amrex::MultiFab* Radiation::m_qheating_rates = nullptr
private

◆ m_rad_freq_in_steps

int Radiation::m_rad_freq_in_steps = 1
private

◆ m_rad_write_fluxes

bool Radiation::m_rad_write_fluxes = false
private

◆ m_step

int Radiation::m_step
private

◆ m_time

amrex::Real Radiation::m_time
private

◆ m_update_rad

bool Radiation::m_update_rad = false
private

Referenced by rad_run_impl().

◆ m_z_phys

amrex::MultiFab* Radiation::m_z_phys = nullptr
private

◆ mu0

real1d Radiation::mu0
private

◆ o3_lay

real1d Radiation::o3_lay
private

◆ p_del

real2d Radiation::p_del
private

◆ p_lay

real2d Radiation::p_lay
private

◆ p_lev

real2d Radiation::p_lev
private

◆ qc_lay

real2d Radiation::qc_lay
private

◆ qi_lay

real2d Radiation::qi_lay
private

◆ qv_lay

real2d Radiation::qv_lay
private

◆ r_lay

real2d Radiation::r_lay
private

◆ rrtmgp_cloud_optics_file_lw

std::string Radiation::rrtmgp_cloud_optics_file_lw
private

◆ rrtmgp_cloud_optics_file_sw

std::string Radiation::rrtmgp_cloud_optics_file_sw
private

◆ rrtmgp_cloud_optics_lw

std::string Radiation::rrtmgp_cloud_optics_lw = "rrtmgp-cloud-optics-coeffs-lw.nc"
private

◆ rrtmgp_cloud_optics_sw

std::string Radiation::rrtmgp_cloud_optics_sw = "rrtmgp-cloud-optics-coeffs-sw.nc"
private

◆ rrtmgp_coeffs_file_lw

std::string Radiation::rrtmgp_coeffs_file_lw
private

◆ rrtmgp_coeffs_file_sw

std::string Radiation::rrtmgp_coeffs_file_sw
private

◆ rrtmgp_coeffs_lw

std::string Radiation::rrtmgp_coeffs_lw = "rrtmgp-data-lw-g224-2018-12-04.nc"
private

◆ rrtmgp_coeffs_sw

std::string Radiation::rrtmgp_coeffs_sw = "rrtmgp-data-sw-g224-2018-12-04.nc"
private

◆ rrtmgp_file_path

std::string Radiation::rrtmgp_file_path = "."
private

◆ sfc_alb_dif

real2d Radiation::sfc_alb_dif
private

◆ sfc_alb_dif_nir

real1d Radiation::sfc_alb_dif_nir
private

◆ sfc_alb_dif_vis

real1d Radiation::sfc_alb_dif_vis
private

◆ sfc_alb_dir

real2d Radiation::sfc_alb_dir
private

◆ sfc_alb_dir_nir

real1d Radiation::sfc_alb_dir_nir
private

◆ sfc_alb_dir_vis

real1d Radiation::sfc_alb_dir_vis
private

◆ sfc_flux_dif_nir

real1d Radiation::sfc_flux_dif_nir
private

◆ sfc_flux_dif_vis

real1d Radiation::sfc_flux_dif_vis
private

◆ sfc_flux_dir_nir

real1d Radiation::sfc_flux_dir_nir
private

◆ sfc_flux_dir_vis

real1d Radiation::sfc_flux_dir_vis
private

◆ sw_bnd_flux_dif

real3d Radiation::sw_bnd_flux_dif
private

◆ sw_bnd_flux_dir

real3d Radiation::sw_bnd_flux_dir
private

◆ sw_bnd_flux_dn

real3d Radiation::sw_bnd_flux_dn
private

◆ sw_bnd_flux_up

real3d Radiation::sw_bnd_flux_up
private

◆ sw_clnclrsky_flux_dn

real2d Radiation::sw_clnclrsky_flux_dn
private

◆ sw_clnclrsky_flux_dn_dir

real2d Radiation::sw_clnclrsky_flux_dn_dir
private

◆ sw_clnclrsky_flux_up

real2d Radiation::sw_clnclrsky_flux_up
private

◆ sw_clnsky_flux_dn

real2d Radiation::sw_clnsky_flux_dn
private

◆ sw_clnsky_flux_dn_dir

real2d Radiation::sw_clnsky_flux_dn_dir
private

◆ sw_clnsky_flux_up

real2d Radiation::sw_clnsky_flux_up
private

◆ sw_clrsky_flux_dn

real2d Radiation::sw_clrsky_flux_dn
private

◆ sw_clrsky_flux_dn_dir

real2d Radiation::sw_clrsky_flux_dn_dir
private

◆ sw_clrsky_flux_up

real2d Radiation::sw_clrsky_flux_up
private

◆ sw_flux_dn

real2d Radiation::sw_flux_dn
private

◆ sw_flux_dn_dir

real2d Radiation::sw_flux_dn_dir
private

◆ sw_flux_up

real2d Radiation::sw_flux_up
private

◆ sw_heating

real2d Radiation::sw_heating
private

◆ t_lay

real2d Radiation::t_lay
private

◆ t_lev

real2d Radiation::t_lev
private

◆ tmp2d

real2d Radiation::tmp2d
private

◆ z_del

real2d Radiation::z_del
private

The documentation for this class was generated from the following files: