Public Member Functions | Private Attributes | List of all members
PDSPHitmonitor_module::PDSPHitMonitorModule Class Reference
Inheritance diagram for PDSPHitmonitor_module::PDSPHitMonitorModule:
art::EDAnalyzer art::detail::Analyzer art::detail::LegacyModule art::Observer art::ModuleBase

Public Member Functions

 PDSPHitMonitorModule (fhicl::ParameterSet const &pset)
 
virtual ~PDSPHitMonitorModule ()
 
void beginJob ()
 
void analyze (const art::Event &evt)
 
void reconfigure (fhicl::ParameterSet const &p)
 
- Public Member Functions inherited from art::EDAnalyzer
 EDAnalyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 EDAnalyzer (Table< Config > const &config)
 
std::string workerType () const
 
- Public Member Functions inherited from art::detail::Analyzer
virtual ~Analyzer () noexcept
 
 Analyzer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 Analyzer (Table< Config > const &config)
 
void doBeginJob (SharedResources const &resources)
 
void doEndJob ()
 
void doRespondToOpenInputFile (FileBlock const &fb)
 
void doRespondToCloseInputFile (FileBlock const &fb)
 
void doRespondToOpenOutputFiles (FileBlock const &fb)
 
void doRespondToCloseOutputFiles (FileBlock const &fb)
 
bool doBeginRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doEndRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doBeginSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEndSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEvent (EventPrincipal &ep, ModuleContext const &mc, std::atomic< std::size_t > &counts_run, std::atomic< std::size_t > &counts_passed, std::atomic< std::size_t > &counts_failed)
 
- Public Member Functions inherited from art::Observer
 ~Observer () noexcept
 
 Observer (Observer const &)=delete
 
 Observer (Observer &&)=delete
 
Observeroperator= (Observer const &)=delete
 
Observeroperator= (Observer &&)=delete
 
void registerProducts (ProductDescriptions &, ModuleDescription const &)
 
void fillDescriptions (ModuleDescription const &)
 
fhicl::ParameterSetID selectorConfig () const
 
- Public Member Functions inherited from art::ModuleBase
virtual ~ModuleBase () noexcept
 
 ModuleBase ()
 
ModuleDescription const & moduleDescription () const
 
void setModuleDescription (ModuleDescription const &)
 
std::array< std::vector< ProductInfo >, NumBranchTypes > const & getConsumables () const
 
void sortConsumables (std::string const &current_process_name)
 
template<typename T , BranchType BT>
ViewToken< T > consumesView (InputTag const &tag)
 
template<typename T , BranchType BT>
ViewToken< T > mayConsumeView (InputTag const &tag)
 

Private Attributes

art::InputTag fTPCHitTag
 
unsigned int fChansPerAPA
 
unsigned int fNofAPA
 
unsigned int fUChanMin
 
unsigned int fVChanMin
 
unsigned int fUChanMax
 
unsigned int fVChanMax
 
unsigned int fZ0ChanMax
 
unsigned int fZ0ChanMin
 
unsigned int fZ1ChanMin
 
unsigned int fZ1ChanMax
 
TH1I * fTotalNHits
 
TH1F * fHitCharge
 
TH1F * fHitRMS
 
TH1F * fHitPeakTime
 
std::vector< TH1I * > fNHitsAPAViewU
 
std::vector< TH1I * > fNHitsAPAViewV
 
std::vector< TH1I * > fNHitsAPAViewZ
 
std::vector< TH1F * > fChargeAPAViewU
 
std::vector< TH1F * > fChargeAPAViewV
 
std::vector< TH1F * > fChargeAPAViewZ
 
std::vector< TH1F * > fRMSAPAViewU
 
std::vector< TH1F * > fRMSAPAViewV
 
std::vector< TH1F * > fRMSAPAViewZ
 
std::vector< TH1F * > fHitPeakTimeAPAViewU
 
std::vector< TH1F * > fHitPeakTimeAPAViewV
 
std::vector< TH1F * > fHitPeakTimeAPAViewZ
 
std::vector< TProfile * > fNHitsAPAViewU_prof
 
std::vector< TProfile * > fNHitsAPAViewV_prof
 
std::vector< TProfile * > fNHitsAPAViewZ_prof
 
std::vector< TProfile * > fChargeAPAViewU_prof
 
std::vector< TProfile * > fChargeAPAViewV_prof
 
std::vector< TProfile * > fChargeAPAViewZ_prof
 
std::vector< TProfile * > fRMSAPAViewU_prof
 
std::vector< TProfile * > fRMSAPAViewV_prof
 
std::vector< TProfile * > fRMSAPAViewZ_prof
 
geo::GeometryCore const * fGeom = &*(art::ServiceHandle<geo::Geometry>())
 
std::vector< unsigned int > fApaLabelNum
 

Additional Inherited Members

- Public Types inherited from art::EDAnalyzer
using WorkerType = WorkerT< EDAnalyzer >
 
using ModuleType = EDAnalyzer
 
- Protected Member Functions inherited from art::Observer
std::string const & processName () const
 
bool wantAllEvents () const noexcept
 
bool wantEvent (ScheduleID id, Event const &e) const
 
Handle< TriggerResultsgetTriggerResults (Event const &e) const
 
 Observer (fhicl::ParameterSet const &config)
 
 Observer (std::vector< std::string > const &select_paths, std::vector< std::string > const &reject_paths, fhicl::ParameterSet const &config)
 
- Protected Member Functions inherited from art::ModuleBase
ConsumesCollectorconsumesCollector ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 

Detailed Description

Definition at line 50 of file PDSPHitMonitor_module.cc.

Constructor & Destructor Documentation

PDSPHitmonitor_module::PDSPHitMonitorModule::PDSPHitMonitorModule ( fhicl::ParameterSet const &  pset)
explicit

Definition at line 124 of file PDSPHitMonitor_module.cc.

124  : EDAnalyzer(pset), fApaLabelNum{3,5,2,6,1,4} {
125  this->reconfigure(pset);
126  }
void reconfigure(fhicl::ParameterSet const &p)
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
PDSPHitmonitor_module::PDSPHitMonitorModule::~PDSPHitMonitorModule ( )
virtual

Definition at line 129 of file PDSPHitMonitor_module.cc.

129 {}

Member Function Documentation

void PDSPHitmonitor_module::PDSPHitMonitorModule::analyze ( const art::Event evt)

Definition at line 265 of file PDSPHitMonitor_module.cc.

265  {
266 
267  // Many options here: gaushit, hitfd, linecluster, lineclusterdc, robusthit, fasthit, etc.
268  auto hitHandle = evt.getHandle<std::vector<recob::Hit> >(fTPCHitTag);
269  if (!hitHandle) {
270  mf::LogWarning("HitMonitor") << " Getting Hits FAIL: " << fTPCHitTag << std::endl;
271  return;
272  }
273 
274  std::vector<art::Ptr<recob::Hit> > hitlist;
275  art::fill_ptr_vector(hitlist, hitHandle);
276 
277  int NHits = hitlist.size();
278 
279  // Arrays ot store number of hits per APA
280  int NHitsPerApa[10][3];
281  for(unsigned int j=0; j<10; j++){
282  for(int k=0; k<3; k++){
283  NHitsPerApa[j][k] = 0;
284  }
285  }
286 
287  int NHitsChannel[20000];
288  for(unsigned int j=0; j<20000; j++){
289  NHitsChannel[j] = 0;
290  }
291 
292  mf::LogVerbatim("HitMonitor") << " Number of hits = " << NHits << std::endl;
293 
294  fTotalNHits->Fill(NHits);
295  for(int i=0; i<NHits; i++){
296  art::Ptr<recob::Hit> hit = hitlist[i];
297 
298  int hit_channel = hit->Channel();
299  unsigned int apa = std::floor(hit->WireID().TPC/2);
300  //unsigned int apa = std::floor(hit_channel/fChansPerAPA);
301 
302  // Protection
303  if(apa > fNofAPA){
304  mf::LogWarning("HitMonitor") << "APA number found (" << apa << ") larger than maximum (" << fNofAPA << "). Skipping hit!" << std::endl;
305  continue;
306  }
307 
308  float hit_charge = hit->Integral();
309  //float hit_chargeSigma = hit->SigmaIntegral();
310  float hit_rms = hit->RMS();
311  float hit_peakT = hit->PeakTime();
312  //float hit_wireID = hit->WireID().Wire;
313  int hit_plane = (int)hit->WireID().Plane;
314 
315  fHitCharge->Fill(hit_charge);
316  fHitRMS->Fill(hit_rms);
317  fHitPeakTime->Fill(hit_peakT);
318 
319  NHitsPerApa[apa][hit_plane]++;
320  NHitsChannel[hit_channel]++;
321 
322  if(hit_plane == 0){
323  fChargeAPAViewU[apa]->Fill(hit_charge);
324  fRMSAPAViewU[apa]->Fill(hit_rms);
325  fHitPeakTimeAPAViewU[apa]->Fill(hit_peakT);
326 
327  fChargeAPAViewU_prof[apa]->Fill(hit_channel, hit_charge, 1);
328  fRMSAPAViewU_prof[apa]->Fill(hit_channel, hit_rms, 1);
329  }
330  else if(hit_plane == 1){
331  fChargeAPAViewV[apa]->Fill(hit_charge);
332  fRMSAPAViewV[apa]->Fill(hit_rms);
333  fHitPeakTimeAPAViewV[apa]->Fill(hit_peakT);
334 
335  fChargeAPAViewV_prof[apa]->Fill(hit_channel, hit_charge, 1);
336  fRMSAPAViewV_prof[apa]->Fill(hit_channel, hit_rms, 1);
337  }
338  else if(hit_plane == 2){
339  fChargeAPAViewZ[apa]->Fill(hit_charge);
340  fRMSAPAViewZ[apa]->Fill(hit_rms);
341  fHitPeakTimeAPAViewZ[apa]->Fill(hit_peakT);
342 
343  fChargeAPAViewZ_prof[apa]->Fill(hit_channel, hit_charge, 1);
344  fRMSAPAViewZ_prof[apa]->Fill(hit_channel, hit_rms, 1);
345  }
346  }
347 
348  // Now fill th number of hits histograms
349  for(unsigned int j=0; j<fNofAPA; j++){
350  fNHitsAPAViewU[j]->Fill(NHitsPerApa[j][0]);
351  fNHitsAPAViewV[j]->Fill(NHitsPerApa[j][1]);
352  fNHitsAPAViewZ[j]->Fill(NHitsPerApa[j][2]);
353 
354  unsigned int UChMin=fUChanMin + j*fChansPerAPA;
355  unsigned int UChMax=fUChanMax + j*fChansPerAPA;
356  unsigned int VChMin=fVChanMin + j*fChansPerAPA;
357  unsigned int VChMax=fVChanMax + j*fChansPerAPA;
358  unsigned int ZChMin=fZ0ChanMin + j*fChansPerAPA;
359  //unsigned int ZChMax=fZ0ChanMax + j*fChansPerAPA;
360  unsigned int ZChMax=fZ1ChanMax + j*fChansPerAPA; //including unused channels
361 
362  for(unsigned int k=UChMin; k<UChMax; k++){
363  fNHitsAPAViewU_prof[j]->Fill(k, NHitsChannel[k], 1);
364  }
365 
366  for(unsigned int k=VChMin; k<VChMax; k++){
367  fNHitsAPAViewV_prof[j]->Fill(k, NHitsChannel[k], 1);
368  }
369 
370  for(unsigned int k=ZChMin; k<ZChMax; k++){
371  fNHitsAPAViewZ_prof[j]->Fill(k, NHitsChannel[k], 1);
372  }
373  }
374 
375  }
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
Handle< PROD > getHandle(SelectorBase const &) const
Definition: DataViewImpl.h:382
geo::WireID WireID() const
Definition: Hit.h:233
float RMS() const
RMS of the hit shape, in tick units.
Definition: Hit.h:220
float Integral() const
Integral under the calibrated signal waveform of the hit, in tick x ADC units.
Definition: Hit.h:224
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:493
Detector simulation of raw signals on wires.
float PeakTime() const
Time of the signal peak, in tick units.
Definition: Hit.h:218
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:297
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:406
raw::ChannelID_t Channel() const
ID of the readout channel the hit was extracted from.
Definition: Hit.h:230
QTextStream & endl(QTextStream &s)
void PDSPHitmonitor_module::PDSPHitMonitorModule::beginJob ( )
virtual

Reimplemented from art::EDAnalyzer.

Definition at line 132 of file PDSPHitMonitor_module.cc.

132  {
133 
135 
136  // Accquiring geometry data
139  mf::LogVerbatim("HitMonitor") << " NTPCs = " << fGeom->NTPC() << " , Ncryostats = " << fGeom->Ncryostats() << " , NofAPA = " << fNofAPA << " , ChansPerAPA = " << fChansPerAPA << std::endl;
140 
141  // loop through channels in the first APA to find the channel boundaries for each view
142  // will adjust for desired APA after
143  fUChanMin = 0;
144  fVChanMin = 0;
145  fUChanMax = 0;
146  fVChanMax = 0;
147  fZ0ChanMax = 0;
148  fZ0ChanMin = 0;
149  fZ1ChanMin = 0;
150  fZ1ChanMax = fChansPerAPA - 1;
151  for(unsigned int c = fUChanMin + 1; c < fZ1ChanMax; c++){
152  if(fGeom->View(c) == geo::kV && fGeom->View(c-1) == geo::kU){
153  fVChanMin = c;
154  fUChanMax = c - 1;
155  }
156  if(fGeom->View(c) == geo::kZ && fGeom->View(c-1) == geo::kV){
157  fZ0ChanMin = c;
158  fVChanMax = c-1;
159  }
160  if(fGeom->View(c) == geo::kZ && fGeom->ChannelToWire(c)[0].TPC == fGeom->ChannelToWire(c-1)[0].TPC + 1){
161  fZ1ChanMin = c;
162  fZ0ChanMax = c-1;
163  }
164  }
165 
166  //unsigned int fNUCh=fUChanMax-fUChanMin+1;
167  //unsigned int fNVCh=fVChanMax-fVChanMin+1;
168  //unsigned int fNZ0Ch=fZ0ChanMax-fZ0ChanMin+1;
169  //unsigned int fNZ1Ch=fZ1ChanMax-fZ1ChanMin+1;
170 
171  //mf::LogVerbatim("HitMonitor") << " U: "<< fNUCh <<" V: "<< fNVCh << " Z0: "<< fNZ0Ch << " Z1: " << fNZ1Ch << std::endl;
172 
173  for(unsigned int i=0;i<fNofAPA;i++){
174 
175  unsigned int j=fApaLabelNum.at(i);
176 
177  unsigned int UChMin=fUChanMin + i*fChansPerAPA;
178  unsigned int UChMax=fUChanMax + i*fChansPerAPA;
179  unsigned int VChMin=fVChanMin + i*fChansPerAPA;
180  unsigned int VChMax=fVChanMax + i*fChansPerAPA;
181  unsigned int ZChMin=fZ0ChanMin + i*fChansPerAPA;
182  //unsigned int ZChMax=fZ0ChanMax + i*fChansPerAPA;
183  unsigned int ZChMax=fZ1ChanMax + i*fChansPerAPA; //including unused channels
184 
185  //std::cout<< "UCh: " << UChMin << " - " << UChMax << std::endl;
186  //std::cout<< "VCh: " << VChMin << " - " << VChMax << std::endl;
187  //std::cout<< "ZCh: " << ZChMin << " - " << ZChMax << std::endl;
188 
189  fNHitsAPAViewU.push_back(tfs->make<TH1I>(Form("NHitsAPA%d_U",j),Form("Number of hits APA%d-U",j),1000,0,20000));
190  fNHitsAPAViewV.push_back(tfs->make<TH1I>(Form("NHitsAPA%d_V",j),Form("Number of hits APA%d-V",j),1000,0,20000));
191  fNHitsAPAViewZ.push_back(tfs->make<TH1I>(Form("NHitsAPA%d_Z",j),Form("Number of hits APA%d-Z",j),1000,0,20000));
192 
193  fChargeAPAViewU.push_back(tfs->make<TH1F>(Form("HitChargeAPA%d_U",j),Form("Hit Charge APA%d-U",j),100,0,1500));
194  fChargeAPAViewV.push_back(tfs->make<TH1F>(Form("HitChargeAPA%d_V",j),Form("Hit Charge APA%d-V",j),100,0,1500));
195  fChargeAPAViewZ.push_back(tfs->make<TH1F>(Form("HitChargeAPA%d_Z",j),Form("Hit Charge APA%d-Z",j),100,0,1500));
196 
197  fRMSAPAViewU.push_back(tfs->make<TH1F>(Form("HitRMSAPA%d_U",j),Form("Hit RMS APA%d-U",j),100,0,10));
198  fRMSAPAViewV.push_back(tfs->make<TH1F>(Form("HitRMSAPA%d_V",j),Form("Hit RMS APA%d-V",j),100,0,10));
199  fRMSAPAViewZ.push_back(tfs->make<TH1F>(Form("HitRMSAPA%d_Z",j),Form("Hit RMS APA%d-Z",j),100,0,10));
200 
201  fHitPeakTimeAPAViewU.push_back(tfs->make<TH1F>(Form("HitPeakTimeAPA%d_U",j),Form("Hit Peak Time APA%d-U",j),100,0,10000));
202  fHitPeakTimeAPAViewV.push_back(tfs->make<TH1F>(Form("HitPeakTimeAPA%d_V",j),Form("Hit Peak Time APA%d-V",j),100,0,10000));
203  fHitPeakTimeAPAViewZ.push_back(tfs->make<TH1F>(Form("HitPeakTimeAPA%d_Z",j),Form("Hit Peak Time APA%d-Z",j),100,0,10000));
204 
205  // U view
206  fNHitsAPAViewU_prof.push_back(tfs->make<TProfile>(Form("fNHitsViewU%d_prof", j),Form("Number of hits distribution vs Channel(Plane U, APA%d)", j), UChMax - UChMin + 1, UChMin-0.5, UChMax+0.5, "s"));
207  fChargeAPAViewU_prof.push_back(tfs->make<TProfile>(Form("fChargeViewU%d_prof", j),Form("Profiled Hit Charge distribution vs Channel(Plane U, APA%d)", j), UChMax - UChMin + 1, UChMin-0.5, UChMax+0.5, "s"));
208  fRMSAPAViewU_prof.push_back(tfs->make<TProfile>(Form("fRMSViewU%d_prof", j),Form("Profiled Hit RMS distribution vs Channel(Plane U, APA%d)", j), UChMax - UChMin + 1, UChMin-0.5, UChMax+0.5, "s"));
209 
210  // V view
211  fNHitsAPAViewV_prof.push_back(tfs->make<TProfile>(Form("fNHitsViewV%d_prof",j),Form("Number of hits distribution vs Channel(Plane V, APA%d)",j), VChMax - VChMin + 1, VChMin-0.5, VChMax+0.5, "s"));
212  fChargeAPAViewV_prof.push_back(tfs->make<TProfile>(Form("fChargeViewV%d_prof",j),Form("Profiled Hit Charge distribution vs Channel(Plane V, APA%d)",j), VChMax - VChMin + 1, VChMin-0.5, VChMax+0.5, "s"));
213  fRMSAPAViewV_prof.push_back(tfs->make<TProfile>(Form("fRMSViewV%d_prof",j),Form("Profiled Hit RMS distribution vs Channel(Plane V, APA%d)",j), VChMax - VChMin + 1, VChMin-0.5, VChMax+0.5, "s"));
214 
215  // Z view
216  fNHitsAPAViewZ_prof.push_back(tfs->make<TProfile>(Form("fNHitsViewZ%d_prof",j),Form("Number of hits distribution vs Channel(Plane Z, APA%d)",j), ZChMax - ZChMin + 1, ZChMin-0.5, ZChMax+0.5, "s"));
217  fChargeAPAViewZ_prof.push_back(tfs->make<TProfile>(Form("fChargeViewZ%d_prof",j),Form("Profiled Hit Charge distribution vs Channel(Plane Z, APA%d)",j), ZChMax - ZChMin + 1, ZChMin-0.5, ZChMax+0.5, "s"));
218  fRMSAPAViewZ_prof.push_back(tfs->make<TProfile>(Form("fRMSViewZ%d_prof",j),Form("Profiled Hit RMS distribution vs Channel(Plane Z, APA%d)",j), ZChMax - ZChMin + 1, ZChMin-0.5, ZChMax+0.5, "s"));
219 
220  // Set titles
221  fNHitsAPAViewU[i]->GetXaxis()->SetTitle("NHits");
222  fNHitsAPAViewV[i]->GetXaxis()->SetTitle("NHits");
223  fNHitsAPAViewZ[i]->GetXaxis()->SetTitle("NHits");
224 
225  fChargeAPAViewU[i]->GetXaxis()->SetTitle("Charge (ADC units)");
226  fChargeAPAViewV[i]->GetXaxis()->SetTitle("Charge (ADC units)");
227  fChargeAPAViewZ[i]->GetXaxis()->SetTitle("Charge (ADC units)");
228 
229  fRMSAPAViewU[i]->GetXaxis()->SetTitle("Hit RMS");
230  fRMSAPAViewV[i]->GetXaxis()->SetTitle("Hit RMS");
231  fRMSAPAViewZ[i]->GetXaxis()->SetTitle("Hit RMS");
232 
233  fHitPeakTimeAPAViewU[i]->GetXaxis()->SetTitle("Hit Peak Time");
234  fHitPeakTimeAPAViewV[i]->GetXaxis()->SetTitle("Hit Peak Time");
235  fHitPeakTimeAPAViewZ[i]->GetXaxis()->SetTitle("Hit Peak Time");
236 
237  fNHitsAPAViewU_prof[i]->GetXaxis()->SetTitle("Channel ID"); fNHitsAPAViewU_prof[i]->GetYaxis()->SetTitle("NHits");
238  fNHitsAPAViewV_prof[i]->GetXaxis()->SetTitle("Channel ID"); fNHitsAPAViewV_prof[i]->GetYaxis()->SetTitle("NHits");
239  fNHitsAPAViewZ_prof[i]->GetXaxis()->SetTitle("Channel ID"); fNHitsAPAViewZ_prof[i]->GetYaxis()->SetTitle("NHits");
240 
241  fChargeAPAViewU_prof[i]->GetXaxis()->SetTitle("Channel ID"); fChargeAPAViewU_prof[i]->GetYaxis()->SetTitle("Charge (ADC units)");
242  fChargeAPAViewV_prof[i]->GetXaxis()->SetTitle("Channel ID"); fChargeAPAViewV_prof[i]->GetYaxis()->SetTitle("Charge (ADC units)");
243  fChargeAPAViewZ_prof[i]->GetXaxis()->SetTitle("Channel ID"); fChargeAPAViewZ_prof[i]->GetYaxis()->SetTitle("Charge (ADC units)");
244 
245  fRMSAPAViewU_prof[i]->GetXaxis()->SetTitle("Channel ID"); fRMSAPAViewU_prof[i]->GetYaxis()->SetTitle("Hit RMS");
246  fRMSAPAViewV_prof[i]->GetXaxis()->SetTitle("Channel ID"); fRMSAPAViewV_prof[i]->GetYaxis()->SetTitle("Hit RMS");
247  fRMSAPAViewZ_prof[i]->GetXaxis()->SetTitle("Channel ID"); fRMSAPAViewZ_prof[i]->GetYaxis()->SetTitle("Hit RMS");
248  }
249 
250  // Summary histograms from all APAs
251  fTotalNHits = tfs->make<TH1I>("fTotalNHits" ,"Total number of hits" ,1000,0,1000000);
252  fHitCharge = tfs->make<TH1F>("fHitCharge" ,"Hit Charge" ,100 ,0,1500);
253  fHitRMS = tfs->make<TH1F>("fHitChargeRMS" ,"Hit RMS" ,100 ,0,10);
254  fHitPeakTime = tfs->make<TH1F>("fHitPeakTime" ,"Hit Peak Time" ,100 ,0,10000);
255 
256  // Set titles
257  fTotalNHits->GetXaxis()->SetTitle("NHits");
258  fHitCharge->GetXaxis()->SetTitle("Charge (ADC units)");
259  fHitRMS->GetXaxis()->SetTitle("Hit RMS");
260  fHitPeakTime->GetXaxis()->SetTitle("Hit Peak Time");
261 
262  }
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
Planes which measure V.
Definition: geo_types.h:130
std::vector< geo::WireID > ChannelToWire(raw::ChannelID_t const channel) const
Returns a list of wires connected to the specified TPC channel.
Planes which measure Z direction.
Definition: geo_types.h:132
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
unsigned int Nchannels() const
Returns the number of TPC readout channels in the detector.
Planes which measure U.
Definition: geo_types.h:129
View_t View(geo::PlaneID const &pid) const
Returns the view (wire orientation) on the channels of specified TPC plane.
unsigned int NTPC(unsigned int cstat=0) const
Returns the total number of TPCs in the specified cryostat.
QTextStream & endl(QTextStream &s)
void PDSPHitmonitor_module::PDSPHitMonitorModule::reconfigure ( fhicl::ParameterSet const &  p)

Definition at line 117 of file PDSPHitMonitor_module.cc.

117  {
118 
119  fTPCHitTag = p.get<art::InputTag>("TPCHitTag", "a:b:c");
120 
121  }
p
Definition: test.py:223

Member Data Documentation

std::vector<unsigned int> PDSPHitmonitor_module::PDSPHitMonitorModule::fApaLabelNum
private

Definition at line 112 of file PDSPHitMonitor_module.cc.

unsigned int PDSPHitmonitor_module::PDSPHitMonitorModule::fChansPerAPA
private

Definition at line 63 of file PDSPHitMonitor_module.cc.

std::vector<TH1F*> PDSPHitmonitor_module::PDSPHitMonitorModule::fChargeAPAViewU
private

Definition at line 85 of file PDSPHitMonitor_module.cc.

std::vector<TProfile*> PDSPHitmonitor_module::PDSPHitMonitorModule::fChargeAPAViewU_prof
private

Definition at line 102 of file PDSPHitMonitor_module.cc.

std::vector<TH1F*> PDSPHitmonitor_module::PDSPHitMonitorModule::fChargeAPAViewV
private

Definition at line 86 of file PDSPHitMonitor_module.cc.

std::vector<TProfile*> PDSPHitmonitor_module::PDSPHitMonitorModule::fChargeAPAViewV_prof
private

Definition at line 103 of file PDSPHitMonitor_module.cc.

std::vector<TH1F*> PDSPHitmonitor_module::PDSPHitMonitorModule::fChargeAPAViewZ
private

Definition at line 87 of file PDSPHitMonitor_module.cc.

std::vector<TProfile*> PDSPHitmonitor_module::PDSPHitMonitorModule::fChargeAPAViewZ_prof
private

Definition at line 104 of file PDSPHitMonitor_module.cc.

geo::GeometryCore const* PDSPHitmonitor_module::PDSPHitMonitorModule::fGeom = &*(art::ServiceHandle<geo::Geometry>())
private

Definition at line 110 of file PDSPHitMonitor_module.cc.

TH1F* PDSPHitmonitor_module::PDSPHitMonitorModule::fHitCharge
private

Definition at line 77 of file PDSPHitMonitor_module.cc.

TH1F* PDSPHitmonitor_module::PDSPHitMonitorModule::fHitPeakTime
private

Definition at line 79 of file PDSPHitMonitor_module.cc.

std::vector<TH1F*> PDSPHitmonitor_module::PDSPHitMonitorModule::fHitPeakTimeAPAViewU
private

Definition at line 93 of file PDSPHitMonitor_module.cc.

std::vector<TH1F*> PDSPHitmonitor_module::PDSPHitMonitorModule::fHitPeakTimeAPAViewV
private

Definition at line 94 of file PDSPHitMonitor_module.cc.

std::vector<TH1F*> PDSPHitmonitor_module::PDSPHitMonitorModule::fHitPeakTimeAPAViewZ
private

Definition at line 95 of file PDSPHitMonitor_module.cc.

TH1F* PDSPHitmonitor_module::PDSPHitMonitorModule::fHitRMS
private

Definition at line 78 of file PDSPHitMonitor_module.cc.

std::vector<TH1I*> PDSPHitmonitor_module::PDSPHitMonitorModule::fNHitsAPAViewU
private

Definition at line 81 of file PDSPHitMonitor_module.cc.

std::vector<TProfile*> PDSPHitmonitor_module::PDSPHitMonitorModule::fNHitsAPAViewU_prof
private

Definition at line 98 of file PDSPHitMonitor_module.cc.

std::vector<TH1I*> PDSPHitmonitor_module::PDSPHitMonitorModule::fNHitsAPAViewV
private

Definition at line 82 of file PDSPHitMonitor_module.cc.

std::vector<TProfile*> PDSPHitmonitor_module::PDSPHitMonitorModule::fNHitsAPAViewV_prof
private

Definition at line 99 of file PDSPHitMonitor_module.cc.

std::vector<TH1I*> PDSPHitmonitor_module::PDSPHitMonitorModule::fNHitsAPAViewZ
private

Definition at line 83 of file PDSPHitMonitor_module.cc.

std::vector<TProfile*> PDSPHitmonitor_module::PDSPHitMonitorModule::fNHitsAPAViewZ_prof
private

Definition at line 100 of file PDSPHitMonitor_module.cc.

unsigned int PDSPHitmonitor_module::PDSPHitMonitorModule::fNofAPA
private

Definition at line 64 of file PDSPHitMonitor_module.cc.

std::vector<TH1F*> PDSPHitmonitor_module::PDSPHitMonitorModule::fRMSAPAViewU
private

Definition at line 89 of file PDSPHitMonitor_module.cc.

std::vector<TProfile*> PDSPHitmonitor_module::PDSPHitMonitorModule::fRMSAPAViewU_prof
private

Definition at line 106 of file PDSPHitMonitor_module.cc.

std::vector<TH1F*> PDSPHitmonitor_module::PDSPHitMonitorModule::fRMSAPAViewV
private

Definition at line 90 of file PDSPHitMonitor_module.cc.

std::vector<TProfile*> PDSPHitmonitor_module::PDSPHitMonitorModule::fRMSAPAViewV_prof
private

Definition at line 107 of file PDSPHitMonitor_module.cc.

std::vector<TH1F*> PDSPHitmonitor_module::PDSPHitMonitorModule::fRMSAPAViewZ
private

Definition at line 91 of file PDSPHitMonitor_module.cc.

std::vector<TProfile*> PDSPHitmonitor_module::PDSPHitMonitorModule::fRMSAPAViewZ_prof
private

Definition at line 108 of file PDSPHitMonitor_module.cc.

TH1I* PDSPHitmonitor_module::PDSPHitMonitorModule::fTotalNHits
private

Definition at line 76 of file PDSPHitMonitor_module.cc.

art::InputTag PDSPHitmonitor_module::PDSPHitMonitorModule::fTPCHitTag
private

Definition at line 61 of file PDSPHitMonitor_module.cc.

unsigned int PDSPHitmonitor_module::PDSPHitMonitorModule::fUChanMax
private

Definition at line 68 of file PDSPHitMonitor_module.cc.

unsigned int PDSPHitmonitor_module::PDSPHitMonitorModule::fUChanMin
private

Definition at line 66 of file PDSPHitMonitor_module.cc.

unsigned int PDSPHitmonitor_module::PDSPHitMonitorModule::fVChanMax
private

Definition at line 69 of file PDSPHitMonitor_module.cc.

unsigned int PDSPHitmonitor_module::PDSPHitMonitorModule::fVChanMin
private

Definition at line 67 of file PDSPHitMonitor_module.cc.

unsigned int PDSPHitmonitor_module::PDSPHitMonitorModule::fZ0ChanMax
private

Definition at line 70 of file PDSPHitMonitor_module.cc.

unsigned int PDSPHitmonitor_module::PDSPHitMonitorModule::fZ0ChanMin
private

Definition at line 71 of file PDSPHitMonitor_module.cc.

unsigned int PDSPHitmonitor_module::PDSPHitMonitorModule::fZ1ChanMax
private

Definition at line 73 of file PDSPHitMonitor_module.cc.

unsigned int PDSPHitmonitor_module::PDSPHitMonitorModule::fZ1ChanMin
private

Definition at line 72 of file PDSPHitMonitor_module.cc.


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