// ------------------------------------------------------------------------- // ----- CbmPsdSimpleDigitizer source file ----- // ----- Created 15/05/12 by Alla & SELIM & FLORIAN ----- // ----- Modified 17/03/18 by Sergey Morozov ----- // ------------------------------------------------------------------------- #include #include #include #include #include "TClonesArray.h" #include "FairRootManager.h" #include "FairLogger.h" #include "CbmPsdDigi.h" #include "CbmPsdSimpleDigitizer.h" #include "CbmPsdPoint.h" #include "TMath.h" #include "TStopwatch.h" using std::cout; using std::endl; using std::setw; using std::right; using std::fixed; using std::setprecision; using std::map; using std::pair; // ----- Default constructor ------------------------------------------- CbmPsdSimpleDigitizer::CbmPsdSimpleDigitizer() : CbmDigitize("PsdDigitizer"), fNofEvents(0), fNofPoints(0.), fNofDigis(0.), fTimeTot(0.), fPointArray(NULL), fDigiArray(NULL) { // Reset(); } // ------------------------------------------------------------------------- // ----- Destructor ---------------------------------------------------- CbmPsdSimpleDigitizer::~CbmPsdSimpleDigitizer() { if ( fDigiArray ) { fDigiArray->Delete(); delete fDigiArray; } } // ------------------------------------------------------------------------- // ----- Public method Init -------------------------------------------- InitStatus CbmPsdSimpleDigitizer::Init() { // Get RootManager FairRootManager* ioman = FairRootManager::Instance(); assert(ioman), // Get input array fPointArray = (TClonesArray*) ioman->GetObject("PsdPoint"); assert(fPointArray); // Create and register output array fDigiArray = new TClonesArray("CbmPsdDigi", 1000); ioman->Register("PsdDigi", "PSD", fDigiArray, IsOutputBranchPersistent("PsdDigi")); // Statistics fNofEvents = 0; fNofPoints = 0; fNofDigis = 0.; fTimeTot = 0.; LOG(info) << fName << ": Initialisation successful " << kSUCCESS ; return kSUCCESS; } // ------------------------------------------------------------------------- // ----- Public method Exec -------------------------------------------- void CbmPsdSimpleDigitizer::Exec(Option_t*) { TStopwatch timer; timer.Start(); LOG(debug) << fName << ": processing event " << fCurrentEvent << " at t = " << fCurrentEventTime << " ns"; // Declare some variables CbmPsdPoint* point = NULL; Int_t modID = -1; // module ID Int_t scinID = -1; // #sciillator Double_t edep[N_PSD_SECT][N_PSD_MODS]; //SELIM: 49 modules, including central & corner modules (rejected in analysis/flow/eventPlane.cxx) memset(edep, 0,(N_PSD_SECT*N_PSD_MODS)*sizeof(Double_t)); map , double> edepmap; TVector3 pos; // Position vector //for (Int_t imod=0; imod<100; imod++) //SELIM: 49 modules, including central & corner modules (rejected in analysis/flow/eventPlane.cxx) for (Int_t imod=0; imodGetEntriesFast(); Int_t sec; for (Int_t iPoint=0; iPointAt(iPoint); if ( ! point ) continue; modID = point->GetModuleID(); //marina 1-44 (45) scinID = point->GetDetectorID();//1-60 Double_t eLoss = point->GetEnergyLoss(); sec = (Int_t)((scinID-1)/6)+1; //marina 1-10 auto insert_result = edepmap.insert( std::make_pair(std::make_pair(modID, sec), point->GetEnergyLoss())); if (!insert_result.second) { // this entry has existed before (*insert_result.first).second += point->GetEnergyLoss(); } //cout <<"PSD modID,scinID,eloss " << modID << ", " << scinID << ", " << eLoss <=0 && (sec-1)=0 && (modID-1)Gaus(eLossMIP * pixPerMIP,sqrt(eLossMIP * pixPerMIP)) / pixPerMIP; Double_t eLossSmeared = eLossMIPSmeared * 0.005; Double_t eNoise = gRandom->Gaus(0,15) / 50. * 0.005; eLossSmeared += eNoise; // V.F. The digi time is set to the event time. This is a workaround only // to integrate PSD in the common digitisation scheme. CbmPsdDigi* digi = new CbmPsdDigi(secID, modID, eLossSmeared, fCurrentEventTime); SendDigi(digi); nDigis++; LOG(debug1) << fName << ": Digi " << nDigis << " Section " << secID << " Module " << modID << " energy " << eLossSmeared; } /* } }// section }//module */ // --- Event log timer.Stop(); LOG(info) << "+ " << setw(15) << GetName() << ": Event " << setw(6) << right << fCurrentEvent << " at " << fixed << setprecision(3) << fCurrentEventTime << " ns, points: " << nPoints << ", digis: " << nDigis << ". Exec time " << setprecision(6) << timer.RealTime() << " s."; // --- Run statistics fNofEvents++; fNofPoints += nPoints; fNofDigis += nDigis; fTimeTot += timer.RealTime(); } // ------------------------------------------------------------------------- // ----- End-of-run ---------------------------------------------------- void CbmPsdSimpleDigitizer::Finish() { std::cout << std::endl; LOG(info) << "====================================="; LOG(info) << GetName() << ": Run summary"; LOG(info) << "Events processed : " << fNofEvents; LOG(info) << "PsdPoint / event : " << setprecision(1) << fNofPoints / Double_t(fNofEvents); LOG(info) << "PsdDigi / event : " << fNofDigis / Double_t(fNofEvents); LOG(info) << "Digis per point : " << setprecision(6) << fNofDigis / fNofPoints; LOG(info) << "Real time per event : " << fTimeTot / Double_t(fNofEvents) << " s"; LOG(info) << "====================================="; } // ------------------------------------------------------------------------- // ----- Reset output arrays ------------------------------------------- void CbmPsdSimpleDigitizer::ResetArrays() { if ( fDigiArray ) fDigiArray->Delete(); } // ------------------------------------------------------------------------- // ----- Write data to arrays ------------------------------------------ void CbmPsdSimpleDigitizer::WriteDigi(CbmDigi* digi, CbmMatch* /*match*/) { CbmPsdDigi* myDigi = dynamic_cast(digi); assert(myDigi); Int_t index = fDigiArray->GetEntriesFast(); new ((*fDigiArray)[index]) CbmPsdDigi(*myDigi); } // ------------------------------------------------------------------------- ClassImp(CbmPsdSimpleDigitizer)