TkN 2.5
Toolkit for Nuclei
Loading...
Searching...
No Matches
tklevel.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 "tklevel.h"
15#include "tkstring.h"
16#include "tklog.h"
17
18#include <iomanip>
19#include <sstream>
20#include <string>
21
22using namespace std;
23
24namespace tkn {
32}
33
34using namespace tkn;
35
37 flevel_id(_id),
38 fspin_parity(make_shared<tkspin_parity>())
39{
40 shared_ptr<tkmeasure> energy = make_shared<tkmeasure>(_struc.value, _struc.unit);
41 if (_struc.value >= 0.) energy->set_info(tkproperty::kKnown);
42 if (_struc.err > 0.)
43 energy->set_error(_struc.err);
44 else if (_struc.err_low > 0. && _struc.err_high > 0.)
45 energy->set_error(_struc.err_low, _struc.err_high);
46 if (!_struc.info.is_empty()) energy->set_info_tag(_struc.info);
47 energy->set_type("Level energy");
48 add_property("energy", energy);
49
50 if (energy->get_info_tag().contains(";OFF=")) {
51 int idx = energy->get_info_tag().index(";OFF=");
52 tkstring offset = energy->get_info_tag().substr(idx + 5, 1);
53 set_energy_offset(offset);
54 }
55}
56
58{
59 shared_ptr<tkmeasure> lifetime = make_shared<tkmeasure>(_struc.value, _struc.unit);
60 lifetime->set_info(tkproperty::kKnown);
61 if (_struc.err > 0.)
62 lifetime->set_error(_struc.err);
63 else if (_struc.err_low > 0. && _struc.err_high > 0.)
64 lifetime->set_error(_struc.err_low, _struc.err_high);
65 if (!_struc.info.is_empty()) lifetime->set_info_tag(_struc.info);
66 if (lifetime->get_info_tag().contains("STABLE")) fis_stable = true;
67 lifetime->set_type("lifetime");
68 add_property("lifetime", lifetime, tkstring::form("%g", lifetime->get_value()));
69
70 if (lifetime->get_info_tag().contains(";ISO=")) {
71 fis_isomer = true;
72 int idx = lifetime->get_info_tag().index(";ISO=");
73 fisomer_level = lifetime->get_info_tag().substr(idx + 6, 1).atoi(); // ex= ;ISO=M1 => 1
74 if (fisomer_level == 0) fisomer_level = 1;
75 }
76}
77
78void tklevel::set_jpi(double _spin, int _parity, const tkstring &_jpi_str)
79{
80 if (_spin >= 0) fspin_parity->set_spin(_spin);
81 if (_parity >= 0) fspin_parity->set_parity(_parity);
82 fspin_parity->set_from_str(_jpi_str);
83 add_property_str("spin_parity", _jpi_str, "");
84}
85
95double tklevel::get_energy(const tkunit_manager::units_keys _unit, bool _with_offset)
96{
97 if (gunits->get_type(_unit) != tkunit_manager::units_type::kEnergy) {
98 glog << error << gunits->get_name(_unit) << " is not of 'energy' type (returns nan(1))" << do_endl;
99 return std::nan("1");
100 }
101 if (get("energy")->get_unit() == "") {
102 cout << fenergy_offset << endl;
103 cout << fbelongs_to_nucleus << endl;
104 print();
105 }
106 if (is_energy_offset() && !_with_offset) {
107 glog << warning << "Level is only known relatively to an offset, use the _with_offset option to get the relative energy" << do_endl;
108 return std::nan("1");
109 }
110 return get("energy")->get_value(_unit);
111}
112
119 */
121{
122 if (!has_property("lifetime")) return std::nan("1");
123 if (gunits->get_type(_unit) != tkunit_manager::units_type::kTime && gunits->get_type(_unit) != tkunit_manager::units_type::kEnergy) {
124 glog << error << gunits->get_name(_unit) << " is not of 'time' or 'energy' type (returns nan(1))" << do_endl;
125 return std::nan("1");
126 }
127 return get("lifetime")->get_value(_unit);
128}
129
135{
136 if (!has_property("lifetime")) return "";
137 tkstring lifetime_str;
138
139 if (is_stable())
140 lifetime_str = "STABLE";
141 else {
142 vector<tkunit_manager::units_keys> units{tkunit_manager::y, tkunit_manager::d, tkunit_manager::h, tkunit_manager::min,
146
147 for (auto &unit : units) {
148 double lifetime = get_lifetime(unit);
149 if (lifetime > 1.) {
150 lifetime_str = tkstring::Form("%.3g %s", lifetime, std::get<0>(tkunit_manager::funits_properties.at(unit)).data());
151 break;
152 }
153 }
154 }
155
156 return lifetime_str;
157}
158
160
161 * @param _with_tentative if enabled, consider also the tentative spin assignments
162 * @details
163 * This method returns true if the level is Yrast: corresponding to the lowest energy level for a given spin value
164 */
165bool tklevel::is_yrast(bool _with_tentative)
166{
167 if (_with_tentative)
168 return fisYrast_uncertain;
169 else
170 return fisYrast_exact;
171}
172
178 */
179void tklevel::print(const tkstring &_option)
180{
181 // in line printout
182 if (_option.contains("tab")) {
183 std::cout << std::setw(15) << left << flevel_id;
184 if (is_energy_offset()) std::cout << get_offset_bandhead() << "+";
185 std::cout << std::setw(15) << left << get_energy_measure()->get_value();
186 if (get_energy_measure()->get_error() > 0.)
187 std::cout << std::setw(15) << left << get_energy_measure()->get_error();
188 else
189 std::cout << std::setw(15) << left << " ";
190 if (fspin_parity)
191 std::cout << std::setw(15) << left << fspin_parity->get_jpi_str();
192 if (get_lifetime_measure() && get_lifetime_measure()->get_value() > 0.) {
193 std::cout << std::setw(15) << left << get_lifetime_measure()->get_value();
194 if (get_lifetime_measure()->get_error() > 0.)
195 std::cout << std::setw(15) << left << get_lifetime_measure()->get_error();
196 else
197 std::cout << std::setw(15) << left << " ";
198 std::cout << std::setw(15) << left << get_lifetime_measure()->get_unit();
199 } else {
200 std::cout << std::setw(15) << left << " "
201 << std::setw(15) << left << " "
202 << std::setw(15) << left << " ";
204 std::cout << std::endl;
205 } else {
206 constexpr size_t print_label_width = 14;
207 std::ostringstream level_message;
208 level_message << std::left << std::setw(print_label_width) << "Nuclear level:";
209
210 tkmeasure energy_print = *get_energy_measure();
211 energy_print.set_type("energy");
212
213 glog << info << level_message.str() << " " << &energy_print;
214 if (fspin_parity)
215 glog << " ; Jpi: " << setw(4) << fspin_parity->get_jpi_str();
216 if (get_lifetime_measure() && get_lifetime_measure()->get_value() > 0.)
217 glog << " ; " << get_lifetime_measure();
218 if (is_isomer()) {
219 if (get_isomer_level() == 1)
220 glog << " [1rst isomer]";
221 else if (get_isomer_level() == 2)
222 glog << " [2nd isomer]";
223 else if (get_isomer_level() == 3)
224 glog << " [3rd isomer]";
225 else
226 glog << " [" << get_isomer_level() << "th isomer]";
227 }
228 if (is_uncertain()) {
229 glog << " uncertain level tag: " << funcertain_level;
230 }
231 glog << do_endl;
232
233 if (has_comment() && _option.contains("com")) {
234 const size_t comment_prefix_size = 13;
235 const size_t line_offset = comment_prefix_size + print_label_width + 1;
236 const std::string first_indent(line_offset - comment_prefix_size, ' ');
237 glog << comment << first_indent << wrap_text(get_comment(), line_offset, line_offset) << do_endl;
238 }
239 }
240}
241
242#ifdef HAS_ROOT
244ClassImp(tklevel);
245#endif
Represents a nuclear excited level (or ground state).
Definition tklevel.h:43
bool is_isomer()
returns true is the level an isomer
Definition tklevel.h:104
bool is_stable()
returns true if the level stable
Definition tklevel.h:102
double get_lifetime(const tkunit_manager::units_keys _unit=tkunit_manager::units_keys::s)
returns the lifetime (in second by default)
Definition tklevel.cpp:203
shared_ptr< tkmeasure > get_lifetime_measure()
returns the lifetime tkmeasure object
Definition tklevel.h:85
const tkstring & get_comment()
get the level comment string
Definition tklevel.h:114
double get_energy(const tkunit_manager::units_keys _unit=tkunit_manager::units_keys::keV, bool _with_offset=false)
returns the energy in keV by default
Definition tklevel.cpp:178
tkstring get_offset_bandhead()
returns the bandhead string offset value (ex: X)
Definition tklevel.h:82
shared_ptr< tkmeasure > get_energy_measure()
returns the energy tkmeasure object
Definition tklevel.h:75
void set_jpi(double _spin, int _parity, const tkstring &_jpi_str)
define the spin parity
Definition tklevel.cpp:161
tklevel(int _id, const tkdb_table::measure_data_struct &_struc)
Definition tklevel.cpp:119
tkstring get_lifetime_str()
returns the lifetime in string (using the best adapted unit)
Definition tklevel.cpp:217
bool is_uncertain()
check if the level is uncertain
Definition tklevel.h:116
bool is_energy_offset()
to know if the energy is known with an energy offset
Definition tklevel.h:80
void print(const tkstring &_option="")
print the level properties
Definition tklevel.cpp:262
int get_isomer_level()
get the isomer level (1 for the ␌1rst (lowest energy) isomer, 2 for the second, etc....
Definition tklevel.h:106
bool is_yrast(bool _with_tentative=false)
returns true if the level is yrast
Definition tklevel.cpp:248
void set_lifetime(tkdb_table::measure_data_struct &_struc)
define the level lifetime
Definition tklevel.cpp:140
bool has_comment()
check if a comment exists
Definition tklevel.h:112
Stores a physical measurement with its value, unit, and uncertainty.
Definition tkmeasure.h:56
void set_type(const tkstring &_type)
set the measured data type (energy, lifetime...)
Definition tkmeasure.h:87
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
Nuclear excited state spin-parity.
std::string with usefull tricks from TString (ROOT) and KVString (KaliVeda) and more....
Definition tkstring.h:32
static const char * form(const char *_format,...)
Definition tkstring.cpp:452
static tkstring Form(const char *_format,...)
Definition tkstring.cpp:366
bool is_empty() const
Definition tkstring.h:141
tkstring substr(size_type __pos=0, size_type __n=npos) const
Inlines.
Definition tkstring.h:157
bool contains(const char *_pat, ECaseCompare _cmp=kExact) const
Definition tkstring.h:175
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
std::string wrap_text(const tkstring &_text, size_t _first_content_col, size_t _continuation_col, size_t _max_line_width=80)
Definition tkstring.cpp:37
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
tklog & comment(tklog &log)
Definition tklog.h:319
data structure used to fill a tkmeasure object from the sqlite database
Definition tkdb_table.h:49