source: trunk/MagicSoft/Mars/mhist/MHMcRate.cc@ 1821

Last change on this file since 1821 was 1800, checked in by moralejo, 22 years ago
*** empty log message ***
File size: 8.1 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 <mailto:tbretz@uni-sw.gwdg.de>
19! Author(s): Harald Kornmayer 1/2001
20! Author(s): Abelardo Moralejo 2/2003
21!
22! Explanations on the rate calculation can be found in
23! chapter 7 of the following diploma thesis:
24! http://www.pd.infn.it/magic/tesi2.ps.gz (in Italian)
25!
26! Copyright: MAGIC Software Development, 2000-2001
27!
28!
29\* ======================================================================== */
30
31#include "MHMcRate.h"
32
33#include "MLog.h"
34#include "MLogManip.h"
35
36ClassImp(MHMcRate);
37
38void MHMcRate::Init(const char *name, const char *title)
39{
40 fName = name ? name : "MMcTriggerRate";
41 fTitle = title ? title : "Task to calc the collection area ";
42
43 fPartId=0; // Type of particle
44
45 fEnergyMax=0.0; // Maximum Energy (TeV)
46 fEnergyMin=1000000.0; // Minimum Energy (TeV)
47
48 fSolidAngle = -1.; // Solid angle within which incident directions
49 // are distributed
50 fThetaMax=0.0; // Maximum theta angle of run
51 fThetaMin=370.0; // Minimum theta angle of run
52 fPhiMax=0.0; // Maximum phi angle of run
53 fPhiMin=370.0; // Minimum phi angle of run
54
55 fImpactMax=0.0; // Maximum impact parameter
56 fImpactMin=100000.0; // Minimum impact parameter
57
58 fBackTrig=-1.0; // Number of triggers from background
59 fBackSim=-1.0; // Number of simulated showers for the background
60
61 fTriggerRate= -1.0; // Trigger rate in Hz
62 fTriggerRateError= -1.0; // Estimated error for the trigger rate in Hz
63
64 fMultiplicity = -1; // Multiplicity of the trigger condition
65 fMeanThreshold = -1.; // Discr. threshold of the trigger condition
66 fTriggerCondNum = 0; // Trigger condition number within camera file
67
68}
69
70// --------------------------------------------------------------------------
71//
72// default constructor
73// fills all member data with initial values
74//
75MHMcRate::MHMcRate(const char *name, const char *title)
76{
77 Init(name, title);
78
79 fSpecIndex=0.0; // dn/dE = k * e^{- fSpecIndex}
80 fFlux0=-1.0; // dn/dE = fFlux0 * E^{-a}
81
82 fShowerRate= -1.0; // Showers rate in Hz
83 fShowerRateError=0.0; // Estimated error of shower rate in Hz
84}
85
86// --------------------------------------------------------------------------
87//
88// overloaded constructor I
89// fills all member data with initial values and sets the rate of
90// incident showers to ShowRate
91//
92MHMcRate::MHMcRate(Float_t showrate,
93 const char *name, const char *title)
94{
95 Init(name, title);
96
97 fSpecIndex=0.0; // dn/dE = k * e^{- fSpecIndex}
98 fFlux0=-1.0; // dn/dE = fFlux0 * E^{-a}
99
100 fShowerRate= showrate; // Showers rate in Hz
101 fShowerRateError=sqrt(showrate); // Estimated error of shower rate in Hz
102}
103
104// --------------------------------------------------------------------------
105//
106// overloaded constructor I
107// fills all member data with initial values and sets the
108// spectral index and the initial flux to SpecIndex and Flux0
109//
110MHMcRate::MHMcRate(Float_t specindex, Float_t flux0,
111 const char *name, const char *title)
112{
113 Init(name, title);
114
115 fSpecIndex=specindex; // dn/dE = k * e^{- fSpecIndex}
116 fFlux0=flux0; // dn/dE = fFlux0 * E^{-a}
117
118 fShowerRate= -1.0;
119 fShowerRateError=0.0;
120}
121
122// --------------------------------------------------------------------------
123//
124// set the particle that produces the showers in the athmosphere
125//
126void MHMcRate:: SetParticle(UShort_t part)
127{
128 fPartId=part;
129}
130
131// --------------------------------------------------------------------------
132//
133// Set the information about trigger due only to the Night Sky Background:
134//
135void MHMcRate::SetBackground (Float_t showers, Float_t triggers)
136{
137 fBackTrig=showers; // Number of triggers from background
138 fBackSim=triggers; // Number of simulated showers for the background
139}
140
141// --------------------------------------------------------------------------
142//
143// set the parameters to compute the incident rate
144//
145void MHMcRate:: SetFlux(Float_t flux0, Float_t specindx)
146{
147 fFlux0=flux0;
148 fSpecIndex=specindx;
149
150}
151
152// --------------------------------------------------------------------------
153//
154// set the incident rate
155//
156void MHMcRate:: SetIncidentRate(Float_t showerrate)
157{
158 fShowerRate=showerrate;
159}
160
161// --------------------------------------------------------------------------
162//
163// update the limits for energy, theta, phi and impact parameter
164//
165void MHMcRate::UpdateBoundaries(Float_t energy, Float_t theta,
166 Float_t phi, Float_t impact)
167{
168 // It updates the limit values
169
170 if (fThetaMax<theta) fThetaMax=theta;
171 if (fThetaMin>theta) fThetaMin=theta;
172
173 if (fPhiMax<phi) fPhiMax=phi;
174 if (fPhiMin>phi) fPhiMin=phi;
175
176 if (fImpactMax<impact) fImpactMax=impact;
177 if (fImpactMin>impact) fImpactMin=impact;
178
179 if (fEnergyMax<energy) fEnergyMax=energy;
180 if (fEnergyMin>energy) fEnergyMin=energy;
181}
182
183// --------------------------------------------------------------------------
184//
185// compute the trigger rate and set the ReadyToSave bit
186//
187void MHMcRate::CalcRate(Float_t trig, Float_t anal, Float_t simu)
188{
189 // It computes the trigger rate
190
191 // First one computes the rate of incident showers.
192 const Double_t specidx = 1.0-fSpecIndex;
193
194 const Double_t epowmax = pow(fEnergyMax, specidx);
195 const Double_t epowmin = pow(fEnergyMin, specidx);
196
197 if (fShowerRate <= 0)
198 fShowerRate = fFlux0/specidx*(epowmax-epowmin);
199
200 if (fSolidAngle < 0.)
201 fSolidAngle = (fPhiMax-fPhiMin)*(cos(fThetaMin)-cos(fThetaMax));
202
203 if (fPartId!=1)
204 fShowerRate *= fSolidAngle;
205
206 fShowerRate *= TMath::Pi()*(fImpactMax/100.0*fImpactMax/100.0 -
207 fImpactMin/100.0*fImpactMin/100.0);
208
209 fShowerRateError = sqrt(fShowerRate);
210
211 // The simulated trigger time in the camera program is 160 ns:
212 // 9/10/2002, AM: Fixed error below in calculation of "anal2"
213 // ( added factor fShowerRate/simu )
214
215 const Double_t anal2 = 1.0-fShowerRate*(anal/simu)*160.0e-9;
216 const Double_t back2 = fBackSim*160.0e-9;
217
218 // Then the trigger rate and its error is evaluated
219 if(fBackTrig<0){
220 fTriggerRateError = sqrt((trig*fShowerRate*fShowerRate/(simu*simu)) +
221 (anal2*anal2*1/(fBackSim*back2*back2)));
222 fBackTrig=0;
223 }
224 else
225 fTriggerRateError = sqrt((trig*fShowerRate*fShowerRate/(simu*simu)) +
226 (anal2*anal2*fBackTrig/(back2*back2)));
227
228 fTriggerRate = trig*fShowerRate/simu + anal2*fBackTrig/back2;
229
230 SetReadyToSave();
231}
232
233// --------------------------------------------------------------------------
234//
235// print the trigger rate
236//
237void MHMcRate::Print(Option_t *) const
238{
239 *fLog << all << "Incident rate " << fShowerRate << " Hz " << endl;
240 *fLog << "Multiplicity: " << fMultiplicity << ", Discr. threshold: " << fMeanThreshold << endl;
241 *fLog << "Trigger Rate " << fTriggerRate << " +- " << fTriggerRateError << " Hz" << endl;
242}
243
244// --------------------------------------------------------------------------
245//
246// draw the trigger rate
247//
248void MHMcRate::Draw(Option_t *)
249{
250 *fLog << all << dbginf << " - MHMcRate::Draw: To be iplemented" << endl;
251}
252
253TObject *MHMcRate::DrawClone(Option_t *) const
254{
255 *fLog << all << dbginf << " - MHMcRate::DrawClone: To be iplemented" << endl;
256 return NULL;
257}
Note: See TracBrowser for help on using the repository browser.