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

Last change on this file since 1522 was 1511, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 12.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#include <math.h> // for Alphas
47
48#include <fstream.h>
49
50#include <TList.h>
51#include <TArrayF.h>
52#include <TArrayD.h>
53#include <TArrayI.h>
54
55#include "MLog.h"
56#include "MLogManip.h"
57
58#include "MParList.h"
59
60#include "MData.h"
61#include "MDataArray.h"
62
63ClassImp(MHMatrix);
64
65static const TString gsDefName = "MHMatrix";
66static const TString gsDefTitle = "Multidimensional Matrix";
67
68// --------------------------------------------------------------------------
69//
70// Default Constructor
71//
72MHMatrix::MHMatrix(const char *name, const char *title)
73 : fNumRow(0), fData(NULL)
74{
75 fName = name ? name : gsDefName.Data();
76 fTitle = title ? title : gsDefTitle.Data();
77}
78
79// --------------------------------------------------------------------------
80//
81// Constructor. Initializes the columns of the matrix with the entries
82// from a MDataArray
83//
84MHMatrix::MHMatrix(MDataArray *mat, const char *name, const char *title)
85 : fNumRow(0), fData(mat)
86{
87 fName = name ? name : gsDefName.Data();
88 fTitle = title ? title : gsDefTitle.Data();
89}
90
91// --------------------------------------------------------------------------
92//
93// Destructor. Does not deleted a user given MDataArray, except IsOwner
94// was called.
95//
96MHMatrix::~MHMatrix()
97{
98 if (TestBit(kIsOwner) && fData)
99 delete fData;
100}
101
102// --------------------------------------------------------------------------
103//
104// Add a new column to the matrix. This can only be done before the first
105// event (row) was filled into the matrix. For the syntax of the rule
106// see MDataChain.
107//
108void MHMatrix::AddColumn(const char *rule)
109{
110 if (fM.IsValid())
111 {
112 *fLog << warn << "Warning - matrix is already in use. Can't add a new column... skipped." << endl;
113 return;
114 }
115
116 if (!fData)
117 {
118 fData = new MDataArray;
119 SetBit(kIsOwner);
120 }
121
122 fData->AddEntry(rule);
123}
124
125void MHMatrix::AddColumns(MDataArray *matrix)
126{
127 if (fM.IsValid())
128 {
129 *fLog << warn << "Warning - matrix is already in use. Can't add new columns... skipped." << endl;
130 return;
131 }
132
133 if (fData)
134 *fLog << warn << "Warning - columns already added... replacing." << endl;
135
136 if (fData && TestBit(kIsOwner))
137 {
138 delete fData;
139 ResetBit(kIsOwner);
140 }
141
142 fData = matrix;
143}
144
145// --------------------------------------------------------------------------
146//
147// Checks whether at least one column is available and PreProcesses all
148// data chains.
149//
150Bool_t MHMatrix::SetupFill(const MParList *plist)
151{
152 if (!fData)
153 {
154 *fLog << err << "Error - No Columns initialized... aborting." << endl;
155 return kFALSE;
156 }
157
158 return fData->PreProcess(plist);
159}
160
161// --------------------------------------------------------------------------
162//
163// If the matrix has not enough rows double the number of available rows.
164//
165void MHMatrix::AddRow()
166{
167 fNumRow++;
168
169 if (fM.GetNrows() > fNumRow)
170 return;
171
172 if (!fM.IsValid())
173 {
174 fM.ResizeTo(1, fData->GetNumEntries());
175 return;
176 }
177
178 TMatrix m(fM);
179
180 fM.ResizeTo(fM.GetNrows()*2, fData->GetNumEntries());
181
182 for (int x=0; x<m.GetNcols(); x++)
183 for (int y=0; y<m.GetNrows(); y++)
184 fM(y, x) = m(y, x);
185}
186
187// --------------------------------------------------------------------------
188//
189// Add the values correspoding to the columns to the new row
190//
191Bool_t MHMatrix::Fill(const MParContainer *par)
192{
193 AddRow();
194
195 for (int col=0; col<fData->GetNumEntries(); col++)
196 fM(fNumRow-1, col) = (*fData)(col);
197
198 return kTRUE;
199}
200
201// --------------------------------------------------------------------------
202//
203// Resize the matrix to a number of rows which corresponds to the number of
204// rows which have really been filled with values.
205//
206Bool_t MHMatrix::Finalize()
207{
208 //
209 // It's not a fatal error so we don't need to stop PostProcessing...
210 //
211 if (fData->GetNumEntries()==0 || fNumRow<1)
212 return kTRUE;
213
214 TMatrix m(fM);
215
216 fM.ResizeTo(fNumRow, fData->GetNumEntries());
217
218 for (int x=0; x<fM.GetNcols(); x++)
219 for (int y=0; y<fM.GetNrows(); y++)
220 fM(y, x) = m(y, x);
221
222 return kTRUE;
223}
224/*
225// --------------------------------------------------------------------------
226//
227// Draw clone of histogram. So that the object can be deleted
228// and the histogram is still visible in the canvas.
229// The cloned object are deleted together with the canvas if the canvas is
230// destroyed. If you want to handle destroying the canvas you can get a
231// pointer to it from this function
232//
233TObject *MHMatrix::DrawClone(Option_t *opt) const
234{
235 TCanvas &c = *MH::MakeDefCanvas(fHist);
236
237 //
238 // This is necessary to get the expected bahviour of DrawClone
239 //
240 gROOT->SetSelectedPad(NULL);
241
242 fHist->DrawCopy(opt);
243
244 TString str(opt);
245 if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2)
246 {
247 TProfile *p = ((TH2*)fHist)->ProfileX();
248 p->Draw("same");
249 p->SetBit(kCanDelete);
250 }
251 if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2)
252 {
253 TProfile *p = ((TH2*)fHist)->ProfileY();
254 p->Draw("same");
255 p->SetBit(kCanDelete);
256 }
257
258 c.Modified();
259 c.Update();
260
261 return &c;
262}
263
264// --------------------------------------------------------------------------
265//
266// Creates a new canvas and draws the histogram into it.
267// Be careful: The histogram belongs to this object and won't get deleted
268// together with the canvas.
269//
270void MHMatrix::Draw(Option_t *opt)
271{
272 if (!gPad)
273 MH::MakeDefCanvas(fHist);
274
275 fHist->Draw(opt);
276
277 TString str(opt);
278 if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2)
279 {
280 TProfile *p = ((TH2*)fHist)->ProfileX();
281 p->Draw("same");
282 p->SetBit(kCanDelete);
283 }
284 if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2)
285 {
286 TProfile *p = ((TH2*)fHist)->ProfileY();
287 p->Draw("same");
288 p->SetBit(kCanDelete);
289 }
290
291 gPad->Modified();
292 gPad->Update();
293}
294*/
295
296// --------------------------------------------------------------------------
297//
298// Prints the meaning of the columns and the contents of the matrix.
299// Becareful, this can take a long time for matrices with many rows.
300//
301void MHMatrix::Print(Option_t *) const
302{
303 if (fData)
304 fData->Print();
305 else
306 *fLog << all << "Sorry, no column information available." << endl;
307
308 fM.Print();
309}
310
311const TMatrix *MHMatrix::InvertPosDef()
312{
313 TMatrix m(fM);
314
315 const Int_t rows = m.GetNrows();
316 const Int_t cols = m.GetNcols();
317
318 for (int x=0; x<cols; x++)
319 {
320 Double_t avg = 0;
321 for (int y=0; y<rows; y++)
322 avg += fM(y, x);
323
324 avg /= rows;
325
326 for (int y=0; y<rows; y++)
327 m(y, x) -= avg;
328 }
329
330 TMatrix *m2 = new TMatrix(m, TMatrix::kTransposeMult, m);
331
332 Double_t det;
333 m2->Invert(&det);
334 if (det==0)
335 {
336 *fLog << err << "ERROR - MHMatrix::InvertPosDef failed (Matrix is sigular)." << endl;
337 delete m2;
338 return NULL;
339 }
340
341 // m2->Print();
342
343 return m2;
344}
345
346// --------------------------------------------------------------------------
347//
348// Calculated the distance of vector evt from the reference sample
349// represented by the covariance metrix m.
350// - If n<0 the kernel method is applied and
351// sum(epx(-d)/n is returned.
352// - For n>=0 the n nearest neighbors are summed and
353// sqrt(sum(d))/n is returned.
354// - if n==0 all distances are summed
355//
356Double_t MHMatrix::CalcDist(const TMatrix &m, const TVector &evt, Int_t num) const
357{
358 if (num==0) // may later be used for another method
359 {
360 TVector d = evt;
361 d *= m;
362 return sqrt(d*evt);
363 }
364
365 const Int_t rows = fM.GetNrows();
366 const Int_t cols = fM.GetNcols();
367
368 TArrayD dists(rows);
369
370 //
371 // Calculate: v^T * M * v
372 //
373 for (int i=0; i<rows; i++)
374 {
375 TVector col(cols);
376 col = TMatrixRow(fM, i);
377
378 TVector d = evt;
379 d -= col;
380
381 TVector d2 = d;
382 d2 *= m;
383
384 dists[i] = d2*d; // square of distance
385
386 //
387 // This corrects for numerical uncertanties in cases of very
388 // small distances...
389 //
390 if (dists[i]<0)
391 dists[i]=0;
392 }
393
394 TArrayI idx(rows);
395 TMath::Sort(dists.GetSize(), dists.GetArray(), idx.GetArray(), kFALSE);
396
397 const Int_t n = abs(num)<rows ? abs(num) : rows;
398
399 if (num<0)
400 {
401 //
402 // Kernel function sum
403 //
404 const Double_t h = 1;
405
406 Double_t res = 0;
407 for (int i=0; i<n; i++)
408 res += exp(-dists[idx[i]]/h);
409
410 return log(res/n);
411 }
412 else
413 {
414 //
415 // Nearest Neighbor sum
416 //
417 Double_t res = 0;
418 for (int i=0; i<n; i++)
419 res += dists[idx[i]];
420
421 return sqrt(res/n);
422 }
423}
424
425// --------------------------------------------------------------------------
426//
427// Calls calc dist. In the case of the first call the covariance matrix
428// fM2 is calculated.
429// - If n<0 it is divided by (nrows-1)/h while h is the kernel factor.
430//
431Double_t MHMatrix::CalcDist(const TVector &evt, Int_t num)
432{
433 if (!fM2.IsValid())
434 {
435 const TMatrix &m = *InvertPosDef();
436 fM2.ResizeTo(m);
437 fM2 = m;
438 fM2 *= fM.GetNrows()-1;
439 delete &m;
440 }
441
442 return CalcDist(fM2, evt, num);
443}
444
445void MHMatrix::Reassign()
446{
447 TMatrix m = fM;
448 fM.ResizeTo(1,1);
449 fM.ResizeTo(m);
450 fM = m;
451}
452
453void MHMatrix::StreamPrimitive(ofstream &out) const
454{
455 Bool_t data = fData && !TestBit(kIsOwner);
456
457 if (data)
458 {
459 fData->SavePrimitive(out);
460 out << endl;
461 }
462
463 out << " MHMatrix " << GetUniqueName();
464
465 if (data || fName!=gsDefName || fTitle!=gsDefTitle)
466 {
467 out << "(";
468 if (data)
469 out << "&" << fData->GetUniqueName();
470 if (fName!=gsDefName || fTitle!=gsDefTitle)
471 {
472 if (data)
473 out << ", ";
474 out << "\"" << fName << "\"";
475 if (fTitle!=gsDefTitle)
476 out << ", \"" << fTitle << "\"";
477 }
478 }
479 out << ");" << endl;
480
481 if (fData && TestBit(kIsOwner))
482 for (int i=0; i<fData->GetNumEntries(); i++)
483 out << " " << GetUniqueName() << ".AddColumn(\"" << (*fData)[i].GetRule() << "\");" << endl;
484}
485
486const TArrayI MHMatrix::GetIndexOfSortedColumn(Int_t ncol, Bool_t desc) const
487{
488 TMatrixColumn col(fM, ncol);
489
490 const Int_t n = fM.GetNrows();
491
492 TArrayF array(n);
493
494 for (int i=0; i<n; i++)
495 array[i] = col(i);
496
497 TArrayI idx(n);
498 TMath::Sort(n, array.GetArray(), idx.GetArray(), desc);
499
500 return idx;
501}
502
503void MHMatrix::SortMatrixByColumn(Int_t ncol, Bool_t desc)
504{
505 TArrayI idx = GetIndexOfSortedColumn(ncol, desc);
506
507 const Int_t n = fM.GetNrows();
508
509 TMatrix m(n, fM.GetNcols());
510 for (int i=0; i<n; i++)
511 {
512 TVector vold(n);
513 vold = TMatrixRow(fM, idx[i]);
514
515 TMatrixRow rownew(m, i);
516 rownew = vold;
517 }
518
519 fM = m;
520}
Note: See TracBrowser for help on using the repository browser.