SND@LHC Software
Loading...
Searching...
No Matches
boxTarget.cxx
Go to the documentation of this file.
1#include "boxTarget.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 "TGeoTube.h"
24#include "TGeoCompositeShape.h"
25#include "TGeoBoolNode.h"
26#include "TGeoMaterial.h"
27#include "TParticle.h"
28#include "TROOT.h"
29#include "TH1D.h"
30#include "TH2D.h"
31#include "TDatabasePDG.h"
32
33#include <iostream>
34using std::cout;
35using std::endl;
36
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{}
49
57
58Bool_t boxTarget::ProcessHits(FairVolume* vol)
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}
98// ----- Private method InitMedium
99Int_t boxTarget::InitMedium(TString name)
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}
115
116
118{
119 FairDetector::Initialize();
120}
121
123{
124
126
127}
128
131
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}
205
206vetoPoint* boxTarget::AddHit(Int_t trackID, Int_t detID,
207 TVector3 pos, TVector3 mom,
208 Double_t time, Double_t length,
209 Double_t eLoss, Int_t pdgCode,TVector3 Lpos, TVector3 Lmom)
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}
216
218{
219
220 FairRootManager::Instance()->Register("vetoPoint", "veto",
222}
223
224TClonesArray* boxTarget::GetCollection(Int_t iColl) const
225{
226 if (iColl == 0) { return fboxTargetPointCollection; }
227 else { return NULL; }
228}
230{
232}
233
Double_t cm
Double_t m
Double_t mm
@ kVETO
virtual TClonesArray * GetCollection(Int_t iColl) const
Int_t InitMedium(TString name)
Definition boxTarget.cxx:99
TClonesArray * fboxTargetPointCollection
Definition boxTarget.h:88
Double_t fLength
time
Definition boxTarget.h:82
virtual void FinishRun()
TString fTargetMaterial
track length
Definition boxTarget.h:83
virtual ~boxTarget()
Definition boxTarget.cxx:50
virtual Bool_t ProcessHits(FairVolume *v=0)
Definition boxTarget.cxx:58
TLorentzVector fPos
volume id
Definition boxTarget.h:79
Int_t fTrackID
Definition boxTarget.h:77
virtual void Reset()
Bool_t fBox
Definition boxTarget.h:86
virtual void Initialize()
Float_t fTargetL
target material
Definition boxTarget.h:84
virtual void Register()
void ConstructGeometry()
virtual void EndOfEvent()
Double_t fTime
momentum at entrance
Definition boxTarget.h:81
TLorentzVector fMom
position at entrance
Definition boxTarget.h:80
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)
ClassImp(ecalContFact) ecalContFact