SND@LHC Software
Loading...
Searching...
No Matches
simpleTarget.cxx
Go to the documentation of this file.
1#include "simpleTarget.h"
2#include <math.h>
3#include "vetoPoint.h"
4
5#include "FairLogger.h" // for FairLogger, MESSAGE_ORIGIN
6#include "FairVolume.h"
7#include "FairGeoVolume.h"
8#include "FairGeoNode.h"
9#include "FairRootManager.h"
10#include "FairGeoLoader.h"
11#include "FairGeoInterface.h"
12#include "FairGeoMedia.h"
13#include "FairGeoBuilder.h"
14#include "FairRun.h"
15#include "FairRuntimeDb.h"
16#include "ShipDetectorList.h"
17#include "ShipStack.h"
18
19#include "TClonesArray.h"
20#include "TVirtualMC.h"
21#include "TGeoManager.h"
22#include "TGeoBBox.h"
23#include "TGeoEltu.h"
24#include "TGeoBoolNode.h"
25#include "TGeoMaterial.h"
26#include "TParticle.h"
27#include "TROOT.h"
28#include "TH1D.h"
29#include "TH2D.h"
30
31#include <iostream>
32using std::cout;
33using std::endl;
34
36 : FairDetector("simpleTarget", kTRUE, kVETO),
37 fTrackID(-1),
38 fVolumeID(-1),
39 fPos(),
40 fMom(),
41 fTime(-1.),
42 fLength(-1.),
43 fThick(-1.),
44 fOnlyMuons(kFALSE),
45 fFastMuon(kFALSE),
46 fzPos(3E8),
47 fTotalEloss(0),
48 fsimpleTargetPointCollection(new TClonesArray("vetoPoint"))
49{}
50
58
59Bool_t simpleTarget::ProcessHits(FairVolume* vol)
60{
62 //Set parameters at entrance of volume. Reset ELoss.
63 if ( gMC->IsTrackEntering() ) {
64 fELoss = 0.;
65 fTime = gMC->TrackTime() * 1.0e09;
66 fLength = gMC->TrackLength();
67 gMC->TrackPosition(fPos);
68 gMC->TrackMomentum(fMom);
69 }
70 // Sum energy loss for all steps in the active volume
71 fELoss += gMC->Edep();
72
73 // Create vetoPoint at exit of active volume
74 if ( gMC->IsTrackExiting() ||
75 gMC->IsTrackStop() ||
76 gMC->IsTrackDisappeared() ) {
77 if (fELoss == 0. ) { return kFALSE; }
78 // if (fOnlyMuons and fPos.Z()<140){ used for LiquidKrypton study
79 if (fOnlyMuons){
81 TClonesArray& clref = *fsimpleTargetPointCollection;
82 new (clref[0]) vetoPoint(0, 0, TVector3(0, 0, 0), TVector3(0, 0, 0),
83 0, 0, fTotalEloss, 0,TVector3(0, 0, 0), TVector3(0, 0, 0));
84 return kTRUE;}
85
86 fTrackID = gMC->GetStack()->GetCurrentTrackNumber();
87
88 Int_t veto_uniqueId;
89 gMC->CurrentVolID(veto_uniqueId);
90
91 TParticle* p=gMC->GetStack()->GetCurrentTrack();
92 Int_t pdgCode = p->GetPdgCode();
93 TLorentzVector Pos;
94 gMC->TrackPosition(Pos);
95 TLorentzVector Mom;
96 gMC->TrackMomentum(Mom);
97 Double_t xmean = (fPos.X()+Pos.X())/2. ;
98 Double_t ymean = (fPos.Y()+Pos.Y())/2. ;
99 Double_t zmean = (fPos.Z()+Pos.Z())/2. ;
100 //cout << veto_uniqueId << " :(" << xmean << ", " << ymean << ", " << zmean << "): " << fELoss << endl;
101 AddHit(fTrackID, veto_uniqueId, TVector3(xmean, ymean, zmean),
102 TVector3(fMom.Px(), fMom.Py(), fMom.Pz()), fTime, fLength,
103 fELoss,pdgCode,TVector3(Pos.X(), Pos.Y(), Pos.Z()),TVector3(Mom.Px(), Mom.Py(), Mom.Pz()) );
104
105 // Increment number of veto det points in TParticle
106 ShipStack* stack = (ShipStack*) gMC->GetStack();
107 stack->AddPoint(kVETO);
108 }
109
110 return kTRUE;
111}
112
113
115{
116 FairDetector::Initialize();
117 TSeqCollection* fileList=gROOT->GetListOfFiles();
118 fout = ((TFile*)fileList->At(0));
119}
120
126
127
129{
130 static FairGeoLoader *geoLoad=FairGeoLoader::Instance();
131 static FairGeoInterface *geoFace=geoLoad->getGeoInterface();
132 static FairGeoMedia *media=geoFace->getMedia();
133 static FairGeoBuilder *geoBuild=geoLoad->getGeoBuilder();
134
135 FairGeoMedium *ShipMedium=media->getMedium(fMaterial);
136 TGeoMedium* mat=gGeoManager->GetMedium(fMaterial);
137 if (mat==NULL)
138 geoBuild->createMedium(ShipMedium);
139 mat =gGeoManager->GetMedium(fMaterial);
140 TGeoVolume *top=gGeoManager->GetTopVolume();
141 TGeoVolume *target = gGeoManager->MakeBox("target",mat,10.*100.,10.*100.,fThick);
142 top->AddNode(target, 1, new TGeoTranslation(0, 0, fzPos));
143 AddSensitiveVolume(target);
144}
145
146vetoPoint* simpleTarget::AddHit(Int_t trackID, Int_t detID,
147 TVector3 pos, TVector3 mom,
148 Double_t time, Double_t length,
149 Double_t eLoss, Int_t pdgCode,TVector3 Lpos, TVector3 Lmom)
150{
151 TClonesArray& clref = *fsimpleTargetPointCollection;
152 Int_t size = clref.GetEntriesFast();
153 return new(clref[size]) vetoPoint(trackID, detID, pos, mom,
154 time, length, eLoss, pdgCode,Lpos,Lmom);
155}
156
158 if (!fFastMuon){return;}
159 if (TMath::Abs(gMC->TrackPid())!=13){
160 gMC->StopTrack();
161 }
162}
163
165{
166
167 FairRootManager::Instance()->Register("vetoPoint", "veto",
169}
170
171TClonesArray* simpleTarget::GetCollection(Int_t iColl) const
172{
173 if (iColl == 0) { return fsimpleTargetPointCollection; }
174 else { return NULL; }
175}
180
@ kVETO
Double_t fLength
time
virtual ~simpleTarget()
virtual void EndOfEvent()
Double_t fELoss
thickness
Bool_t fOnlyMuons
max energy to transport
vetoPoint * AddHit(Int_t trackID, Int_t detID, TVector3 pos, TVector3 mom, Double_t time, Double_t length, Double_t eLoss, Int_t pdgcode, TVector3 Lpos, TVector3 Lmom)
virtual void Reset()
TLorentzVector fMom
position at entrance
Double_t fzPos
length
virtual void Initialize()
Double_t fTime
momentum at entrance
TClonesArray * fsimpleTargetPointCollection
TLorentzVector fPos
volume id
virtual void Register()
Double_t fTotalEloss
void ConstructGeometry()
TFile * fout
for fast processing
virtual void PreTrack()
TString fMaterial
Bool_t fFastMuon
virtual Bool_t ProcessHits(FairVolume *v=0)
Double_t fThick
zPos
virtual TClonesArray * GetCollection(Int_t iColl) const
ClassImp(ecalContFact) ecalContFact