TkN 2.3
Toolkit for Nuclei
Loading...
Searching...
No Matches
tklevel_scheme.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 * This software is governed by the CeCILL-B license under French law and *
11 * abiding by the rules of distribution of free software. You can use, *
12 * modify and/ or redistribute the software under the terms of the *
13 * CeCILL-B license as circulated by CEA, CNRS and INRIA at the following *
14 * URL \"http://www.cecill.info\". *
15 * *
16 * As a counterpart to the access to the source code and rights to copy, *
17 * modify and redistribute granted by the license, users are provided *
18 * only with a limited warranty and the software's author, the holder of *
19 * the economic rights, and the successive licensors have only limited *
20 * liability. *
21 * *
22 * In this respect, the user's attention is drawn to the risks associated *
23 * with loading, using, modifying and/or developing or reproducing the *
24 * software by the user in light of its specific status of free software, *
25 * that may mean that it is complicated to manipulate, and that also *
26 * therefore means that it is reserved for developers and experienced *
27 * professionals having in-depth computer knowledge. Users are therefore *
28 * encouraged to load and test the software's suitability as regards *
29 * their requirements in conditions enabling the security of their *
30 * systems and/or data to be ensured and, more generally, to use and *
31 * operate it in the same conditions as regards security. *
32 * *
33 * The fact that you are presently reading this means that you have had *
34 * knowledge of the CeCILL-B license and that you accept its terms. *
35 ********************************************************************************/
36
37#ifdef HAS_DEBUG
38#include <algorithm>
39#include <chrono>
40#endif
41
42#include "tkdb_table.h"
43#include "tklevel_scheme.h"
44#include "tkstring.h"
45#include "tklog.h"
46#include "tkdatabase.h"
47#include "tkmanager.h"
48
49namespace tkn {
67}
68
69using namespace tkn;
70
71tklevel_scheme::tklevel_scheme(const tkstring &_nuc, int _zz, int _aa):fnucleus(_nuc), fA(_aa), fZ(_zz)
72{
73 if(fA<0||fZ<0) return;
74 init();
75}
76
94std::vector<shared_ptr<tklevel>> tklevel_scheme::get_levels(const std::function<bool (shared_ptr<tklevel>)> &_selection)
95{
96 vector<shared_ptr<tklevel>> res;
97 const auto &levels = get_levels();
98 std::copy_if(levels.begin(), levels.end(), std::back_inserter(res),
99 [&_selection](const auto &lvl) { return _selection(lvl); });
100 return res;
101}
102
117std::vector<shared_ptr<tkdecay>> tklevel_scheme::get_decays(const std::function<bool (shared_ptr<tkdecay>)> &_selection)
119 vector<shared_ptr<tkdecay>> res;
120 const auto &decays = get_decays();
121 std::copy_if(decays.begin(), decays.end(), std::back_inserter(res),
122 [&_selection](const auto &dec) { return _selection(dec); });
123 return res;
124}
125
126
141| | |
142|---|---|
143| tab | print data as sorted in a tabular view |
144| yrast | only confirmed yrast levels will be printed |
145| yrastt | only confirmed and tentative yrast levels will be printed |
146| com | print comments |
147*/
148void tklevel_scheme::print(const tkstring &_data, const tkstring &_option)
149{
150 if(_data.is_empty()) {
151 glog << info << " Level scheme of " << fnucleus << do_endl;
152 glog << info << " Options : 'dataset', 'level', 'decay'" << do_endl;
153 }
154 if(_data.contains("dataset")){
155 glog << info << "Available datasets are : " << do_endl;
156 for(auto &dts : fmap_of_dataset) glog << info << dts.second->get_name() << " (" << dts.first << ")" << do_endl;
157 glog << comment << "Current dataset is '" << get_dataset()->get_name() << "'" << " (" << get_dataset()->get_id() << ")" << do_endl;
158 }
159 if(_data.contains("level")||_data.contains("decay")){
160 fmap_of_dataset[fdatasetid]->print(_data, _option);
161 }
162 if(_data == "*") {
163 int _current_data_set = fdatasetid;
164 for(auto &dataset: fmap_of_dataset) {
165 glog << info << "dataset : " << dataset.second->get_name() << do_endl;
166 select_dataset(dataset.first);
167 dataset.second->print("leveldecay");
168 }
169 // to select back the current dataset
170 select_dataset(_current_data_set);
171 }
172}
173
178bool tklevel_scheme::select_dataset(const tkstring &_dataset_name)
179{
180 if(fmap_of_dataset_name_id.count(_dataset_name)==0) {
181 glog << warning << "dataset '" << _dataset_name << "' do not exist for " << fnucleus << do_endl;
182 print("dataset");
183 return false;
184 }
185 fdatasetid = fmap_of_dataset_name_id[_dataset_name];
186 return select_dataset(fdatasetid);
187}
188
193bool tklevel_scheme::select_dataset(int _dataset_id)
194{
195 if(fmap_of_dataset.count(_dataset_id)==0) {
196 glog << warning << "dataset with id '" << _dataset_id << "' do not exist for " << fnucleus << do_endl;
197 print("dataset");
198 return false;
199 }
200
201 fdatasetid = _dataset_id;
202 fmap_of_dataset[fdatasetid]->load_dataset();
203 return true;
204}
205
218shared_ptr<tklevel> tklevel_scheme::get_level(const tkstring &_name, bool _exact)
219{
220 // get a specific level [spin][parity][rank] (ex: 0+1 for the first 0+ state)
221 string spin;
222 int rank;
223
224 const char* sep = "";
225 if(_name.contains("+")) {sep = "+";}
226 else if(_name.contains("-")) {sep = "-";}
227
228 auto array = _name.tokenize(sep);
229 if(array.size()!=2) {
230 glog << error << " <nuclear_level_scheme::get_level_energy()> invalid argument '" << _name << "'. Syntax: '[spin][parity][rank]'." << do_endl;
231 return nullptr;
232 }
233
234 spin = array.at(0);
235 rank = array.at(1).atoi();
236
237 int idx = 0;
238 for(const auto &lvlp : get_levels([](auto lvl) {return (lvl->get_spin_parity()!=nullptr);})) {
239 auto jpistr = lvlp->get_spin_parity()->get_jpi_str();
240 bool ok=false;
241 if(!_exact) {
242 if(jpistr.equal_to(tkstring::form("%s%s",spin.data(),sep)) ||
243 jpistr.equal_to(tkstring::form("(%s%s)",spin.data(),sep)) ||
244 jpistr.equal_to(tkstring::form("(%s)%s",spin.data(),sep)) ||
245 jpistr.equal_to(tkstring::form("%s(%s)",spin.data(),sep))) {
246 ok = true;
247 }
248 }
249 else if(jpistr.equal_to(tkstring::form("%s%s",spin.data(),sep))) {
250 ok = true;
251 }
252 if(ok) {
253 idx++;
254 if(idx==rank) return lvlp;
255 }
256 }
257 return nullptr;
258}
259
265 * if offset is set to *, all offset values will be considered.
266 * if offset contains "dec", only levels that have known decays will be considered
267 */
268shared_ptr<tklevel> tklevel_scheme::get_level(double _energy, tkstring _offset)
269{
270 shared_ptr<tklevel> clostest_lev;
271 double best_ediff=1e6;
272 for(const auto &lev: get_levels()) {
273 if(!_offset.is_empty() && !_offset.contains("*")) {
274 if(_offset != lev->get_offset_bandhead()) continue;
275 }
276 if(_offset.is_empty() && !lev->get_offset_bandhead().is_empty()) continue;
277 if(_offset.contains("dec") && lev->get_decays_down().size()==0) continue;
278 double ediff = abs(lev->get_energy(tkunit_manager::units_keys::keV,true)-_energy);
279 if(ediff<best_ediff) {
280 best_ediff=ediff;
281 clostest_lev = lev;
282 }
283 }
284 return clostest_lev;
285}
286
287void tklevel_scheme::init()
288{
289 gdatabase->begin("DISTINCT dataset.dataset_id, dataset.dataset_name","level INNER JOIN dataset on level.dataset_id=dataset.dataset_id INNER JOIN isotope on level.isotope_id=isotope.isotope_id INNER JOIN element on isotope.element_id=element.element_id",tkstring::form("element.charge=%d AND isotope.mass=%d",fZ,fA), tkstring::form("read_dataset_%s",fnucleus.data()));
290
291 while(gdatabase->next(tkstring::form("read_dataset_%s",fnucleus.data()))) {
292 tkdb_table& datasets = (*gdatabase)["DATASET"];
293 fmap_of_dataset_name_id[datasets["dataset_name"].get_value()] = datasets["dataset_id"].get_value().atoi();
294 fmap_of_dataset[datasets["dataset_id"].get_value().atoi()] = make_shared<tkdataset>(tkdataset(fnucleus,fZ,fA,datasets["dataset_name"].get_value(),datasets["dataset_id"].get_value().atoi()));
295 }
296 tkstring default_dataset = fnucleus.data();
297 default_dataset += " : ADOPTED LEVELS, GAMMAS";
298 gdebug << "looking for dataset : " << default_dataset << do_endl;
299 bool found = (fmap_of_dataset_name_id.count(default_dataset)>0);
300 gdebug << (found?"":"not ") << "found " << (found?fmap_of_dataset_name_id[default_dataset.data()]:-1) << do_endl;
301
302 if(found) select_dataset(fmap_of_dataset_name_id[default_dataset.data()]);
303 else{
304 for(const auto &dts : fmap_of_dataset) {
305 glog << warning << tkstring::form("%5s",fnucleus.data()) << " -- no 'ADOPTED LEVELS' dataset : '" << dts.second->get_name() << "' used instead." << do_endl;
306 select_dataset(dts.first);
307 break;
308 }
309 }
310 if(fmap_of_dataset.empty()) {
311 tkstring datasetname = "EMPTY";
312 fdatasetid = 0;
313 fmap_of_dataset_name_id[datasetname] = fdatasetid;
314 fmap_of_dataset[fdatasetid] = make_shared<tkdataset>(tkdataset(fnucleus,fZ,fA,datasetname,fdatasetid));
316
317 floaded = true;
318}
319
320void tkdataset::check_yrast(shared_ptr<tklevel> _lvl)
321{
322 // check Yrast
323 if(_lvl->get_spin_parity()->get_spin().is_known()) {
324 int TwoJ = (int)(_lvl->get_spin_parity()->get_spin().get_value()*2+0.5);
325 if(!fyrastmap_exact.count(TwoJ)) {
326 fyrastmap_exact[TwoJ] = _lvl;
327 _lvl->fisYrast_exact = true;
328 // cout<<"exact yrast: 2J:"<<TwoJ<<" "<<lvl->get_spin_parity_str()<<" "<<lvl->get_energy()<<endl;
329 }
330 if(!fyrastmap_uncertain.count(TwoJ)) {
331 fyrastmap_uncertain[TwoJ] = _lvl;
332 _lvl->fisYrast_uncertain = true;
333 // cout<<"uncertain yrast: 2J:"<<TwoJ<<" "<<lvl->get_spin_parity_str()<<" "<<lvl->get_energy()<<endl;
334 }
335 }
336 else if(!_lvl->get_spin_parity_str().is_empty() && !_lvl->get_spin_parity_str().contains(",")) {
337 tkstring tmp = _lvl->get_spin_parity_str().copy().remove_all("(").remove_all(")").remove_all("+").remove_all("-");
338 int TwoJ;
339 if(tmp.contains("/2")) TwoJ = tmp.remove_all("/2").atoi();
340 else TwoJ = tmp.atoi()*2;
341 if(!fyrastmap_uncertain.count(TwoJ)) {
342 fyrastmap_uncertain[TwoJ] = _lvl;
343 _lvl->fisYrast_uncertain = true;
344 // cout<<"uncertain yrast: 2J:"<<TwoJ<<" "<<lvl->get_spin_parity_str()<<" "<<lvl->get_energy()<<endl;
345 }
346 }
347}
348
349void tkdataset::load_dataset()
350{
351 if(floaded) return;
352
353#ifdef HAS_DEBUG
354 auto t_0 = std::chrono::high_resolution_clock::now();
355 gdebug << "loading levels for dataset '" << fname << "' of " << fnucleus << do_endl;
356#endif
357
358 gdatabase->begin("DISTINCT level.*","level",tkstring::form("dataset_id=%d",fid));
359
360 while(gdatabase->next()) {
361 tkdb_table *levels = &(*gdatabase)["LEVEL"];
362 tkdb_column *col = &(*levels)["level_id"];
363 int id = col->get_value().atoi();
364
365 gmanager->set_max_level_id(id);
366
367 tkdb_table::measure_data_struct energy_struc;
368 levels->read_measure(energy_struc, "level_energy");
369
370 shared_ptr<tklevel> lvl = make_shared<tklevel>(id,energy_struc);
371 lvl->fbelongs_to_nucleus = fnucleus;
372
373 double spin = -1.;
374 int parity = -1.;
375 tkstring spin_parity_str = "";
376
377 col = &(*levels)["level_spin_parity"];
378 if(!col->is_empty()) spin_parity_str = col->get_value();
379 col = &(*levels)["level_spin"];
380 if(!col->is_empty()) spin = col->get_value().atof();
381 col = &(*levels)["level_parity"];
382 if(!col->is_empty()) parity = col->get_value().atoi();
383
384 if(!spin_parity_str.is_empty()) lvl->set_jpi(spin,parity,spin_parity_str);
385
386 tkdb_table::measure_data_struct lifetime_struc;
387 levels->read_measure(lifetime_struc, "level_lifetime");
388 if(lifetime_struc.filled) lvl->set_lifetime(lifetime_struc);
389
390 col = &(*levels)["level_comment"];
391 if(!col->is_empty()) lvl->set_comment(col->get_value());
392 col = &(*levels)["level_uncertain"];
393 if(!col->is_empty()) lvl->set_uncertain_level(col->get_value());
394
395 flevels.push_back(lvl);
396
397 check_yrast(lvl);
398 }
399 for(const auto &lvl: flevels) {
400 fmapoflevels[lvl->get_id()] = lvl;
401 }
402
403#ifdef HAS_DEBUG
404 auto t_1 = std::chrono::high_resolution_clock::now();
405 gdebug << "loading decays for dataset '" << fname << "' of " << fnucleus << do_endl;
406#endif
407
408 gdatabase->begin("DISTINCT decay.*","decay INNER JOIN level on decay.level_from_id=level.level_id ",tkstring::form("dataset_id=%d",fid));
409
410 while(gdatabase->next()) {
411 int dectype=-1;
412 int level_from = -1.;
413 int level_to = -1.;
414
415 tkdb_table *decay = &(*gdatabase)["DECAY"];
416 tkdb_column *col = &(*decay)["decay_id"];
417 int id = col->get_value().atoi();
418
419 gmanager->set_max_decay_id(id);
420
421 col = &(*decay)["decay_type"];
422 if(!col->is_empty()) dectype = col->get_value().atoi();
423 col = &(*decay)["level_from_id"];
424 if(!col->is_empty()) level_from = col->get_value().atoi();
425 col = &(*decay)["level_to_id"];
426 if(!col->is_empty()) level_to = col->get_value().atoi();
427
428 shared_ptr<tkdecay> dec = nullptr;
429
430 tkdb_table::measure_data_struct energy_struc;
431 decay->read_measure(energy_struc, "decay_energy");
432
433 if(dectype == tkdecay::kgamma) {
434 dec = make_shared<tkgammadecay>(id, level_from, level_to, energy_struc);
435 shared_ptr<tkgammadecay> gamma = dynamic_pointer_cast<tkgammadecay>(dec);
436
437 tkdb_table::measure_data_struct relative_intensity_struc;
438 decay->read_measure(relative_intensity_struc, "decay_relative_intensity");
439 if(relative_intensity_struc.filled) gamma->set_relative_intensity(relative_intensity_struc);
440
441 col = &(*decay)["decay_multipolarity"];
442 if(!col->is_empty()) gamma->set_multipolarity(col->get_value());
443
444 tkdb_table::measure_data_struct mixing_ratio_struc;
445 decay->read_measure(mixing_ratio_struc, "decay_mixing_ratio");
446 if(mixing_ratio_struc.filled) gamma->set_mixing_ratio(mixing_ratio_struc);
447
448 tkdb_table::measure_data_struct conv_coeff_struc;
449 decay->read_measure(conv_coeff_struc, "decay_conv_coeff");
450 if(conv_coeff_struc.filled) gamma->set_conv_coeff(conv_coeff_struc);
451
452 // transition probabilities
453 tkdb_table::measure_data_struct transprob_struc;
454 decay->read_measure(transprob_struc, "decay_BE");
455 if(transprob_struc.filled) gamma->set_transition_probability(transprob_struc);
456 transprob_struc.clear();
457 decay->read_measure(transprob_struc, "decay_BEW");
458 if(transprob_struc.filled) gamma->set_transition_probability(transprob_struc);
459 transprob_struc.clear();
460 decay->read_measure(transprob_struc, "decay_BM");
461 if(transprob_struc.filled) gamma->set_transition_probability(transprob_struc);
462 transprob_struc.clear();
463 decay->read_measure(transprob_struc, "decay_BMW");
464 if(transprob_struc.filled) gamma->set_transition_probability(transprob_struc);
465 }
466 else continue;
467
468 col = &(*decay)["decay_comment"];
469 if(!col->is_empty()) dec->set_comment(col->get_value());
470 col = &(*decay)["decay_uncertain"];
471 if(!col->is_empty()) dec->set_uncertain_decay(col->get_value());
472
473 fdecays.push_back(dec);
474 fdecays.back()->set_levels(fmapoflevels.at(level_from),fmapoflevels.at(level_to));
475 }
476 for(auto &dec: fdecays) {
477 fmapoflevels.at(dec->get_level_from_id())->add_decay_down(dec);
478 fmapoflevels.at(dec->get_level_to_id())->add_decay_up(dec);
479 }
480
481
482#ifdef HAS_DEBUG
483 auto t_2 = std::chrono::high_resolution_clock::now();
484
485 auto dt01 = std::chrono::duration<double, std::milli>(t_1-t_0).count() * 0.001; // open db time
486 auto dt12 = std::chrono::duration<double, std::milli>(t_2-t_1).count() * 0.001; // first read loop with views
487
488 gdebug << std::fixed << std::setprecision(4) << "time level: "
489 << dt01
490 << " s" << do_endl;
491
492 gdebug << std::fixed << std::setprecision(4) << "time decays: "
493 <<dt12
494 << " s" << do_endl;
495#endif
496
497 floaded = true;
498}
499
504void tkdataset::print(const tkstring &_data, const tkstring &_option)
505{
506 if(_data.contains("level") && _data.contains("decay")) {
507 glog << info << "dataset '" << fname << "' contains " << flevels.size() << " levels" << do_endl;
508 glog << info << "dataset '" << fname << "' contains " << fdecays.size() << " decays" << do_endl;
509 for(auto &lev : flevels) {
510 if(_option.contains("yrastt") && !lev->is_yrast(true)) continue;
511 else if(_option.contains("yrast") && !_option.contains("yrastt") && !lev->is_yrast()) continue;
512 lev->print(_option);
513 for(auto &dec : lev->get_decays_down()) {
514 cout<<setw(13)<<""<<"-> ";dec->print("quiet;levto");
515 }
516 }
517 }
518 else if(_data.contains("level")) {
519 if(_option.contains("tab")){
520 std::cout << std::setw(15) << left << "symbol" << std::setw(15) << "levelid"
521 << std::setw(15) << left << "energy (keV)"
522 << std::setw(15) << left << "energy error"
523 << std::setw(15) << left << "spin parity"
524 << std::setw(15) << left << "t1/2"
525 << std::setw(15) << left << "t1/2 error"
526 << std::setw(15) << left << "t1/2 unit" << std::endl;}
527 else glog << info << "dataset '" << fname << "' contains " << flevels.size() << " levels:" << do_endl;
528 for(auto &lev : flevels) {
529 if(_option.contains("yrastt") && !lev->is_yrast(true)) continue;
530 else if(_option.contains("yrast") && !_option.contains("yrastt") && !lev->is_yrast()) continue;
531 if(_option.contains("tab")) std::cout << std::setw(15)<< left << fnucleus;
532 lev->print(_option);
533 }
534 }
535 else if(_data.contains("decay")) {
536 if(_option.contains("tab")) std::cout << std::setw(15) << left << "symbol" << std::setw(15) << left << "decayid" << std::setw(15) << left << "energy (keV)" << std::setw(15) << left << "energy error (keV)" << std::endl;
537 else glog << info << "dataset '" << fname << "' contains " << fdecays.size() << " decays:" << do_endl;
538 for(auto &dec : fdecays) {
539 if(_option.contains("tab")) std::cout << std::setw(15) << left << fnucleus;
540 dec->print(_option);
541 }
542 }
543}
544
545shared_ptr<tklevel> tkdataset::add_level(double _ener, double _unc, tkstring _unit, tkstring _jpistr)
546{
547 tkspin spin;
548 tkparity parity;
549 double spin_val = -1.;
550 int parity_val = -1.;
551
552 tkdb_table::measure_data_struct lvl_properties{_ener,_unc,0,0,_unit};
553 shared_ptr<tklevel> lvl = make_shared<tklevel>(gmanager->get_new_level_id(),lvl_properties);
554
555 spin.set(_jpistr);parity.set(_jpistr);
556 if(spin.is_known()) spin_val = spin.get_value();
557 if(parity.is_known()) parity_val = (parity.is_parity(tkparity::eparity::kParityMinus));
558
559 lvl->set_jpi(spin_val,parity_val,_jpistr);
560
561 // needs to generate a unique id for the level
562 lvl->fbelongs_to_nucleus = fnucleus;
563
564 // find the good position in the vector to be sorted by energy
565 auto it = std::lower_bound(flevels.begin(), flevels.end(), lvl,
566 [](const std::shared_ptr<tklevel>& a, const std::shared_ptr<tklevel>& b) {
567 return a->get_energy() < b->get_energy();
568 });
569 flevels.insert(it, lvl);
570
571 check_yrast(lvl);
572 fmapoflevels[lvl->get_id()] = lvl;
573
574 return lvl;
575}
576
577shared_ptr<tkgammadecay> tkdataset::add_gamma_decay(shared_ptr<tklevel> _lvlfrom, shared_ptr<tklevel> _lvlto, double _ener, double _unc)
578{
579 if(_ener==0.) _ener = _lvlfrom->get_energy() - _lvlto->get_energy();
580
581 tkdb_table::measure_data_struct energy_struc{_ener,_unc,0.,0.,"keV"};
582
583 shared_ptr<tkgammadecay> decay = make_shared<tkgammadecay>(gmanager->get_new_decay_id(), _lvlfrom->get_id(), _lvlto->get_id(), energy_struc);
584 decay->set_levels(_lvlfrom,_lvlto);
585 fdecays.push_back(decay);
586 _lvlfrom->add_decay_down(decay);
587 _lvlto->add_decay_up(decay);
588
589 return decay;
590}
591
593#ifdef HAS_ROOT
594ClassImp(tklevel_scheme);
595#endif
Stores information on a specific dataset.
shared_ptr< tklevel > add_level(double _ener, double _unc, tkstring _unit, tkstring _jpistr)
manually add a new level to the dataset
shared_ptr< tkgammadecay > add_gamma_decay(shared_ptr< tklevel > _lvlfrom, shared_ptr< tklevel > _lvlto, double _ener=0., double _unc=0.)
add a new decay between two levels to the dataset
friend class tklevel_scheme
tkstring get_value() const
Definition tkdb_column.h:61
Representaiton of a sqlite data table.
Definition tkdb_table.h:52
void read_measure(measure_data_struct &_struct, const tkstring &_col_base_name)
tklevel_scheme(const tkstring &_nuc, int _zz, int _aa)
void print(const tkstring &_data="", const tkstring &_option="")
print the level scheme information
const std::vector< shared_ptr< tklevel > > & get_levels()
get the vector containing all the levels
bool select_dataset(const tkstring &_dataset_name)
select a dataset using its name
const shared_ptr< tkdataset > & get_dataset()
returns the current dataset
shared_ptr< tklevel > get_level(const tkstring &_name, bool _exact=true)
get the level corresponding to the given name
const std::vector< shared_ptr< tkdecay > > & get_decays()
get the vector containing all the decays
virtual bool is_parity(tkparity::eparity _parity)
test the parity value
virtual void set(const tkstring &_st)
define the parity from a string
bool is_known() const
to get some information about this data
Definition tkproperty.h:91
double get_value()
To get the spin as a double.
void set(int n, int d=1)
define the spin value
std::string with usefull tricks from TString (ROOT) and KVString (KaliVeda) and more....
Definition tkstring.h:54
static const char * form(const char *_format,...)
Definition tkstring.cpp:431
std::vector< tkstring > tokenize(const tkstring &_delim=" ") const
Create a vector of string separated by at least one delimiter.
Definition tkstring.cpp:254
bool is_empty() const
Definition tkstring.h:163
tkstring & remove_all(const tkstring &_s1)
Definition tkstring.h:215
int atoi() const
Converts a string to integer value.
Definition tkstring.cpp:175
bool contains(const char *_pat, ECaseCompare _cmp=kExact) const
Definition tkstring.h:197
double atof() const
Converts a string to double value.
Definition tkstring.cpp:191
Definition tklog.cpp:39
tklog & info(tklog &log)
Definition tklog.h:336
tklog & error(tklog &log)
Definition tklog.h:367
tklog & do_endl(tklog &log)
Definition tklog.h:235
tklog & warning(tklog &log)
Definition tklog.h:354
tklog & comment(tklog &log)
Definition tklog.h:342