CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/SimG4Core/HelpfulWatchers/src/G4StepStatistics.h

Go to the documentation of this file.
00001 #ifndef HelpfulWatchers_G4StepStatistics_h
00002 #define HelpfulWatchers_G4StepStatistics_h
00003 // -*- C++ -*-
00004 //
00005 // Package:     HelpfulWatchers
00006 // Class  :     SimTracer
00007 // 
00016 //
00017 // Original Author:  
00018 //         Created:  Tue Nov 22 16:41:33 EST 2005
00019 // $Id: G4StepStatistics.h,v 1.7 2010/01/12 07:14:48 hegner Exp $
00020 //
00021 
00022 // system include files
00023 #include <iostream>
00024 
00025 // user include files
00026 #include "SimG4Core/Notification/interface/Observer.h"
00027 #include "SimG4Core/Watcher/interface/SimWatcher.h"
00028 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00029 #include "G4Step.hh"
00030 #include "G4VProcess.hh"
00031 #include "G4ParticleDefinition.hh"
00032 #include <map>
00033 #include <FWCore/ServiceRegistry/interface/Service.h>
00034 #include <CommonTools/UtilAlgos/interface/TFileService.h>
00035 
00036 #include <TROOT.h>
00037 #include <TTree.h>
00038 #include <TFile.h>
00039 #include <TVector.h>
00040 #include <TString.h>
00041 #include <TClonesArray.h>
00042 //#include<TObjString.h>
00043 
00044 // forward declarations
00045 class DDDWorld;
00046 class BeginOfJob;
00047 class BeginOfRun;
00048 class BeginOfEvent;
00049 class BeginOfTrack;
00050 class G4Step;
00051 
00052 class EndOfRun;
00053 class EndOfEvent;
00054 class EndOfTrack;
00055 
00056 //Define a class StepID
00057 class StepID {
00058 
00059  private:
00060   //G4 Region
00061   G4String theG4RegionName;
00062   //G4 Physical Process
00063   G4String theG4ProcessName;
00064   //Particle PDG ID
00065   G4int theParticlePDGID;
00066   
00067  public:
00068   //Constructor using G4Step
00069   StepID(const G4Step* theG4Step)
00070     :
00071     theG4RegionName("UNDEFINED"),
00072     theG4ProcessName("UNDEFINED"),
00073     theParticlePDGID(theG4Step->GetTrack()->GetDefinition()->GetPDGEncoding())
00074     {
00075       std::cout<<"Start"<<std::endl;
00076       if (theG4Step->GetPreStepPoint()->GetPhysicalVolume()) {
00077           theG4RegionName = theG4Step->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetRegion()->GetName();
00078         }
00079       std::cout<<"Physical Volume"<<std::endl;
00080       if (theG4Step->GetPreStepPoint()->GetProcessDefinedStep()){
00081           theG4ProcessName = theG4Step->GetPreStepPoint()->GetProcessDefinedStep()->GetProcessName();
00082         }
00083       std::cout<<"Process Name"<<std::endl;
00084     }
00085 
00086   //Getters
00087   G4String GetRegionName () const { return theG4RegionName; }
00088   G4String GetProcessName () const { return theG4ProcessName; }
00089   G4int GetParticlePDGID () const { return theParticlePDGID; }
00090 
00091   //Comparison Operators (necessary in order to use StepID as a key in a map)
00092   bool operator==(const StepID& id) const{
00093     return ( strcmp(theG4RegionName,id.GetRegionName())==0 && strcmp(theG4ProcessName,id.GetProcessName())==0 && theParticlePDGID==id.GetParticlePDGID() ) ? true : false;
00094   }
00095   
00096   bool operator<(const StepID& id) const
00097     {
00098       if (theParticlePDGID != id.GetParticlePDGID()){
00099         return (theParticlePDGID > id.GetParticlePDGID());
00100       }
00101       else if (strcmp(theG4RegionName,id.GetRegionName())!=0){
00102         return strcmp(theG4RegionName,id.GetRegionName())>0 ? true : false;
00103       }
00104       else if (strcmp(theG4ProcessName,id.GetProcessName())!=0){
00105         return strcmp(theG4ProcessName,id.GetProcessName())>0 ? true : false;
00106       }
00107       else {//The case in which they are all equal!
00108         return false;
00109       }
00110     }
00111 
00112   bool operator>(const StepID& id) const
00113     {
00114       if(theParticlePDGID != id.GetParticlePDGID()){
00115         return (theParticlePDGID < id.GetParticlePDGID());
00116       }
00117       else if(strcmp(theG4RegionName,id.GetRegionName())!=0){
00118         return strcmp(theG4RegionName,id.GetRegionName())<0 ? true : false;
00119       }
00120       else if (strcmp(theG4ProcessName,id.GetProcessName())!=0){
00121         return strcmp(theG4ProcessName,id.GetProcessName())<0 ? true : false;
00122       }
00123       else {//The case in which they are all equal!
00124         return false;
00125       }
00126     }
00127 };
00128 
00129 #define OBSERVES(type) public Observer<const type*>
00130 #define UPDATE(type) void update(const type*) { std::cout <<"++ signal " #type<<std::endl; }
00131 class G4StepStatistics : public SimWatcher, 
00132 OBSERVES(DDDWorld),
00133 OBSERVES(BeginOfJob),
00134 OBSERVES(BeginOfRun),
00135 OBSERVES(BeginOfEvent),
00136 OBSERVES(BeginOfTrack),
00137 OBSERVES(G4Step),
00138 OBSERVES(EndOfRun),
00139 OBSERVES(EndOfEvent),
00140 OBSERVES(EndOfTrack)
00141 {
00142    public:
00143    G4StepStatistics(const edm::ParameterSet& pSet) : 
00144      m_verbose(pSet.getUntrackedParameter<bool>("verbose",false)),
00145      Event(0)
00146      {
00147        //Adding TFile Service output
00148        G4StepTree = fs->make<TTree>("G4StepTree","G4Step Tree ");
00149        G4StepTree->Branch("Event",&Event,"Event/I");
00150        G4StepTree->Branch("PDGID",&PDGID,"PDGID[100000]/I");
00151        Region = new TClonesArray("TObjString",100000);
00152        G4StepTree->Branch("Region",&Region);
00153        Process = new TClonesArray("TObjString",100000);
00154        G4StepTree->Branch("Process",&Process);
00155        G4StepTree->Branch("G4StepFreq",&G4StepFreq,"G4StepFreq[100000]/I");
00156      }
00157 UPDATE(DDDWorld)
00158 UPDATE(BeginOfJob)
00159 UPDATE(BeginOfRun)
00160   //    void update(const BeginOfRun* iRun) {
00161   //std::cout <<"++ signal BeginOfRun " <<std::endl;
00162   //}
00163 UPDATE(BeginOfEvent)
00164 UPDATE(BeginOfTrack)
00165    void update(const G4Step* iStep) { 
00166    std::cout <<"++ signal G4Step " ;
00167    //Dump the relevant information from the G4 step in the object mysteptest
00168    StepID mysteptest(iStep);
00169    //Add the StepID to the map, or increment if it already exists:
00170    if ( G4StatsMap.find(mysteptest) == G4StatsMap.end() )
00171      {
00172        //Allocating new memory for a pointer to associate with the key mysteptest
00173        //in our map. Initializing it to 1,will be incremented working on the value of the pointer.
00174        unsigned int* MyValue = new unsigned int(1);
00175        //Inserting the new key,value pair
00176        G4StatsMap.insert(std::make_pair(mysteptest, MyValue));
00177      }
00178    else
00179      {
00180        //Incrementing the value of the pointer by 1
00181        *G4StatsMap[mysteptest] = *G4StatsMap[mysteptest] + 1;
00182         }
00183    
00184    //If the verbose flag is set, then dump the information
00185    if (m_verbose) {
00186      std::cout << " StepID RegionName: "<< mysteptest.GetRegionName();
00187      std::cout << " StepID ProcessName: "<< mysteptest.GetProcessName();
00188      std::cout << " StepID ParticlePDGID: "<< mysteptest.GetParticlePDGID();
00189    }
00190    std::cout<<std::endl;
00191 }
00192 //UPDATE(G4Step)
00193 UPDATE(EndOfRun)
00194   //  void update(const EndOfRun* iRun) {
00195   //std::cout <<"++ signal EndOfRun " <<std::endl;
00196   //}
00197   
00198   //UPDATE(EndOfEvent)
00199   void update(const EndOfEvent* iRun) {
00200   std::cout <<"++ signal EndOfEvent " <<std::endl;
00201   Event++;
00202   
00203   //Dumping the map in the log if verbose is chosen:
00204   if(m_verbose) {
00205     std::cout <<" G4StatsMap size is: "<<G4StatsMap.size()<<std::endl;
00206   }
00207   int index(0);
00208   for (std::map<const StepID,unsigned int*>::const_iterator step = G4StatsMap.begin(); step != G4StatsMap.end(); ++step, ++index){
00209     if(m_verbose) {
00210       std::cout <<" G4StatsMap step is: "<<step->first.GetRegionName()<<" "<<step->first.GetProcessName()<<" "<<step->first.GetParticlePDGID();//<<" "<<step->first.GetTrackID() ;
00211       std::cout <<" Number of such steps: "<< *step->second <<std::endl;
00212     }
00213     //Rolling the map into 5 "arrays", containing the StepID information and the G4Step statistics
00214     PDGID[index]=step->first.GetParticlePDGID();
00215     new ((*Region)[index]) TObjString (step->first.GetRegionName());
00216     new ((*Process)[index]) TObjString (step->first.GetProcessName());
00217     G4StepFreq[index]=*step->second;
00218   }
00219   
00220   G4StepTree->Fill();
00221 }
00222 UPDATE(EndOfTrack)
00223 
00224   private:
00225  
00226  bool m_verbose;
00227  
00228 //Adding the G4StatsMap to keep track of statistics in terms of step information... 
00229  std::map<const StepID,unsigned int*> G4StatsMap;
00230  edm::Service<TFileService> fs;
00231  TTree* G4StepTree;
00232  unsigned int Event;
00233  Int_t PDGID[100000];
00234  TClonesArray* Region;
00235  TClonesArray* Process;
00236  Int_t G4StepFreq[100000];
00237 };
00238 
00239 #endif