Manual: The Targets Module

The Targets module implements the abstraction of targets and provides parameters of some commonly-used targets.

eTraj.TargetsModule
module Targets

The Targets module contains abstraction of the targets and provides some pre-defined atoms or molecules.

source

Atoms under SAE approximation

The HydrogenLikeAtom and SAEAtom are both subtypes of the SAEAtomBase base type, which represents an atom under the SAE approximation. A key ingredient of atomic objects lies in the potential function of the residual ion after the electron gets ionized, which is the only difference between the two types.

The HydrogenLikeAtom's potential function is of the form:

\[\begin{equation} V(r) = -\frac{Z}{\sqrt{r^2+a}}, \end{equation}\]

with $a$ the soft-core parameter which avoids singularity of the potential during numerical simulation.

The SAEAtom's potential function is adopted from Tong's model [Tong_2005]:

\[\begin{equation} V(r) = -\frac{Z+a_1\ee^{-b_1 r}+a_2 r\ee^{-b_2 r}+a_3\ee^{-b_3 r}}{\sqrt{r^2+a}}, \end{equation}\]

where $a_i$ and $b_i$ are tunable parameters to fit the effective potential felt by the electron.

eTraj.Targets.HydrogenLikeAtomMethod
HydrogenLikeAtom(Ip, Z [,l=0] [,m=0] [,asymp_coeff=:hartree|<coeff>] [,quan_ax_θ=0.0] [,quan_ax_ϕ=0.0] [,soft_core=1e-10] [,name]) <: SAEAtomBase

Initializes a new HydrogenLikeAtom.

Parameters

  • Ip : Ionization potential of the atom (numerically in a.u. or a Unitful.Quantity).
  • Z : Nuclear charge number.
  • l=0 : Angular quantum number (optional, default 0).
  • m=0 : Magnetic quantum number (optional, default 0).
  • asymp_coeff=:hartree : Asymptotic coefficient related to the wavefunction's behavior when r→∞ (:hartree or a positive number). Passing :hartree (by default) indicates automatic calculation using the Hartree formula.
  • quan_ax_θ=0.0 : Orientation angle θ of the quantization axis relative to the lab frame (numerically in radian or a Unitful.Quantity) (optional, default 0.0).
  • quan_ax_ϕ=0.0 : Orientation angle ϕ of the quantization axis relative to the lab frame (numerically in radian or a Unitful.Quantity) (optional, default 0.0).
  • soft_core=1e-10 : Soft core parameter of the Coulomb potential (optional, default 1e-10).
  • name::String : Name of the atom.

Examples

julia> t = HydrogenLikeAtom(Ip=0.5, Z=1, name="H")
[HydrogenLikeAtom] Atom H, Ip=0.5000 (13.61 eV), Z=1

julia> using eTraj.Units

julia> t = HydrogenLikeAtom(Ip=3.4eV, Z=1, l=1, name="H")
[HydrogenLikeAtom] Atom H (p orbital, m=0), Ip=0.1249 (3.40 eV), Z=1

See Also

The get_atom method provides some atom presets for use.

source
eTraj.Targets.SAEAtomType
struct SAEAtom <: SAEAtomBase <: Targets

Represents an atom under single-active-electron (SAE) approximation.

source
eTraj.Targets.SAEAtomMethod
SAEAtom(Ip, Z [,l=0] [,m=0] [,asymp_coeff=:hartree|<coeff>] [,quan_ax_θ=0.0] [,quan_ax_ϕ=0.0] [,a1,b1,a2,b2,a3,b3] [,soft_core=1e-10] [,name]) <: SAEAtomBase

Initializes a new SAEAtom.

Parameters

  • Ip : Ionization potential of the atom (numerically in a.u. or a Unitful.Quantity).
  • Z : Nuclear charge number.
  • l=0 : Angular quantum number (optional, default 0).
  • m=0 : Magnetic quantum number (optional, default 0).
  • asymp_coeff=:hartree: Asymptotic coefficient related to the wavefunction's behavior when r→∞ (:hartree or a positive number). Passing :hartree (by default) indicates automatic calculation using the Hartree formula.
  • quan_ax_θ=0.0 : Orientation angle θ of the quantization axis relative to the lab frame (numerically in radian or a Unitful.Quantity) (optional, default 0.0).
  • quan_ax_ϕ=0.0 : Orientation angle ϕ of the quantization axis relative to the lab frame (numerically in radian or a Unitful.Quantity) (optional, default 0.0).
  • a1,b1,a2,b2,a3,b3 : Parameters used to fit the atomic potential. See [J. Phys. B 38, 2593 (2005)]
  • soft_core=1e-10 : Soft core parameter to avoid singularity in the potential (optional, default 1e-10).
  • name::String : Name of the atom.

Examples

julia> t = SAEAtom(Ip=0.9035, Z=1, asymp_coeff=:hartree, a1=1.230723, b1=0.6620055, a2=-1.325040, b2=1.236224, a3=-0.2307230, b3=0.4804286, name="He")
[SAEAtom] Atom He, Ip=0.9035 (24.59 eV), Z=1

julia> using eTraj.Units

julia> t = SAEAtom(Ip=12.13eV, Z=1, l=1, a1=51.35554, b1=2.111554, a2=-99.92747, b2=3.737221, a3=1.644457, b3=0.4306465, asymp_coeff=1.3, name="Xe")
[SAEAtom] Atom Xe (p orbital, m=0), Ip=0.4458 (12.13 eV), Z=1

See Also

The get_atom method provides some atom presets for use.

source

For the convenience of the user, there are some presets of commonly-used atoms, which can be accessed through the get_atom method. Available keys are accessed by invoking get_available_atoms.

  • HydrogenLikeAtom : H, He⁺, Li²⁺
  • SAEAtom : He, Ne, Ne⁺, Ne²⁺, Ar, Ar⁺, Ar²⁺, V, Ni, Kr, Kr⁺, Rb, Nb, Pd, Xe, Xe⁺, Ta
eTraj.Targets.get_atomFunction
get_atom(name::String, charge::Integer=0; kwargs...)

Constructs a HydrogenLikeAtom or SAEAtom in the database based on the name and charge given, other initialization parameters can be passed via kwargs. For HydrogenLikeAtom, the init params include Ip,Z and name; for SAEAtom, the init params include Ip,Z,a1,b1,a2,b2,a3,b3,l and name. The available keys can be found by invoking get_available_atoms.

Examples

julia> get_atom("H")
[HydrogenLikeAtom] Atom H, Ip=0.5000 (13.61 eV), Z=1

julia> get_atom("He")
[SAEAtom] Atom He, Ip=0.9036 (24.59 eV), Z=1

julia> get_atom("Ar"; m=1, asymp_coeff=0.950) # asymp_coeff from [Phys.-Uspekhi 47, 855 (2004)] (l=0 case)
[SAEAtom] Atom Ar (p orbital, m=1), Ip=0.5792 (15.76 eV), Z=1

julia> get_atom("Xe", 1)
[SAEAtom] Atom Xe⁺ (p orbital, m=0), Ip=0.7708 (20.97 eV), Z=2

julia> get_atom("Xe", 2)
ERROR: [get_atom] Atom name key `Xe2p` not found in database.
[...]

See also get_available_atoms, HydrogenLikeAtom and SAEAtom.

source

Molecules

For molecule targets, the GenericMolecule type is implemented based on the MoleculeBase supertype, whose structure is more complicated. A GenericMolecule stores information about the atoms that make up the molecule, together with their coordinates, as well as the asymptotic coefficients [$C_{lm}$ in Eq. (38) in Molecular SFA] and WFAT's integral coefficients [Eq. (56) in WFAT], which are obtained using other quantum chemistry packages.

There are two ways to initialize a GenericMolecule: build from zero or from an existing file.

eTraj.Targets.GenericMoleculeMethod
GenericMolecule(atoms, atom_coords [,charge=0] [,spin=0] [,name] [,rot_α=0.0] [,rot_β=0.0] [,rot_γ=0.0])

Initializes a new GenericMolecule with given parameters.

Parameters

  • atoms : Atoms in the molecule, stored as a Vector of String.
  • atom_coords : Atoms' coordinates in the molecule (numerically in Å or a Unitful.Quantity), stored as a N×3 Matrix.
  • charge : Total charge of the molecule (ion) (optional, default 0).
  • spin : Total spin of the molecule (optional, default 0). Note that each unpaired electron contributes 1/2.
  • name : Name of the molecule (optional).
  • rot_α,rot_β,rot_γ : Euler angles (z-y'-z'' convention) specifying the molecule's orientation (numerically in radian or a Unitful.Quantity) (optional, default 0).

Example

julia> m = GenericMolecule(atoms=["H","H"], atom_coords=[0.0 0.0 -0.375; 0.0 0.0 0.375], name="Hydrogen")
[GenericMolecule] Hydrogen

julia> using eTraj.Units

julia> m = GenericMolecule(atoms=["H","H"], atom_coords=[0.0 0.0 -0.375; 0.0 0.0 0.375]*Å, name="Hydrogen", rot_β=90°)
[GenericMolecule] Hydrogen, αβγ=(0.0°,90.0°,0.0°)
source
eTraj.Targets.LoadMoleculeFunction
LoadMolecule(ext_data_path; [rot_α=0.0] [,rot_β=0.0] [,rot_γ=0.0])

Initializes a new GenericMolecule with the data stored in ext_data_path.

Parameters

  • ext_data_path : Path to the molecule's data stored externally.
  • rot_α,rot_β,rot_γ : Euler angles (z-y'-z'' convention) specifying the molecule's orientation (numerically in radian or a Unitful.Quantity) (optional, default 0).
source

We also provide some pre-defined molecules, which can be accessed through the get_mol method.

eTraj.Targets.get_molFunction
get_mol(name::String; [rot_α=0.0] [,rot_β=0.0] [,rot_γ=0.0])

Constructs a GenericMolecule in the database based on the name given. Euler angles which describe the molecule's rotation can be applied via the keyword arguments rot_α, rot_β and rot_γ. The available molecule names can be found by invoking get_available_mols.

Examples


julia> get_mol("Hydrogen")
[GenericMolecule] Hydrogen (H₂)
Asymp coeff of HOMO available
WFAT data of HOMO available
#          E (Ha)  occp
⋮  ⋮         ⋮      ⋮⋮
3  LUMO+1   0.232  ----
2  LUMO     0.148  ----
1  HOMO    -0.594  -↿⇂-

julia> get_mol("Oxygen"; rot_β=π)
[GenericMolecule] Oxygen (O₂), αβγ=(0.0°,180.0°,0.0°)
Asymp coeff of α-HOMO-1 & α-HOMO available
WFAT data of α-HOMO-1 & α-HOMO available
#          Eα(Ha)  occp  Eβ(Ha)
⋮    ⋮        ⋮     ⋮⋮      ⋮     ⋮
11 LUMO+1   0.681  ----   0.741 LUMO+3
10 LUMO     0.377  ----   0.431 LUMO+2
9  HOMO    -0.554  -↿--   0.110 LUMO+1
8  HOMO-1  -0.554  -↿--   0.110 LUMO
7  HOMO-2  -0.763  -↿⇂-  -0.575 HOMO
6  HOMO-3  -0.841  -↿⇂-  -0.575 HOMO-1
5  HOMO-4  -0.841  -↿⇂-  -0.702 HOMO-2
4  HOMO-5  -1.204  -↿⇂-  -0.993 HOMO-3
⋮    ⋮        ⋮     ⋮⋮      ⋮     ⋮

See also get_available_mols and GenericMolecule.

source

Molecule's Orientation

The molecule's orientation is described by a set of Euler angles ($z-y'-z''$ convention), which defines a rotational transformation from the molecular frame (MF) to the lab frame (LF). This property is NOT included in the saved file and thus needs to be specified each time upon initialization of the GenericMolecule object from external files.

Note

Here the three Euler angles (α,β,γ) that describe the GenericMolecule's orientation are completely different from that of the Euler angles (θ,χ) in the WFAT and MO-SFA theory. These theories' "lab frame" is chosen for convenience of theoretical formulation, where the electric field is assumed to be static, pointing towards the $+z$ direction, and has no relation with the lab frame mentioned above.

The orientation of the molecule can be obtained and set via the MolRotation and SetMolRotation! methods.

Quantum Chemistry Calculation

As for the quantum chemistry calculation, we implemented the scheme in PySCFMolecularCalculator using the PySCF, which works on Linux and macOS platforms (as a remedy, Windows users can use the Windows Subsystem of Linux (WSL)). The calculation scheme of WFAT structure factor is adopted from PyStructureFactor. Future extension is possible by implementing the supertype MolecularCalculatorBase. Since there are some presets of molecules available via get_mol, we are not going to detail on the manual of running the calculation in the text.

The example of initializing and calculating the essential data of GenericMolecule in REPL is presented as follows:

julia> using eTraj.Targets, eTraj.Units

julia> mol = GenericMolecule(atoms=["O","C","O"], atom_coords=[0 0 -1.1600; 0 0 0; 0 0 1.1600]*Å, charge=0, name="Carbon Dioxide (CO₂)")
[GenericMolecule] Carbon Dioxide (CO₂)

julia> MolInitCalculator!(mol, basis="cc-pVTZ")
[ Info: [PySCFMolecularCalculator] Running molecular calculation...

julia> MolCalcAsympCoeff!(mol, 0); MolCalcAsympCoeff!(mol, -1)
[ Info: [PySCFMolecularCalculator] Running calculation of asymptotic coefficients... (ionizing orbital HOMO)
[ Info: [PySCFMolecularCalculator] Running calculation of asymptotic coefficients... (ionizing orbital HOMO-1)

julia> MolCalcWFATData!(mol, 0); MolCalcWFATData!(mol, -1)
[ Info: [PySCFMolecularCalculator] Running calculation of WFAT structure factor data... (ionizing orbital HOMO)
[ Info: [PySCFMolecularCalculator] Running calculation of WFAT structure factor data... (ionizing orbital HOMO-1)

julia> MolSaveDataAs!(mol, "Molecule_CO2.jld2")
[ Info: [GenericMolecule] Data saved for molecule Carbon Dioxide (CO₂) at `Molecule_CO2.jld2`.

julia> mol_ = LoadMolecule("Molecule_CO2.jld2")  # load from saved file
[GenericMolecule] Carbon Dioxide (CO₂)
Asymp coeff of HOMO-1 & HOMO available
WFAT data of HOMO-1 & HOMO available
#          E (Ha)  occp
⋮    ⋮       ⋮      ⋮⋮
13 LUMO+1   0.207  ----
12 LUMO     0.175  ----
11 HOMO    -0.542  -↿⇂-
10 HOMO-1  -0.542  -↿⇂-
9  HOMO-2  -0.714  -↿⇂-
8  HOMO-3  -0.714  -↿⇂-
⋮    ⋮       ⋮      ⋮⋮

List of Available Properties & Methods

Generic Target

eTraj.Targets.IonPotentialFunction
IonPotential(t::SAEAtomBase)

Gets the ionization potential of the atom.

IonPotential(mol::GenericMolecule)

Gets the ionization potential of the molecule's HOMO.

IonPotential(mol::GenericMolecule, orbit_ridx)

Gets the ionization potential of the specified MO of molecule.

  • orbit_ridx: Index of selected orbit relative to the HOMO (e.g., 0 indicates HOMO, and -1 indicates HOMO-1). For open-shell molecules, according to α/β spins, should be passed in format (spin, idx) where for α orbitals spin=1 and for β orbitals spin=2.
source
eTraj.Targets.TargetPotentialFunction
TargetPotential(t::SAEAtomBase) -> V(x,y,z)

Gets the potential function of the atom.

TargetPotential(mol::GenericMolecule) -> V(x,y,z)

Gets the asymptotic Coulomb potential function of the molecule.

source
eTraj.Targets.TargetForceFunction
TargetForce(t::SAEAtomBase) -> F(x,y,z) -> (Fx,Fy,Fz)

Gets the Coulomb force exerted on the electron from the atom.

TargetForce(mol::GenericMolecule) -> F(x,y,z) -> (Fx,Fy,Fz)

Gets the asymptotic Coulomb force exerted on the electron from the molecular ion.

source
eTraj.Targets.TrajectoryFunctionFunction
TrajectoryFunction(t::Target, dimension = 2|3, laserFx::Function, laserFy::Function, phase_method = :CTMC|:QTMC|:SCTS)

Gets the trajectory function of the electron according to given parameter.

source

SAEAtomBase


GenericMolecule

eTraj.Targets.MolAtomsFunction
MolAtoms(mol::GenericMolecule) -> ["atom1","atom2",...]

Gets the atoms in the molecule with their coordinates.

source
eTraj.Targets.MolSpinFunction
MolSpin(mol::GenericMolecule)

Gets the total spin of the molecule. Each unpaired electron contributes 1/2.

source
eTraj.Targets.MolEnergyLevelsFunction
MolEnergyLevels(mol::GenericMolecule [,spin=1|2])

Gets the energy levels of the molecule's MOs.

  • spin: For closed-shell molecules, the spin param should be neglected. For open-shell molecules (with non-zero spins), spin=1 indicates α orbitals and spin=2 indicates β orbitals, neglecting spin would return both two sets of orbitals.
source
eTraj.Targets.MolEnergyLevelFunction
MolEnergyLevel(mol::GenericMolecule, orbit_ridx)

Gets the energy level of the molecule's selected MO.

  • orbit_ridx: Index of selected orbit relative to the HOMO (e.g., 0 indicates HOMO, and -1 indicates HOMO-1). For open-shell molecules, according to α/β spins, should be passed in format (spin, idx) where for α orbitals spin=1 and for β orbitals spin=2.
source
eTraj.Targets.MolOrbitalOccupationFunction
MolOrbitalOccupation(mol::GenericMolecule [,spin=1|2])

Gets the occupation of the molecule's MOs.

  • spin: For closed-shell molecules, the spin param should be neglected. For open-shell molecules (with non-zero spins), spin=1 indicates α orbitals and spin=2 indicates β orbitals, neglecting spin would return both two sets of orbitals.
source
eTraj.Targets.MolHOMOEnergyFunction
MolHOMOEnergy(mol::GenericMolecule [,spin=1|2])

Gets the energy of the molecule's HOMO.

  • spin: For closed-shell molecules, the spin param should be neglected. For open-shell molecules (with non-zero spins), spin=1 indicates α orbitals and spin=2 indicates β orbitals, neglecting spin would give both.
source
eTraj.Targets.MolWFATDataFunction
MolWFATData(mol::GenericMolecule, orbit_ridx)

Gets the WFAT data in format (μ, int_data).

  • orbit_ridx: Index of selected orbit relative to the HOMO (e.g., 0 indicates HOMO, and -1 indicates HOMO-1). For open-shell molecules, according to α/β spins, should be passed in format (spin, idx) where for α orbitals spin=1 and for β orbitals spin=2.
source
eTraj.Targets.MolWFATStructureFactor_GFunction
MolWFATStructureFactor_G(mol::GenericMolecule, orbit_ridx, nξ, m, θ, χ)

Gets the WFAT structure factor $G_{n_ξ m}$ according to the given Euler angles θ and χ (z-y'-z'' convention). Note: the rotational Euler angles of the molecule would not be applied.

Parameters

  • orbit_ridx: Index of selected orbit relative to the HOMO (e.g., 0 indicates HOMO, and -1 indicates HOMO-1). For open-shell molecules, according to α/θ spins, should be passed in format (spin, idx) where for α orbitals spin=1 and for β orbitals spin=2.
  • : Parabolic quantum number nξ=0,1,2,⋯ (nξ up to 5 is calculated by default).
  • m : Parabolic quantum number m=⋯,-1,0,1,⋯ (|m| up to 5 is calculated by default).
  • θ : Euler angle θ, passed as a Real value or an AbstractVector of Real.
  • χ : Euler angle χ, passed as a Real value or an AbstractVector of Real.
source
eTraj.Targets.MolWFATMaxChannelsFunction
MolWFATMaxChannels(mol::GenericMolecule, orbit_ridx)

Gets the maximum value of nξ and |m| calculated in the WFAT integral data.

  • orbit_ridx: Index of selected orbit relative to the HOMO (e.g., 0 indicates HOMO, and -1 indicates HOMO-1). For open-shell molecules, according to α/β spins, should be passed in format (spin, idx) where for α orbitals spin=1 and for β orbitals spin=2.
source
eTraj.Targets.MolAsympCoeffFunction
MolAsympCoeff(mol::GenericMolecule, orbit_ridx)

Gets the asymptotic coefficients of the molecule.

  • orbit_ridx: Index of selected orbit relative to the HOMO (e.g., 0 indicates HOMO, and -1 indicates HOMO-1). For open-shell molecules, according to α/β spins, should be passed in format (spin, idx) where for α orbitals spin=1 and for β orbitals spin=2.
source
eTraj.Targets.MolAsympCoeff_lMaxFunction
MolAsympCoeff_lMax(mol::GenericMolecule, orbit_ridx)

Gets the maximum value of l calculated in the asymptotic coefficients.

  • orbit_ridx: Index of selected orbit relative to the HOMO (e.g., 0 indicates HOMO, and -1 indicates HOMO-1). For open-shell molecules, according to α/β spins, should be passed in format (spin, idx) where for α orbitals spin=1 and for β orbitals spin=2.
source
eTraj.Targets.MolRotationFunction
MolRotation(mol::GenericMolecule) -> (α,β,γ)

Gets the Euler angles (z-y'-z'' convention) specifying the molecule's orientation in format (α,β,γ) (in radian).

source
eTraj.Targets.SetMolRotation!Function
SetMolRotation!(mol::GenericMolecule, α,β,γ)
SetMolRotation!(mol::GenericMolecule, (α,β,γ))

Sets the Euler angles (z-y'-z'' convention) specifying the molecule's orientation in format (α,β,γ) (numerically in radian or a Unitful.Quantity).

source
eTraj.Targets.MolInitCalculator!Function
MolInitCalculator!(mol::GenericMolecule, MCType::Type=PySCFMolecularCalculator [;kwargs...])

Initializes the MolecularCalculator of mol with given parameters.

  • MCType : Type of MolecularCalculator if it is not initialized (default is PySCFMolecularCalculator).
  • kwargs... : Keyword arguments to pass to the initializer of MolecularCalculator of MCType, e.g., basis, ...
source
eTraj.Targets.MolCalcWFATData!Function
MolCalcWFATData!(mol::GenericMolecule [,orbit_ridx=0] [;kwargs...])

Calculates the WFAT data of the molecule.

  • orbit_ridx : Index of selected orbit relative to the HOMO (e.g., 0 indicates HOMO, and -1 indicates HOMO-1). For open-shell molecules, according to α/β spins, should be passed in format (spin, idx) where for α orbitals spin=1 and for β orbitals spin=2.
  • kwargs... : Keyword arguments to pass to the calc_WFAT_data method, e.g. grid_rNum, grid_rMax, sf_lMax, ⋯
source
eTraj.Targets.MolCalcAsympCoeff!Function
MolCalcAsympCoeff!(mol::GenericMolecule, orbit_ridx; kwargs...)

Calculates the asymptotic coefficients of the molecule.

  • orbit_ridx : Index of selected orbit relative to the HOMO (e.g., 0 indicates HOMO, and -1 indicates HOMO-1). For open-shell molecules, according to α/β spins, should be passed in format (spin, idx) where for α orbitals spin=1 and for β orbitals spin=2.
  • kwargs... : Keyword arguments to pass to the calc_asymp_coeff, e.g. grid_rNum, l_max.
source
eTraj.Targets.MolSaveDataAs!Function
MolSaveDataAs!(mol::GenericMolecule, data_path [,overwrite=false])

Saves the data of the GenericMolecule to the data_path. To overwrite the existing file, set overwrite=true.

source

MolecularCalculators

eTraj.Targets.PySCFMolecularCalculatorMethod
PySCFMolecularCalculator(; mol::MoleculeBase, basis::String="cc-pVDZ", kwargs...)

Initializes an instance of PySCFMolecularCalculator with given parameter.

Parameters

  • mol::MoleculeBase : The molecule to be calculated.
  • basis="cc-pVDZ" : Basis set used for calculation (default "cc-pVDZ").
source
eTraj.Targets.calc_WFAT_dataFunction

Calculates the data used in WFAT structure factor calculation of the given molecule.

Returns

(μ, int_data) : Orbital dipole momentum and the array which stores the integrals.

Parameters

  • mc : The molecular calculator.
  • orbit_ridx : Index of selected orbit relative to the HOMO (e.g., 0 indicates HOMO, and -1 indicates HOMO-1). For open-shell molecules, according to α/β spins, should be passed in format (spin, idx) where for α orbitals spin=1 and for β orbitals spin=2.
  • grid_rNum : The number of radial grid (default 200).
  • grid_rMax : The maximum radius of the radial grid (default 10.0).
  • grid_θNum : The number of angular grid in the θ direction (default 60).
  • grid_ϕNum : The number of angular grid in the ϕ direction (default 60).
  • sf_nξMax : The maximum number of nξ used in calculation (default 3).
  • sf_mMax : The maximum number of |m| used in calculation (default 3).
  • sf_lMax : The maximum angular quantum number l used in calculation (default 6).
  • swap_HOMO_LUMO=false : If the HOMO & LUMO are degenerate, set swap_HOMO_LUMO=true when calculating LUMO (ridx=1). The program would swap the coefficients of HOMO and LUMO.
source
eTraj.Targets.calc_asymp_coeffFunction

Calculates the asymptotic coefficients (used in ADK, SFA-SPANE, SFA-SPA) of the given molecule.

Parameters

  • mc : The molecular calculator.
  • orbit_ridx : Index of selected orbit relative to the HOMO (e.g., 0 indicates HOMO, and -1 indicates HOMO-1). For open-shell molecules, according to α/β spins, should be passed in format (spin, idx) where for α orbitals spin=1 and for β orbitals spin=2.
  • grid_rNum : The number of radial grid (default 200).
  • grid_rReg : The region of radial distance to fit the wavefunction to obtain the coefficients (default (3,8)).
  • grid_θNum : The number of angular grid in the θ direction (default 60).
  • grid_ϕNum : The number of angular grid in the ϕ direction (default 60).
  • l_max : The maximum number of l calculated (default 6).
source
  • Tong_2005X. M. Tong and C. D. Lin, Empirical formula for static field ionization rates of atoms and molecules by lasers in the barrier-suppression regime, J. Phys. B: At. Mol. Opt. Phys. 38, 2593 (2005). DOI: 10.1088/0953-4075/38/15/001