ampsci
High-precision calculations for one- and two-valence atomic systems
EM_multipole.hpp
1#pragma once
2#include "DiracOperator/Operators/EM_multipole_base.hpp"
3#include "DiracOperator/TensorOperator.hpp"
4#include "IO/InputBlock.hpp"
5#include "Maths/SphericalBessel.hpp"
6#include "Wavefunction/Wavefunction.hpp"
7#include "qip/Maths.hpp"
8
9// include the 'low qr' form here, simply for convenience
10// (EM_multipole_lowqr.hpp already includes EM_multipole_base.hpp)
11#include "EM_multipole_lowqr.hpp"
12
13namespace DiracOperator {
14
15//==============================================================================
16//! @brief Vector electric multipole (transition) operator, Length-form: ( \f$ T^{(+1),{\rm Len}}_K \f$ ).
17/*!
18 @details
19 - nb: q = alpha*omega, so "omega" = c*q
20 - Electric multipole operator in the L-form (transition form
21 convention): t^E_L. The radial dependence is expressed in terms of
22 spherical Bessel functions j_L(q*r) with q = alpha * omega.
23 - The operator supports using either on-the-fly Bessel evaluation or a
24 precomputed `SphericalBessel::JL_table` supplied via the constructor
25 (pointer `jl`). When `jl` is non-null the table is used to look up the
26 closest j_L(q) vector for performance.
27 - Note: The q value for jL(qr) will be the _nearest_ to the requested q
28 - you should ensure the lookup table is close enough, or this can lead to errors.
29*/
30class VEk_Len final : public EM_multipole {
31public:
32 VEk_Len(const Grid &gr, int K, double omega,
33 const SphericalBessel::JL_table *jl = nullptr)
34 : EM_multipole(K, Angular::evenQ(K) ? Parity::even : Parity::odd, 1.0,
35 gr.r(), Realness::imaginary, true, &gr, 'V', 'E', false, jl,
36 'L') {
37 if (omega != 0.0)
38 updateFrequency(omega);
39 }
40 DiracSpinor radial_rhs(const int kappa_a,
41 const DiracSpinor &Fb) const override final;
42
43 double radialIntegral(const DiracSpinor &Fa,
44 const DiracSpinor &Fb) const override final;
45
46 //! nb: q = alpha*omega!
47 void updateFrequency(const double omega) override final;
48
49private:
50 std::vector<double> m_jK{};
51 std::vector<double> m_jKp1{};
52 const std::vector<double> *p_jK{nullptr};
53 const std::vector<double> *p_jKp1{nullptr};
54
55public:
56 VEk_Len(const VEk_Len &other)
57 : EM_multipole(other),
58 m_jK(other.m_jK),
59 m_jKp1(other.m_jKp1),
60 p_jK(other.p_jK == &other.m_jK ? &m_jK : other.p_jK),
61 p_jKp1(other.p_jKp1 == &other.m_jKp1 ? &m_jKp1 : other.p_jKp1) {}
62 VEk_Len &operator=(const VEk_Len &other) {
63 if (this != &other) {
64 EM_multipole::operator=(other);
65 m_jK = other.m_jK;
66 m_jKp1 = other.m_jKp1;
67 p_jK = other.p_jK == &other.m_jK ? &m_jK : other.p_jK;
68 p_jKp1 = other.p_jKp1 == &other.m_jKp1 ? &m_jKp1 : other.p_jKp1;
69 }
70 return *this;
71 }
72};
73
74//==============================================================================
75//! @brief Vector electric multipole (V-form) operator: \f$ V^E_K = T^{(+1)}_K(q) \f$
76/*!
77 @details
78 - nb: q = alpha*omega, so "omega" = c*q
79 - Vector (spatial) electric multipole operator used in the V-form
80 (transition convention). Radial dependence uses spherical
81 Bessel functions j_L(q*r) and derived combinations like j_L(q*r)/(q*r)
82 where needed; q = alpha * omega.
83 - The constructor takes an optional `const SphericalBessel::JL_table *jl` to
84 enable lookup from a precomputed table for improved performance.
85 - Note: The q value for jL(qr) will be the _nearest_ to the requested q
86 - you should ensure the lookup table is close enough, or this can lead to errors.
87*/
88class VEk final : public EM_multipole {
89public:
90 VEk(const Grid &gr, int K, double omega,
91 const SphericalBessel::JL_table *jl = nullptr)
92 : EM_multipole(K, Angular::evenQ(K) ? Parity::even : Parity::odd, 1.0,
93 gr.r(), Realness::imaginary, true, &gr, 'V', 'E', false,
94 jl) {
95 if (omega != 0.0)
96 updateFrequency(omega);
97 }
98 DiracSpinor radial_rhs(const int kappa_a,
99 const DiracSpinor &Fb) const override final;
100
101 double radialIntegral(const DiracSpinor &Fa,
102 const DiracSpinor &Fb) const override final;
103
104 //! nb: q = alpha*omega!
105 void updateFrequency(const double omega) override final;
106
107private:
108 std::vector<double> m_jK_on_qr{};
109 std::vector<double> m_jKp1{};
110 const std::vector<double> *p_jK_on_qr{nullptr};
111 const std::vector<double> *p_jKp1{nullptr};
112
113public:
114 VEk(const VEk &other)
115 : EM_multipole(other),
116 m_jK_on_qr(other.m_jK_on_qr),
117 m_jKp1(other.m_jKp1),
118 p_jK_on_qr(other.p_jK_on_qr == &other.m_jK_on_qr ? &m_jK_on_qr :
119 other.p_jK_on_qr),
120 p_jKp1(other.p_jKp1 == &other.m_jKp1 ? &m_jKp1 : other.p_jKp1) {}
121 VEk &operator=(const VEk &other) {
122 if (this != &other) {
123 EM_multipole::operator=(other);
124 m_jK_on_qr = other.m_jK_on_qr;
125 m_jKp1 = other.m_jKp1;
126 p_jK_on_qr =
127 other.p_jK_on_qr == &other.m_jK_on_qr ? &m_jK_on_qr : other.p_jK_on_qr;
128 p_jKp1 = other.p_jKp1 == &other.m_jKp1 ? &m_jKp1 : other.p_jKp1;
129 }
130 return *this;
131 }
132};
133
134//==============================================================================
135//! @brief Vector longitudinal multipole operator (V-form): \f$ V^L_K = T^{(-1)}_K(q) \f$
136/*!
137 @details
138 - nb: q = alpha*omega, so "omega" = c*q
139 - Implements the longitudinal (scalar-like) component of the vector
140 multipole operator. Radial dependence uses spherical Bessel functions
141 j_L(q*r) and, when required, the combination j_L(q*r)/(q*r).
142 - The constructor takes an optional `const SphericalBessel::JL_table *jl` to
143 enable lookup from a precomputed table for improved performance.
144 - Note: The q value for jL(qr) will be the _nearest_ to the requested q
145 - you should ensure the lookup table is close enough, or this can lead to errors.
146*/
147class VLk final : public EM_multipole {
148public:
149 VLk(const Grid &gr, int K, double omega,
150 const SphericalBessel::JL_table *jl = nullptr)
151 : EM_multipole(K, Angular::evenQ(K) ? Parity::even : Parity::odd, 1.0,
152 gr.r(), Realness::imaginary, true, &gr, 'V', 'L', false,
153 jl) {
154 if (omega != 0.0)
155 updateFrequency(omega);
156 }
157 DiracSpinor radial_rhs(const int kappa_a,
158 const DiracSpinor &Fb) const override final;
159
160 double radialIntegral(const DiracSpinor &Fa,
161 const DiracSpinor &Fb) const override final;
162
163 //! nb: q = alpha*omega!
164 void updateFrequency(const double omega) override final;
165
166private:
167 std::vector<double> m_jK_on_qr{};
168 std::vector<double> m_jKp1{};
169 const std::vector<double> *p_jK_on_qr{nullptr};
170 const std::vector<double> *p_jKp1{nullptr};
171
172public:
173 VLk(const VLk &other)
174 : EM_multipole(other),
175 m_jK_on_qr(other.m_jK_on_qr),
176 m_jKp1(other.m_jKp1),
177 p_jK_on_qr(other.p_jK_on_qr == &other.m_jK_on_qr ? &m_jK_on_qr :
178 other.p_jK_on_qr),
179 p_jKp1(other.p_jKp1 == &other.m_jKp1 ? &m_jKp1 : other.p_jKp1) {}
180 VLk &operator=(const VLk &other) {
181 if (this != &other) {
182 EM_multipole::operator=(other);
183 m_jK_on_qr = other.m_jK_on_qr;
184 m_jKp1 = other.m_jKp1;
185 p_jK_on_qr =
186 other.p_jK_on_qr == &other.m_jK_on_qr ? &m_jK_on_qr : other.p_jK_on_qr;
187 p_jKp1 = other.p_jKp1 == &other.m_jKp1 ? &m_jKp1 : other.p_jKp1;
188 }
189 return *this;
190 }
191};
192
193//==============================================================================
194//! @brief Vector magnetic multipole operator: \f$ V^M_K = T^{(0)}_K(q) \f$
195/*!
196 @details
197 - nb: q = alpha*omega, so "omega" = c*q
198 - Implements the magnetic multipole operator. The radial dependence is
199 expressed via spherical Bessel functions j_L(q*r) with q = alpha * omega.
200 - The constructor takes an optional `const SphericalBessel::JL_table *jl` to
201 enable lookup from a precomputed table for improved performance.
202 - Note: The q value for jL(qr) will be the _nearest_ to the requested q
203 - you should ensure the lookup table is close enough, or this can lead to errors.
204*/
205class VMk final : public EM_multipole {
206public:
207 VMk(const Grid &gr, int K, double omega,
208 const SphericalBessel::JL_table *jl = nullptr)
209 : EM_multipole(K, Angular::evenQ(K) ? Parity::odd : Parity::even, 1.0,
210 gr.r(), Realness::imaginary, true, &gr, 'V', 'M', false,
211 jl) {
212 if (omega != 0.0)
213 updateFrequency(omega);
214 }
215
216 DiracSpinor radial_rhs(const int kappa_a,
217 const DiracSpinor &Fb) const override final;
218
219 double radialIntegral(const DiracSpinor &Fa,
220 const DiracSpinor &Fb) const override final;
221
222 //! nb: q = alpha*omega!
223 void updateFrequency(const double omega) override final;
224
225private:
226 std::vector<double> m_jK{};
227 const std::vector<double> *p_jK{nullptr};
228
229public:
230 VMk(const VMk &other)
231 : EM_multipole(other),
232 m_jK(other.m_jK),
233 p_jK(other.p_jK == &other.m_jK ? &m_jK : other.p_jK) {}
234 VMk &operator=(const VMk &other) {
235 if (this != &other) {
236 EM_multipole::operator=(other);
237 m_jK = other.m_jK;
238 p_jK = other.p_jK == &other.m_jK ? &m_jK : other.p_jK;
239 }
240 return *this;
241 }
242};
243
244//==============================================================================
245//! @brief Temporal component of the vector multipole operator: \f$ \Phi_K = t^K(q) \f$
246/*!
247 @details
248 - nb: q = alpha*omega, so "omega" = c*q
249 - Implements the time-like (temporal) component of the vector multipole
250 operator with explicit frequency dependence via j_L(q*r).
251 - The constructor takes an optional `const SphericalBessel::JL_table *jl` to
252 enable lookup from a precomputed table for improved performance.
253 - Note: The q value for jL(qr) will be the _nearest_ to the requested q
254 - you should ensure the lookup table is close enough, or this can lead to errors.
255*/
256class Phik final : public EM_multipole {
257public:
258 Phik(const Grid &gr, int K, double omega,
259 const SphericalBessel::JL_table *jl = nullptr)
260 : EM_multipole(K, Angular::evenQ(K) ? Parity::even : Parity::odd, 1.0,
261 gr.r(), Realness::real, true, &gr, 'V', 'T', false, jl) {
262 if (omega != 0.0)
263 updateFrequency(omega);
264 }
265 DiracSpinor radial_rhs(const int kappa_a,
266 const DiracSpinor &Fb) const override final;
267
268 double radialIntegral(const DiracSpinor &Fa,
269 const DiracSpinor &Fb) const override final;
270
271 //! nb: q = alpha*omega!
272 void updateFrequency(const double omega) override final;
273
274private:
275 std::vector<double> m_jK{};
276 const std::vector<double> *p_jK{nullptr};
277
278public:
279 Phik(const Phik &other)
280 : EM_multipole(other),
281 m_jK(other.m_jK),
282 p_jK(other.p_jK == &other.m_jK ? &m_jK : other.p_jK) {}
283 Phik &operator=(const Phik &other) {
284 if (this != &other) {
285 EM_multipole::operator=(other);
286 m_jK = other.m_jK;
287 p_jK = other.p_jK == &other.m_jK ? &m_jK : other.p_jK;
288 }
289 return *this;
290 }
291};
292
293//==============================================================================
294//! @brief Scalar multipole operator: \f$ S_K = t^K(q)\gamma^0 \f$
295/*!
296 @details
297 - nb: q = alpha*omega, so "omega" = c*q
298 - Implements the scalar multipole operator whose radial factor is
299 e^{i q r} (represented via spherical Bessel functions j_L(q*r)). The
300 operator includes the gamma^0 Dirac matrix structure.
301 - The constructor takes an optional `const SphericalBessel::JL_table *jl` to
302 enable lookup from a precomputed table for improved performance.
303 - Note: The q value for jL(qr) will be the _nearest_ to the requested q
304 - you should ensure the lookup table is close enough, or this can lead to errors.
305*/
306class Sk final : public EM_multipole {
307public:
308 Sk(const Grid &gr, int K, double omega,
309 const SphericalBessel::JL_table *jl = nullptr)
310 : EM_multipole(K, Angular::evenQ(K) ? Parity::even : Parity::odd, 1.0,
311 gr.r(), Realness::real, true, &gr, 'S', 'T', false, jl) {
312 if (omega != 0.0)
313 updateFrequency(omega);
314 }
315 DiracSpinor radial_rhs(const int kappa_a,
316 const DiracSpinor &Fb) const override final;
317
318 double radialIntegral(const DiracSpinor &Fa,
319 const DiracSpinor &Fb) const override final;
320
321 //! nb: q = alpha*omega!
322 void updateFrequency(const double omega) override final;
323
324private:
325 std::vector<double> m_jK{};
326 const std::vector<double> *p_jK{nullptr};
327
328public:
329 Sk(const Sk &other)
330 : EM_multipole(other),
331 m_jK(other.m_jK),
332 p_jK(other.p_jK == &other.m_jK ? &m_jK : other.p_jK) {}
333 Sk &operator=(const Sk &other) {
334 if (this != &other) {
335 EM_multipole::operator=(other);
336 m_jK = other.m_jK;
337 p_jK = other.p_jK == &other.m_jK ? &m_jK : other.p_jK;
338 }
339 return *this;
340 }
341};
342
343//==============================================================================
344//==============================================================================
345// Gamma^5 versions!
346
347//==============================================================================
348//! @brief Axial electric multipole operator: \f$ A^E_K = T^{(+1)}_K(q)\gamma^5 \f$
349/*!
350 @details
351 - Gamma^5 variant of the electric multipole (vector) operator. Functions
352 analogously to `VEk` but with the gamma^5 Dirac structure applied to the
353 angular part.
354 - Radial functions use spherical Bessel combinations j_L(q*r) or
355 j_L(q*r)/(q*r) with q = alpha * omega.
356 - Accepts an optional `const SphericalBessel::JL_table *jl` to use precomputed
357 Bessel vectors; otherwise computes them on demand.
358*/
359class AEk final : public EM_multipole {
360public:
361 AEk(const Grid &gr, int K, double omega,
362 const SphericalBessel::JL_table *jl = nullptr)
363 : EM_multipole(K, Angular::evenQ(K) ? Parity::odd : Parity::even, 1.0,
364 gr.r(), Realness::real, false, &gr, 'A', 'E', false, jl) {
365 if (omega != 0.0)
366 updateFrequency(omega);
367 }
368 DiracSpinor radial_rhs(const int kappa_a,
369 const DiracSpinor &Fb) const override final;
370
371 double radialIntegral(const DiracSpinor &Fa,
372 const DiracSpinor &Fb) const override final;
373
374 //! nb: q = alpha*omega!
375 void updateFrequency(const double omega) override final;
376
377private:
378 std::vector<double> m_jK_on_qr{};
379 std::vector<double> m_jKp1{};
380 const std::vector<double> *p_jK_on_qr{nullptr};
381 const std::vector<double> *p_jKp1{nullptr};
382
383public:
384 AEk(const AEk &other)
385 : EM_multipole(other),
386 m_jK_on_qr(other.m_jK_on_qr),
387 m_jKp1(other.m_jKp1),
388 p_jK_on_qr(other.p_jK_on_qr == &other.m_jK_on_qr ? &m_jK_on_qr :
389 other.p_jK_on_qr),
390 p_jKp1(other.p_jKp1 == &other.m_jKp1 ? &m_jKp1 : other.p_jKp1) {}
391 AEk &operator=(const AEk &other) {
392 if (this != &other) {
393 EM_multipole::operator=(other);
394 m_jK_on_qr = other.m_jK_on_qr;
395 m_jKp1 = other.m_jKp1;
396 p_jK_on_qr =
397 other.p_jK_on_qr == &other.m_jK_on_qr ? &m_jK_on_qr : other.p_jK_on_qr;
398 p_jKp1 = other.p_jKp1 == &other.m_jKp1 ? &m_jKp1 : other.p_jKp1;
399 }
400 return *this;
401 }
402};
403
404//==============================================================================
405//! @brief Axial longitudinal multipole operator: \f$ A^L_K = T^{(-1)}_K(q)\gamma^5 \f$
406/*!
407 @details
408 - Gamma^5 variant of the longitudinal multipole operator. Works like
409 `VLk` but with the gamma^5 Dirac structure applied to the angular part.
410 - Supports optional `const SphericalBessel::JL_table *jl` for lookup-table
411 acceleration; otherwise uses on-the-fly Bessel evaluation.
412*/
413class ALk final : public EM_multipole {
414public:
415 ALk(const Grid &gr, int K, double omega,
416 const SphericalBessel::JL_table *jl = nullptr)
417 : EM_multipole(K, Angular::evenQ(K) ? Parity::odd : Parity::even, 1.0,
418 gr.r(), Realness::real, true, &gr, 'A', 'L', false, jl) {
419 if (omega != 0.0)
420 updateFrequency(omega);
421 }
422 DiracSpinor radial_rhs(const int kappa_a,
423 const DiracSpinor &Fb) const override final;
424
425 double radialIntegral(const DiracSpinor &Fa,
426 const DiracSpinor &Fb) const override final;
427
428 //! nb: q = alpha*omega!
429 void updateFrequency(const double omega) override final;
430
431private:
432 std::vector<double> m_jK_on_qr{};
433 std::vector<double> m_jKp1{};
434 const std::vector<double> *p_jK_on_qr{nullptr};
435 const std::vector<double> *p_jKp1{nullptr};
436
437public:
438 ALk(const ALk &other)
439 : EM_multipole(other),
440 m_jK_on_qr(other.m_jK_on_qr),
441 m_jKp1(other.m_jKp1),
442 p_jK_on_qr(other.p_jK_on_qr == &other.m_jK_on_qr ? &m_jK_on_qr :
443 other.p_jK_on_qr),
444 p_jKp1(other.p_jKp1 == &other.m_jKp1 ? &m_jKp1 : other.p_jKp1) {}
445 ALk &operator=(const ALk &other) {
446 if (this != &other) {
447 EM_multipole::operator=(other);
448 m_jK_on_qr = other.m_jK_on_qr;
449 m_jKp1 = other.m_jKp1;
450 p_jK_on_qr =
451 other.p_jK_on_qr == &other.m_jK_on_qr ? &m_jK_on_qr : other.p_jK_on_qr;
452 p_jKp1 = other.p_jKp1 == &other.m_jKp1 ? &m_jKp1 : other.p_jKp1;
453 }
454 return *this;
455 }
456};
457
458//==============================================================================
459//! @brief Axial magnetic multipole operator: \f$ A^M_K = T^{(0)}_K(q)\gamma^5 \f$
460/*!
461 @details
462 - Gamma^5 variant of the magnetic multipole operator. Analogous to
463 `VMk` but with the gamma^5 Dirac structure applied where appropriate.
464 - Uses spherical Bessel functions j_L(q*r) for radial dependence and
465 accepts an optional `const SphericalBessel::JL_table *jl`.
466*/
467class AMk final : public EM_multipole {
468public:
469 AMk(const Grid &gr, int K, double omega,
470 const SphericalBessel::JL_table *jl = nullptr)
471 : EM_multipole(K, Angular::evenQ(K) ? Parity::even : Parity::odd, 1.0,
472 gr.r(), Realness::real, true, &gr, 'A', 'M', false, jl) {
473 if (omega != 0.0)
474 updateFrequency(omega);
475 }
476
477 DiracSpinor radial_rhs(const int kappa_a,
478 const DiracSpinor &Fb) const override final;
479
480 double radialIntegral(const DiracSpinor &Fa,
481 const DiracSpinor &Fb) const override final;
482
483 //! nb: q = alpha*omega!
484 void updateFrequency(const double omega) override final;
485
486private:
487 std::vector<double> m_jK{};
488 const std::vector<double> *p_jK{nullptr};
489
490public:
491 AMk(const AMk &other)
492 : EM_multipole(other),
493 m_jK(other.m_jK),
494 p_jK(other.p_jK == &other.m_jK ? &m_jK : other.p_jK) {}
495 AMk &operator=(const AMk &other) {
496 if (this != &other) {
497 EM_multipole::operator=(other);
498 m_jK = other.m_jK;
499 p_jK = other.p_jK == &other.m_jK ? &m_jK : other.p_jK;
500 }
501 return *this;
502 }
503};
504
505//==============================================================================
506//! @brief Temporal component of the axial vector multipole operator
507/*!
508 @details
509 \f$ \Theta_K = \Phi^5_K = t^K(q)\gamma^5 \f$
510
511 - Gamma^5 variant of the temporal (time-like) component of the vector
512 multipole operator. Functions like `Phik` but with gamma^5 applied to
513 the appropriate spin-angular structure.
514 - Radial dependence uses j_L(q*r) and the constructor accepts an
515 optional `const SphericalBessel::JL_table *jl`.
516*/
517class Phi5k final : public EM_multipole {
518public:
519 Phi5k(const Grid &gr, int K, double omega,
520 const SphericalBessel::JL_table *jl = nullptr)
521 : EM_multipole(K, Angular::evenQ(K) ? Parity::odd : Parity::even, 1.0,
522 gr.r(), Realness::real, true, &gr, 'A', 'T', false, jl) {
523 if (omega != 0.0)
524 updateFrequency(omega);
525 }
526 DiracSpinor radial_rhs(const int kappa_a,
527 const DiracSpinor &Fb) const override final;
528
529 double radialIntegral(const DiracSpinor &Fa,
530 const DiracSpinor &Fb) const override final;
531
532 //! nb: q = alpha*omega!
533 void updateFrequency(const double omega) override final;
534
535private:
536 std::vector<double> m_jK{};
537 const std::vector<double> *p_jK{nullptr};
538
539public:
540 Phi5k(const Phi5k &other)
541 : EM_multipole(other),
542 m_jK(other.m_jK),
543 p_jK(other.p_jK == &other.m_jK ? &m_jK : other.p_jK) {}
544 Phi5k &operator=(const Phi5k &other) {
545 if (this != &other) {
546 EM_multipole::operator=(other);
547 m_jK = other.m_jK;
548 p_jK = other.p_jK == &other.m_jK ? &m_jK : other.p_jK;
549 }
550 return *this;
551 }
552};
553
554//==============================================================================
555/*!
556 @brief Pseudoscalar multipole operator: t^k (i g^0 g^5)
557 @details
558
559 \f[ P_K = S^5_K = t^K(q)(i\gamma^0\gamma^5) \f]
560
561 - Implements the pseudoscalar multipole operator ~ \f$ e^{i q r} i \gamma^0 \gamma^5. \f$
562 - Radial dependence is provided via spherical Bessel functions \f$ j_L(q*r). \f$
563 - Supports an optional `const SphericalBessel::JL_table *jl` for precomputed
564 Bessel lookup; otherwise computes on demand.
565*/
566class S5k final : public EM_multipole {
567public:
568 S5k(const Grid &gr, int K, double omega,
569 const SphericalBessel::JL_table *jl = nullptr)
570 : EM_multipole(K, Angular::evenQ(K) ? Parity::odd : Parity::even, 1.0,
571 gr.r(), Realness::real, true, &gr, 'P', 'T', false, jl) {
572 if (omega != 0.0)
573 updateFrequency(omega);
574 }
575 DiracSpinor radial_rhs(const int kappa_a,
576 const DiracSpinor &Fb) const override final;
577
578 double radialIntegral(const DiracSpinor &Fa,
579 const DiracSpinor &Fb) const override final;
580
581 //! nb: q = alpha*omega!
582 void updateFrequency(const double omega) override final;
583
584private:
585 std::vector<double> m_jK{};
586 const std::vector<double> *p_jK{nullptr};
587
588public:
589 S5k(const S5k &other)
590 : EM_multipole(other),
591 m_jK(other.m_jK),
592 p_jK(other.p_jK == &other.m_jK ? &m_jK : other.p_jK) {}
593 S5k &operator=(const S5k &other) {
594 if (this != &other) {
595 EM_multipole::operator=(other);
596 m_jK = other.m_jK;
597 p_jK = other.p_jK == &other.m_jK ? &m_jK : other.p_jK;
598 }
599 return *this;
600 }
601};
602
603//==============================================================================
604//==============================================================================
605
606//! Helper functions for the multipole operators
607namespace multipole {
608
609//! Convert from "transition form" to "moment form"
610inline double moment_factor(int K, double omega) {
611 const auto q = std::abs(PhysConst::alpha * omega);
612 return qip::double_factorial(2 * K + 1) / qip::pow(q, K) *
613 std::sqrt(K / (K + 1.0));
614}
615
616} // namespace multipole
617
618//==============================================================================
619//==============================================================================
620//! @brief Factory class for Multipole operators (never instantiated directly).
621struct Multipole {
622 Multipole() = delete;
623 static std::unique_ptr<TensorOperator> generate(const IO::InputBlock &input,
624 const Wavefunction &wf) {
625 input.check({
626 {"",
627 "Note: This function cannot use the Spherical Bessel looup table. If "
628 "require efficiency for large number of q values, construct directly"},
629 {"k", "Rank: k=1 for E1, =2 for E2 etc. [1]"},
630 {"omega", "Frequency: nb: q := alpha*omega [1.0e-4]"},
631 {"type", "V,A,S,P (Vector, Axial, Scalar, Pseudoscalar) [V]"},
632 {"low_q", "bool. Use low-q formulas (K=0 and 1 only, no L-form) [false]"},
633 {"component", "E,M,L,T (electric, magnetic, longitudanel, temporal). "
634 "Temporal is forced if type = S or P. [E]"},
635 {"form", "L,V (Length, Velocity); only for electric vector [L]"},
636 });
637 if (input.has_option("help")) {
638 return nullptr;
639 }
640 const auto k = input.get("k", 1);
641 const auto omega = input.get("omega", 1.0e-4);
642
643 const auto low_q = input.get("low_q", false);
644
645 using namespace std::string_literals;
646 const auto type = input.get("type", "V"s);
647 const auto component = input.get("component", "E"s);
648 const auto form = input.get("form", "V"s);
649
650 const bool Vector = qip::ci_wc_compare(type, "V*");
651 const bool AxialVector = qip::ci_wc_compare(type, "A*");
652 const bool Scalar = qip::ci_wc_compare(type, "S*");
653 const bool PseudoScalar = qip::ci_wc_compare(type, "P*");
654
655 const bool Electric = qip::ci_wc_compare(component, "E*");
656 const bool Magnetic = qip::ci_wc_compare(component, "M*");
657 const bool Longitudinal = qip::ci_wc_compare(component, "L*");
658 const bool Temporal = qip::ci_wc_compare(component, "T*");
659
660 const bool LengthForm = qip::ci_wc_compare(form, "L*");
661
662 if (LengthForm && !(Electric && Vector)) {
663 std::cout << "Fail; Length form only valid for Electric Vector\n";
664 }
665
666 if (low_q) {
667 if (Electric && Vector)
668 return std::make_unique<VEk_lowq>(wf.grid(), k, omega);
669 if (Electric && AxialVector)
670 return std::make_unique<AEk_lowq>(wf.grid(), k, omega);
671
672 // Longitudinal
673 if (Longitudinal && Vector)
674 return std::make_unique<VLk_lowq>(wf.grid(), k, omega);
675 if (Longitudinal && AxialVector)
676 return std::make_unique<ALk_lowq>(wf.grid(), k, omega);
677
678 // Magnetic
679 if (Magnetic && Vector)
680 return std::make_unique<VMk_lowq>(wf.grid(), k, omega);
681 if (Magnetic && AxialVector)
682 return std::make_unique<AMk_lowq>(wf.grid(), k, omega);
683
684 // Temporal
685 if (Temporal && Vector)
686 return std::make_unique<Phik_lowq>(wf.grid(), k, omega);
687 if (Temporal && AxialVector)
688 return std::make_unique<Phi5k_lowq>(wf.grid(), k, omega);
689
690 if (Scalar)
691 return std::make_unique<Sk_lowq>(wf.grid(), k, omega);
692 if (PseudoScalar)
693 return std::make_unique<S5k_lowq>(wf.grid(), k, omega);
694 }
695
696 // Electric:
697 if (Electric && LengthForm && Vector)
698 return std::make_unique<VEk_Len>(wf.grid(), k, omega);
699 if (Electric && Vector)
700 return std::make_unique<VEk>(wf.grid(), k, omega);
701 if (Electric && AxialVector)
702 return std::make_unique<AEk>(wf.grid(), k, omega);
703
704 // Longitudinal
705 if (Longitudinal && Vector)
706 return std::make_unique<VLk>(wf.grid(), k, omega);
707 if (Longitudinal && AxialVector)
708 return std::make_unique<ALk>(wf.grid(), k, omega);
709
710 // Magnetic
711 if (Magnetic && Vector)
712 return std::make_unique<VMk>(wf.grid(), k, omega);
713 if (Magnetic && AxialVector)
714 return std::make_unique<AMk>(wf.grid(), k, omega);
715
716 // Temporal
717 if (Temporal && Vector)
718 return std::make_unique<Phik>(wf.grid(), k, omega);
719 if (Temporal && AxialVector)
720 return std::make_unique<Phi5k>(wf.grid(), k, omega);
721
722 if (Scalar)
723 return std::make_unique<Sk>(wf.grid(), k, omega);
724 if (PseudoScalar)
725 return std::make_unique<S5k>(wf.grid(), k, omega);
726
727 std::cout << "Fail; Invalid Combination\n";
728 return std::make_unique<NullOperator>();
729 }
730};
731
732//------------------------------------------------------------------------------
733/*!
734 @brief Factory for relativistic multipole operators.
735 @details
736 Constructs and returns a specific multipole operator derived from
737 DiracOperator::TensorOperator, based on the requested Lorentz structure,
738 multipole component, and momentum-transfer regime.
739
740 These are the \f$ t^K_Q \tilde\gamma \f$, \f$ T^{(\sigma)}_{KQ} \tilde\gamma \f$
741 operators from the vector expansion:
742
743 \f[
744 \begin{align}
745 e^{i\vec{q}\cdot\vec{r}}
746 &= \sqrt{4\pi}\sum_{KQ}\sqrt{[K]} \,
747 i^K \, {Y^*_{KQ}}{(\hat q)} \, t^K_Q(q,r),\\
748 \vec{\alpha} \, e^{i\vec{q}\cdot\vec{r}}
749 & = \sqrt{4\pi} \sum_{KQ\sigma} \sqrt{[K]} \, i^{K-\sigma} \,
750 \vec{Y}_{KQ}^{(\sigma)*}(\hat{{q}}) \,
751 T^{(\sigma)}_{KQ}.
752 \end{align}
753 \f]
754
755 The operator corresponds to a spherical multipole of rank @p k with
756 frequency/energy transfer @p omega. The type and component determine
757 the Lorentz structure and spatial character of the interaction.
758
759 These are "frequency"-dependent operators, via momentum transfer, q.
760 For electromagnetic interactions, \f$ \omega = q c\f$.
761 The "updateFrequency()" function expects these units, so even when momentum
762 transfer is not equal to energy exchange, we should pass \f$ qc \f$ to this function.
763
764 @param grid Radial grid on which the operator acts.
765 @param k Multipole rank (total angular momentum of the operator).
766 @param omega Energy (frequency) transfer.
767 @param type Interaction type:
768 - 'V' : Vector
769 - 'A' : Axial-vector
770 - 'S' : Scalar
771 - 'P' : Pseudoscalar
772 @param comp Multipole component (for V/A types):
773 - 'E' : Electric
774 - 'M' : Magnetic
775 - 'L' : Longitudinal
776 - 'T' : Temporal [caution - NOT transverse!]
777 (Ignored for scalar and pseudoscalar operators.)
778 @param low_q If true, construct the low-momentum (long-wavelength)
779 approximation of the operator.
780 @param jl Optional pointer to a precomputed spherical Bessel table.
781 If provided, radial Bessel functions are taken from this
782 table to avoid recomputation. If nullptr, they are generated
783 internally as needed.
784
785 @return A std::unique_ptr to the requested TensorOperator.
786
787 @note Length-form electric operators (if implemented separately) are
788 only valid for vector-electric ('V','E') combinations.
789
790 @warning Invalid combinations of @p type and @p comp may result in
791 a nullptr being returned.
792*/
793std::unique_ptr<DiracOperator::TensorOperator>
794MultipoleOperator(const Grid &grid, int k, double omega, char type, char comp,
795 bool low_q, const SphericalBessel::JL_table *jl = nullptr);
796
797} // namespace DiracOperator
Axial electric multipole operator: .
Definition EM_multipole.hpp:359
DiracSpinor radial_rhs(const int kappa_a, const DiracSpinor &Fb) const override final
Computes the right-hand spinor dF_b for the radial integral.
Definition EM_multipole.cpp:335
void updateFrequency(const double omega) override final
nb: q = alpha*omega!
Definition EM_multipole.cpp:398
double radialIntegral(const DiracSpinor &Fa, const DiracSpinor &Fb) const override final
Radial integral R_ab, defined by RME = angularF(a,b) * radialIntegral(a,b).
Definition EM_multipole.cpp:368
Axial longitudinal multipole operator: .
Definition EM_multipole.hpp:413
DiracSpinor radial_rhs(const int kappa_a, const DiracSpinor &Fb) const override final
Computes the right-hand spinor dF_b for the radial integral.
Definition EM_multipole.cpp:417
void updateFrequency(const double omega) override final
nb: q = alpha*omega!
Definition EM_multipole.cpp:453
double radialIntegral(const DiracSpinor &Fa, const DiracSpinor &Fb) const override final
Radial integral R_ab, defined by RME = angularF(a,b) * radialIntegral(a,b).
Definition EM_multipole.cpp:439
Axial magnetic multipole operator: .
Definition EM_multipole.hpp:467
DiracSpinor radial_rhs(const int kappa_a, const DiracSpinor &Fb) const override final
Computes the right-hand spinor dF_b for the radial integral.
Definition EM_multipole.cpp:475
void updateFrequency(const double omega) override final
nb: q = alpha*omega!
Definition EM_multipole.cpp:514
double radialIntegral(const DiracSpinor &Fa, const DiracSpinor &Fb) const override final
Radial integral R_ab, defined by RME = angularF(a,b) * radialIntegral(a,b).
Definition EM_multipole.cpp:497
Intermediate abstract base class for all EM relativistic multipole operators.
Definition EM_multipole_base.hpp:50
const SphericalBessel::JL_table * jl() const
Returns the precomputed Bessel table pointer (may be nullptr).
Definition EM_multipole_base.hpp:83
Temporal component of the axial vector multipole operator.
Definition EM_multipole.hpp:517
DiracSpinor radial_rhs(const int kappa_a, const DiracSpinor &Fb) const override final
Computes the right-hand spinor dF_b for the radial integral.
Definition EM_multipole.cpp:531
void updateFrequency(const double omega) override final
nb: q = alpha*omega!
Definition EM_multipole.cpp:559
double radialIntegral(const DiracSpinor &Fa, const DiracSpinor &Fb) const override final
Radial integral R_ab, defined by RME = angularF(a,b) * radialIntegral(a,b).
Definition EM_multipole.cpp:548
Temporal component of the vector multipole operator: .
Definition EM_multipole.hpp:256
DiracSpinor radial_rhs(const int kappa_a, const DiracSpinor &Fb) const override final
Computes the right-hand spinor dF_b for the radial integral.
Definition EM_multipole.cpp:248
double radialIntegral(const DiracSpinor &Fa, const DiracSpinor &Fb) const override final
Radial integral R_ab, defined by RME = angularF(a,b) * radialIntegral(a,b).
Definition EM_multipole.cpp:265
void updateFrequency(const double omega) override final
nb: q = alpha*omega!
Definition EM_multipole.cpp:276
Pseudoscalar multipole operator: t^k (i g^0 g^5)
Definition EM_multipole.hpp:566
double radialIntegral(const DiracSpinor &Fa, const DiracSpinor &Fb) const override final
Radial integral R_ab, defined by RME = angularF(a,b) * radialIntegral(a,b).
Definition EM_multipole.cpp:592
DiracSpinor radial_rhs(const int kappa_a, const DiracSpinor &Fb) const override final
Computes the right-hand spinor dF_b for the radial integral.
Definition EM_multipole.cpp:575
void updateFrequency(const double omega) override final
nb: q = alpha*omega!
Definition EM_multipole.cpp:602
Scalar multipole operator: .
Definition EM_multipole.hpp:306
DiracSpinor radial_rhs(const int kappa_a, const DiracSpinor &Fb) const override final
Computes the right-hand spinor dF_b for the radial integral.
Definition EM_multipole.cpp:292
void updateFrequency(const double omega) override final
nb: q = alpha*omega!
Definition EM_multipole.cpp:319
double radialIntegral(const DiracSpinor &Fa, const DiracSpinor &Fb) const override final
Radial integral R_ab, defined by RME = angularF(a,b) * radialIntegral(a,b).
Definition EM_multipole.cpp:309
Vector electric multipole (transition) operator, Length-form: ( ).
Definition EM_multipole.hpp:30
DiracSpinor radial_rhs(const int kappa_a, const DiracSpinor &Fb) const override final
Computes the right-hand spinor dF_b for the radial integral.
Definition EM_multipole.cpp:9
double radialIntegral(const DiracSpinor &Fa, const DiracSpinor &Fb) const override final
Radial integral R_ab, defined by RME = angularF(a,b) * radialIntegral(a,b).
Definition EM_multipole.cpp:31
void updateFrequency(const double omega) override final
nb: q = alpha*omega!
Definition EM_multipole.cpp:47
Vector electric multipole (V-form) operator: .
Definition EM_multipole.hpp:88
void updateFrequency(const double omega) override final
nb: q = alpha*omega!
Definition EM_multipole.cpp:112
double radialIntegral(const DiracSpinor &Fa, const DiracSpinor &Fb) const override final
Radial integral R_ab, defined by RME = angularF(a,b) * radialIntegral(a,b).
Definition EM_multipole.cpp:93
DiracSpinor radial_rhs(const int kappa_a, const DiracSpinor &Fb) const override final
Computes the right-hand spinor dF_b for the radial integral.
Definition EM_multipole.cpp:67
Vector longitudinal multipole operator (V-form): .
Definition EM_multipole.hpp:147
DiracSpinor radial_rhs(const int kappa_a, const DiracSpinor &Fb) const override final
Computes the right-hand spinor dF_b for the radial integral.
Definition EM_multipole.cpp:134
void updateFrequency(const double omega) override final
nb: q = alpha*omega!
Definition EM_multipole.cpp:172
double radialIntegral(const DiracSpinor &Fa, const DiracSpinor &Fb) const override final
Radial integral R_ab, defined by RME = angularF(a,b) * radialIntegral(a,b).
Definition EM_multipole.cpp:158
Vector magnetic multipole operator: .
Definition EM_multipole.hpp:205
void updateFrequency(const double omega) override final
nb: q = alpha*omega!
Definition EM_multipole.cpp:232
DiracSpinor radial_rhs(const int kappa_a, const DiracSpinor &Fb) const override final
Computes the right-hand spinor dF_b for the radial integral.
Definition EM_multipole.cpp:194
double radialIntegral(const DiracSpinor &Fa, const DiracSpinor &Fb) const override final
Radial integral R_ab, defined by RME = angularF(a,b) * radialIntegral(a,b).
Definition EM_multipole.cpp:216
s (spin) operator
Definition jls.hpp:57
Stores radial Dirac spinor: F_nk = (f, g)
Definition DiracSpinor.hpp:42
Non-uniform radial grid with Jacobian, suitable for atomic structure calculations.
Definition Grid.hpp:85
const std::vector< double > & r() const
Full grid vector r.
Definition Grid.hpp:131
Holds a named list of key=value options and nested InputBlocks.
Definition InputBlock.hpp:154
bool check(std::initializer_list< std::string > blocks, const std::vector< std::pair< std::string, std::string > > &list, bool print=false) const
Validates options and sub-blocks against an allowed list.
Definition InputBlock.hpp:649
bool has_option(std::string_view key) const
Returns true if key is present in this block's option list, even if unset.
Definition InputBlock.hpp:247
T get(std::string_view key, T default_value) const
Returns the value of key, or default_value if not found.
Definition InputBlock.hpp:471
Lookup table of spherical Bessel functions: j_L(q*r) = J[L][q][r].
Definition SphericalBessel.hpp:124
Stores Wavefunction (set of valence orbitals, grid, HF etc.)
Definition Wavefunction.hpp:37
const Grid & grid() const
Returns a const reference to the radial grid.
Definition Wavefunction.hpp:82
constexpr bool evenQ(int a)
Returns true if a is even - for integer values.
Definition Wigner369j.hpp:233
double moment_factor(int K, double omega)
Convert from "transition form" to "moment form".
Definition EM_multipole.hpp:610
Dirac operators: TensorOperator base class and derived implementations for single-particle (one-body)...
Definition GenerateOperator.cpp:6
std::unique_ptr< DiracOperator::TensorOperator > MultipoleOperator(const Grid &grid, int k, double omega, char type, char comp, bool low_q, const SphericalBessel::JL_table *jl)
Factory for relativistic multipole operators.
Definition EM_multipole.cpp:629
constexpr double alpha
Fine-structure constant: alpha = 1/137.035 999 177(21) [CODATA 2022].
Definition PhysConst_constants.hpp:24
constexpr double double_factorial(T x)
Double factorial x!! - takes integer, returns double.
Definition Maths.hpp:141
bool ci_wc_compare(std::string_view s1, std::string_view s2)
Case-insensitive version of wildcard_compare.
Definition String.hpp:155
constexpr auto pow(T x)
x^n for compile-time integer n, x any arithmetic type.
Definition Maths.hpp:98
Factory class for Multipole operators (never instantiated directly).
Definition EM_multipole.hpp:621