ClusterCounter5_module.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 // Class: ClusterCounter
3 // Module Type: analyzer
4 // File: ClusterCounter_module.cc
5 //
6 // Same as #4, but use FHiCL parameters validation.
7 // Robert Sulej
8 ////////////////////////////////////////////////////////////////////////
9 
16 #include "art_root_io/TFileService.h"
17 #include "canvas/Persistency/Common/FindManyP.h"
19 #include "fhiclcpp/types/Atom.h"
21 
25 
29 
30 #include "TTree.h"
31 
32 namespace tutorial {
33 
35 public:
36 
37  struct Config {
38  using Name = fhicl::Name;
40 
42  Name("ClusterModuleLabel"),
43  Comment("tag of cluster producer")
44  };
45 
47  Name("MinSize"),
48  Comment("minimum size of clusters")
49  };
50  };
52 
53  explicit ClusterCounter(Parameters const& config);
54 
55  // Plugins should not be copied or assigned.
56  ClusterCounter(ClusterCounter const &) = delete;
57  ClusterCounter(ClusterCounter &&) = delete;
58  ClusterCounter & operator = (ClusterCounter const &) = delete;
60 
61  // Required functions.
62  void analyze(art::Event const & e) override;
63 
64  // Selected optional functions.
65  void beginJob() override;
66  void endJob() override;
67 
68 private:
69  float sumAdc(const std::vector< art::Ptr<recob::Hit> > & hits) const;
70 
72  detinfo::DetectorClocksData const& clockData,
73  const std::vector< art::Ptr<recob::Hit> > & hits,
74  float & fraction, bool & foundEmParent) const;
75 
76  size_t fEvNumber;
77 
78  TTree *fEventTree;
79  size_t fNClusters;
80 
81  TTree *fClusterTree;
82  size_t fNHits;
83  float fAdcSum;
84  float fClean;
85 
86  // ******* fcl parameters *******
88  size_t fMinSize;
89 };
90 
93  fMinSize(config().MinSize())
94 {
95 }
96 
98 {
99  art::ServiceHandle<art::TFileService> tfs; // TTree's are created in the memory managed by ROOT (you don't delete them)
100 
101  fEventTree = tfs->make<TTree>("EventTree", "event by event info");
102  fEventTree->Branch("event", &fEvNumber, "fEvNumber/I");
103  fEventTree->Branch("nclusters", &fNClusters, "fNClusters/I");
104 
105  fClusterTree = tfs->make<TTree>("ClusterTree", "cluster by cluster info");
106  fClusterTree->Branch("event", &fEvNumber, "fEvNumber/I");
107  fClusterTree->Branch("nhits", &fNHits, "fNHits/I");
108  fClusterTree->Branch("adcsum", &fAdcSum, "fAdcSum/F");
109  fClusterTree->Branch("clean", &fClean, "fClean/F");
110 }
111 
113 {
114  fEvNumber = evt.id().event();
115  mf::LogVerbatim("ClusterCounter") << "ClusterCounter module on event #" << fEvNumber;
116 
117  // use auto to make the line shorter when you remember all art types,
118  // the full type of the handle is: art::ValidHandle< std::vector<recob::Cluster> >
119  auto cluHandle = evt.getValidHandle< std::vector<recob::Cluster> >(fClusterModuleLabel);
120 
121  // this will let us look for all hits associated to selected cluster
122  art::FindManyP< recob::Hit > hitsFromClusters(cluHandle, evt, fClusterModuleLabel);
123 
124  fNClusters = 0;
125 
126  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService>()->DataFor(evt);
127  // and here we go by index:
128  for (size_t i = 0; i < cluHandle->size(); ++i)
129  {
130  fNHits = cluHandle->at(i).NHits();
131  fAdcSum = sumAdc(hitsFromClusters.at(i));
132 
133  bool isEM = false;
134  const simb::MCParticle* p = getTruthParticle(clockData, hitsFromClusters.at(i), fClean, isEM);
135  if (p)
136  {
137  if (isEM) { mf::LogVerbatim("ClusterCounter") << "matched mother particle PDG: " << p->PdgCode(); }
138  else { mf::LogVerbatim("ClusterCounter") << "matched particle PDG: " << p->PdgCode(); }
139  }
140  else { mf::LogWarning("ClusterCounter") << "No matcched particle??"; }
141 
142  fClusterTree->Fill();
143 
144  mf::LogVerbatim("ClusterCounter")
145  << "NHits() = " << fNHits << ", assn size = " << hitsFromClusters.at(i).size()
146  << " SummedADC() = " << cluHandle->at(i).SummedADC() << ", sum hits adc = " << fAdcSum;
147 
148  if (fNHits >= fMinSize) { ++fNClusters; }
149  }
150  fEventTree->Fill();
151 }
152 
154  const std::vector< art::Ptr<recob::Hit> > & hits,
155  float & fraction, bool & foundEmParent) const
156 {
157  const simb::MCParticle* mcParticle =0;
158  fraction = 0;
159  foundEmParent = false;
160 
163  std::unordered_map<int, double> trkIDE;
164  for (auto const & h : hits)
165  {
166  for (auto const & ide : bt_serv->HitToTrackIDEs(clockData, h)) // loop over std::vector<sim::TrackIDE>
167  {
168  trkIDE[ide.trackID] += ide.energy; // sum energy contribution by each track ID
169  }
170  }
171 
172  int best_id = 0;
173  double tot_e = 0, max_e = 0;
174  for (auto const & contrib : trkIDE)
175  {
176  tot_e += contrib.second; // sum total energy in these hits
177  if (contrib.second > max_e) // find track ID corresponding to max energy
178  {
179  max_e = contrib.second;
180  best_id = contrib.first;
181  }
182  }
183 
184  if ((max_e > 0) && (tot_e > 0)) // ok, found something reasonable
185  {
186  if (best_id < 0) // NOTE: negative ID means this is EM activity
187  { // caused by track with the same but positive ID
188  best_id = -best_id; // --> we'll find mother MCParticle of these hits
189  foundEmParent = true;
190  }
191  mcParticle = pi_serv->TrackIdToParticle_P(best_id); // MCParticle corresponding to track ID
192  fraction = max_e / tot_e;
193  }
194  else { mf::LogWarning("ClusterCounter") << "No energy deposits??"; }
195 
196  return mcParticle;
197 }
198 
200 {
201  float sum = 0;
202  for (auto const & h : hits)
203  {
204  sum += h->SummedADC();
205  }
206  return sum;
207 }
208 
210 {
211  mf::LogVerbatim("ClusterCounter") << "ClusterCounter finished job";
212 }
213 
214 } // tutorial namespace
215 
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
int PdgCode() const
Definition: MCParticle.h:212
std::vector< sim::TrackIDE > HitToTrackIDEs(detinfo::DetectorClocksData const &clockData, recob::Hit const &hit) const
const simb::MCParticle * TrackIdToParticle_P(int id) const
fhicl::Atom< art::InputTag > ClusterModuleLabel
struct vector vector
ChannelGroupService::Name Name
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
Particle class.
void analyze(art::Event const &e) override
ClusterCounter & operator=(ClusterCounter const &)=delete
const double e
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:67
static Config * config
Definition: config.cpp:1054
float sumAdc(const std::vector< art::Ptr< recob::Hit > > &hits) const
ValidHandle< PROD > getValidHandle(InputTag const &tag) const
Definition: DataViewImpl.h:441
p
Definition: test.py:223
Declaration of signal hit object.
#define Comment
Contains all timing reference information for the detector.
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
EventNumber_t event() const
Definition: EventID.h:116
const simb::MCParticle * getTruthParticle(detinfo::DetectorClocksData const &clockData, const std::vector< art::Ptr< recob::Hit > > &hits, float &fraction, bool &foundEmParent) const
TCEvent evt
Definition: DataStructs.cxx:7
ClusterCounter(Parameters const &config)
EventID id() const
Definition: Event.cc:34