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::MultiFab *lsm_fluxes, amrex::MultiFab *lsm_zenith, amrex::Vector< amrex::MultiFab * > &lsm_input_ptrs, amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs, amrex::MultiFab *qheating_rates, amrex::MultiFab *rad_fluxes, amrex::MultiFab *z_phys, amrex::MultiFab *lat_ptr, amrex::MultiFab *lon_ptr) override
 
void set_grids (int &level, int &step, amrex::Real &time, const amrex::Real &dt, const amrex::BoxArray &ba, amrex::Geometry &geom, amrex::MultiFab *cons_in, amrex::iMultiFab *lmask, amrex::MultiFab *t_surf, amrex::MultiFab *lsm_fluxes, amrex::MultiFab *lsm_zenith, amrex::Vector< amrex::MultiFab * > &lsm_input_ptrs, amrex::MultiFab *qheating_rates, amrex::MultiFab *rad_fluxes, amrex::MultiFab *z_phys, amrex::MultiFab *lat, amrex::MultiFab *lon)
 
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 = {"cos_zenith_angle", "sw_flux_dn", "lw_flux_dn"}
 
amrex::Real m_rad_t_sfc = -1
 
amrex::MultiFab * m_cons_in = nullptr
 
amrex::MultiFab * m_qheating_rates = nullptr
 
amrex::MultiFab * m_rad_fluxes = nullptr
 
amrex::MultiFab * m_z_phys = nullptr
 
amrex::MultiFab * m_lat = nullptr
 
amrex::MultiFab * m_lon = nullptr
 
amrex::Real m_lat_cons = 39.809860
 
amrex::Real m_lon_cons = -98.555183
 
amrex::MultiFab * m_lsm_fluxes = nullptr
 
amrex::MultiFab * m_lsm_zenith = nullptr
 
amrex::MultiFab datalog_mf
 
std::string rrtmgp_file_path = "."
 
std::string rrtmgp_coeffs_sw = "rrtmgp-data-sw-g224-2018-12-04.nc"
 
std::string rrtmgp_coeffs_lw = "rrtmgp-data-lw-g256-2018-12-04.nc"
 
std::string rrtmgp_cloud_optics_sw = "rrtmgp-cloud-optics-coeffs-sw.nc"
 
std::string rrtmgp_cloud_optics_lw = "rrtmgp-cloud-optics-coeffs-lw.nc"
 
std::string rrtmgp_coeffs_file_sw
 
std::string rrtmgp_coeffs_file_lw
 
std::string rrtmgp_cloud_optics_file_sw
 
std::string rrtmgp_cloud_optics_file_lw
 
int m_ngas = 8
 
const std::vector< std::string > m_gas_names
 
const std::vector< amrex::Realm_mol_weight_gas
 
amrex::Real m_co2vmr = 388.717e-6
 
amrex::Vector< amrex::Realm_o3vmr
 
amrex::Real m_n2ovmr = 323.141e-9
 
amrex::Real m_covmr = 1.0e-7
 
amrex::Real m_ch4vmr = 1807.851e-9
 
amrex::Real m_o2vmr = 0.209448
 
amrex::Real m_n2vmr = 0.7906
 
int m_o3_size
 
real1d_k m_gas_mol_weights
 
std::vector< std::string > gas_names_offset
 
GasConcsK< amrex::Real, layout_t, KokkosDefaultDevicem_gas_concs
 
int m_ncol
 
int m_nlay
 
amrex::Vector< int > m_col_offsets
 
bool m_do_aerosol_rad = true
 
bool m_extra_clnsky_diag = false
 
bool m_extra_clnclrsky_diag = false
 
int m_orbital_year = -9999
 
int m_orbital_mon = -9999
 
int m_orbital_day = -9999
 
int m_orbital_sec = -9999
 
bool m_fixed_orbital_year = false
 
amrex::Real m_orbital_eccen = -9999.
 
amrex::Real m_orbital_obliq = -9999.
 
amrex::Real m_orbital_mvelp = -9999.
 
amrex::Real m_fixed_total_solar_irradiance = -9999.
 
amrex::Real m_fixed_solar_zenith_angle = -9999.
 
int m_nswbands
 
int m_nlwbands
 
int m_nswgpts
 
int m_nlwgpts
 
int m_rad_freq_in_steps = 1
 
bool m_do_subcol_sampling = true
 
real1d_k o3_lay
 
real1d_k mu0
 
real1d_k sfc_alb_dir_vis
 
real1d_k sfc_alb_dir_nir
 
real1d_k sfc_alb_dif_vis
 
real1d_k sfc_alb_dif_nir
 
real1d_k sfc_flux_dir_vis
 
real1d_k sfc_flux_dir_nir
 
real1d_k sfc_flux_dif_vis
 
real1d_k sfc_flux_dif_nir
 
real1d_k lat
 
real1d_k lon
 
real1d_k sfc_emis
 
real1d_k t_sfc
 
real1d_k lw_src
 
real2d_k r_lay
 
real2d_k p_lay
 
real2d_k t_lay
 
real2d_k z_del
 
real2d_k qv_lay
 
real2d_k qc_lay
 
real2d_k qi_lay
 
real2d_k cldfrac_tot
 
real2d_k eff_radius_qc
 
real2d_k eff_radius_qi
 
real2d_k lwp
 
real2d_k iwp
 
real2d_k sw_heating
 
real2d_k lw_heating
 
real2d_k sw_clrsky_heating
 
real2d_k lw_clrsky_heating
 
real2d_k d_tint
 
real2d_k p_lev
 
real2d_k t_lev
 
real2d_k sw_flux_up
 
real2d_k sw_flux_dn
 
real2d_k sw_flux_dn_dir
 
real2d_k lw_flux_up
 
real2d_k lw_flux_dn
 
real2d_k sw_clnclrsky_flux_up
 
real2d_k sw_clnclrsky_flux_dn
 
real2d_k sw_clnclrsky_flux_dn_dir
 
real2d_k sw_clrsky_flux_up
 
real2d_k sw_clrsky_flux_dn
 
real2d_k sw_clrsky_flux_dn_dir
 
real2d_k sw_clnsky_flux_up
 
real2d_k sw_clnsky_flux_dn
 
real2d_k sw_clnsky_flux_dn_dir
 
real2d_k lw_clnclrsky_flux_up
 
real2d_k lw_clnclrsky_flux_dn
 
real2d_k lw_clrsky_flux_up
 
real2d_k lw_clrsky_flux_dn
 
real2d_k lw_clnsky_flux_up
 
real2d_k lw_clnsky_flux_dn
 
real3d_k sw_bnd_flux_up
 
real3d_k sw_bnd_flux_dn
 
real3d_k sw_bnd_flux_dir
 
real3d_k sw_bnd_flux_dif
 
real3d_k lw_bnd_flux_up
 
real3d_k lw_bnd_flux_dn
 
real2d_k sfc_alb_dir
 
real2d_k sfc_alb_dif
 
real2d_k emis_sfc
 
real3d_k aero_tau_sw
 
real3d_k aero_ssa_sw
 
real3d_k aero_g_sw
 
real3d_k aero_tau_lw
 
real3d_k cld_tau_sw_bnd
 
real3d_k cld_tau_lw_bnd
 
real3d_k cld_tau_sw_gpt
 
real3d_k cld_tau_lw_gpt
 

Additional Inherited Members

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

Constructor & Destructor Documentation

◆ Radiation()

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

