LArSoft  v10_04_05
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)
 
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 49 of file DrawRawHist_tool.cc.

References configure().

50  {
51  configure(pset);
52  }
void configure(const fhicl::ParameterSet &pset) override

Member Function Documentation

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

Definition at line 142 of file DrawRawHist_tool.cc.

References DEFINE_ART_CLASS_TOOL, fRawDigitHist, evd::ColorDrawingOptions::fRawQHigh, evd::ColorDrawingOptions::fRawQLow, and Get.

Referenced by Fill(), and getMinimum().

143  {
145  auto const& wireReadoutGeom = art::ServiceHandle<geo::WireReadout const>()->Get();
146 
147  // Get rid of the previous histograms
148  if (fRawDigitHist.get()) fRawDigitHist.reset();
149 
150  // figure out the signal type for this plane, assume that
151  // plane n in each TPC/cryostat has the same type
152  geo::SigType_t sigType = wireReadoutGeom.SignalType(channel);
153  int numBins = numTicks;
154 
155  fRawDigitHist = std::make_unique<TH1F>(
156  "fRAWQHisto", ";t [ticks];q [ADC]", numBins, startTick, startTick + numTicks);
157 
158  TH1F* histPtr = fRawDigitHist.get();
159 
160  histPtr->SetMaximum(cst->fRawQHigh[(size_t)sigType]);
161  histPtr->SetMinimum(cst->fRawQLow[(size_t)sigType]);
162 
163  histPtr->SetLineColor(kBlack);
164  histPtr->SetLineWidth(1);
165 
166  histPtr->GetXaxis()->SetLabelSize(0.10); // was 0.15
167  histPtr->GetXaxis()->SetLabelOffset(0.01); // was 0.00
168  histPtr->GetXaxis()->SetTitleSize(0.10); // was 0.15
169  histPtr->GetXaxis()->SetTitleOffset(0.60); // was 0.80
170 
171  histPtr->GetYaxis()->SetLabelSize(0.10); // was 0.15
172  histPtr->GetYaxis()->SetLabelOffset(0.002); // was 0.00
173  histPtr->GetYaxis()->SetTitleSize(0.10); // was 0.15
174  histPtr->GetYaxis()->SetTitleOffset(0.16); // was 0.80
175  }
std::vector< double > fRawQLow
low edge of ADC values for drawing raw digits
cout<< "Opened file "<< fin<< " ixs= "<< ixs<< endl;if(ixs==0) hhh=(TH1F *) fff-> Get("h1")
Definition: AddMC.C:8
enum geo::_plane_sigtype SigType_t
Enumerate the possible plane projections.
std::vector< double > fRawQHigh
high edge of ADC values for drawing raw digits
std::unique_ptr< TH1F > fRawDigitHist
void evdb_tool::DrawRawHist::configure ( const fhicl::ParameterSet pset)
overridevirtual

Implements evdb_tool::IWaveformDrawer.

Definition at line 54 of file DrawRawHist_tool.cc.

Referenced by DrawRawHist().

54 {}
void evdb_tool::DrawRawHist::Draw ( const std::string &  options,
float  maxLowVal,
float  maxHiVal 
)
overridevirtual

Implements evdb_tool::IWaveformDrawer.

Definition at line 130 of file DrawRawHist_tool.cc.

References fRawDigitHist.

131  {
132  TH1F* histPtr = fRawDigitHist.get();
133 
134  // Do we have valid limits to set?
135  histPtr->SetMaximum(maxHiVal);
136  histPtr->SetMinimum(maxLowVal);
137 
138  histPtr->Draw(options.c_str());
139  }
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 56 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().

60  {
62 
63  //grab the singleton with the event
65  if (!event) return;
66 
67  // Handle histograms
68  BookHistogram(channel, lowBin, numTicks);
69 
70  fMinimum = std::numeric_limits<float>::max();
71  fMaximum = std::numeric_limits<float>::lowest();
72 
73  // Loop over the possible producers of RawDigits
74  for (const auto& rawDataLabel : rawOpt->fRawDataLabels) {
75  art::Handle<std::vector<raw::RawDigit>> rawDigitVecHandle;
76  event->getByLabel(rawDataLabel, rawDigitVecHandle);
77 
78  if (!rawDigitVecHandle.isValid()) continue;
79 
80  for (size_t rawDigitIdx = 0; rawDigitIdx < rawDigitVecHandle->size(); rawDigitIdx++) {
81  art::Ptr<raw::RawDigit> rawDigit(rawDigitVecHandle, rawDigitIdx);
82 
83  if (rawDigit->Channel() != channel) continue;
84 
85  // We will need the pedestal service...
86  const lariov::DetPedestalProvider& pedestalRetrievalAlg =
88 
89  // recover the pedestal
90  float pedestal = 0;
91 
92  if (rawOpt->fPedestalOption == 0) { pedestal = pedestalRetrievalAlg.PedMean(channel); }
93  else if (rawOpt->fPedestalOption == 1) {
94  pedestal = rawDigit->GetPedestal();
95  }
96  else if (rawOpt->fPedestalOption == 2) {
97  pedestal = 0;
98  }
99  else {
100  mf::LogWarning("DrawRawHist")
101  << " PedestalOption is not understood: " << rawOpt->fPedestalOption
102  << ". Pedestals not subtracted.";
103  }
104 
105  std::vector<short> uncompressed(rawDigit->Samples());
106  raw::Uncompress(rawDigit->ADCs(), uncompressed, rawDigit->Compression());
107 
108  TH1F* histPtr = fRawDigitHist.get();
109 
110  for (size_t idx = 0; idx < uncompressed.size(); idx++) {
111  float signalVal = float(uncompressed[idx]) - pedestal;
112 
113  histPtr->Fill(float(idx) + 0.5, signalVal);
114  }
115 
116  short minimumVal = *std::min_element(uncompressed.begin(), uncompressed.end());
117  short maximumVal = *std::max_element(uncompressed.begin(), uncompressed.end());
118 
119  fMinimum = float(minimumVal) - pedestal;
120  fMaximum = float(maximumVal) - pedestal;
121 
122  histPtr->SetLineColor(kBlack);
123 
124  // There is only one channel displayed so if here we are done
125  break;
126  }
127  }
128  }
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 35 of file DrawRawHist_tool.cc.

References fMaximum.

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

Implements evdb_tool::IWaveformDrawer.

Definition at line 36 of file DrawRawHist_tool.cc.

References BookHistogram(), and fMinimum.

36 { return fMinimum; };

Member Data Documentation

float evdb_tool::DrawRawHist::fMaximum
private

Definition at line 41 of file DrawRawHist_tool.cc.

Referenced by Fill(), and getMaximum().

float evdb_tool::DrawRawHist::fMinimum
private

Definition at line 42 of file DrawRawHist_tool.cc.

Referenced by Fill(), and getMinimum().

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

Definition at line 44 of file DrawRawHist_tool.cc.

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


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