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

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