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

#include <ERF_RRTMGP.H>

Collaboration diagram for Rrtmgp:

Public Member Functions

 Rrtmgp ()=default
 
 ~Rrtmgp ()=default
 
void initialize (int num_gas, const std::vector< std::string > &active_gas_names, const char *rrtmgp_coefficients_file_sw, const char *rrtmgp_coefficients_file_lw)
 
void finalize ()
 
void run_shortwave_rrtmgp (int ngas, int ncol, int nlay, const real3d &gas_vmr, const real2d &pmid, const real2d &tmid, const real2d &pint, const real1d &coszrs, const real2d &albedo_dir, const real2d &albedo_dif, const real3d &cld_tau_gpt, const real3d &cld_ssa_gpt, const real3d &cld_asm_gpt, const real3d &aer_tau_bnd, const real3d &aer_ssa_bnd, const real3d &aer_asm_bnd, const real2d &allsky_flux_up, const real2d &allsky_flux_dn, const real2d &allsky_flux_net, const real2d &allsky_flux_dn_dir, const real3d &allsky_bnd_flux_up, const real3d &allsky_bnd_flux_dn, const real3d &allsky_bnd_flux_net, const real3d &allsky_bnd_flux_dn_dir, const real2d &clrsky_flux_up, const real2d &clrsky_flux_dn, const real2d &clrsky_flux_net, const real2d &clrsky_flux_dn_dir, const real3d &clrsky_bnd_flux_up, const real3d &clrsky_bnd_flux_dn, const real3d &clrsky_bnd_flux_net, const real3d &clrsky_bnd_flux_dn_dir, double tsi_scaling)
 
void run_longwave_rrtmgp (int ngas, int ncol, int nlay, const real3d &gas_vmr, const real2d &pmid, const real2d &tmid, const real2d &pint, const real2d &tint, const real2d &emis_sfc, const real3d &cld_tau_gpt, const real3d &aer_tau_bnd, const real2d &allsky_flux_up, const real2d &allsky_flux_dn, const real2d &allsky_flux_net, const real3d &allsky_bnd_flux_up, const real3d &allsky_bnd_flux_dn, const real3d &allsky_bnd_flux_net, const real2d &clrsky_flux_up, const real2d &clrsky_flux_dn, const real2d &clrsky_flux_net, const real3d &clrsky_bnd_flux_up, const real3d &clrsky_bnd_flux_dn, const real3d &clrsky_bnd_flux_net)
 
int get_nband_sw ()
 
int get_nband_lw ()
 
int get_ngpt_sw ()
 
int get_ngpt_lw ()
 
double get_min_temperature ()
 
double get_max_temperature ()
 
void get_gpoint_bands_sw (int1d &gpoint_bands)
 
void get_gpoint_bands_lw (int1d &gpoint_bands)
 

Private Attributes

int ngas
 
string1d active_gases
 
char const * coefficients_file_sw
 
char const * coefficients_file_lw
 
GasOpticsRRTMGP k_dist_sw
 
GasOpticsRRTMGP k_dist_lw
 

Constructor & Destructor Documentation

◆ Rrtmgp()

Rrtmgp::Rrtmgp ( )
default

◆ ~Rrtmgp()

Rrtmgp::~Rrtmgp ( )
default

Member Function Documentation

◆ finalize()

void Rrtmgp::finalize ( )
14 {
15  k_dist_sw.finalize();
16  k_dist_lw.finalize();
17  yakl::finalize();
18 }
GasOpticsRRTMGP k_dist_lw
Definition: ERF_RRTMGP.H:122
GasOpticsRRTMGP k_dist_sw
Definition: ERF_RRTMGP.H:121

◆ get_gpoint_bands_lw()

void Rrtmgp::get_gpoint_bands_lw ( int1d &  gpoint_bands)
inline
106  {
107  gpoint_bands = k_dist_lw.get_gpoint_bands();
108  yakl::fence();
109  }

◆ get_gpoint_bands_sw()

void Rrtmgp::get_gpoint_bands_sw ( int1d &  gpoint_bands)
inline
101  {
102  gpoint_bands = k_dist_sw.get_gpoint_bands();
103  yakl::fence();
104  }

◆ get_max_temperature()

