BmnRoot
Loading...
Searching...
No Matches
BmnCSCLayer.cxx
Go to the documentation of this file.
1#include "BmnCSCLayer.h"
2
3Int_t BmnCSCLayer::fUniqueIdL = 0;
4Int_t BmnCSCLayer::fUniqueIdU = 0;
5
7
8 Verbosity = true;
9
10 LayerZoneNumber = 0;
11 LayerType = LowerStripLayer;
12
13 XMinLayer = 0.0;
14 XMaxLayer = 10.0;
15 YMinLayer = 0.0;
16 YMaxLayer = 10.0;
17
18 Pitch = 0.25; //cm
19 AngleDeg = 14.87; //in degrees from a vertical line where a plus value is clockwise
20 AngleRad = AngleDeg*Pi()/180; //in radians
21
22 StripOrder = LeftToRight;
23
24 if(AngleDeg <= 0.0 && AngleDeg >= -90.0) {
26 }
27 else {
28 if(AngleDeg > 0.0 && AngleDeg <= 90.0) {
30 }
31 }
32
33 ClusterFindingThreshold = 0.0;
34}
35
36BmnCSCLayer::BmnCSCLayer(Int_t zone_number, StripLayerType layer_type,
37 Double_t xsize, Double_t ysize,
38 Double_t xorig, Double_t yorig,
39 Double_t pitch, Double_t adeg) {
40
41 Verbosity = true;
42
43 LayerZoneNumber = zone_number;
44 LayerType = layer_type;
45
46 XMinLayer = xorig;
47 XMaxLayer = xorig + xsize;
48 YMinLayer = yorig;
49 YMaxLayer = yorig + ysize;
50
51 Pitch = pitch;
52 AngleDeg = adeg; //in degrees from a vertical line where a plus value is clockwise
53 AngleRad = AngleDeg*Pi()/180; //in radians
54
55 StripOrder = LeftToRight;
56
57 if(AngleDeg <= 0.0 && AngleDeg >= -90.0) {
59 }
60 else {
61 if(AngleDeg > 0.0 && AngleDeg <= 90.0) {
63 }
64 }
65
66 ClusterFindingThreshold = 0.0;
67}
68
72
74 Int_t NStrips = GetNStrips();
75
76 //strips
77 Strips.clear();
78 Strips.resize(NStrips, 0.0);
79
80 //strip matches
83
84 //strip hits
86
87 BmnCSCLayer::fUniqueIdL = 0;
88 BmnCSCLayer::fUniqueIdU = 0;
89}
90
94
96 Double_t xleft_layer = ConvertNormalPointToStripX(XLeftPointOfStripNumbering, YLeftPointOfStripNumbering);
97 Double_t xright_layer = ConvertNormalPointToStripX(XRightPointOfStripNumbering, YRightPointOfStripNumbering);
98 Double_t n_strips = (xright_layer - xleft_layer)/Pitch;
99
100 if( (n_strips - (Int_t)n_strips) < 1E-10 ) {
101 return (Int_t)n_strips;
102 }
103 else {
104 return (Int_t)(n_strips+1);
105 }
106}
107
108void BmnCSCLayer::SetPitch(Double_t pitch) {
109 Pitch = pitch;
110
112}
113
114void BmnCSCLayer::SetLayerSizes(Double_t xsize, Double_t ysize, Double_t xorig, Double_t yorig) {
115 XMinLayer = xorig;
116 XMaxLayer = xorig + xsize;
117 YMinLayer = yorig;
118 YMaxLayer = yorig + ysize;
119
121}
122
123void BmnCSCLayer::SetAngleDeg(Double_t deg) {
124 if(Abs(deg) <= 90.0) {
125 AngleDeg = deg;
126 }
127 else {
128 AngleDeg = 0.0;
129 if(Verbosity) cerr << "WARNING: BmnCSCLayer::SetAngleDeg(): the strip angle is incorrect and has been set to zero\n";
130 }
131
132 AngleRad = AngleDeg*Pi()/180;
133
135}
136
137Bool_t BmnCSCLayer::SetStripNumberingBorders(Double_t x_left, Double_t y_left, Double_t x_right, Double_t y_right) {
138 if(x_right < x_left) return false;
139
140 XLeftPointOfStripNumbering = x_left;
141 YLeftPointOfStripNumbering = y_left;
142 XRightPointOfStripNumbering = x_right;
143 YRightPointOfStripNumbering = y_right;
144
145 ResetLayer();
146
147 return true;
148}
149
151
152 if( left != LeftTop && left != LeftBottom ) {
153 if(Verbosity) cerr << "WARNING: SetStripNumberingBorders: left strip border point is incorrect\n";
154 return false;
155 }
156 if( right != RightTop && right != RightBottom ) {
157 if(Verbosity) cerr << "WARNING: SetStripNumberingBorders: right strip border point is incorrect\n";
158 return false;
159 }
160
161 XLeftPointOfStripNumbering = XMinLayer;
162 XRightPointOfStripNumbering = XMaxLayer;
163
164 switch (left) {
165 case LeftTop:
166 YLeftPointOfStripNumbering = YMaxLayer;
167 break;
168 case LeftBottom:
169 YLeftPointOfStripNumbering = YMinLayer;
170 break;
171 default:
172 return false;
173 }
174
175 switch (right) {
176 case RightTop:
177 YRightPointOfStripNumbering = YMaxLayer;
178 break;
179 case RightBottom:
180 YRightPointOfStripNumbering = YMinLayer;
181 break;
182 default:
183 return false;
184 }
185
186 ResetLayer();
187
188 return true;
189}
190
192 StripOrder = strip_direction;
193
194 ResetLayer();
195
196 return true;
197}
198
199Bool_t BmnCSCLayer::AddDeadZone(Int_t n_points, Double_t *x_points, Double_t *y_points) {
200 DeadZoneOfStripLayer dead_zone;
201 Bool_t status = dead_zone.SetDeadZone(n_points, x_points, y_points);
202 if(status == true) {
203 DeadZones.push_back(dead_zone);
204 return true;
205 }
206 else {
207 return false;
208 }
209}
210
212 if(dead_zone.GetNPoints() > 2) {
213 DeadZones.push_back(dead_zone);
214 return true;
215 }
216 else {
217 return false;
218 }
219}
220
221Bool_t BmnCSCLayer::IsPointInsideDeadZones(Double_t x, Double_t y) {
222 for(size_t izone = 0; izone < DeadZones.size(); ++izone) {
223 if(DeadZones[izone].IsInside(x,y)) return true;
224 }
225
226 return false;
227}
228
229Bool_t BmnCSCLayer::IsPointInsideStripLayer(Double_t x, Double_t y) {
230 if( x >= XMinLayer && x <= XMaxLayer &&
231 y >= YMinLayer && y <= YMaxLayer &&
232 !IsPointInsideDeadZones(x, y) ) return true;
233 else
234 return false;
235}
236
237Bool_t BmnCSCLayer::SetStripSignal(Int_t strip_num, Double_t signal) {
238 if (strip_num >= 0 && strip_num < (Int_t)Strips.size())
239 {
240 Strips[strip_num] = signal;
241 return true;
242 }
243 else return false;
244}
245
246Bool_t BmnCSCLayer::AddStripSignal(Int_t strip_num, Double_t signal) {
247 if(strip_num >= 0 && strip_num < (Int_t)Strips.size()) {
248 Strips[strip_num] += signal;
249 return true;
250 }
251 else return false;
252}
253
254Double_t BmnCSCLayer::GetStripSignal(Int_t strip_num) {
255 if(strip_num >= 0 && strip_num < (Int_t)Strips.size()) {
256 return Strips[strip_num];
257 }
258 else return -1;
259}
260
261Bool_t BmnCSCLayer::SetStripMatch(Int_t strip_num, BmnMatch strip_match) {
262 if(strip_num >= 0 && strip_num < (Int_t)StripMatches.size()) {
263 StripMatches[strip_num] = strip_match;
264 return true;
265 }
266 else return false;
267}
268
269Bool_t BmnCSCLayer::AddLinkToStripMatch(Int_t strip_num, Double_t weight, Int_t refID) {
270 if(strip_num >= 0 && strip_num < (Int_t)StripMatches.size()) {
271 StripMatches[strip_num].AddLink(weight, refID);
272 return true;
273 }
274 else return false;
275}
276
277Bool_t BmnCSCLayer::SetStripDigitNumberMatch(Int_t strip_num, BmnMatch digit_num_match) {
278 if(strip_num >= 0 && strip_num < (Int_t)StripDigitNumberMatches.size()) {
279 StripDigitNumberMatches[strip_num] = digit_num_match;
280 return true;
281 }
282 else return false;
283}
284
285Bool_t BmnCSCLayer::AddLinkToStripDigitNumberMatch(Int_t strip_num, Double_t weight, Int_t digit_num) {
286 if(strip_num >= 0 && strip_num < (Int_t)StripDigitNumberMatches.size()) {
287 StripDigitNumberMatches[strip_num].AddLink(weight, digit_num);
288 return true;
289 }
290 else return false;
291}
292
294 if(strip_num >= 0 && strip_num < (Int_t)StripMatches.size()) {
295 return StripMatches[strip_num];
296 }
297 else return BmnMatch(); //return an empty match
298}
299
301 if(strip_num >= 0 && strip_num < (Int_t)StripDigitNumberMatches.size()) {
302 return StripDigitNumberMatches[strip_num];
303 }
304 else return BmnMatch(); //return an empty match
305}
306
308 Int_t NStrips = GetNStrips();
309
310 StripMatches.clear();
311 StripMatches.resize(NStrips, BmnMatch());
312}
313
315 Int_t NStrips = GetNStrips();
316
317 StripDigitNumberMatches.clear();
318 StripDigitNumberMatches.resize(NStrips, BmnMatch());
319}
320
321Double_t BmnCSCLayer::GetStripHitPos(Int_t num) {
322 if(num >= 0 && num < (Int_t)StripHits.size()) {
323 return StripHits[num];
324 }
325 return -1.0;
326}
327
329 if(num >= 0 && num < (Int_t)StripHitsTotalSignal.size()) {
330 return StripHitsTotalSignal[num];
331 }
332 return -1.0;
333}
334
335Double_t BmnCSCLayer::GetStripHitError(Int_t num) {
336 if(num >= 0 && num < (Int_t)StripHitsErrors.size()) {
337 return StripHitsErrors[num];
338 }
339 return -1.0;
340}
341
343 if(num >= 0 && num < (Int_t)StripHitsClusterSize.size()) {
344 return StripHitsClusterSize[num];
345 }
346 return -1.0;
347}
348
350 StripHits.clear();
351 StripHitsTotalSignal.clear();
352 StripHitsErrors.clear();
353 StripHitsClusterSize.clear();
354 StripClusters.clear();
355}
356
357Double_t BmnCSCLayer::ConvertNormalPointToStripX(Double_t x, Double_t y) {
358 Double_t XRotationCenter;
359 Double_t YRotationCenter;
360
361 if(StripOrder == LeftToRight) {
362 XRotationCenter = XLeftPointOfStripNumbering;
363 YRotationCenter = YLeftPointOfStripNumbering;
364 }
365 else {
366 if(StripOrder == RightToLeft) {
367 XRotationCenter = XRightPointOfStripNumbering;
368 YRotationCenter = YRightPointOfStripNumbering;
369 }
370 }
371 Double_t StripX = (x - XRotationCenter)*Cos(-AngleRad) + (y - YRotationCenter)*Sin(-AngleRad) + XRotationCenter;//see
372 return StripX;
373}
374
375Double_t BmnCSCLayer::ConvertNormalPointToStripY(Double_t x, Double_t y) {
376 Double_t XRotationCenter;
377 Double_t YRotationCenter;
378
379 if(StripOrder == LeftToRight) {
380 XRotationCenter = XLeftPointOfStripNumbering;
381 YRotationCenter = YLeftPointOfStripNumbering;
382 }
383 else {
384 if(StripOrder == RightToLeft) {
385 XRotationCenter = XRightPointOfStripNumbering;
386 YRotationCenter = YRightPointOfStripNumbering;
387 }
388 }
389 Double_t StripY = -(x - XRotationCenter)*Sin(-AngleRad) + (y - YRotationCenter)*Cos(-AngleRad) + YRotationCenter;//see
390 return StripY;
391}
392
393Double_t BmnCSCLayer::ConvertPointToStripPosition(Double_t x, Double_t y) {
394 //This function returns real(double) value,
395 //where integer part - number of a strip, fractional part - position on this strip (as ratio from begin)
396
397 if(StripOrder == LeftToRight) {
398 return (ConvertNormalPointToStripX(x, y)-XLeftPointOfStripNumbering)/Pitch;
399 }
400 else {
401 if(StripOrder == RightToLeft) {
402 return (XRightPointOfStripNumbering-ConvertNormalPointToStripX(x, y))/Pitch;
403 }
404 }
405
406 return -1;
407}
408
409Double_t BmnCSCLayer::CalculateStripEquationB(Double_t strip_pos) {
410
411 Double_t b;
412
413 if(Abs(AngleDeg) != 90.0) { //case: y=a*x+b
414 Double_t x_strip_shift = (strip_pos*Pitch)/Cos(Abs(AngleRad));
415 Double_t xcoord;
416
417 if(StripOrder == LeftToRight) {
418 xcoord = XLeftPointOfStripNumbering + x_strip_shift;
419 b = YLeftPointOfStripNumbering - Tan(PiOver2()-AngleRad)*xcoord;
420 }
421 else {
422 if(StripOrder == RightToLeft) {
423 xcoord = XRightPointOfStripNumbering-x_strip_shift;
424 b = YRightPointOfStripNumbering - Tan(PiOver2()-AngleRad)*xcoord;
425 }
426 }
427 }
428 else { //case: y=b
429 Double_t y_strip_shift = strip_pos*Pitch*Sin(AngleRad);
430
431 if(StripOrder == LeftToRight) {
432 b = YLeftPointOfStripNumbering - y_strip_shift;
433 }
434 else {
435 if(StripOrder == RightToLeft) {
436 b = YRightPointOfStripNumbering + y_strip_shift;
437 }
438 }
439 }
440
441 return b;
442}
443
445
447
448 Double_t threshold = ClusterFindingThreshold;
449 //Double_t threshold = 500.0;
450
451 StripCluster cluster;
452
453 Bool_t ascent = false;
454 Bool_t descent = false;
455
456 //Processing strips
457 vector<Double_t> AnalyzableStrips = Strips;
458
459 //Smooth strip signal ------------------------------------------------------
460 //if(Pitch > 0.079) SmoothStripSignal(AnalyzableStrips, 1, 1, 1.0);
461 //else SmoothStripSignal(AnalyzableStrips, 2, 1, 1.0);
462 //--------------------------------------------------------------------------
463
464 for (Int_t is = 0; is < (Int_t)AnalyzableStrips.size(); is++)
465 {
466 if(AnalyzableStrips.at(is) <= threshold) {
467 if(descent || ascent) {
468 descent = false;
469 ascent = false;
470 //make strip hit
471 MakeStripHit(cluster, AnalyzableStrips, is);
472 }
473 continue;
474 }
475
476 if( cluster.GetClusterSize() > 0 ) {
477 if( AnalyzableStrips.at(is) >= (cluster.Signals.at(cluster.GetClusterSize()-1)) ) {
478 if(descent) {
479 ascent = false;
480 descent = false;
481 //make strip hit
482 MakeStripHit(cluster, AnalyzableStrips, is);
483 //continue;
484 }
485 ascent = true;
486 descent = false;
487 }
488 else {
489 ascent = false;
490 descent = true;
491 }
492 }
493 else {
494 ascent = true;
495 descent = false;
496 }
497
498 cluster.AddStrip(is, AnalyzableStrips.at(is));
499 }
500
501 if(cluster.GetClusterSize() != 0) {
502 //make strip hit
503 Int_t lastnum = AnalyzableStrips.size()-1;
504 MakeStripHit(cluster, AnalyzableStrips, lastnum);
505 }
506}
507
508 void BmnCSCLayer::MakeStripHit(StripCluster &cluster, vector<Double_t> &AnalyzableStrips, Int_t &curcnt) {
509
510 Double_t mean_strip_position = 0.0;
511 Double_t total_cluster_signal = 0.0;
512 Double_t cluster_rms = 0.0;
513
514 Int_t NStripsInCluster = cluster.GetClusterSize();
515
516 //Share the common strip signal of two adjacent clusters in half
517 if(true) {
518 if(AnalyzableStrips.at(curcnt) >= AnalyzableStrips.at(curcnt-1)) {
519 cluster.Signals.at(cluster.Signals.size()-1) *= 0.5;
520 AnalyzableStrips.at(curcnt-1) *= 0.5;
521 }
522 }
523
524 //mean position and cluster signal -----------------------------------------
525 for(Int_t i = 0; i < NStripsInCluster; ++i) {
526 Double_t strip_num = cluster.Strips.at(i);
527 Double_t signal = cluster.Signals.at(i);
528 total_cluster_signal += signal;
529 mean_strip_position += (strip_num+0.5)*signal;
530 }
531 mean_strip_position /= total_cluster_signal;
532
533 if(mean_strip_position < 0.0) mean_strip_position = 0.0;
534 if(mean_strip_position >= AnalyzableStrips.size()) mean_strip_position = AnalyzableStrips.size() - 0.001;
535 //--------------------------------------------------------------------------
536
537 //cluster standard deviation (sigma): RMS ----------------------------------
538// if(NStripsInCluster > 1) {
539// for(Int_t i = 0; i < NStripsInCluster; ++i) {
540// Double_t strip_num = cluster.Strips.at(i);
541// Double_t signal = cluster.Signals.at(i);
542// Double_t residual = (strip_num+0.5) - mean_strip_position;
543// cluster_rms += residual*residual*signal;
544// }
545// cluster_rms /= total_cluster_signal;
546// cluster_rms = TMath::Sqrt(cluster_rms);
547// }
548// else {
549// cluster_rms = 1.0/TMath::Sqrt(12.0);
550// }
551//
552
553 // AZ, STS, method2
554 Double_t sumW = total_cluster_signal;
555 Double_t sumWX = 0.;
556 Double_t sumWX2 = 0.;
557 if(NStripsInCluster > 1) {
558 for(Int_t i = 0; i < NStripsInCluster; ++i) {
559 Double_t strip_num = cluster.Strips.at(i);
560 Double_t signal = cluster.Signals.at(i);
561 sumWX += strip_num * signal;
562 sumWX2 += strip_num * strip_num * signal;
563 }
564 cluster_rms = (sumWX2 - sumWX * sumWX / sumW) / sumW;
565 }
566
567 else {
568 cluster_rms = 1.0 /TMath::Sqrt(12.0);
569 }
570
571// // AZ, STS, Real CF, method3
572// Double_t sumW = total_cluster_signal;
573// Double_t maxStripSignal = 0.;
574//
575// if(NStripsInCluster > 0) {
576// for(Int_t i = 0; i < NStripsInCluster; ++i) {
577// Double_t signal = cluster.Signals.at(i);
578// if (signal > maxStripSignal)
579// maxStripSignal = signal;
580//
581// }
582// cluster_rms = sumW / maxStripSignal;
583// }
584 //--------------------------------------------------------------------------
585
586 StripHits.push_back(mean_strip_position);
587 StripHitsTotalSignal.push_back(total_cluster_signal);
588 StripHitsErrors.push_back(cluster_rms);
589 StripHitsClusterSize.push_back(NStripsInCluster);
590
591 cluster.SetWidth(NStripsInCluster);
592 cluster.MeanPosition = mean_strip_position;
593 cluster.TotalSignal = total_cluster_signal;
594 cluster.Error = cluster_rms;
595 if (LayerType == LowerStripLayer) {
596 cluster.SetType(0);
597 cluster.SetUniqueID(fUniqueIdL++);
598 } else {
599 cluster.SetType(1);
600 cluster.SetUniqueID(fUniqueIdU++);
601 }
602
603 StripClusters.push_back(cluster);
604
605 //return to a previous strip
606 curcnt--;
607 if( curcnt < 0 ) curcnt = 0;
608
609 cluster.Clear();
610 }
611
612 void BmnCSCLayer::SmoothStripSignal(vector<Double_t>& AnalyzableStrips, Int_t NIterations, Int_t SmoothWindow, Double_t Weight) {
613 //It's the Simple Moving Average (SMA) method
614 //AnalyzableStrips - analyzable strip layer (ref)
615 //NIterations - number of smooth iterations (usually 1)
616 //SmoothWindow - number of strips on the left-right of the current strip (usually 1)
617 //Weight - weight of the current strip (usually 1.0 - for simplicity, greater - for weighted value))
618
619 SmoothStrips.clear();
620 Int_t NStrips = AnalyzableStrips.size();
621
622 for(Int_t iteration = 0; iteration < NIterations; ++iteration) {
623 SmoothStrips.clear();
624 for(Int_t istrip = 0; istrip < NStrips; ++istrip) {
625 Double_t mean_value = 0.0;
626 for(Int_t iw = istrip-SmoothWindow; iw <= istrip+SmoothWindow; ++iw) {
627 if(iw >= 0 && iw < NStrips) {
628 if(iw == istrip) mean_value += AnalyzableStrips[iw]*Weight;
629 else mean_value += AnalyzableStrips[iw];
630 }
631 }
632 mean_value /= 2.0*SmoothWindow + Weight;
633 SmoothStrips.push_back(mean_value);
634 }
635 AnalyzableStrips = SmoothStrips;
636 }
637
638 return;
639}
int i
Definition P4_F32vec4.h:22
StripBorderPoint
@ LeftTop
@ RightBottom
@ LeftBottom
@ RightTop
StripNumberingDirection
@ RightToLeft
@ LeftToRight
StripLayerType
@ LowerStripLayer
Double_t CalculateStripEquationB(Double_t strip_pos)
Bool_t SetStripMatch(Int_t strip_num, BmnMatch strip_match)
Int_t GetNStrips()
Double_t GetStripHitPos(Int_t num)
void MakeStripHit(StripCluster &cluster, vector< Double_t > &AnalyzableStrips, Int_t &curcnt)
void InitializeLayer()
Bool_t SetStripDigitNumberMatch(Int_t strip_num, BmnMatch digit_num_match)
void SetAngleDeg(Double_t deg)
Double_t ConvertNormalPointToStripX(Double_t x, Double_t y)
BmnMatch GetStripDigitNumberMatch(Int_t strip_num)
void SmoothStripSignal(vector< Double_t > &AnalyzableStrips, Int_t NIterations, Int_t SmoothWindow, Double_t Weight)
void SetLayerSizes(Double_t xsize, Double_t ysize, Double_t xorig=0.0, Double_t yorig=0.0)
Bool_t AddStripSignal(Int_t strip_num, Double_t signal)
Double_t GetStripHitError(Int_t num)
Double_t GetStripSignal(Int_t strip_num)
Double_t ConvertPointToStripPosition(Double_t x, Double_t y)
void ResetLayer()
Double_t ConvertNormalPointToStripY(Double_t x, Double_t y)
Bool_t AddDeadZone(Int_t n_points, Double_t *x_points, Double_t *y_points)
Bool_t SetStripNumberingOrder(StripNumberingDirection strip_direction)
Bool_t AddLinkToStripMatch(Int_t strip_num, Double_t weight, Int_t refID)
Int_t GetStripHitClusterSize(Int_t num)
Bool_t SetStripNumberingBorders(Double_t x_left, Double_t y_left, Double_t x_right, Double_t y_right)
Bool_t IsPointInsideDeadZones(Double_t x, Double_t y)
Bool_t IsPointInsideStripLayer(Double_t x, Double_t y)
BmnMatch GetStripMatch(Int_t strip_num)
void SetPitch(Double_t pitch)
Bool_t SetStripSignal(Int_t strip_num, Double_t signal)
void ResetStripDigitNumberMatches()
void FindClustersAndStripHits()
Double_t GetStripHitTotalSignal(Int_t num)
void ResetStripHits()
Bool_t AddLinkToStripDigitNumberMatch(Int_t strip_num, Double_t weight, Int_t digit_num)
void ResetStripMatches()
virtual ~BmnCSCLayer()
Bool_t SetDeadZone(Int_t n_points, Double_t *xpoints, Double_t *ypoints)
Double_t TotalSignal
Int_t GetClusterSize()
Double_t MeanPosition
Double_t Error
void SetWidth(Int_t w)
void AddStrip(Int_t strip_num, Double_t strip_signal)
void SetType(Int_t type)
vector< Int_t > Strips
vector< Double_t > Signals