ClusterCounter2_module.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 // Class: ClusterCounter2
3 // Module Type: analyzer
4 // File: ClusterCounter2_module.cc
5 //
6 // Access clusters, fill ROOT tree with a simple info from clusters.
7 // Robert Sulej
8 ////////////////////////////////////////////////////////////////////////
9 
16 #include "art_root_io/TFileService.h"
18 #include "fhiclcpp/ParameterSet.h"
20 
22 
23 #include "TTree.h"
24 
25 namespace tutorial {
26 
28 public:
29  explicit ClusterCounter2(fhicl::ParameterSet const & p);
30 
31  // Plugins should not be copied or assigned.
32  ClusterCounter2(ClusterCounter2 const &) = delete;
33  ClusterCounter2(ClusterCounter2 &&) = delete;
34  ClusterCounter2 & operator = (ClusterCounter2 const &) = delete;
36 
37  // Required functions.
38  void analyze(art::Event const & e) override;
39 
40  // Selected optional functions.
41  void beginJob() override;
42  void endJob() override;
43 
44 private:
45  size_t fEvNumber;
46 
47  TTree *fEventTree;
48  size_t fNClusters;
49 
50  TTree *fClusterTree;
51  size_t fNHits;
52 
53  // ******* fcl parameters *******
55  size_t fMinSize;
56 };
57 
59 {
60  fClusterModuleLabel = p.get< std::string >("ClusterModuleLabel");
61  fMinSize = p.get< size_t >("MinSize");
62 }
63 
65 {
66  art::ServiceHandle<art::TFileService> tfs; // TTree's are created in the memory managed by ROOT (you don't delete them)
67 
68  fEventTree = tfs->make<TTree>("EventTree", "event by event info");
69  fEventTree->Branch("event", &fEvNumber, "fEvNumber/I");
70  fEventTree->Branch("nclusters", &fNClusters, "fNClusters/I");
71 
72  fClusterTree = tfs->make<TTree>("ClusterTree", "cluster by cluster info");
73  fClusterTree->Branch("event", &fEvNumber, "fEvNumber/I");
74  fClusterTree->Branch("nhits", &fNHits, "fNHits/I");
75 }
76 
78 {
79  fEvNumber = evt.id().event();
80  mf::LogVerbatim("ClusterCounter2") << "ClusterCounter2 module on event #" << fEvNumber;
81 
82  // use auto to make the line shorter when you remember all art types,
83  // the full type being de-referenced here is: art::ValidHandle< std::vector<recob::Cluster> >
84  auto const & clusters = *evt.getValidHandle< std::vector<recob::Cluster> >(fClusterModuleLabel);
85 
86  fNClusters = 0;
87 
88  // if you are old c++ granpa (or need cluster index, which may indeed happen!):
89  // for (size_t i = 0; i < clusters.size(); ++i)
90  // or:
91  for (auto const & clu : clusters) // loop over recob::Cluster's stored in the std::vector
92  {
93  fNHits = clu.NHits();
94  fClusterTree->Fill();
95 
96  if (fNHits >= fMinSize) { ++fNClusters; }
97  }
98  fEventTree->Fill();
99 }
100 
102 {
103  mf::LogVerbatim("ClusterCounter2") << "ClusterCounter2 finished job";
104 }
105 
106 } // tutorial namespace
107 
void analyze(art::Event const &e) override
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
std::string string
Definition: nybbler.cc:12
EDAnalyzer(fhicl::ParameterSet const &pset)
Definition: EDAnalyzer.h:25
ClusterCounter2(fhicl::ParameterSet const &p)
const double e
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:67
T get(std::string const &key) const
Definition: ParameterSet.h:271
ValidHandle< PROD > getValidHandle(InputTag const &tag) const
Definition: DataViewImpl.h:441
p
Definition: test.py:223
ClusterCounter2 & operator=(ClusterCounter2 const &)=delete
EventNumber_t event() const
Definition: EventID.h:116
TCEvent evt
Definition: DataStructs.cxx:7
EventID id() const
Definition: Event.cc:34