/scratch0/jsnuveri/BDSIM/BDSIMgit/bdsim/src/BDSRunAction.cc

00001 #include "BDSAcceleratorModel.hh"
00002 #include "BDSBeamline.hh"
00003 #include "BDSAnalysisManager.hh"
00004 #include "BDSDebug.hh"
00005 #include "BDSExecOptions.hh"
00006 #include "BDSGlobalConstants.hh" 
00007 #include "BDSOutputBase.hh" 
00008 #include "BDSRunAction.hh"
00009 
00010 #include "G4Run.hh"
00011 
00012 #include "globals.hh"               // geant4 globals / types
00013 
00014 extern BDSOutputBase* bdsOutput;         // output interface
00015 
00016 BDSRunAction::BDSRunAction()
00017 {}
00018 
00019 BDSRunAction::~BDSRunAction()
00020 {}
00021 
00022 void BDSRunAction::BeginOfRunAction(const G4Run* aRun)
00023 {
00024   //Get the current time
00025   starttime = time(NULL);
00026 
00027   // construct output histograms
00028   // calculate histogram dimensions
00029   G4double smin, smax, binwidth;
00030   G4int nbins;
00031   smin     = 0.0;
00032   smax     = BDSGlobalConstants::Instance()->GetSMax()/CLHEP::m;
00033   binwidth = BDSGlobalConstants::Instance()->GetElossHistoBinWidth(); // CHECK UNITS
00034   nbins    = (int) ceil((smax-smin)/binwidth); // rounding up so last bin definitely covers smax
00035   smax     = smin + (nbins*binwidth);          // redefine smax
00036   // create the histograms
00037   phitsindex = BDSAnalysisManager::Instance()->Create1DHistogram("PhitsHisto","Primary Hits",nbins,smin,smax); //0
00038   plossindex = BDSAnalysisManager::Instance()->Create1DHistogram("PlossHisto","Primary Loss",nbins,smin,smax); //1
00039   elossindex = BDSAnalysisManager::Instance()->Create1DHistogram("ElossHisto","Energy Loss", nbins,smin,smax); //2
00040   // prepare bin edges for a by-element histogram
00041   std::vector<double> binedges;
00042   binedges.push_back(0.0);
00043   BDSBeamline* beamline  = BDSAcceleratorModel::Instance()->GetFlatBeamline();
00044   BDSBeamlineIterator it = beamline->begin();
00045   for(; it != beamline->end(); ++it)
00046     {binedges.push_back((*it)->GetSPositionEnd()/CLHEP::m);}
00047   
00048   // create per element ("pe") bin width histograms
00049   phitspeindex = BDSAnalysisManager::Instance()->Create1DHistogram("PhitsPEHisto","Primary Hits per Element",binedges); //3
00050   plosspeindex = BDSAnalysisManager::Instance()->Create1DHistogram("PlossPEHisto","Primary Loss per Element",binedges); //4
00051   elosspeindex = BDSAnalysisManager::Instance()->Create1DHistogram("ElossPEHisto","Energy Loss per Element" ,binedges); //5
00052   
00053   //Output feedback
00054   G4cout << __METHOD_NAME__ << " Run " << aRun->GetRunID() << " start. Time is " << asctime(localtime(&starttime)) << G4endl;
00055 
00056 }
00057 
00058 void BDSRunAction::EndOfRunAction(const G4Run* aRun)
00059 {
00060   //Get the current time
00061   stoptime = time(NULL);
00062 
00063   //Output feedback
00064   G4cout << __METHOD_NAME__ << "Run " << aRun->GetRunID() << " end. Time is " << asctime(localtime(&stoptime)) << G4endl;
00065   
00066   // Write output
00067   // write histograms to output - do this before potentially closing / opening new files
00068   bdsOutput->WriteHistogram(BDSAnalysisManager::Instance()->GetHistogram(0));
00069   bdsOutput->WriteHistogram(BDSAnalysisManager::Instance()->GetHistogram(1));
00070   bdsOutput->WriteHistogram(BDSAnalysisManager::Instance()->GetHistogram(2));
00071   bdsOutput->WriteHistogram(BDSAnalysisManager::Instance()->GetHistogram(3));
00072   bdsOutput->WriteHistogram(BDSAnalysisManager::Instance()->GetHistogram(4));
00073   bdsOutput->WriteHistogram(BDSAnalysisManager::Instance()->GetHistogram(5));
00074   if(BDSExecOptions::Instance()->GetBatch()) {  // Non-interactive mode
00075     bdsOutput->Write(); // write last file
00076   } else {
00077     bdsOutput->Commit(); // write and open new file
00078   }
00079 
00080   // delete analysis manager
00081   delete BDSAnalysisManager::Instance();
00082   
00083   // note difftime only calculates to the integer second
00084   G4cout << "Run Duration >> " << (int)difftime(stoptime,starttime) << " s" << G4endl;
00085 }

Generated on 28 Jun 2015 for BDSIM by  doxygen 1.4.7