LDMX Software
Public Member Functions | Private Attributes | List of all members
trigscint::TestBeamHitAnalyzer Class Reference

Public Member Functions

 TestBeamHitAnalyzer (const std::string &name, framework::Process &process)
 
void configure (framework::config::Parameters &parameters) override
 Callback for the EventProcessor to configure itself from the given set of parameters.
 
void analyze (const framework::Event &event) override
 Process the event and make histograms or summaries.
 
void onProcessStart () override
 Callback for the EventProcessor to take any necessary action when the processing of events starts, such as creating histograms.
 
void onProcessEnd () override
 Callback for the EventProcessor to take any necessary action when the processing of events finishes, such as calculating job-summary quantities.
 
- Public Member Functions inherited from framework::Analyzer
 Analyzer (const std::string &name, Process &process)
 Class constructor.
 
- 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::vector< std::vector< TH1F * > > vChargeVsTime
 
std::string inputCol_
 
std::string inputPassName_ {""}
 
std::vector< double > peds_
 
int startSample_ {0}
 
int nEv {200}
 
int nChannels {16}
 
TH2F * hEvDisp
 
TH2F * hEvDispPE
 
int fillNb {0}
 
TH1F * hOut [200][16]
 
TH1F * hPE [16]
 
TH1F * hPEinClusters [16]
 
TH2F * hPEVsDelta [16]
 
TH2F * hDeltaPEVsDelta [16]
 
TH2F * hPEmaxVsDelta
 

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::Analyzer
static const int CLASSTYPE {2}
 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 23 of file TestBeamHitAnalyzer.h.

Constructor & Destructor Documentation

◆ TestBeamHitAnalyzer()

trigscint::TestBeamHitAnalyzer::TestBeamHitAnalyzer ( const std::string &  name,
framework::Process process 
)

Definition at line 12 of file TestBeamHitAnalyzer.cxx.

14 : Analyzer(name, process) {}
Analyzer(const std::string &name, Process &process)
Class constructor.

Member Function Documentation

◆ analyze()

void trigscint::TestBeamHitAnalyzer::analyze ( const framework::Event event)
overridevirtual

Process the event and make histograms or summaries.

Parameters
eventThe Event to analyze

Implements framework::Analyzer.

Definition at line 31 of file TestBeamHitAnalyzer.cxx.

31 {
32 const auto channels{
33 event.getCollection<trigscint::TestBeamHit>(inputCol_, inputPassName_)};
34
35 int evNb = event.getEventNumber();
36 // int nChan = channels.size();
37 int leadBar = -1;
38 int subleadBar = -1;
39 float peLead = -1;
40 float peSublead = -1;
41 bool existsIntermediatePE = false;
42
43 for (auto chan : channels) {
44 int bar = chan.getBarID();
45 float PE = chan.getPE();
46 if (evNb < nEv &&
47 bar < nChannels) { // stick within the predefined histogram array
48 hEvDisp->Fill(evNb, bar, PE);
49 } // if within event display range
50 if (PE > peLead) {
51 peSublead = peLead;
52 subleadBar = leadBar;
53 peLead = PE;
54 leadBar = bar;
55 } else if (PE > peSublead) { // need a specific check, bars not sorted in
56 // PE so leadPE might be found before or after
57 peSublead = PE;
58 subleadBar = bar;
59 }
60 hPE[bar]->Fill(PE);
61 if (chan.getQualityFlag() == 0 && bar < 12 && 15 < PE && PE < 40)
62 existsIntermediatePE = true;
63 } // over channels
64
65 if (existsIntermediatePE && fillNb < nEv) {
66 for (auto chan : channels) {
67 int bar = chan.getBarID();
68 if (bar < 12) { // stick within the predefined histogram array
69 float PE = chan.getPE();
70 hEvDispPE->Fill(fillNb, bar, PE);
71 }
72 } // if within event display range
73 fillNb++;
74 } // if PE range triggers writing these event displays
75
76 if (subleadBar == -1) {
77 subleadBar = leadBar;
78 peSublead = peLead;
79 }
80 hPEVsDelta[leadBar]->Fill(leadBar - subleadBar, peLead);
81 hDeltaPEVsDelta[leadBar]->Fill(leadBar - subleadBar, peLead - peSublead);
82
83 // if ( (subleadBar%2) == (leadBar%2) ) // in same layer (or even, hit).
84 // skip if we're not seeing a max across layers
85 // return; -- on the other hand this is evident from the plot: delta is even
86
87 hPEmaxVsDelta->Fill(leadBar - subleadBar, peLead);
88
89 return;
90}
int getBarID() const
Get the bar ID.
This class represents the linearised QIE output from the trigger scintillator, in charge (fC).
Definition TestBeamHit.h:24

References ldmx::TrigScintHit::getBarID().

◆ configure()

void trigscint::TestBeamHitAnalyzer::configure ( framework::config::Parameters parameters)
overridevirtual

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 16 of file TestBeamHitAnalyzer.cxx.

16 {
17 inputCol_ = parameters.getParameter<std::string>("inputCollection");
18 inputPassName_ = parameters.getParameter<std::string>("inputPassName");
19 peds_ = parameters.getParameter<std::vector<double> >("pedestals");
20 startSample_ = parameters.getParameter<int>("startSample");
21
22 std::cout << " [ TestBeamHitAnalyzer ] In configure(), got parameters "
23 << "\n\t inputCollection = " << inputCol_
24 << "\n\t inputPassName = " << inputPassName_
25 << "\n\t startSample = " << startSample_
26 << "\n\t pedestals[0] = " << peds_[0] << "\t." << std::endl;
27
28 return;
29}
T getParameter(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:89

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

◆ onProcessEnd()

void trigscint::TestBeamHitAnalyzer::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 145 of file TestBeamHitAnalyzer.cxx.

145{ return; }

◆ onProcessStart()

void trigscint::TestBeamHitAnalyzer::onProcessStart ( )
overridevirtual

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 92 of file TestBeamHitAnalyzer.cxx.

92 {
93 std::cout << "\n\n Process starts! My analyzer should do something -- like "
94 "print this \n\n"
95 << std::endl;
96
98
99 int nTimeSamp = 40;
100 int PEmax = 400;
101 int nPEbins = 2 * PEmax;
102 // float Qmax = PEmax / (6250. / 4.e6);
103 // float Qmin = -10;
104 // int nQbins = (Qmax - Qmin) / 4;
105
106 for (int iB = 0; iB < nChannels; iB++) {
107 hPE[iB] = new TH1F(Form("hPE_chan%i", iB), Form(";PE, chan%i", iB), nPEbins,
108 0, PEmax);
109 hPEinClusters[iB] = new TH1F(Form("hPEinClusters_chan%i", iB),
110 Form(";PE, chan%i", iB), nPEbins, 0, PEmax);
111 hPEVsDelta[iB] = new TH2F(Form("hPEVsDelta_chan%i", iB),
112 Form(";#Delta_{barID};PE, chan%i has max PE", iB),
113 nChannels + 1, -nChannels / 2 - 0.5,
114 nChannels / 2 + 0.5, nPEbins, 0, PEmax);
115 hDeltaPEVsDelta[iB] = new TH2F(
116 Form("hDeltaPEVsDelta_chan%i", iB),
117 Form(";#Delta_{barID};#Delta_PE, chan%i has max PE", iB), nChannels + 1,
118 -nChannels / 2 - 0.5, nChannels / 2 + 0.5, nPEbins, 0, PEmax);
119 }
120
121 // make event displays for events where there are channels with weird event
122 // hit PE counts
123 for (int iE = 0; iE < nEv; iE++) {
124 for (int iB = 0; iB < nChannels; iB++) {
125 hOut[iE][iB] =
126 new TH1F(Form("hCharge_chan%i_nv%i", iB, iE),
127 Form(";time sample; Q, chan %i, ev %i [fC]", iB, iE),
128 nTimeSamp, -0.5, nTimeSamp - 0.5);
129 }
130 }
131
132 hPEmaxVsDelta =
133 new TH2F("hPEmaxVsDelta", ";#Delta_{barID};PE, max hit", nChannels,
134 -nChannels / 2, nChannels / 2, nPEbins, 0, PEmax);
135 hEvDisp = new TH2F(Form("hEvDisp_ev%i", nEv), ";Event number; Bar ID; PE",
136 nEv, 0.5, nEv + 0.5, nChannels, -0.5, nChannels - 0.5);
137 hEvDispPE = new TH2F("hEvDispPEcut", ";Event number; Bar ID; PE", nEv, 0.5,
138 nEv + 0.5, nChannels, -0.5, nChannels - 0.5);
139
140 fillNb = 0;
141
142 return;
143}
TDirectory * getHistoDirectory()
Access/create a directory in the histogram file for this event processor to create histograms and ana...

Member Data Documentation

◆ fillNb

int trigscint::TestBeamHitAnalyzer::fillNb {0}
private

Definition at line 53 of file TestBeamHitAnalyzer.h.

53{0};

◆ hDeltaPEVsDelta

TH2F* trigscint::TestBeamHitAnalyzer::hDeltaPEVsDelta[16]
private

Definition at line 61 of file TestBeamHitAnalyzer.h.

◆ hEvDisp

TH2F* trigscint::TestBeamHitAnalyzer::hEvDisp
private

Definition at line 50 of file TestBeamHitAnalyzer.h.

◆ hEvDispPE

TH2F* trigscint::TestBeamHitAnalyzer::hEvDispPE
private

Definition at line 51 of file TestBeamHitAnalyzer.h.

◆ hOut

TH1F* trigscint::TestBeamHitAnalyzer::hOut[200][16]
private

Definition at line 56 of file TestBeamHitAnalyzer.h.

◆ hPE

TH1F* trigscint::TestBeamHitAnalyzer::hPE[16]
private

Definition at line 58 of file TestBeamHitAnalyzer.h.

◆ hPEinClusters

TH1F* trigscint::TestBeamHitAnalyzer::hPEinClusters[16]
private

Definition at line 59 of file TestBeamHitAnalyzer.h.

◆ hPEmaxVsDelta

TH2F* trigscint::TestBeamHitAnalyzer::hPEmaxVsDelta
private

Definition at line 62 of file TestBeamHitAnalyzer.h.

◆ hPEVsDelta

TH2F* trigscint::TestBeamHitAnalyzer::hPEVsDelta[16]
private

Definition at line 60 of file TestBeamHitAnalyzer.h.

◆ inputCol_

std::string trigscint::TestBeamHitAnalyzer::inputCol_
private

Definition at line 41 of file TestBeamHitAnalyzer.h.

◆ inputPassName_

std::string trigscint::TestBeamHitAnalyzer::inputPassName_ {""}
private

Definition at line 42 of file TestBeamHitAnalyzer.h.

42{""};

◆ nChannels

int trigscint::TestBeamHitAnalyzer::nChannels {16}
private

Definition at line 48 of file TestBeamHitAnalyzer.h.

48{16};

◆ nEv

int trigscint::TestBeamHitAnalyzer::nEv {200}
private

Definition at line 47 of file TestBeamHitAnalyzer.h.

47{200};

◆ peds_

std::vector<double> trigscint::TestBeamHitAnalyzer::peds_
private

Definition at line 43 of file TestBeamHitAnalyzer.h.

◆ startSample_

int trigscint::TestBeamHitAnalyzer::startSample_ {0}
private

Definition at line 44 of file TestBeamHitAnalyzer.h.

44{0};

◆ vChargeVsTime

std::vector<std::vector<TH1F*> > trigscint::TestBeamHitAnalyzer::vChargeVsTime
private

Definition at line 38 of file TestBeamHitAnalyzer.h.


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