Getting started

... with Julia (in REPL)

Getting started with Julia (in REPL)

Info

Link to the Pluto jl and direct start with example...

Here, we shall not introduce Julia's syntax and concepts for which many tutorials are available on the web. Instead, we just wish to remind and highlight some simple (syntax) features that help to go easier around with JAC, and especially for occasional users from experiment or teaching. This reminder aims to lower the initial threshold for users that have been trained on other languages in the past. Here, we shall pick up some issues whose physics background is explained only later in other tutorial. Obviously, however, Julia is a very rich and powerful language with many features that go well beyond of what is (and will ever) needed for JAC.

In brief, JAC provides tools for performing atomic (structure) calculations of different kind and complexity, and for which further details are given in the tutorials below. To see anything from JAC, we shall first invoke the tools by:

using JenaAtomicCalculator

a line that will appear at the beginning of all subsequent tutorials. – A first powerful and frequently needed feature refers to Julia's help pages or just the "?". By typing, for instance, atom or computation

? atom
search: atomic Atomic AtomicState AtomicCompass AtomicStructure @atomic @atomicswap @atomicreplace

Couldn't find atom
Perhaps you meant atomic, atan, acot, acos, htol, hton, ltoh, ntoh, Atomic, @atomic, pathof, atand, atanh, cat, 
match, catch, stat, acotd, acoth, add, ans, abs, abs2, acosd, acosh, acsc, all, all!, any, any!, asec, asin, axes, 
as, tan or Beam
No documentation found.

Binding atom does not exist.

we see, that atom itself is not a well-defined term in the JAC toolbox but that there exists a number of related terms, such as Atomic, AtomicState (two modules of JAC) and others. We shall not enter here the modular structure of the JAC toolbox but start much simpler with:

? Orbital
JenaAtomicCalculator.Radial.OrbitalType

struct Radial.Orbital ... defines a type for a single-electron radial orbital function with a large and small component, and which can refer to either the standard or an explicitly given grid due to the logical flag useStandardGrid. Bound-state orbitals with energy < 0 are distinguished from free-electron orbitals by the flag isBound.

+ subshell        ::Subshell          ... Relativistic subshell.
+ isBound         ::Bool              ... Logical flag to distinguish between bound (true) and free-electron orbitals (false).
+ useStandardGrid ::Bool              ... Logical flag for using the standard grid (true) or an explicitly given grid (false).
+ energy          ::Float64           ... Single-electron energies of bound orbitals are always negative.
+ P               ::Array{Float64,1}  ... Large and ..
+ Q               ::Array{Float64,1}  ... small component of the radial orbital.
+ Pprime          ::Array{Float64,1}  ... dP/dr.
+ Qprime          ::Array{Float64,1}  ... dQ/dr.
+ grid            ::Array{Float64,1}  ... explic. defined radial grid array for P, Q, if StandardGrid = false.
source

which, apart from its formal meaning, is a particular data structure (struct) of JAC and which represents a relativistic orbital (function) including additional information that appears helpful in the given implementation. There are very many (say, more than 300) of such data struct's specified in the JAC toolbox, and thus quite obvious that nobody will remember the details of all these definitions. Indeed, the "?" is the right and a powerful means to remind yourself and make use of these data structures whenever necessary. Special care has been taken that all data structures and functions/methods comes with a reasonable explanation (docstring) in order to work efficiently with JAC.

For instance, we might ask of what can be added to each other in JAC:

? add
JenaAtomicCalculator.Basics.addFunction

Basics.add(ma::AngularM64, mb::AngularM64) ... adds the projections of the angular momenta ma + mb and returns a mc::AngularM64.

source

Basics.add(pota::Radial.Potential, potb::Radial.Potential) ... to add two radial potentials together that are defined on the same grid. A potential::RadialPotential is returned that inherits its radial size from the potential that is defined in a larger range of r-values.

source

Apart from a short explanation, these docstring always tell the user (i) in which module the method is defined; (ii) which arguments it takes, including Julia's multiple dispatch feature as well as (iii) the type of the return value. All these information are typically relevant to the user, especially if some input or output does not behave as it should. Indeed, the complexity can grow quite rapidly, for instance, if we ask for help of what we can generate:

? generate
JenaAtomicCalculator.Basics.generateFunction

Basics.generate(representation::AtomicState.Representation) ... to generate an atomic representation as specified by the representation.repType::AbstractRepresentationType. All relevant intermediate and final results are printed to screen (stdout). Nothing is returned.

