]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/Upgrade/AliToyMCDrawer.cxx
adding name strings for canvases and histograms, adding z0 resolution
[u/mrichter/AliRoot.git] / TPC / Upgrade / AliToyMCDrawer.cxx
1 #include <iostream>
2 #include <AliMagF.h>
3 #include <TChain.h>
4 #include <TMath.h>
5 #include <TTree.h>
6 #include <TFile.h>
7 #include <TGraph2D.h>
8 #include <TH3F.h>
9 #include <TClonesArray.h>
10 #include <TPolyLine3D.h>
11 #include <TPolyMarker3D.h>
12 #include <TAxis.h>
13 #include <TLegend.h>
14 #include <TPad.h>
15 #include <AliTrackPointArray.h>
16
17 #include "AliToyMCDrawer.h"
18 #include "AliToyMCEvent.h"
19 #include "AliToyMCTrack.h"
20 #include <AliCDBManager.h>
21 #include <AliTPCParam.h>
22 #include <AliGeomManager.h>
23 #include <TGeoGlobalMagField.h>
24 #include <AliTPCcalibDB.h>
25 #include <AliTPCROC.h>
26
27
28 ClassImp(AliToyMCDrawer);
29
30 AliToyMCDrawer::AliToyMCDrawer()
31  : TObject()   
32  ,fInputTree(0x0)
33  ,inFile(0x0)
34  ,fFileName(0x0)
35  ,fEvent(0x0)
36  ,fEventArray(0x0)
37  ,fDispHist(0x0)
38  ,fCenterTime(-1.)
39  ,fDriftVel(-1.)
40    
41  {
42    fEventArray = new TClonesArray("AliToyMCEvent");
43    
44    fTPCParam = AliTPCcalibDB::Instance()->GetParameters();
45    fTPCParam->ReadGeoMatrices();
46    fDriftVel = fTPCParam->GetDriftV();///1000000;
47    fMaxZ0=fTPCParam->GetZLength();
48    fTimeRange = 2*fMaxZ0/fDriftVel;
49    fIFCRadius=  83.5;
50    fOFCRadius= 254.5;   
51  
52    fRoc  = AliTPCROC::Instance();
53    fPoints = new TClonesArray("TPolyMarker3D");
54    fDistPoints = new TClonesArray("TPolyMarker3D");
55    
56  }
57 //________________________________________________________________
58 AliToyMCDrawer::AliToyMCDrawer(const AliToyMCDrawer &drawer)
59   : TObject(drawer)
60   ,fInputTree(drawer.fInputTree)
61   ,inFile(drawer.inFile)
62   ,fFileName(drawer.fFileName)
63   ,fEvent(drawer.fEvent)
64   ,fEventArray(drawer.fEventArray)
65   ,fDispHist(drawer.fDispHist)
66   , fCenterTime(drawer.fCenterTime)
67   ,fDriftVel(drawer.fDriftVel)
68   ,fTPCParam(drawer.fTPCParam)
69   ,fMaxZ0(drawer.fMaxZ0)
70   ,fOFCRadius(drawer.fOFCRadius)
71   ,fIFCRadius(drawer.fIFCRadius)
72   ,fTimeRange(drawer.fTimeRange)
73   ,fRoc(drawer.fRoc)
74   ,fPoints(drawer.fPoints)
75   ,fDistPoints(drawer.fDistPoints)
76 {
77   //
78 }
79 //_____________________________________________________
80 AliToyMCDrawer& AliToyMCDrawer::operator = (const AliToyMCDrawer &drawer)
81 {
82   //assignment operator
83   if (&drawer == this) return *this;
84   new (this) AliToyMCDrawer(drawer);
85
86   return *this;
87 }
88
89 //________________________________________________________________
90 AliToyMCDrawer::~AliToyMCDrawer()
91 {
92   //destructor
93   delete fEvent;
94   delete fEventArray;
95   delete fDispHist;
96   delete fPoints;
97   delete fDistPoints; 
98 }
99 //________________________________________________________________
100 Int_t AliToyMCDrawer::FillEventArray(Double_t snapShotTime)
101 {
102   if(!fFileName) {
103     std::cout << "no input file provided, using default (toyMC.root)" << std::endl;
104     fFileName = "toyMC.root";
105   }
106   
107   inFile = new TFile(fFileName,"read");
108   fInputTree = dynamic_cast<TTree*> (inFile->Get("toyMCtree"));
109   fInputTree->SetBranchAddress("event",&fEvent);
110   fEventArray->Clear();
111   fCenterTime = snapShotTime;
112   Int_t leftSearchIndex  = Int_t(snapShotTime/2e-5);
113   if(leftSearchIndex>fInputTree->GetEntries()) leftSearchIndex = fInputTree->GetEntries()-1;
114   if(leftSearchIndex<0) leftSearchIndex = 0;
115   
116   fInputTree->GetEvent(leftSearchIndex);
117   Double_t leftTime = fEvent ->GetT0();
118   Int_t firstEventIndex;
119   //if fEvent exactly at snapShotTime, check if there are more events at the same time but with lower indices
120   if( TMath::Abs(leftTime - snapShotTime)<5e-9) 
121     {
122       //  printf("close from start, lefttime = %2.2f\n",leftTime);
123       firstEventIndex = leftSearchIndex+1;
124       while( TMath::Abs(leftTime- snapShotTime) < 5e-9){
125         firstEventIndex--;
126         //      printf("close from start, lefttime = %2.2f\n",leftTime);
127         if(fInputTree->GetEvent(firstEventIndex-1))
128           {
129             leftTime=fEvent->GetT0();
130             //     printf("close from start, lefttime = %2.2f\n",leftTime);
131           }
132         else break;
133       }
134
135
136     }
137   else
138     {
139       Int_t rightSearchIndex = leftSearchIndex;
140       Double_t rightTime;
141       
142       Int_t direction = 1; //go right
143       if(leftTime > snapShotTime) {
144         direction = -1;
145       }
146       rightSearchIndex+= direction;
147       fInputTree->GetEvent(rightSearchIndex);
148       rightTime = fEvent->GetT0();
149       //std::cout << "bef search " << leftTime << " " << rightTime<< " " << " " << snapShotTime << " "<< (leftTime-snapShotTime)*(rightTime-snapShotTime)  << std::endl;
150       while ( (leftTime-snapShotTime)*(rightTime-snapShotTime)>0   ){
151         //              printf("searching, leftime %f, righttim %f\n",leftTime,rightTime);
152         rightSearchIndex += direction; 
153         leftSearchIndex +=direction;
154         fInputTree->GetEvent(leftSearchIndex);
155         leftTime = fEvent->GetT0();
156         fInputTree->GetEvent(rightSearchIndex);
157         rightTime = fEvent->GetT0();
158         //      printf("searching, leftime %f, righttim %f\n",leftTime,rightTime);
159       }
160       if (direction==-1) rightSearchIndex = leftSearchIndex;
161       firstEventIndex = rightSearchIndex;
162
163     }
164   fInputTree->GetEvent(firstEventIndex);
165   //std::cout <<"first event after sn time: " << fEvent->GetT0() << std::endl;
166   return FillEventArray(firstEventIndex);
167   
168 }
169 //________________________________________________________________
170 Int_t AliToyMCDrawer::FillEventArray(Int_t middleEventNbr)
171 {
172   if(!fFileName) {
173     std::cout << "no input file provided, using default (toyMC.root)" << std::endl;
174     fFileName = "toyMC.root";
175   }
176   if(!inFile) inFile = new TFile(fFileName,"read");
177   if(!fInputTree) 
178     {
179       fInputTree = dynamic_cast<TTree*> (inFile->Get("toyMCtree"));
180       //   fInputTree->Add(fileName);
181       fInputTree->SetBranchAddress("event",&fEvent);
182       if(fInputTree->GetEvent(middleEventNbr)) fCenterTime = fEvent->GetT0();
183     }
184   Double_t centerEventTime;
185   if(fInputTree->GetEvent(middleEventNbr)) 
186     {
187       centerEventTime = fEvent->GetT0();
188       
189     }
190   else return 0;
191   fEventArray->Clear();
192   if(fInputTree->GetEvent(middleEventNbr)){
193     new((*fEventArray)[fEventArray->GetEntriesFast()]) AliToyMCEvent(*fEvent);
194     
195     //add events after middle event
196     Int_t eventIndex = middleEventNbr + 1;
197     Double_t currentTime = centerEventTime;
198     if(fInputTree->GetEvent(eventIndex)) {
199       currentTime = fEvent->GetT0();
200       while(currentTime - centerEventTime < fTimeRange/2) {
201         new((*fEventArray)[fEventArray->GetEntriesFast()]) AliToyMCEvent(*fEvent);
202         eventIndex++;
203         if(!fInputTree->GetEvent(eventIndex)) break;
204         currentTime = fEvent->GetT0();
205       }
206     }
207     //add events before middle event
208     eventIndex = middleEventNbr - 1;
209     if(fInputTree->GetEvent(eventIndex)){
210       currentTime = fEvent->GetT0();
211       while(centerEventTime-currentTime < fTimeRange/2) {
212       
213         new((*fEventArray)[fEventArray->GetEntriesFast()]) AliToyMCEvent(*fEvent);
214         eventIndex--;
215         if(!fInputTree->GetEvent(eventIndex)) break;
216         currentTime = fEvent->GetT0();
217       }
218     }
219     return 1;
220   }
221   else {
222     printf("Selected event number (%d) out of range!\n",middleEventNbr);
223     return 0;
224   }
225   //return fEventArray;
226   inFile->Close();
227 }
228 //________________________________________________________________
229 void AliToyMCDrawer::DrawEvents(Bool_t both, Bool_t before)
230 {
231   
232   fPoints->Clear();
233   fDistPoints->Clear();
234   DrawGeometry();
235  
236   
237    Double_t phiCut = 1.;
238    TLegend *leg = new TLegend(0.1,0.75,0.6,1,Form("Snapshot of TPC at %f micros after sim start, freq. = , bunchcr.rate = .",1000000*fCenterTime));
239    leg->AddEntry((TObject*)0,Form("%2.2f<#phi<#pi-%2.2f && #pi+%2.2f<#phi<2#pi-%2.2f",phiCut,phiCut,phiCut,phiCut),"");
240
241    
242  
243
244   for(Int_t iEvent = 0; iEvent<fEventArray->GetEntriesFast(); iEvent++){
245  
246  
247      AliToyMCEvent *currentEvent = static_cast<AliToyMCEvent*> (fEventArray->At(iEvent));
248      Double_t currentEventTime = currentEvent->GetT0();
249      
250      if((1000000*currentEventTime)>((1000000*fCenterTime)+0.05   ))  {
251        printf("larger iEvent: %d, current %f, center %f\n",iEvent, currentEventTime,fCenterTime);
252        printf("after\n");
253        if ( !(both || !before)) continue;
254
255      }
256      if( currentEventTime<=fCenterTime)// && !(currentEventTime==fCenterTime &&both   ))
257        { 
258          printf("smaller iEvent: %d, current %f, center %f\n",iEvent, currentEventTime,fCenterTime);
259          printf("before\n");
260          if ( !(both || before)) continue;
261        }
262      
263      //std:: cout << "iev " << iEvent << " " << color << std::endl;
264      Double_t currentEventTimeRelToCentral = fCenterTime - currentEventTime;
265      
266      Int_t color = ((currentEvent->GetEventNumber())%10);
267      if(color == 0||color ==1) {
268        color==0?color=2:color=3;
269        // if(iEvent==0) color = 8;
270        // if(iEvent==1) color = 9;
271        
272      }
273      
274       char bef[] = "before snapshot time";
275       char aft[] = "after snapshot time";
276       char *when; 
277       Int_t sign = 1;
278       if( 1000000*(currentEventTime-fCenterTime) > 0.5) {when = aft;
279         sign = -1;
280           }
281       else when = bef;
282       TGraph *temp = new TGraph();
283       temp->SetName(Form("temp%d",currentEvent->GetEventNumber()));
284        temp->SetLineColor(color);
285        temp->SetLineWidth(10);
286        
287        leg->AddEntry(temp,Form(" %2.2f microseconds %s (global time %2.2f) at z_vertex = %2.2f, centr = ., eventNr: %d", 1000000*sign*currentEventTimeRelToCentral, when, 1000000*currentEvent->GetT0(),  currentEvent->GetZ() ,currentEvent->GetEventNumber() ),"l");
288      
289        
290        DrawEvent(currentEvent, fCenterTime, color);
291        
292
293        
294        
295   
296    }
297    leg->Draw("same");
298
299 }
300 //________________________________________________________________
301
302 void AliToyMCDrawer::DrawEvent(AliToyMCEvent *currentEvent, Double_t centerTime, Int_t color){
303   if(!fDispHist) DrawGeometry();
304   
305   for(Int_t iTrack = 0; iTrack < currentEvent->GetNumberOfTracks(); iTrack++){
306     
307     //Double_t alpha = TMath::ATan2(currentEvent->GetTrack(iTrack)->GetSpacePoint(0)->GetY(),currentEvent->GetTrack(iTrack)->GetSpacePoint(0)->GetX());
308     // if( (     (  !(abs(alpha)>phiCut && abs(alpha) < TMath::Pi() - phiCut)  ) )      )continue;
309     // if(alpha<0) {
310     //  continue;
311     //  delete tempTrack;
312     // }
313     
314     //if(currentEvent->GetTrack(iTrack)->GetSpacePoint(0)->GetY()<0) alpha +=TMath::Pi();
315     //if( (       (alpha>phiCut && alpha < TMath::Pi() - phiCut) || ( alpha>TMath::Pi() + phiCut && alpha < TMath::TwoPi() - phiCut))   )continue;
316     
317     // if(alpha > 0) continue;
318
319     //std::cout << TMath::ASin(tempTrack->GetParameter()[2]) << std::endl;
320     // Double_t x = currentEvent->GetX();
321     // Double_t y = currentEvent->GetY();
322       // Double_t hlxpar[6];
323       // currentEvent->GetTrack(iTrack)->GetHelixParameters(hlxpar,0);
324       // Double_t trackPhi = hlxpar[2];
325       // Double_t phiCut = 1.47;
326       // std::cout << "bef phi cut "<< trackPhi<< std::endl;
327       
328       // std::cout << "aft phi cut " << std::endl;
329       // Double_t trackEta = currentEvent->GetTrack(iTrack)->GetEta();
330       //  Double_t z = currentEvent->GetZ();
331
332      
333      
334       
335       
336       const AliToyMCTrack *tempTrack = currentEvent->GetTrack(iTrack);
337       //AliToyMCTrack *tempTrack = new AliToyMCTrack(*currentEvent->GetTrack(iTrack));
338       DrawTrack(tempTrack, centerTime, currentEvent->GetT0(),color);
339       
340
341
342
343      
344
345       
346
347
348        
349     }
350
351
352
353
354 }
355
356 //________________________________________________________________
357 void AliToyMCDrawer::DrawTrack(const AliToyMCTrack *track,  Double_t centerTime, Double_t currentEventTime, Int_t color){
358   if(!fDispHist) DrawGeometry();
359   
360   TPolyMarker3D *disttrackpoints = new((*fDistPoints)[fDistPoints->GetEntriesFast()]) TPolyMarker3D();
361   TPolyMarker3D *trackpoints = new((*fPoints)[fPoints->GetEntriesFast()]) TPolyMarker3D();
362
363   Double_t currentEventTimeRelToCentral = centerTime - currentEventTime;
364   Int_t nDistPoints = track->GetNumberOfDistSpacePoints();
365   Int_t nPoints = track->GetNumberOfSpacePoints();
366     
367   for(Int_t iPoint = 0; iPoint< (nDistPoints<nPoints?nPoints:nDistPoints);iPoint++){
368     if(iPoint<nPoints) {
369             
370
371       Double_t xp = track->GetSpacePoint(iPoint)->GetX();
372       Double_t yp = track->GetSpacePoint(iPoint)->GetY();
373       Double_t zp = track->GetSpacePoint(iPoint)->GetZ();
374       Double_t zDrifted =  zp+(zp/TMath::Abs(zp))*currentEventTimeRelToCentral * fDriftVel;
375       //Double_t zDrifted = (zp/TMath::Abs(zp))*fMaxZ0  -(zp/TMath::Abs(zp))* fDriftVel*(track->GetSpacePoint(iPoint)->GetTimeBin()      -centerTime   );
376       Float_t xyzp[3] = {xp,yp,zp};
377       AliTrackPoint p;
378       p.SetXYZ(xyzp);
379       Float_t tempcov[6] = {0};
380       p.SetCov(tempcov);
381       Int_t sec = track->GetSpacePoint(iPoint)->GetDetector();
382       Double_t angle=((sec%18)*20.+10.)/TMath::RadToDeg();
383       AliTrackPoint prot = p.Rotate(-angle);
384       xp = prot.GetX();
385       yp = prot.GetY();
386
387
388             
389       if(track->GetSpacePoint(iPoint)->GetRow()!=255) {
390         if(TMath::Abs(zDrifted)<fMaxZ0 && zDrifted*zp >=0 /*&&TMath::Sqrt(xp*xp + yp*yp)<fIFCRadius*/) trackpoints->SetNextPoint(zDrifted,xp,yp);
391               
392       }
393       else std::cout << "row == " << track->GetSpacePoint(iPoint)->GetRow() << std::endl;
394       
395     }
396     if(iPoint<nDistPoints) {
397       Double_t xpdist = track->GetDistortedSpacePoint(iPoint)->GetX();
398       Double_t ypdist = track->GetDistortedSpacePoint(iPoint)->GetY();
399       Double_t zpdist = track->GetDistortedSpacePoint(iPoint)->GetZ();
400       //std::cout << zpdist << std::endl;
401             
402       Float_t xyzpdist[3] = {xpdist,ypdist,zpdist};
403       AliTrackPoint pdist;
404       pdist.SetXYZ(xyzpdist);
405       Float_t tempcovdist[6] = {0};
406       pdist.SetCov(tempcovdist);
407       Int_t secdist = track->GetDistortedSpacePoint(iPoint)->GetDetector();
408       Double_t angledist=((secdist%18)*20.+10.)/TMath::RadToDeg();
409       AliTrackPoint protdist = pdist.Rotate(-angledist);
410       xpdist = protdist.GetX();
411       ypdist = protdist.GetY();
412       
413       
414       UInt_t sector = track->GetDistortedSpacePoint(iPoint)->GetDetector();
415       UInt_t row = track->GetDistortedSpacePoint(iPoint)->GetRow();
416       UInt_t pad  = track->GetDistortedSpacePoint(iPoint)->GetPad();
417       
418       Int_t nPads = fTPCParam->GetNPads(sector,row);
419       Int_t intPad = TMath::Nint(Float_t(pad+nPads/2));
420
421       Float_t xyz[3];
422       //std::cout <<"sector: " << sector << " row: " << row << " pad: " <<pad<< std::endl;
423       fRoc->GetPositionGlobal(sector,row,intPad,xyz);
424       
425       Double_t zDrifteddist = (zpdist/TMath::Abs(zpdist))*fMaxZ0  -(zpdist/TMath::Abs(zpdist))* fDriftVel*(track->GetDistortedSpacePoint(iPoint)->GetTimeBin()      -centerTime   );
426       if(row!=255){
427               
428         if(TMath::Abs(zDrifteddist)<fMaxZ0 && zDrifteddist*zpdist>=0 /*&&TMath::Sqrt(xpdist*xpdist + ypdist*ypdist)<fIFCRadius*/) {
429                 
430           disttrackpoints->SetNextPoint(zDrifteddist,xpdist,ypdist);
431           //if(TMath::Sqrt(xpdist*xpdist + ypdist*ypdist)<fIFCRadius) std::cout << "fMaxZ0 " << fMaxZ0 <<" inside " << xpdist << " " << "zpdist "  << zpdist << " " << "zDrifteddist "<< zDrifteddist << " " << zDrifteddist*zpdist << std::endl;
432         }
433       }
434       else std::cout << "row == " << row << std::endl;
435       //Double_t zDrifteddist =(zpdist/TMath::Abs(zpdist))*fMaxZ0  -(zpdist/TMath::Abs(zpdist))*(currentEvent->GetTrack(iTrack)->GetDistortedSpacePoint(iPoint)->GetTimeBin()- currentEvent->GetT0() )* fDriftVel;
436           
437       
438     } 
439     //  if( (       (trackPhi>phiCut && trackPhi < TMath::Pi() - phiCut) || ( trackPhi>TMath::Pi() + phiCut && trackPhi < TMath::TwoPi() - phiCut))   ) {
440     
441     
442
443
444
445   }
446   if(1){
447     if(trackpoints && trackpoints->GetN()>0) {
448       //   trackpoints->SetMarkerColor(1+currentEvent->GetEventNumber()%9);
449       //trackpoints->SetMarkerStyle(7);
450       trackpoints->Draw("same");
451     }
452     if(disttrackpoints && disttrackpoints->GetN()>0) {
453       //  
454       
455       disttrackpoints->SetMarkerColor(color);
456
457       
458       disttrackpoints->Draw("same");
459       
460     }
461     
462       }
463   
464   
465   
466
467 }
468
469 //________________________________________________________________
470 void AliToyMCDrawer::DrawGeometry() {
471
472   
473   //delete fDispGraph;
474   //fDispGraph = new TGraph2D();
475   delete fDispHist;
476  
477   fDispHist = new TH3F("fDispHist","",100,-fMaxZ0, fMaxZ0, 100,-(fOFCRadius +10), fOFCRadius +10,100,-(fOFCRadius +10), fOFCRadius +10);
478   //if(!fDispGraph) fDispGraph = new TGraph();
479
480   //fDispGraph->Clear();
481   fDispHist->SetStats(0);
482   fDispHist->GetXaxis()->SetTitle("z [cm]");
483   fDispHist->GetYaxis()->SetTitle("x [cm]");
484   fDispHist->GetZaxis()->SetTitle("y [cm]");
485   fDispHist->Draw();
486   gPad->SetPhi(0);
487   gPad->SetTheta(0);
488   TPolyLine3D *endCap1 = new TPolyLine3D();
489   TPolyLine3D *endCap2 = new TPolyLine3D();
490   TPolyLine3D *cage[16] ={0x0};
491   TPolyLine3D *innerCage[16] ={0x0};
492   TPolyLine3D *innerEndCap1 = new TPolyLine3D();
493   TPolyLine3D *innerEndCap2 = new TPolyLine3D();
494   for(Int_t i = 0; i<16; i++){
495     
496     cage[i] = new TPolyLine3D();
497     cage[i]->SetPoint(0,-fMaxZ0,fOFCRadius*TMath::Cos(i*TMath::TwoPi()/16) ,fOFCRadius*TMath::Sin(i*TMath::TwoPi()/16)) ;
498     cage[i]->SetPoint(1,fMaxZ0,fOFCRadius*TMath::Cos(i*TMath::TwoPi()/16) ,fOFCRadius*TMath::Sin(i*TMath::TwoPi()/16)) ;
499     innerCage[i] = new TPolyLine3D();
500     innerCage[i]->SetPoint(0,-fMaxZ0,fIFCRadius*TMath::Cos(i*TMath::TwoPi()/16) ,fIFCRadius*TMath::Sin(i*TMath::TwoPi()/16)) ;
501     innerCage[i]->SetPoint(1,fMaxZ0,fIFCRadius*TMath::Cos(i*TMath::TwoPi()/16) ,fIFCRadius*TMath::Sin(i*TMath::TwoPi()/16)) ;
502
503     endCap1->SetPoint(i,fMaxZ0,fOFCRadius*TMath::Cos(i*TMath::TwoPi()/16) ,fOFCRadius*TMath::Sin(i*TMath::TwoPi()/16)) ;
504     endCap2->SetPoint(i,-fMaxZ0,fOFCRadius*TMath::Cos(i*TMath::TwoPi()/16) ,fOFCRadius*TMath::Sin(i*TMath::TwoPi()/16)) ;
505
506     innerEndCap1->SetPoint(i,fMaxZ0,fIFCRadius*TMath::Cos(i*TMath::TwoPi()/16) ,fIFCRadius*TMath::Sin(i*TMath::TwoPi()/16)) ;
507     innerEndCap2->SetPoint(i,-fMaxZ0,fIFCRadius*TMath::Cos(i*TMath::TwoPi()/16) ,fIFCRadius*TMath::Sin(i*TMath::TwoPi()/16)) ;
508     innerCage[i]->Draw("same");
509     if(!(i%2))  cage[i]->Draw("same");
510     
511   }
512   endCap1->SetPoint(16,fMaxZ0,fOFCRadius*TMath::Cos(16*TMath::TwoPi()/16) ,fOFCRadius*TMath::Sin(16*TMath::TwoPi()/16)) ;
513   endCap2->SetPoint(16,-fMaxZ0,fOFCRadius*TMath::Cos(16*TMath::TwoPi()/16) ,fOFCRadius*TMath::Sin(16*TMath::TwoPi()/16)) ;
514
515   innerEndCap1->SetPoint(16,fMaxZ0,fIFCRadius*TMath::Cos(16*TMath::TwoPi()/16) ,fIFCRadius*TMath::Sin(16*TMath::TwoPi()/16)) ;
516   innerEndCap2->SetPoint(16,-fMaxZ0,fIFCRadius*TMath::Cos(16*TMath::TwoPi()/16) ,fIFCRadius*TMath::Sin(16*TMath::TwoPi()/16)) ;
517     
518
519   //fDispGraph->SetTitle("ToyMC display");
520   
521   endCap1->Draw("same");
522   endCap2->Draw("same");
523
524   innerEndCap2->Draw("same");
525   innerEndCap1->Draw("same");
526   
527
528
529
530 }
531