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

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

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

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

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

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

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

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

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

167{"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.

◆ 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 172 of file EcalVetoProcessor.h.

Referenced by produce().

◆ inverse_skim_

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

Definition at line 164 of file EcalVetoProcessor.h.

164{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 158 of file EcalVetoProcessor.h.

◆ rec_pass_name_

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

Definition at line 157 of file EcalVetoProcessor.h.

◆ recoil_from_tracking_

bool ecal::EcalVetoProcessor::recoil_from_tracking_
private

Definition at line 159 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 169 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 163 of file EcalVetoProcessor.h.

◆ sp_pass_name_

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

Definition at line 156 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};

◆ track_collection_

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

Definition at line 161 of file EcalVetoProcessor.h.

◆ track_pass_name_

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

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