TkN 2.4
Toolkit for Nuclei
Loading...
Searching...
No Matches
tknucleus.cpp
1/********************************************************************************
2 * Copyright (c) : Université de Lyon 1, CNRS/IN2P3, UMR5822, *
3 * IP2I, F-69622 Villeurbanne Cedex, France *
4 * Normandie Université, ENSICAEN, UNICAEN, CNRS/IN2P3, *
5 * LPC Caen, F-14000 Caen, France *
6 * Contibutor(s) : *
7 * Jérémie Dudouet jeremie.dudouet@cnrs.fr [2020] *
8 * Diego Gruyer diego.gruyer@cnrs.fr [2020] *
9 * *
10 * Licensed under the MIT License <http://opensource.org/licenses/MIT>. *
11 * SPDX-License-Identifier: MIT *
12 ********************************************************************************/
13
14#include "tknucleus.h"
15#include "sqlite3.h"
16#include "tklog.h"
17#include "tkstring.h"
18#include "tkmanager.h"
19
20namespace tkn {
29}
30
31using namespace tkn;
32
33vector<int> tknucleus::fmagic_number_list{2,8,20,28,50,82,126};
34
40tknucleus::tknucleus(int _Z, int _A) : tkproperty_list("nucleus")
41{
42 tkstring symbol = gmanager->get_element_symbol(_Z);
43 if(symbol.empty()) {
44 glog << warning << "Unknown nucleus. No properties loaded " << tkstring::form("(%d,%d)",_Z,_A) << do_endl;
45 fZ = _Z;
46 fN = _A-_Z;
47 fA = _A;
48 return;
49 }
50 symbol.prepend(tkstring::form("%d",_A));
51 if(!gmanager->known_nucleus(symbol)){
52 glog << warning << "Unknown nucleus. No properties loaded (" << symbol << ")" << do_endl;
53 return;
54 }
55 const tknucleus & nuc = *gmanager->get_nucleus(symbol).get(); // vérifier que ca fasse une copie et pas une reference
56 (*this) = nuc;
57}
58
64tknucleus::tknucleus(int _Z) : tkproperty_list("nucleus")
65{
66 if(!gmanager->known_element(_Z)) glog << warning << "Unknown element. No properties loaded" << do_endl;
67 else {
68 shared_ptr<tknucleus> bestnuc;
69 for(const auto &nuc: gmanager->get_nuclei_for_z(_Z)) {
70 if(bestnuc==nullptr) bestnuc = nuc;
71 else if(isnan(bestnuc->get_lifetime()) || (nuc->get_lifetime()>bestnuc->get_lifetime())) bestnuc = nuc;
72 else if(nuc->get_lifetime() == bestnuc->get_lifetime() && nuc->get_abundance()>bestnuc->get_abundance()) bestnuc = nuc;
73 }
74 const tknucleus & nuc = *bestnuc.get();
75 (*this) = nuc;
76 }
77}
78
84tknucleus::tknucleus(const char *_symbol) : tkproperty_list("nucleus")
85{
86 if(!gmanager->known_nucleus(tkstring(_symbol))) {
87 int z;
88 if(gmanager->known_element(tkstring(_symbol),z)) {
89 if(!tkstring(_symbol).is_alpha()) {
90 glog << warning << "Unknown nucleus. TkN will choose one for you ;)" << do_endl;
91 }
92 const tknucleus nuc(z);
93 (*this) = nuc;
94 return;
95 }
96 }
97 const tknucleus & nuc = *gmanager->get_nucleus(_symbol).get();
98 (*this) = nuc;
99}
100
101void tknucleus::set_name(const char *_name)
102{
103 fSymbol = _name;
104 felement_symbol = fSymbol.extract_alpha();
105 fA = fSymbol.atoi();
106}
107
109 if(std::find(fmagic_number_list.begin(), fmagic_number_list.end(), fZ) != fmagic_number_list.end()) return true;
110 return false;
111}
112
114 if(std::find(fmagic_number_list.begin(), fmagic_number_list.end(), fN) != fmagic_number_list.end()) return true;
115 return false;
116}
117
118shared_ptr<tklevel_scheme> tknucleus::get_level_scheme()
119{
120 return gmanager->get_level_scheme(get_symbol(), get_z(), get_a());
121}
122
130double tknucleus::get_radius() const
131{
132 if(!has_property("radius")) return -1.;
133 return get("radius")->get_value();
134}
135
144{
145 if(!has_property("magnetic_dipole")) return -1.;
146 return get("magnetic_dipole")->get_value();
147}
148
157{
158 if(!has_property("electric_quadrupole")) return -1.;
159 return get("electric_quadrupole")->get_value();
160}
161
169double tknucleus::get_abundance() const
170{
171 if(!has_property("abundance")) return -1.;
172 return get("abundance")->get_value();
173}
174
183{
184 if(!has_property("mass_excess")) return std::nan("1");
185 if(gunits->get_type(_unit)!=tkunit_manager::units_type::kEnergy){
186 glog << error << gunits->get_name(_unit) << " is not of 'energy' type (returns nan(1))" << do_endl;
187 return std::nan("1");
188 }
189 return get("mass_excess")->get_value(_unit);
190}
191
200{
201 if(!has_property("binding_energy_overA")) return std::nan("1");
202 if(gunits->get_type(_unit)!=tkunit_manager::units_type::kEnergy){
203 glog << error << gunits->get_name(_unit) << " is not of 'energy' type (returns nan(1))" << do_endl;
204 return std::nan("1");
205 }
206 return get("binding_energy_overA")->get_value(_unit);
207}
208
217{
218 if(!has_property("lifetime")) return std::nan("1");
219 if(gunits->get_type(_unit)!=tkunit_manager::units_type::kTime&&gunits->get_type(_unit)!=tkunit_manager::units_type::kEnergy){
220 glog << error << gunits->get_name(_unit) << " is not of 'time' or 'energy' type (returns nan(1))" << do_endl;
221 return std::nan("1");
222 }
223 return get("lifetime")->get_value(_unit);
224}
225
226
227
241double tknucleus::get_fission_yield(tkstring _parent, bool _cumulative) {
242 if(!_parent.begins_with("FY")) _parent.prepend("FY");
243 if(_cumulative) _parent.prepend("c");
244 if(!has_property(_parent)) return 0.;
245 return get(_parent)->get_value();
246}
247
253{
254 if(!has_property("lifetime")) return "";
255 tkstring lifetime_str;
256
257 if(is_stable())
258 lifetime_str = "STABLE";
259 else {
260 vector<tkunit_manager::units_keys> units{tkunit_manager::y, tkunit_manager::d, tkunit_manager::h, tkunit_manager::min,
264
265 for(auto &unit: units) {
266 double lifetime = get_lifetime(unit);
267 if(lifetime>1.) {
268 lifetime_str = tkstring::Form("%.3g %s",lifetime,std::get<0>(tkunit_manager::funits_properties.at(unit)).data());
269 break;
270 }
271 }
272 }
273
274 return lifetime_str;
275}
276
277shared_ptr<tklevel> tknucleus::get_ground_state()
279 if(!get_level_scheme()->get_levels().empty())
280 return get_level_scheme()->get_levels().front();
281 return nullptr;
282}
283
284bool tknucleus::is_bound() const
285{
286 if(!has_property("binding_energy_overA")) return false;
287 return get("binding_energy_overA")->get_value()>0;
288}
289
292
295vector<pair<tkstring, double>> tknucleus::get_decay_modes()
296{
297 vector<pair<tkstring, double>> decay_list;
298 auto decays = get_decay_mode_str();
299 auto tokkens = decays.replace_all("[","]").tokenize("]");
300 for(auto &dec: tokkens) {
301 auto tokk = dec.tokenize(";");
302 decay_list.push_back({tokk.front(),tokk.back().atof()});
303 }
304 return decay_list;
305}
306
307void tknucleus::print() const
308{
309 glog << info << fSymbol << " (Z=" << fZ <<", N=" << fN << ") properties:" << do_endl;
310 glog << info << "Ground state configuration: ";
311 if(get_jpi().is_empty()) glog << "Unkown" << do_endl;
312 else glog << get_jpi() << do_endl;
313 if(is_stable()) glog << info << "Stable nucleus, abundance: " << tkstring::form("%3.2f",get_abundance()) << " %" << do_endl;
314 else {
315 glog << info << "Radioactive nucleus:"<<do_endl;
317 glog << info << "Decay: " << get_property("decay_modes") << do_endl;
318 }
319}
320
321
322#ifdef HAS_ROOT
324ClassImp(tknucleus);
325#endif
A nucleus made of Z protons and N neutrons.
Definition tknucleus.h:31
double get_mass_excess(const tkunit_manager::units_keys _unit=tkunit_manager::units_keys::keV)
return the mass excess in keV by default
const tkstring & get_jpi() const
return the ground state spin and parity as a string
Definition tknucleus.h:146
int get_a()
return the mass number
Definition tknucleus.h:69
vector< pair< tkstring, double > > get_decay_modes()
return the decay modes as a vector of decay
double get_radius() const
returns the radius in fm
bool is_n_magic()
return true in case of n is a magic number
double get_lifetime(const tkunit_manager::units_keys _unit=tkunit_manager::units_keys::s)
returns the lifetime in second
bool is_stable() const
test if the nucleus is stable
Definition tknucleus.h:139
double get_abundance() const
returns the natural abundance in percent
tknucleus()
default constructor
Definition tknucleus.h:53
double get_fission_yield(tkstring _parent, bool _cumulative=false)
returns the fission yield of the current nucleus
tkstring get_lifetime_str()
returns the lifetime in string (using the best adapted unit)
shared_ptr< tklevel > get_ground_state()
return the ground state level (nullptr if no GS known)
bool is_z_magic()
return true in case of z is a magic number
int get_z()
return the proton number
Definition tknucleus.h:65
shared_ptr< tklevel_scheme > get_level_scheme()
return a tklevel_scheme shared pointer to the nucleus level scheme
double get_electric_quadrupole()
return the electric quadrupole moment in barn by default
shared_ptr< tkmeasure > get_lifetime_measure() const
returns the lifetime tkmeasure object
Definition tknucleus.h:115
void print() const
print the main nucleus properties
const tkstring & get_symbol()
return the nucleus symbol
Definition tknucleus.h:72
double get_magnetic_dipole()
return the magnetic dipole moment in mun by default
tkstring get_decay_mode_str()
return the decay modes as a string
Definition tknucleus.h:128
double get_binding_energy_over_a(const tkunit_manager::units_keys _unit=tkunit_manager::units_keys::MeV)
returns the binding energy per mass unit in MeV
bool is_bound() const
test if the nucleus is bound (positive binding energy)
Contains list of properties.
bool has_property(const tkstring &_property) const
to check if the property is available
shared_ptr< tkmeasure > get(const tkstring &_property) const
get the property as tkmeasure
tkstring get_property(const tkstring &_property) const
get the property value as a string
std::string with usefull tricks from TString (ROOT) and KVString (KaliVeda) and more....
Definition tkstring.h:31
static const char * form(const char *_format,...)
Definition tkstring.cpp:408
static tkstring Form(const char *_format,...)
Definition tkstring.cpp:322
tkstring & prepend(const tkstring &_st)
Definition tkstring.h:201
bool begins_with(const char *_s, ECaseCompare _cmp=kExact) const
Definition tkstring.h:162
units_keys
units identifiers
Definition tkunit.h:46
static std::vector< std::tuple< tkstring, units_type, double > > funits_properties
Definition tkunit.h:52
Definition tklog.cpp:16
tklog & info(tklog &log)
Definition tklog.h:313
tklog & error(tklog &log)
Definition tklog.h:344
tklog & do_endl(tklog &log)
Definition tklog.h:212
tklog & warning(tklog &log)
Definition tklog.h:331