BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
ResetEtsAlg.cxx
Go to the documentation of this file.
1#include "GaudiKernel/MsgStream.h"
2#include "GaudiKernel/SmartDataPtr.h"
3#include <iomanip>
4#include <typeinfo>
5
6#include "CalibData/Ets/CorrectedETSCal.h"
7#include "EventModel/EventHeader.h"
8#include "ResetEtsAlg.h"
9#include <algorithm>
10#include <fstream>
11
12using namespace std;
13
15ResetEtsAlg::ResetEtsAlg( const std::string& name, ISvcLocator* pSvcLocator )
16 : Algorithm( name, pSvcLocator )
17 , m_run( -1 )
18 , m_count( 0 )
19 , m_alreadyInDB( false )
20 , m_curIndex( 0 ) {
21 declareProperty( "dump", m_dump = false );
22 declareProperty( "ReadInjSigIntervalFromDB", m_readInjSigIntervalFromDB = true );
23 declareProperty( "ReadInjSigTimeFromDB", m_readInjSigTimeFromDB = true );
24 declareProperty( "InjSigInterval", m_interval = 60 );
25 declareProperty( "RunListFile", m_runListFile );
26 m_printInterval = true;
27}
28
30 MsgStream log( msgSvc(), name() );
31
32 StatusCode sc = service( "CalibDataSvc", m_calibDataSvc );
33 if ( sc.isFailure() )
34 {
35 log << MSG::FATAL << "can not get CalibDataSvc" << endmsg;
36 return StatusCode::FAILURE;
37 }
38
39 sc = service( "InjSigIntervalSvc", m_InjSigIntervalSvc );
40 if ( sc != StatusCode::SUCCESS )
41 { log << MSG::FATAL << "can not use InjSigIntervalSvc" << endmsg; }
42
43 sc = service( "InjSigTimeSvc", m_InjSigTimeSvc );
44 if ( sc != StatusCode::SUCCESS )
45 { log << MSG::FATAL << "can not use InjSigTimeSvc" << endmsg; }
46
47 if ( !m_runListFile.empty() )
48 {
49 int runNo;
50 std::ifstream ifs( m_runListFile.c_str() );
51 ifs >> runNo;
52 while ( ifs.good() )
53 {
54 m_runList.push_back( runNo );
55 ifs >> runNo;
56 }
57 m_fixer = NULL;
58 }
59 else { m_fixer = new EtsFixing(); }
60
61 return StatusCode::SUCCESS;
62}
63
65 MsgStream log( msgSvc(), name() );
66
67 SmartDataPtr<Event::EventHeader> eventHeader( eventSvc(), "/Event/EventHeader" );
68
69 // the corrected ETS is pre-processed
70 if ( m_fixer == NULL )
71 {
72 // first event
73 if ( m_run == -1 )
74 {
75 m_run = eventHeader->runNumber();
76 if ( std::find( m_runList.begin(), m_runList.end(), m_run ) != m_runList.end() )
77 {
78 SmartDataPtr<CalibData::CorrectedETSCal> calConst( m_calibDataSvc,
79 "/Calib/CorrectedETS" );
80 if ( calConst ) { m_alreadyInDB = true; }
81 else
82 {
83 log << MSG::FATAL << "no result in DB for run " << m_run << endmsg;
84 return StatusCode::FAILURE;
85 }
86 }
87 }
88
89 // use corrected ETS in database if available
90 if ( m_alreadyInDB )
91 {
92 if ( m_run != eventHeader->runNumber() )
93 {
94 log << MSG::FATAL << "crossed runs are not supported" << endmsg;
95 return StatusCode::FAILURE;
96 }
97
98 SmartDataPtr<CalibData::CorrectedETSCal> calConst( m_calibDataSvc,
99 "/Calib/CorrectedETS" );
100
101 if ( m_count >= calConst->getNpar() ||
102 eventHeader->eventNumber() != calConst->getEvt( m_count ) )
103 {
104 log << MSG::FATAL << "only 1 input file is supported" << endmsg;
105 return StatusCode::FAILURE;
106 }
107
108 eventHeader->setRawEtsT1( eventHeader->etsT1() );
109 eventHeader->setEtsT1( calConst->getEts1( m_count ) );
110 eventHeader->setEtsT2( calConst->getEts2_pre( m_count ) );
111 eventHeader->setEtsFlag( calConst->getFlag_pre( m_count ) );
112
113 ++m_count;
114
115 return StatusCode::SUCCESS;
116 }
117 /// else: a good run
118 /// no need to correct ETS, and IST will be reset later
119 }
120 else
121 {
122 /// no pre-processed result is available, and then calculate ETS on the fly ...
123 m_fixer->fixT1( eventHeader );
124 }
125
126 /// reset IST and the flag...
127
128 if ( m_readInjSigTimeFromDB )
129 {
130 int Npar = m_InjSigTimeSvc->getNpar();
131 for ( int i = 0; i < Npar; i++ )
132 {
133 int flag = m_InjSigTimeSvc->getFlag( i );
134 ULong64_t time = m_InjSigTimeSvc->getIST( i );
135 m_vecFlag.push_back( flag );
136 m_vecTime.push_back( time );
137 if ( m_dump ) { std::cout << "flag: " << flag << " time: " << time << std::endl; }
138 }
139 m_readInjSigTimeFromDB = false;
140 }
141
142 if ( m_readInjSigIntervalFromDB )
143 {
144 m_interval = double( m_InjSigIntervalSvc->getTInterval() );
145 m_readInjSigIntervalFromDB = false;
146 }
147
148 if ( m_printInterval )
149 {
150 cout << "ResetEts::execute() m_interval = " << m_interval << endl;
151 m_printInterval = false;
152 }
153
154 unsigned long t1 = eventHeader->etsT1(); // 500ns
155 unsigned long t2 = eventHeader->etsT2(); // 500ns
156
157 if ( m_dump )
158 {
159
160 std::cout << "raw data" << std::endl
161 << "Event: " << eventHeader->eventNumber()
162 << " run: " << eventHeader->runNumber() << " time: " << eventHeader->time()
163 << " ETS_old: " << t1 << " IST_old: " << t2 << std::endl;
164 }
165
166 int unit = 2000;
167 double size = m_interval * unit;
168 unsigned long ets2_pre = 0;
169 int flag_pre = 0;
170
171 int vec_size = m_vecTime.size();
172 for ( int i = m_curIndex; i < vec_size; i++ )
173 {
174 if ( t1 < m_vecTime[i] )
175 {
176 if ( i == 0 )
177 {
178 flag_pre = m_vecFlag[i];
179 ets2_pre = m_vecTime[i];
180 m_curIndex = i;
181 break;
182 }
183 else
184 {
185
186 double delta = ( (double)m_vecTime[i - 1] - (double)t1 ) / size;
187
188 if ( delta > 0 && t1 != 0 )
189 {
190 int num = ceil( delta );
191 i = i - num - 1;
192 if ( i < 0 )
193 {
194 i = -1; // for protect
195 }
196 continue;
197 }
198 else
199 {
200 ets2_pre = m_vecTime[i - 1];
201 flag_pre = m_vecFlag[i - 1];
202 m_curIndex = i;
203 break;
204 }
205 }
206 }
207 else if ( t1 >= m_vecTime[vec_size - 1] )
208 {
209 flag_pre = m_vecFlag[vec_size - 1];
210 ets2_pre = m_vecTime[vec_size - 1];
211 break;
212 }
213 }
214
215 // eventHeader->setEtsT1(t1);
216 eventHeader->setEtsT2( ets2_pre );
217 eventHeader->setEtsFlag( flag_pre );
218
219 if ( m_dump )
220 {
221
222 std::cout << "ReSetETSAlg" << std::endl
223 << "Event: " << eventHeader->eventNumber()
224 << " run: " << eventHeader->runNumber() << " time: " << eventHeader->time()
225 << " ETS_new: " << eventHeader->etsT1() << " IST_new: " << eventHeader->etsT2()
226 << " flag_pre: " << eventHeader->etsFlag() << std::endl;
227 }
228
229 return StatusCode::SUCCESS;
230}
231
233 if ( m_fixer != NULL ) { delete m_fixer; }
234
235 return StatusCode::SUCCESS;
236}
DECLARE_COMPONENT(BesBdkRc)
int runNo
Definition DQA_TO_DB.cxx:13
Double_t time
*******INTEGER m_nBinMax INTEGER m_NdiMax !No of bins in histogram for cell exploration division $ !Last vertex $ !Last active cell $ !Last cell in buffer $ !No of sampling when dividing cell $ !No of function total $ !Flag for random ceel for $ !Flag for type of for WtMax $ !Flag which decides whether vertices are included in the sampling $ entire domain is hyp !Maximum effective eevents per saves r n generator level $ !Flag for chat level in !Latex Output unit
Definition FoamA.h:90
IMessageSvc * msgSvc()
StatusCode initialize()
StatusCode finalize()
ResetEtsAlg(const std::string &name, ISvcLocator *pSvcLocator)
StatusCode execute()