LDMX Software
Public Member Functions | Private Attributes | List of all members
trigger::TrigElectronProducer Class Reference

Public Member Functions

 TrigElectronProducer (const std::string &name, framework::Process &process)
 
virtual void configure (framework::config::Parameters &ps)
 Callback for the EventProcessor to configure itself from the given set of parameters.
 
virtual void produce (framework::Event &event)
 Process the event and put new data products into it.
 
virtual void onProcessStart ()
 Callback for the EventProcessor to take any necessary action when the processing of events starts, such as creating histograms.
 
virtual void onProcessEnd ()
 Callback for the EventProcessor to take any necessary action when the processing of events finishes, such as calculating job-summary quantities.
 
void setupMaps (bool isX)
 
float getP (bool isX, float e, float d)
 
float getPx (float e, float d)
 
float getPy (float e, float d)
 
- Public Member Functions inherited from framework::Producer
 Producer (const std::string &name, Process &process)
 Class constructor.
 
virtual void beforeNewRun (ldmx::RunHeader &header)
 Handle allowing producers to modify run headers before the run begins.
 
- Public Member Functions inherited from framework::EventProcessor
 EventProcessor (const std::string &name, Process &process)
 Class constructor.
 
virtual ~EventProcessor ()
 Class destructor.
 
virtual void onNewRun (const ldmx::RunHeader &runHeader)
 Callback for the EventProcessor to take any necessary action when the run being processed changes.
 
virtual void onFileOpen (EventFile &eventFile)
 Callback for the EventProcessor to take any necessary action when a new event input ROOT file is opened.
 
virtual void onFileClose (EventFile &eventFile)
 Callback for the EventProcessor to take any necessary action when a event input ROOT file is closed.
 
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 Attributes

std::string spCollName_
 
std::string clusterCollName_
 
std::string eleCollName_
 
std::string propMapName_
 
TProfile2D * propMapx_ {nullptr}
 
TProfile2D * propMapy_ {nullptr}
 
std::vector< TF1 * > fitsX_ {}
 
std::vector< TF1 * > fitsY_ {}
 

Additional Inherited Members

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

Detailed Description

Definition at line 27 of file TrigElectronProducer.h.

Constructor & Destructor Documentation

◆ TrigElectronProducer()

trigger::TrigElectronProducer::TrigElectronProducer ( const std::string &  name,
framework::Process process 
)
inline

Definition at line 29 of file TrigElectronProducer.h.

30 : framework::Producer(name, process) {}
Base class for a module which produces a data product.

Member Function Documentation

◆ configure()

void trigger::TrigElectronProducer::configure ( framework::config::Parameters parameters)
virtual

Callback for the EventProcessor to configure itself from the given set of parameters.

The parameters a processor has access to are the member variables of the python class in the sequence that has className equal to the EventProcessor class name.

For an example, look at MyProcessor.

Parameters
parametersParameters for configuration.

Reimplemented from framework::EventProcessor.

Definition at line 13 of file TrigElectronProducer.cxx.

13 {
14 spCollName_ = ps.getParameter<std::string>("scoringPlaneCollName");
15 clusterCollName_ = ps.getParameter<std::string>("clusterCollName");
16 eleCollName_ = ps.getParameter<std::string>("eleCollName");
17 propMapName_ = ps.getParameter<std::string>("propMapName");
18}

References framework::config::Parameters::getParameter().

◆ getP()

float trigger::TrigElectronProducer::getP ( bool  isX,
float  e,
float  d 
)

Definition at line 129 of file TrigElectronProducer.cxx.

129 {
130 if (fabs(d) > 300) return 0; // something has gone very wrong
131 const bool debug = false;
132 TProfile2D* prof = isX ? propMapx_ : propMapy_;
133 if (!prof) {
134 if (debug) std::cout << "null pointer" << std::endl;
135 return 0;
136 }
137 int bin1{-9999}, bin2{-9999};
138 bin1 = prof->GetXaxis()->FindBin(e);
139 float frac = e - prof->GetXaxis()->GetBinCenter(bin1);
140 float diff = fabs(prof->GetXaxis()->GetBinCenter(bin1) -
141 prof->GetXaxis()->GetBinCenter(bin2));
142 bin2 = diff > 0 ? bin1 + 1 : bin1 - 1;
143 bin1 = std::max(1, std::min(bin1, prof->GetXaxis()->GetNbins()));
144 bin2 = std::max(1, std::min(bin2, prof->GetXaxis()->GetNbins()));
145 float res1, res2;
146 if (isX) {
147 res1 = fitsX_[bin1 - 1]->GetX(d);
148 res2 = fitsX_[bin2 - 1]->GetX(d);
149 } else {
150 res1 = fitsY_[bin1 - 1]->GetX(d);
151 res2 = fitsY_[bin2 - 1]->GetX(d);
152 }
153
154 if (debug)
155 printf("%f %f %f %f :: %f %f %f \n", d, e,
156 prof->GetXaxis()->GetBinCenter(bin1),
157 prof->GetXaxis()->GetBinCenter(bin2), res1, res2,
158 abs(frac / diff) * res2 + (1 - abs(frac / diff)) * res1);
159 return e * (abs(frac / diff) * res2 + (1 - abs(frac / diff)) * res1);
160}

