/[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.4 - (show annotations) (download)
Tue Jan 16 15:26:26 2007 UTC (17 years, 10 months ago) by mocchiut
Branch: MAIN
CVS Tags: v5r00, v3r03, v4r00, v10RED, v6r00, v9r00, HEAD
Changes since 1.3: +1 -0 lines
File MIME type: text/plain
Added method UpdateRunInfo to synchronize reading of events and runs (see example3.C)

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
55 Int_t nevent = T->GetEntries();
56 //
57 cout << endl<<" Start loop over events ";
58 for (Int_t i=0; i<nevent;i++){
59 //
60 pam_event->Clear();
61
62 T->GetEntry(i);
63 //================================================================================
64 // some general quantities
65 //================================================================================
66 // tracker
67 ntrack->Fill( pam_event->GetTrkLevel2()->GetNTracks() ); //<<
68 // calorimeter
69 qtot->Fill( pam_event->GetCaloLevel2()->qtot ); //<<
70 // ToF
71 Int_t npa=0;
72 for(Int_t ipa=0; ipa<6; ipa++)npa = npa + pam_event->GetToFLevel2()->GetNHitPaddles(ipa); //<<
73 npaddle->Fill(npa);
74 //================================================================================
75 // track related variables
76 //================================================================================
77 for(Int_t it=0; it<pam_event->GetTrkLevel2()->GetNTracks(); it++){ // << loop over the "physical" tracks (no track images)
78 //
79 // << get the it-th physical pamela track
80 // << PamTrack combines the tracker, calorimeter and ToF track-related variables
81 // << (in case of image, choose the best track by means of calorimeter and ToF data)
82 //
83 PamTrack *track = pam_event->GetTrack(it); //<<
84 //
85 // << if the track fit is good, fill some histos
86 //
87 if( track->GetTrkTrack()->chi2 > 0 && track->GetTrkTrack()->chi2 < 100){
88 //
89 // << spatial residuals on the first calorimeter plane for the sorted track
90 //
91 Float_t rxs = track->GetCaloTrack()->tbar[0][0] - pam_event->GetCaloLevel2()->cbar[0][0];
92 Float_t rys = track->GetCaloTrack()->tbar[0][1] - pam_event->GetCaloLevel2()->cbar[0][1];
93 resxs->Fill(rxs);
94 resys->Fill(rys);
95 //
96 rig->Fill( track->GetTrkTrack()->GetRigidity() );
97 //
98 // << spatial residuals on the first calorimeter plane for the image track
99 //
100 if(track->GetTrkTrack()->HasImage()){ // << if the sorted track has an image...
101 //
102 PamTrack *image = pam_event->GetTrackImage(it); // << ...get the image track
103 //
104 Float_t rxi = image->GetCaloTrack()->tbar[0][0] - pam_event->GetCaloLevel2()->cbar[0][0];
105 Float_t ryi = image->GetCaloTrack()->tbar[0][1] - pam_event->GetCaloLevel2()->cbar[0][1];
106 resxi->Fill(rxi);
107 resyi->Fill(ryi);
108 };
109 };
110 }; // end loop over tracks
111 //================================================================================
112
113 }; // end loop over the events
114 cout << endl << endl << " Done "<< endl<<endl;
115 //
116 // close file
117 //
118 f.Close();
119
120 };

  ViewVC Help
Powered by ViewVC 1.1.23