dunezsanalysis_module.cc
Go to the documentation of this file.
1 // dunezsanalysis_module.cc
2 // A basic "skeleton" to read in art::Event records from a file,
3 // access their information, and do something with them.
4 
5 // See
6 // <https://cdcvs.fnal.gov/redmine/projects/larsoftsvn/wiki/Using_the_Framework>
7 // for a description of the ART classes used here.
8 
9 // Almost everything you see in the code below may have to be changed
10 // by you to suit your task. The example task is to make histograms
11 // and n-tuples related to dE/dx of particle tracks in the detector.
12 
13 // As you try to understand why things are done a certain way in this
14 // example ("What's all this stuff about 'auto const&'?"), it will help
15 // to read ADDITIONAL_NOTES.txt in the same directory as this file.
16 
17 #ifndef dunezsanalysis_Module
18 #define dunezsanalysis_Module
19 
20 // LArSoft includes
23 #include "lardataobj/RawData/raw.h"
31 #include "lardataobj/RawData/raw.h"
32 
33 // Framework includes
38 #include "art_root_io/TFileService.h"
40 #include "canvas/Persistency/Common/FindManyP.h"
42 #include "fhiclcpp/ParameterSet.h"
43 
44 // ROOT includes. Note: To look up the properties of the ROOT classes,
45 // use the ROOT web site; e.g.,
46 // <http://root.cern.ch/root/html532/ClassIndex.html>
47 #include "TH1.h"
48 #include "TH2.h"
49 #include "TTree.h"
50 #include "TLorentzVector.h"
51 #include "TVector3.h"
52 
53 // C++ Includes
54 #include <map>
55 #include <vector>
56 #include <algorithm>
57 #include <iostream>
58 #include <string>
59 #include <cmath>
60 
61 namespace dunezsanalysis {
62 
63  //-----------------------------------------------------------------------
64  //-----------------------------------------------------------------------
65  // class definition
66 
68  {
69  public:
70 
71  // Standard constructor and destructor for an ART module.
72  explicit dunezsanalysis(fhicl::ParameterSet const& pset);
73  virtual ~dunezsanalysis();
74 
75  // This method is called once, at the start of the job. In this
76  // example, it will define the histograms and n-tuples we'll write.
77  void beginJob();
78 
79  // This method is called once, at the start of each run. It's a
80  // good place to read databases or files that may have
81  // run-dependent information.
82  void beginRun(const art::Run& run);
83 
84  // This method reads in any parameters from the .fcl files. This
85  // method is called 'reconfigure' because it might be called in the
86  // middle of a job; e.g., if the user changes parameter values in an
87  // interactive event display.
88  void reconfigure(fhicl::ParameterSet const& pset);
89 
90  // The analysis routine, called once per event.
91  void analyze (const art::Event& evt);
92 
93  private:
94 
95  // The stuff below is the part you'll most likely have to change to
96  // go from this custom example to your own task.
97 
98  // The parameters we'll read from the .fcl file.
99  std::string fRawProducerLabel; // The name of the producer that tracked simulated particles through the detector
100 
101  // The n-tuples we'll create.
102  //TTree* fSimulationNtuple; // unused
104 
105  // The variables that will go into the n-tuple.
106  int fEvent;
107  int fRun;
108 
109  // array of sums of charges for different choices of zero-suppression thresholds.
110 
111  double fChargeSum[200];
112 
113  // Other variables that will be shared between different methods.
114  art::ServiceHandle<geo::Geometry> fGeometry; // pointer to Geometry service
115 
116  }; // class dunezsanalysis
117 
118 
119  //-----------------------------------------------------------------------
120  //-----------------------------------------------------------------------
121  // class implementation
122 
123  //-----------------------------------------------------------------------
124  // Constructor
125  dunezsanalysis::dunezsanalysis(fhicl::ParameterSet const& parameterSet) : EDAnalyzer(parameterSet)
126  {
127  // Read in the parameters from the .fcl file.
128  this->reconfigure(parameterSet);
129  }
130 
131  //-----------------------------------------------------------------------
132  // Destructor
133  dunezsanalysis::~dunezsanalysis()
134  {}
135 
136  //-----------------------------------------------------------------------
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  }
162 
163  //-----------------------------------------------------------------------
164  void dunezsanalysis::beginRun(const art::Run& /*run*/)
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  }
172 
173  //-----------------------------------------------------------------------
175  {
176  // no parameters for now. Just read in raw data
177  return;
178  }
179 
180  //-----------------------------------------------------------------------
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  }
244 
245  // This macro has to be defined for this module to be invoked from a
246  // .fcl file; see dunezsanalysis.fcl for more information.
248 
249  } // namespace dunezsanalysis
250 
251 #endif // dunezsanalysis_Module
def analyze(root, level, gtrees, gbranches, doprint)
Definition: rootstat.py:69
Store parameters for running LArG4.
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
std::string string
Definition: nybbler.cc:12
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
Particle class.
art framework interface to geometry description
Definition: Run.h:17
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:67
void beginJob()
Definition: Breakpoints.cc:14
virtual void reconfigure(fhicl::ParameterSet const &pset)
tick_as<> tick
Tick number, represented by std::ptrdiff_t.
Definition: electronics.h:75
Definition of data types for geometry description.
raw::Compress_t Compression() const
Compression algorithm used to store the ADC counts.
Definition: RawDigit.h:216
Declaration of signal hit object.
art::ServiceHandle< geo::Geometry > fGeometry
TCEvent evt
Definition: DataStructs.cxx:7
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
Event finding and building.