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

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