ATLAS Offline Software
Loading...
Searching...
No Matches
TBPlaneTrackingAlgo Class Reference

#include <TBPlaneTrackingAlgo.h>

Inheritance diagram for TBPlaneTrackingAlgo:
Collaboration diagram for TBPlaneTrackingAlgo:

Public Member Functions

 TBPlaneTrackingAlgo (const std::string &name, ISvcLocator *pSvcLocator)
virtual ~TBPlaneTrackingAlgo ()
virtual StatusCode initialize () override
virtual StatusCode execute () override
virtual StatusCode finalize () override
virtual StatusCode sysInitialize () override
 Override sysInitialize.
virtual const DataObjIDColl & extraOutputDeps () const override
 Return the list of extra output dependencies.
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const

Protected Member Functions

void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

void FillRandomHit ()
StatusCode getnewcalib ()
bool fitPlane (const TBHitPlaneCont *hitPlaneCont, double &a1, double &a2, double &chi2, std::vector< double > &residual)
 Fit data to the function u = a1 + a2*w.
bool fitHits (const std::vector< double > &u, const std::vector< double > &w, const std::vector< double > &eu, double &a1, double &a2)
 Fit data to the function u = a1 + a2*w.
double getResidual (const double &u, const double &w, const double &a1, const double &a2)
 Calculates the residual-> r = (u_i - a1 - a2*w_i)
double getChi2 (const std::vector< double > &v_u, const std::vector< double > &v_w, const std::vector< double > &v_eu, const double &a1, const double &a2)
 Calculates the chi2 += ((u_i - a1 - a2*w_i)/eu)^2.
StatusCode buildHits ()
 Build HitPlaneCont from BPC.
Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

bool m_testAlgo
TBHitPlaneCont m_hitPlaneCont_u
TBHitPlaneCont m_hitPlaneCont_v
std::vector< float > m_bpc_posZX
std::vector< float > m_bpc_posZY
std::vector< float > m_bpc_posX
std::vector< float > m_bpc_posY
std::vector< float > m_bpc_errposZX
std::vector< float > m_bpc_errposZY
std::vector< float > m_bpc_errposX
std::vector< float > m_bpc_errposY
std::vector< float > m_bpc_errmeasX
std::vector< float > m_bpc_errmeasY
std::vector< std::string > m_bpc_names
std::string m_calib_filename
unsigned int m_runnumber
DataObjIDColl m_extendedExtraObjects
StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Detailed Description

Definition at line 20 of file TBPlaneTrackingAlgo.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< Algorithm > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ TBPlaneTrackingAlgo()

TBPlaneTrackingAlgo::TBPlaneTrackingAlgo ( const std::string & name,
ISvcLocator * pSvcLocator )

Definition at line 21 of file TBPlaneTrackingAlgo.cxx.

22 : AthAlgorithm(name, pSvcLocator),
24{
25 declareProperty("TestAlgo",m_testAlgo=false);
26 declareProperty("CalibFileName", m_calib_filename="H6HitCalib.txt");
27 declareProperty("BPCnames", m_bpc_names);
28 declareProperty("BPCposX", m_bpc_posX);
29 declareProperty("BPCposY", m_bpc_posY);
30 declareProperty("BPCposZX", m_bpc_posZX);
31 declareProperty("BPCposZY", m_bpc_posZY);
32 declareProperty("BPCerrposX", m_bpc_errposX);
33 declareProperty("BPCerrposY", m_bpc_errposY);
34 declareProperty("BPCerrposZX", m_bpc_errposZX);
35 declareProperty("BPCerrposZY", m_bpc_errposZY);
36 declareProperty("BPCerrmeasX", m_bpc_errmeasX);
37 declareProperty("BPCerrmeasY", m_bpc_errmeasY);
38
39}
AthAlgorithm()
Default constructor:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
std::vector< float > m_bpc_posX
std::vector< float > m_bpc_errmeasX
std::vector< float > m_bpc_errposY
std::vector< float > m_bpc_posZX
std::vector< float > m_bpc_posY
std::vector< std::string > m_bpc_names
std::vector< float > m_bpc_errposZY
std::vector< float > m_bpc_errposZX
std::vector< float > m_bpc_errmeasY
std::vector< float > m_bpc_posZY
std::vector< float > m_bpc_errposX

◆ ~TBPlaneTrackingAlgo()

virtual TBPlaneTrackingAlgo::~TBPlaneTrackingAlgo ( )
inlinevirtual

Definition at line 27 of file TBPlaneTrackingAlgo.h.

27{};

Member Function Documentation

◆ buildHits()

StatusCode TBPlaneTrackingAlgo::buildHits ( )
private

Build HitPlaneCont from BPC.

Definition at line 356 of file TBPlaneTrackingAlgo.cxx.

358{
359 ATH_MSG_DEBUG (" In buildHits() ");
360 m_hitPlaneCont_u.clear();
361 m_hitPlaneCont_v.clear();
362
363 TBBPCCont * bpcCont = nullptr;
364 ATH_CHECK( evtStore()->retrieve(bpcCont, "BPCCont") );
365
366 // Loop over BPC
367 for (const TBBPC* bpc : *bpcCont) {
368 std::string name = bpc->getDetectorName();
369 ATH_MSG_DEBUG (" Hits in BPC "<< name);
370 // Find calibration index for this BPC
371 unsigned int ind=0;
372 while(ind<m_bpc_names.size())
373 {
374 if(name==m_bpc_names[ind]) break;
375 else ind++;
376 }
377 if(ind==m_bpc_names.size()){
378 ATH_MSG_ERROR ("No calibrations for BPC" << name);
379 continue;
380 }
381
382 ATH_MSG_DEBUG (" BPC number "<< ind);
383
384
385 if(!(bpc->isXPosOverflow()||bpc->isXPulseOverflow())){
386 TBHitPlane * hitu = new TBHitPlane();
387 float w = m_bpc_posZX[ind];
388 // float u =m_bpc_posX[ind] + m_bpc_invX[ind]*bpc->getXPos();
389 float u =m_bpc_posX[ind] + bpc->getXPos();
390 ATH_MSG_DEBUG ("BPC" << ind << "X =" << u);
391 hitu->setPosu(u); hitu->setPosw(w);
392 float uerr=sqrt(m_bpc_errposX[ind]*m_bpc_errposX[ind]+m_bpc_errmeasX[ind]*m_bpc_errmeasX[ind]);
393 hitu->setErroru(uerr); hitu->setErrorw(m_bpc_errposZX[ind]);
394 hitu->setValidFlag(bpc->isXPosOverflow()&&bpc->isYPosOverflow());
395 m_hitPlaneCont_u.push_back(hitu);
396 }
397
398 if(!(bpc->isYPosOverflow()||bpc->isYPulseOverflow())){
399 TBHitPlane * hitv = new TBHitPlane();
400 float w = m_bpc_posZY[ind];
401 // float v =m_bpc_posY[ind] + m_bpc_invY[ind]*bpc->getYPos();
402 float v =m_bpc_posY[ind] + bpc->getYPos();
403 ATH_MSG_DEBUG ("BPC" << ind << "Y =" << v);
404 hitv->setPosu(v); hitv->setPosw(w);
405 float verr=sqrt(m_bpc_errposY[ind]*m_bpc_errposY[ind]+m_bpc_errmeasY[ind]*m_bpc_errmeasY[ind]);
406 hitv->setErroru(verr); hitv->setErrorw(m_bpc_errposZY[ind]);
407 hitv->setValidFlag(bpc->isXPosOverflow()&&bpc->isYPosOverflow());
408 m_hitPlaneCont_v.push_back(hitv);
409 }
410 }
411
412 return StatusCode::SUCCESS;
413}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
void setPosw(double w)
Definition TBHitPlane.h:41
void setValidFlag(bool valid)
Definition TBHitPlane.h:45
void setErroru(double eu)
Definition TBHitPlane.h:42
void setPosu(double u)
Definition TBHitPlane.h:40
void setErrorw(double ew)
Definition TBHitPlane.h:43
TBHitPlaneCont m_hitPlaneCont_u
TBHitPlaneCont m_hitPlaneCont_v
@ u
Enums for curvilinear frames.
Definition ParamDefs.h:77
retrieve(aClass, aKey=None)
Definition PyKernel.py:110

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Algorithm > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Algorithm > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< Algorithm > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< Algorithm > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ execute()

StatusCode TBPlaneTrackingAlgo::execute ( )
overridevirtual

Definition at line 49 of file TBPlaneTrackingAlgo.cxx.

51{
52 ATH_MSG_DEBUG ("Executing TBPlaneTracking algorithm");
53 const EventContext& ctx = Gaudi::Hive::currentContext();
54
55 if(m_testAlgo == true){
56 ATH_MSG_WARNING ("TBPlaneTrackingAlgo: " <<
57 "Using FillRandomHit method to get data");
59 }
60
61 // Check if this is not a physics event and exit OK
62 // retrieve Event Info
63 const TBEventInfo* theEventInfo;
64 StatusCode checkOut = evtStore()->retrieve(theEventInfo,"TBEventInfo");
65 if ( checkOut.isFailure() )
66 {
67 ATH_MSG_ERROR ("cannot retrieve TBEventInfo from StoreGate");
68 return StatusCode::FAILURE;
69 }
70 int evtType = theEventInfo->getEventType();
71 ATH_MSG_DEBUG ("Event Type found " << evtType);
72 if (evtType != 1) return StatusCode::SUCCESS;
73
74 StatusCode sc1;
75
76 // Get run number and get new calib constants -----------------------------
77 unsigned int thisrun=ctx.eventID().run_number();
78
79 if(thisrun != m_runnumber)
80 {
81 m_runnumber= thisrun;
83 }
84 // ------------------------------------------------------------------------
85
86
87
88 sc1 = buildHits();
89
90// sc1 = m_StoreGate->retrieve(m_hitPlaneCont_u,"HitPlaneCont_U");
91// sc2 = m_StoreGate->retrieve(m_hitPlaneCont_v,"HitPlaneCont_V");
92
93 if(sc1.isFailure()){
94 ATH_MSG_WARNING ("TBPlaneTrackingAlgo: Retrieval of HitPlane failed");
95 setFilterPassed(false);
96 return StatusCode::SUCCESS;
97 } else {
98
99 int hitPlaneNumU, hitPlaneNumV;
100 hitPlaneNumU = m_hitPlaneCont_u.size();
101 hitPlaneNumV = m_hitPlaneCont_v.size();
102
103 if( (hitPlaneNumU < 2) || (hitPlaneNumV < 2) ){
104 ATH_MSG_WARNING ("TBPlaneTrackingAlgo: "
105 << "Not enough hits in one or both planes, "
106 << "Cannot make track.");
107
108 setFilterPassed(false);
109 return StatusCode::SUCCESS;
110 }
111
112 // Setting the chi2's for the track //
113 double chi2_u = 0;
114 double chi2_v = 0;
115 std::vector<double> residual_u, residual_v;
116 double a1_u = 0;
117 double a2_u = 0;
118 double a1_v = 0;
119 double a2_v = 0;
120 bool check = true;
121
122 check = fitPlane(&m_hitPlaneCont_u, a1_u, a2_u, chi2_u, residual_u);
123 if(check == false){
124 ATH_MSG_ERROR ("TBPlaneTrackingAlgo: " << "Fit failure.");
125 setFilterPassed(false);
126 return StatusCode::SUCCESS;
127 }
128
129 check = fitPlane(&m_hitPlaneCont_v, a1_v, a2_v, chi2_v, residual_v);
130 if(check == false){
131 ATH_MSG_ERROR ("TBPlaneTrackingAlgo: " << "Fit failure.");
132 setFilterPassed(false);
133 return StatusCode::SUCCESS;
134 }
135
136 // Setting the slopes & intercepts for each track //
137 ATH_MSG_DEBUG ("Setting fit parameters of track.");
138
139 TBTrack *track = new TBTrack(hitPlaneNumU, hitPlaneNumV);
140 track->setUintercept(a1_u);
141 track->setVintercept(a1_v);
142 track->setUslope(a2_u);
143 track->setVslope(a2_v);
144
145 // Setting the residual for plane u //
146 for(int i = 0; i < hitPlaneNumU; i++){
147 track->setResidualu(i, residual_u[i]);
148 }
149
150 // Setting the residual for plane v //
151 for(int i = 0; i < hitPlaneNumV; i++){
152 track->setResidualv(i, residual_v[i]);
153 }
154
155 ATH_MSG_DEBUG ("chi2 in u: " << chi2_u);
156 ATH_MSG_DEBUG ("chi2 in v: " << chi2_v);
157 ATH_MSG_DEBUG ("Setting chi2s of track.");
158
159 track->setChi2_u(chi2_u);
160 track->setChi2_v(chi2_v);
161
162 // Record track info. into storeGate //
163 StatusCode sc = evtStore()->record(track,"Track");
164
165 if ( sc.isFailure( ) ) {
166 ATH_MSG_FATAL ("Cannot record Track");
167 }
168 }
169
170 setFilterPassed(true);
171 return StatusCode::SUCCESS;
172}
#define ATH_MSG_FATAL(x)
#define ATH_MSG_WARNING(x)
static Double_t sc
int getEventType() const
Definition TBEventInfo.h:71
StatusCode buildHits()
Build HitPlaneCont from BPC.
bool fitPlane(const TBHitPlaneCont *hitPlaneCont, double &a1, double &a2, double &chi2, std::vector< double > &residual)
Fit data to the function u = a1 + a2*w.
::StatusCode StatusCode
StatusCode definition for legacy code.

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ extraOutputDeps()

