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

Public Member Functions

 FilterAnalyzerRunList (fhicl::ParameterSet const &pset)
 
 FilterAnalyzerRunList (FilterAnalyzerRunList const &)=delete
 
 FilterAnalyzerRunList (FilterAnalyzerRunList &&)=delete
 
FilterAnalyzerRunListoperator= (FilterAnalyzerRunList const &)=delete
 
FilterAnalyzerRunListoperator= (FilterAnalyzerRunList &&)=delete
 
void analyze (art::Event const &evt) override
 
void reconfigure (const fhicl::ParameterSet &pset)
 
- 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 ()
 
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 &)
 
std::string const & processName () const
 
bool wantAllEvents () const
 
bool wantEvent (Event const &e)
 
fhicl::ParameterSetID selectorConfig () const
 
Handle< TriggerResultsgetTriggerResults (Event const &e) 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

std::string fDigitModuleLabel
 
std::string fDigitModuleInstance
 
TF1 * fColFilterFunc
 Parameterized collection filter function. More...
 
TF1 * fIndUFilterFunc
 Parameterized induction filter function. More...
 
TF1 * fIndVFilterFunc
 Parameterized induction filter function. More...
 
const lariov::DetPedestalProviderfPedestalRetrievalAlg = *(lar::providerFrom<lariov::DetPedestalService>())
 
int ThisRunSeq = 0
 
std::vector< int > MyUsefulChans
 
TH2F * RawFFT_Planes
 
TH2F * FixFFT_Planes
 
TH2F * RawFFT_APAs
 
TH2F * FixFFT_APAs
 
TH2F * FFTRaw [36]
 
TH2F * FFTFix [36]
 
int NumRuns = 119
 

Additional Inherited Members

- Public Types inherited from art::EDAnalyzer
using WorkerType = WorkerT< EDAnalyzer >
 
using ModuleType = EDAnalyzer
 
- Protected Member Functions inherited from art::Observer
 Observer (fhicl::ParameterSet const &config)
 
 Observer (std::vector< std::string > const &paths, fhicl::ParameterSet const &config)
 
detail::ProcessAndEventSelectorsprocessAndEventSelectors ()
 
- 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 63 of file FilterAnalyzerRunList_module.cc.

Constructor & Destructor Documentation

DAQToOffline::FilterAnalyzerRunList::FilterAnalyzerRunList ( fhicl::ParameterSet const &  pset)
explicit

Definition at line 104 of file FilterAnalyzerRunList_module.cc.

