Go to the documentation of this file.00001 #include "SimG4Core/CheckSecondary/interface/CheckSecondary.h"
00002
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "FWCore/Framework/interface/EventSetup.h"
00005
00006 #include "SimG4Core/Notification/interface/BeginOfEvent.h"
00007 #include "SimG4Core/Notification/interface/EndOfEvent.h"
00008 #include "SimG4Core/Notification/interface/BeginOfTrack.h"
00009
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011
00012 #include "G4Step.hh"
00013 #include "G4Track.hh"
00014 #include "G4VProcess.hh"
00015 #include "G4HCofThisEvent.hh"
00016 #include "CLHEP/Units/GlobalSystemOfUnits.h"
00017 #include "CLHEP/Units/GlobalPhysicalConstants.h"
00018
00019 #include <cmath>
00020 #include <iostream>
00021 #include <iomanip>
00022
00023 CheckSecondary::CheckSecondary(const edm::ParameterSet &p): treatSecondary(0),
00024 typeEnumerator(0),
00025 nsec(0),procids(0),
00026 px(0),py(0),pz(0),
00027 mass(0),deltae(0),
00028 procs(0),file(0),
00029 tree(0) {
00030
00031 edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("CheckSecondary");
00032 std::string saveFile = m_p.getUntrackedParameter<std::string>("SaveInFile", "None");
00033 treatSecondary = new TreatSecondary(m_p);
00034 typeEnumerator = new G4ProcessTypeEnumerator();
00035
00036 nsec = new std::vector<int>();
00037 px = new std::vector<double>();
00038 py = new std::vector<double>();
00039 pz = new std::vector<double>();
00040 mass = new std::vector<double>();
00041 deltae = new std::vector<double>();
00042 procids = new std::vector<int>();
00043 procs = new std::vector<std::string>();
00044
00045 if (saveFile != "None") {
00046 saveToTree = true;
00047 tree = bookTree (saveFile);
00048 edm::LogInfo("CheckSecondary") << "Instantiate CheckSecondary with first"
00049 << " hadronic interaction information"
00050 << " to be saved in file " << saveFile;
00051 } else {
00052 saveToTree = false;
00053 edm::LogInfo("CheckSecondary") << "Instantiate CheckSecondary with first"
00054 << " hadronic interaction information"
00055 << " not saved";
00056 }
00057 }
00058
00059 CheckSecondary::~CheckSecondary() {
00060 if (saveToTree) endTree();
00061 if (nsec) delete nsec;
00062 if (px) delete px;
00063 if (py) delete py;
00064 if (pz) delete pz;
00065 if (mass) delete mass;
00066 if (deltae) delete deltae;
00067 if (procs) delete procs;
00068 if (procids) delete procids;
00069 if (typeEnumerator) delete typeEnumerator;
00070 if (treatSecondary) delete treatSecondary;
00071 }
00072
00073 TTree* CheckSecondary::bookTree(std::string fileName) {
00074
00075 file = new TFile (fileName.c_str(), "RECREATE");
00076 file->cd();
00077
00078 TTree * t1 = new TTree("T1", "Secondary Particle Information");
00079 t1->Branch("SecondaryPx", "std::vector<double>", &px);
00080 t1->Branch("SecondaryPy", "std::vector<double>", &py);
00081 t1->Branch("SecondaryPz", "std::vector<double>", &pz);
00082 t1->Branch("SecondaryMass", "std::vector<double>", &mass);
00083 t1->Branch("NumberSecondaries", "std::vector<int>", &nsec);
00084 t1->Branch("DeltaEinInteract", "std::vector<double>", &deltae);
00085 t1->Branch("ProcessID", "std::vector<int>", &procids);
00086 t1->Branch("ProcessNames", "std::vector<std::string>", &procs);
00087 return t1;
00088 }
00089
00090 void CheckSecondary::endTree() {
00091
00092 edm::LogInfo("CheckSecondary") << "Save the Secondary Tree "
00093 << tree->GetName() << " (" << tree
00094 << ") in file " << file->GetName() << " ("
00095 << file << ")";
00096 file->cd();
00097 tree->Write();
00098 file->Close();
00099 delete file;
00100 }
00101
00102 void CheckSecondary::update(const BeginOfEvent * evt) {
00103
00104 int iev = (*evt)()->GetEventID();
00105 LogDebug("CheckSecondary") << "CheckSecondary::=====> Begin of event = "
00106 << iev;
00107
00108 (*nsec).clear();
00109 (*procs).clear();
00110 (*procids).clear();
00111 (*deltae).clear();
00112 (*px).clear();
00113 (*py).clear();
00114 (*pz).clear();
00115 (*mass).clear();
00116 }
00117
00118 void CheckSecondary::update(const BeginOfTrack * trk) {
00119
00120 const G4Track * thTk = (*trk)();
00121 treatSecondary->initTrack(thTk);
00122 if (thTk->GetParentID() <= 0) storeIt = true;
00123 else storeIt = false;
00124 nHad = 0;
00125 LogDebug("CheckSecondary") << "CheckSecondary:: Track " << thTk->GetTrackID()
00126 << " Parent " << thTk->GetParentID() << " Flag "
00127 << storeIt;
00128 }
00129
00130 void CheckSecondary::update(const G4Step * aStep) {
00131
00132 std::string name;
00133 int procID;
00134 bool hadrInt;
00135 double deltaE;
00136 std::vector<int> charge;
00137 std::vector<math::XYZTLorentzVector> tracks = treatSecondary->tracks(aStep,
00138 name,
00139 procID,
00140 hadrInt,
00141 deltaE,
00142 charge);
00143 if (storeIt && hadrInt) {
00144 double pInit = (aStep->GetPreStepPoint()->GetMomentum()).mag();
00145 double pEnd = (aStep->GetPostStepPoint()->GetMomentum()).mag();
00146 nHad++;
00147 int sec = (int)(tracks.size());
00148 LogDebug("CheckSecondary") << "CheckSecondary:: Hadronic Interaction "
00149 << nHad << " of type " << name << " ID "
00150 << procID << " with " << sec << " secondaries "
00151 << " and Momentum (MeV/c) at start " << pInit
00152 << " and at end " << pEnd;
00153 (*nsec).push_back(sec);
00154 (*procs).push_back(name);
00155 (*procids).push_back(procID);
00156 (*deltae).push_back(deltaE);
00157 if (nHad == 1) {
00158 for (int i=0; i<sec; i++) {
00159 double m = tracks[i].M();
00160 if (charge[i]<0) m = -m;
00161 (*px).push_back(tracks[i].Px());
00162 (*py).push_back(tracks[i].Py());
00163 (*pz).push_back(tracks[i].Pz());
00164 (*mass).push_back(m);
00165 }
00166 }
00167 }
00168 }
00169
00170 void CheckSecondary::update(const EndOfEvent * evt) {
00171
00172 LogDebug("CheckSecondary") << "CheckSecondary::EndofEvent =====> Event "
00173 << (*evt)()->GetEventID() << " with "
00174 << (*nsec).size() << " hadronic collisions with"
00175 << " secondaries produced in each step";
00176 for (unsigned int i= 0; i < (*nsec).size(); i++)
00177 LogDebug("CheckSecondary") << " " << (*nsec)[i] << " from " << (*procs)[i]
00178 << " ID " << (*procids)[i] << " ("
00179 << typeEnumerator->processG4Name((*procids)[i])
00180 << ") deltaE = " << (*deltae)[i] << " MeV";
00181 LogDebug("CheckSecondary") << "And " << (*mass).size() << " secondaries "
00182 << "produced in the first interactions";
00183 for (unsigned int i= 0; i < (*mass).size(); i++)
00184 LogDebug("CheckSecondary") << "Secondary " << i << " (" << (*px)[i] << ", "
00185 << (*py)[i] << ", " << (*pz)[i] << ", "
00186 << (*mass)[i] << ")";
00187
00188 if (saveToTree) tree->Fill();
00189 }