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

Public Member Functions

 dunezsanalysis (fhicl::ParameterSet const &pset)
 
virtual ~dunezsanalysis ()
 
void beginJob ()
 
void beginRun (const art::Run &run)
 
void reconfigure (fhicl::ParameterSet const &pset)
 
void analyze (const art::Event &evt)
 
- 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

std::string fRawProducerLabel
 
TTree * fReconstructionNtuple
 
int fEvent
 
int fRun
 
double fChargeSum [200]
 
art::ServiceHandle< geo::GeometryfGeometry
 

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 67 of file dunezsanalysis_module.cc.

Constructor & Destructor Documentation

dunezsanalysis::dunezsanalysis::dunezsanalysis ( fhicl::ParameterSet const &  pset)
explicit

Definition at line 125 of file dunezsanalysis_module.cc.

125  : EDAnalyzer(parameterSet)
126  {
127  // Read in the parameters from the .fcl file.
128  this->reconfigure(parameterSet);
129  }
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
void reconfigure(fhicl::ParameterSet const &pset)
dunezsanalysis::dunezsanalysis::~dunezsanalysis ( )
virtual

Definition at line 133 of file dunezsanalysis_module.cc.

134  {}

Member Function Documentation

void dunezsanalysis::dunezsanalysis::analyze ( const art::Event evt)

Definition at line 181 of file dunezsanalysis_module.cc.

182  {
183  // Start by fetching some basic event information for our n-tuple.
184  fEvent = event.id().event();
185  fRun = event.run();
186 
187  auto rawDigitHandle = event.getHandle< std::vector<raw::RawDigit> >("daq");
188 
189  for (int i=0;i<200;i++)
190  {
191  fChargeSum[i] = 0;
192  }
193 
194  unsigned int nNeighbors = 0;
195 
196  // add up all the digits on the collection wires in the entire event, for each value of
197  // the zero-suppression threshold
198 
199  // put it in a more easily usable form
200  std::vector< art::Ptr<raw::RawDigit> > Digits;
201  art::fill_ptr_vector(Digits, rawDigitHandle);
202  //loop through all RawDigits (over entire channels)
203  for(size_t d = 0; d < Digits.size(); d++)
204  {
206  digit=Digits.at(d);
207  uint32_t chan = digit->Channel();
208  std::vector<short> uncompressed(digit->Samples());
209  if (fGeometry->View(chan) == geo::kZ) // for now only do charge sums for collection hits
210  {
211  raw::Uncompress(digit->ADCs(), uncompressed, digit->Compression());
212  for(unsigned int tick=0;tick<uncompressed.size();tick++)
213  {
214  //std::cout << "trjadc: " << fEvent << " " << chan << " " << uncompressed.at(tick) << std::endl;
215  unsigned int tlow = (tick < nNeighbors)? 0: tick - nNeighbors;
216  unsigned int thigh = tick + nNeighbors;
217  if (thigh>=uncompressed.size()) thigh = uncompressed.size()-1;
218  for (short zscut=0;zscut<200;zscut++)
219  {
220  if (uncompressed.at(tick) >= zscut)
221  {
222  fChargeSum[zscut] += uncompressed.at(tick);
223  }
224  else
225  {
226  for (unsigned int k=tlow;k<=thigh;k++)
227  {
228  if (uncompressed.at(k) >= zscut)
229  {
230  fChargeSum[zscut] += uncompressed.at(tick);
231  break;
232  }
233  }
234  }
235  }
236  }
237  }
238  }
239 
240 
241  fReconstructionNtuple->Fill();
242  return;
243  }
const ADCvector_t & ADCs() const
Reference to the compressed ADC count vector.
Definition: RawDigit.h:210
ULong64_t Samples() const
Number of samples in the uncompressed ADC data.
Definition: RawDigit.h:213
ChannelID_t Channel() const
DAQ channel this raw data was read from.
Definition: RawDigit.h:212
Planes which measure Z direction.
Definition: geo_types.h:132
tick_as<> tick
Tick number, represented by std::ptrdiff_t.
Definition: electronics.h:75
View_t View(geo::PlaneID const &pid) const
Returns the view (wire orientation) on the channels of specified TPC plane.
raw::Compress_t Compression() const
Compression algorithm used to store the ADC counts.
Definition: RawDigit.h:216
art::ServiceHandle< geo::Geometry > fGeometry
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:297
void Uncompress(const std::vector< short > &adc, std::vector< short > &uncompressed, raw::Compress_t compress)
Uncompresses a raw data buffer.
Definition: raw.cxx:776
void dunezsanalysis::dunezsanalysis::beginJob ( )
virtual

Reimplemented from art::EDAnalyzer.

Definition at line 137 of file dunezsanalysis_module.cc.

138  {
139 
140  // Access ART's TFileService, which will handle creating and writing
141  // histograms and n-tuples for us.
143 
144  // The arguments to 'make<whatever>' are the same as those passed
145  // to the 'whatever' constructor, provided 'whatever' is a ROOT
146  // class that TFileService recognizes.
147 
148  // Define our n-tuples, which are limited forms of ROOT
149  // TTrees. Start with the TTree itself.
150 
151  fReconstructionNtuple = tfs->make<TTree>("dunezsanalysisReconstruction","dunezsanalysisReconstruction");
152 
153  // Define the branches (columns) of our simulation n-tuple. When
154  // we write a variable, we give the address of the variable to
155  // TTree::Branch.
156  fReconstructionNtuple->Branch("Event", &fEvent, "Event/I");
157  fReconstructionNtuple->Branch("Run", &fRun, "Run/I");
158  // When we write arrays, we give the address of the array to
159  // TTree::Branch; in C++ this is simply the array name.
160  fReconstructionNtuple->Branch("ChargeSum", fChargeSum, "ChargeSum[200]/D");
161  }
void dunezsanalysis::dunezsanalysis::beginRun ( const art::Run run)

Definition at line 164 of file dunezsanalysis_module.cc.

165  {
166  // How to convert from number of electrons to GeV. The ultimate
167  // source of this conversion factor is
168  // ${SRT_PUBLIC_CONTEXT}/SimpleTypesAndConstants/PhysicalConstants.h.
170  //fElectronsToGeV = 1./larParameters->GeVToElectrons();
171  }
void dunezsanalysis::dunezsanalysis::reconfigure ( fhicl::ParameterSet const &  pset)

Definition at line 174 of file dunezsanalysis_module.cc.

175  {
176  // no parameters for now. Just read in raw data
177  return;
178  }

Member Data Documentation

double dunezsanalysis::dunezsanalysis::fChargeSum[200]
private

Definition at line 111 of file dunezsanalysis_module.cc.

int dunezsanalysis::dunezsanalysis::fEvent
private

Definition at line 106 of file dunezsanalysis_module.cc.

art::ServiceHandle<geo::Geometry> dunezsanalysis::dunezsanalysis::fGeometry
private

Definition at line 114 of file dunezsanalysis_module.cc.

std::string dunezsanalysis::dunezsanalysis::fRawProducerLabel
private

Definition at line 99 of file dunezsanalysis_module.cc.

TTree* dunezsanalysis::dunezsanalysis::fReconstructionNtuple
private

Definition at line 103 of file dunezsanalysis_module.cc.

int dunezsanalysis::dunezsanalysis::fRun
private

Definition at line 107 of file dunezsanalysis_module.cc.


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