CGEM BOSS 6.6.5.h
BESIII Offline Software System
Loading...
Searching...
No Matches
BesTuningIO Class Reference

#include <BesTuningIO.hh>

Public Member Functions

 BesTuningIO (std::vector< std::string >)
 
 ~BesTuningIO ()
 
void GetNextEvents (void)
 
void GetMdcHits (void)
 
void GetCgemHits (void)
 
void GetTofHits (void)
 
void GetEmcDigi (void)
 
void GetMucHits (void)
 
void GetRootEvent (int evtID)
 
void GetMdcRootHits ()
 
void GetCgemRootHits ()
 
void GetTofRootHits ()
 
void GetEmcRootDigi ()
 

Public Attributes

TFile * f
 
TTree * HitTree
 
TChain * HitChain
 

Detailed Description

Definition at line 26 of file BesTuningIO.hh.

Constructor & Destructor Documentation

◆ BesTuningIO()

BesTuningIO::BesTuningIO ( std::vector< std::string > name)

Definition at line 41 of file BesTuningIO.cc.

42 :m_tuningFile(name),m_evt(0)
43{
44 m_DigiMan = G4DigiManager::GetDMpointer();
45 m_inputFileStream = new std::ifstream();
46
48 //tuning input root file
49 //TFile *f = new TFile(m_tuningFile);
50 //HitTree = (TTree*)f->Get("HitTree");
51 //m_TMcHitEvent = new TMcHitEvent();
52 //TBranch *branch = HitTree->GetBranch("TMcHitEvent");
53 // begin change from TTree to TChain
54 HitChain = new TChain("HitTree");
55 if (m_tuningFile.size()==0){
56 std::cout << "there is no tuning file" << std::endl;
57 }
58 std::cout << "file number: " << m_tuningFile.size() << std::endl;
59 for (int i = 0 ; i < m_tuningFile.size(); i++){
60
61 //std::cout << "________________________________________________________________________________fileName : " << m_tuningFile[i] << std::endl;
62 //HitChain->Add(&m_tuningFile[i]);
63 HitChain->Add(m_tuningFile[i].c_str());
64 }
65 m_TMcHitEvent = new TMcHitEvent();
66 TBranch *branch = HitChain->GetBranch("TMcHitEvent");
67
68 // end change from TTree to TChain
69 branch->SetAddress(&m_TMcHitEvent);
70 std::cout << "HitChain entries: " << HitChain->GetEntries() << std::endl;
71 }
72 else{// tuning input ascii file
73 //m_inputFileStream->open(m_tuningFile.c_str());
74 //if ((*m_inputFileStream).good()) {
75 // try {
76 // (*m_inputFileStream) >> m_version;
77 // } catch (AsciiDumpException& ) {
78 // std::cerr << "BesTuningIO::Got AsciiDumpException eror while reading VERSION block !!!" << std::endl;
79 // }
80 //}else{
81 // std::cerr << "BesTuningIO::Open tuning input file error!!!" << std::endl;
82 //}
83 }
84}
TChain * HitChain
static G4int GetFormatAR()
char * c_str(Index i)

◆ ~BesTuningIO()

BesTuningIO::~BesTuningIO ( )

Definition at line 86 of file BesTuningIO.cc.

86 {
87 if (m_inputFileStream) delete m_inputFileStream;
88 if (m_evt) delete m_evt;
89}

Member Function Documentation

◆ GetCgemHits()

void BesTuningIO::GetCgemHits ( void )

Definition at line 163 of file BesTuningIO.cc.

164{
165 G4int cgemHitCollID = -1;
166 cgemHitCollID = m_DigiMan->GetHitsCollectionID("BesCgemHitsCollection");
167 if (cgemHitCollID>=0)
168 {
169 BesCgemHitsCollection* cgemDC = (BesCgemHitsCollection*)m_DigiMan->GetHitsCollection(cgemHitCollID);
170 if (cgemDC)
171 {
172 G4int nHit = cgemDC->entries();
173 if (nHit>0)
174 {
175 for (G4int i=0;i<nHit;i++)
176 {
177 delete (*cgemDC)[i];
178 }
179 cgemDC->GetVector()->clear();
180 }
181
182 std::vector<CgemHitType>::iterator iter;
183 iter = (m_evt->cgemHit).hitCol.begin();
184 // Loop over cgem hits
185 for (; iter != (m_evt->cgemHit).hitCol.end(); iter++)
186 {
187 BesCgemHit* newHit = new BesCgemHit();
188
189 newHit->SetTrackID ((*iter).m_ID_track );
190 newHit->SetLayerID ((*iter).m_ID_layer );
191 newHit->SetPDGCode ((*iter).m_pdg_code );
192 newHit->SetGlobalTime ((*iter).m_global_time);
193 newHit->SetTotalEnergyDeposit ((*iter).m_E_deposit );
194 newHit->SetStepLength ((*iter).m_L_step );
195 newHit->SetPositionOfPrePoint (G4ThreeVector((*iter).m_XYZ_pre_x,(*iter).m_XYZ_pre_y,(*iter).m_XYZ_pre_z));
196 newHit->SetPositionOfPostPoint(G4ThreeVector((*iter).m_XYZ_post_x,(*iter).m_XYZ_post_y,(*iter).m_XYZ_post_z));
197 newHit->SetMomentumOfPrePoint (G4ThreeVector((*iter).m_P_pre_x,(*iter).m_P_pre_y,(*iter).m_P_pre_z));
198 newHit->SetMomentumOfPostPoint(G4ThreeVector((*iter).m_P_post_x,(*iter).m_P_post_y,(*iter).m_P_post_z));
199 /*
200 newHit->SetPrePositionInCu (G4ThreeVector((*iter).m_Cu_pre_x,(*iter).m_Cu_pre_y,(*iter).m_Cu_pre_z));
201 newHit->SetPostPositionInCu (G4ThreeVector((*iter).m_Cu_post_x,(*iter).m_Cu_post_y,(*iter).m_Cu_post_z));
202 newHit->SetMomentumOfCuPre (G4ThreeVector((*iter).m_P_Cu_pre_x,(*iter).m_P_Cu_pre_y,(*iter).m_P_Cu_pre_z));
203 newHit->SetMomentumOfCuPost (G4ThreeVector((*iter).m_P_Cu_post_x,(*iter).m_P_Cu_post_y,(*iter).m_P_Cu_post_z));
204 */
205//****************************************************************************************************
206 cgemDC->insert(newHit);
207 }
208 //cgemDC->PrintAllHits();
209 }
210 else
211 {
212 std::cerr << "BesTuningIO::can't get cgemHitsCollection"<<std::endl;
213 }
214 }
215 else
216 {
217 std::cerr << "BesTuningIO::can't get cgemHitCollID"<<std::endl;
218 }
219}
G4THitsCollection< BesCgemHit > BesCgemHitsCollection
EvtStreamInputIterator< typename Generator::result_type > iter(Generator gen, int N=0)
void SetTotalEnergyDeposit(G4double f_E_deposit)
void SetGlobalTime(G4double f_global_time)
void SetMomentumOfPostPoint(G4ThreeVector f_P_post)
void SetPDGCode(G4int f_pdg_code)
void SetLayerID(G4int f_ID_layer)
void SetTrackID(G4int f_ID_track)
void SetPositionOfPrePoint(G4ThreeVector f_XYZ_pre)
void SetMomentumOfPrePoint(G4ThreeVector f_P_pre)
void SetPositionOfPostPoint(G4ThreeVector f_XYZ_post)
void SetStepLength(G4double f_L_step)
CGEMHIT cgemHit
Definition AsciiData.hh:547

◆ GetCgemRootHits()

void BesTuningIO::GetCgemRootHits ( )

Definition at line 343 of file BesTuningIO.cc.

344{
345
346 G4int THCID = -1;
347 THCID = m_DigiMan->GetHitsCollectionID("BesCgemHitsCollection");
348 if (THCID>=0)
349 {
350 BesCgemHitsCollection* cgemDC =(BesCgemHitsCollection*)(m_DigiMan->GetHitsCollection(THCID));
351 if(cgemDC)
352 {
353 int nHit = cgemDC->entries();
354 if(nHit > 0)
355 {
356 for(int i = 0; i < nHit; i++)
357 {
358 delete (*cgemDC)[i];
359 }
360 cgemDC->GetVector()->clear();
361 }
362 }
363
364 int nHits = m_TMcHitEvent->getMcHitCgemCol()->GetEntries();
365 //std::cout << "nHits: " << nHits << std::endl;
366 for(int i = 0; i < nHits; i++)
367 {
368 m_TMcHitCgem = m_TMcHitEvent->getMcHitCgem(i);
369
370 BesCgemHit* cgemHit = new BesCgemHit();
371 cgemHit->SetTrackID (m_TMcHitCgem->GetTrackID ());
372 cgemHit->SetLayerID (m_TMcHitCgem->GetLayerID ());
373 cgemHit->SetPDGCode (m_TMcHitCgem->GetPDGCode ());
374 cgemHit->SetParentID (m_TMcHitCgem->GetParentID ());
375 cgemHit->SetGlobalTime (m_TMcHitCgem->GetGlobalTime ());
376 cgemHit->SetTotalEnergyDeposit (m_TMcHitCgem->GetTotalEnergyDeposit ());
377 cgemHit->SetStepLength (m_TMcHitCgem->GetStepLength ());
378
379 TVector3 XYZ_pre = m_TMcHitCgem->GetPositionOfPrePoint ();
380 TVector3 XYZ_post = m_TMcHitCgem->GetPositionOfPostPoint();
381 TVector3 P_pre = m_TMcHitCgem->GetMomentumOfPrePoint ();
382 TVector3 P_post = m_TMcHitCgem->GetMomentumOfPostPoint();
383
384 G4ThreeVector tmp_XYZ_pre = G4ThreeVector(XYZ_pre.x() , XYZ_pre.y() , XYZ_pre.z() );
385 G4ThreeVector tmp_XYZ_post = G4ThreeVector(XYZ_post.x() , XYZ_post.y() , XYZ_post.z());
386 G4ThreeVector tmp_P_pre = G4ThreeVector(P_pre.x() , P_pre.y() , P_pre.z() );
387 G4ThreeVector tmp_P_post = G4ThreeVector(P_post.x() , P_post.y() , P_post.z() );
388
389 cgemHit->SetPositionOfPrePoint ( tmp_XYZ_pre );
390 cgemHit->SetPositionOfPostPoint ( tmp_XYZ_post );
391 cgemHit->SetMomentumOfPrePoint ( tmp_P_pre );
392 cgemHit->SetMomentumOfPostPoint ( tmp_P_post );
393
394 /*
395// ******************************************************************
396 TVector3 Cu_pre = m_TMcHitCgem->GetPrePositionInCu ();
397 TVector3 Cu_post = m_TMcHitCgem->GetPostPositionInCu();
398 TVector3 P_Cu_pre = m_TMcHitCgem->GetMomentumOfCuPre ();
399 TVector3 P_Cu_post = m_TMcHitCgem->GetMomentumOfCuPost();
400
401 G4ThreeVector tmp_Cu_pre = G4ThreeVector(Cu_pre.x() , Cu_pre.y() , Cu_pre.z() );
402 G4ThreeVector tmp_Cu_post = G4ThreeVector(Cu_post.x() , Cu_post.y() , Cu_post.z() );
403 G4ThreeVector tmp_Cu_P_pre = G4ThreeVector(P_Cu_pre.x() ,P_Cu_pre.y() , P_Cu_pre.z() );
404 G4ThreeVector tmp_Cu_P_post = G4ThreeVector(P_Cu_post.x() ,P_Cu_post.y() , P_Cu_post.z() );
405
406 cgemHit->SetPrePositionInCu ( tmp_Cu_pre );
407 cgemHit->SetPostPositionInCu ( tmp_Cu_post );
408 cgemHit->SetMomentumOfCuPre ( tmp_Cu_P_pre );
409 cgemHit->SetMomentumOfCuPost ( tmp_Cu_P_post );
410// ******************************************************************************************
411 */
412 //cgemHit->Print();
413
414 cgemDC->insert(cgemHit);
415 delete m_TMcHitCgem;
416 }
417 }
418}
void SetParentID(G4int f_ID_parent)
Int_t GetParentID() const
Definition TMcHitCgem.h:39
TVector3 GetPositionOfPostPoint() const
Definition TMcHitCgem.h:44
Double_t GetStepLength() const
Definition TMcHitCgem.h:42
Int_t GetTrackID() const
Definition TMcHitCgem.h:36
Double_t GetGlobalTime() const
Definition TMcHitCgem.h:40
TVector3 GetMomentumOfPrePoint() const
Definition TMcHitCgem.h:45
Int_t GetPDGCode() const
Definition TMcHitCgem.h:38
Double_t GetTotalEnergyDeposit() const
Definition TMcHitCgem.h:41
TVector3 GetPositionOfPrePoint() const
Definition TMcHitCgem.h:43
Int_t GetLayerID() const
Definition TMcHitCgem.h:37
TVector3 GetMomentumOfPostPoint() const
Definition TMcHitCgem.h:46
const TMcHitCgem * getMcHitCgem(Int_t i) const
retrieve a McHitCgem From the collection, using the index into the array
const TObjArray * getMcHitCgemCol() const
retrieve the whole TObjArray of McHitCgem Data
Definition TMcHitEvent.h:42

Referenced by GetRootEvent().

◆ GetEmcDigi()

void BesTuningIO::GetEmcDigi ( void )

Definition at line 221 of file BesTuningIO.cc.

221 {
222 G4int THCID = -1;
223 THCID = m_DigiMan->GetDigiCollectionID("BesEmcDigitsCollection");
224 if (THCID>=0) {
225 BesEmcDigitsCollection* emcDC = new BesEmcDigitsCollection("BesEmcDigitizer","BesEmcDigitsCollection");
226 m_DigiMan->SetDigiCollection(THCID,emcDC);
227 }
228}
G4TDigiCollection< BesEmcDigi > BesEmcDigitsCollection
Definition BesEmcDigi.hh:69

Referenced by GetNextEvents().

◆ GetEmcRootDigi()

void BesTuningIO::GetEmcRootDigi ( )

Definition at line 246 of file BesTuningIO.cc.

246 {
247 G4int THCID = -1;
248 THCID = m_DigiMan->GetDigiCollectionID("BesEmcDigitsCollection");
249 //cout << "THCID " << THCID << endl;
250 if (THCID>=0) {
251 BesEmcDigitsCollection* emcDC = new BesEmcDigitsCollection("BesEmcDigitizer","BesEmcDigitsCollection");
252
253 //std::cout << "GetEmcRootDigi " << emcDC << std::endl;
254 if(emcDC){
255 int nHit = emcDC->entries();
256 //std::cout << "nHit: " << nHit << std::endl;
257 if(nHit > 0){
258 for(int i = 0; i < nHit; i++){
259 delete (*emcDC)[i];
260 }
261 emcDC->GetVector()->clear();
262 }
263 }
264
265 int nHits = m_TMcHitEvent->getMcDigiEmcCol()->GetEntries();
266 //std::cout << "nHits: " << nHits << std::endl;
267 for(int i = 0; i < nHits; i++){
268 m_TMcDigiEmc = m_TMcHitEvent->getMcDigiEmc(i);
269
270 BesEmcDigi* emcDigi = new BesEmcDigi();
271
272 emcDigi->SetPartId(m_TMcDigiEmc->GetPartId());
273 emcDigi->SetThetaNb(m_TMcDigiEmc->GetThetaNb());
274 emcDigi->SetPhiNb(m_TMcDigiEmc->GetPhiNb());
275 emcDigi->SetEnergy(m_TMcDigiEmc->GetEnergy());
276 emcDigi->SetTime(m_TMcDigiEmc->GetTime());
277 emcDigi->SetTrackIndex(m_TMcDigiEmc->GetTrackIndex());
278 //emcDigi->Print();
279 //std::cout << "SetEnergy" << emcDigi->GetEnergy() << std::endl;
280
281 emcDC->insert(emcDigi);
282 delete m_TMcDigiEmc;
283
284 }
285
286 //std::cout << "insert" << std::endl;
287 m_DigiMan->SetDigiCollection(THCID,emcDC);
288
289
290 }
291
292}
void SetTrackIndex(G4int index)
Definition BesEmcDigi.hh:46
void SetTime(G4double time)
Definition BesEmcDigi.hh:45
void SetPartId(G4int id)
Definition BesEmcDigi.hh:41
void SetEnergy(G4double energy)
Definition BesEmcDigi.hh:44
void SetThetaNb(G4int nTheta)
Definition BesEmcDigi.hh:42
void SetPhiNb(G4int nPhi)
Definition BesEmcDigi.hh:43
Int_t GetPartId() const
Definition TMcDigiEmc.h:24
Double_t GetTime() const
Definition TMcDigiEmc.h:28
Int_t GetThetaNb() const
Definition TMcDigiEmc.h:25
Int_t GetPhiNb() const
Definition TMcDigiEmc.h:26
Int_t GetTrackIndex() const
Definition TMcDigiEmc.h:29
Double_t GetEnergy() const
Definition TMcDigiEmc.h:27
const TObjArray * getMcDigiEmcCol() const
retrieve the whole TObjArray of McHitMdc Data
Definition TMcHitEvent.h:51
const TMcDigiEmc * getMcDigiEmc(Int_t i) const
retrieve a McHitMdc From the collection, using the index into the array

Referenced by GetRootEvent().

◆ GetMdcHits()

void BesTuningIO::GetMdcHits ( void )

Definition at line 120 of file BesTuningIO.cc.

120 {
121 G4int mdcHitCollID = -1;
122 mdcHitCollID = m_DigiMan->GetHitsCollectionID("BesMdcHitsCollection");
123 if (mdcHitCollID>=0){
124 BesMdcHitsCollection* mdcDC = (BesMdcHitsCollection*)m_DigiMan->GetHitsCollection(mdcHitCollID);
125 if (mdcDC){
126 G4int nHit = mdcDC->entries();
127 if (nHit>0){
128 for (G4int i=0;i<nHit;i++)
129 {
130 delete (*mdcDC)[i];
131 }
132 mdcDC->GetVector()->clear();
133 }
134
135 std::vector<MdcHitType>::iterator iter;
136 iter = (m_evt->mdcHit).hitCol.begin();
137 // Loop over mdc hits
138 for (; iter != (m_evt->mdcHit).hitCol.end(); iter++) {
139 BesMdcHit* newHit = new BesMdcHit();
140 newHit->SetTrackID ((*iter).trackIndex);
141 newHit->SetLayerNo ((*iter).layerNo);
142 newHit->SetCellNo ((*iter).cellNo);
143 newHit->SetEdep ((*iter).energyDeposit);
144 newHit->SetPos (G4ThreeVector((*iter).posX,(*iter).posY,(*iter).posZ));
145 newHit->SetDriftD ((*iter).driftDistance);
146 newHit->SetTheta((*iter).theta);
147 newHit->SetPosFlag((*iter).posFlag);
148 newHit->SetEnterAngle((*iter).enterAngle);
149 newHit->SetDriftT (0.);
150 newHit->SetGlobalT((*iter).globalT);
151 mdcDC->insert(newHit);
152
153 }
154 //mdcDC->PrintAllHits();
155 }else{
156 std::cerr << "BesTuningIO::can't get mdcHitsCollection"<<std::endl;
157 }
158 }else{
159 std::cerr << "BesTuningIO::can't get mdcHitCollID"<<std::endl;
160 }
161}
G4THitsCollection< BesMdcHit > BesMdcHitsCollection
Definition BesMdcHit.hh:97
void SetEdep(G4double de)
Definition BesMdcHit.hh:41
void SetDriftT(G4double time)
Definition BesMdcHit.hh:44
void SetEnterAngle(G4double angle)
Definition BesMdcHit.hh:47
void SetCellNo(G4int cell)
Definition BesMdcHit.hh:40
void SetPos(G4ThreeVector xyz)
Definition BesMdcHit.hh:42
void SetTrackID(G4int track)
Definition BesMdcHit.hh:38
void SetLayerNo(G4int layer)
Definition BesMdcHit.hh:39
void SetTheta(G4double angle)
Definition BesMdcHit.hh:46
void SetDriftD(G4double distance)
Definition BesMdcHit.hh:43
void SetGlobalT(G4double time)
Definition BesMdcHit.hh:45
void SetPosFlag(G4int flag)
Definition BesMdcHit.hh:48
MDCHIT mdcHit
Definition AsciiData.hh:545

Referenced by GetNextEvents().

◆ GetMdcRootHits()

void BesTuningIO::GetMdcRootHits ( )

Definition at line 294 of file BesTuningIO.cc.

294 {
295
296 G4int THCID = -1;
297 THCID = m_DigiMan->GetHitsCollectionID("BesMdcHitsCollection");
298 if (THCID>=0) {
299 BesMdcHitsCollection* mdcDC = (BesMdcHitsCollection*) (m_DigiMan->GetHitsCollection(THCID));
300 if(mdcDC){
301 int nHit = mdcDC->entries();
302 if(nHit > 0){
303 for(int i = 0; i < nHit; i++){
304 delete (*mdcDC)[i];
305 }
306 mdcDC->GetVector()->clear();
307 }
308 }
309
310 int nHits = m_TMcHitEvent->getMcHitMdcCol()->GetEntries();
311 //std::cout << "nHits: " << nHits << std::endl;
312 for(int i = 0; i < nHits; i++){
313 m_TMcHitMdc = m_TMcHitEvent->getMcHitMdc(i);
314
315 BesMdcHit* mdcHit = new BesMdcHit();
316
317 mdcHit->SetTrackID(m_TMcHitMdc->GetTrackID());
318 mdcHit->SetLayerNo(m_TMcHitMdc->GetLayerNo());
319 mdcHit->SetCellNo(m_TMcHitMdc->GetCellNo());
320 mdcHit->SetEdep(m_TMcHitMdc->GetEdep());
321 mdcHit->SetDriftD(m_TMcHitMdc->GetDriftD());
322 mdcHit->SetDriftT(m_TMcHitMdc->GetDriftT());
323 mdcHit->SetGlobalT(m_TMcHitMdc->GetGlobalT());
324 mdcHit->SetTheta(m_TMcHitMdc->GetTheta());
325 mdcHit->SetEnterAngle(m_TMcHitMdc->GetEnterAngle());
326 mdcHit->SetPosFlag(m_TMcHitMdc->GetPosFlag());
327
328 TVector3 tTemp = m_TMcHitMdc->GetPos();
329 G4ThreeVector gTemp = G4ThreeVector(tTemp.X(), tTemp.Y(), tTemp.Z());
330 mdcHit->SetPos(gTemp);
331 //mdcHit->Print();
332
333 mdcDC->insert(mdcHit);
334 delete m_TMcHitMdc;
335
336 }
337
338
339 }
340
341}
const TMcHitMdc * getMcHitMdc(Int_t i) const
retrieve a McHitMdc From the collection, using the index into the array
const TObjArray * getMcHitMdcCol() const
retrieve the whole TObjArray of McHitMdc Data
Definition TMcHitEvent.h:33
Double_t GetTheta() const
Definition TMcHitMdc.h:36
Double_t GetGlobalT() const
Definition TMcHitMdc.h:35
Int_t GetLayerNo() const
Definition TMcHitMdc.h:29
TVector3 GetPos() const
Definition TMcHitMdc.h:32
Double_t GetEdep() const
Definition TMcHitMdc.h:31
Double_t GetDriftT() const
Definition TMcHitMdc.h:34
Double_t GetEnterAngle() const
Definition TMcHitMdc.h:37
Double_t GetDriftD() const
Definition TMcHitMdc.h:33
Int_t GetPosFlag() const
Definition TMcHitMdc.h:38
Int_t GetCellNo() const
Definition TMcHitMdc.h:30
Int_t GetTrackID() const
Definition TMcHitMdc.h:28

Referenced by GetRootEvent().

◆ GetMucHits()

void BesTuningIO::GetMucHits ( void )
inline

Definition at line 40 of file BesTuningIO.hh.

40{};

Referenced by GetNextEvents(), and GetRootEvent().

◆ GetNextEvents()

void BesTuningIO::GetNextEvents ( void )

Definition at line 91 of file BesTuningIO.cc.

91 {
92 if (m_evt) delete m_evt;
93
94 m_evt = new HitEVENT;
95 try {
96 (*m_inputFileStream) >> *m_evt;
97 } catch (AsciiWrongTag& ex) {
98 std::cerr << "wrong tag, got " << ex.got()
99 << " expected: " << ex.expected()
100 << std::endl;
101 delete m_evt;
102 m_evt=0;
103 return;
104 } catch (AsciiDumpException&) {
105 std::cerr<<"BesTuningIO: Reach file end!"<<std::endl;
106 delete m_evt;
107 m_evt=0;
108 return;
109 }
110
112
114
116
118}
std::string got() const
Definition dmplib.hh:22
std::string expected() const
Definition dmplib.hh:21
void GetMdcHits(void)
void GetEmcDigi(void)
void GetMucHits(void)
void GetTofHits(void)
static G4int GetMdc()
static G4int GetMuc()
static G4int GetTof()
static G4int GetEmc()

Referenced by BesEventAction::EndOfEventAction().

◆ GetRootEvent()

void BesTuningIO::GetRootEvent ( int evtID)

Definition at line 233 of file BesTuningIO.cc.

233 {
234 //std::cout << "evtID: " << evtID << std::endl;
235 //HitTree->GetEntry(evtID);
236 HitChain->GetEntry(evtID);
237 //std::cout << "HitChain" << std::endl;
243}
void GetEmcRootDigi()
void GetTofRootHits()
void GetCgemRootHits()
void GetMdcRootHits()
static G4int GetCgem()

Referenced by BesEventAction::EndOfEventAction().

◆ GetTofHits()

void BesTuningIO::GetTofHits ( void )
inline

Definition at line 38 of file BesTuningIO.hh.

38{};

Referenced by GetNextEvents().

◆ GetTofRootHits()

void BesTuningIO::GetTofRootHits ( )

Definition at line 420 of file BesTuningIO.cc.

420 {
421
422 //retrieve G4Svc
423 ISvcLocator* svcLocator = Gaudi::svcLocator();
424 IG4Svc* tmpSvc;
425 StatusCode sc = svcLocator->service("G4Svc", tmpSvc);
426 G4Svc* m_G4Svc = dynamic_cast<G4Svc *>(tmpSvc);
427
428 const double m_beamTime = m_TMcHitEvent->getBeamTime();
429 m_G4Svc->SetBeamTime(m_beamTime);
430 //std::cout << "beamtime: " << m_beamTime << std::endl;
431
432 G4int THCID = -1;
433 THCID = m_DigiMan->GetHitsCollectionID("BesTofHitsCollection");
434 if (THCID>=0) {
435 BesTofHitsCollection* tofDC = (BesTofHitsCollection*) (m_DigiMan->GetHitsCollection(THCID));
436 if(tofDC){
437 int nHit = tofDC->entries();
438 if(nHit > 0){
439 for(int i = 0; i < nHit; i++){
440 delete (*tofDC)[i];
441 }
442 tofDC->GetVector()->clear();
443 }
444 }
445
446 int nHits = m_TMcHitEvent->getMcHitTofCol()->GetEntries();
447 //std::cout << "nHits: " << nHits << std::endl;
448 for(int i = 0; i < nHits; i++){
449 m_TMcHitTof = m_TMcHitEvent->getMcHitTof(i);
450
451 BesTofHit* tofHit = new BesTofHit();
452
453 tofHit->SetTrackIndex(m_TMcHitTof->GetTrackIndex());
454 tofHit->SetG4Index(m_TMcHitTof->GetG4Index());
455 tofHit->SetPartId(m_TMcHitTof->GetPartId());
456 tofHit->SetScinNb(m_TMcHitTof->GetScinNb());
457 tofHit->SetEdep(m_TMcHitTof->GetEdep());
458 tofHit->SetStepL(m_TMcHitTof->GetStepL());
459 tofHit->SetTrackL(m_TMcHitTof->GetTrackL());
460 tofHit->SetTime(m_TMcHitTof->GetTime());
461 tofHit->SetDeltaT(m_TMcHitTof->GetDeltaT());
462 tofHit->SetCharge(m_TMcHitTof->GetCharge());
463
464 TVector3 tTemp = m_TMcHitTof->GetPos();
465 G4ThreeVector gTemp(tTemp.X(), tTemp.Y(), tTemp.Z());
466 tofHit->SetPos(gTemp);
467
468 tTemp = m_TMcHitTof->GetPDirection();
469 gTemp = G4ThreeVector(tTemp.X(), tTemp.Y(), tTemp.Z());
470 tofHit->SetPDirection(gTemp);
471
472 tTemp = m_TMcHitTof->GetMomentum();
473 gTemp = G4ThreeVector(tTemp.X(), tTemp.Y(), tTemp.Z());
474 tofHit->SetMomentum(gTemp);
475
476 //tofHit->Print();
477
478 tofDC->insert(tofHit);
479 delete m_TMcHitTof;
480
481 }
482
483 //std::cout << "tofDC: " << tofDC->entries() << std::endl;
484
485 }
486
487}
G4THitsCollection< BesTofHit > BesTofHitsCollection
Definition BesTofHit.hh:108
void SetPos(G4ThreeVector pos)
Definition BesTofHit.hh:48
void SetDeltaT(G4double deltaT)
Definition BesTofHit.hh:50
void SetCharge(G4double charge)
Definition BesTofHit.hh:55
void SetTrackIndex(G4int trackIndex)
Definition BesTofHit.hh:41
void SetPDirection(G4ThreeVector pDirection)
Definition BesTofHit.hh:51
void SetPartId(G4int partId)
Definition BesTofHit.hh:43
void SetScinNb(G4int scinNb)
Definition BesTofHit.hh:44
void SetStepL(G4double stepL)
Definition BesTofHit.hh:46
void SetTrackL(G4double length)
Definition BesTofHit.hh:47
void SetTime(G4double time)
Definition BesTofHit.hh:49
void SetEdep(G4double edep)
Definition BesTofHit.hh:45
void SetMomentum(G4ThreeVector momentum)
Definition BesTofHit.hh:52
void SetG4Index(G4int index)
Definition BesTofHit.hh:42
Definition G4Svc.h:32
void SetBeamTime(double value)
Definition G4Svc.h:86
const TObjArray * getMcHitTofCol() const
retrieve the whole TObjArray of McHitTof Data
Definition TMcHitEvent.h:24
const TMcHitTof * getMcHitTof(Int_t i) const
retrieve a McHitTof From the collection, using the index into the array
Double_t getBeamTime() const
Definition TMcHitEvent.h:60
Int_t GetG4Index() const
Definition TMcHitTof.h:32
Int_t GetTrackIndex() const
Definition TMcHitTof.h:31
Double_t GetEdep() const
Definition TMcHitTof.h:35
TVector3 GetPDirection() const
Definition TMcHitTof.h:41
TVector3 GetMomentum() const
Definition TMcHitTof.h:42
Int_t GetPartId() const
Definition TMcHitTof.h:33
Double_t GetStepL() const
Definition TMcHitTof.h:36
Double_t GetTrackL() const
Definition TMcHitTof.h:37
Double_t GetDeltaT() const
Definition TMcHitTof.h:40
Int_t GetScinNb() const
Definition TMcHitTof.h:34
Double_t GetTime() const
Definition TMcHitTof.h:39
Int_t GetCharge() const
Definition TMcHitTof.h:43
TVector3 GetPos() const
Definition TMcHitTof.h:38

Referenced by GetRootEvent().

Member Data Documentation

◆ f

TFile* BesTuningIO::f

Definition at line 43 of file BesTuningIO.hh.

◆ HitChain

TChain* BesTuningIO::HitChain

Definition at line 45 of file BesTuningIO.hh.

Referenced by BesTuningIO(), and GetRootEvent().

◆ HitTree

TTree* BesTuningIO::HitTree

Definition at line 44 of file BesTuningIO.hh.


The documentation for this class was generated from the following files: