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

Public Member Functions

 DrawRawHist (const fhicl::ParameterSet &pset)
 
 ~DrawRawHist ()
 
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::unique_ptr< TH1F > fRawDigitHist
 

Detailed Description

Definition at line 27 of file DrawRawHist_tool.cc.

Constructor & Destructor Documentation

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

Definition at line 51 of file DrawRawHist_tool.cc.

References configure().

52  {
53  configure(pset);
54  }
void configure(const fhicl::ParameterSet &pset) override
evdb_tool::DrawRawHist::~DrawRawHist ( )

Definition at line 56 of file DrawRawHist_tool.cc.

56 {}

Member Function Documentation

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

Definition at line 153 of file DrawRawHist_tool.cc.

References DEFINE_ART_CLASS_TOOL, fRawDigitHist, evd::ColorDrawingOptions::fRawQHigh, evd::ColorDrawingOptions::fRawQLow, and geo::GeometryCore::SignalType().

Referenced by Fill(), and getMinimum().

154  {
157 
158  // Get rid of the previous histograms
159  if (fRawDigitHist.get()) fRawDigitHist.reset();
160 
161  // figure out the signal type for this plane, assume that
162  // plane n in each TPC/cryostat has the same type
163  geo::SigType_t sigType = geo->SignalType(channel);
164  int numBins = numTicks;
165 
166  fRawDigitHist = std::make_unique<TH1F>(
167  "fRAWQHisto", ";t [ticks];q [ADC]", numBins, startTick, startTick + numTicks);
168 
169  TH1F* histPtr = fRawDigitHist.get();
170 
171  histPtr->SetMaximum(cst->fRawQHigh[(size_t)sigType]);
172  histPtr->SetMinimum(cst->fRawQLow[(size_t)sigType]);
173 
174  histPtr->SetLineColor(kBlack);
175  histPtr->SetLineWidth(1);
176 
177  histPtr->GetXaxis()->SetLabelSize(0.10); // was 0.15
178  histPtr->GetXaxis()->SetLabelOffset(0.01); // was 0.00
179  histPtr->GetXaxis()->SetTitleSize(0.10); // was 0.15
180  histPtr->GetXaxis()->SetTitleOffset(0.60); // was 0.80
181 
182  histPtr->GetYaxis()->SetLabelSize(0.10); // was 0.15
183  histPtr->GetYaxis()->SetLabelOffset(0.002); // was 0.00
184  histPtr->GetYaxis()->SetTitleSize(0.10); // was 0.15
185  histPtr->GetYaxis()->SetTitleOffset(0.16); // was 0.80
186  }
std::vector< double > fRawQLow
low edge of ADC values for drawing raw digits
enum geo::_plane_sigtype SigType_t
Enumerate the possible plane projections.
std::vector< double > fRawQHigh
high edge of ADC values for drawing raw digits
SigType_t SignalType(PlaneID const &pid) const
Returns the type of signal on the channels of specified TPC plane.
Namespace collecting geometry-related classes utilities.
std::unique_ptr< TH1F > fRawDigitHist
void evdb_tool::DrawRawHist::configure ( const fhicl::ParameterSet pset)
overridevirtual

Implements evdb_tool::IWaveformDrawer.

Definition at line 58 of file DrawRawHist_tool.cc.

Referenced by DrawRawHist().

59  {
60  return;
61  }
void evdb_tool::DrawRawHist::Draw ( const std::string &  options,
float  maxLowVal,
float  maxHiVal 
)
overridevirtual

Implements evdb_tool::IWaveformDrawer.

Definition at line 139 of file DrawRawHist_tool.cc.

References fRawDigitHist.

140  {
141  TH1F* histPtr = fRawDigitHist.get();
142 
143  // Do we have valid limits to set?
144  histPtr->SetMaximum(maxHiVal);
145  histPtr->SetMinimum(maxLowVal);
146 
147  histPtr->Draw(options.c_str());
148 
149  return;
150  }
std::unique_ptr< TH1F > fRawDigitHist
void evdb_tool::DrawRawHist::Fill ( evdb::View2D view2D,
raw::ChannelID_t channel,
float  lowBin,
float  numTicks 
)
overridevirtual

Implements evdb_tool::IWaveformDrawer.

Definition at line 63 of file DrawRawHist_tool.cc.

References raw::RawDigit::ADCs(), BookHistogram(), raw::RawDigit::Channel(), raw::RawDigit::Compression(), fMaximum, fMinimum, evd::RawDrawingOptions::fPedestalOption, evd::RawDrawingOptions::fRawDataLabels, fRawDigitHist, evdb::EventHolder::GetEvent(), raw::RawDigit::GetPedestal(), evdb::EventHolder::Instance(), art::Handle< T >::isValid(), raw::RawDigit::Samples(), and raw::Uncompress().

67  {
69 
70  //grab the singleton with the event
72  if (!event) return;
73 
74  // Handle histograms
75  BookHistogram(channel, lowBin, numTicks);
76 
77  fMinimum = std::numeric_limits<float>::max();
78  fMaximum = std::numeric_limits<float>::lowest();
79 
80  // Loop over the possible producers of RawDigits
81  for (const auto& rawDataLabel : rawOpt->fRawDataLabels) {
82  art::Handle<std::vector<raw::RawDigit>> rawDigitVecHandle;
83  event->getByLabel(rawDataLabel, rawDigitVecHandle);
84 
85  if (!rawDigitVecHandle.isValid()) continue;
86 
87  for (size_t rawDigitIdx = 0; rawDigitIdx < rawDigitVecHandle->size(); rawDigitIdx++) {
88  art::Ptr<raw::RawDigit> rawDigit(rawDigitVecHandle, rawDigitIdx);
89 
90  if (rawDigit->Channel() != channel) continue;
91 
92  // We will need the pedestal service...
93  const lariov::DetPedestalProvider& pedestalRetrievalAlg =
95 
96  // recover the pedestal
97  float pedestal = 0;
98 
99  if (rawOpt->fPedestalOption == 0) { pedestal = pedestalRetrievalAlg.PedMean(channel); }
100  else if (rawOpt->fPedestalOption == 1) {
101  pedestal = rawDigit->GetPedestal();
102  }
103  else if (rawOpt->fPedestalOption == 2) {
104  pedestal = 0;
105  }
106  else {
107  mf::LogWarning("DrawRawHist")
108  << " PedestalOption is not understood: " << rawOpt->fPedestalOption
109  << ". Pedestals not subtracted.";
110  }
111 
112  std::vector<short> uncompressed(rawDigit->Samples());
113  raw::Uncompress(rawDigit->ADCs(), uncompressed, rawDigit->Compression());
114 
115  TH1F* histPtr = fRawDigitHist.get();
116 
117  for (size_t idx = 0; idx < uncompressed.size(); idx++) {
118  float signalVal = float(uncompressed[idx]) - pedestal;
119 
120  histPtr->Fill(float(idx) + 0.5, signalVal);
121  }
122 
123  short minimumVal = *std::min_element(uncompressed.begin(), uncompressed.end());
124  short maximumVal = *std::max_element(uncompressed.begin(), uncompressed.end());
125 
126  fMinimum = float(minimumVal) - pedestal;
127  fMaximum = float(maximumVal) - pedestal;
128 
129  histPtr->SetLineColor(kBlack);
130 
131  // There is only one channel displayed so if here we are done
132  break;
133  }
134  }
135 
136  return;
137  }
const art::Event * GetEvent() const
Definition: EventHolder.cxx:45
bool isValid() const noexcept
Definition: Handle.h:203
static EventHolder * Instance()
Definition: EventHolder.cxx:15
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
std::vector< art::InputTag > fRawDataLabels
module label that made the raw digits, default is daq
void Uncompress(const std::vector< short > &adc, std::vector< short > &uncompressed, raw::Compress_t compress)
Uncompresses a raw data buffer.
Definition: raw.cxx:744
int fPedestalOption
0: use DetPedestalService; 1: Use pedestal in raw::RawDigt; 2: no ped subtraction ...
std::unique_ptr< TH1F > fRawDigitHist
void BookHistogram(raw::ChannelID_t &, float, float)
Event finding and building.
float evdb_tool::DrawRawHist::getMaximum ( ) const
inlineoverridevirtual

Implements evdb_tool::IWaveformDrawer.

Definition at line 37 of file DrawRawHist_tool.cc.

References fMaximum.

37 { return fMaximum; };
float evdb_tool::DrawRawHist::getMinimum ( ) const
inlineoverridevirtual

Implements evdb_tool::IWaveformDrawer.

Definition at line 38 of file DrawRawHist_tool.cc.

References BookHistogram(), and fMinimum.

38 { return fMinimum; };

Member Data Documentation

float evdb_tool::DrawRawHist::fMaximum
private

Definition at line 43 of file DrawRawHist_tool.cc.

Referenced by Fill(), and getMaximum().

float evdb_tool::DrawRawHist::fMinimum
private

Definition at line 44 of file DrawRawHist_tool.cc.

Referenced by Fill(), and getMinimum().

std::unique_ptr<TH1F> evdb_tool::DrawRawHist::fRawDigitHist
private

Definition at line 46 of file DrawRawHist_tool.cc.

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


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