LDMX Software
ecal::EcalVetoProcessor Class Reference

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

#include <EcalVetoProcessor.h>

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 onNewRun (const ldmx::RunHeader &rh) override
 onNewRun is the first function called for each processor after the conditions are fully configured and accessible.
 
void onProcessEnd () override
 Callback for the EventProcessor to take any necessary action when the processing of events finishes, such as calculating job-summary quantities.
 
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 process (Event &event) final
 Processing an event for a Producer is calling produce.
 
- Public Member Functions inherited from framework::EventProcessor
 DECLARE_FACTORY (EventProcessor, EventProcessor *, const std::string &, Process &)
 declare that we have a factory for this class
 
 EventProcessor (const std::string &name, Process &process)
 Class constructor.
 
virtual ~EventProcessor ()=default
 Class destructor.
 
virtual void beforeNewRun (ldmx::RunHeader &run_header)
 Callback for Producers to add parameters to the run header before conditions are initialized.
 
virtual void onFileOpen (EventFile &event_file)
 Callback for the EventProcessor to take any necessary action when a new event input ROOT file is opened.
 
virtual void onFileClose (EventFile &event_file)
 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.
 
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 > &ecal_rec_hits, float &shower_rms)
 
void fillHitMap (const std::vector< ldmx::EcalHit > &ecal_rec_hits, std::map< ldmx::EcalID, float > &cell_map_)
 Function to load up empty vector of hit maps.
 
void fillIsolatedHitMap (const std::vector< ldmx::EcalHit > &ecal_rec_hits, ldmx::EcalID global_centroid, std::map< ldmx::EcalID, float > &cell_map, std::map< ldmx::EcalID, float > &cell_map_iso, bool doTight=false)
 
std::vector< XYCoords > getTrajectory (std::array< float, 3 > momentum, std::array< float, 3 > position)
 
void buildBDTFeatureVector (const ldmx::EcalVetoResult &result)
 

Private Attributes

int nevents_ {0}
 
float processing_time_ {0.}
 
std::map< std::string, float > profiling_map_
 
std::map< ldmx::EcalID, float > cell_map_
 
std::map< ldmx::EcalID, float > cell_map_tight_iso_
 
std::vector< float > ecal_layer_edep_raw_
 
std::vector< float > ecal_layer_edep_readout_
 
std::vector< float > ecal_layer_time_
 
std::vector< std::vector< float > > roc_range_values_
 
int n_ecal_layers_ {0}
 
int n_readout_hits_ {0}
 
int deepest_layer_hit_ {0}
 
float summed_det_ {0}
 
float summed_tight_iso_ {0}
 
float max_cell_dep_ {0}
 
float shower_rms_ {0}
 
float x_std_ {0}
 
float y_std_ {0}
 
float avg_layer_hit_ {0}
 
float std_layer_hit_ {0}
 
float ecal_back_energy_ {0}
 
int n_tracking_hits_ {0}
 Number of hits outside of the electron roc in the Ecal or if the electron trajectory is missing, all the hits in the Ecal.
 
float ep_ang_ {0}
 Angular separation between the projected photon and electron trajectories as projected at ECAL.
 
float ep_ang_at_target_ {0}
 Angular separation between the projected photon and electron trajectories as at Target.
 
float ep_sep_ {0}
 Distance between the projected photon and electron trajectories at the ECal face.
 
float ep_dot_ {0}
 Dot product of the photon and electron momenta unit vectors at Ecal.
 
float ep_dot_at_target_ {0}
 Dot product of the photon and electron momenta unit vectors at Target.
 
float bdt_cut_val_ {0}
 
float beam_energy_mev_ {0}
 
std::string bdt_file_name_
 
std::string roc_file_name_
 
std::vector< float > bdt_features_
 
std::string feature_list_name_
 
std::string ecal_sp_coll_name_
 
std::string target_sp_coll_name_
 
std::string sp_pass_name_
 
std::string rec_pass_name_
 
std::string rec_coll_name_
 
bool recoil_from_tracking_
 
std::string track_pass_name_
 
std::string track_collection_
 
std::string sim_particles_coll_name_
 
std::string sim_particles_passname_
 
bool inverse_skim_ {false}
 
std::string collection_name_ {"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

- Protected Member Functions inherited from framework::EventProcessor
void abortEvent ()
 Abort the event immediately.
 
- Protected Attributes inherited from framework::EventProcessor
HistogramPool histograms_
 helper object for making and filling histograms
 
NtupleManagerntuple_ {NtupleManager::getInstance()}
 Manager for any ntuples.
 
logging::logger the_log_
 The logger for this EventProcessor.
 

Detailed Description

Determines if event is vetoable using ECAL hit information.

Definition at line 45 of file EcalVetoProcessor.h.

Member Typedef Documentation

◆ CellEnergyPair

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

Definition at line 47 of file EcalVetoProcessor.h.

◆ XYCoords

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

Definition at line 49 of file EcalVetoProcessor.h.

Constructor & Destructor Documentation

◆ EcalVetoProcessor()

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

Definition at line 51 of file EcalVetoProcessor.h.

52 : Producer(name, process) {}
Producer(const std::string &name, Process &process)
Class constructor.
virtual void process(Event &event) final
Processing an event for a Producer is calling produce.

◆ ~EcalVetoProcessor()

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

Definition at line 54 of file EcalVetoProcessor.h.

54{}

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 17 of file EcalVetoProcessor.cxx.

18 {
19 // Base variables
20 bdt_features_.push_back(result.getNReadoutHits());
21 bdt_features_.push_back(result.getSummedDet());
22 bdt_features_.push_back(result.getSummedTightIso());
23 bdt_features_.push_back(result.getMaxCellDep());
24 bdt_features_.push_back(result.getShowerRMS());
25 bdt_features_.push_back(result.getXStd());
26 bdt_features_.push_back(result.getYStd());
27 bdt_features_.push_back(result.getAvgLayerHit());
28 bdt_features_.push_back(result.getStdLayerHit());
29 bdt_features_.push_back(result.getDeepestLayerHit());
30 bdt_features_.push_back(result.getEcalBackEnergy());
31 // MIP tracking
32
33 bdt_features_.push_back(-1.); // NStraight
34 bdt_features_.push_back(-1.); // FirstNearPHLayer
35 bdt_features_.push_back(-1.); // NNearPHHits
36 bdt_features_.push_back(-1.); // PhotonTerritoryHits
37
38 // bdt_features_.push_back(result.getNStraightTracks());
39 // bdt_features_.push_back(result.getFirstNearPhLayer());
40 // bdt_features_.push_back(result.getNNearPhHits());
41 // bdt_features_.push_back(result.getPhotonTerritoryHits());
42 // bdt_features_.push_back(result.getNTrackingHits());
43 bdt_features_.push_back(result.getEPSep());
44 bdt_features_.push_back(result.getEPDot());
45 // Longitudinal segment variables
46 bdt_features_.push_back(result.getEnergySeg()[0]);
47 bdt_features_.push_back(result.getXMeanSeg()[0]);
48 bdt_features_.push_back(result.getYMeanSeg()[0]);
49 bdt_features_.push_back(result.getLayerMeanSeg()[0]);
50 bdt_features_.push_back(result.getEnergySeg()[1]);
51 bdt_features_.push_back(result.getYMeanSeg()[2]);
53 bdt_features_.push_back(result.getEleContEnergy()[0][0]);
54 bdt_features_.push_back(result.getEleContEnergy()[1][0]);
55 bdt_features_.push_back(result.getEleContYMean()[0][0]);
56 bdt_features_.push_back(result.getEleContEnergy()[0][1]);
57 bdt_features_.push_back(result.getEleContEnergy()[1][1]);
58 bdt_features_.push_back(result.getEleContYMean()[0][1]);
60 bdt_features_.push_back(result.getPhContNHits()[0][0]);
61 bdt_features_.push_back(result.getPhContYMean()[0][0]);
62 bdt_features_.push_back(result.getPhContNHits()[0][1]);
64 bdt_features_.push_back(result.getOutContEnergy()[0][0]);
65 bdt_features_.push_back(result.getOutContEnergy()[1][0]);
66 bdt_features_.push_back(result.getOutContEnergy()[2][0]);
67 bdt_features_.push_back(result.getOutContNHits()[0][0]);
68 bdt_features_.push_back(result.getOutContXMean()[0][0]);
69 bdt_features_.push_back(result.getOutContYMean()[0][0]);
70 bdt_features_.push_back(result.getOutContYMean()[1][0]);
71 bdt_features_.push_back(result.getOutContYStd()[0][0]);
72 bdt_features_.push_back(result.getOutContEnergy()[0][1]);
73 bdt_features_.push_back(result.getOutContEnergy()[1][1]);
74 bdt_features_.push_back(result.getOutContEnergy()[2][1]);
75 bdt_features_.push_back(result.getOutContLayerMean()[0][1]);
76 bdt_features_.push_back(result.getOutContLayerStd()[0][1]);
77 bdt_features_.push_back(result.getOutContEnergy()[0][2]);
78 bdt_features_.push_back(result.getOutContLayerMean()[0][2]);
79}

Referenced by produce().

◆ clearProcessor()

void ecal::EcalVetoProcessor::clearProcessor ( )
private

Definition at line 140 of file EcalVetoProcessor.cxx.

140 {
141 cell_map_.clear();
142 cell_map_tight_iso_.clear();
143 bdt_features_.clear();
144
145 n_readout_hits_ = 0;
146 summed_det_ = 0;
147 summed_tight_iso_ = 0;
148 max_cell_dep_ = 0;
149 shower_rms_ = 0;
150 x_std_ = 0;
151 y_std_ = 0;
152 avg_layer_hit_ = 0;
153 std_layer_hit_ = 0;
154 deepest_layer_hit_ = 0;
155 ecal_back_energy_ = 0;
157 ep_ang_ = 0;
159 ep_sep_ = 0;
160 ep_dot_ = 0;
162
163 std::fill(ecal_layer_edep_raw_.begin(), ecal_layer_edep_raw_.end(), 0);
164 std::fill(ecal_layer_edep_readout_.begin(), ecal_layer_edep_readout_.end(),
165 0);
166 std::fill(ecal_layer_time_.begin(), ecal_layer_time_.end(), 0);
167}
float ep_ang_
Angular separation between the projected photon and electron trajectories as projected at ECAL.
float ep_ang_at_target_
Angular separation between the projected photon and electron trajectories as at Target.
float ep_dot_
Dot product of the photon and electron momenta unit vectors at Ecal.
float ep_sep_
Distance between the projected photon and electron trajectories at the ECal face.
int n_tracking_hits_
Number of hits outside of the electron roc in the Ecal or if the electron trajectory is missing,...
float ep_dot_at_target_
Dot product of the photon and electron momenta unit vectors at Target.

◆ 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 81 of file EcalVetoProcessor.cxx.

81 {
82 feature_list_name_ = parameters.get<std::string>("feature_list_name");
83
84 sim_particles_passname_ =
85 parameters.get<std::string>("sim_particles_passname");
86 // Load BDT ONNX file
87 rt_ = std::make_unique<ldmx::ort::ONNXRuntime>(
88 parameters.get<std::string>("bdt_file"));
89
90 // Read in arrays holding 68% containment radius_ per layer_
91 // for different bins in momentum/angle
92 roc_file_name_ = parameters.get<std::string>("roc_file");
93 if (!std::ifstream(roc_file_name_).good()) {
94 EXCEPTION_RAISE(
95 "EcalVetoProcessor",
96 "The specified RoC file '" + roc_file_name_ + "' does not exist!");
97 } else {
98 std::ifstream rocfile(roc_file_name_);
99 std::string line, value;
100
101 // Extract the first line in the file
102 std::getline(rocfile, line);
103 std::vector<float> values;
104
105 // Read data, line by line
106 while (std::getline(rocfile, line)) {
107 std::stringstream ss(line);
108 values.clear();
109 while (std::getline(ss, value, ',')) {
110 float f_value = (value != "") ? std::stof(value) : -1.0;
111 values.push_back(f_value);
112 }
113 roc_range_values_.push_back(values);
114 }
115 }
116
117 n_ecal_layers_ = parameters.get<int>("num_ecal_layers");
118
119 bdt_cut_val_ = parameters.get<double>("disc_cut");
120 ecal_layer_edep_raw_.resize(n_ecal_layers_, 0);
121 ecal_layer_edep_readout_.resize(n_ecal_layers_, 0);
122 ecal_layer_time_.resize(n_ecal_layers_, 0);
123
124 beam_energy_mev_ = parameters.get<double>("beam_energy");
125 // Set the collection name as defined in the configuration
126 ecal_sp_coll_name_ = parameters.get<std::string>("ecal_sp_coll_name");
127 target_sp_coll_name_ = parameters.get<std::string>("target_sp_coll_name");
128 sp_pass_name_ = parameters.get<std::string>("sp_pass_name");
129 sim_particles_coll_name_ =
130 parameters.get<std::string>("sim_particles_coll_name");
131 collection_name_ = parameters.get<std::string>("collection_name");
132 rec_pass_name_ = parameters.get<std::string>("rec_pass_name");
133 rec_coll_name_ = parameters.get<std::string>("rec_coll_name");
134 recoil_from_tracking_ = parameters.get<bool>("recoil_from_tracking");
135 track_collection_ = parameters.get<std::string>("track_collection");
136 track_pass_name_ = parameters.get<std::string>("track_pass_name", "");
137 inverse_skim_ = parameters.get<bool>("inverse_skim");
138}
std::string collection_name_
Name of the collection which will containt the results.
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78

References collection_name_, and framework::config::Parameters::get().

◆ fillHitMap()

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

Function to load up empty vector of hit maps.

Definition at line 1029 of file EcalVetoProcessor.cxx.

1031 {
1032 for (const ldmx::EcalHit &hit : ecal_rec_hits) {
1033 ldmx::EcalID id(hit.getID());
1034 cellMap.emplace(id, hit.getEnergy());
1035 }
1036}
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 > & ecal_rec_hits,
ldmx::EcalID global_centroid,
std::map< ldmx::EcalID, float > & cell_map,
std::map< ldmx::EcalID, float > & cell_map_iso,
bool doTight = false )
private

Definition at line 1038 of file EcalVetoProcessor.cxx.

1041 {
1042 for (const ldmx::EcalHit &hit : ecal_rec_hits) {
1043 auto isolated_hit = std::make_pair(true, ldmx::EcalID());
1044 ldmx::EcalID id(hit.getID());
1045 if (do_tight) {
1046 // Disregard hits_ that are on the centroid.
1047 if (id == global_centroid) continue;
1048
1049 // Skip hits_ that are on centroid inner ring
1050 if (geometry_->isNN(global_centroid, id)) {
1051 continue;
1052 }
1053 }
1054
1055 // Skip hits_ that have a readout neighbor
1056 // Get neighboring cell id's and try to look them up in the full cell map
1057 // (constant speed algo.)
1058 // these ideas are only cell/module_ (must ignore layer_)
1059 std::vector<ldmx::EcalID> cell_nbr_ids = geometry_->getNN(id);
1060
1061 for (int k = 0; k < cell_nbr_ids.size(); k++) {
1062 // update neighbor ID to the current layer_
1063 cell_nbr_ids[k] = ldmx::EcalID(id.layer(), cell_nbr_ids[k].module(),
1064 cell_nbr_ids[k].cell());
1065 // look in cell hit map to see if it is there
1066 if (cellMap.find(cell_nbr_ids[k]) != cellMap.end()) {
1067 isolated_hit = std::make_pair(false, cell_nbr_ids[k]);
1068 break;
1069 }
1070 }
1071 if (!isolated_hit.first) {
1072 continue;
1073 }
1074 // Insert isolated hit
1075 cellMapIso.emplace(id, hit.getEnergy());
1076 }
1077}
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 > & ecal_rec_hits,
float & shower_rms )
private

Definition at line 979 of file EcalVetoProcessor.cxx.

980 {
981 auto wgt_centroid_coords = std::make_pair<float, float>(0., 0.);
982 float sum_edep = 0;
983 ldmx::EcalID return_cell_id;
984
985 // Calculate Energy Weighted Centroid
986 for (const ldmx::EcalHit &hit : ecal_rec_hits) {
987 ldmx::EcalID id(hit.getID());
988 CellEnergyPair cell_energy_pair = std::make_pair(id, hit.getEnergy());
989 auto [rechit_x, rechit_y, rechit_z] = geometry_->getPosition(id);
990 XYCoords centroid_coords = std::make_pair(rechit_x, rechit_y);
991 wgt_centroid_coords.first = wgt_centroid_coords.first +
992 centroid_coords.first * cell_energy_pair.second;
993 wgt_centroid_coords.second =
994 wgt_centroid_coords.second +
995 centroid_coords.second * cell_energy_pair.second;
996 sum_edep += cell_energy_pair.second;
997 } // end loop over rec hits
998 wgt_centroid_coords.first = (sum_edep > 1E-6)
999 ? wgt_centroid_coords.first / sum_edep
1000 : wgt_centroid_coords.first;
1001 wgt_centroid_coords.second = (sum_edep > 1E-6)
1002 ? wgt_centroid_coords.second / sum_edep
1003 : wgt_centroid_coords.second;
1004 // Find Nearest Cell to Centroid
1005 float max_dist = 1e6;
1006 for (const ldmx::EcalHit &hit : ecal_rec_hits) {
1007 auto [rechit_x, rechit_y, rechit_z] = geometry_->getPosition(hit.getID());
1008 XYCoords centroid_coords = std::make_pair(rechit_x, rechit_y);
1009
1010 float delta_r =
1011 sqrt(((centroid_coords.first - wgt_centroid_coords.first) *
1012 (centroid_coords.first - wgt_centroid_coords.first)) +
1013 ((centroid_coords.second - wgt_centroid_coords.second) *
1014 (centroid_coords.second - wgt_centroid_coords.second)));
1015 shower_rms += delta_r * hit.getEnergy();
1016 if (delta_r < max_dist) {
1017 max_dist = delta_r;
1018 return_cell_id = ldmx::EcalID(hit.getID());
1019 }
1020 }
1021 if (sum_edep > 0) shower_rms = shower_rms / sum_edep;
1022 // flatten
1023 return ldmx::EcalID(0, return_cell_id.module(), return_cell_id.cell());
1024}
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::array< float, 3 > momentum,
std::array< float, 3 > position )
private

Definition at line 1081 of file EcalVetoProcessor.cxx.

1082 {
1083 std::vector<XYCoords> positions;
1084 for (int i_layer = 0; i_layer < n_ecal_layers_; i_layer++) {
1085 float pos_x =
1086 position[0] + (momentum[0] / momentum[2]) *
1087 (geometry_->getZPosition(i_layer) - position[2]);
1088 float pos_y =
1089 position[1] + (momentum[1] / momentum[2]) *
1090 (geometry_->getZPosition(i_layer) - position[2]);
1091 positions.push_back(std::make_pair(pos_x, pos_y));
1092 }
1093 return positions;
1094}
double getZPosition(int layer) const
Get the z-coordinate given the layer id.

◆ onNewRun()

void ecal::EcalVetoProcessor::onNewRun ( const ldmx::RunHeader & rh)
overridevirtual

onNewRun is the first function called for each processor after the conditions are fully configured and accessible.

This is where you could create single-processors, multi-event calculation objects.

Reimplemented from framework::EventProcessor.

Definition at line 5 of file EcalVetoProcessor.cxx.

5 {
6 profiling_map_["setup"] = 0.;
7 profiling_map_["recoil_electron"] = 0.;
8 profiling_map_["trajectories"] = 0.;
9 profiling_map_["roc_var"] = 0.;
10 profiling_map_["fill_hitmaps"] = 0.;
11 profiling_map_["containment_var"] = 0.;
12 profiling_map_["mip_tracking_setup"] = 0.;
13 profiling_map_["set_variables"] = 0.;
14 profiling_map_["bdt_variables"] = 0.;
15}

◆ onProcessEnd()

void ecal::EcalVetoProcessor::onProcessEnd ( )
overridevirtual

Callback for the EventProcessor to take any necessary action when the processing of events finishes, such as calculating job-summary quantities.

Reimplemented from framework::EventProcessor.

Definition at line 939 of file EcalVetoProcessor.cxx.

939 {
940 ldmx_log(info) << "AVG Time/Event: " << std::fixed << std::setprecision(2)
941 << processing_time_ / nevents_ << " ms";
942
943 ldmx_log(info) << "Breakdown::";
944
945 ldmx_log(info) << "setup Avg Time/Event = " << std::fixed
946 << std::setprecision(3) << profiling_map_["setup"] / nevents_
947 << " ms";
948
949 ldmx_log(info) << "recoil_electron Avg Time/Event = " << std::fixed
950 << std::setprecision(3)
951 << profiling_map_["recoil_electron"] / nevents_ << " ms";
952 ldmx_log(info) << "trajectories Avg Time/Event = " << std::fixed
953 << std::setprecision(3)
954 << profiling_map_["trajectories"] / nevents_ << " ms";
955
956 ldmx_log(info) << "roc_var Avg Time/Event = " << std::fixed
957 << std::setprecision(3) << profiling_map_["roc_var"] / nevents_
958 << " ms";
959
960 ldmx_log(info) << "fill_hitmaps Avg Time/Event = " << std::fixed
961 << std::setprecision(3)
962 << profiling_map_["fill_hitmaps"] / nevents_ << " ms";
963 ldmx_log(info) << "containment_var Avg Time/Event = " << std::fixed
964 << std::setprecision(3)
965 << profiling_map_["containment_var"] / nevents_ << " ms";
966 ldmx_log(info) << "mip_tracking_setup Avg Time/Event = " << std::fixed
967 << std::setprecision(3)
968 << profiling_map_["mip_tracking_setup"] / nevents_ << " ms";
969
970 ldmx_log(info) << "set_variables Avg Time/Event = " << std::fixed
971 << std::setprecision(3)
972 << profiling_map_["set_variables"] / nevents_ << " ms";
973
974 ldmx_log(info) << "bdt_variables Avg Time/Event = " << std::fixed
975 << std::setprecision(3)
976 << profiling_map_["bdt_variables"] / nevents_ << " ms";
977}

◆ 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 169 of file EcalVetoProcessor.cxx.

169 {
170 auto start = std::chrono::high_resolution_clock::now();
171 nevents_++;
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::array<float, 3> recoil_p = {0., 0., 0.};
185 std::array<float, 3> recoil_pos = {-9999., -9999., -9999.};
186 std::array<float, 3> recoil_p_at_target = {0., 0., 0.};
187 std::array<float, 3> recoil_pos_at_target = {-9999., -9999., -9999.};
188
189 auto setup = std::chrono::high_resolution_clock::now();
190 profiling_map_["setup"] +=
191 std::chrono::duration<float, std::milli>(setup - start).count();
192
193 if (!recoil_from_tracking_ &&
194 event.exists(ecal_sp_coll_name_, sp_pass_name_)) {
195 ldmx_log(trace) << " Loop through all of the sim particles and find the "
196 "recoil electron";
197 //
198 // Loop through all of the sim particles and find the recoil electron.
199 //
200
201 // Get the collection of simulated particles from the event
202 auto particle_map{event.getMap<int, ldmx::SimParticle>(
203 sim_particles_coll_name_, sim_particles_passname_)};
204
205 // Search for the recoil electron
206 auto [recoil_track_id, recoil_electron] = analysis::getRecoil(particle_map);
207
208 // Find ECAL SP hit for recoil electron
209 auto ecal_sp_hits{event.getCollection<ldmx::SimTrackerHit>(
210 ecal_sp_coll_name_, sp_pass_name_)};
211 float pmax = 0;
212 for (ldmx::SimTrackerHit &sp_hit : ecal_sp_hits) {
213 ldmx::SimSpecialID hit_id(sp_hit.getID());
214 auto ecal_sp_momentum = sp_hit.getMomentum();
215 auto ecal_sp_position = sp_hit.getPosition();
216 if (hit_id.plane() != 31 || ecal_sp_momentum[2] <= 0) continue;
217
218 if (sp_hit.getTrackID() == recoil_track_id) {
219 // A*A is faster than pow(A,2)
220 if (sqrt((ecal_sp_momentum[0] * ecal_sp_momentum[0]) +
221 (ecal_sp_momentum[1] * ecal_sp_momentum[1]) +
222 (ecal_sp_momentum[2] * ecal_sp_momentum[2])) > pmax) {
223 recoil_p = {static_cast<float>(ecal_sp_momentum[0]),
224 static_cast<float>(ecal_sp_momentum[1]),
225 static_cast<float>(ecal_sp_momentum[2])};
226 recoil_pos = {(ecal_sp_position[0]), (ecal_sp_position[1]),
227 (ecal_sp_position[2])};
228 pmax = sqrt(recoil_p[0] * recoil_p[0] + recoil_p[1] * recoil_p[1] +
229 recoil_p[2] * recoil_p[2]);
230 }
231 }
232 }
233
234 // Find target SP hit for recoil electron
235 if (event.exists(target_sp_coll_name_, sp_pass_name_)) {
236 std::vector<ldmx::SimTrackerHit> target_sp_hits =
237 event.getCollection<ldmx::SimTrackerHit>(target_sp_coll_name_,
238 sp_pass_name_);
239 pmax = 0;
240 for (ldmx::SimTrackerHit &sp_hit : target_sp_hits) {
241 ldmx::SimSpecialID hit_id(sp_hit.getID());
242 auto target_sp_momentum = sp_hit.getMomentum();
243 auto target_sp_position = sp_hit.getPosition();
244 if (hit_id.plane() != 1 || target_sp_momentum[2] <= 0) continue;
245
246 if (sp_hit.getTrackID() == recoil_track_id) {
247 if (sqrt((target_sp_momentum[0] * target_sp_momentum[0]) +
248 (target_sp_momentum[1] * target_sp_momentum[1]) +
249 (target_sp_momentum[2] * target_sp_momentum[2])) > pmax) {
250 recoil_p_at_target = {static_cast<float>(target_sp_momentum[0]),
251 static_cast<float>(target_sp_momentum[1]),
252 static_cast<float>(target_sp_momentum[2])};
253 recoil_pos_at_target = {target_sp_position[0],
254 target_sp_position[1],
255 target_sp_position[2]};
256 // (A*A) is faster than pow(A,2)
257 pmax = sqrt((recoil_p_at_target[0] * recoil_p_at_target[0]) +
258 (recoil_p_at_target[1] * recoil_p_at_target[1]) +
259 (recoil_p_at_target[2] * recoil_p_at_target[2]));
260 }
261 }
262 } // end loop on target SP hits_
263 } // end condition on target SP
264 } // end condition on ecal SP
265
266 // Get recoil_pos using recoil tracking
267 bool fiducial_in_tracker{false};
268 if (recoil_from_tracking_) {
269 ldmx_log(trace) << " Get recoil tracks collection";
270
271 // Get the recoil track collection
272 auto recoil_tracks{
273 event.getCollection<ldmx::Track>(track_collection_, track_pass_name_)};
274
275 ldmx_log(trace) << " Propagate the recoil ele to the ECAL";
276 ldmx::TrackStateType ts_type = ldmx::AtECAL;
277 auto recoil_track_states_ecal =
278 ecal::trackProp(recoil_tracks, ts_type, "ecal");
279 ldmx_log(trace) << " Propagate the recoil ele to the Target";
280 ldmx::TrackStateType ts_type_target = ldmx::AtTarget;
281 auto recoil_track_states_target =
282 ecal::trackProp(recoil_tracks, ts_type_target, "target");
283
284 ldmx_log(trace) << " Set recoil_pos and recoil_p";
285 // Redefining recoil_pos now to come from the track state
286 // track_state_loc0 is recoil_pos[0] and track_state_loc1 is recoil_pos[1]
287 if (!recoil_track_states_ecal.empty()) {
288 fiducial_in_tracker = true;
289 recoil_pos = {recoil_track_states_ecal[0], recoil_track_states_ecal[1],
290 recoil_track_states_ecal[2]};
291 recoil_p = {(recoil_track_states_ecal[3]), (recoil_track_states_ecal[4]),
292 (recoil_track_states_ecal[5])};
293 } else {
294 ldmx_log(trace) << " No recoil track at ECAL";
295 fiducial_in_tracker = false;
296 }
297 ldmx_log(debug) << " Set recoil_p = (" << recoil_p[0] << ", "
298 << recoil_p[1] << ", " << recoil_p[2]
299 << ") and recoil_pos = (" << recoil_pos[0] << ", "
300 << recoil_pos[1] << ", " << recoil_pos[2] << ")";
301
302 // Repeat the above but now for the taget states
303 if (!recoil_track_states_target.empty()) {
304 recoil_pos_at_target = {(recoil_track_states_target[0]),
305 (recoil_track_states_target[1]),
306 (recoil_track_states_target[2])};
307 recoil_p_at_target = {recoil_track_states_target[3],
308 recoil_track_states_target[4],
309 recoil_track_states_target[5]};
310 }
311 ldmx_log(debug) << " Set recoil_p_at_target = (" << recoil_p_at_target[0]
312 << ", " << recoil_p_at_target[1] << ", "
313 << recoil_p_at_target[2] << ") and recoil_pos_at_target = ("
314 << recoil_pos_at_target[0] << ", "
315 << recoil_pos_at_target[1] << ", "
316 << recoil_pos_at_target[2] << ")";
317 } // condition to do recoil information from tracking
318
319 ldmx_log(trace) << " Get projected trajectories for electron and photon";
320
321 auto recoil_electron = std::chrono::high_resolution_clock::now();
322 profiling_map_["recoil_electron"] +=
323 std::chrono::duration<float, std::milli>(recoil_electron - setup).count();
324
325 // Get projected trajectories for electron and photon
326 std::vector<XYCoords> ele_trajectory, photon_trajectory,
327 ele_trajectory_at_target;
328 // Require that z_-momentum is positive (which will also exclude the default
329 // initializaton) Require that the positions are not the default initializaton
330 if ((recoil_p[2] > 0.) && (recoil_p_at_target[2] > 0.) &&
331 (recoil_pos[0] != -9999.) && (recoil_pos_at_target[0] != -9999.)) {
332 ele_trajectory = getTrajectory(recoil_p, recoil_pos);
333 ele_trajectory_at_target =
334 getTrajectory(recoil_p_at_target, recoil_pos_at_target);
335 // Get the photon projection. This does not require that the photon exists
336 // tho
337 std::array<float, 3> photon_proj_momentum = {
338 -recoil_p_at_target[0], -recoil_p_at_target[1],
339 beam_energy_mev_ - recoil_p_at_target[2]};
340 photon_trajectory =
341 getTrajectory(photon_proj_momentum, recoil_pos_at_target);
342 } else {
343 ldmx_log(trace) << "Ele / photon trajectory cannot be determined, pZ = "
344 << recoil_p[2] << " pZAtTarget = " << recoil_p_at_target[2]
345 << " X = " << recoil_pos[0]
346 << " XAtTarget = " << recoil_pos_at_target[0];
347 }
348
349 float recoil_p_mag = (recoil_p[2] > 0.) ? sqrt((recoil_p[0] * recoil_p[0]) +
350 (recoil_p[1] * recoil_p[1]) +
351 (recoil_p[2] * recoil_p[2]))
352 : -1.0;
353 float recoil_theta =
354 recoil_p_mag > 0 ? acos(recoil_p[2] / recoil_p_mag) * 180.0 / M_PI : -1.0;
355
356 ldmx_log(trace) << " Build Radii of containment (ROC)";
357
358 auto trajectories = std::chrono::high_resolution_clock::now();
359 profiling_map_["trajectories"] +=
360 std::chrono::duration<float, std::milli>(trajectories - recoil_electron)
361 .count();
362
363 // Use the appropriate containment radii for the recoil electron
364 std::vector<float> roc_values_bin_0(roc_range_values_[0].begin() + 4,
365 roc_range_values_[0].end());
366 std::vector<float> ele_radii = roc_values_bin_0;
367 float theta_min, theta_max, p_min, p_max;
368 bool inrange;
369
370 for (int i = 0; i < roc_range_values_.size(); i++) {
371 theta_min = roc_range_values_[i][0];
372 theta_max = roc_range_values_[i][1];
373 p_min = roc_range_values_[i][2];
374 p_max = roc_range_values_[i][3];
375 inrange = true;
376
377 if (theta_min != -1.0) {
378 inrange = inrange && (recoil_theta >= theta_min);
379 }
380 if (theta_max != -1.0) {
381 inrange = inrange && (recoil_theta < theta_max);
382 }
383 if (p_min != -1.0) {
384 inrange = inrange && (recoil_p_mag >= p_min);
385 }
386 if (p_max != -1.0) {
387 inrange = inrange && (recoil_p_mag < p_max);
388 }
389 if (inrange) {
390 std::vector<float> roc_values_bini(roc_range_values_[i].begin() + 4,
391 roc_range_values_[i].end());
392 ele_radii = roc_values_bini;
393 }
394 }
395 // Use default RoC bin for photon
396 std::vector<float> photon_radii = roc_values_bin_0;
397
398 auto roc_var = std::chrono::high_resolution_clock::now();
399 profiling_map_["roc_var"] +=
400 std::chrono::duration<float, std::milli>(roc_var - trajectories).count();
401
402 // Get the collection of digitized Ecal hits_ from the event.
403 const std::vector<ldmx::EcalHit> ecal_rec_hits =
404 event.getCollection<ldmx::EcalHit>(rec_coll_name_, rec_pass_name_);
405
406 ldmx::EcalID global_centroid =
407 getShowerCentroidIdAndRms(ecal_rec_hits, shower_rms_);
408 /* ~~ Fill the hit map ~~ O(n) */
409 fillHitMap(ecal_rec_hits, cell_map_);
410 bool do_tight = true;
411 /* ~~ Fill the isolated hit maps ~~ O(n) */
412 fillIsolatedHitMap(ecal_rec_hits, global_centroid, cell_map_,
413 cell_map_tight_iso_, do_tight);
414
415 auto fill_hitmaps = std::chrono::high_resolution_clock::now();
416 profiling_map_["fill_hitmaps"] +=
417 std::chrono::duration<float, std::milli>(fill_hitmaps - roc_var).count();
418
419 // Loop over the hits_ from the event to calculate the rest of the important
420 // quantities
421
422 float w_avg_layer_hit = 0;
423 float x_mean = 0;
424 float y_mean = 0;
425
426 // Containment variables
427 unsigned int n_regions = 5;
428 std::vector<float> electron_containment_energy(n_regions, 0.0);
429 std::vector<float> photon_containment_energy(n_regions, 0.0);
430 std::vector<float> outside_containment_energy(n_regions, 0.0);
431 std::vector<int> outside_containment_n_hits(n_regions, 0);
432 std::vector<float> outside_containment_x_mean(n_regions, 0.0);
433 std::vector<float> outside_containment_y_mean(n_regions, 0.0);
434 std::vector<float> outside_containment_x_std(n_regions, 0.0);
435 std::vector<float> outside_containment_y_std(n_regions, 0.0);
436 // Longitudinal segmentation
437 std::vector<int> seg_layers = {0, 6, 17, 32};
438 unsigned int n_segments = seg_layers.size() - 1;
439 std::vector<float> energy_seg(n_segments, 0.0);
440 std::vector<float> x_mean_seg(n_segments, 0.0);
441 std::vector<float> x_std_seg(n_segments, 0.0);
442 std::vector<float> y_mean_seg(n_segments, 0.0);
443 std::vector<float> y_std_seg(n_segments, 0.0);
444 std::vector<float> layer_mean_seg(n_segments, 0.0);
445 std::vector<float> layer_std_seg(n_segments, 0.0);
446 std::vector<std::vector<float>> e_cont_energy(
447 n_regions, std::vector<float>(n_segments, 0.0));
448 std::vector<std::vector<float>> e_cont_x_mean(
449 n_regions, std::vector<float>(n_segments, 0.0));
450 std::vector<std::vector<float>> e_cont_y_mean(
451 n_regions, std::vector<float>(n_segments, 0.0));
452 std::vector<std::vector<float>> g_cont_energy(
453 n_regions, std::vector<float>(n_segments, 0.0));
454 std::vector<std::vector<int>> g_cont_n_hits(n_regions,
455 std::vector<int>(n_segments, 0));
456 std::vector<std::vector<float>> g_cont_x_mean(
457 n_regions, std::vector<float>(n_segments, 0.0));
458 std::vector<std::vector<float>> g_cont_y_mean(
459 n_regions, std::vector<float>(n_segments, 0.0));
460 std::vector<std::vector<float>> o_cont_energy(
461 n_regions, std::vector<float>(n_segments, 0.0));
462 std::vector<std::vector<int>> o_cont_n_hits(n_regions,
463 std::vector<int>(n_segments, 0));
464 std::vector<std::vector<float>> o_cont_x_mean(
465 n_regions, std::vector<float>(n_segments, 0.0));
466 std::vector<std::vector<float>> o_cont_y_mean(
467 n_regions, std::vector<float>(n_segments, 0.0));
468 std::vector<std::vector<float>> o_cont_x_std(
469 n_regions, std::vector<float>(n_segments, 0.0));
470 std::vector<std::vector<float>> o_cont_y_std(
471 n_regions, std::vector<float>(n_segments, 0.0));
472 std::vector<std::vector<float>> o_cont_layer_mean(
473 n_regions, std::vector<float>(n_segments, 0.0));
474 std::vector<std::vector<float>> o_cont_layer_std(
475 n_regions, std::vector<float>(n_segments, 0.0));
476
477 auto containment_var = std::chrono::high_resolution_clock::now();
478 profiling_map_["containment_var"] +=
479 std::chrono::duration<float, std::milli>(containment_var - fill_hitmaps)
480 .count();
481
482 // MIP tracking: vector of hits_ to be used in the MIP tracking algorithm.
483 // All hits_ inside the electron ROC (or all hits_ in the ECal if the event is
484 // missing an electron) will be included.
485 std::vector<ldmx::HitData> tracking_hit_list;
486
487 ldmx_log(trace)
488 << " Loop over the hits_ from the event to calculate the BDT features";
489
490 for (const ldmx::EcalHit &hit : ecal_rec_hits) {
491 // Layer-wise quantities
492 ldmx::EcalID id(hit.getID());
493 ecal_layer_edep_raw_[id.layer()] =
494 ecal_layer_edep_raw_[id.layer()] + hit.getEnergy();
495 if (id.layer() >= 20) ecal_back_energy_ += hit.getEnergy();
496 if (max_cell_dep_ < hit.getEnergy()) max_cell_dep_ = hit.getEnergy();
497 if (hit.getEnergy() <= 0) {
498 ldmx_log(fatal)
499 << "ECal hit has negative or zero energy, this should never happen, "
500 "check the threshold settings in HgcrocEmulator";
501 continue;
502 }
503 n_readout_hits_++;
504 ecal_layer_edep_readout_[id.layer()] += hit.getEnergy();
505 ecal_layer_time_[id.layer()] += (hit.getEnergy()) * hit.getTime();
506 auto [rechit_x, rechit_y, rechit_z] = geometry_->getPosition(id);
507 x_mean += rechit_x * hit.getEnergy();
508 y_mean += rechit_y * hit.getEnergy();
509 avg_layer_hit_ += id.layer();
510 w_avg_layer_hit += id.layer() * hit.getEnergy();
511 if (deepest_layer_hit_ < id.layer()) {
512 deepest_layer_hit_ = id.layer();
513 }
514 XYCoords xy_pair = std::make_pair(rechit_x, rechit_y);
515 float distance_ele_trajectory =
516 ele_trajectory.size()
517 ? sqrt((xy_pair.first - ele_trajectory[id.layer()].first) *
518 (xy_pair.first - ele_trajectory[id.layer()].first) +
519 (xy_pair.second - ele_trajectory[id.layer()].second) *
520 (xy_pair.second - ele_trajectory[id.layer()].second))
521 : -1.0;
522 float distance_photon_trajectory =
523 photon_trajectory.size()
524 ? sqrt((xy_pair.first - photon_trajectory[id.layer()].first) *
525 (xy_pair.first - photon_trajectory[id.layer()].first) +
526 (xy_pair.second - photon_trajectory[id.layer()].second) *
527 (xy_pair.second - photon_trajectory[id.layer()].second))
528 : -1.0;
529
530 // Decide which longitudinal segment the hit is in and add to sums
531 for (unsigned int iseg = 0; iseg < n_segments; iseg++) {
532 if (id.layer() >= seg_layers[iseg] &&
533 id.layer() <= seg_layers[iseg + 1] - 1) {
534 energy_seg[iseg] += hit.getEnergy();
535 x_mean_seg[iseg] += xy_pair.first * hit.getEnergy();
536 y_mean_seg[iseg] += xy_pair.second * hit.getEnergy();
537 layer_mean_seg[iseg] += id.layer() * hit.getEnergy();
538
539 // Decide which containment region the hit is in and add to sums
540 for (unsigned int ireg = 0; ireg < n_regions; ireg++) {
541 if (distance_ele_trajectory >= ireg * ele_radii[id.layer()] &&
542 distance_ele_trajectory < (ireg + 1) * ele_radii[id.layer()]) {
543 e_cont_energy[ireg][iseg] += hit.getEnergy();
544 e_cont_x_mean[ireg][iseg] += xy_pair.first * hit.getEnergy();
545 e_cont_y_mean[ireg][iseg] += xy_pair.second * hit.getEnergy();
546 }
547 if (distance_photon_trajectory >= ireg * photon_radii[id.layer()] &&
548 distance_photon_trajectory <
549 (ireg + 1) * photon_radii[id.layer()]) {
550 g_cont_energy[ireg][iseg] += hit.getEnergy();
551 g_cont_n_hits[ireg][iseg] += 1;
552 g_cont_x_mean[ireg][iseg] += xy_pair.first * hit.getEnergy();
553 g_cont_y_mean[ireg][iseg] += xy_pair.second * hit.getEnergy();
554 }
555 if (distance_ele_trajectory > (ireg + 1) * ele_radii[id.layer()] &&
556 distance_photon_trajectory >
557 (ireg + 1) * photon_radii[id.layer()]) {
558 o_cont_energy[ireg][iseg] += hit.getEnergy();
559 o_cont_n_hits[ireg][iseg] += 1;
560 o_cont_x_mean[ireg][iseg] += xy_pair.first * hit.getEnergy();
561 o_cont_y_mean[ireg][iseg] += xy_pair.second * hit.getEnergy();
562 o_cont_layer_mean[ireg][iseg] += id.layer() * hit.getEnergy();
563 }
564 }
565 }
566 }
567
568 // Decide which containment region the hit is in and add to sums
569 for (unsigned int ireg = 0; ireg < n_regions; ireg++) {
570 if (distance_ele_trajectory >= ireg * ele_radii[id.layer()] &&
571 distance_ele_trajectory < (ireg + 1) * ele_radii[id.layer()])
572 electron_containment_energy[ireg] += hit.getEnergy();
573 if (distance_photon_trajectory >= ireg * photon_radii[id.layer()] &&
574 distance_photon_trajectory < (ireg + 1) * photon_radii[id.layer()])
575 photon_containment_energy[ireg] += hit.getEnergy();
576 if (distance_ele_trajectory > (ireg + 1) * ele_radii[id.layer()] &&
577 distance_photon_trajectory > (ireg + 1) * photon_radii[id.layer()]) {
578 outside_containment_energy[ireg] += hit.getEnergy();
579 outside_containment_n_hits[ireg] += 1;
580 outside_containment_x_mean[ireg] += xy_pair.first * hit.getEnergy();
581 outside_containment_y_mean[ireg] += xy_pair.second * hit.getEnergy();
582 }
583 }
584
585 // MIP tracking: Decide whether hit should be added to tracking_hit_list
586 // If inside e- RoC or if etraj is missing, use the hit for tracking:
587 if (distance_ele_trajectory >= ele_radii[id.layer()] ||
588 distance_ele_trajectory == -1.0) {
589 ldmx::HitData hd;
590 hd.pos_ = ROOT::Math::XYZVector(xy_pair.first, xy_pair.second,
591 geometry_->getZPosition(id.layer()));
592 hd.layer_ = id.layer();
593 tracking_hit_list.push_back(hd);
594 }
595 } // end loop over rechits
596
597 n_tracking_hits_ = tracking_hit_list.size();
598
599 for (const auto &[id, energy] : cell_map_tight_iso_) {
600 if (energy > 0) summed_tight_iso_ += energy;
601 }
602
603 for (int i_layer = 0; i_layer < ecal_layer_edep_readout_.size(); i_layer++) {
604 ecal_layer_time_[i_layer] =
605 ecal_layer_time_[i_layer] / ecal_layer_edep_readout_[i_layer];
606 summed_det_ += ecal_layer_edep_readout_[i_layer];
607 }
608
609 if (n_readout_hits_ > 0) {
610 avg_layer_hit_ /= n_readout_hits_;
611 w_avg_layer_hit /= summed_det_;
612 x_mean /= summed_det_;
613 y_mean /= summed_det_;
614 } else {
615 w_avg_layer_hit = 0;
616 avg_layer_hit_ = 0;
617 x_mean = 0;
618 y_mean = 0;
619 }
620
621 // If necessary, quotient out the total energy from the means
622 for (unsigned int iseg = 0; iseg < n_segments; iseg++) {
623 if (energy_seg[iseg] > 0) {
624 x_mean_seg[iseg] /= energy_seg[iseg];
625 y_mean_seg[iseg] /= energy_seg[iseg];
626 layer_mean_seg[iseg] /= energy_seg[iseg];
627 }
628 for (unsigned int ireg = 0; ireg < n_regions; ireg++) {
629 if (e_cont_energy[ireg][iseg] > 0) {
630 e_cont_x_mean[ireg][iseg] /= e_cont_energy[ireg][iseg];
631 e_cont_y_mean[ireg][iseg] /= e_cont_energy[ireg][iseg];
632 }
633 if (g_cont_energy[ireg][iseg] > 0) {
634 g_cont_x_mean[ireg][iseg] /= g_cont_energy[ireg][iseg];
635 g_cont_y_mean[ireg][iseg] /= g_cont_energy[ireg][iseg];
636 }
637 if (o_cont_energy[ireg][iseg] > 0) {
638 o_cont_x_mean[ireg][iseg] /= o_cont_energy[ireg][iseg];
639 o_cont_y_mean[ireg][iseg] /= o_cont_energy[ireg][iseg];
640 o_cont_layer_mean[ireg][iseg] /= o_cont_energy[ireg][iseg];
641 }
642 }
643 }
644
645 for (unsigned int ireg = 0; ireg < n_regions; ireg++) {
646 if (outside_containment_energy[ireg] > 0) {
647 outside_containment_x_mean[ireg] /= outside_containment_energy[ireg];
648 outside_containment_y_mean[ireg] /= outside_containment_energy[ireg];
649 }
650 }
651
652 // Loop over hits_ a second time to find the standard deviations.
653 for (const ldmx::EcalHit &hit : ecal_rec_hits) {
654 ldmx::EcalID id(hit.getID());
655 auto [rechit_x, rechit_y, rechit_z] = geometry_->getPosition(id);
656 if (hit.getEnergy() > 0) {
657 x_std_ += pow((rechit_x - x_mean), 2) * hit.getEnergy();
658 y_std_ += pow((rechit_y - y_mean), 2) * hit.getEnergy();
659 std_layer_hit_ +=
660 pow((id.layer() - w_avg_layer_hit), 2) * hit.getEnergy();
661 }
662 XYCoords xy_pair = std::make_pair(rechit_x, rechit_y);
663 float distance_ele_trajectory =
664 ele_trajectory.size()
665 ? sqrt(pow((xy_pair.first - ele_trajectory[id.layer()].first), 2) +
666 pow((xy_pair.second - ele_trajectory[id.layer()].second), 2))
667 : -1.0;
668 float distance_photon_trajectory =
669 photon_trajectory.size()
670 ? sqrt(pow((xy_pair.first - photon_trajectory[id.layer()].first),
671 2) +
672 pow((xy_pair.second - photon_trajectory[id.layer()].second),
673 2))
674 : -1.0;
675
676 for (unsigned int iseg = 0; iseg < n_segments; iseg++) {
677 if (id.layer() >= seg_layers[iseg] &&
678 id.layer() <= seg_layers[iseg + 1] - 1) {
679 x_std_seg[iseg] +=
680 pow((xy_pair.first - x_mean_seg[iseg]), 2) * hit.getEnergy();
681 y_std_seg[iseg] +=
682 pow((xy_pair.second - y_mean_seg[iseg]), 2) * hit.getEnergy();
683 layer_std_seg[iseg] +=
684 pow((id.layer() - layer_mean_seg[iseg]), 2) * hit.getEnergy();
685
686 for (unsigned int ireg = 0; ireg < n_regions; ireg++) {
687 if (distance_ele_trajectory > (ireg + 1) * ele_radii[id.layer()] &&
688 distance_photon_trajectory >
689 (ireg + 1) * photon_radii[id.layer()]) {
690 o_cont_x_std[ireg][iseg] +=
691 pow((xy_pair.first - o_cont_x_mean[ireg][iseg]), 2) *
692 hit.getEnergy();
693 o_cont_y_std[ireg][iseg] +=
694 pow((xy_pair.second - o_cont_y_mean[ireg][iseg]), 2) *
695 hit.getEnergy();
696 o_cont_layer_std[ireg][iseg] +=
697 pow((id.layer() - o_cont_layer_mean[ireg][iseg]), 2) *
698 hit.getEnergy();
699 }
700 }
701 }
702 }
703
704 for (unsigned int ireg = 0; ireg < n_regions; ireg++) {
705 if (distance_ele_trajectory > (ireg + 1) * ele_radii[id.layer()] &&
706 distance_photon_trajectory > (ireg + 1) * photon_radii[id.layer()]) {
707 outside_containment_x_std[ireg] +=
708 pow((xy_pair.first - outside_containment_x_mean[ireg]), 2) *
709 hit.getEnergy();
710 outside_containment_y_std[ireg] +=
711 pow((xy_pair.second - outside_containment_y_mean[ireg]), 2) *
712 hit.getEnergy();
713 }
714 }
715 } // end loop over rechits (2nd time)
716
717 if (n_readout_hits_ > 0) {
718 x_std_ = sqrt(x_std_ / summed_det_);
719 y_std_ = sqrt(y_std_ / summed_det_);
720 std_layer_hit_ = sqrt(std_layer_hit_ / summed_det_);
721 } else {
722 x_std_ = 0;
723 y_std_ = 0;
724 std_layer_hit_ = 0;
725 }
726
727 // Quotient out the total energies from the standard deviations if possible
728 // and take root
729 for (unsigned int iseg = 0; iseg < n_segments; iseg++) {
730 if (energy_seg[iseg] > 0) {
731 x_std_seg[iseg] = sqrt(x_std_seg[iseg] / energy_seg[iseg]);
732 y_std_seg[iseg] = sqrt(y_std_seg[iseg] / energy_seg[iseg]);
733 layer_std_seg[iseg] = sqrt(layer_std_seg[iseg] / energy_seg[iseg]);
734 }
735 for (unsigned int ireg = 0; ireg < n_regions; ireg++) {
736 if (o_cont_energy[ireg][iseg] > 0) {
737 o_cont_x_std[ireg][iseg] =
738 sqrt(o_cont_x_std[ireg][iseg] / o_cont_energy[ireg][iseg]);
739 o_cont_y_std[ireg][iseg] =
740 sqrt(o_cont_y_std[ireg][iseg] / o_cont_energy[ireg][iseg]);
741 o_cont_layer_std[ireg][iseg] =
742 sqrt(o_cont_layer_std[ireg][iseg] / o_cont_energy[ireg][iseg]);
743 }
744 }
745 }
746
747 for (unsigned int ireg = 0; ireg < n_regions; ireg++) {
748 if (outside_containment_energy[ireg] > 0) {
749 outside_containment_x_std[ireg] = sqrt(outside_containment_x_std[ireg] /
750 outside_containment_energy[ireg]);
751 outside_containment_y_std[ireg] = sqrt(outside_containment_y_std[ireg] /
752 outside_containment_energy[ireg]);
753 }
754 }
755
756 ldmx_log(trace) << " Find out if the recoil electron is fiducial";
757
758 // Find the location of the recoil electron
759 // Ecal face is not where the first layer_ starts,
760 // defined in DetDescr/python/EcalGeometry.py
761 const float dz_from_face =
763 float drifted_recoil_x{-9999.};
764 float drifted_recoil_y{-9999.};
765 if (recoil_p[2] > 0.) {
766 ldmx_log(trace) << " Recoil electron pX = " << recoil_p[0]
767 << " pY = " << recoil_p[1] << " pZ = " << recoil_p[2];
768 ldmx_log(trace) << " Recoil electron PosX = " << recoil_pos[0]
769 << " PosY = " << recoil_pos[1]
770 << " PosZ = " << recoil_pos[2];
771 drifted_recoil_x =
772 (dz_from_face * (recoil_p[0] / recoil_p[2])) + recoil_pos[0];
773 drifted_recoil_y =
774 (dz_from_face * (recoil_p[1] / recoil_p[2])) + recoil_pos[1];
775 }
776 const int recoil_layer_index = 0;
777
778 // Check if it's fiducial
779 bool inside_ecal_cell{false};
780 // At module_ level
781 const auto ecal_id = geometry_->getID(drifted_recoil_x, drifted_recoil_y,
782 recoil_layer_index, true);
783 if (!ecal_id.null()) {
784 // If fiducial at module_ level, check at cell level
785 const auto cell_id =
786 geometry_->getID(drifted_recoil_x, drifted_recoil_y, recoil_layer_index,
787 ecal_id.getModuleID(), true);
788 if (!cell_id.null()) {
789 inside_ecal_cell = true;
790 }
791 }
792
793 ldmx_log(info) << " Is this event is fiducial in ECAL? "
794 << inside_ecal_cell;
795 ROOT::Math::XYZVector e_traj_start;
796 ROOT::Math::XYZVector e_traj_end;
797 ROOT::Math::XYZVector e_traj_target_start;
798 ROOT::Math::XYZVector e_traj_target_end;
799 ROOT::Math::XYZVector p_traj_start;
800 ROOT::Math::XYZVector p_traj_end;
801 if (!ele_trajectory.empty() && !photon_trajectory.empty()) {
802 // Create TVector3s marking the start and endpoints of each projected
803 // trajectory
804 e_traj_start.SetXYZ(ele_trajectory[0].first, ele_trajectory[0].second,
806 e_traj_end.SetXYZ(ele_trajectory[(n_ecal_layers_ - 1)].first,
807 ele_trajectory[(n_ecal_layers_ - 1)].second,
808 geometry_->getZPosition((n_ecal_layers_ - 1)));
809 p_traj_start.SetXYZ(photon_trajectory[0].first, photon_trajectory[0].second,
811 p_traj_end.SetXYZ(photon_trajectory[(n_ecal_layers_ - 1)].first,
812 photon_trajectory[(n_ecal_layers_ - 1)].second,
813 geometry_->getZPosition((n_ecal_layers_ - 1)));
814
815 ROOT::Math::XYZVector evec = e_traj_end - e_traj_start;
816 ROOT::Math::XYZVector e_norm = evec.Unit();
817 ROOT::Math::XYZVector pvec = p_traj_end - p_traj_start;
818 ROOT::Math::XYZVector p_norm = pvec.Unit();
819
820 // Calculate the angle between the projected electron at Ecal and the photon
821 // (at target)
822 ep_dot_ = e_norm.Dot(p_norm);
823 ep_ang_ = acos(ep_dot_) * 180.0 / M_PI;
824 ep_sep_ = sqrt(pow(e_traj_start.X() - p_traj_start.X(), 2) +
825 pow(e_traj_start.Y() - p_traj_start.Y(), 2));
826 // Calculate the electron trajectory with positions and momentum as measured
827 // at the target
828 if (!ele_trajectory_at_target.empty()) {
829 e_traj_target_start.SetXYZ(recoil_pos_at_target[0],
830 recoil_pos_at_target[1],
831 static_cast<float>(0.0));
832 e_traj_target_end.SetXYZ(ele_trajectory_at_target[(0)].first,
833 ele_trajectory_at_target[(0)].second,
835 // Now calculate the ep angle at the target
836 ROOT::Math::XYZVector evec_target =
837 e_traj_target_end - e_traj_target_start;
838 ROOT::Math::XYZVector e_norm_target = evec_target.Unit();
839 ep_dot_at_target_ = e_norm_target.Dot(p_norm);
840 ep_ang_at_target_ = acos(ep_dot_at_target_) * 180.0 / M_PI;
841 }
842 ldmx_log(trace) << " Electron trajectory calculated";
843 } else {
844 // Electron trajectory is missing, so all hits_ in the Ecal are fair game.
845 // Pick e/ptraj so that they won't restrict the tracking algorithm (place
846 // them far outside the ECal).
847 ldmx_log(trace) << " Electron trajectory is missing";
848 e_traj_start = ROOT::Math::XYZVector(999, 999, geometry_->getZPosition(0));
849 e_traj_end = ROOT::Math::XYZVector(
850 999, 999, geometry_->getZPosition((n_ecal_layers_ - 1)));
851 p_traj_start =
852 ROOT::Math::XYZVector(1000, 1000, geometry_->getZPosition(0));
853 p_traj_end = ROOT::Math::XYZVector(
854 1000, 1000, geometry_->getZPosition((n_ecal_layers_ - 1)));
855 /*ensures event will not be vetoed by angle/separation cut */
856 ep_ang_ = 999.;
857 ep_ang_at_target_ = 999.;
858 ep_sep_ = 999.;
859 ep_dot_ = 999.;
860 ep_dot_at_target_ = 999.;
861 }
862 // Took out MIP tracking here (starting at near photon hits_)
863 ldmx::EcalTrajectoryInfo ecal_mip_collection;
864 ldmx_log(trace) << " Set up input info for MIP tracking";
865 ecal_mip_collection.setEleTrajectory(ele_trajectory);
866 ecal_mip_collection.setPhotonTrajectory(photon_trajectory);
867 ecal_mip_collection.setTrackingHitList(tracking_hit_list);
868 event.add("EcalTrajectoryInfo", ecal_mip_collection);
869
870 auto mip_tracking_setup = std::chrono::high_resolution_clock::now();
871 profiling_map_["mip_tracking_setup"] +=
872 std::chrono::duration<double, std::milli>(mip_tracking_setup - start)
873 .count();
874 result.setVariables(
875 n_readout_hits_, deepest_layer_hit_, n_tracking_hits_, summed_det_,
876 summed_tight_iso_, max_cell_dep_, shower_rms_, x_std_, y_std_,
877 avg_layer_hit_, std_layer_hit_, ecal_back_energy_, ep_ang_,
879 electron_containment_energy, photon_containment_energy,
880 outside_containment_energy, outside_containment_n_hits,
881 outside_containment_x_std, outside_containment_y_std, energy_seg,
882 x_mean_seg, y_mean_seg, x_std_seg, y_std_seg, layer_mean_seg,
883 layer_std_seg, e_cont_energy, e_cont_x_mean, e_cont_y_mean, g_cont_energy,
884 g_cont_n_hits, g_cont_x_mean, g_cont_y_mean, o_cont_energy, o_cont_n_hits,
885 o_cont_x_mean, o_cont_y_mean, o_cont_x_std, o_cont_y_std,
886 o_cont_layer_mean, o_cont_layer_std, ecal_layer_edep_readout_, recoil_p,
887 recoil_pos);
888
889 auto set_variables = std::chrono::high_resolution_clock::now();
890 profiling_map_["set_variables"] += std::chrono::duration<double, std::milli>(
891 set_variables - mip_tracking_setup)
892 .count();
893 buildBDTFeatureVector(result);
894 ldmx::ort::FloatArrays inputs({bdt_features_});
895 float pred =
896 rt_->run({feature_list_name_}, inputs, {"probabilities"})[0].at(1);
897 ldmx_log(info) << " BDT was ran, score is " << pred;
898 // Other considerations were (n_linreg_tracks_ == 0) && (first_near_ph_layer_
899 // >= 6)
900 // && (ep_ang_ > 3.0 && ep_ang_ < 900 || ep_sep_ > 10.0 && ep_sep_ < 900)
901 result.setVetoResult(pred > bdt_cut_val_);
902 result.setDiscValue(pred);
903 ldmx_log(info) << " The pred > bdtCutVal = " << (pred > bdt_cut_val_);
904
905 // Persist in the event if the recoil ele is fiducial
906 result.setFiducial(inside_ecal_cell);
907 result.setTrackingFiducial(fiducial_in_tracker);
908
909 // If the event passes the veto, keep it. Otherwise,
910 // drop the event.
911 if (!inverse_skim_) {
912 if (result.passesVeto()) {
914 } else {
916 }
917 } else {
918 // Invert the skimming logic
919 if (result.passesVeto()) {
921 } else {
923 }
924 }
925
926 event.add(collection_name_, result);
927
928 auto bdt_variables = std::chrono::high_resolution_clock::now();
929 profiling_map_["bdt_variables"] +=
930 std::chrono::duration<float, std::milli>(bdt_variables - set_variables)
931 .count();
932
933 auto end = std::chrono::high_resolution_clock::now();
934 auto time_diff = end - start;
935 processing_time_ +=
936 std::chrono::duration<float, std::milli>(time_diff).count();
937}
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.
Definition EcalHelper.cxx:5
void fillHitMap(const std::vector< ldmx::EcalHit > &ecal_rec_hits, std::map< ldmx::EcalID, float > &cell_map_)
Function to load up empty vector of hit maps.
void buildBDTFeatureVector(const ldmx::EcalVetoResult &result)
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:105
EcalID getID(double x, double y, double z, bool fallible=false) const
Get a cell's ID number from its position.
double getEcalFrontZ() const
Get the z-coordinate of the Ecal face.
bool passesVeto() const
Checks if the event passes the Ecal veto.
void setVariables(int n_readout_hits, int deepest_layer_hit, int n_tracking_hits, float summed_det, float summed_tight_iso, float max_cell_dep, float shower_rms, float x_std, float y_std, float avg_layer_hit, float std_layer_hit, float ecal_back_energy, float ep_ang, float ep_ang_at_target, float ep_sep, float ep_dot, float ep_dot_at_target, std::vector< float > electron_containment_energy, std::vector< float > photon_containment_energy, std::vector< float > outside_containment_energy, std::vector< int > outside_containment_n_hits, std::vector< float > outside_containment_x_std, std::vector< float > outside_containment_y_std, std::vector< float > energy_seg, std::vector< float > x_mean_seg, std::vector< float > y_mean_seg, std::vector< float > x_std_seg, std::vector< float > y_std_seg, std::vector< float > layer_mean_seg, std::vector< float > layer_std_seg, std::vector< std::vector< float > > e_cont_energy, std::vector< std::vector< float > > e_cont_x_mean, std::vector< std::vector< float > > e_cont_y_mean, std::vector< std::vector< float > > g_cont_energy, std::vector< std::vector< int > > g_cont_n_hits, std::vector< std::vector< float > > g_cont_x_mean, std::vector< std::vector< float > > g_cont_y_mean, std::vector< std::vector< float > > o_cont_energy, std::vector< std::vector< int > > o_cont_n_hits, std::vector< std::vector< float > > o_cont_x_mean, std::vector< std::vector< float > > o_cont_y_mean, std::vector< std::vector< float > > o_cont_x_std, std::vector< std::vector< float > > o_cont_y_std, std::vector< std::vector< float > > o_cont_layer_mean, std::vector< std::vector< float > > o_cont_layer_std, std::vector< float > ecal_layer_edep_readout, std::array< float, 3 > recoil_p, std::array< float, 3 > recoil_pos)
Set the sim particle and 'is findable' flag.
Class representing a simulated particle.
Definition SimParticle.h:24
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_SHOULD_DROP
storage control hint alias for backwards compatibility
constexpr StorageControl::Hint HINT_SHOULD_KEEP
storage control hint alias for backwards compatibility

References buildBDTFeatureVector(), collection_name_, ep_ang_, ep_ang_at_target_, ep_dot_, ep_dot_at_target_, ep_sep_, framework::Event::exists(), fillHitMap(), geometry_, framework::EventProcessor::getCondition(), ldmx::EcalGeometry::getEcalFrontZ(), ldmx::EcalGeometry::getID(), ldmx::EcalGeometry::getPosition(), ldmx::EcalGeometry::getZPosition(), framework::HINT_SHOULD_DROP, framework::HINT_SHOULD_KEEP, n_tracking_hits_, ldmx::EcalVetoResult::passesVeto(), ldmx::SimSpecialID::plane(), framework::EventProcessor::setStorageHint(), ldmx::EcalVetoResult::setVariables(), and ecal::trackProp().

Member Data Documentation

◆ avg_layer_hit_

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

Definition at line 125 of file EcalVetoProcessor.h.

125{0};

◆ bdt_cut_val_

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

Definition at line 146 of file EcalVetoProcessor.h.

146{0};

◆ bdt_features_

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

Definition at line 152 of file EcalVetoProcessor.h.

◆ bdt_file_name_

std::string ecal::EcalVetoProcessor::bdt_file_name_
private

Definition at line 150 of file EcalVetoProcessor.h.

◆ beam_energy_mev_

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

Definition at line 148 of file EcalVetoProcessor.h.

148{0};

◆ cell_map_

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

Definition at line 106 of file EcalVetoProcessor.h.

◆ cell_map_tight_iso_

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

Definition at line 107 of file EcalVetoProcessor.h.

◆ collection_name_

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

Name of the collection which will containt the results.

Definition at line 170 of file EcalVetoProcessor.h.

170{"EcalVeto"};

Referenced by configure(), and produce().

◆ deepest_layer_hit_

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

Definition at line 117 of file EcalVetoProcessor.h.

117{0};

◆ ecal_back_energy_

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

Definition at line 127 of file EcalVetoProcessor.h.

127{0};

◆ ecal_layer_edep_raw_

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

Definition at line 109 of file EcalVetoProcessor.h.

◆ ecal_layer_edep_readout_

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

Definition at line 110 of file EcalVetoProcessor.h.

◆ ecal_layer_time_

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

Definition at line 111 of file EcalVetoProcessor.h.

◆ ecal_sp_coll_name_

std::string ecal::EcalVetoProcessor::ecal_sp_coll_name_
private

Definition at line 156 of file EcalVetoProcessor.h.

◆ ep_ang_

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

Angular separation between the projected photon and electron trajectories as projected at ECAL.

Definition at line 134 of file EcalVetoProcessor.h.

134{0};

Referenced by produce().

◆ ep_ang_at_target_

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

Angular separation between the projected photon and electron trajectories as at Target.

Definition at line 137 of file EcalVetoProcessor.h.

137{0};

Referenced by produce().

◆ ep_dot_

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

Dot product of the photon and electron momenta unit vectors at Ecal.

Definition at line 142 of file EcalVetoProcessor.h.

142{0};

Referenced by produce().

◆ ep_dot_at_target_

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

Dot product of the photon and electron momenta unit vectors at Target.

Definition at line 144 of file EcalVetoProcessor.h.

144{0};

Referenced by produce().

◆ ep_sep_

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

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

Definition at line 140 of file EcalVetoProcessor.h.

140{0};

Referenced by produce().

◆ feature_list_name_

std::string ecal::EcalVetoProcessor::feature_list_name_
private

Definition at line 153 of file EcalVetoProcessor.h.

◆ geometry_

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

handle to current geometry (to share with member functions)

Definition at line 175 of file EcalVetoProcessor.h.

Referenced by produce().

◆ inverse_skim_

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

Definition at line 167 of file EcalVetoProcessor.h.

167{false};

◆ max_cell_dep_

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

Definition at line 121 of file EcalVetoProcessor.h.

121{0};

◆ n_ecal_layers_

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

Definition at line 115 of file EcalVetoProcessor.h.

115{0};

◆ n_readout_hits_

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

Definition at line 116 of file EcalVetoProcessor.h.

116{0};

◆ n_tracking_hits_

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

Number of hits outside of the electron roc in the Ecal or if the electron trajectory is missing, all the hits in the Ecal.

Definition at line 131 of file EcalVetoProcessor.h.

131{0};

Referenced by produce().

◆ nevents_

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

Definition at line 102 of file EcalVetoProcessor.h.

102{0};

◆ processing_time_

float ecal::EcalVetoProcessor::processing_time_ {0.}
private

Definition at line 103 of file EcalVetoProcessor.h.

103{0.};

◆ profiling_map_

std::map<std::string, float> ecal::EcalVetoProcessor::profiling_map_
private

Definition at line 105 of file EcalVetoProcessor.h.

◆ rec_coll_name_

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

Definition at line 160 of file EcalVetoProcessor.h.

◆ rec_pass_name_

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

Definition at line 159 of file EcalVetoProcessor.h.

◆ recoil_from_tracking_

bool ecal::EcalVetoProcessor::recoil_from_tracking_
private

Definition at line 161 of file EcalVetoProcessor.h.

◆ roc_file_name_

std::string ecal::EcalVetoProcessor::roc_file_name_
private

Definition at line 151 of file EcalVetoProcessor.h.

◆ roc_range_values_

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

Definition at line 113 of file EcalVetoProcessor.h.

◆ rt_

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

Definition at line 172 of file EcalVetoProcessor.h.

◆ shower_rms_

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

Definition at line 122 of file EcalVetoProcessor.h.

122{0};

◆ sim_particles_coll_name_

std::string ecal::EcalVetoProcessor::sim_particles_coll_name_
private

Definition at line 165 of file EcalVetoProcessor.h.

◆ sim_particles_passname_

std::string ecal::EcalVetoProcessor::sim_particles_passname_
private

Definition at line 166 of file EcalVetoProcessor.h.

◆ sp_pass_name_

std::string ecal::EcalVetoProcessor::sp_pass_name_
private

Definition at line 158 of file EcalVetoProcessor.h.

◆ std_layer_hit_

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

Definition at line 126 of file EcalVetoProcessor.h.

126{0};

◆ summed_det_

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

Definition at line 119 of file EcalVetoProcessor.h.

119{0};

◆ summed_tight_iso_

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

Definition at line 120 of file EcalVetoProcessor.h.

120{0};

◆ target_sp_coll_name_

std::string ecal::EcalVetoProcessor::target_sp_coll_name_
private

Definition at line 157 of file EcalVetoProcessor.h.

◆ track_collection_

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

Definition at line 163 of file EcalVetoProcessor.h.

◆ track_pass_name_

std::string ecal::EcalVetoProcessor::track_pass_name_
private

Definition at line 162 of file EcalVetoProcessor.h.

◆ x_std_

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

Definition at line 123 of file EcalVetoProcessor.h.

123{0};

◆ y_std_

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

Definition at line 124 of file EcalVetoProcessor.h.

124{0};

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