const DataObjIDColl & AthAlgorithm::extraOutputDeps ( ) const
overridevirtualinherited

Return the list of extra output dependencies.

This list is extended to include symlinks implied by inheritance relations.

Definition at line 50 of file AthAlgorithm.cxx.

51{
52 // If we didn't find any symlinks to add, just return the collection
53 // from the base class. Otherwise, return the extended collection.
54 if (!m_extendedExtraObjects.empty()) {
56 }
57 return Algorithm::extraOutputDeps();
58}
DataObjIDColl m_extendedExtraObjects

◆ FillRandomHit()

void TBPlaneTrackingAlgo::FillRandomHit ( )
private

Definition at line 315 of file TBPlaneTrackingAlgo.cxx.

317{
318
319 ATH_MSG_DEBUG ("Starting FillRandom");
320
321 TBHitPlaneCont * hitPlaneCont_U = new TBHitPlaneCont();
322 TBHitPlaneCont * hitPlaneCont_V = new TBHitPlaneCont();
323
324 TBHitPlane *hitPaneU1 = new TBHitPlane();
325 hitPaneU1->setPosu(1); hitPaneU1->setPosw(1);
326 hitPlaneCont_U->push_back(hitPaneU1);
327
328
329 TBHitPlane *hitPaneU2 = new TBHitPlane();
330 hitPaneU2->setPosu(2); hitPaneU2->setPosw(2);
331 hitPlaneCont_U->push_back(hitPaneU2);
332
333 StatusCode sc1 = evtStore()->record(hitPlaneCont_U,"HitPlaneCont_U");
334 if ( sc1.isFailure( ) ) {
335 ATH_MSG_FATAL ("Cannot record HitPlaneCont_U");
336 }
337
338 TBHitPlane *hitPaneV1 = new TBHitPlane();
339 hitPaneV1->setPosu(1); hitPaneV1->setPosw(1);
340 hitPlaneCont_V->push_back(hitPaneV1);
341
342 TBHitPlane *hitPaneV2 = new TBHitPlane();
343 hitPaneV2->setPosu(2); hitPaneV2->setPosw(2);
344 hitPlaneCont_V->push_back(hitPaneV2);
345
346 StatusCode sc2 = evtStore()->record(hitPlaneCont_V,"HitPlaneCont_V");
347 if ( sc2.isFailure( ) ) {
348 ATH_MSG_FATAL ("Cannot record HitPlaneCont_V");
349 }
350
351 ATH_MSG_DEBUG ("FillRandom Done");
352
353}
value_type push_back(value_type pElem)
Add an element to the end of the collection.

◆ finalize()

StatusCode TBPlaneTrackingAlgo::finalize ( )
overridevirtual

Definition at line 175 of file TBPlaneTrackingAlgo.cxx.

177{
178 ATH_MSG_INFO ("Finalizing TBTracking algorithm");
179 return StatusCode::SUCCESS;
180}
#define ATH_MSG_INFO(x)

◆ fitHits()

bool TBPlaneTrackingAlgo::fitHits ( const std::vector< double > & u,
const std::vector< double > & w,
const std::vector< double > & eu,
double & a1,
double & a2 )
private

Fit data to the function u = a1 + a2*w.

Definition at line 236 of file TBPlaneTrackingAlgo.cxx.

242{
243 // v_u = vector of u data points
244 // v_w = vector of w data points
245 // v_eu = vector of errors in u data points
246 // v_ew = vector of errors in w data points
247 // a1 and a2 = fit parameters: u = a1 + a2*w
248
249 // *** note that the fit algorithm used (given in 'Numerical Methods
250 // in C' section 15.2) assumes that the w data points are known exactly
251 // i.e. that v_ew[i]=0 for all i
252
253 // 'Numerical Methods' notes that the task of fitting a straight
254 // line model with errors in both coordinates is "considerably harder"
255 // (section 15.3) - but clearly it could be done
256
257 double m_s = 0;
258 double m_su = 0;
259 double m_sww = 0;
260 double m_sw = 0;
261 double m_suw = 0;
262
263 int numHits = v_u.size();
264 for(int i = 0; i < numHits; i++){
265 m_s += 1 / (v_eu[i]*v_eu[i]);
266 m_su += v_u[i] / (v_eu[i]*v_eu[i]);
267 m_sww += v_w[i]*v_w[i] / (v_eu[i]*v_eu[i]);
268 m_sw += v_w[i] / (v_eu[i]*v_eu[i]);
269 m_suw += v_u[i]*v_w[i] / (v_eu[i]*v_eu[i]);
270 }
271
272 const double denom = (m_s*m_sww-m_sw*m_sw);
273 if(denom == 0){
274 return false;
275 }
276
277 // coverity[divide_by_zero]
278 const double inv_denom = 1. / denom;
279 a1 = (m_su*m_sww - m_sw*m_suw) * inv_denom;
280 a2 = (m_s*m_suw - m_su*m_sw) * inv_denom;
281 return true;
282}

◆ fitPlane()

bool TBPlaneTrackingAlgo::fitPlane ( const TBHitPlaneCont * hitPlaneCont,
double & a1,
double & a2,
double & chi2,
std::vector< double > & residual )
private

Fit data to the function u = a1 + a2*w.

and determines intercept, slope, residual for each BPC, and chi2 on fit

Definition at line 183 of file TBPlaneTrackingAlgo.cxx.

186{
187 // Get number of hits in plane //
188 int hitPlaneNum = hitPlaneCont->size();
189
190 ATH_MSG_DEBUG ("The hit plane container size is: " << hitPlaneNum);
191
192 std::vector<double> vec_u;
193 std::vector<double> vec_w;
194 std::vector<double> err_vec_u;
195 std::vector<double> err_vec_w;
196
197 for (const TBHitPlane* hp : *hitPlaneCont) {
198 ATH_MSG_DEBUG ("Position in u: " << hp->getPosu());
199 ATH_MSG_DEBUG ("Position in w: " << hp->getPosw());
200 ATH_MSG_DEBUG ("Error in u: " << hp->getErroru());
201 ATH_MSG_DEBUG ("Error in w: " << hp->getErrorw());
202 vec_u.push_back(hp->getPosu());
203 vec_w.push_back(hp->getPosw());
204 err_vec_u.push_back(hp->getErroru());
205 err_vec_w.push_back(hp->getErrorw());
206 }
207
208 if(vec_u.size() != vec_w.size()){
209 ATH_MSG_ERROR ("TBPlaneTrackingAlgo: Invalid hits");
210 return false;
211 }
212
213 bool check;
214// check = fitHits(vec_u, vec_w, err_vec_u, err_vec_w, a1, a2);
215 check = fitHits(vec_u, vec_w, err_vec_u, a1, a2);
216 if(check == false){
217 ATH_MSG_ERROR ("TBPlaneTrackingAlgo: Invalid denumerator");
218 return false;
219 }
220
221// chi2 = getChi2(vec_u, vec_w, err_vec_u, err_vec_w, a1, a2);
222 chi2 = getChi2(vec_u, vec_w, err_vec_u, a1, a2);
223
224 ATH_MSG_DEBUG ("Fit results:"
225 << " intercept = " << a1
226 << " and slope = " << a2);
227
228 // Must be wrong //
229 for (int i = 0; i < hitPlaneNum; i++) {
230 residual.push_back(getResidual(vec_u[i],vec_w[i],a1, a2));
231 }
232 return true;
233}
size_type size() const noexcept
Returns the number of elements in the collection.
double getChi2(const std::vector< double > &v_u, const std::vector< double > &v_w, const std::vector< double > &v_eu, const double &a1, const double &a2)
Calculates the chi2 += ((u_i - a1 - a2*w_i)/eu)^2.
bool fitHits(const std::vector< double > &u, const std::vector< double > &w, const std::vector< double > &eu, double &a1, double &a2)
Fit data to the function u = a1 + a2*w.
double getResidual(const double &u, const double &w, const double &a1, const double &a2)
Calculates the residual-> r = (u_i - a1 - a2*w_i)
double chi2(TH1 *h0, TH1 *h1)

◆ getChi2()

double TBPlaneTrackingAlgo::getChi2 ( const std::vector< double > & v_u,
const std::vector< double > & v_w,
const std::vector< double > & v_eu,
const double & a1,
const double & a2 )
private

Calculates the chi2 += ((u_i - a1 - a2*w_i)/eu)^2.

Definition at line 292 of file TBPlaneTrackingAlgo.cxx.

298{
299 // *** as in fitHits() above, the chi squared calculation below
300 // (taken from 'Numerical Methods in C' section 15.2)
301 // assumes that the w data points are known exactly i.e. that
302 // v_ew[i]=0 for all i
303
304 int numX = v_u.size();
305 double chi2 = 0;
306
307 for(int i = 0; i < numX; i++){
308 chi2 += (v_u[i] - a1 - a2*v_w[i])*(v_u[i] - a1 - a2*v_w[i])/(v_eu[i]*v_eu[i]);
309 }
310 return chi2;
311}

◆ getnewcalib()

StatusCode TBPlaneTrackingAlgo::getnewcalib ( )
private

Definition at line 416 of file TBPlaneTrackingAlgo.cxx.

417{
418 // Get calib constant from an ASCII file with the following structure :
419 //
420 // runnumber
421 // bpcnumber1 coeff1 coeff2 ... coeff8
422 // bpcnumber2 coeff1 coeff2 ... coeff8
423 // ...
424 // bpcnumber6 coeff1 coeff2 ... coeff8
425 // runnumber
426 // ...
427 //
428 // coeff must have the following order :
429 // bpcnumber posX posY posZX posZY errposX errposY errposZX errposZY errmeasX errmeasY
430
431 ATH_MSG_DEBUG ("Get new calibs for run " << m_runnumber);
432
433 int bpcnumber= m_bpc_names.size();
434
435 m_bpc_posX.clear(); m_bpc_posX.resize(bpcnumber);
436 m_bpc_posY.clear(); m_bpc_posY.resize(bpcnumber);
437 m_bpc_posZX.clear(); m_bpc_posZX.resize(bpcnumber);
438 m_bpc_posZY.clear(); m_bpc_posZY.resize(bpcnumber);
439 m_bpc_errposX.clear(); m_bpc_errposX.resize(bpcnumber);
440 m_bpc_errposY.clear(); m_bpc_errposY.resize(bpcnumber);
441 m_bpc_errposZX.clear(); m_bpc_errposZX.resize(bpcnumber);
442 m_bpc_errposZY.clear(); m_bpc_errposZY.resize(bpcnumber);
443 m_bpc_errmeasX.clear(); m_bpc_errmeasX.resize(bpcnumber);
444 m_bpc_errmeasY.clear(); m_bpc_errmeasY.resize(bpcnumber);
445
446 int pos;
447
448 std::ifstream calibfile;
449 std::string filename = PathResolver::find_file (m_calib_filename, "DATAPATH");
450 calibfile.open(filename.c_str());
451 if(!calibfile.good()){
452 ATH_MSG_WARNING (" Problem with file named "<< m_calib_filename << " in $DATAPATH");
453 return StatusCode::FAILURE;
454 } else {
455 ATH_MSG_DEBUG (" File: "<< filename << " opened");
456 }
457 unsigned int runnumber;
458 calibfile >> runnumber;
459 pos = calibfile.tellg();
460 ATH_MSG_DEBUG (" Run number "<< runnumber);
461 while((runnumber<m_runnumber)&&(!calibfile.eof()))
462 {
463 runnumber=0;
464 pos = calibfile.tellg();
465 // discard next lines
466 for(int j=0;j<bpcnumber+1;j++) calibfile.ignore(5000,'\n');
467 // check next runnumber
468 calibfile >> runnumber;
469 if(runnumber==0) { // reached an empty line : exit.
470 ATH_MSG_DEBUG ("empty line");
471 calibfile.clear();
472 break;
473 }
474 ATH_MSG_DEBUG (" Run number "<< runnumber);
475 }
476
477 // Now we found the good set of constant (the ones following pos)
478 if(runnumber==m_runnumber) pos = calibfile.tellg();
479 ATH_MSG_DEBUG (" Pos = "<< pos);
480 calibfile.seekg(pos);
481 ATH_MSG_DEBUG (" Will use the following constants :");
482 for(int j=0;j<bpcnumber;j++)
483 {
484 int bpcn;
485 calibfile >> bpcn;
486 calibfile >> m_bpc_posX[j];
487 calibfile >> m_bpc_posY[j];
488 calibfile >> m_bpc_posZX[j];
489 calibfile >> m_bpc_posZY[j];
490 calibfile >> m_bpc_errposX[j];
491 calibfile >> m_bpc_errposY[j];
492 calibfile >> m_bpc_errposZX[j];
493 calibfile >> m_bpc_errposZY[j];
494 calibfile >> m_bpc_errmeasX[j];
495 calibfile >> m_bpc_errmeasY[j];
496
497 ATH_MSG_DEBUG (bpcn << " "<<m_bpc_posX[j]
498 << " "<< m_bpc_posY[j]
499 << " "<< m_bpc_posZX[j]
500 << " "<< m_bpc_posZY[j]
501 << " "<< m_bpc_errposX[j]
502 << " "<< m_bpc_errposY[j]
503 << " "<< m_bpc_errposZX[j]
504 << " "<< m_bpc_errposZY[j]
505 << " "<< m_bpc_errmeasX[j]
506 << " "<< m_bpc_errmeasY[j]);
507 }
508
509 calibfile.close();
510
511 return StatusCode::SUCCESS;
512
513}
static std::string find_file(const std::string &logical_file_name, const std::string &search_path)
static std::vector< uint32_t > runnumber
Definition iLumiCalc.h:37

◆ getResidual()

double TBPlaneTrackingAlgo::getResidual ( const double & u,
const double & w,
const double & a1,
const double & a2 )
private

Calculates the residual-> r = (u_i - a1 - a2*w_i)

Definition at line 284 of file TBPlaneTrackingAlgo.cxx.

287{
288 return (u - a1 - a2*w);
289}

◆ initialize()

StatusCode TBPlaneTrackingAlgo::initialize ( )
overridevirtual

Definition at line 42 of file TBPlaneTrackingAlgo.cxx.

44{
45 return StatusCode::SUCCESS;
46}

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Algorithm > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ msg()

MsgStream & AthCommonMsg< Algorithm >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msgLvl()

bool AthCommonMsg< Algorithm >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Algorithm > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< Algorithm > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sysInitialize()

StatusCode AthAlgorithm::sysInitialize ( )
overridevirtualinherited

Override sysInitialize.

Override sysInitialize from the base class.

Loop through all output handles, and if they're WriteCondHandles, automatically register them and this Algorithm with the CondSvc

Scan through all outputHandles, and if they're WriteCondHandles, register them with the CondSvc

Reimplemented from AthCommonDataStore< AthCommonMsg< Algorithm > >.

Reimplemented in AthAnalysisAlgorithm, AthFilterAlgorithm, AthHistogramAlgorithm, and PyAthena::Alg.

Definition at line 66 of file AthAlgorithm.cxx.

66 {
68
69 if (sc.isFailure()) {
70 return sc;
71 }
72 ServiceHandle<ICondSvc> cs("CondSvc",name());
73 for (auto h : outputHandles()) {
74 if (h->isCondition() && h->mode() == Gaudi::DataHandle::Writer) {
75 // do this inside the loop so we don't create the CondSvc until needed
76 if ( cs.retrieve().isFailure() ) {
77 ATH_MSG_WARNING("no CondSvc found: won't autoreg WriteCondHandles");
78 return StatusCode::SUCCESS;
79 }
80 if (cs->regHandle(this,*h).isFailure()) {
81 sc = StatusCode::FAILURE;
82 ATH_MSG_ERROR("unable to register WriteCondHandle " << h->fullKey()
83 << " with CondSvc");
84 }
85 }
86 }
87 return sc;
88}
virtual StatusCode sysInitialize() override
Override sysInitialize.
AthCommonDataStore(const std::string &name, T... args)
virtual std::vector< Gaudi::DataHandle * > outputHandles() const override

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< Algorithm > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< Algorithm > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

Member Data Documentation

◆ m_bpc_errmeasX

std::vector<float> TBPlaneTrackingAlgo::m_bpc_errmeasX
private

Definition at line 83 of file TBPlaneTrackingAlgo.h.

◆ m_bpc_errmeasY

std::vector<float> TBPlaneTrackingAlgo::m_bpc_errmeasY
private

Definition at line 83 of file TBPlaneTrackingAlgo.h.

◆ m_bpc_errposX

std::vector<float> TBPlaneTrackingAlgo::m_bpc_errposX
private

Definition at line 82 of file TBPlaneTrackingAlgo.h.

◆ m_bpc_errposY

std::vector<float> TBPlaneTrackingAlgo::m_bpc_errposY
private

Definition at line 82 of file TBPlaneTrackingAlgo.h.

◆ m_bpc_errposZX

std::vector<float> TBPlaneTrackingAlgo::m_bpc_errposZX
private

Definition at line 82 of file TBPlaneTrackingAlgo.h.

◆ m_bpc_errposZY

std::vector<float> TBPlaneTrackingAlgo::m_bpc_errposZY
private

Definition at line 82 of file TBPlaneTrackingAlgo.h.

◆ m_bpc_names

std::vector<std::string> TBPlaneTrackingAlgo::m_bpc_names
private

Definition at line 85 of file TBPlaneTrackingAlgo.h.

◆ m_bpc_posX

std::vector<float> TBPlaneTrackingAlgo::m_bpc_posX
private

Definition at line 81 of file TBPlaneTrackingAlgo.h.

◆ m_bpc_posY

std::vector<float> TBPlaneTrackingAlgo::m_bpc_posY
private

Definition at line 81 of file TBPlaneTrackingAlgo.h.

◆ m_bpc_posZX

std::vector<float> TBPlaneTrackingAlgo::m_bpc_posZX
private

Definition at line 81 of file TBPlaneTrackingAlgo.h.

◆ m_bpc_posZY

std::vector<float> TBPlaneTrackingAlgo::m_bpc_posZY
private

Definition at line 81 of file TBPlaneTrackingAlgo.h.

◆ m_calib_filename

std::string TBPlaneTrackingAlgo::m_calib_filename
private

Definition at line 88 of file TBPlaneTrackingAlgo.h.

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< Algorithm > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< Algorithm > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_extendedExtraObjects

DataObjIDColl AthAlgorithm::m_extendedExtraObjects
privateinherited

Definition at line 79 of file AthAlgorithm.h.

◆ m_hitPlaneCont_u

TBHitPlaneCont TBPlaneTrackingAlgo::m_hitPlaneCont_u
private

Definition at line 77 of file TBPlaneTrackingAlgo.h.

◆ m_hitPlaneCont_v

TBHitPlaneCont TBPlaneTrackingAlgo::m_hitPlaneCont_v
private

Definition at line 78 of file TBPlaneTrackingAlgo.h.

◆ m_runnumber

unsigned int TBPlaneTrackingAlgo::m_runnumber
private

Definition at line 89 of file TBPlaneTrackingAlgo.h.

◆ m_testAlgo

bool TBPlaneTrackingAlgo::m_testAlgo
private

Definition at line 76 of file TBPlaneTrackingAlgo.h.

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< Algorithm > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< Algorithm > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.


The documentation for this class was generated from the following files: