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

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