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 = 39.809860
 
amrex::Real m_lon_cons = -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 = 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 }
ParmParse pp("prob")
std::string rrtmgp_coeffs_file_sw
Definition: ERF_Radiation.H:286
std::string rrtmgp_coeffs_sw
Definition: ERF_Radiation.H:282
int m_rad_freq_in_steps
Definition: ERF_Radiation.H:364
amrex::Real m_lon_cons
Definition: ERF_Radiation.H:275
int m_nswbands
Definition: ERF_Radiation.H:358
bool m_do_aerosol_rad
Definition: ERF_Radiation.H:326
std::string rrtmgp_coeffs_lw
Definition: ERF_Radiation.H:283
bool m_moist
Definition: ERF_Radiation.H:237
std::string rrtmgp_cloud_optics_file_lw
Definition: ERF_Radiation.H:289
bool m_lsm
Definition: ERF_Radiation.H:241
bool m_do_subcol_sampling
Definition: ERF_Radiation.H:367
bool m_rad_write_fluxes
Definition: ERF_Radiation.H:232
std::string rrtmgp_cloud_optics_file_sw
Definition: ERF_Radiation.H:288
amrex::Real m_orbital_mvelp
Definition: ERF_Radiation.H:346
amrex::Vector< amrex::Real > m_o3vmr
Definition: ERF_Radiation.H:300
std::string rrtmgp_cloud_optics_sw
Definition: ERF_Radiation.H:284
int m_nlwgpts
Definition: ERF_Radiation.H:361
std::string rrtmgp_cloud_optics_lw
Definition: ERF_Radiation.H:285
amrex::Real m_co2vmr
Definition: ERF_Radiation.H:299
bool m_extra_clnsky_diag
Definition: ERF_Radiation.H:329
amrex::Real m_lat_cons
Definition: ERF_Radiation.H:274
amrex::Real m_fixed_total_solar_irradiance
Definition: ERF_Radiation.H:351
amrex::Real m_o2vmr
Definition: ERF_Radiation.H:304
amrex::Real m_n2ovmr
Definition: ERF_Radiation.H:301
amrex::Real m_rad_t_sfc
Definition: ERF_Radiation.H:255
std::string rrtmgp_file_path
Definition: ERF_Radiation.H:281
int m_nlwbands
Definition: ERF_Radiation.H:359
amrex::Real m_orbital_eccen
Definition: ERF_Radiation.H:344
amrex::Real m_covmr
Definition: ERF_Radiation.H:302
std::string rrtmgp_coeffs_file_lw
Definition: ERF_Radiation.H:287
amrex::Real m_ch4vmr
Definition: ERF_Radiation.H:303
bool m_ice
Definition: ERF_Radiation.H:238
amrex::Real m_fixed_solar_zenith_angle
Definition: ERF_Radiation.H:355
amrex::Real m_orbital_obliq
Definition: ERF_Radiation.H:345
amrex::Real m_n2vmr
Definition: ERF_Radiation.H:305
bool m_fixed_orbital_year
Definition: ERF_Radiation.H:343
int m_nswgpts
Definition: ERF_Radiation.H:360
bool m_extra_clnclrsky_diag
Definition: ERF_Radiation.H:330
int m_orbital_year
Definition: ERF_Radiation.H:335
static NCFile open(const std::string &name, const int cmode=NC_NOWRITE)
Definition: ERF_NCInterface.cpp:707
LandSurfaceType lsm_type
Definition: ERF_DataStruct.H:1197
MoistureType moisture_type
Definition: ERF_DataStruct.H:1194
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 ( )
202 {
203  // 1d size (m_ngas)
204  const Real* mol_weight_gas_p = m_mol_weight_gas.data();
205  const std::string* gas_names_p = m_gas_names.data();
206  m_gas_mol_weights = real1d_k("m_gas_mol_weights", m_ngas);
207  realHost1d_k m_gas_mol_weights_h("m_gas_mol_weights_h", m_ngas);
208  gas_names_offset.clear(); gas_names_offset.resize(m_ngas);
209  std::string* gas_names_offset_p = gas_names_offset.data();
210  Kokkos::parallel_for(Kokkos::RangePolicy<Kokkos::Serial>(0, m_ngas),
211  [&] (int igas)
212  {
213  m_gas_mol_weights_h(igas) = mol_weight_gas_p[igas];
214  gas_names_offset_p[igas] = gas_names_p[igas];
215  });
216  Kokkos::deep_copy(m_gas_mol_weights, m_gas_mol_weights_h);
217 
218  // 1d size (1 or nlay)
219  m_o3_size = m_o3vmr.size();
220  AMREX_ALWAYS_ASSERT_WITH_MESSAGE(((m_o3_size==1) || (m_o3_size==m_nlay)),
221  "O3 VMR array must be length 1 or nlay");
222  Real* o3vmr_p = m_o3vmr.data();
223  o3_lay = real1d_k("o3_lay", m_o3_size);
224  realHost1d_k o3_lay_h("o3_lay_h", m_o3_size);
225  Kokkos::parallel_for(Kokkos::RangePolicy<Kokkos::Serial>(0, m_o3_size),
226  [&] (int io3)
227  {
228  o3_lay_h(io3) = o3vmr_p[io3];
229  });
230  Kokkos::deep_copy(o3_lay, o3_lay_h);
231 
232  // 1d size (ncol)
233  mu0 = real1d_k("mu0" , m_ncol);
234  sfc_alb_dir_vis = real1d_k("sfc_alb_dir_vis" , m_ncol);
235  sfc_alb_dir_nir = real1d_k("sfc_alb_dir_nir" , m_ncol);
236  sfc_alb_dif_vis = real1d_k("sfc_alb_dif_vis" , m_ncol);
237  sfc_alb_dif_nir = real1d_k("sfc_alb_dif_nir" , m_ncol);
238  sfc_flux_dir_vis = real1d_k("sfc_flux_dir_vis", m_ncol);
239  sfc_flux_dir_nir = real1d_k("sfc_flux_dir_nir", m_ncol);
240  sfc_flux_dif_vis = real1d_k("sfc_flux_dif_vis", m_ncol);
241  sfc_flux_dif_nir = real1d_k("sfc_flux_dif_nir", m_ncol);
242  lat = real1d_k("lat" , m_ncol);
243  lon = real1d_k("lon" , m_ncol);
244  sfc_emis = real1d_k("sfc_emis" , m_ncol);
245  t_sfc = real1d_k("t_sfc" , m_ncol);
246  lw_src = real1d_k("lw_src" , m_ncol);
247 
248  // 2d size (ncol, nlay)
249  r_lay = real2d_k("r_lay" , m_ncol, m_nlay);
250  p_lay = real2d_k("p_lay" , m_ncol, m_nlay);
251  t_lay = real2d_k("t_lay" , m_ncol, m_nlay);
252  z_del = real2d_k("z_del" , m_ncol, m_nlay);
253  qv_lay = real2d_k("qv" , m_ncol, m_nlay);
254  qc_lay = real2d_k("qc" , m_ncol, m_nlay);
255  qi_lay = real2d_k("qi" , m_ncol, m_nlay);
256  cldfrac_tot = real2d_k("cldfrac_tot" , m_ncol, m_nlay);
257  eff_radius_qc = real2d_k("eff_radius_qc", m_ncol, m_nlay);
258  eff_radius_qi = real2d_k("eff_radius_qi", m_ncol, m_nlay);
259  lwp = real2d_k("lwp" , m_ncol, m_nlay);
260  iwp = real2d_k("iwp" , m_ncol, m_nlay);
261  sw_heating = real2d_k("sw_heating" , m_ncol, m_nlay);
262  lw_heating = real2d_k("lw_heating" , m_ncol, m_nlay);
263  sw_clrsky_heating = real2d_k("sw_clrsky_heating", m_ncol, m_nlay);
264  lw_clrsky_heating = real2d_k("lw_clrsky_heating", m_ncol, m_nlay);
265 
266  // 2d size (ncol, nlay+1)
267  d_tint = real2d_k("d_tint" , m_ncol, m_nlay+1);
268  p_lev = real2d_k("p_lev" , m_ncol, m_nlay+1);
269  t_lev = real2d_k("t_lev" , m_ncol, m_nlay+1);
270 
271  sw_flux_up = real2d_k("sw_flux_up" , m_ncol, m_nlay+1);
272  sw_flux_dn = real2d_k("sw_flux_dn" , m_ncol, m_nlay+1);
273  sw_flux_dn_dir = real2d_k("sw_flux_dn_dir" , m_ncol, m_nlay+1);
274 
275  lw_flux_up = real2d_k("lw_flux_up" , m_ncol, m_nlay+1);
276  lw_flux_dn = real2d_k("lw_flux_dn" , m_ncol, m_nlay+1);
277 
278  sw_clnclrsky_flux_up = real2d_k("sw_clnclrsky_flux_up" , m_ncol, m_nlay+1);
279  sw_clnclrsky_flux_dn = real2d_k("sw_clnclrsky_flux_dn" , m_ncol, m_nlay+1);
280  sw_clnclrsky_flux_dn_dir = real2d_k("sw_clnclrsky_flux_dn_dir", m_ncol, m_nlay+1);
281  sw_clrsky_flux_up = real2d_k("sw_clrsky_flux_up" , m_ncol, m_nlay+1);
282  sw_clrsky_flux_dn = real2d_k("sw_clrsky_flux_dn" , m_ncol, m_nlay+1);
283  sw_clrsky_flux_dn_dir = real2d_k("sw_clrsky_flux_dn_dir" , m_ncol, m_nlay+1);
284  sw_clnsky_flux_up = real2d_k("sw_clnsky_flux_up" , m_ncol, m_nlay+1);
285  sw_clnsky_flux_dn = real2d_k("sw_clnsky_flux_dn" , m_ncol, m_nlay+1);
286  sw_clnsky_flux_dn_dir = real2d_k("sw_clnsky_flux_dn_dir" , m_ncol, m_nlay+1);
287 
288  lw_clnclrsky_flux_up = real2d_k("lw_clnclrsky_flux_up" , m_ncol, m_nlay+1);
289  lw_clnclrsky_flux_dn = real2d_k("lw_clnclrsky_flux_dn" , m_ncol, m_nlay+1);
290  lw_clrsky_flux_up = real2d_k("lw_clrsky_flux_up" , m_ncol, m_nlay+1);
291  lw_clrsky_flux_dn = real2d_k("lw_clrsky_flux_dn" , m_ncol, m_nlay+1);
292  lw_clnsky_flux_up = real2d_k("lw_clnsky_flux_up" , m_ncol, m_nlay+1);
293  lw_clnsky_flux_dn = real2d_k("lw_clnsky_flux_dn" , m_ncol, m_nlay+1);
294 
295  // 3d size (ncol, nlay+1, nswbands)
296  sw_bnd_flux_up = real3d_k("sw_bnd_flux_up" , m_ncol, m_nlay+1, m_nswbands);
297  sw_bnd_flux_dn = real3d_k("sw_bnd_flux_dn" , m_ncol, m_nlay+1, m_nswbands);
298  sw_bnd_flux_dir = real3d_k("sw_bnd_flux_dir", m_ncol, m_nlay+1, m_nswbands);
299  sw_bnd_flux_dif = real3d_k("sw_bnd_flux_dif", m_ncol, m_nlay+1, m_nswbands);
300 
301  // 3d size (ncol, nlay+1, nlwbands)
302  lw_bnd_flux_up = real3d_k("lw_bnd_flux_up" , m_ncol, m_nlay+1, m_nlwbands);
303  lw_bnd_flux_dn = real3d_k("lw_bnd_flux_dn" , m_ncol, m_nlay+1, m_nlwbands);
304 
305  // 2d size (ncol, nswbands)
306  sfc_alb_dir = real2d_k("sfc_alb_dir", m_ncol, m_nswbands);
307  sfc_alb_dif = real2d_k("sfc_alb_dif", m_ncol, m_nswbands);
308 
309  // 2d size (ncol, nlwbands)
310  //emis_sfc = real2d_k("emis_sfc", m_ncol, m_nlwbands);
311 
312  /*
313  // 3d size (ncol, nlay, n[sw,lw]bands)
314  aero_tau_sw = real3d_k("aero_tau_sw", m_ncol, m_nlay, m_nswbands);
315  aero_ssa_sw = real3d_k("aero_ssa_sw", m_ncol, m_nlay, m_nswbands);
316  aero_g_sw = real3d_k("aero_g_sw" , m_ncol, m_nlay, m_nswbands);
317  aero_tau_lw = real3d_k("aero_tau_lw", m_ncol, m_nlay, m_nlwbands);
318 
319  // 3d size (ncol, nlay, n[sw,lw]bnds)
320  cld_tau_sw_bnd = real3d_k("cld_tau_sw_bnd", m_ncol, m_nlay, m_nswbands);
321  cld_tau_lw_bnd = real3d_k("cld_tau_lw_bnd", m_ncol, m_nlay, m_nlwbands);
322 
323  // 3d size (ncol, nlay, n[sw,lw]gpts)
324  cld_tau_sw_gpt = real3d_k("cld_tau_sw_gpt", m_ncol, m_nlay, m_nswgpts);
325  cld_tau_lw_gpt = real3d_k("cld_tau_lw_gpt", m_ncol, m_nlay, m_nlwgpts);
326  */
327 }
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:433
real2d_k lw_flux_up
Definition: ERF_Radiation.H:413
int m_o3_size
Definition: ERF_Radiation.H:309
real3d_k sw_bnd_flux_dir
Definition: ERF_Radiation.H:434
real2d_k sw_clnsky_flux_dn
Definition: ERF_Radiation.H:422
real2d_k d_tint
Definition: ERF_Radiation.H:407
real2d_k lw_clnclrsky_flux_dn
Definition: ERF_Radiation.H:425
real2d_k lwp
Definition: ERF_Radiation.H:399
real1d_k lw_src
Definition: ERF_Radiation.H:386
real2d_k eff_radius_qi
Definition: ERF_Radiation.H:398
real1d_k m_gas_mol_weights
Definition: ERF_Radiation.H:310
real2d_k sw_heating
Definition: ERF_Radiation.H:401
real3d_k lw_bnd_flux_dn
Definition: ERF_Radiation.H:439
real3d_k sw_bnd_flux_up
Definition: ERF_Radiation.H:432
real2d_k qv_lay
Definition: ERF_Radiation.H:393
real2d_k sw_clnclrsky_flux_dn_dir
Definition: ERF_Radiation.H:417
real1d_k sfc_flux_dif_vis
Definition: ERF_Radiation.H:380
real2d_k lw_clnclrsky_flux_up
Definition: ERF_Radiation.H:424
real1d_k lat
Definition: ERF_Radiation.H:382
real2d_k qi_lay
Definition: ERF_Radiation.H:395
real2d_k sw_clrsky_flux_up
Definition: ERF_Radiation.H:418
real2d_k t_lev
Definition: ERF_Radiation.H:409
real1d_k o3_lay
Definition: ERF_Radiation.H:370
real1d_k sfc_alb_dif_vis
Definition: ERF_Radiation.H:376
real1d_k mu0
Definition: ERF_Radiation.H:373
real2d_k cldfrac_tot
Definition: ERF_Radiation.H:396
real1d_k sfc_flux_dir_nir
Definition: ERF_Radiation.H:379
real2d_k sw_clnclrsky_flux_up
Definition: ERF_Radiation.H:415
real2d_k sw_flux_up
Definition: ERF_Radiation.H:410
real3d_k sw_bnd_flux_dif
Definition: ERF_Radiation.H:435
real2d_k sfc_alb_dif
Definition: ERF_Radiation.H:443
real1d_k sfc_alb_dif_nir
Definition: ERF_Radiation.H:377
real2d_k lw_clnsky_flux_dn
Definition: ERF_Radiation.H:429
real2d_k p_lay
Definition: ERF_Radiation.H:390
real2d_k r_lay
Definition: ERF_Radiation.H:389
int m_ncol
Definition: ERF_Radiation.H:319
real2d_k sw_flux_dn_dir
Definition: ERF_Radiation.H:412
const std::vector< amrex::Real > m_mol_weight_gas
Definition: ERF_Radiation.H:295
std::vector< std::string > gas_names_offset
Definition: ERF_Radiation.H:311
real2d_k sw_clrsky_flux_dn
Definition: ERF_Radiation.H:419
real2d_k sw_clnsky_flux_dn_dir
Definition: ERF_Radiation.H:423
real2d_k qc_lay
Definition: ERF_Radiation.H:394
real2d_k lw_clrsky_flux_up
Definition: ERF_Radiation.H:426
real2d_k sw_flux_dn
Definition: ERF_Radiation.H:411
real2d_k lw_clrsky_flux_dn
Definition: ERF_Radiation.H:427
real2d_k z_del
Definition: ERF_Radiation.H:392
real2d_k lw_flux_dn
Definition: ERF_Radiation.H:414
real2d_k sw_clrsky_flux_dn_dir
Definition: ERF_Radiation.H:420
real1d_k sfc_flux_dif_nir
Definition: ERF_Radiation.H:381
int m_ngas
Definition: ERF_Radiation.H:292
real2d_k lw_heating
Definition: ERF_Radiation.H:402
real1d_k sfc_alb_dir_nir
Definition: ERF_Radiation.H:375
real2d_k lw_clnsky_flux_up
Definition: ERF_Radiation.H:428
real2d_k sfc_alb_dir
Definition: ERF_Radiation.H:442
real1d_k lon
Definition: ERF_Radiation.H:383
real1d_k sfc_emis
Definition: ERF_Radiation.H:384
real2d_k sw_clnsky_flux_up
Definition: ERF_Radiation.H:421
int m_nlay
Definition: ERF_Radiation.H:320
real3d_k lw_bnd_flux_up
Definition: ERF_Radiation.H:438
real1d_k sfc_alb_dir_vis
Definition: ERF_Radiation.H:374
const std::vector< std::string > m_gas_names
Definition: ERF_Radiation.H:293
real2d_k sw_clrsky_heating
Definition: ERF_Radiation.H:403
real2d_k t_lay
Definition: ERF_Radiation.H:391
real1d_k sfc_flux_dir_vis
Definition: ERF_Radiation.H:378
real2d_k iwp
Definition: ERF_Radiation.H:400
real2d_k p_lev
Definition: ERF_Radiation.H:408
real2d_k lw_clrsky_heating
Definition: ERF_Radiation.H:404
real2d_k eff_radius_qc
Definition: ERF_Radiation.H:397
real1d_k t_sfc
Definition: ERF_Radiation.H:385
real2d_k sw_clnclrsky_flux_dn
Definition: ERF_Radiation.H:416

◆ dealloc_buffers()

void Radiation::dealloc_buffers ( )
331 {
332  // 1d size (m_ngas)
334 
335  // 1d size (1 or nlay)
336  o3_lay = real1d_k();
337 
338  // 1d size (ncol)
339  mu0 = real1d_k();
348  lat = real1d_k();
349  lon = real1d_k();
350  sfc_emis = real1d_k();
351  t_sfc = real1d_k();
352  lw_src = real1d_k();
353 
354  // 2d size (ncol, nlay)
355  r_lay = real2d_k();
356  p_lay = real2d_k();
357  t_lay = real2d_k();
358  z_del = real2d_k();
359  qv_lay = real2d_k();
360  qc_lay = real2d_k();
361  qi_lay = real2d_k();
362  cldfrac_tot = real2d_k();
365  lwp = real2d_k();
366  iwp = real2d_k();
367  sw_heating = real2d_k();
368  lw_heating = real2d_k();
371  sw_heating = real2d_k();
372  lw_heating = real2d_k();
375 
376  // 2d size (ncol, nlay+1)
377  d_tint = real2d_k();
378  p_lev = real2d_k();
379  t_lev = real2d_k();
380  sw_flux_up = real2d_k();
381  sw_flux_dn = real2d_k();
383  lw_flux_up = real2d_k();
384  lw_flux_dn = real2d_k();
400 
401  // 3d size (ncol, nlay+1, nswbands)
406 
407  // 3d size (ncol, nlay+1, nlwbands)
410 
411  // 2d size (ncol, nswbands)
412  sfc_alb_dir = real2d_k();
413  sfc_alb_dif = real2d_k();
414 
415  // 2d size (ncol, nlwbands)
416  //emis_sfc = real2d_k();
417 
418  /*
419  // 3d size (ncol, nlay, n[sw,lw]bands)
420  aero_tau_sw = real3d_k();
421  aero_ssa_sw = real3d_k();
422  aero_g_sw = real3d_k();
423  aero_tau_lw = real3d_k();
424 
425  // 3d size (ncol, nlay, n[sw,lw]bnds)
426  cld_tau_sw_bnd = real3d_k();
427  cld_tau_lw_bnd = real3d_k();
428 
429  // 3d size (ncol, nlay, n[sw,lw]gpts)
430  cld_tau_sw_gpt = real3d_k();
431  cld_tau_lw_gpt = real3d_k();
432  */
433 }

◆ finalize_impl()

void Radiation::finalize_impl ( amrex::Vector< amrex::MultiFab * > &  lsm_output_ptrs)
1281 {
1282  // Finish rrtmgp
1283  m_gas_concs.reset();
1285 
1286  // Fill the AMReX MFs from Kokkos Views
1287  kokkos_buffers_to_mf(lsm_output_ptrs);
1288 
1289  // Write fluxes if requested
1291 
1292  // Fill output data for datalog before deallocating
1293  if (datalog_int > 0) {
1296 
1298  }
1299 
1300  // Deallocate the buffer arrays
1301  dealloc_buffers();
1302 }
int datalog_int
Definition: ERF_RadiationInterface.H:88
void dealloc_buffers()
Definition: ERF_Radiation.cpp:330
void populateDatalogMF()
Definition: ERF_Radiation.cpp:768
void kokkos_buffers_to_mf(amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
Definition: ERF_Radiation.cpp:649
GasConcsK< amrex::Real, layout_t, KokkosDefaultDevice > m_gas_concs
Definition: ERF_Radiation.H:313
void write_rrtmgp_fluxes()
Definition: ERF_Radiation.cpp:728
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.

184  {
185  return m_lsm_input_names;
186  }
amrex::Vector< std::string > m_lsm_input_names
Definition: ERF_Radiation.H:244

◆ get_lsm_output_varnames()

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

Reimplemented from IRadiation.

192  {
193  return m_lsm_output_names;
194  }
amrex::Vector< std::string > m_lsm_output_names
Definition: ERF_Radiation.H:249

◆ 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  };
const int klo
Definition: ERF_InitCustomPertVels_ParticleTests.H:4
const int khi
Definition: ERF_InitCustomPert_Bubble.H:21
amrex::Vector< int > m_col_offsets
Definition: ERF_Radiation.H:323

◆ initialize_impl()

void Radiation::initialize_impl ( )
1010 {
1011  // Call API to initialize
1016 }
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)
650 {
651  // Heating rate, fluxes, zenith, lsm ptrs
652 
653  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))});
654  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))});
655  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))});
656  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))});
657  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))});
658  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))});
659  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))});
660 
661  TableData<Real,1> sfc_flux_sw_dn; sfc_flux_sw_dn.resize({0}, {static_cast<int>(sw_flux_dn.extent(0))});
662  TableData<Real,1> sfc_flux_lw_dn; sfc_flux_lw_dn.resize({0}, {static_cast<int>(lw_flux_dn.extent(0))});
663  Table1D<Real> sfc_flux_sw_dn_tab = sfc_flux_sw_dn.table();
664  Table1D<Real> sfc_flux_lw_dn_tab = sfc_flux_lw_dn.table();
665  Table1D<Real> sfc_flux_sw_dir_vis_tab(sfc_flux_dir_vis.data(), {0}, {static_cast<int>(sfc_flux_dir_vis.extent(0))});
666  Table1D<Real> sfc_flux_sw_dir_nir_tab(sfc_flux_dir_nir.data(), {0}, {static_cast<int>(sfc_flux_dir_nir.extent(0))});
667  Table1D<Real> sfc_flux_sw_dif_vis_tab(sfc_flux_dif_vis.data(), {0}, {static_cast<int>(sfc_flux_dif_vis.extent(0))});
668  Table1D<Real> sfc_flux_sw_dif_nir_tab(sfc_flux_dif_nir.data(), {0}, {static_cast<int>(sfc_flux_dif_nir.extent(0))});
669  Table1D<Real> mu0_tab(mu0.data(), {0}, {static_cast<int>(mu0.extent(0))});
670  Vector<Table1D<Real>> rrtmgp_out_vars = {mu0_tab , sfc_flux_sw_dn_tab ,
671  sfc_flux_sw_dir_vis_tab, sfc_flux_sw_dir_nir_tab,
672  sfc_flux_sw_dif_vis_tab, sfc_flux_sw_dif_nir_tab,
673  sfc_flux_lw_dn_tab };
674 
675  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
676  const auto& vbx = mfi.validbox();
677  const auto& sbx = makeSlab(vbx,2,vbx.smallEnd(2));
678  const int nx = vbx.length(0);
679  const int imin = vbx.smallEnd(0);
680  const int jmin = vbx.smallEnd(1);
681  const int offset = m_col_offsets[mfi.index()];
682  const Array4<Real>& q_arr = m_qheating_rates->array(mfi);
683  const Array4<Real>& f_arr = m_rad_fluxes->array(mfi);
684  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
685  {
686  // map [i,j,k] 0-based to [icol, ilay] 0-based
687  const int icol = (j-jmin)*nx + (i-imin) + offset;
688  const int ilay = k;
689 
690  // Temperature heating rate for SW and LW
691  q_arr(i,j,k,0) = sw_heating_tab(icol,ilay);
692  q_arr(i,j,k,1) = lw_heating_tab(icol,ilay);
693 
694  // Convert the dT/dz to dTheta/dz
695  Real iexner = 1./getExnergivenP(Real(p_lay_tab(icol,ilay)), R_d/Cp_d);
696  q_arr(i,j,k,0) *= iexner;
697  q_arr(i,j,k,1) *= iexner;
698 
699  // Populate the fluxes
700  f_arr(i,j,k,0) = sw_flux_up_tab(icol,ilay);
701  f_arr(i,j,k,1) = sw_flux_dn_tab(icol,ilay);
702  f_arr(i,j,k,2) = lw_flux_up_tab(icol,ilay);
703  f_arr(i,j,k,3) = lw_flux_dn_tab(icol,ilay);
704 
705  if (k==0) {
706  sfc_flux_sw_dn_tab(icol) = sw_flux_dn_tab(icol,ilay);
707  sfc_flux_lw_dn_tab(icol) = lw_flux_dn_tab(icol,ilay);
708  }
709  });
710  for (int ivar(0); ivar<lsm_output_ptrs.size(); ivar++) {
711  if (lsm_output_ptrs[ivar]) {
712  auto rrtmgp_for_fill = rrtmgp_out_vars[ivar];
713  const Array4<Real>& lsm_out_arr = lsm_output_ptrs[ivar]->array(mfi);
714  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
715  {
716  // map [i,j,k] 0-based to [icol, ilay] 0-based
717  const int icol = (j-jmin)*nx + (i-imin) + offset;
718 
719  // export the desired variable at surface
720  lsm_out_arr(i,j,k) = rrtmgp_for_fill(icol);
721  });
722  } // valid ptr
723  } // ivar
724  }// mfi
725 }
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
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const auto prob_lo=geomdata.ProbLo();const auto dx=geomdata.CellSize();const Real x=(prob_lo[0]+(i+0.5) *dx[0])/mf_m(i, j, 0);const Real z=z_cc(i, j, k);Real L=std::sqrt(std::pow((x - x_c)/x_r, 2)+std::pow((z - z_c)/z_r, 2));if(L<=1.0) { Real dT=T_pert *(std::cos(PI *L)+1.0)/2.0;Real Tbar_hse=p_hse(i, j, k)/(R_d *r_hse(i, j, k));Real theta_perturbed=(Tbar_hse+dT) *std::pow(p_0/p_hse(i, j, k), rdOcp);Real theta_0=(Tbar_hse) *std::pow(p_0/p_hse(i, j, k), rdOcp);if(const_rho) { state_pert(i, j, k, RhoTheta_comp)=r_hse(i, j, k) *(theta_perturbed - theta_0);} else { state_pert(i, j, k, Rho_comp)=getRhoThetagivenP(p_hse(i, j, k))/theta_perturbed - r_hse(i, j, k);} } })
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:258
amrex::MultiFab * m_rad_fluxes
Definition: ERF_Radiation.H:264
amrex::MultiFab * m_qheating_rates
Definition: ERF_Radiation.H:261
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 
)
440 {
441  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))});
442  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))});
443  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))});
444  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))});
445  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))});
446  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))});
447  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))});
448  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))});
449 
450  Table2D<Real,Order::C> lwp_tab(lwp.data(), {0,0}, {static_cast<int>(lwp.extent(0)),static_cast<int>(lwp.extent(1))});
451  Table2D<Real,Order::C> iwp_tab(iwp.data(), {0,0}, {static_cast<int>(iwp.extent(0)),static_cast<int>(iwp.extent(1))});
452  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))});
453  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))});
454 
455  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))});
456  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))});
457 
458  Table1D<Real> lat_tab(lat.data(), {0}, {static_cast<int>(lat.extent(0))});
459  Table1D<Real> lon_tab(lon.data(), {0}, {static_cast<int>(lon.extent(0))});
460  Table1D<Real> t_sfc_tab(t_sfc.data(), {0}, {static_cast<int>(t_sfc.extent(0))});
461 
462  bool moist = m_moist;
463  bool ice = m_ice;
464  const bool has_lsm = m_lsm;
465  const bool has_lat = m_lat;
466  const bool has_lon = m_lon;
467  const bool has_surflayer = (t_surf);
468  int ncol = m_ncol;
469  int nlay = m_nlay;
470  Real dz = m_geom.CellSize(2);
471  Real cons_lat = m_lat_cons;
472  Real cons_lon = m_lon_cons;
473  Real rad_t_sfc = m_rad_t_sfc;
474 
475  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
476  const auto& vbx = mfi.validbox();
477  const int nx = vbx.length(0);
478  const int imin = vbx.smallEnd(0);
479  const int jmin = vbx.smallEnd(1);
480  const int offset = m_col_offsets[mfi.index()];
481  const Array4<const Real>& cons_arr = m_cons_in->const_array(mfi);
482  const Array4<const Real>& z_arr = (m_z_phys) ? m_z_phys->const_array(mfi) :
483  Array4<const Real>{};
484  const Array4<const Real>& lat_arr = (m_lat) ? m_lat->const_array(mfi) :
485  Array4<const Real>{};
486  const Array4<const Real>& lon_arr = (m_lon) ? m_lon->const_array(mfi) :
487  Array4<const Real>{};
488  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
489  {
490  // map [i,j,k] 0-based to [icol, ilay] 0-based
491  const int icol = (j-jmin)*nx + (i-imin) + offset;
492  const int ilay = k;
493 
494  // EOS input (at CC)
495  Real r = cons_arr(i,j,k,Rho_comp);
496  Real rt = cons_arr(i,j,k,RhoTheta_comp);
497  Real qv = (moist) ? std::max(cons_arr(i,j,k,RhoQ1_comp)/r,0.0) : 0.0;
498  Real qc = (moist) ? std::max(cons_arr(i,j,k,RhoQ2_comp)/r,0.0) : 0.0;
499  Real qi = (ice) ? std::max(cons_arr(i,j,k,RhoQ3_comp)/r,0.0) : 0.0;
500 
501  // EOS avg to z-face
502  Real r_lo = cons_arr(i,j,k-1,Rho_comp);
503  Real rt_lo = cons_arr(i,j,k-1,RhoTheta_comp);
504  Real qv_lo = (moist) ? cons_arr(i,j,k-1,RhoQ1_comp)/r_lo : 0.0;
505  Real dz_k = (z_arr) ? 0.125 * ( (z_arr(i ,j ,k+1) - z_arr(i ,j ,k))
506  + (z_arr(i+1,j ,k+1) - z_arr(i+1,j ,k))
507  + (z_arr(i ,j+1,k+1) - z_arr(i ,j+1,k))
508  + (z_arr(i+1,j+1,k+1) - z_arr(i+1,j+1,k)) ) : 0.5*dz; // Dist from w-face to CC at k
509  Real dz_km1 = (z_arr) ? 0.125 * ( (z_arr(i ,j ,k ) - z_arr(i ,j ,k-1))
510  + (z_arr(i+1,j ,k ) - z_arr(i+1,j ,k-1))
511  + (z_arr(i ,j+1,k ) - z_arr(i ,j+1,k-1))
512  + (z_arr(i+1,j+1,k ) - z_arr(i+1,j+1,k-1)) ) : 0.5*dz; // Dist from w-face to CC at k-1
513  Real r_avg = (dz_k*r + dz_km1*r_lo ) / (dz_k + dz_km1);
514  Real rt_avg = (dz_k*rt + dz_km1*rt_lo) / (dz_k + dz_km1);
515  Real qv_avg = (dz_k*qv + dz_km1*qv_lo) / (dz_k + dz_km1);
516 
517  // Views at CC
518  r_lay_tab(icol,ilay) = r;
519 
520  p_lay_tab(icol,ilay) = getPgivenRTh(rt, qv);
521  t_lay_tab(icol,ilay) = getTgivenRandRTh(r, rt, qv);
522  z_del_tab(icol,ilay) = (z_arr) ? 0.25 * ( (z_arr(i ,j ,k+1) - z_arr(i ,j ,k))
523  + (z_arr(i+1,j ,k+1) - z_arr(i+1,j ,k))
524  + (z_arr(i ,j+1,k+1) - z_arr(i ,j+1,k))
525  + (z_arr(i+1,j+1,k+1) - z_arr(i+1,j+1,k)) ) : dz;
526  qv_lay_tab(icol,ilay) = qv;
527  qc_lay_tab(icol,ilay) = qc;
528  qi_lay_tab(icol,ilay) = qi;
529  cldfrac_tot_tab(icol,ilay) = ((qc+qi)>0.0) ? 1. : 0.;
530 
531  // NOTE: These are populated in 'mixing_ratio_to_cloud_mass'
532  lwp_tab(icol,ilay) = 0.0;
533  iwp_tab(icol,ilay) = 0.0;
534 
535  // NOTE: These would be populated from P3 (we use the constants in p3_main_impl.hpp)
536  // NOTE: These are in units of micron!
537  eff_radius_qc_tab(icol,ilay) = (qc>0.0) ? 10.0 : 0.0;
538  eff_radius_qi_tab(icol,ilay) = (qi>0.0) ? 25.0 : 0.0;
539 
540  // Buffers on z-faces (nlay+1)
541  p_lev_tab(icol,ilay) = getPgivenRTh(rt_avg, qv_avg);
542  t_lev_tab(icol,ilay) = getTgivenRandRTh(r_avg, rt_avg, qv_avg);
543  if (ilay==(nlay-1)) {
544  Real r_hi = cons_arr(i,j,k+1,Rho_comp);
545  Real rt_hi = cons_arr(i,j,k+1,RhoTheta_comp);
546  Real qv_hi = (moist) ? std::max(cons_arr(i,j,k+1,RhoQ1_comp)/r_hi,0.0) : 0.0;
547  Real dz_kp1 = (z_arr) ? 0.125 * ( (z_arr(i ,j ,k+2) - z_arr(i ,j ,k+1))
548  + (z_arr(i+1,j ,k+2) - z_arr(i+1,j ,k+1))
549  + (z_arr(i ,j+1,k+2) - z_arr(i ,j+1,k+1))
550  + (z_arr(i+1,j+1,k+2) - z_arr(i+1,j+1,k+1)) ) : 0.5*dz; // Dist from w-face to CC at k+1
551  r_avg = (dz_k*r + dz_kp1*r_hi ) / (dz_k + dz_kp1);
552  rt_avg = (dz_k*rt + dz_kp1*rt_hi) / (dz_k + dz_kp1);
553  qv_avg = (dz_k*qv + dz_kp1*qv_hi) / (dz_k + dz_kp1);
554  p_lev_tab(icol,ilay+1) = getPgivenRTh(rt_avg, qv_avg);
555  t_lev_tab(icol,ilay+1) = getTgivenRandRTh(r_avg, rt_avg, qv_avg);
556  }
557 
558  // 1D data structures
559  if (k==0) {
560  lat_tab(icol) = (has_lat) ? lat_arr(i,j,0) : cons_lat;
561  lon_tab(icol) = (has_lon) ? lon_arr(i,j,0) : cons_lon;
562  }
563 
564  });
565  } // mfi
566 
567  // Populate vars LSM would provide
568  if (!has_lsm && !has_surflayer) {
569  // Parsed surface temp
570  Kokkos::deep_copy(t_sfc, rad_t_sfc);
571 
572  // EAMXX dummy atmos constants
573  Kokkos::deep_copy(sfc_alb_dir_vis, 0.06);
574  Kokkos::deep_copy(sfc_alb_dir_nir, 0.06);
575  Kokkos::deep_copy(sfc_alb_dif_vis, 0.06);
576  Kokkos::deep_copy(sfc_alb_dif_nir, 0.06);
577 
578  // AML NOTE: These are not used in current EAMXX, I've left
579  // the code to plug into these if we need it.
580  //
581  // Current EAMXX constants
582  Kokkos::deep_copy(sfc_emis, 0.98);
583  Kokkos::deep_copy(lw_src , 0.0 );
584  } else {
585  Vector<real1d_k> rrtmgp_in_vars = {t_sfc, sfc_emis,
588  Vector<Real> rrtmgp_default_vals = {rad_t_sfc, 0.98,
589  0.06, 0.06,
590  0.06, 0.06};
591  for (int ivar(0); ivar<lsm_input_ptrs.size(); ivar++) {
592  auto rrtmgp_default_val = rrtmgp_default_vals[ivar];
593  auto rrtmgp_to_fill_k = rrtmgp_in_vars[ivar];
594  amrex::Table1D<amrex::Real> rrtmgp_to_fill(rrtmgp_to_fill_k.data(),
595  0, rrtmgp_to_fill_k.extent(0));
596  for (MFIter mfi(*m_cons_in); mfi.isValid(); ++mfi) {
597  const auto& vbx = mfi.validbox();
598  const auto& sbx = makeSlab(vbx,2,vbx.smallEnd(2));
599  const int nx = vbx.length(0);
600  const int imin = vbx.smallEnd(0);
601  const int jmin = vbx.smallEnd(1);
602  const int offset = m_col_offsets[mfi.index()];
603  const Array4<const int>& lmask_arr = (lmask) ? lmask->const_array(mfi) :
604  Array4<const int> {};
605  const Array4<const Real>& tsurf_arr = (t_surf) ? t_surf->const_array(mfi) :
606  Array4<const Real> {};
607  const Array4<const Real>& lsm_in_arr = (lsm_input_ptrs[ivar]) ? lsm_input_ptrs[ivar]->const_array(mfi) :
608  Array4<const Real> {};
609  ParallelFor(sbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
610  {
611  // map [i,j,k] 0-based to [icol, ilay] 0-based
612  const int icol = (j-jmin)*nx + (i-imin) + offset;
613 
614  // Check if over land
615  bool is_land = (lmask_arr) ? lmask_arr(i,j,k) : 1;
616 
617  // Check if valid LSM data
618  bool valid_lsm_data{false};
619  if (lsm_in_arr) { valid_lsm_data = (lsm_in_arr(i,j,k) >= 0.); }
620 
621  // Have LSM and are over land
622  if (is_land && valid_lsm_data) {
623  rrtmgp_to_fill(icol) = lsm_in_arr(i,j,k);
624  }
625  // We have a SurfLayer (enforce consistency with temperature)
626  else if (tsurf_arr && (ivar==0)) {
627  rrtmgp_to_fill(icol) = tsurf_arr(i,j,k);
628  }
629  // Use the default value
630  else {
631  rrtmgp_to_fill(icol) = rrtmgp_default_val;
632  }
633  });
634  } //mfi
635  } // ivar
636  Kokkos::deep_copy(lw_src, 0.0 );
637  } // have lsm
638 
639  // Enforce consistency between t_sfc and t_lev at bottom surface
640  Kokkos::parallel_for(Kokkos::RangePolicy(0, ncol),
641  KOKKOS_LAMBDA (int icol)
642  {
643  t_lev_tab(icol,0) = t_sfc_tab(icol);
644  });
645 }
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
const Real dz
Definition: ERF_InitCustomPert_Bubble.H:25
amrex::MultiFab * m_lon
Definition: ERF_Radiation.H:271
amrex::MultiFab * m_z_phys
Definition: ERF_Radiation.H:267
amrex::Geometry m_geom
Definition: ERF_Radiation.H:223
amrex::MultiFab * m_lat
Definition: ERF_Radiation.H:270
@ qv
Definition: ERF_Kessler.H:28
@ qc
Definition: ERF_SatAdj.H:36
Here is the call graph for this function:

◆ populateDatalogMF()

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

◆ rad_run_impl()

void Radiation::rad_run_impl ( amrex::Vector< amrex::MultiFab * > &  lsm_output_ptrs)
inline
169  {
170  if (m_update_rad) {
171  amrex::Print() << "Radiation advancing level " << m_lev << " at (YY-MM-DD SS) " << m_orbital_year << '-'
172  << m_orbital_mon << '-' << m_orbital_day << ' ' << m_orbital_sec << " ...";
173  this->initialize_impl();
174  this->run_impl();
175  this->finalize_impl(lsm_output_ptrs);
176  amrex::Print() << "DONE\n";
177  }
178  }
void initialize_impl()
Definition: ERF_Radiation.cpp:1009
void run_impl()
Definition: ERF_Radiation.cpp:1020
int m_orbital_mon
Definition: ERF_Radiation.H:336
void finalize_impl(amrex::Vector< amrex::MultiFab * > &lsm_output_ptrs)
Definition: ERF_Radiation.cpp:1280
int m_orbital_sec
Definition: ERF_Radiation.H:338
int m_lev
Definition: ERF_Radiation.H:211
bool m_update_rad
Definition: ERF_Radiation.H:229
int m_orbital_day
Definition: ERF_Radiation.H:337

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.

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

◆ run_impl()

void Radiation::run_impl ( )
1021 {
1022  // Local copies
1023  const auto ncol = m_ncol;
1024  const auto nlay = m_nlay;
1025  const auto nswbands = m_nswbands;
1026 
1027  // Compute orbital parameters; these are used both for computing
1028  // the solar zenith angle and also for computing total solar
1029  // irradiance scaling (tsi_scaling).
1030  double obliqr, lambm0, mvelpp;
1031  int orbital_year = m_orbital_year;
1032  double eccen = m_orbital_eccen;
1033  double obliq = m_orbital_obliq;
1034  double mvelp = m_orbital_mvelp;
1035  if (eccen >= 0 && obliq >= 0 && mvelp >= 0) {
1036  // fixed orbital parameters forced with orbital_year == ORB_UNDEF_INT
1037  orbital_year = ORB_UNDEF_INT;
1038  }
1039  orbital_params(orbital_year, eccen, obliq,
1040  mvelp, obliqr, lambm0, mvelpp);
1041 
1042  // Use the orbital parameters to calculate the solar declination and eccentricity factor
1043  double delta, eccf;
1044  // Want day + fraction; calday 1 == Jan 1 0Z
1045  static constexpr double dpy[] = {0.0 , 31.0, 59.0, 90.0, 120.0, 151.0,
1046  181.0, 212.0, 243.0, 273.0, 304.0, 334.0};
1047  bool leap = (m_orbital_year % 4 == 0 && (!(m_orbital_year % 100 == 0) || (m_orbital_year % 400 == 0))) ? true : false;
1048  double calday = 1.0 + dpy[m_orbital_mon-1] + (m_orbital_day-1.0) + m_orbital_sec/86400.0;
1049  // add extra day if leap year and past February
1050  if (leap && m_orbital_mon>2) { calday += 1.0; }
1051  orbital_decl(calday, eccen, mvelpp, lambm0, obliqr, delta, eccf);
1052 
1053  // Overwrite eccf if using a fixed solar constant.
1054  auto fixed_total_solar_irradiance = m_fixed_total_solar_irradiance;
1055  if (fixed_total_solar_irradiance >= 0){
1056  eccf = fixed_total_solar_irradiance/1360.9;
1057  }
1058 
1059  // Precompute volume mixing ratio (VMR) for all gases
1060  //
1061  // H2O is obtained from qv.
1062  // O3 may be a constant or a 1D vector
1063  // All other comps are set to constants for now
1064  for (int igas(0); igas < m_ngas; ++igas) {
1065  auto tmp2d = Kokkos::View<RealT**,layout_t,KokkosDefaultMem>("tmp2d", ncol, nlay);
1066  auto name = m_gas_names[igas];
1067  auto gas_mol_weight = m_mol_weight_gas[igas];
1068  if (name == "H2O") {
1069  auto qv_lay_d = qv_lay;
1070  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1071  KOKKOS_LAMBDA (int icol, int ilay)
1072  {
1073  tmp2d(icol,ilay) = qv_lay_d(icol,ilay) * mwdair/gas_mol_weight;
1074  });
1075  } else if (name == "CO2") {
1076  Kokkos::deep_copy(tmp2d, m_co2vmr);
1077  } else if (name == "O3") {
1078  if (m_o3_size==1) {
1079  Kokkos::deep_copy(tmp2d, m_o3vmr[0] );
1080  } else {
1081  auto o3_lay_d = o3_lay;
1082  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1083  KOKKOS_LAMBDA (int icol, int ilay)
1084  {
1085  tmp2d(icol,ilay) = o3_lay_d(ilay);
1086  });
1087  }
1088  } else if (name == "N2O") {
1089  Kokkos::deep_copy(tmp2d, m_n2ovmr);
1090  } else if (name == "CO") {
1091  Kokkos::deep_copy(tmp2d, m_covmr );
1092  } else if (name == "CH4") {
1093  Kokkos::deep_copy(tmp2d, m_ch4vmr);
1094  } else if (name == "O2") {
1095  Kokkos::deep_copy(tmp2d, m_o2vmr );
1096  } else if (name == "N2") {
1097  Kokkos::deep_copy(tmp2d, m_n2vmr );
1098  } else {
1099  Abort("Radiation: Unknown gas component.");
1100  }
1101 
1102  // Populate GasConcs object
1103  m_gas_concs.set_vmr(name, tmp2d);
1104  Kokkos::fence();
1105  }
1106 
1107  // Populate mu0 1D array
1108  // This must be done on HOST and copied to device.
1109  auto h_mu0 = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), mu0);
1110  if (m_fixed_solar_zenith_angle > 0) {
1111  Kokkos::deep_copy(h_mu0, m_fixed_solar_zenith_angle);
1112  } else {
1113  auto h_lat = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), lat);
1114  auto h_lon = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), lon);
1115  double dt = double(m_dt);
1116  auto rad_freq_in_steps = m_rad_freq_in_steps;
1117  Kokkos::parallel_for(Kokkos::RangePolicy<Kokkos::Serial>(0, ncol),
1118  [&] (int icol)
1119  {
1120  // Convert lat/lon to radians
1121  double lat_col = h_lat(icol)*PI/180.0;
1122  double lon_col = h_lon(icol)*PI/180.0;
1123  double lcalday = calday;
1124  double ldelta = delta;
1125  double dt_avg = static_cast<double>(rad_freq_in_steps) * dt;
1126  h_mu0(icol) = Real(orbital_cos_zenith(lcalday, lat_col, lon_col, ldelta, dt_avg));
1127  });
1128  }
1129  Kokkos::deep_copy(mu0, h_mu0);
1130 
1131  // Compute layer cloud mass per unit area (populates lwp/iwp)
1134 
1135  // Convert to g/m2 (needed by RRTMGP)
1136  Table2D<Real,Order::C> lwp_tab(lwp.data(), {0,0}, {static_cast<int>(lwp.extent(0)),static_cast<int>(lwp.extent(1))});
1137  Table2D<Real,Order::C> iwp_tab(iwp.data(), {0,0}, {static_cast<int>(iwp.extent(0)),static_cast<int>(iwp.extent(1))});
1138  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1139  KOKKOS_LAMBDA (int icol, int ilay)
1140  {
1141  lwp_tab(icol,ilay) *= 1.e3;
1142  iwp_tab(icol,ilay) *= 1.e3;
1143  });
1144 
1145  // Expand surface_albedos along nswbands.
1146  // This is needed since rrtmgp require band-by-band.
1151  // Run RRTMGP driver
1153  p_lay, t_lay,
1154  p_lev, t_lev,
1155  m_gas_concs,
1157  t_sfc, sfc_emis, lw_src,
1173 
1174 #if 0
1175  // UNIT TEST
1176  //================================================================================
1177  Kokkos::deep_copy(mu0, 0.86);
1178  Kokkos::deep_copy(sfc_alb_dir_vis, 0.06);
1179  Kokkos::deep_copy(sfc_alb_dir_nir, 0.06);
1180  Kokkos::deep_copy(sfc_alb_dif_vis, 0.06);
1181  Kokkos::deep_copy(sfc_alb_dif_nir, 0.06);
1182 
1183  Kokkos::deep_copy(aero_tau_sw, 0.0);
1184  Kokkos::deep_copy(aero_ssa_sw, 0.0);
1185  Kokkos::deep_copy(aero_g_sw , 0.0);
1186  Kokkos::deep_copy(aero_tau_lw, 0.0);
1187 
1188  // Generate some fake liquid and ice water data. We pick values to be midway between
1189  // the min and max of the valid lookup table values for effective radii
1190  real rel_val = 0.5 * (rrtmgp::cloud_optics_sw_k->get_min_radius_liq()
1191  + rrtmgp::cloud_optics_sw_k->get_max_radius_liq());
1192  real rei_val = 0.5 * (rrtmgp::cloud_optics_sw_k->get_min_radius_ice()
1193  + rrtmgp::cloud_optics_sw_k->get_max_radius_ice());
1194 
1195  // Restrict clouds to troposphere (> 100 hPa = 100*100 Pa) and not very close to the ground (< 900 hPa), and
1196  // put them in 2/3 of the columns since that's roughly the total cloudiness of earth.
1197  // Set sane values for liquid and ice water path.
1198  // NOTE: these "sane" values are in g/m2!
1199  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<2>>({0, 0}, {ncol, nlay}),
1200  KOKKOS_LAMBDA (int icol, int ilay)
1201  {
1202  cldfrac_tot(icol,ilay) = (p_lay(icol,ilay) > 100. * 100.) &&
1203  (p_lay(icol,ilay) < 900. * 100.) &&
1204  (icol%3 != 0);
1205  // Ice and liquid will overlap in a few layers
1206  lwp(icol,ilay) = (cldfrac_tot(icol,ilay) && t_lay(icol,ilay) > 263.) ? 10. : 0.;
1207  iwp(icol,ilay) = (cldfrac_tot(icol,ilay) && t_lay(icol,ilay) < 273.) ? 10. : 0.;
1208  eff_radius_qc(icol,ilay) = (lwp(icol,ilay) > 0.) ? rel_val : 0.;
1209  eff_radius_qi(icol,ilay) = (iwp(icol,ilay) > 0.) ? rei_val : 0.;
1210  });
1211 
1216 
1218  p_lay, t_lay,
1219  p_lev, t_lev,
1220  m_gas_concs,
1222  t_sfc, sfc_emis, lw_src,
1237  1.0, false, false);
1238  //================================================================================
1239 #endif
1240 
1241  // Update heating tendency
1244 
1245  /*
1246  // AML DEBUG
1247  Kokkos::parallel_for(nlay+1, KOKKOS_LAMBDA (int ilay)
1248  {
1249  printf("Fluxes: %i %e %e %e %e %e\n",ilay,
1250  sw_flux_up(5,ilay), sw_flux_dn(5,ilay), sw_flux_dn_dir(5,ilay),
1251  lw_flux_up(5,ilay), lw_flux_dn(5,ilay));
1252  });
1253  Kokkos::parallel_for(nlay, KOKKOS_LAMBDA (int ilay)
1254  {
1255  printf("Heating Rate: %i %e %e\n",ilay,sw_heating(5,ilay),lw_heating(5,ilay));
1256  });
1257  */
1258 
1259  // Compute surface fluxes
1260  const int kbot = 0;
1261  Table3D<Real,Order::C> sw_bnd_flux_dif_tab(sw_bnd_flux_dif.data(), {0,0,0},
1262  {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))});
1263  Table3D<Real,Order::C> sw_bnd_flux_dn_tab(sw_bnd_flux_dn.data(), {0,0,0},
1264  {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))});
1265  Table3D<Real,Order::C> sw_bnd_flux_dir_tab(sw_bnd_flux_dir.data(), {0,0,0},
1266  {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))});
1267  Kokkos::parallel_for(Kokkos::MDRangePolicy<Kokkos::Rank<3>>({0, 0, 0}, {ncol, nlay+1, nswbands}),
1268  KOKKOS_LAMBDA (int icol, int ilay, int ibnd)
1269  {
1270  sw_bnd_flux_dif_tab(icol,ilay,ibnd) = sw_bnd_flux_dn_tab(icol,ilay,ibnd) - sw_bnd_flux_dir_tab(icol,ilay,ibnd);
1271  });
1272  rrtmgp::compute_broadband_surface_fluxes(ncol, kbot, nswbands,
1276 }
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: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
real3d_k aero_tau_sw
Definition: ERF_Radiation.H:449
real3d_k cld_tau_lw_gpt
Definition: ERF_Radiation.H:460
real3d_k cld_tau_sw_gpt
Definition: ERF_Radiation.H:459
real3d_k cld_tau_lw_bnd
Definition: ERF_Radiation.H:456
real3d_k aero_g_sw
Definition: ERF_Radiation.H:451
amrex::Real m_dt
Definition: ERF_Radiation.H:220
real3d_k aero_ssa_sw
Definition: ERF_Radiation.H:450
real3d_k aero_tau_lw
Definition: ERF_Radiation.H:452
real3d_k cld_tau_sw_bnd
Definition: ERF_Radiation.H:455
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::Vector< amrex::MultiFab * > &  lsm_input_ptrs,
amrex::MultiFab *  qheating_rates,
amrex::MultiFab *  rad_fluxes,
amrex::MultiFab *  z_phys,
amrex::MultiFab *  lat,
amrex::MultiFab *  lon 
)
150 {
151  // Set data members that may change
152  m_lev = level;
153  m_step = step;
154  m_time = time;
155  m_dt = dt;
156  m_geom = geom;
157  m_cons_in = cons_in;
158  m_qheating_rates = qheating_rates;
159  m_rad_fluxes = rad_fluxes;
160  m_z_phys = z_phys;
161  m_lat = lat;
162  m_lon = lon;
163 
164  // Update the day and month
165  time_t timestamp = time_t(time);
166  struct tm *timeinfo = gmtime(&timestamp);
167  if (m_fixed_orbital_year) {
168  m_orbital_mon = timeinfo->tm_mon + 1;
169  m_orbital_day = timeinfo->tm_mday;
170  m_orbital_sec = timeinfo->tm_hour*3600 + timeinfo->tm_min*60 + timeinfo->tm_sec;
171  } else {
172  m_orbital_year = timeinfo->tm_year + 1900;
173  m_orbital_mon = timeinfo->tm_mon + 1;
174  m_orbital_day = timeinfo->tm_mday;
175  m_orbital_sec = timeinfo->tm_hour*3600 + timeinfo->tm_min*60 + timeinfo->tm_sec;
176  }
177 
178  // Only allocate and proceed if we are going to update radiation
179  m_update_rad = false;
180  if (m_rad_freq_in_steps > 0) { m_update_rad = ( (m_step == 0) || (m_step % m_rad_freq_in_steps == 0) ); }
181 
182  if (m_update_rad) {
183  // Call to Init() has set the dimensions: ncol & nlay
184 
185  // Allocate the buffer arrays
186  alloc_buffers();
187 
188  // Fill the KOKKOS Views from AMReX MFs
189  mf_to_kokkos_buffers(lmask, t_surf, lsm_input_ptrs);
190 
191  // Initialize datalog MF on first step
192  if (m_first_step) {
193  m_first_step = false;
194  datalog_mf.define(cons_in->boxArray(), cons_in->DistributionMap(), 25, 0);
195  datalog_mf.setVal(0.0);
196  }
197  }
198 }
int m_step
Definition: ERF_Radiation.H:214
void mf_to_kokkos_buffers(amrex::iMultiFab *lmask, amrex::MultiFab *t_surf, amrex::Vector< amrex::MultiFab * > &lsm_input_ptrs)
Definition: ERF_Radiation.cpp:437
bool m_first_step
Definition: ERF_Radiation.H:234
void alloc_buffers()
Definition: ERF_Radiation.cpp:201
amrex::Real m_time
Definition: ERF_Radiation.H:217

Referenced by Run().

Here is the caller graph for this function:

◆ write_rrtmgp_fluxes()

void Radiation::write_rrtmgp_fluxes ( )
729 {
730  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))});
731  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))});
732  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))});
733  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))});
734  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))});
735 
736  int n_fluxes = 5;
737  MultiFab mf_flux(m_cons_in->boxArray(), m_cons_in->DistributionMap(), n_fluxes, 0);
738 
739  for (MFIter mfi(mf_flux); mfi.isValid(); ++mfi) {
740  const auto& vbx = mfi.validbox();
741  const int nx = vbx.length(0);
742  const int imin = vbx.smallEnd(0);
743  const int jmin = vbx.smallEnd(1);
744  const int offset = m_col_offsets[mfi.index()];
745  const Array4<Real>& dst_arr = mf_flux.array(mfi);
746  ParallelFor(vbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
747  {
748  // map [i,j,k] 0-based to [icol, ilay] 0-based
749  const int icol = (j-jmin)*nx + (i-imin) + offset;
750  const int ilay = k;
751 
752  // SW and LW fluxes
753  dst_arr(i,j,k,0) = sw_flux_up_tab(icol,ilay);
754  dst_arr(i,j,k,1) = sw_flux_dn_tab(icol,ilay);
755  dst_arr(i,j,k,2) = sw_flux_dn_dir_tab(icol,ilay);
756  dst_arr(i,j,k,3) = lw_flux_up_tab(icol,ilay);
757  dst_arr(i,j,k,4) = lw_flux_dn_tab(icol,ilay);
758  });
759  }
760 
761 
762  std::string plotfilename = amrex::Concatenate("plt_rad", m_step, 5);
763  Vector<std::string> flux_names = {"sw_flux_up", "sw_flux_dn", "sw_flux_dir",
764  "lw_flux_up", "lw_flux_dn"};
765  WriteSingleLevelPlotfile(plotfilename, mf_flux, flux_names, m_geom, m_time, m_step);
766 }
Here is the call graph for this function:

◆ WriteDataLog()

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

Implements IRadiation.

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

Member Data Documentation

◆ aero_g_sw

real3d_k Radiation::aero_g_sw
private

◆ aero_ssa_sw

real3d_k Radiation::aero_ssa_sw
private

◆ aero_tau_lw

real3d_k Radiation::aero_tau_lw
private

◆ aero_tau_sw

real3d_k Radiation::aero_tau_sw
private

◆ cld_tau_lw_bnd

real3d_k Radiation::cld_tau_lw_bnd
private

◆ cld_tau_lw_gpt

real3d_k Radiation::cld_tau_lw_gpt
private

◆ cld_tau_sw_bnd

real3d_k Radiation::cld_tau_sw_bnd
private

◆ cld_tau_sw_gpt

real3d_k Radiation::cld_tau_sw_gpt
private

◆ cldfrac_tot

real2d_k Radiation::cldfrac_tot
private

◆ d_tint

real2d_k Radiation::d_tint
private

◆ datalog_mf

amrex::MultiFab Radiation::datalog_mf
private

◆ eff_radius_qc

real2d_k Radiation::eff_radius_qc
private

◆ eff_radius_qi

real2d_k Radiation::eff_radius_qi
private

◆ emis_sfc

real2d_k Radiation::emis_sfc
private

◆ gas_names_offset

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

◆ iwp

real2d_k Radiation::iwp
private

◆ lat

real1d_k Radiation::lat
private

◆ lon

real1d_k Radiation::lon
private

◆ lw_bnd_flux_dn

real3d_k Radiation::lw_bnd_flux_dn
private

◆ lw_bnd_flux_up

real3d_k Radiation::lw_bnd_flux_up
private

◆ lw_clnclrsky_flux_dn

real2d_k Radiation::lw_clnclrsky_flux_dn
private

◆ lw_clnclrsky_flux_up

real2d_k Radiation::lw_clnclrsky_flux_up
private

◆ lw_clnsky_flux_dn

real2d_k Radiation::lw_clnsky_flux_dn
private

◆ lw_clnsky_flux_up

real2d_k Radiation::lw_clnsky_flux_up
private

◆ lw_clrsky_flux_dn

real2d_k Radiation::lw_clrsky_flux_dn
private

◆ lw_clrsky_flux_up

real2d_k Radiation::lw_clrsky_flux_up
private

◆ lw_clrsky_heating

real2d_k Radiation::lw_clrsky_heating
private

◆ lw_flux_dn

real2d_k Radiation::lw_flux_dn
private

◆ lw_flux_up

real2d_k Radiation::lw_flux_up
private

◆ lw_heating

real2d_k Radiation::lw_heating
private

◆ lw_src

real1d_k Radiation::lw_src
private

◆ lwp

real2d_k Radiation::lwp
private

◆ m_ba

amrex::BoxArray Radiation::m_ba
private

◆ m_ch4vmr

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

◆ m_co2vmr

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

◆ m_col_offsets

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

Referenced by Init().

◆ m_cons_in

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

◆ m_covmr

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

◆ m_do_aerosol_rad

bool Radiation::m_do_aerosol_rad = true
private

◆ m_do_subcol_sampling

bool Radiation::m_do_subcol_sampling = true
private

◆ m_dt

amrex::Real Radiation::m_dt
private

◆ m_extra_clnclrsky_diag

bool Radiation::m_extra_clnclrsky_diag = false
private

◆ m_extra_clnsky_diag

bool Radiation::m_extra_clnsky_diag = false
private

◆ m_first_step

bool Radiation::m_first_step = true
private

◆ m_fixed_orbital_year

bool Radiation::m_fixed_orbital_year = false
private

◆ m_fixed_solar_zenith_angle

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

◆ m_fixed_total_solar_irradiance

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

◆ m_gas_concs

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

◆ m_gas_mol_weights

real1d_k Radiation::m_gas_mol_weights
private

◆ m_gas_names

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

◆ m_geom

amrex::Geometry Radiation::m_geom
private

◆ m_ice

bool Radiation::m_ice = false
private

◆ m_lat

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

◆ m_lat_cons

amrex::Real Radiation::m_lat_cons = 39.809860
private

◆ m_lev

int Radiation::m_lev
private

Referenced by rad_run_impl().

◆ m_lon

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

◆ m_lon_cons

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

◆ m_lsm

bool Radiation::m_lsm = false
private

◆ m_lsm_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:
= {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: