LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
SpacePointAna_module.cc
Go to the documentation of this file.
1 //
2 // Name: SpacePointAna_module.cc
3 //
4 // Purpose: Module SpacePointAna.
5 //
6 //
7 // Configuration parameters.
8 //
9 // HitModuleLabel: Hit module label.
10 // UseClusterHits: If true, use clustered hits (otherwise all hits).
11 // ClusterModuleLabel: Cluster module label.
12 // UseMC: Use MC truth information.
13 // SpacePointAlgTime: SpacePointAlg configuration (loose time cuts).
14 // SpacePointAlgSep: SpacePointAlg configuration (loose separation cuts).
15 // SpacePointAlgDefault: SpacePointAlg configuration (default cuts).
16 //
17 // Created: 2-Aug-2011 H. Greenlee
18 //
19 
20 #include <iostream>
21 #include <map>
22 #include <vector>
23 #include <algorithm>
24 
32 
44 
45 #include "TH1F.h"
46 #include "TH2F.h"
47 
48 
49 namespace trkf {
50 
52  {
53  public:
54 
55  // Constructors, destructor
56 
57  explicit SpacePointAna(fhicl::ParameterSet const& pset);
58  virtual ~SpacePointAna();
59 
60  // Book histograms.
61 
62  void bookHistograms(bool mc);
63 
64  // Overrides.
65 
66  void analyze(const art::Event& evt);
67 
68  private:
69 
70  // Fcl Attributes.
71 
72  const SpacePointAlg fSptalgTime; // Algorithm object (increased time cut).
73  const SpacePointAlg fSptalgSep; // Algorithm object (increased sepataion cut).
74  const SpacePointAlg fSptalgDefault; // Algorithm object (default cuts).
75  std::string fHitModuleLabel;
77  std::string fClusterModuleLabel;
78  bool fUseMC;
79  double fMinX; // Minimum x.
80  double fMaxX; // Maximum x.
81  double fMinY; // Minimum y.
82  double fMaxY; // Maximum y.
83  double fMinZ; // Minimum z.
84  double fMaxZ; // Maximum z.
85 
86  // Histograms.
87 
88  bool fBooked; // Have histograms been booked yet?
89  TH1F* fHDTUE; // U-drift electrons time difference.
90  TH1F* fHDTVE; // V-drift electrons time difference.
91  TH1F* fHDTWE; // W-drift electrons time difference.
92  TH1F* fHDTUPull; // U-drift electrons time pull.
93  TH1F* fHDTVPull; // V-drift electrons time pull.
94  TH1F* fHDTWPull; // W-drift electrons time pull.
95  TH1F* fHDTUV; // U-V time difference.
96  TH1F* fHDTVW; // V-W time difference.
97  TH1F* fHDTWU; // W-U time difference.
98  TH2F* fHDTUVU; // U-V time difference vs. U.
99  TH2F* fHDTUVV; // U-V time difference vs. V.
100  TH2F* fHDTVWV; // V-W time difference vs. V.
101  TH2F* fHDTVWW; // V-W time difference vs. W.
102  TH2F* fHDTWUW; // W-U time difference vs. W.
103  TH2F* fHDTWUU; // W-U time difference vs. U.
104  TH1F* fHS; // Spatial separation.
105  TH1F* fHchisq; // Space point chisquare.
106  TH1F* fHx; // X position.
107  TH1F* fHy; // Y position.
108  TH1F* fHz; // Z position.
109  TH1F* fHAmpU; // U hit amplitude.
110  TH1F* fHAmpV; // V hit amplitude.
111  TH1F* fHAmpW; // W hit amplitude.
112  TH1F* fHAreaU; // U hit area.
113  TH1F* fHAreaV; // V hit area.
114  TH1F* fHAreaW; // W hit area.
115  TH1F* fHSumU; // U hit sum ADC.
116  TH1F* fHSumV; // V hit sum ADC.
117  TH1F* fHSumW; // W hit sum ADC.
118  TH1F* fHMCdx; // X residual (reco vs. mc truth).
119  TH1F* fHMCdy; // Y residual (reco vs. mc truth).
120  TH1F* fHMCdz; // Z residual (reco vs. mc truth).
121  TH1F* fHMCxpull; // X pull (reco vs. mc truth).
122  TH1F* fHMCypull; // Y pull (reco vs. mc truth).
123  TH1F* fHMCzpull; // Z pull (reco vs. mc truth).
124 
125  // Statistics.
126 
128  };
129 
131 
132  SpacePointAna::SpacePointAna(const fhicl::ParameterSet& pset)
133  //
134  // Purpose: Constructor.
135  //
136  // Arguments: pset - Module parameters.
137  //
138  : EDAnalyzer(pset)
139  , fSptalgTime(pset.get<fhicl::ParameterSet>("SpacePointAlgTime"))
140  , fSptalgSep(pset.get<fhicl::ParameterSet>("SpacePointAlgSep"))
141  , fSptalgDefault(pset.get<fhicl::ParameterSet>("SpacePointAlgDefault"))
142  , fHitModuleLabel(pset.get<std::string>("HitModuleLabel"))
143  , fUseClusterHits(pset.get<bool>("UseClusterHits"))
144  , fClusterModuleLabel(pset.get<std::string>("ClusterModuleLabel"))
145  , fUseMC(pset.get<bool>("UseMC"))
146  , fMinX(pset.get<double>("MinX", -1.e10))
147  , fMaxX(pset.get<double>("MaxX", 1.e10))
148  , fMinY(pset.get<double>("MinY", -1.e10))
149  , fMaxY(pset.get<double>("MaxY", 1.e10))
150  , fMinZ(pset.get<double>("MinZ", -1.e10))
151  , fMaxZ(pset.get<double>("MaxZ", 1.e10))
152  , fBooked(false)
153  , fHDTUE(0)
154  , fHDTVE(0)
155  , fHDTWE(0)
156  , fHDTUPull(0)
157  , fHDTVPull(0)
158  , fHDTWPull(0)
159  , fHDTUV(0)
160  , fHDTVW(0)
161  , fHDTWU(0)
162  , fHDTUVU(0)
163  , fHDTUVV(0)
164  , fHDTVWV(0)
165  , fHDTVWW(0)
166  , fHDTWUW(0)
167  , fHDTWUU(0)
168  , fHS(0)
169  , fHchisq(0)
170  , fHx(0)
171  , fHy(0)
172  , fHz(0)
173  , fHAmpU(0)
174  , fHAmpV(0)
175  , fHAmpW(0)
176  , fHAreaU(0)
177  , fHAreaV(0)
178  , fHAreaW(0)
179  , fHSumU(0)
180  , fHSumV(0)
181  , fHSumW(0)
182  , fHMCdx(0)
183  , fHMCdy(0)
184  , fHMCdz(0)
185  , fHMCxpull(0)
186  , fHMCypull(0)
187  , fHMCzpull(0)
188  , fNumEvent(0)
189  {
190 
191  // Report.
192 
193  mf::LogInfo("SpacePointAna")
194  << "SpacePointAna configured with the following parameters:\n"
195  << " HitModuleLabel = " << fHitModuleLabel << "\n"
196  << " UseClusterHits = " << fUseClusterHits << "\n"
197  << " ClusterModuleLabel = " << fClusterModuleLabel << "\n"
198  << " UseMC = " << fUseMC;
199  }
200 
202  //
203  // Purpose: Destructor.
204  //
205  {}
206 
208  //
209  // Purpose: Book histograms.
210  //
211  {
212  if(!fBooked) {
213  fBooked = true;
214 
217  art::TFileDirectory dir = tfs->mkdir("sptana", "SpacePointAna histograms");
218 
219  unsigned int nwiresU=0, nwiresV=0, nwiresW=0;
220 
221  // Figure out the number of wires in U, V, and W planes.
222 
223  // Loop over cryostats, tpcs, and planes.
224 
225  for(unsigned int cstat = 0; cstat < geom->Ncryostats(); ++cstat){
226 
227  const geo::CryostatGeo& cryogeom = geom->Cryostat(cstat);
228  unsigned int const ntpc = cryogeom.NTPC();
229 
230  for(unsigned int tpc = 0; tpc < ntpc; ++tpc) {
231 
232  const geo::TPCGeo& tpcgeom = cryogeom.TPC(tpc);
233  unsigned int const nplane = tpcgeom.Nplanes();
234 
235  for(unsigned int plane = 0; plane < nplane; ++plane) {
236 
237  const geo::PlaneGeo& pgeom = tpcgeom.Plane(plane);
238  unsigned int nwires = pgeom.Nwires();
239  geo::View_t view = pgeom.View();
240  if(view == geo::kU)
241  nwiresU = nwires;
242  else if(view == geo::kV)
243  nwiresV = nwires;
244  else if(view == geo::kZ)
245  nwiresW = nwires;
246  }
247  }
248  }
249 
250 
251 
252 
253  if(mc && fUseMC) {
254  fHDTUE = dir.make<TH1F>("MCDTUE", "U-Drift Electrons Time Difference", 100, -5., 5.);
255  fHDTVE = dir.make<TH1F>("MCDTVE", "V-Drift Electrons Time Difference", 100, -5., 5.);
256  fHDTWE = dir.make<TH1F>("MCDTWE", "W-Drift Electrons Time Difference", 100, -5., 5.);
257  fHDTUPull = dir.make<TH1F>("MCDTUPull", "U-Drift Electrons Time Pull", 100, -50., 50.);
258  fHDTVPull = dir.make<TH1F>("MCDTVPull", "V-Drift Electrons Time Pull", 100, -50., 50.);
259  fHDTWPull = dir.make<TH1F>("MCDTWPull", "W-Drift Electrons Time Pull", 100, -50., 50.);
260  }
261  if(!fSptalgTime.merge()) {
262  fHDTUV = dir.make<TH1F>("DTUV", "U-V time difference", 100, -20., 20.);
263  fHDTVW = dir.make<TH1F>("DTVW", "V-W time difference", 100, -20., 20.);
264  fHDTWU = dir.make<TH1F>("DTWU", "W-U time difference", 100, -20., 20.);
265  fHDTUVU = dir.make<TH2F>("DTUVU", "U-V time difference vs. U",
266  100, 0., double(nwiresU), 100, -20., 20.);
267  fHDTUVV = dir.make<TH2F>("DTUVV", "U-V time difference vs. V",
268  100, 0., double(nwiresV), 100, -20., 20.);
269  fHDTVWV = dir.make<TH2F>("DTVWV", "V-W time difference vs. V",
270  100, 0., double(nwiresV), 100, -20., 20.);
271  fHDTVWW = dir.make<TH2F>("DTVWW", "V-W time difference vs. W",
272  100, 0., double(nwiresW), 100, -20., 20.);
273  fHDTWUW = dir.make<TH2F>("DTWUW", "W-U time difference vs. W",
274  100, 0., double(nwiresW), 100, -20., 20.);
275  fHDTWUU = dir.make<TH2F>("DTWUU", "W-U time difference vs. U",
276  100, 0., double(nwiresU), 100, -20., 20.);
277  fHS = dir.make<TH1F>("DS", "Spatial Separation", 100, -2., 2.);
278  }
279  fHchisq = dir.make<TH1F>("chisq", "Chisquare", 100, 0., 20.);
280 
281  fHx = dir.make<TH1F>("xpos", "X Position",
282  100, 0., 2.*geom->DetHalfWidth());
283  fHy = dir.make<TH1F>("ypos", "Y Position",
284  100, -geom->DetHalfHeight(), geom->DetHalfHeight());
285  fHz = dir.make<TH1F>("zpos", "Z Position",
286  100, 0., geom->DetLength());
287  fHAmpU = dir.make<TH1F>("ampU", "U Hit Amplitude", 50, 0., 50.);
288  fHAmpV = dir.make<TH1F>("ampV", "V Hit Amplitude", 50, 0., 50.);
289  fHAmpW = dir.make<TH1F>("ampW", "W Hit Amplitude", 50, 0., 50.);
290  fHAreaU = dir.make<TH1F>("areaU", "U Hit Area", 100, 0., 500.);
291  fHAreaV = dir.make<TH1F>("areaV", "V Hit Area", 100, 0., 500.);
292  fHAreaW = dir.make<TH1F>("areaW", "W Hit Area", 100, 0., 500.);
293  fHSumU = dir.make<TH1F>("sumU", "U Hit Sum ADC", 100, 0., 500.);
294  fHSumV = dir.make<TH1F>("sumV", "V Hit Sum ADC", 100, 0., 500.);
295  fHSumW = dir.make<TH1F>("sumW", "W Hit Sum ADC", 100, 0., 500.);
296  if(mc && fUseMC) {
297  fHMCdx = dir.make<TH1F>("MCdx", "X MC Residual", 100, -1., 1.);
298  fHMCdy = dir.make<TH1F>("MCdy", "Y MC Residual", 100, -1., 1.);
299  fHMCdz = dir.make<TH1F>("MCdz", "Z MC Residual", 100, -1., 1.);
300  fHMCxpull = dir.make<TH1F>("MCxpull", "X MC Pull", 100, -50., 50.);
301  fHMCypull = dir.make<TH1F>("MCypull", "Y MC Pull", 100, -50., 50.);
302  fHMCzpull = dir.make<TH1F>("MCzpull", "Z MC Pull", 100, -50., 50.);
303  }
304  }
305  }
306 
308  //
309  // Purpose: Analyze method.
310  //
311  // Arguments: event - Art event.
312  //
313  {
314  ++fNumEvent;
315 
316  // Make sure histograms are booked.
317 
318  bool mc = !evt.isRealData();
319  bookHistograms(mc);
320 
321  // Get Services.
322 
324  const detinfo::DetectorProperties* detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
325 
326  // Get Hits.
327 
329 
330  if(fUseClusterHits) {
331 
332  // Get clusters.
333 
335  evt.getByLabel(fClusterModuleLabel, clusterh);
336 
337  // Get hits from all clusters.
339 
340  if(clusterh.isValid()) {
341  int nclus = clusterh->size();
342 
343  for(int i = 0; i < nclus; ++i) {
344  art::Ptr<recob::Cluster> pclus(clusterh, i);
345  std::vector< art::Ptr<recob::Hit> > clushits = fm.at(i);
346  int nhits = clushits.size();
347  hits.reserve(hits.size() + nhits);
348 
349  for(std::vector< art::Ptr<recob::Hit> >::const_iterator ihit = clushits.begin(); ihit != clushits.end(); ++ihit) {
350  hits.push_back(*ihit);
351  }
352  }
353  }
354  }
355  else {
356 
357  // Get unclustered hits.
358 
360  evt.getByLabel(fHitModuleLabel, hith);
361  if(hith.isValid()) {
362  int nhits = hith->size();
363  hits.reserve(nhits);
364 
365  for(int i = 0; i < nhits; ++i)
366  hits.push_back(art::Ptr<recob::Hit>(hith, i));
367  }
368  }
369 
370  // Fill histograms that don't depend on space points.
371 
372  if(mc && fUseMC) {
373 
375 
376  // Loop over hits and fill hit-electron time difference histogram.
377 
378  for(art::PtrVector<recob::Hit>::const_iterator ihit = hits.begin();
379  ihit != hits.end(); ++ihit) {
380  const recob::Hit& hit = **ihit;
381 
382  //unsigned int channel = hit.Channel();
383  //geo::View_t geo_view = geom->View(channel);
384  //geo::View_t hit_view = hit.View();
385  //assert(geo_view == hit_view);
386  double tpeak = hit.PeakTime();
387  double terr = hit.SigmaPeakTime();
388 
389  //assert(channel == hit.Channel());
390 
391  // Loop over electrons associated with this hit/channel and fill
392  // hit-electron time difference histograms.
393 
394  // loop over the map of TDC to sim::IDE to get the TDC for each energy dep
395  // Find the average time in ticks for this hit.
396 
397  //double sumw = 0.;
398  //double sumt = 0.;
399 
400  //const std::map<unsigned short, std::vector<sim::IDE> > &idemap = simchan.TDCIDEMap();
401  //std::map<unsigned short, std::vector<sim::IDE> >::const_iterator mitr = idemap.begin();
402  //for(mitr = idemap.begin(); mitr != idemap.end(); mitr++) {
403  // double tdc = double((*mitr).first);
404  // if(tdc >= tstart && tdc <= tend) {
405  // const std::vector<sim::IDE>& idevec = (*mitr).second;
406  // for(std::vector<sim::IDE>::const_iterator iide=idevec.begin();
407  // iide != idevec.end(); ++iide) {
408  // const sim::IDE& ide = *iide;
409  // double w = ide.numElectrons;
410  // sumw += w;
411  // sumt += w*tdc;
412  // }
413  // }
414  //}
415  //double tav = 0.;
416  //if(sumw != 0.)
417  // tav = sumt / sumw;
418 
419  bool tav_ok = true;
420  double tav = 0.;
421  try {
422  std::vector<double> hitxyz = bt_serv->HitToXYZ(*ihit);
423  tav = detprop->ConvertXToTicks(hitxyz[0], (*ihit)->WireID().Plane, (*ihit)->WireID().TPC, (*ihit)->WireID().Cryostat);
424  }
425  catch(cet::exception& x) {
426  tav_ok = false;
427  }
428  if(tav_ok) {
429  if((*ihit)->View() == geo::kU) {
430  fHDTUE->Fill(tpeak - tav);
431  fHDTUPull->Fill((tpeak - tav) / terr);
432  }
433  else if((*ihit)->View() == geo::kV) {
434  fHDTVE->Fill(tpeak - tav);
435  fHDTVPull->Fill((tpeak - tav) / terr);
436  }
437  else if((*ihit)->View() == geo::kZ) {
438  fHDTWE->Fill(tpeak - tav);
439  fHDTWPull->Fill((tpeak - tav) / terr);
440  }
441  else
442  throw cet::exception("SpacePointAna") << "Bad view = " << (*ihit)->View() << "\n";
443  }
444  }
445  }
446 
447  std::vector<recob::SpacePoint> spts1; // For time histograms.
448  std::vector<recob::SpacePoint> spts2; // For separation histogram.
449  std::vector<recob::SpacePoint> spts3; // Default cuts.
450 
451  // If nonzero time cut is specified, make space points using that
452  // time cut (for time histograms).
453 
454  if(!fSptalgTime.merge()) {
455  if(mc && fUseMC)
456  fSptalgTime.makeMCTruthSpacePoints(hits, spts1);
457  else
458  fSptalgTime.makeSpacePoints(hits, spts1);
459 
460  // Report number of space points.
461 
462  LOG_DEBUG("SpacePointAna") << "Found " << spts1.size()
463  << " space points using special time cut.";
464  }
465 
466  // If nonzero separation cut is specified, make space points using that
467  // separation cut (for separation histogram).
468 
469  if(!fSptalgSep.merge()) {
470  if(mc && fUseMC)
471  fSptalgSep.makeMCTruthSpacePoints(hits, spts2);
472  else
473  fSptalgSep.makeSpacePoints(hits, spts2);
474 
475  // Report number of space points.
476 
477  LOG_DEBUG("SpacePointAna") << "Found " << spts2.size()
478  << " space points using special seperation cut.";
479  }
480 
481  // Make space points using default cuts.
482 
483  if(mc && fUseMC)
485  else
486  fSptalgDefault.makeSpacePoints(hits, spts3);
487 
488  // Report number of space points.
489 
490  LOG_DEBUG("SpacePointAna") << "Found " << spts3.size()
491  << " space points using default cuts.";
492 
493  if(!fSptalgTime.merge()) {
494 
495  // Loop over space points and fill time histograms.
496 
497  for(std::vector<recob::SpacePoint>::const_iterator i = spts1.begin();
498  i != spts1.end(); ++i) {
499  const recob::SpacePoint& spt = *i;
500  if(spt.XYZ()[0] < fMinX || spt.XYZ()[0] > fMaxX ||
501  spt.XYZ()[1] < fMinY || spt.XYZ()[1] > fMaxY ||
502  spt.XYZ()[2] < fMinZ || spt.XYZ()[2] > fMaxZ)
503  continue;
504 
505  // Get hits associated with this SpacePoint.
506 
508 
509  // Make a double loop over hits and fill hit time difference histograms.
510 
512  ihit != spthits.end(); ++ihit) {
513  const recob::Hit& hit1 = **ihit;
514 
515  geo::WireID hit1WireID = hit1.WireID();
516  unsigned int tpc1, plane1, wire1;
517  tpc1 = hit1WireID.TPC;
518  plane1 = hit1WireID.Plane;
519  wire1 = hit1WireID.Wire;
520 
521  geo::View_t view1 = hit1.View();
522  double t1 = fSptalgTime.correctedTime(hit1);
523 
525  jhit != spthits.end(); ++jhit) {
526  const recob::Hit& hit2 = **jhit;
527 
528  geo::WireID hit2WireID = hit2.WireID();
529  unsigned int tpc2, plane2, wire2;
530  tpc2 = hit2WireID.TPC;
531  plane2 = hit2WireID.Plane;
532  wire2 = hit2WireID.Wire;
533 
534  // Require same tpc, different view.
535 
536  if(tpc1 == tpc2 && plane1 != plane2) {
537 
538  geo::View_t view2 = hit2.View();
539  double t2 = fSptalgTime.correctedTime(hit2);
540 
541  if(view1 == geo::kU) {
542  if(view2 == geo::kV) {
543  fHDTUV->Fill(t1-t2);
544  fHDTUVU->Fill(double(wire1), t1-t2);
545  fHDTUVV->Fill(double(wire2), t1-t2);
546  }
547  if(view2 == geo::kZ) {
548  fHDTWU->Fill(t2-t1);
549  fHDTWUW->Fill(double(wire2), t2-t1);
550  fHDTWUU->Fill(double(wire1), t2-t1);
551  }
552  }
553  if(view1 == geo::kV) {
554  if(view2 == geo::kZ) {
555  fHDTVW->Fill(t1-t2);
556  fHDTVWV->Fill(double(wire1), t1-t2);
557  fHDTVWW->Fill(double(wire2), t1-t2);
558  }
559  if(view2 == geo::kU) {
560  fHDTUV->Fill(t2-t1);
561  fHDTUVU->Fill(double(wire2), t2-t1);
562  fHDTUVV->Fill(double(wire1), t2-t1);
563  }
564  }
565  if(view1 == geo::kZ) {
566  if(view2 == geo::kU) {
567  fHDTWU->Fill(t1-t2);
568  fHDTWUW->Fill(double(wire1), t1-t2);
569  fHDTWUU->Fill(double(wire2), t1-t2);
570  }
571  if(view2 == geo::kV) {
572  fHDTVW->Fill(t2-t1);
573  fHDTVWV->Fill(double(wire2), t2-t1);
574  fHDTVWW->Fill(double(wire1), t2-t1);
575  }
576  }
577  }
578  }
579  }
580  }
581 
582  // Loop over space points and fill seperation histograms.
583 
584  for(std::vector<recob::SpacePoint>::const_iterator i = spts2.begin();
585  i != spts2.end(); ++i) {
586  const recob::SpacePoint& spt = *i;
587  if(spt.XYZ()[0] < fMinX || spt.XYZ()[0] > fMaxX ||
588  spt.XYZ()[1] < fMinY || spt.XYZ()[1] > fMaxY ||
589  spt.XYZ()[2] < fMinZ || spt.XYZ()[2] > fMaxZ)
590  continue;
591 
592  // Get hits associated with this SpacePoint.
593 
595 
596  // Fill separation histogram.
597 
598  double sep = fSptalgSep.separation(spthits);
599  fHS->Fill(sep);
600  }
601  }
602 
603  // Loop over default space points and fill histograms.
604 
606  i != spts3.end(); ++i) {
607  const recob::SpacePoint& spt = *i;
608  if(spt.XYZ()[0] < fMinX || spt.XYZ()[0] > fMaxX ||
609  spt.XYZ()[1] < fMinY || spt.XYZ()[1] > fMaxY ||
610  spt.XYZ()[2] < fMinZ || spt.XYZ()[2] > fMaxZ)
611  continue;
612 
613  fHchisq->Fill(spt.Chisq());
614  fHx->Fill(spt.XYZ()[0]);
615  fHy->Fill(spt.XYZ()[1]);
616  fHz->Fill(spt.XYZ()[2]);
617 
618  // Get hits associated with this SpacePoint.
619 
620  std::vector<art::Ptr<recob::Hit> > spthits;
622  for( auto const& ptr : av_spthits ){ spthits.push_back(ptr);}
623 
624  // Fill single hit histograms.
625 
626  for(art::PtrVector<recob::Hit>::const_iterator ihit = spthits.begin();
627  ihit != spthits.end(); ++ihit) {
628  const recob::Hit& hit = **ihit;
629 
630  geo::View_t view = hit.View();
631 
632  if(view == geo::kU) {
633  fHAmpU->Fill(hit.PeakAmplitude());
634  fHAreaU->Fill(hit.Integral());
635  fHSumU->Fill(hit.SummedADC());
636  }
637  if(view == geo::kV) {
638  fHAmpV->Fill(hit.PeakAmplitude());
639  fHAreaV->Fill(hit.Integral());
640  fHSumV->Fill(hit.SummedADC());
641  }
642  if(view == geo::kZ) {
643  fHAmpW->Fill(hit.PeakAmplitude());
644  fHAreaW->Fill(hit.Integral());
645  fHSumW->Fill(hit.SummedADC());
646  }
647  }
648 
649  if(mc && fUseMC) {
650 
652 
653  try {
654  std::vector<double> mcxyz = bt_serv->SpacePointHitsToWeightedXYZ(spthits);
655  fHMCdx->Fill(spt.XYZ()[0] - mcxyz[0]);
656  fHMCdy->Fill(spt.XYZ()[1] - mcxyz[1]);
657  fHMCdz->Fill(spt.XYZ()[2] - mcxyz[2]);
658  if(spt.ErrXYZ()[0] > 0.)
659  fHMCxpull->Fill((spt.XYZ()[0] - mcxyz[0]) / std::sqrt(spt.ErrXYZ()[0]));
660  if(spt.ErrXYZ()[2] > 0.)
661  fHMCypull->Fill((spt.XYZ()[1] - mcxyz[1]) / std::sqrt(spt.ErrXYZ()[2]));
662  if(spt.ErrXYZ()[5] > 0.)
663  fHMCzpull->Fill((spt.XYZ()[2] - mcxyz[2]) / std::sqrt(spt.ErrXYZ()[5]));
664  }
665  catch(cet::exception& x) {}
666  }
667  }
668  }
669 }
Float_t x
Definition: compare.C:6
double correctedTime(const recob::Hit &hit) const
geo::Length_t DetHalfWidth(geo::TPCID const &tpcid) const
Returns the half width of the active volume of the specified TPC.
TTree * t1
Definition: plottest35.C:26
const std::vector< double > HitToXYZ(const recob::Hit &hit)
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
Planes which measure V.
Definition: geo_types.h:77
geo::WireID WireID() const
Initial tdc tick for hit.
Definition: Hit.h:234
unsigned int Nplanes() const
Number of planes in this tpc.
Definition: TPCGeo.h:145
iterator begin()
Definition: PtrVector.h:223
Declaration of signal hit object.
Geometry information for a single TPC.
Definition: TPCGeo.h:37
Double32_t Chisq() const
Definition: SpacePoint.h:67
STL namespace.
Planes which measure Z direction.
Definition: geo_types.h:79
float Integral() const
Integral under the calibrated signal waveform of the hit, in tick x ADC units.
Definition: Hit.h:225
geo::View_t View() const
View for the plane of the hit.
Definition: Hit.h:233
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:313
bool isRealData() const
Definition: Event.h:83
Geometry information for a single cryostat.
Definition: CryostatGeo.h:36
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
const art::PtrVector< recob::Hit > & getAssociatedHits(const recob::SpacePoint &spt) const
TFileDirectory mkdir(std::string const &dir, std::string const &descr="")
const SpacePointAlg fSptalgSep
void makeMCTruthSpacePoints(const art::PtrVector< recob::Hit > &hits, std::vector< recob::SpacePoint > &spts) const
virtual double ConvertXToTicks(double X, int p, int t, int c) const =0
float PeakAmplitude() const
The estimated amplitude of the hit at its peak, in ADC units.
Definition: Hit.h:222
double separation(const art::PtrVector< recob::Hit > &hits) const
bool isValid() const
Definition: Handle.h:190
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
Planes which measure U.
Definition: geo_types.h:76
View_t View() const
Which coordinate does this plane measure.
Definition: PlaneGeo.h:171
void hits()
Definition: readHits.C:15
geo::Length_t DetHalfHeight(geo::TPCID const &tpcid) const
Returns the half height of the active volume of the specified TPC.
intermediate_table::const_iterator const_iterator
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
parameter set interface
Geometry information for a single wire plane.The plane is represented in the geometry by a solid whic...
Definition: PlaneGeo.h:78
iterator end()
Definition: PtrVector.h:237
unsigned int NTPC() const
Number of TPCs in this cryostat.
Definition: CryostatGeo.h:155
geo::Length_t DetLength(geo::TPCID const &tpcid) const
Returns the length of the active volume of the specified TPC.
CryostatGeo const & Cryostat(geo::CryostatID const &cryoid) const
Returns the specified cryostat.
const Double32_t * XYZ() const
Definition: SpacePoint.h:65
EDAnalyzer(Table< Config > const &config)
Definition: EDAnalyzer.h:100
const Double32_t * ErrXYZ() const
Definition: SpacePoint.h:66
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:258
Declaration of cluster object.
const std::vector< double > SpacePointHitsToWeightedXYZ(std::vector< art::Ptr< recob::Hit >> const &hits)
Detector simulation of raw signals on wires.
TTree * t2
Definition: plottest35.C:36
data_t::const_iterator const_iterator
Definition: PtrVector.h:61
void makeSpacePoints(const art::PtrVector< recob::Hit > &hits, std::vector< recob::SpacePoint > &spts) const
float PeakTime() const
Time of the signal peak, in tick units.
Definition: Hit.h:219
T * make(ARGS...args) const
SpacePointAna(fhicl::ParameterSet const &pset)
Utility object to perform functions of association.
Encapsulate the construction of a single detector plane.
const SpacePointAlg fSptalgTime
TDirectory * dir
Definition: macro.C:5
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
const TPCGeo & TPC(unsigned int itpc) const
Return the itpc&#39;th TPC in the cryostat.
void analyze(const art::Event &evt)
unsigned int Nwires() const
Number of wires in this plane.
Definition: PlaneGeo.h:250
#define LOG_DEBUG(id)
float SummedADC() const
The sum of calibrated ADC counts of the hit (0. by default)
Definition: Hit.h:224
const SpacePointAlg fSptalgDefault
float SigmaPeakTime() const
Uncertainty for the signal peak, in tick units.
Definition: Hit.h:220
bool merge() const
Definition: SpacePointAlg.h:89
PlaneGeo const & Plane(geo::View_t view) const
Return the plane in the tpc with View_t view.
Definition: TPCGeo.cxx:299
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:49
TCEvent evt
Definition: DataStructs.cxx:5
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:203
Algorithm for generating space points from hits.
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
Encapsulate the construction of a single detector plane.