104  : art::EDAnalyzer(pset) {
105 
106  this->reconfigure(pset);
107  gStyle->SetOptStat(0);
108 
109  MyUsefulChans.push_back(10); MyUsefulChans.push_back(70); MyUsefulChans.push_back(110);
110  MyUsefulChans.push_back(160); MyUsefulChans.push_back(205); MyUsefulChans.push_back(260);
111  MyUsefulChans.push_back(300); MyUsefulChans.push_back(410); MyUsefulChans.push_back(500);
112 
113  MyUsefulChans.push_back(525); MyUsefulChans.push_back(585); MyUsefulChans.push_back(625);
114  MyUsefulChans.push_back(675); MyUsefulChans.push_back(720); MyUsefulChans.push_back(775);
115  MyUsefulChans.push_back(815); MyUsefulChans.push_back(915); MyUsefulChans.push_back(1015);
116 
117  MyUsefulChans.push_back(1040); MyUsefulChans.push_back(1095); MyUsefulChans.push_back(1135);
118  MyUsefulChans.push_back(1185); MyUsefulChans.push_back(1230); MyUsefulChans.push_back(1285);
119  MyUsefulChans.push_back(1325); MyUsefulChans.push_back(1435); MyUsefulChans.push_back(1525);
120 
121  MyUsefulChans.push_back(1545); MyUsefulChans.push_back(1605); MyUsefulChans.push_back(1645);
122  MyUsefulChans.push_back(1695); MyUsefulChans.push_back(1740); MyUsefulChans.push_back(1795);
123  MyUsefulChans.push_back(1835); MyUsefulChans.push_back(1945); MyUsefulChans.push_back(2035);
124 
125  std::cout << "I have made my useful channel vector. It has size " << MyUsefulChans.size() << ". It's elements are as follows." << std::endl;
126  for (size_t xx=0; xx<MyUsefulChans.size(); ++xx) {
127  std::cout << "MyUsefulChans["<<xx<<"] " << MyUsefulChans[xx] << std::endl;
128  }
129 
131 
132  RawFFT_Planes = tfs->make<TH2F>("RawFFT_Planes", "Raw FFT for all channels grouped by plane; Channel Number; Frequency (KHz)", 36, 0, 36, 1000, 0, 1000 );
133  FixFFT_Planes = tfs->make<TH2F>("FixFFT_Planes", "Raw FFT for all channels grouped by plane; Channel Number; Frequency (KHz)", 36, 0, 36, 1000, 0, 1000 );
134 
135  RawFFT_APAs = tfs->make<TH2F>("RawFFT_APAs", "Raw FFT for all channels grouped by plane; Channel Number; Frequency (KHz)", 36, 0, 36, 1000, 0, 1000 );
136  FixFFT_APAs = tfs->make<TH2F>("FixFFT_APAs", "Raw FFT for all channels grouped by plane; Channel Number; Frequency (KHz)", 36, 0, 36, 1000, 0, 1000 );
137 
138  for (int XBin=0; XBin<36; ++XBin) {
139  std::stringstream oss;
140  oss<<"Chan_"<<MyUsefulChans[XBin];
141  int WhAPA = XBin/9;
142  int WhPla = (XBin/3)%3;
143  int PlaIn = XBin%3;
144  RawFFT_Planes->GetXaxis()->SetBinLabel((WhAPA*3)+(WhPla*12)+PlaIn+1, oss.str().c_str());
145  FixFFT_Planes->GetXaxis()->SetBinLabel((WhAPA*3)+(WhPla*12)+PlaIn+1, oss.str().c_str());
146 
147  RawFFT_APAs->GetXaxis()->SetBinLabel(XBin+1, oss.str().c_str());
148  FixFFT_APAs->GetXaxis()->SetBinLabel(XBin+1, oss.str().c_str());
149  }
150 
151  for (size_t HistoChan=0; HistoChan<MyUsefulChans.size(); HistoChan++) {
152  std::stringstream oss1a, oss1b;
153  oss1a << "RawFFT_"<<HistoChan;
154  oss1b << "Raw FFT for Channel "<<MyUsefulChans[HistoChan]<<"; Run Number; Frequency (KHz)";
155  FFTRaw[HistoChan] = tfs->make<TH2F>(oss1a.str().c_str(), oss1b.str().c_str(), NumRuns, 0, NumRuns, 1000, 0, 1000);
156 
157  std::stringstream oss2a, oss2b;
158  oss2a << "FixFFT_"<<HistoChan;
159  oss2b << "Fix FFT for Channel "<<MyUsefulChans[HistoChan]<<"; Run Number; Frequency (KHz)";
160  FFTFix[HistoChan] = tfs->make<TH2F>(oss2a.str().c_str(), oss2b.str().c_str(), NumRuns, 0, NumRuns, 1000, 0, 1000);
161  }
162 
163 
164 }
void reconfigure(const fhicl::ParameterSet &pset)
QTextStream & endl(QTextStream &s)
DAQToOffline::FilterAnalyzerRunList::FilterAnalyzerRunList ( FilterAnalyzerRunList const &  )
delete
DAQToOffline::FilterAnalyzerRunList::FilterAnalyzerRunList ( FilterAnalyzerRunList &&  )
delete

Member Function Documentation

void DAQToOffline::FilterAnalyzerRunList::analyze ( art::Event const &  evt)
overridevirtual

Implements art::EDAnalyzer.

Definition at line 190 of file FilterAnalyzerRunList_module.cc.

190  {
191  if (evt.event() != 1) return;
192 
193  ++ThisRunSeq;
194  std::cout << "\n\n\n\nThis Run sequence is " << ThisRunSeq << "\n\n\n" << std::endl;
195  for (size_t HistoChan=0; HistoChan<MyUsefulChans.size(); HistoChan++) {
196  std::stringstream oss;
197  oss<<"Run "<<(int)evt.run();
198  FFTRaw[HistoChan]->GetXaxis()->SetBinLabel(ThisRunSeq, oss.str().c_str());
199  FFTFix[HistoChan]->GetXaxis()->SetBinLabel(ThisRunSeq, oss.str().c_str());
200  }
201 
203 
205  evt.getByLabel(fDigitModuleLabel, fDigitModuleInstance, rawDigitHandle);
206  std::vector<raw::RawDigit> const& rawDigitVector(*rawDigitHandle);
207 
208 
209  std::vector< std::pair<int,int> > ZeroFreq;
210  //ZeroFreq.push_back( std::make_pair(276 , 285 ) );
211  //ZeroFreq.push_back( std::make_pair(558 , 568 ) );
212  //ZeroFreq.push_back( std::make_pair(837 , 849 ) );
213  //ZeroFreq.push_back( std::make_pair(1116, 1127) );
214  //ZeroFreq.push_back( std::make_pair(4340, 5205) );
215 
216  for (size_t DigLoop=0; DigLoop < rawDigitVector.size(); ++DigLoop) {
217  int Channel = rawDigitVector[DigLoop].Channel();
218  size_t NADC = rawDigitVector[DigLoop].NADC();
219  double Pedestal = rawDigitVector[DigLoop].GetPedestal();
220  const geo::View_t view = geo->View(Channel);
221  // Check if this channel is one of the ones I want.
222  for (size_t GotChan=0; GotChan<MyUsefulChans.size(); ++GotChan) {
223  if (MyUsefulChans[GotChan] != Channel) continue;
224 
225  int WhAPA = GotChan/9;
226  int WhPla = (GotChan/3)%3;
227  int PlaIn = GotChan%3;
228  //std::cout << "Looking at Chan " << Channel << ", GotChan " << GotChan << ". I am going to fill bin " << (WhAPA*3)+(WhPla*12)+PlaIn+0.5 << std::endl;
229 
230  //std::cout << "Looking at rawDigitVector["<<DigLoop<<"] it was on channel " << rawDigitVector[DigLoop].Channel() << "("<<Channel<<") it is in View " << view
231  // << ", NADC is " << rawDigitVector[DigLoop].NADC() << " ("<<NADC<<")"
232  // << ", pedestal is " << rawDigitVector[DigLoop].GetPedestal() << " ("<<Pedestal<<")"
233  // << std::endl;
234 
235  // Fill the RawDigit histogram for this histogram.
236  TH1F* hRawDigit = new TH1F("hRawDigit","",NADC,0,NADC/2);
237  TH1F* hRawFFT = new TH1F("hRawFFT" ,"",NADC,0,NADC);
238  for (size_t ADCs=0; ADCs < NADC; ++ADCs) {
239  hRawDigit -> SetBinContent( ADCs+1, rawDigitVector[DigLoop].ADC(ADCs)-Pedestal );
240  }
241  for (size_t ww=NADC; ww<NADC; ++ww)
242  hRawFFT -> SetBinContent( ww, 0 );
243  // Make the FFT for this channel.
244  hRawDigit -> FFT( hRawFFT ,"MAG");
245  for (size_t bin = 0; bin < NADC; ++bin) {
246  double BinVal = hRawFFT->GetBinContent(bin+1);
247  double freq = 2000. * bin / (double)NADC;
248  if (freq < 1000 && BinVal < 1e5) {
249  FFTRaw[GotChan] -> Fill( (ThisRunSeq-0.5) , freq, BinVal );
250  RawFFT_Planes -> Fill( (WhAPA*3)+(WhPla*12)+PlaIn+0.5, freq, BinVal );
251  RawFFT_APAs -> Fill( GotChan+0.5, freq, BinVal );
252  }
253  }
254 
255  // I want to do an inverse FFT, so need to convert the tranformed FFT into an array....
256  //double Re[NADC], Im[NADC];
257  std::unique_ptr<double[]> Re( new double[NADC]);
258  std::unique_ptr<double[]> Im( new double[NADC]);
259  TVirtualFFT *fft = TVirtualFFT::GetCurrentTransform();
260  fft->GetPointsComplex(Re.get(),Im.get());
261 
262  // Set the noisy frequency range bins to an average value.
263  for (size_t aa=0; aa<ZeroFreq.size(); ++aa) {
264  for (int bb=ZeroFreq[aa].first; bb<ZeroFreq[aa].second; ++bb) {
265  double ReMeanVal=0;
266  double ImMeanVal=0;
267  int Range = 50;
268  for (int cc=0; cc<Range; ++cc) {
269  ReMeanVal += Re[ZeroFreq[aa].first-cc] + Re[ZeroFreq[aa].second+cc];
270  ImMeanVal += Im[ZeroFreq[aa].first-cc] + Im[ZeroFreq[aa].second+cc];
271  }
272  ReMeanVal = ReMeanVal / Range;
273  Re[bb] = Re[1500-bb] = ReMeanVal;
274  ImMeanVal = ImMeanVal / Range;
275  Im[bb] = Im[1500-bb] = ImMeanVal;
276  }
277  }
278  Re[0] = Re[1];
279  Im[0] = Im[1];
280 
281  // Apply the filter...
282  for (size_t bin = 0; bin < NADC; ++bin) {
283  double freq = 2000. * bin / NADC;
284  if (view == geo::kU) { // U plane
285  Re[bin] = Re[bin]*fIndUFilterFunc->Eval(freq);
286  Im[bin] = Im[bin]*fIndUFilterFunc->Eval(freq);
287  } else if ( view == geo::kV) { // V plane
288  Re[bin] = Re[bin]*fIndVFilterFunc->Eval(freq);
289  Im[bin] = Im[bin]*fIndVFilterFunc->Eval(freq);
290  } else if ( view == geo::kZ) { // Collection plane
291  Re[bin] = Re[bin]*fColFilterFunc->Eval(freq);
292  Im[bin] = Im[bin]*fColFilterFunc->Eval(freq);
293  }
294  double MagVal = pow ( Re[bin]*Re[bin] + Im[bin]*Im[bin], 0.5);
295  if (TMath::IsNaN(MagVal)) MagVal = 0;
296 
297  // Now do the big histograms...
298  if (freq < 1000 && MagVal < 1e5) {
299  FFTFix[GotChan] -> Fill( (ThisRunSeq-0.5) , freq, MagVal );
300  FixFFT_Planes -> Fill( (WhAPA*3)+(WhPla*12)+PlaIn+0.5, freq, MagVal );
301  FixFFT_APAs -> Fill( GotChan+0.5, freq, MagVal );
302  }
303  }
304  } // Checking if one of my chosen useful channels
305  } // RawDigit Loop
306  return;
307 }
unsigned int event
Definition: DataStructs.h:574
unsigned int run
Definition: DataStructs.h:575
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
Planes which measure V.
Definition: geo_types.h:126
constexpr T pow(T x)
Definition: pow.h:72
Planes which measure Z direction.
Definition: geo_types.h:128
TF1 * fIndUFilterFunc
Parameterized induction filter function.
Planes which measure U.
Definition: geo_types.h:125
TF1 * fIndVFilterFunc
Parameterized induction filter function.
TF1 * fColFilterFunc
Parameterized collection filter function.
View_t View(geo::PlaneID const &pid) const
Returns the view (wire orientation) on the channels of specified TPC plane.
ChannelMappingService::Channel Channel
QTextStream & bin(QTextStream &s)
TCEvent evt
Definition: DataStructs.cxx:7
Namespace collecting geometry-related classes utilities.
QTextStream & endl(QTextStream &s)
FilterAnalyzerRunList& DAQToOffline::FilterAnalyzerRunList::operator= ( FilterAnalyzerRunList const &  )
delete
FilterAnalyzerRunList& DAQToOffline::FilterAnalyzerRunList::operator= ( FilterAnalyzerRunList &&  )
delete
void DAQToOffline::FilterAnalyzerRunList::reconfigure ( const fhicl::ParameterSet pset)

Definition at line 166 of file FilterAnalyzerRunList_module.cc.

166  {
167  fDigitModuleLabel = pset.get<std::string>("DigitModuleLabel");
168  fDigitModuleInstance = pset.get<std::string>("DigitModuleInstance");
169 
170  // Make the filter functions.
171  std::string colFilt = pset.get<std::string>("ColFilter");
172  std::vector<double> colFiltParams = pset.get<std::vector<double> >("ColFilterParams");
173  fColFilterFunc = new TF1("colFilter", colFilt.c_str());
174  for(unsigned int i=0; i<colFiltParams.size(); ++i)
175  fColFilterFunc->SetParameter(i, colFiltParams[i]);
176 
177  std::string indUFilt = pset.get<std::string>("IndUFilter");
178  std::vector<double> indUFiltParams = pset.get<std::vector<double> >("IndUFilterParams");
179  fIndUFilterFunc = new TF1("indUFilter", indUFilt.c_str());
180  for(unsigned int i=0; i<indUFiltParams.size(); ++i)
181  fIndUFilterFunc->SetParameter(i, indUFiltParams[i]);
182 
183  std::string indVFilt = pset.get<std::string>("IndVFilter");
184  std::vector<double> indVFiltParams = pset.get<std::vector<double> >("IndVFilterParams");
185  fIndVFilterFunc = new TF1("indVFilter", indVFilt.c_str());
186  for(unsigned int i=0; i<indVFiltParams.size(); ++i)
187  fIndVFilterFunc->SetParameter(i, indVFiltParams[i]);
188 }
std::string string
Definition: nybbler.cc:12
TF1 * fIndUFilterFunc
Parameterized induction filter function.
T get(std::string const &key) const
Definition: ParameterSet.h:231
TF1 * fIndVFilterFunc
Parameterized induction filter function.
TF1 * fColFilterFunc
Parameterized collection filter function.

Member Data Documentation

TF1* DAQToOffline::FilterAnalyzerRunList::fColFilterFunc
private

Parameterized collection filter function.

Definition at line 84 of file FilterAnalyzerRunList_module.cc.

std::string DAQToOffline::FilterAnalyzerRunList::fDigitModuleInstance
private

Definition at line 82 of file FilterAnalyzerRunList_module.cc.

std::string DAQToOffline::FilterAnalyzerRunList::fDigitModuleLabel
private

Definition at line 81 of file FilterAnalyzerRunList_module.cc.

TH2F* DAQToOffline::FilterAnalyzerRunList::FFTFix[36]
private

Definition at line 99 of file FilterAnalyzerRunList_module.cc.

TH2F* DAQToOffline::FilterAnalyzerRunList::FFTRaw[36]
private

Definition at line 98 of file FilterAnalyzerRunList_module.cc.

TF1* DAQToOffline::FilterAnalyzerRunList::fIndUFilterFunc
private

Parameterized induction filter function.

Definition at line 85 of file FilterAnalyzerRunList_module.cc.

TF1* DAQToOffline::FilterAnalyzerRunList::fIndVFilterFunc
private

Parameterized induction filter function.

Definition at line 86 of file FilterAnalyzerRunList_module.cc.

TH2F* DAQToOffline::FilterAnalyzerRunList::FixFFT_APAs
private

Definition at line 96 of file FilterAnalyzerRunList_module.cc.

TH2F* DAQToOffline::FilterAnalyzerRunList::FixFFT_Planes
private

Definition at line 94 of file FilterAnalyzerRunList_module.cc.

const lariov::DetPedestalProvider& DAQToOffline::FilterAnalyzerRunList::fPedestalRetrievalAlg = *(lar::providerFrom<lariov::DetPedestalService>())
private

Definition at line 88 of file FilterAnalyzerRunList_module.cc.

std::vector<int> DAQToOffline::FilterAnalyzerRunList::MyUsefulChans
private

Definition at line 91 of file FilterAnalyzerRunList_module.cc.

int DAQToOffline::FilterAnalyzerRunList::NumRuns = 119
private

Definition at line 101 of file FilterAnalyzerRunList_module.cc.

TH2F* DAQToOffline::FilterAnalyzerRunList::RawFFT_APAs
private

Definition at line 95 of file FilterAnalyzerRunList_module.cc.

TH2F* DAQToOffline::FilterAnalyzerRunList::RawFFT_Planes
private

Definition at line 93 of file FilterAnalyzerRunList_module.cc.

int DAQToOffline::FilterAnalyzerRunList::ThisRunSeq = 0
private

Definition at line 90 of file FilterAnalyzerRunList_module.cc.


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