double Rrtmgp::get_max_temperature ( )
inline
97  {
98  return std::max(k_dist_sw.temp_ref_max, k_dist_lw.temp_ref_max);
99  }

◆ get_min_temperature()

double Rrtmgp::get_min_temperature ( )
inline
93  {
94  return std::min(k_dist_sw.temp_ref_min, k_dist_lw.temp_ref_min);
95  }

◆ get_nband_lw()

int Rrtmgp::get_nband_lw ( )
inline
81  {
82  return k_dist_lw.get_nband();
83  }

◆ get_nband_sw()

int Rrtmgp::get_nband_sw ( )
inline
77  {
78  return k_dist_sw.get_nband();
79  }

◆ get_ngpt_lw()

int Rrtmgp::get_ngpt_lw ( )
inline
89  {
90  return k_dist_lw.get_ngpt();
91  }

◆ get_ngpt_sw()

int Rrtmgp::get_ngpt_sw ( )
inline
85  {
86  return k_dist_sw.get_ngpt();
87  }

◆ initialize()

void Rrtmgp::initialize ( int  num_gas,
const std::vector< std::string > &  active_gas_names,
const char *  rrtmgp_coefficients_file_sw,
const char *  rrtmgp_coefficients_file_lw 
)
16 {
17  // Read gas optics coefficients from file
18  // Need to initialize available_gases here! The only field of the
19  // available_gases type that is used in the kdist initialize is
20  // available_gases%gas_name, which gives the name of each gas that would be
21  // present in the ty_gas_concs object. So, we can just set this here, rather
22  // than trying to fully populate the ty_gas_concs object here, which would be
23  // impossible from this initialization routine because I do not think the
24  // rad_cnst objects are setup yet.
25  // the other tasks!
26  ngas = num_gas;
27  coefficients_file_sw = rrtmgp_coefficients_file_sw;
28  coefficients_file_lw = rrtmgp_coefficients_file_lw;
29 
30  active_gases = string1d("active_gases", ngas);
31  for (int igas=0; igas<ngas; igas++)
32  active_gases(igas+1) = active_gas_names[igas];
33 
34  GasConcs available_gases;
35  available_gases.init(active_gases, 1, 1);
36  load_and_init(k_dist_sw, coefficients_file_sw, available_gases);
37  load_and_init(k_dist_lw, coefficients_file_lw, available_gases);
38 }
char const * coefficients_file_lw
Definition: ERF_RRTMGP.H:118
string1d active_gases
Definition: ERF_RRTMGP.H:114
int ngas
Definition: ERF_RRTMGP.H:113
char const * coefficients_file_sw
Definition: ERF_RRTMGP.H:117

◆ run_longwave_rrtmgp()

