2#include "DiracOperator/TensorOperator.hpp"
3#include "IO/InputBlock.hpp"
4#include "Physics/NuclearData.hpp"
5#include "Potentials/NuclearPotentials.hpp"
6#include "Wavefunction/Wavefunction.hpp"
7#include "fmt/color.hpp"
8#include "qip/Maths.hpp"
50std::vector<double>
tk_radial(
int k,
double rN,
const std::vector<double> &r,
128 double R,
bool u_option,
bool print =
true);
166 double l1,
int gl1,
double I2,
double l2,
271 using RadialFunction = std::function<double(
double,
double)>;
289 hfs(
int in_k,
double GQ,
double rN_au,
const Grid &rgrid,
292 Hyperfine::tk_radial(in_k, rN_au, rgrid.r(), hfs_F)),
294 magnetic(k % 2 != 0),
295 cfg(magnetic ? 1.0 : 0.0),
296 cff(magnetic ? 0.0 : 1.0),
299 const auto power = magnetic ? (k - 1) / 2 : k / 2;
305 std::pow(PhysConst::barn_au, power);
306 m_unit = mMHzQ ? unit_au * PhysConst::Hartree_MHz : 1.0;
309 std::string
name() const override final {
310 return "hfs" + std::to_string(k) +
"";
312 std::string
units() const override final {
return mMHzQ ?
"MHz" :
"au"; }
314 double angularF(
const int ka,
const int kb)
const override final {
315 return magnetic ? -double(ka + kb) / double(k) *
320 double angularCff(
int,
int)
const override final {
return cff; }
321 double angularCgg(
int,
int)
const override final {
return cff; }
322 double angularCfg(
int,
int)
const override final {
return cfg; }
323 double angularCgf(
int,
int)
const override final {
return cfg; }
325 static std::unique_ptr<TensorOperator> generate(
const IO::InputBlock &input,
329 {{
"",
"Most following will be taken from the default nucleus if "
330 "not explicitely given"},
331 {
"mu",
"Magnetic moment in mu_N"},
332 {
"Q",
"Nuclear quadrupole moment, in barns. Also used as overall "
333 "constant for any higher-order moments [1.0]"},
334 {
"k",
"Multipolarity. 1=mag. dipole, 2=elec. quad, etc. [1]"},
336 "nuclear (magnetic) rms radius, in Fermi (fm) (defult is charge rms)"},
338 "Units for output (only for k=1,k=2). MHz or au. For MHz, "
339 "assumes nuclear moments g/Q are in mu_N*barns^p. For atomic "
340 "units, best to set g=Q=1 to get raw t^k matrix elements [MHz]"},
341 {
"F",
"F(r): Nuclear moment distribution: ball, point, shell, "
342 "SingleParticle, doublyOddSP, uSP, Fermi, or Gaussian [ball]"},
343 {
"",
"The following options are all for F(r) details. See "
344 "https://ampsci.dev for full details (search `hyperfine`)"},
345 {
"printF",
"Writes F(r) to a text file [false]"},
346 {
"print",
"Write F(r) info to screen [true]"},
347 {
"",
"The following are only for F=SingleParticle or doublyOddSP"},
348 {
"I",
"Nuclear spin. Taken from nucleus"},
349 {
"parity",
"Nulcear parity: +/-1"},
350 {
"l",
"l for unpaired nucleon (automatically derived from I and "
351 "parity; best to leave as default)"},
352 {
"gl",
"=1 for proton, =0 for neutron"},
353 {
"",
"The following are only used if F=doublyOddSP"},
354 {
"mu1",
"mag moment of 'first' unpaired nucleon"},
355 {
"gl1",
"gl of 'first' unpaired nucleon"},
356 {
"l1",
"l of 'first' unpaired nucleon"},
357 {
"l2",
"l of 'second' unpaired nucleon"},
358 {
"I1",
"total spin (J) of 'first' unpaired nucleon"},
359 {
"I2",
"total spin (J) of 'second' unpaired nucleon"},
360 {
"",
"The following are only for u(r) function (F=uSP)"},
361 {
"u",
"u1 or u2 : u1(r) = (R-r)^n, u2(r) = r^n [u1]"},
362 {
"n",
"n that appears above. Should be between 0 and 2 [0]"},
363 {
"",
"The following are only for Fermi function (F=Fermi)"},
365 "skin thickness (normally 2.3). note: c taken from wavefunction, or "
366 "from rrms above [2.3]"}});
374 auto mu = input.
get(
"mu", isotope.mu ? *isotope.mu : 1.0);
375 auto I_nuc = input.
get(
"I", isotope.I_N ? *isotope.I_N : 1.0);
376 const auto print = input.
get(
"print",
true);
377 const auto k = input.
get(
"k", 1);
383 fmt2::styled_print(fg(fmt::color::red),
"\nError 246:\n");
384 std::cout <<
"In hyperfine: invalid K=" << k <<
"! meaningless results\n";
386 if (I_nuc <= 0 && k == 1) {
387 fmt2::styled_print(fg(fmt::color::orange),
"\nWarning 253:\n");
388 std::cout <<
"In hyperfine: invalid I_nuc=" << I_nuc
389 <<
"! meaningless results\nSetting I=1\n";
392 if (mu == 0.0 && k == 1) {
393 fmt2::styled_print(fg(fmt::color::orange),
"\nWarning 352:\n");
394 std::cout <<
"Setting mu=1\n";
399 (k == 1) ? (mu / I_nuc) :
400 (k == 2) ? input.
get(
"Q", isotope.q ? *isotope.q : 1.0) :
403 enum class DistroType {
415 const std::string default_distribution =
"ball";
420 input.
get<std::string>(
"F", default_distribution) :
421 input.has_option(
"nuc_mag") ?
422 input.get<std::string>(
"nuc_mag", default_distribution) :
423 input.get<std::string>(
"F(r)", default_distribution);
425 const auto distro_type =
437 if (distro_type == DistroType::Error) {
438 fmt2::styled_print(fg(fmt::color::red),
"\nError 271:\n");
439 std::cout <<
"\nIn hyperfine. Unkown F(r) - " << Fr_str <<
"\n";
440 std::cout <<
"Defaulting to pointlike!\n";
444 distro_type == DistroType::point ? 0.0 : input.
get(
"rrms", nuc.r_rms());
445 const auto r_nucfm = std::sqrt(5.0 / 3) * r_rmsfm;
446 const auto r_nucau = r_nucfm / PhysConst::aB_fm;
449 std::cout <<
"\nHyperfine structure: " << wf.
atom() <<
"\n";
450 std::cout <<
"K=" << k <<
" ("
451 << (k == 1 ?
"magnetic dipole" :
452 k == 2 ?
"electric quadrupole" :
453 k % 2 == 0 ?
"electric multipole" :
454 "magnetic multipole")
456 std::cout <<
"Using " << Fr_str <<
" nuclear distro for F(r)\n"
457 <<
"w/ r_N = " << r_nucfm <<
"fm = " << r_nucau
458 <<
"au (r_rms=" << r_rmsfm <<
"fm)\n";
459 std::cout <<
"Points inside nucleus: " << wf.
grid().
getIndex(r_nucau)
462 std::cout <<
"mu = " << mu <<
", I = " << I_nuc <<
", g = " << g_or_Q
465 std::cout <<
"Q = " << g_or_Q <<
"\n";
472 if (distro_type == DistroType::ball) {
475 }
else if (distro_type == DistroType::shell) {
478 }
else if (distro_type == DistroType::SingleParticle) {
479 const auto pi = input.
get(
"parity", isotope.parity ? *isotope.parity : 0);
480 const auto l_tmp = int(I_nuc + 0.5 + 0.0001);
481 auto l = ((l_tmp % 2 == 0) == (pi == 1)) ? l_tmp : l_tmp - 1;
482 l = input.
get(
"l", l);
483 const auto gl_default = wf.
Znuc() % 2 == 0 ? 0 : 1;
484 const auto gl = input.
get<
int>(
"gl", gl_default);
486 std::cout <<
"Single-Particle (Volotka formula) for unpaired";
488 std::cout <<
" proton ";
490 std::cout <<
" neturon ";
492 std::cout <<
" gl=" << gl
493 <<
"??? program will run, but prob wrong!\n";
494 std::cout <<
"with l=" << l <<
" (pi=" << pi <<
")\n";
498 }
else if (distro_type == DistroType::uSP) {
499 const auto pi = input.
get(
"parity", isotope.parity ? *isotope.parity : 0);
501 input.
get(
"u", std::string{
"u1"});
502 const bool u_option = u_func == std::string{
"u1"};
503 const auto n = input.
get(
"n", 1.0);
504 const auto l_tmp = int(I_nuc + 0.5 + 0.0001);
505 auto l = ((l_tmp % 2 == 0) == (pi == 1)) ? l_tmp : l_tmp - 1;
506 l = input.
get(
"l", l);
507 const auto gl_default = wf.
Znuc() % 2 == 0 ? 0 : 1;
508 const auto gl = input.
get<
int>(
"gl", gl_default);
510 std::cout <<
"Single-Particle (Volotka formula) with u(r) for unpaired";
512 std::cout <<
" proton ";
514 std::cout <<
" neturon ";
516 std::cout <<
" gl=" << gl
517 <<
"??? program will run, but prob wrong!\n";
518 std::cout <<
"with l=" << l <<
" (pi=" << pi <<
")\n";
521 std::cout <<
"u(r) = " << u_func <<
": "
522 << (u_option ?
"(R-r)^n" :
"r^n") <<
", n=" << n <<
"\n";
524 Fr =
Hyperfine::uSP(mu, I_nuc, l, gl, n, r_nucau, u_option, print);
526 }
else if (distro_type == DistroType::doublyOddSP) {
527 const auto mu1 = input.
get<
double>(
"mu1", 1.0);
528 const auto gl1 = input.
get<
int>(
"gl1", -1);
529 if (gl1 != 0 && gl1 != 1) {
530 fmt2::styled_print(fg(fmt::color::red),
"\nError 324:\n");
531 std::cout <<
"In " << input.
name() <<
" " << Fr_str
532 <<
"; have gl1=" << gl1 <<
" but need 1 or 0\n";
533 return std::make_unique<NullOperator>(NullOperator());
535 const auto l1 = input.
get<
double>(
"l1", -1.0);
536 const auto l2 = input.
get<
double>(
"l2", -1.0);
537 const auto I1 = input.
get<
double>(
"I1", -1.0);
538 const auto I2 = input.
get<
double>(
"I2", -1.0);
542 }
else if (distro_type == DistroType::Fermi) {
543 const auto t_fm = input.
get(
"t", nuc.t());
546 std::cout <<
"Fermi distro: c=" << c_fm <<
" fm, t=" << t_fm <<
" fm\n";
551 }
else if (distro_type == DistroType::Gaussian) {
552 const auto rN_au = r_rmsfm / PhysConst::aB_fm;
554 std::cout <<
"Gaussian distro: r_rms=" << r_rmsfm <<
" fm\n";
555 std::vector<double> rho;
557 for (
const auto ri : wf.grid().r()) {
558 rho.push_back(std::exp(-1.5 * ri * ri / (rN_au * rN_au)));
561 }
else if (distro_type == DistroType::Error) {
562 std::cout <<
"Using pointlike F(r):\n";
567 if (input.
get<
bool>(
"printF",
false)) {
568 std::ofstream of(wf.
identity() +
"_" + Fr_str +
".txt");
570 for (
auto r : wf.grid()) {
571 of << r * PhysConst::aB_fm <<
" "
572 << Fr(r * PhysConst::aB_fm, r_nucau * PhysConst::aB_fm) <<
"\n";
576 return std::make_unique<hfs>(k, g_or_Q, r_nucau, wf.
grid(), Fr, use_MHz);
General tensor operator (virtual base class); all single-particle (one-body) tenosor operators derive...
Definition TensorOperator.hpp:198
Generalised hyperfine-structure operator, including relevant nuclear moment.
Definition hfs.hpp:270
double angularCgg(int, int) const override final
Angular coefficient C_gg for the g_a*g_b term of the radial integral.
Definition hfs.hpp:321
double angularCgf(int, int) const override final
Angular coefficient C_gf for the g_a*f_b term of the radial integral.
Definition hfs.hpp:323
double angularCfg(int, int) const override final
Angular coefficient C_fg for the f_a*g_b term of the radial integral.
Definition hfs.hpp:322
std::string units() const override final
Returns units of operator as a string (usually au, may be MHz, etc.)
Definition hfs.hpp:312
std::string name() const override final
Returns "name" of operator (e.g., 'E1')
Definition hfs.hpp:309
hfs(int in_k, double GQ, double rN_au, const Grid &rgrid, const RadialFunction &hfs_F=Hyperfine::pointlike_F(), bool MHzQ=true)
Constructs hyperfine operator of multipolarity k.
Definition hfs.hpp:289
double angularF(const int ka, const int kb) const override final
Angular factor A_ab linking the radial integral to the RME.
Definition hfs.hpp:314
double angularCff(int, int) const override final
Angular coefficient C_ff for the f_a*f_b term of the radial integral.
Definition hfs.hpp:320
Non-uniform radial grid with Jacobian, suitable for atomic structure calculations.
Definition Grid.hpp:85
auto num_points() const
Number of grid points.
Definition Grid.hpp:120
std::size_t getIndex(double x, bool require_nearest=false) const
Returns the grid index closest to x. If require_nearest is false, returns the index of the largest gr...
Definition Grid.cpp:76
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
int Znuc() const
Nuclear charge, Z.
Definition Wavefunction.hpp:99
std::string identity() const
Atomic symbol, including core ionisation degree and run_label.
Definition Wavefunction.cpp:690
std::string atom() const
String of atom info (e.g., "Cs, Z=55, A=133")
Definition Wavefunction.hpp:190
const Nuclear::Nucleus & nucleus() const
Returns Nuclear::nucleus object (contains nuc. parameters)
Definition Wavefunction.hpp:97
double Ck_kk(int k, int ka, int kb)
Reduced relativistic angular ME: <ka||C^k||kb>. Takes kappa values.
Definition Wigner369j.hpp:486
double convert_RME_to_HFSconstant_2J(int k, int tja, int tjb)
Converts reduced matrix element to A/B coeficients (takes k, 2J, 2J)
Definition hfs.cpp:206
RadialFunction doublyOddSP_F(double mut, double It, double mu1, double I1, double l1, int gl1, double I2, double l2, bool print)
Volotka single-particle model for doubly-odd nuclei.
Definition hfs.cpp:150
std::vector< double > tk_radial(int k, double rN, const std::vector< double > &r, const RadialFunction &hfs_F)
Forms the hyperfine radial operator: F(r, r_N) / r^{k+1}.
Definition hfs.cpp:14
double convert_RME_to_HFSconstant(int k, int ka, int kb)
Converts reduced matrix element to A/B coeficients (takes k, kappa, kappa)
Definition hfs.cpp:240
RadialFunction pointlike_F()
Pointlike F(r): F = 1 everywhere.
Definition hfs.cpp:38
RadialFunction sphericalBall_F(int k)
Spherical ball model for F(r, r_N) [default]. Uniformly distributed point k-poles.
Definition hfs.cpp:26
RadialFunction sphericalShell_F()
Spherical shell F(r, r_N): 0 for r < r_N, 1 for r > r_N.
Definition hfs.cpp:33
RadialFunction generic_F(const Grid &grid, const std::vector< double > &rho, int k)
Constructs F(r) from a numerical density rho(r) given on the radial grid.
Definition hfs.cpp:168
RadialFunction uSP(double mu, double I_nuc, double l_pn, int gl, double n, double R, bool u_option, bool print)
Elizarov single-particle magnetisation model [extended Volotka].
Definition hfs.cpp:81
std::function< double(double r, double r_Nuc)> RadialFunction
Type for radial function F(r, r_Nuc) – nuclear magnetisation/charge distribution.
Definition hfs.hpp:30
RadialFunction VolotkaSP_F(double mu, double I_nuc, double l_pn, int gl, bool print)
Volotka single-particle nuclear model: F(r, rN)
Definition hfs.cpp:44
Dirac operators: TensorOperator base class and derived implementations for single-particle (one-body)...
Definition GenerateOperator.cpp:6
Parity
Parity of operator.
Definition TensorOperator.hpp:57
double c_hdr_formula_rrms_t(double rrms, double t)
Calculates c from rrms and t.
Definition NuclearData.cpp:73
Isotope findIsotopeData(int z, int a)
Looks up + returns an isotope from the list. If not in list, partially blank.
Definition NuclearData.cpp:6
std::vector< double > fermiNuclearDensity_tcN(double t, double c, double Z_norm, const Grid &grid)
Fermi charge distribution, rho(r) - normalised to Z_norm.
Definition NuclearPotentials.cpp:319
constexpr double muN_CGS
Nulcear magneton (in Gaussian CGS-derived atomic units):
Definition PhysConst_constants.hpp:101
General-purpose utility library.
Definition Array.hpp:23
bool ci_wc_compare(std::string_view s1, std::string_view s2)
Case-insensitive version of wildcard_compare.
Definition String.hpp:155
bool ci_compare(std::string_view s1, std::string_view s2)
Case-insensitive string comparison; equivalent to tolower(s1) == tolower(s2).
Definition String.hpp:143