BOSS 7.0.8
BESIII Offline Software System
Loading...
Searching...
No Matches
EtsFixing Class Reference

#include <EtsFixing.h>

Public Member Functions

 EtsFixing (int preRefCount=30, double factor=1.0000115)
 
virtual ~EtsFixing ()
 
void fixT1 (Event::EventHeader *header)
 

Detailed Description

Definition at line 6 of file EtsFixing.h.

Constructor & Destructor Documentation

◆ EtsFixing()

EtsFixing::EtsFixing ( int  preRefCount = 30,
double  factor = 1.0000115 
)

Definition at line 21 of file EtsFixing.cxx.

22 : m_nPre(preRefCount),
23 m_factor(factor),
24 m_t0Sec(0),
25 m_t0SecLocal(0),
26 m_count(0),
27 m_refSlope(0.0),
28 m_critical(OneSecVar+999),
29 m_pileup(0)
30{
31 m_preEvt = new int[m_nPre];
32 m_preSec = new int[m_nPre];
33 m_preT1 = new long[m_nPre];
34
35 m_preT1Sec = new int[m_nPre];
36 m_preT1Old = new long[m_nPre];
37 m_prePileup = new long[m_nPre];
38
39 for ( int i = 0; i < m_nPre; ++i ) {
40 m_preEvt[i] = -999;
41 m_preSec[i] = -999;
42 m_preT1[i] = -999;
43
44 m_preT1Sec[i] = -999;
45 m_preT1Old[i] = -999;
46 m_prePileup[i] = -999;
47 }
48
49 m_lastEvt = new int[NLastBackup];
50 m_lastT1 = new long[NLastBackup];
51
52 for ( int i = 0; i < NLastBackup; ++i ) {
53 m_lastEvt[i] = -999;
54 m_lastT1[i] = -999;
55 }
56
57 //if (_debug) {
58 // _file = new TFile("fff.root", "RECREATE");
59 // _tree = new TTree("ff", "for test");
60 // _tree->Branch("count", &_count);
61 // _tree->Branch("event", &_event);
62 // _tree->Branch("time", &_time);
63 // _tree->Branch("pileup",&_pileup);
64 // _tree->Branch("t1", &_t1, "t1/L");
65 // _tree->Branch("rawt1", &_rawt1, "rawt1/L");
66 // _tree->Branch("tdiff", &_tdiff, "tdiff/L");
67 // _tree->Branch("slope", &_slope, "slope/F");
68 //}
69}

◆ ~EtsFixing()

EtsFixing::~EtsFixing ( )
virtual

Definition at line 71 of file EtsFixing.cxx.

72{
73 delete [] m_preEvt;
74 delete [] m_preSec;
75 delete [] m_preT1;
76
77 delete [] m_lastEvt;
78 delete [] m_lastT1;
79
80 delete [] m_preT1Sec;
81 delete [] m_preT1Old;
82 delete [] m_prePileup;
83
84 //if (_debug) {
85 // _file->Write();
86 // _file->Close();
87 // delete _file;
88 //}
89}

Member Function Documentation

◆ fixT1()

void EtsFixing::fixT1 ( Event::EventHeader header)

fix with ideal slope

adjust with event build time

fix with neighboring event

the time difference should not be too large for very near (event id) events

use nearest event if no neighboring event is found

in case of negative or abnormal slope

update slope with reliable events

update parameters for future event Note: DO NOT use the finalT1 here. t1New and finalT1 are not consistent

Refined fixing of ETS with the factor (1.0000115 by default)

Definition at line 91 of file EtsFixing.cxx.

92{
93 long t1Old = header->etsT1();
94
95 if ( t1Old == 0 ) {
96 return;
97 }
98
99 int evtNo = header->eventNumber();
100 int tnSec = header->time() - m_t0Sec;
101 long t1NanoSec = t1Old % OneSecVar;
102 long t1New = t1NanoSec; //tnSec*OneSecVar + t1NanoSec;
103
104 header->setRawEtsT1(t1Old);
105
106 if ( m_count != 0 ) {
107
108 int curIdx = tnSec % m_nPre;
109
110 int preIdx = curIdx;
111 int secShift = 0;
112 if ( tnSec != m_preSec[preIdx] ) {
113 preIdx = (preIdx+m_nPre-1) % m_nPre;
114 secShift = 1;
115 while ( m_preEvt[preIdx] == -999 || m_preSec[preIdx] < m_preSec[(preIdx+1)%m_nPre] ) {
116 preIdx = (preIdx+m_nPre-1) % m_nPre;
117 ++secShift;
118 }
119 if ( (m_preSec[preIdx]+secShift) <= (tnSec-m_nPre) ) {
120 secShift += m_nPre * ((tnSec-secShift-m_preSec[preIdx]) / m_nPre);
121 }
122 //if (_debug > 2) std::cout << tnSec << ", " << m_preSec[preIdx] << ", " << secShift << std::endl;
123 }
124 t1New += (m_preT1[preIdx]/OneSecVar + secShift) * OneSecVar;
125 //if (_debug > 1) std::cout << __LINE__ << ": " << t1New << std::endl;
126
127 /// fix with ideal slope
128 int eDiff = evtNo - m_preEvt[preIdx];
129 long tDiff = t1New - m_preT1[preIdx];
130 long t1Expect = m_preT1[preIdx] + m_refSlope * eDiff;
131 long t1Shift = t1New - t1Expect;
132 if ( t1Shift > ShiftThreshold ) {
133 t1New -= OneSecVar;
134 tDiff -= OneSecVar;
135 }
136 else if ( t1Shift < -ShiftThreshold ) {
137 t1New += OneSecVar;
138 tDiff += OneSecVar;
139 }
140 //if (_debug) _tdiff = t1New - t1Expect;
141 //if (_debug > 1) std::cout << __LINE__ << ": " << t1New << std::endl;
142
143 /// adjust with event build time
144 if (tnSec - m_t0SecLocal > 3) {
145 long tXX = long(tnSec)*OneSecVar - t1New - long(tnSec-t1Old/OneSecVar)*OneSecVar*(m_factor-1);
146 //if (_debug > 2) std::cout << "xxx: " << tnSec << ", " << m_critical << ", " << tXX << std::endl;
147 //expected range: (-m_critical, OneSecVar-m_critical], with a toleration of 0.04s
148 while ( tXX < -m_critical-OneSecVar*0.04 ) {
149 t1New -= OneSecVar;
150 tDiff -= OneSecVar;
151 tXX += OneSecVar;
152 }
153 while ( tXX > OneSecVar*1.04-m_critical ) {
154 t1New += OneSecVar;
155 tDiff += OneSecVar;
156 tXX -= OneSecVar;
157 }
158 }
159 else {
160 if ( secShift == 1 ) {
161 m_critical = t1New-long(tnSec-1)*OneSecVar;
162 //if (_debug > 2) std::cout << "zzz: " << tnSec << ", " << t1New << ", " << t1Old << std::endl;
163 }
164 }
165 //if (_debug > 1) std::cout << __LINE__ << ": " << t1New << std::endl;
166
167 /// fix with neighboring event
168 bool noNeighbor = true;
169 int nearestIdx = NLastBackup-1;
170 int nearestDE = 10000000;
171 for (int i = nearestIdx; i >= 0; --i) {
172 int lidx = (m_count+i)%NLastBackup;
173 int lediff = evtNo-m_lastEvt[lidx];
174 if (abs(lediff) < 20) {
175 long ltdiff = t1New - m_lastT1[lidx];
176 /// the time difference should not be too large for very near (event id) events
177 while ( ltdiff > ShiftThreshold ) {
178 t1New -= OneSecVar;
179 tDiff -= OneSecVar;
180 ltdiff -= OneSecVar;
181 }
182 while ( ltdiff < -ShiftThreshold ) {
183 t1New += OneSecVar;
184 tDiff += OneSecVar;
185 ltdiff += OneSecVar;
186 }
187 if ( lediff > 0 ) {
188 while ( ltdiff < 0 ) {
189 t1New += OneSecVar;
190 tDiff += OneSecVar;
191 ltdiff += OneSecVar;
192 }
193 }
194 else {
195 while ( ltdiff > 0 ) {
196 t1New -= OneSecVar;
197 tDiff -= OneSecVar;
198 ltdiff -= OneSecVar;
199 }
200 }
201 noNeighbor = false;
202 break;
203 }
204 else {
205 if (abs(lediff) < nearestDE) {
206 nearestIdx = i;
207 }
208 }
209 }
210 //if (_debug > 1) std::cout << __LINE__ << ": " << t1New << std::endl;
211
212 /// use nearest event if no neighboring event is found
213 if ( noNeighbor ) {
214 int lidx = (m_count+nearestIdx)%NLastBackup;
215 int lediff = evtNo-m_lastEvt[lidx];
216 long ltdiff = t1New - m_lastT1[lidx];
217 /// in case of negative or abnormal slope
218 if ( lediff > 0 ) {
219 while ( ltdiff*1.0/lediff < m_refSlope*0.2 ) {
220 t1New += OneSecVar;
221 tDiff += OneSecVar;
222 ltdiff += OneSecVar;
223 if ( ltdiff > OneSecVar && ltdiff*1.0/lediff > m_refSlope*5 ) {
224 t1New -= OneSecVar;
225 tDiff -= OneSecVar;
226 ltdiff -= OneSecVar;
227 break;
228 }
229 }
230 }
231 else {
232 while ( ltdiff*1.0/lediff > m_refSlope*5 ) {
233 t1New += OneSecVar;
234 tDiff += OneSecVar;
235 ltdiff += OneSecVar;
236 }
237 while ( ltdiff > 0 ) {
238 t1New -= OneSecVar;
239 tDiff -= OneSecVar;
240 ltdiff -= OneSecVar;
241 }
242 while ( ltdiff*1.0/lediff < m_refSlope*0.2 ) {
243 t1New -= OneSecVar;
244 tDiff -= OneSecVar;
245 ltdiff -= OneSecVar;
246 }
247 }
248 }
249 //if (_debug > 1) std::cout << __LINE__ << ": " << t1New << std::endl;
250
251 /// update slope with reliable events
252 float curSlope = (1.0 * tDiff) / eDiff;
253 if ( curSlope > 0 ) {
254 if ( m_count > 30) {
255 float guard = curSlope / m_refSlope;
256 if ( guard > 0.75 && guard < 1.5 ) {
257 m_refSlope = m_refSlope*0.9 + curSlope*0.1;
258 }
259 }
260 else {
261 m_refSlope = (m_refSlope*(m_count-1) + curSlope) / m_count;
262 }
263 }
264
265 /// update parameters for future event
266 /// Note: DO NOT use the finalT1 here. t1New and finalT1 are not consistent
267 if ( secShift == 1 ) {
268 //do not update these parameters when an event comes earlier more than 1 second
269 m_preEvt[curIdx] = evtNo;
270 m_preSec[curIdx] = tnSec;
271 m_preT1[curIdx] = t1New;
272 }
273 m_lastEvt[m_count%NLastBackup] = evtNo;
274 m_lastT1[m_count%NLastBackup] = t1New;
275
276 /// Refined fixing of ETS with the factor (1.0000115 by default)
277 if ( t1New != t1Old ) {
278 int t1Sec = t1New/OneSecVar;
279 int fIdx = t1Sec % m_nPre;
280 if ( t1Sec == m_preT1Sec[fIdx] ) {
281 m_pileup = m_prePileup[fIdx];
282 }
283 else {
284 for (int i = 1; i < m_nPre; ++i) {
285 int pT1Sec = t1Sec - i;
286 int pFIdx = (pT1Sec + m_nPre) % m_nPre;
287 if ( pT1Sec == m_preT1Sec[pFIdx] ) {
288 int fTotalPileup = t1Sec - t1Old/OneSecVar;
289 int pTotalPileup = m_preT1Sec[pFIdx] - m_preT1Old[pFIdx]/OneSecVar;
290 if ( fTotalPileup != pTotalPileup ) {
291 m_pileup = m_prePileup[pFIdx] + long(fTotalPileup - pTotalPileup)*OneSecVar;
292 }
293 else { //the GPS clock is recovered
294 m_pileup = 0;
295 }
296 break;
297 }
298 }
299 m_preT1Sec[fIdx] = t1Sec;
300 m_preT1Old[fIdx] = t1Old;
301 m_prePileup[fIdx] = m_pileup;
302 }
303 long finalT1 = t1New + (m_factor-1.0)*m_pileup;
304 if ( finalT1 > 0 ) {
305 header->setEtsT1(finalT1);
306 }
307 else {
308 header->setEtsT1(0);
309 }
310 }
311
312 //if (_debug > 1) std::cout << "event " << evtNo << " : "
313 // << header->etsT1() << " - " << header->rawEtsT1()
314 // << " = " << long(header->etsT1()) - long(header->rawEtsT1())
315 // << ", " << m_refSlope
316 // << std::endl;
317 }
318 else { //first event
319 int firstSec = t1Old / OneSecVar;
320 int firstIdx = firstSec % m_nPre;
321 m_t0Sec = header->time() - firstSec;
322 m_t0SecLocal = firstSec;
323 m_preEvt[firstIdx] = evtNo;
324 m_preSec[firstIdx] = firstSec;
325 m_preT1[firstIdx] = t1Old;
326
327 m_lastEvt[0] = evtNo;
328 m_lastT1[0] = t1Old;
329
330 m_pileup = 0;
331 m_preT1Sec[firstIdx] = firstSec;
332 m_preT1Old[firstIdx] = t1Old;
333 m_prePileup[firstIdx] = 0;
334 }
335
336 ++m_count;
337
338 //if (_debug) {
339 // _count = m_count;
340 // _event = header->eventNumber();
341 // _time = header->time();
342 // _pileup = m_pileup/OneSecVar;
343 // _t1 = header->etsT1();
344 // _rawt1 = header->rawEtsT1();
345 // _slope = m_refSlope;
346 // _tree->Fill();
347 //}
348}
void setRawEtsT1(unsigned long value)
Definition: EventHeader.h:99
int eventNumber() const
Retrieve event number.
Definition: EventHeader.h:37
unsigned int time() const
Definition: EventHeader.h:46
void setEtsT1(unsigned long value)
Update ETS.
Definition: EventHeader.h:69
unsigned long etsT1()
Retrieve ETS.
Definition: EventHeader.h:63

Referenced by ResetEtsAlg::execute().


The documentation for this class was generated from the following files: