LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
CheckCalorimeterHits_module.cc
Go to the documentation of this file.
1 //
2 // __ __ __ __ __
3 // ____ ______/ /_____ _/ // / / /_/ /__
4 // / __ `/ ___/ __/ __ `/ // /_/ __/ //_/
5 // / /_/ / / / /_/ /_/ /__ __/ /_/ ,<
6 // \__,_/_/ \__/\__, / /_/ \__/_/|_|
7 // /____/
8 //
9 // artg4tk: art based Geant 4 Toolkit
10 //
11 //=============================================================================
12 // CheckCalorimeterHits_module.cc: Analyzer module that demonstrates access to
13 // Calorimeter hits
14 // and makes some histograms
15 //
16 // Author: Hans Wenzel (Fermilab)
17 //=============================================================================
18 
25 #include "art_root_io/TFileService.h"
26 #include "artg4tk/pluginDetectors/gdml/CalorimeterHit.hh"
27 
28 #include "CLHEP/Units/SystemOfUnits.h"
29 
30 #include "TH1F.h"
31 #include "TNtuple.h"
32 
33 using namespace std;
34 namespace artg4tk {
35  class CheckCalorimeterHits;
36 }
37 
39 public:
40  explicit CheckCalorimeterHits(fhicl::ParameterSet const& p);
41  void beginJob() override;
42  void analyze(const art::Event& event) override;
43 
44 private:
45  bool _DumpGDML; // enable/disable dumping of GDML geometry information
46  TH1F* _hnHits; // number of CaloHits
47  TH1F* _hEdep; // total energy deposition in CaloHits
48  TH1F* _haEdep; // average energy deposition in CaloHits
49  TNtuple* _ntuple;
50 };
51 
53  : art::EDAnalyzer(p)
54  , _DumpGDML(p.get<bool>("DumpGDML"))
55  , _hnHits(0)
56  , _hEdep(0)
57  , _haEdep(0)
58  , _ntuple(0)
59 {}
60 
61 void
63 {
65  _hnHits = tfs->make<TH1F>("hnHits", "Number of CaloArtHits", 300, 0, 0);
66  _hEdep = tfs->make<TH1F>("hEdep", "total Energy deposition in CaloArtHits", 2000, 0, 0);
67  _haEdep = tfs->make<TH1F>("haEdep", "z of Energy deposition in CaloArtHits", 200, -500., 500.);
68  _ntuple = tfs->make<TNtuple>(
69  "ntuple", "Demo ntuple", "Event:ID:Edep:em_Edep:nonem_Edep:xpos:ypos:zpos:time");
70 } // end beginJob
71 
72 void
74 {
75  typedef std::vector<art::Handle<CalorimeterHitCollection>> HandleVector;
76  auto allSims = event.getMany<CalorimeterHitCollection>();
77  double sumE = 0.0;
78  for (HandleVector::const_iterator i = allSims.begin(); i != allSims.end(); ++i) {
79  const CalorimeterHitCollection& sims(**i);
80  sumE = 0.0;
81  _hnHits->Fill(sims.size());
82  for (CalorimeterHitCollection::const_iterator j = sims.begin(); j != sims.end(); ++j) {
83  const CalorimeterHit& hit = *j;
84  sumE = sumE + hit.GetEdep();
85  _haEdep->Fill(hit.GetZpos(), hit.GetEdep());
86  _ntuple->Fill(event.event(),
87  hit.GetID(),
88  hit.GetEdep(),
89  hit.GetEdepEM(),
90  hit.GetEdepnonEM(),
91  hit.GetXpos(),
92  hit.GetYpos(),
93  hit.GetZpos(),
94  hit.GetTime());
95  }
96  _hEdep->Fill(sumE / CLHEP::GeV);
97  }
98 } // end analyze
99 
STL namespace.
intermediate_table::const_iterator const_iterator
CheckCalorimeterHits(fhicl::ParameterSet const &p)
void analyze(const art::Event &event) override
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
void beginJob()
Definition: Breakpoints.cc:14
EventNumber_t event() const
Definition: Event.cc:41
Detector simulation of raw signals on wires.
decltype(auto) get(T &&obj)
ADL-aware version of std::to_string.
Definition: StdUtils.h:120
Definition: MVAAlg.h:12
Event finding and building.