def photod_rate(time, time_scale, target, id_type, observatory, time_format, mol_tag): # imported here to avoid circular dependency with activity.gas from ..activity.gas import photo_timescale epoch = Time(time, scale=time_scale, format=time_format) obj = Horizons(id=target, epochs=epoch.jd, location=observatory, id_type=id_type) try: orb = obj.ephemerides() except ValueError: raise delta = (orb['delta'].data * u.au).to('m') r = (orb['r'].data) cat = JPLSpec.get_species_table() mol = cat[cat['TAG'] == mol_tag] name = mol['NAME'].data[0] timescale = photo_timescale(name) if timescale.ndim != 0: # array timescale = timescale[0] beta = (timescale) * r**2 return beta, delta
def test_einstein(): transition_freq_list = [(1611.7935180 * u.GHz).to('MHz'), (177.26111120 * u.GHz).to('MHz')] mol_tag_list = [28001, 27001] temp_estimate = 300. * u.K result = [] catalog_result = [] cat = JPLSpec.get_species_table() for i in range(2): transition_freq = transition_freq_list[i] mol_tag = mol_tag_list[i] mol_data = Phys.from_jplspec(temp_estimate, transition_freq, mol_tag) intl = intensity_conversion(mol_data) mol_data.apply([intl.value] * intl.unit, name='intl') au = einstein_coeff(mol_data) result.append(au.value) mol = cat[cat['TAG'] == mol_tag] mol_name = mol['NAME'].data[0] lam_search = Lamda.query(mol=mol_name.lower()) lam_result = lam_search[1] tran = transition_freq.to('GHz').value lam_found = lam_result[lam_result['Frequency'] == tran] au_cat = lam_found['EinsteinA'] au_cat = au_cat.data[0] catalog_result.append(au_cat) err = (abs((np.array(catalog_result) - np.array(result)) / np.array(catalog_result) * 100)) assert np.all(err < 23.5)
def get_molecular_parameters(molecule_name, tex=50, fmin=1 * u.GHz, fmax=1 * u.THz, catalog='JPL', **kwargs): """ Get the molecular parameters for a molecule from the JPL or CDMS catalog (this version should, in principle, be entirely self-consistent) Parameters ---------- molecule_name : string The string name of the molecule (normal name, like CH3OH or CH3CH2OH, but it has to match the JPL catalog spec) tex : float Optional excitation temperature (basically checks if the partition function calculator works) catalog : 'JPL' or 'CDMS' Which catalog to pull from fmin : quantity with frequency units fmax : quantity with frequency units The minimum and maximum frequency to search over Examples -------- >>> from pyspeckit.spectrum.models.lte_molecule import get_molecular_parameters >>> freqs, aij, deg, EU, partfunc = get_molecular_parameters('CH2CHCN', ... fmin=220*u.GHz, ... fmax=222*u.GHz, ) >>> freqs, aij, deg, EU, partfunc = get_molecular_parameters('CH3OH', ... fmin=90*u.GHz, ... fmax=100*u.GHz) """ if catalog == 'JPL': from astroquery.jplspec import JPLSpec as QueryTool elif catalog == 'CDMS': from astroquery.linelists.cdms import CDMS as QueryTool else: raise ValueError("Invalid catalog specification") speciestab = QueryTool.get_species_table() jpltable = speciestab[speciestab['NAME'] == molecule_name] if len(jpltable) != 1: raise ValueError(f"Too many or too few matches to {molecule_name}") jpltbl = QueryTool.query_lines(fmin, fmax, molecule=molecule_name, parse_name_locally=True) def partfunc(tem): """ interpolate the partition function WARNING: this can be very wrong """ tem = u.Quantity(tem, u.K).value tems = np.array(jpltable.meta['Temperature (K)']) keys = [k for k in jpltable.keys() if 'q' in k.lower()] logQs = jpltable[keys] logQs = np.array(list(logQs[0])) inds = np.argsort(tems) #logQ = np.interp(tem, tems[inds], logQs[inds]) # linear interpolation is appropriate; Q is linear with T... for some cases... # it's a safer interpolation, anyway. # to get a better solution, you can fit a functional form as shown in the # JPLSpec docs, but that is... left as an exercise. # (we can test the degree of deviation there) linQ = np.interp(tem, tems[inds], 10**logQs[inds]) return linQ freqs = jpltbl['FREQ'].quantity freq_MHz = freqs.to(u.MHz).value deg = np.array(jpltbl['GUP']) EL = jpltbl['ELO'].quantity.to(u.erg, u.spectral()) dE = freqs.to(u.erg, u.spectral()) EU = EL + dE # need elower, eupper in inverse centimeter units elower_icm = jpltbl['ELO'].quantity.to(u.cm**-1).value eupper_icm = elower_icm + (freqs.to(u.cm**-1, u.spectral()).value) # from Brett McGuire https://github.com/bmcguir2/simulate_lte/blob/1f3f7c666946bc88c8d83c53389556a4c75c2bbd/simulate_lte.py#L2580-L2587 # LGINT: Base 10 logarithm of the integrated intensity in units of nm2 ·MHz at 300 K. # (See Section 3 for conversions to other units.) # see also https://cdms.astro.uni-koeln.de/classic/predictions/description.html#description CT = 300 logint = np.array(jpltbl['LGINT']) # this should just be a number #from CDMS website sijmu = (np.exp(np.float64(-(elower_icm / 0.695) / CT)) - np.exp(np.float64(-(eupper_icm / 0.695) / CT)))**(-1) * ( (10**logint) / freq_MHz) * (24025.120666) * partfunc(CT) #aij formula from CDMS. Verfied it matched spalatalogue's values aij = 1.16395 * 10**(-20) * freq_MHz**3 * sijmu / deg # we want logA for consistency with use in generate_model below aij = np.log10(aij) EU = EU.to(u.erg).value ok = np.isfinite(aij) & np.isfinite(EU) & np.isfinite(deg) & np.isfinite( freqs) return freqs[ok], aij[ok], deg[ok], EU[ok], partfunc
def from_jplspec(cls, temp_estimate, transition_freq, mol_tag): """Returns relevant constants from JPLSpec catalog and energy calculations Parameters ---------- temp_estimate : `~astropy.units.Quantity` Estimated temperature in Kelvins transition_freq : `~astropy.units.Quantity` Transition frequency in MHz mol_tag : int or str Molecule identifier. Make sure it is an exclusive identifier, although this function can take a regex as your molecule tag, it will return an error if there is ambiguity on what the molecule of interest is. The function `~astroquery.jplspec.JPLSpec.query_lines_async` with the option `parse_name_locally=True` can be used to parse for the exclusive identifier of a molecule you might be interested in. For more information, visit `astroquery.jplspec` documentation. Returns ------- Molecular data : `~sbpy.data.Phys` instance Quantities in the following order from JPL Spectral Molecular Catalog: | Transition frequency | Temperature | Integrated line intensity at 300 K | Partition function at 300 K | Partition function at designated temperature | Upper state degeneracy | Upper level energy in Joules | Lower level energy in Joules | Degrees of freedom """ if isinstance(mol_tag, str): query = JPLSpec.query_lines_async( min_frequency=(transition_freq - (1 * u.GHz)), max_frequency=(transition_freq + (1 * u.GHz)), molecule=mol_tag, parse_name_locally=True, get_query_payload=True) res = dict(query) # python request payloads aren't stable (could be # dictionary or list) # depending on the version, so make # sure to check back from time to time if len(res['Mol']) > 1: raise JPLSpecQueryFailed(("Ambiguous choice for molecule,\ more than one molecule was found for \ the given mol_tag. Please refine \ your search to one of the following tags\ {} by using JPLSpec.get_species_table()\ (as shown in JPLSpec documentation)\ to parse their names and choose your \ molecule of interest, or refine your\ regex to be more specific (hint '^name$'\ will match 'name' exactly with no\ ambiguity).").format(res['Mol'])) else: mol_tag = res['Mol'][0] query = JPLSpec.query_lines( min_frequency=(transition_freq - (1 * u.GHz)), max_frequency=(transition_freq + (1 * u.GHz)), molecule=mol_tag) freq_list = query['FREQ'] if freq_list[0] == 'Zero lines we': raise JPLSpecQueryFailed( ("Zero lines were found by JPLSpec in a +/- 1 GHz " "range from your provided transition frequency for " "molecule tag {}.").format(mol_tag)) t_freq = min(list(freq_list.quantity), key=lambda x: abs(x - transition_freq)) data = query[query['FREQ'] == t_freq.value] df = int(data['DR'].data) lgint = float(data['LGINT'].data) lgint = 10**lgint * u.nm * u.nm * u.MHz elo = float(data['ELO'].data) / u.cm gu = float(data['GUP'].data) cat = JPLSpec.get_species_table() mol = cat[cat['TAG'] == mol_tag] temp_list = cat.meta['Temperature (K)'] * u.K part = list(mol['QLOG1', 'QLOG2', 'QLOG3', 'QLOG4', 'QLOG5', 'QLOG6', 'QLOG7'][0]) temp = temp_estimate f = interp(log(temp.value), log(temp_list.value[::-1]), log(part[::-1])) f = exp(f) partition = 10**(f) part300 = 10**(float(mol['QLOG1'].data)) # yields in 1/cm energy = elo + (t_freq.to(1 / u.cm, equivalencies=u.spectral())) energy_J = energy.to(u.J, equivalencies=u.spectral()) elo_J = elo.to(u.J, equivalencies=u.spectral()) quantities = [ t_freq, temp, lgint, part300, partition, gu, energy_J, elo_J, df, mol_tag ] names = [ 't_freq', 'temp', 'lgint300', 'partfn300', 'partfn', 'dgup', 'eup_J', 'elo_J', 'degfreedom', 'mol_tag' ] # names = ('Transition frequency', # 'Temperature', # 'Integrated line intensity at 300 K', # 'Partition function at 300 K', # 'Partition function at designated temperature', # 'Upper state degeneracy', # 'Upper level energy in Joules', # 'Lower level energy in Joules', # 'Degrees of freedom', 'Molecule Identifier') result = cls.from_dict(dict(zip(names, quantities))) return result
def molecular_data(temp_estimate, transition_freq, mol_tag): """ Returns relevant constants from JPLSpec catalog and energy calculations Parameters ---------- transition_freq : `~astropy.units.Quantity` Transition frequency in MHz temp_estimate : `~astropy.units.Quantity` Estimated temperature in Kelvins mol_tag : int or str Molecule identifier. Make sure it is an exclusive identifier. Returns ------- Molecular data : list List of constants in the following order: | Transtion frequency | Temperature | Integrated line intensity at 300 K | Partition function at 300 K | Partition function at designated temperature | Upper state degeneracy | Upper level energy in Joules | Lower level energy in Joules | Degrees of freedom """ query = JPLSpec.query_lines(min_frequency=(transition_freq - (1 * u.MHz)), max_frequency=(transition_freq + (1 * u.MHz)), molecule=mol_tag) freq_list = query['FREQ'] t_freq = min(list(freq_list.quantity), key=lambda x: abs(x-transition_freq)) data = query[query['FREQ'] == t_freq.value] df = int(data['DR'].data) lgint = float(data['LGINT'].data) lgint = 10**lgint * u.nm * u.nm * u.MHz elo = float(data['ELO'].data) / u.cm gu = float(data['GUP'].data) cat = JPLSpec.get_species_table() mol = cat[cat['TAG'] == mol_tag] temp_list = cat.meta['Temperature (K)'] * u.K part = list(mol['QLOG1', 'QLOG2', 'QLOG3', 'QLOG4', 'QLOG5', 'QLOG6', 'QLOG7'][0]) temp = temp_estimate f = interpolate.interp1d(temp_list, part, 'linear') partition = 10**(f(temp_estimate.value)) part300 = 10 ** (float(mol['QLOG1'].data)) # yields in 1/cm energy = elo + (t_freq.to(1/u.cm, equivalencies=u.spectral())) energy_J = energy.to(u.J, equivalencies=u.spectral()) elo_J = elo.to(u.J, equivalencies=u.spectral()) result = [] result.append(t_freq) result.extend((temp, lgint, part300, partition, gu, energy_J, elo_J, df)) return result
def get_molecular_parameters(molecule_name, molecule_name_jpl=None, tex=50, fmin=1 * u.GHz, fmax=1 * u.THz, line_lists=['SLAIM'], chem_re_flags=0, **kwargs): """ Get the molecular parameters for a molecule from the CDMS database using vamdclib Parameters ---------- molecule_name : string The string name of the molecule (normal name, like CH3OH or CH3CH2OH) molecule_name_jpl : string or None Default to molecule_name, but if jplspec doesn't have your molecule by the right name, specify it here tex : float Optional excitation temperature (basically checks if the partition function calculator works) fmin : quantity with frequency units fmax : quantity with frequency units The minimum and maximum frequency to search over line_lists : list A list of Splatalogue line list catalogs to search. Valid options include SLAIM, CDMS, JPL. Only a single catalog should be used to avoid repetition of transitions and species chem_re_flags : int An integer flag to be passed to splatalogue's chemical name matching tool Examples -------- >>> freqs, aij, deg, EU, partfunc = get_molecular_parameters(molecule_name='CH2CHCN', ... fmin=220*u.GHz, ... fmax=222*u.GHz, ) >>> freqs, aij, deg, EU, partfunc = get_molecular_parameters('CH3OH', ... fmin=90*u.GHz, ... fmax=100*u.GHz) """ from astroquery.splatalogue import Splatalogue from astroquery.jplspec import JPLSpec # do this here, before trying to compute the partition function, because # this query could fail tbl = Splatalogue.query_lines(fmin, fmax, chemical_name=molecule_name, line_lists=line_lists, show_upper_degeneracy=True, **kwargs) if molecule_name_jpl is None: molecule_name_jpl = molecule_name jpltable = JPLSpec.get_species_table()[JPLSpec.get_species_table()['NAME'] == molecule_name_jpl] if len(jpltable) != 1: raise ValueError(f"Too many or too few matches to {molecule_name_jpl}") def partfunc(tem): """ interpolate the partition function WARNING: this can be very wrong """ tem = u.Quantity(tem, u.K).value tems = np.array(jpltable.meta['Temperature (K)']) logQs = jpltable['QLOG1 QLOG2 QLOG3 QLOG4 QLOG5 QLOG6 QLOG7'.split()] logQs = np.array(list(logQs[0])) inds = np.argsort(tems) logQ = np.interp(tem, tems[inds], logQs[inds]) return 10**logQ freqs = (np.array(tbl['Freq-GHz']) * u.GHz if 'Freq-GHz' in tbl.colnames else np.array(tbl['Freq-GHz(rest frame,redshifted)']) * u.GHz) aij = tbl['Log<sub>10</sub> (A<sub>ij</sub>)'] deg = tbl['Upper State Degeneracy'] EU = (np.array(tbl['E_U (K)']) * u.K * constants.k_B).to(u.erg).value return freqs, aij, deg, EU, partfunc
def from_pyradex(self, integrated_flux, mol_data, line_width=1.0 * u.km / u.s, escapeProbGeom='lvg', iter=100, collider_density={'H2': 900*2.2}): """ Calculate production rate from the Non-LTE iterative code pyradex Presently, only the LAMDA catalog is supported by this function. In the future a function will be provided by sbpy to build your own molecular data file from JPLSpec for use in this function. Collider is assumed to be H2O for the default setting since we are only considering comet production rates. See documentation for pyradex installation tips Parameters ---------- integrated_flux : `~astropy.units.Quantity` The integrated flux in K * km/s mol_data : `sbpy.data.phys` `sbpy.data.phys` object that contains AT LEAST the following data: | mol_tag: molecule of interest as string or int JPLSpec identifier | temp: kinetic temperature in gas coma (unit K) | cdensity : cdensity estimate (can be calculated from cdensity_Bockelee) (unit 1/cm^2) | temp_back: (optional) background temperature in K (default is 2.730 K) | lamda_name: (optional) LAMDA molecule identifier to avoid ambiguity. `Lamda.molecule_dict` provides list Keywords that can be used for these values are found under `~sbpy.data.conf.fieldnames` documentation. Make sure to inform yourself on the values needed for each function, their units, and their interchangeable keywords as part of the Phys data class. line_width : `~astropy.units.Quantity` The FWHM line width (really, the single-zone velocity width to scale the column density by: this is most sensibly interpreted as a velocity gradient (dv/length)) in km/s (default is 1.0 km/s) escapeProbGeom : str Which escape probability method to use, available choices are 'sphere', 'lvg', and 'slab' iter : int Number of iterations you wish to perform. Default is 100, more iterations will take more time to do but will offer a better range of results to compare against. The range of guesses is built by having the column density guess and subtracting/adding an order of magnitude for the start and end values of the loop, respectively. i.e. a guess of 1e15 will create a range between 1e14 and 1e16 collider_density : dict Dictionary of colliders and their densities in cm^-3. Allowed entries are any of the following : h2,oh2,ph2,e,He,H,H+ See `~Pyradex` documentation for more information. Default dictionary is {'H2' : 900*2.2} where 900 is the collider density of H2 and 2.2 is the value for the square root of the ratio of reduced masses of H2O/H2 as follows: (mu_H2O/mu_H2)**0.5 = ((18**2/18*2)/((18*2)/(18+2)))**0.5 = 2.2 in order to scale the collisional excitation to H2O as the main collisional partner. (Walker, et al. 2014; de val Borro, et al. 2017; & Schoier, et al. 2004) Returns ------- column density : `~astropy.units.Quantity` column density to use for the Haser model ratio calculation Note: it is normal for pyradex/RADEX to output warnings depending on the setup the user gives it (such as not having found a molecular data file so it searches for it on LAMDA catalog) Examples -------- >>> from sbpy.activity import NonLTE >>> from sbpy.data import Phys >>> import astropy.units as u >>> transition_freq = (177.196 * u.GHz).to(u.MHz) >>> mol_tag = 29002 >>> cdensity_guess = (1.89*10.**(14) / (u.cm * u.cm)) >>> temp_estimate = 20. * u.K >>> mol_data = Phys.from_jplspec(temp_estimate, transition_freq, ... mol_tag) # doctest: +REMOTE_DATA >>> mol_data.apply([cdensity_guess.value] * cdensity_guess.unit, ... name= 'cdensity') # doctest: +REMOTE_DATA >>> mol_data.apply(['HCO+@xpol'], ... name='lamda_name') # doctest: +REMOTE_DATA >>> nonLTE = NonLTE() >>> try: # doctest: +SKIP ... cdensity = nonLTE.from_pyradex(1.234 * u.K * u.km / u.s, ... mol_data, iter=600, ... collider_density={'H2': 900}) # doctest: +REMOTE_DATA ... print(cdensity) # doctest: +REMOTE_DATA ... except ImportError: ... pass Closest Integrated Flux:[1.24925956] K km / s Given Integrated Flux: 1.234 K km / s [1.06363773e+14] 1 / cm2 References ---------- Haser 1957, Bulletin de la Societe Royale des Sciences de Liege 43, 740. Walker, et al., On the Validity of Collider-mass Scaling for Molecular Rotational Excitation, APJ, August 2014. van der Tak, et al., A computer program for fast non-LTE analysis of interstellar line spectra. With diagnostic plots to interpret observed line intensity ratios. A&A, February 12 2013. de Val Borro, et al., Measuring molecular abundances in comet C/2014 Q2 (Lovejoy) using the APEX telescope. Monthly Notices of the Royal Astronomical Society, October 27 2017. Schoier, et al., An atomic and molecular database for analysis of submillimetre line observations. A&A, November 4 2004. """ try: import pyradex except ImportError: raise ImportError('Pyradex not installed. Please see \ https://github.com/keflavich/pyradex/blob/master/INSTALL.rst') if not isinstance(mol_data, Phys): raise ValueError('mol_data must be a `sbpy.data.phys` instance.') register('Production Rates', {'Radex': '2007A&A...468..627V'}) # convert mol_tag JPLSpec identifier to verbose name if needed try: mol_data['lamda_name'] name = mol_data['lamda_name'][0] name = name.lower() except KeyError: if not isinstance(mol_data['mol_tag'][0], str): cat = JPLSpec.get_species_table() mol = cat[cat['TAG'] == mol_data['mol_tag'][0]] name = mol['NAME'].data[0] name = name.lower() else: name = mol_data['mol_tag'][0] name = name.lower() # try various common instances of molecule names and check them against LAMDA before complaining try: Lamda.molecule_dict[name] except KeyError: try_name = "{}@xpol".format(name) try: Lamda.molecule_dict[try_name] name = try_name except KeyError: print('Molecule name {} not found in LAMDA, module tried {} and also\ found no molecule with this identifier within LAMDA. Please\ enter LAMDA identifiable name using mol_data["lamda_name"]\ . Use Lamda.molecule_dict to see all available options.'.format(name, try_name)) raise # define Temperature temp = mol_data['temp'] # check for optional values within mol_data if 'temp_back' in mol_data: tbackground = mol_data['temp_back'] else: tbackground = 2.730 * u.K # define cdensity and iteration parameters cdensity = mol_data['cdensity'].to(1 / (u.cm * u.cm)) cdensity_low = cdensity - (cdensity*0.9) cdensity_high = cdensity + (cdensity*9) # range for 400 iterations cdensity_range = np.linspace(cdensity_low, cdensity_high, iter) fluxes = [] column_density = [] with tempfile.TemporaryDirectory() as datapath: for i in cdensity_range: R = pyradex.Radex(column=i, deltav=line_width, tbackground=tbackground, species=name, temperature=temp, datapath=datapath, escapeProbGeom=escapeProbGeom, collider_densities=collider_density) table = R() # find closest matching frequency to user defined indx = (np.abs(table['frequency']-mol_data['t_freq'])).argmin() radexfreq = table['frequency'][indx] # get table for that frequency values = table[table['frequency'] == radexfreq] # use eq in io.f from Pyradex to get integrated flux in K * km/s int_flux_pyradex = 1.0645 * values['T_B'] * line_width fluxes.append(int_flux_pyradex) column_density.append(i) # closest matching integrated flux from pyradex fluxes = np.array(fluxes) index_flux = (np.abs(fluxes-integrated_flux.to(u.K * u.km / u.s).value)).argmin() # corresponding column density in 1/cm^2 column_density = column_density[index_flux] print('Closest Integrated Flux:{}'.format(fluxes[index_flux] * u.K * u.km / u.s)) print('Given Integrated Flux: {}'.format(integrated_flux)) return column_density
def beta_factor(mol_data, ephemobj): """ Returns beta factor based on timescales from `~sbpy.activity.gas` and distance from the Sun using an `~sbpy.data.ephem` object. The calculation is: parent photodissociation timescale * (distance from comet to Sun)**2 If you wish to provide your own beta factor, you can calculate the equation expressed in units of AU**2 * s , all that is needed is the timescale of the molecule and the distance of the comet from the Sun. Once you have the beta factor you can append it to your mol_data phys object with the name 'beta' or any of its alternative names. Parameters ---------- mol_data : `~sbpy.data.phys` `sbpy.data.phys` object that contains AT LEAST the following data: | mol_tag: Molecular identifier (`int` or `str`) This field can be given by the user directly or found using `~sbpy.data.phys.from_jplspec`. If the mol_tag is an integer, the program will assume it is the JPL Spectral Molecular Catalog identifier of the molecule and will treat it as such. If mol_tag is a string, then it will be assumed to be the human-readable name of the molecule. The molecule MUST be defined in `sbpy.activity.gas.timescale`, otherwise this function cannot be used and the beta factor will have to be provided by the user directly for calculations. The user can obtain the beta factor from the formula provided above. Keywords that can be used for these values are found under `~sbpy.data.conf.fieldnames` documentation. We recommend the use of the JPL Molecular Spectral Catalog and the use of `~sbpy.data.Phys.from_jplspec` to obtain these values in order to maintain consistency. Yet, if you wish to use your own molecular data, it is possible. Make sure to inform yourself on the values needed for each function, their units, and their interchangeable keywords as part of the Phys data class. ephemobj : `~sbpy.data.ephem` `sbpy.data.ephem` object holding ephemeride information including distance from comet to Sun ['r'] and from comet to observer ['delta'] Returns ------- q : `~astropy.units.Quantity` Beta factor 'beta', which can be appended to the original `sbpy.phys` object for future calculations """ # imported here to avoid circular dependency with activity.gas from .core import photo_timescale from ...data import Ephem if not isinstance(ephemobj, Ephem): raise ValueError('ephemobj must be a `sbpy.data.ephem` instance.') if not isinstance(mol_data, Phys): raise ValueError('mol_data must be a `sbpy.data.phys` instance.') orb = ephemobj delta = (orb['delta']).to('m') r = (orb['r']) if not isinstance(mol_data['mol_tag'][0], str): cat = JPLSpec.get_species_table() mol = cat[cat['TAG'] == mol_data['mol_tag'][0]] name = mol['NAME'].data[0] else: name = mol_data['mol_tag'][0] timescale = photo_timescale(name) if timescale.ndim != 0: # array timescale = timescale[0] beta = (timescale) * r**2 return beta