ATLAS Offline Software
Loading...
Searching...
No Matches
PropResultRootWriterSvc.h
Go to the documentation of this file.
1
2/*
3 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
4*/
5
6#ifndef TRKEXALGS_PROPRESULTROOTWRITERSVC_H
7#define TRKEXALGS_PROPRESULTROOTWRITERSVC_H
8
9#include "Gaudi/Property.h" /*no forward decl: typedef*/
10#include "GaudiKernel/ITHistSvc.h"
12
13#include "TTree.h"
14
15#include <mutex>
16
17class TFile;
18
19namespace Trk {
20
22 public:
24 PropResultRootWriterSvc(const std::string& name, ISvcLocator* svcloc);
27
28 virtual StatusCode initialize() override;
29
30 template <typename T>
31 void write(const T* initialPerigee,
32 const T* fwdParameters=nullptr, double fwdtime=std::numeric_limits<float>::quiet_NaN(),
33 const T* bkwParameters=nullptr, double bkwtime=std::numeric_limits<float>::quiet_NaN() );
34
35private:
36 std::mutex m_writeMutex;
37
38 // jobOptions properties
40 TTree* m_tree;
41 Gaudi::Property<std::string> m_ntupleDirName{this, "DirName", "/ExtrapolationStudies/", ""};
42 Gaudi::Property<std::string> m_treeName{this, "TreeName", "ATLAS", ""};
43
44 int m_eventNum = 0;
45
46 float m_start_d0 = 0.0F ;
47 float m_start_z0 = 0.0F ;
48 float m_start_phi = 0.0F ;
49 float m_start_theta = 0.0F ;
50 float m_start_qop = 0.0F ;
51
52 float m_fwd_success = 0.0F ;
53 float m_fwd_time = 0.0F ;
54 float m_fwd_final_l0 = 0.0F ;
55 float m_fwd_final_l1 = 0.0F ;
56 float m_fwd_final_x = 0.0F ;
57 float m_fwd_final_y = 0.0F ;
58 float m_fwd_final_z = 0.0F ;
59 float m_fwd_final_phi = 0.0F ;
60 float m_fwd_final_theta = 0.0F ;
61 float m_fwd_final_qop = 0.0F ;
62 float m_fwd_final_sigma_l0 = 0.0F ;
63 float m_fwd_final_sigma_l1 = 0.0F ;
64 float m_fwd_final_sigma_phi = 0.0F ;
66 float m_fwd_final_sigma_qop = 0.0F ;
67
68 float m_bkw_success = 0.0F ;
69 float m_bkw_time = 0.0F ;
70 float m_bkw_final_d0 = 0.0F ;
71 float m_bkw_final_z0 = 0.0F ;
72 float m_bkw_final_phi = 0.0F ;
73 float m_bkw_final_theta = 0.0F ;
74 float m_bkw_final_qop = 0.0F ;
75 float m_bkw_final_sigma_d0 = 0.0F ;
76 float m_bkw_final_sigma_z0 = 0.0F ;
77 float m_bkw_final_sigma_phi = 0.0F ;
79 float m_bkw_final_sigma_qop = 0.0F ;
80
81 };
82
83 template <typename T>
84 void PropResultRootWriterSvc::write(const T* initialPerigee,
85 const T* fwdParameters, double fwdtime,
86 const T* bkwParameters, double bkwtime ) {
87
88 const auto& ctx = Gaudi::Hive::currentContext();
89
90 m_eventNum = ctx.eventID().event_number();
91
92 if (initialPerigee) {
93 m_start_d0 = initialPerigee->parameters()[0];
94 m_start_z0 = initialPerigee->parameters()[1];
95 m_start_phi = initialPerigee->parameters()[2];
96 m_start_theta= initialPerigee->parameters()[3];
97 m_start_qop = initialPerigee->parameters()[4];
98 } else {
99 m_start_d0 = std::numeric_limits<float>::quiet_NaN();
100 m_start_z0 = std::numeric_limits<float>::quiet_NaN();
101 m_start_phi = std::numeric_limits<float>::quiet_NaN();
102 m_start_theta= std::numeric_limits<float>::quiet_NaN();
103 m_start_qop = std::numeric_limits<float>::quiet_NaN();
104 }
105
106 if (fwdParameters) {
107 m_fwd_success = 1;
108 m_fwd_time = fwdtime;
109 m_fwd_final_l0 = fwdParameters->parameters()[0];
110 m_fwd_final_l1 = fwdParameters->parameters()[1];
111 m_fwd_final_phi = fwdParameters->parameters()[2];
112 m_fwd_final_theta = fwdParameters->parameters()[3];
113 m_fwd_final_qop = fwdParameters->parameters()[4];
114 m_fwd_final_x = fwdParameters->position()[0];
115 m_fwd_final_y = fwdParameters->position()[1];
116 m_fwd_final_z = fwdParameters->position()[2];
117 if (fwdParameters->covariance()) {
118 m_fwd_final_sigma_l0 = (*fwdParameters->covariance())(0,0);
119 m_fwd_final_sigma_l1 = (*fwdParameters->covariance())(1,1);
120 m_fwd_final_sigma_phi = (*fwdParameters->covariance())(2,2);
121 m_fwd_final_sigma_theta = (*fwdParameters->covariance())(3,3);
122 m_fwd_final_sigma_qop = (*fwdParameters->covariance())(4,4);
123 } else {
124 m_fwd_final_sigma_l0 = std::numeric_limits<float>::quiet_NaN();
125 m_fwd_final_sigma_l1 = std::numeric_limits<float>::quiet_NaN();
126 m_fwd_final_sigma_phi = std::numeric_limits<float>::quiet_NaN();
127 m_fwd_final_sigma_theta = std::numeric_limits<float>::quiet_NaN();
128 m_fwd_final_sigma_qop = std::numeric_limits<float>::quiet_NaN();
129 }
130 } else {
131 m_fwd_success = std::numeric_limits<float>::quiet_NaN();
132 m_fwd_time = std::numeric_limits<float>::quiet_NaN();
133 m_fwd_final_l0 = std::numeric_limits<float>::quiet_NaN();
134 m_fwd_final_l1 = std::numeric_limits<float>::quiet_NaN();
135 m_fwd_final_phi = std::numeric_limits<float>::quiet_NaN();
136 m_fwd_final_theta = std::numeric_limits<float>::quiet_NaN();
137 m_fwd_final_qop = std::numeric_limits<float>::quiet_NaN();
138 m_fwd_final_x = std::numeric_limits<float>::quiet_NaN();
139 m_fwd_final_y = std::numeric_limits<float>::quiet_NaN();
140 m_fwd_final_z = std::numeric_limits<float>::quiet_NaN();
141 m_fwd_final_sigma_l0 = std::numeric_limits<float>::quiet_NaN();
142 m_fwd_final_sigma_l1 = std::numeric_limits<float>::quiet_NaN();
143 m_fwd_final_sigma_phi = std::numeric_limits<float>::quiet_NaN();
144 m_fwd_final_sigma_theta = std::numeric_limits<float>::quiet_NaN();
145 m_fwd_final_sigma_qop = std::numeric_limits<float>::quiet_NaN();
146 }
147
148 if (bkwParameters) {
149 m_bkw_success = 1;
150 m_bkw_time = bkwtime;
151 m_bkw_final_d0 = bkwParameters->parameters()[0];
152 m_bkw_final_z0 = bkwParameters->parameters()[1];
153 m_bkw_final_phi = bkwParameters->parameters()[2];
154 m_bkw_final_theta = bkwParameters->parameters()[3];
155 m_bkw_final_qop = bkwParameters->parameters()[4];
156 if (bkwParameters->covariance()) {
157 m_bkw_final_sigma_d0 = (*bkwParameters->covariance())(0,0);
158 m_bkw_final_sigma_z0 = (*bkwParameters->covariance())(1,1);
159 m_bkw_final_sigma_phi = (*bkwParameters->covariance())(2,2);
160 m_bkw_final_sigma_theta = (*bkwParameters->covariance())(3,3);
161 m_bkw_final_sigma_qop = (*bkwParameters->covariance())(4,4);
162 } else {
163 m_bkw_final_sigma_d0 = std::numeric_limits<float>::quiet_NaN();
164 m_bkw_final_sigma_z0 = std::numeric_limits<float>::quiet_NaN();
165 m_bkw_final_sigma_phi = std::numeric_limits<float>::quiet_NaN();
166 m_bkw_final_sigma_theta = std::numeric_limits<float>::quiet_NaN();
167 m_bkw_final_sigma_qop = std::numeric_limits<float>::quiet_NaN();
168 }
169 } else {
170 m_bkw_success = std::numeric_limits<float>::quiet_NaN();
171 m_bkw_time = std::numeric_limits<float>::quiet_NaN();
172 m_bkw_final_d0 = std::numeric_limits<float>::quiet_NaN();
173 m_bkw_final_z0 = std::numeric_limits<float>::quiet_NaN();
174 m_bkw_final_phi = std::numeric_limits<float>::quiet_NaN();
175 m_bkw_final_theta = std::numeric_limits<float>::quiet_NaN();
176 m_bkw_final_qop = std::numeric_limits<float>::quiet_NaN();
177 m_bkw_final_sigma_d0 = std::numeric_limits<float>::quiet_NaN();
178 m_bkw_final_sigma_z0 = std::numeric_limits<float>::quiet_NaN();
179 m_bkw_final_sigma_phi = std::numeric_limits<float>::quiet_NaN();
180 m_bkw_final_sigma_theta = std::numeric_limits<float>::quiet_NaN();
181 m_bkw_final_sigma_qop = std::numeric_limits<float>::quiet_NaN();
182 }
183 m_tree->Fill();
184 }
185
186
187
188}
189
190#endif
Gaudi::Property< std::string > m_ntupleDirName
virtual ~PropResultRootWriterSvc()
destructor
Gaudi::Property< std::string > m_treeName
ServiceHandle< ITHistSvc > m_thistSvc
virtual StatusCode initialize() override
PropResultRootWriterSvc(const std::string &name, ISvcLocator *svcloc)
constructor
void write(const T *initialPerigee, const T *fwdParameters=nullptr, double fwdtime=std::numeric_limits< float >::quiet_NaN(), const T *bkwParameters=nullptr, double bkwtime=std::numeric_limits< float >::quiet_NaN())
Ensure that the ATLAS eigen extensions are properly loaded.