Belle II Software development
Segment.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9//-----------------------------------------------------------------------------
10// Description : A class to represent a wire in CDC.
11//-----------------------------------------------------------------------------
12
13#define TRG_SHORT_NAMES
14#define TRGCDC_SHORT_NAMES
15
16#include <iostream>
17#include "trg/trg/Utilities.h"
18#include "trg/trg/Debug.h"
19#include "trg/cdc/Wire.h"
20#include "trg/cdc/WireHit.h"
21#include "trg/cdc/Segment.h"
22#include "trg/cdc/SegmentHit.h"
23#include "trg/cdc/LUT.h"
24
25#include <framework/datastore/StoreArray.h>
26#include <cdc/dataobjects/CDCHit.h>
27#include <trg/cdc/dataobjects/CDCTriggerSegmentHit.h>
28#include <mdst/dataobjects/MCParticle.h>
29
30using namespace std;
31
32#define P3D HepGeom::Point3D<double>
33
34namespace Belle2 {
41 const TCLayer& layer,
42 const TCWire& w,
43 const TRGClock& clock,
44 const std::string& TSLUTFile,
45 const std::vector<const TCWire*>& cells)
46 : TCCell(id,
47 layer.size(),
48 layer,
49 w.forwardPosition(),
50 w.backwardPosition()),
51 _wires(cells),
52 _center(), // 2019/07/31 by ytlai
53 _signal(std::string("TS_") + TRGUtil::itostring(id), clock),
54 _storeHits{},
55 m_TSLUTFileName(TSLUTFile)
56 {
57 }
58
59
61 {
62 }
63
64 void
66 {
68 }
69
70 void
71 TRGCDCSegment::dump(const string& msg,
72 const string& pre) const
73 {
74 cout << pre << name() << " (ptn=" << hitPattern() << ")" << endl;
75 if ((msg.find("geometry") != string::npos) ||
76 (msg.find("detail") != string::npos)) {
77 cout << pre << "id " << id();
78 cout << ",local " << localId();
79 cout << ",layer " << layerId();
80 cout << ",super layer " << superLayerId();
81 cout << ",local layer " << localLayerId();
82 cout << endl;
83 }
84 if ((msg.find("hit") != string::npos) ||
85 (msg.find("detail") != string::npos)) {
86 cout << pre << "Wires ";
87 for (unsigned i = 0; i < _wires.size(); i++) {
88 cout << _wires[i]->name();
89 if (i < _wires.size() - 1)
90 cout << ",";
91 else
92 cout << endl;
93 }
94 if (_hits.size() == 0) {
95 cout << pre << "no wire hit" << endl;
96 } else {
97 cout << pre << "WHit dump : ";
98 for (unsigned i = 0; i < _hits.size(); i++) {
99 cout << _hits[i]->cell().name();
100 if (i < _hits.size() - 1)
101 cout << ",";
102 else
103 cout << endl;
104 }
105 for (unsigned i = 0; i < _hits.size(); i++) {
106 _hits[i]->dump(msg, pre + " ");
107 }
108 }
109 if (hit()) {
110 cout << pre << "SHit dump" << endl;
111 hit()->dump(msg, pre + " ");
112 } else {
113 cout << pre << "no TSHit" << endl;
114 }
115 }
116// if (msg.find("neighbor") != string::npos ||
117// msg.find("detail") != string::npos) {
118// for (unsigned i = 0; i < 7; i++)
119// if (neighbor(i))
120// neighbor(i)->dump("", pre + TRGCDC::itostring(i) + " ");
121// }
122 if ((msg.find("trigger") != string::npos) ||
123 (msg.find("detail") != string::npos)) {
124 if (_signal.active())
125 _signal.dump(msg, pre + " ");
126 else
127 cout << pre << "no trigger signal" << endl;
128 }
129 }
130
131 void
133 {
134 TCCell::clear();
135 _signal.clear();
136 _hits.clear();
137 _storeHits.clear();
138 }
139
140 string
142 {
143 string t;
144 if (axial())
145 t = "-";
146 else
147 t = "=";
148 string n0 = string("TS") + TRGUtil::itostring(layerId());
149 string n1 = TRGUtil::itostring(localId());
150 return n0 + t + n1;
151 }
152
153 void
154 TRGCDCSegment::simulate(bool clockSimulation, bool logicLUTFlag,
155 const std::string& cdcCollectionName, const std::string& tsCollectionName)
156 {
157 //...Get wire informtion for speed-up...
158 unsigned nHits = 0;
159 for (unsigned i = 0, n = _wires.size(); i < n; i++) {
160 if (_wires[i]->signal().active())
161 ++nHits;
162 }
163
164 //..No wire hit case...
165 if (nHits == 0)
166 return;
167
168 if (clockSimulation) {
169 simulateWithClock(cdcCollectionName, tsCollectionName);
170 } else {
171 simulateWithoutClock(logicLUTFlag);
172 }
173 }
174
175 void
177 {
178 TRGDebug::enterStage("TS sim");
179
180 //...Get wire informtion...
181 const unsigned n = _wires.size();
182 unsigned nHits = 0;
183 vector<TRGSignal> signals;
184 for (unsigned i = 0; i < n; i++) {
185
186 //...Store wire hit information...
187 const TCWHit* h = _wires[i]->hit();
188 if (h)
189 _hits.push_back(h);
190
191 //...Copy signal from a wire...
192 const TRGSignal& s = _wires[i]->signal();
193 signals.push_back(s);
194
195 //...Widen it...
196 const unsigned width = signals.back().clock().unit(1000);
197 signals.back().widen(width);
198
199 if (s.active())
200 ++nHits;
201 }
202
203 if (logicLUTFlag == 0) {
205 //...Check number of hit wires...
206 //cout<<"TSF: nHits is "<<nHits<<endl;
207 if (nHits < 4) {
208 TRGDebug::leaveStage("TS sim");
209 return;
210 }
211
212 //...Signal simulation...
213 TRGSignal l0, l1, l2, l3, l4;
214 TRGSignal wo1, wo2, wo3, wo4;
215 TRGSignal all;
216 if (n == 11) {
217
218 //...Simple simulation assuming 3:2:1:2:3 shape...
219 l0 = signals[0] | signals[1] | signals[2];
220 l1 = signals[3] | signals[4];
221 l2 = signals[5];
222 l3 = signals[6] | signals[7];
223 l4 = signals[8] | signals[9] | signals[10];
224 //l0.dump();
225 //l1.dump();
226 //l2.dump();
227 //l3.dump();
228 //l4.dump();
229 wo1 = l1 & l3 & l4;
230 wo2 = l0 & l3 & l4;
231 wo3 = l0 & l1 & l4;
232 wo4 = l0 & l1 & l3;
233 all = l2 & (wo1 | wo2 | wo3 | wo4);
234
235 } else if (n == 15) {
236
237 //...Simple simulation assuming 1:2:3:4:5 shape...
238 l0 = signals[0];
239 l1 = signals[1] | signals[2];
240 l2 = signals[3] | signals[4] | signals[5];
241 l3 = signals[6] | signals[7] | signals[8] | signals[9];
242 l4 = signals[10] | signals[11] | signals[12] | signals[13] | signals[14];
243 wo1 = l2 & l3 & l4;
244 wo2 = l1 & l3 & l4;
245 wo3 = l1 & l2 & l4;
246 wo4 = l1 & l2 & l3;
247 all = l0 & (wo1 | wo2 | wo3 | wo4);
248 }
249
250 //...Coincidence of all layers...
251// TRGSignal all = l0 & l1 & l2 & l3 & l4;
252
253 if (all.nEdges()) {
254 //cout<<"TSF is found"<<endl;
255 all.name(name());
256 _signal = all;
257 //cout<<all.name()<<":#signals="<<all.nSignals()<<endl;;
258 //all.dump();
259 }
261 }
262
263 if (logicLUTFlag == 1) {
265 //... Find hit wires ...
266 vector<TRGSignal> hitSignals;
267 for (unsigned iWire = 0; iWire < signals.size(); iWire++) {
268 if (signals[iWire].active()) hitSignals.push_back(signals[iWire]);
269 }
270 //... Coincidence all hit wires ...
271 TRGSignal allSignals;
272 if (hitSignals.size() != 0) {
273 allSignals = hitSignals[0];
274 for (unsigned iHitWire = 1; iHitWire < hitSignals.size(); iHitWire++) {
275 allSignals = allSignals & hitSignals[iHitWire];
276 }
277 }
278
279 int lutValue = LUT()->getValue(lutPattern());
280 if ((lutValue != 0) && (priority().signal().active() != 0)) {
281 allSignals.name(name());
282 _signal = allSignals;
283 }
285 }
286
287 TRGDebug::leaveStage("TS sim");
288 }
289
290 void
291 TRGCDCSegment::simulateWithClock(string cdcCollectionName, string tsCollectionName)
292 {
293 // check LUT pattern without clock -> if there is no hit, skip clock simulation
294 if (m_TSLUT->getValue(lutPattern()) == 0) return;
295
296 TRGDebug::enterStage("TS sim with clock");
297
298 StoreArray<CDCHit> cdcHits(cdcCollectionName);
299 StoreArray<CDCTriggerSegmentHit> segmentHits(tsCollectionName);
300
301 // get data clock of first and last hit
302 const TRGClock& wireClock = _wires[0]->signal().clock();
303 int clkMin = 1000;
304 int clkMax = -1000;
305 for (unsigned i = 0, n = _wires.size(); i < n; ++i) {
306 const TRGSignal& s = _wires[i]->signal();
307 if (s.active()) {
308 int clk0 = s[0]->time();
309 int clk1 = s[s.nEdges() - 2]->time();
310 if (clk0 < clkMin) clkMin = clk0;
311 if (clk1 > clkMax) clkMax = clk1;
312 }
313 }
314 // loop over data clock cycles
315 //const int step = wireClock.frequency() / TRGCDC::getTRGCDC()->dataClock().frequency();
316 const int step = wireClock.frequency() / signal().clock().frequency();
317 const int width = 16 * step;
318 clkMin -= clkMin % step;
319 clkMax -= clkMax % step;
320 int lastLutValue = 0;
321 int lastPriority = 0;
322 int lastFastest = 0;
323 for (int iclk = clkMin; iclk <= clkMax; iclk += step) {
324 // check pattern in the last width clock cycles
325 unsigned pattern = lutPattern(iclk - width, iclk + step);
326 int lutValue = m_TSLUT->getValue(pattern);
327 if (lutValue) {
328 int priorityPos = priorityPosition(iclk - width, iclk + step);
329 int fastest = fastestTime(iclk - width);
330 // make a new hit if L/R changes to known, if priority changes to first
331 // or if the fastest hit changes
332 if ((lastLutValue == 3 && lutValue != 3) ||
333 (lastPriority != 3 && priorityPos == 3) ||
334 fastest != lastFastest) {
335 // add new edge to signal
336 TRGTime rise = TRGTime(wireClock.absoluteTime(iclk), true, _signal.clock());
337 TRGTime fall = rise;
338 fall.shift(1).reverse();
339 _signal |= TRGSignal(rise & fall);
340 // get priority wire from position flag
341 int ipr = (priorityPos == 3) ? 0 : priorityPos;
342 const TRGCDCWire* priorityWire = (_wires.size() == 15) ? _wires[ipr] : _wires[ipr + 5];
343 // get priority time (first hit on priority wire in time window)
344 int tdc = priorityWire->signal()[0]->time();
345 if (tdc < iclk - width) {
346 for (unsigned itdc = 2, edges = priorityWire->signal().nEdges(); itdc < edges; itdc += 2) {
347 tdc = priorityWire->signal()[itdc]->time();
348 if (tdc >= iclk - width) break;
349 }
350 }
351 // create hit
352 const CDCHit* priorityHit = cdcHits[priorityWire->hit()->iCDCHit()];
353 const CDCTriggerSegmentHit* storeHit =
354 segmentHits.appendNew(*priorityHit,
355 id(),
356 priorityPos,
357 lutValue,
358 tdc,
359 fastest,
360 iclk + step);
361 addStoreHit(storeHit);
362 // relation to all CDCHits in segment
363 for (unsigned iw = 0; iw < _wires.size(); ++iw) {
364 if (_wires[iw]->signal().active(iclk - width, iclk + step)) {
365 // priority wire has relation weight 2
366 double weight = (_wires[iw] == priorityWire) ? 2. : 1.;
367 storeHit->addRelationTo(cdcHits[_wires[iw]->hit()->iCDCHit()], weight);
368 }
369 }
370 // relation to MCParticles (same as priority hit)
372 for (unsigned imc = 0; imc < mcrel.size(); ++imc) {
373 mcrel[imc]->addRelationTo(storeHit, mcrel.weight(imc));
374 }
375 // store values of this hit to compare with the next hit
376 lastLutValue = lutValue;
377 lastPriority = priorityPos;
378 lastFastest = fastest;
379 }
380 }
381 }
382
383 TRGDebug::leaveStage("TS sim with clock");
384 }
385
386 float
388 {
389 if ((LUT()->getValue(lutPattern()))) {
390 float tmpFastTime = 9999;
391 for (unsigned i = 0; i < _wires.size(); i++) {
392 if (_wires[i]->signal().active()) {
393 float dt = _wires[i]->signal()[0]->time();
394 if (dt < tmpFastTime) {
395 tmpFastTime = dt;
396 }
397 }
398 }
399 return tmpFastTime;
400 } else
401 return -1;
402 }
403
404 float
406 {
407 int fastest = 9999;
408 for (unsigned iw = 0; iw < _wires.size(); ++iw) {
409 if (_wires[iw]->signal().active()) {
410 for (unsigned itdc = 0, edges = _wires[iw]->signal().nEdges(); itdc < edges; itdc += 2) {
411 float dt = _wires[iw]->signal()[itdc]->time();
412 if (dt >= clk0) {
413 if (dt < fastest) fastest = dt;
414 break;
415 }
416 }
417 }
418 }
419 return fastest;
420 }
421
422 float
424 {
425 if ((LUT()->getValue(lutPattern()))) {
426 float tmpFoundTime[5] = {9999, 9999, 9999, 9999, 9999};
427 for (unsigned i = 0; i < _wires.size(); i++) {
428 if (!_wires[i]->signal().active()) continue;
429 float dt = _wires[i]->signal()[0]->time();
430 if (_wires.size() == 11) {
431 if (i < 3) {
432 if (tmpFoundTime[0] > dt) tmpFoundTime[0] = dt;
433 } else if (i < 5) {
434 if (tmpFoundTime[1] > dt) tmpFoundTime[1] = dt;
435 } else if (i == 5) {
436 if (tmpFoundTime[2] > dt) tmpFoundTime[2] = dt;
437 } else if (i < 8) {
438 if (tmpFoundTime[3] > dt) tmpFoundTime[3] = dt;
439 } else {
440 if (tmpFoundTime[4] > dt) tmpFoundTime[4] = dt;
441 }
442 } else {
443 if (i == 0) {
444 if (tmpFoundTime[0] > dt) tmpFoundTime[0] = dt;
445 } else if (i < 3) {
446 if (tmpFoundTime[1] > dt) tmpFoundTime[1] = dt;
447 } else if (i < 6) {
448 if (tmpFoundTime[2] > dt) tmpFoundTime[2] = dt;
449 } else if (i < 10) {
450 if (tmpFoundTime[3] > dt) tmpFoundTime[3] = dt;
451 } else {
452 if (tmpFoundTime[4] > dt) tmpFoundTime[4] = dt;
453 }
454 }
455 }
456 sort(tmpFoundTime, tmpFoundTime + 5);
457 return tmpFoundTime[3];
458 } else
459 return -1;
460 }
461
462 float
464 {
465 const TRGSignal& prioritySignal = priority().signal();
466 if (prioritySignal.active()) {
467 return prioritySignal[0]->time();
468 }
469 return -1;
470 }
471
472 int
474 {
475 if (center().signal().active()) {
476 return 3;
477 } else {
478 const TRGCDCWire* priorityL;
479 const TRGCDCWire* priorityR;
480 if (_wires.size() == 15) {
481 priorityL = _wires[2];
482 priorityR = _wires[1];
483 } else {
484 priorityL = _wires[7];
485 priorityR = _wires[6];
486 }
487 if (priorityL->signal().active()) {
488 if (priorityR->signal().active()) {
489 if ((priorityL->signal()[0]->time()) >= (priorityR->signal()[0]->time())) return 1;
490 else return 2;
491 } else return 2;
492 } else if (priorityR->signal().active()) {
493 return 1;
494 } else return 0;
495 }
496 }
497
498 int
499 TRGCDCSegment::priorityPosition(int clk0, int clk1) const
500 {
501 if (center().signal().active(clk0, clk1)) {
502 return 3;
503 } else {
504 const TRGCDCWire* priorityL;
505 const TRGCDCWire* priorityR;
506 if (_wires.size() == 15) {
507 priorityL = _wires[2];
508 priorityR = _wires[1];
509 } else {
510 priorityL = _wires[7];
511 priorityR = _wires[6];
512 }
513 if (priorityL->signal().active(clk0, clk1)) {
514 if (priorityR->signal().active(clk0, clk1)) {
515 if ((priorityL->signal()[0]->time()) >= (priorityR->signal()[0]->time())) return 1;
516 else return 2;
517 } else return 2;
518 } else if (priorityR->signal().active(clk0, clk1)) {
519 return 1;
520 } else return 0;
521 }
522 }
523
524 const TRGCDCWire&
526 {
528 int offset = (_wires.size() == 15) ? 0 : 5;
529 if (priority == 1 || priority == 2)
530 return *_wires[offset + priority];
531 return *_wires[offset];
532 }
533
534 unsigned
536 {
537 unsigned ptn = 0;
538 for (unsigned i = 0; i < _wires.size(); i++) {
539 const TRGSignal& s = _wires[i]->signal();
540 if (s.active())
541 ptn |= (1 << i);
542 }
543 return ptn;
544 }
545
546 unsigned
548 {
549 unsigned ptn = 0;
550 for (unsigned i = 0; i < _wires.size(); i++) {
551 const TRGSignal& s = _wires[i]->signal_adc();
552 if (s.active())
553 ptn |= (1 << i);
554 }
555 return ptn;
556 }
557
558 unsigned
559 TRGCDCSegment::hitPattern(int clk0, int clk1) const
560 {
561 unsigned ptn = 0;
562 for (unsigned i = 0; i < _wires.size(); i++) {
563 const TRGSignal& s = _wires[i]->signal();
564 if (s.active(clk0, clk1))
565 ptn |= (1 << i);
566 }
567 return ptn;
568 }
569
570 unsigned
572 {
573 unsigned outValue = (hitPattern()) * 2;
574 if (priorityPosition() == 2) {
575 outValue += 1;
576 }
577 return outValue;
578 }
579
580 unsigned
581 TRGCDCSegment::lutPattern(int clk0, int clk1) const
582 {
583 unsigned outValue = (hitPattern(clk0, clk1)) * 2;
584 if (priorityPosition(clk0, clk1) == 2) {
585 outValue += 1;
586 }
587 return outValue;
588 }
589
590 bool
591 TRGCDCSegment::hasMember(const std::string& a) const
592 {
593 const unsigned n = _wires.size();
594 for (unsigned i = 0; i < n; i++) {
595 if (_wires[i]->hasMember(a))
596 return true;
597 }
598 return false;
599 }
600
602} // namespace Belle2
Class containing the result of the unpacker in raw data and the result of the digitizer in simulation...
Definition: CDCHit.h:40
Combination of several CDCHits to a track segment hit for the trigger.
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
float weight(int index) const
Get weight with index.
void addRelationTo(const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
Add a relation from this object to another object (with caching).
RelationVector< FROM > getRelationsFrom(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from another store array to this object.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
T * appendNew()
Construct a new T object at the end of the array.
Definition: StoreArray.h:246
std::vector< const TRGCDCWire * > _wires
LookUp Table.
Definition: Segment.h:170
std::vector< const CDCTriggerSegmentHit * > _storeHits
list of DataStore hits.
Definition: Segment.h:185
TRGSignal _signal
Trigger signal.
Definition: Segment.h:176
std::vector< const TRGCDCWireHit * > _hits
Wire hits.
Definition: Segment.h:182
std::string m_TSLUTFileName
TS LUT file name.
Definition: Segment.h:188
TRGCDCLUT * m_TSLUT
LookUp Table. 0: no hit, 1: right, 2: left, 3: not determined.
Definition: Segment.h:164
A class to represent a wire in CDC.
Definition: Wire.h:56
A class to represent a digitized signal. Unit is nano second.
Definition: Clock.h:38
A class to represent a digitized signal. Unit is nano second.
Definition: Signal.h:23
A class to represent a signal timing in the trigger system.
Definition: Time.h:25
virtual bool hasMember(const std::string &a) const override
returns true this has member named a.
Definition: Segment.cc:591
const TRGCDCWire & center(void) const
returns a center wire.
Definition: Segment.h:265
unsigned nEdges(void) const
returns # of edges.
Definition: Signal.h:258
double absoluteTime(int clockPosition) const
returns absolute time of clock position
Definition: Clock.cc:128
int priorityPosition(void) const
return priority cell position in TSHit. 0: no hit, 3: 1st priority, 1: 2nd right, 2: 2nd left
Definition: Segment.cc:473
float fastestTime(void) const
return fastest time in TSHit.
Definition: Segment.cc:387
unsigned id(void) const
returns id.
Definition: Cell.h:200
void addStoreHit(const CDCTriggerSegmentHit *)
sets a pointer to a CDCTriggerSegmentHit.
Definition: Segment.h:244
virtual ~TRGCDCSegment()
Destructor.
Definition: Segment.cc:60
void simulateWithoutClock(bool logicLUTFlag)
simulates time-indegrated TF hit
Definition: Segment.cc:176
const TRGClock & clock(void) const
returns clock.
Definition: Signal.h:331
unsigned layerId(void) const
returns layer id.
Definition: Cell.h:214
void simulateWithClock(std::string cdcCollectionName, std::string tsCollectionName)
simulates TF hit time-dependently
Definition: Segment.cc:291
unsigned localLayerId(void) const
returns local layer id in a super layer.
Definition: Cell.h:228
unsigned hitPattern(void) const
returns hit pattern.
Definition: Segment.cc:535
TRGTime & reverse(void)
reverse edge.
Definition: Time.h:141
int getValue(unsigned) const
get LUT Values
Definition: LUT.cc:47
const std::string & name(void) const
returns name.
Definition: Signal.h:188
static void enterStage(const std::string &stageName)
Declare that you enter new stage.
Definition: Debug.cc:24
unsigned superLayerId(void) const
returns super layer id.
Definition: Cell.h:221
const TRGSignal & signal(void) const override
returns trigger output. Null will returned if no signal.
Definition: Segment.h:214
const TRGCDCWire & priority(void) const
returns priority wire.
Definition: Segment.cc:525
static TRGCDCLUT * getLUT(const std::string &filename, int)
get LUT from dictionary, load new LUT if it doesn't exist
Definition: LUT.cc:83
TRGCDCSegment(unsigned id, const TRGCDCLayer &layer, const TRGCDCWire &w, const TRGClock &clock, const std::string &TSLUTFile, const std::vector< const TRGCDCWire * > &wires)
Constructor.
Definition: Segment.cc:40
TRGTime & shift(int unit)
delays by clock unit.
Definition: Time.h:163
const TRGCDCWireHit * hit(void) const
returns a pointer to a TRGCDCWireHit.
Definition: Wire.h:153
void initialize(void)
initilize variables.
Definition: Segment.cc:65
std::string name(void) const override
returns name.
Definition: Segment.cc:141
unsigned lutPattern(void) const
hit pattern containing bit for priority position
Definition: Segment.cc:571
void simulate(bool clockSimulation, bool logicLUTFlag, const std::string &cdcCollectionName=std::string(""), const std::string &tsCollectionName=std::string(""))
simulates TF hit using wire information.
Definition: Segment.cc:154
const TRGCDCLUT * LUT(void) const
returns LUT
Definition: Segment.h:258
void dump(const std::string &message=std::string(""), const std::string &prefix=std::string("")) const override
dumps debug information.
Definition: Segment.cc:71
void clear(void) override
clears information.
Definition: Segment.cc:132
unsigned iCDCHit(void) const
returns an index to CDCHit.
Definition: CellHit.h:360
bool active(void) const
returns true if there is a signal.
Definition: Signal.h:277
float priorityTime(void) const
return priority time in TSHit.
Definition: Segment.cc:463
static void leaveStage(const std::string &stageName)
Declare that you leave a stage.
Definition: Debug.cc:34
unsigned hitPattern_adc(void) const
returns hit pattern with adc cut.
Definition: Segment.cc:547
float foundTime(void) const
return found time in TSHit.
Definition: Segment.cc:423
void clear(void)
clears contents.
Definition: Signal.h:265
const TRGCDCSegmentHit * hit(void) const
returns a pointer to a TRGCDCSegmentHit.
Definition: Segment.h:236
bool axial(void) const
returns true if this wire is in an axial layer.
Definition: Cell.h:256
double frequency(void) const
returns frequency in MHz.
Definition: Clock.h:187
unsigned localId(void) const
returns local id in a layer.
Definition: Cell.h:207
void dump(const std::string &message="", const std::string &pre="") const
dumps contents.
Definition: Signal.cc:139
Abstract base class for different kinds of events.
STL namespace.