source: trunk/MagicSoft/Mars/mhbase/MHMatrix.cc@ 3077

Last change on this file since 3077 was 2738, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 30.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! Rudy Boeck 2003 <mailto:
20! Wolfgang Wittek2003 <mailto:wittek@mppmu.mpg.de>
21!
22! Copyright: MAGIC Software Development, 2000-2003
23!
24!
25\* ======================================================================== */
26
27/////////////////////////////////////////////////////////////////////////////
28//
29// MHMatrix
30//
31// This is a histogram container which holds a matrix with one column per
32// data variable. The data variable can be a complex rule (MDataChain).
33// Each event for wich Fill is called (by MFillH) is added as a new
34// row to the matrix.
35//
36// For example:
37// MHMatrix m;
38// m.AddColumn("MHillas.fSize");
39// m.AddColumn("MMcEvt.fImpact/100");
40// m.AddColumn("HillasSource.fDist*MGeomCam.fConvMm2Deg");
41// MFillH fillm(&m);
42// taskliost.AddToList(&fillm);
43// [...]
44// m.Print();
45//
46/////////////////////////////////////////////////////////////////////////////
47#include "MHMatrix.h"
48
49#include <fstream>
50
51#include <TList.h>
52#include <TArrayF.h>
53#include <TArrayD.h>
54#include <TArrayI.h>
55
56#include <TH1.h>
57#include <TCanvas.h>
58#include <TRandom3.h>
59
60#include "MLog.h"
61#include "MLogManip.h"
62
63#include "MFillH.h"
64#include "MEvtLoop.h"
65#include "MParList.h"
66#include "MTaskList.h"
67#include "MProgressBar.h"
68
69#include "MData.h"
70#include "MDataArray.h"
71#include "MFilter.h"
72
73ClassImp(MHMatrix);
74
75using namespace std;
76
77const TString MHMatrix::gsDefName = "MHMatrix";
78const TString MHMatrix::gsDefTitle = "Multidimensional Matrix";
79
80// --------------------------------------------------------------------------
81//
82// Default Constructor
83//
84MHMatrix::MHMatrix(const char *name, const char *title)
85 : fNumRows(0), fData(NULL)
86{
87 fName = name ? name : gsDefName.Data();
88 fTitle = title ? title : gsDefTitle.Data();
89}
90
91// --------------------------------------------------------------------------
92//
93// Default Constructor
94//
95MHMatrix::MHMatrix(const TMatrix &m, const char *name, const char *title)
96 : fNumRows(m.GetNrows()), fM(m), fData(NULL)
97{
98 fName = name ? name : gsDefName.Data();
99 fTitle = title ? title : gsDefTitle.Data();
100}
101
102// --------------------------------------------------------------------------
103//
104// Constructor. Initializes the columns of the matrix with the entries
105// from a MDataArray
106//
107MHMatrix::MHMatrix(MDataArray *mat, const char *name, const char *title)
108 : fNumRows(0), fData(mat)
109{
110 fName = name ? name : gsDefName.Data();
111 fTitle = title ? title : gsDefTitle.Data();
112}
113
114// --------------------------------------------------------------------------
115//
116// Destructor. Does not deleted a user given MDataArray, except IsOwner
117// was called.
118//
119MHMatrix::~MHMatrix()
120{
121 if (TestBit(kIsOwner) && fData)
122 delete fData;
123}
124
125// --------------------------------------------------------------------------
126//
127Bool_t MHMatrix::SetNumRow(Int_t row)
128{
129 if (row>=fM.GetNrows() || row<0) return kFALSE;
130 fRow = row;
131 return kTRUE;
132}
133
134// --------------------------------------------------------------------------
135//
136// Add a new column to the matrix. This can only be done before the first
137// event (row) was filled into the matrix. For the syntax of the rule
138// see MDataChain.
139// Returns the index of the new column, -1 in case of failure.
140// (0, 1, 2, ... for the 1st, 2nd, 3rd, ...)
141//
142Int_t MHMatrix::AddColumn(const char *rule)
143{
144 if (fM.IsValid())
145 {
146 *fLog << warn << "Warning - matrix is already in use. Can't add a new column... skipped." << endl;
147 return -1;
148 }
149
150 if (TestBit(kIsLocked))
151 {
152 *fLog << warn << "Warning - matrix is locked. Can't add new column... skipped." << endl;
153 return -1;
154 }
155
156 if (!fData)
157 {
158 fData = new MDataArray;
159 SetBit(kIsOwner);
160 }
161
162 fData->AddEntry(rule);
163 return fData->GetNumEntries()-1;
164}
165
166// --------------------------------------------------------------------------
167//
168void MHMatrix::AddColumns(MDataArray *matrix)
169{
170 if (fM.IsValid())
171 {
172 *fLog << warn << "Warning - matrix is already in use. Can't add new columns... skipped." << endl;
173 return;
174 }
175
176 if (TestBit(kIsLocked))
177 {
178 *fLog << warn << "Warning - matrix is locked. Can't add new columns... skipped." << endl;
179 return;
180 }
181
182 if (fData)
183 *fLog << warn << "Warning - columns already added... replacing." << endl;
184
185 if (fData && TestBit(kIsOwner))
186 {
187 delete fData;
188 ResetBit(kIsOwner);
189 }
190
191 fData = matrix;
192}
193
194// --------------------------------------------------------------------------
195//
196// Checks whether at least one column is available and PreProcesses all
197// data chains.
198//
199Bool_t MHMatrix::SetupFill(const MParList *plist)
200{
201 if (!fData)
202 {
203 *fLog << err << "Error - No Columns initialized... aborting." << endl;
204 return kFALSE;
205 }
206
207 return fData->PreProcess(plist);
208}
209
210// --------------------------------------------------------------------------
211//
212// If the matrix has not enough rows double the number of available rows.
213//
214void MHMatrix::AddRow()
215{
216 fNumRows++;
217
218 if (fM.GetNrows() > fNumRows)
219 return;
220
221 if (!fM.IsValid())
222 {
223 fM.ResizeTo(1, fData->GetNumEntries());
224 return;
225 }
226
227#if ROOT_VERSION_CODE < ROOT_VERSION(3,05,07)
228 TMatrix m(fM);
229#endif
230 fM.ResizeTo(fM.GetNrows()*2, fData->GetNumEntries());
231
232#if ROOT_VERSION_CODE < ROOT_VERSION(3,05,07)
233 TVector vold(fM.GetNcols());
234 for (int x=0; x<m.GetNrows(); x++)
235 TMatrixRow(fM, x) = vold = TMatrixRow(m, x);
236#endif
237}
238
239// --------------------------------------------------------------------------
240//
241// Add the values correspoding to the columns to the new row
242//
243Bool_t MHMatrix::Fill(const MParContainer *par, const Stat_t w)
244{
245 AddRow();
246
247 for (int col=0; col<fData->GetNumEntries(); col++)
248 fM(fNumRows-1, col) = (*fData)(col);
249
250 return kTRUE;
251}
252
253// --------------------------------------------------------------------------
254//
255// Resize the matrix to a number of rows which corresponds to the number of
256// rows which have really been filled with values.
257//
258Bool_t MHMatrix::Finalize()
259{
260 //
261 // It's not a fatal error so we don't need to stop PostProcessing...
262 //
263 if (fData->GetNumEntries()==0 || fNumRows<1)
264 return kTRUE;
265
266 if (fNumRows != fM.GetNrows())
267 {
268 TMatrix m(fM);
269 CopyCrop(fM, m, fNumRows);
270 }
271
272 return kTRUE;
273}
274
275/*
276// --------------------------------------------------------------------------
277//
278// Draw clone of histogram. So that the object can be deleted
279// and the histogram is still visible in the canvas.
280// The cloned object are deleted together with the canvas if the canvas is
281// destroyed. If you want to handle destroying the canvas you can get a
282// pointer to it from this function
283//
284TObject *MHMatrix::DrawClone(Option_t *opt) const
285{
286 TCanvas &c = *MH::MakeDefCanvas(fHist);
287
288 //
289 // This is necessary to get the expected bahviour of DrawClone
290 //
291 gROOT->SetSelectedPad(NULL);
292
293 fHist->DrawCopy(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 c.Modified();
310 c.Update();
311
312 return &c;
313}
314
315// --------------------------------------------------------------------------
316//
317// Creates a new canvas and draws the histogram into it.
318// Be careful: The histogram belongs to this object and won't get deleted
319// together with the canvas.
320//
321void MHMatrix::Draw(Option_t *opt)
322{
323 if (!gPad)
324 MH::MakeDefCanvas(fHist);
325
326 fHist->Draw(opt);
327
328 TString str(opt);
329 if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2)
330 {
331 TProfile *p = ((TH2*)fHist)->ProfileX();
332 p->Draw("same");
333 p->SetBit(kCanDelete);
334 }
335 if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2)
336 {
337 TProfile *p = ((TH2*)fHist)->ProfileY();
338 p->Draw("same");
339 p->SetBit(kCanDelete);
340 }
341
342 gPad->Modified();
343 gPad->Update();
344}
345*/
346
347// --------------------------------------------------------------------------
348//
349// Prints the meaning of the columns and the contents of the matrix.
350// Becareful, this can take a long time for matrices with many rows.
351// Use the option 'size' to print the size of the matrix.
352// Use the option 'cols' to print the culumns
353// Use the option 'data' to print the contents
354//
355void MHMatrix::Print(Option_t *o) const
356{
357 TString str(o);
358
359 *fLog << all << flush;
360
361 if (str.Contains("size", TString::kIgnoreCase))
362 {
363 *fLog << GetDescriptor() << ": NumColumns=" << fM.GetNcols();
364 *fLog << " NumRows=" << fM.GetNrows() << endl;
365 }
366
367 if (!fData && str.Contains("cols", TString::kIgnoreCase))
368 *fLog << "Sorry, no column information available." << endl;
369
370 if (fData && str.Contains("cols", TString::kIgnoreCase))
371 fData->Print();
372
373 if (str.Contains("data", TString::kIgnoreCase))
374 fM.Print();
375}
376
377// --------------------------------------------------------------------------
378//
379const TMatrix *MHMatrix::InvertPosDef()
380{
381 TMatrix m(fM);
382
383 const Int_t rows = m.GetNrows();
384 const Int_t cols = m.GetNcols();
385
386 for (int x=0; x<cols; x++)
387 {
388 Double_t avg = 0;
389 for (int y=0; y<rows; y++)
390 avg += fM(y, x);
391
392 avg /= rows;
393
394 TMatrixColumn(m, x) += -avg;
395 }
396
397 TMatrix *m2 = new TMatrix(m, TMatrix::kTransposeMult, m);
398
399 Double_t det;
400 m2->Invert(&det);
401 if (det==0)
402 {
403 *fLog << err << "ERROR - MHMatrix::InvertPosDef failed (Matrix is singular)." << endl;
404 delete m2;
405 return NULL;
406 }
407
408 // m2->Print();
409
410 return m2;
411}
412
413// --------------------------------------------------------------------------
414//
415// Calculated the distance of vector evt from the reference sample
416// represented by the covariance metrix m.
417// - If n<0 the kernel method is applied and
418// -log(sum(epx(-d/h))/n) is returned.
419// - For n>0 the n nearest neighbors are summed and
420// sqrt(sum(d)/n) is returned.
421// - if n==0 all distances are summed
422//
423Double_t MHMatrix::CalcDist(const TMatrix &m, const TVector &evt, Int_t num) const
424{
425 if (num==0) // may later be used for another method
426 {
427 TVector d = evt;
428 d *= m;
429 return TMath::Sqrt(d*evt);
430 }
431
432 const Int_t rows = fM.GetNrows();
433 const Int_t cols = fM.GetNcols();
434
435 TArrayD dists(rows);
436
437 //
438 // Calculate: v^T * M * v
439 //
440 for (int i=0; i<rows; i++)
441 {
442 TVector col(cols);
443 col = TMatrixRow(fM, i);
444
445 TVector d = evt;
446 d -= col;
447
448 TVector d2 = d;
449 d2 *= m;
450
451 dists[i] = d2*d; // square of distance
452
453 //
454 // This corrects for numerical uncertanties in cases of very
455 // small distances...
456 //
457 if (dists[i]<0)
458 dists[i]=0;
459 }
460
461 TArrayI idx(rows);
462 TMath::Sort(dists.GetSize(), dists.GetArray(), idx.GetArray(), kFALSE);
463
464 Int_t from = 0;
465 Int_t to = TMath::Abs(num)<rows ? TMath::Abs(num) : rows;
466 //
467 // This is a zero-suppression for the case a test- and trainings
468 // sample is identical. This would result in an unwanted leading
469 // zero in the array. To suppress also numerical uncertanties of
470 // zero we cut at 1e-5. Due to Rudy this should be enough. If
471 // you encounter problems we can also use (eg) 1e-25
472 //
473 if (dists[idx[0]]<1e-5)
474 {
475 from++;
476 to ++;
477 if (to>rows)
478 to = rows;
479 }
480
481 if (num<0)
482 {
483 //
484 // Kernel function sum (window size h set according to literature)
485 //
486 const Double_t h = TMath::Power(rows, -1./(cols+4));
487 const Double_t hwin = h*h*2;
488
489 Double_t res = 0;
490 for (int i=from; i<to; i++)
491 res += TMath::Exp(-dists[idx[i]]/hwin);
492
493 return -TMath::Log(res/(to-from));
494 }
495 else
496 {
497 //
498 // Nearest Neighbor sum
499 //
500 Double_t res = 0;
501 for (int i=from; i<to; i++)
502 res += dists[idx[i]];
503
504 return TMath::Sqrt(res/(to-from));
505 }
506}
507
508// --------------------------------------------------------------------------
509//
510// Calls calc dist. In the case of the first call the covariance matrix
511// fM2 is calculated.
512// - If n<0 it is divided by (nrows-1)/h while h is the kernel factor.
513//
514Double_t MHMatrix::CalcDist(const TVector &evt, Int_t num)
515{
516 if (!fM2.IsValid())
517 {
518 if (!fM.IsValid())
519 {
520 *fLog << err << "MHMatrix::CalcDist - ERROR: fM not valid." << endl;
521 return -1;
522 }
523
524 const TMatrix *m = InvertPosDef();
525 if (!m)
526 return -1;
527
528 fM2.ResizeTo(*m);
529 fM2 = *m;
530 fM2 *= fM.GetNrows()-1;
531 delete m;
532 }
533
534 return CalcDist(fM2, evt, num);
535}
536
537// --------------------------------------------------------------------------
538//
539void MHMatrix::Reassign()
540{
541 TMatrix m = fM;
542 fM.ResizeTo(1,1);
543 fM.ResizeTo(m);
544 fM = m;
545}
546
547// --------------------------------------------------------------------------
548//
549// Implementation of SavePrimitive. Used to write the call to a constructor
550// to a macro. In the original root implementation it is used to write
551// gui elements to a macro-file.
552//
553void MHMatrix::StreamPrimitive(ofstream &out) const
554{
555 Bool_t data = fData && !TestBit(kIsOwner);
556
557 if (data)
558 {
559 fData->SavePrimitive(out);
560 out << endl;
561 }
562
563 out << " MHMatrix " << GetUniqueName();
564
565 if (data || fName!=gsDefName || fTitle!=gsDefTitle)
566 {
567 out << "(";
568 if (data)
569 out << "&" << fData->GetUniqueName();
570 if (fName!=gsDefName || fTitle!=gsDefTitle)
571 {
572 if (data)
573 out << ", ";
574 out << "\"" << fName << "\"";
575 if (fTitle!=gsDefTitle)
576 out << ", \"" << fTitle << "\"";
577 }
578 }
579 out << ");" << endl;
580
581 if (fData && TestBit(kIsOwner))
582 for (int i=0; i<fData->GetNumEntries(); i++)
583 out << " " << GetUniqueName() << ".AddColumn(\"" << (*fData)[i].GetRule() << "\");" << endl;
584}
585
586// --------------------------------------------------------------------------
587//
588const TArrayI MHMatrix::GetIndexOfSortedColumn(Int_t ncol, Bool_t desc) const
589{
590 TMatrixColumn col(fM, ncol);
591
592 const Int_t n = fM.GetNrows();
593
594 TArrayF array(n);
595
596 for (int i=0; i<n; i++)
597 array[i] = col(i);
598
599 TArrayI idx(n);
600 TMath::Sort(n, array.GetArray(), idx.GetArray(), desc);
601
602 return idx;
603}
604
605// --------------------------------------------------------------------------
606//
607void MHMatrix::SortMatrixByColumn(Int_t ncol, Bool_t desc)
608{
609 TArrayI idx = GetIndexOfSortedColumn(ncol, desc);
610
611 const Int_t n = fM.GetNrows();
612
613 TMatrix m(n, fM.GetNcols());
614 TVector vold(fM.GetNcols());
615 for (int i=0; i<n; i++)
616 TMatrixRow(m, i) = vold = TMatrixRow(fM, idx[i]);
617
618 fM = m;
619}
620
621// --------------------------------------------------------------------------
622//
623Bool_t MHMatrix::Fill(MParList *plist, MTask *read, MFilter *filter)
624{
625 //
626 // Read data into Matrix
627 //
628 const Bool_t is = plist->IsOwner();
629 plist->SetOwner(kFALSE);
630
631 MTaskList tlist;
632 plist->Replace(&tlist);
633
634 MFillH fillh(this);
635
636 tlist.AddToList(read);
637
638 if (filter)
639 {
640 tlist.AddToList(filter);
641 fillh.SetFilter(filter);
642 }
643
644 tlist.AddToList(&fillh);
645
646 //MProgressBar bar;
647 MEvtLoop evtloop("MHMatrix::Fill-EvtLoop");
648 evtloop.SetParList(plist);
649 //evtloop.SetProgressBar(&bar);
650
651 if (!evtloop.Eventloop())
652 return kFALSE;
653
654 tlist.PrintStatistics();
655
656 plist->Remove(&tlist);
657 plist->SetOwner(is);
658
659 return kTRUE;
660}
661
662// --------------------------------------------------------------------------
663//
664// Return a comma seperated list of all data members used in the matrix.
665// This is mainly used in MTask::AddToBranchList
666//
667TString MHMatrix::GetDataMember() const
668{
669 return fData ? fData->GetDataMember() : TString("");
670}
671
672// --------------------------------------------------------------------------
673//
674//
675void MHMatrix::ReduceNumberOfRows(UInt_t numrows, const TString opt)
676{
677 UInt_t rows = fM.GetNrows();
678
679 if (rows==numrows)
680 {
681 *fLog << warn << "Matrix has already the correct number of rows..." << endl;
682 return;
683 }
684
685 Float_t ratio = (Float_t)numrows/fM.GetNrows();
686
687 if (ratio>=1)
688 {
689 *fLog << warn << "Matrix cannot be enlarged..." << endl;
690 return;
691 }
692
693 Double_t sum = 0;
694
695 UInt_t oldrow = 0;
696 UInt_t newrow = 0;
697
698 TVector vold(fM.GetNcols());
699 while (oldrow<rows)
700 {
701 sum += ratio;
702
703 if (newrow<=(unsigned int)sum)
704 TMatrixRow(fM, newrow++) = vold = TMatrixRow(fM, oldrow);
705
706 oldrow++;
707 }
708}
709
710// ------------------------------------------------------------------------
711//
712// Used in DefRefMatrix to display the result graphically
713//
714void MHMatrix::DrawDefRefInfo(const TH1 &hth, const TH1 &hthd, const TH1 &thsh, Int_t refcolumn)
715{
716 //
717 // Fill a histogram with the distribution after raduction
718 //
719 TH1F hta;
720 hta.SetDirectory(NULL);
721 hta.SetName("hta");
722 hta.SetTitle("Distribution after reduction");
723 SetBinning(&hta, &hth);
724
725 for (Int_t i=0; i<fM.GetNrows(); i++)
726 hta.Fill(fM(i, refcolumn));
727
728 TCanvas *th1 = MakeDefCanvas(this);
729 th1->Divide(2,2);
730
731 th1->cd(1);
732 ((TH1&)hth).DrawCopy(); // real histogram before
733
734 th1->cd(2);
735 ((TH1&)hta).DrawCopy(); // histogram after
736
737 th1->cd(3);
738 ((TH1&)hthd).DrawCopy(); // correction factors
739
740 th1->cd(4);
741 ((TH1&)thsh).DrawCopy(); // target
742}
743
744// ------------------------------------------------------------------------
745//
746// Resizes th etarget matrix to rows*source.GetNcol() and copies
747// the data from the first (n)rows or the source into the target matrix.
748//
749void MHMatrix::CopyCrop(TMatrix &target, const TMatrix &source, Int_t rows)
750{
751 TVector v(source.GetNcols());
752
753 target.ResizeTo(rows, source.GetNcols());
754 for (Int_t ir=0; ir<rows; ir++)
755 TMatrixRow(target, ir) = v = TMatrixRow(source, ir);
756}
757
758// ------------------------------------------------------------------------
759//
760// Define the reference matrix
761// refcolumn number of the column (starting at 0) containing the variable,
762// for which a target distribution may be given;
763// thsh histogram containing the target distribution of the variable
764// nmaxevts the number of events the reference matrix should have after
765// the renormalization
766// rest a TMatrix conatining the resulting (not choosen)
767// columns of the primary matrix. Maybe NULL if you
768// are not interested in this
769//
770Bool_t MHMatrix::DefRefMatrix(const UInt_t refcolumn, const TH1F &thsh,
771 Int_t nmaxevts, TMatrix *rest)
772{
773 if (!fM.IsValid())
774 {
775 *fLog << err << dbginf << "Matrix not initialized" << endl;
776 return kFALSE;
777 }
778
779 if (thsh.GetMinimum()<0)
780 {
781 *fLog << err << dbginf << "Renormalization not possible: ";
782 *fLog << "Target Distribution has values < 0" << endl;
783 return kFALSE;
784 }
785
786
787 if (nmaxevts>fM.GetNrows())
788 {
789 *fLog << warn << dbginf << "No.requested (" << nmaxevts;
790 *fLog << ") > available events (" << fM.GetNrows() << ")... ";
791 *fLog << "setting equal." << endl;
792 nmaxevts = fM.GetNrows();
793 }
794
795
796 if (nmaxevts<0)
797 {
798 *fLog << err << dbginf << "Number of requested events < 0" << endl;
799 return kFALSE;
800 }
801
802 if (nmaxevts==0)
803 nmaxevts = fM.GetNrows();
804
805 //
806 // refcol is the column number starting at 0; it is >= 0
807 //
808 // number of the column (count from 0) containing
809 // the variable for which the target distribution is given
810 //
811
812 //
813 // Calculate normalization factors
814 //
815 //const int nbins = thsh.GetNbinsX();
816 //const double frombin = thsh.GetBinLowEdge(1);
817 //const double tobin = thsh.GetBinLowEdge(nbins+1);
818 //const double dbin = thsh.GetBinWidth(1);
819
820 const Int_t nrows = fM.GetNrows();
821 const Int_t ncols = fM.GetNcols();
822
823 //
824 // set up the real histogram (distribution before)
825 //
826 //TH1F hth("th", "Distribution before reduction", nbins, frombin, tobin);
827 TH1F hth;
828 hth.SetNameTitle("th", "Distribution before reduction");
829 SetBinning(&hth, &thsh);
830 hth.SetDirectory(NULL);
831 for (Int_t j=0; j<nrows; j++)
832 hth.Fill(fM(j, refcolumn));
833
834 //TH1F hthd("thd", "Correction factors", nbins, frombin, tobin);
835 TH1F hthd;
836 hthd.SetNameTitle("thd", "Correction factors");
837 SetBinning(&hthd, &thsh);
838 hthd.SetDirectory(NULL);
839 hthd.Divide((TH1F*)&thsh, &hth, 1, 1);
840
841 if (hthd.GetMaximum() <= 0)
842 {
843 *fLog << err << dbginf << "Maximum correction factor <= 0... abort." << endl;
844 return kFALSE;
845 }
846
847 //
848 // ===== obtain correction factors (normalization factors)
849 //
850 hthd.Scale(1/hthd.GetMaximum());
851
852 //
853 // get random access
854 //
855 TArrayI ind(nrows);
856 GetRandomArrayI(ind);
857
858 //
859 // define new matrix
860 //
861 Int_t evtcount1 = -1;
862 Int_t evtcount2 = 0;
863
864 TMatrix mnewtmp(nrows, ncols);
865 TMatrix mrest(nrows, ncols);
866
867 TArrayF cumulweight(nrows); // keep track for each bin how many events
868
869 //
870 // Project values in reference column into [0,1]
871 //
872 TVector v(fM.GetNrows());
873 v = TMatrixColumn(fM, refcolumn);
874 //v += -frombin;
875 //v *= 1/dbin;
876
877 //
878 // select events (distribution after renormalization)
879 //
880 Int_t ir;
881 TVector vold(fM.GetNcols());
882 for (ir=0; ir<nrows; ir++)
883 {
884 // const Int_t indref = (Int_t)v(ind[ir]);
885 const Int_t indref = hthd.FindBin(v(ind[ir])) - 1;
886 cumulweight[indref] += hthd.GetBinContent(indref+1);
887 if (cumulweight[indref]<=0.5)
888 {
889 TMatrixRow(mrest, evtcount2++) = vold = TMatrixRow(fM, ind[ir]);
890 continue;
891 }
892
893 cumulweight[indref] -= 1.;
894 if (++evtcount1 >= nmaxevts)
895 break;
896
897 TMatrixRow(mnewtmp, evtcount1) = vold = TMatrixRow(fM, ind[ir]);
898 }
899
900 for (/*empty*/; ir<nrows; ir++)
901 TMatrixRow(mrest, evtcount2++) = vold = TMatrixRow(fM, ind[ir]);
902
903 //
904 // reduce size
905 //
906 // matrix fM having the requested distribution
907 // and the requested number of rows;
908 // this is the matrix to be used in the g/h separation
909 //
910 CopyCrop(fM, mnewtmp, evtcount1);
911 fNumRows = evtcount1;
912
913 if (evtcount1 < nmaxevts)
914 *fLog << warn << "Reference sample contains less events (" << evtcount1 << ") than requested (" << nmaxevts << ")" << endl;
915
916 if (TestBit(kEnableGraphicalOutput))
917 DrawDefRefInfo(hth, hthd, thsh, refcolumn);
918
919 if (rest)
920 CopyCrop(*rest, mrest, evtcount2);
921
922 return kTRUE;
923}
924
925// ------------------------------------------------------------------------
926//
927// Returns a array containing randomly sorted indices
928//
929void MHMatrix::GetRandomArrayI(TArrayI &ind) const
930{
931 const Int_t rows = ind.GetSize();
932
933 TArrayF ranx(rows);
934
935 TRandom3 rnd(0);
936 for (Int_t i=0; i<rows; i++)
937 ranx[i] = rnd.Rndm(i);
938
939 TMath::Sort(rows, ranx.GetArray(), ind.GetArray(), kTRUE);
940}
941
942// ------------------------------------------------------------------------
943//
944// Define the reference matrix
945// nmaxevts maximum number of events in the reference matrix
946// rest a TMatrix conatining the resulting (not choosen)
947// columns of the primary matrix. Maybe NULL if you
948// are not interested in this
949//
950// the target distribution will be set
951// equal to the real distribution; the events in the reference
952// matrix will then be simply a random selection of the events
953// in the original matrix.
954//
955Bool_t MHMatrix::DefRefMatrix(Int_t nmaxevts, TMatrix *rest)
956{
957 if (!fM.IsValid())
958 {
959 *fLog << err << dbginf << "Matrix not initialized" << endl;
960 return kFALSE;
961 }
962
963 if (nmaxevts>fM.GetNrows())
964 {
965 *fLog << dbginf << "No.of requested events (" << nmaxevts
966 << ") exceeds no.of available events (" << fM.GetNrows()
967 << ")" << endl;
968 *fLog << dbginf
969 << " set no.of requested events = no.of available events"
970 << endl;
971 nmaxevts = fM.GetNrows();
972 }
973
974 if (nmaxevts<0)
975 {
976 *fLog << err << dbginf << "Number of requested events < 0" << endl;
977 return kFALSE;
978 }
979
980 if (nmaxevts==0)
981 nmaxevts = fM.GetNrows();
982
983 const Int_t nrows = fM.GetNrows();
984 const Int_t ncols = fM.GetNcols();
985
986 //
987 // get random access
988 //
989 TArrayI ind(nrows);
990 GetRandomArrayI(ind);
991
992 //
993 // define new matrix
994 //
995 Int_t evtcount1 = 0;
996 Int_t evtcount2 = 0;
997
998 TMatrix mnewtmp(nrows, ncols);
999 TMatrix mrest(nrows, ncols);
1000
1001 //
1002 // select events (distribution after renormalization)
1003 //
1004 TVector vold(fM.GetNcols());
1005 for (Int_t ir=0; ir<nmaxevts; ir++)
1006 TMatrixRow(mnewtmp, evtcount1++) = vold = TMatrixRow(fM, ind[ir]);
1007
1008 for (Int_t ir=nmaxevts; ir<nrows; ir++)
1009 TMatrixRow(mrest, evtcount2++) = vold = TMatrixRow(fM, ind[ir]);
1010
1011 //
1012 // reduce size
1013 //
1014 // matrix fM having the requested distribution
1015 // and the requested number of rows;
1016 // this is the matrix to be used in the g/h separation
1017 //
1018 CopyCrop(fM, mnewtmp, evtcount1);
1019 fNumRows = evtcount1;
1020
1021 if (evtcount1 < nmaxevts)
1022 *fLog << warn << "The reference sample contains less events (" << evtcount1 << ") than requested (" << nmaxevts << ")" << endl;
1023
1024 if (!rest)
1025 return kTRUE;
1026
1027 CopyCrop(*rest, mrest, evtcount2);
1028
1029 return kTRUE;
1030}
1031
1032// --------------------------------------------------------------------------
1033//
1034// overload TOject member function read
1035// in order to reset the name of the object read
1036//
1037Int_t MHMatrix::Read(const char *name)
1038{
1039 Int_t ret = TObject::Read(name);
1040 SetName(name);
1041
1042 return ret;
1043}
1044
1045// --------------------------------------------------------------------------
1046//
1047// Read the setup from a TEnv:
1048// Column0, Column1, Column2, ..., Column10, ..., Column100, ...
1049//
1050// Searching stops if the first key isn't found in the TEnv. Empty
1051// columns are not allowed
1052//
1053// eg.
1054// MHMatrix.Column0: cos(MMcEvt.fTelescopeTheta)
1055// MHMatrix.Column1: MHillasSrc.fAlpha
1056//
1057Bool_t MHMatrix::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
1058{
1059 if (fM.IsValid())
1060 {
1061 *fLog << err << "ERROR - matrix is already in use. Can't add a new column from TEnv... skipped." << endl;
1062 return kERROR;
1063 }
1064
1065 if (TestBit(kIsLocked))
1066 {
1067 *fLog << err << "ERROR - matrix is locked. Can't add new column from TEnv... skipped." << endl;
1068 return kERROR;
1069 }
1070
1071 //
1072 // Search (beginning with 0) all keys
1073 //
1074 int i=0;
1075 while (1)
1076 {
1077 TString idx = "Column";
1078 idx += i;
1079
1080 // Output if print set to kTRUE
1081 if (!IsEnvDefined(env, prefix, idx, print))
1082 break;
1083
1084 // Try to get the file name
1085 TString name = GetEnvValue(env, prefix, idx, "");
1086 if (name.IsNull())
1087 {
1088 *fLog << warn << prefix+"."+idx << " empty." << endl;
1089 continue;
1090 }
1091
1092 if (i==0)
1093 if (fData)
1094 {
1095 *fLog << inf << "Removing all existing columns in " << GetDescriptor() << endl;
1096 fData->Delete();
1097 }
1098 else
1099 {
1100 fData = new MDataArray;
1101 SetBit(kIsOwner);
1102 }
1103
1104 fData->AddEntry(name);
1105 i++;
1106 }
1107
1108 return i!=0;
1109}
1110
1111// --------------------------------------------------------------------------
1112//
1113// ShuffleEvents. Shuffles the order of the matrix rows.
1114//
1115//
1116void MHMatrix::ShuffleRows(UInt_t seed)
1117{
1118 TRandom rnd(seed);
1119
1120 TVector v(fM.GetNcols());
1121 TVector tmp(fM.GetNcols());
1122 for (Int_t irow = 0; irow<fNumRows; irow++)
1123 {
1124 const Int_t jrow = rnd.Integer(fNumRows);
1125
1126 v = TMatrixRow(fM, irow);
1127 TMatrixRow(fM, irow) = tmp = TMatrixRow(fM, jrow);
1128 TMatrixRow(fM, jrow) = v;
1129 }
1130
1131 *fLog << warn << GetDescriptor() << ": Attention! Matrix rows have been shuffled." << endl;
1132}
1133
1134// --------------------------------------------------------------------------
1135//
1136// Reduces the number of rows to the given number num by cutting out the
1137// last rows.
1138//
1139void MHMatrix::ReduceRows(UInt_t num)
1140{
1141 if ((Int_t)num>=fM.GetNrows())
1142 {
1143 *fLog << warn << GetDescriptor() << ": Warning - " << num << " >= rows=" << fM.GetNrows() << endl;
1144 return;
1145 }
1146
1147#if ROOT_VERSION_CODE < ROOT_VERSION(3,05,07)
1148 const TMatrix m(fM);
1149#endif
1150 fM.ResizeTo(num, fM.GetNcols());
1151
1152#if ROOT_VERSION_CODE < ROOT_VERSION(3,05,07)
1153 TVector tmp(fM.GetNcols());
1154 for (UInt_t irow=0; irow<num; irow++)
1155 TMatrixRow(fM, irow) = tmp = TMatrixRow(m, irow);
1156#endif
1157}
1158
1159// --------------------------------------------------------------------------
1160//
1161// Remove rows which contains numbers not fullfilling TMath::Finite
1162//
1163Bool_t MHMatrix::RemoveInvalidRows()
1164{
1165 TMatrix m(fM);
1166
1167 const Int_t ncol=fM.GetNcols();
1168 TVector vold(ncol);
1169 int irow=0;
1170
1171 for (int i=0; i<m.GetNrows(); i++)
1172 {
1173 const TMatrixRow &row = TMatrixRow(m, i);
1174
1175 // finite (-> math.h) checks for NaN as well as inf
1176 int jcol;
1177 for (jcol=0; jcol<ncol; jcol++)
1178 if (!TMath::Finite(vold(jcol)))
1179 break;
1180
1181 if (jcol==ncol)
1182 TMatrixRow(fM, irow++) = vold = row;
1183 else
1184 *fLog << warn << "Warning - MHMatrix::RemoveInvalidRows: row #" << i<< " removed." << endl;
1185 }
1186
1187 // Do not use ResizeTo (in older root versions this doesn't save the contents
1188 ReduceRows(irow);
1189
1190 return kTRUE;
1191}
Note: See TracBrowser for help on using the repository browser.