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

#include <ERF_Radiation.H>

Inheritance diagram for Radiation:
Collaboration diagram for Radiation:

Public Member Functions

 Radiation (const int &lev, SolverChoice &sc)
 
 ~Radiation ()
 
virtual void Init (const amrex::Geometry &geom, const amrex::BoxArray &ba, amrex::MultiFab *cons_in) override
 
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::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
 
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::Vector< amrex::MultiFab * > &lsm_input_ptrs, amrex::MultiFab *qheating_rates, amrex::MultiFab *rad_fluxes, amrex::MultiFab *z_phys, amrex::MultiFab *lat, amrex::MultiFab *lon)
 
void alloc_buffers ()
 
void dealloc_buffers ()
 
void mf_to_kokkos_buffers (amrex::Vector< amrex::MultiFab * > &lsm_input_ptrs)
 
void kokkos_buffers_to_mf (amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
 
void write_rrtmgp_fluxes ()
 
void initialize_impl ()
 
void run_impl ()
 
void finalize_impl (amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
 
void rad_run_impl (amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
 
virtual amrex::Vector< std::string > get_lsm_input_varnames () override
 
virtual amrex::Vector< std::string > get_lsm_output_varnames () override
 
void populateDatalogMF ()
 
virtual void WriteDataLog (const amrex::Real &time) override
 
- Public Member Functions inherited from IRadiation
virtual ~IRadiation ()=default
 
void setupDataLog ()
 
void setDataLogFrequency (const int nstep)
 
bool hasDatalog ()
 

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_first_step = true
 
bool m_moist = false
 
bool m_ice = false
 
bool m_lsm = false
 
amrex::Vector< std::string > m_lsm_input_names
 
amrex::Vector< std::string > m_lsm_output_names = {"cos_zenith_angle", "sw_flux_dn", "lw_flux_dn"}
 
amrex::Real m_rad_t_sfc = -1
 
amrex::MultiFab * m_cons_in = nullptr
 
amrex::MultiFab * m_qheating_rates = nullptr
 
amrex::MultiFab * m_rad_fluxes = 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
 
amrex::MultiFab datalog_mf
 
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-g256-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::Realm_mol_weight_gas
 
amrex::Real m_co2vmr = 388.717e-6
 
amrex::Vector< amrex::Realm_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_k m_gas_mol_weights
 
std::vector< std::string > gas_names_offset
 
GasConcsK< amrex::Real, layout_t, KokkosDefaultDevicem_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
 
int m_nlwbands
 
int m_nswgpts
 
int m_nlwgpts
 
int m_rad_freq_in_steps = 1
 
bool m_do_subcol_sampling = true
 
real1d_k o3_lay
 
real1d_k mu0
 
real1d_k sfc_alb_dir_vis
 
real1d_k sfc_alb_dir_nir
 
real1d_k sfc_alb_dif_vis
 
real1d_k sfc_alb_dif_nir
 
real1d_k sfc_flux_dir_vis
 
real1d_k sfc_flux_dir_nir
 
real1d_k sfc_flux_dif_vis
 
real1d_k sfc_flux_dif_nir
 
real1d_k lat
 
real1d_k lon
 
real1d_k sfc_emis
 
real1d_k t_sfc
 
real1d_k lw_src
 
real2d_k r_lay
 
real2d_k p_lay
 
real2d_k t_lay
 
real2d_k z_del
 
real2d_k qv_lay
 
real2d_k qc_lay
 
real2d_k qi_lay
 
real2d_k cldfrac_tot
 
real2d_k eff_radius_qc
 
real2d_k eff_radius_qi
 
real2d_k lwp
 
real2d_k iwp
 
real2d_k sw_heating
 
real2d_k lw_heating
 
real2d_k sw_clrsky_heating
 
real2d_k lw_clrsky_heating
 
real2d_k d_tint
 
real2d_k p_lev
 
real2d_k t_lev
 
real2d_k sw_flux_up
 
real2d_k sw_flux_dn
 
real2d_k sw_flux_dn_dir
 
real2d_k lw_flux_up
 
real2d_k lw_flux_dn
 
real2d_k sw_clnclrsky_flux_up
 
real2d_k sw_clnclrsky_flux_dn
 
real2d_k sw_clnclrsky_flux_dn_dir
 
real2d_k sw_clrsky_flux_up
 
real2d_k sw_clrsky_flux_dn
 
real2d_k sw_clrsky_flux_dn_dir
 
real2d_k sw_clnsky_flux_up
 
real2d_k sw_clnsky_flux_dn
 
real2d_k sw_clnsky_flux_dn_dir
 
real2d_k lw_clnclrsky_flux_up
 
real2d_k lw_clnclrsky_flux_dn
 
real2d_k lw_clrsky_flux_up
 
real2d_k lw_clrsky_flux_dn
 
real2d_k lw_clnsky_flux_up
 
real2d_k lw_clnsky_flux_dn
 
real3d_k sw_bnd_flux_up
 
real3d_k sw_bnd_flux_dn
 
real3d_k sw_bnd_flux_dir
 
real3d_k sw_bnd_flux_dif
 
real3d_k lw_bnd_flux_up
 
real3d_k lw_bnd_flux_dn
 
real2d_k sfc_alb_dir
 
real2d_k sfc_alb_dif
 
real2d_k emis_sfc
 
real3d_k aero_tau_sw
 
real3d_k aero_ssa_sw
 
real3d_k aero_g_sw
 
real3d_k aero_tau_lw
 
real3d_k cld_tau_sw_bnd
 
real3d_k cld_tau_lw_bnd
 
real3d_k cld_tau_sw_gpt
 
real3d_k cld_tau_lw_gpt
 

Additional Inherited Members

- Protected Attributes inherited from IRadiation
std::unique_ptr< std::fstream > datalog = nullptr
 
std::string datalogname
 
int datalog_int = -1
 

Constructor & Destructor Documentation

◆ Radiation()

Radiation::Radiation ( const int &  lev,
SolverChoice sc 
)
22 {
23  // Note that Kokkos is now initialized in main.cpp
24 
25  // Check if we have a valid moisture model
26  if (sc.moisture_type != MoistureType::None) { m_moist = true; }
27 
28  // Check if we have a moisture model with ice
29  if (sc.moisture_type == MoistureType::SAM) { m_ice = true; }
30 
31  // Check if we have a land surface model enabled
32  if (sc.lsm_type != LandSurfaceType::None) { m_lsm = true; }
33 
34  // Construct parser object for following reads
35  ParmParse pp("erf");
36 
37  // Must specify a surface temp without a LSM
38  if (!m_lsm) { pp.get("rad_t_sfc",m_rad_t_sfc); }
39 
40  // Radiation timestep, as a number of atm steps
41  pp.query("rad_freq_in_steps", m_rad_freq_in_steps);
42 
43  // Flag to write fluxes to plt file
44  pp.query("rad_write_fluxes", m_rad_write_fluxes);
45 
46  // Do MCICA subcolumn sampling
47  pp.query("rad_do_subcol_sampling", m_do_subcol_sampling);
48 
49  // Determine orbital year. If orbital_year is negative, use current year
50  // from timestamp for orbital year; if positive, use provided orbital year
51  // for duration of simulation.
52  m_fixed_orbital_year = pp.query("rad_orbital_year", m_orbital_year);
53 
54  // Get orbital parameters from inputs file
55  pp.query("rad_orbital_eccentricity", m_orbital_eccen);
56  pp.query("rad_orbital_obliquity" , m_orbital_obliq);
57  pp.query("rad_orbital_mvelp" , m_orbital_mvelp);
58 
59  // Get a constant lat/lon for idealized simulations
60  pp.query("rad_cons_lat", m_lat_cons);
61  pp.query("rad_cons_lon", m_lon_cons);
62 
63  // Value for prescribing an invariant solar constant (i.e. total solar irradiance at
64  // TOA). Used for idealized experiments such as RCE. Disabled when value is less than 0.
65  pp.query("fixed_total_solar_irradiance", m_fixed_total_solar_irradiance);
66 
67  // Determine whether or not we are using a fixed solar zenith angle (positive value)
68  pp.query("fixed_solar_zenith_angle", m_fixed_solar_zenith_angle);
69 
70  // Get prescribed surface values of greenhouse gases
71  pp.query("co2vmr", m_co2vmr);
72  pp.queryarr("o3vmr" , m_o3vmr );
73  pp.query("n2ovmr", m_n2ovmr);
74  pp.query("covmr" , m_covmr );
75  pp.query("ch4vmr", m_ch4vmr);
76  pp.query("o2vmr" , m_o2vmr );
77  pp.query("n2vmr" , m_n2vmr );
78 
79  // Required aerosol optical properties from SPA
80  pp.query("rad_do_aerosol", m_do_aerosol_rad);
81 
82  // Whether we do extra clean/clear sky calculations
83  pp.query("rad_extra_clnclrsky_diag", m_extra_clnclrsky_diag);
84  pp.query("rad_extra_clnsky_diag" , m_extra_clnsky_diag);
85 
86  // Parse path and file names
87  pp.query("rrtmgp_file_path" , rrtmgp_file_path);
88  pp.query("rrtmgp_coeffs_sw" , rrtmgp_coeffs_sw );
89  pp.query("rrtmgp_coeffs_lw" , rrtmgp_coeffs_lw );
90  pp.query("rrtmgp_cloud_optics_sw", rrtmgp_cloud_optics_sw);
91  pp.query("rrtmgp_cloud_optics_lw", rrtmgp_cloud_optics_lw);
92 
93  // Append file names to path
98 
99  // Get dimensions from lookup data
100  if (ParallelDescriptor::IOProcessor()) {
101  auto ncf_sw = ncutils::NCFile::open(rrtmgp_coeffs_file_sw, NC_CLOBBER | NC_NETCDF4);
102  m_nswbands = ncf_sw.dim("bnd").len();
103  m_nswgpts = ncf_sw.dim("gpt").len();
104  ncf_sw.close();
105 
106  auto ncf_lw = ncutils::NCFile::open(rrtmgp_coeffs_file_lw, NC_CLOBBER | NC_NETCDF4);
107  m_nlwbands = ncf_lw.dim("bnd").len();
108  m_nlwgpts = ncf_lw.dim("gpt").len();
109  ncf_lw.close();
110  }
111  int ioproc = ParallelDescriptor::IOProcessorNumber(); // I/O rank
112  ParallelDescriptor::Bcast(&m_nswbands, 1, ioproc);
113  ParallelDescriptor::Bcast(&m_nlwbands, 1, ioproc);
114  ParallelDescriptor::Bcast(&m_nswgpts, 1, ioproc);
115  ParallelDescriptor::Bcast(&m_nlwgpts, 1, ioproc);
116 
117  // Output for user
118  if (lev == 0) {
119  Print() << "Radiation interface constructed:\n";
120  Print() << "========================================================\n";
121  Print() << "Coeff SW file: " << rrtmgp_coeffs_file_sw << "\n";
122  Print() << "Coeff LW file: " << rrtmgp_coeffs_file_lw << "\n";
123  Print() << "Cloud SW file: " << rrtmgp_cloud_optics_file_sw << "\n";
124  Print() << "Cloud LW file: " << rrtmgp_cloud_optics_file_lw << "\n";
125  Print() << "Number of short/longwave bands: "
126  << m_nswbands << " " << m_nlwbands << "\n";
127  Print() << "Number of short/longwave gauss points: "
128  << m_nswgpts << " " << m_nlwgpts << "\n";
129  Print() << "========================================================\n";
130  }
131 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real pp(amrex::Real y)
Definition: ERF_MicrophysicsUtils.H:233
std::string rrtmgp_coeffs_file_sw
Definition: ERF_Radiation.H:285
std::string rrtmgp_coeffs_sw
Definition: ERF_Radiation.H:281
int m_rad_freq_in_steps
Definition: ERF_Radiation.H:363
amrex::Real m_lon_cons
Definition: ERF_Radiation.H:270
int m_nswbands
Definition: ERF_Radiation.H:357
bool m_do_aerosol_rad
Definition: ERF_Radiation.H:325
std::string rrtmgp_coeffs_lw
Definition: ERF_Radiation.H:282
bool m_moist
Definition: ERF_Radiation.H:235
std::string rrtmgp_cloud_optics_file_lw
Definition: ERF_Radiation.H:288
bool m_lsm
Definition: ERF_Radiation.H:239
bool m_do_subcol_sampling
Definition: ERF_Radiation.H:366
bool m_rad_write_fluxes
Definition: ERF_Radiation.H:230
std::string rrtmgp_cloud_optics_file_sw
Definition: ERF_Radiation.H:287
amrex::Real m_orbital_mvelp
Definition: ERF_Radiation.H:345
amrex::Vector< amrex::Real > m_o3vmr
Definition: ERF_Radiation.H:299
std::string rrtmgp_cloud_optics_sw
Definition: ERF_Radiation.H:283
int m_nlwgpts
Definition: ERF_Radiation.H:360
std::string rrtmgp_cloud_optics_lw
Definition: ERF_Radiation.H:284
amrex::Real m_co2vmr
Definition: ERF_Radiation.H:298
bool m_extra_clnsky_diag
Definition: ERF_Radiation.H:328
amrex::Real m_lat_cons
Definition: ERF_Radiation.H:269
amrex::Real m_fixed_total_solar_irradiance
Definition: ERF_Radiation.H:350
amrex::Real m_o2vmr
Definition: ERF_Radiation.H:303
amrex::Real m_n2ovmr
Definition: ERF_Radiation.H:300
amrex::Real m_rad_t_sfc
Definition: ERF_Radiation.H:250
std::string rrtmgp_file_path
Definition: ERF_Radiation.H:280
int m_nlwbands
Definition: ERF_Radiation.H:358
amrex::Real m_orbital_eccen
Definition: ERF_Radiation.H:343
amrex::Real m_covmr
Definition: ERF_Radiation.H:301
std::string rrtmgp_coeffs_file_lw
Definition: ERF_Radiation.H:286
amrex::Real m_ch4vmr
Definition: ERF_Radiation.H:302
bool m_ice
Definition: ERF_Radiation.H:236
amrex::Real m_fixed_solar_zenith_angle
Definition: ERF_Radiation.H:354
amrex::Real m_orbital_obliq
Definition: ERF_Radiation.H:344
amrex::Real m_n2vmr
Definition: ERF_Radiation.H:304
bool m_fixed_orbital_year
Definition: ERF_Radiation.H:342
int m_nswgpts
Definition: ERF_Radiation.H:359
bool m_extra_clnclrsky_diag
Definition: ERF_Radiation.H:329
int m_orbital_year
Definition: ERF_Radiation.H:334
static NCFile open(const std::string &name, const int cmode=NC_NOWRITE)
Definition: ERF_NCInterface.cpp:707
LandSurfaceType lsm_type
Definition: ERF_DataStruct.H:1023
MoistureType moisture_type
Definition: ERF_DataStruct.H:1020
Here is the call graph for this function:

◆ ~Radiation()

Radiation::~Radiation ( )
inline
56  {
57  // Note that Kokkos is now finalized in main.cpp
58  }

Member Function Documentation

◆ alloc_buffers()

void Radiation::alloc_buffers ( )
204 {
205  // 1d size (m_ngas)
206  const Real* mol_weight_gas_p = m_mol_weight_gas.data();
207  const std::string* gas_names_p = m_gas_names.data();
208  m_gas_mol_weights = real1d_k("m_gas_mol_weights", m_ngas);
209  realHost1d_k m_gas_mol_weights_h("m_gas_mol_weights_h", m_ngas);
210  gas_names_offset.clear(); gas_names_offset.resize(m_ngas);
211  std::string* gas_names_offset_p = gas_names_offset.data();
212  Kokkos::parallel_for(Kokkos::RangePolicy<Kokkos::Serial>(0, m_ngas),
213  [&] (int igas)
214  {
215  m_gas_mol_weights_h(igas) = mol_weight_gas_p[igas];
216  gas_names_offset_p[igas] = gas_names_p[igas];
217  });
218  Kokkos::deep_copy(m_gas_mol_weights, m_gas_mol_weights_h);
219 
220  // 1d size (1 or nlay)
221  m_o3_size = m_o3vmr.size();
222  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(((m_o3_size==1) || (m_o3_size==m_nlay)),
223  "O3 VMR array must be length 1 or nlay");
224  Real* o3vmr_p = m_o3vmr.data();
225  o3_lay = real1d_k("o3_lay", m_o3_size);
226  realHost1d_k o3_lay_h("o3_lay_h", m_o3_size);
227  Kokkos::parallel_for(Kokkos::RangePolicy<Kokkos::Serial>(0, m_o3_size),
228  [&] (int io3)
229  {
230  o3_lay_h(io3) = o3vmr_p[io3];
231  });
232  Kokkos::deep_copy(o3_lay, o3_lay_h);
233 
234  // 1d size (ncol)
235  mu0 = real1d_k("mu0" , m_ncol);
236  sfc_alb_dir_vis = real1d_k("sfc_alb_dir_vis" , m_ncol);
237  sfc_alb_dir_nir = real1d_k("sfc_alb_dir_nir" , m_ncol);
238  sfc_alb_dif_vis = real1d_k("sfc_alb_dif_vis" , m_ncol);
239  sfc_alb_dif_nir = real1d_k("sfc_alb_dif_nir" , m_ncol);
240  sfc_flux_dir_vis = real1d_k("sfc_flux_dir_vis", m_ncol);
241  sfc_flux_dir_nir = real1d_k("sfc_flux_dir_nir", m_ncol);
242  sfc_flux_dif_vis = real1d_k("sfc_flux_dif_vis", m_ncol);
243  sfc_flux_dif_nir = real1d_k("sfc_flux_dif_nir", m_ncol);
244  lat = real1d_k("lat" , m_ncol);
245  lon = real1d_k("lon" , m_ncol);
246  sfc_emis = real1d_k("sfc_emis" , m_ncol);
247  t_sfc = real1d_k("t_sfc" , m_ncol);
248  lw_src = real1d_k("lw_src" , m_ncol);
249 
250  // 2d size (ncol, nlay)
251  r_lay = real2d_k("r_lay" , m_ncol, m_nlay);
252  p_lay = real2d_k("p_lay" , m_ncol, m_nlay);
253  t_lay = real2d_k("t_lay" , m_ncol, m_nlay);
254  z_del = real2d_k("z_del" , m_ncol, m_nlay);
255  qv_lay = real2d_k("qv" , m_ncol, m_nlay);
256  qc_lay = real2d_k("qc" , m_ncol, m_nlay);
257  qi_lay = real2d_k("qi" , m_ncol, m_nlay);
258  cldfrac_tot = real2d_k("cldfrac_tot" , m_ncol, m_nlay);
259  eff_radius_qc = real2d_k("eff_radius_qc", m_ncol, m_nlay);
260  eff_radius_qi = real2d_k("eff_radius_qi", m_ncol, m_nlay);
261  lwp = real2d_k("lwp" , m_ncol, m_nlay);
262  iwp = real2d_k("iwp" , m_ncol, m_nlay);
263  sw_heating = real2d_k("sw_heating" , m_ncol, m_nlay);
264  lw_heating = real2d_k("lw_heating" , m_ncol, m_nlay);
265  sw_clrsky_heating = real2d_k("sw_clrsky_heating", m_ncol, m_nlay);
266  lw_clrsky_heating = real2d_k("lw_clrsky_heating", m_ncol, m_nlay);
267 
268  // 2d size (ncol, nlay+1)
269  d_tint = real2d_k("d_tint" , m_ncol, m_nlay+1);
270  p_lev = real2d_k("p_lev" , m_ncol, m_nlay+1);
271  t_lev = real2d_k("t_lev" , m_ncol, m_nlay+1);
272  sw_flux_up = real2d_k("sw_flux_up" , m_ncol, m_nlay+1);
273  sw_flux_dn = real2d_k("sw_flux_dn" , m_ncol, m_nlay+1);
274  sw_flux_dn_dir = real2d_k("sw_flux_dn_dir" , m_ncol, m_nlay+1);
275  lw_flux_up = real2d_k("sw_flux_up" , m_ncol, m_nlay+1);
276  lw_flux_dn = real2d_k("sw_flux_dn" , m_ncol, m_nlay+1);
277  sw_clnclrsky_flux_up = real2d_k("sw_clnclrsky_flux_up" , m_ncol, m_nlay+1);
278  sw_clnclrsky_flux_dn = real2d_k("sw_clnclrsky_flux_dn" , m_ncol, m_nlay+1);
279  sw_clnclrsky_flux_dn_dir = real2d_k("sw_clnclrsky_flux_dn_dir", m_ncol, m_nlay+1);
280  sw_clrsky_flux_up = real2d_k("sw_clrsky_flux_up" , m_ncol, m_nlay+1);
281  sw_clrsky_flux_dn = real2d_k("sw_clrsky_flux_dn" , m_ncol, m_nlay+1);
282  sw_clrsky_flux_dn_dir = real2d_k("sw_clrsky_flux_dn_dir" , m_ncol, m_nlay+1);
283  sw_clnsky_flux_up = real2d_k("sw_clnsky_flux_up" , m_ncol, m_nlay+1);
284  sw_clnsky_flux_dn = real2d_k("sw_clnsky_flux_dn" , m_ncol, m_nlay+1);
285  sw_clnsky_flux_dn_dir = real2d_k("sw_clnsky_flux_dn_dir" , m_ncol, m_nlay+1);
286  lw_clnclrsky_flux_up = real2d_k("lw_clnclrsky_flux_up" , m_ncol, m_nlay+1);
287  lw_clnclrsky_flux_dn = real2d_k("lw_clnclrsky_flux_dn" , m_ncol, m_nlay+1);
288  lw_clrsky_flux_up = real2d_k("lw_clrsky_flux_up" , m_ncol, m_nlay+1);
289  lw_clrsky_flux_dn = real2d_k("lw_clrsky_flux_dn" , m_ncol, m_nlay+1);
290  lw_clnsky_flux_up = real2d_k("lw_clnsky_flux_up" , m_ncol, m_nlay+1);
291  lw_clnsky_flux_dn = real2d_k("lw_clnsky_flux_dn" , m_ncol, m_nlay+1);
292 
293  // 3d size (ncol, nlay+1, nswbands)
294  sw_bnd_flux_up = real3d_k("sw_bnd_flux_up" , m_ncol, m_nlay+1, m_nswbands);
295  sw_bnd_flux_dn = real3d_k("sw_bnd_flux_dn" , m_ncol, m_nlay+1, m_nswbands);
296  sw_bnd_flux_dir = real3d_k("sw_bnd_flux_dir", m_ncol, m_nlay+1, m_nswbands);
297  sw_bnd_flux_dif = real3d_k("sw_bnd_flux_dif", m_ncol, m_nlay+1, m_nswbands);
298 
299  // 3d size (ncol, nlay+1, nlwbands)
300  lw_bnd_flux_up = real3d_k("lw_bnd_flux_up" , m_ncol, m_nlay+1, m_nlwbands);
301  lw_bnd_flux_dn = real3d_k("lw_bnd_flux_up" , m_ncol, m_nlay+1, m_nlwbands);
302 
303  // 2d size (ncol, nswbands)
304  sfc_alb_dir = real2d_k("sfc_alb_dir", m_ncol, m_nswbands);
305  sfc_alb_dif = real2d_k("sfc_alb_dif", m_ncol, m_nswbands);
306 
307  // 2d size (ncol, nlwbands)
308  emis_sfc = real2d_k("emis_sfc", m_ncol, m_nlwbands);
309 
310  /*
311  // 3d size (ncol, nlay, n[sw,lw]bands)
312  aero_tau_sw = real3d_k("aero_tau_sw", m_ncol, m_nlay, m_nswbands);
313  aero_ssa_sw = real3d_k("aero_ssa_sw", m_ncol, m_nlay, m_nswbands);
314  aero_g_sw = real3d_k("aero_g_sw" , m_ncol, m_nlay, m_nswbands);
315  aero_tau_lw = real3d_k("aero_tau_lw", m_ncol, m_nlay, m_nlwbands);
316 
317  // 3d size (ncol, nlay, n[sw,lw]bnds)
318  cld_tau_sw_bnd = real3d_k("cld_tau_sw_bnd", m_ncol, m_nlay, m_nswbands);
319  cld_tau_lw_bnd = real3d_k("cld_tau_lw_bnd", m_ncol, m_nlay, m_nlwbands);
320 
321  // 3d size (ncol, nlay, n[sw,lw]gpts)
322  cld_tau_sw_gpt = real3d_k("cld_tau_sw_gpt", m_ncol, m_nlay, m_nswgpts);
323  cld_tau_lw_gpt = real3d_k("cld_tau_lw_gpt", m_ncol, m_nlay, m_nlwgpts);
324  */
325 }
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
Kokkos::View< RealT *, KokkosHostDevice > realHost1d_k
Definition: ERF_Kokkos.H:16
amrex::Real Real
Definition: ERF_ShocInterface.H:19
real3d_k sw_bnd_flux_dn
Definition: ERF_Radiation.H:432
real2d_k lw_flux_up
Definition: ERF_Radiation.H:412
int m_o3_size
Definition: ERF_Radiation.H:308
real3d_k sw_bnd_flux_dir
Definition: ERF_Radiation.H:433
real2d_k sw_clnsky_flux_dn
Definition: ERF_Radiation.H:421
real2d_k d_tint
Definition: ERF_Radiation.H:406
real2d_k lw_clnclrsky_flux_dn
Definition: ERF_Radiation.H:424
real2d_k lwp
Definition: ERF_Radiation.H:398
real1d_k lw_src
Definition: ERF_Radiation.H:385
real2d_k eff_radius_qi
Definition: ERF_Radiation.H:397
real1d_k m_gas_mol_weights
Definition: ERF_Radiation.H:309
real2d_k sw_heating
Definition: ERF_Radiation.H:400
real3d_k lw_bnd_flux_dn
Definition: ERF_Radiation.H:438
real3d_k sw_bnd_flux_up
Definition: ERF_Radiation.H:431
real2d_k qv_lay
Definition: ERF_Radiation.H:392
real2d_k sw_clnclrsky_flux_dn_dir
Definition: ERF_Radiation.H:416
real1d_k sfc_flux_dif_vis
Definition: ERF_Radiation.H:379
real2d_k lw_clnclrsky_flux_up
Definition: ERF_Radiation.H:423
real1d_k lat
Definition: ERF_Radiation.H:381
real2d_k qi_lay
Definition: ERF_Radiation.H:394
real2d_k sw_clrsky_flux_up
Definition: ERF_Radiation.H:417
real2d_k t_lev
Definition: ERF_Radiation.H:408
real1d_k o3_lay
Definition: ERF_Radiation.H:369
real1d_k sfc_alb_dif_vis
Definition: ERF_Radiation.H:375
real1d_k mu0
Definition: ERF_Radiation.H:372
real2d_k cldfrac_tot
Definition: ERF_Radiation.H:395
real1d_k sfc_flux_dir_nir
Definition: ERF_Radiation.H:378
real2d_k sw_clnclrsky_flux_up
Definition: ERF_Radiation.H:414
real2d_k sw_flux_up
Definition: ERF_Radiation.H:409
real3d_k sw_bnd_flux_dif
Definition: ERF_Radiation.H:434
real2d_k sfc_alb_dif
Definition: ERF_Radiation.H:442
real2d_k emis_sfc
Definition: ERF_Radiation.H:445
real1d_k sfc_alb_dif_nir
Definition: ERF_Radiation.H:376
real2d_k lw_clnsky_flux_dn
Definition: ERF_Radiation.H:428
real2d_k p_lay
Definition: ERF_Radiation.H:389
real2d_k r_lay
Definition: ERF_Radiation.H:388
int m_ncol
Definition: ERF_Radiation.H:318
real2d_k sw_flux_dn_dir
Definition: ERF_Radiation.H:411
const std::vector< amrex::Real > m_mol_weight_gas
Definition: ERF_Radiation.H:294
std::vector< std::string > gas_names_offset
Definition: ERF_Radiation.H:310
real2d_k sw_clrsky_flux_dn
Definition: ERF_Radiation.H:418
real2d_k sw_clnsky_flux_dn_dir
Definition: ERF_Radiation.H:422
real2d_k qc_lay
Definition: ERF_Radiation.H:393
real2d_k lw_clrsky_flux_up
Definition: ERF_Radiation.H:425
real2d_k sw_flux_dn
Definition: ERF_Radiation.H:410
real2d_k lw_clrsky_flux_dn
Definition: ERF_Radiation.H:426
real2d_k z_del
Definition: ERF_Radiation.H:391
real2d_k lw_flux_dn
Definition: ERF_Radiation.H:413
real2d_k sw_clrsky_flux_dn_dir
Definition: ERF_Radiation.H:419
real1d_k sfc_flux_dif_nir
Definition: ERF_Radiation.H:380
int m_ngas
Definition: ERF_Radiation.H:291
real2d_k lw_heating
Definition: ERF_Radiation.H:401
real1d_k sfc_alb_dir_nir
Definition: ERF_Radiation.H:374
real2d_k lw_clnsky_flux_up
Definition: ERF_Radiation.H:427
real2d_k sfc_alb_dir
Definition: ERF_Radiation.H:441
real1d_k lon
Definition: ERF_Radiation.H:382
real1d_k sfc_emis
Definition: ERF_Radiation.H:383
real2d_k sw_clnsky_flux_up
Definition: ERF_Radiation.H:420
int m_nlay
Definition: ERF_Radiation.H:319
real3d_k lw_bnd_flux_up
Definition: ERF_Radiation.H:437
real1d_k sfc_alb_dir_vis
Definition: ERF_Radiation.H:373
const std::vector< std::string > m_gas_names
Definition: ERF_Radiation.H:292
real2d_k sw_clrsky_heating
Definition: ERF_Radiation.H:402
real2d_k t_lay
Definition: ERF_Radiation.H:390
real1d_k sfc_flux_dir_vis
Definition: ERF_Radiation.H:377
real2d_k iwp
Definition: ERF_Radiation.H:399
real2d_k p_lev
Definition: ERF_Radiation.H:407
real2d_k lw_clrsky_heating
Definition: ERF_Radiation.H:403
real2d_k eff_radius_qc
Definition: ERF_Radiation.H:396
real1d_k t_sfc
Definition: ERF_Radiation.H:384
real2d_k sw_clnclrsky_flux_dn
Definition: ERF_Radiation.H:415

◆ dealloc_buffers()

void Radiation::dealloc_buffers ( )
329 {
330  // 1d size (m_ngas)
332 
333  // 1d size (1 or nlay)
334  o3_lay = real1d_k();
335 
336  // 1d size (ncol)
337  mu0 = real1d_k();
346  lat = real1d_k();
347  lon = real1d_k();
348  sfc_emis = real1d_k();
349  t_sfc = real1d_k();
350  lw_src = real1d_k();
351 
352  // 2d size (ncol, nlay)
353  r_lay = real2d_k();
354  p_lay = real2d_k();
355  t_lay = real2d_k();
356  z_del = real2d_k();
357  qv_lay = real2d_k();
358  qc_lay = real2d_k();
359  qi_lay = real2d_k();
360  cldfrac_tot = real2d_k();
363  lwp = real2d_k();
364  iwp = real2d_k();
365  sw_heating = real2d_k();
366  lw_heating = real2d_k();
369  sw_heating = real2d_k();
370  lw_heating = real2d_k();
373 
374  // 2d size (ncol, nlay+1)
375  d_tint = real2d_k();
376  p_lev = real2d_k();
377  t_lev = real2d_k();
378  sw_flux_up = real2d_k();
379  sw_flux_dn = real2d_k();
381  lw_flux_up = real2d_k();
382  lw_flux_dn = real2d_k();
398 
399  // 3d size (ncol, nlay+1, nswbands)
404 
405  // 3d size (ncol, nlay+1, nlwbands)
408 
409  // 2d size (ncol, nswbands)
410  sfc_alb_dir = real2d_k();
411  sfc_alb_dif = real2d_k();
412 
413  // 2d size (ncol, nlwbands)
414  emis_sfc = real2d_k();
415 
416  /*
417  // 3d size (ncol, nlay, n[sw,lw]bands)
418  aero_tau_sw = real3d_k();
419  aero_ssa_sw = real3d_k();
420  aero_g_sw = real3d_k();
421  aero_tau_lw = real3d_k();
422 
423  // 3d size (ncol, nlay, n[sw,lw]bnds)
424  cld_tau_sw_bnd = real3d_k();
425  cld_tau_lw_bnd = real3d_k();
426 
427  // 3d size (ncol, nlay, n[sw,lw]gpts)
428  cld_tau_sw_gpt = real3d_k();
429  cld_tau_lw_gpt = real3d_k();
430  */
431 }

◆ finalize_impl()

void Radiation::finalize_impl ( amrex::Vector< amrex::MultiFab * > &  lsm_output_ptrs)
1249 {
1250  // Finish rrtmgp
1251  m_gas_concs.reset();
1253 
1254  // Fill the AMReX MFs from Kokkos Views
1255  kokkos_buffers_to_mf(lsm_output_ptrs);
1256 
1257  // Write fluxes if requested
1259 
1260  // Fill output data for datalog before deallocating
1261  if (datalog_int > 0) {
1264 
1266  }
1267 
1268  // Deallocate the buffer arrays
1269  dealloc_buffers();
1270 }
int datalog_int
Definition: ERF_RadiationInterface.H:88
void dealloc_buffers()
Definition: ERF_Radiation.cpp:328
void populateDatalogMF()
Definition: ERF_Radiation.cpp:757
void kokkos_buffers_to_mf(amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
Definition: ERF_Radiation.cpp:606
GasConcsK< amrex::Real, layout_t, KokkosDefaultDevice > m_gas_concs
Definition: ERF_Radiation.H:312
void write_rrtmgp_fluxes()
Definition: ERF_Radiation.cpp:716
void rrtmgp_finalize()
Definition: ERF_RRTMGP_Interface.cpp:266
void compute_heating_rate(View1 const &flux_up, View2 const &flux_dn, View3 const &rho, View4 const &dz, View5 &heating_rate)
Definition: ERF_RRTMGP_Utils.H:81

Referenced by rad_run_impl().

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

◆ get_lsm_input_varnames()

virtual amrex::Vector<std::string> Radiation::get_lsm_input_varnames ( )
inlineoverridevirtual

Reimplemented from IRadiation.

182  {
183  return m_lsm_input_names;
184  }
amrex::Vector< std::string > m_lsm_input_names
Definition: ERF_Radiation.H:242

◆ get_lsm_output_varnames()

virtual amrex::Vector<std::string> Radiation::get_lsm_output_varnames ( )
inlineoverridevirtual

Reimplemented from IRadiation.

190  {
191  return m_lsm_output_names;
192  }
amrex::Vector< std::string > m_lsm_output_names
Definition: ERF_Radiation.H:247

◆ Init()

virtual void Radiation::Init ( const amrex::Geometry &  geom,
const amrex::BoxArray &  ba,
amrex::MultiFab *  cons_in 
)
inlineoverridevirtual

Implements IRadiation.

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  };
amrex::Vector< int > m_col_offsets
Definition: ERF_Radiation.H:322

◆ initialize_impl()

void Radiation::initialize_impl ( )
982 {
983  // Call API to initialize
988 }
void rrtmgp_initialize(gas_concs_t &gas_concs_k, 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:229

Referenced by rad_run_impl().

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

◆ kokkos_buffers_to_mf()

void Radiation::kokkos_buffers_to_mf ( amrex::Vector< amrex::MultiFab * > &  lsm_output_ptrs)
607 {
608  // Heating rate, fluxes, zenith, lsm ptrs
609  Vector<real2d_k> rrtmgp_out_vars = {sw_flux_dn, lw_flux_dn};
610 
611  // Expose for device
612  auto sw_heating_d = sw_heating;
613  auto lw_heating_d = lw_heating;
614  auto p_lay_d = p_lay;
615  auto sfc_flux_dir_vis_d = sfc_flux_dir_vis;
616  auto sfc_flux_dir_nir_d = sfc_flux_dir_nir;
617  auto sfc_flux_dif_vis_d = sfc_flux_dif_vis;
618  auto sfc_flux_dif_nir_d = sfc_flux_dif_nir;
619  auto sw_flux_up_d = sw_flux_up;
620  auto sw_flux_dn_d = sw_flux_dn;
621  auto lw_flux_up_d = lw_flux_up;
622  auto lw_flux_dn_d = lw_flux_dn;
623  auto mu0_d = mu0;
624 
625  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
626  const auto& vbx = mfi.validbox();
627  const auto& sbx = makeSlab(vbx,2,vbx.smallEnd(2));
628  const int nx = vbx.length(0);
629  const int imin = vbx.smallEnd(0);
630  const int jmin = vbx.smallEnd(1);
631  const int offset = m_col_offsets[mfi.index()];
632  const Array4<Real>& q_arr = m_qheating_rates->array(mfi);
633  const Array4<Real>& f_arr = m_rad_fluxes->array(mfi);
634  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
635  {
636  // map [i,j,k] 0-based to [icol, ilay] 0-based
637  const int icol = (j-jmin)*nx + (i-imin) + offset;
638  const int ilay = k;
639 
640  // Temperature heating rate for SW and LW
641  q_arr(i,j,k,0) = sw_heating_d(icol,ilay);
642  q_arr(i,j,k,1) = lw_heating_d(icol,ilay);
643 
644  // Convert the dT/dz to dTheta/dz
645  Real iexner = 1./getExnergivenP(Real(p_lay_d(icol,ilay)), R_d/Cp_d);
646  q_arr(i,j,k,0) *= iexner;
647  q_arr(i,j,k,1) *= iexner;
648 
649  // Populate the fluxes
650  f_arr(i,j,k,0) = sw_flux_up_d(icol,ilay);
651  f_arr(i,j,k,1) = sw_flux_dn_d(icol,ilay);
652  f_arr(i,j,k,2) = lw_flux_up_d(icol,ilay);
653  f_arr(i,j,k,3) = lw_flux_dn_d(icol,ilay);
654  });
655  if (m_lsm_fluxes) {
656  const Array4<Real>& lsm_arr = m_lsm_fluxes->array(mfi);
657  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
658  {
659  // map [i,j,k] 0-based to [icol, ilay] 0-based
660  const int icol = (j-jmin)*nx + (i-imin) + offset;
661 
662  // SW fluxes for LSM
663  lsm_arr(i,j,k,0) = sfc_flux_dir_vis_d(icol);
664  lsm_arr(i,j,k,1) = sfc_flux_dir_nir_d(icol);
665  lsm_arr(i,j,k,2) = sfc_flux_dif_vis_d(icol);
666  lsm_arr(i,j,k,3) = sfc_flux_dif_nir_d(icol);
667 
668  // Net SW flux for LSM
669  lsm_arr(i,j,k,4) = sfc_flux_dir_vis_d(icol) + sfc_flux_dir_nir_d(icol)
670  + sfc_flux_dif_vis_d(icol) + sfc_flux_dif_nir_d(icol);
671 
672  // LW flux for LSM (at bottom surface)
673  lsm_arr(i,j,k,5) = lw_flux_dn_d(icol,0);
674  });
675  }
676  if (m_lsm_zenith) {
677  const Array4<Real>& lsm_zenith_arr = m_lsm_zenith->array(mfi);
678  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
679  {
680  // map [i,j,k] 0-based to [icol, ilay] 0-based
681  const int icol = (j-jmin)*nx + (i-imin) + offset;
682 
683  // export cosine zenith angle for LSM
684  lsm_zenith_arr(i,j,k) = mu0_d(icol);
685  });
686  }
687  for (int ivar(0); ivar<lsm_output_ptrs.size(); ivar++) {
688  if (lsm_output_ptrs[ivar]) {
689  const Array4<Real>& lsm_out_arr = lsm_output_ptrs[ivar]->array(mfi);
690  if (ivar==0) {
691  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
692  {
693  // map [i,j,k] 0-based to [icol, ilay] 0-based
694  const int icol = (j-jmin)*nx + (i-imin) + offset;
695 
696  // export the desired variable at surface
697  lsm_out_arr(i,j,k) = mu0_d(icol);
698  });
699  } else {
700  auto rrtmgp_for_fill = rrtmgp_out_vars[ivar-1];
701  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
702  {
703  // map [i,j,k] 0-based to [icol, ilay] 0-based
704  const int icol = (j-jmin)*nx + (i-imin) + offset;
705 
706  // export the desired variable at surface
707  lsm_out_arr(i,j,k) = rrtmgp_for_fill(icol,0);
708  });
709  } // ivar
710  } // valid ptr
711  } // ivar
712  }// mfi
713 }
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:141
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
amrex::MultiFab * m_cons_in
Definition: ERF_Radiation.H:253
amrex::MultiFab * m_rad_fluxes
Definition: ERF_Radiation.H:259
amrex::MultiFab * m_lsm_fluxes
Definition: ERF_Radiation.H:273
amrex::MultiFab * m_qheating_rates
Definition: ERF_Radiation.H:256
amrex::MultiFab * m_lsm_zenith
Definition: ERF_Radiation.H:274
Here is the call graph for this function:

◆ mf_to_kokkos_buffers()

void Radiation::mf_to_kokkos_buffers ( amrex::Vector< amrex::MultiFab * > &  lsm_input_ptrs)
436 {
437  // Expose for device
438  auto r_lay_d = r_lay;
439  auto p_lay_d = p_lay;
440  auto t_lay_d = t_lay;
441  auto z_del_d = z_del;
442  auto qv_lay_d = qv_lay;
443  auto qc_lay_d = qc_lay;
444  auto qi_lay_d = qi_lay;
445  auto cldfrac_tot_d = cldfrac_tot;
446  auto lwp_d = lwp;
447  auto iwp_d = iwp;
448  auto eff_radius_qc_d = eff_radius_qc;
449  auto eff_radius_qi_d = eff_radius_qi;
450  auto p_lev_d = p_lev;
451  auto t_lev_d = t_lev;
452  auto lat_d = lat;
453  auto lon_d = lon;
454  auto t_sfc_d = t_sfc;
455 
456  bool moist = m_moist;
457  bool ice = m_ice;
458  const bool has_lsm = m_lsm;
459  const bool has_lat = m_lat;
460  const bool has_lon = m_lon;
461  int ncol = m_ncol;
462  int nlay = m_nlay;
463  Real dz = m_geom.CellSize(2);
464  Real cons_lat = m_lat_cons;
465  Real cons_lon = m_lon_cons;
466  Real rad_t_sfc = m_rad_t_sfc;
467  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
468  const auto& vbx = mfi.validbox();
469  const int nx = vbx.length(0);
470  const int imin = vbx.smallEnd(0);
471  const int jmin = vbx.smallEnd(1);
472  const int offset = m_col_offsets[mfi.index()];
473  const Array4<const Real>& cons_arr = m_cons_in->const_array(mfi);
474  const Array4<const Real>& z_arr = (m_z_phys) ? m_z_phys->const_array(mfi) :
475  Array4<const Real>{};
476  const Array4<const Real>& lat_arr = (m_lat) ? m_lat->const_array(mfi) :
477  Array4<const Real>{};
478  const Array4<const Real>& lon_arr = (m_lon) ? m_lon->const_array(mfi) :
479  Array4<const Real>{};
480  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
481  {
482  // map [i,j,k] 0-based to [icol, ilay] 0-based
483  const int icol = (j-jmin)*nx + (i-imin) + offset;
484  const int ilay = k;
485 
486  // EOS input (at CC)
487  Real r = cons_arr(i,j,k,Rho_comp);
488  Real rt = cons_arr(i,j,k,RhoTheta_comp);
489  Real qv = (moist) ? std::max(cons_arr(i,j,k,RhoQ1_comp)/r,0.0) : 0.0;
490  Real qc = (moist) ? std::max(cons_arr(i,j,k,RhoQ2_comp)/r,0.0) : 0.0;
491  Real qi = (ice) ? std::max(cons_arr(i,j,k,RhoQ3_comp)/r,0.0) : 0.0;
492 
493  // EOS avg to z-face
494  Real r_lo = cons_arr(i,j,k-1,Rho_comp);
495  Real rt_lo = cons_arr(i,j,k-1,RhoTheta_comp);
496  Real qv_lo = (moist) ? cons_arr(i,j,k-1,RhoQ1_comp)/r_lo : 0.0;
497  Real r_avg = 0.5 * (r + r_lo);
498  Real rt_avg = 0.5 * (rt + rt_lo);
499  Real qv_avg = 0.5 * (qv + qv_lo);
500 
501  // Views at CC
502  r_lay_d(icol,ilay) = r;
503  p_lay_d(icol,ilay) = getPgivenRTh(rt, qv);
504  t_lay_d(icol,ilay) = getTgivenRandRTh(r, rt, qv);
505  z_del_d(icol,ilay) = (z_arr) ? 0.25 * ( (z_arr(i ,j ,k+1) - z_arr(i ,j ,k))
506  + (z_arr(i+1,j ,k+1) - z_arr(i+1,j ,k))
507  + (z_arr(i ,j+1,k+1) - z_arr(i ,j+1,k))
508  + (z_arr(i+1,j+1,k+1) - z_arr(i+1,j+1,k)) ) : dz;
509  qv_lay_d(icol,ilay) = qv;
510  qc_lay_d(icol,ilay) = qc;
511  qi_lay_d(icol,ilay) = qi;
512  cldfrac_tot_d(icol,ilay) = ((qc+qi)>0.0) ? 1. : 0.;
513 
514  // NOTE: These are populated in 'mixing_ratio_to_cloud_mass'
515  lwp_d(icol,ilay) = 0.0;
516  iwp_d(icol,ilay) = 0.0;
517 
518  // NOTE: These would be populated from P3 (we use the constants in p3_main_impl.hpp)
519  eff_radius_qc_d(icol,ilay) = (qc>0.0) ? 10.0e-6 : 0.0;
520  eff_radius_qi_d(icol,ilay) = (qi>0.0) ? 25.0e-6 : 0.0;
521 
522  // Buffers on z-faces (nlay+1)
523  p_lev_d(icol,ilay) = getPgivenRTh(rt_avg, qv_avg);
524  t_lev_d(icol,ilay) = getTgivenRandRTh(r_avg, rt_avg, qv_avg);
525  if (ilay==(nlay-1)) {
526  Real r_hi = cons_arr(i,j,k+1,Rho_comp);
527  Real rt_hi = cons_arr(i,j,k+1,RhoTheta_comp);
528  Real qv_hi = (moist) ? std::max(cons_arr(i,j,k+1,RhoQ1_comp)/r_hi,0.0) : 0.0;
529  r_avg = 0.5 * (r + r_hi);
530  rt_avg = 0.5 * (rt + rt_hi);
531  qv_avg = 0.5 * (qv + qv_hi);
532  p_lev_d(icol,ilay+1) = getPgivenRTh(rt_avg, qv_avg);
533  t_lev_d(icol,ilay+1) = getTgivenRandRTh(r_avg, rt_avg, qv_avg);
534  }
535 
536  // 1D data structures
537  if (k==0) {
538  lat_d(icol) = (has_lat) ? lat_arr(i,j,0) : cons_lat;
539  lon_d(icol) = (has_lon) ? lon_arr(i,j,0) : cons_lon;
540  }
541 
542  });
543  } // mfi
544 
545  // Populate vars LSM would provide
546  if (!has_lsm) {
547  // Parsed surface temp
548  Kokkos::deep_copy(t_sfc, rad_t_sfc);
549 
550  // EAMXX dummy atmos constants
551  Kokkos::deep_copy(sfc_alb_dir_vis, 0.06);
552  Kokkos::deep_copy(sfc_alb_dir_nir, 0.06);
553  Kokkos::deep_copy(sfc_alb_dif_vis, 0.06);
554  Kokkos::deep_copy(sfc_alb_dif_nir, 0.06);
555 
556  // AML NOTE: These are not used in current EAMXX, I've left
557  // the code to plug into these if we need it.
558  //
559  // Current EAMXX constants
560  Kokkos::deep_copy(emis_sfc, 0.98);
561  Kokkos::deep_copy(lw_src , 0.0 );
562  } else {
563  Vector<real1d_k> rrtmgp_in_vars = {t_sfc, sfc_emis,
566  Vector<Real> rrtmgp_default_vals = {rad_t_sfc, 0.98,
567  0.06, 0.06,
568  0.06, 0.06};
569  for (int ivar(0); ivar<lsm_input_ptrs.size(); ivar++) {
570  auto rrtmgp_to_fill = rrtmgp_in_vars[ivar];
571  if (!lsm_input_ptrs[ivar]) {
572  Kokkos::deep_copy(rrtmgp_to_fill, rrtmgp_default_vals[ivar]);
573  } else {
574  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
575  const auto& vbx = mfi.validbox();
576  const auto& sbx = makeSlab(vbx,2,vbx.smallEnd(2));
577  const int nx = vbx.length(0);
578  const int imin = vbx.smallEnd(0);
579  const int jmin = vbx.smallEnd(1);
580  const int offset = m_col_offsets[mfi.index()];
581  const Array4<const Real>& lsm_in_arr = lsm_input_ptrs[ivar]->const_array(mfi);
582  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
583  {
584  // map [i,j,k] 0-based to [icol, ilay] 0-based
585  const int icol = (j-jmin)*nx + (i-imin) + offset;
586 
587  // 2D mf and 1D view
588  rrtmgp_to_fill(icol) = lsm_in_arr(i,j,k);
589  });
590  } //mfi
591  } // valid lsm ptr
592  } // ivar
593  Kokkos::deep_copy(lw_src, 0.0 );
594  } // have lsm
595 
596  // Enforce consistency between t_sfc and t_lev at bottom surface
597  Kokkos::parallel_for(Kokkos::RangePolicy(0, ncol),
598  KOKKOS_LAMBDA (int icol)
599  {
600  t_lev_d(icol,0) = t_sfc_d(icol);
601  });
602 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=0.)
Definition: ERF_EOS.H:81
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::MultiFab * m_lon
Definition: ERF_Radiation.H:266
amrex::MultiFab * m_z_phys
Definition: ERF_Radiation.H:262
amrex::Geometry m_geom
Definition: ERF_Radiation.H:221
amrex::MultiFab * m_lat
Definition: ERF_Radiation.H:265
@ qv
Definition: ERF_Kessler.H:28
@ qc
Definition: ERF_SatAdj.H:36
Here is the call graph for this function:

◆ populateDatalogMF()

void Radiation::populateDatalogMF ( )
758 {
759  // Expose for device
760  auto sw_flux_up_d = sw_flux_up;
761  auto sw_flux_dn_d = sw_flux_dn;
762  auto sw_flux_dn_dir_d = sw_flux_dn_dir;
763  auto lw_flux_up_d = lw_flux_up;
764  auto lw_flux_dn_d = lw_flux_dn;
765  auto mu0_d = mu0;
766  auto sw_clrsky_heating_d = sw_clrsky_heating;
767  auto lw_clrsky_heating_d = lw_clrsky_heating;
768 
769  auto sw_clrsky_flux_up_d = sw_clrsky_flux_up;
770  auto sw_clrsky_flux_dn_d = sw_clrsky_flux_dn;
771  auto sw_clrsky_flux_dn_dir_d = sw_clrsky_flux_dn_dir;
772  auto lw_clrsky_flux_up_d = lw_clrsky_flux_up;
773  auto lw_clrsky_flux_dn_d = lw_clrsky_flux_dn;
774  auto sw_clnsky_flux_up_d = sw_clnsky_flux_up;
775  auto sw_clnsky_flux_dn_d = sw_clnsky_flux_dn;
776  auto sw_clnsky_flux_dn_dir_d = sw_clnsky_flux_dn_dir;
777  auto lw_clnsky_flux_up_d = lw_clnsky_flux_up;
778  auto lw_clnsky_flux_dn_d = lw_clnsky_flux_dn;
779  auto sw_clnclrsky_flux_up_d = sw_clnclrsky_flux_up;
780  auto sw_clnclrsky_flux_dn_d = sw_clnclrsky_flux_dn;
781  auto sw_clnclrsky_flux_dn_dir_d = sw_clnclrsky_flux_dn_dir;
782  auto lw_clnclrsky_flux_up_d = lw_clnclrsky_flux_up;
783  auto lw_clnclrsky_flux_dn_d = lw_clnclrsky_flux_dn;
784 
785  auto extra_clnsky_diag = m_extra_clnsky_diag;
786  auto extra_clnclrsky_diag = m_extra_clnclrsky_diag;
787 
788  for (MFIter mfi(datalog_mf); mfi.isValid(); ++mfi) {
789  const auto& vbx = mfi.validbox();
790  const int nx = vbx.length(0);
791  const int imin = vbx.smallEnd(0);
792  const int jmin = vbx.smallEnd(1);
793  const int offset = m_col_offsets[mfi.index()];
794  const Array4<Real>& dst_arr = datalog_mf.array(mfi);
795  const Array4<Real>& q_arr = m_qheating_rates->array(mfi);
796  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
797  {
798  // map [i,j,k] 0-based to [icol, ilay] 0-based
799  const int icol = (j-jmin)*nx + (i-imin) + offset;
800  const int ilay = k;
801 
802  dst_arr(i,j,k,0) = q_arr(i, j, k, 0);
803  dst_arr(i,j,k,1) = q_arr(i, j, k, 1);
804 
805  // SW and LW fluxes
806  dst_arr(i,j,k,2) = sw_flux_up_d(icol,ilay);
807  dst_arr(i,j,k,3) = sw_flux_dn_d(icol,ilay);
808  dst_arr(i,j,k,4) = sw_flux_dn_dir_d(icol,ilay);
809  dst_arr(i,j,k,5) = lw_flux_up_d(icol,ilay);
810  dst_arr(i,j,k,6) = lw_flux_dn_d(icol,ilay);
811 
812  // Cosine zenith angle
813  dst_arr(i,j,k,7) = mu0_d(icol);
814 
815  // Clear sky heating rates and fluxes:
816  dst_arr(i,j,k,8) = sw_clrsky_heating_d(icol, ilay);
817  dst_arr(i,j,k,9) = lw_clrsky_heating_d(icol, ilay);
818 
819  dst_arr(i,j,k,10) = sw_clrsky_flux_up_d(icol,ilay);
820  dst_arr(i,j,k,11) = sw_clrsky_flux_dn_d(icol,ilay);
821  dst_arr(i,j,k,12) = sw_clrsky_flux_dn_dir_d(icol,ilay);
822  dst_arr(i,j,k,13) = lw_clrsky_flux_up_d(icol,ilay);
823  dst_arr(i,j,k,14) = lw_clrsky_flux_dn_d(icol,ilay);
824 
825  // Clean sky fluxes:
826  if (extra_clnsky_diag) {
827  dst_arr(i,j,k,15) = sw_clnsky_flux_up_d(icol,ilay);
828  dst_arr(i,j,k,16) = sw_clnsky_flux_dn_d(icol,ilay);
829  dst_arr(i,j,k,17) = sw_clnsky_flux_dn_dir_d(icol,ilay);
830  dst_arr(i,j,k,18) = lw_clnsky_flux_up_d(icol,ilay);
831  dst_arr(i,j,k,19) = lw_clnsky_flux_dn_d(icol,ilay);
832  }
833 
834  // Clean-clear sky fluxes:
835  if (extra_clnclrsky_diag) {
836  dst_arr(i,j,k,20) = sw_clnclrsky_flux_up_d(icol,ilay);
837  dst_arr(i,j,k,21) = sw_clnclrsky_flux_dn_d(icol,ilay);
838  dst_arr(i,j,k,22) = sw_clnclrsky_flux_dn_dir_d(icol,ilay);
839  dst_arr(i,j,k,23) = lw_clnclrsky_flux_up_d(icol,ilay);
840  dst_arr(i,j,k,24) = lw_clnclrsky_flux_dn_d(icol,ilay);
841  }
842  });
843  }
844 }
amrex::MultiFab datalog_mf
Definition: ERF_Radiation.H:277
Here is the call graph for this function:

◆ rad_run_impl()

void Radiation::rad_run_impl ( amrex::Vector< amrex::MultiFab * > &  lsm_output_ptrs)
inline
167  {
168  if (m_update_rad) {
169  amrex::Print() << "Radiation advancing level " << m_lev << " at (YY-MM-DD SS) " << m_orbital_year << '-'
170  << m_orbital_mon << '-' << m_orbital_day << ' ' << m_orbital_sec << " ...";
171  this->initialize_impl();
172  this->run_impl();
173  this->finalize_impl(lsm_output_ptrs);
174  amrex::Print() << "DONE\n";
175  }
176  }
void initialize_impl()
Definition: ERF_Radiation.cpp:981
void run_impl()
Definition: ERF_Radiation.cpp:992
int m_orbital_mon
Definition: ERF_Radiation.H:335
void finalize_impl(amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
Definition: ERF_Radiation.cpp:1248
int m_orbital_sec
Definition: ERF_Radiation.H:337
int m_lev
Definition: ERF_Radiation.H:209
bool m_update_rad
Definition: ERF_Radiation.H:227
int m_orbital_day
Definition: ERF_Radiation.H:336

Referenced by Run().

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

◆ Run()

virtual void Radiation::Run ( 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::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 
)
inlineoverridevirtual

Implements IRadiation.

106  {
107  set_grids(level, step, time, dt, ba, geom,
108  cons_in, lsm_fluxes, lsm_zenith,
109  lsm_input_ptrs, qheating_rates,
110  rad_fluxes, z_phys, lat_ptr, lon_ptr);
111  rad_run_impl(lsm_output_ptrs);
112  }
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::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 rad_run_impl(amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
Definition: ERF_Radiation.H:166
Here is the call graph for this function:

◆ run_impl()

void Radiation::run_impl ( )
993 {
994  // Local copies
995  const auto ncol = m_ncol;
996  const auto nlay = m_nlay;
997  const auto nswbands = m_nswbands;
998 
999  // Compute orbital parameters; these are used both for computing
1000  // the solar zenith angle and also for computing total solar
1001  // irradiance scaling (tsi_scaling).
1002  double obliqr, lambm0, mvelpp;
1003  int orbital_year = m_orbital_year;
1004  double eccen = m_orbital_eccen;
1005  double obliq = m_orbital_obliq;
1006  double mvelp = m_orbital_mvelp;
1007  if (eccen >= 0 && obliq >= 0 && mvelp >= 0) {
1008  // fixed orbital parameters forced with orbital_year == ORB_UNDEF_INT
1009  orbital_year = ORB_UNDEF_INT;
1010  }
1011  orbital_params(orbital_year, eccen, obliq,
1012  mvelp, obliqr, lambm0, mvelpp);
1013 
1014  // Use the orbital parameters to calculate the solar declination and eccentricity factor
1015  double delta, eccf;
1016  // Want day + fraction; calday 1 == Jan 1 0Z
1017  static constexpr double dpy[] = {0.0 , 31.0, 59.0, 90.0, 120.0, 151.0,
1018  181.0, 212.0, 243.0, 273.0, 304.0, 334.0};
1019  bool leap = (m_orbital_year % 4 == 0 && (!(m_orbital_year % 100 == 0) || (m_orbital_year % 400 == 0))) ? true : false;
1020  double calday = dpy[m_orbital_mon-1] + (m_orbital_day-1.0) + m_orbital_sec/86400.0;
1021  // add extra day if leap year
1022  if (leap) { calday += 1.0; }
1023  orbital_decl(calday, eccen, mvelpp, lambm0, obliqr, delta, eccf);
1024 
1025  // Overwrite eccf if using a fixed solar constant.
1026  auto fixed_total_solar_irradiance = m_fixed_total_solar_irradiance;
1027  if (fixed_total_solar_irradiance >= 0){
1028  eccf = fixed_total_solar_irradiance/1360.9;
1029  }
1030 
1031  // Precompute volume mixing ratio (VMR) for all gases
1032  //
1033  // H2O is obtained from qv.
1034  // O3 may be a constant or a 1D vector
1035  // All other comps are set to constants for now
1036  for (int igas(0); igas < m_ngas; ++igas) {
1037  auto tmp2d = Kokkos::View<RealT**,layout_t,KokkosDefaultMem>("tmp2d", ncol, nlay);
1038  auto name = m_gas_names[igas];
1039  auto gas_mol_weight = m_mol_weight_gas[igas];
1040  if (name == "H2O") {
1041  auto qv_lay_d = qv_lay;
1042  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1043  KOKKOS_LAMBDA (int icol, int ilay)
1044  {
1045  tmp2d(icol,ilay) = qv_lay_d(icol,ilay) * mwdair/gas_mol_weight;
1046  });
1047  } else if (name == "CO2") {
1048  Kokkos::deep_copy(tmp2d, m_co2vmr);
1049  } else if (name == "O3") {
1050  if (m_o3_size==1) {
1051  Kokkos::deep_copy(tmp2d, m_o3vmr[0] );
1052  } else {
1053  auto o3_lay_d = o3_lay;
1054  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1055  KOKKOS_LAMBDA (int icol, int ilay)
1056  {
1057  tmp2d(icol,ilay) = o3_lay_d(ilay);
1058  });
1059  }
1060  } else if (name == "N2O") {
1061  Kokkos::deep_copy(tmp2d, m_n2ovmr);
1062  } else if (name == "CO") {
1063  Kokkos::deep_copy(tmp2d, m_covmr );
1064  } else if (name == "CH4") {
1065  Kokkos::deep_copy(tmp2d, m_ch4vmr);
1066  } else if (name == "O2") {
1067  Kokkos::deep_copy(tmp2d, m_o2vmr );
1068  } else if (name == "N2") {
1069  Kokkos::deep_copy(tmp2d, m_n2vmr );
1070  } else {
1071  Abort("Radiation: Unknown gas component.");
1072  }
1073 
1074  // Populate GasConcs object
1075  m_gas_concs.set_vmr(name, tmp2d);
1076  Kokkos::fence();
1077  }
1078 
1079  // Populate mu0 1D array
1080  // This must be done on HOST and copied to device.
1081  auto h_mu0 = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), mu0);
1082  if (m_fixed_solar_zenith_angle > 0) {
1083  Kokkos::deep_copy(h_mu0, m_fixed_solar_zenith_angle);
1084  } else {
1085  auto h_lat = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), lat);
1086  auto h_lon = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), lon);
1087  double dt = double(m_dt);
1088  auto rad_freq_in_steps = m_rad_freq_in_steps;
1089  Kokkos::parallel_for(Kokkos::RangePolicy<Kokkos::Serial>(0, ncol),
1090  [&] (int icol)
1091  {
1092  // Convert lat/lon to radians
1093  double lat_col = h_lat(icol)*PI/180.0;
1094  double lon_col = h_lon(icol)*PI/180.0;
1095  double lcalday = calday;
1096  double ldelta = delta;
1097  h_mu0(icol) = Real(orbital_cos_zenith(lcalday, lat_col, lon_col, ldelta, rad_freq_in_steps * dt));
1098  });
1099  }
1100  Kokkos::deep_copy(mu0, h_mu0);
1101 
1102  // Compute layer cloud mass per unit area (populates lwp/iwp)
1105 
1106  // Convert to g/m2 (needed by RRTMGP)
1107  auto lwp_d = lwp;
1108  auto iwp_d = iwp;
1109  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1110  KOKKOS_LAMBDA (int icol, int ilay)
1111  {
1112  lwp_d(icol,ilay) *= 1.e3;
1113  iwp_d(icol,ilay) *= 1.e3;
1114  });
1115 
1116  // Expand surface_albedos along nswbands.
1117  // This is needed since rrtmgp require band-by-band.
1122  // Run RRTMGP driver
1124  p_lay, t_lay,
1125  p_lev, t_lev,
1126  m_gas_concs,
1128  t_sfc, emis_sfc, lw_src,
1144 
1145 #if 0
1146  // UNIT TEST
1147  //================================================================================
1148  Kokkos::deep_copy(mu0, 0.86);
1149  Kokkos::deep_copy(sfc_alb_dir_vis, 0.06);
1150  Kokkos::deep_copy(sfc_alb_dir_nir, 0.06);
1151  Kokkos::deep_copy(sfc_alb_dif_vis, 0.06);
1152  Kokkos::deep_copy(sfc_alb_dif_nir, 0.06);
1153 
1154  Kokkos::deep_copy(aero_tau_sw, 0.0);
1155  Kokkos::deep_copy(aero_ssa_sw, 0.0);
1156  Kokkos::deep_copy(aero_g_sw , 0.0);
1157  Kokkos::deep_copy(aero_tau_lw, 0.0);
1158 
1159  // Generate some fake liquid and ice water data. We pick values to be midway between
1160  // the min and max of the valid lookup table values for effective radii
1161  real rel_val = 0.5 * (rrtmgp::cloud_optics_sw_k->get_min_radius_liq()
1162  + rrtmgp::cloud_optics_sw_k->get_max_radius_liq());
1163  real rei_val = 0.5 * (rrtmgp::cloud_optics_sw_k->get_min_radius_ice()
1164  + rrtmgp::cloud_optics_sw_k->get_max_radius_ice());
1165 
1166  // Restrict clouds to troposphere (> 100 hPa = 100*100 Pa) and not very close to the ground (< 900 hPa), and
1167  // put them in 2/3 of the columns since that's roughly the total cloudiness of earth.
1168  // Set sane values for liquid and ice water path.
1169  // NOTE: these "sane" values are in g/m2!
1170  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1171  KOKKOS_LAMBDA (int icol, int ilay)
1172  {
1173  cldfrac_tot(icol,ilay) = (p_lay(icol,ilay) > 100. * 100.) &&
1174  (p_lay(icol,ilay) < 900. * 100.) &&
1175  (icol%3 != 0);
1176  // Ice and liquid will overlap in a few layers
1177  lwp(icol,ilay) = (cldfrac_tot(icol,ilay) && t_lay(icol,ilay) > 263.) ? 10. : 0.;
1178  iwp(icol,ilay) = (cldfrac_tot(icol,ilay) && t_lay(icol,ilay) < 273.) ? 10. : 0.;
1179  eff_radius_qc(icol,ilay) = (lwp(icol,ilay) > 0.) ? rel_val : 0.;
1180  eff_radius_qi(icol,ilay) = (iwp(icol,ilay) > 0.) ? rei_val : 0.;
1181  });
1182 
1187 
1189  p_lay, t_lay,
1190  p_lev, t_lev,
1191  m_gas_concs,
1193  t_sfc, emis_sfc, lw_src,
1208  1.0, false, false);
1209  //================================================================================
1210 #endif
1211 
1212  // Update heating tendency
1215 
1216  /*
1217  // AML DEBUG
1218  Kokkos::parallel_for(nlay+1, KOKKOS_LAMBDA (int ilay)
1219  {
1220  printf("Fluxes: %i %e %e %e %e %e\n",ilay,
1221  sw_flux_up(5,ilay), sw_flux_dn(5,ilay), sw_flux_dn_dir(5,ilay),
1222  lw_flux_up(5,ilay), lw_flux_dn(5,ilay));
1223  });
1224  Kokkos::parallel_for(nlay, KOKKOS_LAMBDA (int ilay)
1225  {
1226  printf("Heating Rate: %i %e %e\n",ilay,sw_heating(5,ilay),lw_heating(5,ilay));
1227  });
1228  */
1229 
1230  // Compute surface fluxes
1231  const int kbot = 0;
1232  auto sw_bnd_flux_dif_d = sw_bnd_flux_dif;
1233  auto sw_bnd_flux_dn_d = sw_bnd_flux_dn;
1234  auto sw_bnd_flux_dir_d = sw_bnd_flux_dir;
1235  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<3>>({0, 0, 0}, {ncol, nlay+1, nswbands}),
1236  KOKKOS_LAMBDA (int icol, int ilay, int ibnd)
1237  {
1238  sw_bnd_flux_dif_d(icol,ilay,ibnd) = sw_bnd_flux_dn_d(icol,ilay,ibnd) - sw_bnd_flux_dir_d(icol,ilay,ibnd);
1239  });
1240  rrtmgp::compute_broadband_surface_fluxes(ncol, kbot, nswbands,
1244 }
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
AMREX_GPU_HOST AMREX_FORCE_INLINE real orbital_cos_zenith(real &jday, real &lat, real &lon, real &declin, real dt_avg=-1., real uniform_angle=-1., real constant_zenith_angle_deg=-1.)
Definition: ERF_OrbCosZenith.H:559
AMREX_GPU_HOST AMREX_FORCE_INLINE void orbital_decl(real &calday, real &eccen, real &mvelpp, real &lambm0, real &obliqr, real &delta, real &eccf)
Definition: ERF_OrbCosZenith.H:15
AMREX_GPU_HOST AMREX_FORCE_INLINE void orbital_params(int &iyear_AD, real &eccen, real &obliq, real &mvelp, real &obliqr, real &lambm0, real &mvelpp)
Definition: ERF_OrbCosZenith.H:81
real3d_k aero_tau_sw
Definition: ERF_Radiation.H:448
real3d_k cld_tau_lw_gpt
Definition: ERF_Radiation.H:459
real3d_k cld_tau_sw_gpt
Definition: ERF_Radiation.H:458
real3d_k cld_tau_lw_bnd
Definition: ERF_Radiation.H:455
real3d_k aero_g_sw
Definition: ERF_Radiation.H:450
amrex::Real m_dt
Definition: ERF_Radiation.H:218
real3d_k aero_ssa_sw
Definition: ERF_Radiation.H:449
real3d_k aero_tau_lw
Definition: ERF_Radiation.H:451
real3d_k cld_tau_sw_bnd
Definition: ERF_Radiation.H:454
void rrtmgp_main(const int ncol, const int nlay, real2d_k &p_lay, real2d_k &t_lay, real2d_k &p_lev, real2d_k &t_lev, gas_concs_t &gas_concs, real2d_k &sfc_alb_dir, real2d_k &sfc_alb_dif, real1d_k &mu0, real1d_k &t_sfc, real2d_k &emis_sfc, real1d_k &lw_src, real2d_k &lwp, real2d_k &iwp, real2d_k &rel, real2d_k &rei, real2d_k &cldfrac, real3d_k &, real3d_k &, real3d_k &, real3d_k &, real3d_k &, real3d_k &, real3d_k &, real3d_k &, real2d_k &sw_flux_up, real2d_k &sw_flux_dn, real2d_k &sw_flux_dn_dir, real2d_k &lw_flux_up, real2d_k &lw_flux_dn, real2d_k &sw_clnclrsky_flux_up, real2d_k &sw_clnclrsky_flux_dn, real2d_k &sw_clnclrsky_flux_dn_dir, real2d_k &sw_clrsky_flux_up, real2d_k &sw_clrsky_flux_dn, real2d_k &sw_clrsky_flux_dn_dir, real2d_k &sw_clnsky_flux_up, real2d_k &sw_clnsky_flux_dn, real2d_k &sw_clnsky_flux_dn_dir, real2d_k &lw_clnclrsky_flux_up, real2d_k &lw_clnclrsky_flux_dn, real2d_k &lw_clrsky_flux_up, real2d_k &lw_clrsky_flux_dn, real2d_k &lw_clnsky_flux_up, real2d_k &lw_clnsky_flux_dn, real3d_k &sw_bnd_flux_up, real3d_k &sw_bnd_flux_dn, real3d_k &sw_bnd_flux_dn_dir, real3d_k &lw_bnd_flux_up, real3d_k &lw_bnd_flux_dn, const RealT tsi_scaling, const bool extra_clnclrsky_diag, const bool extra_clnsky_diag)
Definition: ERF_RRTMGP_Interface.cpp:384
void compute_band_by_band_surface_albedos(const int ncol, const int nswbands, real1d_k &sfc_alb_dir_vis, real1d_k &sfc_alb_dir_nir, real1d_k &sfc_alb_dif_vis, real1d_k &sfc_alb_dif_nir, real2d_k &sfc_alb_dir, real2d_k &sfc_alb_dif)
Definition: ERF_RRTMGP_Interface.cpp:282
void compute_broadband_surface_fluxes(const int ncol, const int kbot, const int nswbands, real3d_k &sw_bnd_flux_dir, real3d_k &sw_bnd_flux_dif, real1d_k &sfc_flux_dir_vis, real1d_k &sfc_flux_dir_nir, real1d_k &sfc_flux_dif_vis, real1d_k &sfc_flux_dif_nir)
Definition: ERF_RRTMGP_Interface.cpp:326
void mixing_ratio_to_cloud_mass(View1 const &mixing_ratio, View2 const &cloud_fraction, View3 const &rho, View4 const &dz, View5 const &cloud_mass)
Definition: ERF_RRTMGP_Utils.H:12
std::unique_ptr< cloud_optics_t > cloud_optics_sw_k
Definition: ERF_RRTMGP_Interface.cpp:18

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::Vector< amrex::MultiFab * > &  lsm_input_ptrs,
amrex::MultiFab *  qheating_rates,
amrex::MultiFab *  rad_fluxes,
amrex::MultiFab *  z_phys,
amrex::MultiFab *  lat,
amrex::MultiFab *  lon 
)
150 {
151  // Set data members that may change
152  m_lev = level;
153  m_step = step;
154  m_time = time;
155  m_dt = dt;
156  m_geom = geom;
157  m_cons_in = cons_in;
158  m_lsm_fluxes = lsm_fluxes;
159  m_lsm_zenith = lsm_zenith;
160  m_qheating_rates = qheating_rates;
161  m_rad_fluxes = rad_fluxes;
162  m_z_phys = z_phys;
163  m_lat = lat;
164  m_lon = lon;
165 
166  // Update the day and month
167  time_t timestamp = time_t(time);
168  struct tm *timeinfo = gmtime(&timestamp);
169  if (m_fixed_orbital_year) {
170  m_orbital_mon = timeinfo->tm_mon + 1;
171  m_orbital_day = timeinfo->tm_mday;
172  m_orbital_sec = timeinfo->tm_hour*3600 + timeinfo->tm_min*60 + timeinfo->tm_sec;
173  } else {
174  m_orbital_year = timeinfo->tm_year + 1900;
175  m_orbital_mon = timeinfo->tm_mon + 1;
176  m_orbital_day = timeinfo->tm_mday;
177  m_orbital_sec = timeinfo->tm_hour*3600 + timeinfo->tm_min*60 + timeinfo->tm_sec;
178  }
179 
180  // Only allocate and proceed if we are going to update radiation
181  m_update_rad = false;
182  if (m_rad_freq_in_steps > 0) { m_update_rad = ( (m_step == 0) || (m_step % m_rad_freq_in_steps == 0) ); }
183 
184  if (m_update_rad) {
185  // Call to Init() has set the dimensions: ncol & nlay
186 
187  // Allocate the buffer arrays
188  alloc_buffers();
189 
190  // Fill the KOKKOS Views from AMReX MFs
191  mf_to_kokkos_buffers(lsm_input_ptrs);
192 
193  // Initialize datalog MF on first step
194  if (m_first_step) {
195  m_first_step = false;
196  datalog_mf.define(cons_in->boxArray(), cons_in->DistributionMap(), 25, 0);
197  datalog_mf.setVal(0.0);
198  }
199  }
200 }
int m_step
Definition: ERF_Radiation.H:212
void mf_to_kokkos_buffers(amrex::Vector< amrex::MultiFab * > &lsm_input_ptrs)
Definition: ERF_Radiation.cpp:435
bool m_first_step
Definition: ERF_Radiation.H:232
void alloc_buffers()
Definition: ERF_Radiation.cpp:203
amrex::Real m_time
Definition: ERF_Radiation.H:215

Referenced by Run().

Here is the caller graph for this function:

◆ write_rrtmgp_fluxes()

void Radiation::write_rrtmgp_fluxes ( )
717 {
718  // Expose for device
719  auto sw_flux_up_d = sw_flux_up;
720  auto sw_flux_dn_d = sw_flux_dn;
721  auto sw_flux_dn_dir_d = sw_flux_dn_dir;
722  auto lw_flux_up_d = lw_flux_up;
723  auto lw_flux_dn_d = lw_flux_dn;
724 
725  int n_fluxes = 5;
726  MultiFab mf_flux(m_cons_in->boxArray(), m_cons_in->DistributionMap(), n_fluxes, 0);
727 
728  for (MFIter mfi(mf_flux); mfi.isValid(); ++mfi) {
729  const auto& vbx = mfi.validbox();
730  const int nx = vbx.length(0);
731  const int imin = vbx.smallEnd(0);
732  const int jmin = vbx.smallEnd(1);
733  const int offset = m_col_offsets[mfi.index()];
734  const Array4<Real>& dst_arr = mf_flux.array(mfi);
735  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
736  {
737  // map [i,j,k] 0-based to [icol, ilay] 0-based
738  const int icol = (j-jmin)*nx + (i-imin) + offset;
739  const int ilay = k;
740 
741  // SW and LW fluxes
742  dst_arr(i,j,k,0) = sw_flux_up_d(icol,ilay);
743  dst_arr(i,j,k,1) = sw_flux_dn_d(icol,ilay);
744  dst_arr(i,j,k,2) = sw_flux_dn_dir_d(icol,ilay);
745  dst_arr(i,j,k,3) = lw_flux_up_d(icol,ilay);
746  dst_arr(i,j,k,4) = lw_flux_dn_d(icol,ilay);
747  });
748  }
749 
750 
751  std::string plotfilename = amrex::Concatenate("plt_rad", m_step, 5);
752  Vector<std::string> flux_names = {"sw_flux_up", "sw_flux_dn", "sw_flux_dir",
753  "lw_flux_up", "lw_flux_dn"};
754  WriteSingleLevelPlotfile(plotfilename, mf_flux, flux_names, m_geom, m_time, m_step);
755 }
Here is the call graph for this function:

◆ WriteDataLog()

void Radiation::WriteDataLog ( const amrex::Real time)
overridevirtual

Implements IRadiation.

847 {
848  constexpr int datwidth = 14;
849  constexpr int datprecision = 9;
850  constexpr int timeprecision = 13;
851 
852  Gpu::HostVector<Real> h_avg_radqrsw, h_avg_radqrlw, h_avg_sw_up, h_avg_sw_dn, h_avg_sw_dn_dir, h_avg_lw_up, h_avg_lw_dn, h_avg_zenith;
853  // Clear sky
854  Gpu::HostVector<Real> h_avg_radqrcsw, h_avg_radqrclw, h_avg_sw_clr_up, h_avg_sw_clr_dn, h_avg_sw_clr_dn_dir, h_avg_lw_clr_up, h_avg_lw_clr_dn;
855  // Clean sky
856  Gpu::HostVector<Real> h_avg_sw_cln_up, h_avg_sw_cln_dn, h_avg_sw_cln_dn_dir, h_avg_lw_cln_up, h_avg_lw_cln_dn;
857  // Clean clear sky
858  Gpu::HostVector<Real> h_avg_sw_clnclr_up, h_avg_sw_clnclr_dn, h_avg_sw_clnclr_dn_dir, h_avg_lw_clnclr_up, h_avg_lw_clnclr_dn;
859 
860 
861  auto domain = m_geom.Domain();
862  h_avg_radqrsw = sumToLine(datalog_mf, 0, 1, domain, 2);
863  h_avg_radqrlw = sumToLine(datalog_mf, 1, 1, domain, 2);
864  h_avg_sw_up = sumToLine(datalog_mf, 2, 1, domain, 2);
865  h_avg_sw_dn = sumToLine(datalog_mf, 3, 1, domain, 2);
866  h_avg_sw_dn_dir = sumToLine(datalog_mf, 4, 1, domain, 2);
867  h_avg_lw_up = sumToLine(datalog_mf, 5, 1, domain, 2);
868  h_avg_lw_dn = sumToLine(datalog_mf, 6, 1, domain, 2);
869  h_avg_zenith = sumToLine(datalog_mf, 7, 1, domain, 2);
870 
871  h_avg_radqrcsw = sumToLine(datalog_mf, 8, 1, domain, 2);
872  h_avg_radqrclw = sumToLine(datalog_mf, 9, 1, domain, 2);
873  h_avg_sw_clr_up = sumToLine(datalog_mf, 10, 1, domain, 2);
874  h_avg_sw_clr_dn = sumToLine(datalog_mf, 11, 1, domain, 2);
875  h_avg_sw_clr_dn_dir = sumToLine(datalog_mf, 12, 1, domain, 2);
876  h_avg_lw_clr_up = sumToLine(datalog_mf, 13, 1, domain, 2);
877  h_avg_lw_clr_dn = sumToLine(datalog_mf, 14, 1, domain, 2);
878 
879  if (m_extra_clnsky_diag) {
880  h_avg_sw_cln_up = sumToLine(datalog_mf, 15, 1, domain, 2);
881  h_avg_sw_cln_dn = sumToLine(datalog_mf, 16, 1, domain, 2);
882  h_avg_sw_cln_dn_dir = sumToLine(datalog_mf, 17, 1, domain, 2);
883  h_avg_lw_cln_up = sumToLine(datalog_mf, 18, 1, domain, 2);
884  h_avg_lw_cln_dn = sumToLine(datalog_mf, 19, 1, domain, 2);
885  }
886 
888  h_avg_sw_clnclr_up = sumToLine(datalog_mf, 20, 1, domain, 2);
889  h_avg_sw_clnclr_dn = sumToLine(datalog_mf, 21, 1, domain, 2);
890  h_avg_sw_clnclr_dn_dir = sumToLine(datalog_mf, 22, 1, domain, 2);
891  h_avg_lw_clnclr_up = sumToLine(datalog_mf, 23, 1, domain, 2);
892  h_avg_lw_clnclr_dn = sumToLine(datalog_mf, 24, 1, domain, 2);
893  }
894 
895  Real area_z = static_cast<Real>(domain.length(0)*domain.length(1));
896  int nz = domain.length(2);
897  for (int k = 0; k < nz; k++) {
898  h_avg_radqrsw[k] /= area_z;
899  h_avg_radqrlw[k] /= area_z;
900  h_avg_sw_up[k] /= area_z;
901  h_avg_sw_dn[k] /= area_z;
902  h_avg_sw_dn_dir[k] /= area_z;
903  h_avg_lw_up[k] /= area_z;
904  h_avg_lw_dn[k] /= area_z;
905  h_avg_zenith[k] /= area_z;
906 
907  h_avg_radqrcsw[k] /= area_z;
908  h_avg_radqrclw[k] /= area_z;
909  h_avg_sw_clr_up[k] /= area_z;
910  h_avg_sw_clr_dn[k] /= area_z;
911  h_avg_sw_clr_dn_dir[k] /= area_z;
912  h_avg_lw_clr_up[k] /= area_z;
913  h_avg_lw_clr_dn[k] /= area_z;
914  }
915 
916  if (m_extra_clnsky_diag) {
917  for (int k = 0; k < nz; k++) {
918  h_avg_sw_cln_up[k] /= area_z;
919  h_avg_sw_cln_dn[k] /= area_z;
920  h_avg_sw_cln_dn_dir[k] /= area_z;
921  h_avg_lw_cln_up[k] /= area_z;
922  h_avg_lw_cln_dn[k] /= area_z;
923  }
924  }
925 
927  for (int k = 0; k < nz; k++) {
928  h_avg_sw_clnclr_up[k] /= area_z;
929  h_avg_sw_clnclr_dn[k] /= area_z;
930  h_avg_sw_clnclr_dn_dir[k] /= area_z;
931  h_avg_lw_clnclr_up[k] /= area_z;
932  h_avg_lw_clnclr_dn[k] /= area_z;
933  }
934  }
935 
936  if (ParallelDescriptor::IOProcessor()) {
937  std::ostream& log = *datalog;
938  if (log.good()) {
939 
940  for (int k = 0; k < nz; k++)
941  {
942  Real z = k * m_geom.CellSize(2);
943  log << std::setw(datwidth) << std::setprecision(timeprecision) << time << " "
944  << std::setw(datwidth) << std::setprecision(datprecision) << z << " "
945  << h_avg_radqrsw[k] << " " << h_avg_radqrlw[k] << " " << h_avg_sw_up[k] << " "
946  << h_avg_sw_dn[k] << " " << h_avg_sw_dn_dir[k] << " " << h_avg_lw_up[k] << " "
947  << h_avg_lw_dn[k] << " " << h_avg_zenith[k] << " "
948  << h_avg_radqrcsw[k] << " " << h_avg_radqrclw[k] << " " << h_avg_sw_clr_up[k] << " "
949  << h_avg_sw_clr_dn[k] << " " << h_avg_sw_clr_dn_dir[k] << " " << h_avg_lw_clr_up[k] << " "
950  << h_avg_lw_clr_dn[k] << " ";
951  if (m_extra_clnsky_diag) {
952  log << h_avg_sw_cln_up[k] << " " << h_avg_sw_cln_dn[k] << " " << h_avg_sw_cln_dn_dir[k] << " "
953  << h_avg_lw_cln_up[k] << " " << h_avg_lw_cln_dn[k] << " ";
954  } else {
955  log << 0.0 << " " << 0.0 << " " << 0.0 << " " << 0.0 << " " << 0.0 << " ";
956  }
957 
959  log << h_avg_sw_clnclr_up[k] << " " << h_avg_sw_clnclr_dn[k] << " " << h_avg_sw_clnclr_dn_dir[k] << " "
960  << h_avg_lw_clnclr_up[k] << " " << h_avg_lw_clnclr_dn[k] << std::endl;
961  } else {
962  log << 0.0 << " " << 0.0 << " " << 0.0 << " " << 0.0 << " " << 0.0 << std::endl;
963  }
964  }
965  // Write top face values
966  Real z = nz * m_geom.CellSize(2);
967  log << std::setw(datwidth) << std::setprecision(timeprecision) << time << " "
968  << std::setw(datwidth) << std::setprecision(datprecision) << z << " "
969  << 0.0 << " " << 0.0 << " " << 0.0 << " " << 0.0 << " " << 0.0 << " " << 0.0 << " "
970  << 0.0 << " " << 0.0 << " "
971  << 0.0 << " " << 0.0 << " " << 0.0 << " " << 0.0 << " " << 0.0 << " " << 0.0 << " "
972  << 0.0 << " "
973  << 0.0 << " " << 0.0 << " " << 0.0 << " " << 0.0 << " " << 0.0 << " "
974  << 0.0 << " " << 0.0 << " " << 0.0 << " " << 0.0 << " " << 0.0
975  << std::endl;
976  }
977  }
978 }
std::unique_ptr< std::fstream > datalog
Definition: ERF_RadiationInterface.H:85

Member Data Documentation

◆ aero_g_sw

real3d_k Radiation::aero_g_sw
private

◆ aero_ssa_sw

real3d_k Radiation::aero_ssa_sw
private

◆ aero_tau_lw

real3d_k Radiation::aero_tau_lw
private

◆ aero_tau_sw

real3d_k Radiation::aero_tau_sw
private

◆ cld_tau_lw_bnd

real3d_k Radiation::cld_tau_lw_bnd
private

◆ cld_tau_lw_gpt

real3d_k Radiation::cld_tau_lw_gpt
private

◆ cld_tau_sw_bnd

real3d_k Radiation::cld_tau_sw_bnd
private

◆ cld_tau_sw_gpt

real3d_k Radiation::cld_tau_sw_gpt
private

◆ cldfrac_tot

real2d_k Radiation::cldfrac_tot
private

◆ d_tint

real2d_k Radiation::d_tint
private

◆ datalog_mf

amrex::MultiFab Radiation::datalog_mf
private

◆ eff_radius_qc

real2d_k Radiation::eff_radius_qc
private

◆ eff_radius_qi

real2d_k Radiation::eff_radius_qi
private

◆ emis_sfc

real2d_k Radiation::emis_sfc
private

◆ gas_names_offset

std::vector<std::string> Radiation::gas_names_offset
private

◆ iwp

real2d_k Radiation::iwp
private

◆ lat

real1d_k Radiation::lat
private

◆ lon

real1d_k Radiation::lon
private

◆ lw_bnd_flux_dn

real3d_k Radiation::lw_bnd_flux_dn
private

◆ lw_bnd_flux_up

real3d_k Radiation::lw_bnd_flux_up
private

◆ lw_clnclrsky_flux_dn

real2d_k Radiation::lw_clnclrsky_flux_dn
private

◆ lw_clnclrsky_flux_up

real2d_k Radiation::lw_clnclrsky_flux_up
private

◆ lw_clnsky_flux_dn

real2d_k Radiation::lw_clnsky_flux_dn
private

◆ lw_clnsky_flux_up

real2d_k Radiation::lw_clnsky_flux_up
private

◆ lw_clrsky_flux_dn

real2d_k Radiation::lw_clrsky_flux_dn
private

◆ lw_clrsky_flux_up

real2d_k Radiation::lw_clrsky_flux_up
private

◆ lw_clrsky_heating

real2d_k Radiation::lw_clrsky_heating
private

◆ lw_flux_dn

real2d_k Radiation::lw_flux_dn
private

◆ lw_flux_up

real2d_k Radiation::lw_flux_up
private

◆ lw_heating

real2d_k Radiation::lw_heating
private

◆ lw_src

real1d_k Radiation::lw_src
private

◆ lwp

real2d_k 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

Referenced by Init().

◆ 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_first_step

bool Radiation::m_first_step = true
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

GasConcsK<amrex::Real, layout_t, KokkosDefaultDevice> Radiation::m_gas_concs
private

◆ m_gas_mol_weights

real1d_k 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

bool Radiation::m_lsm = false
private

◆ m_lsm_fluxes

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

◆ m_lsm_input_names

amrex::Vector<std::string> Radiation::m_lsm_input_names
private
Initial value:
= {"t_sfc" , "sfc_emis" ,
"sfc_alb_dir_vis", "sfc_alb_dir_nir",
"sfc_alb_dif_vis", "sfc_alb_dif_nir"}

Referenced by get_lsm_input_varnames().

◆ m_lsm_output_names

amrex::Vector<std::string> Radiation::m_lsm_output_names = {"cos_zenith_angle", "sw_flux_dn", "lw_flux_dn"}
private

Referenced by get_lsm_output_varnames().

◆ 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

Referenced by Init().

◆ m_ngas

int Radiation::m_ngas = 8
private

◆ m_nlay

int Radiation::m_nlay
private

Referenced by Init().

◆ m_nlwbands

int Radiation::m_nlwbands
private

◆ m_nlwgpts

int Radiation::m_nlwgpts
private

◆ m_nswbands

int Radiation::m_nswbands
private

◆ m_nswgpts

int Radiation::m_nswgpts
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

Referenced by rad_run_impl().

◆ m_orbital_eccen

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

◆ m_orbital_mon

int Radiation::m_orbital_mon = -9999
private

Referenced by rad_run_impl().

◆ 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

Referenced by rad_run_impl().

◆ m_orbital_year

int Radiation::m_orbital_year = -9999
private

Referenced by rad_run_impl().

◆ m_qheating_rates

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

◆ m_rad_fluxes

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

◆ m_rad_freq_in_steps

int Radiation::m_rad_freq_in_steps = 1
private

◆ m_rad_t_sfc

amrex::Real Radiation::m_rad_t_sfc = -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_k Radiation::mu0
private

◆ o3_lay

real1d_k Radiation::o3_lay
private

◆ p_lay

real2d_k Radiation::p_lay
private

◆ p_lev

real2d_k Radiation::p_lev
private

◆ qc_lay

real2d_k Radiation::qc_lay
private

◆ qi_lay

real2d_k Radiation::qi_lay
private

◆ qv_lay

real2d_k Radiation::qv_lay
private

◆ r_lay

real2d_k 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-g256-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_k Radiation::sfc_alb_dif
private

◆ sfc_alb_dif_nir

real1d_k Radiation::sfc_alb_dif_nir
private

◆ sfc_alb_dif_vis

real1d_k Radiation::sfc_alb_dif_vis
private

◆ sfc_alb_dir

real2d_k Radiation::sfc_alb_dir
private

◆ sfc_alb_dir_nir

real1d_k Radiation::sfc_alb_dir_nir
private

◆ sfc_alb_dir_vis

real1d_k Radiation::sfc_alb_dir_vis
private

◆ sfc_emis

real1d_k Radiation::sfc_emis
private

◆ sfc_flux_dif_nir

real1d_k Radiation::sfc_flux_dif_nir
private

◆ sfc_flux_dif_vis

real1d_k Radiation::sfc_flux_dif_vis
private

◆ sfc_flux_dir_nir

real1d_k Radiation::sfc_flux_dir_nir
private

◆ sfc_flux_dir_vis

real1d_k Radiation::sfc_flux_dir_vis
private

◆ sw_bnd_flux_dif

real3d_k Radiation::sw_bnd_flux_dif
private

◆ sw_bnd_flux_dir

real3d_k Radiation::sw_bnd_flux_dir
private

◆ sw_bnd_flux_dn

real3d_k Radiation::sw_bnd_flux_dn
private

◆ sw_bnd_flux_up

real3d_k Radiation::sw_bnd_flux_up
private

◆ sw_clnclrsky_flux_dn

real2d_k Radiation::sw_clnclrsky_flux_dn
private

◆ sw_clnclrsky_flux_dn_dir

real2d_k Radiation::sw_clnclrsky_flux_dn_dir
private

◆ sw_clnclrsky_flux_up

real2d_k Radiation::sw_clnclrsky_flux_up
private

◆ sw_clnsky_flux_dn

real2d_k Radiation::sw_clnsky_flux_dn
private

◆ sw_clnsky_flux_dn_dir

real2d_k Radiation::sw_clnsky_flux_dn_dir
private

◆ sw_clnsky_flux_up

real2d_k Radiation::sw_clnsky_flux_up
private

◆ sw_clrsky_flux_dn

real2d_k Radiation::sw_clrsky_flux_dn
private

◆ sw_clrsky_flux_dn_dir

real2d_k Radiation::sw_clrsky_flux_dn_dir
private

◆ sw_clrsky_flux_up

real2d_k Radiation::sw_clrsky_flux_up
private

◆ sw_clrsky_heating

real2d_k Radiation::sw_clrsky_heating
private

◆ sw_flux_dn

real2d_k Radiation::sw_flux_dn
private

◆ sw_flux_dn_dir

real2d_k Radiation::sw_flux_dn_dir
private

◆ sw_flux_up

real2d_k Radiation::sw_flux_up
private

◆ sw_heating

real2d_k Radiation::sw_heating
private

◆ t_lay

real2d_k Radiation::t_lay
private

◆ t_lev

real2d_k Radiation::t_lev
private

◆ t_sfc

real1d_k Radiation::t_sfc
private

◆ z_del

real2d_k Radiation::z_del
private

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