GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TSiLiHit.cxx
Go to the documentation of this file.
1#include "TSiLi.h"
2#include "TSiLiHit.h"
3
8
10 : fFitCharge(frag.GetCharge())
11{
12 frag.Copy(*this);
13 SetWavefit(frag);
14}
15
17{
18 Clear();
19 rhs.Copy(*this);
20}
21
22void TSiLiHit::Copy(TObject& rhs, bool suppress) const
23{
25
26 static_cast<TSiLiHit&>(rhs).fTimeFit = fTimeFit;
27 static_cast<TSiLiHit&>(rhs).fSig2Noise = fSig2Noise;
28 static_cast<TSiLiHit&>(rhs).fSmirnov = fSmirnov;
29 static_cast<TSiLiHit&>(rhs).fFitCharge = fFitCharge;
30 static_cast<TSiLiHit&>(rhs).fFitBase = fFitBase;
31 static_cast<TSiLiHit&>(rhs).fSiLiHitBits = 0;
32 if(!suppress) {
33 static_cast<TSiLiHit&>(rhs).fAddBackSegments = fAddBackSegments;
34 static_cast<TSiLiHit&>(rhs).fAddBackEnergy = fAddBackEnergy;
35 }
36}
37
38void TSiLiHit::Clear(Option_t* opt)
39{
41 // fSegment = -1;
42 fTimeFit = -1;
43 fFitCharge = -1;
44 fFitBase = -1;
45 fSig2Noise = -1;
46 fSmirnov = -1;
47
48 fAddBackSegments.clear();
49 fAddBackEnergy.clear();
52}
53
55{
57 if(pulse != nullptr) {
58 fTimeFit = pulse->Get_wpar_T0();
61 fSig2Noise = pulse->get_sig2noise();
62 fSmirnov = pulse->GetsiliSmirnov();
63 delete pulse;
64 }
65}
66
67// Broken up for external analysis script use
68TPulseAnalyzer* TSiLiHit::FitFrag(const TFragment& frag, int ShapeFit, int segment)
69{
70 return FitFrag(frag, ShapeFit, GetSiLiHitChannel(segment));
71}
72
73TPulseAnalyzer* TSiLiHit::FitFrag(const TFragment& frag, int ShapeFit, TChannel* channel)
74{
75 auto* pulse = new TPulseAnalyzer(frag, TSiLi::fSiLiNoiseFac);
76 if(FitPulseAnalyzer(pulse, ShapeFit, channel) != 0) {
77 return pulse;
78 }
79 delete pulse;
80 return nullptr;
81}
82
84{
85 std::stringstream ss;
86 ss << "SPI00XN" << std::uppercase << std::uppercase << std::setfill('0') << std::setw(2) << std::hex << segment;
87 // std::cout<<std::endl<<ss.str().c_str();
88 return TChannel::FindChannelByName(ss.str().c_str());
89}
90
91int TSiLiHit::FitPulseAnalyzer(TPulseAnalyzer* pulse, int ShapeFit, int segment)
92{
93 return FitPulseAnalyzer(pulse, ShapeFit, GetSiLiHitChannel(segment));
94}
95
96int TSiLiHit::FitPulseAnalyzer(TPulseAnalyzer* pulse, int ShapeFit, TChannel* channel)
97{
98 if(pulse == nullptr) { return 0; }
99 if(pulse->IsSet()) {
100 double Decay = 0.;
101 double Rise = 0.;
102 double Base = 0.;
103
104 if(channel != nullptr) {
105 if(channel->UseWaveParam()) {
106 Rise = channel->GetWaveRise();
107 Decay = channel->GetWaveDecay();
108 Base = channel->GetWaveBaseLine();
109 }
110 }
111 // std::cout<<std::endl<<Decay<<" "<<Rise<<" "<<Base;
112
113 if(Decay == 0.) { Decay = TSiLi::fSiLiDefaultDecay; }
114 if(Rise == 0.) { Rise = TSiLi::fSiLiDefaultRise; }
115 if(Base == 0.) { Base = TSiLi::fSiLiDefaultBaseline; }
116
117 bool goodfit = false;
118 if(ShapeFit < 2) { goodfit = pulse->GetSiliShape(Decay, Rise); }
119 if(ShapeFit == 1 && !goodfit) { ShapeFit++; } //So currently it does a TF1 fit if initial fit fails, this might be a bad idea
120 if(ShapeFit == 2) { goodfit = pulse->GetSiliShapeTF1(Decay, Rise, Base); }
121 if(ShapeFit == 3) { goodfit = pulse->GetSiliShapeTF1(Decay, Rise, Base, TSiLi::fBaseFreq); }
122 if(goodfit) { return 1 + ShapeFit; }
123 }
124 return 0;
125}
126
127TVector3 TSiLiHit::GetPosition(Double_t, bool smear) const
128{
129 return TSiLi::GetPosition(GetRing(), GetSector(), smear);
130}
131
132TVector3 TSiLiHit::GetPosition(bool smear) const
133{
134 return GetPosition(GetDefaultDistance(), smear);
135}
136
137void TSiLiHit::Print(Option_t*) const
138{
139 Print(std::cout);
140}
141
142void TSiLiHit::Print(std::ostream& out) const
143{
144 std::ostringstream str;
145 str << "===============" << std::endl;
146 str << "not yet written" << std::endl;
147 str << "===============" << std::endl;
148 out << str.str();
149}
150
152{
153 if(this == hit) {
154 return;
155 }
156
157 if(fAddBackSegments.empty()) {
158 hit->Copy(*this, true); // suppresses copying of addback
159 fAddBackSegments.clear();
160 fAddBackEnergy.clear();
161 SetEnergy(0);
163 }
164
165 SetEnergy(GetEnergy() + hit->GetEnergy());
166 fAddBackSegments.push_back(hit->GetSegment());
167 fAddBackEnergy.push_back(hit->GetEnergy());
168}
169
170Int_t TSiLiHit::GetRing() const
171{
172 return TSiLi::GetRing(GetSegment());
173}
175{
177}
179{
181}
182Int_t TSiLiHit::GetPin() const
183{
184 return TSiLi::GetPin(GetSegment());
185}
187{
189}
190
192{
195 }
196 TChannel* chan = GetChannel();
197 if(chan == nullptr) {
198 return fFitCharge;
199 }
200 return chan->CalibrateENG(fFitCharge, 0);
201}
202
203double TSiLiHit::GetEnergy(Option_t*) const
204{
206 return TDetectorHit::GetEnergy(); // If not fitting waveforms, be normal.
207 }
208 TChannel* chan = GetChannel();
209 if(chan == nullptr) {
210 return SetEnergy(fFitCharge);
211 }
212
213 return SetEnergy(chan->CalibrateENG(fFitCharge, 0)); // this will use the integration value
214} // in the TChannel if it exists.
double CalibrateENG(double) const
Definition TChannel.cxx:646
double GetWaveDecay() const
Definition TChannel.h:277
static TChannel * FindChannelByName(const char *ccName)
Definition TChannel.cxx:494
double GetWaveRise() const
Definition TChannel.h:276
double GetWaveBaseLine() const
Definition TChannel.h:278
bool UseWaveParam() const
Definition TChannel.h:279
Double_t SetEnergy(const double &energy) const
virtual double GetEnergy(Option_t *opt="") const
virtual void ClearTransients() const
TChannel * GetChannel() const
!
bool TestHitBit(EBitFlag flag) const
void Clear(Option_t *opt="") override
!
virtual Int_t GetSegment() const
!
void Copy(TObject &) const override
!
void SetHitBit(EBitFlag, Bool_t set=true) const
bool GetSiliShapeTF1(double tauDecay, double tauRise, double baseline, double basefreq=0)
bool IsSet() const
bool GetSiliShape(double tauDecay, double tauRise)
double Get_wpar_amplitude()
double Get_wpar_baselinefin()
std::vector< int16_t > fAddBackSegments
!
Definition TSiLiHit.h:133
Double_t fTimeFit
Definition TSiLiHit.h:138
static TChannel * GetSiLiHitChannel(int segment)
Definition TSiLiHit.cxx:83
double GetEnergy(Option_t *opt=nullptr) const override
Definition TSiLiHit.cxx:203
void Copy(TObject &, bool=false) const override
!
Definition TSiLiHit.cxx:22
static TPulseAnalyzer * FitFrag(const TFragment &frag, int ShapeFit, int segment)
Definition TSiLiHit.cxx:68
Int_t GetPreamp() const
Definition TSiLiHit.cxx:178
double GetFitEnergy() const
Definition TSiLiHit.cxx:191
Double_t fFitCharge
Definition TSiLiHit.h:141
Double_t fFitBase
Definition TSiLiHit.h:142
Int_t GetSector() const
Definition TSiLiHit.cxx:174
Int_t GetRing() const
Definition TSiLiHit.cxx:170
static int FitPulseAnalyzer(TPulseAnalyzer *pulse, int ShapeFit, int segment)
Definition TSiLiHit.cxx:91
void Clear(Option_t *opt="") override
Definition TSiLiHit.cxx:38
Double_t fSmirnov
Definition TSiLiHit.h:140
TTransientBits< UChar_t > fSiLiHitBits
Definition TSiLiHit.h:136
void SetWavefit(const TFragment &)
Definition TSiLiHit.cxx:54
bool MagnetShadow() const
Definition TSiLiHit.cxx:186
Int_t GetPin() const
Definition TSiLiHit.cxx:182
std::vector< double > fAddBackEnergy
!
Definition TSiLiHit.h:134
void SumHit(TSiLiHit *)
Definition TSiLiHit.cxx:151
Double_t fSig2Noise
Definition TSiLiHit.h:139
void Print(Option_t *opt="") const override
Definition TSiLiHit.cxx:137
Double_t GetDefaultDistance() const
Definition TSiLiHit.h:131
TVector3 GetPosition() const override
!
Definition TSiLiHit.h:77
static Int_t GetPreamp(Int_t seg)
Definition TSiLi.h:88
static int fFitSiLiShape
!
Definition TSiLi.h:115
static double fSiLiDefaultDecay
Definition TSiLi.h:82
static double fBaseFreq
!
Definition TSiLi.h:116
static Int_t GetPin(Int_t seg)
Definition TSiLi.h:89
static double fSiLiDefaultBaseline
Definition TSiLi.h:84
static double fSiLiNoiseFac
Definition TSiLi.h:81
static double fSiLiDefaultRise
Definition TSiLi.h:83
static Int_t GetRing(Int_t seg)
Definition TSiLi.h:86
static bool MagnetShadow(Int_t seg)
Definition TSiLi.h:102
static Int_t GetSector(Int_t seg)
Definition TSiLi.h:87
static TVector3 GetPosition(int ring, int sector, bool smear=false)
Definition TSiLi.cxx:77
Bool_t TestBit(T bit) const