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::iMultiFab *lmask, amrex::MultiFab *t_surf, 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::iMultiFab *lmask, amrex::MultiFab *t_surf, 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::iMultiFab *lmask, amrex::MultiFab *t_surf, 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
 
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 = amrex::Real(39.809860)
 
amrex::Real m_lon_cons = -amrex::Real(98.555183)
 
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 = amrex::Real(388.717e-6)
 
amrex::Vector< amrex::Realm_o3vmr
 
amrex::Real m_n2ovmr = amrex::Real(323.141e-9)
 
amrex::Real m_covmr = amrex::Real(1.0e-7)
 
amrex::Real m_ch4vmr = amrex::Real(1807.851e-9)
 
amrex::Real m_o2vmr = amrex::Real(0.209448)
 
amrex::Real m_n2vmr = amrex::Real(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 = false
 
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 = -amrex::Real(9999.)
 
amrex::Real m_orbital_obliq = -amrex::Real(9999.)
 
amrex::Real m_orbital_mvelp = -amrex::Real(9999.)
 
amrex::Real m_fixed_total_solar_irradiance = -amrex::Real(9999.)
 
amrex::Real m_fixed_solar_zenith_angle = -amrex::Real(9999.)
 
int m_nswbands
 
int m_nlwbands
 
int m_nswgpts
 
int m_nlwgpts
 
int m_rad_freq_in_steps = 1
 
int m_ncol_chunk = 1024
 
int m_rad_nvar = 12
 
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
 
real3d_k aero_tau_sw
 
real3d_k aero_ssa_sw
 
real3d_k aero_g_sw
 
real3d_k aero_tau_lw
 

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 (LSM can overwrite)
38  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  // Get nvar if specified
44  pp.query("rad_nvar", m_rad_nvar);
45  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(m_rad_nvar >= 0,
46  "erf.rad_nvar must be greater than 0. "
47  "It controls the amount of memory allocated for temporaries with RRTMGP; "
48  "a value of 0 would allocate no memory.");
49 
50  // Number of columns per RRTMGP chunk (controls peak GPU memory)
51  pp.query("rad_ncol_chunk", m_ncol_chunk);
52  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(m_ncol_chunk > 0,
53  "erf.rad_ncol_chunk must be a positive integer (default 5000). "
54  "It controls the number of columns processed per RRTMGP kernel launch; "
55  "a value of 0 or negative would produce an infinite loop.");
56 
57  // Flag to write fluxes to plt file
58  pp.query("rad_write_fluxes", m_rad_write_fluxes);
59 
60  // Do MCICA subcolumn sampling
61  pp.query("rad_do_subcol_sampling", m_do_subcol_sampling);
62 
63  // Determine orbital year. If orbital_year is negative, use current year
64  // from timestamp for orbital year; if positive, use provided orbital year
65  // for duration of simulation.
66  m_fixed_orbital_year = pp.query("rad_orbital_year", m_orbital_year);
67 
68  // Get orbital parameters from inputs file
69  pp.query("rad_orbital_eccentricity", m_orbital_eccen);
70  pp.query("rad_orbital_obliquity" , m_orbital_obliq);
71  pp.query("rad_orbital_mvelp" , m_orbital_mvelp);
72 
73  // Get a constant lat/lon for idealized simulations
74  pp.query("rad_cons_lat", m_lat_cons);
75  pp.query("rad_cons_lon", m_lon_cons);
76 
77  // Value for prescribing an invariant solar constant (i.e. total solar irradiance at
78  // TOA). Used for idealized experiments such as RCE. Disabled when value is less than zero
79  pp.query("fixed_total_solar_irradiance", m_fixed_total_solar_irradiance);
80 
81  // Determine whether or not we are using a fixed solar zenith angle (positive value)
82  pp.query("fixed_solar_zenith_angle", m_fixed_solar_zenith_angle);
83 
84  // Get prescribed surface values of greenhouse gases
85  pp.query("co2vmr", m_co2vmr);
86  pp.queryarr("o3vmr" , m_o3vmr );
87  pp.query("n2ovmr", m_n2ovmr);
88  pp.query("covmr" , m_covmr );
89  pp.query("ch4vmr", m_ch4vmr);
90  pp.query("o2vmr" , m_o2vmr );
91  pp.query("n2vmr" , m_n2vmr );
92 
93  // Aerosol forcing hook (not implemented). The aerosol arrays that used to be
94  // passed through rrtmgp_main were never populated with real data, so enabling
95  // this flag only ever multiplied radiation by zero aerosol optics. The hook is
96  // kept so a future SPA/prescribed-aerosol scheme can wire in without touching
97  // the ParmParse surface.
98  pp.query("rad_do_aerosol", m_do_aerosol_rad);
99  if (m_do_aerosol_rad) {
100  amrex::Abort("erf.rad_do_aerosol = true is not supported: aerosol forcing is "
101  "currently not implemented in the ERF RRTMGP interface. The hook "
102  "is retained for a future aerosol coupling; set rad_do_aerosol = "
103  "false (or remove it) to continue.");
104  }
105 
106  // Whether we do extra clean/clear sky calculations
107  pp.query("rad_extra_clnclrsky_diag", m_extra_clnclrsky_diag);
108  pp.query("rad_extra_clnsky_diag" , m_extra_clnsky_diag);
109 
110  // Parse path and file names
111  pp.query("rrtmgp_file_path" , rrtmgp_file_path);
112  pp.query("rrtmgp_coeffs_sw" , rrtmgp_coeffs_sw );
113  pp.query("rrtmgp_coeffs_lw" , rrtmgp_coeffs_lw );
114  pp.query("rrtmgp_cloud_optics_sw", rrtmgp_cloud_optics_sw);
115  pp.query("rrtmgp_cloud_optics_lw", rrtmgp_cloud_optics_lw);
116 
117  // Append file names to path
122 
123  // Get dimensions from lookup data
124  if (ParallelDescriptor::IOProcessor()) {
125  auto ncf_sw = ncutils::NCFile::open(rrtmgp_coeffs_file_sw, NC_CLOBBER | NC_NETCDF4);
126  m_nswbands = ncf_sw.dim("bnd").len();
127  m_nswgpts = ncf_sw.dim("gpt").len();
128  ncf_sw.close();
129 
130  auto ncf_lw = ncutils::NCFile::open(rrtmgp_coeffs_file_lw, NC_CLOBBER | NC_NETCDF4);
131  m_nlwbands = ncf_lw.dim("bnd").len();
132  m_nlwgpts = ncf_lw.dim("gpt").len();
133  ncf_lw.close();
134  }
135  int ioproc = ParallelDescriptor::IOProcessorNumber(); // I/O rank
136  ParallelDescriptor::Bcast(&m_nswbands, 1, ioproc);
137  ParallelDescriptor::Bcast(&m_nlwbands, 1, ioproc);
138  ParallelDescriptor::Bcast(&m_nswgpts, 1, ioproc);
139  ParallelDescriptor::Bcast(&m_nlwgpts, 1, ioproc);
140 
141  // Output for user
142  if (lev == 0) {
143  Print() << "Radiation interface constructed:\n";
144  Print() << "========================================================\n";
145  Print() << "Coeff SW file: " << rrtmgp_coeffs_file_sw << "\n";
146  Print() << "Coeff LW file: " << rrtmgp_coeffs_file_lw << "\n";
147  Print() << "Cloud SW file: " << rrtmgp_cloud_optics_file_sw << "\n";
148  Print() << "Cloud LW file: " << rrtmgp_cloud_optics_file_lw << "\n";
149  Print() << "Number of short/longwave bands: "
150  << m_nswbands << " " << m_nlwbands << "\n";
151  Print() << "Number of short/longwave gauss points: "
152  << m_nswgpts << " " << m_nlwgpts << "\n";
153  Print() << "========================================================\n";
154  }
155 }
ParmParse pp("prob")
std::string rrtmgp_coeffs_file_sw
Definition: ERF_Radiation.H:292
std::string rrtmgp_coeffs_sw
Definition: ERF_Radiation.H:288
int m_rad_freq_in_steps
Definition: ERF_Radiation.H:373
amrex::Real m_lon_cons
Definition: ERF_Radiation.H:281
int m_nswbands
Definition: ERF_Radiation.H:367
bool m_do_aerosol_rad
Definition: ERF_Radiation.H:335
std::string rrtmgp_coeffs_lw
Definition: ERF_Radiation.H:289
bool m_moist
Definition: ERF_Radiation.H:243
std::string rrtmgp_cloud_optics_file_lw
Definition: ERF_Radiation.H:295
bool m_lsm
Definition: ERF_Radiation.H:247
bool m_do_subcol_sampling
Definition: ERF_Radiation.H:382
bool m_rad_write_fluxes
Definition: ERF_Radiation.H:238
std::string rrtmgp_cloud_optics_file_sw
Definition: ERF_Radiation.H:294
amrex::Real m_orbital_mvelp
Definition: ERF_Radiation.H:355
amrex::Vector< amrex::Real > m_o3vmr
Definition: ERF_Radiation.H:306
std::string rrtmgp_cloud_optics_sw
Definition: ERF_Radiation.H:290
int m_nlwgpts
Definition: ERF_Radiation.H:370
std::string rrtmgp_cloud_optics_lw
Definition: ERF_Radiation.H:291
amrex::Real m_co2vmr
Definition: ERF_Radiation.H:305
bool m_extra_clnsky_diag
Definition: ERF_Radiation.H:338
amrex::Real m_lat_cons
Definition: ERF_Radiation.H:280
amrex::Real m_fixed_total_solar_irradiance
Definition: ERF_Radiation.H:360
amrex::Real m_o2vmr
Definition: ERF_Radiation.H:310
amrex::Real m_n2ovmr
Definition: ERF_Radiation.H:307
amrex::Real m_rad_t_sfc
Definition: ERF_Radiation.H:261
std::string rrtmgp_file_path
Definition: ERF_Radiation.H:287
int m_nlwbands
Definition: ERF_Radiation.H:368
amrex::Real m_orbital_eccen
Definition: ERF_Radiation.H:353
amrex::Real m_covmr
Definition: ERF_Radiation.H:308
int m_ncol_chunk
Definition: ERF_Radiation.H:376
std::string rrtmgp_coeffs_file_lw
Definition: ERF_Radiation.H:293
amrex::Real m_ch4vmr
Definition: ERF_Radiation.H:309
bool m_ice
Definition: ERF_Radiation.H:244
amrex::Real m_fixed_solar_zenith_angle
Definition: ERF_Radiation.H:364
amrex::Real m_orbital_obliq
Definition: ERF_Radiation.H:354
amrex::Real m_n2vmr
Definition: ERF_Radiation.H:311
bool m_fixed_orbital_year
Definition: ERF_Radiation.H:352
int m_nswgpts
Definition: ERF_Radiation.H:369
bool m_extra_clnclrsky_diag
Definition: ERF_Radiation.H:339
int m_rad_nvar
Definition: ERF_Radiation.H:379
int m_orbital_year
Definition: ERF_Radiation.H:344
static NCFile open(const std::string &name, const int cmode=NC_NOWRITE)
Definition: ERF_NCInterface.cpp:707
LandSurfaceType lsm_type
Definition: ERF_DataStruct.H:1302
MoistureType moisture_type
Definition: ERF_DataStruct.H:1299
Here is the call graph for this function:

◆ ~Radiation()

Radiation::~Radiation ( )
inline
56  {
57  // Release k-distribution data and memory pool
58  if (rrtmgp::initialized) {
60  }
61  // Note that Kokkos is now finalized in main.cpp
62  }
void rrtmgp_finalize()
Definition: ERF_RRTMGP_Interface.cpp:266
bool initialized
Definition: ERF_RRTMGP_Interface.cpp:24
Here is the call graph for this function:

Member Function Documentation

◆ alloc_buffers()