Basics.perform(representation::AtomicState.Representation; output=true) ... to generate the same but to return the complete output in a dictionary; the particular output depends on the type and specifications of the representation but can easily accessed by the keys of this dictionary.

source

Basics.generate(repType::AtomicState.MeanFieldBasis, representation::AtomicState.Representation) ... to generate a mean-field basis (representation) for a set of reference configurations; all relevant intermediate and final results are printed to screen (stdout). Nothing is returned.

Basics.generate(repType::AtomicState.MeanFieldBasis, representation::AtomicState.Representation; output=true) ... to generate the same but to return the complete output in a dictionary; the particular output depends on the type and specifications of the representation but can easily accessed by the keys of this dictionary.

source

Basics.generate(repType::AtomicState.MeanFieldMultiplet, representation::AtomicState.Representation) ... to generate a mean-field basis (representation) for a set of reference configurations; all relevant intermediate and final results are printed to screen (stdout). Nothing is returned.

Basics.generate(repType::AtomicState.MeanFieldMultiplet, representation::AtomicState.Representation; output=true) ... to generate the same but to return the complete output in a dictionary; the particular output depends on the type and specifications of the representation but can easily accessed by the keys of this dictionary.

source

Basics.generate(repType::AtomicState.OneElectronSpectrum, representation::AtomicState.Representation) ... to generate a one-electron spectrum for the atomic potential from the (given) levels, based on a set of reference configurations as well as for given settings. Relevant intermediate and final results are printed to screen (stdout). Nothing is returned in this case.

+ `(repType::AtomicState.OneElectronSpectrum, representation::AtomicState.Representation; output=true)`  
... to generate the same but to return the complete output in a orbitals::Dict{Subshell, Orbital}.
source

Basics.generate(repType::AtomicState.CiExpansion, representation::AtomicState.Representation) ... to generate a configuration-interaction expansion for a single level symmetry and based on a set of reference configurations and a number of pre-specified steps. All relevant intermediate and final results are printed to screen (stdout). Nothing is returned.

Basics.generate(repType::AtomicState.CiExpansion, representation::AtomicState.Representation; output=true) ... to generate the same but to return the complete output in a dictionary; the particular output depends on the type and specifications of the computations but can easily accessed by the keys of this dictionary.

source

Basics.generate(repType::AtomicState.RasExpansion, representation::AtomicState.Representation) ... to generate a restricted active-space expansion for a single level symmetry and based on a set of reference configurations and a number of pre-specified steps. All relevant intermediate and final results are printed to screen (stdout). Nothing is returned.

Basics.generate(repType::AtomicState.RasExpansion, representation::AtomicState.Representation; output=true) ... to generate the same but to return the complete output in a dictionary; the particular output depends on the type and specifications of the computations but can easily accessed by the keys of this dictionary.

source

Basics.generate(repType::AtomicState.GreenExpansion, representation::AtomicState.Representation) ... to generate a Green (function) expansion for a given approach and excitation scheme of the electron, based on a set of reference configurations, a list of level symmetries as well as for given settings. All relevant intermediate and final results are printed to screen (stdout). Nothing is returned.

Basics.generate(repType::AtomicState.GreenExpansion, representation::AtomicState.Representation; output=true) ... to generate the same but to return the complete output in a dictionary; the particular output depends on the type and specifications of the representation but can easily accessed by the keys of this dictionary.

source

Basics.generate("condensed multiplet: by single weight", multiplet::Multiplet) ... to condense/reduce the number of CSF in the basis of the given multiplet due to a single 'weight'; a multiplet::Multiplet is returned. Not yet implemented !

source

Basics.generate("configuration list: NR, from basis", basis::Basis) ... to (re-) generate the list of NR configurations from the given basis; a confList::Array{Configuration,1} is returned.

source

Basics.generate("configuration list: NR, single-configuration", refConf::Configuration, NoExcitations::Int64, fromShells::Array{Shell,1}, toShells::Array{Shell,1}) ... to generate a non-relativistic configuration list, including the given reference configuration (refConf) and with all configurations that differ by NoExcitations from the fromShells into the toShells; an Array{Configuration,1} is returned.

source

Basics.generate("shells: ordered list for NR configurations", confs::Array{Configuration,1}) ... to generate for confs, i.e. all the given (non-relativistic) configurations, a common and ordered shell list; a list::Array{Shell,1} is returned.

source

