def __init__(self, **kwargs): e = constants('elementary charge') k = constants('Boltzmann constant') eps0 = epsilon_0 #constants('vacuum electric permittivity') amu = constants('atomic mass constant') if 'Z' in kwargs: self.q = kwargs['Z'] * e elif 'q' in kwargs: self.q = kwargs['q'] else: self.q = self.def_Z * e if 'amu' in kwargs: self.m = kwargs['amu'] * amu elif 'm' in kwargs: self.m = kwargs['m'] else: self.m = self.def_amu * amu if 'eV' in kwargs: self.T = kwargs['eV'] * e / k elif 'vth' in kwargs: self.T = self.m * kwargs['vth']**2 / k elif 'T' in kwargs: self.T = kwargs['T'] else: self.T = 1000 if self.T < 0: self.T = 0 logger.warning('Negative temperature interpreted as zero') self.n = kwargs.pop('n', 1e11) self.alpha = kwargs.pop('alpha', 0) self.kappa = kwargs.pop('kappa', float('inf')) self.Z = self.q / e self.amu = self.m / amu self.eV = self.T * k / e self.vth = np.sqrt(k * self.T / self.m) self.omega_p = np.sqrt(self.q**2 * self.n / (eps0 * self.m)) self.freq_p = self.omega_p / (2 * np.pi) self.period_p = 1 / self.freq_p if self.kappa == float('inf'): self.debye = np.sqrt(eps0 * k * self.T / (self.q**2 * self.n)) * np.sqrt( (1.0 + 15.0 * self.alpha) / (1.0 + 3.0 * self.alpha)) else: self.debye = np.sqrt(eps0 * k * self.T / (self.q**2 * self.n)) *\ np.sqrt(((self.kappa - 1.5) / (self.kappa - 0.5)) *\ ((1.0 + 15.0 * self.alpha * ((self.kappa - 1.5) / (self.kappa - 2.5))) / (1.0 + 3.0 * self.alpha * ((self.kappa - 1.5) / (self.kappa - 0.5)))))
def test_species_defaults(): T = 1000 n = 1e11 e = constants('elementary charge') m = constants('electron mass') k = constants('Boltzmann constant') s = Species(n=n, T=T) assert(s.m == approx(m)) assert(s.q == approx(-e)) assert(s.n == approx(n)) assert(s.T == approx(T)) assert(s.vth == approx(np.sqrt(k*1000/m))) assert(s.alpha == 0) assert(s.kappa == float('inf'))
def test_eta(current, electron): args = getargspec(current).args kwargs = {} if 'V' in args or 'eta' in args: # Should have both or none assert 'V' in args assert 'eta' in args e = constants('elementary charge') k = constants('Boltzmann constant') T = 1000 eta = 10 geometry = Cylinder(r=electron.debye, l=10 * electron.debye) I_eta = current(geometry, electron, eta=eta) I_V = current(geometry, electron, V=eta * k * T / e) assert I_eta == approx(I_V)
def jacobsen_density(geometry, biases, currents, species=Electron()): """ Density computed according to the slope in current squared versus voltage (the Jacobsen-Bekkeng method). This assumes that OML theory for a cylindrical probe is exact. Parameters ---------- geometry: Cylinder Probe geometry biases: array-like of floats Probe bias voltages [V] with respect to some common reference currents: 2D array-like of floats Current measurements. currents[i,j] is sample i, probe current j species: Species The attracted species. The density and temperature in the species object is disregarded by this function. Returns ------- Array of computed densities, one element for each row in currents. """ currents = make_array(currents) biases = make_array(biases) if not isinstance(geometry, Cylinder): raise ValueError('Geometry not supported: {}'.format(geometry)) if len(biases) != currents.shape[1]: raise ValueError('The number of columns in currents must equal the'\ 'number of elements in biases') m = species.m q = species.q k = constants('Boltzmann constant') C = 2 / np.sqrt(np.pi) S = 2 * np.pi * geometry.r * geometry.l # surface area slope_factor = -2 * np.pi * m / ((C * S)**2 * q**3) m = len(biases) sv = np.sum(biases) svv = np.sum(biases**2) curr_sq = currents**2 si = np.sum(curr_sq, axis=1) svi = np.sum(biases * curr_sq, axis=1) slope = (m * svi - sv * si) / (m * svv - sv**2) density = np.sqrt(slope_factor * slope) return density
Plots a moving exponential average of "y" versus "x" in a matplotlib plot while showing the raw values of "y" in the background. tau is the relaxation time in the same unit as the value on the x-axis. """ dx = x[1] - x[0] if tau != 0.0: plt.plot(x, y, '#CCCCCC', linewidth=1, zorder=0) p = plt.plot(x, expAvg(y, dx, tau), linewidth=1, label=label) # returning this allows color to be extracted return p eps0 = constants('electric constant') e = constants('elementary charge') me = constants('electron mass') mp = constants('proton mass') kB = constants('Boltzmann constant') ne = 1e10 debye = 1 wpe = np.sqrt(e**2 * ne / (eps0 * me)) vthe = wpe * debye vthi = vthe * np.sqrt(me / mp) Rp = 1 * debye Te = me * vthe**2 / kB V0 = kB * Te / e I0_sp = -e * ne * Rp**2 * np.sqrt(8 * np.pi * kB * Te / me)
pat = re.compile('\d+') files.sort(key=lambda x: int(pat.findall(x)[-1])) files = files[args.f:args.t:args.s] vars = parse_xoopic_input(pjoin(folder, 'input.inp')) lx = vars['Grid'][0]['x1f'] dt = vars['Control'][0]['dt'] mI = vars['Species'][1]['m'] mE = vars['Species'][0]['m'] ######## Normalized Units ############### ne = 1E13 #vars['Load'][0]['density'] ni = ne eps0 = constants('electric constant') kb = constants('Boltzmann constant') e = constants('elementary charge') gamma_e = 5. / 3 if 'BeamEmitter' in vars: units_0 = vars['BeamEmitter'][0]['units'] if units_0 == 'MKS': vthE = vars['BeamEmitter'][0]['temperature'] tEeV = 0.5 * mE * (vthE * vthE) / e tEK = tEeV * 11604.525 units_1 = vars['BeamEmitter'][1]['units'] if units_1 == 'MKS': vthI = vars['BeamEmitter'][1]['temperature'] tIeV = 0.5 * mI * (vthI * vthI) / e
k = 2 * np.pi * np.arange(Mx) / (Mx * dx) print('Length of k: ', len(k)) print('Max of k: ', np.max(k)) Omega, K = np.meshgrid(omega, k, indexing='ij') print('Shape of Omega: ', Omega.shape) F = np.fft.fftn(f, norm='ortho') halflen = np.array(F.shape, dtype=int) // 2 Omega = Omega[:halflen[0], :halflen[1]] K = K[:halflen[0], :halflen[1]] F = F[:halflen[0], :halflen[1]] # Analytical ion-acoustic dispresion relation ne = vars['Region']['Load'][0]['density'] ni = ne eps0 = constants('electric constant') kb = constants('Boltzmann constant') me = constants('electron mass') e = constants('elementary charge') c0 = constants('speed of light in vacuum') mi = vars['Region']['Species'][1][ 'm'] #40*constants('atomic mass constant') nK = vars['Region']['Grid'][0]['J'] gamma_e = 5. / 3 # print(vars['Region']['Load'][0]['units']) if 'BeamEmitter' in vars['Region']: units_0 = vars['Region']['BeamEmitter'][0]['units'] if units_0 == 'MKS': vthE = vars['Region']['BeamEmitter'][0]['temperature']
def OML_current(geometry, species, V=None, eta=None, normalization=None): """ Current collected by a probe according to the Orbital Motion Limited (OML) theory. The model assumes a probe of infinitely small radius compared to the Debye length, and for a cylindrical probe, that it is infinitely long. Probes with radii up to 0.2 Debye lengths (for spherical probes) or 1.0 Debye lengths (for cylindrical probes) are very well approximated by this theory, although the literature is diverse as to how long cylindrical probes must be for this theory to be a good approximation. Parameters ---------- geometry: Plane, Cylinder or Sphere Probe geometry species: Species or array-like of Species Species constituting the background plasma V: float or array-like of floats Probe voltage(s) in [V]. Overrides eta. eta: float or array-like of floats Probe voltage(s) normalized by k*T/q, where q and T are the species' charge and temperature and k is Boltzmann's constant. normalization: 'th', 'thmax', None Wether to normalize the output current by, respectively, the thermal current, the Maxwellian thermal current, or not at all, i.e., current in [A/m]. Returns ------- float if voltage is float. array of floats corresponding to voltage if voltage is array-like. """ if isinstance(species, list): if normalization is not None: logger.error('Cannot normalize current to more than one species') return None if eta is not None: logger.error('Cannot normalize voltage to more than one species') return None I = 0 for s in species: I += OML_current(geometry, s, V, eta) return I q, T = species.q, species.T kappa, alpha = species.kappa, species.alpha k = constants('Boltzmann constant') if V is not None: eta_isarray = isinstance(V, (np.ndarray, list, tuple)) V = make_array(V) eta = -q * V / (k * T) else: eta_isarray = isinstance(eta, (np.ndarray, list, tuple)) eta = make_array(eta) eta = deepcopy(eta) # Prevents this function from changing eta in caller I = np.zeros_like(eta, dtype=float) eta[np.where( eta == 0)] = np.finfo(float).eps # replace zeros with machine eps indices_n = np.where(eta < 0)[0] # indices for repelled particles indices_p = np.where(eta >= 0)[0] # indices for attracted particles if kappa == float('inf'): C = 1.0 D = (1. + 24 * alpha) / (1. + 15 * alpha) E = 4. * alpha / (1. + 24 * alpha) F = (1. + 8 * alpha) / (1. + 24 * alpha) else: C = np.sqrt(kappa - 1.5) * gamma(kappa - 1.) / gamma(kappa - 0.5) D = (1. + 24 * alpha * ((kappa - 1.5)**2 / ((kappa - 2.) * (kappa - 3.)))) / (1. + 15 * alpha * ((kappa - 1.5) / (kappa - 2.5))) E = 4. * alpha * kappa * (kappa - 1.) / ((kappa - 2.) * (kappa - 3.) + 24 * alpha * (kappa - 1.5)**2) F = ((kappa - 1.) * (kappa - 2.) * (kappa - 3.) + 8 * alpha * (kappa - 3.) * (kappa - 1.5)**2) / ((kappa - 2.) * (kappa - 3.) * (kappa - 1.5) + 24 * alpha * (kappa - 1.5)**3) if normalization is None: I0 = normalization_current(geometry, species) elif normalization.lower() == 'thmax': I0 = 1 elif normalization.lower() == 'th': I0 = normalization_current(geometry, species)/\ thermal_current(geometry, species) else: raise ValueError( 'Normalization not supported: {}'.format(normalization)) if isinstance(geometry, Sphere): # repelled particles: if species.kappa == float('inf'): I[indices_n] = I0 * C * D * np.exp( eta[indices_n]) * (1. + E * -eta[indices_n] * (-eta[indices_n] + 4.)) else: I[indices_n] = I0 * C * D * (1. - eta[indices_n] / (kappa - 1.5))**(1. - kappa) * ( 1. + E * (-eta[indices_n]) * (-eta[indices_n] + 4. * ((kappa - 1.5) / (kappa - 1.)))) # attracted particles: eta[indices_p] = np.abs(eta[indices_p]) I[indices_p] = I0 * C * D * (1. + F * eta[indices_p]) elif isinstance(geometry, Cylinder): # repelled particles: if species.kappa == float('inf'): I[indices_n] = I0 * C * D * \ np.exp(eta[indices_n]) * (1. + E * (-eta[indices_n]) * (-eta[indices_n] + 4.)) else: I[indices_n] = I0 * C * D * (1. - eta[indices_n] / (kappa - 1.5))**(1. - kappa) * ( 1. + E * (-eta[indices_n]) * (-eta[indices_n] + 4. * ((kappa - 1.5) / (kappa - 1.)))) # attracted particles: eta[indices_p] = np.abs(eta[indices_p]) if species.kappa == float('inf'): I[indices_p] = I0 * C * D * \ ((2. / np.sqrt(np.pi)) * ( 1 - 0.5 * E * eta[indices_p]) * np.sqrt(eta[indices_p]) + \ np.exp(eta[indices_p]) * (1. + E * eta[indices_p] * (eta[indices_p] - 4.)) * erfc(np.sqrt(eta[indices_p]))) else: C = np.sqrt(kappa - 1.5) * (kappa - .5) / (kappa - 1.0) D = (1. + 3 * alpha * ((kappa - 1.5) / (kappa - 0.5))) / \ (1. + 15 * alpha * ((kappa - 1.5) / (kappa - 2.5))) E = 4. * alpha * kappa * \ (kappa - 1.) / ((kappa - .5) * (kappa - 1.5) + 3. * alpha * (kappa - 1.5)**2) I[indices_p] = (2./np.sqrt(np.pi))*I0 * C * D * (eta[indices_p]/(kappa-1.5))**(1.-kappa) * \ (((kappa - 1.) / (kappa - 3.)) * E * (eta[indices_p]**2) * hyp2f1(kappa - 3, kappa + .5, kappa - 2., 1. - (kappa - 1.5) / (eta[indices_p])) + \ ((kappa - 1.5 - 2. * (kappa - 1.) * eta[indices_p]) / (kappa - 2.)) * E * eta[indices_p] * hyp2f1(kappa - 2, kappa + .5, kappa - 1., 1. - (kappa - 1.5) / (eta[indices_p])) + (1. + E * eta[indices_p] * (eta[indices_p]-((kappa-1.5)/(kappa-1.)))) * hyp2f1(kappa - 1., kappa + .5, kappa, 1. - (kappa - 1.5) / (eta[indices_p]))) else: raise ValueError('Geometry not supported: {}'.format(geometry)) return I if eta_isarray else I[0]
def finite_radius_current(geometry, species, V=None, eta=None, normalization=None, table='laframboise-darian-marholm'): """ A current model taking into account the effects of finite radius by interpolating between tabulated normalized currents. The model only accounts for the attracted-species currents (for which eta<0). It does not extrapolate, but returns ``nan`` when the input parameters are outside the domain of the model. This happens when the normalized potential for any given species is less than -25, when kappa is less than 4, when alpha is more than 0.2 or when the radius is more than 10 or sometimes all the way up towards 100 (as the distribution approaches Maxwellian). Normally finite radius effects are negligible for radii less than 0.2 Debye lengths (spheres) or 1.0 Debye lengths (cylinders). The model can be based on the following tables, as decided by the ``table`` parameter: - ``'laframboise'``. The de-facto standard tables for finite radius currents, tables 5c and 6c in Laframboise, "Theory of Spherical and Cylindrical Langmuir Probes in a Collisionless, Maxwellian Plasma at Rest", PhD Thesis. Covers Maxwellian plasmas only, probe radii ranging from 0 to 100 Debye lengths. - ``'darian-marholm uncomplete'``. These tables covers Maxwellian, Kappa, Cairns and Kappa-Cairns distributions for radii ranging from 0.2 Debye lengths (spheres) or 1.0 Debye length (cylinders) up to 10 Debye lengths. They are not as accurate as ``'laframboise'`` for pure the Maxwellian, but usually within one or two percent. - ``'darian-marholm'``. Same as above, but this is complemented by adding analytical values from OML theory, thereby extending the range of valid radii down to zero Debye lengths. In addition, the values for zero potential are replaced by analytical values (i.e. the thermal current), since these are amongst the most inaccurate in the above, and more accurate values can be analytically computed. - ``'laframboise-darian-marholm'``. This replaces the tabulated values for the Maxwellian distribution in ``'darian-marholm'`` with those of Laframboise. Accordingly this table produces the most accurate result available in any situation, and has the widest available parameter domain, with the probe radius gradually increasing from 10 towards 100 Debye lengths as the distribution approaches the Maxwellian. Parameters ---------- geometry: Plane, Cylinder or Sphere Probe geometry species: Species or array-like of Species Species constituting the background plasma V: float or array-like of floats Probe voltage(s) in [V]. Overrides eta. eta: float or array-like of floats Probe voltage(s) normalized by k*T/q, where q and T are the species' charge and temperature and k is Boltzmann's constant. normalization: 'th', 'thmax', 'oml', None Wether to normalize the output current by, respectively, the thermal current, the Maxwellian thermal current, the OML current, or not at all, i.e., current in [A/m]. table: string Which table to use for interpolation. See detailed description above. Returns ------- float if voltage is float. array of floats corresponding to voltage if voltage is array-like. """ if isinstance(species, list): if normalization is not None: logger.error('Cannot normalize current to more than one species') return None if eta is not None: logger.error('Cannot normalize voltage to more than one species') return None I = 0 for s in species: I += finite_radius_current(geometry, s, V, eta, table=table) return I q, m, n, T = species.q, species.m, species.n, species.T kappa, alpha = species.kappa, species.alpha k = constants('Boltzmann constant') if V is not None: V = make_array(V) eta = -q * V / (k * T) else: eta = make_array(eta) eta = deepcopy(eta) I = np.zeros_like(eta) indices_n = np.where(eta < 0)[0] # indices for repelled particles indices_p = np.where(eta >= 0)[0] # indices for attracted particles if normalization is None: I0 = normalization_current(geometry, species) elif normalization.lower() == 'thmax': I0 = 1 elif normalization.lower() == 'th': I0 = normalization_current(geometry, species)/\ thermal_current(geometry, species) elif normalization.lower() == 'oml': I0 = normalization_current(geometry, species)/\ OML_current(geometry, species, eta=eta) else: raise ValueError( 'Normalization not supported: {}'.format(normalization)) if isinstance(geometry, Sphere): table += ' sphere' elif isinstance(geometry, Cylinder): table += ' cylinder' else: raise ValueError('Geometry not supported: {}'.format(geometry)) R = geometry.r / species.debye if "darian-marholm" in table: table = get_table(table) pts = table['points'] vals = table['values'].reshape(-1) I[indices_p] = I0 * griddata(pts, vals, (1 / kappa, alpha, R, eta[indices_p])) else: table = get_table(table) pts = table['points'] vals = table['values'].reshape(-1) I[indices_p] = I0 * griddata(pts, vals, (R, eta[indices_p])) if (kappa != float('inf') or alpha != 0): logger.warning( "Using pure Laframboise tables discards spectral indices kappa and alpha" ) I[indices_n] = I0 * OML_current( geometry, species, eta=eta[indices_n], normalization='thmax') if any(np.isnan(I)): logger.warning( "Data points occurred outside the domain of tabulated values resulting in nan" ) return I[0] if len(I) == 1 else I
def test_species_eV(): s = Species(n=1e11, eV=0.2) e = constants('elementary charge') k = constants('Boltzmann constant') assert(s.T == approx(0.2*e/k))
def test_species_vth(): s = Species(n=1e11, vth=1e6) assert(s.T == approx(1e6**2*s.m/constants('Boltzmann constant'))) assert(s.vth == approx(1e6))
def test_species_amu(): s = Species(n=1e11, T=1000, amu=16) assert(s.m == approx(16*constants('atomic mass constant')))
def test_species_Z(): s = Species(n=1e11, T=1000, Z=2) assert(s.q == approx(2*constants('elementary charge')))
def test_species_positron(): s = Positron() assert(s.m == approx(constants('electron mass'))) assert(s.q == approx(constants('elementary charge')))
def finite_length_current(geometry, species, V=None, eta=None, normalization=None, interpolate='I'): """ Current collected by a cylindrical probe according to the Marholm-Marchand finite-length model. Assumes small radius compared to the Debye length. Parameters ---------- geometry: Cylinder Probe geometry species: Species or array-like of Species Species constituting the background plasma V: float or array-like of floats Probe voltage(s) in [V]. Overrides eta. eta: float or array-like of floats Probe voltage(s) normalized by k*T/q, where q and T are the species' charge and temperature and k is Boltzmann's constant. normalization: 'th', 'thmax', 'oml', None Whether to normalize the output current by, respectively, the thermal current, the Maxwellian thermal current, the OML current, or not at all, i.e., current in [A/m]. interpolate: 'I', 'g' Whether to interpolate the coefficients of the profile function g and then integrate g to get the current (faster), or if g is integrated it its present grid to get a grid of currents which can be interpolated from. This makes the interpolation linear in I, and avoids the irregular behaviour sometimes experienced for shorter probes due to irregularities in the coefficients. Returns ------- float if voltage is float. array of floats corresponding to voltage if voltage is array-like. """ if isinstance(species, list): if normalization is not None: logger.error('Cannot normalize current to more than one species') return None if eta is not None: logger.error('Cannot normalize voltage to more than one species') return None I = 0 for s in species: I += finite_length_current(geometry, s, V, eta) return I q, m, n, T = species.q, species.m, species.n, species.T kappa, alpha = species.kappa, species.alpha k = constants('Boltzmann constant') if kappa != float('inf') or alpha != 0: logger.error("Finite length effect data only available for Maxwellian") if V is not None: eta_isarray = isinstance(V, (np.ndarray, list, tuple)) V = make_array(V) eta = -q * V / (k * T) else: eta_isarray = isinstance(eta, (np.ndarray, list, tuple)) eta = make_array(eta) if (np.any(eta > 100.)): logger.warning( 'Finite-length theory yields erroneous results for eta>100') if not isinstance(geometry, Cylinder): raise ValueError('Geometry not supported: {}'.format(geometry)) lambd_p = geometry.l / species.debye # Normalized probe length lambd_l = geometry.lguard / species.debye # Normalized left guard length lambd_r = geometry.rguard / species.debye # Normalized right guard length lambd_t = lambd_l + lambd_p + lambd_r # Normalized total length if normalization is None: # I0 = I_OML => I = I_OML * integral of g = actual current geonorm = deepcopy(geometry) geonorm.l = 1 I0 = OML_current(geonorm, species, eta=eta) elif normalization.lower() in ['th', 'thmax']: # I = actual current / I_th geonorm = deepcopy(geometry) geonorm.l = 1 I0 = OML_current(geonorm, species, eta=eta, normalization='th') / geometry.l # Notice the difference between i_th [A/m] and I_th [A] # geonorm.l = 1 makes it per unit length, i.e. [A/m]. # Hence OML_current returns division by i_th and not I_th. # To correct this we need to divide by geometry.l elif normalization.lower() == 'oml': # I = integral of g I0 = (1 / geometry.l) * np.ones_like(eta) else: raise ValueError( 'Normalization not supported: {}'.format(normalization)) if interpolate.lower() == 'g': C, A, alpha, delta = get_lerped_coeffs(lambd_t, eta) def H(zeta): if zeta == float('inf'): return np.zeros_like(alpha) return A * (alpha * (delta - zeta) - 2) * np.exp(-alpha * zeta) / alpha**2 int_g = C * (lambd_p + H(lambd_p + lambd_l) + H(lambd_p + lambd_r) - H(lambd_l) - H(lambd_r)) I = I0 * species.debye * int_g elif interpolate.lower() in ['i', 'i2']: fname = os.path.join(os.path.dirname(os.path.abspath(__file__)), 'params_structured.npz') file = np.load(fname) lambds = file['lambds'] etas = file['etas'] ind = np.where(lambds < lambd_t)[0][-1] if ind == len(lambds) - 1: # extrapolate from highest lambda lambds = lambds[ind] As = file['As'][ind] Cs = file['Cs'][ind] alphas = file['alphas'][ind] deltas = file['deltas'][ind] def H(zeta): res = np.zeros_like(alphas) if (zeta == float('inf')): return np.zeros_like(alphas) else: return As * (alphas * (deltas - zeta) - 2) * np.exp( -alphas * zeta) / alphas**2 int_gs = Cs*(lambd_p+H(lambd_p+lambd_l)+H(lambd_p+lambd_r) \ -H(lambd_l)-H(lambd_r)) else: # interpolate between lambdas As = file['As'][ind:ind + 2] Cs = file['Cs'][ind:ind + 2] alphas = file['alphas'][ind:ind + 2] deltas = file['deltas'][ind:ind + 2] # Stretch whole probe, including guards lambd_ts = lambds[ind:ind + 2] lambd_ls = lambd_l * lambd_ts / lambd_t lambd_ps = lambd_p * lambd_ts / lambd_t lambd_rs = lambd_r * lambd_ts / lambd_t # Stretch only probe segment, constant guards # lambd_ts = lambds[ind:ind+2] # lambd_ps = lambd_ts-lambd_l-lambd_r # lambd_ls = np.array([lambd_l, lambd_l]) # lambd_rs = np.array([lambd_r, lambd_r]) # if lambd_ps[0]<0: # # Probe segment can't be shorter than zero. # # Let it be zero, and (rightly) let it collect zero current. # lambd_ps[0] = 0 # lambd_ts[0] = lambd_l+lambd_r # Cs[0] = 0 # print("WARNING") def H(zetas): return As*(alphas*(deltas-zetas[:,None])-2) \ *np.exp(-alphas*zetas[:,None])/alphas**2 int_gs = Cs*( lambd_ps[:,None] \ +H(lambd_ps+lambd_ls) \ +H(lambd_ps+lambd_rs) \ -H(lambd_ls) \ -H(lambd_rs)) weight = (lambd_ps[1] - lambd_p) / (lambd_ps[1] - lambd_ps[0]) int_gs = weight * int_gs[0] + (1 - weight) * int_gs[1] attracted = np.where(eta >= 0.)[0] repelled = np.where(eta < 0.)[0] over = np.where(eta > 100.)[0] eta[over] = 100 I = I0 * species.debye * np.ones_like(eta) func = interp1d(etas, int_gs) I[attracted] *= func(eta[attracted]) I[repelled] *= lambd_p else: raise ValueError("interpolate must be either 'g' or 'I'") return I if eta_isarray else I[0]
class Species(object): """ Defines a species using a set of flags and keyword parameters. Maxwellian electrons are the default but Kappa, Cairns, and Kappa-Cairns distributions can be specified by setting the spectral indices kappa and alpha. Other parameters can be specified/overridden using the keyword parameters. Examples:: >>> # Maxwellian electrons with density 1e11 m^(-3) and 1000 K >>> species = Electron(n=1e11, T=1000) >>> # Kappa-distributed electrons with kappa=2 >>> species = Electron(n=1e11, T=1000, kappa=2) >>> # Doubly charged Maxwellain Oxygen ions with 0.26 eV >>> species = Species(Z=2, amu=16, n=1e11, eV=0.26) >>> # Kappa-Cairns-distributed protons with kappa = 5 and alpha=0.2 >>> species = Proton(n=1e11, T=1000, kappa=5, alpha=0.2) A plasma is fully specified by a list of species. E.g. for a Maxwellian electron-proton plasma:: >>> plasma = [] >>> plasma.append(Electron(n=1e11, T=1000)) >>> plasma.append(Proton(n=1e11, T=1000)) Keyword parameters: ------------------- m : mass [kg] amu : mass [amu] q : charge [C] Z : charge [elementary charges] n : density [m^(-3)] T : temperature [K] eV : temperature [eV] vth : thermal velocity [m/s] kappa : spectral index kappa for Kappa and Kappa-Cairns distribution alpha : spectral index alpha for Cairns and Kappa-Cairns distribution """ def_Z = 1 def_amu = constants('electron mass') / constants('atomic mass constant') def __init__(self, **kwargs): e = constants('elementary charge') k = constants('Boltzmann constant') eps0 = epsilon_0 #constants('vacuum electric permittivity') amu = constants('atomic mass constant') if 'Z' in kwargs: self.q = kwargs['Z'] * e elif 'q' in kwargs: self.q = kwargs['q'] else: self.q = self.def_Z * e if 'amu' in kwargs: self.m = kwargs['amu'] * amu elif 'm' in kwargs: self.m = kwargs['m'] else: self.m = self.def_amu * amu if 'eV' in kwargs: self.T = kwargs['eV'] * e / k elif 'vth' in kwargs: self.T = self.m * kwargs['vth']**2 / k elif 'T' in kwargs: self.T = kwargs['T'] else: self.T = 1000 if self.T < 0: self.T = 0 logger.warning('Negative temperature interpreted as zero') self.n = kwargs.pop('n', 1e11) self.alpha = kwargs.pop('alpha', 0) self.kappa = kwargs.pop('kappa', float('inf')) self.Z = self.q / e self.amu = self.m / amu self.eV = self.T * k / e self.vth = np.sqrt(k * self.T / self.m) self.omega_p = np.sqrt(self.q**2 * self.n / (eps0 * self.m)) self.freq_p = self.omega_p / (2 * np.pi) self.period_p = 1 / self.freq_p if self.kappa == float('inf'): self.debye = np.sqrt(eps0 * k * self.T / (self.q**2 * self.n)) * np.sqrt( (1.0 + 15.0 * self.alpha) / (1.0 + 3.0 * self.alpha)) else: self.debye = np.sqrt(eps0 * k * self.T / (self.q**2 * self.n)) *\ np.sqrt(((self.kappa - 1.5) / (self.kappa - 0.5)) *\ ((1.0 + 15.0 * self.alpha * ((self.kappa - 1.5) / (self.kappa - 2.5))) / (1.0 + 3.0 * self.alpha * ((self.kappa - 1.5) / (self.kappa - 0.5))))) def omega_c(self, B): """ Returns the angular cyclotron frequency of the species for a given B-field """ return np.abs(B * self.q / self.m) def freq_c(self, B): """ Returns the cyclotron frequency of the species for a given B-field """ return self.omega_c(B) / (2 * np.pi) def period_c(self, B): """ Returns the cyclotron period of the species for a given B-field """ return 1 / self.freq_c(B) def larmor(self, B): """ Returns the Larmor radius of the species for a given B-field """ return self.vth / self.omega_c(B) def __repr__(self): s = "Species(q={}, m={}, n={}, T={}".format(self.q, self.m, self.n, self.T) if (self.kappa != float('inf')): s += ", kappa={}".format(self.kappa) if (self.alpha != 0): s += ", alpha={}".format(self.alpha) s += ")" return s
def finite_length_current_density(geometry, species, V=None, eta=None, z=None, zeta=None, normalization=None): """ Current collected per unit length on a cylindrical probe according to the Marholm-Marchand finite-length model. Assumes small radius compared to the Debye length. Parameters ---------- geometry: Cylinder Probe geometry species: Species or array-like of Species Species constituting the background plasma V: float or array-like of floats Probe voltage(s) in [V]. Overrides eta. eta: float or array-like of floats Probe voltage(s) normalized by k*T/q, where q and T are the species' charge and temperature and k is Boltzmann's constant. z: float or array-like of floats Position(s) along the probe in [m]. Overrides zeta. zeta: float or array-like of floats Position(s) along the probe normalized by the Debye length. normalization: 'th', 'thmax', 'oml', None Wether to normalize the output current per unit length by, respectively, the thermal current per unit length, the Maxwellian thermal current per unit length, the OML current per unit length, or not at all, i.e., current in [A/m]. Returns ------- float if voltage and position are floats. 1D array of floats corresponding to voltage or position if one of them is array-like. 2D array of floats if voltage and position are both array-like, one row per voltage. """ if isinstance(species, list): if normalization is not None: logger.error('Cannot normalize current to more than one species') return None if eta is not None: logger.error('Cannot normalize voltage to more than one species') return None if zeta is not None: logger.error('Cannot normalize position to more than one species') return None i = 0 for s in species: i += finite_length_current_density(geometry, s, V, eta, z, zeta) return i q, m, n, T = species.q, species.m, species.n, species.T kappa, alpha = species.kappa, species.alpha k = constants('Boltzmann constant') if kappa != float('inf') or alpha != 0: logger.error("Finite length effect data only available for Maxwellian") if V is not None: eta_isarray = isinstance(V, (np.ndarray, list, tuple)) V = make_array(V) eta = -q * V / (k * T) # V must be array (not list) to allow this else: eta_isarray = isinstance(eta, (np.ndarray, list, tuple)) eta = make_array(eta) # eta = eta[:, None] # Make eta rows if not isinstance(geometry, Cylinder): raise ValueError('Geometry not supported: {}'.format(geometry)) if zeta is None: zeta = 0.5 * geometry.l / species.debye if z is not None: zeta_isarray = isinstance(z, (np.ndarray, list, tuple)) z = make_array(z) zeta = z / species.debye # z must be array (not list) to allow this else: zeta_isarray = isinstance(zeta, (np.ndarray, list, tuple)) zeta = make_array(zeta) lambd_p = geometry.l / species.debye # Normalized probe length lambd_l = geometry.lguard / species.debye # Normalized left guard length lambd_r = geometry.rguard / species.debye # Normalized right guard length lambd_t = lambd_l + lambd_p + lambd_r # Normalized total length C, A, alpha, delta = get_lerped_coeffs(lambd_t, eta) # Make these column vectors, i.e. one eta per row C = C[:, None] A = A[:, None] alpha = alpha[:, None] delta = delta[:, None] eta = eta[:, None] if normalization is None: # i0 = i_OML => i = i_OML * g geonorm = deepcopy(geometry) geonorm.l = 1 i0 = OML_current(geonorm, species, eta=eta) elif normalization.lower() in ['th', 'thmax']: # i0 = i_th => i = i_th * g geonorm = deepcopy(geometry) geonorm.l = 1 i0 = OML_current(geonorm, species, eta=eta, normalization='th') elif normalization.lower() == 'oml': # i0 = 1 => i = g i0 = 1 else: raise ValueError( 'Normalization not supported: {}'.format(normalization)) def h(zeta): res = np.zeros((len(alpha), len(zeta))) ind = np.where(zeta != np.inf)[0] # res=0 where zeta=inf zeta = zeta[ind] res[:, ind] = A * (zeta - delta + alpha**(-1)) * np.exp(-alpha * zeta) return res g = C * (1 + h(lambd_l + zeta) + h(lambd_p + lambd_r - zeta)) i = i0 * g if zeta_isarray: if eta_isarray: return i else: return i.ravel() else: if eta_isarray: return i.ravel() else: return i[0][0]
from langmuir import * from scipy.constants import value as constants from scipy.optimize import root_scalar n = 1e11 T = 1000 e = constants('elementary charge') k = constants('Boltzmann constant') plasma = [Electron(n=n, T=T), Hydrogen(n=n, T=T)] geometry = Sphere(r=0.2 * debye(plasma)) def res(V): return OML_current(geometry, plasma, V) sol = root_scalar(res, x0=-3, x1=0) print(sol.root) print(-2.5 * k * T / e)