Plasma computations
JenaAtomicCalculator.Plasma.AbstractPlasmaScheme
— Typeabstract type Plasma.AbstractPlasmaScheme
... defines an abstract type to distinguish different kinds of plasma computations; see also:
+ struct AverageAtomScheme
... to perform an average-atom computation.
+ struct LineShiftScheme
... to compute the energy shifts and properties of atomic/ionic lines in some selected plasma model.
+ struct RateSchemePhotoionization
... to compute (empirical) photoionization plasma rates and rate coefficients.
+ struct RateSchemePhotorecombination
... to compute (empirical) photorecombination plasma rates and rate coefficients.
+ struct SahaBoltzmannScheme
... to compute thermodynamic properties of a Saha-Boltzmann LTE mixture.
JenaAtomicCalculator.Plasma.AverageAtomScheme
— Typestruct Plasma.AverageAtomScheme <: Plasma.AbstractPlasmaScheme
... a struct to perform an average-atom computation.
+ nMax ::Int64 ... maximum principal quantum number for the subshells of the AA model.
+ lMax ::Int64 ... maximum orbital quantum number for the subshells of the AA model.
+ scField ::AbstractScField ... maximum orbital quantum number for the subshells of the AA model.
+ calcPhotoionizationCs ::Bool ... True, if photoionization cross sections are to be calculated.
+ calcFormFactor ::Bool ... True, if the form factor need to be calculated.
+ calcScatteringFactor ::Bool ... True, if scattering factor need to be calculated.
+ piSubshells ::Array{Subshell,1} ... Bound subshells to be included into the photoionization cross sections.
+ omegas ::Array{Float64,1} ... energies for the photoionization & scattering factors [in a.u.].
+ qValues ::Array{Float64,1} ... q-values for calculating the form factor [in a.u.].
JenaAtomicCalculator.Plasma.AverageAtomScheme
— MethodPlasma.AverageAtomScheme()
... constructor for an 'default' instance of a Plasma.AverageAtomScheme.
JenaAtomicCalculator.Plasma.Computation
— Typestruct Computation
... defines a type for defining (simple) plasma computation for atoms and ions in a given set of reference configurations. It also support different atomic processes under plasma conditions. The plasma environment is typical described in terms of its temperature, density, etc.
+ scheme ::AbstractPlasmaScheme ... Scheme (kind) of plasma computation.
+ nuclearModel ::Nuclear.Model ... Model, charge and parameters of the nucleus.
+ grid ::Radial.Grid ... The radial grid to be used for the computation.
+ refConfigs ::Array{Configuration,1} ... A list of non-relativistic configurations.
+ asfSettings ::AsfSettings
... Provides the settings for the SCF process (under plasma conditions) and the associated CI calculations.
+ settings ::Plasma.Settings ... communicates the properties of the plasma
JenaAtomicCalculator.Plasma.Computation
— MethodPlasma.Computation( ... example for plasma SCF computations)
grid = Radial.Grid(true)
nuclearM = Nuclear.Model(18., "Fermi")
...
refConfigs = [Configuration("[Ne] 3s^2 3p^5")]
Plasma.Computation(Plasma.Computation(), grid=grid, nuclearModel=nuclearM, refConfigs=refConfigs, asfSettings=... )
... These simple examples can be further improved by overwriting the corresponding parameters.
JenaAtomicCalculator.Plasma.Computation
— MethodPlasma.Computation(comp::Plasma.Computation;
scheme=.., nuclearModel=.., grid=.., refConfigs=.., asfSettings=..,
settings=..,
printout::Bool=false)
... constructor for modifying the given Plasma.Computation by 'overwriting' the previously selected parameters.
JenaAtomicCalculator.Plasma.Computation
— MethodPlasma.Computation()
... constructor for an 'empty' instance::Plasma.Computation.
JenaAtomicCalculator.Plasma.LineShiftScheme
— Typestruct Plasma.LineShiftScheme <: Plasma.AbstractPlasmaScheme
... defines a type for the details and parameters of computing level energies with plasma interactions.
+ plasmaModel ::AbstractPlasmaModel ... Specify a particular plasma model, e.g. ion-sphere, Debye.
+ initialConfigs ::Array{Configuration,1} ... List of one or several configurations that define the initial-state multiplet.
+ finalConfigs ::Array{Configuration,1} ... List of one or several configurations that define the final-state multiplet.
+ settings ::AbstractLineShiftSettings ... Specify the process and settings for which line-shifts need to be computed.
## + NoBoundElectrons ::Int64 ... Effective number of bound electrons.
JenaAtomicCalculator.Plasma.LineShiftScheme
— MethodPlasma.LineShiftScheme()
... constructor for a standard instance of Plasma.LineShiftScheme.
JenaAtomicCalculator.Plasma.PartialWaveData
— Typestruct Plasma.PartialWaveData
... defines a type to collect for partial-waves (kappa) cross sections, rates, etc. at different energies
+ kappa ::Int64 ... kappa
+ pairs ::Vector{Tuple{Float64, Float64}}
... vector of pairs, for instance (energy, cs), to later combine data in a proper manner.
JenaAtomicCalculator.Plasma.PartialWaveData
— MethodPlasma.PartialWaveData()
... constructor for the default values of Plasma.PartialWaveData set
JenaAtomicCalculator.Plasma.SahaBoltzmannScheme
— Typestruct Plasma.SahaBoltzmannScheme <: Plasma.AbstractPlasmaScheme
... a struct to thermodynamic properties of a Saha-Boltzmann LTE mixture..
+ plasmaModel ::Basics.AbstractPlasmaModel
... A plasma model that "restricts" the Saha-Boltzmann equilibrium densities by some additional parameters, for instance,
due to ionization-potential-depression (IPD) or others.
+ calcLTE ::Bool ... True, if the Saha-Boltzmann equilibrium densities should be calculated.
+ printIonLevels ::Bool ... True, for printing detailed information about all ionic levels.
+ qRange ::UnitRange{Int64} ... Range of charge states q for which ions are taken into account.
+ maxNoIonLevels ::Int64 ... (maximum) No of ionic levels for any charge state of the ions in the mixture.
+ NoExcitations ::Int64
... No of excitations (S, D, T) that are taken into account with regard to the reference (ground) configuration of the ions.
This number is taken as a second qualifier to characterize the quality of the ionic-level data. Usually, NoExcitations = 1,2.
+ upperShellNo ::Int64
... upper-most princicpal quantum number n for which orbitals are taken into account into the ionic-level computations;
this is often chosen upperShellNo= 4..10; if ionic-level data are given, this number decides which of the levels
are taken into account.
+ isotopicMixture ::Array{Basics.IsotopicFraction,1}
... List of (non-normlized) isotopic fractions that form the requested mixture; the fractions will first be renormalized
to 1 in course of the Saha-Boltzmann LTE computations.
+ isotopeFilenames ::Array{String,1}
... set of files names from which the ionic-level data can be read in for the different isotopes (Z,A) in the mixture.
JenaAtomicCalculator.Plasma.SahaBoltzmannScheme
— MethodPlasma.SahaBoltzmannScheme()
... constructor for an 'default' instance of a Plasma.SahaBoltzmannScheme.
JenaAtomicCalculator.Plasma.Settings
— Typestruct Plasma.Settings
... defines a type for the details and parameters of computing photoionization lines.
+ temperature ::Float64 ... Plasma temperature in [K].
+ density ::Float64 ... Plasma density in [g/cm^3].
+ useNumberDensity ::Bool
... true, if the density above is taken as (total ion) number density ni, and false otherwise.
JenaAtomicCalculator.Plasma.Settings
— MethodPlasma.Settings(set::Plasma.Settings;
temperature=.., density=.., useNumberDensity =..)
... constructor for modifying the given Plasma.Settings by 'overwriting' the previously selected parameters.
JenaAtomicCalculator.Plasma.Settings
— MethodPlasma.Settings()
... constructor for the default values of plasma computations
JenaAtomicCalculator.Basics.perform
— MethodBasics.perform(comp::Plasma.Computation)
... to set-up and perform a plasma computation that starts from a given set of reference configurations and support both, an atomic-average SCF procedure and the computation of various plasma properties and processe. The results of all individual steps are printed to screen but nothing is returned otherwise.
Basics.perform(comp::Plasma.Computation; output::Bool=true)
... to perform the same but to return the complete output in a dictionary; the particular output depends on the type and specifications of the plasma computation but can easily accessed by the keys of this dictionary.
Empirical plasma rates
Line-shift computations
JenaAtomicCalculator.Plasma.computeOutcomes
— MethodPlasma.computeOutcomes(multiplet::Multiplet, nm::Nuclear.Model, grid::Radial.Grid, asfSettings::AsfSettings, settings::Plasma.Settings; output=true)
... to compute a new multiplet from the plasma-modified Hamiltonian matrix as specified by the given settings. The diagonalization of the (plasma-modified) Hamiltonian matrix follows the asfSettings that were applied before for the computation of the multiplet; warnings are issued if features were selected that are not supported for plasma calculations, such as the Breit interaction, QED shifts and others. The energies and plasma shifts (with regard to the unperturbed multiplet) are printed to screen and into the summary file. The generated (plasma) multiplet is returned if output=true and nothing otherwise,
JenaAtomicCalculator.Plasma.displayResults
— MethodPlasma.displayResults(stream::IO, multiplet::Multiplet, newMultiplet::Multiplet, settings::Plasma.Settings)
... to display the energies, M_ms and F-parameters, etc. for the selected levels. A neat table is printed but nothing is returned otherwise.
JenaAtomicCalculator.Plasma.perform
— MethodPlasma.perform(scheme::Plasma.LineShiftScheme, computation::Plasma.Computation; output::Bool=true)
... to perform a Saha-Boltzmann equilibrium computation for a given ion mixture. For output=true, a dictionary is returned from which the relevant results can be can easily accessed by proper keys.
Saha-Boltzmann computations of ionic mixtures
JenaAtomicCalculator.Plasma.IonicClass
— Typestruct Plasma.IonicClass
... a struct to comprise the (necessary) information about a single ionic charge state in Saha-Boltzmann computations.
+ q ::Int64 ... charge state q of the ions.
+ groundEnergy ::Float64 ... ground level energy of the ion.
+ ionLevels ::Array{Plasma.IonicLevel,1} ... ionic levels of this charge state
JenaAtomicCalculator.Plasma.IonicLevel
— Typestruct Plasma.IonicLevel
... a struct to comprise the (necessary) information about a single ionic level in Saha-Boltzmann computations.
+ energy ::Float64 ... total energy E (Z, A, q) of the ionic level.
+ g ::Float64 ... degeneracy g (Z, A, q) of the ionic level.
+ nDensity ::Float64 ... number density of ions I(Z, A, q) in this level.
+ uppermostShell ::Shell ... uppermost shell in the underlying configuration.
JenaAtomicCalculator.Plasma.IsotopeClass
— Typestruct Plasma.IsotopeClass
... a struct to comprise the (necessary) information about a single isotopic fraction in Saha-Boltzmann computations.
+ isotopicDensity ::Float64 ... (number) density n (Z,A) of all ions of this isotope class.
+ Lambda ::Float64 ... thermal length Lambda of the given isotope.
+ dominantEnergy ::Float64
... energy of the "dominant" ion level of this isotope; this depends on energy and help improve convergence
+ isotopicFraction ::Basics.IsotopicFraction ... isotopic fraction in the Saha-Boltzmann mixture.
+ ionClasses ::Array{Plasma.IonicClass,1} ... list of all ions of the isotope I(Z,A).
JenaAtomicCalculator.Plasma.computeDominantIsotopeEnergy
— MethodPlasma.computeDominantIsotopeEnergy(isoClass::IsotopeClass)
... computes the (dominant) energy of all ions of the given isotope class, i.e. an energy at which ionic levels contributes most to the overall mixture; this (dominant) energy is utilized to make the computation more stable; an energy::Float64 is returned. ... Not used in the present code.
JenaAtomicCalculator.Plasma.computeElectronChemicalPotential
— MethodPlasma.computeElectronChemicalPotential(temp::Float64, ne::Float64)
... computes the electron chemical potential for the given electron (number) density; a chemMu::Float64 is returned.
JenaAtomicCalculator.Plasma.computeElectronNumberDensity
— MethodPlasma.computeElectronNumberDensity(temp::Float64, isoClasses::Array{IsotopeClass,1})
... computes the electron (number) density n_e (T) from the density of all ionic levels, charge states and isotopic mixture at temperature T; an nDensity::Float64 is returned.
JenaAtomicCalculator.Plasma.computeIonLevelNumberDensity
— MethodPlasma.computeIonLevelNumberDensity(temp::Float64, q::Int64, ionLevel::IonicLevel, pfIsoClass::Float64, isoDensity::Float64, chemMuE::Float64, dominantEnergy::Float64)
... computes the (number) density of an individual ionLevel of ions of charge state q at temperature T, the electron chemical potential chemMuE and for a given partion function of the associated IsotopeClass pfIonClass; an nDensity::Float64 is returned.
JenaAtomicCalculator.Plasma.computeIonLevelNumberDensityTotal
— MethodPlasma.computeIonLevelNumberDensityTotal(isoClass::IsotopeClass)
... computes the total (number) density of all charge states and individual ionLevel's just by summation. an nDensity::Float64 is returned.
JenaAtomicCalculator.Plasma.computeIonicLevelChemicalPotential
— MethodPlasma.computeIonicLevelChemicalPotential(temp::Float64, ionLevel::IonicLevel, Lambda::Float64)
... computes the ionic-level chemical potential for the given level (number) density; a chemMu::Float64 is returned. ... Not used in the present code.
JenaAtomicCalculator.Plasma.computeIsotopicPartitionFunction
— MethodPlasma.computeIsotopicPartitionFunction(temp::Float64, isoClass::IsotopeClass, chemMuE::Float64, dominantEnergy::Float64)
... computes the (unnormalized) partition function of the given isotopic class in terms of a summation over all charge states and ionic levels of these class. A pf::Float64 is returned.
JenaAtomicCalculator.Plasma.computeMeanIsotopeChargeState
— MethodPlasma.computeMeanIsotopeChargeState(isoClass::IsotopeClass)
... computes the (mean) charge states <q> of all ions of the given isotope; a qbar::Float64 is returned.
JenaAtomicCalculator.Plasma.determineInitialIonDensitiesPropterties
— MethodPlasma.determineInitialIonDensitiesPropterties(scheme::Plasma.SahaBoltzmannScheme, temp::Float64, ionDensity::Float64, isoClass::IsotopeClass)
... determines the initial level densities and all temperature-dependent properties which are not yet defined for solving the Saha-Boltzmann equations. The procedure assumes however that the ionic level information has been provided before. A list of isoClasses::Array{Plasma.IsotopeClass,1} is returned.
JenaAtomicCalculator.Plasma.determineIonNumberDensity
— MethodPlasma.determineIonNumberDensity(rho::Float64, mixture::Array{IsotopicFraction,1})
... determines the ion number density for a given matter density [g/cm^3] and mixture of isotopic fractions. It re-normalizes the isotopic fractions and performs the conversion by taking the mass of ^12C as the basis. A ion number density ni [ions/a_o^3] is returned.
JenaAtomicCalculator.Plasma.determineIonicClasses
— MethodPlasma.determineIonicClasses(scheme::Plasma.SahaBoltzmannScheme, temp::Float64, isoFraction::Basics.IsotopicFraction)
... determines all relevant ionic classes of the given isotope with fraction x at temperature temp which (sufficiently) contribute to the ionic mixture. The procedure selects the relevant charge states but does neither determine the number nor the properties of the ionic levels of these ion classes. These levels need to be later initialized by either explicit computations or reading them from a file. A list of ionClasses::Array{Plasma.IonicClass,1} is returned.
JenaAtomicCalculator.Plasma.determineIpShifts
— MethodPlasma.determineIpShifts(pm::Basics.AbstractPlasmaModel, temp::Float64, ni::Float64, isotopeClass::IsotopeClass)
... determines the ground-state energie as well as the shifts Delta I_p (q –> q+1) for all charge states q = 0 ... qmax of the ions defined by isotopeClass. The procedure assumes that the plasma model pm brings in the current plasma parameters. In some models, the ion number density ni is required as well. Two dictionaries groundE[q]::Dict{Int64,Float64}, deltaIp[q]::Dict{Int64,Float64} are returned with all the energies (shifts) in [Hartree].
JenaAtomicCalculator.Plasma.determineIsotopeClasses
— MethodPlasma.determineIsotopeClasses(scheme::Plasma.SahaBoltzmannScheme, temp::Float64)
... determines which isotope classes which need to be considered for the given Saha-Boltzmann mixture; it first normalizes all isotopic fractions to 1 and then selects the relevant charge states in the mixture due to the given temperature and the number of charge states involved the self-consistent treatment of the mixture. However, this procedure does neither determine the number nor the detailed properties of the ionic levels of these ion classes. These ionic levels and all further properties are first left out (empty) for later specification. A list of isoClasses::Array{Plasma.IsotopicClass,1}
JenaAtomicCalculator.Plasma.determineReferenceConfiguration
— MethodPlasma.determineReferenceConfiguration(ne::Int64)
... determines the reference configuration for an ion with ne electrons. A refConfig::Configuration is returned.
JenaAtomicCalculator.Plasma.determineValenceShells
— MethodPlasma.determineValenceShells(ne::Int64)
... determines the valance shells for an ion with ne electrons. A shellList::Array{Shell,1} is returned.
JenaAtomicCalculator.Plasma.displayIsotopeClasses
— MethodPlasma.displayIsotopeClasses(stream::IO, isoClasses::Array{IsotopeClass,1})
... list all basic information about the ionic mixture, i.e. the isotopic fractions, the thermal length, the isotopic density as well as the number of charge states and ionic levels involved into the computations. A neat table is printed but nothing is returned otherwise.
JenaAtomicCalculator.Plasma.displaySahaBoltzmannEquilibrium
— MethodPlasma.displaySahaBoltzmannEquilibrium(stream::IO, isoClasses::Array{IsotopeClass,1}, pm::Basics.AbstractPlasmaModel; printLevelDensities::Bool=false)
... list some basic information about the ionic mixture, along with the charge states and the ionic densities; it also prints the ionic level number densities if this is required. A neat table is printed but nothing is returned otherwise.
JenaAtomicCalculator.Plasma.displayThermodynamicProperties
— MethodPlasma.displayThermodynamicProperties(stream::IO, temp::Float64, totalIonDensity::Float64, isoClasses::Array{IsotopeClass,1})
... evaluates, converts units and displays the thermodynamic properties of the ionic mixture, based on the given isotopic-class data. A neat table is printed but nothing is returned otherwise.
JenaAtomicCalculator.Plasma.freeEnergyDensity
— MethodPlasma.freeEnergyDensity(temp::Float64, totalIonDensity::Float64, isoClasses::Array{IsotopeClass,1})
... computes the internal energy density u(T) of an isotopic mixture at temperature T; a u::Float64 is returned.
JenaAtomicCalculator.Plasma.generateIonLevelData
— MethodPlasma.generateIonLevelData(scheme::Plasma.SahaBoltzmannScheme, isoClass::IsotopeClass, q::Int64, grid::Radial.Grid)
... generates the ionic-level data for the given isotope class and charge state q of the ionic mixture. These ionic-level data are generated for the requested number in scheme.NoIonLevels and excitations in scheme.NoExcitations of levels. The procedure presently assumes that the electron shells in the reference configuration are filled in standard order with the Z-q electrons. A special treatment of the reference configuration need to be introduced if this is not the case. The energy levels are generated by including scheme.NoExcitations excitations w.r.t. the reference configuration, and all levels with an excitation energy < maxEn are taken into account, up to scheme.NoIonLevels. An ionLevels::Array{IonicLevels,1} is returned where all ionic-level data are properly places but all (number) densities are still set to 0.
JenaAtomicCalculator.Plasma.internalEnergyDensity
— MethodPlasma.internalEnergyDensity(temp::Float64, totalIonDensity::Float64, isoClasses::Array{IsotopeClass,1})
... computes the internal energy density u(T) of an isotopic mixture at temperature T; a u::Float64 is returned.
JenaAtomicCalculator.Plasma.perform
— MethodPlasma.perform(scheme::Plasma.SahaBoltzmannScheme, computation::Plasma.Computation; output::Bool=true)
... to perform a Saha-Boltzmann equilibrium computation for a given ion mixture. For output=true, a dictionary is returned from which the relevant results can be can easily accessed by proper keys.
JenaAtomicCalculator.Plasma.pressure
— MethodPlasma.pressure(temp::Float64, totalIonDensity::Float64, isoClasses::Array{IsotopeClass,1})
... computes the pressure of an isotopic mixture at temperature T; a pressure::Float64 is returned.
JenaAtomicCalculator.Plasma.readEvaluateIonLevelData
— MethodPlasma.readEvaluateIonLevelData(scheme::Plasma.SahaBoltzmannScheme, isoClass::IsotopeClass, grid::Radial.Grid)
... reads in or evaluates the ionic-level data for the given isotope class of the ionic mixture. This steps basically generates/provides all atomic level data upon which Saha-Boltzmann ionic mixture is based on. To generate these data it analyses of whether ion-level data are provided by one of the files in scheme.isotopeFilenames::{String,1}; if this is the case, it further analyzes of whether the provided data fullfill the requested number scheme.NoIonLevels and excitations scheme.NoExcitations of levels for each charge state of the isotope. If no suitable ion-level data are found for an isotope (Z,A) in isoClasses, the routine generates a (new) set of ion-level data for the current computation. In addition, a (new) file with name NewIonLevelsZxxxAyyyy.... is written to disk and can be utilized in any subsequent computation. It is typically assumed that these file are renamed with the proper physical settings for their generation in mind.
A newIsoClass::IsotopeClass is returned where all ionic-level data are properly places
but all (number) densities are still set to 0. All other fields remain unchanged.
JenaAtomicCalculator.Plasma.readEvaluateIonLevelData
— MethodPlasma.readEvaluateIonLevelData(filename::String, isoClass::IsotopeClass, NoExcitations::Int64, upperShellNo::Int64)
... reads in, if available, the ionic-level data for the given isotope class from filename. The ionic-level data are accepted for return, if (1) they belong to isotope (Z,A) in isoClasses and if they fullfill (2) upperShellNo <= filename:upperShellNo, (3) NoExcitations <= filename:NoExcitations. If proper data are found, an updated newIsoClass::IsotopeClass returned and missing otherwise.
JenaAtomicCalculator.Plasma.readUpdateIonLevelDataObsolete
— MethodPlasma.readUpdateIonLevelDataObsolete(filename::String, isoClass::IsotopeClass, NoExcitations::Int64, upperShellNo::Int64, newIonClasses::Array{Plasma.IonicClass,1})
... reads in, if available, the ionic-level data for the given isotope class from filename. The ionic-level data are accepted for return, if (1) they belong to isotope (Z,A) in isoClasses and if they fullfill (2) upperShellNo <= filename:upperShellNo, (3) NoExcitations <= filename:NoExcitations. The newIonClasses are updated/appended in the given ionClasses data and a new ionic-level data file is written out. The procedure terminates with an error if no proper ionic level data are found in filename. Nothing is returned. ... This procedures is not used a present.
JenaAtomicCalculator.Plasma.restrictIonLevelData
— MethodPlasma.restrictIonLevelData(scheme::Plasma.SahaBoltzmannScheme, isoClass::IsotopeClass)
... restricts the compiled ion-level-data tp those levels that need to be included into the Saha-Boltzmann equilibrium iteration. a newIsoClass::IsotopeClass is returned.
JenaAtomicCalculator.Plasma.solveSahaBoltzmannEquilibrium
— MethodPlasma.solveSahaBoltzmannEquilibrium(temp::Float64, isoClasses::Array{IsotopeClass,1}, pm::Basics.AbstractPlasmaModel)
... solves the Saha-Boltzmann equations for a given ionic mixture and temperature; it iterates the ion-level number densities until no relevant changes occur. A newIsoClasses::Array{IsotopeClass,1} with the proper number densities is returned.
JenaAtomicCalculator.Plasma.thermalLength
— MethodPlasma.thermalLength(temp::Float64, isoFraction::Basics.IsotopicFraction)
... returns the thermal length Lambda(T, M) for a particle with mass M at temperature T. A length Lambda::Float64 in {a_o] is returned.
JenaAtomicCalculator.Plasma.updateIpdPlasmaModel
— MethodPlasma.updateIpdPlasmaModel(pm::Basics.AbstractPlasmaModel, temp::Float64, ne::Float64, isotopeClasses::Array{IsotopeClass,1})
... updates the (parameters) of the plasma model pm due to the given temperature temp, electron density ne and the level population of the mixture. A new plasma model npm::type(pm) is returned.
JenaAtomicCalculator.Plasma.writeIonLevelData
— MethodPlasma.writeIonLevelData(filename::String, isoClass::IsotopeClass, NoExcitations::Int64, upperShellNo::Int64)
... writes-out the ionic-level data for the given isotope class to filename. Nothing is returned.
JenaAtomicCalculator.Plasma.writeIonLevelDataRobin
— MethodPlasma.writeIonLevelDataRobin(filenameRobin::String, filenameJac::String)
... reads-in ion-level data from filenameRobin (in tabular form) and writes a proper isoData file for the given isotope class to filenameJac. Nothing is returned.