LDMX Software
ecal::EcalVetoProcessor Class Reference

Determines if event is vetoable using ECAL hit information. More...

#include <EcalVetoProcessor.h>

Classes

struct  HitData
 

Public Types

typedef std::pair< ldmx::EcalID, float > CellEnergyPair
 
typedef std::pair< float, float > XYCoords
 

Public Member Functions

 EcalVetoProcessor (const std::string &name, framework::Process &process)
 
void configure (framework::config::Parameters &parameters) override
 Configure the processor using the given user specified parameters.
 
void produce (framework::Event &event) override
 Process the event and put new data products into it.
 
- Public Member Functions inherited from framework::Producer
 Producer (const std::string &name, Process &process)
 Class constructor.
 
virtual void beforeNewRun (ldmx::RunHeader &header)
 Handle allowing producers to modify run headers before the run begins.
 
- Public Member Functions inherited from framework::EventProcessor
 EventProcessor (const std::string &name, Process &process)
 Class constructor.
 
virtual ~EventProcessor ()
 Class destructor.
 
virtual void onNewRun (const ldmx::RunHeader &runHeader)
 Callback for the EventProcessor to take any necessary action when the run being processed changes.
 
virtual void onFileOpen (EventFile &eventFile)
 Callback for the EventProcessor to take any necessary action when a new event input ROOT file is opened.
 
virtual void onFileClose (EventFile &eventFile)
 Callback for the EventProcessor to take any necessary action when a event input ROOT file is closed.
 
virtual void onProcessStart ()
 Callback for the EventProcessor to take any necessary action when the processing of events starts, such as creating histograms.
 
virtual void onProcessEnd ()
 Callback for the EventProcessor to take any necessary action when the processing of events finishes, such as calculating job-summary quantities.
 
template<class T >
const T & getCondition (const std::string &condition_name)
 Access a conditions object for the current event.
 
TDirectory * getHistoDirectory ()
 Access/create a directory in the histogram file for this event processor to create histograms and analysis tuples.
 
void setStorageHint (framework::StorageControl::Hint hint)
 Mark the current event as having the given storage control hint from this module.
 
void setStorageHint (framework::StorageControl::Hint hint, const std::string &purposeString)
 Mark the current event as having the given storage control hint from this module and the given purpose string.
 
int getLogFrequency () const
 Get the current logging frequency from the process.
 
int getRunNumber () const
 Get the run number from the process.
 
std::string getName () const
 Get the processor name.
 
void createHistograms (const std::vector< framework::config::Parameters > &histos)
 Internal function which is used to create histograms passed from the python configuration @parma histos vector of Parameters that configure histograms to create.
 

Private Member Functions

void clearProcessor ()
 
ldmx::EcalID GetShowerCentroidIDAndRMS (const std::vector< ldmx::EcalHit > &ecalRecHits, double &showerRMS)
 
void fillHitMap (const std::vector< ldmx::EcalHit > &ecalRecHits, std::map< ldmx::EcalID, float > &cellMap_)
 Function to load up empty vector of hit maps.
 
void fillIsolatedHitMap (const std::vector< ldmx::EcalHit > &ecalRecHits, ldmx::EcalID globalCentroid, std::map< ldmx::EcalID, float > &cellMap, std::map< ldmx::EcalID, float > &cellMapIso, bool doTight=false)
 
std::vector< XYCoords > getTrajectory (std::vector< double > momentum, std::vector< float > position)
 
void buildBDTFeatureVector (const ldmx::EcalVetoResult &result)
 
float distTwoLines (TVector3 v1, TVector3 v2, TVector3 w1, TVector3 w2)
 Returns the distance between the lines v and w, with v defined to pass through the points (v1,v2) (and similarly for w).
 
float distPtToLine (TVector3 h1, TVector3 p1, TVector3 p2)
 Return the minimum distance between the point h1 and the line passing through points p1 and p2.
 
std::vector< float > trackProp (const ldmx::Tracks &tracks, ldmx::TrackStateType ts_type, const std::string &ts_title)
 Return a vector of parameters for a propagated recoil track.
 

Private Attributes

std::map< ldmx::EcalID, float > cellMap_
 
std::map< ldmx::EcalID, float > cellMapTightIso_
 
std::vector< float > ecalLayerEdepRaw_
 
std::vector< float > ecalLayerEdepReadout_
 
std::vector< float > ecalLayerTime_
 
std::vector< std::vector< double > > roc_range_values_
 
int nEcalLayers_ {0}
 
int nReadoutHits_ {0}
 
int deepestLayerHit_ {0}
 
double summedDet_ {0}
 
double summedTightIso_ {0}
 
double maxCellDep_ {0}
 
double showerRMS_ {0}
 
double xStd_ {0}
 
double yStd_ {0}
 
double avgLayerHit_ {0}
 
double stdLayerHit_ {0}
 
double ecalBackEnergy_ {0}
 
int nStraightTracks_ {0}
 Number of "straight" tracks found in the event.
 
int nLinregTracks_ {0}
 Number of "linreg" tracks found in the event.
 
int firstNearPhLayer_ {0}
 First ECal layer in which a hit is found near the photon.
 
int nNearPhHits_ {0}
 Number of hits near the photon trajectory.
 
float epAng_ {0}
 Angular separation between the projected photon and electron trajectories (currently unused)
 
float epSep_ {0}
 Distance between the projected photon and electron trajectories at the ECal face.
 
float epDot_ {0}
 Dot product of the photon and electron momenta unit vectors.
 
int photonTerritoryHits_ {0}
 Number of hits in the photon territory.
 
double bdtCutVal_ {0}
 
double beamEnergyMeV_ {0}
 
double linreg_radius_ {0}
 
bool verbose_ {false}
 
std::string bdtFileName_
 
std::string rocFileName_
 
std::vector< float > bdtFeatures_
 
std::string featureListName_
 
std::string rec_pass_name_
 
std::string rec_coll_name_
 
bool recoil_from_tracking_
 
std::string track_collection_
 
bool inverse_skim_ {false}
 
std::string collectionName_ {"EcalVeto"}
 Name of the collection which will containt the results.
 
std::unique_ptr< ldmx::Ort::ONNXRuntimert_
 
const ldmx::EcalGeometrygeometry_
 handle to current geometry (to share with member functions)
 

Additional Inherited Members

- Static Public Member Functions inherited from framework::EventProcessor
static void declare (const std::string &classname, int classtype, EventProcessorMaker *)
 Internal function which is part of the PluginFactory machinery.
 
- Static Public Attributes inherited from framework::Producer
static const int CLASSTYPE {1}
 Constant used to track EventProcessor types by the PluginFactory.
 
- Protected Member Functions inherited from framework::EventProcessor
void abortEvent ()
 Abort the event immediately.
 
- Protected Attributes inherited from framework::EventProcessor
HistogramHelper histograms_
 Interface class for making and filling histograms.
 
NtupleManagerntuple_ {NtupleManager::getInstance()}
 Manager for any ntuples.
 
logging::logger theLog_
 The logger for this EventProcessor.
 

Detailed Description

Determines if event is vetoable using ECAL hit information.

Definition at line 35 of file EcalVetoProcessor.h.

Member Typedef Documentation

◆ CellEnergyPair

std::pair<ldmx::EcalID, float> ecal::EcalVetoProcessor::CellEnergyPair

Definition at line 37 of file EcalVetoProcessor.h.

◆ XYCoords

std::pair<float, float> ecal::EcalVetoProcessor::XYCoords

Definition at line 39 of file EcalVetoProcessor.h.

Constructor & Destructor Documentation

◆ EcalVetoProcessor()

ecal::EcalVetoProcessor::EcalVetoProcessor ( const std::string & name,
framework::Process & process )
inline

Definition at line 41 of file EcalVetoProcessor.h.

42 : Producer(name, process) {}
Producer(const std::string &name, Process &process)
Class constructor.

◆ ~EcalVetoProcessor()

virtual ecal::EcalVetoProcessor::~EcalVetoProcessor ( )
inlinevirtual

Definition at line 44 of file EcalVetoProcessor.h.

44{}

Member Function Documentation

◆ buildBDTFeatureVector()

void ecal::EcalVetoProcessor::buildBDTFeatureVector ( const ldmx::EcalVetoResult & result)
private

Electron RoC variables

Photon RoC variables

Outside RoC variables

Definition at line 30 of file EcalVetoProcessor.cxx.

31 {
32 // Base variables
33 bdtFeatures_.push_back(result.getNReadoutHits());
34 bdtFeatures_.push_back(result.getSummedDet());
35 bdtFeatures_.push_back(result.getSummedTightIso());
36 bdtFeatures_.push_back(result.getMaxCellDep());
37 bdtFeatures_.push_back(result.getShowerRMS());
38 bdtFeatures_.push_back(result.getXStd());
39 bdtFeatures_.push_back(result.getYStd());
40 bdtFeatures_.push_back(result.getAvgLayerHit());
41 bdtFeatures_.push_back(result.getStdLayerHit());
42 bdtFeatures_.push_back(result.getDeepestLayerHit());
43 bdtFeatures_.push_back(result.getEcalBackEnergy());
44 // MIP tracking
45 bdtFeatures_.push_back(result.getNStraightTracks());
46 // bdtFeatures_.push_back(result.getNLinregTracks());
47 bdtFeatures_.push_back(result.getFirstNearPhLayer());
48 bdtFeatures_.push_back(result.getNNearPhHits());
49 bdtFeatures_.push_back(result.getPhotonTerritoryHits());
50 bdtFeatures_.push_back(result.getEPSep());
51 bdtFeatures_.push_back(result.getEPDot());
52 // Longitudinal segment variables
53 bdtFeatures_.push_back(result.getEnergySeg()[0]);
54 bdtFeatures_.push_back(result.getXMeanSeg()[0]);
55 bdtFeatures_.push_back(result.getYMeanSeg()[0]);
56 bdtFeatures_.push_back(result.getLayerMeanSeg()[0]);
57 bdtFeatures_.push_back(result.getEnergySeg()[1]);
58 bdtFeatures_.push_back(result.getYMeanSeg()[2]);
60 bdtFeatures_.push_back(result.getEleContEnergy()[0][0]);
61 bdtFeatures_.push_back(result.getEleContEnergy()[1][0]);
62 bdtFeatures_.push_back(result.getEleContYMean()[0][0]);
63 bdtFeatures_.push_back(result.getEleContEnergy()[0][1]);
64 bdtFeatures_.push_back(result.getEleContEnergy()[1][1]);
65 bdtFeatures_.push_back(result.getEleContYMean()[0][1]);
67 bdtFeatures_.push_back(result.getPhContNHits()[0][0]);
68 bdtFeatures_.push_back(result.getPhContYMean()[0][0]);
69 bdtFeatures_.push_back(result.getPhContNHits()[0][1]);
71 bdtFeatures_.push_back(result.getOutContEnergy()[0][0]);
72 bdtFeatures_.push_back(result.getOutContEnergy()[1][0]);
73 bdtFeatures_.push_back(result.getOutContEnergy()[2][0]);
74 bdtFeatures_.push_back(result.getOutContNHits()[0][0]);
75 bdtFeatures_.push_back(result.getOutContXMean()[0][0]);
76 bdtFeatures_.push_back(result.getOutContYMean()[0][0]);
77 bdtFeatures_.push_back(result.getOutContYMean()[1][0]);
78 bdtFeatures_.push_back(result.getOutContYStd()[0][0]);
79 bdtFeatures_.push_back(result.getOutContEnergy()[0][1]);
80 bdtFeatures_.push_back(result.getOutContEnergy()[1][1]);
81 bdtFeatures_.push_back(result.getOutContEnergy()[2][1]);
82 bdtFeatures_.push_back(result.getOutContLayerMean()[0][1]);
83 bdtFeatures_.push_back(result.getOutContLayerStd()[0][1]);
84 bdtFeatures_.push_back(result.getOutContEnergy()[0][2]);
85 bdtFeatures_.push_back(result.getOutContLayerMean()[0][2]);
86}
int getNStraightTracks() const
Number of straight tracks found.

References ldmx::EcalVetoResult::getNStraightTracks().

Referenced by produce().

◆ clearProcessor()

void ecal::EcalVetoProcessor::clearProcessor ( )
private

Definition at line 141 of file EcalVetoProcessor.cxx.

141 {
142 cellMap_.clear();
143 cellMapTightIso_.clear();
144 bdtFeatures_.clear();
145
146 nReadoutHits_ = 0;
147 summedDet_ = 0;
148 summedTightIso_ = 0;
149 maxCellDep_ = 0;
150 showerRMS_ = 0;
151 xStd_ = 0;
152 yStd_ = 0;
153 avgLayerHit_ = 0;
154 stdLayerHit_ = 0;
155 deepestLayerHit_ = 0;
156 ecalBackEnergy_ = 0;
157 // MIP tracking
159 nLinregTracks_ = 0;
161 nNearPhHits_ = 0;
162 epAng_ = 0;
163 epSep_ = 0;
164 epDot_ = 0;
166
167 std::fill(ecalLayerEdepRaw_.begin(), ecalLayerEdepRaw_.end(), 0);
168 std::fill(ecalLayerEdepReadout_.begin(), ecalLayerEdepReadout_.end(), 0);
169 std::fill(ecalLayerTime_.begin(), ecalLayerTime_.end(), 0);
170}
float epSep_
Distance between the projected photon and electron trajectories at the ECal face.
int photonTerritoryHits_
Number of hits in the photon territory.
float epAng_
Angular separation between the projected photon and electron trajectories (currently unused)
int nStraightTracks_
Number of "straight" tracks found in the event.
int nLinregTracks_
Number of "linreg" tracks found in the event.
int firstNearPhLayer_
First ECal layer in which a hit is found near the photon.
float epDot_
Dot product of the photon and electron momenta unit vectors.
int nNearPhHits_
Number of hits near the photon trajectory.

◆ configure()

void ecal::EcalVetoProcessor::configure ( framework::config::Parameters & parameters)
overridevirtual

Configure the processor using the given user specified parameters.

Parameters
parametersSet of parameters used to configure this processor.

Reimplemented from framework::EventProcessor.

Definition at line 88 of file EcalVetoProcessor.cxx.

88 {
89 verbose_ = parameters.getParameter<bool>("verbose");
90 featureListName_ = parameters.getParameter<std::string>("feature_list_name");
91 // Load BDT ONNX file
92 rt_ = std::make_unique<ldmx::Ort::ONNXRuntime>(
93 parameters.getParameter<std::string>("bdt_file"));
94
95 // Read in arrays holding 68% containment radius per layer
96 // for different bins in momentum/angle
97 rocFileName_ = parameters.getParameter<std::string>("roc_file");
98 if (!std::ifstream(rocFileName_).good()) {
99 EXCEPTION_RAISE(
100 "EcalVetoProcessor",
101 "The specified RoC file '" + rocFileName_ + "' does not exist!");
102 } else {
103 std::ifstream rocfile(rocFileName_);
104 std::string line, value;
105
106 // Extract the first line in the file
107 std::getline(rocfile, line);
108 std::vector<double> values;
109
110 // Read data, line by line
111 while (std::getline(rocfile, line)) {
112 std::stringstream ss(line);
113 values.clear();
114 while (std::getline(ss, value, ',')) {
115 double f_value = (value != "") ? std::stof(value) : -1.0;
116 values.push_back(f_value);
117 }
118 roc_range_values_.push_back(values);
119 }
120 }
121
122 nEcalLayers_ = parameters.getParameter<int>("num_ecal_layers");
123
124 bdtCutVal_ = parameters.getParameter<double>("disc_cut");
125 ecalLayerEdepRaw_.resize(nEcalLayers_, 0);
126 ecalLayerEdepReadout_.resize(nEcalLayers_, 0);
127 ecalLayerTime_.resize(nEcalLayers_, 0);
128
129 beamEnergyMeV_ = parameters.getParameter<double>("beam_energy");
130 linreg_radius_ = parameters.getParameter<double>("linreg_radius");
131
132 // Set the collection name as defined in the configuration
133 collectionName_ = parameters.getParameter<std::string>("collection_name");
134 rec_pass_name_ = parameters.getParameter<std::string>("rec_pass_name");
135 rec_coll_name_ = parameters.getParameter<std::string>("rec_coll_name");
136 recoil_from_tracking_ = parameters.getParameter<bool>("recoil_from_tracking");
137 track_collection_ = parameters.getParameter<std::string>("track_collection");
138 inverse_skim_ = parameters.getParameter<bool>("inverse_skim");
139}
std::string collectionName_
Name of the collection which will containt the results.
T getParameter(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:89

References collectionName_, and framework::config::Parameters::getParameter().

◆ distPtToLine()

float ecal::EcalVetoProcessor::distPtToLine ( TVector3 h1,
TVector3 p1,
TVector3 p2 )
private

Return the minimum distance between the point h1 and the line passing through points p1 and p2.

Parameters
[in]h1Point to find the distance to
[in]p1An arbitrary point on the line
[in]p2A second, distinct point on the line
Returns
Minimum distance between h1 and the line

Definition at line 1302 of file EcalVetoProcessor.cxx.

1302 {
1303 return ((h1 - p1).Cross(h1 - p2)).Mag() / (p1 - p2).Mag();
1304}

Referenced by produce().

◆ distTwoLines()

float ecal::EcalVetoProcessor::distTwoLines ( TVector3 v1,
TVector3 v2,
TVector3 w1,
TVector3 w2 )
private

Returns the distance between the lines v and w, with v defined to pass through the points (v1,v2) (and similarly for w).

Parameters
[in]v1An arbitrary point on line v
[in]v2A second, distinct point on line v
[in]w1An arbitrary point on line w
[in]w2A second, distinct point on line w
Returns
Closest distance of approach of lines u and v

Definition at line 1289 of file EcalVetoProcessor.cxx.

1290 {
1291 TVector3 e1 = v1 - v2;
1292 TVector3 e2 = w1 - w2;
1293 TVector3 crs = e1.Cross(e2);
1294 if (crs.Mag() == 0) {
1295 return 100.0; // arbitrary large number; edge case that shouldn't cause
1296 // problems.
1297 } else {
1298 return std::abs(crs.Dot(v1 - w1) / crs.Mag());
1299 }
1300}

Referenced by produce().

◆ fillHitMap()

void ecal::EcalVetoProcessor::fillHitMap ( const std::vector< ldmx::EcalHit > & ecalRecHits,
std::map< ldmx::EcalID, float > & cellMap_ )
private

Function to load up empty vector of hit maps.

Definition at line 1220 of file EcalVetoProcessor.cxx.

1222 {
1223 for (const ldmx::EcalHit &hit : ecalRecHits) {
1224 ldmx::EcalID id(hit.getID());
1225 cellMap.emplace(id, hit.getEnergy());
1226 }
1227}
Stores reconstructed hit information from the ECAL.
Definition EcalHit.h:19
Extension of DetectorID providing access to ECal layers and cell numbers in a hex grid.
Definition EcalID.h:20

Referenced by produce().

◆ fillIsolatedHitMap()

void ecal::EcalVetoProcessor::fillIsolatedHitMap ( const std::vector< ldmx::EcalHit > & ecalRecHits,
ldmx::EcalID globalCentroid,
std::map< ldmx::EcalID, float > & cellMap,
std::map< ldmx::EcalID, float > & cellMapIso,
bool doTight = false )
private

Definition at line 1229 of file EcalVetoProcessor.cxx.

1232 {
1233 for (const ldmx::EcalHit &hit : ecalRecHits) {
1234 auto isolatedHit = std::make_pair(true, ldmx::EcalID());
1235 ldmx::EcalID id(hit.getID());
1236 if (doTight) {
1237 // Disregard hits that are on the centroid.
1238 if (id == globalCentroid) continue;
1239
1240 // Skip hits that are on centroid inner ring
1241 if (geometry_->isNN(globalCentroid, id)) {
1242 continue;
1243 }
1244 }
1245
1246 // Skip hits that have a readout neighbor
1247 // Get neighboring cell id's and try to look them up in the full cell map
1248 // (constant speed algo.)
1249 // these ideas are only cell/module (must ignore layer)
1250 std::vector<ldmx::EcalID> cellNbrIds = geometry_->getNN(id);
1251
1252 for (int k = 0; k < cellNbrIds.size(); k++) {
1253 // update neighbor ID to the current layer
1254 cellNbrIds[k] = ldmx::EcalID(id.layer(), cellNbrIds[k].module(),
1255 cellNbrIds[k].cell());
1256 // look in cell hit map to see if it is there
1257 if (cellMap.find(cellNbrIds[k]) != cellMap.end()) {
1258 isolatedHit = std::make_pair(false, cellNbrIds[k]);
1259 break;
1260 }
1261 }
1262 if (!isolatedHit.first) {
1263 continue;
1264 }
1265 // Insert isolated hit
1266 cellMapIso.emplace(id, hit.getEnergy());
1267 }
1268}
const ldmx::EcalGeometry * geometry_
handle to current geometry (to share with member functions)
std::vector< EcalID > getNN(EcalID id) const
Get the Nearest Neighbors of the input ID.
bool isNN(EcalID centroid, EcalID probe) const
Check if the probe id is one of the nearest neightbors of the centroid id.

◆ GetShowerCentroidIDAndRMS()

ldmx::EcalID ecal::EcalVetoProcessor::GetShowerCentroidIDAndRMS ( const std::vector< ldmx::EcalHit > & ecalRecHits,
double & showerRMS )
private

Definition at line 1173 of file EcalVetoProcessor.cxx.

1174 {
1175 auto wgtCentroidCoords = std::make_pair<float, float>(0., 0.);
1176 float sumEdep = 0;
1177 ldmx::EcalID returnCellId;
1178
1179 // Calculate Energy Weighted Centroid
1180 for (const ldmx::EcalHit &hit : ecalRecHits) {
1181 ldmx::EcalID id(hit.getID());
1182 CellEnergyPair cell_energy_pair = std::make_pair(id, hit.getEnergy());
1183 auto [x, y, z] = geometry_->getPosition(id);
1184 XYCoords centroidCoords = std::make_pair(x, y);
1185 wgtCentroidCoords.first = wgtCentroidCoords.first +
1186 centroidCoords.first * cell_energy_pair.second;
1187 wgtCentroidCoords.second = wgtCentroidCoords.second +
1188 centroidCoords.second * cell_energy_pair.second;
1189 sumEdep += cell_energy_pair.second;
1190 }
1191 wgtCentroidCoords.first = (sumEdep > 1E-6) ? wgtCentroidCoords.first / sumEdep
1192 : wgtCentroidCoords.first;
1193 wgtCentroidCoords.second = (sumEdep > 1E-6)
1194 ? wgtCentroidCoords.second / sumEdep
1195 : wgtCentroidCoords.second;
1196 // Find Nearest Cell to Centroid
1197 float maxDist = 1e6;
1198 for (const ldmx::EcalHit &hit : ecalRecHits) {
1199 auto [x, y, z] = geometry_->getPosition(hit.getID());
1200 XYCoords centroidCoords = std::make_pair(x, y);
1201
1202 float deltaR =
1203 pow(pow((centroidCoords.first - wgtCentroidCoords.first), 2) +
1204 pow((centroidCoords.second - wgtCentroidCoords.second), 2),
1205 .5);
1206 showerRMS += deltaR * hit.getEnergy();
1207 if (deltaR < maxDist) {
1208 maxDist = deltaR;
1209 returnCellId = ldmx::EcalID(hit.getID());
1210 }
1211 }
1212 if (sumEdep > 0) showerRMS = showerRMS / sumEdep;
1213 // flatten
1214 return ldmx::EcalID(0, returnCellId.module(), returnCellId.cell());
1215}
std::tuple< double, double, double > getPosition(EcalID id) const
Get a cell's position from its ID number.
int cell() const
Get the value of the cell field from the ID.
Definition EcalID.h:111
int module() const
Get the value of the module field from the ID.
Definition EcalID.h:87

◆ getTrajectory()

std::vector< std::pair< float, float > > ecal::EcalVetoProcessor::getTrajectory ( std::vector< double > momentum,
std::vector< float > position )
private

Definition at line 1272 of file EcalVetoProcessor.cxx.

1273 {
1274 std::vector<XYCoords> positions;
1275 for (int iLayer = 0; iLayer < nEcalLayers_; iLayer++) {
1276 float posX =
1277 position[0] + (momentum[0] / momentum[2]) *
1278 (geometry_->getZPosition(iLayer) - position[2]);
1279 float posY =
1280 position[1] + (momentum[1] / momentum[2]) *
1281 (geometry_->getZPosition(iLayer) - position[2]);
1282 positions.push_back(std::make_pair(posX, posY));
1283 }
1284 return positions;
1285}
double getZPosition(int layer) const
Get the z-coordinate given the layer id.

◆ produce()

void ecal::EcalVetoProcessor::produce ( framework::Event & event)
overridevirtual

Process the event and put new data products into it.

Parameters
eventThe Event to process.

Implements framework::Producer.

Definition at line 172 of file EcalVetoProcessor.cxx.

172 {
173 // Get the Ecal Geometry
175 ldmx::EcalGeometry::CONDITIONS_OBJECT_NAME);
176
178
179 clearProcessor();
180
181 // Get the collection of Ecal scoring plane hits. If it doesn't exist,
182 // don't bother adding any truth tracking information.
183
184 std::vector<double> recoilP;
185 std::vector<float> recoilPos;
186 std::vector<double> recoilPAtTarget;
187 std::vector<float> recoilPosAtTarget;
188
189 if (verbose_) {
190 ldmx_log(debug) << " Loop through all of the sim particles and find the "
191 "recoil electron";
192 }
193
194 if (event.exists("EcalScoringPlaneHits")) {
195 //
196 // Loop through all of the sim particles and find the recoil electron.
197 //
198
199 // Get the collection of simulated particles from the event
200 auto particleMap{event.getMap<int, ldmx::SimParticle>("SimParticles")};
201
202 // Search for the recoil electron
203 auto [recoilTrackID, recoilElectron] = Analysis::getRecoil(particleMap);
204
205 // Find ECAL SP hit for recoil electron
206 auto ecalSpHits{
207 event.getCollection<ldmx::SimTrackerHit>("EcalScoringPlaneHits")};
208 float pmax = 0;
209 for (ldmx::SimTrackerHit &spHit : ecalSpHits) {
210 ldmx::SimSpecialID hit_id(spHit.getID());
211 if (hit_id.plane() != 31 || spHit.getMomentum()[2] <= 0) continue;
212
213 if (spHit.getTrackID() == recoilTrackID) {
214 if (sqrt(pow(spHit.getMomentum()[0], 2) +
215 pow(spHit.getMomentum()[1], 2) +
216 pow(spHit.getMomentum()[2], 2)) > pmax) {
217 recoilP = spHit.getMomentum();
218 recoilPos = spHit.getPosition();
219 pmax = sqrt(pow(recoilP[0], 2) + pow(recoilP[1], 2) +
220 pow(recoilP[2], 2));
221 }
222 }
223 }
224
225 // Find target SP hit for recoil electron
226 if (event.exists("TargetScoringPlaneHits")) {
227 std::vector<ldmx::SimTrackerHit> targetSpHits =
228 event.getCollection<ldmx::SimTrackerHit>("TargetScoringPlaneHits");
229 pmax = 0;
230 for (ldmx::SimTrackerHit &spHit : targetSpHits) {
231 ldmx::SimSpecialID hit_id(spHit.getID());
232 if (hit_id.plane() != 1 || spHit.getMomentum()[2] <= 0) continue;
233
234 if (spHit.getTrackID() == recoilTrackID) {
235 if (sqrt(pow(spHit.getMomentum()[0], 2) +
236 pow(spHit.getMomentum()[1], 2) +
237 pow(spHit.getMomentum()[2], 2)) > pmax) {
238 recoilPAtTarget = spHit.getMomentum();
239 recoilPosAtTarget = spHit.getPosition();
240 pmax =
241 sqrt(pow(recoilPAtTarget[0], 2) + pow(recoilPAtTarget[1], 2) +
242 pow(recoilPAtTarget[2], 2));
243 }
244 }
245 }
246 }
247 }
248 // Get recoilPos using recoil tracking
249 if (recoil_from_tracking_) {
250 std::vector<float> recoil_track_states;
251 if (verbose_) {
252 ldmx_log(debug) << " Propagate recoil tracks to ECAL face";
253 }
254 // Get the recoil track collection
255 auto recoil_tracks{event.getCollection<ldmx::Track>(track_collection_)};
256
257 ldmx::TrackStateType ts_type = ldmx::TrackStateType::AtECAL;
258 recoil_track_states = trackProp(recoil_tracks, ts_type, "ecal");
259 // Redefining recoilPos now to come from the track state
260 // track_state_loc0 is recoilPos[0] and track_state_loc1 is recoilPos[1]
261 recoilPos = recoil_track_states;
262 }
263
264 if (verbose_) {
265 ldmx_log(debug) << " Get projected trajectories for electron and photon";
266 }
267 // Get projected trajectories for electron and photon
268 std::vector<XYCoords> ele_trajectory, photon_trajectory;
269 if (!recoilP.empty() && !recoilPos.empty()) {
270 ele_trajectory = getTrajectory(recoilP, recoilPos);
271 std::vector<double> pvec = recoilPAtTarget.size()
272 ? recoilPAtTarget
273 : std::vector<double>{0.0, 0.0, 0.0};
274 std::vector<float> posvec = recoilPosAtTarget.size()
275 ? recoilPosAtTarget
276 : std::vector<float>{0.0, 0.0, 0.0};
277 photon_trajectory =
278 getTrajectory({-pvec[0], -pvec[1], beamEnergyMeV_ - pvec[2]}, posvec);
279 }
280
281 float recoilPMag =
282 recoilP.size()
283 ? sqrt(pow(recoilP[0], 2) + pow(recoilP[1], 2) + pow(recoilP[2], 2))
284 : -1.0;
285 float recoilTheta =
286 recoilPMag > 0 ? acos(recoilP[2] / recoilPMag) * 180.0 / M_PI : -1.0;
287
288 if (verbose_) {
289 ldmx_log(debug) << " Build Radii of containment (ROC)";
290 }
291 // Use the appropriate containment radii for the recoil electron
292 std::vector<double> roc_values_bin0(roc_range_values_[0].begin() + 4,
293 roc_range_values_[0].end());
294 std::vector<double> ele_radii = roc_values_bin0;
295 double theta_min, theta_max, p_min, p_max;
296 bool inrange;
297
298 for (int i = 0; i < roc_range_values_.size(); i++) {
299 theta_min = roc_range_values_[i][0];
300 theta_max = roc_range_values_[i][1];
301 p_min = roc_range_values_[i][2];
302 p_max = roc_range_values_[i][3];
303 inrange = true;
304
305 if (theta_min != -1.0) {
306 inrange = inrange && (recoilTheta >= theta_min);
307 }
308 if (theta_max != -1.0) {
309 inrange = inrange && (recoilTheta < theta_max);
310 }
311 if (p_min != -1.0) {
312 inrange = inrange && (recoilPMag >= p_min);
313 }
314 if (p_max != -1.0) {
315 inrange = inrange && (recoilPMag < p_max);
316 }
317 if (inrange) {
318 std::vector<double> roc_values_bini(roc_range_values_[i].begin() + 4,
319 roc_range_values_[i].end());
320 ele_radii = roc_values_bini;
321 }
322 }
323 // Use default RoC bin for photon
324 std::vector<double> photon_radii = roc_values_bin0;
325
326 // Get the collection of digitized Ecal hits from the event.
327 const std::vector<ldmx::EcalHit> ecalRecHits =
328 event.getCollection<ldmx::EcalHit>(rec_coll_name_, rec_pass_name_);
329
330 ldmx::EcalID globalCentroid =
331 GetShowerCentroidIDAndRMS(ecalRecHits, showerRMS_);
332 /* ~~ Fill the hit map ~~ O(n) */
333 fillHitMap(ecalRecHits, cellMap_);
334 bool doTight = true;
335 /* ~~ Fill the isolated hit maps ~~ O(n) */
336 fillIsolatedHitMap(ecalRecHits, globalCentroid, cellMap_, cellMapTightIso_,
337 doTight);
338
339 // Loop over the hits from the event to calculate the rest of the important
340 // quantities
341
342 float wavgLayerHit = 0;
343 float xMean = 0;
344 float yMean = 0;
345
346 // Containment variables
347 unsigned int nregions = 5;
348 std::vector<float> electronContainmentEnergy(nregions, 0.0);
349 std::vector<float> photonContainmentEnergy(nregions, 0.0);
350 std::vector<float> outsideContainmentEnergy(nregions, 0.0);
351 std::vector<int> outsideContainmentNHits(nregions, 0);
352 std::vector<float> outsideContainmentXmean(nregions, 0.0);
353 std::vector<float> outsideContainmentYmean(nregions, 0.0);
354 std::vector<float> outsideContainmentXstd(nregions, 0.0);
355 std::vector<float> outsideContainmentYstd(nregions, 0.0);
356 // Longitudinal segmentation
357 std::vector<int> segLayers = {0, 6, 17, 34};
358 unsigned int nsegments = segLayers.size() - 1;
359 std::vector<float> energySeg(nsegments, 0.0);
360 std::vector<float> xMeanSeg(nsegments, 0.0);
361 std::vector<float> xStdSeg(nsegments, 0.0);
362 std::vector<float> yMeanSeg(nsegments, 0.0);
363 std::vector<float> yStdSeg(nsegments, 0.0);
364 std::vector<float> layerMeanSeg(nsegments, 0.0);
365 std::vector<float> layerStdSeg(nsegments, 0.0);
366 std::vector<std::vector<float>> eContEnergy(
367 nregions, std::vector<float>(nsegments, 0.0));
368 std::vector<std::vector<float>> eContXMean(
369 nregions, std::vector<float>(nsegments, 0.0));
370 std::vector<std::vector<float>> eContYMean(
371 nregions, std::vector<float>(nsegments, 0.0));
372 std::vector<std::vector<float>> gContEnergy(
373 nregions, std::vector<float>(nsegments, 0.0));
374 std::vector<std::vector<int>> gContNHits(nregions,
375 std::vector<int>(nsegments, 0));
376 std::vector<std::vector<float>> gContXMean(
377 nregions, std::vector<float>(nsegments, 0.0));
378 std::vector<std::vector<float>> gContYMean(
379 nregions, std::vector<float>(nsegments, 0.0));
380 std::vector<std::vector<float>> oContEnergy(
381 nregions, std::vector<float>(nsegments, 0.0));
382 std::vector<std::vector<int>> oContNHits(nregions,
383 std::vector<int>(nsegments, 0));
384 std::vector<std::vector<float>> oContXMean(
385 nregions, std::vector<float>(nsegments, 0.0));
386 std::vector<std::vector<float>> oContYMean(
387 nregions, std::vector<float>(nsegments, 0.0));
388 std::vector<std::vector<float>> oContXStd(nregions,
389 std::vector<float>(nsegments, 0.0));
390 std::vector<std::vector<float>> oContYStd(nregions,
391 std::vector<float>(nsegments, 0.0));
392 std::vector<std::vector<float>> oContLayerMean(
393 nregions, std::vector<float>(nsegments, 0.0));
394 std::vector<std::vector<float>> oContLayerStd(
395 nregions, std::vector<float>(nsegments, 0.0));
396
397 // MIP tracking: vector of hits to be used in the MIP tracking algorithm. All
398 // hits inside the electron ROC (or all hits in the ECal if the event is
399 // missing an electron) will be included.
400 std::vector<HitData> trackingHitList;
401
402 if (verbose_) {
403 ldmx_log(debug)
404 << " Loop over the hits from the event to calculate the BDT features";
405 }
406
407 for (const ldmx::EcalHit &hit : ecalRecHits) {
408 // Layer-wise quantities
409 ldmx::EcalID id(hit.getID());
410 ecalLayerEdepRaw_[id.layer()] =
411 ecalLayerEdepRaw_[id.layer()] + hit.getEnergy();
412 if (id.layer() >= 20) ecalBackEnergy_ += hit.getEnergy();
413 if (maxCellDep_ < hit.getEnergy()) maxCellDep_ = hit.getEnergy();
414 if (hit.getEnergy() <= 0) {
415 ldmx_log(fatal)
416 << "ECal hit has negative or zero energy, this should never happen, "
417 "check the threshold settings in HgcrocEmulator";
418 continue;
419 }
420 nReadoutHits_++;
421 ecalLayerEdepReadout_[id.layer()] += hit.getEnergy();
422 ecalLayerTime_[id.layer()] += (hit.getEnergy()) * hit.getTime();
423 auto [x, y, z] = geometry_->getPosition(id);
424 xMean += x * hit.getEnergy();
425 yMean += y * hit.getEnergy();
426 avgLayerHit_ += id.layer();
427 wavgLayerHit += id.layer() * hit.getEnergy();
428 if (deepestLayerHit_ < id.layer()) {
429 deepestLayerHit_ = id.layer();
430 }
431 XYCoords xy_pair = std::make_pair(x, y);
432 float distance_ele_trajectory =
433 ele_trajectory.size()
434 ? sqrt(pow((xy_pair.first - ele_trajectory[id.layer()].first), 2) +
435 pow((xy_pair.second - ele_trajectory[id.layer()].second), 2))
436 : -1.0;
437 float distance_photon_trajectory =
438 photon_trajectory.size()
439 ? sqrt(pow((xy_pair.first - photon_trajectory[id.layer()].first),
440 2) +
441 pow((xy_pair.second - photon_trajectory[id.layer()].second),
442 2))
443 : -1.0;
444
445 // Decide which longitudinal segment the hit is in and add to sums
446 for (unsigned int iseg = 0; iseg < nsegments; iseg++) {
447 if (id.layer() >= segLayers[iseg] &&
448 id.layer() <= segLayers[iseg + 1] - 1) {
449 energySeg[iseg] += hit.getEnergy();
450 xMeanSeg[iseg] += xy_pair.first * hit.getEnergy();
451 yMeanSeg[iseg] += xy_pair.second * hit.getEnergy();
452 layerMeanSeg[iseg] += id.layer() * hit.getEnergy();
453
454 // Decide which containment region the hit is in and add to sums
455 for (unsigned int ireg = 0; ireg < nregions; ireg++) {
456 if (distance_ele_trajectory >= ireg * ele_radii[id.layer()] &&
457 distance_ele_trajectory < (ireg + 1) * ele_radii[id.layer()]) {
458 eContEnergy[ireg][iseg] += hit.getEnergy();
459 eContXMean[ireg][iseg] += xy_pair.first * hit.getEnergy();
460 eContYMean[ireg][iseg] += xy_pair.second * hit.getEnergy();
461 }
462 if (distance_photon_trajectory >= ireg * photon_radii[id.layer()] &&
463 distance_photon_trajectory <
464 (ireg + 1) * photon_radii[id.layer()]) {
465 gContEnergy[ireg][iseg] += hit.getEnergy();
466 gContNHits[ireg][iseg] += 1;
467 gContXMean[ireg][iseg] += xy_pair.first * hit.getEnergy();
468 gContYMean[ireg][iseg] += xy_pair.second * hit.getEnergy();
469 }
470 if (distance_ele_trajectory > (ireg + 1) * ele_radii[id.layer()] &&
471 distance_photon_trajectory >
472 (ireg + 1) * photon_radii[id.layer()]) {
473 oContEnergy[ireg][iseg] += hit.getEnergy();
474 oContNHits[ireg][iseg] += 1;
475 oContXMean[ireg][iseg] += xy_pair.first * hit.getEnergy();
476 oContYMean[ireg][iseg] += xy_pair.second * hit.getEnergy();
477 oContLayerMean[ireg][iseg] += id.layer() * hit.getEnergy();
478 }
479 }
480 }
481 }
482
483 // Decide which containment region the hit is in and add to sums
484 for (unsigned int ireg = 0; ireg < nregions; ireg++) {
485 if (distance_ele_trajectory >= ireg * ele_radii[id.layer()] &&
486 distance_ele_trajectory < (ireg + 1) * ele_radii[id.layer()])
487 electronContainmentEnergy[ireg] += hit.getEnergy();
488 if (distance_photon_trajectory >= ireg * photon_radii[id.layer()] &&
489 distance_photon_trajectory < (ireg + 1) * photon_radii[id.layer()])
490 photonContainmentEnergy[ireg] += hit.getEnergy();
491 if (distance_ele_trajectory > (ireg + 1) * ele_radii[id.layer()] &&
492 distance_photon_trajectory > (ireg + 1) * photon_radii[id.layer()]) {
493 outsideContainmentEnergy[ireg] += hit.getEnergy();
494 outsideContainmentNHits[ireg] += 1;
495 outsideContainmentXmean[ireg] += xy_pair.first * hit.getEnergy();
496 outsideContainmentYmean[ireg] += xy_pair.second * hit.getEnergy();
497 }
498 }
499
500 // MIP tracking: Decide whether hit should be added to trackingHitList
501 // If inside e- RoC or if etraj is missing, use the hit for tracking:
502 if (distance_ele_trajectory >= ele_radii[id.layer()] ||
503 distance_ele_trajectory == -1.0) {
504 HitData hd;
505 hd.pos = TVector3(xy_pair.first, xy_pair.second,
506 geometry_->getZPosition(id.layer()));
507 hd.layer = id.layer();
508 trackingHitList.push_back(hd);
509 }
510 } // end loop over rechits
511
512 for (const auto &[id, energy] : cellMapTightIso_) {
513 if (energy > 0) summedTightIso_ += energy;
514 }
515
516 for (int iLayer = 0; iLayer < ecalLayerEdepReadout_.size(); iLayer++) {
517 ecalLayerTime_[iLayer] =
518 ecalLayerTime_[iLayer] / ecalLayerEdepReadout_[iLayer];
519 summedDet_ += ecalLayerEdepReadout_[iLayer];
520 }
521
522 if (nReadoutHits_ > 0) {
523 avgLayerHit_ /= nReadoutHits_;
524 wavgLayerHit /= summedDet_;
525 xMean /= summedDet_;
526 yMean /= summedDet_;
527 } else {
528 wavgLayerHit = 0;
529 avgLayerHit_ = 0;
530 xMean = 0;
531 yMean = 0;
532 }
533
534 // If necessary, quotient out the total energy from the means
535 for (unsigned int iseg = 0; iseg < nsegments; iseg++) {
536 if (energySeg[iseg] > 0) {
537 xMeanSeg[iseg] /= energySeg[iseg];
538 yMeanSeg[iseg] /= energySeg[iseg];
539 layerMeanSeg[iseg] /= energySeg[iseg];
540 }
541 for (unsigned int ireg = 0; ireg < nregions; ireg++) {
542 if (eContEnergy[ireg][iseg] > 0) {
543 eContXMean[ireg][iseg] /= eContEnergy[ireg][iseg];
544 eContYMean[ireg][iseg] /= eContEnergy[ireg][iseg];
545 }
546 if (gContEnergy[ireg][iseg] > 0) {
547 gContXMean[ireg][iseg] /= gContEnergy[ireg][iseg];
548 gContYMean[ireg][iseg] /= gContEnergy[ireg][iseg];
549 }
550 if (oContEnergy[ireg][iseg] > 0) {
551 oContXMean[ireg][iseg] /= oContEnergy[ireg][iseg];
552 oContYMean[ireg][iseg] /= oContEnergy[ireg][iseg];
553 oContLayerMean[ireg][iseg] /= oContEnergy[ireg][iseg];
554 }
555 }
556 }
557
558 for (unsigned int ireg = 0; ireg < nregions; ireg++) {
559 if (outsideContainmentEnergy[ireg] > 0) {
560 outsideContainmentXmean[ireg] /= outsideContainmentEnergy[ireg];
561 outsideContainmentYmean[ireg] /= outsideContainmentEnergy[ireg];
562 }
563 }
564
565 // Loop over hits a second time to find the standard deviations.
566 for (const ldmx::EcalHit &hit : ecalRecHits) {
567 ldmx::EcalID id(hit.getID());
568 auto [x, y, z] = geometry_->getPosition(id);
569 if (hit.getEnergy() > 0) {
570 xStd_ += pow((x - xMean), 2) * hit.getEnergy();
571 yStd_ += pow((y - yMean), 2) * hit.getEnergy();
572 stdLayerHit_ += pow((id.layer() - wavgLayerHit), 2) * hit.getEnergy();
573 }
574 XYCoords xy_pair = std::make_pair(x, y);
575 float distance_ele_trajectory =
576 ele_trajectory.size()
577 ? sqrt(pow((xy_pair.first - ele_trajectory[id.layer()].first), 2) +
578 pow((xy_pair.second - ele_trajectory[id.layer()].second), 2))
579 : -1.0;
580 float distance_photon_trajectory =
581 photon_trajectory.size()
582 ? sqrt(pow((xy_pair.first - photon_trajectory[id.layer()].first),
583 2) +
584 pow((xy_pair.second - photon_trajectory[id.layer()].second),
585 2))
586 : -1.0;
587
588 for (unsigned int iseg = 0; iseg < nsegments; iseg++) {
589 if (id.layer() >= segLayers[iseg] &&
590 id.layer() <= segLayers[iseg + 1] - 1) {
591 xStdSeg[iseg] +=
592 pow((xy_pair.first - xMeanSeg[iseg]), 2) * hit.getEnergy();
593 yStdSeg[iseg] +=
594 pow((xy_pair.second - yMeanSeg[iseg]), 2) * hit.getEnergy();
595 layerStdSeg[iseg] +=
596 pow((id.layer() - layerMeanSeg[iseg]), 2) * hit.getEnergy();
597
598 for (unsigned int ireg = 0; ireg < nregions; ireg++) {
599 if (distance_ele_trajectory > (ireg + 1) * ele_radii[id.layer()] &&
600 distance_photon_trajectory >
601 (ireg + 1) * photon_radii[id.layer()]) {
602 oContXStd[ireg][iseg] +=
603 pow((xy_pair.first - oContXMean[ireg][iseg]), 2) *
604 hit.getEnergy();
605 oContYStd[ireg][iseg] +=
606 pow((xy_pair.second - oContYMean[ireg][iseg]), 2) *
607 hit.getEnergy();
608 oContLayerStd[ireg][iseg] +=
609 pow((id.layer() - oContLayerMean[ireg][iseg]), 2) *
610 hit.getEnergy();
611 }
612 }
613 }
614 }
615
616 for (unsigned int ireg = 0; ireg < nregions; ireg++) {
617 if (distance_ele_trajectory > (ireg + 1) * ele_radii[id.layer()] &&
618 distance_photon_trajectory > (ireg + 1) * photon_radii[id.layer()]) {
619 outsideContainmentXstd[ireg] +=
620 pow((xy_pair.first - outsideContainmentXmean[ireg]), 2) *
621 hit.getEnergy();
622 outsideContainmentYstd[ireg] +=
623 pow((xy_pair.second - outsideContainmentYmean[ireg]), 2) *
624 hit.getEnergy();
625 }
626 }
627 } // end loop over rechits (2nd time)
628
629 if (nReadoutHits_ > 0) {
630 xStd_ = sqrt(xStd_ / summedDet_);
631 yStd_ = sqrt(yStd_ / summedDet_);
632 stdLayerHit_ = sqrt(stdLayerHit_ / summedDet_);
633 } else {
634 xStd_ = 0;
635 yStd_ = 0;
636 stdLayerHit_ = 0;
637 }
638
639 // Quotient out the total energies from the standard deviations if possible
640 // and take root
641 for (unsigned int iseg = 0; iseg < nsegments; iseg++) {
642 if (energySeg[iseg] > 0) {
643 xStdSeg[iseg] = sqrt(xStdSeg[iseg] / energySeg[iseg]);
644 yStdSeg[iseg] = sqrt(yStdSeg[iseg] / energySeg[iseg]);
645 layerStdSeg[iseg] = sqrt(layerStdSeg[iseg] / energySeg[iseg]);
646 }
647 for (unsigned int ireg = 0; ireg < nregions; ireg++) {
648 if (oContEnergy[ireg][iseg] > 0) {
649 oContXStd[ireg][iseg] =
650 sqrt(oContXStd[ireg][iseg] / oContEnergy[ireg][iseg]);
651 oContYStd[ireg][iseg] =
652 sqrt(oContYStd[ireg][iseg] / oContEnergy[ireg][iseg]);
653 oContLayerStd[ireg][iseg] =
654 sqrt(oContLayerStd[ireg][iseg] / oContEnergy[ireg][iseg]);
655 }
656 }
657 }
658
659 for (unsigned int ireg = 0; ireg < nregions; ireg++) {
660 if (outsideContainmentEnergy[ireg] > 0) {
661 outsideContainmentXstd[ireg] =
662 sqrt(outsideContainmentXstd[ireg] / outsideContainmentEnergy[ireg]);
663 outsideContainmentYstd[ireg] =
664 sqrt(outsideContainmentYstd[ireg] / outsideContainmentEnergy[ireg]);
665 }
666 }
667
668 if (verbose_) {
669 ldmx_log(debug) << " Find out if the recoil electron is fiducial";
670 }
671
672 // Find the location of the recoil electron
673 // Ecal face is not where the first layer starts,
674 // defined in DetDescr/python/EcalGeometry.py
675 const float dz_from_face{7.932};
676 float drifted_recoil_x{-9999.};
677 float drifted_recoil_y{-9999.};
678 if (!recoilP.empty()) {
679 drifted_recoil_x =
680 (dz_from_face * (recoilP[0] / recoilP[2])) + recoilPos[0];
681 drifted_recoil_y =
682 (dz_from_face * (recoilP[1] / recoilP[2])) + recoilPos[1];
683 }
684 const int recoil_layer_index = 0;
685
686 // Check if it's fiducial
687 bool inside{false};
688 // At module level
689 const auto ecalID = geometry_->getID(drifted_recoil_x, drifted_recoil_y,
690 recoil_layer_index, true);
691 if (!ecalID.null()) {
692 // If fiducial at module level, check at cell level
693 const auto cellID =
694 geometry_->getID(drifted_recoil_x, drifted_recoil_y, recoil_layer_index,
695 ecalID.getModuleID(), true);
696 if (!cellID.null()) {
697 inside = true;
698 }
699 }
700
701 if (!inside) {
702 ldmx_log(info) << "This event is non-fiducial in ECAL";
703 }
704
705 // ------------------------------------------------------
706 // MIP tracking starts here
707
708 /* Goal: Calculate
709 * nStraightTracks (self-explanatory),
710 * nLinregTracks (tracks found by linreg algorithm),
711 */
712
713 // Find epAng and epSep, and prepare EP trajectory vectors:
714 TVector3 e_traj_start;
715 TVector3 e_traj_end;
716 TVector3 p_traj_start;
717 TVector3 p_traj_end;
718 if (!ele_trajectory.empty() && !photon_trajectory.empty()) {
719 // Create TVector3s marking the start and endpoints of each projected
720 // trajectory
721 e_traj_start.SetXYZ(ele_trajectory[0].first, ele_trajectory[0].second,
723 e_traj_end.SetXYZ(ele_trajectory[(nEcalLayers_ - 1)].first,
724 ele_trajectory[(nEcalLayers_ - 1)].second,
725 geometry_->getZPosition((nEcalLayers_ - 1)));
726 p_traj_start.SetXYZ(photon_trajectory[0].first, photon_trajectory[0].second,
728 p_traj_end.SetXYZ(photon_trajectory[(nEcalLayers_ - 1)].first,
729 photon_trajectory[(nEcalLayers_ - 1)].second,
730 geometry_->getZPosition((nEcalLayers_ - 1)));
731
732 TVector3 evec = e_traj_end - e_traj_start;
733 TVector3 e_norm = evec.Unit();
734 TVector3 pvec = p_traj_end - p_traj_start;
735 TVector3 p_norm = pvec.Unit();
736 epDot_ = e_norm.Dot(p_norm);
737 epAng_ = acos(epDot_) * 180.0 / M_PI;
738 epSep_ = sqrt(pow(e_traj_start.X() - p_traj_start.X(), 2) +
739 pow(e_traj_start.Y() - p_traj_start.Y(), 2));
740 } else {
741 // Electron trajectory is missing, so all hits in the Ecal are fair game.
742 // Pick e/ptraj so that they won't restrict the tracking algorithm (place
743 // them far outside the ECal).
744 e_traj_start = TVector3(999, 999, geometry_->getZPosition(0)); // 0);
745 e_traj_end = TVector3(
746 999, 999, geometry_->getZPosition((nEcalLayers_ - 1))); // 999);
747 p_traj_start = TVector3(1000, 1000, geometry_->getZPosition(0)); // 0);
748 p_traj_end = TVector3(
749 1000, 1000, geometry_->getZPosition((nEcalLayers_ - 1))); // 1000);
750 /*ensures event will not be vetoed by angle/separation cut */
751 epAng_ = 999.;
752 epSep_ = 999.;
753 epDot_ = 999.;
754 }
755
756 // Near photon step: Find the first layer of the ECal where a hit near the
757 // projected photon trajectory is found Currently unusued pending further
758 // study; performance has dropped between v9 and v12. Currently used in
759 // segmipBDT
760 firstNearPhLayer_ = nEcalLayers_ - 1;
761
762 // If no photon trajectory, leave this at the default (ECal back)
763 if (!photon_trajectory.empty()) {
764 for (std::vector<HitData>::iterator it = trackingHitList.begin();
765 it != trackingHitList.end(); ++it) {
766 float ehDist =
767 sqrt(pow((*it).pos.X() - photon_trajectory[(*it).layer].first, 2) +
768 pow((*it).pos.Y() - photon_trajectory[(*it).layer].second, 2));
769 if (ehDist < 8.7) {
770 nNearPhHits_++;
771 if ((*it).layer < firstNearPhLayer_) {
772 firstNearPhLayer_ = (*it).layer;
773 }
774 }
775 }
776 }
777
778 // Territories limited to trackingHitList
779 TVector3 gToe = (e_traj_start - p_traj_start).Unit();
780 TVector3 origin = p_traj_start + 0.5 * 8.7 * gToe;
781 if (!ele_trajectory.empty()) {
782 for (auto &hitData : trackingHitList) {
783 TVector3 hitPos = hitData.pos;
784 TVector3 hitPrime = hitPos - origin;
785 if (hitPrime.Dot(gToe) <= 0) {
787 }
788 }
789 } else {
790 photonTerritoryHits_ = nReadoutHits_;
791 }
792
793 // ------------------------------------------------------
794 // Find straight MIP tracks:
795
796 std::sort(trackingHitList.begin(), trackingHitList.end(),
797 [](HitData ha, HitData hb) { return ha.layer > hb.layer; });
798 // For merging tracks: Need to keep track of existing tracks
799 // Candidate tracks to merge in will always be in front of the current track
800 // (lower z), so only store the last hit 3-layer vector: each track = vector
801 // of 3-tuples (xy+layer).
802 std::vector<std::vector<HitData>> track_list;
803
804 // print trackingHitList
805 if (verbose_) {
806 ldmx_log(debug) << "====== Tracking hit list (original) length "
807 << trackingHitList.size() << " ======";
808 for (int i = 0; i < trackingHitList.size(); i++) {
809 std::cout << "[" << trackingHitList[i].pos.X() << ", "
810 << trackingHitList[i].pos.Y() << ", "
811 << trackingHitList[i].layer << "], ";
812 }
813 std::cout << std::endl;
814 ldmx_log(debug) << "====== END OF Tracking hit list ======";
815 }
816
817 // in v14 minR is 4.17 mm
818 // while maxR is 4.81 mm
819 float cellWidth = 2 * geometry_->getCellMaxR();
820 for (int iHit = 0; iHit < trackingHitList.size(); iHit++) {
821 // list of hit numbers in track (34 = maximum theoretical length)
822 int track[34];
823 int currenthit{iHit};
824 int trackLen{1};
825
826 track[0] = iHit;
827
828 // Search for hits to add to the track:
829 // repeatedly find hits in the front two layers with same x & y positions
830 // but since v14 the odd layers are offset, so we allow half a cellWidth
831 // deviation and then add to track until no more hits are found
832 int jHit = iHit;
833 while (jHit < trackingHitList.size()) {
834 if ((trackingHitList[jHit].layer ==
835 trackingHitList[currenthit].layer - 1 ||
836 trackingHitList[jHit].layer ==
837 trackingHitList[currenthit].layer - 2) &&
838 abs(trackingHitList[jHit].pos.X() -
839 trackingHitList[currenthit].pos.X()) <= 0.5 * cellWidth &&
840 abs(trackingHitList[jHit].pos.Y() -
841 trackingHitList[currenthit].pos.Y()) <= 0.5 * cellWidth) {
842 track[trackLen] = jHit;
843 trackLen++;
844 currenthit = jHit;
845 }
846 jHit++;
847 }
848
849 // Confirm that the track is valid:
850 if (trackLen < 2) continue; // Track must contain at least 2 hits
851 float closest_e = distTwoLines(trackingHitList[track[0]].pos,
852 trackingHitList[track[trackLen - 1]].pos,
853 e_traj_start, e_traj_end);
854 float closest_p = distTwoLines(trackingHitList[track[0]].pos,
855 trackingHitList[track[trackLen - 1]].pos,
856 p_traj_start, p_traj_end);
857 // Make sure that the track is near the photon trajectory and away from the
858 // electron trajectory Details of these constraints may be revised
859 if (closest_p > cellWidth and closest_e < 2 * cellWidth) continue;
860 if (trackLen < 4 and closest_e > closest_p) continue;
861 if (verbose_) {
862 ldmx_log(debug) << "====== After rejection for MIP tracking ======";
863 ldmx_log(debug) << "current hit: [" << trackingHitList[iHit].pos.X()
864 << ", " << trackingHitList[iHit].pos.Y() << ", "
865 << trackingHitList[iHit].layer << "]";
866
867 for (int k = 0; k < trackLen; k++) {
868 ldmx_log(debug) << "track[" << k << "] position = ["
869 << trackingHitList[track[k]].pos.X() << ", "
870 << trackingHitList[track[k]].pos.Y() << ", "
871 << trackingHitList[track[k]].layer << "]";
872 }
873 }
874
875 // if track found, increment nStraightTracks and remove all hits in track
876 // from future consideration
877 if (trackLen >= 2) {
878 std::vector<HitData> temp_track_list;
879 int n_remove = 0;
880 for (int kHit = 0; kHit < trackLen; kHit++) {
881 temp_track_list.push_back(trackingHitList[track[kHit] - n_remove]);
882 trackingHitList.erase(trackingHitList.begin() + track[kHit] - n_remove);
883 n_remove++;
884 }
885 // print trackingHitList
886 if (verbose_) {
887 ldmx_log(debug) << "====== Tracking hit list (after erase) length "
888 << trackingHitList.size() << " ======";
889 for (int i = 0; i < trackingHitList.size(); i++) {
890 std::cout << "[" << trackingHitList[i].pos.X() << ", "
891 << trackingHitList[i].pos.Y() << ", "
892 << trackingHitList[i].layer << "] ";
893 }
894 std::cout << std::endl;
895 ldmx_log(debug) << "====== END OF Tracking hit list ======";
896 }
897
898 track_list.push_back(temp_track_list);
899 // The *current* hit will have been removed, so iHit is currently pointing
900 // to the next hit. Decrement iHit so no hits will get skipped by iHit++
901 iHit--;
902 }
903 }
904
905 ldmx_log(debug) << "Straight tracks found (before merge): "
906 << track_list.size();
907 if (verbose_) {
908 for (int iTrack = 0; iTrack < track_list.size(); iTrack++) {
909 ldmx_log(debug) << "Track " << iTrack << ":";
910 for (int iHit = 0; iHit < track_list[iTrack].size(); iHit++) {
911 std::cout << " Hit " << iHit << ": ["
912 << track_list[iTrack][iHit].pos.X() << ", "
913 << track_list[iTrack][iHit].pos.Y() << ", "
914 << track_list[iTrack][iHit].layer << "]" << std::endl;
915 }
916 std::cout << std::endl;
917 }
918 }
919
920 // Optional addition: Merge nearby straight tracks. Not necessary for veto.
921 // Criteria: consider tail of track. Merge if head of next track is 1/2
922 // layers behind, within 1 cell of xy position.
923 ldmx_log(debug) << "Beginning track merging using " << track_list.size()
924 << " tracks";
925
926 for (int track_i = 0; track_i < track_list.size(); track_i++) {
927 // for each track, check the remainder of the track list for compatible
928 // tracks
929 std::vector<HitData> base_track = track_list[track_i];
930 HitData tail_hitdata = base_track.back(); // xylayer of last hit in track
931 if (verbose_) ldmx_log(debug) << " Considering track " << track_i;
932 for (int track_j = track_i + 1; track_j < track_list.size(); track_j++) {
933 std::vector<HitData> checking_track = track_list[track_j];
934 HitData head_hitdata = checking_track.front();
935 // if 1-2 layers behind, and xy within one cell...
936 if ((head_hitdata.layer == tail_hitdata.layer + 1 ||
937 head_hitdata.layer == tail_hitdata.layer + 2) &&
938 pow(pow(head_hitdata.pos.X() - tail_hitdata.pos.X(), 2) +
939 pow(head_hitdata.pos.Y() - tail_hitdata.pos.Y(), 2),
940 0.5) <= cellWidth) {
941 // ...then append the second track to the first one and delete it
942 // NOTE: TO ADD: (trackingHitList[iHit].pos -
943 // trackingHitList[jHit].pos).Mag()
944 if (verbose_) {
945 ldmx_log(debug) << " ** Compatible track found at index "
946 << track_j;
947 ldmx_log(debug) << " Tail xylayer: " << head_hitdata.pos.X()
948 << "," << head_hitdata.pos.Y() << ","
949 << head_hitdata.layer;
950 ldmx_log(debug) << " Head xylayer: " << tail_hitdata.pos.X()
951 << "," << tail_hitdata.pos.Y() << ","
952 << tail_hitdata.layer;
953 }
954 for (int hit_k = 0; hit_k < checking_track.size(); hit_k++) {
955 base_track.push_back(track_list[track_j][hit_k]);
956 }
957 track_list[track_i] = base_track;
958 track_list.erase(track_list.begin() + track_j);
959 break;
960 }
961 }
962 }
963 nStraightTracks_ = track_list.size();
964 // print the track list
965 ldmx_log(debug) << "Straight tracks found (after merge): "
967 for (int track_i = 0; track_i < track_list.size(); track_i++) {
968 ldmx_log(debug) << "Track " << track_i << ":";
969 for (int hit_i = 0; hit_i < track_list[track_i].size(); hit_i++) {
970 ldmx_log(debug) << " Hit " << hit_i << ": ["
971 << track_list[track_i][hit_i].pos.X() << ", "
972 << track_list[track_i][hit_i].pos.Y() << ", "
973 << track_list[track_i][hit_i].layer << "]";
974 }
975 }
976
977 // ------------------------------------------------------
978 // Linreg tracking:
979 ldmx_log(info) << "Finding linreg tracks from a total of "
980 << trackingHitList.size() << " hits using a radius of "
981 << linreg_radius_ << " mm";
982
983 for (int iHit = 0; iHit < trackingHitList.size(); iHit++) {
984 // Hits being considered at a given time
985 std::vector<int> hitsInRegion;
986 TMatrixD Vm(3, 3);
987 TMatrixD hdt(3, 3);
988 TVector3 slopeVec;
989 TVector3 hmean;
990 TVector3 hpoint;
991 float r_corr_best{0.0};
992 // Temp array having 3 potential hits
993 int hitNums[3];
994 // From the above which are passing the correlation reqs
995 int bestHitNums[3];
996
997 hitsInRegion.push_back(iHit);
998 // Find all hits within 2 cells of the primary hit:
999 for (int jHit = 0; jHit < trackingHitList.size(); jHit++) {
1000 // Dont try to put hits on the same layer to the lin-reg track
1001 if (trackingHitList[iHit].pos(2) == trackingHitList[jHit].pos(2)) {
1002 continue;
1003 }
1004 float dstToHit =
1005 (trackingHitList[iHit].pos - trackingHitList[jHit].pos).Mag();
1006 // This distance optimized to give the best significance
1007 // it used to be 2*cellWidth, i.e. 4.81 mm
1008 // note, the layers in the back have a separation of 22.3
1009 if (dstToHit <= 2 * linreg_radius_) {
1010 hitsInRegion.push_back(jHit);
1011 }
1012 }
1013 // Found a track that passed the lin-reg reqs
1014 bool bestLinRegFound{false};
1015 if (verbose_) {
1016 ldmx_log(debug) << "There are " << hitsInRegion.size()
1017 << " hits within a radius of " << linreg_radius_ << " mm";
1018 }
1019 // Look at combinations of hits within the region (do not consider the same
1020 // combination twice):
1021 hitNums[0] = iHit;
1022 for (int jHitInReg = 1; jHitInReg < hitsInRegion.size() - 1; jHitInReg++) {
1023 // We require (exactly) 3 hits for the lin-reg track building
1024 if (hitsInRegion.size() < 3) break;
1025 hitNums[1] = hitsInRegion[jHitInReg];
1026 for (int kHitReg = jHitInReg + 1; kHitReg < hitsInRegion.size();
1027 kHitReg++) {
1028 hitNums[2] = hitsInRegion[kHitReg];
1029 for (int hInd = 0; hInd < 3; hInd++) {
1030 // hmean = geometric mean, subtract off from hits to improve SVD
1031 // performance
1032 hmean(hInd) = (trackingHitList[hitNums[0]].pos(hInd) +
1033 trackingHitList[hitNums[1]].pos(hInd) +
1034 trackingHitList[hitNums[2]].pos(hInd)) /
1035 3.0;
1036 }
1037 for (int hInd = 0; hInd < 3; hInd++) {
1038 for (int lInd = 0; lInd < 3; lInd++) {
1039 hdt(hInd, lInd) =
1040 trackingHitList[hitNums[hInd]].pos(lInd) - hmean(lInd);
1041 }
1042 }
1043
1044 // Perform "linreg" on selected points
1045 // Calculate the determinant of the matrix
1046 double determinant =
1047 hdt(0, 0) * (hdt(1, 1) * hdt(2, 2) - hdt(1, 2) * hdt(2, 1)) -
1048 hdt(0, 1) * (hdt(1, 0) * hdt(2, 2) - hdt(1, 2) * hdt(2, 0)) +
1049 hdt(0, 2) * (hdt(1, 0) * hdt(2, 1) - hdt(1, 1) * hdt(2, 0));
1050 // Exit early if the matrix is singular (i.e. det = 0)
1051 if (determinant == 0) continue;
1052 // Perform matrix decomposition with SVD
1053 TDecompSVD svdObj(hdt);
1054 bool decomposed = svdObj.Decompose();
1055 if (!decomposed) continue;
1056
1057 // First col of V matrix is the slope of the best-fit line
1058 Vm = svdObj.GetV();
1059 for (int hInd = 0; hInd < 3; hInd++) {
1060 slopeVec(hInd) = Vm[0][hInd];
1061 }
1062 // hmean, hpoint are points on the best-fit line
1063 hpoint = slopeVec + hmean;
1064 // linreg complete: Now have best-fit line for 3 hits under
1065 // consideration Check whether the track is valid: r^2 must be high,
1066 // and the track must plausibly originate from the photon
1067 float closest_e = distTwoLines(hmean, hpoint, e_traj_start, e_traj_end);
1068 float closest_p = distTwoLines(hmean, hpoint, p_traj_start, p_traj_end);
1069 // Projected track must be close to the photon; details may change after
1070 // future study.
1071 if (closest_p > cellWidth or closest_e < 1.5 * cellWidth) continue;
1072 // find r^2
1073 // ~variance
1074 float vrnc = (trackingHitList[hitNums[0]].pos - hmean).Mag() +
1075 (trackingHitList[hitNums[1]].pos - hmean).Mag() +
1076 (trackingHitList[hitNums[2]].pos - hmean).Mag();
1077 // sum of |errors|
1078 float sumerr =
1079 distPtToLine(trackingHitList[hitNums[0]].pos, hmean, hpoint) +
1080 distPtToLine(trackingHitList[hitNums[1]].pos, hmean, hpoint) +
1081 distPtToLine(trackingHitList[hitNums[2]].pos, hmean, hpoint);
1082 float r_corr = 1 - sumerr / vrnc;
1083 // Check whether r^2 exceeds a low minimum r_corr: "Fake" tracks are
1084 // still much more common in background, so making the algorithm
1085 // oversensitive doesn't lower performance significantly
1086 if (r_corr > r_corr_best and r_corr > .6) {
1087 r_corr_best = r_corr;
1088 // Only looking for 3-hit tracks currently
1089 bestLinRegFound = true;
1090 for (int k = 0; k < 3; k++) {
1091 bestHitNums[k] = hitNums[k];
1092 }
1093 }
1094 } // end loop on hits in the region
1095 } // end 2nd loop on hits in the region
1096
1097 // Continue early if not hits on track
1098 if (!bestLinRegFound) continue;
1099 // Otherwise increase the number of lin-reg tracks
1101 ldmx_log(debug) << " Lin-reg track " << nLinregTracks_;
1102 for (int finalHitIndx = 0; finalHitIndx < 3; finalHitIndx++) {
1103 ldmx_log(debug) << " Hit " << finalHitIndx << " ["
1104 << trackingHitList[bestHitNums[finalHitIndx]].pos(0)
1105 << ", "
1106 << trackingHitList[bestHitNums[finalHitIndx]].pos(1)
1107 << ", "
1108 << trackingHitList[bestHitNums[finalHitIndx]].pos(2)
1109 << "] ";
1110 }
1111
1112 // Exclude all hits in a found track from further consideration:
1113 for (int lHit = 0; lHit < 3; lHit++) {
1114 trackingHitList.erase(trackingHitList.begin() + bestHitNums[lHit]);
1115 }
1116 iHit--;
1117 } // end loop on all hits
1118
1119 ldmx_log(info) << " MIP tracking completed; found " << nStraightTracks_
1120 << " straight tracks and " << nLinregTracks_
1121 << " lin-reg tracks";
1122
1123 result.setVariables(
1124 nReadoutHits_, deepestLayerHit_, summedDet_, summedTightIso_, maxCellDep_,
1125 showerRMS_, xStd_, yStd_, avgLayerHit_, stdLayerHit_, ecalBackEnergy_,
1127 photonTerritoryHits_, epAng_, epSep_, epDot_, electronContainmentEnergy,
1128 photonContainmentEnergy, outsideContainmentEnergy,
1129 outsideContainmentNHits, outsideContainmentXstd, outsideContainmentYstd,
1130 energySeg, xMeanSeg, yMeanSeg, xStdSeg, yStdSeg, layerMeanSeg,
1131 layerStdSeg, eContEnergy, eContXMean, eContYMean, gContEnergy, gContNHits,
1132 gContXMean, gContYMean, oContEnergy, oContNHits, oContXMean, oContYMean,
1133 oContXStd, oContYStd, oContLayerMean, oContLayerStd,
1134 ecalLayerEdepReadout_, recoilP, recoilPos);
1135
1136 buildBDTFeatureVector(result);
1137 ldmx::Ort::FloatArrays inputs({bdtFeatures_});
1138 float pred = rt_->run({featureListName_}, inputs, {"probabilities"})[0].at(1);
1139 ldmx_log(info) << " BDT was ran, score is " << pred;
1140 // Other considerations were (nLinregTracks_ == 0) && (firstNearPhLayer_ >=
1141 // 6)
1142 // && (epAng_ > 3.0 && epAng_ < 900 || epSep_ > 10.0 && epSep_ < 900)
1143 bool passesTrackingVeto = (nStraightTracks_ < 3);
1144 result.setVetoResult(pred > bdtCutVal_ && passesTrackingVeto);
1145 result.setDiscValue(pred);
1146 ldmx_log(info) << " The pred > bdtCutVal = " << (pred > bdtCutVal_)
1147 << " and MIP tracking passed = " << passesTrackingVeto;
1148
1149 // Persist in the event if the recoil ele is fiducial
1150 result.setFiducial(inside);
1151
1152 // If the event passes the veto, keep it. Otherwise,
1153 // drop the event.
1154 if (!inverse_skim_) {
1155 if (result.passesVeto()) {
1157 } else {
1159 }
1160 } else {
1161 // Invert the skimming logic
1162 if (result.passesVeto()) {
1164 } else {
1166 }
1167 }
1168
1169 event.add(collectionName_, result);
1170}
std::tuple< int, const ldmx::SimParticle * > getRecoil(const std::map< int, ldmx::SimParticle > &particleMap)
Find and return the sim particle associated with the recoil electron.
std::vector< float > trackProp(const ldmx::Tracks &tracks, ldmx::TrackStateType ts_type, const std::string &ts_title)
Return a vector of parameters for a propagated recoil track.
float distTwoLines(TVector3 v1, TVector3 v2, TVector3 w1, TVector3 w2)
Returns the distance between the lines v and w, with v defined to pass through the points (v1,...
void buildBDTFeatureVector(const ldmx::EcalVetoResult &result)
float distPtToLine(TVector3 h1, TVector3 p1, TVector3 p2)
Return the minimum distance between the point h1 and the line passing through points p1 and p2.
void fillHitMap(const std::vector< ldmx::EcalHit > &ecalRecHits, std::map< ldmx::EcalID, float > &cellMap_)
Function to load up empty vector of hit maps.
const T & getCondition(const std::string &condition_name)
Access a conditions object for the current event.
void setStorageHint(framework::StorageControl::Hint hint)
Mark the current event as having the given storage control hint from this module.
bool exists(const std::string &name, const std::string &passName="", bool unique=true) const
Check for the existence of an object or collection with the given name and pass name in the event.
Definition Event.cxx:92
EcalID getID(double x, double y, double z, bool fallible=false) const
Get a cell's ID number from its position.
double getCellMaxR() const
Get the center-to-corner radius of the cell hexagons.
void setVariables(int nReadoutHits, int deepestLayerHit, float summedDet, float summedTightIso, float maxCellDep, float showerRMS, float xStd, float yStd, float avgLayerHit, float stdLayerHit, float ecalBackEnergy, int nStraightTracks, int nLinregTracks, int firstNearPhLayer, int nNearPhHits, int photonTerritoryHits, float epAng, float epSep, float epDot, std::vector< float > electronContainmentEnergy, std::vector< float > photonContainmentEnergy, std::vector< float > outsideContainmentEnergy, std::vector< int > outsideContainmentNHits, std::vector< float > outsideContainmentXStd, std::vector< float > outsideContainmentYStd, std::vector< float > energySeg, std::vector< float > xMeanSeg, std::vector< float > yMeanSeg, std::vector< float > xStdSeg, std::vector< float > yStdSeg, std::vector< float > layerMeanSeg, std::vector< float > layerStdSeg, std::vector< std::vector< float > > eContEnergy, std::vector< std::vector< float > > eContXMean, std::vector< std::vector< float > > eContYMean, std::vector< std::vector< float > > gContEnergy, std::vector< std::vector< int > > gContNHits, std::vector< std::vector< float > > gContXMean, std::vector< std::vector< float > > gContYMean, std::vector< std::vector< float > > oContEnergy, std::vector< std::vector< int > > oContNHits, std::vector< std::vector< float > > oContXMean, std::vector< std::vector< float > > oContYMean, std::vector< std::vector< float > > oContXStd, std::vector< std::vector< float > > oContYStd, std::vector< std::vector< float > > oContLayerMean, std::vector< std::vector< float > > oContLayerStd, std::vector< float > EcalLayerEdepReadout, std::vector< double > recoilP, std::vector< float > recoilPos)
Set the sim particle and 'is findable' flag.
bool passesVeto() const
Checks if the event passes the Ecal veto.
Class representing a simulated particle.
Definition SimParticle.h:23
Implements detector ids for special simulation-derived hits like scoring planes.
Represents a simulated tracker hit in the simulation.
Implementation of a track object.
Definition Track.h:53
constexpr StorageControl::Hint hint_shouldKeep
storage control hint alias for backwards compatibility
constexpr StorageControl::Hint hint_shouldDrop
storage control hint alias for backwards compatibility

References buildBDTFeatureVector(), collectionName_, distPtToLine(), distTwoLines(), epAng_, epDot_, epSep_, framework::Event::exists(), fillHitMap(), firstNearPhLayer_, geometry_, ldmx::EcalGeometry::getCellMaxR(), framework::EventProcessor::getCondition(), ldmx::EcalGeometry::getID(), ldmx::EcalGeometry::getPosition(), ldmx::EcalGeometry::getZPosition(), framework::hint_shouldDrop, framework::hint_shouldKeep, nLinregTracks_, nNearPhHits_, nStraightTracks_, ldmx::EcalVetoResult::passesVeto(), photonTerritoryHits_, ldmx::SimSpecialID::plane(), framework::EventProcessor::setStorageHint(), ldmx::EcalVetoResult::setVariables(), and trackProp().

◆ trackProp()

std::vector< float > ecal::EcalVetoProcessor::trackProp ( const ldmx::Tracks & tracks,
ldmx::TrackStateType ts_type,
const std::string & ts_title )
private

Return a vector of parameters for a propagated recoil track.

Parameters
[in]tracksThe track collection
[in]ts_typeThe track state type, i.e. tracks state at the ECAL face
[in]ts_titleThe track state title, most likely "ecal"
Returns
Vector of parameters for a propagated recoil track

Definition at line 1306 of file EcalVetoProcessor.cxx.

1308 {
1309 // Vector to hold the new track state variables
1310 std::vector<float> new_track_states;
1311
1312 // Return if no tracks
1313 if (tracks.empty()) return new_track_states;
1314
1315 // Otherwise loop on the tracks
1316 for (auto &track : tracks) {
1317 // Get track state for ts_type
1318 auto trk_ts = track.getTrackState(ts_type);
1319 // Continue if there's no value
1320 if (!trk_ts.has_value()) continue;
1321 ldmx::Track::TrackState &ecal_track_state = trk_ts.value();
1322
1323 // Check that the track state is filled
1324 if (ecal_track_state.params.size() < 5) continue;
1325
1326 float track_state_loc0 = static_cast<float>(ecal_track_state.params[0]);
1327 float track_state_loc1 = static_cast<float>(ecal_track_state.params[1]);
1328
1329 // Store the new track state variables
1330 new_track_states.push_back(track_state_loc0);
1331 new_track_states.push_back(track_state_loc1);
1332 // z-position at the ECAL scoring plane
1333 new_track_states.push_back(239.999);
1334
1335 // Break after getting the first valid track state
1336 // TODO: interface this with CLUE to make sure the propageted track
1337 // has an associated cluster in the ECAL
1338 break;
1339 }
1340
1341 return new_track_states;
1342}

Referenced by produce().

Member Data Documentation

◆ avgLayerHit_

double ecal::EcalVetoProcessor::avgLayerHit_ {0}
private

Definition at line 139 of file EcalVetoProcessor.h.

139{0};

◆ bdtCutVal_

double ecal::EcalVetoProcessor::bdtCutVal_ {0}
private

Definition at line 162 of file EcalVetoProcessor.h.

162{0};

◆ bdtFeatures_

std::vector<float> ecal::EcalVetoProcessor::bdtFeatures_
private

Definition at line 171 of file EcalVetoProcessor.h.

◆ bdtFileName_

std::string ecal::EcalVetoProcessor::bdtFileName_
private

Definition at line 169 of file EcalVetoProcessor.h.

◆ beamEnergyMeV_

double ecal::EcalVetoProcessor::beamEnergyMeV_ {0}
private

Definition at line 164 of file EcalVetoProcessor.h.

164{0};

◆ cellMap_

std::map<ldmx::EcalID, float> ecal::EcalVetoProcessor::cellMap_
private

Definition at line 120 of file EcalVetoProcessor.h.

◆ cellMapTightIso_

std::map<ldmx::EcalID, float> ecal::EcalVetoProcessor::cellMapTightIso_
private

Definition at line 121 of file EcalVetoProcessor.h.

◆ collectionName_

std::string ecal::EcalVetoProcessor::collectionName_ {"EcalVeto"}
private

Name of the collection which will containt the results.

Definition at line 181 of file EcalVetoProcessor.h.

181{"EcalVeto"};

Referenced by configure(), and produce().

◆ deepestLayerHit_

int ecal::EcalVetoProcessor::deepestLayerHit_ {0}
private

Definition at line 131 of file EcalVetoProcessor.h.

131{0};

◆ ecalBackEnergy_

double ecal::EcalVetoProcessor::ecalBackEnergy_ {0}
private

Definition at line 141 of file EcalVetoProcessor.h.

141{0};

◆ ecalLayerEdepRaw_

std::vector<float> ecal::EcalVetoProcessor::ecalLayerEdepRaw_
private

Definition at line 123 of file EcalVetoProcessor.h.

◆ ecalLayerEdepReadout_

std::vector<float> ecal::EcalVetoProcessor::ecalLayerEdepReadout_
private

Definition at line 124 of file EcalVetoProcessor.h.

◆ ecalLayerTime_

std::vector<float> ecal::EcalVetoProcessor::ecalLayerTime_
private

Definition at line 125 of file EcalVetoProcessor.h.

◆ epAng_

float ecal::EcalVetoProcessor::epAng_ {0}
private

Angular separation between the projected photon and electron trajectories (currently unused)

Definition at line 153 of file EcalVetoProcessor.h.

153{0};

Referenced by produce().

◆ epDot_

float ecal::EcalVetoProcessor::epDot_ {0}
private

Dot product of the photon and electron momenta unit vectors.

Definition at line 158 of file EcalVetoProcessor.h.

158{0};

Referenced by produce().

◆ epSep_

float ecal::EcalVetoProcessor::epSep_ {0}
private

Distance between the projected photon and electron trajectories at the ECal face.

Definition at line 156 of file EcalVetoProcessor.h.

156{0};

Referenced by produce().

◆ featureListName_

std::string ecal::EcalVetoProcessor::featureListName_
private

Definition at line 172 of file EcalVetoProcessor.h.

◆ firstNearPhLayer_

int ecal::EcalVetoProcessor::firstNearPhLayer_ {0}
private

First ECal layer in which a hit is found near the photon.

Definition at line 148 of file EcalVetoProcessor.h.

148{0};

Referenced by produce().

◆ geometry_

const ldmx::EcalGeometry* ecal::EcalVetoProcessor::geometry_
private

handle to current geometry (to share with member functions)

Definition at line 186 of file EcalVetoProcessor.h.

Referenced by produce().

◆ inverse_skim_

bool ecal::EcalVetoProcessor::inverse_skim_ {false}
private

Definition at line 178 of file EcalVetoProcessor.h.

178{false};

◆ linreg_radius_

double ecal::EcalVetoProcessor::linreg_radius_ {0}
private

Definition at line 165 of file EcalVetoProcessor.h.

165{0};

◆ maxCellDep_

double ecal::EcalVetoProcessor::maxCellDep_ {0}
private

Definition at line 135 of file EcalVetoProcessor.h.

135{0};

◆ nEcalLayers_

int ecal::EcalVetoProcessor::nEcalLayers_ {0}
private

Definition at line 129 of file EcalVetoProcessor.h.

129{0};

◆ nLinregTracks_

int ecal::EcalVetoProcessor::nLinregTracks_ {0}
private

Number of "linreg" tracks found in the event.

Definition at line 146 of file EcalVetoProcessor.h.

146{0};

Referenced by produce().

◆ nNearPhHits_

int ecal::EcalVetoProcessor::nNearPhHits_ {0}
private

Number of hits near the photon trajectory.

Definition at line 150 of file EcalVetoProcessor.h.

150{0};

Referenced by produce().

◆ nReadoutHits_

int ecal::EcalVetoProcessor::nReadoutHits_ {0}
private

Definition at line 130 of file EcalVetoProcessor.h.

130{0};

◆ nStraightTracks_

int ecal::EcalVetoProcessor::nStraightTracks_ {0}
private

Number of "straight" tracks found in the event.

Definition at line 144 of file EcalVetoProcessor.h.

144{0};

Referenced by produce().

◆ photonTerritoryHits_

int ecal::EcalVetoProcessor::photonTerritoryHits_ {0}
private

Number of hits in the photon territory.

Definition at line 160 of file EcalVetoProcessor.h.

160{0};

Referenced by produce().

◆ rec_coll_name_

std::string ecal::EcalVetoProcessor::rec_coll_name_
private

Definition at line 175 of file EcalVetoProcessor.h.

◆ rec_pass_name_

std::string ecal::EcalVetoProcessor::rec_pass_name_
private

Definition at line 174 of file EcalVetoProcessor.h.

◆ recoil_from_tracking_

bool ecal::EcalVetoProcessor::recoil_from_tracking_
private

Definition at line 176 of file EcalVetoProcessor.h.

◆ roc_range_values_

std::vector<std::vector<double> > ecal::EcalVetoProcessor::roc_range_values_
private

Definition at line 127 of file EcalVetoProcessor.h.

◆ rocFileName_

std::string ecal::EcalVetoProcessor::rocFileName_
private

Definition at line 170 of file EcalVetoProcessor.h.

◆ rt_

std::unique_ptr<ldmx::Ort::ONNXRuntime> ecal::EcalVetoProcessor::rt_
private

Definition at line 183 of file EcalVetoProcessor.h.

◆ showerRMS_

double ecal::EcalVetoProcessor::showerRMS_ {0}
private

Definition at line 136 of file EcalVetoProcessor.h.

136{0};

◆ stdLayerHit_

double ecal::EcalVetoProcessor::stdLayerHit_ {0}
private

Definition at line 140 of file EcalVetoProcessor.h.

140{0};

◆ summedDet_

double ecal::EcalVetoProcessor::summedDet_ {0}
private

Definition at line 133 of file EcalVetoProcessor.h.

133{0};

◆ summedTightIso_

double ecal::EcalVetoProcessor::summedTightIso_ {0}
private

Definition at line 134 of file EcalVetoProcessor.h.

134{0};

◆ track_collection_

std::string ecal::EcalVetoProcessor::track_collection_
private

Definition at line 177 of file EcalVetoProcessor.h.

◆ verbose_

bool ecal::EcalVetoProcessor::verbose_ {false}
private

Definition at line 167 of file EcalVetoProcessor.h.

167{false};

◆ xStd_

double ecal::EcalVetoProcessor::xStd_ {0}
private

Definition at line 137 of file EcalVetoProcessor.h.

137{0};

◆ yStd_

double ecal::EcalVetoProcessor::yStd_ {0}
private

Definition at line 138 of file EcalVetoProcessor.h.

138{0};

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