BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
BesEmcSD Class Reference

#include <BesEmcSD.hh>

Inheritance diagram for BesEmcSD:

Public Member Functions

 BesEmcSD (G4String, BesEmcConstruction *, BesEmcGeometry *)
 ~BesEmcSD ()
void Initialize (G4HCofThisEvent *)
G4bool ProcessHits (G4Step *, G4TouchableHistory *)
void EndOfEvent (G4HCofThisEvent *)
void clear ()
void DrawAll ()
void PrintAll ()
void BeginOfTruthEvent (const G4Event *)
void EndOfTruthEvent (const G4Event *)
void ComputeThetaPhi (G4int, G4int, G4int)
G4int EndPhiNumberForGDML (G4int)
G4int ComputeEndCopyNb (G4int)
 BesEmcSD (G4String, BesEmcConstruction *, BesEmcGeometry *)
 ~BesEmcSD ()
void Initialize (G4HCofThisEvent *)
G4bool ProcessHits (G4Step *, G4TouchableHistory *)
void EndOfEvent (G4HCofThisEvent *)
void clear ()
void DrawAll ()
void PrintAll ()
void BeginOfTruthEvent (const G4Event *)
void EndOfTruthEvent (const G4Event *)
void ComputeThetaPhi (G4int, G4int, G4int)
G4int EndPhiNumberForGDML (G4int)
G4int ComputeEndCopyNb (G4int)
 BesEmcSD (G4String, BesEmcConstruction *, BesEmcGeometry *)
 ~BesEmcSD ()
void Initialize (G4HCofThisEvent *)
G4bool ProcessHits (G4Step *, G4TouchableHistory *)
void EndOfEvent (G4HCofThisEvent *)
void clear ()
void DrawAll ()
void PrintAll ()
void BeginOfTruthEvent (const G4Event *)
void EndOfTruthEvent (const G4Event *)
void ComputeThetaPhi (G4int, G4int, G4int)
G4int EndPhiNumberForGDML (G4int)
G4int ComputeEndCopyNb (G4int)
Public Member Functions inherited from BesSensitiveDetector
 BesSensitiveDetector (const G4String name)
virtual ~BesSensitiveDetector ()
virtual void BeginOfTrack (const G4Track *)
virtual void EndOfTrack (const G4Track *)
 BesSensitiveDetector (const G4String name)
virtual ~BesSensitiveDetector ()
virtual void BeginOfTrack (const G4Track *)
virtual void EndOfTrack (const G4Track *)
 BesSensitiveDetector (const G4String name)
virtual ~BesSensitiveDetector ()
virtual void BeginOfTrack (const G4Track *)
virtual void EndOfTrack (const G4Track *)

Additional Inherited Members

Protected Member Functions inherited from BesSensitiveDetector
void GetCurrentTrackIndex (G4int &trackIndex, G4int &g4TrackId) const
void GetCurrentTrackIndex (G4int &trackIndex, G4int &g4TrackId) const
void GetCurrentTrackIndex (G4int &trackIndex, G4int &g4TrackId) const

Detailed Description

Constructor & Destructor Documentation

◆ BesEmcSD() [1/3]

BesEmcSD::BesEmcSD ( G4String name,
BesEmcConstruction * det,
BesEmcGeometry * besEMCGeometry )

Definition at line 31 of file BesEmcSD.cc.

32 : BesSensitiveDetector( name ), Detector( det ), fBesEmcGeometry( besEMCGeometry ) {
33 collectionName.insert( "BesEmcHitsCollection" ); // for digitization
34 collectionName.insert( "BesEmcHitsList" ); // for MC truth
35 collectionName.insert( "BesEmcTruthHitsList" );
36 HitID = new G4int[40000];
37}
BesSensitiveDetector(const G4String name)

◆ ~BesEmcSD() [1/3]

BesEmcSD::~BesEmcSD ( )

Definition at line 41 of file BesEmcSD.cc.

41{ delete[] HitID; }

◆ BesEmcSD() [2/3]

BesEmcSD::BesEmcSD ( G4String ,
BesEmcConstruction * ,
BesEmcGeometry *  )

◆ ~BesEmcSD() [2/3]

BesEmcSD::~BesEmcSD ( )

◆ BesEmcSD() [3/3]

BesEmcSD::BesEmcSD ( G4String ,
BesEmcConstruction * ,
BesEmcGeometry *  )

◆ ~BesEmcSD() [3/3]

BesEmcSD::~BesEmcSD ( )

Member Function Documentation

◆ BeginOfTruthEvent() [1/3]

void BesEmcSD::BeginOfTruthEvent ( const G4Event * )
virtual

Reimplemented from BesSensitiveDetector.

Definition at line 51 of file BesEmcSD.cc.

51 {
52 CalList = new BesEmcHitsCollection( SensitiveDetectorName, collectionName[1] );
53 CalTruthList = new BesEmcTruthHitsCollection( SensitiveDetectorName, collectionName[2] );
54 m_trackIndex = -99;
55}
G4THitsCollection< BesEmcTruthHit > BesEmcTruthHitsCollection
G4THitsCollection< BesEmcHit > BesEmcHitsCollection

◆ BeginOfTruthEvent() [2/3]

void BesEmcSD::BeginOfTruthEvent ( const G4Event * )
virtual

Reimplemented from BesSensitiveDetector.

◆ BeginOfTruthEvent() [3/3]

void BesEmcSD::BeginOfTruthEvent ( const G4Event * )
virtual

Reimplemented from BesSensitiveDetector.

◆ clear() [1/3]

void BesEmcSD::clear ( )

Definition at line 441 of file BesEmcSD.cc.

441{}

◆ clear() [2/3]

void BesEmcSD::clear ( )

◆ clear() [3/3]

void BesEmcSD::clear ( )

◆ ComputeEndCopyNb() [1/3]

G4int BesEmcSD::ComputeEndCopyNb ( G4int num)

Definition at line 416 of file BesEmcSD.cc.

416 {
417 G4int copyNb;
418 switch ( num )
419 {
420 case 30: copyNb = 5; break;
421 case 31: copyNb = 6; break;
422 case 32: copyNb = 14; break;
423 case 33: copyNb = 15; break;
424 case 34: copyNb = 16; break;
425 default: copyNb = num; break;
426 }
427 return copyNb;
428}

Referenced by ProcessHits().

◆ ComputeEndCopyNb() [2/3]

G4int BesEmcSD::ComputeEndCopyNb ( G4int )

◆ ComputeEndCopyNb() [3/3]

G4int BesEmcSD::ComputeEndCopyNb ( G4int )

◆ ComputeThetaPhi() [1/3]

void BesEmcSD::ComputeThetaPhi ( G4int id,
G4int sector,
G4int nb )

Definition at line 301 of file BesEmcSD.cc.

301 {
302 if ( ( sector >= 0 ) && ( sector < 3 ) ) sector += 16;
303 // if((sector!=7)&&(sector!=15))
304 {
305 if ( ( nb >= 0 ) && ( nb < 4 ) )
306 {
307 CryNumberTheta = 0;
308 CryNumberPhi = ( sector - 3 ) * 4 + nb;
309 }
310 else if ( ( nb >= 4 ) && ( nb < 8 ) )
311 {
312 CryNumberTheta = 1;
313 CryNumberPhi = ( sector - 3 ) * 4 + nb - 4;
314 }
315 else if ( ( nb >= 8 ) && ( nb < 13 ) )
316 {
317 CryNumberTheta = 2;
318 CryNumberPhi = ( sector - 3 ) * 5 + nb - 8;
319 }
320 else if ( ( nb >= 13 ) && ( nb < 18 ) )
321 {
322 CryNumberTheta = 3;
323 CryNumberPhi = ( sector - 3 ) * 5 + nb - 13;
324 }
325 else if ( ( nb >= 18 ) && ( nb < 24 ) )
326 {
327 CryNumberTheta = 4;
328 CryNumberPhi = ( sector - 3 ) * 6 + nb - 18;
329 }
330 else if ( ( nb >= 24 ) && ( nb < 30 ) )
331 {
332 CryNumberTheta = 5;
333 CryNumberPhi = ( sector - 3 ) * 6 + nb - 24;
334 }
335 }
336 /*else// if((sector=7)||(sector==15))
337 {
338 if((nb>=0)&&(nb<4))
339 {
340 CryNumberTheta = 0;
341 CryNumberPhi = (sector-3)*4+3-nb;
342 }
343 else if((nb>=4)&&(nb<8))
344 {
345 CryNumberTheta = 1;
346 CryNumberPhi = (sector-3)*4+7-nb;
347 }
348 else if((nb>=8)&&(nb<13))
349 {
350 CryNumberTheta = 2;
351 CryNumberPhi = (sector-3)*5+12-nb;
352 }
353 else if((nb>=13)&&(nb<18))
354 {
355 CryNumberTheta = 3;
356 CryNumberPhi = (sector-3)*5+17-nb;
357 }
358 else if((nb>=18)&&(nb<24))
359 {
360 CryNumberTheta = 4;
361 CryNumberPhi = (sector-3)*6+23-nb;
362 }
363 else if((nb>=24)&&(nb<30))
364 {
365 CryNumberTheta = 5;
366 CryNumberPhi = (sector-3)*6+29-nb;
367 }
368 }*/
369
370 if ( id == 2 )
371 {
372 switch ( CryNumberTheta )
373 {
374 case 0:
375 if ( CryNumberPhi < 32 ) CryNumberPhi = 31 - CryNumberPhi;
376 else CryNumberPhi = 95 - CryNumberPhi;
377 break;
378 case 1:
379 if ( CryNumberPhi < 32 ) CryNumberPhi = 31 - CryNumberPhi;
380 else CryNumberPhi = 95 - CryNumberPhi;
381 break;
382 case 2:
383 if ( CryNumberPhi < 40 ) CryNumberPhi = 39 - CryNumberPhi;
384 else CryNumberPhi = 119 - CryNumberPhi;
385 break;
386 case 3:
387 if ( CryNumberPhi < 40 ) CryNumberPhi = 39 - CryNumberPhi;
388 else CryNumberPhi = 119 - CryNumberPhi;
389 break;
390 case 4:
391 if ( CryNumberPhi < 48 ) CryNumberPhi = 47 - CryNumberPhi;
392 else CryNumberPhi = 143 - CryNumberPhi;
393 break;
394 case 5:
395 if ( CryNumberPhi < 48 ) CryNumberPhi = 47 - CryNumberPhi;
396 else CryNumberPhi = 143 - CryNumberPhi;
397 default:;
398 }
399 }
400}

Referenced by ProcessHits().

◆ ComputeThetaPhi() [2/3]

void BesEmcSD::ComputeThetaPhi ( G4int ,
G4int ,
G4int  )

◆ ComputeThetaPhi() [3/3]

void BesEmcSD::ComputeThetaPhi ( G4int ,
G4int ,
G4int  )

◆ DrawAll() [1/3]

void BesEmcSD::DrawAll ( )

Definition at line 445 of file BesEmcSD.cc.

445{}

◆ DrawAll() [2/3]

void BesEmcSD::DrawAll ( )

◆ DrawAll() [3/3]

void BesEmcSD::DrawAll ( )

◆ EndOfEvent() [1/3]

void BesEmcSD::EndOfEvent ( G4HCofThisEvent * HCE)

Definition at line 432 of file BesEmcSD.cc.

432 {
433 static G4int HCID = -1;
434 if ( HCID < 0 )
435 { HCID = G4SDManager::GetSDMpointer()->GetCollectionID( collectionName[0] ); }
436 HCE->AddHitsCollection( HCID, CalCollection );
437}

◆ EndOfEvent() [2/3]

void BesEmcSD::EndOfEvent ( G4HCofThisEvent * )

◆ EndOfEvent() [3/3]

void BesEmcSD::EndOfEvent ( G4HCofThisEvent * )

◆ EndOfTruthEvent() [1/3]

void BesEmcSD::EndOfTruthEvent ( const G4Event * evt)
virtual

Reimplemented from BesSensitiveDetector.

Definition at line 57 of file BesEmcSD.cc.

57 {
58 static G4int HLID = -1;
59 if ( HLID < 0 ) HLID = G4SDManager::GetSDMpointer()->GetCollectionID( collectionName[1] );
60 G4HCofThisEvent* HCE = evt->GetHCofThisEvent();
61 HCE->AddHitsCollection( HLID, CalList );
62
63 static G4int HTID = -1;
64 if ( HTID < 0 ) HTID = G4SDManager::GetSDMpointer()->GetCollectionID( collectionName[2] );
65 HCE->AddHitsCollection( HTID, CalTruthList );
66}

◆ EndOfTruthEvent() [2/3]

void BesEmcSD::EndOfTruthEvent ( const G4Event * )
virtual

Reimplemented from BesSensitiveDetector.

◆ EndOfTruthEvent() [3/3]

void BesEmcSD::EndOfTruthEvent ( const G4Event * )
virtual

Reimplemented from BesSensitiveDetector.

◆ EndPhiNumberForGDML() [1/3]

G4int BesEmcSD::EndPhiNumberForGDML ( G4int copyNb)

Definition at line 402 of file BesEmcSD.cc.

402 {
403 if ( copyNb < 0 || copyNb > 15 )
404 {
405 // G4Exception("Wrong copy number of EMC Endcap Phi Volume!");
406 G4cout << "Wrong copy number of EMC Endcap Phi Volume!" << G4endl;
407 exit( -1 );
408 }
409
410 if ( copyNb == 0 || copyNb == 1 ) { return 15 - copyNb * 8; }
411 else if ( copyNb == 2 || copyNb == 3 ) { return 30 - copyNb * 8; }
412 else if ( copyNb <= 9 ) { return 17 - copyNb; }
413 else { return 15 - copyNb; }
414}

Referenced by ProcessHits().

◆ EndPhiNumberForGDML() [2/3]

G4int BesEmcSD::EndPhiNumberForGDML ( G4int )

◆ EndPhiNumberForGDML() [3/3]

G4int BesEmcSD::EndPhiNumberForGDML ( G4int )

◆ Initialize() [1/3]

void BesEmcSD::Initialize ( G4HCofThisEvent * )

Definition at line 45 of file BesEmcSD.cc.

45 {
46 CalCollection = new BesEmcHitsCollection( SensitiveDetectorName, collectionName[0] );
47 for ( G4int j = 0; j < 40000; j++ ) { HitID[j] = -1; }
48 nHit = 0;
49}

◆ Initialize() [2/3]

void BesEmcSD::Initialize ( G4HCofThisEvent * )

◆ Initialize() [3/3]

void BesEmcSD::Initialize ( G4HCofThisEvent * )

◆ PrintAll() [1/3]

void BesEmcSD::PrintAll ( )

Definition at line 449 of file BesEmcSD.cc.

449{}

◆ PrintAll() [2/3]

void BesEmcSD::PrintAll ( )

◆ PrintAll() [3/3]

void BesEmcSD::PrintAll ( )

◆ ProcessHits() [1/3]

G4bool BesEmcSD::ProcessHits ( G4Step * aStep,
G4TouchableHistory *  )

Definition at line 69 of file BesEmcSD.cc.

69 {
70
71 // Ma Qiumei Add for removing Warning Infomations
72 gErrorIgnoreLevel = kFatal;
73
74 G4double edep = aStep->GetTotalEnergyDeposit();
75 G4double stepl = aStep->GetStepLength();
76 if ( ( edep == 0. ) && ( stepl == 0. ) ) return false;
77
78 G4TouchableHistory* theTouchable =
79 (G4TouchableHistory*)( aStep->GetPreStepPoint()->GetTouchable() );
80 G4VPhysicalVolume* physVol = theTouchable->GetVolume();
81
82 if ( physVol->GetName().contains( "physicalCrystal" ) ) // barrel code
83 {
84 PartId = 1;
85 CryNumberPhi = theTouchable->GetReplicaNumber( 2 );
86 CryNumberTheta = theTouchable->GetReplicaNumber( 1 );
87 }
88 else if ( physVol->GetName().contains( "logicalCrystal" ) ) // barrel gdml
89 {
90 PartId = 1;
91 std::istringstream thetaBuf( ( theTouchable->GetVolume( 1 )->GetName() ).substr( 16, 2 ) );
92 thetaBuf >> CryNumberTheta;
93 CryNumberPhi = 308 - theTouchable->GetCopyNumber( 2 );
94 }
95 else if ( physVol->GetName().contains( "physicalEndCrystal" ) ) // endcap code
96 {
97 PartId = theTouchable->GetReplicaNumber( 3 );
98 G4int endSector = theTouchable->GetReplicaNumber( 2 );
99 G4int endNb = theTouchable->GetReplicaNumber( 0 );
100 ComputeThetaPhi( PartId, endSector, ComputeEndCopyNb( endNb ) );
101 }
102 else if ( physVol->GetName().contains( "logicalEndCrystal" ) ) // endcap gdml
103 {
104 PartId = 2 - 2 * theTouchable->GetCopyNumber( 3 );
105 G4int endSector = theTouchable->GetCopyNumber( 2 );
106 G4int endNb, endNbGDML;
107 std::istringstream thetaBuf( ( theTouchable->GetVolume( 0 )->GetName() ).substr( 20, 2 ) );
108 thetaBuf >> endNb;
109 endNbGDML = ComputeEndCopyNb( endNb );
110 G4int endSectorGDML = EndPhiNumberForGDML( endSector );
111 ComputeThetaPhi( PartId, endSectorGDML, endNbGDML );
112 }
113 else if ( physVol->GetName().contains( "physicalPD" ) ) // barrel PD code
114 {
115 edep *= 50.;
116 PartId = 1;
117 CryNumberPhi = theTouchable->GetReplicaNumber( 2 );
118 CryNumberTheta = theTouchable->GetReplicaNumber( 1 );
119 }
120 else if ( physVol->GetName().contains( "logicalPD" ) ) // barrel PD gdml
121 {
122 edep *= 50.;
123 PartId = 1;
124 G4int nb = theTouchable->GetCopyNumber( 1 );
125 if ( nb == 216 ) { CryNumberTheta = 0; }
126 else { CryNumberTheta = 43 - ( nb - 2 ) / 5; }
127 CryNumberPhi = 308 - theTouchable->GetCopyNumber( 2 );
128 }
129
130 if ( verboseLevel > 1 )
131 G4cout << "(Check ID)New EMC Hit on crystal: " << CryNumberPhi << "(phi)" << CryNumberTheta
132 << "(theta)"
133 << " and the volume is " << physVol->GetName() << G4endl;
134
135 //-----------------------------------------------------
136 G4int trackId = aStep->GetTrack()->GetTrackID();
137
138 BesEmcHit* calHit = new BesEmcHit();
139 calHit->SetTrackIndex( -99 );
140 calHit->SetG4Index( trackId );
141 calHit->SetEdepCrystal( edep );
142 calHit->SetTrakCrystal( stepl );
143 calHit->SetPosCrystal( aStep->GetPreStepPoint()->GetPosition() );
144 calHit->SetTimeCrystal( aStep->GetPreStepPoint()->GetGlobalTime() );
145 calHit->SetMomentum( aStep->GetPreStepPoint()->GetMomentum() );
146 calHit->SetNumCrystal( PartId, CryNumberTheta, CryNumberPhi );
147 // calHit->Print(1);
148 if ( edep > 0 )
149 {
150 G4int idhit = CalCollection->insert( calHit ) - 1;
151 if ( nHit < 40000 ) HitID[nHit] = idhit;
152 nHit++;
153 }
154 else
155 {
156 delete calHit;
157 // if(nHit==40000)
158 // G4cout << "Hits number exceed store space!!!" << G4endl;
159 }
160
161 // for MC Truth
162 if ( CalList )
163 {
164 G4int trackIndex, g4TrackId;
165 GetCurrentTrackIndex( trackIndex, g4TrackId );
166 calHit->SetTrackIndex( trackIndex );
167
168 if ( m_trackIndex != trackIndex ) // find a new track
169 {
170 m_trackIndex = trackIndex;
171
172 G4int flag = 1;
173 G4int pdg = abs( aStep->GetTrack()->GetDefinition()->GetPDGEncoding() );
174 if ( pdg == 12 || pdg == 14 || pdg == 16 )
175 { // neutrino
176 flag = 0;
177 }
178
179 if ( flag && aStep->GetTrack()->GetTrackID() == g4TrackId ) // is primary particle
180 {
181 BesEmcHit* truHit = new BesEmcHit();
182 truHit->SetTrackIndex( trackIndex );
183 truHit->SetG4Index( trackId );
184 truHit->SetEdepCrystal( edep );
185 truHit->SetTrakCrystal( stepl );
186 truHit->SetPosCrystal( aStep->GetPreStepPoint()->GetPosition() );
187 truHit->SetTimeCrystal( aStep->GetPreStepPoint()->GetGlobalTime() );
188 truHit->SetMomentum( aStep->GetPreStepPoint()->GetMomentum() );
189 truHit->SetNumCrystal( PartId, CryNumberTheta, CryNumberPhi );
190 CalList->insert( truHit );
191 }
192 }
193
194 else if ( m_trackIndex == trackIndex )
195 {
196 G4int length = CalList->entries();
197 if ( length >= 1 )
198 {
199 for ( G4int i = 0; i < length; i++ )
200 {
201 BesEmcHit* oldHit = ( *CalList )[i];
202 if ( oldHit->GetTrackIndex() == trackIndex ) // find the same trackIndex in CalList
203 {
204 oldHit->SetEdepCrystal( oldHit->GetEdepCrystal() + edep );
205 break;
206 }
207 }
208 }
209 }
210 }
211
212 // for full Mc Truth
213 if ( CalTruthList )
214 {
215 G4int trackIndex, g4TrackId;
216 GetCurrentTrackIndex( trackIndex, g4TrackId );
217 Identifier id = EmcID::crystal_id( PartId, CryNumberTheta, CryNumberPhi );
218
219 G4int flag = 1;
220 if ( CalTruthList->entries() > 0 )
221 {
222 for ( G4int i = 0; i < CalTruthList->entries(); i++ )
223 {
224 BesEmcTruthHit* oldHit = ( *CalTruthList )[i];
225 if ( oldHit->GetTrackIndex() == trackIndex )
226 { // find the same trackIndex in CalList
227 flag = 0;
228 oldHit->SetEDep( oldHit->GetEDep() + edep ); // total energy
229 if ( oldHit->Find( id ) != oldHit->End() ) { oldHit->AddEHit( id, edep ); }
230 else { oldHit->Insert( id, edep ); }
231 break;
232 }
233 }
234 }
235
236 if ( flag == 1 )
237 { // new track
238 G4int pdg = abs( aStep->GetTrack()->GetDefinition()->GetPDGEncoding() );
239 if ( !( pdg == 12 || pdg == 14 || pdg == 16 ) )
240 { // not neutrino
241 BesEmcTruthHit* truHit = new BesEmcTruthHit;
242 truHit->SetTrackIndex( trackIndex );
243 truHit->SetG4TrackId( g4TrackId );
244 truHit->SetEDep( edep );
245 truHit->SetIdentify( id );
246 truHit->Insert( id, edep );
247
248 if ( aStep->GetTrack()->GetTrackID() == g4TrackId )
249 { // is primary particle
250 truHit->SetHitEmc( 1 );
251 truHit->SetPDGCode( aStep->GetTrack()->GetDefinition()->GetPDGEncoding() );
252 truHit->SetPDGCharge( aStep->GetTrack()->GetDefinition()->GetPDGCharge() );
253 truHit->SetParticleName( aStep->GetTrack()->GetDefinition()->GetParticleName() );
254 truHit->SetTime( aStep->GetPreStepPoint()->GetGlobalTime() );
255 truHit->SetPosition( aStep->GetPreStepPoint()->GetPosition() );
256 truHit->SetMomentum( aStep->GetPreStepPoint()->GetMomentum() );
257 }
258 else
259 { // the primary particle doesn't hit emc
260
261 BesSensitiveManager* sensitiveManager = BesSensitiveManager::GetSensitiveManager();
262 std::vector<BesTruthTrack*>* trackList = sensitiveManager->GetTrackList();
263
264 for ( unsigned i = 0; i < trackList->size(); i++ )
265 {
266 BesTruthTrack* truthTrack = ( *trackList )[i];
267
268 if ( trackIndex == truthTrack->GetIndex() )
269 { // find primary particle
270 truHit->SetHitEmc( 0 );
271 truHit->SetPDGCode( truthTrack->GetPDGCode() );
272 truHit->SetPDGCharge( truthTrack->GetPDGCharge() );
273 truHit->SetParticleName( truthTrack->GetParticleName() );
274
275 if ( truthTrack->GetTerminalVertex() )
276 {
277 truHit->SetTime( truthTrack->GetTerminalVertex()->GetTime() );
278 truHit->SetPosition( truthTrack->GetTerminalVertex()->GetPosition() );
279 }
280 else
281 { // have not found terminal vertex
282 truHit->SetTime( -99 );
283 truHit->SetPosition( G4ThreeVector( -99, -99, -99 ) );
284 }
285
286 break;
287 } // end if find primary particle
288 }
289 }
290
291 CalTruthList->insert( truHit );
292 }
293 }
294 }
295
296 return true;
297}
void SetNumCrystal(G4int id, G4int numTheta, G4int numPhi)
G4int ComputeEndCopyNb(G4int)
Definition BesEmcSD.cc:416
void ComputeThetaPhi(G4int, G4int, G4int)
Definition BesEmcSD.cc:301
G4int EndPhiNumberForGDML(G4int)
Definition BesEmcSD.cc:402
void AddEHit(Identifier, G4double)
Definition BesEmcHit.cc:173
std::map< Identifier, G4double >::const_iterator End() const
Definition BesEmcHit.cc:163
std::map< Identifier, G4double >::const_iterator Find(Identifier) const
Definition BesEmcHit.cc:167
void Insert(Identifier, G4double)
Definition BesEmcHit.cc:177
void GetCurrentTrackIndex(G4int &trackIndex, G4int &g4TrackId) const
static Identifier crystal_id(const unsigned int barrel_ec, const unsigned int theta_module, const unsigned int phi_module)
For a single crystal.
Definition EmcID.cxx:63

◆ ProcessHits() [2/3]

G4bool BesEmcSD::ProcessHits ( G4Step * ,
G4TouchableHistory *  )

◆ ProcessHits() [3/3]

G4bool BesEmcSD::ProcessHits ( G4Step * ,
G4TouchableHistory *  )

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