ampsci
High-precision calculations for one- and two-valence atomic systems
DiracOperator::Skfinal

Scalar multipole operator: \( S_K = t^K(q)\gamma^0 \). More...

#include <EM_multipole.hpp>

+ Inheritance diagram for DiracOperator::Sk:

Public Member Functions

 Sk (const Grid &gr, int K, double omega, const SphericalBessel::JL_table *jl=nullptr)
 
DiracSpinor radial_rhs (const int kappa_a, const DiracSpinor &Fb) const override final
 radial_int = Fa * radial_rhs(a, Fb) (a needed for angular factor)
 
double radialIntegral (const DiracSpinor &Fa, const DiracSpinor &Fb) const override final
 Radial part of integral R_ab = (Fa|t|Fb).
 
void updateFrequency (const double omega) override final
 nb: q = alpha*omega!
 
 Sk (const Sk &other)
 
Skoperator= (const Sk &other)
 
- Public Member Functions inherited from DiracOperator::EM_multipole
const SphericalBessel::JL_tablejl () const
 Returns the precomputed Bessel table pointer (may be nullptr).
 
std::string name () const override
 Returns a human-readable label, e.g. "T^E_1", "T^M5_2", "t_1", "P_1".
 
double angularF (const int ka, const int kb) const override
 Angular factor linking the radial integral to the reduced matrix element: \( \langle a \| h \| b \rangle = F(\kappa_a,\kappa_b) \cdot \int \! dr \).
 
void updateRank (int new_K) override
 Updates the tensor rank and adjusts parity accordingly.
 
std::unique_ptr< TensorOperatorclone () const override
 Creates a fully independent copy of this operator at its current (rank, frequency) state via the MultipoleOperator factory.
 
 EM_multipole (const EM_multipole &)=default
 
EM_multipoleoperator= (const EM_multipole &)=default
 
 EM_multipole (EM_multipole &&)=default
 
EM_multipoleoperator= (EM_multipole &&)=default
 
- Public Member Functions inherited from DiracOperator::TensorOperator
 TensorOperator (const TensorOperator &)=default
 
TensorOperatoroperator= (const TensorOperator &)=default
 
 TensorOperator (TensorOperator &&)=default
 
TensorOperatoroperator= (TensorOperator &&)=default
 
bool freqDependantQ () const
 
bool isZero (int ka, int kb) const
 If matrix element <a|h|b> is zero, returns true.
 
bool isZero (const DiracSpinor &Fa, const DiracSpinor &Fb) const
 
bool selectrion_rule (int twoJA, int piA, int twoJB, int piB) const
 
void scale (double lambda)
 Permanently re-scales the operator by constant, lambda.
 
const std::vector< double > & getv () const
 Returns a const ref to vector v.
 
double getc () const
 Returns a const ref to constant c.
 
int get_d_order () const
 
bool imaginaryQ () const
 returns true if operator is imaginary (has imag MEs)
 
int rank () const
 Rank k of operator.
 
int parity () const
 returns parity, as integer (+1 or -1)
 
int symm_sign (const DiracSpinor &Fa, const DiracSpinor &Fb) const
 returns relative sign between <a||x||b> and <b||x||a>
 
virtual std::string units () const
 Returns units of operator (usually au, may be MHz, etc.)
 
virtual double angularCff (int, int) const
 Angular factor for f_a*f_b part of radial integral.
 
virtual double angularCgg (int, int) const
 Angular factor for g_a*g_b part of radial integral.
 
virtual double angularCfg (int, int) const
 Angular factor for f_a*g_b part of radial integral.
 
virtual double angularCgf (int, int) const
 Angular factor for g_a*f_b part of radial integral.
 
double rme3js (int twoja, int twojb, int two_mb=1, int two_q=0) const
 ME = rme3js * RME.
 
double rme3js (const DiracSpinor &Fa, const DiracSpinor &Fb, int two_mb=1, int two_q=0) const
 ME = rme3js * RME.
 
DiracSpinor reduced_rhs (const int ka, const DiracSpinor &Fb) const
 <a||h||b> = Fa * reduced_rhs(a, Fb) (a needed for angular factor)
 
DiracSpinor reduced_lhs (const int ka, const DiracSpinor &Fb) const
 <b||h||a> = Fa * reduced_lhs(a, Fb) (a needed for angular factor)
 
double reducedME (const DiracSpinor &Fa, const DiracSpinor &Fb) const
 The reduced matrix element.
 
double fullME (const DiracSpinor &Fa, const DiracSpinor &Fb, std::optional< int > two_ma=std::nullopt, std::optional< int > two_mb=std::nullopt, std::optional< int > two_q=std::nullopt) const
 Returns "full" matrix element, for optional (ma, mb, q) [taken as int 2*]. If not specified, returns z-component (q=0), with ma=mb=min(ja,jb)
 
double matel_factor (MatrixElementType type, int twoJa, int twoJb) const
 Converts reduced matrix element to different "type" (MatrixElementType)
 
double matel_factor (MatrixElementType type, const DiracSpinor &Fa, const DiracSpinor &Fb) const
 

Additional Inherited Members

- Protected Member Functions inherited from DiracOperator::EM_multipole
 EM_multipole (int rank_k, Parity pi, double constant, const std::vector< double > &vec, Realness RorI, bool freq_dep, const Grid *grid, char type, char comp, bool low_q, const SphericalBessel::JL_table *jl=nullptr, char form='V')
 Initialise the EM_multipole base layer.
 
- Protected Member Functions inherited from DiracOperator::TensorOperator
 TensorOperator (int rank_k, Parity pi, double constant=1.0, const std::vector< double > &vec={}, int diff_order=0, Realness RorI=Realness::real, bool freq_dep=false)
 Constructs a tensor operator description.
 
- Protected Attributes inherited from DiracOperator::EM_multipole
const Gridm_grid {nullptr}
 
double m_omega
 Current frequency; cached by each derived updateFrequency().
 
char m_type {}
 
char m_comp {}
 
bool m_low_q {}
 
char m_form {}
 
const SphericalBessel::JL_tablem_jl {nullptr}
 
- Protected Attributes inherited from DiracOperator::TensorOperator
int m_rank
 
Parity m_parity
 
int m_diff_order
 
Realness opC
 
bool m_freqDependantQ {false}
 
double m_constant
 
std::vector< double > m_vec
 

Detailed Description

Scalar multipole operator: \( S_K = t^K(q)\gamma^0 \).

  • nb: q = alpha*omega, so "omega" = c*q
  • Implements the scalar multipole operator whose radial factor is e^{i q r} (represented via spherical Bessel functions j_L(q*r)). The operator includes the gamma^0 Dirac matrix structure.
  • The constructor takes an optional const SphericalBessel::JL_table *jl to enable lookup from a precomputed table for improved performance.
  • Note: The q value for jL(qr) will be the nearest to the requested q
    • you should ensure the lookup table is close enough, or this can lead to errors.

Member Function Documentation

◆ radial_rhs()

DiracSpinor DiracOperator::Sk::radial_rhs ( const int  kappa_a,
const DiracSpinor Fb 
) const
finaloverridevirtual

radial_int = Fa * radial_rhs(a, Fb) (a needed for angular factor)

By default, is defined as:

\[ \delta F_{\rm b} = C_{\rm overall}\,v(r) \begin{pmatrix} C_{ff}f_b(r) + C_{fg}g_b(r) \\ C_{gf}f_b(r) + C_{gg}g_b(r) \end{pmatrix} \]

\( C_{f/g} \) are angular coeficients - see TensorOperator::angularCff

See also: TensorOperator::radial_rhs

Reimplemented from DiracOperator::TensorOperator.

◆ radialIntegral()

double DiracOperator::Sk::radialIntegral ( const DiracSpinor Fa,
const DiracSpinor Fb 
) const
finaloverridevirtual

Radial part of integral R_ab = (Fa|t|Fb).

<a||t||b> = angularF(a,b) * radialIntegral(a,b)

The 'radial' integral \(R_{ab}\) for operator \(t\) is defined as

\[ \langle a ||t|| b \rangle \equiv R_{ab} * A_{ab} \]

where \(A_{ab}\) is TensorOperator::angularF()

By default, is implemented as:

\[ R = \int_0^\infty F_a\cdot\delta F_{\rm b} \,{\rm d}r \]

where \( \delta F_{\rm b} \) is defined in TensorOperator::radial_rhs

So, if TensorOperator::radial_rhs is defined normally, we have:

\[ C_{\rm overall}\int_0^\infty v(r)\left( C_{ff}f_a(r)f_b(r) + C_{fg}f_a(r)g_b(r) + C_{gf}g_a(r)f_b(r) + C_{gg}g_a(r)g_b(r) \right)\,{\rm d}r \]

\( C_{f/g} \) are angular coeficients - see TensorOperator::angularCff

Warning
  • This is defined for the "standard" case, which doesn't always suit
  • TensorOperator::radial_rhs may be overridden for special cases
  • If radial_rhs is overridden, then radialIntegral must also be_

Reimplemented from DiracOperator::TensorOperator.

◆ updateFrequency()

void DiracOperator::Sk::updateFrequency ( const double  omega)
finaloverridevirtual

nb: q = alpha*omega!

Reimplemented from DiracOperator::TensorOperator.


The documentation for this class was generated from the following files: