LArSoft  v10_04_05
Liquid Argon Software toolkit - https://larsoft.org/
evdb_tool::DrawWireHist Class Reference
Inheritance diagram for evdb_tool::DrawWireHist:
evdb_tool::IWaveformDrawer

Public Member Functions

 DrawWireHist (const fhicl::ParameterSet &pset)
 
void configure (const fhicl::ParameterSet &pset) override
 
void Fill (evdb::View2D &, raw::ChannelID_t &, float, float) override
 
void Draw (const std::string &, float, float) override
 
float getMaximum () const override
 
float getMinimum () const override
 

Private Member Functions

void BookHistogram (raw::ChannelID_t &, float, float)
 

Private Attributes

float fMaximum
 
float fMinimum
 
std::vector< int > fColorMap
 
std::unordered_map< std::string, std::unique_ptr< TH1F > > fRecoHistMap
 

Detailed Description

Definition at line 24 of file DrawWireHist_tool.cc.

Constructor & Destructor Documentation

evdb_tool::DrawWireHist::DrawWireHist ( const fhicl::ParameterSet pset)
explicit

Definition at line 47 of file DrawWireHist_tool.cc.

References configure().

48  {
49  configure(pset);
50  }
void configure(const fhicl::ParameterSet &pset) override

Member Function Documentation

void evdb_tool::DrawWireHist::BookHistogram ( raw::ChannelID_t channel,
float  startTick,
float  numTicks 
)
private

Definition at line 130 of file DrawWireHist_tool.cc.

References DEFINE_ART_CLASS_TOOL, fRecoHistMap, evd::ColorDrawingOptions::fRecoQHigh, evd::ColorDrawingOptions::fRecoQLow, evd::RecoDrawingOptions::fWireLabels, and Get.

Referenced by Fill(), and getMinimum().

131  {
135  auto const& wireReadoutGeom = art::ServiceHandle<geo::WireReadout const>()->Get();
136 
137  // Get rid of the previous histograms
138  fRecoHistMap.clear();
139 
140  // Now add a histogram for each of the wire labels
141  for (auto& tag : recoOpt->fWireLabels) {
142  // figure out the signal type for this plane, assume that
143  // plane n in each TPC/cryostat has the same type
144  geo::SigType_t sigType = wireReadoutGeom.SignalType(channel);
145  std::string tagString(tag.encode());
146  int numBins = numTicks;
147 
148  fRecoHistMap[tagString] = std::make_unique<TH1F>(
149  "fCALTQHisto", ";t [ticks];q [ADC]", numBins, startTick, startTick + numTicks);
150 
151  TH1F* histPtr = fRecoHistMap.at(tagString).get();
152 
153  histPtr->SetMaximum(cst->fRecoQHigh[(size_t)sigType]);
154  histPtr->SetMinimum(cst->fRecoQLow[(size_t)sigType]);
155 
156  histPtr->SetLineColor(kBlue);
157  histPtr->SetLineWidth(1);
158 
159  histPtr->GetXaxis()->SetLabelSize(0.10); // was 0.15
160  histPtr->GetXaxis()->SetLabelOffset(0.01); // was 0.00
161  histPtr->GetXaxis()->SetTitleSize(0.10); // was 0.15
162  histPtr->GetXaxis()->SetTitleOffset(0.60); // was 0.80
163 
164  histPtr->GetYaxis()->SetLabelSize(0.10); // was 0.15
165  histPtr->GetYaxis()->SetLabelOffset(0.002); // was 0.00
166  histPtr->GetYaxis()->SetTitleSize(0.10); // was 0.15
167  histPtr->GetYaxis()->SetTitleOffset(0.16); // was 0.80
168  }
169  }
std::unordered_map< std::string, std::unique_ptr< TH1F > > fRecoHistMap
cout<< "Opened file "<< fin<< " ixs= "<< ixs<< endl;if(ixs==0) hhh=(TH1F *) fff-> Get("h1")
Definition: AddMC.C:8
std::vector< double > fRecoQHigh
high edge of ADC values for drawing raw digits
std::vector< double > fRecoQLow
low edge of ADC values for drawing raw digits
std::vector< art::InputTag > fWireLabels
module labels that produced wires
enum geo::_plane_sigtype SigType_t
Enumerate the possible plane projections.
void evdb_tool::DrawWireHist::configure ( const fhicl::ParameterSet pset)
overridevirtual

Implements evdb_tool::IWaveformDrawer.

Definition at line 52 of file DrawWireHist_tool.cc.

References fColorMap, and fRecoHistMap.

Referenced by DrawWireHist().

53  {
54  fColorMap.push_back(kBlue);
55  fColorMap.push_back(kMagenta);
56  fColorMap.push_back(kBlack);
57  fColorMap.push_back(kRed);
58 
59  fRecoHistMap.clear();
60  }
std::unordered_map< std::string, std::unique_ptr< TH1F > > fRecoHistMap
std::vector< int > fColorMap
void evdb_tool::DrawWireHist::Draw ( const std::string &  options,
float  maxLowVal,
float  maxHiVal 
)
overridevirtual

Implements evdb_tool::IWaveformDrawer.

Definition at line 116 of file DrawWireHist_tool.cc.

References fRecoHistMap.

117  {
118  for (const auto& histMap : fRecoHistMap) {
119  TH1F* histPtr = histMap.second.get();
120 
121  // Set the limits
122  histPtr->SetMaximum(maxHiVal);
123  histPtr->SetMinimum(maxLowVal);
124 
125  histPtr->Draw(options.c_str());
126  }
127  }
std::unordered_map< std::string, std::unique_ptr< TH1F > > fRecoHistMap
void evdb_tool::DrawWireHist::Fill ( evdb::View2D view2D,
raw::ChannelID_t channel,
float  lowBin,
float  numTicks 
)
overridevirtual

Implements evdb_tool::IWaveformDrawer.

Definition at line 62 of file DrawWireHist_tool.cc.

References BookHistogram(), recob::Wire::Channel(), art::InputTag::encode(), fColorMap, evd::RawDrawingOptions::fDrawRawDataOrCalibWires, fMaximum, fMinimum, fRecoHistMap, evd::RecoDrawingOptions::fWireLabels, evdb::EventHolder::GetEvent(), evdb::EventHolder::Instance(), and recob::Wire::Signal().

66  {
69 
70  // Check if we're supposed to draw raw hits at all
71  if (rawOpt->fDrawRawDataOrCalibWires == 0) return;
72 
73  //grab the singleton with the event
75  if (!event) return;
76 
77  // Handle histograms
78  BookHistogram(channel, lowBin, numTicks);
79 
80  fMinimum = std::numeric_limits<float>::max();
81  fMaximum = std::numeric_limits<float>::lowest();
82 
83  int nWireLabels = 0;
84  for (size_t imod = 0; imod < recoOpt->fWireLabels.size(); ++imod) {
85  // Step one is to recover the hits for this label that match the input channel
86  art::InputTag const which = recoOpt->fWireLabels[imod];
87 
89  if (!event->getByLabel(which, wireVecHandle)) continue;
90  ++nWireLabels;
91 
92  for (size_t wireIdx = 0; wireIdx < wireVecHandle->size(); wireIdx++) {
93  art::Ptr<recob::Wire> wire(wireVecHandle, wireIdx);
94 
95  if (wire->Channel() != channel) continue;
96 
97  const std::vector<float>& signalVec = wire->Signal();
98 
99  TH1F* histPtr = fRecoHistMap.at(which.encode()).get();
100 
101  for (size_t idx = 0; idx < signalVec.size(); idx++) {
102  histPtr->Fill(float(idx) + 0.5, signalVec[idx]);
103 
104  fMinimum = std::min(fMinimum, signalVec[idx]);
105  fMaximum = std::max(fMaximum, signalVec[idx]);
106  }
107 
108  histPtr->SetLineColor(fColorMap.at((nWireLabels - 1) % recoOpt->fWireLabels.size()));
109 
110  // There is only one channel displayed so if here we are done
111  break;
112  }
113  } //end loop over HitFinding modules
114  }
std::unordered_map< std::string, std::unique_ptr< TH1F > > fRecoHistMap
const art::Event * GetEvent() const
Definition: EventHolder.cxx:45
int fDrawRawDataOrCalibWires
0 for raw
void BookHistogram(raw::ChannelID_t &, float, float)
std::string encode() const
Definition: InputTag.cc:97
std::vector< art::InputTag > fWireLabels
module labels that produced wires
static EventHolder * Instance()
Definition: EventHolder.cxx:15
std::vector< int > fColorMap
Definition: fwd.h:26
Event finding and building.
float evdb_tool::DrawWireHist::getMaximum ( ) const
inlineoverridevirtual

Implements evdb_tool::IWaveformDrawer.

Definition at line 32 of file DrawWireHist_tool.cc.

References fMaximum.

32 { return fMaximum; };
float evdb_tool::DrawWireHist::getMinimum ( ) const
inlineoverridevirtual

Implements evdb_tool::IWaveformDrawer.

Definition at line 33 of file DrawWireHist_tool.cc.

References BookHistogram(), and fMinimum.

33 { return fMinimum; };

Member Data Documentation

std::vector<int> evdb_tool::DrawWireHist::fColorMap
private

Definition at line 41 of file DrawWireHist_tool.cc.

Referenced by configure(), and Fill().

float evdb_tool::DrawWireHist::fMaximum
private

Definition at line 38 of file DrawWireHist_tool.cc.

Referenced by Fill(), and getMaximum().

float evdb_tool::DrawWireHist::fMinimum
private

Definition at line 39 of file DrawWireHist_tool.cc.

Referenced by Fill(), and getMinimum().

std::unordered_map<std::string, std::unique_ptr<TH1F> > evdb_tool::DrawWireHist::fRecoHistMap
private

Definition at line 42 of file DrawWireHist_tool.cc.

Referenced by BookHistogram(), configure(), Draw(), and Fill().


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