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

#include <FPGATrackSimMatrixIO.h>

Collaboration diagram for FPGATrackSimMatrixWriter:

Public Member Functions

 FPGATrackSimMatrixWriter (TTree *tree, int nLayers, int nCoords)
size_t fill (std::vector< module_t > modules, FPGATrackSimMatrixAccumulator &acc)

Private Attributes

TTree * m_tree
size_t m_nEntries
int m_nLayers
int m_nCoords
int m_nCoords2
float m_coverage = 0.0F
std::vector< short > m_bins_QoP
std::vector< short > m_bins_phi
std::vector< short > m_bins_d0
std::vector< short > m_bins_z0
std::vector< short > m_bins_eta
std::vector< short > * m_pQoP
std::vector< short > * m_pphi
std::vector< short > * m_pd0
std::vector< short > * m_pz0
std::vector< short > * m_peta

Detailed Description

Definition at line 86 of file FPGATrackSimMatrixIO.h.

Constructor & Destructor Documentation

◆ FPGATrackSimMatrixWriter()

FPGATrackSimMatrixWriter::FPGATrackSimMatrixWriter ( TTree * tree,
int nLayers,
int nCoords )

Definition at line 125 of file FPGATrackSimMatrixIO.cxx.

125 :
126 m_tree(tree),
127 m_nEntries(0),
128 m_nLayers(nLayers),
129 m_nCoords(nCoords),
130 m_nCoords2(nCoords * nCoords),
136{
137 // Dummy variables for typing (the branch address will be set for each sector)
138 int anInt;
139 double aDouble;
140 float aFloat;
141
142 m_tree->Branch("ndim", &m_nCoords, "ndim/I");
143 m_tree->Branch("ndim2", &m_nCoords2, "ndim2/I");
144 m_tree->Branch("nplanes", &m_nLayers, "nplanes/I");
145
146 m_tree->Branch("sectorID", &anInt, "sectorID[nplanes]/I");
147 m_tree->Branch("hashID", &anInt, "hashID[nplanes]/I");
148
149 m_tree->Branch("tmpC", &aDouble, "tmpC/D");
150 m_tree->Branch("tmpD", &aDouble, "tmpD/D");
151 m_tree->Branch("tmpPhi", &aDouble, "tmpPhi/D");
152 m_tree->Branch("tmpCoto", &aDouble, "tmpCoto/D");
153 m_tree->Branch("tmpZ", &aDouble, "tmpZ/D");
154
155 m_tree->Branch("Vec", &aFloat, "Vec[ndim]/F");
156 m_tree->Branch("VecG", &aFloat, "VecG[ndim]/F");
157 m_tree->Branch("tmpxC", &aDouble, "tmpxC[ndim]/D");
158 m_tree->Branch("tmpxGC", &aDouble, "tmpxGC[ndim]/D");
159 m_tree->Branch("tmpxD", &aDouble, "tmpxD[ndim]/D");
160 m_tree->Branch("tmpxPhi", &aDouble, "tmpxPhi[ndim]/D");
161 m_tree->Branch("tmpxCoto", &aDouble, "tmpxCoto[ndim]/D");
162 m_tree->Branch("tmpxGCoto", &aDouble, "tmpxGCoto[ndim]/D");
163 m_tree->Branch("tmpxZ", &aDouble, "tmpxZ[ndim]/D");
164 m_tree->Branch("tmpcovx", &aDouble, "tmpcovx[ndim2]/D");
165 m_tree->Branch("tmpcovxG", &aDouble, "tmpcovxG[ndim2]/D");
166
167 m_tree->Branch("nhit", &m_coverage, "nhit/F");
168 m_tree->Branch("tmpintc", &m_pQoP);
169 m_tree->Branch("tmpintphi", &m_pphi);
170 m_tree->Branch("tmpintd0", &m_pd0);
171 m_tree->Branch("tmpintz0", &m_pz0);
172 m_tree->Branch("tmpinteta", &m_peta);
173}
std::vector< short > m_bins_d0
std::vector< short > m_bins_eta
std::vector< short > m_bins_z0
std::vector< short > * m_peta
std::vector< short > m_bins_phi
std::vector< short > * m_pz0
std::vector< short > m_bins_QoP
std::vector< short > * m_pd0
std::vector< short > * m_pphi
std::vector< short > * m_pQoP
TChain * tree

Member Function Documentation

◆ fill()

size_t FPGATrackSimMatrixWriter::fill ( std::vector< module_t > modules,
FPGATrackSimMatrixAccumulator & acc )

Definition at line 175 of file FPGATrackSimMatrixIO.cxx.

176{
177
178 m_tree->SetBranchAddress("sectorID", acc.FTK_modules.data());
179 m_tree->SetBranchAddress("hashID", modules.data());
180
181 m_tree->SetBranchAddress("tmpC", &acc.pars.qOverPt);
182 m_tree->SetBranchAddress("tmpD", &acc.pars.d0);
183 m_tree->SetBranchAddress("tmpPhi", &acc.pars.phi);
184 m_tree->SetBranchAddress("tmpCoto", &acc.pars.eta);
185 m_tree->SetBranchAddress("tmpZ", &acc.pars.z0);
186
187 m_tree->SetBranchAddress("Vec", acc.hit_coords.data());
188 m_tree->SetBranchAddress("VecG", acc.hit_coordsG.data());
189 m_tree->SetBranchAddress("tmpxC", acc.hit_x_QoP.data());
190 m_tree->SetBranchAddress("tmpxGC", acc.hit_xG_HIP.data());
191 m_tree->SetBranchAddress("tmpxD", acc.hit_x_d0.data());
192 m_tree->SetBranchAddress("tmpxPhi", acc.hit_x_phi.data());
193 m_tree->SetBranchAddress("tmpxCoto", acc.hit_x_eta.data());
194 m_tree->SetBranchAddress("tmpxGCoto", acc.hit_xG_eta.data());
195 m_tree->SetBranchAddress("tmpxZ", acc.hit_x_z0.data());
196 m_tree->SetBranchAddress("tmpcovx", acc.covariance.data());
197 m_tree->SetBranchAddress("tmpcovxG", acc.covarianceG.data());
198
199 m_bins_QoP.clear();
200 m_bins_phi.clear();
201 m_bins_z0.clear();
202 m_bins_d0.clear();
203 m_bins_eta.clear();
204 for (FPGATrackSimTrackParsI const & pars : acc.track_bins)
205 {
206 m_bins_QoP.push_back(pars.qOverPt);
207 m_bins_phi.push_back(pars.phi);
208 m_bins_d0.push_back(pars.d0);
209 m_bins_z0.push_back(pars.z0);
210 m_bins_eta.push_back(pars.eta);
211 }
212 m_coverage = static_cast<float>(acc.track_bins.size());
213 m_tree->Fill();
214 return m_nEntries++;
215}

Member Data Documentation

◆ m_bins_d0

std::vector<short> FPGATrackSimMatrixWriter::m_bins_d0
private

Definition at line 107 of file FPGATrackSimMatrixIO.h.

◆ m_bins_eta

std::vector<short> FPGATrackSimMatrixWriter::m_bins_eta
private

Definition at line 109 of file FPGATrackSimMatrixIO.h.

◆ m_bins_phi

std::vector<short> FPGATrackSimMatrixWriter::m_bins_phi
private

Definition at line 106 of file FPGATrackSimMatrixIO.h.

◆ m_bins_QoP

std::vector<short> FPGATrackSimMatrixWriter::m_bins_QoP
private

Definition at line 105 of file FPGATrackSimMatrixIO.h.

◆ m_bins_z0

std::vector<short> FPGATrackSimMatrixWriter::m_bins_z0
private

Definition at line 108 of file FPGATrackSimMatrixIO.h.

◆ m_coverage

float FPGATrackSimMatrixWriter::m_coverage = 0.0F
private

Definition at line 103 of file FPGATrackSimMatrixIO.h.

◆ m_nCoords

int FPGATrackSimMatrixWriter::m_nCoords
private

Definition at line 101 of file FPGATrackSimMatrixIO.h.

◆ m_nCoords2

int FPGATrackSimMatrixWriter::m_nCoords2
private

Definition at line 102 of file FPGATrackSimMatrixIO.h.

◆ m_nEntries

size_t FPGATrackSimMatrixWriter::m_nEntries
private

Definition at line 99 of file FPGATrackSimMatrixIO.h.

◆ m_nLayers

int FPGATrackSimMatrixWriter::m_nLayers
private

Definition at line 100 of file FPGATrackSimMatrixIO.h.

◆ m_pd0

std::vector<short>* FPGATrackSimMatrixWriter::m_pd0
private

Definition at line 114 of file FPGATrackSimMatrixIO.h.

◆ m_peta

std::vector<short>* FPGATrackSimMatrixWriter::m_peta
private

Definition at line 116 of file FPGATrackSimMatrixIO.h.

◆ m_pphi

std::vector<short>* FPGATrackSimMatrixWriter::m_pphi
private

Definition at line 113 of file FPGATrackSimMatrixIO.h.

◆ m_pQoP

std::vector<short>* FPGATrackSimMatrixWriter::m_pQoP
private

Definition at line 112 of file FPGATrackSimMatrixIO.h.

◆ m_pz0

std::vector<short>* FPGATrackSimMatrixWriter::m_pz0
private

Definition at line 115 of file FPGATrackSimMatrixIO.h.

◆ m_tree

TTree* FPGATrackSimMatrixWriter::m_tree
private

Definition at line 97 of file FPGATrackSimMatrixIO.h.


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