/[PAMELA software]/PamelaLevel2/doc/examples/example1.C
ViewVC logotype

Contents of /PamelaLevel2/doc/examples/example1.C

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.3 - (show annotations) (download)
Mon Jan 15 11:51:39 2007 UTC (17 years, 10 months ago) by pam-fi
Branch: MAIN
Changes since 1.2: +81 -51 lines
File MIME type: text/plain
v3r00 **NEW**

1 //
2 // Example to integrate all detector info by means of the PamLevel2 class
3 //
4 // An object of this class has access to all the public members of all the detector classes.
5 // Some methods allows to sort tracks and solve the tracker y-view ambiguity (reject track images)
6 //
7 // --- Modified on January 2007 ---
8 //
9 // In order to access the members of the detector classes you have now to call the
10 // methods
11 //
12 // TrkLevel2* PamLevel2::GetTrkLevel2()
13 // CaloLevel2* PamLevel2::GetCaloLevel2()
14 // OrbitalInfo* PamLevel2::GetOrbitalInfo()
15 // ecc...
16 //
17 // For example:
18 //
19 // PamLevel2* pam_event = new PamLevel2();
20 // cout << pam_event->GetCaloLevel2()->qtot << endl;
21 //
22 // The same is for the track. The class PamTrack combines the track-related info from
23 // tracker, calorimeter and TOF. To access the members of these detector track-related members
24 // you have to call the methods:
25 //
26 // TrkTrack* PamTrack::GetTrkTrack()
27 // CaloTrkVar* PamTrack::GetCaloTrack()
28 // ToFTrkVar* PamTrack::GetToFTrack()
29 //
30 // For example:
31 //
32 // PamTrack *track = pam_event->GetTrack(0); //<< retrieve first track, solving the track-image ambiguity
33 // if(track)cout<< track->GetTrkTrack()->GetRigidity()<<endl;
34 //
35 void example1(TString file){
36 //
37 // histos
38 //
39 TH1F* qtot = new TH1F("qtot","Total energy in the calorimeter",1000,0.,2000.);
40 TH1F* npaddle = new TH1F("npaddle","Total number of hit ToF paddles",16,-0.5,15.5);
41 TH1F* ntrack = new TH1F("ntr","Number of fitted tracks",6,-0.5,5.5);
42 //
43 TH1F* rig = new TH1F("rig","Particle rigidity",100,0.,100.);
44 TH1F* resxs = new TH1F("resxs","Spatial residual (X) on the 1^ calo plane (sorted tracks)",100,-10.,10.);
45 TH1F* resys = new TH1F("resys","Spatial residual (Y) on the 1^ calo plane (sorted tracks)",100,-10.,10.);
46 TH1F* resxi = new TH1F("resxi","Spatial residual (X) on the 1^ calo plane (image tracks)",100,-10.,10.);
47 TH1F* resyi = new TH1F("resyi","Spatial residual (Y) on the 1^ calo plane (image tracks)",100,-10.,10.);
48 //
49 PamLevel2* pam_event = new PamLevel2(); // << create pamela event
50 //
51 TFile f(file);
52 //
53 TTree *T = pam_event->GetPamTree(&f); // << get Pamela trees from file f
54 Int_t nevent = T->GetEntries();
55 //
56 cout << endl<<" Start loop over events ";
57 for (Int_t i=0; i<nevent;i++){
58 //
59 pam_event->Clear();
60
61 T->GetEntry(i);
62 //================================================================================
63 // some general quantities
64 //================================================================================
65 // tracker
66 ntrack->Fill( pam_event->GetTrkLevel2()->GetNTracks() ); //<<
67 // calorimeter
68 qtot->Fill( pam_event->GetCaloLevel2()->qtot ); //<<
69 // ToF
70 Int_t npa=0;
71 for(Int_t ipa=0; ipa<6; ipa++)npa = npa + pam_event->GetToFLevel2()->GetNHitPaddles(ipa); //<<
72 npaddle->Fill(npa);
73 //================================================================================
74 // track related variables
75 //================================================================================
76 for(Int_t it=0; it<pam_event->GetTrkLevel2()->GetNTracks(); it++){ // << loop over the "physical" tracks (no track images)
77 //
78 // << get the it-th physical pamela track
79 // << PamTrack combines the tracker, calorimeter and ToF track-related variables
80 // << (in case of image, choose the best track by means of calorimeter and ToF data)
81 //
82 PamTrack *track = pam_event->GetTrack(it); //<<
83 //
84 // << if the track fit is good, fill some histos
85 //
86 if( track->GetTrkTrack()->chi2 > 0 && track->GetTrkTrack()->chi2 < 100){
87 //
88 // << spatial residuals on the first calorimeter plane for the sorted track
89 //
90 Float_t rxs = track->GetCaloTrack()->tbar[0][0] - pam_event->GetCaloLevel2()->cbar[0][0];
91 Float_t rys = track->GetCaloTrack()->tbar[0][1] - pam_event->GetCaloLevel2()->cbar[0][1];
92 resxs->Fill(rxs);
93 resys->Fill(rys);
94 //
95 rig->Fill( track->GetTrkTrack()->GetRigidity() );
96 //
97 // << spatial residuals on the first calorimeter plane for the image track
98 //
99 if(track->GetTrkTrack()->HasImage()){ // << if the sorted track has an image...
100 //
101 PamTrack *image = pam_event->GetTrackImage(it); // << ...get the image track
102 //
103 Float_t rxi = image->GetCaloTrack()->tbar[0][0] - pam_event->GetCaloLevel2()->cbar[0][0];
104 Float_t ryi = image->GetCaloTrack()->tbar[0][1] - pam_event->GetCaloLevel2()->cbar[0][1];
105 resxi->Fill(rxi);
106 resyi->Fill(ryi);
107 };
108 };
109 }; // end loop over tracks
110 //================================================================================
111
112 }; // end loop over the events
113 cout << endl << endl << " Done "<< endl<<endl;
114 //
115 // close file
116 //
117 f.Close();
118
119 };

  ViewVC Help
Powered by ViewVC 1.1.23