CbmRoot
Loading...
Searching...
No Matches
CbmPsdIdealDigitizer.cxx
Go to the documentation of this file.
1/* Copyright (C) 2012-2019 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Alla Maevskaya, Selim Seddiki, Volker Friese [committer], Evgeny Kashirin */
4
5// -------------------------------------------------------------------------
6// ----- CbmPsdIdealDigitizer source file -----
7// ----- Created 15/05/12 by Alla & SELIM & FLORIAN -----
8// -------------------------------------------------------------------------
10
11#include "CbmPsdDigi.h"
12#include "CbmPsdPoint.h"
13
14#include "FairRootManager.h"
15#include <Logger.h>
16
17#include "TClonesArray.h"
18#include "TMath.h"
19
20#include <iostream>
21#include <map>
22
23using std::cout;
24using std::endl;
25using std::map;
26using std::pair;
27
28
29// ----- Default constructor -------------------------------------------
31 : FairTask("Ideal Psd Digitizer", 1)
32 , fNDigis(0)
33 , fPointArray(NULL)
34 , fDigiArray(NULL)
35{
36 // Reset();
37}
38// -------------------------------------------------------------------------
39
40
41// ----- Destructor ----------------------------------------------------
43{
44 if (fDigiArray) {
45 fDigiArray->Delete();
46 delete fDigiArray;
47 }
48}
49// -------------------------------------------------------------------------
50
51
52// ----- Public method Init --------------------------------------------
54{
55
56 // Get RootManager
57 FairRootManager* ioman = FairRootManager::Instance();
58 if (!ioman) {
59 LOG(fatal) << "CbmPsdIdealDigitizer::Init: RootManager not instantised!"; //FLORIAN
60 return kFATAL;
61 }
62
63 // Get input array
64 fPointArray = (TClonesArray*) ioman->GetObject("PsdPoint");
65 if (!fPointArray) {
66 LOG(fatal) << "CbmPsdIdealDigitizer::Init: No PSDPoint array!"; //FLORIAN
67 return kERROR;
68 }
69
70 // Create and register output array
71 fDigiArray = new TClonesArray("CbmPsdDigi", 1000);
72 ioman->Register("PsdDigi", "PSD", fDigiArray, IsOutputBranchPersistent("PsdDigi"));
73
74 cout << "-I- CbmPsdIdealDigitizer: Intialisation successfull " << kSUCCESS << endl;
75 return kSUCCESS;
76}
77
78
79// -------------------------------------------------------------------------
80
81
82// ----- Public method Exec --------------------------------------------
83void CbmPsdIdealDigitizer::Exec(Option_t* /*opt*/)
84{
85
86 cout << " CbmPsdIdealDigitizer::Exec begin " << endl;
87 // Reset output array
88 if (!fDigiArray) Fatal("Exec", "No PsdDigi array");
89 Reset(); // SELIM: reset!!!
90
91 // Declare some variables
92 CbmPsdPoint* point = NULL;
93 Int_t modID = -1; // module ID
94 Int_t scinID = -1; // #sciillator
95
96 Double_t edep
98 [NB_PSD_MODS]; //SELIM: 49 modules, including central & corner modules (rejected in analysis/flow/eventPlane.cxx)
99 memset(edep, 0, (NB_PSD_SECT * NB_PSD_MODS) * sizeof(Double_t));
100
101 TVector3 pos; // Position vector
102 fNDigis = 0;
103
104 //for (Int_t imod=0; imod<100; imod++) //SELIM: 49 modules, including central & corner modules (rejected in analysis/flow/eventPlane.cxx)
105 for (Int_t imod = 0; imod < NB_PSD_MODS; imod++) //marina
106 {
107 for (Int_t isec = 0; isec < NB_PSD_SECT; isec++) {
108 edep[isec][imod] = 0.;
109 }
110 }
111
112 map<pair<int, int>, double> edepmap;
113
114 // Loop over PsdPoints
115 Int_t nPoints = fPointArray->GetEntriesFast();
116 cout << " nPoints " << nPoints << endl;
117
118 Int_t sec;
119
120 for (Int_t iPoint = 0; iPoint < nPoints; iPoint++) {
121 point = (CbmPsdPoint*) fPointArray->At(iPoint);
122 if (!point) continue;
123
124 // Detector ID
125 //modID = point->GetModuleID(); // !!!!! correction SELIM: scintID <-> modID !!!!!!
126 //scinID = point->GetDetectorID();
127 //edep[sec][modID] += point->GetEnergyLoss();
128
129 modID = point->GetModuleID(); //marina 1-44 (45)
130 scinID = point->GetDetectorID(); //1-60
131
132 //sec = (Int_t)((scinID-1)/6); //marina 0-9
133 sec = (Int_t)((scinID - 1) / 6) + 1; //marina 1-10
134
135 auto insert_result = edepmap.insert(std::make_pair(std::make_pair(modID, sec), point->GetEnergyLoss()));
136
137 if (!insert_result.second) { // this entry has existed before
138 (*insert_result.first).second += point->GetEnergyLoss();
139 }
140
141 //edep[sec-1][modID-1] += (Double_t) point->GetEnergyLoss();
142 //cout <<"MARINA modID,scinID " <<modID <<" " <<scinID <<" " <<sec <<endl;
143
144 //if ( sec > modID_max) modID_max = sec;
145 //if ( sec < modID_min) modID_min = sec;
146 } // Loop over MCPoints
147
148 //cout << "modID in: " << modID_min << ", " << modID_max << endl;
149
150 //for (Int_t imod=1; imod<50; imod++) //SELIM: 49 modules, including central & corner modules (rejected in analysis/flow/eventPlane.cxx)
151
152 /*
153 for (Int_t imod=0; imod<NB_PSD_MODS; imod++)//marina
154 {
155 for (Int_t isec=0; isec<NB_PSD_SECT; isec++)
156 {
157 //if (edep[isec][imod]<=0.) cout << "!! edep !! : " << edep[isec][imod] << endl;
158 if ( edep[isec][imod] <= 0. ) continue;
159 else {
160 new ((*fDigiArray)[fNDigis]) CbmPsdDigi(isec+1, imod+1, edep[isec][imod]);
161 fNDigis++;
162 //cout <<"MARINA CbmPsdIdealDigitizer " <<fNDigis <<" " <<isec+1 <<" " <<imod+1 <<" " <<edep[isec][imod] <<endl;
163 }
164 } // section
165 }//module
166 */
167
168 fNDigis = 0;
169 for (auto edep_entry : edepmap) {
170 modID = edep_entry.first.first;
171 int secID = edep_entry.first.second;
172 double edep1 = edep_entry.second;
173 new ((*fDigiArray)[fNDigis]) CbmPsdDigi(secID, modID, edep1);
174 fNDigis++;
175 }
176
177 // Event summary
178 cout << "-I- CbmPsdIdealDigitizer: " << fNDigis << " CbmPsdDigi created." << endl;
179}
180// -------------------------------------------------------------------------
181
182// ----- Private method Reset ------------------------------------------
184{
185 fNDigis = 0;
186 if (fDigiArray) fDigiArray->Delete();
187}
188
189
ClassImp(CbmConverterManager)
const Int_t NB_PSD_MODS
const Int_t NB_PSD_SECT
Data class for PSD digital information.
Definition CbmPsdDigi.h:36
virtual InitStatus Init()
virtual void Exec(Option_t *opt)
int32_t GetModuleID() const
Definition CbmPsdPoint.h:61