source: trunk/MagicSoft/Mars/macros/trigrate.C@ 1788

Last change on this file since 1788 was 1788, checked in by moralejo, 22 years ago
*** empty log message ***
File size: 5.6 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz 12/2000 (tbretz@uni-sw.gwdg.de)
19!
20! Copyright: MAGIC Software Development, 2000-2001
21!
22! Modified 4/7/2002, Abelardo Moralejo:
23! Added one optional input parameter: a camera .root file containing
24! pure NSB events. One such file is generated running the camera over an
25! "empty" reflector file, with the NSB option on, and asking the camera
26! program (see camera program manual) to do the simulation even if no
27! photoelectron from the shower arrives at the camera. One also needs to
28! write to the output file all the events, both triggered and untriggered
29! (see again camera manual). These nsb camera files must contain the same
30! trigger conditions as the proton file.
31!
32! If no nsb file is supplied, the macro will assume no triggers from
33! pure NSB fluctuations.
34!
35\* ======================================================================== */
36
37Float_t GetNSBEvents(TString name, Float_t *BgR, int dim)
38{
39 Int_t numnsbevents;
40
41 TFile bgfile(name);
42 TTree *events = (TTree*) bgfile.Get("Events");
43
44 TH1F h("h","",5,.5,5.5);
45
46 UInt_t from = dim>0 ? 1 : -dim;
47 UInt_t to = dim>0 ? dim : -dim;
48
49 cout << endl << "Calculating NSB triggers from file " << name << "..." << endl << endl;
50
51 for (UInt_t i = from; i <= to; i++)
52 {
53 TString cond("MMcTrig;");
54 cond += i;
55 cond += ".";
56 TBranch *b = events->GetBranch(cond);
57 MMcTrig *mctrig;
58 b->SetAddress(&mctrig);
59
60 Int_t tottrig = 0.;
61
62 UInt_t imax = b->GetEntries();
63
64 for (UInt_t iev = 0; iev < imax; iev++)
65 {
66 b->GetEvent(iev);
67 tottrig += mctrig->GetFirstLevel();
68 }
69 // Set total number of L1 triggers:
70 BgR[dim>0? i-1: 0] = (Float_t) tottrig;
71
72 numnsbevents = (Float_t) imax;
73 }
74
75 return numnsbevents;
76}
77
78
79void trigrate(int dim=0, char *filename = "data/camera.root",
80 char *nsbfile = NULL)
81{
82 // The dim parameter has three possible values:
83 // dim : = 0 -> root file with 1 trigger condition.
84 // > 0 -> number of trigger condition to be analised
85 // in multi conditon file.
86 // < 0 -> selects the -dim trigger condition.
87 //
88 // first we have to create our empty lists
89 //
90 MParList parlist;
91 MTaskList tasklist;
92
93 //
94 // Setup the parameter list.
95 // - we do not need to create any other container. All of them
96 // are created automatically without loss - we don't have to
97 // access them-
98 //
99 // - we need to create MHMcRate only. The other containers
100 // are created automatically without loss - we don't have to
101 // access them-
102 // - MHMcRate must be created by us because we need the pointer
103 // to it and if it would get created automatically it would also be
104 // deleted automatically
105 // - Actually, depending on using a single trigger option MonteCarlo
106 // file or a multyple trigger option, a MHMcRate or an array of
107 // MHMcRate are needed.
108 //
109 parlist.AddToList(&tasklist);
110
111 //
112 // You don't have to add the MHMcRate container here by hand.
113 // But if you want to print or display these containers later on
114 // it is necessary (Rem: No printing or displaying is done in this
115 // macro yet)
116 //
117 UInt_t from = dim>0 ? 1 : -dim;
118 UInt_t to = dim>0 ? dim : -dim;
119
120 Int_t num = to-from+1;
121
122 TObjArray hists(MParList::CreateObjList("MHMcRate", from, to));
123 hists.SetOwner();
124
125 //
126 // Check if the list really contains the right number of histograms
127 //
128 if (hists.GetEntriesFast() != num)
129 return;
130
131 //
132 // Add the histograms to the paramater list.
133 //
134 parlist.AddToList(&hists);
135
136 //
137 // Setup out tasks:
138 // - First we have to read the events
139 // - Then we can calculate rates, for what the number of
140 // triggered showers from a empty reflector file for the
141 // analised trigger conditions should be set (BgR[])
142 //
143
144 MReadMarsFile reader("Events", filename);
145 tasklist.AddToList(&reader);
146
147 // Now we have to build the BgR array, containing the number
148 // of triggers (may be more than 1 trigger/event!) from the
149 // numnsbevents simulated in the nsbfile (3rd input parameter).
150 // If no nsbfile is supplied, we assume no triggers from NSB
151
152 Float_t* BgR = new Float_t[num];
153 memset(BgR, 0, num*sizeof(Float_t));
154
155 Float_t numnsbevents = 5.e4; // some default value.
156 if (nsbfile)
157 numnsbevents = GetNSBEvents(nsbfile, BgR, dim);
158
159 cout << "Number of Trigger conditions: " << num << endl;
160
161 MMcTriggerRateCalc rate(dim, BgR, numnsbevents);
162
163 tasklist.AddToList(&rate);
164
165 //
166 // set up the loop for the processing
167 //
168 MEvtLoop magic;
169 magic.SetParList(&parlist);
170
171 //
172 // Start to loop over all events
173 //
174 MProgressBar bar;
175 magic.SetProgressBar(&bar);
176 if (!magic.Eventloop())
177 return;
178
179 delete BgR;
180
181 hists.Print();
182}
Note: See TracBrowser for help on using the repository browser.