| 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 | #include <TRandom.h>
|
|---|
| 52 |
|
|---|
| 53 | #include "MParList.h"
|
|---|
| 54 |
|
|---|
| 55 | #include "MLog.h"
|
|---|
| 56 | #include "MLogManip.h"
|
|---|
| 57 |
|
|---|
| 58 | #include "MMcEvt.hxx"
|
|---|
| 59 | #include "MMcCorsikaRunHeader.h"
|
|---|
| 60 |
|
|---|
| 61 | ClassImp(MFEnergySlope);
|
|---|
| 62 |
|
|---|
| 63 | using namespace std;
|
|---|
| 64 |
|
|---|
| 65 | // --------------------------------------------------------------------------
|
|---|
| 66 | //
|
|---|
| 67 | // Default Constructor
|
|---|
| 68 | //
|
|---|
| 69 | MFEnergySlope::MFEnergySlope(const char *name, const char *title):
|
|---|
| 70 | fNewSlope(-1), fMcMinEnergy(-1.), fMcMaxEnergy(-1.)
|
|---|
| 71 | {
|
|---|
| 72 | fName = name ? name : "MFEnergySlope";
|
|---|
| 73 | fTitle = title ? title : "Filter to select energy with given slope";
|
|---|
| 74 | }
|
|---|
| 75 |
|
|---|
| 76 | // --------------------------------------------------------------------------
|
|---|
| 77 | //
|
|---|
| 78 | // Constructor
|
|---|
| 79 | //
|
|---|
| 80 | MFEnergySlope::MFEnergySlope(Float_t slope, const char *name, const char *title):
|
|---|
| 81 | fNewSlope(TMath::Abs(slope)), fMcMinEnergy(-1.), fMcMaxEnergy(-1.)
|
|---|
| 82 | {
|
|---|
| 83 | fName = name ? name : "MFEnergySlope";
|
|---|
| 84 | fTitle = title ? title : "Filter to select energy with given slope";
|
|---|
| 85 | }
|
|---|
| 86 |
|
|---|
| 87 | // --------------------------------------------------------------------------
|
|---|
| 88 | //
|
|---|
| 89 | // Constructor
|
|---|
| 90 | //
|
|---|
| 91 | MFEnergySlope::MFEnergySlope(Float_t slope, Float_t emin, const char *name, const char *title):
|
|---|
| 92 | fNewSlope(TMath::Abs(slope)), fMcMinEnergy(emin), fMcMaxEnergy(-1.)
|
|---|
| 93 | {
|
|---|
| 94 | fName = name ? name : "MFEnergySlope";
|
|---|
| 95 | fTitle = title ? title : "Filter to select energy with given slope";
|
|---|
| 96 | }
|
|---|
| 97 |
|
|---|
| 98 | // --------------------------------------------------------------------------
|
|---|
| 99 | //
|
|---|
| 100 | // Preprocess
|
|---|
| 101 | //
|
|---|
| 102 | // MC slope and min/max energy are read
|
|---|
| 103 | // Normalization factor is computed
|
|---|
| 104 | //
|
|---|
| 105 | Int_t MFEnergySlope::PreProcess(MParList *pList)
|
|---|
| 106 | {
|
|---|
| 107 | if (fNewSlope<0)
|
|---|
| 108 | {
|
|---|
| 109 | *fLog << err << "New slope still undefined... aborting." << endl;
|
|---|
| 110 | return kFALSE;
|
|---|
| 111 | }
|
|---|
| 112 |
|
|---|
| 113 | fEvt = (MMcEvt*)pList->FindObject("MMcEvt");
|
|---|
| 114 | if (!fEvt)
|
|---|
| 115 | {
|
|---|
| 116 | *fLog << err << "MMcEvt not found... aborting." << endl;
|
|---|
| 117 | return kFALSE;
|
|---|
| 118 | }
|
|---|
| 119 |
|
|---|
| 120 | return kTRUE;
|
|---|
| 121 | }
|
|---|
| 122 |
|
|---|
| 123 | // --------------------------------------------------------------------------
|
|---|
| 124 | //
|
|---|
| 125 | // Preprocess
|
|---|
| 126 | //
|
|---|
| 127 | // MC slope and min/max energy are read
|
|---|
| 128 | // Normalization factor is computed
|
|---|
| 129 | //
|
|---|
| 130 | Bool_t MFEnergySlope::ReInit(MParList *pList)
|
|---|
| 131 | {
|
|---|
| 132 | MMcCorsikaRunHeader *runheader = (MMcCorsikaRunHeader*)pList->FindObject("MMcCorsikaRunHeader");
|
|---|
| 133 | if (!runheader)
|
|---|
| 134 | {
|
|---|
| 135 | *fLog << err << "MMcCorsikaRunHeader not found... aborting." << endl;
|
|---|
| 136 | return kFALSE;
|
|---|
| 137 | }
|
|---|
| 138 |
|
|---|
| 139 | //
|
|---|
| 140 | // Read info from0 the MC sample (it must be generated with
|
|---|
| 141 | // reflector ver.0.6 and camera ver. 0.6)
|
|---|
| 142 | //
|
|---|
| 143 | fMcSlope = runheader->GetSlopeSpec();
|
|---|
| 144 | if (fMcMinEnergy<0)
|
|---|
| 145 | fMcMinEnergy = runheader->GetELowLim();
|
|---|
| 146 | if (fMcMinEnergy<0)
|
|---|
| 147 | fMcMaxEnergy = runheader->GetEUppLim();
|
|---|
| 148 |
|
|---|
| 149 | // Slope is used with positive values in the code
|
|---|
| 150 | if (fNewSlope < 0)
|
|---|
| 151 | fNewSlope *= -1;
|
|---|
| 152 | if (fMcSlope < 0)
|
|---|
| 153 | fMcSlope *= -1;
|
|---|
| 154 |
|
|---|
| 155 | // Calculate normalization factor
|
|---|
| 156 | fN0 = TMath::Power(fNewSlope>fMcSlope ? fMcMinEnergy : fMcMaxEnergy, fNewSlope-fMcSlope);
|
|---|
| 157 |
|
|---|
| 158 | // Print some infos
|
|---|
| 159 | *fLog << inf;
|
|---|
| 160 | *fLog << "Fetched MC info:" << endl;
|
|---|
| 161 | *fLog << " Change E Slope from " << -fMcSlope << " (" << fMcMinEnergy << " < E < ";
|
|---|
| 162 | *fLog << fMcMaxEnergy << ") to " << -fNewSlope << endl;
|
|---|
| 163 | *fLog << " Norm factor: " << fN0 << endl;
|
|---|
| 164 |
|
|---|
| 165 | return kTRUE;
|
|---|
| 166 | }
|
|---|
| 167 |
|
|---|
| 168 | // --------------------------------------------------------------------------
|
|---|
| 169 | //
|
|---|
| 170 | // Select events randomly according to the MC ("old") and required ("new")
|
|---|
| 171 | // energy slope.
|
|---|
| 172 | // Old E slope is fMcSlope
|
|---|
| 173 | // New E slope is set by the user (fval; fNewSlope)
|
|---|
| 174 | // If old and new energy slope are the same skip the selection.
|
|---|
| 175 | // The MC energy slope and lower and upper limits are taken from the
|
|---|
| 176 | // run header (requires reflector ver.>0.6 and camera ver.>0.6)
|
|---|
| 177 | //
|
|---|
| 178 | Int_t MFEnergySlope::Process()
|
|---|
| 179 | {
|
|---|
| 180 | fResult = kTRUE;
|
|---|
| 181 |
|
|---|
| 182 | // Energy slopes are the same: skip it
|
|---|
| 183 | if (fNewSlope == fMcSlope)
|
|---|
| 184 | return kTRUE;
|
|---|
| 185 |
|
|---|
| 186 | // The value of the normalized spectrum is compared with a
|
|---|
| 187 | // random value in [0,1];
|
|---|
| 188 | // if the latter is higher the event is accepted
|
|---|
| 189 | const Float_t energy = fEvt->GetEnergy();
|
|---|
| 190 |
|
|---|
| 191 | /*
|
|---|
| 192 | //
|
|---|
| 193 | // If energy bounds different from MC ones are selected, then
|
|---|
| 194 | // events outside these bounds are rejected, as their slope has
|
|---|
| 195 | // not been changed.
|
|---|
| 196 | //
|
|---|
| 197 | if (energy > fMcMaxEnergy || energy < fMcMinEnergy)
|
|---|
| 198 | {
|
|---|
| 199 | fResult = kFALSE;
|
|---|
| 200 | return kTRUE;
|
|---|
| 201 | }
|
|---|
| 202 | */
|
|---|
| 203 |
|
|---|
| 204 | const Float_t Nexp = fN0 * pow(energy, fMcSlope-fNewSlope);
|
|---|
| 205 | const Float_t Nrnd = gRandom->Uniform();
|
|---|
| 206 |
|
|---|
| 207 | fResult = Nexp > Nrnd;
|
|---|
| 208 |
|
|---|
| 209 | return kTRUE;
|
|---|
| 210 | }
|
|---|
| 211 |
|
|---|
| 212 |
|
|---|
| 213 | // --------------------------------------------------------------------------
|
|---|
| 214 | //
|
|---|
| 215 | // MFEnergySlope.NewSlope: -2.8
|
|---|
| 216 | //
|
|---|
| 217 | Int_t MFEnergySlope::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
|
|---|
| 218 | {
|
|---|
| 219 | Bool_t rc = kFALSE;
|
|---|
| 220 | if (IsEnvDefined(env, prefix, "NewSlope", print))
|
|---|
| 221 | {
|
|---|
| 222 | rc = kTRUE;
|
|---|
| 223 | SetNewSlope(GetEnvValue(env, prefix, "NewSlope", fNewSlope));
|
|---|
| 224 | }
|
|---|
| 225 | return rc;
|
|---|
| 226 | }
|
|---|