Hex  1.0
Hydrogen-electron collision solver
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
amplitudes.cpp File Reference
#include <algorithm>
#include <cmath>
#include <complex>
#include <cstring>
#include <cstdlib>
#include <regex>
#include <vector>
#include <fftw3.h>
#include "arrays.h"
#include "bspline.h"
#include "chebyshev.h"
#include "clenshawcurtis.h"
#include "radial.h"
#include "special.h"
#include "matrix.h"
Include dependency graph for amplitudes.cpp:

Functions

cArray computeLambda (Bspline const &bspline, rArray const &kf, rArray const &ki, int maxell, int L, int Spin, int Pi, int ni, int li, int mi, rArray const &Ei, int lf, cArray const &Pf_overlaps, std::vector< std::pair< int, int >> const &coupled_states)
 Extract radial part of scattering amplitude. More...
 
Chebyshev< double, Complexfcheb (Bspline const &bspline, cArrayView const &PsiSc, double kmax, int l1, int l2)
 
cArrays computeXi (Bspline const &bspline, int maxell, int L, int Spin, int Pi, int ni, int li, int mi, rArray const &Ei, rArray &ics, std::vector< std::pair< int, int >> const &coupled_states)
 Extract radial part of ionization amplitude. More...
 

Function Documentation

cArray computeLambda ( Bspline const &  bspline,
rArray const &  kf,
rArray const &  ki,
int  maxell,
int  L,
int  Spin,
int  Pi,
int  ni,
int  li,
int  mi,
rArray const &  Ei,
int  lf,
cArray const &  Pf_overlaps,
std::vector< std::pair< int, int >> const &  coupled_states 
)

Compute radial integrals for evaluation of the discrete T-matrices,

\[ \Lambda_l^{(1)LMS} = \int_0^{R_0} P_{n_f l_f}(r_1) \mathcal{W}\left[ \psi_{\ell_1 \ell_2}^{LMS}(r_1,\bullet), \hat{j}_l(k_f\bullet) \right]_{R_0} \mathrm{d}r_1 \ . \]

Parameters
bsplineThe B-spline environment to use when evaluating the solutions.
kfOutgoing projectile momenta.
kiincoming projectile momenta.
maxellMaximal angular momentum of the projectile.
LTotal angular momentum (partial wave).
SpinConserved total spin (partial wave).
PiTotal parity (partial wave).
niInitial atomic state - principal quantum number.
liInitial atomic state - orbital quantum number.
miInitial atomic state - magnetic quantum number.
EiInitial projectile energies.
lfFinal atomic orbital momentum.
Pf_overlapsHydrogenic function B-spline overlap integrals.
coupled_statesList of all coupled-state angular momenta pairs.
Returns
Vector of radial integrals.
cArrays computeXi ( Bspline const &  bspline,
int  maxell,
int  L,
int  Spin,
int  Pi,
int  ni,
int  li,
int  mi,
rArray const &  Ei,
rArray ics,
std::vector< std::pair< int, int >> const &  coupled_states 
)

Compute hyperangular integrals for evaluation of the ionization T-matrices,

\[ \Xi_{\ell_1 \ell_2}^{LS}(k_1,k_2) = \int_0^{\pi/2} \left( F_{\ell_1}(k_1,r_1) F_{\ell_2}(k_2,r_2) \frac{\partial}{\partial\rho} \psi^{LS}_{\ell_1\ell_2}(r_1,r_2) - \psi^{LS}_{\ell_1\ell_2}(r_1,r_2) \frac{\partial}{\partial\rho} F_{\ell_1}(k_1,r_1) F_{\ell_2}(k_2,r_2) \right) \rho\mathrm{d}\alpha \ , \]

where \( r_1 = \rho \cos\alpha \) and \( r_2 = \rho \sin\alpha \).

Parameters
bsplineThe B-spline environment to use when evaluating the solutions.
maxellMaximal angular momentum of the free electrons.
LTotal angular momentum (partial wave).
SpinConserved total spin (partial wave).
PiTotal parity (partial wave).
niInitial atomic state - principal quantum number.
liInitial atomic state - orbital quantum number.
miInitial atomic state - magnetic quantum number.
EiInitial projectile energies.
icsIonization cross section (on return).
coupled_statesList of all coupled-state angular momenta pairs.
Returns
Vector of radial integrals.
Chebyshev<double,Complex> fcheb ( Bspline const &  bspline,
cArrayView const &  PsiSc,
double  kmax,
int  l1,
int  l2 
)

DEBUG

DEBUG