source: trunk/MagicSoft/Mars/mhist/MHMatrix.cc@ 1388

Last change on this file since 1388 was 1355, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 9.4 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 2002 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2002
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MHMatrix
28//
29// This is a histogram container which holds a matrix with one column per
30// data variable. The data variable can be a complex rule (MDataChain).
31// Each event for wich Fill is called (by MFillH) is added as a new
32// row to the matrix.
33//
34// For example:
35// MHMatrix m;
36// m.AddColumn("MHillas.fSize");
37// m.AddColumn("MMcEvt.fImpact/100");
38// m.AddColumn("HillasSource.fDist*MGeomCam.fConvMm2Deg");
39// MFillH fillm(&m);
40// taskliost.AddToList(&fillm);
41// [...]
42// m.Print();
43//
44/////////////////////////////////////////////////////////////////////////////
45#include "MHMatrix.h"
46
47#include <TList.h>
48#include <TArrayD.h>
49#include <TArrayI.h>
50
51#include "MLog.h"
52#include "MLogManip.h"
53
54#include "MParList.h"
55
56#include "MDataChain.h"
57
58ClassImp(MHMatrix);
59
60// --------------------------------------------------------------------------
61//
62// Default Constructor
63//
64MHMatrix::MHMatrix(const char *name, const char *title)
65 : fNumRow(0)
66{
67 fName = name ? name : "MHMatrix";
68 fTitle = title ? title : "Multidimensional Matrix";
69
70 fData = new TList;
71 fData->SetOwner();
72
73 fRules = new TList;
74 fRules->SetOwner();
75}
76
77// --------------------------------------------------------------------------
78//
79// Destructor
80//
81MHMatrix::~MHMatrix()
82{
83 delete fData;
84 delete fRules;
85}
86
87// --------------------------------------------------------------------------
88//
89// Add a new column to the matrix. This can only be done before the first
90// event (row) was filled into the matrix. For the syntax of the rule
91// see MDataChain.
92//
93void MHMatrix::AddColumn(const char *rule)
94{
95 if (fM.IsValid())
96 {
97 *fLog << warn << "Warning - matrix is already in use. Can't add a new column... skipped." << endl;
98 return;
99 }
100
101 MDataChain &chain = *new MDataChain(rule);
102 if (!chain.IsValid())
103 {
104 *fLog << err << "Error - Rule cannot be translated... ignored." << endl;
105 delete &chain;
106 return;
107 }
108 fData->Add(&chain);
109
110 TNamed *name = new TNamed(rule, rule); // Fimxe, in 3.02/07 the title can't be "", why?
111 fRules->Add(name);
112}
113
114void MHMatrix::AddColumns(const MHMatrix *matrix)
115{
116 TIter Next(matrix->fRules);
117 TObject *rule=NULL;
118 while ((rule=Next()))
119 AddColumn(rule->GetName());
120}
121
122// --------------------------------------------------------------------------
123//
124// Checks whether at least one column is available and PreProcesses all
125// data chains.
126//
127Bool_t MHMatrix::SetupFill(const MParList *plist)
128{
129 if (fData->GetSize()==0)
130 {
131 *fLog << err << "Error - No Column specified... aborting." << endl;
132 return kFALSE;
133 }
134
135 TIter Next(fData);
136 MData *data = NULL;
137 while ((data=(MData*)Next()))
138 if (!data->PreProcess(plist))
139 return kFALSE;
140
141 return kTRUE;
142}
143
144// --------------------------------------------------------------------------
145//
146// If the matrix has not enough rows double the number of available rows.
147//
148void MHMatrix::AddRow()
149{
150 fNumRow++;
151
152 if (fM.GetNrows() > fNumRow)
153 return;
154
155 if (!fM.IsValid())
156 {
157 fM.ResizeTo(1, fData->GetSize());
158 return;
159 }
160
161 TMatrix m(fM);
162
163 fM.ResizeTo(fM.GetNrows()*2, fData->GetSize());
164
165 for (int x=0; x<m.GetNcols(); x++)
166 for (int y=0; y<m.GetNrows(); y++)
167 fM(y, x) = m(y, x);
168}
169
170// --------------------------------------------------------------------------
171//
172// Add the values correspoding to the columns to the new row
173//
174Bool_t MHMatrix::Fill(const MParContainer *par)
175{
176 AddRow();
177
178 int col=0;
179 TIter Next(fData);
180 MData *data = NULL;
181 while ((data=(MData*)Next()))
182 fM(fNumRow-1, col++) = data->GetValue();
183
184 return kTRUE;
185}
186
187// --------------------------------------------------------------------------
188//
189// Resize the matrix to a number of rows which corresponds to the number of
190// rows which have really been filled with values.
191//
192Bool_t MHMatrix::Finalize()
193{
194 //
195 // It's not a fatal error so we don't need to stop PostProcessing...
196 //
197 if (fData->GetSize()<1 || fNumRow<1)
198 return kTRUE;
199
200 TMatrix m(fM);
201
202 fM.ResizeTo(fNumRow, fData->GetSize());
203
204 for (int x=0; x<fM.GetNcols(); x++)
205 for (int y=0; y<fM.GetNrows(); y++)
206 fM(y, x) = m(y, x);
207
208 return kTRUE;
209}
210/*
211// --------------------------------------------------------------------------
212//
213// Draw clone of histogram. So that the object can be deleted
214// and the histogram is still visible in the canvas.
215// The cloned object are deleted together with the canvas if the canvas is
216// destroyed. If you want to handle destroying the canvas you can get a
217// pointer to it from this function
218//
219TObject *MHMatrix::DrawClone(Option_t *opt) const
220{
221 TCanvas &c = *MH::MakeDefCanvas(fHist);
222
223 //
224 // This is necessary to get the expected bahviour of DrawClone
225 //
226 gROOT->SetSelectedPad(NULL);
227
228 fHist->DrawCopy(opt);
229
230 TString str(opt);
231 if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2)
232 {
233 TProfile *p = ((TH2*)fHist)->ProfileX();
234 p->Draw("same");
235 p->SetBit(kCanDelete);
236 }
237 if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2)
238 {
239 TProfile *p = ((TH2*)fHist)->ProfileY();
240 p->Draw("same");
241 p->SetBit(kCanDelete);
242 }
243
244 c.Modified();
245 c.Update();
246
247 return &c;
248}
249
250// --------------------------------------------------------------------------
251//
252// Creates a new canvas and draws the histogram into it.
253// Be careful: The histogram belongs to this object and won't get deleted
254// together with the canvas.
255//
256void MHMatrix::Draw(Option_t *opt)
257{
258 if (!gPad)
259 MH::MakeDefCanvas(fHist);
260
261 fHist->Draw(opt);
262
263 TString str(opt);
264 if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2)
265 {
266 TProfile *p = ((TH2*)fHist)->ProfileX();
267 p->Draw("same");
268 p->SetBit(kCanDelete);
269 }
270 if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2)
271 {
272 TProfile *p = ((TH2*)fHist)->ProfileY();
273 p->Draw("same");
274 p->SetBit(kCanDelete);
275 }
276
277 gPad->Modified();
278 gPad->Update();
279}
280*/
281
282// --------------------------------------------------------------------------
283//
284// Prints the meaning of the columns and the contents of the matrix.
285// Becareful, this can take a long time for matrices with many rows.
286//
287void MHMatrix::Print(Option_t *) const
288{
289 int n=0;
290
291 TIter Next(fData->GetSize() ? fData : fRules);
292 MData *data = NULL;
293 while ((data=(MData*)Next()))
294 {
295 *fLog << all << " Column " << setw(3) << n++ << ": " << flush;
296 data->Print();
297 *fLog << endl;
298 }
299
300 fM.Print();
301}
302
303const TMatrix *MHMatrix::InvertPosDef()
304{
305 TMatrix m(fM);
306
307 const Int_t rows = m.GetNrows();
308 const Int_t cols = m.GetNcols();
309
310 for (int x=0; x<cols; x++)
311 {
312 Double_t avg = 0;
313 for (int y=0; y<rows; y++)
314 avg += fM(y, x);
315
316 avg /= rows;
317
318 for (int y=0; y<rows; y++)
319 m(y, x) -= avg;
320 }
321
322 TMatrix *m2 = new TMatrix(m, TMatrix::kTransposeMult, m);
323
324 Double_t det;
325 m2->Invert(&det);
326 if (det==0)
327 {
328 *fLog << err << "ERROR - MHMatrix::InvertPosDef failed (Matrix is sigular)." << endl;
329 delete m2;
330 return NULL;
331 }
332
333 // m2->Print();
334
335 return m2;
336}
337
338Double_t MHMatrix::CalcDist(const TMatrix &m, const TVector &evt, Int_t num) const
339{
340 const Int_t rows = fM.GetNrows();
341 const Int_t cols = fM.GetNcols();
342
343 TArrayD dists(rows);
344
345 //
346 // Calculate: v^T * M * v
347 //
348 for (int i=0; i<rows; i++)
349 {
350 TVector col(cols);
351 col = TMatrixRow(fM, i);
352
353 TVector d = evt;
354 d -= col;
355
356 TVector d2 = d;
357 d2 *= m;
358
359 dists[i] = d2*d; // square of distance
360 }
361
362 TArrayI idx(rows);
363 TMath::Sort(dists.GetSize(), dists.GetArray(), idx.GetArray(), kFALSE);
364
365 const Int_t n = num<rows ? num : rows;
366
367 Double_t res = 0;
368 for (int i=0; i<n; i++)
369 res += dists[idx[i]];
370
371 return sqrt(res)/n;
372}
373
374Double_t MHMatrix::CalcDist(const TVector &evt, Int_t num)
375{
376 if (!fM2.IsValid())
377 {
378 const TMatrix &m = *InvertPosDef();
379 fM2.ResizeTo(m);
380 fM2 = m;
381 delete &m;
382 }
383
384 return CalcDist(fM2, evt, num);
385}
386
387void MHMatrix::Reassign()
388{
389 TMatrix m = fM;
390 fM.ResizeTo(1,1);
391 fM.ResizeTo(m);
392 fM = m;
393
394 TIter Next(fRules);
395 TObject *obj = NULL;
396 while ((obj=Next()))
397 fData->Add(new MDataChain(obj->GetName()));
398}
Note: See TracBrowser for help on using the repository browser.