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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.1.1.1 - (hide annotations) (download) (vendor branch)
Fri Jun 16 16:43:55 2006 UTC (18 years, 5 months ago) by pam-fi
Branch: PamelaLevel2
CVS Tags: v0r01
Changes since 1.1: +0 -0 lines
File MIME type: text/plain
First CVS release of PamelaLevel2 class library ( compatible with DarthVader from v0r02 )

1 pam-fi 1.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 local methods allows to sort tracks and solve the tracker y-view ambiguity (reject track images)
6     //
7     void example1(TString file){
8     //
9     // histos
10     //
11     TH1F* qtot = new TH1F("qtot","Total energy in the calorimeter",1000,0.,2000.);
12     TH1F* npaddle = new TH1F("npaddle","Total number of hit ToF paddles",16,-0.5,15.5);
13     TH1F* ntrack = new TH1F("ntr","Number of fitted tracks",6,-0.5,5.5);
14     //
15     TH1F* rig = new TH1F("rig","Particle rigidity",100,0.,100.);
16     TH1F* resxs = new TH1F("resxs","Spatial residual (X) on the 1^ calo plane (sorted tracks)",100,-10.,10.);
17     TH1F* resys = new TH1F("resys","Spatial residual (Y) on the 1^ calo plane (sorted tracks)",100,-10.,10.);
18     TH1F* resxi = new TH1F("resxi","Spatial residual (X) on the 1^ calo plane (image tracks)",100,-10.,10.);
19     TH1F* resyi = new TH1F("resyi","Spatial residual (Y) on the 1^ calo plane (image tracks)",100,-10.,10.);
20     //
21     PamLevel2* pam_event = new PamLevel2(); // << create pamela event
22     //
23     TFile f(file);
24     //
25     TTree *T = pam_event->LoadPamTrees(&f); // << load Pamela trees from file f
26     Int_t nevent = T->GetEntries();
27     //
28     cout << endl<<" Start loop over events ";
29     for (Int_t i=0; i<nevent;i++){
30     //
31     T->GetEntry(i);
32     //================================================================================
33     // some general quantities
34     //================================================================================
35     // tracker
36     ntrack->Fill( pam_event->GetNTracks());
37     // calorimeter
38     qtot->Fill(pam_event->qtot);
39     // ToF
40     Int_t npa=0;
41     for(Int_t ipa=0; ipa<6; ipa++)npa = npa + pam_event->GetNHitPaddles(ipa);
42     npaddle->Fill(npa);
43     //================================================================================
44     // track related variables
45     //================================================================================
46     for(Int_t it=0; it<pam_event->GetNTracks(); it++){ // << loop over the "physical" tracks (no track images)
47     //
48     // << get the it-th physical pamela track
49     // << PamTrack combines the tracker, calorimeter and ToF track-related variables
50     // << (in case of image, choose the best track by means of calorimeter and ToF data)
51     //
52     PamTrack *track = pam_event->GetTrack(it);
53     //
54     // << if the track fit is good, fill some histos
55     //
56     if( track->chi2 > 0 && track->chi2 < 100){
57     //
58     // << spatial residuals on the first calorimeter plane for the sorted track
59     //
60     Float_t rxs = track->tbar[0][0] - pam_event->cbar[0][0];
61     Float_t rys = track->tbar[0][1] - pam_event->cbar[0][1];
62     resxs->Fill(rxs);
63     resys->Fill(rys);
64     //
65     rig->Fill( track->GetRigidity() );
66     //
67     // << spatial residuals on the first calorimeter plane for the image track
68     //
69     if(track->HasImage()){ // << if the sorted track has an image...
70     //
71     PamTrack *image = pam_event->GetTrackImage(it); // << ...get the image track
72     //
73     Float_t rxi = image->tbar[0][0] - pam_event->cbar[0][0];
74     Float_t ryi = image->tbar[0][1] - pam_event->cbar[0][1];
75     resxi->Fill(rxi);
76     resyi->Fill(ryi);
77     };
78     };
79     }; // end loop over tracks
80     //================================================================================
81    
82     }; // end loop over the events
83     cout << endl << endl << " Done "<< endl<<endl;
84     //
85     // close file
86     //
87     f.Close();
88    
89     };

  ViewVC Help
Powered by ViewVC 1.1.23