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