File indexing completed on 2025-01-30 09:19:32
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045
0046
0047
0048
0049
0050 #include "TFile.h"
0051 #include "TSystem.h"
0052 #include "TTree.h"
0053 #include "TH1.h"
0054 #include "TH2.h"
0055
0056 #include "Event.hh"
0057 #include "lArTPCHit.hh"
0058 #include "PhotonHit.hh"
0059
0060 int main(int argc, char** argv)
0061 {
0062
0063 TSystem ts;
0064 gSystem->Load("libCaTSClassesDict");
0065 if(argc < 4)
0066 {
0067 G4cout << "Program requires 3 arguments: name of input file, name of "
0068 "output file, Volume that sensitive detector is attached to"
0069 << G4endl;
0070 exit(1);
0071 }
0072 TFile* outfile = new TFile(argv[2], "RECREATE");
0073 outfile->cd();
0074 TH2F* pos2 = new TH2F("position", "position of Photon Hits", 400, -1000.,
0075 1000., 400, -500, 500);
0076 TH1F* time = new TH1F("time", "timing of photon hits", 1000, 0., 250.);
0077 TH1F* time0 = new TH1F("time0", "timing of photon hits", 1000, 0., 250.);
0078 TH1F* time1 = new TH1F("time1", "timing of photon hits", 1000, 0., 250.);
0079 TH1F* time2 = new TH1F("time2", "timing of photon hits", 1000, 0., 250.);
0080 TH1F* wl = new TH1F("wl", "wavelength of detected photons", 1000, 0., 250.);
0081 TH1F* np = new TH1F("np", "number of detected photons", 100, 0., 50.);
0082 TFile fo(argv[1]);
0083 fo.GetListOfKeys()->Print();
0084 Event* event = new Event();
0085 TTree* Tevt = (TTree*) fo.Get("Events");
0086 Tevt->SetBranchAddress("event.", &event);
0087 TBranch* fevtbranch = Tevt->GetBranch("event.");
0088 Int_t nevent = fevtbranch->GetEntries();
0089 G4cout << "Nr. of Events: " << nevent << G4endl;
0090 std::string CollectionName = argv[3];
0091 CollectionName = CollectionName + "_Photondetector_HC";
0092 for(Int_t i = 0; i < nevent; i++)
0093 {
0094 fevtbranch->GetEntry(i);
0095 auto* hcmap = event->GetHCMap();
0096 for(const auto& ele : *hcmap)
0097 {
0098 auto hits = ele.second;
0099 if(ele.first.compare(CollectionName) == 0)
0100 {
0101 auto hits = ele.second;
0102 G4int NbHits = hits.size();
0103 G4cout << "Event: " << i << " Number of Hits: " << NbHits << G4endl;
0104 np->Fill(NbHits);
0105 for(G4int ii = 0; ii < NbHits; ii++)
0106 {
0107 PhotonHit* photonHit = dynamic_cast<PhotonHit*>(hits.at(ii));
0108 time->Fill(photonHit->GetTime());
0109 wl->Fill(photonHit->GetWavelength());
0110 if(photonHit->GetPosition().getZ() < -100.)
0111 time0->Fill(photonHit->GetTime());
0112 if(photonHit->GetPosition().getZ() > -100. &&
0113 photonHit->GetPosition().getZ() < 100)
0114 time1->Fill(photonHit->GetTime());
0115 if(photonHit->GetPosition().getZ() < 100.)
0116 time2->Fill(photonHit->GetTime());
0117 pos2->Fill(photonHit->GetPosition().getZ(),
0118 photonHit->GetPosition().getY());
0119 }
0120 }
0121 }
0122 }
0123 outfile->cd();
0124 outfile->Write();
0125 }