◆ getPx()

float trigger::TrigElectronProducer::getPx ( float  e,
float  d 
)
inline

Definition at line 42 of file TrigElectronProducer.h.

42{ return getP(true, e, d); }

◆ getPy()

float trigger::TrigElectronProducer::getPy ( float  e,
float  d 
)
inline

Definition at line 43 of file TrigElectronProducer.h.

43{ return getP(false, e, d); }

◆ onProcessEnd()

void trigger::TrigElectronProducer::onProcessEnd ( )
virtual

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 178 of file TrigElectronProducer.cxx.

178 {
179 ldmx_log(debug) << "Process ends!";
180
181 // free maps
182 for (int i = 0; i < fitsX_.size(); i++) delete fitsX_[i];
183 for (int i = 0; i < fitsY_.size(); i++) delete fitsY_[i];
184
185 return;
186}

◆ onProcessStart()

void trigger::TrigElectronProducer::onProcessStart ( )
virtual

Callback for the EventProcessor to take any necessary action when the processing of events starts, such as creating histograms.

Reimplemented from framework::EventProcessor.

Definition at line 162 of file TrigElectronProducer.cxx.

162 {
163 ldmx_log(debug) << "Process starts!";
164
165 TFile* f = new TFile(propMapName_.c_str(), "read");
166 propMapx_ = (TProfile2D*)f->Get("profx");
167 propMapx_->SetDirectory(0);
168 propMapy_ = (TProfile2D*)f->Get("profy");
169 propMapy_->SetDirectory(0);
170 f->Close();
171
172 setupMaps(0); // X
173 setupMaps(1); // Y
174
175 return;
176}

◆ produce()

void trigger::TrigElectronProducer::produce ( framework::Event event)
virtual

Process the event and put new data products into it.

Parameters
eventThe Event to process.

Implements framework::Producer.

Definition at line 20 of file TrigElectronProducer.cxx.

20 {
21 if (!event.exists(clusterCollName_)) return;
22 auto ecalClusters{
23 event.getObject<TrigCaloClusterCollection>(clusterCollName_)};
24
25 if (!event.exists(spCollName_)) return;
26 const std::vector<ldmx::SimTrackerHit> TargetSPHit =
27 event.getCollection<ldmx::SimTrackerHit>(spCollName_);
28 // ldmx::SimTrackerHit targetPrimary;
29 // std::map<int,int> tk_to_iTargetSPHit;
30 float xT = 0, yT = 0;
31 for (const auto& hit : TargetSPHit) {
32 if (hit.getTrackID() != 1) continue;
33 if (!(abs(hit.getPdgID()) == 11)) continue;
34 auto xyz = hit.getPosition();
35 if (xyz[2] < 0 || xyz[2] > 1) continue; // select one sp
36 xT = xyz[0];
37 yT = xyz[1];
38 // more details:
39 // https://github.com/LDMX-Software/ldmx-analysis/blob/ch/dev/src/TriggerAnalyzer.cxx
40 }
41
42 // std::vector<ldmx::SimParticle> eles;
43 TrigParticleCollection eles;
44 for (const auto& clus : ecalClusters) {
45 TrigParticle el;
46 // ldmx::SimParticle el;
47
48 // We fit for the response in an eGun sample like this:
49 // Events->Draw("Electron_e[0]/Truth_e[0]:Truth_e[0]","Truth_e[0]>500 &&
50 // Truth_e[0]<3500","prof")
51 // --> update fit to use the electron energy at the ecal face (not target)
52 // -->
53 // Events->Draw("Electron_eClus[0]/TruthEcal_e[0]:TruthEcal_e[0]","TruthEcal_e[0]>500
54 // && TruthEcal_e[0]<3500","prof") TF1* func = new
55 // TF1("func","[0]/x+[1]",500,3500); htemp->Fit(func) float calib_e =
56 // clus.e() / (98.099700/clus.e() + 0.77202700);
57 float calib_e =
58 clus.e() / (184.103 / clus.e() + 0.753756); // divide by response
59
60 // const float dX = clus.x() - xT;
61 // float R = calib_e*(11500/4000.); // convert 11.5m/4GeV (in mm/MeV)
62 // float zd = 240.; // 240mm z detector
63 // float a = dX/zd;
64 // float b = (dX*dX+zd*zd)/(2*R*zd);
65 // float pred_px = clus.e() * (-b+a*sqrt(1+a*a-b*b))/(1+a*a);
66 // float pred_py = clus.e()/4e3 * (clus.y() - yT) * 13.3; //mev/mm
67
68 // Events->Draw("TruthEcal_y-Truth_y -
69 // Electron_dy*240/Electron_zClus","TruthEcal_e>1e3")
70
71 float pred_px = clus.z() > 0
72 ? getPx(calib_e, (240. / clus.z()) * (clus.x() - xT))
73 : 0; // adjust, tracing back to ecal face
74 float pred_py =
75 clus.z() > 0 ? getPy(calib_e, (240. / clus.z()) * (clus.y() - yT)) : 0;
76
77 // std::cout << "expectX " << pred_px << " versus " << fit_px << std::endl;
78 // std::cout << "expectY " << pred_py << " versus " << fit_py << std::endl;
79 float pred_pz = sqrt(
80 std::max(pow(calib_e, 2) - (pow(pred_px, 2) + pow(pred_py, 2)), 0.));
81
82 // produce el
83 Point targ(xT, yT, 0);
84 LorentzVector p4(pred_px, pred_py, pred_pz, calib_e);
85 eles.emplace_back(p4, targ, 11);
86 Point calo(clus.x(), clus.y(),
87 clus.z()); // n.b. "Z" is not at the ECal Face!
88 eles.back().setEndPoint(
89 calo); // clus.x(), clus.y(), clus.z()); how to handle?
90 eles.back().setClusEnergy(clus.e());
91 eles.back().setClusTP(clus.nTP());
92 eles.back().setClusDepth(clus.depth());
93
94 // el.setEnergy(clus.e());
95 // el.setPdgID(11);
96 // el.setCharge(-1);
97 // el.setMass(0.000511);
98 // el.setVertex(0, xT, yT);
99 // el.setMomentum(pred_px, pred_py, pred_pz);
100 // eles.push_back(el);
101 }
102 event.add(eleCollName_, eles);
103}
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
Represents a simulated tracker hit in the simulation.
std::vector< float > getPosition() const
Get the XYZ position of the hit [mm].

References framework::Event::exists(), and ldmx::SimTrackerHit::getPosition().

◆ setupMaps()

void trigger::TrigElectronProducer::setupMaps ( bool  isX)

Definition at line 105 of file TrigElectronProducer.cxx.

105 {
106 TProfile2D* prof = isX ? propMapx_ : propMapy_;
107 const int N = prof->GetXaxis()->GetNbins();
108 for (int i = 1; i <= N; i++) {
109 TF1* func = new TF1("func", "pol1", -20, 20); // going to fit for px/e
110 TProfile* proj = prof->ProfileY("h", i, i);
111 proj->Fit("func", "q", "", -1, 1);
112 bool debug = false;
113 if (debug) {
114 TCanvas c("c", "");
115 proj->Draw();
116 func->Draw("same");
117 c.SaveAs(TString::Format("debugFit_%d_%d.pdf", int(isX), i));
118 }
119 delete proj;
120 if (isX) {
121 fitsX_.push_back(func);
122 } else {
123 fitsY_.push_back(func);
124 }
125 }
126 return;
127}

Member Data Documentation

◆ clusterCollName_

std::string trigger::TrigElectronProducer::clusterCollName_
private

Definition at line 49 of file TrigElectronProducer.h.

◆ eleCollName_

std::string trigger::TrigElectronProducer::eleCollName_
private

Definition at line 51 of file TrigElectronProducer.h.

◆ fitsX_

std::vector<TF1*> trigger::TrigElectronProducer::fitsX_ {}
private

Definition at line 56 of file TrigElectronProducer.h.

56{};

◆ fitsY_

std::vector<TF1*> trigger::TrigElectronProducer::fitsY_ {}
private

Definition at line 57 of file TrigElectronProducer.h.

57{};

◆ propMapName_

std::string trigger::TrigElectronProducer::propMapName_
private

Definition at line 53 of file TrigElectronProducer.h.

◆ propMapx_

TProfile2D* trigger::TrigElectronProducer::propMapx_ {nullptr}
private

Definition at line 54 of file TrigElectronProducer.h.

54{nullptr};

◆ propMapy_

TProfile2D* trigger::TrigElectronProducer::propMapy_ {nullptr}
private

Definition at line 55 of file TrigElectronProducer.h.

55{nullptr};

◆ spCollName_

std::string trigger::TrigElectronProducer::spCollName_
private

Definition at line 47 of file TrigElectronProducer.h.


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