ATLAS Offline Software
Loading...
Searching...
No Matches
TauShot.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include "tauEvent/TauShot.h"
7
8
9namespace Analysis
10{
11 // Default Constructor
13 : P4EEtaPhiM()
15 , m_cluster()
16 , m_seedCell()
17 , m_nCellsInEta(TauDetails::DEFAULT_INT)
18 , m_pt1(TauDetails::DEFAULT)
19 , m_pt3(TauDetails::DEFAULT)
20 , m_pt5(TauDetails::DEFAULT)
21 , m_ws5(TauDetails::DEFAULT)
24 , m_sdevPt5(TauDetails::DEFAULT)
25 , m_deltaPt12_min(TauDetails::DEFAULT)
26 , m_Fside_3not1(TauDetails::DEFAULT)
27 , m_Fside_5not1(TauDetails::DEFAULT)
28 , m_Fside_5not3(TauDetails::DEFAULT)
32 , m_pt1OverPt3(TauDetails::DEFAULT)
33 , m_pt3OverPt5(TauDetails::DEFAULT)
34 , m_mergedScore(TauDetails::DEFAULT)
35 , m_signalScore(TauDetails::DEFAULT)
36 , m_nPhotons(TauDetails::DEFAULT_INT)
37 {}
38
39
40 // Copy Constructor
73
74 // Assignment operator
76 {
77 if (this!=&rhs){
78 P4EEtaPhiM::operator= (*this);
79 m_cluster = rhs.m_cluster;
82 m_pt1 = rhs.m_pt1;
83 m_pt3 = rhs.m_pt3;
84 m_pt5 = rhs.m_pt5;
85 m_ws5 = rhs.m_ws5;
88 m_sdevPt5 = rhs.m_sdevPt5;
101 }
102 return *this;
103 }
104
105 // Comparison Operator
106 bool TauShot::operator== ( const TauShot& rhs ) const
107 {
108 if(m_cluster != rhs.m_cluster) return false;
109 if(m_seedCell != rhs.m_seedCell) return false;
110 if(m_nCellsInEta != rhs.m_nCellsInEta) return false;
111 if(m_pt1 != rhs.m_pt1) return false;
112 if(m_pt3 != rhs.m_pt3) return false;
113 if(m_pt5 != rhs.m_pt5) return false;
114 if(m_ws5 != rhs.m_ws5) return false;
115 if(m_sdevEta5_WRTmean != rhs.m_sdevEta5_WRTmean) return false;
116 if(m_sdevEta5_WRTmode != rhs.m_sdevEta5_WRTmode) return false;
117 if(m_sdevPt5 != rhs.m_sdevPt5) return false;
118 if(m_deltaPt12_min != rhs.m_deltaPt12_min) return false;
119 if(m_Fside_3not1 != rhs.m_Fside_3not1) return false;
120 if(m_Fside_5not1 != rhs.m_Fside_5not1) return false;
121 if(m_Fside_5not3 != rhs.m_Fside_5not3) return false;
122 if(m_fracSide_3not1 != rhs.m_fracSide_3not1) return false;
123 if(m_fracSide_5not1 != rhs.m_fracSide_5not1) return false;
124 if(m_fracSide_5not3 != rhs.m_fracSide_5not3) return false;
125 if(m_pt1OverPt3 != rhs.m_pt1OverPt3) return false;
126 if(m_pt3OverPt5 != rhs.m_pt3OverPt5) return false;
127 if(m_mergedScore != rhs.m_mergedScore) return false;
128 if(m_signalScore != rhs.m_signalScore) return false;
129 if(m_nPhotons != rhs.m_nPhotons) return false;
130 return true;
131 }
132
133 // Destructor
135 {
136 }
137
138
139 std::vector<std::vector<const CaloCell*> > TauShot::getCellVector(const CaloCell_ID* calo_id) const
140 {
141 std::vector<std::vector<const CaloCell*> > cellVector;
142
143
144 // have two layers in phi
145 cellVector.emplace_back(m_nCellsInEta, nullptr);
146 cellVector.emplace_back(m_nCellsInEta, nullptr);
147 const CaloCell* seedCell = this->seedCell();
148 const CaloCell* mergedCell = nullptr;
149 // get merged cell in phi. Keep nullptr if shot is not merged across phi
150 std::vector<IdentifierHash> nextInPhi;
151 std::vector<IdentifierHash> prevInPhi;
152 calo_id->get_neighbours(seedCell->caloDDE()->calo_hash(),LArNeighbours::nextInPhi,nextInPhi);
153 calo_id->get_neighbours(seedCell->caloDDE()->calo_hash(),LArNeighbours::prevInPhi,prevInPhi);
154 CaloCluster::cell_iterator cellItr = this->cluster()->cell_begin();
155 CaloCluster::cell_iterator cellItrE = this->cluster()->cell_end();
156 for(;cellItr!=cellItrE;++cellItr){
157 std::vector<IdentifierHash>::iterator itr = nextInPhi.begin();
158 for( ; itr!=nextInPhi.end(); ++itr ){
159 if((*cellItr)->caloDDE()->calo_hash() != (*itr)) continue;
160 mergedCell = (*cellItr);
161 break;
162 }
163 if(mergedCell!=nullptr) break;
164 itr = prevInPhi.begin();
165 for( ; itr!=prevInPhi.end(); ++itr ){
166 if((*cellItr)->caloDDE()->calo_hash() != (*itr)) continue;
167 mergedCell = (*cellItr);
168 break;
169 }
170 if(mergedCell!=nullptr) break;
171 }
172 // store cells in the eta layer, which contains the seed cell
173 int nCellsFromSeed = 1;
174 const CaloCell* lastCell = seedCell;
175 cellVector.at(0).at(m_nCellsInEta/2) = seedCell; // store seed cell
176 std::vector<IdentifierHash> next;
177 while(lastCell!=nullptr && nCellsFromSeed<m_nCellsInEta/2+1){
178 calo_id->get_neighbours(lastCell->caloDDE()->calo_hash(),LArNeighbours::nextInEta,next);
179 lastCell = nullptr;
180 for(cellItr=this->cluster()->cell_begin();cellItr!=cellItrE;++cellItr){
181 std::vector<IdentifierHash>::iterator itr = next.begin();
182 for( ; itr!=next.end(); ++itr ){
183 if((*cellItr)->caloDDE()->calo_hash() != (*itr)) continue;
184 cellVector.at(0).at(m_nCellsInEta/2+nCellsFromSeed) = (*cellItr);
185 lastCell = (*cellItr);
186 }
187 }
188 nCellsFromSeed++;
189 }
190 nCellsFromSeed = 1;
191 lastCell = seedCell;
192 while(lastCell!=nullptr && nCellsFromSeed<m_nCellsInEta/2+1){
193 calo_id->get_neighbours(lastCell->caloDDE()->calo_hash(),LArNeighbours::prevInEta,next);
194 lastCell = nullptr;
195 for(cellItr=this->cluster()->cell_begin();cellItr!=cellItrE;++cellItr){
196 std::vector<IdentifierHash>::iterator itr = next.begin();
197 for( ; itr!=next.end(); ++itr ){
198 if((*cellItr)->caloDDE()->calo_hash() != (*itr)) continue;
199 cellVector.at(0).at(m_nCellsInEta/2-nCellsFromSeed) = (*cellItr);
200 lastCell = (*cellItr);
201 }
202 }
203 nCellsFromSeed++;
204 }
205 // store cells in the eta layer, which contains the merged cell
206 int nCellsFromMerged = 1;
207 lastCell = mergedCell; // is nullptr if shot is not merged
208 cellVector.at(1).at(m_nCellsInEta/2) = mergedCell; // store merged cell
209 while(lastCell!=nullptr && nCellsFromMerged<m_nCellsInEta/2+1){
210 calo_id->get_neighbours(lastCell->caloDDE()->calo_hash(),LArNeighbours::nextInEta,next);
211 lastCell = nullptr;
212 for(cellItr=this->cluster()->cell_begin();cellItr!=cellItrE;++cellItr){
213 std::vector<IdentifierHash>::iterator itr = next.begin();
214 for( ; itr!=next.end(); ++itr ){
215 if((*cellItr)->caloDDE()->calo_hash() != (*itr)) continue;
216 cellVector.at(1).at(m_nCellsInEta/2+nCellsFromMerged) = (*cellItr);
217 lastCell = (*cellItr);
218 }
219 }
220 nCellsFromMerged++;
221 }
222 nCellsFromMerged = 1;
223 lastCell = mergedCell;
224 while(lastCell!=nullptr && nCellsFromMerged<m_nCellsInEta/2+1){
225 calo_id->get_neighbours(lastCell->caloDDE()->calo_hash(),LArNeighbours::prevInEta,next);
226 lastCell = nullptr;
227 for(cellItr=this->cluster()->cell_begin();cellItr!=cellItrE;++cellItr){
228 std::vector<IdentifierHash>::iterator itr = next.begin();
229 for( ; itr!=next.end(); ++itr ){
230 if((*cellItr)->caloDDE()->calo_hash() != (*itr)) continue;
231 cellVector.at(1).at(m_nCellsInEta/2-nCellsFromMerged) = (*cellItr);
232 lastCell = (*cellItr);
233 }
234 }
235 nCellsFromMerged++;
236 }
237 return cellVector;
238 }
239
240
241 // just for testing
242 void TauShot::print() const{
243 int oldpr = std::cout.precision(5);
244 // can probably do this better through IMomentum::dump(), but whatevs
245 std::cout << "in TauShot::dump()" << std::endl;
246 std::cout << "pt: " << this->pt() << std::endl;
247 std::cout << "eta: " << this->eta() << std::endl;
248 std::cout << "phi: " << this->phi() << std::endl;
249 std::cout << "m: " << this->m() << std::endl;
250 std::cout.precision(oldpr);
251 }
252
253}//end of namespace definitions
254
Declaration of tau details base class.
void print() const
print method
Definition TauShot.cxx:242
float m_deltaPt12_min
Definition TauShot.h:159
const CaloCluster * cluster() const
pointer to CaloCluster
float m_fracSide_5not1
Definition TauShot.h:164
ElementLink< CaloClusterContainer > m_cluster
element link to cluster
Definition TauShot.h:145
int m_nCellsInEta
other shot variables
Definition TauShot.h:151
float m_fracSide_3not1
Definition TauShot.h:163
TauShot & operator=(const TauShot &)
Assignment operator.
Definition TauShot.cxx:75
float m_Fside_5not3
Definition TauShot.h:162
const CaloCell * seedCell() const
pointer to seed cell
std::vector< std::vector< const CaloCell * > > getCellVector(const CaloCell_ID *) const
get vector<vector<CaloCell*> > used for calculation of shape vars
Definition TauShot.cxx:139
TauShot()
Default Constructor standard constructor which sets everything to 0, needed for persistency.
Definition TauShot.cxx:12
float m_signalScore
Definition TauShot.h:169
float m_Fside_3not1
Definition TauShot.h:160
float m_sdevEta5_WRTmean
Definition TauShot.h:156
ElementLink< CaloCellContainer > m_seedCell
pointer to seed cell
Definition TauShot.h:148
float m_Fside_5not1
Definition TauShot.h:161
float m_mergedScore
Definition TauShot.h:168
virtual ~TauShot()
Destructor.
Definition TauShot.cxx:134
float m_sdevEta5_WRTmode
Definition TauShot.h:157
float m_fracSide_5not3
Definition TauShot.h:165
bool operator==(const TauShot &) const
Comparison operator.
Definition TauShot.cxx:106
int get_neighbours(const IdentifierHash caloHash, const LArNeighbours::neighbourOption &option, std::vector< IdentifierHash > &neighbourList) const
access to hashes for neighbours return == 0 for neighbours found
Helper class for offline cell identifiers.
Definition CaloCell_ID.h:34
Data object for each calorimeter readout cell.
Definition CaloCell.h:57
const CaloDetDescrElement * caloDDE() const
get pointer to CaloDetDescrElement (data member)
Definition CaloCell.h:321
CaloCompositeCellBase< CaloClusterNavigable >::cell_iterator cell_iterator
Iterator on CaloCell s.
cell_iterator cell_end() const
Retrieve a STL-type end() iterator for the cell store.
cell_iterator cell_begin() const
Retrieve a STL-type begin() iterator for the cell store.
I4Momentum is an abstract base class providing 4-momentum behavior.
Definition I4Momentum.h:31
virtual double pt() const =0
transverse momentum
virtual double phi() const
get phi data member
Definition P4EEtaPhiM.h:108
P4EEtaPhiM(const double e, const double eta, const double phi, const double m)
constructor with all data members
Definition P4EEtaPhiM.cxx:7
virtual double eta() const
get eta data member
Definition P4EEtaPhiM.h:105
virtual double m() const
get mass data member
Definition P4EEtaPhiM.h:111
The namespace of all packages in PhysicsAnalysis/JetTagging.