LDMX Software
TrigElectronProducer.cxx
2
4// #include "SimCore/Event/SimParticle.h"
6#include "TCanvas.h"
7#include "TString.h"
8#include "Trigger/Event/TrigCaloCluster.h"
9#include "Trigger/Event/TrigParticle.h"
10
11namespace trigger {
12
14 sp_coll_name_ = ps.get<std::string>("scoringPlaneCollName");
15 cluster_coll_name_ = ps.get<std::string>("clusterCollName");
16 ele_coll_name_ = ps.get<std::string>("eleCollName");
17 prop_map_name_ = ps.get<std::string>("propMapName");
18 target_sp_passname_ = ps.get<std::string>("target_sp_passname");
19 cluster_coll_passname_ = ps.get<std::string>("cluster_coll_passname");
20 cluster_collname_events_passname_ =
21 ps.get<std::string>("cluster_collname_events_passname");
22 sp_collname_events_passname_ =
23 ps.get<std::string>("sp_collname_events_passname_");
24}
25
27 if (!event.exists(cluster_coll_name_, cluster_collname_events_passname_))
28 return;
29 auto ecal_clusters{event.getObject<TrigCaloClusterCollection>(
30 cluster_coll_name_, cluster_coll_passname_)};
31
32 if (!event.exists(sp_coll_name_, sp_collname_events_passname_)) return;
33 const std::vector<ldmx::SimTrackerHit> target_sp_hit =
34 event.getCollection<ldmx::SimTrackerHit>(sp_coll_name_,
35 target_sp_passname_);
36 // ldmx::SimTrackerHit targetPrimary;
37 // std::map<int,int> tk_to_iTargetSPHit;
38 float x_t = 0, y_t = 0;
39 for (const auto& hit : target_sp_hit) {
40 if (hit.getTrackID() != 1) continue;
41 if (!(abs(hit.getPdgID()) == 11)) continue;
42 auto xyz = hit.getPosition();
43 if (xyz[2] < 0 || xyz[2] > 1) continue; // select one sp
44 x_t = xyz[0];
45 y_t = xyz[1];
46 // more details:
47 // https://github.com/LDMX-Software/ldmx-analysis/blob/ch/dev/src/TriggerAnalyzer.cxx
48 }
49
50 // std::vector<ldmx::SimParticle> eles;
51 TrigParticleCollection eles;
52 for (const auto& clus : ecal_clusters) {
53 TrigParticle el;
54 // ldmx::SimParticle el;
55
56 // We fit for the response in an eGun sample like this:
57 // Events->Draw("Electron_e[0]/Truth_e[0]:Truth_e[0]","Truth_e[0]>500 &&
58 // Truth_e[0]<3500","prof")
59 // --> update fit to use the electron energy at the ecal face (not target)
60 // -->
61 // Events->Draw("Electron_eClus[0]/TruthEcal_e[0]:TruthEcal_e[0]","TruthEcal_e[0]>500
62 // && TruthEcal_e[0]<3500","prof") TF1* func = new
63 // TF1("func","[0]/x+[1]",500,3500); htemp->Fit(func) float calib_e =
64 // clus.e() / (98.099700/clus.e() + 0.77202700);
65 float calib_e =
66 clus.e() / (184.103 / clus.e() + 0.753756); // divide by response
67
68 // const float dX = clus.x() - xT;
69 // float R = calib_e*(11500/4000.); // convert 11.5m/4GeV (in mm/MeV)
70 // float zd = 240.; // 240mm z detector
71 // float a = dX/zd;
72 // float b = (dX*dX+zd*zd)/(2*R*zd);
73 // float pred_px = clus.e() * (-b+a*sqrt(1+a*a-b*b))/(1+a*a);
74 // float pred_py = clus.e()/4e3 * (clus.y() - yT) * 13.3; //mev/mm
75
76 // Events->Draw("TruthEcal_y-Truth_y -
77 // Electron_dy*240/Electron_zClus","TruthEcal_e>1e3")
78
79 float pred_px = clus.z() > 0
80 ? getPx(calib_e, (240. / clus.z()) * (clus.x() - x_t))
81 : 0; // adjust, tracing back to ecal face
82 float pred_py =
83 clus.z() > 0 ? getPy(calib_e, (240. / clus.z()) * (clus.y() - y_t)) : 0;
84
85 // std::cout << "expectX " << pred_px << " versus " << fit_px << std::endl;
86 // std::cout << "expectY " << pred_py << " versus " << fit_py << std::endl;
87 float pred_pz = sqrt(
88 std::max(pow(calib_e, 2) - (pow(pred_px, 2) + pow(pred_py, 2)), 0.));
89
90 // produce el
91 Point targ(x_t, y_t, 0);
92 LorentzVector p4(pred_px, pred_py, pred_pz, calib_e);
93 eles.emplace_back(p4, targ, 11);
94 Point calo(clus.x(), clus.y(),
95 clus.z()); // n.b. "Z" is not at the ECal Face!
96 eles.back().setEndPoint(
97 calo); // clus.x(), clus.y(), clus.z()); how to handle?
98 eles.back().setClusEnergy(clus.e());
99 eles.back().setClusTP(clus.nTP());
100 eles.back().setClusDepth(clus.depth());
101
102 // el.setEnergy(clus.e());
103 // el.setPdgID(11);
104 // el.setCharge(-1);
105 // el.setMass(0.000511);
106 // el.setVertex(0, xT, yT);
107 // el.setMomentum(pred_px, pred_py, pred_pz);
108 // eles.push_back(el);
109 }
110 event.add(ele_coll_name_, eles);
111}
112
113void TrigElectronProducer::setupMaps(bool isX) {
114 TProfile2D* prof = isX ? prop_mapx_ : prop_mapy_;
115 const int n = prof->GetXaxis()->GetNbins();
116 for (int i = 1; i <= n; i++) {
117 TF1* func = new TF1("func", "pol1", -20, 20); // going to fit for px/e
118 TProfile* proj = prof->ProfileY("h", i, i);
119 proj->Fit("func", "q", "", -1, 1);
120 bool debug = false;
121 if (debug) {
122 TCanvas c("c", "");
123 proj->Draw();
124 func->Draw("same");
125 c.SaveAs(TString::Format("debug_fit_%d_%d.pdf", int(isX), i));
126 }
127 delete proj;
128 if (isX) {
129 fits_x_.push_back(func);
130 } else {
131 fits_y_.push_back(func);
132 }
133 }
134 return;
135}
136
137float TrigElectronProducer::getP(bool isX, float e, float d) {
138 if (fabs(d) > 300) return 0; // something has gone very wrong
139 const bool debug = false;
140 TProfile2D* prof = isX ? prop_mapx_ : prop_mapy_;
141 if (!prof) {
142 if (debug) std::cout << "null pointer" << std::endl;
143 return 0;
144 }
145 int bin1{-9999}, bin2{-9999};
146 bin1 = prof->GetXaxis()->FindBin(e);
147 float frac = e - prof->GetXaxis()->GetBinCenter(bin1);
148 float diff = fabs(prof->GetXaxis()->GetBinCenter(bin1) -
149 prof->GetXaxis()->GetBinCenter(bin2));
150 bin2 = diff > 0 ? bin1 + 1 : bin1 - 1;
151 bin1 = std::max(1, std::min(bin1, prof->GetXaxis()->GetNbins()));
152 bin2 = std::max(1, std::min(bin2, prof->GetXaxis()->GetNbins()));
153 float res1, res2;
154 if (isX) {
155 res1 = fits_x_[bin1 - 1]->GetX(d);
156 res2 = fits_x_[bin2 - 1]->GetX(d);
157 } else {
158 res1 = fits_y_[bin1 - 1]->GetX(d);
159 res2 = fits_y_[bin2 - 1]->GetX(d);
160 }
161
162 if (debug)
163 printf("%f %f %f %f :: %f %f %f \n", d, e,
164 prof->GetXaxis()->GetBinCenter(bin1),
165 prof->GetXaxis()->GetBinCenter(bin2), res1, res2,
166 std::abs(frac / diff) * res2 + (1 - std::abs(frac / diff)) * res1);
167 return e *
168 (std::abs(frac / diff) * res2 + (1 - std::abs(frac / diff)) * res1);
169}
170
172 ldmx_log(debug) << "Process starts!";
173
174 TFile* f = new TFile(prop_map_name_.c_str(), "read");
175 prop_mapx_ = (TProfile2D*)f->Get("profx");
176 prop_mapx_->SetDirectory(0);
177 prop_mapy_ = (TProfile2D*)f->Get("profy");
178 prop_mapy_->SetDirectory(0);
179 f->Close();
180
181 setupMaps(0); // X
182 setupMaps(1); // Y
183
184 return;
185}
186
188 ldmx_log(debug) << "Process ends!";
189
190 // free maps
191 for (int i = 0; i < fits_x_.size(); i++) delete fits_x_[i];
192 for (int i = 0; i < fits_y_.size(); i++) delete fits_y_[i];
193
194 return;
195}
196
197} // namespace trigger
198
Class that translates raw positions of ECal module hits into cells in a hexagonal readout.
#define DECLARE_PRODUCER(CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
Class which encapsulates information from a hit in a simulated tracking detector.
Trigger electron reconstruction algorithm.
Implements an event buffer system for storing event data.
Definition Event.h:42
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
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:29
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78
Represents a simulated tracker hit in the simulation.
virtual void onProcessStart()
Callback for the EventProcessor to take any necessary action when the processing of events starts,...
virtual void produce(framework::Event &event)
Process the event and put new data products into it.
virtual void onProcessEnd()
Callback for the EventProcessor to take any necessary action when the processing of events finishes,...
virtual void configure(framework::config::Parameters &ps)
Callback for the EventProcessor to configure itself from the given set of parameters.
Class for particles reconstructed by the trigger system.