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

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