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

#include <ERF_Slingo.H>

Static Public Member Functions

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)
 
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)
 

Member Function Documentation

◆ slingo_liq_optics_lw()

static void Slingo::slingo_liq_optics_lw ( int  ncol,
int  nlev,
int  nlwbands,
const real2d &  cldn,
const real2d &  iclwpth,
const real2d &  iciwpth,
const real3d &  abs_od 
)
inlinestatic
128  {
129  real2d ficemr("ficemr", ncol, nlev);
130  real2d cwp("cwp", ncol, nlev);
131  real2d cldtau("cldtau", ncol, nlev);
132 
133  parallel_for(SimpleBounds<2>(nlev, ncol), YAKL_LAMBDA (int k, int i)
134  {
135  cwp (i,k) = 1000.0 * iclwpth(i,k) + 1000.0 * iciwpth(i, k);
136  ficemr(i,k) = 1000.0 * iciwpth(i,k)/(std::max(1.e-18, cwp(i,k)));
137  });
138 
139  parallel_for(SimpleBounds<2>(nlev, ncol), YAKL_LAMBDA (int k, int i)
140  {
141  // Note from Andrew Conley:
142  // Optics for RK no longer supported, This is constructed to get
143  // close to bit for bit. Otherwise we could simply use liquid water path
144  //note that optical properties for ice valid only
145  //in range of 13 > rei > 130 micron (Ebert and Curry 92)
146  real kabsl = 0.090361;
147  auto kabs = kabsl*(1.-ficemr(i,k));
148  cldtau(i,k) = kabs*cwp(i,k);
149  });
150 
151  parallel_for(SimpleBounds<3>(nlwbands, ncol, nlev), YAKL_LAMBDA (int lwband, int icol, int ilev)
152  {
153  abs_od(lwband,icol,ilev) = cldtau(icol,ilev);
154  });
155  }

Referenced by Optics::get_cloud_optics_lw().

Here is the caller graph for this function:

◆ slingo_liq_optics_sw()

static void Slingo::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 
)
inlinestatic
16  {
17  real1d wavmin("wavmin", nswbands);
18  real1d wavmax("wavmax", nswbands);
19 
20  // Minimum cloud amount (as a fraction of the grid-box area) to
21  // distinguish from clear sky
22  const real cldmin = 1.0e-80;
23 
24  // Decimal precision of cloud amount (0 -> preserve full resolution;
25  // 10^-n -> preserve n digits of cloud amount)
26  const real cldeps = 0.0;
27 
28  // A. Slingo's data for cloud particle radiative properties (from 'A GCM
29  // Parameterization for the Shortwave Properties of Water Clouds' JAS
30  // vol. 46 may 1989 pp 1419-1427)
31  real1d abarl("abarl", 4); // A coefficient for extinction optical depth
32  real1d bbarl("bbarl", 4); // B coefficient for extinction optical depth
33  real1d cbarl("cbarl", 4); // C coefficient for extinction optical depth
34  real1d dbarl("dbarl", 4); // D coefficient for extinction optical depth
35  real1d ebarl("ebarl", 4); // E coefficient for extinction optical depth
36  real1d fbarl("fbarl", 4); // F coefficient for extinction optical depth
37  parallel_for(SimpleBounds<1>(1), YAKL_LAMBDA (int i)
38  {
39  abarl(1) = 2.817e-02;
40  abarl(2) = 2.682e-02;
41  abarl(3) = 2.264e-02;
42  abarl(4) = 1.281e-02;
43  bbarl(1) = 1.305;
44  bbarl(2) = 1.346;
45  bbarl(3) = 1.454;
46  bbarl(4) = 1.641;
47  cbarl(1) = -5.62e-08;
48  cbarl(2) = -6.94e-06;
49  cbarl(3) = 4.64e-04;
50  abarl(4) = 0.201;
51  dbarl(1) = 1.63e-07;
52  dbarl(2) = 2.35e-05;
53  dbarl(3) = 1.24e-03;
54  dbarl(4) = 7.56e-03;
55  ebarl(1) = 0.829;
56  ebarl(2) = 0.794;
57  ebarl(3) = 0.754;
58  ebarl(4) = 0.826;
59  fbarl(1) = 2.482e-03;
60  fbarl(2) = 4.226e-03;
61  fbarl(3) = 6.560e-03;
62  fbarl(4) = 4.353e-03;
63  });
64 
65  // Caution... A. Slingo recommends no less than 4.0 micro-meters nor
66  // greater than 20 micro-meters. Here we set effective radius limits
67  // for liquid to the range 4.2 < rel < 16 micron (Slingo 89)
68  const real rel_min = 4.2;
69  const real rel_max = 16.;
70 
71  int indxsl;
72 
74 
75  for (auto ns=0; ns<nswbands; ++ns) {
76  // Set index for cloud particle properties based on the wavelength,
77  // according to A. Slingo (1989) equations 1-3:
78  // Use index 1 (0.25 to 0.69 micrometers) for visible
79  // Use index 2 (0.69 - 1.19 micrometers) for near-infrared
80  // Use index 3 (1.19 to 2.38 micrometers) for near-infrared
81  // Use index 4 (2.38 to 4.00 micrometers) for near-infrared
82  if(wavmax(ns) <= 0.7) {
83  indxsl = 1;
84  } else if(wavmax(ns) <= 1.25) {
85  indxsl = 2;
86  } else if(wavmax(ns) <= 2.38) {
87  indxsl = 3;
88  } else if(wavmax(ns) > 2.38) {
89  indxsl = 4;
90  }
91 
92  // Set cloud extinction optical depth, single scatter albedo,
93  // asymmetry parameter, and forward scattered fraction:
94  parallel_for(SimpleBounds<2>(nlev, ncol), YAKL_LAMBDA (int k, int i)
95  {
96  auto abarli = abarl(indxsl);
97  auto bbarli = bbarl(indxsl);
98  auto cbarli = cbarl(indxsl);
99  auto dbarli = dbarl(indxsl);
100  auto ebarli = ebarl(indxsl);
101  auto fbarli = fbarl(indxsl);
102  real tmp1l, tmp2l, tmp3l, g;
103 
104  // note that optical properties for liquid valid only
105  // in range of 4.2 > rel > 16 micron (Slingo 89)
106  if (cldn(i,k) >= cldmin && cldn(i,k) >= cldeps) {
107  tmp1l = abarli + bbarli/std::min(std::max(rel_min,rel(i,k)),rel_max);
108  liq_tau(ns,i,k) = 1000.*cliqwp(i,k)*tmp1l;
109  } else {
110  liq_tau(ns,i,k) = 0.0;
111  }
112 
113  tmp2l = 1. - cbarli - dbarli*std::min(std::max(rel_min,rel(i,k)),rel_max);
114  tmp3l = fbarli*std::min(std::max(rel_min,rel(i,k)),rel_max);
115  // Do not let single scatter albedo be 1. Delta-eddington solution
116  // for non-conservative case has different analytic form from solution
117  // for conservative case, and raddedmx is written for non-conservative case.
118  liq_tau_w(ns,i,k) = liq_tau(ns,i,k) * std::min(tmp2l,.999999);
119  g = ebarli + tmp3l;
120  liq_tau_w_g(ns,i,k) = liq_tau_w(ns,i,k) * g;
121  liq_tau_w_f(ns,i,k) = liq_tau_w(ns,i,k) * g * g;
122  });
123  } // nswbands
124  }
@ micrometer
Definition: ERF_Rad_constants.H:109
static void get_sw_spectral_boundaries(real1d &low_boundaries, real1d &high_boundaries, Units units)
Definition: ERF_Rad_constants.H:205

Referenced by Optics::get_cloud_optics_sw().

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

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