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