Autoionization

JenaAtomicCalculator.AutoIonization.ChannelType

struct Channel ... defines a type for a AutoIonization channel to help characterize a scattering (continuum) state of many electron-states with a single free electron.

+ kappa          ::Int64                ... partial-wave of the free electron
+ symmetry       ::LevelSymmetry        ... total angular momentum and parity of the scattering state
+ phase          ::Float64              ... phase of the partial wave
+ amplitude      ::Complex{Float64}     ... Auger amplitude associated with the given channel.
source
JenaAtomicCalculator.AutoIonization.LineType

struct Line ... defines a type for a AutoIonization line that may include the definition of sublines and their corresponding amplitudes.

+ initialLevel   ::Level           ... initial-(state) level
+ finalLevel     ::Level           ... final-(state) level
+ electronEnergy ::Float64         ... Energy of the (incoming free) electron.
+ totalRate      ::Float64         ... Total rate of this line.
+ angularAlpha   ::Float64         ... Angular alpha_2 coefficient.
+ channels       ::Array{AutoIonization.Channel,1}  ... List of AutoIonization channels of this line.
source
JenaAtomicCalculator.AutoIonization.PlasmaSettingsType

struct AutoIonization.PlasmaSettings <: Basics.AbstractLineShiftSettings ... defines a type for the details and parameters of computing Auger rates with plasma interactions.

+ printBefore         ::Bool             ... True, if all energies and lines are printed before their evaluation.
+ lineSelection       ::LineSelection    ... Specifies the selected levels, if any.
source
JenaAtomicCalculator.AutoIonization.SettingsType

struct Settings <: AbstractProcessSettings ... defines a type for the details and parameters of computing Auger lines.

+ calcAnisotropy      ::Bool               ... True, if the intrinsic alpha_2,4 angular parameters are to be 
                                                calculated, and false otherwise.
+ calcTeAuger         ::Bool               
    ... True, if contributions of the two-electron Auger transitions are to be calculated, and false otherwise;
        this flag requires a proper (resonant) Green function that supports the TEA transitions.
+ printBefore         ::Bool               ... True, if all energies and lines are printed before their evaluation.
+ lineSelection       ::LineSelection      ... Specifies the selected levels, if any.
+ augerEnergyShift    ::Float64            ... An overall energy shift for all Auger (free-electron) energies.
+ minAugerEnergy      ::Float64            ... Minimum energy of free (Auger) electrons to be included.
+ maxAugerEnergy      ::Float64            ... Maximum energy of free (Auger) electrons to be included.
+ maxKappa            ::Int64              ... Maximum kappa value of partial waves to be included.
+ operator            ::AbstractEeInteraction   
    ... Auger operator that is to be used for evaluating the Auger amplitudes; allowed values are: 
        CoulombInteraction(), BreitInteraction(), ...
+ gMultiplet          ::Multiplet      
    ... Mean-field multiplet of intermediate levels in the computations, sometimes referred to as
        (resonant) Green function.
source
JenaAtomicCalculator.AutoIonization.SettingsMethod

AutoIonization.Settings(set::AutoIonization.Settings;

    calcAnisotropy=..,      calcTeAuger..,              printBefore=..,         augerEnergyShift=.., 
    minAugerEnergy=..,      maxAugerEnergy=..,          maxKappa=..,            operator=..,
    gMultiplet=.. )
                
... constructor for modifying the given AutoIonization.Settings by 'overwriting' the previously selected parameters.
source
JenaAtomicCalculator.AutoIonization.amplitudeMethod

AutoIonization.amplitude(kind::AbstractEeInteraction, channel::AutoIonization.Channel, continuumLevel::Level, initialLevel::Level, grid::Radial.Grid; printout::Bool=true) ... to compute the kind in CoulombInteraction(), BreitInteraction(), CoulombBreit(), CoulombGaunt() Auger amplitude <(alphaf Jf, kappa) Ji || O^(Auger, kind) || alphai J_i> due to the interelectronic interaction for the given final and initial level. A value::ComplexF64 is returned.

source
JenaAtomicCalculator.AutoIonization.channelAmplitudeMethod

AutoIonization.channelAmplitude(kind::String, channel::AutoIonization.Channel, energy::Float64, finalLevel::Level, initialLevel::Level, grid::Radial.Grid) ... to compute the kind = (CoulombInteraction(), BreitInteraction(), CoulombBreit(), CoulombGaunt()) Auger amplitude <(alphaf Jf, kappa) Ji || O^(Auger, kind) || alphai J_i> due to the interelectronic interaction for the given final and initial level. A newChannel::AutoIonization.Channel is returned.

source
JenaAtomicCalculator.AutoIonization.computeAmplitudesPropertiesMethod

AutoIonization.computeAmplitudesProperties(line::AutoIonization.Line, nm::Nuclear.Model, grid::Radial.Grid, nrContinuum::Int64, settings::AutoIonization.Settings; printout::Bool=true) ... to compute all amplitudes and properties of the given line; a line::AutoIonization.Line is returned for which the amplitudes and properties are now evaluated.

source
JenaAtomicCalculator.AutoIonization.computeAmplitudesPropertiesPlasmaMethod

AutoIonization.computeAmplitudesPropertiesPlasma(line::AutoIonization.Line, nm::Nuclear.Model, grid::Radial.Grid, nrContinuum::Int64, settings::AutoIonization.PlasmaSettings; printout::Bool=true) ... to compute all amplitudes and properties of the given line but for the given plasma model; a line::AutoIonization.Line is returned for which the amplitudes and properties are now evaluated.

source
JenaAtomicCalculator.AutoIonization.computeLinesMethod

AutoIonization.computeLines(finalMultiplet::Multiplet, initialMultiplet::Multiplet, nm::Nuclear.Model, grid::Radial.Grid, settings::AutoIonization.Settings; output=true, printout::Bool=true) ... to compute the Auger transition amplitudes and all properties as requested by the given settings. A list of lines::Array{AutoIonization.Lines} is returned.

source
JenaAtomicCalculator.AutoIonization.computeLinesCascadeMethod

AutoIonization.computeLinesCascade(finalMultiplet::Multiplet, initialMultiplet::Multiplet, nm::Nuclear.Model, grid::Radial.Grid, settings::AutoIonization.Settings; output::Bool=true, printout::Bool=true) ... to compute the Auger transition amplitudes and all properties as requested by the given settings. The computations and printout is adapted for large cascade computations by including only lines with at least one channel and by sending all printout to a summary file only. A list of lines::Array{AutoIonization.Lines} is returned.

source
JenaAtomicCalculator.AutoIonization.computeLinesFromOrbitalsMethod

AutoIonization.computeLinesFromOrbitals(finalMultiplet::Multiplet, initialMultiplet::Multiplet, nm::Nuclear.Model, grid::Radial.Grid, settings::AutoIonization.Settings, contOrbitals::Dict{Subshell, Orbital}; output::Bool=true, printout::Bool=true) ... to compute the Auger transition amplitudes and all properties as requested by the given settings but by using the given set of continuum orbitals. The computations and printout is adapted for large cascade computations by including only lines with at least one channel and by sending all printout to a summary file only. A list of lines::Array{AutoIonization.Lines} is returned.

source
JenaAtomicCalculator.AutoIonization.computeLinesPlasmaMethod

AutoIonization.computeLinesPlasma(finalMultiplet::Multiplet, initialMultiplet::Multiplet, nm::Nuclear.Model, grid::Radial.Grid, settings::AutoIonization.PlasmaSettings; output=true) ... to compute the Auger transition amplitudes and all properties as requested by the given settings. A list of lines::Array{AutoIonization.Lines} is returned.

source
JenaAtomicCalculator.AutoIonization.computeTeaAmplitudeMethod

AutoIonization.computeTeaAmplitude(kind::AbstractEeInteraction, channel::AutoIonization.Channel, continuumLevel::Level, gMultiplet::Multiplet, initialLevel::Level, grid::Radial.Grid; printout::Bool=true) ... to compute the kind in CoulombInteraction(), BreitInteraction(), CoulombBreit() total Auger amplitude for the two-electron Auger transitions via the gMultiplet as the resonant Green function:

    <(alpha_f J_f, kappa) J_i || O^(TEA, kind) || alpha_i J_i> 
    
                <(alpha_f J_f, kappa) J_i || O^(Auger, kind) || alpha_in J_n> <J_n || O^(e-e, kind) || alpha_i J_i>
            = ---------------------------------------------------------------------------------------------------
                                                                E_i  -  E_n
                                                                
    due to the interelectronic interaction as well as the given initial, intermediate (gMultiplet), and final 
    (continuum) levels and the given kind of interaction. A value::ComplexF64 is returned.
source
JenaAtomicCalculator.AutoIonization.determineChannelsMethod

AutoIonization.determineChannels(finalLevel::Level, initialLevel::Level, settings::AutoIonization.Settings) ... to determine a list of Auger Channel for a transitions from the initial to final level and by taking into account the particular settings of for this computation; an Array{AutoIonization.Channel,1} is returned.

source
JenaAtomicCalculator.AutoIonization.determineLinesMethod

AutoIonization.determineLines(finalMultiplet::Multiplet, initialMultiplet::Multiplet, settings::AutoIonization.Settings) ... to determine a list of AutoIonization.Line's for transitions between levels from the initial- and final-state multiplets, and by taking into account the particular selections and settings for this computation; an Array{AutoIonization.Line,1} is returned. Apart from the level specification, all physical properties are set to zero during the initialization process.

source
JenaAtomicCalculator.AutoIonization.displayLinesMethod

AutoIonization.displayLines(stream::IO, lines::Array{AutoIonization.Line,1}) ... to display a list of lines and channels that have been selected due to the prior settings. A neat table of all selected transitions and energies is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.AutoIonization.displayRatesMethod

AutoIonization.displayRates(stream::IO, lines::Array{AutoIonization.Line,1}, settings::AutoIonization.Settings) ... to list all results, energies, rates, etc. of the selected lines. A neat table is printed but nothing is returned otherwise.

source

Dielectronic Recombination

JenaAtomicCalculator.DielectronicRecombination.AbstractCorrectionsType

abstract type DielectronicRecombination.AbstractCorrections ... defines an abstract type to distinguish different types of corrections to the decay rates and strength. These corrections are based on the classification of shell:

    n^(core)  <   n^(final)  <  n^(hydrogenic)  <  n^(lowest-captured)  <  n^(lower-empirical)    
              <=  n^(upper-empirical)              ... where
              
    n^(core)            ... refers to the (maximum) principal quantum number to which initial core electrons are excited;
    n^(final)           ... the maximum number for which shells are treated explicitly in the representation of the final levels f;
    n^(hydrogenic)      ... to the maximum n-shell, to which the radiative decay is modeled by scaled-hydrogenic rates, 
                        ... and which can be omitted also from the list. 
    n^(lowest-captured) ... is the lowest, high-n shell, into which the additional electron is captured and which must
                            (of course) occur explicitly in the basis of the intermediate and final levels. 
    [n^(lower-empirical)  <=  n^(upper-empirical)]  
                        ... designates additional (empirical) high-n shells for which the contributions to the DR resonances 
                            are still estimated empirically by using arguments from quantum-defect theory. All shells with 
                            n > n^(upper-empirical) are neglected completely for their contributions to the DR spectra; see also:

+ struct DielectronicRecombination.EmpiricalCorrections  
    ... to estimate empirically the contributions of additional resonances for the capture of an electron into
        shells with [n^(lower-empirical)  <=  n^(upper-empirical)]. A simple scaling of the rates, calculated initially
        for n^(lowest-captured), ... only, is utilized for estimating the associated strength for these additional
        resonances.
+ struct DielectronicRecombination.HydrogenicCorrections  
    ... to add for missing final decay levels to the (total) photon decay rates by scaling the corresponding rates
        of non-relativistic hydrogenic ions with a suitable effective charge (Zeff); these hydrogenic corrections improve
        goth, the total photon rate as well as the resonance strength.
+ struct DielectronicRecombination.MaximumlCorrection  
    ... to exclude all subshells with l > l_max in the hydrogenic corrections; this restriction does not apply to the 
        given resonance levels, which can be controlled (and are specified) by the list of intermediate configurations.
source
JenaAtomicCalculator.DielectronicRecombination.EmpiricalCorrectionsType

struct DielectronicRecombination.EmpiricalCorrections <: DielectronicRecombination.AbstractCorrections ... to include empirical corrections for the shells with [n^(lower-empirical) <= n^(upper-empirical)]. A rather rude model is used so far.

+ nUpperEmpirical ::Union{Int64,Missing}   
    ... The upper-empirical shell for which rate contributions are estimated; the lower-empirical shell = n^(captured-max + 1) 
        is derived from the given configuration lists. No corrections are made for nUpperEmpirical <= n^(captured-max + 1).
+ effectiveZ      ::Union{Float64,Missing}  ... effective charge Z_eff for the hydrogenic correction (inactive).
+ rateScaling     ::Union{Float64,Missing}  ... scaling factor to modify the estimated rates.
source
JenaAtomicCalculator.DielectronicRecombination.EmpiricalTreatmentType

struct DielectronicRecombination.EmpiricalTreatment ... defines an (internal) type to communicate and distribute the physical (and technical) parameters that are utilized to make the requested empirical corrections or just nothing. This data type should not be applied by the user but is initialized by the given (set of) corrections. Otherwise, it is treated like any other type in JAC. All parameters are made physically "explicit", even if they were "missing" originally, and can be directly applied in the empirical treatment of the DR process. The following hierarchy of shells is used:

    n^(core)  <   n^(final)  <  n^(hydrogenic)  <  n^(lowest-captured)  <  n^(lower-empirical)    
              <=  n^(upper-empirical) 
    
+ doEmpiricalCorrections      ::Bool    ... True, if empirical corrections are needed, false o/w.
+ doHydrogenicCorrections     ::Bool    ... True, if hydrogenic corrections are needed, false o/w.
+ doMaximumlCorrection        ::Bool    ... True, if a maximum l values is used, false o/w.
+ doResonanceWindowCorrection ::Bool    ... True, if a window of resonances is specified, false o/w.
+ nCore                       ::Int64   
    ... (maximum) principal quantum number to which initial core electrons are excited;
+ nFinal                      ::Int64   
    ... the maximum number for which shells are treated explicitly in the representation of the final levels f;
+ nHydrogenic                 ::Int64   
    ... maximum n-shell, to which the radiative decay is modeled by scaled-hydrogenic rates. 
+ nLowestCaptured             ::Int64   
    ... lowest, high-n shell, into which the additional electron is captured and which must (of course) occur 
        explicitly in the basis of the intermediate and final levels.
+ nLowerEmpirical             ::Int64   
    ... maximum n-shell, to which the radiative decay is modeled by scaled-hydrogenic rates. 
+ nUpperEmpirical             ::Int64   
    ... additional (empirical) high-n shells for which the contributions to the DR resonances are still 
        estimated empirically by using arguments from quantum-defect theory.
+ maximum_l                   ::Int64    ... maximum l value; is set to a large value if not specified by the user.
+ hydrogenicEffectiveZ        ::Float64  ... effective charge Z_eff for the hydrogenic correction (inactive).
+ hydrogenicRateScaling       ::Float64  ... scaling factor to modify the estimated hydrogenic rates.
+ empiricalEffectiveZ         ::Float64  ... effective Z for empirical estimates
+ empiricalRateScaling        ::Float64  ... scaling factor to modify the empirical rates.
+ resonanceEnergyMin:         ::Float64  ... minimum energy [Hartree] of the resonances to be considered.
+ resonanceEnergyMax:         ::Float64  ... maximum energy [Hartree] of the resonances to be considered.
source
JenaAtomicCalculator.DielectronicRecombination.HydrogenicCorrectionsType

struct DielectronicRecombination.HydrogenicCorrections <: DielectronicRecombination.AbstractCorrections ... to add for missing final decay levels the photon decay rates for non-relativistic hydrogenic ions; this improves the total photon rate as well as the resonance strength. These corrections are taken into account for all shells with n^{final}+1 <= n <= nHydrogenic

+ nHydrogenic       ::Union{Int64,Missing}   
    ... upper principal quantum number nHydrogenic for which hydrogenic correctios to the radiative photon rates are 
        calculated explicitly; the photon rates are further scaled if some proper effectiveZ and/or rateScaling
        is provided.
+ effectiveZ      ::Union{Float64,Missing}   ... effective charge Z_eff for the hydrogenic correction.
+ rateScaling     ::Union{Float64,Missing}   ... scaling factor to scale the photon rates
source
JenaAtomicCalculator.DielectronicRecombination.MaximumlCorrectionType

struct DielectronicRecombination.MaximumlCorrection <: DielectronicRecombination.AbstractCorrections ... to exclude all subshells with l > l_max, both in the treatment of the corrections shells.

+ maximum_l    ::Union{Int64,Missing}   
    ... maximum orbital angular momentum quantum number for which contributions to the DR strengths are 
        taken into account. This number applies for all subshells for which other corrections are 
        requested, whereas the "physical subshells" are defined by the configuration lists.
source
JenaAtomicCalculator.DielectronicRecombination.PassageType

struct DielectronicRecombination.Passage ... defines a type for a dielectronic recombination passage, i.e. a (reduced) pathways, that include the definition of channels and their corresponding amplitudes for the individual i –> m resonances, whereas the subsequent radiative stabilization is considered only later.

+ initialLevel      ::Level                   ... initial-(state) level
+ intermediateLevel ::Level                   ... intermediate-(state) level
+ electronEnergy    ::Float64                 ... energy of the (incoming, captured) electron
+ captureRate       ::Float64                 ... rate for the electron capture (Auger rate)
+ photonRate        ::EmProperty              ... rate for the photon emission
+ reducedStrength   ::EmProperty              
    ... reduced resonance strength Sum_f S(i -> d -> f) * Gamma_d of this passage; this reduced strength does 
        not require the knowledge of Gamma_d for the individual passage.
+ captureChannels   ::Array{AutoIonization.Channel,1}   ... List of |i> -->  |n>   dielectronic (Auger) capture channels.
source
JenaAtomicCalculator.DielectronicRecombination.PathwayType

struct DielectronicRecombination.Pathway ... defines a type for a dielectronic recombination pathways that may include the definition of channels and their corresponding amplitudes.

+ initialLevel      ::Level                   ... initial-(state) level
+ intermediateLevel ::Level                   ... intermediate-(state) level
+ finalLevel        ::Level                   ... final-(state) level
+ electronEnergy    ::Float64                 ... energy of the (incoming, captured) electron
+ photonEnergy      ::Float64                 ... energy of the (emitted) photon
+ captureRate       ::Float64                 ... rate for the electron capture (Auger rate)
+ photonRate        ::EmProperty              ... rate for the photon emission
+ angularBeta       ::EmProperty              ... beta parameter of the photon emission
+ reducedStrength   ::EmProperty              ... reduced resonance strength S(i -> d -> f) * Gamma_d of this pathway;
                                                    this reduced strength does not require the knowledge of Gamma_d for each pathway.
+ captureChannels   ::Array{AutoIonization.Channel,1}   ... List of |i> -->  |n>   dielectronic (Auger) capture channels.
+ photonChannels    ::Array{PhotoEmission.Channel,1}    ... List of |n> -->  |f>   radiative stabilization channels.
source
JenaAtomicCalculator.DielectronicRecombination.ResonanceType

struct DielectronicRecombination.Resonance ... defines a type for a dielectronic resonance as defined by a given initial and resonance level but by summing over all final levels

+ initialLevel      ::Level             ... initial-(state) level
+ intermediateLevel ::Level             ... intermediate-(state) level
+ resonanceEnergy   ::Float64           ... energy of the resonance w.r.t. the inital-state
+ resonanceStrength ::EmProperty        ... strength of this resonance due to the stabilization into any of the allowed final levels.
+ captureRate       ::Float64           ... capture (Auger) rate to form the intermediate resonance, starting from the initial level.
+ augerRate         ::Float64           ... total (Auger) rate for an electron emission of the intermediate resonance
+ photonRate        ::EmProperty        ... total photon rate for a photon emission, i.e. for stabilization.
source
JenaAtomicCalculator.DielectronicRecombination.ResonanceSelectionType

struct DielectronicRecombination.ResonanceSelection ... defines a type for selecting classes of resonances in terms of leading configurations.

+ active          ::Bool              ... initial-(state) level
+ fromShells      ::Array{Shell,1}    ... List of shells from which excitations are to be considered.
+ toShells        ::Array{Shell,1}    ... List of shells to which (core-shell) excitations are to be considered.
+ intoShells      ::Array{Shell,1}    ... List of shells into which electrons are initially placed (captured).
source
JenaAtomicCalculator.DielectronicRecombination.ResonanceWindowCorrectionType

struct DielectronicRecombination.ResonanceWindowCorrection <: DielectronicRecombination.AbstractCorrections ... to exclude all DR resonances outside of a given "window [Emin, Emax]" of resonance energies with regard to the initial level.

+ energyMin  ::Float64   ... minimum energy [Hartree] of the resonances to be considered.  
+ energyMax  ::Float64   ... maximum energy [Hartree] of the resonances to be considered.
source
JenaAtomicCalculator.DielectronicRecombination.SettingsType

struct DielectronicRecombination.Settings <: AbstractProcessSettings ... defines a type for the details and parameters of computing dielectronic recombination pathways.

+ multipoles            ::Array{EmMultipoles}  ... Multipoles of the radiation field that are to be included.
+ gauges                ::Array{UseGauge}      ... Specifies the gauges to be included into the computations.
+ calcOnlyPassages      ::Bool                 
    ... Only compute resonance strength but without making all the pathways explicit. This option is useful
        for the capture into high-n shells or if the photons are not considered explicit. It also treats the 
        shells differently due to the given core shells < final-state shells < hydrogenically-scaled shells <
        capture-shells < asymptotic-shells. Various correction and multi-threading techiques can be applied
        to deal with or omit different classes of these shells.
+ calcRateAlpha         ::Bool                 
    ... True, if the DR rate coefficients are to be calculated, and false o/w.
+ printBefore           ::Bool                 
    ... True, if all energies and pathways are printed before their evaluation.
+ pathwaySelection      ::PathwaySelection     ... Specifies the selected levels/pathways, if any.
+ electronEnergyShift   ::Float64              
    ... An overall energy shift for all electron energies (i.e. from the initial to the resonance levels [Hartree].
+ photonEnergyShift     ::Float64              
    ... An overall energy shift for all photon energies (i.e. from the resonance to the final levels.
+ mimimumPhotonEnergy   ::Float64              
    ... minimum transition energy for which photon transitions are  included into the evaluation.
+ temperatures          ::Array{Float64,1}     
    ... list of temperatures for which plasma rate coefficients are displayed; however, these rate coefficients
        only include the contributions from those pathsways that are calculated here explicitly.
+ corrections           ::Array{DielectronicRecombination.AbstractCorrections,1}
    ... Specify, if appropriate, the inclusion of additional corrections to the rates and DR strengths.
+ augerOperator         ::AbstractEeInteraction 
    ... Auger operator that is to be used for evaluating the Auger amplitude's; the allowed values are: 
        CoulombInteraction(), BreitInteration(), CoulombBreit(), CoulombGaunt().
source
JenaAtomicCalculator.DielectronicRecombination.SettingsMethod

(set::DielectronicRecombination.Settings;

    multipoles=..,             gauges=..,                  
    calcOnlyPassages=..,       calcRateAlpha=..,         printBefore=..,           pathwaySelection=..,     
    electronEnergyShift=..,    photonEnergyShift=..,       
    mimimumPhotonEnergy=..,    temperatures=..,          corrections=..,           augerOperator=..)
                
... constructor for modifying the given DielectronicRecombination.Settings by 'overwriting' the previously selected parameters.
source
JenaAtomicCalculator.DielectronicRecombination.addEmpiricalPassages!Method

DielectronicRecombination.addEmpiricalPassages!(passages::Array{DielectronicRecombination.Passage,1}, empTreatment::EmpiricalTreatment) ... to add further (empirical) passages due to the capture of electrons into higher n-shells; information about the last "captured" shell as well as the interval [nLowerEmpirical, nUpperEmpirical] are taken from empTreatment. A simple scaling rules of energies and rates are presently applied bu could be improved if needed. The array passages::Array{DielectronicRecombination.Passage,1} is modified and nothing is returned otherwise.

source
JenaAtomicCalculator.DielectronicRecombination.checkConsistentMultipletsMethod

DielectronicRecombination.checkConsistentMultiplets(finalMultiplet::Multiplet, intermediateMultiplet::Multiplet, initialMultiplet::Multiplet) ... to check that the given initial-, intermediate- and final-state levels and multiplets are consistent to each other and to avoid later problems with the computations. An error message is issued if an inconsistency occurs, and nothing is returned otherwise.

source
JenaAtomicCalculator.DielectronicRecombination.checkOrbitalRepresentationMethod

DielectronicRecombination.checkOrbitalRepresentation(finalMultiplet::Multiplet, intermediateMultiplet::Multiplet, initialMultiplet::Multiplet) ... to check (and analyze) that all high nl orbitals in these multiplets are properly represented on the given grid. The function prints for each symmetry block kappa the high-nl orbitals and checks that they are all bound. An error message is issued if this is not the case, and nothing is returned otherwise.

source
JenaAtomicCalculator.DielectronicRecombination.computeAmplitudesPropertiesMethod

DielectronicRecombination.computeAmplitudesProperties(passage::DielectronicRecombination.Passage, finalMultiplet::Multiplet, nm::Nuclear.Model, grid::Radial.Grid, nrContinuum::Int64, empTreatment::EmpiricalTreatment, settings::DielectronicRecombination.Settings) ... to compute all amplitudes and properties of the given line; a line::DielectronicRecombination.Pathway is returned for which the amplitudes and properties have now been evaluated.

source
JenaAtomicCalculator.DielectronicRecombination.computeAmplitudesPropertiesMethod

(pathway::DielectronicRecombination.Pathway, nm::Nuclear.Model, grid::Radial.Grid, nrContinuum::Int64, settings::DielectronicRecombination.Settings, hasCaptureChannels::Bool, lastCaptureChannels::Array{AutoIonization.Channel,1}) ... to compute all amplitudes and properties of the given line; a line::DielectronicRecombination.Pathway is returned for which the amplitudes and properties have now been evaluated.

source
JenaAtomicCalculator.DielectronicRecombination.computeHydrogenicRateMethod

DielectronicRecombination.computeHydrogenicRate(ni::Int64, li::Int64, nf::Int64, lf::Int64, Zeff::Float64) ... to compute the non-relativistic electric-dipole rate for the transition from shell ni,li –> nf,lf of a hydrogenic ion with effective charge Zeff. The recursion formulas by Infeld and Hull (1951) are used together with the absorption oscillator strength. This makes the overall formulation/computation rather obscure, unfortunately. Uses SpecialFunctions.logfactorial. A rate::Float64 [a.u.] is returned. This procedure has been worked out by Stefan Schippers (2023).

source
JenaAtomicCalculator.DielectronicRecombination.computePassagesMethod

DielectronicRecombination.computePassages(finalMultiplet::Multiplet, intermediateMultiplet::Multiplet, initialMultiplet::Multiplet, nm::Nuclear.Model, grid::Radial.Grid, empTreatment::EmpiricalTreatment, settings::DielectronicRecombination.Settings) ... to compute the data for all resonances (resonance lines) directly from the given multiplets of the initial-, intermediate- and final states. It also enables one to (successively) include a set of corrections to the resonance strength to incorporate the contributions of shells that were not considered explicitly. A list of resonances::Array{DielectronicRecombination.Resonance,1} is returned.

source
JenaAtomicCalculator.DielectronicRecombination.computePathwaysMethod

DielectronicRecombination.computePathways(finalMultiplet::Multiplet, intermediateMultiplet::Multiplet, initialMultiplet::Multiplet, nm::Nuclear.Model, grid::Radial.Grid, settings::DielectronicRecombination.Settings; output=true) ... to compute the dielectronic recombination amplitudes and all properties as requested by the given settings. A list of pathways::Array{DielectronicRecombination.Pathway,1} is returned.

source
JenaAtomicCalculator.DielectronicRecombination.computeRateCoefficientMethod

DielectronicRecombination.computeRateCoefficient(resonance::DielectronicRecombination.Resonance, temp::Float64) ... computes for a delta-like resonance the DR rate coefficient alpha_d (i, Te) from the given resonance strength and temperature [K], and for both, Coulomb and Babushkin gauge. All values are directly returned in [cm^3/s]. An alphaDR::EmProperty is returned.

source
JenaAtomicCalculator.DielectronicRecombination.computeResonancesMethod

DielectronicRecombination.computeResonances(passages::Array{DielectronicRecombination.Passage,1}, settings::DielectronicRecombination.Settings) ... to compute the data for all resonances (resonance lines) as defined by the given passages and and settings. In fact, passages and resonances are treated rather similar to each other in JAC and can be readily extracted from each other. A list of resonances::Array{DielectronicRecombination.Resonance,1} is returned.

source
JenaAtomicCalculator.DielectronicRecombination.computeResonancesMethod

(pathways::Array{DielectronicRecombination.Pathway,1}, settings::DielectronicRecombination.Settings) ... to compute the data for all resonances (resonance lines) as defined by the given pathways and and settings. A list of resonances::Array{DielectronicRecombination.Resonance,1} is returned.

source
JenaAtomicCalculator.DielectronicRecombination.determineCaptureChannelsMethod

DielectronicRecombination.determineCaptureChannels(intermediateLevel::Level, initialLevel::Level, settings::DielectronicRecombination.Settings) ... to determine a list of AutoIonization.Channel for a (Auger) capture transitions from the initial to an intermediate level, and by taking into account the particular settings of for this computation; an Array{AutoIonization.Channel,1} is returned.

source
JenaAtomicCalculator.DielectronicRecombination.determineEmpiricalTreatmentMethod

DielectronicRecombination.determineEmpiricalTreatment(finalMultiplet::Multiplet, intermediateMultiplet::Multiplet, nm::Nuclear.Model, initialMultiplet::Multiplet, settings::DielectronicRecombination.Settings) ... to determine an instance of empiricalTreatment::EmpiricalTreatment that is (internally) applied to simplify the use of "corrections" to the DR strenghts. This data structure summarizes all parameters that help introduce several empirical corrections. The procedure is simple but slightly sophisticated as we wish to support "missing" parameters in the individual corrections as well as the knowledge that can be derived internally. The definition of EmpiricalTreatment() can readily be extended as the need arises from the user side.

source
JenaAtomicCalculator.DielectronicRecombination.determinePassagesMethod

DielectronicRecombination.determinePassages(intermediateMultiplet::Multiplet, initialMultiplet::Multiplet, empTreatment::EmpiricalTreatment, settings::DielectronicRecombination.Settings) ... to determine a list of dielectronic-recombination resonances between the levels from the given initial- and intermediate- states, whereas the final states are considered "on-fly"; the particular selections and settings for this computation are taken into account; an Array{DielectronicRecombination.Passsage,1} is returned. Apart from the level specification, all physical properties are set to zero during the initialization process.

source
JenaAtomicCalculator.DielectronicRecombination.determinePathwaysMethod

DielectronicRecombination.determinePathways(finalMultiplet::Multiplet, intermediateMultiplet::Multiplet, initialMultiplet::Multiplet, settings::DielectronicRecombination.Settings) ... to determine a list of dielectronic-recombination pathways between the levels from the given initial-, intermediate- and final-state multiplets and by taking into account the particular selections and settings for this computation; an Array{DielectronicRecombination.Pathway,1} is returned. Apart from the level specification, all physical properties are set to zero during the initialization process.

source
JenaAtomicCalculator.DielectronicRecombination.determinePhotonChannelsMethod

DielectronicRecombination.determinePhotonChannels(finalLevel::Level, intermediateLevel::Level, settings::DielectronicRecombination.Settings) ... to determine a list of PhotoEmission.Channel for the photon transitions from the intermediate and to a final level, and by taking into account the particular settings of for this computation; an Array{PhotoEmission.Channel,1} is returned.

source
JenaAtomicCalculator.DielectronicRecombination.displayPassagesMethod

DielectronicRecombination.displayPassages(stream::IO, passages::Array{DielectronicRecombination.Passage,1}) ... to display a list of passages and channels that have been selected due to the prior settings. A neat table of all selected transitions and energies is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.DielectronicRecombination.displayPathwaysMethod

DielectronicRecombination.displayPathways(stream::IO, pathways::Array{DielectronicRecombination.Pathway,1}) ... to display a list of pathways and channels that have been selected due to the prior settings. A neat table of all selected transitions and energies is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.DielectronicRecombination.displayRateCoefficientsMethod

DielectronicRecombination.displayRateCoefficients(stream::IO, resonances::Array{DielectronicRecombination.Resonance,1}, settings::DielectronicRecombination.Settings) ... to list, if settings.calcRateAlpha, all rate coefficients for the selected temperatures. Both, the individual as well as the total DR plasma rate coefficients are printed in neat tables, though nothing is returned otherwise.

source
JenaAtomicCalculator.DielectronicRecombination.displayResultsMethod

DielectronicRecombination.displayResults(stream::IO, pathways::Array{DielectronicRecombination.Pathway,1}, settings::DielectronicRecombination.Settings) ... to list all results, energies, cross sections, etc. of the selected lines. A neat table is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.DielectronicRecombination.extractRateCoefficientsMethod

DielectronicRecombination.extractRateCoefficients(resonances::Array{DielectronicRecombination.Resonance,1}, settings::DielectronicRecombination.Settings) ... to extract, if settings.calcRateAlpha, the total DR rate coefficients for all temperatures. A list of total rate coefficients [cm^3/s] alphaDR::Array{EmProperty,1} is returned.

source
JenaAtomicCalculator.DielectronicRecombination.isResonanceToBeExcludedMethod

DielectronicRecombination.isResonanceToBeExcluded(level::Level, refLevel, rSelection::ResonanceSelection) returns true, if level is to be excluded from the valid resonances, and false otherwise. It returns false if the ResonanceSelection() is inactive or if level belongs to the selected resoances. It is true only of ResonanceSelection() is active but the level does not belong to the selected resonances.

source

Electron Impact Excitation

JenaAtomicCalculator.ImpactExcitation.ChannelType

struct ImpactExcitation.Channel ... defines a type for a electron-impact excitaiton channel to help characterize the incoming and outgoing (continuum) states of many electron-states with a single free electron

+ initialKappa     ::Int64              ... partial-wave of the incoming free electron
+ finalKappa       ::Int64              ... partial-wave of the outgoing free electron
+ symmetry         ::LevelSymmetry      ... total angular momentum and parity of the scattering state
+ initialPhase     ::Float64            ... phase of the incoming partial wave
+ finalPhase       ::Float64            ... phase of the outgoing partial wave
+ amplitude        ::Complex{Float64}   ... Collision amplitude associated with the given channel.
source
JenaAtomicCalculator.ImpactExcitation.LineType

struct ImpactExcitation.Line ... defines a type for a electron-impact excitation line that may include the definition of channels and their corresponding amplitudes.

+ initialLevel           ::Level         ... initial- (bound-state) level
+ finalLevel             ::Level         ... final- (bound-state) level
+ initialElectronEnergy  ::Float64       ... energy of the incoming (initial-state) free-electron
+ finalElectronEnergy    ::Float64       ... energy of the outgoing (final-state) free-electron
+ crossSection           ::Float64       ... total cross section of this line
+ collisionStrength      ::Float64       ... total collision strength of this line
+ channels               ::Array{ImpactExcitation.Channel,1}  ... List of ImpactExcitation channels of this line.
+ convergence            ::Float64       ... convergence of calculation
source
JenaAtomicCalculator.ImpactExcitation.LineMethod

ImpactExcitation.Line(initialLevel::Level, finalLevel::Level, crossSection::Float64) ... constructor for an electron-impact excitation line between a specified initial and final level.

source
JenaAtomicCalculator.ImpactExcitation.RateCoefficientsType

struct ImpactExcitation.RateCoefficients ... Defines a type for the output results from excitation rate or effective collision strengths calculations

+ initialLevel        ::Level               ... initial- (bound-state) level
+ finalLevel          ::Level               ... final- (bound-state) level
+ temperatures        ::Array{Float64,1}    ... Temperatures in [K] to calculate excitation rates and effective collision strengths
+ alphas              ::Array{Float64,1}    ... Excitation rate coefficients in [cm^3/s] for the input temperatures
+ effOmegas           ::Array{Float64,1}    ... Effective collision strengths for the input temperatures
source
JenaAtomicCalculator.ImpactExcitation.SettingsType

struct ImpactExcitation.Settings <: AbstractProcessSettings ... defines a type for the details and parameters of computing electron-impact excitation lines.

+ lineSelection            ::LineSelection      ... Specifies the selected levels, if any.
+ electronEnergies         ::Array{Float64,1}   ... List of impact-energies of the incoming elecgtrons (in user-defined units).
+ energyShift              ::Float64            ... An overall energy shift for all transitions |i> --> |f>.
+ maxKappa                 ::Int64              ... Maximum kappa value of partial waves to be included.
+ calcRateCoefficient      ::Bool               ... True, if the plasma rate coefficients to be calculated, false otherwise.
+ maxEnergyMultiplier      ::Float64            ... Maximum initial electron energy for eff. collision strength integration.
                                 (maxEnergyMultiplier * Excitation threshold energy), after this assymptotic limit is applied.
+ numElectronEnergies      ::Int64              ... No. of different electron energy points at which collision strengths to compute
                                                    in electron energy range [0, (maxEnergyMultiplier * Excitation threshold energy)]
+ temperatures             ::Array{Float64, 1}  ... Electron temperatures [K] for eff. collision strengths and rate coefficients.
+ printBefore              ::Bool               ... True, if all energies and lines are printed before their evaluation.
+ operator                 ::AbstractEeInteraction   
    ... Interaction operator that is to be used for evaluating the e-e interaction amplitudes; allowed values are: 
        CoulombInteraction(), BreitInteraction(), ...
source
JenaAtomicCalculator.ImpactExcitation.SettingsMethod

ImpactExcitation.Settings(set::ImpactExcitation.Settings;

lineSelection...,  electronEnergies..., energyShift..., maxKappa...,
calcRateCoefficient..., maxEnergyMultiplier..., numElectronEnergies..., temperatures...,
printBefore..., operator...)

... constructor for modifying the given ImpactExcitation.Settings by 'overwriting' the previously selected parameters.
source
JenaAtomicCalculator.ImpactExcitation.amplitudeMethod

ImpactExcitation.amplitude(kind::AbstractEeInteraction, channel::ImpactExcitation.Channel, cFinalLevel::Level, cInitialLevel::Level, grid::Radial.Grid; printout::Bool=true) ... to compute the kind in CoulombInteraction(), BreitInteraction(), CoulombBreit() electron-impact interaction amplitude <(alphaf Jf, kappaf) Jt || O^(e-e, kind) || (alphai Ji, kappai) Jt> due to the interelectronic interaction for the given final and initial (continuum) level. A value::ComplexF64 is returned.

source
JenaAtomicCalculator.ImpactExcitation.computeAmplitudesPropertiesMethod

ImpactExcitation.computeAmplitudesProperties(line::ImpactExcitation.Line, nm::Nuclear.Model, grid::Radial.Grid, settings::ImpactExcitation.Settings; printout::Bool=true) ... to compute all amplitudes and properties of the given line; a line::ImpactExcitation.Line is returned for which the amplitudes and properties are now evaluated.

source
JenaAtomicCalculator.ImpactExcitation.computeEffStrengthsMethod

ImpactExcitation.computeEffStrengths(lines::Array{ImpactExcitation.Line, 1}, settings::ImpactExcitation.Settings) ... computes Effective collision strengths from the calculated line collision strengths at temperature(s) [K]. Returns an Array{ImpactExcitation.RateCoefficients,1}.

source
JenaAtomicCalculator.ImpactExcitation.computeLinesMethod

ImpactExcitation.computeLines(finalMultiplet::Multiplet, initialMultiplet::Multiplet, nm::Nuclear.Model, grid::Radial.Grid, settings::ImpactExcitation.Settings; output=true) ... to compute the electron-impact excitation transition amplitudes and all properties as requested by the given settings. A list of lines::Array{ImpactExcitation.Lines} is returned.

source
JenaAtomicCalculator.ImpactExcitation.computeRateCoefficientsMethod

ImpactExcitation.computeRateCoefficients(effStrengths::Vector{RateCoefficients}, settings::ImpactExcitation.Settings) ... computes Exitation rate coefficients from the calculated collsion strengths at a temperature(s) [K]. The rate coefficients are returned in [cm^3/s]. Returns an Array{ImpactExcitation.RateCoefficients,1}.

source
JenaAtomicCalculator.ImpactExcitation.determineChannelsMethod

ImpactExcitation.determineChannels(finalLevel::Level, initialLevel::Level, settings::ImpactExcitation.Settings) ... to determine a list of electron-impact excitation Channels for a transitions from the initial to the final level and by taking into account the particular settings of for this computation; an Array{ImpactExcitation.Channel,1} is returned.

source
JenaAtomicCalculator.ImpactExcitation.determineLinesMethod

ImpactExcitation.determineLines(finalMultiplet::Multiplet, initialMultiplet::Multiplet, settings::ImpactExcitation.Settings) ... to determine a list of ImpactExcitation.Line's for transitions between levels from the initial- and final-state multiplets, and by taking into account the particular selections and settings for this computation; an Array{ImpactExcitation.Line,1} is returned. Apart from the level specification, all physical properties are set to zero during the initialization process.

source
JenaAtomicCalculator.ImpactExcitation.displayLinesMethod

ImpactExcitation.displayLines(lines::Array{ImpactExcitation.Line,1}) ... to display a list of lines and channels that have been selected due to the prior settings. A neat table of all selected transitions and energies is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.ImpactExcitation.displayResultsMethod

ImpactExcitation.displayResults(allRates::Array{RateCoefficients,1}) ... to list the excitation rate coefficients and effective collision strengths for the selected lines at the selected temperatures. A neat table is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.ImpactExcitation.groupLinesMethod

ImpactExcitation.groupLines(lines::Array{ImpactExcitation.Line,1}, settings::ImpactExcitation.Settings) ... groups lines having the same initial and final level but different energies returns an Array{Array{ImpactExcitation.Line,1},1}

source
JenaAtomicCalculator.ImpactExcitation.interpolateCSMethod

ImpactExcitation.interpolateCS(x::Float64, xa::Vector{Float64}, ya::Vector{Float64}, isE1Allowed::Bool) ... interpolates the cross section or collision stregths for a given initial electron Energy if the energy is beyond the upper bound then uses asymptotic approximation for extrapolation, depending on the transition is E1 allowed or not. Returns a Float64.

source

Photoemission

JenaAtomicCalculator.PhotoEmission.ChannelType

struct PhotoEmission.Channel ... defines a type for a single radiative emission/absorption channel that specifies the multipole, gauge and amplitude.

+ multipole         ::EmMultipole        ... Multipole of the photon emission/absorption.
+ gauge             ::EmGauge            ... Gauge for dealing with the (coupled) radiation field.
+ amplitude         ::Complex{Float64}   ... Amplitude of this multiple channel.
source
JenaAtomicCalculator.PhotoEmission.LineType

struct PhotoEmission.Line ... defines a type for a radiative line that may include the definition of sublines and their corresponding amplitudes.

+ initialLevel   ::Level               ... initial-(state) level
+ finalLevel     ::Level               ... final-(state) level
+ omega          ::Float64             ... Transition frequency of this line; can be shifted w.r.t. the level energies.
+ photonRate     ::EmProperty          ... Total rate of this line.
+ angularBeta    ::EmProperty          ... Angular beta_2 coefficient.
+ hasSublines    ::Bool                ... Determines whether the sublines are defined in terms of their multipolarity, amplitude, or not.
+ channels       ::Array{PhotoEmission.Channel,1}  ... List of radiative (photon) channels
source
JenaAtomicCalculator.PhotoEmission.SettingsType

struct PhotoEmission.Settings <: AbstractProcessSettings ... defines a type for the details and parameters of computing radiative lines.

+ multipoles              ::Array{EmMultipoles}     ... Specifies the (radiat. field) multipoles to be included.
+ gauges                  ::Array{UseGauge}         ... Gauges to be included into the computations.
+ calcAnisotropy          ::Bool                    ... True, if the anisotropy (structure) functions are to be 
                                                        calculated and false otherwise 
+ printBefore             ::Bool                    ... True, if all energies and lines are printed before comput.
+ corePolarization        ::CorePolarization        ... Parametrization of the core-polarization potential/contribution.
+ lineSelection           ::LineSelection           ... Specifies the selected levels, if any.
+ photonEnergyShift       ::Float64                 ... An overall energy shift for all photon energies.
+ mimimumPhotonEnergy     ::Float64                 ... minimum transition energy for which (photon) transitions 
                                                        are included into the computation.
+ maximumPhotonEnergy     ::Float64                 ... maximum transition energy for which (photon) transitions 
                                                        are included.
source
JenaAtomicCalculator.PhotoEmission.SettingsMethod

PhotoEmission.Settings(set::PhotoEmission.Settings;

    multipoles::=..,        gauges=..,                calcAnisotropy=..,          printBefore=..,
    corePolarization=..,    lineSelection=..,         photonEnergyShift=..,       
    mimimumPhotonEnergy=.., maximumPhotonEnergy=..) 
                
... constructor for modifying the given PhotoEmission.Settings by 'overwriting' the previously selected parameters.
source
JenaAtomicCalculator.PhotoEmission.amplitudeMethod

+ (kind::String, cp::CorePolarization, omega::Float64, finalLevel::Level, initialLevel::Level, grid::Radial.Grid; display::Bool=false, printout::Bool=false) ... to compute the kind = E1 with core-polarization emission amplitude <alphaf Jf || O^(E1, emission with core-polarization) || alphai Ji> in length gauge and for the given transition energy. A value::ComplexF64 is returned. The amplitude value is printed to screen if display=true.

source
JenaAtomicCalculator.PhotoEmission.amplitudeMethod

PhotoEmission.amplitude(kind::String, Mp::EmMultipole, gauge::EmGauge, omega::Float64, finalLevel::Level, initialLevel::Level, grid::Radial.Grid; display::Bool=false, printout::Bool=false) ... to compute the kind = (absorption or emission) amplitude <alphaf Jf || O^(Mp, kind) || alphai Ji> for the interaction with photon of multipolarity Mp and for the given transition energy and gauge. A value::ComplexF64 is returned. The amplitude value is printed to screen if display=true.

source
JenaAtomicCalculator.PhotoEmission.amplitude_WuMethod

+ amplitude_Wu(kind::String, Mp::EmMultipole, gauge::EmGauge, omega::Float64, finalLevel::Level, initialLevel::Level, grid::Radial.Grid; display::Bool=false, printout::Bool=false) ... to compute the kind = (absorption or emission) amplitude <alphaf Jf || O^(Mp, kind) || alphai Ji> for the interaction with photon of multipolarity Mp and for the given transition energy and gauge. A value::ComplexF64 is returned. The radial function is calculated by InteractionStrength.MabEmissionJohnsony_Wu. The amplitude value is printed to screen if display=true.

source
JenaAtomicCalculator.PhotoEmission.computeAmplitudesPropertiesMethod

PhotoEmission.computeAmplitudesProperties(line::PhotoEmission.Line, grid::Radial.Grid, settings::Einstein.Settings; printout::Bool=true) ... to compute all amplitudes and properties of the given line; a line::Einstein.Line is returned for which the amplitudes and properties are now evaluated.

source
JenaAtomicCalculator.PhotoEmission.computeLinesMethod

PhotoEmission.computeLines(finalMultiplet::Multiplet, initialMultiplet::Multiplet, grid::Radial.Grid, settings::PhotoEmission.Settings; output=true) ... to compute the radiative transition amplitudes and all properties as requested by the given settings. A list of lines::Array{PhotoEmission.Lines} is returned.

source
JenaAtomicCalculator.PhotoEmission.computeLinesCascadeMethod

PhotoEmission.computeLinesCascade(finalMultiplet::Multiplet, initialMultiplet::Multiplet, grid::Radial.Grid, settings::PhotoEmission.Settings; output::Bool=true, printout::Bool=true) ... to compute the radiative transition amplitudes and all properties as requested by the given settings. The computations and printout is adapted for larger cascade computations by including only lines with at least one channel and by sending all printout to a summary file only. A list of lines::Array{PhotoEmission.Lines} is returned.

source
JenaAtomicCalculator.PhotoEmission.determineChannelsMethod

PhotoEmission.determineChannels(finalLevel::Level, initialLevel::Level, settings::PhotoEmission.Settings) ... to determine a list of PhotoEmission.Channel for a transitions from the initial to final level and by taking into account the particular settings of for this computation; an Array{PhotoEmission.Channel,1} is returned.

source
JenaAtomicCalculator.PhotoEmission.determineLinesMethod

PhotoEmission.determineLines(finalMultiplet::Multiplet, initialMultiplet::Multiplet, settings::PhotoEmission.Settings) ... to determine a list of PhotoEmission Line's for transitions between the levels from the given initial- and final-state multiplets and by taking into account the particular selections and settings for this computation; an Array{PhotoEmission.Line,1} is returned. Apart from the level specification, all physical properties are set to zero during the initialization process.

source
JenaAtomicCalculator.PhotoEmission.displayAnisotropiesMethod

PhotoEmission.displayAnisotropies(stream::IO, lines::Array{PhotoEmission.Line,1}, settings::PhotoEmission.Settings) ... to list all energies and anisotropy parameters of the selected lines. A neat table is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.PhotoEmission.displayLifetimesMethod

PhotoEmission.displayLifetimes(stream::IO, lines::Array{PhotoEmission.Line,1}, settings::PhotoEmission.Settings) ... to list all lifetimes as derived from the selected lines. A neat table is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.PhotoEmission.displayLinesMethod

PhotoEmission.displayLines(stream::IO, lines::Array{PhotoEmission.Line,1}) ... to display a list of lines and channels that have been selected due to the prior settings. A neat table of all selected transitions and energies is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.PhotoEmission.displayRatesMethod

PhotoEmission.displayRates(stream::IO, lines::Array{PhotoEmission.Line,1}, settings::PhotoEmission.Settings) ... to list all results, energies, rates, etc. of the selected lines. A neat table is printed but nothing is returned otherwise.

source

Photoexcitation

JenaAtomicCalculator.PhotoExcitation.LineType

struct PhotoExcitation.Line ... defines a type for a photo-excitation line that may include the definition of sublines and their corresponding amplitudes.

+ initialLevel   ::Level                       ... initial-(state) level
+ finalLevel     ::Level                       ... final-(state) level
+ omega          ::Float64                     ... Transition frequency of this line; can be shifted w.r.t. the level energies.
+ oscStrength    ::EmProperty                  ... Absorption oscillator strength
+ crossSection   ::EmProperty                  ... Total cross section of this line.
+ staTensor      ::Array{TensorComp,1}         ... Array of statistical tensor components rho_kq
+ hasSublines    ::Bool                        ... Determines whether the individual sublines are defined in terms of their 
                                                    multipolarity, amplitude, or not; cf. PhotoEmission.Channel
+ channels       ::Array{PhotoEmission.Channel,1}  ... List of radiative (photon) channels
source
JenaAtomicCalculator.PhotoExcitation.LineMethod

PhotoExcitation.Line(initialLevel::Level, finalLevel::Level, omega::Float64, crossSection::EmProperty) ... constructor an photo-excitation line between a specified initial and final level.

source
JenaAtomicCalculator.PhotoExcitation.SettingsType

struct PhotoExcitation.Settings <: AbstractProcessSettings ... defines a type for the details and parameters of computing photo-excitation lines.

+ multipoles              ::Array{EmMultipole,1}    ... Specifies the multipoles of the radiation field that are to be included.
+ gauges                  ::Array{UseGauge,1}       ... Specifies the gauges to be included into the computations.
+ calcForStokes           ::Bool                    ... True, if the excitation cross sections are to be calculated (and false otherwise)
                                                        for given Stokes parameter of the incident plane-wave photons.
+ calcPhotonDm            ::Bool                    ... True, if the photon density matrix of a subsequently emitted fluorescence photon 
                                                        is to be calculated and false otherwise. 
+ calcTensors             ::Bool                    ... True, if statistical tensors of the excited atom are to be calculated, false otherwise. 
+ printBefore             ::Bool                    ... True, if all energies and lines are printed before their evaluation.
+ lineSelection           ::LineSelection           ... Specifies the selected levels, if any.
+ photonEnergyShift       ::Float64                 ... An overall energy shift for all photon energies.
+ mimimumPhotonEnergy     ::Float64                 ... minimum transition energy for which (photon) transitions are included into the
                                                        computation.
+ maximumPhotonEnergy     ::Float64                 ... maximum transition energy for which (photon) transitions are included.
+ stokes                  ::ExpStokes               ... Stokes parameters of the incident radiation.
source
JenaAtomicCalculator.PhotoExcitation.SettingsMethod

PhotoExcitation.Settings(set::PhotoExcitation.Settings;

    multipoles=..,          gauges=..,                  calcForStokes=..,           calcPhotonDm=..,    
    calcTensors=..,         printBefore=..,             lineSelection=..,    
    photonEnergyShift=..,   mimimumPhotonEnergy=..,     maximumPhotonEnergy=..,     stokes=..)
                
... constructor for modifying the given PhotoExcitation.Settings by 'overwriting' the previously selected parameters.
source
JenaAtomicCalculator.PhotoExcitation.computeAmplitudesPropertiesMethod

PhotoExcitation.computeAmplitudesProperties(line::PhotoExcitation.Line, grid::Radial.Grid, settings::PhotoExcitation.Settings; printout::Bool=true) ... to compute all amplitudes and properties of the given line; a line::PhotoExcitation.Line is returned for which the amplitudes and properties have now been evaluated.

source
JenaAtomicCalculator.PhotoExcitation.computeCrossSectionMethod

PhotoExcitation.computeCrossSection(line::PhotoExcitation.Line, stokes::ExpStokes) ... to compute the excitation cross section for the excitation of unpolarized atoms by plane-wave photons, whose polarization is described by the given (experimental) Stokes parameters. A cs::EmProperty is returned.

source
JenaAtomicCalculator.PhotoExcitation.computeLinesMethod

PhotoExcitation.computeLines(finalMultiplet::Multiplet, initialMultiplet::Multiplet, grid::Radial.Grid, settings::PhotoExcitation.Settings; output=true) ... to compute the photo-excitation amplitudes and all properties as requested by the given settings. A list of lines::Array{PhotoExcitation.Lines} is returned.

source
JenaAtomicCalculator.PhotoExcitation.computeLinesCascadeMethod

PhotoExcitation.computeLinesCascade(finalMultiplet::Multiplet, initialMultiplet::Multiplet, grid::Radial.Grid, settings::PhotoExcitation.Settings, initialLevelSelection::LevelSelection; output::Bool=true, printout::Bool=true) ... to compute the excitation (absorption) transition amplitudes and all properties as requested by the given settings. The computations and printout is adapted for larger cascade computations by including only lines with at least one channel and by sending all printout to a summary file only. A list of lines::Array{PhotoExcitation.Lines} is returned.

source
JenaAtomicCalculator.PhotoExcitation.computeStatisticalTensorMethod

PhotoExcitation.computeStatisticalTensor(k::Int64, q::Int64, line::PhotoExcitation.Line, stokes::ExpStokes) ... to compute the statistical tensor (component) rho{k,q} of the final level for the excitation of unpolarized atoms by plane-wave photons, whose polarization is described by the given (experimental) Stokes parameters. A rhokq::EmProperty is returned.

source
JenaAtomicCalculator.PhotoExcitation.determineChannelsMethod

PhotoExcitation.determineChannels(finalLevel::Level, initialLevel::Level, settings::PhotoExcitation.Settings) ... to determine a list of PhotoExcitation.Channel for a transitions from the initial to final level and by taking into account the particular settings of for this computation; an Array{PhotoExcitation.Channel,1} is returned.

source
JenaAtomicCalculator.PhotoExcitation.determineLinesMethod

PhotoExcitation.determineLines(finalMultiplet::Multiplet, initialMultiplet::Multiplet, settings::PhotoExcitation.Settings) ... to determine a list of photo-excitation Line's for transitions between the levels from the given initial- and final-state multiplets and by taking into account the particular selections and settings for this computation; an Array{PhotoExcitation.Line,1} is returned. Apart from the level specification, all physical properties are set to zero during the initialization process.

source
JenaAtomicCalculator.PhotoExcitation.displayCrossSectionsMethod

PhotoExcitation.displayCrossSections(stream::IO, lines::Array{PhotoExcitation.Line,1}, settings::PhotoExcitation.Settings) ... to list all results, energies, cross sections, etc. of the selected lines. A neat table is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.PhotoExcitation.displayLineDataMethod

PhotoExcitation.displayLineData(stream::IO, lines::Array{PhotoExcitation.Line,1}) ... to display the calculated data, ordered by the initial levels and the photon energies involved. Neat tables of all initial levels and photon energies as well as all associated total cross sections are printed but nothing is returned otherwise.

source
JenaAtomicCalculator.PhotoExcitation.displayLinesMethod

PhotoExcitation.displayLines(stream::IO, lines::Array{PhotoExcitation.Line,1}) ... to display a list of lines and channels that have been selected due to the prior settings. A neat table of all selected transitions and energies is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.PhotoExcitation.estimateCrossSectionMethod

PhotoExcitation.estimateCrossSection(lines::Array{PhotoExcitation.Line,1}, omega::Float64, gamma::Float64, initialLevel) ... to estimate from lines the total PE cross section for any omega and for initial level. The procedure assumes a Gaussian line shape for each PhotoExcitation.Line with widths gamma and distributes the (total cross section according to the line shape within [Er - 10*gamma, Er + 10*gamma]. Other line shapes can be readily implemented; a cross section cs::EmProperty is returned.

source

Photoionization

JenaAtomicCalculator.PhotoIonization.ChannelType

struct PhotoIonization.Channel ... defines a type for a photoionization channel to help characterize a single multipole and scattering (continuum) state of many electron-states with a single free electron.

+ multipole      ::EmMultipole          ... Multipole of the photon absorption.
+ gauge          ::EmGauge              ... Gauge for dealing with the (coupled) radiation field.
+ kappa          ::Int64                ... partial-wave of the free electron
+ symmetry       ::LevelSymmetry        ... total angular momentum and parity of the scattering state
+ phase          ::Float64              ... phase of the partial wave
+ amplitude      ::Complex{Float64}     ... Photoionization amplitude associated with the given channel.
source
JenaAtomicCalculator.PhotoIonization.LineType

struct Line ... defines a type for a photoionization line that may include the definition of channels.

+ initialLevel   ::Level                  ... initial-(state) level
+ finalLevel     ::Level                  ... final-(state) level
+ electronEnergy ::Float64                ... Energy of the (outgoing free) electron.
+ photonEnergy   ::Float64                ... Energy of the absorbed photon.
+ crossSection   ::EmProperty             ... Cross section for this photoionization.
+ angularBeta    ::EmProperty             ... beta -parameter for unpolarized targets with J=0, 1/2, 1
+ coherentDelay  ::EmProperty             ... coherent time-delay due to the selected averaging of phases.
+ incoherentDelay::EmProperty             ... incoherent time-delay due to the selected averaging of phases.
+ channels       ::Array{PhotoIonization.Channel,1}  ... List of PhotoIonization.Channels of this line.
source
JenaAtomicCalculator.PhotoIonization.PlasmaSettingsType

struct PhotoIonization.PlasmaSettings <: Basics.AbstractLineShiftSettings ... defines a type for the details and parameters of computing photoionization rates with plasma interactions.

+ multipoles             ::Array{Basics.EmMultipole}     ... Specifies the multipoles of the radiation field that are to be included.
+ gauges                 ::Array{Basics.UseGauge}        ... Specifies the gauges to be included into the computations.
+ photonEnergies         ::Array{Float64,1}              ... List of photon energies.  
+ printBefore            ::Bool                          ... True, if all energies and lines are printed before their evaluation.
+ lineSelection          ::LineSelection                 ... Specifies the selected levels, if any.
source
JenaAtomicCalculator.PhotoIonization.SettingsType

struct PhotoIonization.Settings <: AbstractProcessSettings ... defines a type for the details and parameters of computing photoionization lines.

+ multipoles                    ::Array{EmMultipole}  ... Specifies the multipoles of the radiation field that are to be included.
+ gauges                        ::Array{UseGauge}     ... Specifies the gauges to be included into the computations.
+ photonEnergies                ::Array{Float64,1}    ... List of photon energies [in user-selected units].  
+ electronEnergies              ::Array{Float64,1}    ... List of electron energies; usually only one of these lists are utilized. 
+ thetas                        ::Array{Float64,1}    ... List of theta-values if angle-differential CS are calculated explicitly. 
+ calcAnisotropy                ::Bool                ... True, if the beta anisotropy parameters are to be calculated and false otherwise (o/w). 
+ calcPartialCs                 ::Bool                ... True, if partial cross sections are to be calculated and false otherwise.  
+ calcTimeDelay                 ::Bool                ... True, if time-delays are to be calculated and false otherwise.  
+ calcNonE1AngleDifferentialCS  ::Bool                ... True, if non-E1 angle-differential CS are be calculated and false otherwise.  
+ calcTensors                   ::Bool                ... True, if statistical tensors of the excited atom are to be calculated and false o/w. 
+ printBefore                   ::Bool                ... True, if all energies and lines are printed before their evaluation.
+ lineSelection                 ::LineSelection       ... Specifies the selected levels, if any.
+ stokes                        ::ExpStokes           ... Stokes parameters of the incident radiation.
+ freeElectronShift             ::Float64             ... An overall energy shift of all free-electron energies [user-specified units].
+ lValues                       ::Array{Int64,1}      ... Orbital angular momentum of free-electrons, for which partial waves are considered.
source
JenaAtomicCalculator.PhotoIonization.SettingsMethod

PhotoIonization.Settings(set::PhotoIonization.Settings;

    multipoles=..,                      gauges=..,                  photonEnergies=..,          electronEnergies=..,     
    thetas=..,                          calcAnisotropy=..,          calcPartialCs..,            calcTimeDelay=..,           
    calcNonE1AngleDifferentialCS=..,    calcTensors=..,             printBefore=..,             lineSelection=..,           
    stokes=..,                          freeElectronShift=..,       lValues=.. )
                
... constructor for modifying the given PhotoIonization.Settings by 'overwriting' the previously selected parameters.
source
JenaAtomicCalculator.PhotoIonization.amplitudeMethod

PhotoIonization.amplitude(kind::String, channel::PhotoIonization.Channel, omega::Float64, continuumLevel::Level, initialLevel::Level, grid::Radial.Grid) ... to compute the kind = (photoionization) amplitude <(alphaf Jf, epsilon kappa) Jt || O^(photoionization) || alphai J_i> due to the electron-photon interaction for the given final and initial level, the partial wave of the outgoing electron as well as the given multipole and gauge. A value::ComplexF64 is returned.

source
JenaAtomicCalculator.PhotoIonization.angularFunctionKMethod

PhotoIonization.angularFunctionK(L1::Int64, L2::Int64, X::Int64, Ji::AngularJ64, Jf::AngularJ64, kappa1::Int64, J1::AngularJ64, kappa2::Int64, J2::AngularJ64) ... to compute angular function K(...) as defined for the non-E1 angle-differential cross sections by Nishita Hosea (2025). No tests are made that the triangular conditions of the quantum numbers are fulfilled. A wa::Float64 is returned.

source
JenaAtomicCalculator.PhotoIonization.angularFunctionWMethod

PhotoIonization.angularFunctionW(theta::Float64, L1::Int64, L2::Int64, X::Int64, lambda1::Int64, lambda2::Int64, kappa1::Int64, mu1::Rational{Int64}, kappa2::Int64, mu2::Rational{Int64}) ... to compute angular function W(theta; ...) as defined for the non-E1 angle-differential cross sections by Nishita Hosea (2025). No tests are made that the triangular conditions of the quantum numbers are fulfilled. A wa::Float64 is returned.

source
JenaAtomicCalculator.PhotoIonization.computeAmplitudesPropertiesMethod

PhotoIonization.computeAmplitudesProperties(line::PhotoIonization.Line, nm::Nuclear.Model, grid::Radial.Grid, nrContinuum::Int64, settings::PhotoIonization.Settings; printout::Bool=false) ... to compute all amplitudes and properties of the given line; a line::PhotoIonization.Line is returned for which the amplitudes and properties are now evaluated.

source
JenaAtomicCalculator.PhotoIonization.computeAmplitudesPropertiesPlasmaMethod

PhotoIonization.computeAmplitudesPropertiesPlasma(line::PhotoIonization.Line, nm::Nuclear.Model, grid::Radial.Grid, settings::PhotoIonization.PlasmaSettings) ... to compute all amplitudes and properties of the given line but for the given plasma model; a line::PhotoIonization.Line is returned for which the amplitudes and properties are now evaluated.

source
JenaAtomicCalculator.PhotoIonization.computeAngularBetaMethod

PhotoIonization.computeAngularBeta(iLevel::Level, fLevel::Level, channels::Array{PhotoIonization.Channel,1}) ... to compute the beta anisotropy parameter for the photoionization transition i -> f with the given channels; here, the formula from Balashov (1994, Eq. 2.135) has been utilized. A beta::EmProperty parameter is returned. These (gauge-dependent) beta parameters are set to -9., if no amplitudes are calculated for the given gauge.

source
JenaAtomicCalculator.PhotoIonization.computeDisplayNonE1AngleDifferentialCSMethod

PhotoIonization.computeDisplayNonE1AngleDifferentialCS(stream::IO, lines::Array{PhotoIonization.Line,1}, settings::PhotoIonization.Settings) ... to compute & display the non-E1 angle-differential photoionization cross sections for all PhotoIonization.Line's and at all angles theta as defined in the settings. The general formula by Nishita Hosea (2025) is applied here. A neat table is printed for each line but nothing is returned otherwise.

source
JenaAtomicCalculator.PhotoIonization.computeLinesMethod

PhotoIonization.computeLines(finalMultiplet::Multiplet, initialMultiplet::Multiplet, nm::Nuclear.Model, grid::Radial.Grid, settings::PhotoIonization.Settings; output::Bool=true) ... to compute the photoIonization transition amplitudes and all properties as requested by the given settings. A list of lines::Array{PhotoIonization.Lines} is returned.

source
JenaAtomicCalculator.PhotoIonization.computeLinesCascadeMethod

PhotoIonization.computeLinesCascade(finalMultiplet::Multiplet, initialMultiplet::Multiplet, nm::Nuclear.Model, grid::Radial.Grid, settings::PhotoIonization.Settings, initialLevelSelection::LevelSelection; output=true, printout::Bool=true) ... to compute the photoionization transition amplitudes and all properties as requested by the given settings. The computations and printout is adapted for large cascade computations by including only lines with at least one channel and by sending all printout to a summary file only. A list of lines::Array{PhotoIonization.Lines} is returned.

source
JenaAtomicCalculator.PhotoIonization.computeLinesPlasmaMethod

PhotoIonization.computeLinesPlasma(finalMultiplet::Multiplet, initialMultiplet::Multiplet, nm::Nuclear.Model, grid::Radial.Grid, settings::PhotoIonization.PlasmaSettings; output::Bool=true) ... to compute the photoIonization transition amplitudes and all properties as requested by the given settings. A list of lines::Array{PhotoIonization.Lines} is returned.

source
JenaAtomicCalculator.PhotoIonization.computeStatisticalTensorUnpolarizedMethod

PhotoIonization.computeStatisticalTensorUnpolarized(k::Int64, q::Int64, gauge::EmGauge, line::PhotoIonization.Line, settings::PhotoIonization.Settings) ... to compute the statistical tensor of the photoion in its final level after the photoionization of initially unpolarized atoms by plane-wave photons with given Stokes parameters (density matrix). A value::ComplexF64 is returned.

source
JenaAtomicCalculator.PhotoIonization.computeTimeDelaysMethod

PhotoIonization.computeTimeDelays(channels::Array{PhotoIonization.Channel,1}, xchannels::Array{PhotoIonization.Channel,1}, deltaE::Float64, Jf::AngularJ64) ... to compute the – coherent and incoherent – time delay from the channels as calculated for two neighboured photon energies (deltaE = xE - E). Two tuple of two time delays (coherentDelay::EmProperty, incoherentDelay::EmProperty) is returned.

source
JenaAtomicCalculator.PhotoIonization.determineChannelsMethod

PhotoIonization.determineChannels(finalLevel::Level, initialLevel::Level, settings::PhotoIonization.Settings) ... to determine a list of photoionization Channel for a transitions from the initial to final level and by taking into account the particular settings of for this computation; an Array{PhotoIonization.Channel,1} is returned.

source
JenaAtomicCalculator.PhotoIonization.determineLinesMethod

PhotoIonization.determineLines(finalMultiplet::Multiplet, initialMultiplet::Multiplet, settings::PhotoIonization.Settings) ... to determine a list of PhotoIonization.Line's for transitions between levels from the initial- and final-state multiplets, and by taking into account the particular selections and settings for this computation; an Array{PhotoIonization.Line,1} is returned. Apart from the level specification, all physical properties are set to zero during the initialization process.

source
JenaAtomicCalculator.PhotoIonization.displayLineDataMethod

PhotoIonization.displayLineData(stream::IO, lines::Array{PhotoIonization.Line,1}) ... to display the calculated data, ordered by the initial levels and the photon energies involved. Neat tables of all initial levels and photon energies as well as all associated cross sections are printed but nothing is returned otherwise.

source
JenaAtomicCalculator.PhotoIonization.displayLinesMethod

PhotoIonization.displayLines(stream::IO, lines::Array{PhotoIonization.Line,1}) ... to display a list of lines and channels that have been selected due to the prior settings. A neat table of all selected transitions and energies is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.PhotoIonization.displayPhasesMethod

PhotoIonization.displayPhases(lines::Array{PhotoIonization.Line,1}) ... to display a list of lines, channels and phases of the continuum wave that have been selected due to the prior settings. A neat table of all selected transitions and energies is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.PhotoIonization.displayResultsMethod

PhotoIonization.displayResults(stream::IO, lines::Array{PhotoIonization.Line,1}, settings::PhotoIonization.Settings) ... to list all results, energies, cross sections, etc. of the selected lines. A neat table is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.PhotoIonization.extractCrossSectionMethod

PhotoIonization.extractCrossSection(lines::Array{PhotoIonization.Line,1}, omega::Float64, initialLevel) ... to extract from lines the total PI cross section that refer to the given omega and initial level; a cross section cs::EmProperty is returned.

source
JenaAtomicCalculator.PhotoIonization.extractCrossSectionMethod

PhotoIonization.extractCrossSection(lines::Array{PhotoIonization.Line,1}, omega::Float64, shell::Shell, initialLevel) ... to extract from lines the total PI cross section that refer to the given omega and initial level and to to the ionization of an electron from shell; a cross section cs::EmProperty is returned.

source
JenaAtomicCalculator.PhotoIonization.extractLinesMethod

PhotoIonization.extractLines(lines::Array{PhotoIonization.Line,1}, omega::Float64) ... to extract from lines all those that refer to the given omega; a reduced list rLines::Array{PhotoIonization.Line,1} is returned.

source
JenaAtomicCalculator.PhotoIonization.getLineKappasMethod

PhotoIonization.getLineKappas(line::PhotoIonization.Line) ... returns a list of kappa-values (partial waves) which contribute to the given line, to which one or several channels are assigned. An kappaList::Array{Int64,1} is returned.

source
JenaAtomicCalculator.PhotoIonization.interpolateCrossSectionMethod

PhotoIonization.interpolateCrossSection(lines::Array{PhotoIonization.Line,1}, omega::Float64, initialLevel) ... to interpolate (or extrapolate) from lines the total PI cross section for any given omega and initial level. The procedure applies a linear interpolation/extrapolation by just using the cross sections from the two nearest (given) omega points; a cross section cs::EmProperty is returned.

source

Photorecombination

JenaAtomicCalculator.PhotoRecombination.ChannelType

struct PhotoRecombination.Channel ... defines a type for a photorecombination channel to help characterize a single multipole and scattering (continuum) state of many electron-states with a single free electron.

+ multipole      ::EmMultipole          ... Multipole of the photon emission/absorption.
+ gauge          ::EmGauge              ... Gauge for dealing with the (coupled) radiation field.
+ kappa          ::Int64                ... partial-wave of the free electron
+ symmetry       ::LevelSymmetry        ... total angular momentum and parity of the scattering state
+ phase          ::Float64              ... phase of the partial wave
+ amplitude      ::Complex{Float64}     ... Rec amplitude associated with the given channel.
source
JenaAtomicCalculator.PhotoRecombination.LineType

struct PhotoRecombination.Line ... defines a type for a Photorecombination line that may include the definition of channels.

+ initialLevel   ::Level                  ... initial-(state) level
+ finalLevel     ::Level                  ... final-(state) level
+ electronEnergy ::Float64                ... Energy of the (incoming free) electron.
+ photonEnergy   ::Float64                ... Energy of the emitted photon.
+ betaGamma2     ::Float64                ... beta^2 * gamma^2.
+ weight         ::Float64                ... weight of line in the integration over electron energies.
+ crossSection   ::EmProperty             ... Cross section for this electron capture.
+ channels       ::Array{PhotoRecombination.Channel,1}    ... List of photorecombination channels of this line.
source
JenaAtomicCalculator.PhotoRecombination.SettingsType

struct PhotoRecombination.Settings <: AbstractProcessSettings ... defines a type for the details and parameters of computing photo recombination lines.

+ multipoles          ::Array{EmMultipole}  ... Multipoles of the radiation field that are to be included.
+ gauges              ::Array{UseGauge}     ... Gauges to be included into the computations.
+ electronEnergies    ::Array{Float64,1}    ... List of electron energies [in default units].
+ ionEnergies         ::Array{Float64,1}    ... List of ion energies [in MeV/u].
+ useIonEnergies      ::Bool                ... Make use of ion energies in [MeV/u] to obtain the electron energies.
+ calcTotalCs         ::Bool                ... True, if the total cross sections is to be calculated/displayed for all initial levels.
+ calcAnisotropy      ::Bool                ... True, if the overall anisotropy is to be calculated.
+ calcTensors         ::Bool                ... True, if the statistical tensors are to be calculated and 
                                                false otherwise.
+ printBefore         ::Bool                ... True, if all energies and lines are printed before their evaluation.
+ maxKappa            ::Int64               ... Maximum kappa value of partial waves to be included.
+ lineSelection       ::LineSelection       ... Specifies the selected levels, if any.
source
JenaAtomicCalculator.PhotoRecombination.SettingsMethod

PhotoRecombination.Settings(set::PhotoRecombination..Settings;

    multipoles=..,          gauges=..,              electronEnergies=..,          ionEnergies=..,     
    useIonEnergies=..,      calcTotalCs..,          calcAnisotropy=..,            calcTensors=..,             
    printBefore=..,         maxKappa=..,            lineSelection=..)
                
... constructor for modifying the given PhotoRecombination..Settings by 'overwriting' the previously selected parameters.
source
JenaAtomicCalculator.PhotoRecombination.amplitudeMethod

PhotoRecombination.amplitude(kind::String, channel::PhotoRecombination.Channel, energy::Float64, finalLevel::Level, continuumLevel::Level, grid::Radial.Grid) ... to compute the kind = (photorecombination) amplitude < alphaf Jf || O^(photorecombination) || (alphai Ji, epsilon kappa) J_t> due to the electron-photon interaction for the given final and continuum level, the partial wave of the outgoing electron as well as the given multipole and gauge. A value::ComplexF64 is returned.

source
JenaAtomicCalculator.PhotoRecombination.checkConsistentMultipletsMethod

PhotoRecombination.checkConsistentMultiplets(finalMultiplet::Multiplet, initialMultiplet::Multiplet) ... to check that the given initial- and final-state levels and multiplets are consistent to each other and to avoid later problems with the computations. An error message is issued if an inconsistency occurs, and nothing is returned otherwise.

source
JenaAtomicCalculator.PhotoRecombination.compareCrossSectionEmpiricalMethod

PhotoRecombination.compareCrossSectionEmpirical(energies::Array{Float64,1}, Z::Float64) ... to evaluate and compare different (non-relativistic) shell-resolved and total cross section for the RR of a free electron with energy into initially bare ions with nuclear charge Z; an cs::Float64 is returned.

source
JenaAtomicCalculator.PhotoRecombination.comparePlasmaRateEmpiricalMethod

PhotoRecombination.comparePlasmaRateEmpirical(temps::Array{Float64,1}, Z::Float64) ... to evaluate and compare different plasma rate coefficients for the capture of a Maxwellian-distributed free electron with temperature Te in temps into initially bare ions with nuclear charge Z; nothing is returned.

source
JenaAtomicCalculator.PhotoRecombination.computeAmplitudesPropertiesMethod

PhotoRecombination.computeAmplitudesProperties(line::PhotoRecombination.Line, nm::Nuclear.Model, grid::Radial.Grid, nrContinuum::Int64, settings::PhotoRecombination.Settings) ... to compute all amplitudes and properties of the given line; a line::PhotoRecombination.Line is returned for which the amplitudes and properties are now evaluated.

source
JenaAtomicCalculator.PhotoRecombination.computeCrossSectionBareIonMethod

PhotoRecombination.computeCrossSectionBareIon(energy_eV::Float64, subshell::Subshell, multipoles::Array{EmMultipole,1}, gauge::EmGauge, nm::Nuclear.Model, grid::Radial.Grid, nrContinuum::Int64) ... to compute the (hydrogenic) RR cross section for the capture of an electron into a single subshell; only the amplitudes for the given multipoles are taken into account; an cs::Float64 [a.u.] is returned.

source
JenaAtomicCalculator.PhotoRecombination.computeLinesMethod

PhotoRecombination.computeLines(finalMultiplet::Multiplet, initialMultiplet::Multiplet, nm::Nuclear.Model, grid::Radial.Grid, settings::PhotoRecombination.Settings; output::Bool=true) ... to compute the photo recombination transition amplitudes and all properties as requested by the given settings. A list of lines::Array{PhotoRecombination.Lines} is returned.

source
JenaAtomicCalculator.PhotoRecombination.computeLinesWithContinuumOrbitalMethod

PhotoRecombination.computeLinesWithContinuumOrbital(finalMultiplet::Multiplet, initialMultiplet::Multiplet, nm::Nuclear.Model, grid::Radial.Grid, cOrbitals::Dict{Subshell, Orbital}, energyGrid::Radial.GridGL, settings::PhotoRecombination.Settings; output::Bool=true) ... to compute the photo recombination transition amplitudes and all properties as requested by the given settings but with the given continuum orbitals cOrbitals and for the (free-) electronEnergies. A error message is issued if these energies are not the same a given by the settings. The continuum orbital with electronEnergies[i] has the principal quantum number 100+i. A list of lines::Array{PhotoRecombination.Lines} is returned.

source
JenaAtomicCalculator.PhotoRecombination.crossSectionBellTotalMethod

PhotoRecombination.crossSectionBellTotal(energy::Float64, Z::Float64) ... to evaluate the (non-relativistic) total Bell & Bell cross section for the RR of a free electron with energy into any shell of initially bare ions with nuclear charge Z; an cs::Float64 is returned.

source
JenaAtomicCalculator.PhotoRecombination.crossSectionKramersMethod

PhotoRecombination.crossSectionKramers(energy::Float64, Z::Float64, nLowUp::Tuple{Int64,Int64}) ... to evaluate the (non-relativistic) Kramers cross section for the RR of a free electron with energy into all shells with n = nLow ... nUp of initially bare ions with nuclear charge Z; an cs::Float64 is returned.

source
JenaAtomicCalculator.PhotoRecombination.crossSectionStobbeMethod

PhotoRecombination.crossSectionStobbe(energy::Float64, Z::Float64) ... to evaluate the (non-relativistic) Stobbe cross section for the RR of a free electron with energy into the 1s state of initially bare ions with nuclear charge Z; an cs::Float64 is returned.

source
JenaAtomicCalculator.PhotoRecombination.determineChannelsMethod

PhotoRecombination.determineChannels(finalLevel::Level, initialLevel::Level, settings::PhotoRecombination.Settings) ... to determine a list of RecChannel for a transitions from the initial to final level and by taking into account the particular settings of for this computation; an Array{PhotoRecombination.Channel,1} is returned.

source
JenaAtomicCalculator.PhotoRecombination.determineLinesMethod

PhotoRecombination.determineLines(finalMultiplet::Multiplet, initialMultiplet::Multiplet, settings::PhotoRecombination.Settings) ... to determine a list of PhotoRecombination.Line's for transitions between levels from the initial- and final-state multiplets, and by taking into account the particular selections and settings for this computation; an Array{PhotoRecombination.Line,1} is returned. Apart from the level specification, all physical properties are set to zero during the initialization process.

source
JenaAtomicCalculator.PhotoRecombination.displayLinesMethod

PhotoRecombination.displayLines(stream::IO, lines::Array{PhotoRecombination.Line,1}) ... to display a list of lines and channels that have been selected due to the prior settings. A neat table of all selected transitions and energies is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.PhotoRecombination.displayResultsMethod

PhotoRecombination.displayResults(stream::IO, lines::Array{PhotoRecombination.Line,1}, settings::PhotoRecombination.Settings) ... to list all results, energies, cross sections, etc. of the selected lines. A neat table is printed but nothing is returned otherwise.

source
JenaAtomicCalculator.PhotoRecombination.plasmaRateKotelnikovMethod

PhotoRecombination.plasmaRateKotelnikov(Te::Float64, Z::Float64) ... to evaluate the RR plasma rate coefficient at temperature Te and for a Maxwellian distribution of electron energies and an initially bare ions with nuclear charge Z; an alpha::Float64 is returned.

source
JenaAtomicCalculator.PhotoRecombination.plasmaRateKotelnikov_1sMethod

PhotoRecombination.plasmaRateKotelnikov_1s(Te::Float64, Z::Float64) ... to evaluate the RR plasma rate coefficient at temperature Te and for a Maxwellian distribution of electron energies and an initially bare ions with nuclear charge Z; here only the capture into 1s is taken into account. An alpha::Float64 is returned.

source
JenaAtomicCalculator.PhotoRecombination.plasmaRatePartialSeatonMethod

PhotoRecombination.plasmaRatePartialSeaton(Te::Float64, Z::Float64, n::Int64) ... to evaluate the RR plasma rate coefficient at temperature Te and for a Maxwellian distribution of electron energies and an initially bare ions with nuclear charge Z. An alpha::Float64 is returned.

source
JenaAtomicCalculator.PhotoRecombination.plasmaRateSeatonMethod

PhotoRecombination.plasmaRateSeaton(Te::Float64, Z::Float64) ... to evaluate the RR plasma rate coefficient at temperature Te and for a Maxwellian distribution of electron energies and an initially bare ions with nuclear charge Z; this formula is valid for 2*Te / Z^2 << 1. An alpha::Float64 is returned.

source