SND@LHC Software
Loading...
Searching...
No Matches
boxTarget Class Reference

#include <boxTarget.h>

Inheritance diagram for boxTarget:
Collaboration diagram for boxTarget:

Public Member Functions

 boxTarget (const char *Name, Bool_t Active)
 
 boxTarget ()
 
virtual ~boxTarget ()
 
virtual void Initialize ()
 
virtual Bool_t ProcessHits (FairVolume *v=0)
 
virtual void Register ()
 
virtual TClonesArray * GetCollection (Int_t iColl) const
 
virtual void Reset ()
 
void ConstructGeometry ()
 
virtual void CopyClones (TClonesArray *cl1, TClonesArray *cl2, Int_t offset)
 
virtual void SetSpecialPhysicsCuts ()
 
virtual void EndOfEvent ()
 
virtual void FinishPrimary ()
 
virtual void FinishRun ()
 
virtual void BeginPrimary ()
 
virtual void PostTrack ()
 
virtual void PreTrack ()
 
virtual void BeginEvent ()
 
vetoPointAddHit (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)
 
void SetTarget (TString material, Float_t L, Bool_t choice)
 

Private Member Functions

Int_t InitMedium (TString name)
 

Private Attributes

Int_t fTrackID
 
Int_t fVolumeID
 track index
 
TLorentzVector fPos
 volume id
 
TLorentzVector fMom
 position at entrance
 
Double_t fTime
 momentum at entrance
 
Double_t fLength
 time
 
TString fTargetMaterial
 track length
 
Float_t fTargetL
 target material
 
Int_t index
 target length
 
Bool_t fBox
 
TClonesArray * fboxTargetPointCollection
 

Detailed Description

Definition at line 15 of file boxTarget.h.

Constructor & Destructor Documentation

◆ boxTarget() [1/2]

boxTarget::boxTarget ( const char *  Name,
Bool_t  Active 
)

Name : Detector Name Active: kTRUE for active detectors (ProcessHits() will be called) kFALSE for inactive detectors

◆ boxTarget() [2/2]

boxTarget::boxTarget ( )

default constructor

Definition at line 37 of file boxTarget.cxx.

38 : FairDetector("boxTarget", kTRUE, kVETO),
39 fTrackID(-1),
40 fVolumeID(-1),
41 fPos(),
42 fMom(),
43 fTime(-1.),
44 fLength(-1.),
45 fTargetMaterial("tungsten"),
46 fTargetL(1.),
47 fboxTargetPointCollection(new TClonesArray("vetoPoint"))
48{}
@ kVETO
TClonesArray * fboxTargetPointCollection
Definition boxTarget.h:88
Double_t fLength
time
Definition boxTarget.h:82
TString fTargetMaterial
track length
Definition boxTarget.h:83
TLorentzVector fPos
volume id
Definition boxTarget.h:79
Int_t fTrackID
Definition boxTarget.h:77
Int_t fVolumeID
track index
Definition boxTarget.h:78
Float_t fTargetL
target material
Definition boxTarget.h:84
Double_t fTime
momentum at entrance
Definition boxTarget.h:81
TLorentzVector fMom
position at entrance
Definition boxTarget.h:80

◆ ~boxTarget()

boxTarget::~boxTarget ( )
virtual

destructor

Definition at line 50 of file boxTarget.cxx.

51{
55 }
56}

Member Function Documentation

◆ AddHit()

vetoPoint * boxTarget::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 
)

Definition at line 206 of file boxTarget.cxx.

210{
211 TClonesArray& clref = *fboxTargetPointCollection;
212 Int_t size = clref.GetEntriesFast();
213 return new(clref[size]) vetoPoint(trackID, detID, pos, mom,
214 time, length, eLoss, pdgCode,Lpos,Lmom);
215}

◆ BeginEvent()

virtual void boxTarget::BeginEvent ( )
inlinevirtual

Definition at line 65 of file boxTarget.h.

65{;}

◆ BeginPrimary()

virtual void boxTarget::BeginPrimary ( )
inlinevirtual

Definition at line 62 of file boxTarget.h.

62{;}

◆ ConstructGeometry()

void boxTarget::ConstructGeometry ( )

Create the detector geometry

Definition at line 132 of file boxTarget.cxx.

133{
134 Double_t cm = 1; // cm
135 Double_t m = 100*cm; // m
136 Double_t mm = 0.1*cm; // mm
137
138 TGeoVolume *top=gGeoManager->FindVolumeFast("Detector");
139 if(!top) LOG(ERROR) << "no Detector volume found " ;
140
142 TGeoMedium *TargetMaterial = gGeoManager->GetMedium(fTargetMaterial);
143 InitMedium("vacuums");
144 TGeoMedium *vac = gGeoManager->GetMedium("vacuums");
145 InitMedium("Concrete");
146 TGeoMedium *concrete = gGeoManager->GetMedium("Concrete");
147 if (fBox){
148 // for studying absorption length
149 TGeoVolume* target = gGeoManager->MakeBox("Target",TargetMaterial,199.*cm,199.*cm,(fTargetL/2)*cm);
150 top->AddNode(target, 1, new TGeoTranslation(0, 0, (fTargetL/2)*cm));
151
152 TGeoVolume *sensPlane = gGeoManager->MakeBox("sensPlane",vac,199.*cm,199.*cm,1.*mm);
153 sensPlane->SetLineColor(kGreen);
154 top->AddNode(sensPlane, 13, new TGeoTranslation(0, 0, (fTargetL+ 0.11)*cm));
155 AddSensitiveVolume(sensPlane);
156 AddSensitiveVolume(target);
157 }else{
158 // more realistic setup of coldBox
159 // x1=-1140 mm, x2= 400 mm
160 // y1=100 mm, y2=1114 mm
161 // z1=50 mm, z2=1500mm
162 // 20cm x 2cm hole for SciFi cables
163 // 20cm x 20cm hole for cooling pipes
164 float x1=-1140;
165 float x2= 400;
166 float y1=100;
167 float y2=1114;
168 float z1=50;
169 float z2=1500;
170 float dx = (x2-x1)/10./2.;
171 float dy = (y2-y1)/10./2.;
172 float dz = (z2-z1)/10./2.;
173 float d = fTargetL;
174 float rAir = 10.;
175 float cablesDx = 10;
176 float cablesDy = 1;
177
178 TGeoBBox* box_I = new TGeoBBox("box_I",dx,dy,dz);
179 double origin[3] = {0,d/2,0};
180 TGeoBBox* box_O = new TGeoBBox("box_O",dx+d,dy+d/2.,dz+d,origin);
181 TGeoTube* hole_air = new TGeoTube("hole_air",0.,rAir,dx);
182 TGeoRotation* RF = new TGeoRotation("R_airTube");
183 RF->SetAngles(0.,0.,90.);
184 TGeoCombiTrans* CombiTrans1 = new TGeoCombiTrans("T_AirTube",-dx/2,0,10.,RF);
185 CombiTrans1->RegisterYourself();
186 TGeoCombiTrans* CombiTrans2 = new TGeoCombiTrans("T_CablesDuct",dx/2,0,10.,RF);
187 CombiTrans2->RegisterYourself();
188 TGeoBBox* hole_cables = new TGeoBBox("hole_cables",cablesDx,cablesDy,dx);
189 double floor_thickness = 25.;
190 double originF[3] = {0,-dy-floor_thickness,0};
191 TGeoBBox* floor = new TGeoBBox("floor",1.5*dx,floor_thickness,1.5*dz,originF);
192 TGeoCompositeShape* box = new TGeoCompositeShape("box","box_O-box_I-hole_air:T_AirTube-hole_cables:T_CablesDuct");
193 TGeoVolume *volBox = new TGeoVolume("vbox",box,TargetMaterial);
194 TGeoVolume *volFloor = new TGeoVolume("vfloor",floor,concrete);
195 TGeoVolume *sensBox = new TGeoVolume("sensBox",box_I,vac);
196 sensBox->SetLineColor(kGreen);
197 top->AddNode(sensBox, 13);
198 AddSensitiveVolume(sensBox);
199 top->AddNode(volBox, 1);
200 AddSensitiveVolume(volBox);
201 top->AddNode(volFloor, 2);
202 }
203
204}
Double_t cm
Double_t m
Double_t mm
Int_t InitMedium(TString name)
Definition boxTarget.cxx:99
Bool_t fBox
Definition boxTarget.h:86
origin(sTree, it)

◆ CopyClones()

virtual void boxTarget::CopyClones ( TClonesArray *  cl1,
TClonesArray *  cl2,
Int_t  offset 
)
inlinevirtual

The following methods can be implemented if you need to make any optional action in your detector during the transport.

Definition at line 56 of file boxTarget.h.

57 {;}

◆ EndOfEvent()

void boxTarget::EndOfEvent ( )
virtual

Definition at line 122 of file boxTarget.cxx.

123{
124
126
127}

◆ FinishPrimary()

virtual void boxTarget::FinishPrimary ( )
inlinevirtual

Definition at line 60 of file boxTarget.h.

60{;}

◆ FinishRun()

void boxTarget::FinishRun ( )
virtual

Definition at line 129 of file boxTarget.cxx.

129 {
130}

◆ GetCollection()

TClonesArray * boxTarget::GetCollection ( Int_t  iColl) const
virtual

Gets the produced collections

Definition at line 224 of file boxTarget.cxx.

225{
226 if (iColl == 0) { return fboxTargetPointCollection; }
227 else { return NULL; }
228}

◆ Initialize()

void boxTarget::Initialize ( )
virtual

Initialization of the detector is done here

Definition at line 117 of file boxTarget.cxx.

118{
119 FairDetector::Initialize();
120}

◆ InitMedium()

Int_t boxTarget::InitMedium ( TString  name)
private

Definition at line 99 of file boxTarget.cxx.

100{
101 static FairGeoLoader *geoLoad=FairGeoLoader::Instance();
102 static FairGeoInterface *geoFace=geoLoad->getGeoInterface();
103 static FairGeoMedia *media=geoFace->getMedia();
104 static FairGeoBuilder *geoBuild=geoLoad->getGeoBuilder();
105
106 FairGeoMedium *ShipMedium=media->getMedium(name);
107
108 if (!ShipMedium)
109 Fatal("InitMedium","Material %s not defined in media file.", name.Data());
110 TGeoMedium* medium=gGeoManager->GetMedium(name);
111 if (medium)
112 return ShipMedium->getMediumIndex();
113 return geoBuild->createMedium(ShipMedium);
114}

◆ PostTrack()

virtual void boxTarget::PostTrack ( )
inlinevirtual

Definition at line 63 of file boxTarget.h.

63{;}

◆ PreTrack()

virtual void boxTarget::PreTrack ( )
inlinevirtual

Definition at line 64 of file boxTarget.h.

64{;}

◆ ProcessHits()

Bool_t boxTarget::ProcessHits ( FairVolume *  v = 0)
virtual

this method is called for each step during simulation (see FairMCApplication::Stepping())

This method is called from the MC stepping

Definition at line 58 of file boxTarget.cxx.

59{
61 if ( gMC->IsTrackEntering() ) {
62 fTime = gMC->TrackTime() * 1.0e09;
63 fLength = gMC->TrackLength();
64 gMC->TrackPosition(fPos);
65 gMC->TrackMomentum(fMom);
66 LOG(DEBUG) <<"track enters"<<gMC->GetStack()->GetCurrentTrack()->GetPdgCode();
67 }
68 if (!gMC->IsTrackEntering() ) {
69 LOG(DEBUG) <<"track is not entering"<<gMC->GetStack()->GetCurrentTrack()->GetPdgCode();
70 }
71 if ( gMC->IsTrackExiting() ||
72 gMC->IsTrackStop() ||
73 gMC->IsTrackDisappeared() ) {
74 LOG(DEBUG) <<"track stopped"<<gMC->GetStack()->GetCurrentTrack()->GetPdgCode();
75
76 fTrackID = gMC->GetStack()->GetCurrentTrackNumber();
77 Int_t veto_uniqueId;
78 gMC->CurrentVolID(veto_uniqueId);
79
80 TParticle* p=gMC->GetStack()->GetCurrentTrack();
81 Int_t pdgCode = p->GetPdgCode();
82 TLorentzVector Pos;
83 gMC->TrackPosition(Pos);
84 TLorentzVector Mom;
85 gMC->TrackMomentum(Mom);
86 Double_t xmean = (fPos.X()+Pos.X())/2. ;
87 Double_t ymean = (fPos.Y()+Pos.Y())/2. ;
88 Double_t zmean = (fPos.Z()+Pos.Z())/2. ;
89 AddHit(fTrackID, veto_uniqueId, TVector3(xmean, ymean, zmean),
90 TVector3(fMom.Px(), fMom.Py(), fMom.Pz()), fTime, fLength,
91 0.,pdgCode,TVector3(Pos.X(), Pos.Y(), Pos.Z()),TVector3(Mom.Px(), Mom.Py(), Mom.Pz()) );
92 // Increment number of veto det points in TParticle
93 ShipStack* stack = (ShipStack*) gMC->GetStack();
94 stack->AddPoint(kVETO);
95 }
96 return kTRUE;
97}
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)

◆ Register()

void boxTarget::Register ( )
virtual

Registers the produced collections in FAIRRootManager.

Definition at line 217 of file boxTarget.cxx.

218{
219
220 FairRootManager::Instance()->Register("vetoPoint", "veto",
222}

◆ Reset()

void boxTarget::Reset ( )
virtual

has to be called after each event to reset the containers

Definition at line 229 of file boxTarget.cxx.

230{
232}

◆ SetSpecialPhysicsCuts()

virtual void boxTarget::SetSpecialPhysicsCuts ( )
inlinevirtual

Definition at line 58 of file boxTarget.h.

58{;}

◆ SetTarget()

void boxTarget::SetTarget ( TString  material,
Float_t  L,
Bool_t  choice 
)
inline

Member Data Documentation

◆ fBox

Bool_t boxTarget::fBox
private

Definition at line 86 of file boxTarget.h.

◆ fboxTargetPointCollection

TClonesArray* boxTarget::fboxTargetPointCollection
private

container for data points

Definition at line 88 of file boxTarget.h.

◆ fLength

Double_t boxTarget::fLength
private

time

Definition at line 82 of file boxTarget.h.

◆ fMom

TLorentzVector boxTarget::fMom
private

position at entrance

Definition at line 80 of file boxTarget.h.

◆ fPos

TLorentzVector boxTarget::fPos
private

volume id

Definition at line 79 of file boxTarget.h.

◆ fTargetL

Float_t boxTarget::fTargetL
private

target material

Definition at line 84 of file boxTarget.h.

◆ fTargetMaterial

TString boxTarget::fTargetMaterial
private

track length

Definition at line 83 of file boxTarget.h.

◆ fTime

Double_t boxTarget::fTime
private

momentum at entrance

Definition at line 81 of file boxTarget.h.

◆ fTrackID

Int_t boxTarget::fTrackID
private

Track information to be stored until the track leaves the active volume.

Definition at line 77 of file boxTarget.h.

◆ fVolumeID

Int_t boxTarget::fVolumeID
private

track index

Definition at line 78 of file boxTarget.h.

◆ index

Int_t boxTarget::index
private

target length

Definition at line 85 of file boxTarget.h.


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