All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
CheckSimEnergyDeposit_module.cc
Go to the documentation of this file.
1 // art Framework includes.
2 #include "art/Framework/Core/EDAnalyzer.h"
3 #include "art/Framework/Principal/Event.h"
4 #include "art/Framework/Principal/Handle.h"
5 #include "art/Framework/Core/ModuleMacros.h"
6 #include "art_root_io/TFileService.h"
7 #include "art_root_io/TFileDirectory.h"
8 #include "art/Framework/Services/Registry/ServiceHandle.h"
9 #include "artg4tk/services/DetectorHolder_service.hh"
11 
12 // Root includes.
13 #include "TH1F.h"
14 #include "TNtuple.h"
15 
16 // STL includes.
17 #include <cmath>
18 
19 // Other includes.
20 #include "CLHEP/Units/SystemOfUnits.h"
21 
22 using namespace std;
23 namespace larg4 {
25 }
26 
27 class larg4::CheckSimEnergyDeposit : public art::EDAnalyzer {
28 public:
29 
30  explicit CheckSimEnergyDeposit(fhicl::ParameterSet const& p);
31 
32 private:
33  void beginJob() override;
34  void analyze(const art::Event& event) override;
35 
36  TH1F* _hnHits{nullptr}; // number of SimEnergyDepositHits
37  TH1F* _hEdep{nullptr}; // average energy deposition in SimEnergyDepositHits
38  TH1F* _hnumPhotons{nullptr}; // number of Photons per SimEnergyDepositHits
39  TH1F* _hLandauPhotons{nullptr}; // Edep/cm SimEnergyDepositHits
40  TH1F* _hLandauEdep{nullptr}; // number of Photons/cm SimEnergyDepositHits
41  TH1F* _hSteplength{nullptr}; // Geant 4 step length
42  TNtuple* _ntuple{nullptr};
43 };
44 
46  art::EDAnalyzer(p)
47 {}
48 
50 {
51  art::ServiceHandle<art::TFileService const> tfs;
52  _hnHits = tfs->make<TH1F>("hnHits", "Number of SimEnergyDeposits", 300, 0, 0);
53  _hEdep = tfs->make<TH1F>("hEdep", "Energy deposition in SimEnergyDeposits", 100,0.,0.02);
54  _hnumPhotons = tfs->make<TH1F>("hnumPhotons", "number of photons per SimEnergyDeposit", 100,0.,500.);
55  _hLandauPhotons= tfs->make<TH1F>("hLandauPhotons", "number of photons/cm", 100,0.,2000000.);
56  _hLandauEdep= tfs->make<TH1F>("hLandauEdep", "Edep/cm", 100,0.,10.);
57  _hSteplength= tfs->make<TH1F>("hSteplength", "geant 4 step length", 100,0.,0.05);
58  _ntuple = tfs->make<TNtuple>("ntuple","Demo ntuple",
59  "Event:Edep:em_Edep:nonem_Edep:xpos:ypos:zpos:time");
60 } // end beginJob
61 
62 void larg4::CheckSimEnergyDeposit::analyze(const art::Event& event)
63 {
64  auto allSims = event.getMany<sim::SimEnergyDepositCollection>();
65  for (auto const& sims : allSims) {
66  double sumPhotons=0.0;
67  double sumE = 0.0;
68  _hnHits->Fill(sims->size());
69  for (auto const& hit : *sims) {
70  // sum up energy deposit in a 1cm slice of liquid Argon.
71  if (std::abs(hit.EndZ())<0.5) {
72  sumPhotons= sumPhotons + hit.NumPhotons();
73  sumE= sumE +hit.Energy();
74  }
75  _hnumPhotons->Fill( hit.NumPhotons());
76  _hEdep->Fill( hit.Energy()); // energy deposit in MeV
77  _hSteplength->Fill( hit.StepLength()); // step length in cm
78  }
79  _hLandauPhotons->Fill(sumPhotons);
80  _hLandauEdep->Fill(sumE);
81  }
82 } // end analyze
83 
84 DEFINE_ART_MODULE(larg4::CheckSimEnergyDeposit)
process_name opflashCryo1 flashfilter analyze
double std(const std::vector< short > &wf, const double ped_mean, size_t start, size_t nsample)
Definition: UtilFunc.cxx:42
pdgs p
Definition: selectors.fcl:22
process_name hit
Definition: cheaterreco.fcl:51
CheckSimEnergyDeposit(fhicl::ParameterSet const &p)
T abs(T value)
standard_singlep gaussian distribution ie ie ie gaussian gaussian gaussian simWire CheckSimEnergyDeposit
Definition: zengen.fcl:23
std::vector< SimEnergyDeposit > SimEnergyDepositCollection
contains information for a single step in the detector simulation
art::ServiceHandle< art::TFileService > tfs
void analyze(const art::Event &event) override