void Radiation::alloc_buffers ( )
226 {
227  // 1d size (m_ngas)
228  const Real* mol_weight_gas_p = m_mol_weight_gas.data();
229  const std::string* gas_names_p = m_gas_names.data();
230  m_gas_mol_weights = real1d_k("m_gas_mol_weights", m_ngas);
231  realHost1d_k m_gas_mol_weights_h("m_gas_mol_weights_h", m_ngas);
232  gas_names_offset.clear(); gas_names_offset.resize(m_ngas);
233  std::string* gas_names_offset_p = gas_names_offset.data();
234  Kokkos::parallel_for(Kokkos::RangePolicy<Kokkos::Serial>(0, m_ngas),
235  [&] (int igas)
236  {
237  m_gas_mol_weights_h(igas) = mol_weight_gas_p[igas];
238  gas_names_offset_p[igas] = gas_names_p[igas];
239  });
240  Kokkos::deep_copy(m_gas_mol_weights, m_gas_mol_weights_h);
241 
242  // 1d size (1 or nlay)
243  m_o3_size = m_o3vmr.size();
244  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(((m_o3_size==1) || (m_o3_size==m_nlay)),
245  "O3 VMR array must be length 1 or nlay");
246  Real* o3vmr_p = m_o3vmr.data();
247  o3_lay = real1d_k("o3_lay", m_o3_size);
248  realHost1d_k o3_lay_h("o3_lay_h", m_o3_size);
249  Kokkos::parallel_for(Kokkos::RangePolicy<Kokkos::Serial>(0, m_o3_size),
250  [&] (int io3)
251  {
252  o3_lay_h(io3) = o3vmr_p[io3];
253  });
254  Kokkos::deep_copy(o3_lay, o3_lay_h);
255 
256  // 1d size (ncol)
257  mu0 = real1d_k("mu0" , m_ncol);
258  sfc_alb_dir_vis = real1d_k("sfc_alb_dir_vis" , m_ncol);
259  sfc_alb_dir_nir = real1d_k("sfc_alb_dir_nir" , m_ncol);
260  sfc_alb_dif_vis = real1d_k("sfc_alb_dif_vis" , m_ncol);
261  sfc_alb_dif_nir = real1d_k("sfc_alb_dif_nir" , m_ncol);
262  sfc_flux_dir_vis = real1d_k("sfc_flux_dir_vis", m_ncol);
263  sfc_flux_dir_nir = real1d_k("sfc_flux_dir_nir", m_ncol);
264  sfc_flux_dif_vis = real1d_k("sfc_flux_dif_vis", m_ncol);
265  sfc_flux_dif_nir = real1d_k("sfc_flux_dif_nir", m_ncol);
266  lat = real1d_k("lat" , m_ncol);
267  lon = real1d_k("lon" , m_ncol);
268  sfc_emis = real1d_k("sfc_emis" , m_ncol);
269  t_sfc = real1d_k("t_sfc" , m_ncol);
270  lw_src = real1d_k("lw_src" , m_ncol);
271 
272  // 2d size (ncol, nlay)
273  r_lay = real2d_k("r_lay" , m_ncol, m_nlay);
274  p_lay = real2d_k("p_lay" , m_ncol, m_nlay);
275  t_lay = real2d_k("t_lay" , m_ncol, m_nlay);
276  z_del = real2d_k("z_del" , m_ncol, m_nlay);
277  qv_lay = real2d_k("qv" , m_ncol, m_nlay);
278  qc_lay = real2d_k("qc" , m_ncol, m_nlay);
279  qi_lay = real2d_k("qi" , m_ncol, m_nlay);
280  cldfrac_tot = real2d_k("cldfrac_tot" , m_ncol, m_nlay);
281  eff_radius_qc = real2d_k("eff_radius_qc", m_ncol, m_nlay);
282  eff_radius_qi = real2d_k("eff_radius_qi", m_ncol, m_nlay);
283  lwp = real2d_k("lwp" , m_ncol, m_nlay);
284  iwp = real2d_k("iwp" , m_ncol, m_nlay);
285  sw_heating = real2d_k("sw_heating" , m_ncol, m_nlay);
286  lw_heating = real2d_k("lw_heating" , m_ncol, m_nlay);
287  sw_clrsky_heating = real2d_k("sw_clrsky_heating", m_ncol, m_nlay);
288  lw_clrsky_heating = real2d_k("lw_clrsky_heating", m_ncol, m_nlay);
289 
290  // 2d size (ncol, nlay+1)
291  d_tint = real2d_k("d_tint" , m_ncol, m_nlay+1);
292  p_lev = real2d_k("p_lev" , m_ncol, m_nlay+1);
293  t_lev = real2d_k("t_lev" , m_ncol, m_nlay+1);
294 
295  sw_flux_up = real2d_k("sw_flux_up" , m_ncol, m_nlay+1);
296  sw_flux_dn = real2d_k("sw_flux_dn" , m_ncol, m_nlay+1);
297  sw_flux_dn_dir = real2d_k("sw_flux_dn_dir" , m_ncol, m_nlay+1);
298 
299  lw_flux_up = real2d_k("lw_flux_up" , m_ncol, m_nlay+1);
300  lw_flux_dn = real2d_k("lw_flux_dn" , m_ncol, m_nlay+1);
301 
302  // Clear-sky flux arrays are always needed
303  sw_clrsky_flux_up = real2d_k("sw_clrsky_flux_up" , m_ncol, m_nlay+1);
304  sw_clrsky_flux_dn = real2d_k("sw_clrsky_flux_dn" , m_ncol, m_nlay+1);
305  sw_clrsky_flux_dn_dir = real2d_k("sw_clrsky_flux_dn_dir" , m_ncol, m_nlay+1);
306  lw_clrsky_flux_up = real2d_k("lw_clrsky_flux_up" , m_ncol, m_nlay+1);
307  lw_clrsky_flux_dn = real2d_k("lw_clrsky_flux_dn" , m_ncol, m_nlay+1);
308 
309  // Clean-clear-sky diagnostic fluxes (only when enabled)
311  sw_clnclrsky_flux_up = real2d_k("sw_clnclrsky_flux_up" , m_ncol, m_nlay+1);
312  sw_clnclrsky_flux_dn = real2d_k("sw_clnclrsky_flux_dn" , m_ncol, m_nlay+1);
313  sw_clnclrsky_flux_dn_dir = real2d_k("sw_clnclrsky_flux_dn_dir", m_ncol, m_nlay+1);
314  lw_clnclrsky_flux_up = real2d_k("lw_clnclrsky_flux_up" , m_ncol, m_nlay+1);
315  lw_clnclrsky_flux_dn = real2d_k("lw_clnclrsky_flux_dn" , m_ncol, m_nlay+1);
316  } else {
317  sw_clnclrsky_flux_up = real2d_k("sw_clnclrsky_flux_up" , 1, 1);
318  sw_clnclrsky_flux_dn = real2d_k("sw_clnclrsky_flux_dn" , 1, 1);
319  sw_clnclrsky_flux_dn_dir = real2d_k("sw_clnclrsky_flux_dn_dir", 1, 1);
320  lw_clnclrsky_flux_up = real2d_k("lw_clnclrsky_flux_up" , 1, 1);
321  lw_clnclrsky_flux_dn = real2d_k("lw_clnclrsky_flux_dn" , 1, 1);
322  }
323 
324  // Clean-sky diagnostic fluxes (only when enabled)
325  if (m_extra_clnsky_diag) {
326  sw_clnsky_flux_up = real2d_k("sw_clnsky_flux_up" , m_ncol, m_nlay+1);
327  sw_clnsky_flux_dn = real2d_k("sw_clnsky_flux_dn" , m_ncol, m_nlay+1);
328  sw_clnsky_flux_dn_dir = real2d_k("sw_clnsky_flux_dn_dir" , m_ncol, m_nlay+1);
329  lw_clnsky_flux_up = real2d_k("lw_clnsky_flux_up" , m_ncol, m_nlay+1);
330  lw_clnsky_flux_dn = real2d_k("lw_clnsky_flux_dn" , m_ncol, m_nlay+1);
331  } else {
332  sw_clnsky_flux_up = real2d_k("sw_clnsky_flux_up" , 1, 1);
333  sw_clnsky_flux_dn = real2d_k("sw_clnsky_flux_dn" , 1, 1);
334  sw_clnsky_flux_dn_dir = real2d_k("sw_clnsky_flux_dn_dir" , 1, 1);
335  lw_clnsky_flux_up = real2d_k("lw_clnsky_flux_up" , 1, 1);
336  lw_clnsky_flux_dn = real2d_k("lw_clnsky_flux_dn" , 1, 1);
337  }
338 
339  // 3d size (ncol, nlay+1, nswbands)
340  sw_bnd_flux_up = real3d_k("sw_bnd_flux_up" , m_ncol, m_nlay+1, m_nswbands);
341  sw_bnd_flux_dn = real3d_k("sw_bnd_flux_dn" , m_ncol, m_nlay+1, m_nswbands);
342  sw_bnd_flux_dir = real3d_k("sw_bnd_flux_dir", m_ncol, m_nlay+1, m_nswbands);
343  sw_bnd_flux_dif = real3d_k("sw_bnd_flux_dif", m_ncol, m_nlay+1, m_nswbands);
344 
345  // 3d size (ncol, nlay+1, nlwbands)
346  lw_bnd_flux_up = real3d_k("lw_bnd_flux_up" , m_ncol, m_nlay+1, m_nlwbands);
347  lw_bnd_flux_dn = real3d_k("lw_bnd_flux_dn" , m_ncol, m_nlay+1, m_nlwbands);
348 
349  // 2d size (ncol, nswbands)
350  sfc_alb_dir = real2d_k("sfc_alb_dir", m_ncol, m_nswbands);
351  sfc_alb_dif = real2d_k("sfc_alb_dif", m_ncol, m_nswbands);
352 
353  // Aerosol optical properties — allocated only when aerosol coupling is on.
354  // The flag gates allocation so today (coupling not implemented, abort fires
355  // in the constructor) these stay as empty Views and cost nothing. When a
356  // future aerosol scheme populates them, hook up the plumbing into
357  // rrtmgp_main as well.
358  if (m_do_aerosol_rad) {
359  aero_tau_sw = real3d_k("aero_tau_sw", m_ncol, m_nlay, m_nswbands);
360  aero_ssa_sw = real3d_k("aero_ssa_sw", m_ncol, m_nlay, m_nswbands);
361  aero_g_sw = real3d_k("aero_g_sw", m_ncol, m_nlay, m_nswbands);
362  aero_tau_lw = real3d_k("aero_tau_lw", m_ncol, m_nlay, m_nlwbands);
363  }
364 }
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:448
real2d_k lw_flux_up
Definition: ERF_Radiation.H:428
real3d_k aero_tau_sw
Definition: ERF_Radiation.H:467
int m_o3_size
Definition: ERF_Radiation.H:315
real3d_k sw_bnd_flux_dir
Definition: ERF_Radiation.H:449
real2d_k sw_clnsky_flux_dn
Definition: ERF_Radiation.H:437
real2d_k d_tint
Definition: ERF_Radiation.H:422
real2d_k lw_clnclrsky_flux_dn
Definition: ERF_Radiation.H:440
real2d_k lwp
Definition: ERF_Radiation.H:414
real1d_k lw_src
Definition: ERF_Radiation.H:401
real2d_k eff_radius_qi
Definition: ERF_Radiation.H:413
real1d_k m_gas_mol_weights
Definition: ERF_Radiation.H:316
real2d_k sw_heating
Definition: ERF_Radiation.H:416
real3d_k lw_bnd_flux_dn
Definition: ERF_Radiation.H:454
real3d_k sw_bnd_flux_up
Definition: ERF_Radiation.H:447
real2d_k qv_lay
Definition: ERF_Radiation.H:408
real2d_k sw_clnclrsky_flux_dn_dir
Definition: ERF_Radiation.H:432
real1d_k sfc_flux_dif_vis
Definition: ERF_Radiation.H:395
real2d_k lw_clnclrsky_flux_up
Definition: ERF_Radiation.H:439
real1d_k lat
Definition: ERF_Radiation.H:397
real2d_k qi_lay
Definition: ERF_Radiation.H:410
real2d_k sw_clrsky_flux_up
Definition: ERF_Radiation.H:433
real2d_k t_lev
Definition: ERF_Radiation.H:424
real1d_k o3_lay
Definition: ERF_Radiation.H:385
real1d_k sfc_alb_dif_vis
Definition: ERF_Radiation.H:391
real1d_k mu0
Definition: ERF_Radiation.H:388
real2d_k cldfrac_tot
Definition: ERF_Radiation.H:411
real1d_k sfc_flux_dir_nir
Definition: ERF_Radiation.H:394
real2d_k sw_clnclrsky_flux_up
Definition: ERF_Radiation.H:430
real3d_k aero_g_sw
Definition: ERF_Radiation.H:469
real2d_k sw_flux_up
Definition: ERF_Radiation.H:425
real3d_k sw_bnd_flux_dif
Definition: ERF_Radiation.H:450
real2d_k sfc_alb_dif
Definition: ERF_Radiation.H:458
real1d_k sfc_alb_dif_nir
Definition: ERF_Radiation.H:392
real2d_k lw_clnsky_flux_dn
Definition: ERF_Radiation.H:444
real2d_k p_lay
Definition: ERF_Radiation.H:405
real2d_k r_lay
Definition: ERF_Radiation.H:404
int m_ncol
Definition: ERF_Radiation.H:325
real2d_k sw_flux_dn_dir
Definition: ERF_Radiation.H:427
const std::vector< amrex::Real > m_mol_weight_gas
Definition: ERF_Radiation.H:301
std::vector< std::string > gas_names_offset
Definition: ERF_Radiation.H:317
real2d_k sw_clrsky_flux_dn
Definition: ERF_Radiation.H:434
real3d_k aero_ssa_sw
Definition: ERF_Radiation.H:468
real2d_k sw_clnsky_flux_dn_dir
Definition: ERF_Radiation.H:438
real2d_k qc_lay
Definition: ERF_Radiation.H:409
real2d_k lw_clrsky_flux_up
Definition: ERF_Radiation.H:441
real2d_k sw_flux_dn
Definition: ERF_Radiation.H:426
real2d_k lw_clrsky_flux_dn
Definition: ERF_Radiation.H:442
real2d_k z_del
Definition: ERF_Radiation.H:407
real3d_k aero_tau_lw
Definition: ERF_Radiation.H:470
real2d_k lw_flux_dn
Definition: ERF_Radiation.H:429
real2d_k sw_clrsky_flux_dn_dir
Definition: ERF_Radiation.H:435
real1d_k sfc_flux_dif_nir
Definition: ERF_Radiation.H:396
int m_ngas
Definition: ERF_Radiation.H:298
real2d_k lw_heating
Definition: ERF_Radiation.H:417
real1d_k sfc_alb_dir_nir
Definition: ERF_Radiation.H:390
real2d_k lw_clnsky_flux_up
Definition: ERF_Radiation.H:443
real2d_k sfc_alb_dir
Definition: ERF_Radiation.H:457
real1d_k lon
Definition: ERF_Radiation.H:398
real1d_k sfc_emis
Definition: ERF_Radiation.H:399
real2d_k sw_clnsky_flux_up
Definition: ERF_Radiation.H:436
int m_nlay
Definition: ERF_Radiation.H:326
real3d_k lw_bnd_flux_up
Definition: ERF_Radiation.H:453
real1d_k sfc_alb_dir_vis
Definition: ERF_Radiation.H:389
const std::vector< std::string > m_gas_names
Definition: ERF_Radiation.H:299
real2d_k sw_clrsky_heating
Definition: ERF_Radiation.H:418
real2d_k t_lay
Definition: ERF_Radiation.H:406
real1d_k sfc_flux_dir_vis
Definition: ERF_Radiation.H:393
real2d_k iwp
Definition: ERF_Radiation.H:415
real2d_k p_lev
Definition: ERF_Radiation.H:423
real2d_k lw_clrsky_heating
Definition: ERF_Radiation.H:419
real2d_k eff_radius_qc
Definition: ERF_Radiation.H:412
real1d_k t_sfc
Definition: ERF_Radiation.H:400
real2d_k sw_clnclrsky_flux_dn
Definition: ERF_Radiation.H:431

◆ dealloc_buffers()

void Radiation::dealloc_buffers ( )
368 {
369  // 1d size (m_ngas)
371 
372  // 1d size (1 or nlay)
373  o3_lay = real1d_k();
374 
375  // 1d size (ncol)
376  mu0 = real1d_k();
385  lat = real1d_k();
386  lon = real1d_k();
387  sfc_emis = real1d_k();
388  t_sfc = real1d_k();
389  lw_src = real1d_k();
390 
391  // 2d size (ncol, nlay)
392  r_lay = real2d_k();
393  p_lay = real2d_k();
394  t_lay = real2d_k();
395  z_del = real2d_k();
396  qv_lay = real2d_k();
397  qc_lay = real2d_k();
398  qi_lay = real2d_k();
399  cldfrac_tot = real2d_k();
402  lwp = real2d_k();
403  iwp = real2d_k();
404  sw_heating = real2d_k();
405  lw_heating = real2d_k();
408  sw_heating = real2d_k();
409  lw_heating = real2d_k();
412 
413  // 2d size (ncol, nlay+1)
414  d_tint = real2d_k();
415  p_lev = real2d_k();
416  t_lev = real2d_k();
417  sw_flux_up = real2d_k();
418  sw_flux_dn = real2d_k();
420  lw_flux_up = real2d_k();
421  lw_flux_dn = real2d_k();
437 
438  // 3d size (ncol, nlay+1, nswbands)
443 
444  // 3d size (ncol, nlay+1, nlwbands)
447 
448  // 2d size (ncol, nswbands)
449  sfc_alb_dir = real2d_k();
450  sfc_alb_dif = real2d_k();
451 
452  // Aerosol scaffolding (no-op unless m_do_aerosol_rad enabled allocation above)
453  aero_tau_sw = real3d_k();
454  aero_ssa_sw = real3d_k();
455  aero_g_sw = real3d_k();
456  aero_tau_lw = real3d_k();
457 }

◆ finalize_impl()

void Radiation::finalize_impl ( amrex::Vector< amrex::MultiFab * > &  lsm_output_ptrs)
1348 {
1349  // Reset gas concentrations (k-dist data persists across steps)
1350  m_gas_concs.reset();
1351 
1352  // Fill the AMReX MFs from Kokkos Views
1353  kokkos_buffers_to_mf(lsm_output_ptrs);
1354 
1355  // Write fluxes if requested
1357 
1358  // Fill output data for datalog before deallocating
1359  if (datalog_int > 0) {
1362  Kokkos::fence();
1364  }
1365 
1366  // Deallocate the buffer arrays
1367  dealloc_buffers();
1368 }
int datalog_int
Definition: ERF_RadiationInterface.H:88
void dealloc_buffers()
Definition: ERF_Radiation.cpp:367
void populateDatalogMF()
Definition: ERF_Radiation.cpp:792
void kokkos_buffers_to_mf(amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
Definition: ERF_Radiation.cpp:673
GasConcsK< amrex::Real, layout_t, KokkosDefaultDevice > m_gas_concs
Definition: ERF_Radiation.H:319
void write_rrtmgp_fluxes()
Definition: ERF_Radiation.cpp:752
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.

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

◆ get_lsm_output_varnames()

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

Reimplemented from IRadiation.

198  {
199  return m_lsm_output_names;
200  }
amrex::Vector< std::string > m_lsm_output_names
Definition: ERF_Radiation.H:255

◆ Init()

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

Implements IRadiation.

69  {
70  // Ensure the boxes span klo -> khi
71  int klo = geom.Domain().smallEnd(2);
72  int khi = geom.Domain().bigEnd(2);
73 
74  // Reset vector of offsets for columnar data
75  m_nlay = geom.Domain().length(2);
76 
77  m_ncol = 0;
78  m_col_offsets.clear();
79  m_col_offsets.resize(int(ba.size()));
80  for (amrex::MFIter mfi(*cons_in); mfi.isValid(); ++mfi) {
81  const amrex::Box& vbx = mfi.validbox();
82  AMREX_ALWAYS_ASSERT_WITH_MESSAGE((klo == vbx.smallEnd(2)) &&
83  (khi == vbx.bigEnd(2)),
84  "Vertical decomposition with radiation is not allowed.");
85  int nx = vbx.length(0);
86  int ny = vbx.length(1);
87  m_col_offsets[mfi.index()] = m_ncol;
88  m_ncol += nx * ny;
89  }
90 
91  m_ncol_chunk = amrex::min(m_ncol, m_ncol_chunk);
92  };
const int khi
Definition: ERF_InitCustomPert_Bubble.H:21
amrex::Vector< int > m_col_offsets
Definition: ERF_Radiation.H:329

◆ initialize_impl()

void Radiation::initialize_impl ( )
1034 {
1035  // Initialize gas concentrations for this step
1037 
1038  // Load k-distribution and cloud optics data only once.
1039  // These are static lookup tables that never change.
1040  // Size the memory pool for m_ncol_chunk (not min with current m_ncol) so that
1041  // the pool remains valid even if m_ncol grows after regridding/load balancing.
1042  if (!rrtmgp::initialized) {
1043  gas_concs_t gas_concs_pool;
1044  gas_concs_pool.init(gas_names_offset, m_ncol_chunk, m_nlay);
1045  rrtmgp::rrtmgp_initialize(gas_concs_pool,
1048  m_rad_nvar);
1049  gas_concs_pool.reset();
1050  }
1051 }
GasConcsK< RealT, layout_t, KokkosDefaultDevice > gas_concs_t
Definition: ERF_RRTMGP_Interface.H:32
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, const int &nvar)
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)
674 {
675  // Heating rate, fluxes, zenith, lsm ptrs
676 
677  Table2D<Real,Order::C> p_lay_tab(p_lay.data(), {0,0}, {static_cast<int>(p_lay.extent(0)),static_cast<int>(p_lay.extent(1))});
678  Table2D<Real,Order::C> sw_heating_tab(sw_heating.data(), {0,0}, {static_cast<int>(sw_heating.extent(0)),static_cast<int>(sw_heating.extent(1))});
679  Table2D<Real,Order::C> lw_heating_tab(lw_heating.data(), {0,0}, {static_cast<int>(lw_heating.extent(0)),static_cast<int>(lw_heating.extent(1))});
680  Table2D<Real,Order::C> sw_flux_up_tab(sw_flux_up.data(), {0,0}, {static_cast<int>(sw_flux_up.extent(0)),static_cast<int>(sw_flux_up.extent(1))});
681  Table2D<Real,Order::C> sw_flux_dn_tab(sw_flux_dn.data(), {0,0}, {static_cast<int>(sw_flux_dn.extent(0)),static_cast<int>(sw_flux_dn.extent(1))});
682  Table2D<Real,Order::C> lw_flux_up_tab(lw_flux_up.data(), {0,0}, {static_cast<int>(lw_flux_up.extent(0)),static_cast<int>(lw_flux_up.extent(1))});
683  Table2D<Real,Order::C> lw_flux_dn_tab(lw_flux_dn.data(), {0,0}, {static_cast<int>(lw_flux_dn.extent(0)),static_cast<int>(lw_flux_dn.extent(1))});
684 
685  TableData<Real,1> sfc_flux_sw_dn; sfc_flux_sw_dn.resize({0}, {static_cast<int>(sw_flux_dn.extent(0))});
686  TableData<Real,1> sfc_flux_lw_dn; sfc_flux_lw_dn.resize({0}, {static_cast<int>(lw_flux_dn.extent(0))});
687  Table1D<Real> sfc_flux_sw_dn_tab = sfc_flux_sw_dn.table();
688  Table1D<Real> sfc_flux_lw_dn_tab = sfc_flux_lw_dn.table();
689  Table1D<Real> sfc_flux_sw_dir_vis_tab(sfc_flux_dir_vis.data(), {0}, {static_cast<int>(sfc_flux_dir_vis.extent(0))});
690  Table1D<Real> sfc_flux_sw_dir_nir_tab(sfc_flux_dir_nir.data(), {0}, {static_cast<int>(sfc_flux_dir_nir.extent(0))});
691  Table1D<Real> sfc_flux_sw_dif_vis_tab(sfc_flux_dif_vis.data(), {0}, {static_cast<int>(sfc_flux_dif_vis.extent(0))});
692  Table1D<Real> sfc_flux_sw_dif_nir_tab(sfc_flux_dif_nir.data(), {0}, {static_cast<int>(sfc_flux_dif_nir.extent(0))});
693  Table1D<Real> mu0_tab(mu0.data(), {0}, {static_cast<int>(mu0.extent(0))});
694  Vector<Table1D<Real>> rrtmgp_out_vars = {mu0_tab , sfc_flux_sw_dn_tab ,
695  sfc_flux_sw_dir_vis_tab, sfc_flux_sw_dir_nir_tab,
696  sfc_flux_sw_dif_vis_tab, sfc_flux_sw_dif_nir_tab,
697  sfc_flux_lw_dn_tab };
698 
699  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
700  const auto& vbx = mfi.validbox();
701  const auto& sbx = makeSlab(vbx,2,vbx.smallEnd(2));
702  const int nx = vbx.length(0);
703  const int imin = vbx.smallEnd(0);
704  const int jmin = vbx.smallEnd(1);
705  const int offset = m_col_offsets[mfi.index()];
706  const Array4<Real>& q_arr = m_qheating_rates->array(mfi);
707  const Array4<Real>& f_arr = m_rad_fluxes->array(mfi);
708  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
709  {
710  // map [i,j,k] 0-based to [icol, ilay] 0-based
711  const int icol = (j-jmin)*nx + (i-imin) + offset;
712  const int ilay = k;
713 
714  // Temperature heating rate for SW and LW
715  q_arr(i,j,k,0) = sw_heating_tab(icol,ilay);
716  q_arr(i,j,k,1) = lw_heating_tab(icol,ilay);
717 
718  // Convert the dT/dz to dTheta/dz
719  Real iexner = one/getExnergivenP(Real(p_lay_tab(icol,ilay)), R_d/Cp_d);
720  q_arr(i,j,k,0) *= iexner;
721  q_arr(i,j,k,1) *= iexner;
722 
723  // Populate the fluxes
724  f_arr(i,j,k,0) = sw_flux_up_tab(icol,ilay);
725  f_arr(i,j,k,1) = sw_flux_dn_tab(icol,ilay);
726  f_arr(i,j,k,2) = lw_flux_up_tab(icol,ilay);
727  f_arr(i,j,k,3) = lw_flux_dn_tab(icol,ilay);
728 
729  if (k==0) {
730  sfc_flux_sw_dn_tab(icol) = sw_flux_dn_tab(icol,ilay);
731  sfc_flux_lw_dn_tab(icol) = lw_flux_dn_tab(icol,ilay);
732  }
733  });
734  for (int ivar(0); ivar<lsm_output_ptrs.size(); ivar++) {
735  if (lsm_output_ptrs[ivar]) {
736  auto rrtmgp_for_fill = rrtmgp_out_vars[ivar];
737  const Array4<Real>& lsm_out_arr = lsm_output_ptrs[ivar]->array(mfi);
738  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
739  {
740  // map [i,j,k] 0-based to [icol, ilay] 0-based
741  const int icol = (j-jmin)*nx + (i-imin) + offset;
742 
743  // export the desired variable at surface
744  lsm_out_arr(i,j,k) = rrtmgp_for_fill(icol);
745  });
746  } // valid ptr
747  } // ivar
748  }// mfi
749 }
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:31
constexpr amrex::Real one
Definition: ERF_Constants.H:7
constexpr amrex::Real R_d
Definition: ERF_Constants.H:29
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getExnergivenP(const amrex::Real P, const amrex::Real rdOcp)
Definition: ERF_EOS.H:141
ParallelFor(grown_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) { qrcuten_arr(i, j, k)=Real(0);qscuten_arr(i, j, k)=Real(0);qicuten_arr(i, j, k)=Real(0);})
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:264
amrex::MultiFab * m_rad_fluxes
Definition: ERF_Radiation.H:270
amrex::MultiFab * m_qheating_rates
Definition: ERF_Radiation.H:267
Here is the call graph for this function:

◆ mf_to_kokkos_buffers()

void Radiation::mf_to_kokkos_buffers ( amrex::iMultiFab *  lmask,
amrex::MultiFab *  t_surf,
amrex::Vector< amrex::MultiFab * > &  lsm_input_ptrs 
)
464 {
465  Table2D<Real,Order::C> r_lay_tab(r_lay.data(), {0,0}, {static_cast<int>(r_lay.extent(0)),static_cast<int>(r_lay.extent(1))});
466  Table2D<Real,Order::C> p_lay_tab(p_lay.data(), {0,0}, {static_cast<int>(p_lay.extent(0)),static_cast<int>(p_lay.extent(1))});
467  Table2D<Real,Order::C> t_lay_tab(t_lay.data(), {0,0}, {static_cast<int>(t_lay.extent(0)),static_cast<int>(t_lay.extent(1))});
468  Table2D<Real,Order::C> z_del_tab(z_del.data(), {0,0}, {static_cast<int>(z_del.extent(0)),static_cast<int>(z_del.extent(1))});
469  Table2D<Real,Order::C> qv_lay_tab(qv_lay.data(), {0,0}, {static_cast<int>(qv_lay.extent(0)),static_cast<int>(qv_lay.extent(1))});
470  Table2D<Real,Order::C> qc_lay_tab(qc_lay.data(), {0,0}, {static_cast<int>(qc_lay.extent(0)),static_cast<int>(qc_lay.extent(1))});
471  Table2D<Real,Order::C> qi_lay_tab(qi_lay.data(), {0,0}, {static_cast<int>(qi_lay.extent(0)),static_cast<int>(qi_lay.extent(1))});
472  Table2D<Real,Order::C> cldfrac_tot_tab(cldfrac_tot.data(), {0,0}, {static_cast<int>(cldfrac_tot.extent(0)),static_cast<int>(cldfrac_tot.extent(1))});
473 
474  Table2D<Real,Order::C> lwp_tab(lwp.data(), {0,0}, {static_cast<int>(lwp.extent(0)),static_cast<int>(lwp.extent(1))});
475  Table2D<Real,Order::C> iwp_tab(iwp.data(), {0,0}, {static_cast<int>(iwp.extent(0)),static_cast<int>(iwp.extent(1))});
476  Table2D<Real,Order::C> eff_radius_qc_tab(eff_radius_qc.data(), {0,0}, {static_cast<int>(eff_radius_qc.extent(0)),static_cast<int>(eff_radius_qc.extent(1))});
477  Table2D<Real,Order::C> eff_radius_qi_tab(eff_radius_qi.data(), {0,0}, {static_cast<int>(eff_radius_qi.extent(0)),static_cast<int>(eff_radius_qi.extent(1))});
478 
479  Table2D<Real,Order::C> p_lev_tab(p_lev.data(), {0,0}, {static_cast<int>(p_lev.extent(0)),static_cast<int>(p_lev.extent(1))});
480  Table2D<Real,Order::C> t_lev_tab(t_lev.data(), {0,0}, {static_cast<int>(t_lev.extent(0)),static_cast<int>(t_lev.extent(1))});
481 
482  Table1D<Real> lat_tab(lat.data(), {0}, {static_cast<int>(lat.extent(0))});
483  Table1D<Real> lon_tab(lon.data(), {0}, {static_cast<int>(lon.extent(0))});
484  Table1D<Real> t_sfc_tab(t_sfc.data(), {0}, {static_cast<int>(t_sfc.extent(0))});
485 
486  bool moist = m_moist;
487  bool ice = m_ice;
488  const bool has_lsm = m_lsm;
489  const bool has_lat = m_lat;
490  const bool has_lon = m_lon;
491  const bool has_surflayer = (t_surf);
492  int ncol = m_ncol;
493  int nlay = m_nlay;
494  Real dz = m_geom.CellSize(2);
495  Real cons_lat = m_lat_cons;
496  Real cons_lon = m_lon_cons;
497  Real rad_t_sfc = m_rad_t_sfc;
498 
499  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
500  const auto& vbx = mfi.validbox();
501  const int nx = vbx.length(0);
502  const int imin = vbx.smallEnd(0);
503  const int jmin = vbx.smallEnd(1);
504  const int offset = m_col_offsets[mfi.index()];
505  const Array4<const Real>& cons_arr = m_cons_in->const_array(mfi);
506  const Array4<const Real>& z_arr = (m_z_phys) ? m_z_phys->const_array(mfi) :
507  Array4<const Real>{};
508  const Array4<const Real>& lat_arr = (m_lat) ? m_lat->const_array(mfi) :
509  Array4<const Real>{};
510  const Array4<const Real>& lon_arr = (m_lon) ? m_lon->const_array(mfi) :
511  Array4<const Real>{};
512  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
513  {
514  // map [i,j,k] 0-based to [icol, ilay] 0-based
515  const int icol = (j-jmin)*nx + (i-imin) + offset;
516  const int ilay = k;
517 
518  // EOS input (at CC)
519  Real r = cons_arr(i,j,k,Rho_comp);
520  Real rt = cons_arr(i,j,k,RhoTheta_comp);
521  Real qv = (moist) ? std::max(cons_arr(i,j,k,RhoQ1_comp)/r,Real(0.)) : Real(0.);
522  Real qc = (moist) ? std::max(cons_arr(i,j,k,RhoQ2_comp)/r,Real(0.)) : Real(0.);
523  Real qi = (ice) ? std::max(cons_arr(i,j,k,RhoQ3_comp)/r,Real(0.)) : Real(0.);
524 
525  // EOS avg to z-face
526  Real r_lo = cons_arr(i,j,k-1,Rho_comp);
527  Real rt_lo = cons_arr(i,j,k-1,RhoTheta_comp);
528  Real qv_lo = (moist) ? cons_arr(i,j,k-1,RhoQ1_comp)/r_lo : Real(0.);
529  Real dz_k = (z_arr) ? Real(0.125) * ( (z_arr(i ,j ,k+1) - z_arr(i ,j ,k))
530  + (z_arr(i+1,j ,k+1) - z_arr(i+1,j ,k))
531  + (z_arr(i ,j+1,k+1) - z_arr(i ,j+1,k))
532  + (z_arr(i+1,j+1,k+1) - z_arr(i+1,j+1,k)) ) : Real(0.5)*dz; // Dist from w-face to CC at k
533  Real dz_km1 = (z_arr) ? Real(0.125) * ( (z_arr(i ,j ,k ) - z_arr(i ,j ,k-1))
534  + (z_arr(i+1,j ,k ) - z_arr(i+1,j ,k-1))
535  + (z_arr(i ,j+1,k ) - z_arr(i ,j+1,k-1))
536  + (z_arr(i+1,j+1,k ) - z_arr(i+1,j+1,k-1)) ) : Real(0.5)*dz; // Dist from w-face to CC at k-1
537  Real r_avg = (dz_k*r + dz_km1*r_lo ) / (dz_k + dz_km1);
538  Real rt_avg = (dz_k*rt + dz_km1*rt_lo) / (dz_k + dz_km1);
539  Real qv_avg = (dz_k*qv + dz_km1*qv_lo) / (dz_k + dz_km1);
540 
541  // Views at CC
542  r_lay_tab(icol,ilay) = r;
543 
544  p_lay_tab(icol,ilay) = getPgivenRTh(rt, qv);
545  t_lay_tab(icol,ilay) = getTgivenRandRTh(r, rt, qv);
546  z_del_tab(icol,ilay) = (z_arr) ? Real(0.25) * ( (z_arr(i ,j ,k+1) - z_arr(i ,j ,k))
547  + (z_arr(i+1,j ,k+1) - z_arr(i+1,j ,k))
548  + (z_arr(i ,j+1,k+1) - z_arr(i ,j+1,k))
549  + (z_arr(i+1,j+1,k+1) - z_arr(i+1,j+1,k)) ) : dz;
550  qv_lay_tab(icol,ilay) = qv;
551  qc_lay_tab(icol,ilay) = qc;
552  qi_lay_tab(icol,ilay) = qi;
553  cldfrac_tot_tab(icol,ilay) = ((qc+qi)>Real(0.)) ? Real(1.) : Real(0.);
554 
555  // NOTE: These are populated in 'mixing_ratio_to_cloud_mass'
556  lwp_tab(icol,ilay) = Real(0.);
557  iwp_tab(icol,ilay) = Real(0.);
558 
559  // NOTE: These would be populated from P3 (we use the constants in p3_main_impl.hpp)
560  // NOTE: These are in units of micron!
561  eff_radius_qc_tab(icol,ilay) = (qc>Real(0.)) ? Real(10.0) : Real(0.);
562  eff_radius_qi_tab(icol,ilay) = (qi>Real(0.)) ? Real(25.0) : Real(0.);
563 
564  // Buffers on z-faces (nlay+1)
565  p_lev_tab(icol,ilay) = getPgivenRTh(rt_avg, qv_avg);
566  t_lev_tab(icol,ilay) = getTgivenRandRTh(r_avg, rt_avg, qv_avg);
567  if (ilay==(nlay-1)) {
568  Real r_hi = cons_arr(i,j,k+1,Rho_comp);
569  Real rt_hi = cons_arr(i,j,k+1,RhoTheta_comp);
570  Real qv_hi = (moist) ? std::max(cons_arr(i,j,k+1,RhoQ1_comp)/r_hi,Real(0.)) : Real(0.);
571  Real dz_kp1 = (z_arr) ? Real(0.125) * ( (z_arr(i ,j ,k+2) - z_arr(i ,j ,k+1))
572  + (z_arr(i+1,j ,k+2) - z_arr(i+1,j ,k+1))
573  + (z_arr(i ,j+1,k+2) - z_arr(i ,j+1,k+1))
574  + (z_arr(i+1,j+1,k+2) - z_arr(i+1,j+1,k+1)) ) : Real(0.5)*dz; // Dist from w-face to CC at k+1
575  r_avg = (dz_k*r + dz_kp1*r_hi ) / (dz_k + dz_kp1);
576  rt_avg = (dz_k*rt + dz_kp1*rt_hi) / (dz_k + dz_kp1);
577  qv_avg = (dz_k*qv + dz_kp1*qv_hi) / (dz_k + dz_kp1);
578  p_lev_tab(icol,ilay+1) = getPgivenRTh(rt_avg, qv_avg);
579  t_lev_tab(icol,ilay+1) = getTgivenRandRTh(r_avg, rt_avg, qv_avg);
580  }
581 
582  // 1D data structures
583  if (k==0) {
584  lat_tab(icol) = (has_lat) ? lat_arr(i,j,0) : cons_lat;
585  lon_tab(icol) = (has_lon) ? lon_arr(i,j,0) : cons_lon;
586  }
587 
588  });
589  } // mfi
590 
591  // Populate vars LSM would provide
592  if (!has_lsm && !has_surflayer) {
593  // Parsed surface temp
594  Kokkos::deep_copy(t_sfc, rad_t_sfc);
595 
596  // EAMXX dummy atmos constants
597  Kokkos::deep_copy(sfc_alb_dir_vis, Real(0.06));
598  Kokkos::deep_copy(sfc_alb_dir_nir, Real(0.06));
599  Kokkos::deep_copy(sfc_alb_dif_vis, Real(0.06));
600  Kokkos::deep_copy(sfc_alb_dif_nir, Real(0.06));
601 
602  // AML NOTE: These are not used in current EAMXX, I've left
603  // the code to plug into these if we need it.
604  //
605  // Current EAMXX constants
606  Kokkos::deep_copy(sfc_emis, Real(0.98));
607  Kokkos::deep_copy(lw_src , zero );
608  } else {
609  Vector<real1d_k> rrtmgp_in_vars = {t_sfc, sfc_emis,
612  Vector<Real> rrtmgp_default_vals = {rad_t_sfc, Real(0.98),
613  Real(0.06), Real(0.06),
614  Real(0.06), Real(0.06)};
615  for (int ivar(0); ivar<lsm_input_ptrs.size(); ivar++) {
616  auto rrtmgp_default_val = rrtmgp_default_vals[ivar];
617  auto rrtmgp_to_fill_k = rrtmgp_in_vars[ivar];
618  amrex::Table1D<amrex::Real> rrtmgp_to_fill(rrtmgp_to_fill_k.data(),
619  0, rrtmgp_to_fill_k.extent(0));
620  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
621  const auto& vbx = mfi.validbox();
622  const auto& sbx = makeSlab(vbx,2,vbx.smallEnd(2));
623  const int nx = vbx.length(0);
624  const int imin = vbx.smallEnd(0);
625  const int jmin = vbx.smallEnd(1);
626  const int offset = m_col_offsets[mfi.index()];
627  const Array4<const int>& lmask_arr = (lmask) ? lmask->const_array(mfi) :
628  Array4<const int> {};
629  const Array4<const Real>& tsurf_arr = (t_surf) ? t_surf->const_array(mfi) :
630  Array4<const Real> {};
631  const Array4<const Real>& lsm_in_arr = (lsm_input_ptrs[ivar]) ? lsm_input_ptrs[ivar]->const_array(mfi) :
632  Array4<const Real> {};
633  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
634  {
635  // map [i,j,k] 0-based to [icol, ilay] 0-based
636  const int icol = (j-jmin)*nx + (i-imin) + offset;
637 
638  // Check if over land
639  bool is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
640 
641  // Check if valid LSM data
642  bool valid_lsm_data{false};
643  if (lsm_in_arr) { valid_lsm_data = (lsm_in_arr(i,j,k) >= Real(0.)); }
644 
645  // Have LSM and are over land
646  if (is_land && valid_lsm_data) {
647  rrtmgp_to_fill(icol) = lsm_in_arr(i,j,k);
648  }
649  // We have a SurfLayer (enforce consistency with temperature)
650  else if (tsurf_arr && (ivar==0)) {
651  rrtmgp_to_fill(icol) = tsurf_arr(i,j,k);
652  }
653  // Use the default value
654  else {
655  rrtmgp_to_fill(icol) = rrtmgp_default_val;
656  }
657  });
658  } //mfi
659  } // ivar
660  Kokkos::deep_copy(lw_src, zero );
661  } // have lsm
662 
663  // Enforce consistency between t_sfc and t_lev at bottom surface
664  Kokkos::parallel_for(Kokkos::RangePolicy(0, ncol),
665  KOKKOS_LAMBDA (int icol)
666  {
667  t_lev_tab(icol,0) = t_sfc_tab(icol);
668  });
669 }
constexpr amrex::Real zero
Definition: ERF_Constants.H:6
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:46
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=amrex::Real(0))
Definition: ERF_EOS.H:81
#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:277
amrex::MultiFab * m_z_phys
Definition: ERF_Radiation.H:273
amrex::Geometry m_geom
Definition: ERF_Radiation.H:229
amrex::MultiFab * m_lat
Definition: ERF_Radiation.H:276
@ t_surf
Definition: ERF_OceanSurf.H:14
@ qv
Definition: ERF_Kessler.H:29
@ qc
Definition: ERF_SatAdj.H:40
@ qi
Definition: ERF_WSM6.H:26
@ dz
Definition: ERF_AdvanceWSM6.cpp:104
Here is the call graph for this function:

◆ populateDatalogMF()

void Radiation::populateDatalogMF ( )
793 {
794  Table2D<Real,Order::C> sw_flux_up_tab(sw_flux_up.data(), {0,0}, {static_cast<int>(sw_flux_up.extent(0)),static_cast<int>(sw_flux_up.extent(1))});
795  Table2D<Real,Order::C> sw_flux_dn_tab(sw_flux_dn.data(), {0,0}, {static_cast<int>(sw_flux_dn.extent(0)),static_cast<int>(sw_flux_dn.extent(1))});
796  Table2D<Real,Order::C> sw_flux_dn_dir_tab(sw_flux_dn_dir.data(), {0,0}, {static_cast<int>(sw_flux_dn_dir.extent(0)),static_cast<int>(sw_flux_dn_dir.extent(1))});
797  Table2D<Real,Order::C> lw_flux_up_tab(lw_flux_up.data(), {0,0}, {static_cast<int>(lw_flux_up.extent(0)),static_cast<int>(lw_flux_up.extent(1))});
798  Table2D<Real,Order::C> lw_flux_dn_tab(lw_flux_dn.data(), {0,0}, {static_cast<int>(lw_flux_dn.extent(0)),static_cast<int>(lw_flux_dn.extent(1))});
799 
800  Table2D<Real,Order::C> sw_clrsky_flux_up_tab(sw_clrsky_flux_up.data(), {0,0},
801  {static_cast<int>(sw_clrsky_flux_up.extent(0)),static_cast<int>(sw_clrsky_flux_up.extent(1))});
802  Table2D<Real,Order::C> sw_clrsky_flux_dn_tab(sw_clrsky_flux_dn.data(), {0,0},
803  {static_cast<int>(sw_clrsky_flux_dn.extent(0)),static_cast<int>(sw_clrsky_flux_dn.extent(1))});
804  Table2D<Real,Order::C> sw_clrsky_flux_dn_dir_tab(sw_clrsky_flux_dn_dir.data(), {0,0},
805  {static_cast<int>(sw_clrsky_flux_dn_dir.extent(0)),static_cast<int>(sw_clrsky_flux_dn_dir.extent(1))});
806  Table2D<Real,Order::C> lw_clrsky_flux_up_tab(lw_clrsky_flux_up.data(), {0,0},
807  {static_cast<int>(lw_clrsky_flux_up.extent(0)),static_cast<int>(lw_clrsky_flux_up.extent(1))});
808  Table2D<Real,Order::C> lw_clrsky_flux_dn_tab(lw_clrsky_flux_dn.data(), {0,0},
809  {static_cast<int>(lw_clrsky_flux_dn.extent(0)),static_cast<int>(lw_clrsky_flux_dn.extent(1))});
810  Table2D<Real,Order::C> sw_clrsky_heating_tab(sw_clrsky_heating.data(), {0,0},
811  {static_cast<int>(sw_clrsky_heating.extent(0)),static_cast<int>(sw_clrsky_heating.extent(1))});
812  Table2D<Real,Order::C> lw_clrsky_heating_tab(lw_clrsky_heating.data(), {0,0},
813  {static_cast<int>(lw_clrsky_heating.extent(0)),static_cast<int>(lw_clrsky_heating.extent(1))});
814  Table2D<Real,Order::C> sw_clnsky_flux_up_tab(sw_clnsky_flux_up.data(), {0,0},
815  {static_cast<int>(sw_clnsky_flux_up.extent(0)),static_cast<int>(sw_clnsky_flux_up.extent(1))});
816  Table2D<Real,Order::C> sw_clnsky_flux_dn_tab(sw_clnsky_flux_dn.data(), {0,0},
817  {static_cast<int>(sw_clnsky_flux_dn.extent(0)),static_cast<int>(sw_clnsky_flux_dn.extent(1))});
818  Table2D<Real,Order::C> sw_clnsky_flux_dn_dir_tab(sw_clnsky_flux_dn_dir.data(), {0,0},
819  {static_cast<int>(sw_clnsky_flux_dn_dir.extent(0)),static_cast<int>(sw_clnsky_flux_dn_dir.extent(1))});
820  Table2D<Real,Order::C> lw_clnsky_flux_up_tab(lw_clnsky_flux_up.data(), {0,0},
821  {static_cast<int>(lw_clnsky_flux_up.extent(0)),static_cast<int>(lw_clnsky_flux_up.extent(1))});
822  Table2D<Real,Order::C> lw_clnsky_flux_dn_tab(lw_clnsky_flux_dn.data(), {0,0},
823  {static_cast<int>(lw_clnsky_flux_dn.extent(0)),static_cast<int>(lw_clnsky_flux_dn.extent(1))});
824  Table2D<Real,Order::C> sw_clnclrsky_flux_up_tab(sw_clnclrsky_flux_up.data(), {0,0},
825  {static_cast<int>(sw_clnclrsky_flux_up.extent(0)),static_cast<int>(sw_clnclrsky_flux_up.extent(1))});
826  Table2D<Real,Order::C> sw_clnclrsky_flux_dn_tab(sw_clnclrsky_flux_dn.data(), {0,0},
827  {static_cast<int>(sw_clnclrsky_flux_dn.extent(0)),static_cast<int>(sw_clnclrsky_flux_dn.extent(1))});
828  Table2D<Real,Order::C> sw_clnclrsky_flux_dn_dir_tab(sw_clnclrsky_flux_dn_dir.data(), {0,0},
829  {static_cast<int>(sw_clnclrsky_flux_dn_dir.extent(0)),static_cast<int>(sw_clnclrsky_flux_dn_dir.extent(1))});
830  Table2D<Real,Order::C> lw_clnclrsky_flux_up_tab(lw_clnclrsky_flux_up.data(), {0,0},
831  {static_cast<int>(lw_clnclrsky_flux_up.extent(0)),static_cast<int>(lw_clnclrsky_flux_up.extent(1))});
832  Table2D<Real,Order::C> lw_clnclrsky_flux_dn_tab(lw_clnclrsky_flux_dn.data(), {0,0},
833  {static_cast<int>(lw_clnclrsky_flux_dn.extent(0)),static_cast<int>(lw_clnclrsky_flux_dn.extent(1))});
834 
835  Table1D<Real> mu0_tab(mu0.data(), {0}, {static_cast<int>(mu0.extent(0))});
836 
837  auto extra_clnsky_diag = m_extra_clnsky_diag;
838  auto extra_clnclrsky_diag = m_extra_clnclrsky_diag;
839 
840  for (MFIter mfi(datalog_mf); mfi.isValid(); ++mfi) {
841  const auto& vbx = mfi.validbox();
842  const int nx = vbx.length(0);
843  const int imin = vbx.smallEnd(0);
844  const int jmin = vbx.smallEnd(1);
845  const int offset = m_col_offsets[mfi.index()];
846  const Array4<Real>& dst_arr = datalog_mf.array(mfi);
847  const Array4<Real>& q_arr = m_qheating_rates->array(mfi);
848  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
849  {
850  // map [i,j,k] 0-based to [icol, ilay] 0-based
851  const int icol = (j-jmin)*nx + (i-imin) + offset;
852  const int ilay = k;
853 
854  dst_arr(i,j,k,0) = q_arr(i, j, k, 0);
855  dst_arr(i,j,k,1) = q_arr(i, j, k, 1);
856 
857  // SW and LW fluxes
858  dst_arr(i,j,k,2) = sw_flux_up_tab(icol,ilay);
859  dst_arr(i,j,k,3) = sw_flux_dn_tab(icol,ilay);
860  dst_arr(i,j,k,4) = sw_flux_dn_dir_tab(icol,ilay);
861  dst_arr(i,j,k,5) = lw_flux_up_tab(icol,ilay);
862  dst_arr(i,j,k,6) = lw_flux_dn_tab(icol,ilay);
863 
864  // Cosine zenith angle
865  dst_arr(i,j,k,7) = mu0_tab(icol);
866 
867  // Clear sky heating rates and fluxes:
868  dst_arr(i,j,k,8) = sw_clrsky_heating_tab(icol, ilay);
869  dst_arr(i,j,k,9) = lw_clrsky_heating_tab(icol, ilay);
870 
871  dst_arr(i,j,k,10) = sw_clrsky_flux_up_tab(icol,ilay);
872  dst_arr(i,j,k,11) = sw_clrsky_flux_dn_tab(icol,ilay);
873  dst_arr(i,j,k,12) = sw_clrsky_flux_dn_dir_tab(icol,ilay);
874  dst_arr(i,j,k,13) = lw_clrsky_flux_up_tab(icol,ilay);
875  dst_arr(i,j,k,14) = lw_clrsky_flux_dn_tab(icol,ilay);
876 
877  // Clean sky fluxes:
878  if (extra_clnsky_diag) {
879  dst_arr(i,j,k,15) = sw_clnsky_flux_up_tab(icol,ilay);
880  dst_arr(i,j,k,16) = sw_clnsky_flux_dn_tab(icol,ilay);
881  dst_arr(i,j,k,17) = sw_clnsky_flux_dn_dir_tab(icol,ilay);
882  dst_arr(i,j,k,18) = lw_clnsky_flux_up_tab(icol,ilay);
883  dst_arr(i,j,k,19) = lw_clnsky_flux_dn_tab(icol,ilay);
884  }
885 
886  // Clean-clear sky fluxes:
887  if (extra_clnclrsky_diag) {
888  dst_arr(i,j,k,20) = sw_clnclrsky_flux_up_tab(icol,ilay);
889  dst_arr(i,j,k,21) = sw_clnclrsky_flux_dn_tab(icol,ilay);
890  dst_arr(i,j,k,22) = sw_clnclrsky_flux_dn_dir_tab(icol,ilay);
891  dst_arr(i,j,k,23) = lw_clnclrsky_flux_up_tab(icol,ilay);
892  dst_arr(i,j,k,24) = lw_clnclrsky_flux_dn_tab(icol,ilay);
893  }
894  });
895  }
896 }
amrex::MultiFab datalog_mf
Definition: ERF_Radiation.H:284
Here is the call graph for this function:

◆ rad_run_impl()

void Radiation::rad_run_impl ( amrex::Vector< amrex::MultiFab * > &  lsm_output_ptrs)
inline
175  {
176  if (m_update_rad) {
177  amrex::Print() << "Radiation advancing level " << m_lev << " at (YY-MM-DD SS) " << m_orbital_year << '-'
178  << m_orbital_mon << '-' << m_orbital_day << ' ' << m_orbital_sec << " ...";
179  this->initialize_impl();
180  this->run_impl();
181  this->finalize_impl(lsm_output_ptrs);
182  amrex::Print() << "DONE\n";
183  }
184  }
void initialize_impl()
Definition: ERF_Radiation.cpp:1033
void run_impl()
Definition: ERF_Radiation.cpp:1055
int m_orbital_mon
Definition: ERF_Radiation.H:345
void finalize_impl(amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
Definition: ERF_Radiation.cpp:1347
int m_orbital_sec
Definition: ERF_Radiation.H:347
int m_lev
Definition: ERF_Radiation.H:217
bool m_update_rad
Definition: ERF_Radiation.H:235
int m_orbital_day
Definition: ERF_Radiation.H:346

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::iMultiFab *  lmask,
amrex::MultiFab *  t_surf,
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.

112  {
113  set_grids(level, step, time, dt, ba, geom,
114  cons_in, lmask, t_surf,
115  lsm_input_ptrs, qheating_rates,
116  rad_fluxes, z_phys, lat_ptr, lon_ptr);
117  rad_run_impl(lsm_output_ptrs);
118  }
void set_grids(int &level, int &step, amrex::Real &time, const amrex::Real &dt, const amrex::BoxArray &ba, amrex::Geometry &geom, amrex::MultiFab *cons_in, amrex::iMultiFab *lmask, amrex::MultiFab *t_surf, amrex::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:158
void rad_run_impl(amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
Definition: ERF_Radiation.H:174
Here is the call graph for this function:

◆ run_impl()

void Radiation::run_impl ( )
1056 {
1057  // Local copies
1058  const auto ncol = m_ncol;
1059  const auto nlay = m_nlay;
1060  const auto nswbands = m_nswbands;
1061 
1062  // Compute orbital parameters; these are used both for computing
1063  // the solar zenith angle and also for computing total solar
1064  // irradiance scaling (tsi_scaling).
1065  double obliqr, lambm0, mvelpp;
1066  int orbital_year = m_orbital_year;
1067  double eccen = m_orbital_eccen;
1068  double obliq = m_orbital_obliq;
1069  double mvelp = m_orbital_mvelp;
1070  if (eccen >= 0 && obliq >= 0 && mvelp >= 0) {
1071  // fixed orbital parameters forced with orbital_year == ORB_UNDEF_INT
1072  orbital_year = ORB_UNDEF_INT;
1073  }
1074  orbital_params(orbital_year, eccen, obliq,
1075  mvelp, obliqr, lambm0, mvelpp);
1076 
1077  // Use the orbital parameters to calculate the solar declination and eccentricity factor
1078  double delta, eccf;
1079  // Want day + fraction; calday 1 == Jan 1 0Z
1080  static constexpr double dpy[] = {zero , Real(31.0), Real(59.0), Real(90.0), Real(120.0), Real(151.0),
1081  Real(181.0), Real(212.0), Real(243.0), Real(273.0), Real(304.0), Real(334.0)};
1082  bool leap = (m_orbital_year % 4 == 0 && (!(m_orbital_year % 100 == 0) || (m_orbital_year % 400 == 0))) ? true : false;
1083  double calday = one + dpy[m_orbital_mon-1] + (m_orbital_day-one) + m_orbital_sec/Real(86400.0);
1084  // add extra day if leap year and past February
1085  if (leap && m_orbital_mon>2) { calday += one; }
1086  orbital_decl(calday, eccen, mvelpp, lambm0, obliqr, delta, eccf);
1087 
1088  // Overwrite eccf if using a fixed solar constant.
1089  auto fixed_total_solar_irradiance = m_fixed_total_solar_irradiance;
1090  if (fixed_total_solar_irradiance >= 0){
1091  eccf = fixed_total_solar_irradiance/Real(1360.9);
1092  }
1093 
1094  // Precompute volume mixing ratio (VMR) for all gases
1095  //
1096  // H2O is obtained from qv.
1097  // O3 may be a constant or a 1D vector
1098  // All other comps are set to constants for now
1099  Vector<real2d_k> vmr_full_vec(m_ngas);
1100  for (int igas(0); igas < m_ngas; ++igas) {
1101  auto name = m_gas_names[igas];
1102  vmr_full_vec[igas] = real2d_k("vmr_full_" + name, ncol, nlay);
1103  auto tmp2d = vmr_full_vec[igas];
1104  auto gas_mol_weight = m_mol_weight_gas[igas];
1105  if (name == "H2O") {
1106  auto qv_lay_d = qv_lay;
1107  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1108  KOKKOS_LAMBDA (int icol, int ilay)
1109  {
1110  tmp2d(icol,ilay) = qv_lay_d(icol,ilay) * mwdair/gas_mol_weight;
1111  });
1112  } else if (name == "CO2") {
1113  Kokkos::deep_copy(tmp2d, m_co2vmr);
1114  } else if (name == "O3") {
1115  if (m_o3_size==1) {
1116  Kokkos::deep_copy(tmp2d, m_o3vmr[0] );
1117  } else {
1118  auto o3_lay_d = o3_lay;
1119  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1120  KOKKOS_LAMBDA (int icol, int ilay)
1121  {
1122  tmp2d(icol,ilay) = o3_lay_d(ilay);
1123  });
1124  }
1125  } else if (name == "N2O") {
1126  Kokkos::deep_copy(tmp2d, m_n2ovmr);
1127  } else if (name == "CO") {
1128  Kokkos::deep_copy(tmp2d, m_covmr );
1129  } else if (name == "CH4") {
1130  Kokkos::deep_copy(tmp2d, m_ch4vmr);
1131  } else if (name == "O2") {
1132  Kokkos::deep_copy(tmp2d, m_o2vmr );
1133  } else if (name == "N2") {
1134  Kokkos::deep_copy(tmp2d, m_n2vmr );
1135  } else {
1136  Abort("Radiation: Unknown gas component.");
1137  }
1138 
1139  // Populate GasConcs object
1140  m_gas_concs.set_vmr(name, tmp2d);
1141  Kokkos::fence();
1142  }
1143 
1144  // Populate mu0 1D array
1145  // This must be done on HOST and copied to device.
1146  auto h_mu0 = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), mu0);
1147  if (m_fixed_solar_zenith_angle > 0) {
1148  Kokkos::deep_copy(h_mu0, m_fixed_solar_zenith_angle);
1149  } else {
1150  auto h_lat = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), lat);
1151  auto h_lon = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), lon);
1152  double dt = double(m_dt);
1153  auto rad_freq_in_steps = m_rad_freq_in_steps;
1154  Kokkos::parallel_for(Kokkos::RangePolicy<Kokkos::Serial>(0, ncol),
1155  [&] (int icol)
1156  {
1157  // Convert lat/lon to radians
1158  double lat_col = h_lat(icol)*PI/Real(180.0);
1159  double lon_col = h_lon(icol)*PI/Real(180.0);
1160  double lcalday = calday;
1161  double ldelta = delta;
1162  double dt_avg = static_cast<double>(rad_freq_in_steps) * dt;
1163  h_mu0(icol) = Real(orbital_cos_zenith(lcalday, lat_col, lon_col, ldelta, dt_avg));
1164  });
1165  }
1166  Kokkos::deep_copy(mu0, h_mu0);
1167 
1168  // Compute layer cloud mass per unit area (populates lwp/iwp)
1171 
1172  // Convert to g/m2 (needed by RRTMGP)
1173  Table2D<Real,Order::C> lwp_tab(lwp.data(), {0,0}, {static_cast<int>(lwp.extent(0)),static_cast<int>(lwp.extent(1))});
1174  Table2D<Real,Order::C> iwp_tab(iwp.data(), {0,0}, {static_cast<int>(iwp.extent(0)),static_cast<int>(iwp.extent(1))});
1175  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1176  KOKKOS_LAMBDA (int icol, int ilay)
1177  {
1178  lwp_tab(icol,ilay) *= Real(1.e3);
1179  iwp_tab(icol,ilay) *= Real(1.e3);
1180  });
1181 
1182  // -----------------------------------------------------------------------
1183  // Process radiation in column chunks to limit peak GPU memory.
1184  // Radiation columns are independent (no horizontal coupling), so
1185  // chunking produces bit-identical results.
1186  // -----------------------------------------------------------------------
1187  const int ncol_chunk = std::min(m_ncol_chunk, ncol);
1188  const int kbot = 0;
1189 
1190  for (int col_s = 0; col_s < ncol; col_s += ncol_chunk) {
1191  const int ncol_c = std::min(ncol_chunk, ncol - col_s);
1192  const int col_e = col_s + ncol_c;
1193  auto cr = std::make_pair(col_s, col_e);
1194 
1195  // --- Chunk subviews: 1D (ncol) ---
1196  real1d_k mu0_c (mu0.data() + col_s, ncol_c);
1197  real1d_k sfc_alb_dir_vis_c (sfc_alb_dir_vis.data() + col_s, ncol_c);
1198  real1d_k sfc_alb_dir_nir_c (sfc_alb_dir_nir.data() + col_s, ncol_c);
1199  real1d_k sfc_alb_dif_vis_c (sfc_alb_dif_vis.data() + col_s, ncol_c);
1200  real1d_k sfc_alb_dif_nir_c (sfc_alb_dif_nir.data() + col_s, ncol_c);
1201  real1d_k sfc_flux_dir_vis_c (sfc_flux_dir_vis.data() + col_s, ncol_c);
1202  real1d_k sfc_flux_dir_nir_c (sfc_flux_dir_nir.data() + col_s, ncol_c);
1203  real1d_k sfc_flux_dif_vis_c (sfc_flux_dif_vis.data() + col_s, ncol_c);
1204  real1d_k sfc_flux_dif_nir_c (sfc_flux_dif_nir.data() + col_s, ncol_c);
1205  real1d_k t_sfc_c (t_sfc.data() + col_s, ncol_c);
1206  real1d_k sfc_emis_c (sfc_emis.data() + col_s, ncol_c);
1207  real1d_k lw_src_c (lw_src.data() + col_s, ncol_c);
1208 
1209  // --- Chunk subviews: 2D (ncol, nlay) via LayoutRight pointer offset ---
1210  const int stride2_nlay = nlay;
1211  const int stride2_nlayp1 = nlay + 1;
1212  real2d_k p_lay_c (p_lay.data() + col_s*stride2_nlay, ncol_c, nlay);
1213  real2d_k t_lay_c (t_lay.data() + col_s*stride2_nlay, ncol_c, nlay);
1214  real2d_k r_lay_c (r_lay.data() + col_s*stride2_nlay, ncol_c, nlay);
1215  real2d_k z_del_c (z_del.data() + col_s*stride2_nlay, ncol_c, nlay);
1216  real2d_k lwp_c (lwp.data() + col_s*stride2_nlay, ncol_c, nlay);
1217  real2d_k iwp_c (iwp.data() + col_s*stride2_nlay, ncol_c, nlay);
1218  real2d_k eff_radius_qc_c(eff_radius_qc.data()+ col_s*stride2_nlay, ncol_c, nlay);
1219  real2d_k eff_radius_qi_c(eff_radius_qi.data()+ col_s*stride2_nlay, ncol_c, nlay);
1220  real2d_k cldfrac_tot_c (cldfrac_tot.data() + col_s*stride2_nlay, ncol_c, nlay);
1221  real2d_k sw_heating_c (sw_heating.data() + col_s*stride2_nlay, ncol_c, nlay);
1222  real2d_k lw_heating_c (lw_heating.data() + col_s*stride2_nlay, ncol_c, nlay);
1223 
1224  // --- Chunk subviews: 2D (ncol, nlay+1) ---
1225  real2d_k p_lev_c (p_lev.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1226  real2d_k t_lev_c (t_lev.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1227  real2d_k sw_flux_up_c (sw_flux_up.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1228  real2d_k sw_flux_dn_c (sw_flux_dn.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1229  real2d_k sw_flux_dn_dir_c (sw_flux_dn_dir.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1230  real2d_k lw_flux_up_c (lw_flux_up.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1231  real2d_k lw_flux_dn_c (lw_flux_dn.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1232  // Clear-sky flux subviews (always active)
1233  real2d_k sw_clrsky_flux_up_c (sw_clrsky_flux_up.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1234  real2d_k sw_clrsky_flux_dn_c (sw_clrsky_flux_dn.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1235  real2d_k sw_clrsky_flux_dn_dir_c (sw_clrsky_flux_dn_dir.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1236  real2d_k lw_clrsky_flux_up_c (lw_clrsky_flux_up.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1237  real2d_k lw_clrsky_flux_dn_c (lw_clrsky_flux_dn.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1238 
1239  // Diagnostic flux subviews (placeholder when disabled)
1240  real2d_k sw_clnclrsky_flux_up_c, sw_clnclrsky_flux_dn_c, sw_clnclrsky_flux_dn_dir_c;
1241  real2d_k lw_clnclrsky_flux_up_c, lw_clnclrsky_flux_dn_c;
1242  if (m_extra_clnclrsky_diag) {
1243  sw_clnclrsky_flux_up_c = real2d_k(sw_clnclrsky_flux_up.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1244  sw_clnclrsky_flux_dn_c = real2d_k(sw_clnclrsky_flux_dn.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1245  sw_clnclrsky_flux_dn_dir_c = real2d_k(sw_clnclrsky_flux_dn_dir.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1246  lw_clnclrsky_flux_up_c = real2d_k(lw_clnclrsky_flux_up.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1247  lw_clnclrsky_flux_dn_c = real2d_k(lw_clnclrsky_flux_dn.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1248  } else {
1249  sw_clnclrsky_flux_up_c = real2d_k("sw_clnclrsky_flux_up_c" , 1, 1);
1250  sw_clnclrsky_flux_dn_c = real2d_k("sw_clnclrsky_flux_dn_c" , 1, 1);
1251  sw_clnclrsky_flux_dn_dir_c = real2d_k("sw_clnclrsky_flux_dn_dir_c", 1, 1);
1252  lw_clnclrsky_flux_up_c = real2d_k("lw_clnclrsky_flux_up_c" , 1, 1);
1253  lw_clnclrsky_flux_dn_c = real2d_k("lw_clnclrsky_flux_dn_c" , 1, 1);
1254  }
1255 
1256  real2d_k sw_clnsky_flux_up_c, sw_clnsky_flux_dn_c, sw_clnsky_flux_dn_dir_c;
1257  real2d_k lw_clnsky_flux_up_c, lw_clnsky_flux_dn_c;
1258  if (m_extra_clnsky_diag) {
1259  sw_clnsky_flux_up_c = real2d_k(sw_clnsky_flux_up.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1260  sw_clnsky_flux_dn_c = real2d_k(sw_clnsky_flux_dn.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1261  sw_clnsky_flux_dn_dir_c = real2d_k(sw_clnsky_flux_dn_dir.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1262  lw_clnsky_flux_up_c = real2d_k(lw_clnsky_flux_up.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1263  lw_clnsky_flux_dn_c = real2d_k(lw_clnsky_flux_dn.data() + col_s*stride2_nlayp1, ncol_c, nlay+1);
1264  } else {
1265  sw_clnsky_flux_up_c = real2d_k("sw_clnsky_flux_up_c" , 1, 1);
1266  sw_clnsky_flux_dn_c = real2d_k("sw_clnsky_flux_dn_c" , 1, 1);
1267  sw_clnsky_flux_dn_dir_c = real2d_k("sw_clnsky_flux_dn_dir_c", 1, 1);
1268  lw_clnsky_flux_up_c = real2d_k("lw_clnsky_flux_up_c" , 1, 1);
1269  lw_clnsky_flux_dn_c = real2d_k("lw_clnsky_flux_dn_c" , 1, 1);
1270  }
1271 
1272  // --- Chunk subviews: 2D (ncol, nswbands) ---
1273  real2d_k sfc_alb_dir_c(sfc_alb_dir.data() + col_s*nswbands, ncol_c, nswbands);
1274  real2d_k sfc_alb_dif_c(sfc_alb_dif.data() + col_s*nswbands, ncol_c, nswbands);
1275 
1276  // --- Chunk subviews: 3D (ncol, nlay+1, nbands) ---
1277  const int stride3_sw = stride2_nlayp1 * nswbands;
1278  const int stride3_lw = stride2_nlayp1 * m_nlwbands;
1279  real3d_k sw_bnd_flux_up_c (sw_bnd_flux_up.data() + col_s*stride3_sw, ncol_c, nlay+1, nswbands);
1280  real3d_k sw_bnd_flux_dn_c (sw_bnd_flux_dn.data() + col_s*stride3_sw, ncol_c, nlay+1, nswbands);
1281  real3d_k sw_bnd_flux_dir_c(sw_bnd_flux_dir.data() + col_s*stride3_sw, ncol_c, nlay+1, nswbands);
1282  real3d_k sw_bnd_flux_dif_c(sw_bnd_flux_dif.data() + col_s*stride3_sw, ncol_c, nlay+1, nswbands);
1283  real3d_k lw_bnd_flux_up_c (lw_bnd_flux_up.data() + col_s*stride3_lw, ncol_c, nlay+1, m_nlwbands);
1284  real3d_k lw_bnd_flux_dn_c (lw_bnd_flux_dn.data() + col_s*stride3_lw, ncol_c, nlay+1, m_nlwbands);
1285 
1286  // --- Create chunk gas concentrations by subsetting from pre-fetched VMR ---
1287  gas_concs_t gas_concs_c;
1288  gas_concs_c.init(gas_names_offset, ncol_c, nlay);
1289  for (int igas = 0; igas < m_ngas; ++igas) {
1290  real2d_k vmr_c("vmr_c", ncol_c, nlay);
1291  auto vmr_full = vmr_full_vec[igas];
1292  auto cs = col_s;
1293  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol_c, nlay}),
1294  KOKKOS_LAMBDA (int i, int j) {
1295  vmr_c(i, j) = vmr_full(cs + i, j);
1296  });
1297  gas_concs_c.set_vmr(m_gas_names[igas], vmr_c);
1298  }
1299 
1300  // Expand surface albedos along nswbands for this chunk
1302  sfc_alb_dir_vis_c, sfc_alb_dir_nir_c,
1303  sfc_alb_dif_vis_c, sfc_alb_dif_nir_c,
1304  sfc_alb_dir_c , sfc_alb_dif_c);
1305 
1306  // Run RRTMGP driver for this column chunk
1307  rrtmgp::rrtmgp_main(ncol_c, m_nlay,
1308  p_lay_c, t_lay_c,
1309  p_lev_c, t_lev_c,
1310  gas_concs_c,
1311  sfc_alb_dir_c, sfc_alb_dif_c, mu0_c,
1312  t_sfc_c, sfc_emis_c, lw_src_c,
1313  lwp_c, iwp_c, eff_radius_qc_c, eff_radius_qi_c, cldfrac_tot_c,
1314  sw_flux_up_c, sw_flux_dn_c, sw_flux_dn_dir_c,
1315  lw_flux_up_c, lw_flux_dn_c,
1316  sw_clnclrsky_flux_up_c, sw_clnclrsky_flux_dn_c, sw_clnclrsky_flux_dn_dir_c,
1317  sw_clrsky_flux_up_c, sw_clrsky_flux_dn_c, sw_clrsky_flux_dn_dir_c,
1318  sw_clnsky_flux_up_c, sw_clnsky_flux_dn_c, sw_clnsky_flux_dn_dir_c,
1319  lw_clnclrsky_flux_up_c, lw_clnclrsky_flux_dn_c,
1320  lw_clrsky_flux_up_c, lw_clrsky_flux_dn_c,
1321  lw_clnsky_flux_up_c, lw_clnsky_flux_dn_c,
1322  sw_bnd_flux_up_c, sw_bnd_flux_dn_c, sw_bnd_flux_dir_c,
1323  lw_bnd_flux_up_c, lw_bnd_flux_dn_c,
1325 
1326  // Compute heating rates for this chunk
1327  rrtmgp::compute_heating_rate(sw_flux_up_c, sw_flux_dn_c, r_lay_c, z_del_c, sw_heating_c);
1328  rrtmgp::compute_heating_rate(lw_flux_up_c, lw_flux_dn_c, r_lay_c, z_del_c, lw_heating_c);
1329 
1330  // Compute diffuse band fluxes and broadband surface fluxes for this chunk
1331  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<3>>({0, 0, 0}, {ncol_c, nlay+1, nswbands}),
1332  KOKKOS_LAMBDA (int icol, int ilay, int ibnd)
1333  {
1334  sw_bnd_flux_dif_c(icol,ilay,ibnd) = sw_bnd_flux_dn_c(icol,ilay,ibnd) - sw_bnd_flux_dir_c(icol,ilay,ibnd);
1335  });
1336  rrtmgp::compute_broadband_surface_fluxes(ncol_c, kbot, nswbands,
1337  sw_bnd_flux_dir_c , sw_bnd_flux_dif_c ,
1338  sfc_flux_dir_vis_c, sfc_flux_dir_nir_c,
1339  sfc_flux_dif_vis_c, sfc_flux_dif_nir_c);
1340 
1341  gas_concs_c.reset();
1342  } // end column chunk loop
1343 }
static constexpr int ORB_UNDEF_INT
Definition: ERF_Constants.H:115
constexpr amrex::Real mwdair
Definition: ERF_Constants.H:83
constexpr amrex::Real PI
Definition: ERF_Constants.H:24
AMREX_GPU_HOST AMREX_FORCE_INLINE real orbital_cos_zenith(real &jday, real &lat, real &lon, real &declin, real dt_avg=-one, real uniform_angle=-one, real constant_zenith_angle_deg=-one)
Definition: ERF_OrbCosZenith.H:562
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
amrex::Real m_dt
Definition: ERF_Radiation.H:226
real(c_double), private cs
Definition: ERF_module_mp_morr_two_moment.F90:203
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, real1d_k &sfc_emis, real1d_k &lw_src, real2d_k &lwp, real2d_k &iwp, real2d_k &rel, real2d_k &rei, real2d_k &cldfrac, 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

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::iMultiFab *  lmask,
amrex::MultiFab *  t_surf,
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 
)
174 {
175  // Set data members that may change
176  m_lev = level;
177  m_step = step;
178  m_time = time;
179  m_dt = dt;
180  m_geom = geom;
181  m_cons_in = cons_in;
182  m_qheating_rates = qheating_rates;
183  m_rad_fluxes = rad_fluxes;
184  m_z_phys = z_phys;
185  m_lat = lat;
186  m_lon = lon;
187 
188  // Update the day and month
189  time_t timestamp = time_t(time);
190  struct tm *timeinfo = gmtime(&timestamp);
191  if (m_fixed_orbital_year) {
192  m_orbital_mon = timeinfo->tm_mon + 1;
193  m_orbital_day = timeinfo->tm_mday;
194  m_orbital_sec = timeinfo->tm_hour*3600 + timeinfo->tm_min*60 + timeinfo->tm_sec;
195  } else {
196  m_orbital_year = timeinfo->tm_year + 1900;
197  m_orbital_mon = timeinfo->tm_mon + 1;
198  m_orbital_day = timeinfo->tm_mday;
199  m_orbital_sec = timeinfo->tm_hour*3600 + timeinfo->tm_min*60 + timeinfo->tm_sec;
200  }
201 
202  // Only allocate and proceed if we are going to update radiation
203  m_update_rad = false;
204  if (m_rad_freq_in_steps > 0) { m_update_rad = ( (m_step == 0) || (m_step % m_rad_freq_in_steps == 0) ); }
205 
206  if (m_update_rad) {
207  // Call to Init() has set the dimensions: ncol & nlay
208 
209  // Allocate the buffer arrays
210  alloc_buffers();
211 
212  // Fill the KOKKOS Views from AMReX MFs
213  mf_to_kokkos_buffers(lmask, t_surf, lsm_input_ptrs);
214 
215  // Initialize datalog MF on first step
216  if (m_first_step) {
217  m_first_step = false;
218  datalog_mf.define(cons_in->boxArray(), cons_in->DistributionMap(), 25, 0);
219  datalog_mf.setVal(0.0);
220  }
221  }
222 }
int m_step
Definition: ERF_Radiation.H:220
void mf_to_kokkos_buffers(amrex::iMultiFab *lmask, amrex::MultiFab *t_surf, amrex::Vector< amrex::MultiFab * > &lsm_input_ptrs)
Definition: ERF_Radiation.cpp:461
bool m_first_step
Definition: ERF_Radiation.H:240
void alloc_buffers()
Definition: ERF_Radiation.cpp:225
amrex::Real m_time
Definition: ERF_Radiation.H:223

Referenced by Run().

Here is the caller graph for this function:

◆ write_rrtmgp_fluxes()

void Radiation::write_rrtmgp_fluxes ( )
753 {
754  Table2D<Real,Order::C> sw_flux_up_tab(sw_flux_up.data(), {0,0}, {static_cast<int>(sw_flux_up.extent(0)),static_cast<int>(sw_flux_up.extent(1))});
755  Table2D<Real,Order::C> sw_flux_dn_tab(sw_flux_dn.data(), {0,0}, {static_cast<int>(sw_flux_dn.extent(0)),static_cast<int>(sw_flux_dn.extent(1))});
756  Table2D<Real,Order::C> sw_flux_dn_dir_tab(sw_flux_dn_dir.data(), {0,0}, {static_cast<int>(sw_flux_dn_dir.extent(0)),static_cast<int>(sw_flux_dn_dir.extent(1))});
757  Table2D<Real,Order::C> lw_flux_up_tab(lw_flux_up.data(), {0,0}, {static_cast<int>(lw_flux_up.extent(0)),static_cast<int>(lw_flux_up.extent(1))});
758  Table2D<Real,Order::C> lw_flux_dn_tab(lw_flux_dn.data(), {0,0}, {static_cast<int>(lw_flux_dn.extent(0)),static_cast<int>(lw_flux_dn.extent(1))});
759 
760  int n_fluxes = 5;
761  MultiFab mf_flux(m_cons_in->boxArray(), m_cons_in->DistributionMap(), n_fluxes, 0);
762 
763  for (MFIter mfi(mf_flux); mfi.isValid(); ++mfi) {
764  const auto& vbx = mfi.validbox();
765  const int nx = vbx.length(0);
766  const int imin = vbx.smallEnd(0);
767  const int jmin = vbx.smallEnd(1);
768  const int offset = m_col_offsets[mfi.index()];
769  const Array4<Real>& dst_arr = mf_flux.array(mfi);
770  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
771  {
772  // map [i,j,k] 0-based to [icol, ilay] 0-based
773  const int icol = (j-jmin)*nx + (i-imin) + offset;
774  const int ilay = k;
775 
776  // SW and LW fluxes
777  dst_arr(i,j,k,0) = sw_flux_up_tab(icol,ilay);
778  dst_arr(i,j,k,1) = sw_flux_dn_tab(icol,ilay);
779  dst_arr(i,j,k,2) = sw_flux_dn_dir_tab(icol,ilay);
780  dst_arr(i,j,k,3) = lw_flux_up_tab(icol,ilay);
781  dst_arr(i,j,k,4) = lw_flux_dn_tab(icol,ilay);
782  });
783  }
784 
785 
786  std::string plotfilename = amrex::Concatenate("plt_rad", m_step, 5);
787  Vector<std::string> flux_names = {"sw_flux_up", "sw_flux_dn", "sw_flux_dir",
788  "lw_flux_up", "lw_flux_dn"};
789  WriteSingleLevelPlotfile(plotfilename, mf_flux, flux_names, m_geom, m_time, m_step);
790 }
Here is the call graph for this function:

◆ WriteDataLog()

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

Implements IRadiation.

899 {
900  constexpr int datwidth = 14;
901  constexpr int datprecision = 9;
902  constexpr int timeprecision = 13;
903 
904  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;
905  // Clear sky
906  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;
907  // Clean sky
908  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;
909  // Clean clear sky
910  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;
911 
912 
913  auto domain = m_geom.Domain();
914  h_avg_radqrsw = sumToLine(datalog_mf, 0, 1, domain, 2);
915  h_avg_radqrlw = sumToLine(datalog_mf, 1, 1, domain, 2);
916  h_avg_sw_up = sumToLine(datalog_mf, 2, 1, domain, 2);
917  h_avg_sw_dn = sumToLine(datalog_mf, 3, 1, domain, 2);
918  h_avg_sw_dn_dir = sumToLine(datalog_mf, 4, 1, domain, 2);
919  h_avg_lw_up = sumToLine(datalog_mf, 5, 1, domain, 2);
920  h_avg_lw_dn = sumToLine(datalog_mf, 6, 1, domain, 2);
921  h_avg_zenith = sumToLine(datalog_mf, 7, 1, domain, 2);
922 
923  h_avg_radqrcsw = sumToLine(datalog_mf, 8, 1, domain, 2);
924  h_avg_radqrclw = sumToLine(datalog_mf, 9, 1, domain, 2);
925  h_avg_sw_clr_up = sumToLine(datalog_mf, 10, 1, domain, 2);
926  h_avg_sw_clr_dn = sumToLine(datalog_mf, 11, 1, domain, 2);
927  h_avg_sw_clr_dn_dir = sumToLine(datalog_mf, 12, 1, domain, 2);
928  h_avg_lw_clr_up = sumToLine(datalog_mf, 13, 1, domain, 2);
929  h_avg_lw_clr_dn = sumToLine(datalog_mf, 14, 1, domain, 2);
930 
931  if (m_extra_clnsky_diag) {
932  h_avg_sw_cln_up = sumToLine(datalog_mf, 15, 1, domain, 2);
933  h_avg_sw_cln_dn = sumToLine(datalog_mf, 16, 1, domain, 2);
934  h_avg_sw_cln_dn_dir = sumToLine(datalog_mf, 17, 1, domain, 2);
935  h_avg_lw_cln_up = sumToLine(datalog_mf, 18, 1, domain, 2);
936  h_avg_lw_cln_dn = sumToLine(datalog_mf, 19, 1, domain, 2);
937  }
938 
940  h_avg_sw_clnclr_up = sumToLine(datalog_mf, 20, 1, domain, 2);
941  h_avg_sw_clnclr_dn = sumToLine(datalog_mf, 21, 1, domain, 2);
942  h_avg_sw_clnclr_dn_dir = sumToLine(datalog_mf, 22, 1, domain, 2);
943  h_avg_lw_clnclr_up = sumToLine(datalog_mf, 23, 1, domain, 2);
944  h_avg_lw_clnclr_dn = sumToLine(datalog_mf, 24, 1, domain, 2);
945  }
946 
947  Real area_z = static_cast<Real>(domain.length(0)*domain.length(1));
948  int nz = domain.length(2);
949  for (int k = 0; k < nz; k++) {
950  h_avg_radqrsw[k] /= area_z;
951  h_avg_radqrlw[k] /= area_z;
952  h_avg_sw_up[k] /= area_z;
953  h_avg_sw_dn[k] /= area_z;
954  h_avg_sw_dn_dir[k] /= area_z;
955  h_avg_lw_up[k] /= area_z;
956  h_avg_lw_dn[k] /= area_z;
957  h_avg_zenith[k] /= area_z;
958 
959  h_avg_radqrcsw[k] /= area_z;
960  h_avg_radqrclw[k] /= area_z;
961  h_avg_sw_clr_up[k] /= area_z;
962  h_avg_sw_clr_dn[k] /= area_z;
963  h_avg_sw_clr_dn_dir[k] /= area_z;
964  h_avg_lw_clr_up[k] /= area_z;
965  h_avg_lw_clr_dn[k] /= area_z;
966  }
967 
968  if (m_extra_clnsky_diag) {
969  for (int k = 0; k < nz; k++) {
970  h_avg_sw_cln_up[k] /= area_z;
971  h_avg_sw_cln_dn[k] /= area_z;
972  h_avg_sw_cln_dn_dir[k] /= area_z;
973  h_avg_lw_cln_up[k] /= area_z;
974  h_avg_lw_cln_dn[k] /= area_z;
975  }
976  }
977 
979  for (int k = 0; k < nz; k++) {
980  h_avg_sw_clnclr_up[k] /= area_z;
981  h_avg_sw_clnclr_dn[k] /= area_z;
982  h_avg_sw_clnclr_dn_dir[k] /= area_z;
983  h_avg_lw_clnclr_up[k] /= area_z;
984  h_avg_lw_clnclr_dn[k] /= area_z;
985  }
986  }
987 
988  if (ParallelDescriptor::IOProcessor()) {
989  std::ostream& log = *datalog;
990  if (log.good()) {
991 
992  for (int k = 0; k < nz; k++)
993  {
994  Real z = k * m_geom.CellSize(2);
995  log << std::setw(datwidth) << std::setprecision(timeprecision) << time << " "
996  << std::setw(datwidth) << std::setprecision(datprecision) << z << " "
997  << h_avg_radqrsw[k] << " " << h_avg_radqrlw[k] << " " << h_avg_sw_up[k] << " "
998  << h_avg_sw_dn[k] << " " << h_avg_sw_dn_dir[k] << " " << h_avg_lw_up[k] << " "
999  << h_avg_lw_dn[k] << " " << h_avg_zenith[k] << " "
1000  << h_avg_radqrcsw[k] << " " << h_avg_radqrclw[k] << " " << h_avg_sw_clr_up[k] << " "
1001  << h_avg_sw_clr_dn[k] << " " << h_avg_sw_clr_dn_dir[k] << " " << h_avg_lw_clr_up[k] << " "
1002  << h_avg_lw_clr_dn[k] << " ";
1003  if (m_extra_clnsky_diag) {
1004  log << h_avg_sw_cln_up[k] << " " << h_avg_sw_cln_dn[k] << " " << h_avg_sw_cln_dn_dir[k] << " "
1005  << h_avg_lw_cln_up[k] << " " << h_avg_lw_cln_dn[k] << " ";
1006  } else {
1007  log << zero << " " << zero << " " << zero << " " << zero << " " << zero << " ";
1008  }
1009 
1010  if (m_extra_clnclrsky_diag) {
1011  log << h_avg_sw_clnclr_up[k] << " " << h_avg_sw_clnclr_dn[k] << " " << h_avg_sw_clnclr_dn_dir[k] << " "
1012  << h_avg_lw_clnclr_up[k] << " " << h_avg_lw_clnclr_dn[k] << std::endl;
1013  } else {
1014  log << zero << " " << zero << " " << zero << " " << zero << " " << zero << std::endl;
1015  }
1016  }
1017  // Write top face values
1018  Real z = nz * m_geom.CellSize(2);
1019  log << std::setw(datwidth) << std::setprecision(timeprecision) << time << " "
1020  << std::setw(datwidth) << std::setprecision(datprecision) << z << " "
1021  << zero << " " << zero << " " << zero << " " << zero << " " << zero << " " << zero << " "
1022  << zero << " " << zero << " "
1023  << zero << " " << zero << " " << zero << " " << zero << " " << zero << " " << zero << " "
1024  << zero << " "
1025  << zero << " " << zero << " " << zero << " " << zero << " " << zero << " "
1026  << zero << " " << zero << " " << zero << " " << zero << " " << zero
1027  << std::endl;
1028  }
1029  }
1030 }
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

◆ 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

◆ 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 = amrex::Real(1807.851e-9)
private

◆ m_co2vmr

amrex::Real Radiation::m_co2vmr = amrex::Real(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 = amrex::Real(1.0e-7)
private

◆ m_do_aerosol_rad

bool Radiation::m_do_aerosol_rad = false
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 = -amrex::Real(9999.)
private

◆ m_fixed_total_solar_irradiance

amrex::Real Radiation::m_fixed_total_solar_irradiance = -amrex::Real(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 = amrex::Real(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 = -amrex::Real(98.555183)
private

◆ m_lsm

bool Radiation::m_lsm = false
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
private
Initial value:
= {"cos_zenith_angle" , "sw_flux_dn" ,
"sw_flux_dn_dir_vis", "sw_flux_dn_dir_nir",
"sw_flux_dn_dif_vis", "sw_flux_dn_dif_nir",
"lw_flux_dn"}

Referenced by get_lsm_output_varnames().

◆ 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:
= {amrex::Real(18.01528), amrex::Real(44.00950), amrex::Real(47.9982), amrex::Real(44.0128),
amrex::Real(28.01010), amrex::Real(16.04246), amrex::Real(31.9980), amrex::Real(28.0134)}

◆ m_n2ovmr

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

◆ m_n2vmr

amrex::Real Radiation::m_n2vmr = amrex::Real(0.7906)
private

◆ m_ncol

int Radiation::m_ncol
private

Referenced by Init().

◆ m_ncol_chunk

int Radiation::m_ncol_chunk = 1024
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 = amrex::Real(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 = -amrex::Real(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 = -amrex::Real(9999.)
private

◆ m_orbital_obliq

amrex::Real Radiation::m_orbital_obliq = -amrex::Real(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_nvar

int Radiation::m_rad_nvar = 12
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: