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

138 {
139 cell_map_.clear();
140 cell_map_tight_iso_.clear();
141 bdt_features_.clear();
142
143 n_readout_hits_ = 0;
144 summed_det_ = 0;
145 summed_tight_iso_ = 0;
146 max_cell_dep_ = 0;
147 shower_rms_ = 0;
148 x_std_ = 0;
149 y_std_ = 0;
150 avg_layer_hit_ = 0;
151 std_layer_hit_ = 0;
152 deepest_layer_hit_ = 0;
153 ecal_back_energy_ = 0;
155 ep_ang_ = 0;
157 ep_sep_ = 0;
158 ep_dot_ = 0;
160
161 std::fill(ecal_layer_edep_raw_.begin(), ecal_layer_edep_raw_.end(), 0);
162 std::fill(ecal_layer_edep_readout_.begin(), ecal_layer_edep_readout_.end(),
163 0);
164 std::fill(ecal_layer_time_.begin(), ecal_layer_time_.end(), 0);
165}
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 collection_name_ = parameters.get<std::string>("collection_name");
130 rec_pass_name_ = parameters.get<std::string>("rec_pass_name");
131 rec_coll_name_ = parameters.get<std::string>("rec_coll_name");
132 recoil_from_tracking_ = parameters.get<bool>("recoil_from_tracking");
133 track_collection_ = parameters.get<std::string>("track_collection");
134 track_pass_name_ = parameters.get<std::string>("track_pass_name", "");
135 inverse_skim_ = parameters.get<bool>("inverse_skim");
136}
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 1027 of file EcalVetoProcessor.cxx.

1029 {
1030 for (const ldmx::EcalHit &hit : ecal_rec_hits) {
1031 ldmx::EcalID id(hit.getID());
1032 cellMap.emplace(id, hit.getEnergy());
1033 }
1034}
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 1036 of file EcalVetoProcessor.cxx.

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

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

1080 {
1081 std::vector<XYCoords> positions;
1082 for (int i_layer = 0; i_layer < n_ecal_layers_; i_layer++) {
1083 float pos_x =
1084 position[0] + (momentum[0] / momentum[2]) *
1085 (geometry_->getZPosition(i_layer) - position[2]);
1086 float pos_y =
1087 position[1] + (momentum[1] / momentum[2]) *
1088 (geometry_->getZPosition(i_layer) - position[2]);
1089 positions.push_back(std::make_pair(pos_x, pos_y));
1090 }
1091 return positions;
1092}
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 937 of file EcalVetoProcessor.cxx.

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

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

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

169{"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 174 of file EcalVetoProcessor.h.

Referenced by produce().

◆ inverse_skim_

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

Definition at line 166 of file EcalVetoProcessor.h.

166{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 171 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_passname_

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

Definition at line 165 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: