source: trunk/examples/advanced/Rich/src/RichTbSD.cc@ 1187

Last change on this file since 1187 was 807, checked in by garnier, 17 years ago

update

File size: 5.4 KB
RevLine 
[807]1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// Rich advanced example for Geant4
27// RichTbSD.cc for Rich of LHCb
28// History:
29// Created: Sajan Easo (Sajan.Easo@cern.ch)
30// Revision and changes: Patricia Mendez (Patricia.Mendez@cern.ch)
31/////////////////////////////////////////////////////////////////////////////
32#include "RichTbSD.hh"
33#include "RichTbHit.hh"
34#include "G4Step.hh"
35#include "G4HCofThisEvent.hh"
36#include "G4Track.hh"
37#include "G4SDManager.hh"
38#include "G4ios.hh"
39#include "G4VProcess.hh"
40
41RichTbSD::RichTbSD(G4String DetName, G4int Numhpd , G4int NumSect,
42 G4int NumPixel,
43 G4String ) // CollectionName)
44 :G4VSensitiveDetector(DetName),NumberOfSensitiveHpds(Numhpd),
45 NumberOfSensitiveSectorsInHpd(NumSect),
46 NumberOfSensitivePixelsInSector(NumPixel),
47 HpdSDID(Numhpd*NumSect*NumPixel),HCID(-1)
48{
49 G4String HCname;
50 collectionName.insert(HCname="RichTbHitsCollection");
51
52}
53
54RichTbSD::~RichTbSD(){
55
56}
57
58void RichTbSD::Initialize(G4HCofThisEvent*)
59{
60 RichTbHitCollection = new RichTbHitsCollection
61 (SensitiveDetectorName,collectionName[0]);
62 for (G4int jhpd =0; jhpd < NumberOfHpds; jhpd++){
63 for(G4int jsect=0; jsect < NumberOfSiDetSectors; jsect++){
64 for (G4int jpixel = 0; jpixel < NumberOfPadHpdSiPixels ; jpixel++ ){
65
66 HpdSDID[NumberOfSensitiveSectorsInHpd*
67 NumberOfSensitivePixelsInSector*jhpd+
68 NumberOfSensitivePixelsInSector*jsect+jpixel]=-1;
69 }
70 }
71 }
72
73}
74
75
76G4bool RichTbSD::ProcessHits(G4Step*aStep,G4TouchableHistory*ROhist)
77{
78
79
80
81 G4Track* aTrack = aStep->GetTrack();
82 if(aTrack ->GetDefinition()->GetPDGCharge() == 0.0 ) return false;
83 G4double edep = aStep->GetTotalEnergyDeposit();
84
85 if(edep < 0.0001) return false;
86 G4StepPoint* PrePosition = aStep->GetPreStepPoint();
87
88
89 G4int CurrentHpdNumber = PrePosition->GetTouchableHandle()
90 -> GetReplicaNumber(1);
91 G4int CurrentSectorNumber= PrePosition->GetTouchableHandle()
92 -> GetReplicaNumber();
93 G4VPhysicalVolume* ROphysVol = ROhist -> GetVolume();
94 G4int CurrentPixelNumber = ROphysVol->GetCopyNo();
95
96 G4ThreeVector HitAtPhotoCathode;
97 if(aTrack->GetCreatorProcess() -> GetProcessName() == "PadHpdPhot") {
98 HitAtPhotoCathode=aTrack->GetVertexPosition();
99 }
100
101 G4int CopyId= NumberOfSensitiveSectorsInHpd*
102 NumberOfSensitivePixelsInSector*CurrentHpdNumber+
103 NumberOfSensitivePixelsInSector*CurrentSectorNumber+
104 CurrentPixelNumber;
105
106
107 if( HpdSDID[CopyId ] == -1 ) {
108 RichTbHit* newHit = new RichTbHit();
109 newHit->SetEdep( edep );
110 newHit->SetPos( aStep->GetPreStepPoint()->GetPosition() );
111 newHit->SetPosPC(HitAtPhotoCathode);
112 newHit->SetCurHpdNum ( CurrentHpdNumber );
113 newHit->SetCurSectNum ( CurrentSectorNumber );
114 newHit->SetCurPixNum ( CurrentPixelNumber );
115 //For now the same SD detector cell is allowed to have
116 //multiple hits from the same event. This is ok for
117 // analogue readout for photoelectron counting.
118 //This is to be changed in the future.
119 G4int NumHits = RichTbHitCollection->insert( newHit );
120
121 HpdSDID[CopyId]= NumHits -1 ;
122
123 } else {
124 // the current pixel is already hit in this event.
125 // here we can add extra energy (adc counts) to the
126 // existing hit. But this is not relevant to the
127 // current Tb analysis.
128 (*RichTbHitCollection)[HpdSDID[CopyId]]->AddEdep( edep );
129
130 G4cout << " Multiple hits in Hpd sector pixel " << CurrentHpdNumber
131 <<" "<< CurrentSectorNumber<<" "<< CurrentPixelNumber<< G4endl;
132 }
133
134 return true;
135}
136
137void RichTbSD::EndOfEvent(G4HCofThisEvent*HCE){
138 if( HCID < 0 ){
139 HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]); }
140 HCE->AddHitsCollection( HCID, RichTbHitCollection );
141
142}
143
144void RichTbSD::clear(){}
145
146void RichTbSD::DrawAll(){ }
147
148void RichTbSD::PrintAll(){ }
149
150
151
152
Note: See TracBrowser for help on using the repository browser.