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

#include <ERF_Optics.H>

Collaboration diagram for Optics:

Public Member Functions

 Optics ()=default
 
 Optics (int ngases, char *gas_names[])
 
 ~Optics ()=default
 
void get_cloud_optics_sw (int ncol, int nlev, int nbnd, bool do_snow, const real2d &cld, const real2d &cldfsnow, const real2d &iclwp, const real2d &iciwp, const real2d &icswp, const real2d &lambdac, const real2d &mu, const real2d &dei, const real2d &des, const real2d &rel, const real2d &rei, const real3d &tau_out, const real3d &ssa_out, const real3d &asm_out, const real3d &liq_tau_out, const real3d &ice_tau_out, const real3d &snw_tau_out)
 
void get_cloud_optics_lw (int ncol, int nlev, int nbnd, bool do_snow, const real2d &cld, const real2d &cldfsnow, const real2d &iclwp, const real2d &iciwp, const real2d &icswp, const real2d &lambdac, const real2d &mu, const real2d &dei, const real2d &des, const real2d &rei, const real3d &tau_out, const real3d &liq_tau_out, const real3d &ice_tau_out, const real3d &snw_tau_out)
 
void sample_cloud_optics_sw (int ncol, int nlev, int ngpt, const int1d &gpt2bnd, const real2d &pmid, const real2d &cld, const real2d &cldfsnow, const real3d &tau_bnd, const real3d &ssa_bnd, const real3d &asm_bnd, const real3d &tau_gpt, const real3d &ssa_gpt, const real3d &asm_gpt)
 
void sample_cloud_optics_lw (int ncol, int nlev, int ngpt, const int1d &gpt2bnd, const real2d &pmid, const real2d &cld, const real2d &cldfsnow, const real3d &tau_bnd, const real3d &tau_gpt)
 
void set_aerosol_optics_sw (int icall, int ncol, int nlev, int nswbands, real dt, const int1d &night_indices, bool is_cmip6_volc, const real3d &tau_out, const real3d &ssa_out, const real3d &asm_out, const real2d &clear_rh)
 
void set_aerosol_optics_lw (int icall, real dt, bool is_cmip6_volc, const real2d &zi, const real3d &tau, const real2d &clear_rh)
 
void mcica_subcol_mask (int ngpt, int ncol, int nlev, const real2d &cldfrac, const bool3d &iscloudy)
 
void combine_properties (int nbands, int ncols, int nlevs, const real2d &fraction1, const real3d &property1, const real2d &fraction2, const real3d &property2, const real3d &combined_property)
 
void initialize (int ngas, int nmodes, int num_aeros, int nswbands, int nlwbands, int ncol, int nlev, int nrh, int top_lev, const std::vector< std::string > &aero_names, const real2d &zi, const real2d &pmid, const real2d &pdel, const real2d &temp, const real2d &qi, const real2d &geom_radius)
 
void finalize ()
 

Private Attributes

int ngas
 
char ** gas_names
 
std::string icecldoptics
 
std::string liqcldoptics
 
CloudRadProps cloud_optics
 
AerRadProps aero_optics
 

Constructor & Destructor Documentation

◆ Optics() [1/2]

Optics::Optics ( )
default

◆ Optics() [2/2]

Optics::Optics ( int  ngases,
char *  gas_names[] 
)
explicit

◆ ~Optics()

Optics::~Optics ( )
default

Member Function Documentation

◆ combine_properties()

void Optics::combine_properties ( int  nbands,
int  ncols,
int  nlevs,
const real2d &  fraction1,
const real3d &  property1,
const real2d &  fraction2,
const real3d &  property2,
const real3d &  combined_property 
)
260 {
261  // Combined fraction (max of property 1 and 2)
262  real2d combined_fraction("combined_fraction", ncols,nlevs);
263 
264  // Combined fraction
265  parallel_for(SimpleBounds<2>(ncols, nlevs), YAKL_LAMBDA (int icol, int ilev)
266  {
267  combined_fraction(icol, ilev) = std::max(fraction1(icol, ilev), fraction2(icol, ilev));
268  });
269 
270  // Combine optical properties by weighting by amount of cloud and snow
271  yakl::memset(combined_property, 0.);
272  parallel_for(SimpleBounds<3>(nlevs, ncols, nbands), YAKL_LAMBDA (int ilev, int icol, int iband)
273  {
274  if (combined_fraction(icol,ilev) > 0) {
275  combined_property(iband,icol,ilev) = (
276  fraction1(icol,ilev) * property1(iband,icol,ilev)
277  + fraction2(icol,ilev) * property2(iband,icol,ilev)
278  ) / combined_fraction(icol,ilev);
279  }
280  });
281 }

Referenced by get_cloud_optics_lw(), and get_cloud_optics_sw().

Here is the caller graph for this function:

◆ finalize()

void Optics::finalize ( )
27 {
28  // nothing needs to be done here
29 }

◆ get_cloud_optics_lw()

void Optics::get_cloud_optics_lw ( int  ncol,
int  nlev,
int  nbnd,
bool  do_snow,
const real2d &  cld,
const real2d &  cldfsnow,
const real2d &  iclwp,
const real2d &  iciwp,
const real2d &  icswp,
const real2d &  lambdac,
const real2d &  mu,
const real2d &  dei,
const real2d &  des,
const real2d &  rei,
const real3d &  tau_out,
const real3d &  liq_tau_out,
const real3d &  ice_tau_out,
const real3d &  snw_tau_out 
)
185 {
186  // Temporary variables to hold absorption optical depth
187  real3d ice_tau("ice_tau", nbnd, ncol, nlev);
188  real3d liq_tau("liq_tau", nbnd, ncol, nlev);
189  real3d snow_tau("snow_tau", nbnd, ncol, nlev);
190  real3d cld_tau("cld_tau", nbnd, ncol, nlev);
191  real3d combined_tau("combined_tau", nbnd, ncol, nlev);
192 
193  // Initialize outputs
194  yakl::memset(tau_out, 0.);
195  yakl::memset(liq_tau_out, 0.);
196  yakl::memset(ice_tau_out, 0.);
197  yakl::memset(snw_tau_out, 0.);
198 
199  // Initialize local variables
200  yakl::memset(ice_tau, 0.);
201  yakl::memset(liq_tau, 0.);
202  yakl::memset(snow_tau, 0.);
203  yakl::memset(cld_tau, 0.);
204  yakl::memset(combined_tau, 0.);
205 
206  // Get ice optics
207  if (icecldoptics == "mitchell") {
208  cloud_optics.mitchell_ice_optics_lw(ncol, nlev, iciwp, dei, ice_tau);
209  }
210  else if (icecldoptics == "ebertcurry") {
211  EbertCurry::ec_ice_optics_lw(ncol, nlev, nbnd, cld, iclwp, iciwp, rei, ice_tau);
212  }
213 
214  // Get liquid optics
215  if (liqcldoptics == "gammadist") {
216  cloud_optics.gammadist_liq_optics_lw(ncol, nlev, iclwp, lambdac, mu, liq_tau);
217  }
218  else if (liqcldoptics == "slingo") {
219  Slingo::slingo_liq_optics_lw(ncol, nlev, nbnd, cld, iclwp, iciwp, liq_tau);
220  }
221 
222  // Combined cloud optics
223  parallel_for(SimpleBounds<3>(nbnd, ncol, nlev), YAKL_LAMBDA (int ibnd, int icol, int ilev)
224  {
225  cld_tau(ibnd, icol, ilev) = liq_tau(ibnd, icol, ilev) + ice_tau(ibnd, icol, ilev);
226  });
227 
228  // Get snow optics?
229  if (do_snow) {
230  cloud_optics.mitchell_ice_optics_lw(ncol, nlev, icswp, des, snow_tau);
231  combine_properties(nbnd, ncol, nlev,
232  cld, cld_tau, cldfsnow, snow_tau, combined_tau);
233  }
234  else {
235  parallel_for(SimpleBounds<3>(nbnd, ncol, nlev), YAKL_LAMBDA (int ibnd, int icol, int ilev)
236  {
237  combined_tau(ibnd,icol,ilev) = cld_tau(ibnd,icol,ilev);
238  });
239  }
240 
241  // Set output optics
242  parallel_for(SimpleBounds<3>(nbnd, ncol, nlev), YAKL_LAMBDA (int ibnd, int icol, int ilev)
243  {
244  tau_out(icol,ilev,ibnd) = combined_tau(ibnd,icol,ilev);
245  liq_tau_out(icol,ilev,ibnd) = liq_tau(ibnd,icol,ilev);
246  ice_tau_out(icol,ilev,ibnd) = ice_tau(ibnd,icol,ilev);
247  snw_tau_out(icol,ilev,ibnd) = snow_tau(ibnd,icol,ilev);
248  });
249 }
void mitchell_ice_optics_lw(const int &ncol, const int &nlev, const real2d &iciwpth, const real2d &dei, real3d &abs_od)
Definition: ERF_CloudRadProps.cpp:215
void gammadist_liq_optics_lw(const int &ncol, const int &nlev, const real2d &iclwpth, const real2d &lamc, const real2d &pgam, real3d &abs_od)
Definition: ERF_CloudRadProps.cpp:132
static void ec_ice_optics_lw(int ncol, int nlev, int nlwbands, const real2d &cldn, const real2d &iclwpth, const real2d &iciwpth, const real2d &rei, const real3d &abs_od)
Definition: ERF_EbertCurry.H:116
std::string icecldoptics
Definition: ERF_Optics.H:91
void combine_properties(int nbands, int ncols, int nlevs, const real2d &fraction1, const real3d &property1, const real2d &fraction2, const real3d &property2, const real3d &combined_property)
Definition: ERF_Optics.cpp:256
std::string liqcldoptics
Definition: ERF_Optics.H:92
CloudRadProps cloud_optics
Definition: ERF_Optics.H:94
static void slingo_liq_optics_lw(int ncol, int nlev, int nlwbands, const real2d &cldn, const real2d &iclwpth, const real2d &iciwpth, const real3d &abs_od)
Definition: ERF_Slingo.H:126
Here is the call graph for this function:

◆ get_cloud_optics_sw()

void Optics::get_cloud_optics_sw ( int  ncol,
int  nlev,
int  nbnd,
bool  do_snow,
const real2d &  cld,
const real2d &  cldfsnow,
const real2d &  iclwp,
const real2d &  iciwp,
const real2d &  icswp,
const real2d &  lambdac,
const real2d &  mu,
const real2d &  dei,
const real2d &  des,
const real2d &  rel,
const real2d &  rei,
const real3d &  tau_out,
const real3d &  ssa_out,
const real3d &  asm_out,
const real3d &  liq_tau_out,
const real3d &  ice_tau_out,
const real3d &  snw_tau_out 
)
38 {
39  //Temporary variables to hold cloud optical properties before combining into
40  // output arrays. Same shape as output arrays, so get shapes from output.
41  real3d liq_tau("liq_tau",nbnd, ncol, nlev);
42  real3d liq_tau_ssa("liq_tau_ssa",nbnd, ncol, nlev);
43  real3d liq_tau_ssa_g("liq_tau_ssa_g", nbnd, ncol, nlev);
44  real3d liq_tau_ssa_f("liq_tau_ssa_f", nbnd, ncol, nlev);
45  real3d ice_tau("ice_tau", nbnd, ncol, nlev);
46  real3d ice_tau_ssa("ice_tau_ssa", nbnd, ncol, nlev);
47  real3d ice_tau_ssa_g("ice_tau_ssa_g", nbnd, ncol, nlev);
48  real3d ice_tau_ssa_f("ice_tau_ssa_f", nbnd, ncol, nlev);
49  real3d cld_tau("cld_tau", nbnd, ncol, nlev);
50  real3d cld_tau_ssa("cld_tau_ssa", nbnd, ncol, nlev);
51  real3d cld_tau_ssa_g("cld_tau_ssa_g", nbnd, ncol, nlev);
52  real3d cld_tau_ssa_f("cld_tau_ssa_f", nbnd, ncol, nlev);
53  real3d snow_tau("snow_tau", nbnd, ncol, nlev );
54  real3d snow_tau_ssa("snow_tau_ssa", nbnd, ncol, nlev);
55  real3d snow_tau_ssa_g("snow_tau_ssa_g", nbnd, ncol, nlev);
56  real3d snow_tau_ssa_f("snow_tau_ssa_f", nbnd, ncol, nlev);
57  real3d combined_tau("combined_tau", nbnd, ncol, nlev);
58  real3d combined_tau_ssa("combined_tau_ssa", nbnd, ncol, nlev);
59  real3d combined_tau_ssa_g("combined_tau_ssa_g", nbnd, ncol, nlev);
60  real3d combined_tau_ssa_f("combined_tau_ssa_f", nbnd, ncol, nlev);
61 
62  // Initialize outputs
63  yakl::memset(tau_out, 0.);
64  yakl::memset(ssa_out, 0.);
65  yakl::memset(asm_out, 0.);
66  yakl::memset(liq_tau_out, 0.);
67  yakl::memset(ice_tau_out, 0.);
68  yakl::memset(snw_tau_out, 0.);
69 
70  // Initialize local variables
71  yakl::memset(ice_tau, 0.);
72  yakl::memset(ice_tau_ssa, 0.);
73  yakl::memset(ice_tau_ssa_g, 0.);
74  yakl::memset(ice_tau_ssa_f, 0.);
75  yakl::memset(liq_tau, 0.);
76  yakl::memset(liq_tau_ssa, 0.);
77  yakl::memset(liq_tau_ssa_g, 0.);
78  yakl::memset(liq_tau_ssa_f, 0.);
79  yakl::memset(snow_tau, 0.);
80  yakl::memset(snow_tau_ssa, 0.);
81  yakl::memset(snow_tau_ssa_g, 0.);
82  yakl::memset(snow_tau_ssa_f, 0.);
83  yakl::memset(combined_tau, 0.);
84  yakl::memset(combined_tau_ssa, 0.);
85  yakl::memset(combined_tau_ssa_g, 0.);
86  yakl::memset(combined_tau_ssa_f, 0.);
87 
88  // Get ice cloud optics
89  if (icecldoptics == "mitchell") {
90  cloud_optics.mitchell_ice_optics_sw(ncol, nlev, iciwp, dei,
91  ice_tau, ice_tau_ssa,
92  ice_tau_ssa_g, ice_tau_ssa_f);
93  }
94  else if (icecldoptics == "ebertcurry") {
95  EbertCurry::ec_ice_optics_sw(ncol, nlev, nbnd, cld, iciwp,
96  rei, ice_tau, ice_tau_ssa,
97  ice_tau_ssa_g, ice_tau_ssa_f);
98  }
99 
100  // Get liquid cloud optics
101  if (liqcldoptics == "gammadist") {
102  cloud_optics.gammadist_liq_optics_sw(ncol, nlev, iclwp, lambdac, mu,
103  liq_tau, liq_tau_ssa,
104  liq_tau_ssa_g, liq_tau_ssa_f);
105  }
106  else if (liqcldoptics == "slingo") {
107  Slingo::slingo_liq_optics_sw(ncol, nlev, nbnd, cld, iclwp, rel,
108  liq_tau, liq_tau_ssa,
109  liq_tau_ssa_g, liq_tau_ssa_f);
110  }
111 
112  // Get snow cloud optics
113  if (do_snow) {
114  cloud_optics.mitchell_ice_optics_sw(ncol, nlev, icswp, des,
115  snow_tau, snow_tau_ssa,
116  snow_tau_ssa_g, snow_tau_ssa_f);
117  }
118  else {
119  // We are not doing snow optics, so set these to zero so we can still use
120  // the arrays without additional logic
121  yakl::memset(snow_tau, 0.);
122  yakl::memset(snow_tau_ssa, 0.);
123  yakl::memset(snow_tau_ssa_g, 0.);
124  yakl::memset(snow_tau_ssa_f, 0.);
125  }
126 
127  // Combine all cloud optics from CAM routines
128  parallel_for(SimpleBounds<3>(nbnd, ncol, nlev), YAKL_LAMBDA (int ibnd, int icol, int ilev)
129  {
130  cld_tau(ibnd, icol, ilev) = ice_tau(ibnd, icol, ilev) + liq_tau(ibnd, icol, ilev);
131  cld_tau_ssa(ibnd, icol, ilev) = ice_tau_ssa(ibnd, icol, ilev) + liq_tau_ssa(ibnd, icol, ilev);
132  cld_tau_ssa_g(ibnd, icol, ilev) = ice_tau_ssa_g(ibnd, icol, ilev) + liq_tau_ssa_g(ibnd, icol, ilev);
133  });
134 
135  if (do_snow) {
136  combine_properties( nbnd, ncol, nlev,
137  cld, cld_tau, cldfsnow, snow_tau, combined_tau);
138 
139  combine_properties( nbnd, ncol, nlev,
140  cld, cld_tau_ssa, cldfsnow, snow_tau_ssa, combined_tau_ssa);
141 
142  combine_properties( nbnd, ncol, nlev,
143  cld, cld_tau_ssa_g, cldfsnow, snow_tau_ssa_g, combined_tau_ssa_g);
144  }
145  else {
146  parallel_for(SimpleBounds<3>(nbnd, ncol, nlev), YAKL_LAMBDA (int ibnd, int icol, int ilev)
147  {
148  combined_tau(ibnd, icol, ilev) = cld_tau(ibnd, icol, ilev);
149  combined_tau_ssa(ibnd, icol, ilev) = cld_tau_ssa(ibnd, icol, ilev);
150  combined_tau_ssa_g(ibnd, icol, ilev) = cld_tau_ssa_g(ibnd, icol, ilev);
151  });
152  }
153 
154  // Copy to output arrays, converting to optical depth, single scattering
155  // albedo, and asymmetry parameter from the products that the CAM routines
156  // return. Make sure we do not try to divide by zero...
157  parallel_for(SimpleBounds<3>(nbnd, ncol, nlev), YAKL_LAMBDA (int iband, int icol, int ilev)
158  {
159  tau_out(icol,ilev,iband) = combined_tau(iband,icol,ilev);
160  if (combined_tau(iband,icol,ilev) > 0) {
161  ssa_out(icol,ilev,iband)
162  = combined_tau_ssa(iband,icol,ilev) / combined_tau(iband,icol,ilev);
163  } else {
164  ssa_out(icol,ilev,iband) = 1.;
165  }
166 
167  if (combined_tau_ssa(iband,icol,ilev) > 0) {
168  asm_out(icol,ilev,iband)
169  = combined_tau_ssa_g(iband,icol,ilev) / combined_tau_ssa(iband,icol,ilev);
170  } else {
171  asm_out(icol,ilev,iband) = 0.;
172  }
173 
174  // Re-order diagnostics outputs
175  liq_tau_out(icol,ilev,iband) = liq_tau(iband,icol,ilev);
176  ice_tau_out(icol,ilev,iband) = ice_tau(iband,icol,ilev);
177  snw_tau_out(icol,ilev,iband) = snow_tau(iband,icol,ilev);
178  });
179 }
void mitchell_ice_optics_sw(const int &ncol, const int &nlev, const real2d &iciwpth, const real2d &dei, real3d &tau, real3d &tau_w, real3d &tau_w_g, real3d &tau_w_f)
Definition: ERF_CloudRadProps.cpp:153
void gammadist_liq_optics_sw(const int &ncol, const int &nlev, const real2d &iclwpth, const real2d &lamc, const real2d &pgam, real3d &tau, real3d &tau_w, real3d &tau_w_g, real3d &tau_w_f)
Definition: ERF_CloudRadProps.cpp:96
static void ec_ice_optics_sw(int ncol, int nlev, int nswbands, const real2d &cldn, const real2d &cicewp, const real2d &rei, const real3d &ice_tau, const real3d &ice_tau_w, const real3d &ice_tau_w_g, const real3d &ice_tau_w_f)
Definition: ERF_EbertCurry.H:11
static void slingo_liq_optics_sw(int ncol, int nlev, int nswbands, const real2d &cldn, const real2d &cliqwp, const real2d &rel, const real3d &liq_tau, const real3d &liq_tau_w, const real3d &liq_tau_w_g, const real3d &liq_tau_w_f)
Definition: ERF_Slingo.H:12
Here is the call graph for this function:

◆ initialize()

void Optics::initialize ( int  ngas,
int  nmodes,
int  num_aeros,
int  nswbands,
int  nlwbands,
int  ncol,
int  nlev,
int  nrh,
int  top_lev,
const std::vector< std::string > &  aero_names,
const real2d &  zi,
const real2d &  pmid,
const real2d &  pdel,
const real2d &  temp,
const real2d &  qi,
const real2d &  geom_radius 
)
19 {
21  aero_optics.initialize(ngas, nmodes, num_aeros,
22  nswbands, nlwbands, ncol, nlev, nrh, top_lev,
23  aero_names, zi, pmid, pdel, temp, qi, geom_radius);
24 }
void initialize(int num_gas, int num_modes, int naeroes, int nswbands_, int nlwbands_, int ncoloum, int nlevel, int num_rh, int top_levels, const std::vector< std::string > &aerosol_names, const real2d &zint, const real2d &pmiddle, const real2d &pdel, const real2d &temperature, const real2d &qtotal, const real2d &geom_rad)
Definition: ERF_AeroRadProps.cpp:12
void initialize()
Definition: ERF_CloudRadProps.cpp:11
AerRadProps aero_optics
Definition: ERF_Optics.H:95
int ngas
Definition: ERF_Optics.H:88
Here is the call graph for this function:

◆ mcica_subcol_mask()

void Optics::mcica_subcol_mask ( int  ngpt,
int  ncol,
int  nlev,
const real2d &  cldfrac,
const bool3d &  iscloudy 
)
430 {
431  // Local vars
432  const real cldmin = 1.0e-80; // min cloud fraction
433  real2d cldf("cldf",ncol,nlev); // cloud fraction clipped to cldmin
434  real3d cdf("cdf", ngpt, ncol, nlev); // random numbers
435 
436  // clip cloud fraction
437  parallel_for(SimpleBounds<2>(ncol, nlev), YAKL_LAMBDA (int icol, int ilev)
438  {
439  cldf(icol,ilev) = cldfrac(icol,ilev);
440  if (cldf(icol,ilev) < cldmin) cldf(icol,ilev) = 0.;
441  });
442 
443  amrex::RandomEngine engine;
444  // Generate random numbers in each subcolumn at every level
445  parallel_for(SimpleBounds<3>(ngpt, ncol, nlev), YAKL_LAMBDA (int isubcol, int icol, int ilev)
446  {
447  cdf(isubcol,icol,ilev) = amrex::Random(engine);
448  });
449 
450  // Maximum-Random overlap
451  // i) pick a random number for top layer.
452  // ii) walk down the column:
453  // - if the layer above is cloudy, use the same random number as in the layer above
454  // - if the layer above is clear, use a new random number
455  parallel_for(SimpleBounds<3>(nlev, ncol, ngpt), YAKL_LAMBDA (int k, int i, int isubcol)
456  {
457  if (k > 1) {
458  if (cdf(isubcol,i,k-1) > 1. - cldf(i,k-1) ) {
459  cdf(isubcol,i,k) = cdf(isubcol,i,k-1);
460  } else {
461  cdf(isubcol,i,k) = cdf(isubcol,i,k) * (1. - cldf(i,k-1));
462  }
463  iscloudy(isubcol,i,k) = cdf(isubcol,i,k) >= 1.-cldf(i,k) ? true : false;
464  }
465  });
466 }

Referenced by sample_cloud_optics_lw(), and sample_cloud_optics_sw().

Here is the caller graph for this function:

◆ sample_cloud_optics_lw()

void Optics::sample_cloud_optics_lw ( int  ncol,
int  nlev,
int  ngpt,
const int1d &  gpt2bnd,
const real2d &  pmid,
const real2d &  cld,
const real2d &  cldfsnow,
const real3d &  tau_bnd,
const real3d &  tau_gpt 
)
327 {
328  //real(r8), dimension(ncol,nlev) :: combined_cld
329  real2d combined_cld("combined_cld", ncol, nlev);
330 
331  //logical, dimension(ngpt,ncol,nlev) :: iscloudy
332  bool3d iscloudy("iscloudy", ngpt, ncol, nlev);
333 
334  // Combine cloud and snow fractions for MCICA sampling
335  parallel_for(SimpleBounds<2>(ncol, nlev), YAKL_LAMBDA (int icol, int ilev)
336  {
337  combined_cld(icol,ilev) = std::max(cld(icol,ilev), cldfsnow(icol,ilev));
338  });
339 
340  // Get the stochastic subcolumn cloudy mask
341  mcica_subcol_mask(ngpt, ncol, nlev, combined_cld, iscloudy);
342 
343  // Map optics to g-points, selecting a single subcolumn for each
344  // g-point. This implementation generates homogeneous clouds, but it would be
345  // straightforward to extend this to handle horizontally heterogeneous clouds
346  // as well.
347  parallel_for(SimpleBounds<3>(ngpt, nlev, ncol), YAKL_LAMBDA (int igpt, int ilev, int icol)
348  {
349  if (iscloudy(igpt,icol,ilev) && combined_cld(icol,ilev) > 0.) {
350  tau_gpt(icol,ilev,igpt) = tau_bnd(icol,ilev,gpt2bnd(igpt));
351  } else {
352  tau_gpt(icol,ilev,igpt) = 0.;
353  }
354  });
355  }
void mcica_subcol_mask(int ngpt, int ncol, int nlev, const real2d &cldfrac, const bool3d &iscloudy)
Definition: ERF_Optics.cpp:428
Here is the call graph for this function:

◆ sample_cloud_optics_sw()

void Optics::sample_cloud_optics_sw ( int  ncol,
int  nlev,
int  ngpt,
const int1d &  gpt2bnd,
const real2d &  pmid,
const real2d &  cld,
const real2d &  cldfsnow,
const real3d &  tau_bnd,
const real3d &  ssa_bnd,
const real3d &  asm_bnd,
const real3d &  tau_gpt,
const real3d &  ssa_gpt,
const real3d &  asm_gpt 
)
290 {
291  //real(r8), dimension(ncol,nlev) :: combined_cld
292  real2d combined_cld("combined_cld", ncol, nlev);
293 
294  //logical, dimension(ngpt,ncol,nlev) :: iscloudy
295  bool3d iscloudy("iscloudy", ngpt, ncol, nlev);
296 
297  // Combined snow and cloud fraction
298  parallel_for(SimpleBounds<2>(ncol, nlev), YAKL_LAMBDA (int icol, int ilev)
299  {
300  combined_cld(icol,ilev) = std::max(cld(icol,ilev), cldfsnow(icol,ilev));
301  });
302 
303  // Get stochastic subcolumn cloud mask
304  mcica_subcol_mask(ngpt, ncol, nlev, combined_cld, iscloudy);
305 
306  // Generate subcolumns for homogeneous clouds
307  parallel_for(SimpleBounds<3>(ngpt, nlev, ncol), YAKL_LAMBDA (int igpt, int ilev, int icol)
308  {
309  if (iscloudy(igpt,icol,ilev) && combined_cld(icol,ilev) > 0.) {
310  tau_gpt(icol,ilev,igpt) = tau_bnd(icol,ilev,gpt2bnd(igpt));
311  ssa_gpt(icol,ilev,igpt) = ssa_bnd(icol,ilev,gpt2bnd(igpt));
312  asm_gpt(icol,ilev,igpt) = asm_bnd(icol,ilev,gpt2bnd(igpt));
313  } else {
314  tau_gpt(icol,ilev,igpt) = 0.;
315  ssa_gpt(icol,ilev,igpt) = 1.;
316  asm_gpt(icol,ilev,igpt) = 0.;
317  }
318  });
319  }
Here is the call graph for this function:

◆ set_aerosol_optics_lw()

void Optics::set_aerosol_optics_lw ( int  icall,
real  dt,
bool  is_cmip6_volc,
const real2d &  zi,
const real3d &  tau,
const real2d &  clear_rh 
)
416 {
417  // Get aerosol absorption optical depth from CAM routine
418  yakl::memset(tau, 0.);
419  aero_optics.aer_rad_props_lw(is_cmip6_volc, icall, dt, zi, tau, clear_rh);
420 }
void aer_rad_props_lw(const bool &is_cmip6_volc, const int &list_idx, const real &dt, const real2d &zi, const real3d &odap_aer, const real2d &clear_rh)
Definition: ERF_AeroRadProps.cpp:304
Here is the call graph for this function:

◆ set_aerosol_optics_sw()

void Optics::set_aerosol_optics_sw ( int  icall,
int  ncol,
int  nlev,
int  nswbands,
real  dt,
const int1d &  night_indices,
bool  is_cmip6_volc,
const real3d &  tau_out,
const real3d &  ssa_out,
const real3d &  asm_out,
const real2d &  clear_rh 
)
361 {
362  // NOTE: aer_rad_props expects 0:pver indexing on these! It appears this is to
363  // account for the extra layer added above model top, but it is not entirely
364  // clear. This is not done for the longwave, and it is not really documented
365  // anywhere that I can find. Regardless, optical properties for the zero index
366  // are set to zero in aer_rad_props_sw as far as I can tell.
367  //
368  // NOTE: dimension ordering is different than for cloud optics!
369  real3d tau("tau", ncol, nlev+1, nswbands);
370  real3d tau_w("tau_w", ncol, nlev+1, nswbands);
371  real3d tau_w_g("tau_w_g", ncol, nlev+1, nswbands);
372  real3d tau_w_f("tau_w_f", ncol, nlev+1, nswbands);
373 
374  // Get aerosol absorption optical depth from CAM routine
375  yakl::memset(tau, 0.);
376  yakl::memset(tau_w, 0.);
377  yakl::memset(tau_w_g, 0.);
378  yakl::memset(tau_w_f, 0.);
379 
380  int1d ic("icount",1);
381  intHost1d ic_host("ic_host",1);
382  parallel_for(SimpleBounds<1>(ncol), YAKL_LAMBDA (int i)
383  {
384  if (night_indices(i) > 0) ++ic(1);
385  });
386  ic.deep_copy_to(ic_host);
387 
388  aero_optics.aer_rad_props_sw(icall, dt,
389  ic_host(1), night_indices, is_cmip6_volc,
390  tau, tau_w, tau_w_g, tau_w_f, clear_rh);
391 
392  // Extract quantities from products
393  parallel_for(SimpleBounds<3>(ncol, nlev, nswbands), YAKL_LAMBDA (int icol, int ilev, int iband)
394  {
395  // Copy cloud optical depth over directly
396  tau_out(icol,ilev,iband) = tau(icol,ilev,iband);
397  // Extract single scattering albedo from the product-defined fields
398  if (tau(icol,ilev,iband) > 0) {
399  ssa_out(icol,ilev,iband) = tau_w(icol,ilev,iband) / tau(icol,ilev,iband);
400  } else {
401  ssa_out(icol,ilev,iband) = 1.;
402  }
403 
404  // Extract asymmetry parameter from the product-defined fields
405  if (tau_w(icol,ilev,iband) > 0) {
406  asm_out(icol,ilev,iband) = tau_w_g(icol,ilev,iband) / tau_w(icol,ilev,iband);
407  } else {
408  asm_out(icol,ilev,iband) = 0.;
409  }
410  });
411  }
void aer_rad_props_sw(const int &list_idx, const real &dt, const int &nnite, const int1d &idxnite, const bool is_cmip6_volc, const real3d &tau, const real3d &tau_w, const real3d &tau_w_g, const real3d &tau_w_f, const real2d &clear_rh)
Definition: ERF_AeroRadProps.cpp:56
Here is the call graph for this function:

Member Data Documentation

◆ aero_optics

AerRadProps Optics::aero_optics
private

◆ cloud_optics

CloudRadProps Optics::cloud_optics
private

◆ gas_names

char** Optics::gas_names
private

◆ icecldoptics

std::string Optics::icecldoptics
private

◆ liqcldoptics

std::string Optics::liqcldoptics
private

◆ ngas

int Optics::ngas
private

Referenced by initialize().


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