◆ ~Radiation()

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

Member Function Documentation

◆ alloc_buffers()

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

◆ dealloc_buffers()

void Radiation::dealloc_buffers ( )
335 {
336  // 1d size (m_ngas)
338 
339  // 1d size (1 or nlay)
340  o3_lay = real1d_k();
341 
342  // 1d size (ncol)
343  mu0 = real1d_k();
352  lat = real1d_k();
353  lon = real1d_k();
354  sfc_emis = real1d_k();
355  t_sfc = real1d_k();
356  lw_src = real1d_k();
357 
358  // 2d size (ncol, nlay)
359  r_lay = real2d_k();
360  p_lay = real2d_k();
361  t_lay = real2d_k();
362  z_del = real2d_k();
363  qv_lay = real2d_k();
364  qc_lay = real2d_k();
365  qi_lay = real2d_k();
366  cldfrac_tot = real2d_k();
369  lwp = real2d_k();
370  iwp = real2d_k();
371  sw_heating = real2d_k();
372  lw_heating = real2d_k();
375  sw_heating = real2d_k();
376  lw_heating = real2d_k();
379 
380  // 2d size (ncol, nlay+1)
381  d_tint = real2d_k();
382  p_lev = real2d_k();
383  t_lev = real2d_k();
384  sw_flux_up = real2d_k();
385  sw_flux_dn = real2d_k();
387  lw_flux_up = real2d_k();
388  lw_flux_dn = real2d_k();
404 
405  // 3d size (ncol, nlay+1, nswbands)
410 
411  // 3d size (ncol, nlay+1, nlwbands)
414 
415  // 2d size (ncol, nswbands)
416  sfc_alb_dir = real2d_k();
417  sfc_alb_dif = real2d_k();
418 
419  // 2d size (ncol, nlwbands)
420  //emis_sfc = real2d_k();
421 
422  /*
423  // 3d size (ncol, nlay, n[sw,lw]bands)
424  aero_tau_sw = real3d_k();
425  aero_ssa_sw = real3d_k();
426  aero_g_sw = real3d_k();
427  aero_tau_lw = real3d_k();
428 
429  // 3d size (ncol, nlay, n[sw,lw]bnds)
430  cld_tau_sw_bnd = real3d_k();
431  cld_tau_lw_bnd = real3d_k();
432 
433  // 3d size (ncol, nlay, n[sw,lw]gpts)
434  cld_tau_sw_gpt = real3d_k();
435  cld_tau_lw_gpt = real3d_k();
436  */
437 }

◆ finalize_impl()

void Radiation::finalize_impl ( amrex::Vector< amrex::MultiFab * > &  lsm_output_ptrs)
1302 {
1303  // Finish rrtmgp
1304  m_gas_concs.reset();
1306 
1307  // Fill the AMReX MFs from Kokkos Views
1308  kokkos_buffers_to_mf(lsm_output_ptrs);
1309 
1310  // Write fluxes if requested
1312 
1313  // Fill output data for datalog before deallocating
1314  if (datalog_int > 0) {
1317 
1319  }
1320 
1321  // Deallocate the buffer arrays
1322  dealloc_buffers();
1323 }
int datalog_int
Definition: ERF_RadiationInterface.H:90
void dealloc_buffers()
Definition: ERF_Radiation.cpp:334
void populateDatalogMF()
Definition: ERF_Radiation.cpp:790
void kokkos_buffers_to_mf(amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
Definition: ERF_Radiation.cpp:638
GasConcsK< amrex::Real, layout_t, KokkosDefaultDevice > m_gas_concs
Definition: ERF_Radiation.H:318
void write_rrtmgp_fluxes()
Definition: ERF_Radiation.cpp:750
void rrtmgp_finalize()
Definition: ERF_RRTMGP_Interface.cpp:266
void compute_heating_rate(View1 const &flux_up, View2 const &flux_dn, View3 const &rho, View4 const &dz, View5 &heating_rate)
Definition: ERF_RRTMGP_Utils.H:81

Referenced by rad_run_impl().

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

◆ get_lsm_input_varnames()

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

Reimplemented from IRadiation.

188  {
189  return m_lsm_input_names;
190  }
amrex::Vector< std::string > m_lsm_input_names
Definition: ERF_Radiation.H:248

◆ get_lsm_output_varnames()

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

Reimplemented from IRadiation.

196  {
197  return m_lsm_output_names;
198  }
amrex::Vector< std::string > m_lsm_output_names
Definition: ERF_Radiation.H:253

◆ Init()

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

Implements IRadiation.

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

◆ initialize_impl()

void Radiation::initialize_impl ( )
1032 {
1033  // Call API to initialize
1038 }
void rrtmgp_initialize(gas_concs_t &gas_concs_k, const std::string &coefficients_file_sw, const std::string &coefficients_file_lw, const std::string &cloud_optics_file_sw, const std::string &cloud_optics_file_lw)
Definition: ERF_RRTMGP_Interface.cpp:229

Referenced by rad_run_impl().

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

◆ kokkos_buffers_to_mf()

void Radiation::kokkos_buffers_to_mf ( amrex::Vector< amrex::MultiFab * > &  lsm_output_ptrs)
639 {
640  // Heating rate, fluxes, zenith, lsm ptrs
641  Vector<real2d_k> rrtmgp_out_vars = {sw_flux_dn, lw_flux_dn};
642 
643  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))});
644  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))});
645  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))});
646  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))});
647  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))});
648  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))});
649  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))});
650 
651  Table1D<Real> sfc_flux_dir_vis_tab(sfc_flux_dir_vis.data(), {0}, {static_cast<int>(sfc_flux_dir_vis.extent(0))});
652  Table1D<Real> sfc_flux_dir_nir_tab(sfc_flux_dir_nir.data(), {0}, {static_cast<int>(sfc_flux_dir_nir.extent(0))});
653  Table1D<Real> sfc_flux_dif_vis_tab(sfc_flux_dif_vis.data(), {0}, {static_cast<int>(sfc_flux_dif_vis.extent(0))});
654  Table1D<Real> sfc_flux_dif_nir_tab(sfc_flux_dif_nir.data(), {0}, {static_cast<int>(sfc_flux_dif_nir.extent(0))});
655  Table1D<Real> mu0_tab(mu0.data(), {0}, {static_cast<int>(mu0.extent(0))});
656 
657  // Expose for device
658 
659  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
660  const auto& vbx = mfi.validbox();
661  const auto& sbx = makeSlab(vbx,2,vbx.smallEnd(2));
662  const int nx = vbx.length(0);
663  const int imin = vbx.smallEnd(0);
664  const int jmin = vbx.smallEnd(1);
665  const int offset = m_col_offsets[mfi.index()];
666  const Array4<Real>& q_arr = m_qheating_rates->array(mfi);
667  const Array4<Real>& f_arr = m_rad_fluxes->array(mfi);
668  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
669  {
670  // map [i,j,k] 0-based to [icol, ilay] 0-based
671  const int icol = (j-jmin)*nx + (i-imin) + offset;
672  const int ilay = k;
673 
674  // Temperature heating rate for SW and LW
675  q_arr(i,j,k,0) = sw_heating_tab(icol,ilay);
676  q_arr(i,j,k,1) = lw_heating_tab(icol,ilay);
677 
678  // Convert the dT/dz to dTheta/dz
679  Real iexner = 1./getExnergivenP(Real(p_lay_tab(icol,ilay)), R_d/Cp_d);
680  q_arr(i,j,k,0) *= iexner;
681  q_arr(i,j,k,1) *= iexner;
682 
683  // Populate the fluxes
684  f_arr(i,j,k,0) = sw_flux_up_tab(icol,ilay);
685  f_arr(i,j,k,1) = sw_flux_dn_tab(icol,ilay);
686  f_arr(i,j,k,2) = lw_flux_up_tab(icol,ilay);
687  f_arr(i,j,k,3) = lw_flux_dn_tab(icol,ilay);
688  });
689  if (m_lsm_fluxes) {
690  const Array4<Real>& lsm_arr = m_lsm_fluxes->array(mfi);
691  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
692  {
693  // map [i,j,k] 0-based to [icol, ilay] 0-based
694  const int icol = (j-jmin)*nx + (i-imin) + offset;
695 
696  // SW fluxes for LSM
697  lsm_arr(i,j,k,0) = sfc_flux_dir_vis_tab(icol);
698  lsm_arr(i,j,k,1) = sfc_flux_dir_nir_tab(icol);
699  lsm_arr(i,j,k,2) = sfc_flux_dif_vis_tab(icol);
700  lsm_arr(i,j,k,3) = sfc_flux_dif_nir_tab(icol);
701 
702  // Net SW flux for LSM
703  lsm_arr(i,j,k,4) = sfc_flux_dir_vis_tab(icol) + sfc_flux_dir_nir_tab(icol)
704  + sfc_flux_dif_vis_tab(icol) + sfc_flux_dif_nir_tab(icol);
705 
706  // LW flux for LSM (at bottom surface)
707  lsm_arr(i,j,k,5) = lw_flux_dn_tab(icol,0);
708  });
709  }
710  if (m_lsm_zenith) {
711  const Array4<Real>& lsm_zenith_arr = m_lsm_zenith->array(mfi);
712  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
713  {
714  // map [i,j,k] 0-based to [icol, ilay] 0-based
715  const int icol = (j-jmin)*nx + (i-imin) + offset;
716 
717  // export cosine zenith angle for LSM
718  lsm_zenith_arr(i,j,k) = mu0_tab(icol);
719  });
720  }
721  for (int ivar(0); ivar<lsm_output_ptrs.size(); ivar++) {
722  if (lsm_output_ptrs[ivar]) {
723  const Array4<Real>& lsm_out_arr = lsm_output_ptrs[ivar]->array(mfi);
724  if (ivar==0) {
725  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
726  {
727  // map [i,j,k] 0-based to [icol, ilay] 0-based
728  const int icol = (j-jmin)*nx + (i-imin) + offset;
729 
730  // export the desired variable at surface
731  lsm_out_arr(i,j,k) = mu0_tab(icol);
732  });
733  } else {
734  auto rrtmgp_for_fill = rrtmgp_out_vars[ivar-1];
735  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
736  {
737  // map [i,j,k] 0-based to [icol, ilay] 0-based
738  const int icol = (j-jmin)*nx + (i-imin) + offset;
739 
740  // export the desired variable at surface
741  lsm_out_arr(i,j,k) = rrtmgp_for_fill(icol,0);
742  });
743  } // ivar
744  } // valid ptr
745  } // ivar
746  }// mfi
747 }
constexpr amrex::Real Cp_d
Definition: ERF_Constants.H:12
constexpr amrex::Real R_d
Definition: ERF_Constants.H:10
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getExnergivenP(const amrex::Real P, const amrex::Real rdOcp)
Definition: ERF_EOS.H:141
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
amrex::MultiFab * m_cons_in
Definition: ERF_Radiation.H:259
amrex::MultiFab * m_rad_fluxes
Definition: ERF_Radiation.H:265
amrex::MultiFab * m_lsm_fluxes
Definition: ERF_Radiation.H:279
amrex::MultiFab * m_qheating_rates
Definition: ERF_Radiation.H:262
amrex::MultiFab * m_lsm_zenith
Definition: ERF_Radiation.H:280
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 
)
444 {
445  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))});
446  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))});
447  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))});
448  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))});
449  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))});
450  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))});
451  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))});
452  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))});
453 
454  Table2D<Real,Order::C> lwp_tab(lwp.data(), {0,0}, {static_cast<int>(lwp.extent(0)),static_cast<int>(lwp.extent(1))});
455  Table2D<Real,Order::C> iwp_tab(iwp.data(), {0,0}, {static_cast<int>(iwp.extent(0)),static_cast<int>(iwp.extent(1))});
456  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))});
457  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))});
458 
459  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))});
460  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))});
461 
462  Table1D<Real> lat_tab(lat.data(), {0}, {static_cast<int>(lat.extent(0))});
463  Table1D<Real> lon_tab(lon.data(), {0}, {static_cast<int>(lon.extent(0))});
464  Table1D<Real> t_sfc_tab(t_sfc.data(), {0}, {static_cast<int>(t_sfc.extent(0))});
465 
466  bool moist = m_moist;
467  bool ice = m_ice;
468  const bool has_lsm = m_lsm;
469  const bool has_lat = m_lat;
470  const bool has_lon = m_lon;
471  const bool has_surflayer = (t_surf);
472  int ncol = m_ncol;
473  int nlay = m_nlay;
474  Real dz = m_geom.CellSize(2);
475  Real cons_lat = m_lat_cons;
476  Real cons_lon = m_lon_cons;
477  Real rad_t_sfc = m_rad_t_sfc;
478 
479  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
480  const auto& vbx = mfi.validbox();
481  const int nx = vbx.length(0);
482  const int imin = vbx.smallEnd(0);
483  const int jmin = vbx.smallEnd(1);
484  const int offset = m_col_offsets[mfi.index()];
485  const Array4<const Real>& cons_arr = m_cons_in->const_array(mfi);
486  const Array4<const Real>& z_arr = (m_z_phys) ? m_z_phys->const_array(mfi) :
487  Array4<const Real>{};
488  const Array4<const Real>& lat_arr = (m_lat) ? m_lat->const_array(mfi) :
489  Array4<const Real>{};
490  const Array4<const Real>& lon_arr = (m_lon) ? m_lon->const_array(mfi) :
491  Array4<const Real>{};
492  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
493  {
494  // map [i,j,k] 0-based to [icol, ilay] 0-based
495  const int icol = (j-jmin)*nx + (i-imin) + offset;
496  const int ilay = k;
497 
498  // EOS input (at CC)
499  Real r = cons_arr(i,j,k,Rho_comp);
500  Real rt = cons_arr(i,j,k,RhoTheta_comp);
501  Real qv = (moist) ? std::max(cons_arr(i,j,k,RhoQ1_comp)/r,0.0) : 0.0;
502  Real qc = (moist) ? std::max(cons_arr(i,j,k,RhoQ2_comp)/r,0.0) : 0.0;
503  Real qi = (ice) ? std::max(cons_arr(i,j,k,RhoQ3_comp)/r,0.0) : 0.0;
504 
505  // EOS avg to z-face
506  Real r_lo = cons_arr(i,j,k-1,Rho_comp);
507  Real rt_lo = cons_arr(i,j,k-1,RhoTheta_comp);
508  Real qv_lo = (moist) ? cons_arr(i,j,k-1,RhoQ1_comp)/r_lo : 0.0;
509  Real r_avg = 0.5 * (r + r_lo);
510  Real rt_avg = 0.5 * (rt + rt_lo);
511  Real qv_avg = 0.5 * (qv + qv_lo);
512 
513  // Views at CC
514  r_lay_tab(icol,ilay) = r;
515 
516  p_lay_tab(icol,ilay) = getPgivenRTh(rt, qv);
517  t_lay_tab(icol,ilay) = getTgivenRandRTh(r, rt, qv);
518  z_del_tab(icol,ilay) = (z_arr) ? 0.25 * ( (z_arr(i ,j ,k+1) - z_arr(i ,j ,k))
519  + (z_arr(i+1,j ,k+1) - z_arr(i+1,j ,k))
520  + (z_arr(i ,j+1,k+1) - z_arr(i ,j+1,k))
521  + (z_arr(i+1,j+1,k+1) - z_arr(i+1,j+1,k)) ) : dz;
522  qv_lay_tab(icol,ilay) = qv;
523  qc_lay_tab(icol,ilay) = qc;
524  qi_lay_tab(icol,ilay) = qi;
525  cldfrac_tot_tab(icol,ilay) = ((qc+qi)>0.0) ? 1. : 0.;
526 
527  // NOTE: These are populated in 'mixing_ratio_to_cloud_mass'
528  lwp_tab(icol,ilay) = 0.0;
529  iwp_tab(icol,ilay) = 0.0;
530 
531  // NOTE: These would be populated from P3 (we use the constants in p3_main_impl.hpp)
532  eff_radius_qc_tab(icol,ilay) = (qc>0.0) ? 10.0e-6 : 0.0;
533  eff_radius_qi_tab(icol,ilay) = (qi>0.0) ? 25.0e-6 : 0.0;
534 
535  // Buffers on z-faces (nlay+1)
536  p_lev_tab(icol,ilay) = getPgivenRTh(rt_avg, qv_avg);
537  t_lev_tab(icol,ilay) = getTgivenRandRTh(r_avg, rt_avg, qv_avg);
538  if (ilay==(nlay-1)) {
539  Real r_hi = cons_arr(i,j,k+1,Rho_comp);
540  Real rt_hi = cons_arr(i,j,k+1,RhoTheta_comp);
541  Real qv_hi = (moist) ? std::max(cons_arr(i,j,k+1,RhoQ1_comp)/r_hi,0.0) : 0.0;
542  r_avg = 0.5 * (r + r_hi);
543  rt_avg = 0.5 * (rt + rt_hi);
544  qv_avg = 0.5 * (qv + qv_hi);
545  p_lev_tab(icol,ilay+1) = getPgivenRTh(rt_avg, qv_avg);
546  t_lev_tab(icol,ilay+1) = getTgivenRandRTh(r_avg, rt_avg, qv_avg);
547  }
548 
549  // 1D data structures
550  if (k==0) {
551  lat_tab(icol) = (has_lat) ? lat_arr(i,j,0) : cons_lat;
552  lon_tab(icol) = (has_lon) ? lon_arr(i,j,0) : cons_lon;
553  }
554 
555  });
556  } // mfi
557 
558  // Populate vars LSM would provide
559  if (!has_lsm && !has_surflayer) {
560  // Parsed surface temp
561  Kokkos::deep_copy(t_sfc, rad_t_sfc);
562 
563  // EAMXX dummy atmos constants
564  Kokkos::deep_copy(sfc_alb_dir_vis, 0.06);
565  Kokkos::deep_copy(sfc_alb_dir_nir, 0.06);
566  Kokkos::deep_copy(sfc_alb_dif_vis, 0.06);
567  Kokkos::deep_copy(sfc_alb_dif_nir, 0.06);
568 
569  // AML NOTE: These are not used in current EAMXX, I've left
570  // the code to plug into these if we need it.
571  //
572  // Current EAMXX constants
573  Kokkos::deep_copy(sfc_emis, 0.98);
574  Kokkos::deep_copy(lw_src , 0.0 );
575  } else {
576  Vector<real1d_k> rrtmgp_in_vars = {t_sfc, sfc_emis,
579  Vector<Real> rrtmgp_default_vals = {rad_t_sfc, 0.98,
580  0.06, 0.06,
581  0.06, 0.06};
582  for (int ivar(0); ivar<lsm_input_ptrs.size(); ivar++) {
583  auto rrtmgp_default_val = rrtmgp_default_vals[ivar];
584  auto rrtmgp_to_fill = rrtmgp_in_vars[ivar];
585  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
586  const auto& vbx = mfi.validbox();
587  const auto& sbx = makeSlab(vbx,2,vbx.smallEnd(2));
588  const int nx = vbx.length(0);
589  const int imin = vbx.smallEnd(0);
590  const int jmin = vbx.smallEnd(1);
591  const int offset = m_col_offsets[mfi.index()];
592  const Array4<const int>& lmask_arr = (lmask) ? lmask->const_array(mfi) :
593  Array4<const int> {};
594  const Array4<const Real>& tsurf_arr = (t_surf) ? t_surf->const_array(mfi) :
595  Array4<const Real> {};
596  const Array4<const Real>& lsm_in_arr = (lsm_input_ptrs[ivar]) ? lsm_input_ptrs[ivar]->const_array(mfi) :
597  Array4<const Real> {};
598  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
599  {
600  // map [i,j,k] 0-based to [icol, ilay] 0-based
601  const int icol = (j-jmin)*nx + (i-imin) + offset;
602 
603  // Check if over land
604  bool is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
605 
606  // Check if valid LSM data
607  bool valid_lsm_data{false};
608  if (lsm_in_arr) { valid_lsm_data = (lsm_in_arr(i,j,k) > 0.); }
609 
610  // Have LSM and are over land
611  if (is_land && valid_lsm_data) {
612  rrtmgp_to_fill(icol) = lsm_in_arr(i,j,k);
613  }
614  // We have a SurfLayer (enforce consistency with temperature)
615  else if (tsurf_arr && (ivar==0)) {
616  rrtmgp_to_fill(icol) = tsurf_arr(i,j,k);
617  }
618  // Use the default value
619  else {
620  rrtmgp_to_fill(icol) = rrtmgp_default_val;
621  }
622  });
623  } //mfi
624  } // ivar
625  Kokkos::deep_copy(lw_src, 0.0 );
626  } // have lsm
627 
628  // Enforce consistency between t_sfc and t_lev at bottom surface
629  Kokkos::parallel_for(Kokkos::RangePolicy(0, ncol),
630  KOKKOS_LAMBDA (int icol)
631  {
632  t_lev_tab(icol,0) = t_sfc_tab(icol);
633  });
634 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getPgivenRTh(const amrex::Real rhotheta, const amrex::Real qv=0.)
Definition: ERF_EOS.H:81
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real getTgivenRandRTh(const amrex::Real rho, const amrex::Real rhotheta, const amrex::Real qv=0.0)
Definition: ERF_EOS.H:46
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
#define RhoQ2_comp
Definition: ERF_IndexDefines.H:43
#define RhoQ3_comp
Definition: ERF_IndexDefines.H:44
#define RhoQ1_comp
Definition: ERF_IndexDefines.H:42
amrex::MultiFab * m_lon
Definition: ERF_Radiation.H:272
amrex::MultiFab * m_z_phys
Definition: ERF_Radiation.H:268
amrex::Geometry m_geom
Definition: ERF_Radiation.H:227
amrex::MultiFab * m_lat
Definition: ERF_Radiation.H:271
@ qv
Definition: ERF_Kessler.H:28
@ qc
Definition: ERF_SatAdj.H:36
Here is the call graph for this function:

◆ populateDatalogMF()

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

◆ rad_run_impl()

void Radiation::rad_run_impl ( amrex::Vector< amrex::MultiFab * > &  lsm_output_ptrs)
inline
173  {
174  if (m_update_rad) {
175  amrex::Print() << "Radiation advancing level " << m_lev << " at (YY-MM-DD SS) " << m_orbital_year << '-'
176  << m_orbital_mon << '-' << m_orbital_day << ' ' << m_orbital_sec << " ...";
177  this->initialize_impl();
178  this->run_impl();
179  this->finalize_impl(lsm_output_ptrs);
180  amrex::Print() << "DONE\n";
181  }
182  }
void initialize_impl()
Definition: ERF_Radiation.cpp:1031
void run_impl()
Definition: ERF_Radiation.cpp:1042
int m_orbital_mon
Definition: ERF_Radiation.H:341
void finalize_impl(amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
Definition: ERF_Radiation.cpp:1301
int m_orbital_sec
Definition: ERF_Radiation.H:343
int m_lev
Definition: ERF_Radiation.H:215
bool m_update_rad
Definition: ERF_Radiation.H:233
int m_orbital_day
Definition: ERF_Radiation.H:342

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::MultiFab *  lsm_fluxes,
amrex::MultiFab *  lsm_zenith,
amrex::Vector< amrex::MultiFab * > &  lsm_input_ptrs,
amrex::Vector< amrex::MultiFab * > &  lsm_output_ptrs,
amrex::MultiFab *  qheating_rates,
amrex::MultiFab *  rad_fluxes,
amrex::MultiFab *  z_phys,
amrex::MultiFab *  lat_ptr,
amrex::MultiFab *  lon_ptr 
)
inlineoverridevirtual