void Rrtmgp::run_longwave_rrtmgp ( int  ngas,
int  ncol,
int  nlay,
const real3d &  gas_vmr,
const real2d &  pmid,
const real2d &  tmid,
const real2d &  pint,
const real2d &  tint,
const real2d &  emis_sfc,
const real3d &  cld_tau_gpt,
const real3d &  aer_tau_bnd,
const real2d &  allsky_flux_up,
const real2d &  allsky_flux_dn,
const real2d &  allsky_flux_net,
const real3d &  allsky_bnd_flux_up,
const real3d &  allsky_bnd_flux_dn,
const real3d &  allsky_bnd_flux_net,
const real2d &  clrsky_flux_up,
const real2d &  clrsky_flux_dn,
const real2d &  clrsky_flux_net,
const real3d &  clrsky_bnd_flux_up,
const real3d &  clrsky_bnd_flux_dn,
const real3d &  clrsky_bnd_flux_net 
)
19 {
20  // Wrap pointers in YAKL arrays
21  int nlwbands = k_dist_lw.get_nband();
22  int nlwgpts = k_dist_lw.get_ngpt();
23 
24  // Populate gas concentrations
25  GasConcs gas_concs;
26  gas_concs.init(active_gases, ncol, nlay);
27  real2d tmp2d;
28  tmp2d = real2d("tmp", ncol, nlay);
29  for (int igas = 1; igas <= ngas; igas++) {
30  parallel_for(SimpleBounds<2>(nlay,ncol), YAKL_LAMBDA(int ilay, int icol)
31  {
32  tmp2d(icol,ilay) = gas_vmr(igas,icol,ilay);
33  });
34  gas_concs.set_vmr(active_gases(igas), tmp2d);
35  }
36 
37  // Boundary conditions
38  SourceFuncLW lw_sources;
39  lw_sources.alloc(ncol, nlay, k_dist_lw);
40 
41  // Weights and angle secants for first order (k=1) Gaussian quadrature.
42  // Values from Table 2, Clough et al, 1992, doi:10.1029/92JD01419
43  // after Abramowitz & Stegun 1972, page 921
44  int constexpr max_gauss_pts = 4;
45  realHost2d gauss_Ds_host ("gauss_Ds" ,max_gauss_pts,max_gauss_pts);
46  gauss_Ds_host(1,1) = 1.66_wp ; gauss_Ds_host(2,1) = 0._wp; gauss_Ds_host(3,1) = 0._wp; gauss_Ds_host(4,1) = 0._wp;
47  gauss_Ds_host(1,2) = 1.18350343_wp; gauss_Ds_host(2,2) = 2.81649655_wp; gauss_Ds_host(3,2) = 0._wp; gauss_Ds_host(4,2) = 0._wp;
48  gauss_Ds_host(1,3) = 1.09719858_wp; gauss_Ds_host(2,3) = 1.69338507_wp; gauss_Ds_host(3,3) = 4.70941630_wp; gauss_Ds_host(4,3) = 0._wp;
49  gauss_Ds_host(1,4) = 1.06056257_wp; gauss_Ds_host(2,4) = 1.38282560_wp; gauss_Ds_host(3,4) = 2.40148179_wp; gauss_Ds_host(4,4) = 7.15513024_wp;
50 
51  realHost2d gauss_wts_host("gauss_wts",max_gauss_pts,max_gauss_pts);
52  gauss_wts_host(1,1) = 0.5_wp ; gauss_wts_host(2,1) = 0._wp ; gauss_wts_host(3,1) = 0._wp ; gauss_wts_host(4,1) = 0._wp ;
53  gauss_wts_host(1,2) = 0.3180413817_wp; gauss_wts_host(2,2) = 0.1819586183_wp; gauss_wts_host(3,2) = 0._wp ; gauss_wts_host(4,2) = 0._wp ;
54  gauss_wts_host(1,3) = 0.2009319137_wp; gauss_wts_host(2,3) = 0.2292411064_wp; gauss_wts_host(3,3) = 0.0698269799_wp; gauss_wts_host(4,3) = 0._wp ;
55  gauss_wts_host(1,4) = 0.1355069134_wp; gauss_wts_host(2,4) = 0.2034645680_wp; gauss_wts_host(3,4) = 0.1298475476_wp; gauss_wts_host(4,4) = 0.0311809710_wp;
56 
57  real2d gauss_Ds ("gauss_Ds" ,max_gauss_pts,max_gauss_pts);
58  real2d gauss_wts("gauss_wts",max_gauss_pts,max_gauss_pts);
59  gauss_Ds_host .deep_copy_to(gauss_Ds );
60  gauss_wts_host.deep_copy_to(gauss_wts);
61 
62  // Populate optical property objects
63  OpticalProps1scl combined_optics;
64  combined_optics.alloc_1scl(ncol, nlay, k_dist_lw);
65  bool1d top_at_1_g("top_at_1_g",1);
66  boolHost1d top_at_1_h("top_at_1_h",1);
67  bool top_at_1;
68  real1d t_sfc("t_sfc", ncol);
69  parallel_for(SimpleBounds<1>(ncol), YAKL_LAMBDA (int icol)
70  {
71  t_sfc(icol) = tint(icol,nlay+1);
72  top_at_1_g(1) = pmid(1, 1) < pmid (1, nlay);
73  });
74  top_at_1_g.deep_copy_to(top_at_1_h);
75  top_at_1 = top_at_1_h(1);
76  k_dist_lw.gas_optics(ncol, nlay, top_at_1, pmid, pint, tmid, t_sfc, gas_concs, combined_optics, lw_sources, real2d(), tint);
77 
78  // Add in aerosol; we can define this by bands or gpoints. If we define by
79  // bands, then internally when increment() is called it will map these to
80  // gpoints. Not sure if there is a beneift one way or another.
81  OpticalProps1scl aerosol_optics;
82  auto &aerosol_optics_tau = aerosol_optics.tau;
83  if (false) {
84  aerosol_optics.alloc_1scl(ncol, nlay, k_dist_lw);
85  auto gpt_bnd = aerosol_optics.get_gpoint_bands();
86  parallel_for(SimpleBounds<3>(nlwgpts,nlay,ncol) , YAKL_LAMBDA (int igpt, int ilay, int icol)
87  {
88  aerosol_optics_tau(icol,ilay,igpt) = aer_tau_bnd(icol,ilay,gpt_bnd(igpt));
89  });
90  } else {
91  aerosol_optics.alloc_1scl(ncol, nlay, k_dist_lw.get_band_lims_wavenumber());
92  parallel_for(SimpleBounds<3>(nlwbands,nlay,ncol), YAKL_LAMBDA (int ibnd, int ilay, int icol)
93  {
94  aerosol_optics_tau(icol,ilay,ibnd) = aer_tau_bnd(icol,ilay,ibnd);
95  });
96  }
97  aerosol_optics.increment(combined_optics);
98 
99  // Do the clearsky calculation before adding in clouds
100  FluxesByband fluxes_clrsky;
101  fluxes_clrsky.flux_up = real2d("clrsky_flux_up" , ncol, nlay+1); // clrsky_flux_up;
102  fluxes_clrsky.flux_dn = real2d("clrsky_flux_dn" , ncol, nlay+1); //clrsky_flux_dn;
103  fluxes_clrsky.flux_net = real2d("clrsky_flux_net", ncol, nlay+1); //clrsky_flux_net;
104  fluxes_clrsky.bnd_flux_up = real3d("clrsky_bnd_flux_up" , ncol, nlay+1, nlwbands); //clrsky_bnd_flux_up;
105  fluxes_clrsky.bnd_flux_dn = real3d("clrsky_bnd_flux_dn" , ncol, nlay+1, nlwbands); //clrsky_bnd_flux_dn;
106  fluxes_clrsky.bnd_flux_net = real3d("clrsky_bnd_flux_net", ncol, nlay+1, nlwbands); //clrsky_bnd_flux_net;
107 
108  rte_lw(max_gauss_pts, gauss_Ds, gauss_wts, combined_optics, top_at_1, lw_sources, emis_sfc, fluxes_clrsky);
109 
110  // Copy fluxes back out of FluxesByband object
111  fluxes_clrsky.flux_up.deep_copy_to(clrsky_flux_up);
112  fluxes_clrsky.flux_dn.deep_copy_to(clrsky_flux_dn);
113  fluxes_clrsky.flux_net.deep_copy_to(clrsky_flux_net);
114  fluxes_clrsky.bnd_flux_up.deep_copy_to(clrsky_bnd_flux_up);
115  fluxes_clrsky.bnd_flux_dn.deep_copy_to(clrsky_bnd_flux_dn);
116  fluxes_clrsky.bnd_flux_net.deep_copy_to(clrsky_bnd_flux_net);
117 
118  // Add in clouds
119  OpticalProps1scl cloud_optics;
120  cloud_optics.alloc_1scl(ncol, nlay, k_dist_lw);
121  auto &cloud_optics_tau = cloud_optics.tau;
122  parallel_for(SimpleBounds<3>(nlwgpts,nlay,ncol) , YAKL_LAMBDA (int igpt, int ilay, int icol) {
123  cloud_optics_tau(icol,ilay,igpt) = cld_tau_gpt(icol,ilay,igpt);
124  });
125  cloud_optics.increment(combined_optics);
126 
127  // Call LW flux driver
128  FluxesByband fluxes_allsky;
129  fluxes_allsky.flux_up = real2d("allsky_flux_up" , ncol, nlay+1); //allsky_flux_up;
130  fluxes_allsky.flux_dn = real2d("allsky_flux_dn" , ncol, nlay+1); //allsky_flux_dn;
131  fluxes_allsky.flux_net = real2d("allsky_flux_net", ncol, nlay+1); //allsky_flux_net;
132  fluxes_allsky.bnd_flux_up = real3d("allsky_bnd_flux_up" , ncol, nlay+1, nlwbands); //allsky_bnd_flux_up;
133  fluxes_allsky.bnd_flux_dn = real3d("allsky_bnd_flux_dn" , ncol, nlay+1, nlwbands); //allsky_bnd_flux_dn;
134  fluxes_allsky.bnd_flux_net = real3d("allsky_bnd_flux_net", ncol, nlay+1, nlwbands); //allsky_bnd_flux_net;
135 
136  rte_lw(max_gauss_pts, gauss_Ds, gauss_wts, combined_optics, top_at_1, lw_sources, emis_sfc, fluxes_allsky);
137 
138  // Copy fluxes back out of FluxesByband object
139  fluxes_allsky.flux_up.deep_copy_to(allsky_flux_up);
140  fluxes_allsky.flux_dn.deep_copy_to(allsky_flux_dn);
141  fluxes_allsky.flux_net.deep_copy_to(allsky_flux_net);
142  fluxes_allsky.bnd_flux_up.deep_copy_to(allsky_bnd_flux_up);
143  fluxes_allsky.bnd_flux_dn.deep_copy_to(allsky_bnd_flux_dn);
144  fluxes_allsky.bnd_flux_net.deep_copy_to(allsky_bnd_flux_net);
145  yakl::fence();
146 }

◆ run_shortwave_rrtmgp()

void Rrtmgp::run_shortwave_rrtmgp ( int  ngas,
int  ncol,
int  nlay,
const real3d &  gas_vmr,
const real2d &  pmid,
const real2d &  tmid,
const real2d &  pint,
const real1d &  coszrs,
const real2d &  albedo_dir,
const real2d &  albedo_dif,
const real3d &  cld_tau_gpt,
const real3d &  cld_ssa_gpt,
const real3d &  cld_asm_gpt,
const real3d &  aer_tau_bnd,
const real3d &  aer_ssa_bnd,
const real3d &  aer_asm_bnd,
const real2d &  allsky_flux_up,
const real2d &  allsky_flux_dn,
const real2d &  allsky_flux_net,
const real2d &  allsky_flux_dn_dir,
const real3d &  allsky_bnd_flux_up,
const real3d &  allsky_bnd_flux_dn,
const real3d &  allsky_bnd_flux_net,
const real3d &  allsky_bnd_flux_dn_dir,
const real2d &  clrsky_flux_up,
const real2d &  clrsky_flux_dn,
const real2d &  clrsky_flux_net,
const real2d &  clrsky_flux_dn_dir,
const real3d &  clrsky_bnd_flux_up,
const real3d &  clrsky_bnd_flux_dn,
const real3d &  clrsky_bnd_flux_net,
const real3d &  clrsky_bnd_flux_dn_dir,
double  tsi_scaling 
)
14 {
15  // Wrap pointers in YAKL arrays
16  int nswbands = k_dist_sw.get_nband();
17  int nswgpts = k_dist_sw.get_ngpt();
18 
19  // Populate gas concentrations object
20  GasConcs gas_concs;
21  gas_concs.init(active_gases, ncol, nlay);
22  real2d tmp2d;
23  tmp2d = real2d("tmp", ncol, nlay);
24  for (int igas = 1; igas <= ngas; igas++) {
25  parallel_for(SimpleBounds<2>(nlay,ncol), YAKL_LAMBDA(int ilay, int icol) {
26  tmp2d(icol,ilay) = gas_vmr(igas,icol,ilay);
27  });
28  gas_concs.set_vmr(active_gases(igas), tmp2d);
29  }
30 
31  // Do gas optics
32  // TODO: should we avoid allocating here?
33  OpticalProps2str combined_optics;
34  combined_optics.alloc_2str(ncol, nlay, k_dist_sw);
35  bool1d top_at_1_g("top_at_1_g",1);
36  boolHost1d top_at_1_h("top_at_1_h",1);
37  bool top_at_1;
38  parallel_for(SimpleBounds<1>(1), YAKL_LAMBDA (int ilay) { // HACK: Single loop kernel is not efficient
39  top_at_1_g(1) = pmid(1, 1) < pmid (1, nlay);
40  });
41  top_at_1_g.deep_copy_to(top_at_1_h);
42  top_at_1 = top_at_1_h(1);
43  real2d toa_flux("toa_flux", ncol, nswgpts);
44 
45 #ifdef AMREX_DEBUG
46  // print extra info from RRTMGP
47  k_dist_sw.print_norms();
48 
49  // check pressure is within bounds
50  real pint_min = yakl::intrinsics::minval<real>(pint);
51  real pint_max = yakl::intrinsics::maxval<real>(pint);
52  AMREX_ASSERT(pint_max < k_dist_sw.get_press_max());
53  AMREX_ASSERT(pint_min > k_dist_sw.get_press_min());
54 #endif
55 
56  k_dist_sw.gas_optics(ncol, nlay, top_at_1, pmid, pint, tmid, gas_concs, combined_optics, toa_flux);
57 
58  // Apply TOA flux scaling
59  parallel_for(SimpleBounds<2>(nswgpts,ncol), YAKL_LAMBDA (int igpt, int icol) {
60  toa_flux(icol, igpt) = tsi_scaling * toa_flux(icol,igpt);
61  });
62 
63  // Add in aerosol, allocate on gpts and map aer_*_bnd from bands
64  OpticalProps2str aerosol_optics;
65  aerosol_optics.alloc_2str(ncol, nlay, k_dist_sw);
66  auto gpt_bnd = aerosol_optics.get_gpoint_bands();
67  parallel_for(SimpleBounds<3>(nswgpts,nlay,ncol) , YAKL_LAMBDA (int igpt, int ilay, int icol)
68  {
69  aerosol_optics.tau(icol,ilay,igpt) = aer_tau_bnd(icol,ilay,gpt_bnd(igpt));
70  aerosol_optics.ssa(icol,ilay,igpt) = aer_ssa_bnd(icol,ilay,gpt_bnd(igpt));
71  aerosol_optics.g (icol,ilay,igpt) = aer_asm_bnd(icol,ilay,gpt_bnd(igpt));
72  });
73 
74  aerosol_optics.delta_scale();
75 
76  // NOTE: aero_optics is allocated on nswgpts and combined_optics is on nswgpts
77  // The `increment` call below can handle matching or differing (nswgpts/nswbnds)
78  // sizes; see calls in `mo_optical_props_kernels.cpp`. Since we have matching
79  // gpt sizes, we will call `increment_2stream_by_2stream`
80  aerosol_optics.increment(combined_optics);
81 
82  // Do the clearsky calculation before adding in clouds
83  FluxesByband fluxes_clrsky;
84  fluxes_clrsky.flux_up = real2d("clrsky_flux_up" , ncol, nlay+1); // clrsky_flux_up;
85  fluxes_clrsky.flux_dn = real2d("clrsky_flux_nd" , ncol, nlay+1); //clrsky_flux_dn;
86  fluxes_clrsky.flux_net = real2d("clrsky_flux_net", ncol, nlay+1); //clrsky_flux_net;
87  fluxes_clrsky.flux_dn_dir = real2d("clrsky_flux_dn_dir", ncol, nlay+1); //clrsky_flux_dn_dir;
88  fluxes_clrsky.bnd_flux_up = real3d("clrsky_bnd_flux_up" , ncol, nlay+1, nswbands); //clrsky_bnd_flux_up;
89  fluxes_clrsky.bnd_flux_dn = real3d("clrsky_bnd_flux_dn" , ncol, nlay+1, nswbands); //clrsky_bnd_flux_dn;
90  fluxes_clrsky.bnd_flux_net = real3d("clrsky_bnd_flux_net", ncol, nlay+1, nswbands); //clrsky_bnd_flux_net;
91  fluxes_clrsky.bnd_flux_dn_dir = real3d("clrsky_bnd_flux_dn_dir", ncol, nlay+1, nswbands); //clrsky_bnd_flux_dn_dir;
92 
93  rte_sw(combined_optics, top_at_1, coszrs, toa_flux, albedo_dir, albedo_dif, fluxes_clrsky);
94 
95  // Copy fluxes back out of FluxesByband object
96  fluxes_clrsky.flux_up.deep_copy_to (clrsky_flux_up);
97  fluxes_clrsky.flux_dn.deep_copy_to (clrsky_flux_dn);
98  fluxes_clrsky.flux_net.deep_copy_to(clrsky_flux_net);
99  fluxes_clrsky.flux_dn_dir.deep_copy_to(clrsky_flux_dn_dir);
100  fluxes_clrsky.bnd_flux_up.deep_copy_to (clrsky_bnd_flux_up);
101  fluxes_clrsky.bnd_flux_dn.deep_copy_to (clrsky_bnd_flux_dn);
102  fluxes_clrsky.bnd_flux_net.deep_copy_to(clrsky_bnd_flux_net);
103  fluxes_clrsky.bnd_flux_dn_dir.deep_copy_to(clrsky_bnd_flux_dn_dir);
104 
105  // Add in clouds, which are already on gpts
106  OpticalProps2str cloud_optics;
107  cloud_optics.alloc_2str(ncol, nlay, k_dist_sw);
108  auto &cloud_optics_tau = cloud_optics.tau;
109  auto &cloud_optics_ssa = cloud_optics.ssa;
110  auto &cloud_optics_g = cloud_optics.g ;
111  parallel_for(SimpleBounds<3>(nswgpts,nlay,ncol) , YAKL_LAMBDA (int igpt, int ilay, int icol) {
112  cloud_optics_tau(icol,ilay,igpt) = cld_tau_gpt(icol,ilay,igpt);
113  cloud_optics_ssa(icol,ilay,igpt) = cld_ssa_gpt(icol,ilay,igpt);
114  cloud_optics_g (icol,ilay,igpt) = cld_asm_gpt(icol,ilay,igpt);
115  });
116 
117  cloud_optics.delta_scale();
118 
119  // TODO: Sort through cloud optics, this increment does not
120  // modify the optical properties (tau/ssa/g) so we get
121  // the same fluxes and heating rate.
122 
123  // NOTE: See above and here we again have matching gpt sizes
124  cloud_optics.increment(combined_optics);
125 
126  // Call SW flux driver
127  FluxesByband fluxes_allsky;
128  fluxes_allsky.flux_up = real2d("allsky_flux_up" , ncol, nlay+1); //allsky_flux_up;
129  fluxes_allsky.flux_dn = real2d("allsky_flux_dn" , ncol, nlay+1); //allsky_flux_dn;
130  fluxes_allsky.flux_net = real2d("allsky_flux_net", ncol, nlay+1); //allsky_flux_net;
131  fluxes_allsky.flux_dn_dir = real2d("allsky_flux_dn_dir", ncol, nlay+1); //allsky_flux_dn_dir;
132  fluxes_allsky.bnd_flux_up = real3d("allsky_bnd_flux_up" , ncol, nlay+1, nswbands); //allsky_bnd_flux_up;
133  fluxes_allsky.bnd_flux_dn = real3d("allsky_bnd_flux_dn" , ncol, nlay+1, nswbands); //allsky_bnd_flux_dn;
134  fluxes_allsky.bnd_flux_net = real3d("allsky_bnd_flux_net", ncol, nlay+1, nswbands); //allsky_bnd_flux_net;
135  fluxes_allsky.bnd_flux_dn_dir = real3d("allsky_bnd_flux_dn_dir", ncol, nlay+1, nswbands); //allsky_bnd_flux_dn_dir;
136 
137  rte_sw(combined_optics, top_at_1, coszrs, toa_flux, albedo_dir, albedo_dif, fluxes_allsky);
138 
139  // Copy fluxes back out of FluxesByband object
140  fluxes_allsky.flux_up.deep_copy_to(allsky_flux_up);
141  fluxes_allsky.flux_dn.deep_copy_to(allsky_flux_dn);
142  fluxes_allsky.flux_net.deep_copy_to(allsky_flux_net);
143  fluxes_allsky.flux_dn_dir.deep_copy_to(allsky_flux_dn_dir);
144  fluxes_allsky.bnd_flux_up.deep_copy_to(allsky_bnd_flux_up);
145  fluxes_allsky.bnd_flux_dn.deep_copy_to(allsky_bnd_flux_dn);
146  fluxes_allsky.bnd_flux_net.deep_copy_to(allsky_bnd_flux_net);
147  fluxes_allsky.bnd_flux_dn_dir.deep_copy_to(allsky_bnd_flux_dn_dir);
148  yakl::fence();
149 }

Member Data Documentation

◆ active_gases

string1d Rrtmgp::active_gases
private

◆ coefficients_file_lw

char const* Rrtmgp::coefficients_file_lw
private

Referenced by initialize().

◆ coefficients_file_sw

char const* Rrtmgp::coefficients_file_sw
private

Referenced by initialize().

◆ k_dist_lw

◆ k_dist_sw

◆ ngas

int Rrtmgp::ngas
private

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