/[PAMELA software]/PamVMC_update/trk/include/PamVMCTrkSD.h
ViewVC logotype

Contents of /PamVMC_update/trk/include/PamVMCTrkSD.h

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.1.1.1 - (show annotations) (download) (vendor branch)
Tue Oct 15 15:51:27 2013 UTC (11 years, 1 month ago) by formato
Branch: MAIN, rel
CVS Tags: reltag, HEAD
Changes since 1.1: +0 -0 lines
File MIME type: text/plain
PamVMC update

1 #ifndef PAMVMCTRKSD_H
2 #define PAMVMCTRKSD_H
3
4 #ifdef __cplusplus
5 extern "C" {
6 #endif
7
8 #define GPINIT gpinit_
9 #define GPUDIFFUSION gpudiffusion_
10 #define GPBEGIN gpbegin_
11 #define GPDSPE gpdspe_
12 #define GPDCSPE gpdcspe_
13
14 extern void GPINIT(float*, float*);
15
16 extern void GPBEGIN();
17
18 extern void GPUDIFFUSION(int* ,float* ,int*, float* ,float* ,int*);
19
20 extern void GPDSPE();
21
22 extern void GPDCSPE(float*, float*, float*, float*, float*, float* );
23
24 #ifdef __cplusplus
25 }
26 #endif
27
28
29
30 #include <iostream>
31
32 #include "PamVMCDetectorSD.h"
33 #include "PamVMCTrkHit.h"
34 #include "PamVMCTrkF77GpsSpe.h"
35 #include "PamVMCTrkF77GpcSpe.h"
36 #include "PamVMCTrkDig.h"
37 #include <TDatabasePDG.h>
38 #include <TGeoManager.h>
39 #include <TSystem.h>
40 #include <TFile.h>
41 #include <TH1F.h>
42
43
44 class pGPSSPEData;
45
46
47 class PamVMCTrkSD: public PamVMCDetectorSD{
48
49 private:
50 cGPCSPE * fGPCSPE;
51 cGPSSPE * fGPSSPE;
52 pGPSSPEHits * fgps;
53
54 public:
55 PamVMCTrkSD():PamVMCDetectorSD("PamVMCTrkHit","TSPA",1000)
56 {
57 fGPCSPE = 0;
58 fGPSSPE = 0;
59 fgps = new pGPSSPEHits();
60 Init();
61 };
62
63 virtual ~PamVMCTrkSD(){ delete fgps;}
64
65 //Loading SPE hit resolution data from resxy_new.root
66 //which should present in the same directory then application library
67
68 void Init(){
69 TString PAM_VMC=gSystem->Getenv("PAM_VMC");
70 TString PLATFORM=gSystem->Getenv("PLATFORM");
71 TFile * resspe =
72 new TFile(PAM_VMC+"/lib/tgt_"+PLATFORM+"/resxy_new.root");
73 if(resspe->IsOpen()){
74 TH1F * hsx;
75 TH1F * hsy;
76 TString num;
77 float cwx[11][1000];
78 float cwy[11][1000];
79 for(Int_t j=0; j<11; j++){
80 num=(::Form("%4i",1000+2*j));
81 resspe->GetObject("h"+num+";1",hsx);
82 num=(::Form("%4i",2000+2*j));
83 resspe->GetObject("h"+num+";1",hsy);
84 if ((hsx)&&(hsy)){
85 for(Int_t i = 0; i<1000; i++){
86 cwx[j][i]=hsx->GetBinContent(i);
87 cwy[j][i]=hsy->GetBinContent(i);
88 }
89 }else{
90 cout<<"!!!WARNING: SPE hit resolution data file corrupted!!!"<<endl
91 <<"!!!Application terminated..."<<endl;
92 // exit(-1);
93 }
94 }
95
96 gpinit_((*cwx), (*cwy));
97
98 delete hsx;
99 delete hsy;
100 } else cout<<"!!!WARNING: SPE hit resolution data file resxy_new.root"
101 <<" is not present in: " <<PAM_VMC+"/lib/tgt_"+PLATFORM
102 <<" or corrupted!!!"<<endl;
103
104 delete resspe;
105 }
106
107 /*This method for zeroing PROXTANTI, PROYTANTI, GLOBSTRIPX, GLOBSTRIPY
108 arrays in fortran common. Calls before each event */
109 virtual void PreEvent(){ gpbegin_(); }
110
111 void Register(){
112
113 PamRootManager::Instance()->
114 Register(fdname.Data(),"TClonesArray", &fHitColl);
115 TString a=fHitColl?"exist":"not exist";
116 PamRootManager::Instance()->
117 Register("GPSSPE","pGPSSPEHits",&fgps);
118
119 if (fdig) ((PamVMCTrkDig*)fdig)->RegisterTrkData(fgps);
120
121 }
122
123 virtual void ClearHitColl(){
124 fnohit=0;
125 (*fHitColl).Clear("C");
126 fgps->Clear();
127 }
128
129
130 virtual void Compress(){fgps->Compress(); (*fHitColl).Compress();}
131
132 virtual void SaveTrkHit(const char * dname){
133 PamVMCDetectorHit * t=(PamVMCDetectorHit*) CreateTrkHit(dname);
134 *t=*fhit;
135 }
136
137 virtual PamVMCTrkHit * CreateTrkHit(const char * dname)
138 {
139 fdetID=(PamVMCGeoID*)fdmap(dname);
140 return CreateTrkHit();
141 }
142
143 virtual PamVMCTrkHit * CreateTrkHit()
144 {
145 PamVMCTrkHit *hit= (PamVMCTrkHit *) fHitColl->New(fnohit++);
146 InitHit();
147 SetTrkHit(hit);
148 return hit;
149 }
150
151 void FillHit(fin f,TVirtualMC *g, Bool_t is_prim){
152 ffi=f;
153 //g->SetMaxStep(0.001); // 10 um
154 switch(f) {
155 case ENTERING:
156 CleanHit();
157 FillVolID();
158 InitHit();
159 //cout<<">>>ENTERING in PLANE: "<<GetPlaneID()<<endl;
160 default:
161 UpdateHit(g, is_prim);
162 //TParticlePDG* particlePDG = TDatabasePDG::Instance()->GetParticle(gMC->TrackPid());
163 //Double_t mass = particlePDG->Mass();
164 //Double_t x[3];
165 //Double_t B[3];
166 //Double_t * xp = &x[0];
167 //Double_t *Bp = &B[0];
168 //g->TrackPosition(x[0],x[1],x[2]);
169 //g->TrackMomentum(x,y,z,e);
170 //Double_t P0 = TMath::Sqrt(x*x+y*y+z*z);
171 //Double_t energy = TMath::Sqrt(P0*P0+mass*mass)-mass;
172 //cout<<"stepping... Stepsize"<<gMC->TrackStep()*10000<<"um "<<" particle PDG is:"<<gMC->TrackPid()<<" Kin.energy:"<<energy*1000000<<" keV"<<endl;
173 //PamVMCFieldMgr::Instance()->Field(xp,Bp);
174 //cout<<"X,Y,Z: "<<x[0]<<" "<<x[1]<<" "<<x[2]<<" "<<"Bx, By, Bz: "<<B[0]<<" "<<B[1]<<" "<<B[2]<<endl;
175 //fhit->Print();
176 //cout<<"Edep"<<g->Edep()*1e6<<" kev"<<endl;
177 if (g->Edep()) PreDigit();
178
179 break;
180 }
181
182 switch(f){
183 case EXITING:
184 //cout<<"<<<EXITING out"<<gMC->TrackPid()<<endl;
185 // Save hit if energy release is greater than zero
186 //cout<<"PDG:"<<fhit->GetPDG()<<"TOF="<<fhit->GetTOF()<<endl;
187
188 if(fhit->GetEREL() ){
189 //cout<<"Saving.."<<endl;
190 SaveTrkHit(fdname.Data());
191 }
192 break;
193 default:
194 break;
195 }
196 }
197
198
199
200 /*This method for simulation spatial tracker resolution, use data from
201 resxy_new.root file. Calls on track exiting when we know XOUT,YOUT,ZOUT */
202 void SetTrkHit(PamVMCTrkHit *t){
203 float xin = (float)fhit->GetXIN();
204 float yin = (float)fhit->GetYIN();
205 float zin = (float)fhit->GetZIN();
206 float xout = (float)fhit->GetXOUT();
207 float yout = (float)fhit->GetYOUT();
208 float zout = (float)fhit->GetZOUT();
209
210 gpdcspe_(&xin, &yin, &zin, &xout, &yout, &zout);
211
212 fGPCSPE=&gpcspe_;
213
214 // cout<<"GPCSPE"<<fGPCSPE->xave<<" "<<fGPCSPE->yave<<" "<<
215 // fGPCSPE->zave<<" " <<fGPCSPE->nxmult<<" "<<fGPCSPE->nymult<<endl;
216
217 fGPCSPE->SetHit(t);
218 }
219
220
221 /*This method designed for calculaton energy release on trakcer strips
222 is calls for each step*/
223
224 virtual void PreDigit(){
225
226 Double_t x,y,z,px,py,pz,e,p;
227
228 //cout<<"Calling pre-digit.."<<endl;
229
230 gMC->TrackMomentum(px,py,pz,e);
231
232
233 p=sqrt(px*px+py*py+pz*pz);
234
235
236 float *trpar = new float[6];
237
238 if(p){ //sometimes FLUKA vives u p=0
239 trpar[3]=(float)px/p;
240 trpar[4]=(float)py/p;
241 trpar[5]=(float)pz/p;
242 } else {
243 trpar[3]=0.;
244 trpar[4]=0.;
245 trpar[5]=0.;
246 }
247
248 gGeoManager->cd(gMC->CurrentVolPath());
249
250 int iact=ffi;
251
252 gMC->TrackPosition(x,y,z);
253
254 trpar[0]=(float)x;
255 trpar[1]=(float)y;
256 trpar[2]=(float)z;
257
258
259
260 int *numvol = new int[2];
261 numvol[0]=GetPlaneID();
262 //if (numvol[0]==6) numvol[0]=0;
263 numvol[1]=GetPadID();
264
265 float deloss=(float)gMC->Edep();
266
267 float step=(float)gMC->TrackStep();
268
269 int itypar=TDatabasePDG::Instance()->ConvertPdgToGeant3(gMC->TrackPid());
270
271 gpudiffusion_(&iact, trpar, numvol, &deloss, &step, &itypar);
272
273 delete numvol;
274 delete trpar;
275
276
277 //cout<<"Exiting pre-digit"<<endl;
278 }
279
280
281
282 /*This method calculate ADC amplitudes from strips. It use data arrays
283 PROXTANTI, PROYTANTI, GLOBSTRIPX, GLOBSTRIPY in common. Calls after
284 each event from application */
285
286 virtual void Digitize(){
287
288 gpdspe_();
289
290 fGPSSPE = &gpsspe_;
291
292 fgps->SetData(fGPSSPE);
293 }
294
295
296
297 Int_t GetPlaneID(){
298 Int_t pl=Int_t((fhit->GetPOS()-1)/6.)+1;
299 return pl;
300 }
301
302 Int_t GetPadID(){
303 return (fhit->GetPOS()-(GetPlaneID()-1)*6);
304 }
305
306 ClassDef(PamVMCTrkSD,1)
307
308 };
309
310 #endif // PAMVMCTRKSD_H
311

  ViewVC Help
Powered by ViewVC 1.1.23