source: trunk/MagicSoft/Mars/mfilter/MFEnergySlope.cc@ 9620

Last change on this file since 9620 was 8907, checked in by tbretz, 17 years ago
*** empty log message ***
File size: 7.5 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): Antonio Stamerra 02/2003 <mailto:antonio.stamerra@pi.infn.it>
19!
20! Copyright: MAGIC Software Development, 2000-2005
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26// //
27// MFEnergySlope //
28// //
29// A filter to select MC events (generated with a energy slope MCSlope) //
30// with a different energy slope NewSlope set by the user. //
31// //
32// The new slope is set by the user with SetSlope(). //
33// Only negative slopes are admitted; positive ones are internally //
34// converted. //
35// Events are selected following the new energy power slope, and the //
36// sample is normalized to the number of events at: //
37// 1. McMaxEnergy, if abs(NewSlope) < abs(McSlope); //
38// 2. McMinEnergy, if abs(NewSlope) > abs(McSlope); //
39// Mc{Max,Min}Energy are set with SetMcMinEnergy() and SetMcMaxEnergy(); //
40// with GeV values. //
41// Default values are the min/max energy of the MC sample. //
42// //
43// With this filter the statistics of the MC sample is reduced. //
44// camera ver.0.6 and reflector ver.0.6 are required to fetch //
45// the correct MC information //
46// //
47/////////////////////////////////////////////////////////////////////////////
48#include "MFEnergySlope.h"
49
50#include <fstream>
51
52#include <TMath.h>
53#include <TRandom.h>
54
55#include "MParList.h"
56
57#include "MLog.h"
58#include "MLogManip.h"
59
60#include "MMcEvt.hxx"
61#include "MMcCorsikaRunHeader.h"
62
63ClassImp(MFEnergySlope);
64
65using namespace std;
66
67// --------------------------------------------------------------------------
68//
69// Default Constructor
70//
71MFEnergySlope::MFEnergySlope(const char *name, const char *title):
72 fNewSlope(-1), fMcMinEnergy(-1.), fMcMaxEnergy(-1.)
73{
74 fName = name ? name : "MFEnergySlope";
75 fTitle = title ? title : "Filter to select energy with given slope";
76}
77
78// --------------------------------------------------------------------------
79//
80// Constructor
81//
82MFEnergySlope::MFEnergySlope(Float_t slope, const char *name, const char *title):
83 fNewSlope(TMath::Abs(slope)), fMcMinEnergy(-1.), fMcMaxEnergy(-1.)
84{
85 fName = name ? name : "MFEnergySlope";
86 fTitle = title ? title : "Filter to select energy with given slope";
87}
88
89// --------------------------------------------------------------------------
90//
91// Constructor
92//
93MFEnergySlope::MFEnergySlope(Float_t slope, Float_t emin, const char *name, const char *title):
94 fNewSlope(TMath::Abs(slope)), fMcMinEnergy(emin), fMcMaxEnergy(-1.)
95{
96 fName = name ? name : "MFEnergySlope";
97 fTitle = title ? title : "Filter to select energy with given slope";
98}
99
100// --------------------------------------------------------------------------
101//
102// Preprocess
103//
104// MC slope and min/max energy are read
105// Normalization factor is computed
106//
107Int_t MFEnergySlope::PreProcess(MParList *pList)
108{
109 if (fNewSlope<0)
110 {
111 *fLog << err << "New slope still undefined... aborting." << endl;
112 return kFALSE;
113 }
114
115 fEvt = (MMcEvt*)pList->FindObject("MMcEvt");
116 if (!fEvt)
117 {
118 *fLog << err << "MMcEvt not found... aborting." << endl;
119 return kFALSE;
120 }
121
122 return kTRUE;
123}
124
125// --------------------------------------------------------------------------
126//
127// Preprocess
128//
129// MC slope and min/max energy are read
130// Normalization factor is computed
131//
132Bool_t MFEnergySlope::ReInit(MParList *pList)
133{
134 MMcCorsikaRunHeader *runheader = (MMcCorsikaRunHeader*)pList->FindObject("MMcCorsikaRunHeader");
135 if (!runheader)
136 {
137 *fLog << err << "MMcCorsikaRunHeader not found... aborting." << endl;
138 return kFALSE;
139 }
140
141 //
142 // Read info from0 the MC sample (it must be generated with
143 // reflector ver.0.6 and camera ver. 0.6)
144 //
145 fMcSlope = runheader->GetSlopeSpec();
146 if (fMcMinEnergy<0)
147 fMcMinEnergy = runheader->GetELowLim();
148 if (fMcMinEnergy<0)
149 fMcMaxEnergy = runheader->GetEUppLim();
150
151 // Slope is used with positive values in the code
152 if (fNewSlope < 0)
153 fNewSlope *= -1;
154 if (fMcSlope < 0)
155 fMcSlope *= -1;
156
157 // Calculate normalization factor
158 fN0 = TMath::Power(fNewSlope>fMcSlope ? fMcMinEnergy : fMcMaxEnergy, fNewSlope-fMcSlope);
159
160 // Print some infos
161 *fLog << inf;
162 *fLog << "Fetched MC info:" << endl;
163 *fLog << " Change E Slope from " << -fMcSlope << " (" << fMcMinEnergy << " < E < ";
164 *fLog << fMcMaxEnergy << ") to " << -fNewSlope << endl;
165 *fLog << " Norm factor: " << fN0 << endl;
166
167 return kTRUE;
168}
169
170// --------------------------------------------------------------------------
171//
172// Select events randomly according to the MC ("old") and required ("new")
173// energy slope.
174// Old E slope is fMcSlope
175// New E slope is set by the user (fval; fNewSlope)
176// If old and new energy slope are the same skip the selection.
177// The MC energy slope and lower and upper limits are taken from the
178// run header (requires reflector ver.>0.6 and camera ver.>0.6)
179//
180Int_t MFEnergySlope::Process()
181{
182 fResult = kTRUE;
183
184 // Energy slopes are the same: skip it
185 if (fNewSlope == fMcSlope)
186 return kTRUE;
187
188 // The value of the normalized spectrum is compared with a
189 // random value in [0,1];
190 // if the latter is higher the event is accepted
191 const Float_t energy = fEvt->GetEnergy();
192
193 /*
194 //
195 // If energy bounds different from MC ones are selected, then
196 // events outside these bounds are rejected, as their slope has
197 // not been changed.
198 //
199 if (energy > fMcMaxEnergy || energy < fMcMinEnergy)
200 {
201 fResult = kFALSE;
202 return kTRUE;
203 }
204 */
205
206 const Float_t Nexp = fN0 * pow(energy, fMcSlope-fNewSlope);
207 const Float_t Nrnd = gRandom->Uniform();
208
209 fResult = Nexp > Nrnd;
210
211 return kTRUE;
212}
213
214
215// --------------------------------------------------------------------------
216//
217// MFEnergySlope.NewSlope: -2.8
218//
219Int_t MFEnergySlope::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
220{
221 Bool_t rc = kFALSE;
222 if (IsEnvDefined(env, prefix, "NewSlope", print))
223 {
224 rc = kTRUE;
225 SetNewSlope(GetEnvValue(env, prefix, "NewSlope", fNewSlope));
226 }
227 return rc;
228}
Note: See TracBrowser for help on using the repository browser.