CarsSpectrum

Synthesize coherent anti-stokes Raman spectra for common gas species.

class carspy.cars_synth.CarsSpectrum(pressure=1, init_comp=None, chi_set='SET 1', **kwargs)[source]

Basic class for synthesizing CARS spectrum.

Input sample volume conditions.

Parameters
pressurefloat, optional

Pressure in bars, by default 1.

init_compdict, optional

Initial composition in the measurement volume (prior to chemical reaction). If not given standard air composition is assumed as:

{
'N2': 0.79,
'Ar': 0.0,
'CO2': 0,
'CO': 0,
'H2': 0,
'O2': 0.21,
'H2O': 0
}
chi_setstr, optional

Choose from the available set of non-resonant susceptibilities, by default ‘SET 1’:

  • ‘SET 1’: based on CARSFT [Pal89].

  • ‘SET 2’: based on NRC [PS88].

  • ‘SET 3’: based on Eckbreth [Eck96].

Other Parameters
speciesstr, optional

Specify the molecule, by default ‘N2’. Currently only supports ‘N2’.

custom_dictdict, optional

Specify custom molecular constants, by default None. This can be used to modify default molecular constants and/or add custom species. The dictionary should contain all the necessary keys as in the default:

['we', 'wx', 'wy', 'wz', 'Be', 'De', 'alpha_e', 'beta_e',
'gamma_e', 'H0', 'He', 'mu', 'MW', 'Const_Raman', 'G/A']
chi_nrs_est(temperature, local_comp=None)[source]

Calculate effective local nonresonant susceptibility.

Parameters
temperaturefloat

Temperature in [K].

local_compdict, optional

A dictionary of local gas composition (percentile) for the purpose of calculating nonresonant background and collisional narrowing, by default None. If not given, initial composition will be used.

chi_rs_gmat(nu_s, temperature, vs=3, js=70, branches=(0), del_Tv=0.0)[source]

Resonant susceptibility based on G-matrix.

This is implemented following [KFP85].

Parameters
nu_s1d array of floats

Stokes frequencies (i.e., the calculation spectral domain).

temperaturefloat

Temperature in [K].

vsint, optional

Total number of vibrational levels to be considered, by default 3.

jsint, optional

Total number of rotational levels to be considered, by default 70.

brancheslist of int, optional

Branches to be considered, by default (0,) (i.e., only Q-branch).

del_Tvfloat, optional

Absolute difference between vibrational and rotational temperature, by default 0.

Returns
1d array (complex)

Theoretical resonant contributions (complex) based on G-matrix.

chi_rs_isolated(nu_s, temperature, vs=3, js=70, branches=(0, 2, - 2), del_Tv=0)[source]

Resonant susceptibility based on isolated line approximation.

Parameters
nu_s1d array of floats

Stokes frequencies (i.e., the calculation spectral domain).

temperaturefloat

Temperature in the probe volume.

vsint, optional

Total number of vibrational levels to be considered, by default 3.

jsint, optional

Total number of rotational levels to be considered, by default 70.

brancheslist of int, optional

Branches to be considered, by default (2, -2, 0).

del_Tvfloat, optional

Absolute difference between vibrational and rotational temperature, by default 0.

Returns
1d array, 1d array

Theoretical resonant contributions (complex) based on isolated line assumption.

num_dens(temperature)[source]

Calculate number density in the probe volume.

Parameters
temperaturefloat

Temperature in [K].

Returns
float

Number density in the probe volume in [\(\mathrm{molecules}/\mathrm{cm}^3\)].

peak_check(x_mol, temperature, v, j, branch=0, Gamma_j=None)[source]

Calculate the theoretical peak intensity.

This is done for the specified transition, assuming isolated lines.

Note

The purpose of this method is to check against the values tabulated in the NRC report [PS88] on D4-D9. The agreement is excellent if x_mol is set to 1 and if the Gamma_j from NRC is used.

Parameters
x_molfloat

Mole fraction of probed molecule within [0, 1].

temperaturefloat

Temperature in the probe volume.

vint

Vibrational quantum number.

jint

Rotational quantum number.

branchint, optional

Q, S or O-branch with a value of 0, 2 or -2, by default 0.

Gamma_jfloat, optional

Raman linewidth, by default None. When not given, the value will be calculated from carspy.line_strength.linewidth_isolated.

Returns
float

Peak transition amplitude based on specified v, j and branch.

relax_mat(temperature, js=70, mode='MEG')[source]

Calculate relaxation rate matrix.

The matrix elements are calculated using different exponential gap laws. The calculations are conducted for a fixed number of j levels, independent of v or branch.

Note

The (absolute) diagonal values are the HWHM of Raman linewidth. They are cross checked against NRC report [PS88] on D4-D9 at 1 atm and 1550 K.

Parameters
temperaturefloat

Temperature in [K].

jsint, optional

Total number of rotational levels to be considered, by default 70.

modestr, optional

By default ‘MEG’ is used. Possible options are:

  • ‘MEG’: Modified exponential gap law considering only N2-N2 collisions [RP86]. The fitting parameters are taken from [WS90].

  • ‘EMEG’: Extended MEG which weighs the contributions of N2, CO2 and H2O [WS90] (to be implemented).

  • ‘XMEG’: Extended MEG plus vibrational dephasing [PGS+90] (to be implemented).

Returns
2d array with the shape of (js, js)

Relaxation matrix with the shape of js x js.

relax_rate(temperature, j_i, j_j, fit_param=None)[source]

Calculate relaxation rates.

This is carried out using the modified exponential gap law (MEG) [RP86].

Note

Only considers Q-branch. O- and S-branches are taken as the same as Q-branch. It is also independent of vibrational state. Only valid for N2 at the moment.

Attention

Due to nuclear-spin selection rules for N2, gamma_ji is zero when (j_j-j_i) is odd. This would need to be adjusted for other types of molecules.

Parameters
temperaturefloat

Temperature in [K].

j_iint

Rotational quantum number of the initial state.

j_jint

Rotational quantum number of the final state.

fit_paramlist

The four fitting parameters in a list: alpha, beta, sigma, m, default None. Needs to be provided.

Returns
gamma_ji, gamma_ijfloat

Upward and downward relaxation rate.

signal_as(temperature, x_mol=None, nu_s=None, pump_lw=None, synth_mode=None, eq_func=<function _ensureCantera.<locals>.no_op>, **kwargs)[source]

Anti-Stokes signal convoluted with a chosen laser lineshape.

Note

The Kataoka (or T-K) convolution [KMH82, Tee84] is a simplified form derived by [FR85]. It is essentially an average of convolve(chi, pump)**2 and convolve(chi**2, pump), which takes care of the cross coherence effect. Excellent agreements are found when compared with Figs.13-15 in the NRC report [PS88].

Parameters
temperaturefloat

Temperature in the probe volume.

x_molfloat, optional

Mole fraction of probed molecule, by default None.

nu_s1d array of float, optional

Stokes frequencies (i.e., the calculation spectral domain), by default None.

pump_lwfloat, optional

Pump laser linewdith (FWHM) in [\(\mathrm{cm}^{-1}\)], by default None (i.e., no laser convolution is performed).

synth_modedict, optional

A dictionary containing the control parameters for creating the CARS spectrum, by default:

pump_ls‘Gaussian’

Choose between pump laser lineshape between ‘Gaussian’ and ‘Lorentzian’.

chi_rs‘G-matrix’

The method to compute resonant susceptibility: ‘isolated’ or ‘G-matrix’ (collisional narrowing).

convol‘Kataoka’ or ‘K’

Two ways to convolve the laser line with the resonant CARS susceptibilities:

  • ‘Yuratich’ or ‘Y’: One convolution, no cross coherence effects accounted for [Yur79].

  • ‘Kataoka’ or ‘K’: Cross-coherence convolution, following [KMH82, Tee84]. This is necessary if pump linewidth is comparable to the Raman linewidth (as is often the case at low T) and if nonresonant contribution competes with resonant signal.

doppler_effectTrue

Whether or not to take into account additional Doppler broadening.

chem_eqFalse

Whether or not to assume chemical equilibrium. If True, an eq_func needs to be provided.

eq_funcfunc, optional

A function used to calculate local gas composition based on temperature and initial gas composition. Default is a simple template employing cantera.

Returns
1d array of floats, 1d array of floats

Spectral locations in wavenumbers and theoretical CARS spectrum. If pump_lw is not specified, I_as needs to be multiplied by 1e-30 to retain its actual physical magnitude in [\(\mathrm{cm}^3/\mathrm{erg}\)].

Other Parameters
vsint, optional

Total number of vibrational levels to be considered, by default 3.

jsint, optional

Total number of rotational levels to be considered, by default 70.

brancheslist of int, optional

Branches to be considered, by default (2, -2, 0).

del_Tvfloat, optional

Absolute difference between vibrational and rotational temperature, by default 0.

trans_amp(v, j, branch=0)[source]

Transition amplitude d or polarizability matrix element.

It’s related to Raman differential scattering cross-section and number density. Note that it is molecule based (divided number density N).

Parameters
vint

Vibrational quantum number.

jint

Rotational quantum number.

branchint, optional

Q, S or O-branch with a value of 0, 2 or -2, by default 0.

Returns
d_squaredfloat

Transition amplitude squared based on specified v, j and branch.