Basics.generate("subshells: ordered list for two bases", basisA::Basis, basisB::Basis) ... to generate common and ordered subshell list for the two basis A and B; a list::Array{Subshell,1} is returned.

source

Basics.generate("single-electron spectrum: STO", N::Int64, potential::Radial.Potential, grid::Radial.Grid; N_0::Int64=30, alpha_0::Float64=1.0, beta_0::Float64=1.1) ... to generate a complete one-electron spectrum with N positive and N negative states, and by using even-tempered Slater-type orbitals (STO) with parameters $lpha_i = lpha_0 eta_0^i$; a spectrum::SingleElecSpectrum is returned where just N0 positive and N_0 negative are kept for later use. Not yet implemented !

Basics.generate("single-electron spectrum: STO, positive", N::Int64, potential::Radial.Potential, grid::Radial.Grid; N_0::Int64=30, alpha_0::Float64=1.0, beta_0::Float64=1.1) ... to generate the same but to return only the N_0 positive states. Not yet implemented !

source

Well, this is quite a lot, and we shall explain some of these methods below; a similar or even larger output, you can generate by ? perform as well as few other terms that are central to the implementation of JAC.

Constructors & program control:$\quad$ Another frequent use of the (help) "?" concerns the data flow and control of almost all computations. In JAC, we often make use of (so-called) Settings that enable the user to overwrite default values or to control the computation to the extent, he or she wishes to have control. These Settings are context dependent and are different for each atomic property or process that can be computed by the JAC toolbox. They are defined in the various modules and need to be specified accordingly. For instance, to control the computation of transition probabilities for the (fine-structure) levels between given initial- and final-state configuration, one has to overwrite the (defaults) settings:

? PhotoEmission.Settings
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

We shall meet these and (many) other settings quite often in the tutorials below. –- Beside of Julia's help features (?), however, it is sometimes difficult to remember the right term or function name. In this case, it easy to make a <double-tab> after the dot (notation) or to make use of the (Unix/Linux) grep command within the JenaAtomicCalculator/src directly. Similar line-search commands will exist also at other platforms. In particular, for those of you who wishes to support and extend the JAC toolbox, the dot expansion and the grep command will be found very helpful, perhaps more than other search tools.

Info

Users can also refer to the API Reference section that provides a comprehensive list of JAC declared Types and Functions for selected atomic processes.

Use of constructors:$\quad$ Another Julia feature, that is frequently applied in JAC, is the successive definition of constructors in order to set-up complex data structures. This features is applied, for instance, in order to define an Atomic.Computation or a Cascade.Computation as a whole. We shall explain these rather complex data types below in different tutorials. The same issue appears however already at a much simpler level. For example, if we wish to select (specify) a number of levels from a multiplet prior to some particular – configuration interaction – computation, we can make use of a

? LevelSelection
JenaAtomicCalculator.Basics.LevelSelectionType

struct Basics.LevelSelection < Basics.AbstractSelection ... defines a struct to specify a list of levels by means of their (level) indices or level symmetries.

+ active       ::Bool                     ... true, if some selection has been made.
+ indices      ::Array{Int64,1}           ... List of selected indices.
+ symmetries   ::Array{LevelSymmetry,1}   ... List of selected symmetries
source

Apart from the logical flag active, such a level selection requires to either specify a list of level numbers (indices) or level symmetries

LevelSelection(true, indices= [i for i in 1:11])
LevelSelection:  indices = [1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11];    symmetries = LevelSymmetry[].
LevelSelection(true, symmetries= [LevelSymmetry(1//2, Basics.plus)])
LevelSelection:  indices = Int64[];    symmetries = LevelSymmetry[1/2 +].

Here, we made use of a LevelSymmetry to specify the overall rotational $J^P$ symmetry of atomic levels.

? LevelSymmetry
JenaAtomicCalculator.Basics.LevelSymmetryType

struct Basics.LevelSymmetry < AbstractAngularMomentum ... defines a struct for defining the overall J^P symmetry of a level.

+ J          ::AngularJ64  ... total angular momentum of a level
+ parity     ::Parity      ... total parity of the level
source

As seen from this definition, the level symmetry just comprises the total angular momentum (of type AngularJ64) and the parity of the level (of type Parity). Therefore, the specification of a list of level symmetries in LevelSelection already requires to nest four constructors in order make the level selection explicit: (i) For the angular momentum, (ii) the parity, (iii) the level symmetry and (iv) to create a list (array) of such level symmetries. All the constructors can be specified and built together also in subsequent steps, such as:

J1    = AngularJ64(1//2);           J2 = AngularJ64(5//2)
pl    = Basics.plus;                mn = Basics.minus
lsym1 = LevelSymmetry(J1, pl);      lsym2 = LevelSymmetry(J2, mn)
levelsyms = [lsym1, lsym2]
2-element Vector{LevelSymmetry}:
 1/2 +
 5/2 -

or simply by nesting all the information within a single step

levelsyms = [LevelSymmetry(AngularJ64(1//2), Basics.plus), LevelSymmetry(AngularJ64(5//2), Basics.minus)]
2-element Vector{LevelSymmetry}:
 1/2 +
 5/2 -

Both way have their pros and cons, and often some mixture is applied where complex constructors are first assigned to some variables, and which are later utilized to built up constructors of higher complexity. –- To finally specify aǹ instance of a LevelSelection, we use (onc more) its second constructor above:

LevelSelection(true, indices=[1,2,3], symmetries=levelsyms)
LevelSelection:  indices = [1, 2, 3];    symmetries = LevelSymmetry[1/2 +, 5/2 -].

and which will tell the JAC program to compute the lowest three levels (1, 2, 3) as well as all levels with 1/2+ and 5/2- symmetry. Apart from the selection of individual levels, it is often helpful for the computation of atomic processes to make a prior LineSelection and in some cases even a PathwaySelection as, for instance, for dielectronic recombination processes. Some of these features will be explained below in subsequent tutorials of JAC:

Functions & methods:$\quad$ Like most other languages, Julia is based on the successive work through functions and methods; a function is first of all specified by its name and it maps a tuple of argument values upon a return value. For instance, the function

function addSomething(a, b)
    c = a + b
end
addSomething(3.1, 5//2)
5.6

can be used to add two numbers, rather independent of their particular type, and which are infered here automatically. However, additional type declarations might help to specialize a function and to ensure type stability:

function addSomething(a::Int64, b::Int64, c::Int64)
    d = a + b + c
end
addSomething (generic function with 2 methods)

While the function name is the same in both of these examples above, Julia carefully distinguishes between these two methods of the function addSomething that may differ by the type and/or the number of arguments. This multiple use (dispatch) of function name enables the user to write highly specialized code. Although a proper (and specialized) definition of functions is often very important for the performance of the program, we shall not discuss such technical issues here. Let us just mention, that a function/method may also return nothing:

typeof(nothing)
Nothing

In JAC, the value nothing is usually returned by all display functions that print some selected data or tabulation to screen or elsewhere but does not return a value otherwise.

Code failures:$\quad$ Beside of its large flexibility and user-friendliness, JAC might terminate from time to time for non-obvious reasons. Since JAC is first of all a physics code, no attempt has been made that all possible errors are fully captured and recovered by the program. Wrong input parameters or an inappropriate use of contructors will often lead to errors that cannot be resolved by the program. While some of this input can be readily recognized as wrong, and then lead to a proper error message, other wrong data may appear dynamically and cannot be captured with a reasonable overhead of the code. In JAC, therefore, several conditional if ... elseif ... else ... end blocks include an additional clause error("stop a") or similar; these are clauses, which due to a first design of a function should never be entered, but this appears not to be true in all cases. The use of these (fully) non-instructive error message have still a great advantage due to Julia: If not switched-off explicitly, Julia always reports for all program failures the hierarchy of call's, that are made before the error occurs, and lists these calls together with the file and line number of source code. For this reason, an error("stop a") readily shows the position where something unexpected occurs. A short inspection of the corresponding (line of the) source code often help to understand of what went wrong internally.

Julia macros:$\quad$ What can one do, if the (source) code itself does not tell so much about the problem ? –- In this case, it is often useful to include some additional printouts near to the line in question into the code and to re-run it again. There are different ways (@show, print(), println()) to place printout in the code; cf. https://julialang.org/learning/ A particular quick and useful way makes use of the Julia macro:

@show levelsyms
2-element Vector{LevelSymmetry}:
 1/2 +
 5/2 -

which simply repeats the names of the variables together with their values. Of course, the values of several such variables can be shown within the same call:

wa = 5;   wb = [2.0, pi];   wc = ones(3)
@show wa, wb, wc
(5, [2.0, 3.141592653589793], [1.0, 1.0, 1.0])

Indeed, this @show macro makes printout very easy. There are many macros (all starting with @) in Julia which need not to be considered here. We just mention that @time in front of a Julia command (block) will take and display the CPU time that is necessary to run this line(s):

@time rand(50000)   ;
50000-element Vector{Float64}:
 0.7322902666689579
 0.665010952196498
 0.17016319549006076
 0.1282412830728945
 0.374715612495395
 0.625735078543251
 0.8837019201191297
 0.7320370587800145
 0.6387568001243382
 0.19593637125649577
 ⋮
 0.8890106906820144
 0.2565097041012029
 0.1149205542635019
 0.5410890914990949
 0.8406938214185473
 0.5465005720732459
 0.5453958996342793
 0.674136907182943
 0.21976626807970978

... with JAC (in REPL)

Getting started with JAC (in REPL)

Info

JAC user guide pdf .... link

using JenaAtomicCalculator

Welcome to JAC, the Jena Atomic Calculator

... that provides various tools for performing atomic (structure) calculations of different kinds and complexities. Apart from the computation of atomic (many-electron) amplitudes, properties and processes, JAC supports interactive, restricted-active space (RAS) and cascade computations. It also help perform a few simple hydrogenic and semi-empirical estimates as well as simplify symbolic expressions from Racah's algebra. –- Let's first use ? JenaAtomicCalculator in order to obtain more information about this toolbox:

? JenaAtomicCalculator
JenaAtomicCalculatorModule

module JenaAtomicCalculator ... Jena Atomic Calculator (JAC) provides tools for performing atomic (structure) calculations at various degrees of complexity and sophistication. It has been designed to not only calculate atomic level structures and properties [such as g-factors or hyperfine and isotope-shift parameters] but also transition amplitudes between bound-state levels [for the dipole operator, etc.] and, in particular, (atomic) transition probabilities, Auger rates, photoionization cross sections, radiative and dielectronic recombination rates as well as cross sections for several other – elementary or composed – processes.

source

H'm, this tells us a lot of details which we still need to better understand. To quickly list the atomic properties, that have been (partly) considered in JAC, we can use ? Details.properties or some other of the listed calls:

? Details.properties
Missing docstring.

Missing docstring for Details.properties. Check Documenter's build log for details.

In the design of JAC, we first of all aim for a precise language that (i) is simple enough for both, seldom and a more frequent use of this package, (ii) highlights the underlying physics and (iii) avoids most technical slang that is often unnecessary but quite common to many other codes. An intuitive picture about the level or hyperfine structure of an atom, its properties as well as possible excitation and/or decay processes should (always) come first in order to generate the desired data: By making use of suitable data types (struct), we indeed wish to introduce a language close to the underlying formalism. –- While JAC is overall based on a rather large number $(> 300)$ of such types, a few simple examples are:

  • (atomic) Shell: $\quad$1s, 2s, 2p, ...
  • Subshell: $\quad$1s1/2, 2s1/2, 2p1/2, 2p3/2, ...
  • (electron) Configuration: $\quad$1s^2 2s^2 2p^6 3s $\quad$ or $\quad$ [Ne] 3s, ...
  • Level: $\quad$1s^2 2s^2 ^1S_0, ...

and many other terms (types) that we shall explain later.

Let us simply start, for instance, with specifying and assigning the $1s$ and $2p$ shells:

w1s = Shell("1s")
w2p = Shell("2p")
2p

Similarly, we can readily specify and assign any (relativistic) subshell:

Subshell("2p_1/2"),   Subshell("2p_3/2")
(2p_1/2, 2p_3/2)

In JAC, we make use of these Shell's and Subshell's whenever they will naturally occur in describing the level structure or the excitation, decay or occupation of an atom, and this both at input and output. If you have forgotten how to specify such a subshell (constructor), simply ask:

? Subshell
JenaAtomicCalculator.Basics.SubshellType

struct Basics.Subshell ... defines a type for the allowed values of a relativistic subhell.

+ n        ::Int64  ... principal quantum number 
+ kappa    ::Int64  ... relativistic angular quantum number
source

Of course, we can interactively also specify any electron configuration:

? Configuration
wc1 = Configuration("1s^2 2s^2 2p^5")
wc2 = Configuration("[Ar] 4s^2 3d^5")
Configuration: 1s^2 2s^2 2p^6 3s^2 3p^6 3d^5 4s^2 
Info

For specific processes users can find the list of types and functions in the API Reference

This input just shows three (very) simple examples and how the details of some computation can be readily specified in line with our basic understanding of the atomic shell model. One can use ? Details.datatypes in order to see a more complete list of most data structures that are speficic to the JAC module ... and which will give you a very first impression about the size of the JAC program.

? Details.datatypes
Missing docstring.

Missing docstring for Details.datatypes. Check Documenter's build log for details.

This list gives further details why Julia (and JAC) is a very suitable and powerful framework for running – many-electron – atomic computations.

Of course, there are many other features that make Julia & JAC as powerful as it is: For example, the user may pre-define and overwrite the units in which he wishes to communicate with JAC. These units determine how (most of) the input data are interpreted as well as output data are displayed in tabulations or to screen. The current defaults settings for the units can be seen by typing:

Basics.display("settings")
Current settings of the JAC module:
-----------------------------------

  + Framework:                              relativistic
  + Energy unit:                            eV
  + Rate and transition probability unit:   1/s
  + Cross section unit:                     barn
  + Time unit:                              sec

  + A standard grid has been defined; cf. Defaults.getDefaults()

which show that energies are taken/printed in eV, rates in 1/s, etc. Apart from modifying these defaults directly in the source code, the can be overwritten by the user at any time of the program executation. This is done by means of the function

? Defaults.setDefaults
JenaAtomicCalculator.Defaults.setDefaultsFunction

Defaults.setDefaults() ... (re-) defines some 'standard' settings which are common to all the computations with the JAC module, and which can be 'overwritten' by the user. –- An improper setting of some variable may lead to an error message, if recognized immediately. The following defaults apply if not specified otherwise by the user: the framework is 'relativistic', energies are given in eV and cross sections in barn. Note that, internally, atomic units are used throughout for all the computations within the program. nothing is returned if not indicated otherwise.

+ ("framework: relativistic") or ("framework: non-relativistic")

... to define a relativistic or non-relativistic framework for all subsequent computations.

  • ("method: continuum, spherical Bessel") or ("method: continuum, pure sine") or ("method: continuum, asymptotic Coulomb") or ("method: continuum, nonrelativistic Coulomb") or ("method: continuum, Galerkin") ... to define a a method for the generation of the continuum orbitals as (pure) spherical Bessel, pure sine, asymptotic Coulomb, nonrelativistic Coulomb orbital or by means of the B-spline-Galerkin method.

  • ("method: normalization, pure sine") or ("method: normalization, pure Coulomb") or ("method: normalization, Ong-Russek") ... to define a method for the normalization of the continuum orbitals as asymptotically (pure) sine or Coulomb functions, or following the procedure by Ong & Russek (1978).

  • ("QED model: Petersburg") or ("QED model: Sydney") ... to define a model for the computation of the QED corrections following the work by Shabaev et al. (2011; Petersburg) or Flambaum and Ginges (2004; Syney).

  • ("unit: energy", "eV") or ("unit: energy", "Kayser") or ("unit: energy", "Hartree") or ("unit: energy", "Hz") or ("unit: energy", "Hz") ... to (pre-) define the energy units for all further printouts and communications with the JAC module.

  • ("unit: cross section", "a.u.") or ("unit: cross section", "barn") or ("unit: cross section", "Mbarn") ... to (pre-) define the unit for the printout of cross sections.

  • ("unit: rate", "a.u.") or ("unit: rate", "1/s") ... to (pre-) define the unit for the printout of rates.

  • ("unit: resonance strength", "a.u.") or ("unit: resonance strength", "barn eV") or ("unit: resonance strength", "cm^2 eV") ... to (pre-) define the unit for the printout of resonance strengths.

  • ("unit: time", "a.u.") or ("unit: time", "sec") or ("unit: time", "fs") or ("unit: time", "as") ... to (pre-) define the unit for the printout and communications of times with the JAC module.

source
  • ("relativistic subshell list", subshells::Array{Subshell,1}; printout::Bool=true) ... to (pre-) define internally the standard relativistic subshell list on which the standard order of orbitals is based.
source
  • ("standard grid", grid::Radial.Grid; printout::Bool=true) ... to (pre-) define internally the standard radial grid which is used to represent most orbitals.
source
  • ("continuum: potential", scField::Basics.AbstractScField) ... to (re-) define the potential that is applied for the generation of the continuum orbitals.
source
  • ("QED: damped-hydrogenic", Znuc::Float64, wa::Array{Float64,1}) ... to (re-) define the lambda-C damped overlap integrals of the lowest kappa-orbitals [ wa1s1/2, wa2p1/2, wa2p3/2, wa3d3/2, wa3d5/2 ] for the (new) nuclear charge Znuc; nothing is returned.
source

which enables one to re-define various global values of JAC. If we wish to enter/display energies in Kaysers or cross sections in atomic units, we can simply type:

Defaults.setDefaults("unit: energy", "Kayser")
Defaults.setDefaults("unit: cross section", "a.u.")

Here, again nothing is returned but the corresponding global constants are now changed.

Basics.display("settings")
Current settings of the JAC module:
-----------------------------------

  + Framework:                              relativistic
  + Energy unit:                            Kayser
  + Rate and transition probability unit:   1/s
  + Cross section unit:                     a.u.
  + Time unit:                              sec

  + A standard grid has been defined; cf. Defaults.getDefaults()

Apart from the default units, one can similarly overwrite the method that is use for the generation and normalization of continuum orbitals and several others. Although called global, the corresponding values can be accesses just in two ways. (i) The global constants, such as the electron mass, the speed of light, the fine-structure constant $\alpha$, etc., are accessed via the function:

? Defaults.getDefaults
JenaAtomicCalculator.Defaults.getDefaultsFunction

Defaults.getDefaults() ... gives/supplies different information about the (present) framework of the computation or about some given data; cf. Defaults.setDefaults().

  • ("alpha") or ("fine-structure constant alpha") ... to get the (current) value::Float64 of the fine-structure constant alpha.

  • ("electron mass: kg") or ("electron mass: amu") ... to get the (current) value::Float64 of the electron mass in the specified unit.

  • ("framework") ... to give the (current) setting::String of the overall framework.

  • ("electron rest energy") or ("mc^2") ... to get the electron rest energy.

  • ("electron g-factor") ... to give the electron g-factor g_s = 2.00232.

  • ("unit: energy") or ("unit: cross section") or ("unit: rate") or ("unit: strength") or ("unit: time") ... to get the corresponding (user-defined) unit::String for the current computations.

  • ("standard grid") ... to get the (current standard) grid::Array{Float64,1} to which all radial orbital functions usually refer.

  • ("speed of light: c") ... to get the speed of light in atomic units.

  • ("summary flag/stream") ... to get the logical flag and stream for printing a summary file; a tupel (flag, iostream) is returned.

source
  • ("ordered shell list: non-relativistic", n_max::Int64) ... to give an ordered list of non-relativistic shells::Array{Shell,1} up to the (maximum) principal number n_max.
  • ("ordered subshell list: relativistic", n_max::Int64) ... to give an ordered list of relativistic subshells::Array{Subshell,1} up to the (maximum) principal number n_max.
source
Defaults.getDefaults("alpha")
Defaults.getDefaults("electron rest energy")
Defaults.getDefaults("unit: energy")
"Kayser"

(ii) These global values are frequently applied in order to – internally or externally – convert physical numbers into units of the same dimension. This is done by the function:

? Defaults.convertUnits
JenaAtomicCalculator.Defaults.convertUnitsFunction

Defaults.convertUnits() ... converts some data from one format/unit into another one; cf. the supported keystrings and return values.

  • ("cross section: from atomic to predefined unit", value::Float64) or ("cross section: from atomic", value::Float64) ... to convert an cross section value from atomic to the predefined cross section unit; a Float64 is returned.

  • ("cross section: from barn to atomic unit", value::Float64) ... to convert an cross section value from barn atomic section unit; a Float64 is returned.

  • ("cross section: from atomic to barn", value::Float64) or ("cross section: from atomic to Mbarn", value::Float64) or ("cross section: from atomic to cm^2", value::Float64) ... to convert an energy value from atomic to the speficied cross section unit; a Float64 is returned.

  • ("cross section: from predefined to atomic unit", value::Float64) or ("cross section: to atomic", value::Float64) ... to convert a cross section value from the predefined to the atomic cross section unit; a Float64 is returned.

  • ("einstein B: from atomic", value::Float64) ... to convert a Einstein B coefficient from atomic to the speficied energy units; a Float64 is returned.

  • ("density: from [g/cm^3] to atomic", value::Float64) ... to convert a mass density from [g/cm^3] to atomic units [u/a_o^3]; a Float64 is returned.

  • ("energy-diff. cross section: from atomic to predefined unit", value::Float64) or ("energy-diff. cross section: from atomic", value::Float64) ... to convert an energy-diff. cross section value from atomic to the predefined energy-diff. cross section unit; a Float64 is returned.

  • ("energy: from atomic to eV", value::Float64) or ("energy: from atomic to Kayser", value::Float64) or ("energy: from atomic to Hz", value::Float64) or ("energy: from atomic to Angstrom", value::Float64) or ("energy: from atomic to Ws", value::Float64) ... to convert an energy value from atomic to the speficied energy unit; a Float64 is returned.

  • ("energy: from predefined to atomic unit", value::Float64) or ("energy: to atomic", value::Float64)... to convert an energy value from the predefined to the atomic energy unit; a Float64 is returned.

  • ("energy: from eV to atomic", value::Float64) ... to convert an energy value from eV to the atomic energy unit; a Float64 is returned.

  • ("energy: from wavelength [nm] to atomic", value::Float64) ... to convert a wavelength [nm] to the atomic energy unit; a Float64 is returned.

  • ("intensity: from W/cm^2 to atomic", value::Float64) ... to convert the intensity [in W/cm^2] to the atomic intensity unit; a Float64 is returned.

  • ("kinetic energy to wave number: atomic units", value::Float64) ... to convert a kinetic energy value (in a.u.) into a wave number k (a.u.); a Float64 is returned.

  • ("kinetic energy to wavelength: atomic units", value::Float64) ... to convert a kinetic energy value (in a.u.) into a wavelength (a.u.); a Float64 is returned.

  • ("length: from fm to atomic", value::Float64) ... to convert a length value (in fm) into a.u.; a Float64 is returned.

  • ("length: from atomic to fm", value::Float64) or ("energy: from atomic to Kayser", value::Float64) ... to convert a length value (in Bohr's a.u.) to the speficied length unit; a Float64 is returned.

  • ("rate: from atomic to predefined unit", value::Float64) or ("rate: from atomic", value::Float64) ... to convert a rate value from atomic to the predefined rate unit; a Float64 is returned.

  • ("rate: from atomic to 1/s", value::Float64) ... to convert an rate value from atomic to the speficied rate unit; a Float64 is returned.

  • ("rate: from predefined to atomic unit", value::Float64) or `("rate: to atomic", value::Float64)'... to convert a rate value from the predefined to the atomic rate unit; a Float64 is returned.

  • ("strength: from atomic to predefined unit", value::Float64) or ("strength: from atomic", value::Float64) ... to convert a (resonance) strength value from atomic to the predefined rate unit; a Float64 is returned.

  • ("time: from atomic to predefined unit", value::Float64) or ("time: from atomic", value::Float64) ... to convert an time value from atomic to the predefined time unit; a Float64 is returned.

  • ("time: from atomic to sec", value::Float64) or ("time: from atomic to fs", value::Float64)' or("time: from atomic to as", value::Float64)` ... to convert a time value from atomic to the speficied time unit; a Float64 is returned.

  • ("time: from predefined to atomic unit", value::Float64) or `("time: to atomic", value::Float64)' ... to convert a time value from the predefined to the atomic time unit; a Float64 is returned.

  • ("temperature: from Kelvin to (Hartree) units", value::Float64) ... to convert a temperature in Kelvin into atomic (Hartree) units; a Float64 is returned.

  • ("temperature: from atomic to Kelvin", value::Float64) ... to convert an atomic (energy) unit into Kelvin; 1 Hartree = 315774.64 K; a Float64 is returned.

  • ("wave number to total electron energy: atomic units", value::Float64) ... to convert a wavenumber (a.u.) into the total electron energy, including the rest energy; a Float64 is returned.

  • ("wave number to kinetic energy: atomic units", value::Float64) ... to convert a wavenumber (a.u.) into the kinetic energy; a Float64 is returned.

source
  • (string, values::Array{Float64,1}) ... to convert for the same strings as above but for an list of values; a corresponding Array{Float64,1} is returned.
source

This function is called at many places within JAC to generate tables where all physical data are printed out in the pre-specified units:

Defaults.convertUnits("energy: from atomic", 1.0)
219474.63068

With the given user-selection, this is equivalent to:

Defaults.convertUnits("energy: from atomic to Kayser", 1.0)
219474.63068

In JAC, the call of this function is often combined with some proper formatting of the results, such as:

using Printf
@sprintf("%.4e", Defaults.convertUnits("energy: from atomic", 1.0))
"2.1947e+05"