LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
GausHitFinderAna_module.cc
Go to the documentation of this file.
1 //
3 // Gaus(s)HitFinder class designed to analyze signal on a wire in the TPC
4 //
5 // jaasaadi@syr.edu
6 //
7 // Note: This is a rework of the original hit finder ana module
8 // The only histograms that are saved are ones that can be used
9 // to make sure the hit finder is functioning...the rest is
10 // outputted to a TTree for offline analysis.
12 
13 // LArSoft includes
21 
22 // ROOT includes
23 #include "TH1.h"
24 #include "TH2.h"
25 #include "TTree.h"
26 
27 // C++ includes
28 #include <string>
29 
30 // Framework includes
36 #include "art_root_io/TFileService.h"
38 #include "fhiclcpp/ParameterSet.h"
40 
41 constexpr int kMaxHits = 20000;
42 
43 namespace hit {
44 
47  public:
48  explicit GausHitFinderAna(fhicl::ParameterSet const& pset);
49 
50  private:
51  void analyze(const art::Event& evt) override;
52  void beginJob() override;
53 
54  std::string fHitFinderModuleLabel;
55  std::string fLArG4ModuleLabel;
56  std::string fCalDataModuleLabel;
57 
62 
63  // ### TTree for offline analysis ###
64  TTree* fHTree;
65 
66  // === Event Information ===
67  Int_t fRun; // Run Number
68  Int_t fEvt; // Event Number
69 
70  // === Wire Information ====
71  Float_t fWireTotalCharge; // Charge on all wires
72 
73  // === Hit Information ===
74  Int_t fnHits; // Number of Hits in the Event
75  Int_t fWire[kMaxHits]; // Wire Number
76  Float_t fStartTime[kMaxHits]; // Start Time
77  Float_t fEndTime[kMaxHits]; // End Time
78  Float_t fPeakTime[kMaxHits]; // Peak Time
79  Float_t fPeakTimeUncert[kMaxHits]; // Peak Time Uncertainty
80  Float_t fCharge[kMaxHits]; // Charge of this hit
81  Float_t fChargeUncert[kMaxHits]; // Charge Uncertainty of this hit
82  Int_t fMultiplicity[kMaxHits]; // Hit pulse multiplicity
83  Float_t fGOF[kMaxHits]; // Goodness of Fit (Chi2/NDF)
84 
85  // === Total Hit Information ===
86  Float_t fTotalHitChargePerEvent; //Total charge recorded in each event
87 
88  // === Truth Hit Info from BackTrackerService ===
89  Float_t fTruePeakPos[kMaxHits]; // Truth Time Tick info from BackTrackerService
90 
91  }; // class GausHitFinderAna
92 
93  //-------------------------------------------------
95  {
96  fHitFinderModuleLabel = pset.get<std::string>("HitsModuleLabel");
97  fLArG4ModuleLabel = pset.get<std::string>("LArGeantModuleLabel");
98  fCalDataModuleLabel = pset.get<std::string>("CalDataModuleLabel");
99  }
100 
101  //-------------------------------------------------
103  {
105  fHitResidualAll = tfs->make<TH1F>("fHitResidualAll", "Hit Residual All", 1600, -400, 400);
106  fHitResidualAllAlt = tfs->make<TH1F>("fHitResidualAllAlt", "Hit Residual All", 1600, -400, 400);
108  tfs->make<TH1F>("fNumberOfHitsPerEvent", "Number of Hits in Each Event", 10000, 0, 10000);
110  tfs->make<TH2F>("fPeakTimeVsWire", "Peak Time vs Wire Number", 3200, 0, 3200, 9500, 0, 9500);
111 
112  fHTree = tfs->make<TTree>("HTree", "HTree");
113  fHTree->Branch("Evt", &fEvt, "Evt/I");
114  fHTree->Branch("Run", &fRun, "Run/I");
115  fHTree->Branch("WireTotalCharge", &fWireTotalCharge, "WireTotalCharge/F");
116 
117  // === Hit Info ===
118  fHTree->Branch("nHits", &fnHits, "nHits/I");
119  fHTree->Branch("Wire", &fWire, "Wire[nHits]/I");
120  fHTree->Branch("StartTime", &fStartTime, "fStartTime[nHits]/F");
121  fHTree->Branch("EndTime", &fEndTime, "fEndTime[nHits]/F");
122  fHTree->Branch("PeakTime", &fPeakTime, "fPeakTime[nHits]/F");
123  fHTree->Branch("PeakTimeUncert", &fPeakTimeUncert, "fPeakTimeUncert[nHits]/F");
124  fHTree->Branch("Charge", &fCharge, "fCharge[nHits]/F");
125  fHTree->Branch("ChargeUncert", &fChargeUncert, "fChargeUncert[nHits]/F");
126  fHTree->Branch("Multiplicity", &fMultiplicity, "fMultiplicity[nHits]/I");
127  fHTree->Branch("GOF", &fGOF, "fGOF[nHits]/F");
128 
129  // === Total Hit Information ===
130  fHTree->Branch("TotalHitChargePerEvent", &fTotalHitChargePerEvent, "TotalHitChargePerEvent/F");
131 
132  // === Truth Hit Information from BackTrackerService ===
133  fHTree->Branch("TruePeakPos", &fTruePeakPos, "fTruePeakPos[nHits]/F");
134  }
135 
136  //-------------------------------------------------
138  {
139  // ### TTree Run/Event ###
140  fEvt = evt.id().event();
141  fRun = evt.run();
142 
144  auto const clock_data =
146  auto const det_prop =
148 
150  evt.getByLabel(fCalDataModuleLabel, wireVecHandle);
151 
152  // Charge directly from wire info
153  float TotWireCharge = 0;
154 
155  for (size_t wireIter = 0; wireIter < wireVecHandle->size(); wireIter++) {
156  art::Ptr<recob::Wire> wire(wireVecHandle, wireIter);
157  std::vector<float> signal(wire->Signal());
158 
159  for (auto timeIter = signal.begin(); timeIter + 2 < signal.end(); timeIter++) {
160 
161  if (*timeIter < 2) { continue; }
162 
163  TotWireCharge += *timeIter;
164  }
165  }
166 
167  fWireTotalCharge = TotWireCharge;
168 
169  // Reconstructed hit information
171  evt.getByLabel(fHitFinderModuleLabel, hitHandle);
172 
173  std::vector<art::Ptr<recob::Hit>> hits;
174  art::fill_ptr_vector(hits, hitHandle);
175 
176  float TotCharge = 0;
177  int hitCount = 0;
178  fnHits = hitHandle->size();
179  fNumberOfHitsPerEvent->Fill(hitHandle->size());
180 
181  for (size_t numHit = 0; numHit < hitHandle->size(); ++numHit) {
182  // === Finding Channel associated with the hit ===
183  art::Ptr<recob::Hit> hit(hitHandle, numHit);
184 
185  fWire[hitCount] = hit->WireID().Wire;
186  fStartTime[hitCount] = hit->PeakTimeMinusRMS();
187  fEndTime[hitCount] = hit->PeakTimePlusRMS();
188  fPeakTime[hitCount] = hit->PeakTime();
189  fPeakTimeUncert[hitCount] = hit->SigmaPeakTime();
190  fCharge[hitCount] = hit->Integral();
191  fChargeUncert[hitCount] = hit->SigmaIntegral();
192  fMultiplicity[hitCount] = hit->Multiplicity();
193  fGOF[hitCount] = hit->GoodnessOfFit();
194 
195  hitCount++;
196  TotCharge += hit->Integral();
197 
198  fPeakTimeVsWire->Fill(hit->WireID().Wire, hit->PeakTime());
199  } //<---End numHit
200  fTotalHitChargePerEvent = TotCharge;
201 
202  // Truth hit info from BackTracker
203 
204  Float_t TruthHitTime = 0, TruthHitCalculated = 0;
205  int count = 0;
206 
207  double time_tick = sampling_rate(clock_data) / 1000.;
208  double drift_velocity = det_prop.DriftVelocity(det_prop.Efield(), det_prop.Temperature());
209 
210  for (size_t nh = 0; nh < hitHandle->size(); nh++) {
211  // === Finding Channel associated with the hit ===
212  art::Ptr<recob::Hit> hitPoint(hitHandle, nh);
213  auto const& planeID = hitPoint->WireID().asPlaneID();
214 
215  // ===================================================================
216  // Using Track IDE's to locate the XYZ location from truth information
217  // ===================================================================
218  std::vector<sim::TrackIDE> trackides;
219  std::vector<double> xyz;
220  try {
222  trackides = bt_serv->HitToTrackIDEs(clock_data, hitPoint);
223  xyz = bt_serv->HitToXYZ(clock_data, hitPoint);
224  }
225  catch (cet::exception const&) {
226  mf::LogWarning("GausHitFinderAna") << "BackTrackerService Failed";
227  continue;
228  }
229 
230  // ==============================================================
231  // Calculating the truth tick position of the hit using 2 methods
232  // Method 1: ConvertXtoTicks from the detector properties package
233  // Method 2: Actually do the calculation myself to double check things
234  // ==============================================================
235 
236  // ### Method 1 ###
237  TruthHitTime = det_prop.ConvertXToTicks(xyz[0], planeID);
238 
239  // ### Method 2 ###
240  // ================================================
241  // Establishing the x-position of the current plane
242  // ================================================
243  auto const pos = geom->Plane(planeID).GetBoxCenter();
244  double planePos_timeCorr = (pos.X() / drift_velocity) * (1. / time_tick) + 60;
245  //<---x position of plane / drift velocity + 60 (Trigger offset)
246 
247  TruthHitCalculated = ((xyz[0]) / (drift_velocity * time_tick)) + planePos_timeCorr;
248 
249  fTruePeakPos[count] = TruthHitTime;
250  count++;
251  double hitresid = ((TruthHitTime - hitPoint->PeakTime()) / hitPoint->SigmaPeakTime());
252  fHitResidualAll->Fill(hitresid);
253 
254  double hitresidAlt =
255  ((TruthHitCalculated - hitPoint->PeakTime()) / hitPoint->SigmaPeakTime());
256  fHitResidualAllAlt->Fill(hitresidAlt);
257 
258  } //<---End nh loop
259 
260  fHTree->Fill();
261  } // end analyze method
262 
264 
265 } // end of hit namespace
constexpr int kMaxHits
std::vector< sim::TrackIDE > HitToTrackIDEs(detinfo::DetectorClocksData const &clockData, recob::Hit const &hit) const
Declaration of signal hit object.
float SigmaIntegral() const
Initial tdc tick for hit.
Definition: Hit.h:248
float Integral() const
Integral under the calibrated signal waveform of the hit, in tick x ADC units.
Definition: Hit.h:244
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:563
geo::WireID const & WireID() const
Initial tdc tick for hit.
Definition: Hit.h:280
std::vector< double > HitToXYZ(detinfo::DetectorClocksData const &clockData, const recob::Hit &hit) const
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.cc:6
float GoodnessOfFit() const
Degrees of freedom in the determination of the hit signal shape (-1 by default)
Definition: Hit.h:260
short int Multiplicity() const
How many hits could this one be shared with.
Definition: Hit.h:252
Float_t fPeakTimeUncert[kMaxHits]
PlaneGeo const & Plane(PlaneID const &planeid) const
Returns the specified wire.
void hits()
Definition: readHits.C:15
Point_t GetBoxCenter() const
Returns the centre of the box representing the plane.
Definition: PlaneGeo.h:449
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
T get(std::string const &key) const
Definition: ParameterSet.h:314
GausHitFinderAna(fhicl::ParameterSet const &pset)
constexpr PlaneID const & asPlaneID() const
Conversion to PlaneID (for convenience of notation).
Definition: geo_types.h:520
float PeakTimeMinusRMS(float sigmas=+1.) const
Returns a time sigmas RMS away from the peak time.
Definition: Hit.h:290
Detector simulation of raw signals on wires.
std::vector< float > Signal() const
Return a zero-padded full length vector filled with RoI signal.
Definition: Wire.cxx:30
float PeakTime() const
Time of the signal peak, in tick units.
Definition: Hit.h:220
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Encapsulate the construction of a single detector plane.
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
EventNumber_t event() const
Definition: EventID.h:116
Declaration of basic channel signal object.
Base class for creation of raw signals on wires.
float SigmaPeakTime() const
Uncertainty for the signal peak, in tick units.
Definition: Hit.h:224
float PeakTimePlusRMS(float sigmas=+1.) const
Returns a time sigmas RMS away from the peak time.
Definition: Hit.h:285
TCEvent evt
Definition: DataStructs.cxx:8
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:306
RunNumber_t run() const
Definition: Event.cc:29
double sampling_rate(DetectorClocksData const &data)
Returns the period of the TPC readout electronics clock.
Definition: fwd.h:26
EventID id() const
Definition: Event.cc:23
void analyze(const art::Event &evt) override
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33