source: trunk/MagicSoft/Mars/manalysis/MApplyPadding.cc@ 1951

Last change on this file since 1951 was 1951, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 9.3 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): Robert Wagner <magicsoft@rwagner.de> 10/2002
19!
20! Copyright: MAGIC Software Development, 2000-2002
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26// //
27// MApplyPadding //
28// //
29// This task applies padding to a given Sigmabar target value. //
30// The task checks whether the data stream it is applied to has to be //
31// padded or not and if so, Gaussian noise with the difference in Sigma //
32// is produced and added to the particular event. Number of photons and //
33// error on this number is altered. //
34// //
35// There are three ways to define the sigmabar value to which the events //
36// are padded: //
37// //
38// 1) Set a fixed level with SetTargetLevel //
39// //
40// 2) Give a TH1D* which defines the Sigmabar as function of Theta //
41// with SetDefiningHistogram //
42// The given histogram should have the same binning in Theta as //
43// is used in the analysis //
44// //
45// 3) Do nothing, then PreProcess will try to read in a (workaround!) //
46// propriety format ASCII database for the CT1 test. //
47// the name of this file is set by SetDatabaseFile //
48// Better know what you are doing or use methods 1 or 2. //
49// //
50// This implementation is still PRELIMINARY and requires some workarounds //
51// put in SPECIFICALLY FOR THE CT1 TESTS, since a database to access is //
52// missing. It is not the FINAL MAGIC VERSION. //
53// //
54/////////////////////////////////////////////////////////////////////////////
55#include "MApplyPadding.h"
56
57#include <stdio.h>
58
59#include "TH1.h"
60#include "TH2.h"
61#include "TRandom.h"
62
63#include "MBinning.h"
64#include "MSigmabar.h"
65#include "MMcEvt.hxx"
66#include "MLog.h"
67#include "MLogManip.h"
68#include "MParList.h"
69#include "MGeomCam.h"
70#include "MCerPhotPix.h"
71#include "MCerPhotEvt.h"
72#include "MPedestalPix.h"
73
74ClassImp(MApplyPadding);
75
76// --------------------------------------------------------------------------
77//
78// Default constructor.
79//
80MApplyPadding::MApplyPadding(const char *name, const char *title) : fRunType(0), fGroup(0), fUseHistogram(kTRUE), fFixedSigmabar(0.0), fRnd(0)
81{
82 fName = name ? name : "MApplyPadding";
83 fTitle = title ? title : "Task to apply padding";
84 Print();
85}
86
87// --------------------------------------------------------------------------
88//
89// You can provide a TH1D* histogram containing the target Sigmabar in
90// bins of theta. Be sure to use the same binning as for the analysis
91//
92void MApplyPadding::SetDefiningHistogram(TH1D *histo)
93{
94 fHSigmabarMax = histo;
95}
96
97
98// --------------------------------------------------------------------------
99//
100// check if MEvtHeader exists in the Parameter list already.
101// if not create one and add them to the list
102//
103Bool_t MApplyPadding::PreProcess(MParList *pList)
104{
105 fMcEvt = (MMcEvt*)pList->FindObject("MMcEvt");
106 if (!fMcEvt)
107 {
108 *fLog << err << dbginf << "MMcEvt not found... aborting." << endl;
109 return kFALSE;
110 }
111
112 fPed = (MPedestalCam*)pList->FindObject("MPedestalCam");
113 if (!fPed)
114 {
115 *fLog << dbginf << "MPedestalCam not found... aborting." << endl;
116 return kFALSE;
117 }
118
119 fCam = (MGeomCam*)pList->FindObject("MGeomCam");
120 if (!fCam)
121 {
122 *fLog << dbginf << "MGeomCam not found (no geometry information available)... aborting." << endl;
123 return kFALSE;
124 }
125
126 fEvt = (MCerPhotEvt*)pList->FindObject("MCerPhotEvt");
127 if (!fEvt)
128 {
129 *fLog << dbginf << "MCerPhotEvt not found... aborting." << endl;
130 return kFALSE;
131 }
132
133 fSigmabar = (MSigmabar*)pList->FindCreateObj("MSigmabar");
134 if (!fSigmabar)
135 {
136 *fLog << dbginf << "MSigmabar not found... aborting." << endl;
137 return kFALSE;
138 }
139
140 // Get Theta Binning
141 MBinning* binstheta = (MBinning*)pList->FindObject("BinningTheta");
142 if (!binstheta)
143 {
144 *fLog << err << dbginf << "BinningTheta not found... aborting." << endl;
145 return kFALSE;
146 }
147
148 // Create fSigmabarMax histogram
149 // (only if no fixed Sigmabar target value or a histogram have already been
150 // provided)
151 if ((!fUseHistogram) && (fHSigmabarMax==NULL)) {
152
153 fHSigmabarMax = new TH1D;
154 fHSigmabarMax->SetNameTitle("fHSigmabarMax","Sigmabarmax for this analysis");
155
156 MH::SetBinning(fHSigmabarMax, binstheta);
157
158 // -------------------------------------------------
159 // read in SigmabarParams
160 // workaround--proprietary file format--CT1test only BEGIN
161 // -------------------------------------------------
162
163 FILE *f;
164 if( !(f =fopen(fDatabaseFilename, "r")) ) {
165 *fLog << err << dbginf << "Database file " << fDatabaseFilename << "was not found... (specify with MApplyPadding::SetDatabaseFile) aborting." << endl;
166 return kFALSE;
167 }
168 char line[80];
169 Float_t sigmabarMin, sigmabarMax, thetaMin, thetaMax, ra, dec2;
170 Int_t type, group, mjd, nr;
171 while ( fgets(line, sizeof(line), f) != NULL) {
172 if ((line[0]!='#')) {
173 sscanf(line,"%d %d %f %f %d %d %f %f %f %f",&type, &group, &ra, &dec2, &mjd, &nr, &sigmabarMin,&sigmabarMax,&thetaMin,&thetaMax);
174 if ((group==fGroup)||(type==1)) //selected ON group or OFF
175 {
176 // find out which bin(s) we have to look at
177 for (Int_t i=fHSigmabarMax->GetXaxis()->FindBin(thetaMin);
178 i<fHSigmabarMax->GetXaxis()->FindBin(thetaMax)+1; i++)
179 if (sigmabarMax > fHSigmabarMax->GetBinContent(i))
180 fHSigmabarMax->SetBinContent(i, sigmabarMax);
181 }
182 }
183 }//while
184
185 // workaround--proprietary file format--CT1test only END
186
187 // fHSigmabarMax->DrawClone();
188 // fTest = new TH2D("fTest", "Test if padding works", 201, -0.05, 2.05, 201, -0.05, 2.05);
189
190 } //!fUseHistogram
191 return kTRUE;
192}
193
194// --------------------------------------------------------------------------
195//
196// Calculate Sigmabar for current event
197// Then apply padding
198//
199// 1) have current event
200// 2) get sigmabar(theta)
201// 3) pad event
202//
203Bool_t MApplyPadding::Process()
204{
205 // Calculate sigmabar of event
206 fSigmabar->Calc(*fCam, *fPed);
207 Double_t mySig = fSigmabar->GetSigmabar();
208
209 // Get sigmabar which we have to pad to
210 Double_t otherSig = fFixedSigmabar;
211 if (fUseHistogram) {
212 Int_t binNumber = fHSigmabarMax->GetXaxis()->FindBin(fMcEvt->GetTelescopeTheta()*kRad2Deg);
213 otherSig = fHSigmabarMax->GetBinContent(binNumber);
214 }
215
216 // Determine quadratic difference other-mine
217 const Double_t quadraticDiff = otherSig*otherSig - mySig*mySig;
218
219 // crosscheck, should never happen
220 if (quadraticDiff < 0) {
221 *fLog << err << dbginf << "Event has higher Sigmabar="<<mySig<<" than Sigmabarmax="<<otherSig << " ...Skipping this event" <<endl;
222 return kCONTINUE; //skip
223 }
224
225 if (quadraticDiff == 0) return kTRUE; //no padding necessary.
226
227 // Pad if quadratic difference > 0
228 MPedestalCam newPed;
229 newPed.InitSize(fPed->GetSize());
230
231 // const UInt_t npix = fPed->GetSize();
232 const UInt_t npix = fEvt->GetNumPixels(); // Total number of pixels
233
234 for (UInt_t i=0; i<npix; i++) {
235 MCerPhotPix &pix = (*fEvt)[i];
236 if (!pix.IsPixelUsed())
237 continue;
238 pix.SetNumPhotons(pix.GetNumPhotons() +
239 sqrt(quadraticDiff)*
240 fRnd.Gaus(0.0, 1.0)/
241 fCam->GetPixRatio(pix.GetPixId())
242 );
243 // error: add sigma of padded noise quadratically
244 Double_t error = pix.GetErrorPhot();
245 pix.SetErrorPhot(sqrt(error*error + quadraticDiff));
246 MPedestalPix &ppix = (*fPed)[i];
247 newPed[i].SetSigma(sqrt(ppix.GetSigma()*ppix.GetSigma() + quadraticDiff));
248 } //for
249 // Calculate Sigmabar again and crosscheck
250 fSigmabar->Calc(*fCam, newPed);
251 //mySig = fSigmabar->GetSigmabar();
252 // fTest->Fill(otherSig,mySig);
253 return kTRUE;
254}
255
Note: See TracBrowser for help on using the repository browser.