TkN 2.4
Toolkit for Nuclei
Loading...
Searching...
No Matches
tkmanager.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 <algorithm>
15#include <fstream>
16#include <iterator>
17
18#include "tkmanager.h"
19#include "tklog.h"
20#include "tkdatabase.h"
21#include "json.hpp"
22
23namespace tkn {
31
32}
33
34using namespace tkn;
35using json = nlohmann::json;
36
37tkmanager* tkmanager::g_data_manager = nullptr;
38
40{
41 if(g_data_manager) return g_data_manager;
42
43 std::lock_guard<std::recursive_mutex> lock(get_mutex());
44
45 if(g_data_manager == nullptr) {
46 g_data_manager = new tkmanager();
47 }
48
49 return g_data_manager;
50}
51
53{
54 std::lock_guard<std::recursive_mutex> lock(get_mutex());
55
56 preload_nuclei();
57
58 g_data_manager = this;
59}
60
62{
63 g_data_manager = nullptr;
64}
65
66shared_ptr<tklevel_scheme> tkmanager::get_level_scheme(const tkn::tkstring &_nuc, int _zz, int _aa)
67{
68 if(fmap_of_level_scheme.count(_nuc)) return fmap_of_level_scheme[_nuc];
69
70 std::lock_guard<std::recursive_mutex> lock(get_mutex());
71
72 auto [it, inserted] = fmap_of_level_scheme.try_emplace(
73 _nuc, make_shared<tklevel_scheme>(_nuc,_zz,_aa));
74 return it->second;
75}
76
77void tkmanager::preload_nuclei()
78{
79 // try with isotope.*, element.*
80 gdatabase->begin("isotope.*, element.*","isotope INNER JOIN element on isotope.element_id=element.element_id");
81
82 while(gdatabase->next()) {
83 tkdb_table& element = (*gdatabase)["ELEMENT"];
84 tkdb_table& iso_table = (*gdatabase)["ISOTOPE"];
85
86 shared_ptr<tknucleus> anucleus = make_shared<tknucleus>();
87
88 tkstring name = element["symbol"].get_value();
89 anucleus->felement_symbol = name;
90 anucleus->felement_name = element["name"].get_value();
91
92 int zz = element["charge"].get_value().atoi();
93 int aa = iso_table["mass"].get_value().atoi();
94
95 gdebug << aa << name << " loaded" << do_endl;
97 name.prepend(iso_table["mass"].get_value());
98
99 anucleus->set_name(name.data());
100 anucleus->fA = aa;
101 anucleus->fZ = zz;
102 anucleus->fN = aa - zz;
103 anucleus->fis_known = true;
104
105 // load element properties
106 for (auto &col : element.get_columns()) {
107 tkstring colname = col.first;
108 if(colname=="element_id") continue;
109 if(colname.ends_with("_unit")) continue;
110
111 tkstring colname_unit = colname + "_unit";
112 tkstring unit("");
113 tkstring value = col.second.get_value();
114 if(element.get_columns().count(colname_unit)) {
115 unit = element.get_columns().at(colname_unit).get_value();
116 }
117 if( colname == "atomic_mass" ||
118 colname == "atomic_radius_van_der_Waals" ||
119 colname == "boiling_point" ||
120 colname == "density" ||
121 colname == "ionization_energy" ||
122 colname == "melting_point"
123 ) {
124 shared_ptr<tkmeasure> data = make_shared<tkmeasure>(value.atof(),unit);
125 data->set_info(tkproperty::kKnown);
126 data->set_type(colname);
127 anucleus->add_property(colname,data);
128 }
129 else if(colname.begins_with("XRay")) {
130 unit = "keV";
131 shared_ptr<tkmeasure> data = make_shared<tkmeasure>(value.atof(),unit);
132 data->set_info(tkproperty::kKnown);
133 data->set_type(colname);
134 if(value.atof()>0) anucleus->add_property(colname,data);
135 }
136 else {
137 anucleus->add_property_str(colname,value,unit);
138 }
139 }
140
141 // load isotope properties
142 vector<tkstring> fNames = {
143 "lifetime", // From ENSDF
144 "abundance", "quadrupoleDeformation", "FY235U", "FY238U", "FY239Pu", "FY241Pu", "FY252Cf", "cFY235U", "cFY238U", "cFY239Pu", "cFY241Pu", "cFY252Cf", // From NUDAT
145 "mass_excess", "radius", "magnetic_dipole", "electric_quadrupole" // From IAEA-NDS API
146 };
147
148 for(auto &property_name: fNames) {
150 iso_table.read_measure(data_struct,property_name);
151 if(data_struct.filled) {
152 shared_ptr<tkmeasure> data = make_shared<tkmeasure>(data_struct.value,data_struct.unit,data_struct.err);
153 if(data_struct.info.contains("ERR=SY")) data->set_info(tkproperty::kSystematic);
154 else data->set_info(tkproperty::kKnown);
155 data->set_type(property_name);
156
157 if(property_name == "lifetime") {
158 if(!data_struct.info.is_empty()) data->set_info_tag(data_struct.info);
159 if(data->get_info_tag().contains("STABLE")) {
160 anucleus->fis_stable = true;
161 data->set_error(0);
162 anucleus->add_property(property_name,data,"STABLE","");
163 }
164 else anucleus->add_property(property_name,data, tkstring::form("%g",data->get_value()));
165 }
166 else {
167 anucleus->add_property(property_name,data);
168 }
169 }
170 }
171
172 anucleus->add_property_str("isotope_year_discovered",iso_table["isotope_year_discovered"].get_value(),"");
173
174 tkstring decay_modes = iso_table["decay_modes"].get_value();
175 anucleus->add_property_str("decay_modes",decay_modes,"");
176
177 tkstring spin_parity_str = iso_table["spin_parity"].get_value();
178 double spin = iso_table["spin"].get_value().atof();
179 int parity = iso_table["parity"].get_value().atoi();
180
181 anucleus->fspin_parity.set_spin(spin);
182 anucleus->fspin_parity.set_parity(parity);
183 anucleus->fspin_parity.set_from_str(spin_parity_str);
184
185 anucleus->add_property_str("spin_parity",spin_parity_str,"");
186
187 fmap_of_nuclei[name] = anucleus;
188 fmap_of_nuclei_per_z[anucleus->fZ].push_back(anucleus);
189 fmap_of_nuclei_per_a[anucleus->fA].push_back(anucleus);
190 fmap_of_nuclei_per_n[anucleus->fN].push_back(anucleus);
191 fmap_of_nuclei_per_z_and_a[anucleus->fZ][anucleus->fA] = anucleus;
192 fvector_of_nuclei.push_back(anucleus);
193
194 if(!fmap_of_symbols.count(anucleus->get_z()))
195 fmap_of_symbols[anucleus->get_z()] = anucleus->get_element_symbol();
196 }
197
198 // separation energy calculation
199 if(!get_nucleus(0,1)||!get_nucleus(1,1)||!get_nucleus(2,4)) return;
200 auto neutron_ME = fmap_of_nuclei_per_z_and_a.at(0).at(1)->get("mass_excess");
201 auto proton_ME = fmap_of_nuclei_per_z_and_a.at(1).at(1)->get("mass_excess");
202 auto alpha_ME = fmap_of_nuclei_per_z_and_a.at(2).at(4)->get("mass_excess");
203 tkmeasure electron_mass(0.51099895000,"MeV",0.00000000015);
204 for(auto &nuc: fmap_of_nuclei) {
205 auto ME = nuc.second->get("mass_excess");
206 if(!ME) continue;
207 // Sn calculation: Sn(A,Z) = − ME(A,Z) + ME(A−1,Z) + ME(n)
208 auto daughter = get_nucleus(nuc.second->get_z(),nuc.second->get_a()-1);
209 if(daughter && daughter->get("mass_excess")) {
210 shared_ptr<tkmeasure> Sn = make_shared<tkmeasure>(-*ME + *daughter->get("mass_excess") + *neutron_ME);
211 Sn->set_type("neutronSeparationEnergy");
212 nuc.second->add_property("neutronSeparationEnergy",Sn);
213 }
214 // S2n calculation: S2n(A,Z) = − ME(A,Z) + ME(A−2,Z) + 2*ME(n)
215 daughter = get_nucleus(nuc.second->get_z(),nuc.second->get_a()-2);
216 if(daughter && daughter->get("mass_excess")) {
217 shared_ptr<tkmeasure> S2n = make_shared<tkmeasure>(-*ME + *daughter->get("mass_excess") + 2.*(*neutron_ME));
218 S2n->set_type("twoNeutronSeparationEnergy");
219 nuc.second->add_property("twoNeutronSeparationEnergy",S2n);
220 }
221 // Sp calculation: Sp(A,Z) = − ME(A,Z) + ME(A−1,Z-1) + ME(p)
222 daughter = get_nucleus(nuc.second->get_z()-1,nuc.second->get_a()-1);
223 if(daughter && daughter->get("mass_excess")) {
224 shared_ptr<tkmeasure> Sp = make_shared<tkmeasure>(-*ME + *daughter->get("mass_excess") + *proton_ME);
225 Sp->set_type("protonSeparationEnergy");
226 nuc.second->add_property("protonSeparationEnergy",Sp);
227 }
228 // S2p calculation: S2p(A,Z) = − ME(A,Z) + ME(A−2,Z-2) + 2*ME(p)
229 daughter = get_nucleus(nuc.second->get_z()-2,nuc.second->get_a()-2);
230 if(daughter && daughter->get("mass_excess")) {
231 shared_ptr<tkmeasure> S2p = make_shared<tkmeasure>(-*ME + *daughter->get("mass_excess") + 2.*(*proton_ME));
232 S2p->set_type("twoProtonSeparationEnergy");
233 nuc.second->add_property("twoProtonSeparationEnergy",S2p);
234 }
235
236 // Qalpha
237 daughter = get_nucleus(nuc.second->get_z()-2,nuc.second->get_a()-4);
238 if(daughter && daughter->get("mass_excess")) {
239 shared_ptr<tkmeasure> Q = make_shared<tkmeasure>( *ME - *daughter->get("mass_excess") - *alpha_ME);
240 Q->set_type("Qalpha");
241 nuc.second->add_property("Qalpha",Q);
242 }
243 // QbetaMinus
244 daughter = get_nucleus(nuc.second->get_z()+1,nuc.second->get_a());
245 if(daughter && daughter->get("mass_excess")) {
246 shared_ptr<tkmeasure> Q = make_shared<tkmeasure>( *ME - *daughter->get("mass_excess"));
247 Q->set_type("QbetaMinus");
248 nuc.second->add_property("QbetaMinus",Q);
249 }
250 // QbetaMinusOneNeutronEmission
251 daughter = get_nucleus(nuc.second->get_z()+1,nuc.second->get_a()-1);
252 if(daughter && daughter->get("mass_excess")) {
253 shared_ptr<tkmeasure> Q = make_shared<tkmeasure>( *ME - *daughter->get("mass_excess") - *neutron_ME);
254 Q->set_type("QbetaMinusOneNeutronEmission");
255 nuc.second->add_property("QbetaMinusOneNeutronEmission",Q);
256 }
257 // QbetaMinusTwoNeutronEmission
258 daughter = get_nucleus(nuc.second->get_z()+1,nuc.second->get_a()-2);
259 if(daughter && daughter->get("mass_excess")) {
260 shared_ptr<tkmeasure> Q = make_shared<tkmeasure>( *ME - *daughter->get("mass_excess") - 2.*(*neutron_ME) );
261 Q->set_type("QbetaMinusTwoNeutronEmission");
262 nuc.second->add_property("QbetaMinusTwoNeutronEmission",Q);
263 }
264 // QdeltaAlpha: 0.5*(Qa(Z+2,N+2)-Qa(Z,N))
265 daughter = get_nucleus(nuc.second->get_z()+2,nuc.second->get_a()+4);
266 if(daughter && daughter->get("mass_excess") && nuc.second->get("Qalpha")) {
267 auto Qa = nuc.second->get("Qalpha");
268 auto Qa2 = *daughter->get("mass_excess") - *ME - *alpha_ME;
269 shared_ptr<tkmeasure> Q = make_shared<tkmeasure>( 0.5*(Qa2 - *Qa ));
270 Q->set_type("QdeltaAlpha");
271 nuc.second->add_property("QdeltaAlpha",Q);
272 }
273 // QdoubleBetaMinus
274 daughter = get_nucleus(nuc.second->get_z()+2,nuc.second->get_a());
275 if(daughter && daughter->get("mass_excess")) {
276 shared_ptr<tkmeasure> Q = make_shared<tkmeasure>( *ME - *daughter->get("mass_excess"));
277 Q->set_type("QdoubleBetaMinus");
278 nuc.second->add_property("QdoubleBetaMinus",Q);
279 }
280 // QelectronCapture
281 daughter = get_nucleus(nuc.second->get_z()-1,nuc.second->get_a());
282 if(daughter && daughter->get("mass_excess")) {
283 shared_ptr<tkmeasure> Q = make_shared<tkmeasure>( *ME - *daughter->get("mass_excess"));
284 Q->set_type("QelectronCapture");
285 nuc.second->add_property("QelectronCapture",Q);
286 }
287 // QdoubleElectronCapture
288 daughter = get_nucleus(nuc.second->get_z()-2,nuc.second->get_a());
289 if(daughter && daughter->get("mass_excess")) {
290 shared_ptr<tkmeasure> Q = make_shared<tkmeasure>( *ME - *daughter->get("mass_excess"));
291 Q->set_type("QdoubleElectronCapture");
292 nuc.second->add_property("QdoubleElectronCapture",Q);
293 }
294 // QelectronCaptureOneProtonEmission
295 daughter = get_nucleus(nuc.second->get_z()-2,nuc.second->get_a()-1);
296 if(daughter && daughter->get("mass_excess")) {
297 shared_ptr<tkmeasure> Q = make_shared<tkmeasure>( *ME - *daughter->get("mass_excess") - *proton_ME);
298 Q->set_type("QelectronCaptureOneProtonEmission");
299 nuc.second->add_property("QelectronCaptureOneProtonEmission",Q);
300 }
301 // QpositronEmission
302 daughter = get_nucleus(nuc.second->get_z()-1,nuc.second->get_a());
303 if(daughter && daughter->get("mass_excess")) {
304 auto q_positron = make_shared<tkmeasure>( *ME - *daughter->get("mass_excess") - 2.*electron_mass);
305 q_positron->set_type("QpositronEmission");
306 nuc.second->add_property("QpositronEmission",q_positron);
307 }
308
309 // binding energy
310 auto q_binding = make_shared<tkmeasure>( 1./nuc.second->get_a()*(nuc.second->get_z()*(*proton_ME) + nuc.second->get_n()*(*neutron_ME) - *ME));
311 q_binding->set_type("binding_energy_overA");
312 nuc.second->add_property("binding_energy_overA",q_binding);
313
314 // pairingGap PG=0.5 x (-1)N x ( 2BE(Z,N) - BE(Z,N-1) -BE(Z,N+1) )
315 auto nuc_minus = get_nucleus(nuc.second->get_z(),nuc.second->get_a()-1);
316 auto nuc_plus = get_nucleus(nuc.second->get_z(),nuc.second->get_a()+1);
317 if(nuc_minus && nuc_plus && nuc_minus->get("mass_excess") && nuc_plus->get("mass_excess")) {
318 auto q_pairing = make_shared<tkmeasure>( 0.5 * pow(-1,nuc.second->get_n()) * ( -2* (*ME) + *nuc_minus->get("mass_excess") + *nuc_plus->get("mass_excess")));
319 q_pairing->set_type("pairingGap");
320 nuc.second->add_property("pairingGap",q_pairing);
321 }
322 }
323}
324
325void tkmanager::preload_level_schemes(bool _verbose)
326{
327 size_t nnuc=fvector_of_nuclei.size(),inuc=0;
328 for(const auto &nuc : fvector_of_nuclei){
329 if(_verbose) glog.progress_bar(nnuc,inuc,"... pre-loading");
330 get_level_scheme(nuc->get_symbol(),nuc->get_z(),nuc->get_a());
331 inuc++;
332 }
333}
334
350vector<shared_ptr<tknucleus> > tkmanager::get_nuclei(const std::function<bool (shared_ptr<tknucleus>)> &_selection)
351{
352 vector<shared_ptr<tknucleus>> res;
353 std::copy_if(fvector_of_nuclei.begin(), fvector_of_nuclei.end(), std::back_inserter(res),
354 [&_selection](const auto &nuc) { return _selection(nuc); });
355 return res;
356}
357
358bool tkmanager::known_element(tkstring _nuc, int &_z) {
359 tkstring symbol=_nuc.extract_alpha();
360 for(auto &nuc_test: fmap_of_nuclei) {
361 if(nuc_test.first.copy().extract_alpha() == symbol) {
362 _z = nuc_test.second->get_z();
363 return true;
364 }
365 }
366 return false;
367}
368
370{
371 std::lock_guard<std::recursive_mutex> lock(get_mutex());
372 fmax_level_id++;
373 return fmax_level_id;
374}
375
377{
378 std::lock_guard<std::recursive_mutex> lock(get_mutex());
379 fmax_decay_id++;
380 return fmax_decay_id;
381}
382
394 * - `S1p`, `S1p_sigma`
395 * - `S2p`, `S2p_sigma`
396 *
397 * @param json_path Path to the JSON file (default: "%s/databases/driplines.json").
398 *
399 * @note After loading, all points are sorted by (Z,N) for deterministic order.
400 *
401 * @note Data extracted from the supplementary material of:
402 * L. Neufcourt, Y. Cao, S. A. Giuliani, W. Nazarewicz, E. Olsen, O. B. Tarasov,
403 * "Quantified limits of the nuclear landscape", Phys. Rev. C 101, 044307 (2020).
404 *
405 * @code
406 * // Example usage:
407 * gmanager->load_drip_lines("path/to/nuclei_quantities.json");
408 * auto s2n = gmanager->get_drip_line("S2n");
409 * @endcode
410 */
411void tkmanager::load_drip_lines(const tkstring& json_path)
412{
413 std::ifstream f(json_path);
414 if (!f.is_open()) {
415 glog << error << "tkmanager::load_drip_lines: cannot open " << json_path << do_endl;
416 return;
417 }
418
419 json j;
420 try {
421 f >> j;
422 } catch (const std::exception& e) {
423 glog << error << "tkmanager::load_drip_lines: JSON parse error in "
424 << json_path << " : " << e.what() << do_endl;
425 return;
426 }
427
428 if (!j.is_array()) {
429 glog << error << "tkmanager::load_drip_lines: root is not an array in "
430 << json_path << do_endl;
431 return;
432 }
433
434 fdrip_lines.clear();
435
436 auto push_if_present = [&](const json& o, const char* qname, const char* sname){
437 if (!o.contains(qname)) return;
438 if (o[qname].is_null()) return;
439
440 tkn_drip_point p;
441 p.Z = o.value("Z", 0);
442 p.N = o.value("N", 0);
443 p.value_mev = o.value(qname, std::numeric_limits<double>::quiet_NaN());
444 // sigma may be absent → keep NaN
445 if (o.contains(sname) && !o[sname].is_null()) {
446 p.sigma_mev = o.value(sname, std::numeric_limits<double>::quiet_NaN());
447 }
448 p.type = qname;
449
450 fdrip_lines[p.type].push_back(std::move(p));
451 };
452
453 try {
454 for (const auto& o : j) {
455 if (!o.is_object()) continue;
456 // Z and N must exist; skip otherwise
457 if (!o.contains("Z") || !o.contains("N")) continue;
458
459 push_if_present(o, "S1n", "S1n_sigma");
460 push_if_present(o, "S2n", "S2n_sigma");
461 push_if_present(o, "S1p", "S1p_sigma");
462 push_if_present(o, "S2p", "S2p_sigma");
463 }
464 } catch (const std::exception& e) {
465 glog << error << "tkmanager::load_drip_lines: error while iterating JSON: "
466 << e.what() << do_endl;
467 fdrip_lines.clear();
468 return;
469 }
470
471 // Sort deterministically by Z, then N
472 for (auto& kv : fdrip_lines) {
473 auto& vec = kv.second;
474 std::sort(vec.begin(), vec.end(), [](const tkn_drip_point& a, const tkn_drip_point& b){
475 if (a.Z != b.Z) return a.Z < b.Z;
476 return a.N < b.N;
477 });
478 }
479
480 glog << info << "load_drip_lines: loaded " << j.size()
481 << " nuclei objects from " << json_path << do_endl;
482}
483
484
509std::vector<tkn_drip_point> tkmanager::get_drip_line(const std::string& _type, bool reduce_per_Z)
510{
511 if (fdrip_lines.empty()) {
512 // charge par défaut si pas déjà fait
513 load_drip_lines(tkstring::Form("%s/databases/driplines.json", TKN_SYS));
514 }
515
516 auto it = fdrip_lines.find(_type);
517 if (it == fdrip_lines.end()) {
518 glog << warning << "get_drip_line: unknown type '" << _type << "'" << do_endl;
519 return {};
520 }
521
522 const auto& vec = it->second;
523
524 if (!reduce_per_Z) {
525 return vec;
526 }
527
528 // --- Réduction : garder seulement le point avec min |value| par Z
529 std::unordered_map<int, tkn_drip_point> best_by_Z;
530 for (const auto& p : vec) {
531 auto itZ = best_by_Z.find(p.Z);
532 if (itZ == best_by_Z.end()) {
533 best_by_Z[p.Z] = p;
534 } else {
535 const auto& cur = itZ->second;
536 if (std::isnan(cur.value_mev) ||
537 (!std::isnan(p.value_mev) && std::abs(p.value_mev) < std::abs(cur.value_mev))) {
538 best_by_Z[p.Z] = p;
539 }
540 }
541 }
542
543 std::vector<tkn_drip_point> reduced;
544 reduced.reserve(best_by_Z.size());
545 std::transform(best_by_Z.begin(), best_by_Z.end(), std::back_inserter(reduced),
546 [](const auto &kv) { return kv.second; });
547
548 std::sort(reduced.begin(), reduced.end(), [](const tkn_drip_point& a, const tkn_drip_point& b){
549 if (a.Z != b.Z) return a.Z < b.Z;
550 return a.N < b.N;
551 });
552
553 return reduced;
554}
555
556#ifdef HAS_ROOT
558ClassImp(tkmanager)
559#endif
Representaiton of a sqlite data table.
Definition tkdb_table.h:29
void read_measure(measure_data_struct &_struct, const tkstring &_col_base_name)
row & get_columns()
Definition tkdb_table.h:116
void clear()
Definition tklog.h:156
Manages the database loading and provides access to the physics properties.
Definition tkmanager.h:51
static tkmanager * the_data_manager()
Definition tkmanager.cpp:83
int get_new_level_id()
define a new unique level id
std::vector< tkn_drip_point > get_drip_line(const std::string &_type, bool reduce_per_Z=true)
Get the drip line for a given type (optionally reduced per Z)
shared_ptr< tknucleus > get_nucleus(const tkstring &_nuc)
return a shared pointer to a nucleus from its name
Definition tkmanager.h:125
const vector< shared_ptr< tknucleus > > & get_nuclei()
return a vector containing all the known nuclei
Definition tkmanager.h:106
int get_new_decay_id()
define a new unique decay id
virtual ~tkmanager()
bool known_element(tkstring _nuc, int &_z)
is the element symbol is known (ex: "C")
void preload_level_schemes(bool _verbose=false)
preload all the level schemes from the database
std::string with usefull tricks from TString (ROOT) and KVString (KaliVeda) and more....
Definition tkstring.h:31
tkstring extract_alpha()
Returns a tkstring composed only of the alphabetic letters of the original tkstring.
Definition tkstring.cpp:364
static const char * form(const char *_format,...)
Definition tkstring.cpp:408
static tkstring Form(const char *_format,...)
Definition tkstring.cpp:322
bool is_empty() const
Definition tkstring.h:140
bool ends_with(const char *_s, ECaseCompare _cmp=kExact) const
Definition tkstring.cpp:218
bool contains(const char *_pat, ECaseCompare _cmp=kExact) const
Definition tkstring.h:174
tkstring & prepend(const tkstring &_st)
Definition tkstring.h:201
bool begins_with(const char *_s, ECaseCompare _cmp=kExact) const
Definition tkstring.h:162
double atof() const
Converts a string to double value.
Definition tkstring.cpp:168
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
data structure used to fill a tkmeasure object from the sqlite database
Definition tkdb_table.h:49
int Z
Proton number.
Definition tkmanager.h:41
double value_mev
Quantity value (MeV)
Definition tkmanager.h:43
std::string type
Quantity type: "S1n", "S2n", "S1p", or "S2p".
Definition tkmanager.h:45
int N
Neutron number.
Definition tkmanager.h:42
double sigma_mev
Quantity uncertainty (MeV)
Definition tkmanager.h:44