def eval_albedo_mono(self, w: pint.Quantity) -> pint.Quantity: with xr.open_dataset(self.dataset) as ds: wavelengths = w.m_as(ds.w.attrs["units"]) interpolated_albedo = ds.albedo.interp(w=wavelengths) albedo = to_quantity(interpolated_albedo) albedo_array = albedo * np.ones(self.n_layers) return albedo_array
def eval_mono(self, w: pint.Quantity) -> pint.Quantity: w_units = ureg(self.data.ssi.w.attrs["units"]) irradiance = to_quantity( self.data.ssi.interp(w=w.m_as(w_units), method="linear")) # Raise if out of bounds or ill-formed dataset if np.any(np.isnan(irradiance.magnitude)): raise ValueError("dataset interpolation returned nan") return irradiance * self.scale * self._scale_earth_sun_distance()
def integral(self, wmin: pint.Quantity, wmax: pint.Quantity) -> pint.Quantity: # Collect spectral coordinates wavelength_units = self.wavelengths.units s_w = self.wavelengths.m s_wmin = s_w.min() s_wmax = s_w.max() # Select all spectral mesh vertices between wmin and wmax, as well as # the bounds themselves wmin = wmin.m_as(wavelength_units) wmax = wmax.m_as(wavelength_units) w = [ wmin, *s_w[np.where(np.logical_and(wmin < s_w, s_w < wmax))[0]], wmax ] # Make explicit the fact that the function underlying this spectrum has # a finite support by padding the s_wmin and s_wmax values with a very # small margin eps = 1e-12 # nm try: w.insert(w.index(s_wmin), s_wmin - eps) except ValueError: pass try: w.insert(w.index(s_wmax) + 1, s_wmax + eps) except ValueError: pass # Evaluate spectrum at wavelengths w.sort() w = w * wavelength_units interp = self.eval_mono(w) # Compute integral return np.trapz(interp, w)
def air_refractive_index( wavelength: pint.Quantity = ureg.Quantity(550.0, "nm"), number_density: pint.Quantity = _STANDARD_AIR_NUMBER_DENSITY, ) -> np.ndarray: """ Computes the air refractive index. The wavelength dependence of the refractive index is computed using equation 2 from :cite:`Peck1972DispersionAir`. This formula is a fit of measurements of the air refractive index in the range of wavelength from :math:`\\lambda = 240` nm to :math:`1690` nm. The number density dependence is computed using a simple proportionality rule. Parameters ---------- wavelength : quantity Wavelength. number_density : quantity Number density. Default: Air number density at 101325 Pa and 288.15 K. Returns ------- float or ndarray Air refractive index value(s). """ # wavenumber in inverse micrometer sigma = 1 / wavelength.m_as("micrometer") sigma2 = np.square(sigma) # refractivity in parts per 1e8 x = (5791817.0 / (238.0183 - sigma2)) + 167909.0 / (57.362 - sigma2) if isinstance(x, np.ndarray) and isinstance(number_density.magnitude, np.ndarray): x = x[:, np.newaxis] number_density = number_density[np.newaxis, :] # number density scaling x_scaled = x * (number_density / _STANDARD_AIR_NUMBER_DENSITY).m_as("dimensionless") # refractive index index = 1 + x_scaled * 1e-8 return index
def interpolate_radprops(radprops: xr.Dataset, new_z_layer: pint.Quantity) -> xr.Dataset: """ Interpolate a radiative property data set onto a new level altitude grid. Out of bounds values are replaced with zeros. Parameters ---------- radprops : :class:`~xarray.Dataset` Radiative property data set. new_z_layer : :class:`~pint.Quantity`) Layer altitude grid to interpolate onto. Returns ------- interpolated : Dataset Interpolated radiative property data set. """ mask = (new_z_layer >= to_quantity(radprops.z_level).min()) & ( new_z_layer <= to_quantity(radprops.z_level).max()) masked = new_z_layer[ mask] # altitudes that fall within the bounds of radprops # interpolate within radprops altitude bounds (safe) with xr.set_options(keep_attrs=True): interpolated_safe = radprops.interp( z_layer=masked.m_as(radprops.z_layer.units), kwargs=dict(fill_value="extrapolate" ), # handle floating point arithmetic issue method= "nearest", # radiative properties are assumed uniform in altitude layers ) # interpolate over the full range with xr.set_options(keep_attrs=True): interpolated = interpolated_safe.interp( z_layer=new_z_layer.m_as(radprops.z_layer.units), kwargs=dict(fill_value=0.0), method="nearest", ) return interpolated
def eval_mono(self, w: pint.Quantity) -> np.ndarray: """ Evaluate phase function in monochromatic modes. Parameters ---------- w : :class:`pint.Quantity` Wavelength values at which the phase function is to be evaluated. Returns ------- ndarray Evaluated phase function as a 1D or 2D array depending on the shape of `w` (angle dimension comes last). """ return (self.data.sel(i=0, j=0).interp( w=w.m_as(self.data.w.units), mu=np.linspace(-1, 1, self._n_mu), kwargs=dict(bounds_error=True), ).data)
def eval_sigma_t_mono(self, w: pint.Quantity) -> pint.Quantity: with xr.open_dataset(self.dataset) as ds: ds_w_units = ureg(ds.w.attrs["units"]) # find the extinction data variable for dv in ds.data_vars: standard_name = ds[dv].standard_name if "extinction" in standard_name: extinction = ds[dv] wavelength = w.m_as(ds_w_units) xs_t = to_quantity(extinction.interp(w=wavelength)) xs_t_550 = to_quantity( extinction.interp(w=ureg.convert(550.0, ureg.nm, ds_w_units))) fractions = self.eval_fractions() sigma_t_array = xs_t_550 * fractions dz = (self.top - self.bottom) / self.n_layers normalized_sigma_t_array = self._normalize_to_tau( ki=sigma_t_array.magnitude, dz=dz, tau=self.tau_550, ) return normalized_sigma_t_array * xs_t / xs_t_550
def compute_sigma_a( ds: xr.Dataset, wl: pint.Quantity = ureg.Quantity(550.0, "nm"), p: pint.Quantity = ureg.Quantity(101325.0, "Pa"), t: pint.Quantity = ureg.Quantity(288.15, "K"), n: t.Optional[pint.Quantity] = None, fill_values: t.Optional[float] = None, methods: t.Optional[t.Dict[str, str]] = None, ) -> pint.Quantity: """ Compute monochromatic absorption coefficient at given wavelength, pressure and temperature values. Parameters ---------- ds : Dataset Absorption cross section data set. wl : quantity Wavelength [nm]. p : quantity Pressure [Pa]. .. note:: If ``p``, ``t`` and ``n`` are arrays, their lengths must be the same. t : quantity Temperature [K]. .. note:: If the coordinate ``t`` is not in the input dataset ``ds``, the interpolation on temperature is not performed. n : quantity Number density [m^-3]. .. note:: If ``n`` is ``None``, the values of ``t`` and ``p`` are then used only to compute the corresponding number density. fill_values : dict, optional Mapping of coordinates (in ``["w", "pt"]``) and fill values (either ``None`` or float). If not ``None``, out of bounds values are assigned the fill value during interpolation along the wavelength or pressure and temperature coordinates. If ``None``, out of bounds values will trigger the raise of a ``ValueError``. Only one fill value can be provided for both pressure and temperature coordinates. methods : dict, optional Mapping of coordinates (in ``["w", "pt"]``) and interpolation methods. Default interpolation method is linear. Only one interpolation method can be specified for both pressure and temperature coordinates. Returns ------- quantity Absorption coefficient values. Raises ------ ValueError When wavelength, pressure, or temperature values are out of the range of the data set and the corresponding fill value in ``fill_values`` is ``None``. Warnings -------- The values of the absorption cross section at the desired wavelength, pressure and temperature values, :math:`\\sigma_{a\\lambda} (p, T)`, are obtained by interpolating the input absorption cross section data set along the corresponding dimensions. Notes ----- The absorption coefficient is given by: .. math:: k_{a\\lambda} = n \\, \\sigma_{a\\lambda} (p, T) where * :math:`k_{a\\lambda}` is the absorption coefficient [:math:`L^{-1}`], * :math:`\\lambda` is the wavelength [:math:`L`], * :math:`n` is the number density [:math:`L^{-3}`], * :math:`\\sigma_a` is the absorption cross section [:math:`L^2`], * :math:`p` is the pressure [:math:`ML^{-1}T^{-2}`] and * :math:`t` is the temperature [:math:`\\Theta`]. """ if fill_values is None: fill_values = dict(w=None, pt=None) if methods is None: methods = dict(w="linear", pt="linear") for name in ["w", "pt"]: if name not in fill_values: fill_values[name] = None if name not in methods: methods[name] = "linear" # Interpolate along wavenumber dimension xsw = ds.interp( w=(1.0 / wl).m_as(ds.w.units), # wavenumber in cm^-1 method=methods["w"], kwargs=dict( bounds_error=(fill_values["w"] is None), fill_value=fill_values["w"], ), ) # If the data set includes a temperature coordinate, we interpolate along # both pressure and temperature dimensions. # Else, we interpolate only along the pressure dimension. p_m = p.m_as(ds.p.units) p_values = np.array([p_m]) if isinstance(p_m, float) else p_m pz = xr.DataArray(p_values, dims="pt") if "t" in ds.coords: t_m = t.m_as(ds.t.units) t_values = np.array([t_m] * len(p_values)) if isinstance(t_m, float) else t_m tz = xr.DataArray(t_values, dims="pt") interpolated = xsw.interp( p=pz, t=tz, method=methods["pt"], kwargs=dict( bounds_error=(fill_values["pt"] is None), fill_value=fill_values["pt"], ), ) else: interpolated = xsw.interp( p=pz, method=methods["pt"], kwargs=dict( bounds_error=(fill_values["pt"] is None), fill_value=fill_values["pt"], ), ) xs = to_quantity(interpolated.xs) n = p / (_BOLTZMANN * t) if n is None else n return (n * xs).to("km^-1")