BOSS 7.1.0
BESIII Offline Software System
Loading...
Searching...
No Matches
BesMdcSD.hh
Go to the documentation of this file.
1//---------------------------------------------------------------------------//
2// BOOST --- BESIII Object_Oriented Simulation Tool //
3//---------------------------------------------------------------------------
4//Description: Sensitive detector definition for axial layers of MDC
5//Author: Yuan Ye([email protected])
6//Created: 4 Dec, 2003
7//Modified:
8//Comment: Used in "BesMdc"
9//---------------------------------------------------------------------------//
10
11#ifndef BesMdcSD_h
12#define BesMdcSD_h 1
13class TFile;
14class TH1F;
15class G4Svc;
16
17//#include "G4VSensitiveDetector.hh"
19#include "BesMdcHit.hh"
20#include "BesMdcGeoParameter.hh"
21#include "BesMdcCalTransfer.hh"
26#include "GaudiKernel/NTuple.h"
27#include "TF1.h"
28
29#include <vector>
30
31class G4Step;
32class G4HCofThisEvent;
34class G4Svc;
35
37{
38public:
39 BesMdcSD(G4String);
40 ~BesMdcSD();
41
42 void Initialize(G4HCofThisEvent*);
43 G4bool ProcessHits(G4Step*, G4TouchableHistory*);
44 void EndOfEvent(G4HCofThisEvent*);
45
46 void BeginOfTruthEvent(const G4Event*);
47 void EndOfTruthEvent(const G4Event*);
48
49 G4double Distance(G4int, G4int,G4ThreeVector,G4ThreeVector,G4ThreeVector&,G4double&);
50
51 void dedxFuncInti(void);
52
53private:
54 G4int hitPointer[43][288],truthPointer[43][288];
55 BesMdcHitsCollection* hitsCollection;
56 BesMdcHitsCollection* truthCollection;
57 BesMdcGeoParameter* mdcGeoPointer;
58 BesMdcCalTransfer* mdcCalPointer;
59 MdcGeomSvc* mdcGeomSvc;
60 G4Svc* m_G4Svc;
61 TF1 *dEdE_mylanfunc;
62
63 ///dedx sim ---------------------------
64 CalibDataSvc* m_calibDataSvc;
65 DedxSimSvc* m_pDedxSimSvc;
66 IDedxCurSvc* m_pDedxCurSvc;
67 std::vector<TH1F>* m_dedx_hists;
68 G4int m_version;
69 G4int m_numDedxHists;
70 G4int m_numBg;
71 G4int m_numTheta;
72 std::vector<double>* m_bgRange;
73 G4int GetBetagammaIndex(G4double bg);
74 G4int GetAngleIndex(G4double);
75 G4int GetChargeIndex(G4int);
76 G4double GetValDedxCurve(G4double bg, G4double charge);
77 G4double dedxSample(G4double betagamma, G4double length, G4double theta);
78
79 //dedx ntuple
80 NTuple::Tuple* m_tupleMdc;
81 NTuple::Item<double> m_betaGamma;
82 NTuple::Item<double> m_fitval;
83 NTuple::Item<double> m_random;
84 NTuple::Item<double> m_dedx;
85 NTuple::Item<double> m_de;
86 //NTuple::Item<double> m_length;
87 NTuple::Item<double> m_charge;
88 NTuple::Item<double> m_costheta;
89 ///------------------------------------
90};
91
92#endif
93
G4THitsCollection< BesMdcHit > BesMdcHitsCollection
Definition: BesMdcHit.hh:78
void dedxFuncInti(void)
Definition: BesMdcSD.cc:412
void BeginOfTruthEvent(const G4Event *)
Definition: BesMdcSD.cc:107
void EndOfTruthEvent(const G4Event *)
Definition: BesMdcSD.cc:114
G4double Distance(G4int, G4int, G4ThreeVector, G4ThreeVector, G4ThreeVector &, G4double &)
Definition: BesMdcSD.cc:351
G4bool ProcessHits(G4Step *, G4TouchableHistory *)
Definition: BesMdcSD.cc:125
void EndOfEvent(G4HCofThisEvent *)
Definition: BesMdcSD.cc:338
~BesMdcSD()
Definition: BesMdcSD.cc:86
Definition: G4Svc.h:33
float charge
float bg