BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
GetRawETS.cxx
Go to the documentation of this file.
1#include "GaudiKernel/MsgStream.h"
2#include "GaudiKernel/SmartDataPtr.h"
3#include <iomanip>
4
5#include "EventModel/EventHeader.h"
6
7#include "GetRawETS.h"
8
9using namespace std;
11GetRawETS::GetRawETS( const std::string& name, ISvcLocator* pSvcLocator )
12 : Algorithm( name, pSvcLocator ), m_flag( 0 ), m_ets2( 0 ) {
13 declareProperty( "dump", m_dump = false );
14 declareProperty( "InjectionInterval", m_interval = 60 ); /*in ms*/
15 declareProperty( "NumFill", m_nfill = 1 );
16 declareProperty( "FirstEvtInTopUpFlag", m_1stFlag = 0 );
17 declareProperty( "ETS_FILE", m_etsfile = "ets2.root" );
18 declareProperty( "ReadFromDB", m_readFromDB = true );
19
20 m_printInterval = true;
21}
22
23// static TTree* _TestTree;
24// static unsigned long _Event;
25// static unsigned int _EvtTime;
26// static unsigned long _ETS1;
27// static unsigned long _ETS2;
28
30 MsgStream log( msgSvc(), name() );
31
32 m_tree = new TTree( "ist", "The IST time for BES3 raw data" );
33 m_tree->Branch( "flag", &m_flag, "flag/I" );
34 m_tree->Branch( "IST", &m_ets2, "IST/l" );
35
36 m_root = TFile::Open( m_etsfile.c_str(), "RECREATE" );
37 m_tree->SetDirectory( m_root );
38
39 //_TestTree = new TTree("test", "The test TTree");
40 //_TestTree->Branch("event", &_Event, "event/I");
41 //_TestTree->Branch("time", &_EvtTime, "time/I");
42 //_TestTree->Branch("ets1", &_ETS1, "ets1/l");
43 //_TestTree->Branch("ets2", &_ETS2, "ets2/l");
44 //_TestTree->SetDirectory(m_root);
45
46 StatusCode sc = service( "InjSigIntervalSvc", m_InjSigIntervalSvc );
47 if ( sc != StatusCode::SUCCESS )
48 { log << MSG::FATAL << "can not use InjSigIntervalSvc" << endmsg; }
49 return StatusCode::SUCCESS;
50}
51
52StatusCode GetRawETS::execute() {
53 MsgStream log( msgSvc(), name() );
54
55 SmartDataPtr<Event::EventHeader> eventHeader( eventSvc(), "/Event/EventHeader" );
56 unsigned long t1 = eventHeader->etsT1(); // 500ns
57 unsigned long t2 = eventHeader->etsT2(); // 500ns
58
59 if ( m_readFromDB )
60 {
61 m_interval = m_InjSigIntervalSvc->getTInterval();
62 m_readFromDB = false;
63 }
64
65 if ( m_printInterval )
66 {
67 cout << "GetRawETS::execute() m_interval = " << m_interval << endl;
68 m_printInterval = false;
69 }
70 if ( m_dump )
71 {
72 double dt1 = ( double( t1 ) / double( 2000000. ) ); // in seconds
73 double dt2 = ( double( t2 ) / double( 2000000. ) ); // in seconds
74
75 std::cout << "Event: " << eventHeader->eventNumber()
76 << " run: " << eventHeader->runNumber() << " time: " << eventHeader->time()
77 << " EtsT1: " << t1 << " EtsT2: " << t2 << std::endl;
78 }
79
80 if ( t2 != 0 && t2 != m_ets2 )
81 {
82 m_ets2 = t2;
83 m_vec.push_back( m_ets2 );
84 }
85
86 //{
87 // _Event = eventHeader->eventNumber();
88 // _EvtTime = eventHeader->time();
89 // _ETS1 = t1;
90 // _ETS2 = t2;
91 // _TestTree->Fill();
92 //}
93
94 return StatusCode::SUCCESS;
95}
96
97StatusCode GetRawETS::finalize() {
98 if ( m_vec.empty() )
99 {
100 // DO NOT CRASH when m_vec is empty
101 m_root->Write();
102 return StatusCode::SUCCESS;
103 }
104
105 std::sort( m_vec.begin(), m_vec.end() );
106 m_vec.erase( std::unique( m_vec.begin(), m_vec.end() ), m_vec.end() );
107
108 unsigned long _interval = m_interval * 2000; // -> 500ns
109
110 // Fill the Tree
111 {
112 unsigned long tmp = m_vec.front();
113 for ( int j = m_nfill; j > 0; --j )
114 {
115 if ( tmp <= _interval * j ) continue;
116 m_flag = 10 + j;
117 m_ets2 = tmp - _interval * j;
118 m_tree->Fill(); // case 13,12,11: prepend N addtional ets to current injection period
119 }
120 m_flag = m_1stFlag;
121 m_ets2 = tmp;
122 m_tree->Fill(); // the first real ets
123 }
124 for ( unsigned int i = 1; i < m_vec.size(); ++i )
125 {
126 unsigned long diff = m_vec[i] - m_vec[i - 1];
127 bool _1stEvtInTopUp = false;
128
129 if ( diff > _interval * ( m_nfill + 1 ) * 2 + 10000 )
130 { // long gap: +5ms for safe
131 _1stEvtInTopUp = true;
132 for ( int j = 1; j < m_nfill + 1; ++j )
133 {
134 m_flag = 20 + j;
135 m_ets2 = m_vec[i - 1] + _interval * j;
136 m_tree->Fill(); // case 21,22,23: append N addtional ets to previous injection period
137 }
138 for ( int j = m_nfill; j > 0; --j )
139 {
140 m_flag = 10 + j;
141 m_ets2 = m_vec[i] - _interval * j;
142 m_tree->Fill(); // case 13,12,11: prepend N addtional ets to current injection period
143 }
144 }
145 else if ( diff > _interval + 10000 )
146 { // missing ets2: +5ms for safe
147 for ( int j = 1; j < 10; ++j )
148 {
149 m_flag = 1;
150 m_ets2 = m_vec[i - 1] + _interval * j; // m_interval * ms * j
151 if ( m_ets2 >= m_vec[i] || m_vec[i] - m_ets2 < 10000 )
152 { // 5ms
153 break;
154 }
155 m_tree->Fill(); // case 1: fill the missing ets2 explicitly
156 }
157 }
158
159 m_flag = _1stEvtInTopUp ? m_1stFlag : 0;
160 m_ets2 = m_vec[i];
161 m_tree->Fill(); // case 0: good ets2
162 }
163 {
164 unsigned long tmp = m_vec.back();
165 for ( int j = 1; j < m_nfill + 1; ++j )
166 {
167 m_flag = 20 + j;
168 m_ets2 = tmp + _interval * j;
169 m_tree->Fill(); // case 21,22,23: append N addtional ets to the last injection period
170 }
171 }
172
173 // Write ROOT file
174 m_root->Write();
175
176 return StatusCode::SUCCESS;
177}
DECLARE_COMPONENT(BesBdkRc)
IMessageSvc * msgSvc()
GetRawETS(const std::string &name, ISvcLocator *pSvcLocator)
Definition GetRawETS.cxx:11
StatusCode initialize()
Definition GetRawETS.cxx:29
StatusCode execute()
Definition GetRawETS.cxx:52
StatusCode finalize()
Definition GetRawETS.cxx:97