Implements IRadiation.

108  {
109  set_grids(level, step, time, dt, ba, geom,
110  cons_in, lmask, t_surf, lsm_fluxes,
111  lsm_zenith, lsm_input_ptrs, qheating_rates,
112  rad_fluxes, z_phys, lat_ptr, lon_ptr);
113  rad_run_impl(lsm_output_ptrs);
114  }
void set_grids(int &level, int &step, amrex::Real &time, const amrex::Real &dt, const amrex::BoxArray &ba, amrex::Geometry &geom, amrex::MultiFab *cons_in, amrex::iMultiFab *lmask, amrex::MultiFab *t_surf, amrex::MultiFab *lsm_fluxes, amrex::MultiFab *lsm_zenith, amrex::Vector< amrex::MultiFab * > &lsm_input_ptrs, amrex::MultiFab *qheating_rates, amrex::MultiFab *rad_fluxes, amrex::MultiFab *z_phys, amrex::MultiFab *lat, amrex::MultiFab *lon)
Definition: ERF_Radiation.cpp:134
void rad_run_impl(amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
Definition: ERF_Radiation.H:172
Here is the call graph for this function:

◆ run_impl()

void Radiation::run_impl ( )
1043 {
1044  // Local copies
1045  const auto ncol = m_ncol;
1046  const auto nlay = m_nlay;
1047  const auto nswbands = m_nswbands;
1048 
1049  // Compute orbital parameters; these are used both for computing
1050  // the solar zenith angle and also for computing total solar
1051  // irradiance scaling (tsi_scaling).
1052  double obliqr, lambm0, mvelpp;
1053  int orbital_year = m_orbital_year;
1054  double eccen = m_orbital_eccen;
1055  double obliq = m_orbital_obliq;
1056  double mvelp = m_orbital_mvelp;
1057  if (eccen >= 0 && obliq >= 0 && mvelp >= 0) {
1058  // fixed orbital parameters forced with orbital_year == ORB_UNDEF_INT
1059  orbital_year = ORB_UNDEF_INT;
1060  }
1061  orbital_params(orbital_year, eccen, obliq,
1062  mvelp, obliqr, lambm0, mvelpp);
1063 
1064  // Use the orbital parameters to calculate the solar declination and eccentricity factor
1065  double delta, eccf;
1066  // Want day + fraction; calday 1 == Jan 1 0Z
1067  static constexpr double dpy[] = {0.0 , 31.0, 59.0, 90.0, 120.0, 151.0,
1068  181.0, 212.0, 243.0, 273.0, 304.0, 334.0};
1069  bool leap = (m_orbital_year % 4 == 0 && (!(m_orbital_year % 100 == 0) || (m_orbital_year % 400 == 0))) ? true : false;
1070  double calday = dpy[m_orbital_mon-1] + (m_orbital_day-1.0) + m_orbital_sec/86400.0;
1071  // add extra day if leap year
1072  if (leap) { calday += 1.0; }
1073  orbital_decl(calday, eccen, mvelpp, lambm0, obliqr, delta, eccf);
1074 
1075  // Overwrite eccf if using a fixed solar constant.
1076  auto fixed_total_solar_irradiance = m_fixed_total_solar_irradiance;
1077  if (fixed_total_solar_irradiance >= 0){
1078  eccf = fixed_total_solar_irradiance/1360.9;
1079  }
1080 
1081  // Precompute volume mixing ratio (VMR) for all gases
1082  //
1083  // H2O is obtained from qv.
1084  // O3 may be a constant or a 1D vector
1085  // All other comps are set to constants for now
1086  for (int igas(0); igas < m_ngas; ++igas) {
1087  auto tmp2d = Kokkos::View<RealT**,layout_t,KokkosDefaultMem>("tmp2d", ncol, nlay);
1088  auto name = m_gas_names[igas];
1089  auto gas_mol_weight = m_mol_weight_gas[igas];
1090  if (name == "H2O") {
1091  auto qv_lay_d = qv_lay;
1092  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1093  KOKKOS_LAMBDA (int icol, int ilay)
1094  {
1095  tmp2d(icol,ilay) = qv_lay_d(icol,ilay) * mwdair/gas_mol_weight;
1096  });
1097  } else if (name == "CO2") {
1098  Kokkos::deep_copy(tmp2d, m_co2vmr);
1099  } else if (name == "O3") {
1100  if (m_o3_size==1) {
1101  Kokkos::deep_copy(tmp2d, m_o3vmr[0] );
1102  } else {
1103  auto o3_lay_d = o3_lay;
1104  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1105  KOKKOS_LAMBDA (int icol, int ilay)
1106  {
1107  tmp2d(icol,ilay) = o3_lay_d(ilay);
1108  });
1109  }
1110  } else if (name == "N2O") {
1111  Kokkos::deep_copy(tmp2d, m_n2ovmr);
1112  } else if (name == "CO") {
1113  Kokkos::deep_copy(tmp2d, m_covmr );
1114  } else if (name == "CH4") {
1115  Kokkos::deep_copy(tmp2d, m_ch4vmr);
1116  } else if (name == "O2") {
1117  Kokkos::deep_copy(tmp2d, m_o2vmr );
1118  } else if (name == "N2") {
1119  Kokkos::deep_copy(tmp2d, m_n2vmr );
1120  } else {
1121  Abort("Radiation: Unknown gas component.");
1122  }
1123 
1124  // Populate GasConcs object
1125  m_gas_concs.set_vmr(name, tmp2d);
1126  Kokkos::fence();
1127  }
1128 
1129  // Populate mu0 1D array
1130  // This must be done on HOST and copied to device.
1131  auto h_mu0 = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), mu0);
1132  if (m_fixed_solar_zenith_angle > 0) {
1133  Kokkos::deep_copy(h_mu0, m_fixed_solar_zenith_angle);
1134  } else {
1135  auto h_lat = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), lat);
1136  auto h_lon = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), lon);
1137  double dt = double(m_dt);
1138  auto rad_freq_in_steps = m_rad_freq_in_steps;
1139  Kokkos::parallel_for(Kokkos::RangePolicy<Kokkos::Serial>(0, ncol),
1140  [&] (int icol)
1141  {
1142  // Convert lat/lon to radians
1143  double lat_col = h_lat(icol)*PI/180.0;
1144  double lon_col = h_lon(icol)*PI/180.0;
1145  double lcalday = calday;
1146  double ldelta = delta;
1147  h_mu0(icol) = Real(orbital_cos_zenith(lcalday, lat_col, lon_col, ldelta, rad_freq_in_steps * dt));
1148  });
1149  }
1150  Kokkos::deep_copy(mu0, h_mu0);
1151 
1152  // Compute layer cloud mass per unit area (populates lwp/iwp)
1155 
1156  // Convert to g/m2 (needed by RRTMGP)
1157  Table2D<Real,Order::C> lwp_tab(lwp.data(), {0,0}, {static_cast<int>(lwp.extent(0)),static_cast<int>(lwp.extent(1))});
1158  Table2D<Real,Order::C> iwp_tab(iwp.data(), {0,0}, {static_cast<int>(iwp.extent(0)),static_cast<int>(iwp.extent(1))});
1159  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1160  KOKKOS_LAMBDA (int icol, int ilay)
1161  {
1162  lwp_tab(icol,ilay) *= 1.e3;
1163  iwp_tab(icol,ilay) *= 1.e3;
1164  });
1165 
1166  // Expand surface_albedos along nswbands.
1167  // This is needed since rrtmgp require band-by-band.
1172  // Run RRTMGP driver
1174  p_lay, t_lay,
1175  p_lev, t_lev,
1176  m_gas_concs,
1178  t_sfc, sfc_emis, lw_src,
1194 
1195 #if 0
1196  // UNIT TEST
1197  //================================================================================
1198  Kokkos::deep_copy(mu0, 0.86);
1199  Kokkos::deep_copy(sfc_alb_dir_vis, 0.06);
1200  Kokkos::deep_copy(sfc_alb_dir_nir, 0.06);
1201  Kokkos::deep_copy(sfc_alb_dif_vis, 0.06);
1202  Kokkos::deep_copy(sfc_alb_dif_nir, 0.06);
1203 
1204  Kokkos::deep_copy(aero_tau_sw, 0.0);
1205  Kokkos::deep_copy(aero_ssa_sw, 0.0);
1206  Kokkos::deep_copy(aero_g_sw , 0.0);
1207  Kokkos::deep_copy(aero_tau_lw, 0.0);
1208 
1209  // Generate some fake liquid and ice water data. We pick values to be midway between
1210  // the min and max of the valid lookup table values for effective radii
1211  real rel_val = 0.5 * (rrtmgp::cloud_optics_sw_k->get_min_radius_liq()
1212  + rrtmgp::cloud_optics_sw_k->get_max_radius_liq());
1213  real rei_val = 0.5 * (rrtmgp::cloud_optics_sw_k->get_min_radius_ice()
1214  + rrtmgp::cloud_optics_sw_k->get_max_radius_ice());
1215 
1216  // Restrict clouds to troposphere (> 100 hPa = 100*100 Pa) and not very close to the ground (< 900 hPa), and
1217  // put them in 2/3 of the columns since that's roughly the total cloudiness of earth.
1218  // Set sane values for liquid and ice water path.
1219  // NOTE: these "sane" values are in g/m2!
1220  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1221  KOKKOS_LAMBDA (int icol, int ilay)
1222  {
1223  cldfrac_tot(icol,ilay) = (p_lay(icol,ilay) > 100. * 100.) &&
1224  (p_lay(icol,ilay) < 900. * 100.) &&
1225  (icol%3 != 0);
1226  // Ice and liquid will overlap in a few layers
1227  lwp(icol,ilay) = (cldfrac_tot(icol,ilay) && t_lay(icol,ilay) > 263.) ? 10. : 0.;
1228  iwp(icol,ilay) = (cldfrac_tot(icol,ilay) && t_lay(icol,ilay) < 273.) ? 10. : 0.;
1229  eff_radius_qc(icol,ilay) = (lwp(icol,ilay) > 0.) ? rel_val : 0.;
1230  eff_radius_qi(icol,ilay) = (iwp(icol,ilay) > 0.) ? rei_val : 0.;
1231  });
1232 
1237 
1239  p_lay, t_lay,
1240  p_lev, t_lev,
1241  m_gas_concs,
1243  t_sfc, sfc_emis, lw_src,
1258  1.0, false, false);
1259  //================================================================================
1260 #endif
1261 
1262  // Update heating tendency
1265 
1266  /*
1267  // AML DEBUG
1268  Kokkos::parallel_for(nlay+1, KOKKOS_LAMBDA (int ilay)
1269  {
1270  printf("Fluxes: %i %e %e %e %e %e\n",ilay,
1271  sw_flux_up(5,ilay), sw_flux_dn(5,ilay), sw_flux_dn_dir(5,ilay),
1272  lw_flux_up(5,ilay), lw_flux_dn(5,ilay));
1273  });
1274  Kokkos::parallel_for(nlay, KOKKOS_LAMBDA (int ilay)
1275  {
1276  printf("Heating Rate: %i %e %e\n",ilay,sw_heating(5,ilay),lw_heating(5,ilay));
1277  });
1278  */
1279 
1280  // Compute surface fluxes
1281  const int kbot = 0;
1282  Table3D<Real,Order::C> sw_bnd_flux_dif_tab(sw_bnd_flux_dif.data(), {0,0,0},
1283  {static_cast<int>(sw_bnd_flux_dif.extent(0)),static_cast<int>(sw_bnd_flux_dif.extent(1)),static_cast<int>(sw_bnd_flux_dif.extent(2))});
1284  Table3D<Real,Order::C> sw_bnd_flux_dn_tab(sw_bnd_flux_dn.data(), {0,0,0},
1285  {static_cast<int>(sw_bnd_flux_dn.extent(0)),static_cast<int>(sw_bnd_flux_dn.extent(1)),static_cast<int>(sw_bnd_flux_dn.extent(2))});
1286  Table3D<Real,Order::C> sw_bnd_flux_dir_tab(sw_bnd_flux_dir.data(), {0,0,0},
1287  {static_cast<int>(sw_bnd_flux_dir.extent(0)),static_cast<int>(sw_bnd_flux_dir.extent(1)),static_cast<int>(sw_bnd_flux_dir.extent(2))});
1288  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<3>>({0, 0, 0}, {ncol, nlay+1, nswbands}),
1289  KOKKOS_LAMBDA (int icol, int ilay, int ibnd)
1290  {
1291  sw_bnd_flux_dif_tab(icol,ilay,ibnd) = sw_bnd_flux_dn_tab(icol,ilay,ibnd) - sw_bnd_flux_dir_tab(icol,ilay,ibnd);
1292  });
1293  rrtmgp::compute_broadband_surface_fluxes(ncol, kbot, nswbands,
1297 }
static constexpr int ORB_UNDEF_INT
Definition: ERF_Constants.H:96
constexpr amrex::Real mwdair
Definition: ERF_Constants.H:64
constexpr amrex::Real PI
Definition: ERF_Constants.H:6
double real
Definition: ERF_OrbCosZenith.H:9
AMREX_GPU_HOST AMREX_FORCE_INLINE real orbital_cos_zenith(real &jday, real &lat, real &lon, real &declin, real dt_avg=-1., real uniform_angle=-1., real constant_zenith_angle_deg=-1.)
Definition: ERF_OrbCosZenith.H:559
AMREX_GPU_HOST AMREX_FORCE_INLINE void orbital_decl(real &calday, real &eccen, real &mvelpp, real &lambm0, real &obliqr, real &delta, real &eccf)
Definition: ERF_OrbCosZenith.H:15
AMREX_GPU_HOST AMREX_FORCE_INLINE void orbital_params(int &iyear_AD, real &eccen, real &obliq, real &mvelp, real &obliqr, real &lambm0, real &mvelpp)
Definition: ERF_OrbCosZenith.H:81
real3d_k aero_tau_sw
Definition: ERF_Radiation.H:454
real3d_k cld_tau_lw_gpt
Definition: ERF_Radiation.H:465
real3d_k cld_tau_sw_gpt
Definition: ERF_Radiation.H:464
real3d_k cld_tau_lw_bnd
Definition: ERF_Radiation.H:461
real3d_k aero_g_sw
Definition: ERF_Radiation.H:456
amrex::Real m_dt
Definition: ERF_Radiation.H:224
real3d_k aero_ssa_sw
Definition: ERF_Radiation.H:455
real3d_k aero_tau_lw
Definition: ERF_Radiation.H:457
real3d_k cld_tau_sw_bnd
Definition: ERF_Radiation.H:460
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 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, real3d_k &, real3d_k &, real3d_k &, real3d_k &, real3d_k &, real3d_k &, real3d_k &, real3d_k &, real2d_k &sw_flux_up, real2d_k &sw_flux_dn, real2d_k &sw_flux_dn_dir, real2d_k &lw_flux_up, real2d_k &lw_flux_dn, real2d_k &sw_clnclrsky_flux_up, real2d_k &sw_clnclrsky_flux_dn, real2d_k &sw_clnclrsky_flux_dn_dir, real2d_k &sw_clrsky_flux_up, real2d_k &sw_clrsky_flux_dn, real2d_k &sw_clrsky_flux_dn_dir, real2d_k &sw_clnsky_flux_up, real2d_k &sw_clnsky_flux_dn, real2d_k &sw_clnsky_flux_dn_dir, real2d_k &lw_clnclrsky_flux_up, real2d_k &lw_clnclrsky_flux_dn, real2d_k &lw_clrsky_flux_up, real2d_k &lw_clrsky_flux_dn, real2d_k &lw_clnsky_flux_up, real2d_k &lw_clnsky_flux_dn, real3d_k &sw_bnd_flux_up, real3d_k &sw_bnd_flux_dn, real3d_k &sw_bnd_flux_dn_dir, real3d_k &lw_bnd_flux_up, real3d_k &lw_bnd_flux_dn, const RealT tsi_scaling, const bool extra_clnclrsky_diag, const bool extra_clnsky_diag)
Definition: ERF_RRTMGP_Interface.cpp:384
void mixing_ratio_to_cloud_mass(View1 const &mixing_ratio, View2 const &cloud_fraction, View3 const &rho, View4 const &dz, View5 const &cloud_mass)
Definition: ERF_RRTMGP_Utils.H:12
std::unique_ptr< cloud_optics_t > cloud_optics_sw_k
Definition: ERF_RRTMGP_Interface.cpp:18

Referenced by rad_run_impl().

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

◆ set_grids()

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

Referenced by Run().

Here is the caller graph for this function:

◆ write_rrtmgp_fluxes()

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

◆ WriteDataLog()

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

Implements IRadiation.

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

Member Data Documentation

◆ aero_g_sw

real3d_k Radiation::aero_g_sw
private

◆ aero_ssa_sw

real3d_k Radiation::aero_ssa_sw
private

◆ aero_tau_lw

real3d_k Radiation::aero_tau_lw
private

◆ aero_tau_sw

real3d_k Radiation::aero_tau_sw
private

◆ cld_tau_lw_bnd

real3d_k Radiation::cld_tau_lw_bnd
private

◆ cld_tau_lw_gpt

real3d_k Radiation::cld_tau_lw_gpt
private

◆ cld_tau_sw_bnd

real3d_k Radiation::cld_tau_sw_bnd
private

◆ cld_tau_sw_gpt

real3d_k Radiation::cld_tau_sw_gpt
private

◆ cldfrac_tot

real2d_k Radiation::cldfrac_tot
private

◆ d_tint

real2d_k Radiation::d_tint
private

◆ datalog_mf

amrex::MultiFab Radiation::datalog_mf
private

◆ eff_radius_qc

real2d_k Radiation::eff_radius_qc
private

◆ eff_radius_qi

real2d_k Radiation::eff_radius_qi
private

◆ emis_sfc

real2d_k Radiation::emis_sfc
private

◆ gas_names_offset

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

◆ iwp

real2d_k Radiation::iwp
private

◆ lat

real1d_k Radiation::lat
private

◆ lon

real1d_k Radiation::lon
private

◆ lw_bnd_flux_dn

real3d_k Radiation::lw_bnd_flux_dn
private

◆ lw_bnd_flux_up

real3d_k Radiation::lw_bnd_flux_up
private

◆ lw_clnclrsky_flux_dn

real2d_k Radiation::lw_clnclrsky_flux_dn
private

◆ lw_clnclrsky_flux_up

real2d_k Radiation::lw_clnclrsky_flux_up
private

◆ lw_clnsky_flux_dn

real2d_k Radiation::lw_clnsky_flux_dn
private

◆ lw_clnsky_flux_up

real2d_k Radiation::lw_clnsky_flux_up
private

◆ lw_clrsky_flux_dn

real2d_k Radiation::lw_clrsky_flux_dn
private

◆ lw_clrsky_flux_up

real2d_k Radiation::lw_clrsky_flux_up
private

◆ lw_clrsky_heating

real2d_k Radiation::lw_clrsky_heating
private

◆ lw_flux_dn

real2d_k Radiation::lw_flux_dn
private

◆ lw_flux_up

real2d_k Radiation::lw_flux_up
private

◆ lw_heating

real2d_k Radiation::lw_heating
private

◆ lw_src

real1d_k Radiation::lw_src
private

◆ lwp

real2d_k Radiation::lwp
private

◆ m_ba

amrex::BoxArray Radiation::m_ba
private

◆ m_ch4vmr

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

◆ m_co2vmr

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

◆ m_col_offsets

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

Referenced by Init().

◆ m_cons_in

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

◆ m_covmr

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

◆ m_do_aerosol_rad

bool Radiation::m_do_aerosol_rad = true
private

◆ m_do_subcol_sampling

bool Radiation::m_do_subcol_sampling = true
private

◆ m_dt

amrex::Real Radiation::m_dt
private

◆ m_extra_clnclrsky_diag

bool Radiation::m_extra_clnclrsky_diag = false
private

◆ m_extra_clnsky_diag

bool Radiation::m_extra_clnsky_diag = false
private

◆ m_first_step

bool Radiation::m_first_step = true
private

◆ m_fixed_orbital_year

bool Radiation::m_fixed_orbital_year = false
private

◆ m_fixed_solar_zenith_angle

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

◆ m_fixed_total_solar_irradiance

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

◆ m_gas_concs

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

◆ m_gas_mol_weights

real1d_k Radiation::m_gas_mol_weights
private

◆ m_gas_names

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

◆ m_geom

amrex::Geometry Radiation::m_geom
private

◆ m_ice

bool Radiation::m_ice = false
private

◆ m_lat

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

◆ m_lat_cons

amrex::Real Radiation::m_lat_cons = 39.809860
private

◆ m_lev

int Radiation::m_lev
private

Referenced by rad_run_impl().

◆ m_lon

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

◆ m_lon_cons

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

◆ m_lsm

bool Radiation::m_lsm = false
private

◆ m_lsm_fluxes

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

◆ m_lsm_input_names

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

Referenced by get_lsm_input_varnames().

◆ m_lsm_output_names

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

Referenced by get_lsm_output_varnames().

◆ m_lsm_zenith

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

◆ m_moist

bool Radiation::m_moist = false
private

◆ m_mol_weight_gas

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

◆ m_n2ovmr

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

◆ m_n2vmr

amrex::Real Radiation::m_n2vmr = 0.7906
private

◆ m_ncol

int Radiation::m_ncol
private

Referenced by Init().

◆ m_ngas

int Radiation::m_ngas = 8
private

◆ m_nlay

int Radiation::m_nlay
private

Referenced by Init().

◆ m_nlwbands

int Radiation::m_nlwbands
private

◆ m_nlwgpts

int Radiation::m_nlwgpts
private

◆ m_nswbands

int Radiation::m_nswbands
private

◆ m_nswgpts

int Radiation::m_nswgpts
private

◆ m_o2vmr

amrex::Real Radiation::m_o2vmr = 0.209448
private

◆ m_o3_size

int Radiation::m_o3_size
private

◆ m_o3vmr

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

◆ m_orbital_day

int Radiation::m_orbital_day = -9999
private

Referenced by rad_run_impl().

◆ m_orbital_eccen

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

◆ m_orbital_mon

int Radiation::m_orbital_mon = -9999
private

Referenced by rad_run_impl().

◆ m_orbital_mvelp

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

◆ m_orbital_obliq

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

◆ m_orbital_sec

int Radiation::m_orbital_sec = -9999
private

Referenced by rad_run_impl().

◆ m_orbital_year

int Radiation::m_orbital_year = -9999
private

Referenced by rad_run_impl().

◆ m_qheating_rates

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

◆ m_rad_fluxes

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

◆ m_rad_freq_in_steps

int Radiation::m_rad_freq_in_steps = 1
private

◆ m_rad_t_sfc

amrex::Real Radiation::m_rad_t_sfc = -1
private

◆ m_rad_write_fluxes

bool Radiation::m_rad_write_fluxes = false
private

◆ m_step

int Radiation::m_step
private

◆ m_time

amrex::Real Radiation::m_time
private

◆ m_update_rad

bool Radiation::m_update_rad = false
private

Referenced by rad_run_impl().

◆ m_z_phys

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

◆ mu0

real1d_k Radiation::mu0
private

◆ o3_lay

real1d_k Radiation::o3_lay
private

◆ p_lay

real2d_k Radiation::p_lay
private

◆ p_lev

real2d_k Radiation::p_lev
private

◆ qc_lay

real2d_k Radiation::qc_lay
private

◆ qi_lay

real2d_k Radiation::qi_lay
private

◆ qv_lay

real2d_k Radiation::qv_lay
private

◆ r_lay

real2d_k Radiation::r_lay
private

◆ rrtmgp_cloud_optics_file_lw

std::string Radiation::rrtmgp_cloud_optics_file_lw
private

◆ rrtmgp_cloud_optics_file_sw

std::string Radiation::rrtmgp_cloud_optics_file_sw
private

◆ rrtmgp_cloud_optics_lw

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

◆ rrtmgp_cloud_optics_sw

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

◆ rrtmgp_coeffs_file_lw

std::string Radiation::rrtmgp_coeffs_file_lw
private

◆ rrtmgp_coeffs_file_sw

std::string Radiation::rrtmgp_coeffs_file_sw
private

◆ rrtmgp_coeffs_lw

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

◆ rrtmgp_coeffs_sw

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

◆ rrtmgp_file_path

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

◆ sfc_alb_dif

real2d_k Radiation::sfc_alb_dif
private

◆ sfc_alb_dif_nir

real1d_k Radiation::sfc_alb_dif_nir
private

◆ sfc_alb_dif_vis

real1d_k Radiation::sfc_alb_dif_vis
private

◆ sfc_alb_dir

real2d_k Radiation::sfc_alb_dir
private

◆ sfc_alb_dir_nir

real1d_k Radiation::sfc_alb_dir_nir
private

◆ sfc_alb_dir_vis

real1d_k Radiation::sfc_alb_dir_vis
private

◆ sfc_emis

real1d_k Radiation::sfc_emis
private

◆ sfc_flux_dif_nir

real1d_k Radiation::sfc_flux_dif_nir
private

◆ sfc_flux_dif_vis

real1d_k Radiation::sfc_flux_dif_vis
private

◆ sfc_flux_dir_nir

real1d_k Radiation::sfc_flux_dir_nir
private

◆ sfc_flux_dir_vis

real1d_k Radiation::sfc_flux_dir_vis
private

◆ sw_bnd_flux_dif

real3d_k Radiation::sw_bnd_flux_dif
private

◆ sw_bnd_flux_dir

real3d_k Radiation::sw_bnd_flux_dir
private

◆ sw_bnd_flux_dn

real3d_k Radiation::sw_bnd_flux_dn
private

◆ sw_bnd_flux_up

real3d_k Radiation::sw_bnd_flux_up
private

◆ sw_clnclrsky_flux_dn

real2d_k Radiation::sw_clnclrsky_flux_dn
private

◆ sw_clnclrsky_flux_dn_dir

real2d_k Radiation::sw_clnclrsky_flux_dn_dir
private

◆ sw_clnclrsky_flux_up

real2d_k Radiation::sw_clnclrsky_flux_up
private

◆ sw_clnsky_flux_dn

real2d_k Radiation::sw_clnsky_flux_dn
private

◆ sw_clnsky_flux_dn_dir

real2d_k Radiation::sw_clnsky_flux_dn_dir
private

◆ sw_clnsky_flux_up

real2d_k Radiation::sw_clnsky_flux_up
private

◆ sw_clrsky_flux_dn

real2d_k Radiation::sw_clrsky_flux_dn
private

◆ sw_clrsky_flux_dn_dir

real2d_k Radiation::sw_clrsky_flux_dn_dir
private

◆ sw_clrsky_flux_up

real2d_k Radiation::sw_clrsky_flux_up
private

◆ sw_clrsky_heating

real2d_k Radiation::sw_clrsky_heating
private

◆ sw_flux_dn

real2d_k Radiation::sw_flux_dn
private

◆ sw_flux_dn_dir

real2d_k Radiation::sw_flux_dn_dir
private

◆ sw_flux_up

real2d_k Radiation::sw_flux_up
private

◆ sw_heating

real2d_k Radiation::sw_heating
private

◆ t_lay

real2d_k Radiation::t_lay
private

◆ t_lev

real2d_k Radiation::t_lev
private

◆ t_sfc

real1d_k Radiation::t_sfc
private

◆ z_del

real2d_k Radiation::z_del
private

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