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

Last change on this file since 5303 was 5100, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 32.3 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 {
372 fData->SetName(fName);
373 fData->Print();
374 }
375
376 if (str.Contains("data", TString::kIgnoreCase))
377 fM.Print();
378}
379
380// --------------------------------------------------------------------------
381//
382const TMatrix *MHMatrix::InvertPosDef()
383{
384 TMatrix m(fM);
385
386 const Int_t rows = m.GetNrows();
387 const Int_t cols = m.GetNcols();
388
389 for (int x=0; x<cols; x++)
390 {
391 Double_t avg = 0;
392 for (int y=0; y<rows; y++)
393 avg += fM(y, x);
394
395 avg /= rows;
396
397 TMatrixColumn(m, x) += -avg;
398 }
399
400 TMatrix *m2 = new TMatrix(m, TMatrix::kTransposeMult, m);
401
402 Double_t det;
403 m2->Invert(&det);
404 if (det==0)
405 {
406 *fLog << err << "ERROR - MHMatrix::InvertPosDef failed (Matrix is singular)." << endl;
407 delete m2;
408 return NULL;
409 }
410
411 // m2->Print();
412
413 return m2;
414}
415
416// --------------------------------------------------------------------------
417//
418// Calculated the distance of vector evt from the reference sample
419// represented by the covariance metrix m.
420// - If n<0 the kernel method is applied and
421// -log(sum(epx(-d/h))/n) is returned.
422// - For n>0 the n nearest neighbors are summed and
423// sqrt(sum(d)/n) is returned.
424// - if n==0 all distances are summed
425//
426Double_t MHMatrix::CalcDist(const TMatrix &m, const TVector &evt, Int_t num) const
427{
428 if (num==0) // may later be used for another method
429 {
430 TVector d = evt;
431 d *= m;
432 return TMath::Sqrt(d*evt);
433 }
434
435 const Int_t rows = fM.GetNrows();
436 const Int_t cols = fM.GetNcols();
437
438 TArrayD dists(rows);
439
440 //
441 // Calculate: v^T * M * v
442 //
443 for (int i=0; i<rows; i++)
444 {
445 TVector col(cols);
446#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
447 col = TMatrixRow(fM, i);
448#else
449 col = TMatrixFRow_const(fM, i);
450#endif
451 TVector d = evt;
452 d -= col;
453
454 TVector d2 = d;
455 d2 *= m;
456
457 dists[i] = d2*d; // square of distance
458
459 //
460 // This corrects for numerical uncertanties in cases of very
461 // small distances...
462 //
463 if (dists[i]<0)
464 dists[i]=0;
465 }
466
467 TArrayI idx(rows);
468 TMath::Sort(dists.GetSize(), dists.GetArray(), idx.GetArray(), kFALSE);
469
470 Int_t from = 0;
471 Int_t to = TMath::Abs(num)<rows ? TMath::Abs(num) : rows;
472 //
473 // This is a zero-suppression for the case a test- and trainings
474 // sample is identical. This would result in an unwanted leading
475 // zero in the array. To suppress also numerical uncertanties of
476 // zero we cut at 1e-5. Due to Rudy this should be enough. If
477 // you encounter problems we can also use (eg) 1e-25
478 //
479 if (dists[idx[0]]<1e-5)
480 {
481 from++;
482 to ++;
483 if (to>rows)
484 to = rows;
485 }
486
487 if (num<0)
488 {
489 //
490 // Kernel function sum (window size h set according to literature)
491 //
492 const Double_t h = TMath::Power(rows, -1./(cols+4));
493 const Double_t hwin = h*h*2;
494
495 Double_t res = 0;
496 for (int i=from; i<to; i++)
497 res += TMath::Exp(-dists[idx[i]]/hwin);
498
499 return -TMath::Log(res/(to-from));
500 }
501 else
502 {
503 //
504 // Nearest Neighbor sum
505 //
506 Double_t res = 0;
507 for (int i=from; i<to; i++)
508 res += dists[idx[i]];
509
510 return TMath::Sqrt(res/(to-from));
511 }
512}
513
514// --------------------------------------------------------------------------
515//
516// Calls calc dist. In the case of the first call the covariance matrix
517// fM2 is calculated.
518// - If n<0 it is divided by (nrows-1)/h while h is the kernel factor.
519//
520Double_t MHMatrix::CalcDist(const TVector &evt, Int_t num)
521{
522 if (!fM2.IsValid())
523 {
524 if (!fM.IsValid())
525 {
526 *fLog << err << "MHMatrix::CalcDist - ERROR: fM not valid." << endl;
527 return -1;
528 }
529
530 const TMatrix *m = InvertPosDef();
531 if (!m)
532 return -1;
533
534 fM2.ResizeTo(*m);
535 fM2 = *m;
536 fM2 *= fM.GetNrows()-1;
537 delete m;
538 }
539
540 return CalcDist(fM2, evt, num);
541}
542
543// --------------------------------------------------------------------------
544//
545void MHMatrix::Reassign()
546{
547 TMatrix m = fM;
548 fM.ResizeTo(1,1);
549 fM.ResizeTo(m);
550 fM = m;
551}
552
553// --------------------------------------------------------------------------
554//
555// Implementation of SavePrimitive. Used to write the call to a constructor
556// to a macro. In the original root implementation it is used to write
557// gui elements to a macro-file.
558//
559void MHMatrix::StreamPrimitive(ofstream &out) const
560{
561 Bool_t data = fData && !TestBit(kIsOwner);
562
563 if (data)
564 {
565 fData->SavePrimitive(out);
566 out << endl;
567 }
568
569 out << " MHMatrix " << GetUniqueName();
570
571 if (data || fName!=gsDefName || fTitle!=gsDefTitle)
572 {
573 out << "(";
574 if (data)
575 out << "&" << fData->GetUniqueName();
576 if (fName!=gsDefName || fTitle!=gsDefTitle)
577 {
578 if (data)
579 out << ", ";
580 out << "\"" << fName << "\"";
581 if (fTitle!=gsDefTitle)
582 out << ", \"" << fTitle << "\"";
583 }
584 }
585 out << ");" << endl;
586
587 if (fData && TestBit(kIsOwner))
588 for (int i=0; i<fData->GetNumEntries(); i++)
589 out << " " << GetUniqueName() << ".AddColumn(\"" << (*fData)[i].GetRule() << "\");" << endl;
590}
591
592// --------------------------------------------------------------------------
593//
594const TArrayI MHMatrix::GetIndexOfSortedColumn(Int_t ncol, Bool_t desc) const
595{
596#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
597 TMatrixColumn col(fM, ncol);
598#else
599 TMatrixFColumn_const col(fM, ncol);
600#endif
601
602 const Int_t n = fM.GetNrows();
603
604 TArrayF array(n);
605
606 for (int i=0; i<n; i++)
607 array[i] = col(i);
608
609 TArrayI idx(n);
610 TMath::Sort(n, array.GetArray(), idx.GetArray(), desc);
611
612 return idx;
613}
614
615// --------------------------------------------------------------------------
616//
617void MHMatrix::SortMatrixByColumn(Int_t ncol, Bool_t desc)
618{
619 TArrayI idx = GetIndexOfSortedColumn(ncol, desc);
620
621 const Int_t n = fM.GetNrows();
622
623#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
624 TVector vold(fM.GetNcols());
625#endif
626
627 TMatrix m(n, fM.GetNcols());
628 for (int i=0; i<n; i++)
629#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
630 TMatrixRow(m, i) = vold = TMatrixRow(fM, idx[i]);
631#else
632 TMatrixFRow(m, i) = TMatrixFRow_const(fM, idx[i]);
633#endif
634 fM = m;
635}
636
637// --------------------------------------------------------------------------
638//
639Bool_t MHMatrix::Fill(MParList *plist, MTask *read, MFilter *filter)
640{
641 //
642 // Read data into Matrix
643 //
644 const Bool_t is = plist->IsOwner();
645 plist->SetOwner(kFALSE);
646
647 MTaskList tlist;
648 plist->Replace(&tlist);
649
650 MFillH fillh(this);
651
652 tlist.AddToList(read);
653
654 if (filter)
655 {
656 tlist.AddToList(filter);
657 fillh.SetFilter(filter);
658 }
659
660 tlist.AddToList(&fillh);
661
662 //MProgressBar bar;
663 MEvtLoop evtloop("MHMatrix::Fill-EvtLoop");
664 evtloop.SetParList(plist);
665 //evtloop.SetProgressBar(&bar);
666
667 if (!evtloop.Eventloop())
668 return kFALSE;
669
670 tlist.PrintStatistics();
671
672 plist->Remove(&tlist);
673 plist->SetOwner(is);
674
675 return kTRUE;
676}
677
678// --------------------------------------------------------------------------
679//
680// Return a comma seperated list of all data members used in the matrix.
681// This is mainly used in MTask::AddToBranchList
682//
683TString MHMatrix::GetDataMember() const
684{
685 return fData ? fData->GetDataMember() : TString("");
686}
687
688// --------------------------------------------------------------------------
689//
690//
691void MHMatrix::ReduceNumberOfRows(UInt_t numrows, const TString opt)
692{
693 UInt_t rows = fM.GetNrows();
694
695 if (rows==numrows)
696 {
697 *fLog << warn << "Matrix has already the correct number of rows..." << endl;
698 return;
699 }
700
701 Float_t ratio = (Float_t)numrows/fM.GetNrows();
702
703 if (ratio>=1)
704 {
705 *fLog << warn << "Matrix cannot be enlarged..." << endl;
706 return;
707 }
708
709 Double_t sum = 0;
710
711 UInt_t oldrow = 0;
712 UInt_t newrow = 0;
713
714#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
715 TVector vold(fM.GetNcols());
716#endif
717 while (oldrow<rows)
718 {
719 sum += ratio;
720
721 if (newrow<=(unsigned int)sum)
722#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
723 TMatrixRow(fM, newrow++) = vold = TMatrixRow(fM, oldrow);
724#else
725 TMatrixFRow(fM, newrow++) = TMatrixFRow_const(fM, oldrow);
726#endif
727
728 oldrow++;
729 }
730}
731
732// ------------------------------------------------------------------------
733//
734// Used in DefRefMatrix to display the result graphically
735//
736void MHMatrix::DrawDefRefInfo(const TH1 &hth, const TH1 &hthd, const TH1 &thsh, Int_t refcolumn)
737{
738 //
739 // Fill a histogram with the distribution after raduction
740 //
741 TH1F hta;
742 hta.SetDirectory(NULL);
743 hta.SetName("hta");
744 hta.SetTitle("Distribution after reduction");
745 SetBinning(&hta, &hth);
746
747 for (Int_t i=0; i<fM.GetNrows(); i++)
748 hta.Fill(fM(i, refcolumn));
749
750 TCanvas *th1 = MakeDefCanvas(this);
751 th1->Divide(2,2);
752
753 th1->cd(1);
754 ((TH1&)hth).DrawCopy(); // real histogram before
755
756 th1->cd(2);
757 ((TH1&)hta).DrawCopy(); // histogram after
758
759 th1->cd(3);
760 ((TH1&)hthd).DrawCopy(); // correction factors
761
762 th1->cd(4);
763 ((TH1&)thsh).DrawCopy(); // target
764}
765
766// ------------------------------------------------------------------------
767//
768// Resizes th etarget matrix to rows*source.GetNcol() and copies
769// the data from the first (n)rows or the source into the target matrix.
770//
771void MHMatrix::CopyCrop(TMatrix &target, const TMatrix &source, Int_t rows)
772{
773#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
774 TVector v(source.GetNcols());
775#endif
776 target.ResizeTo(rows, source.GetNcols());
777 for (Int_t ir=0; ir<rows; ir++)
778#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
779 TMatrixRow(target, ir) = v = TMatrixRow(source, ir);
780#else
781 TMatrixFRow(target, ir) = TMatrixFRow_const(source, ir);
782#endif
783}
784
785// ------------------------------------------------------------------------
786//
787// Define the reference matrix
788// refcolumn number of the column (starting at 0) containing the variable,
789// for which a target distribution may be given;
790// thsh histogram containing the target distribution of the variable
791// nmaxevts the number of events the reference matrix should have after
792// the renormalization
793// rest a TMatrix conatining the resulting (not choosen)
794// columns of the primary matrix. Maybe NULL if you
795// are not interested in this
796//
797Bool_t MHMatrix::DefRefMatrix(const UInt_t refcolumn, const TH1F &thsh,
798 Int_t nmaxevts, TMatrix *rest)
799{
800 if (!fM.IsValid())
801 {
802 *fLog << err << dbginf << "Matrix not initialized" << endl;
803 return kFALSE;
804 }
805
806 if (thsh.GetMinimum()<0)
807 {
808 *fLog << err << dbginf << "Renormalization not possible: ";
809 *fLog << "Target Distribution has values < 0" << endl;
810 return kFALSE;
811 }
812
813
814 if (nmaxevts>fM.GetNrows())
815 {
816 *fLog << warn << dbginf << "No.requested (" << nmaxevts;
817 *fLog << ") > available events (" << fM.GetNrows() << ")... ";
818 *fLog << "setting equal." << endl;
819 nmaxevts = fM.GetNrows();
820 }
821
822
823 if (nmaxevts<0)
824 {
825 *fLog << err << dbginf << "Number of requested events < 0" << endl;
826 return kFALSE;
827 }
828
829 if (nmaxevts==0)
830 nmaxevts = fM.GetNrows();
831
832 //
833 // refcol is the column number starting at 0; it is >= 0
834 //
835 // number of the column (count from 0) containing
836 // the variable for which the target distribution is given
837 //
838
839 //
840 // Calculate normalization factors
841 //
842 //const int nbins = thsh.GetNbinsX();
843 //const double frombin = thsh.GetBinLowEdge(1);
844 //const double tobin = thsh.GetBinLowEdge(nbins+1);
845 //const double dbin = thsh.GetBinWidth(1);
846
847 const Int_t nrows = fM.GetNrows();
848 const Int_t ncols = fM.GetNcols();
849
850 //
851 // set up the real histogram (distribution before)
852 //
853 //TH1F hth("th", "Distribution before reduction", nbins, frombin, tobin);
854 TH1F hth;
855 hth.SetNameTitle("th", "Distribution before reduction");
856 SetBinning(&hth, &thsh);
857 hth.SetDirectory(NULL);
858 for (Int_t j=0; j<nrows; j++)
859 hth.Fill(fM(j, refcolumn));
860
861 //TH1F hthd("thd", "Correction factors", nbins, frombin, tobin);
862 TH1F hthd;
863 hthd.SetNameTitle("thd", "Correction factors");
864 SetBinning(&hthd, &thsh);
865 hthd.SetDirectory(NULL);
866 hthd.Divide((TH1F*)&thsh, &hth, 1, 1);
867
868 if (hthd.GetMaximum() <= 0)
869 {
870 *fLog << err << dbginf << "Maximum correction factor <= 0... abort." << endl;
871 return kFALSE;
872 }
873
874 //
875 // ===== obtain correction factors (normalization factors)
876 //
877 hthd.Scale(1/hthd.GetMaximum());
878
879 //
880 // get random access
881 //
882 TArrayI ind(nrows);
883 GetRandomArrayI(ind);
884
885 //
886 // define new matrix
887 //
888 Int_t evtcount1 = -1;
889 Int_t evtcount2 = 0;
890
891 TMatrix mnewtmp(nrows, ncols);
892 TMatrix mrest(nrows, ncols);
893
894 TArrayF cumulweight(nrows); // keep track for each bin how many events
895
896 //
897 // Project values in reference column into [0,1]
898 //
899 TVector v(fM.GetNrows());
900 v = TMatrixColumn(fM, refcolumn);
901 //v += -frombin;
902 //v *= 1/dbin;
903
904 //
905 // select events (distribution after renormalization)
906 //
907 Int_t ir;
908#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
909 TVector vold(fM.GetNcols());
910#endif
911 for (ir=0; ir<nrows; ir++)
912 {
913 // const Int_t indref = (Int_t)v(ind[ir]);
914 const Int_t indref = hthd.FindBin(v(ind[ir])) - 1;
915 cumulweight[indref] += hthd.GetBinContent(indref+1);
916 if (cumulweight[indref]<=0.5)
917 {
918#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
919 TMatrixRow(mrest, evtcount2++) = vold = TMatrixRow(fM, ind[ir]);
920#else
921 TMatrixFRow(mrest, evtcount2++) = TMatrixFRow_const(fM, ind[ir]);
922#endif
923 continue;
924 }
925
926 cumulweight[indref] -= 1.;
927 if (++evtcount1 >= nmaxevts)
928 break;
929
930#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
931 TMatrixRow(mnewtmp, evtcount1) = vold = TMatrixRow(fM, ind[ir]);
932#else
933 TMatrixFRow(mnewtmp, evtcount1) = TMatrixFRow_const(fM, ind[ir]);
934#endif
935 }
936
937 for (/*empty*/; ir<nrows; ir++)
938#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
939 TMatrixRow(mrest, evtcount2++) = vold = TMatrixRow(fM, ind[ir]);
940#else
941 TMatrixFRow(mrest, evtcount2++) = TMatrixFRow_const(fM, ind[ir]);
942#endif
943
944 //
945 // reduce size
946 //
947 // matrix fM having the requested distribution
948 // and the requested number of rows;
949 // this is the matrix to be used in the g/h separation
950 //
951 CopyCrop(fM, mnewtmp, evtcount1);
952 fNumRows = evtcount1;
953
954 if (evtcount1 < nmaxevts)
955 *fLog << warn << "Reference sample contains less events (" << evtcount1 << ") than requested (" << nmaxevts << ")" << endl;
956
957 if (TestBit(kEnableGraphicalOutput))
958 DrawDefRefInfo(hth, hthd, thsh, refcolumn);
959
960 if (rest)
961 CopyCrop(*rest, mrest, evtcount2);
962
963 return kTRUE;
964}
965
966// ------------------------------------------------------------------------
967//
968// Returns a array containing randomly sorted indices
969//
970void MHMatrix::GetRandomArrayI(TArrayI &ind) const
971{
972 const Int_t rows = ind.GetSize();
973
974 TArrayF ranx(rows);
975
976 TRandom3 rnd(0);
977 for (Int_t i=0; i<rows; i++)
978 ranx[i] = rnd.Rndm(i);
979
980 TMath::Sort(rows, ranx.GetArray(), ind.GetArray(), kTRUE);
981}
982
983// ------------------------------------------------------------------------
984//
985// Define the reference matrix
986// nmaxevts maximum number of events in the reference matrix
987// rest a TMatrix conatining the resulting (not choosen)
988// columns of the primary matrix. Maybe NULL if you
989// are not interested in this
990//
991// the target distribution will be set
992// equal to the real distribution; the events in the reference
993// matrix will then be simply a random selection of the events
994// in the original matrix.
995//
996Bool_t MHMatrix::DefRefMatrix(Int_t nmaxevts, TMatrix *rest)
997{
998 if (!fM.IsValid())
999 {
1000 *fLog << err << dbginf << "Matrix not initialized" << endl;
1001 return kFALSE;
1002 }
1003
1004 if (nmaxevts>fM.GetNrows())
1005 {
1006 *fLog << dbginf << "No.of requested events (" << nmaxevts
1007 << ") exceeds no.of available events (" << fM.GetNrows()
1008 << ")" << endl;
1009 *fLog << dbginf
1010 << " set no.of requested events = no.of available events"
1011 << endl;
1012 nmaxevts = fM.GetNrows();
1013 }
1014
1015 if (nmaxevts<0)
1016 {
1017 *fLog << err << dbginf << "Number of requested events < 0" << endl;
1018 return kFALSE;
1019 }
1020
1021 if (nmaxevts==0)
1022 nmaxevts = fM.GetNrows();
1023
1024 const Int_t nrows = fM.GetNrows();
1025 const Int_t ncols = fM.GetNcols();
1026
1027 //
1028 // get random access
1029 //
1030 TArrayI ind(nrows);
1031 GetRandomArrayI(ind);
1032
1033 //
1034 // define new matrix
1035 //
1036 Int_t evtcount1 = 0;
1037 Int_t evtcount2 = 0;
1038
1039 TMatrix mnewtmp(nrows, ncols);
1040 TMatrix mrest(nrows, ncols);
1041
1042 //
1043 // select events (distribution after renormalization)
1044 //
1045#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
1046 TVector vold(fM.GetNcols());
1047#endif
1048 for (Int_t ir=0; ir<nmaxevts; ir++)
1049#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
1050 TMatrixRow(mnewtmp, evtcount1++) = vold = TMatrixRow(fM, ind[ir]);
1051#else
1052 TMatrixFRow(mnewtmp, evtcount1++) = TMatrixFRow_const(fM, ind[ir]);
1053#endif
1054
1055 for (Int_t ir=nmaxevts; ir<nrows; ir++)
1056#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
1057 TMatrixRow(mrest, evtcount2++) = vold = TMatrixRow(fM, ind[ir]);
1058#else
1059 TMatrixFRow(mrest, evtcount2++) = TMatrixFRow_const(fM, ind[ir]);
1060#endif
1061
1062 //
1063 // reduce size
1064 //
1065 // matrix fM having the requested distribution
1066 // and the requested number of rows;
1067 // this is the matrix to be used in the g/h separation
1068 //
1069 CopyCrop(fM, mnewtmp, evtcount1);
1070 fNumRows = evtcount1;
1071
1072 if (evtcount1 < nmaxevts)
1073 *fLog << warn << "The reference sample contains less events (" << evtcount1 << ") than requested (" << nmaxevts << ")" << endl;
1074
1075 if (!rest)
1076 return kTRUE;
1077
1078 CopyCrop(*rest, mrest, evtcount2);
1079
1080 return kTRUE;
1081}
1082
1083// --------------------------------------------------------------------------
1084//
1085// overload TOject member function read
1086// in order to reset the name of the object read
1087//
1088Int_t MHMatrix::Read(const char *name)
1089{
1090 Int_t ret = TObject::Read(name);
1091 SetName(name);
1092
1093 return ret;
1094}
1095
1096// --------------------------------------------------------------------------
1097//
1098// Read the setup from a TEnv:
1099// Column0, Column1, Column2, ..., Column10, ..., Column100, ...
1100//
1101// Searching stops if the first key isn't found in the TEnv. Empty
1102// columns are not allowed
1103//
1104// eg.
1105// MHMatrix.Column0: cos(MMcEvt.fTelescopeTheta)
1106// MHMatrix.Column1: MHillasSrc.fAlpha
1107//
1108Int_t MHMatrix::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
1109{
1110 if (fM.IsValid())
1111 {
1112 *fLog << err << "ERROR - matrix is already in use. Can't add a new column from TEnv... skipped." << endl;
1113 return kERROR;
1114 }
1115
1116 if (TestBit(kIsLocked))
1117 {
1118 *fLog << err << "ERROR - matrix is locked. Can't add new column from TEnv... skipped." << endl;
1119 return kERROR;
1120 }
1121
1122 //
1123 // Search (beginning with 0) all keys
1124 //
1125 int i=0;
1126 while (1)
1127 {
1128 TString idx = "Column";
1129 idx += i;
1130
1131 // Output if print set to kTRUE
1132 if (!IsEnvDefined(env, prefix, idx, print))
1133 break;
1134
1135 // Try to get the file name
1136 TString name = GetEnvValue(env, prefix, idx, "");
1137 if (name.IsNull())
1138 {
1139 *fLog << warn << prefix+"."+idx << " empty." << endl;
1140 continue;
1141 }
1142
1143 if (i==0)
1144 if (fData)
1145 {
1146 *fLog << inf << "Removing all existing columns in " << GetDescriptor() << endl;
1147 fData->Delete();
1148 }
1149 else
1150 {
1151 fData = new MDataArray;
1152 SetBit(kIsOwner);
1153 }
1154
1155 fData->AddEntry(name);
1156 i++;
1157 }
1158
1159 return i!=0;
1160}
1161
1162// --------------------------------------------------------------------------
1163//
1164// ShuffleEvents. Shuffles the order of the matrix rows.
1165//
1166//
1167void MHMatrix::ShuffleRows(UInt_t seed)
1168{
1169 TRandom rnd(seed);
1170
1171 TVector v(fM.GetNcols());
1172#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
1173 TVector tmp(fM.GetNcols());
1174#endif
1175 for (Int_t irow = 0; irow<fNumRows; irow++)
1176 {
1177 const Int_t jrow = rnd.Integer(fNumRows);
1178
1179#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
1180 v = TMatrixRow(fM, irow);
1181 TMatrixRow(fM, irow) = tmp = TMatrixRow(fM, jrow);
1182 TMatrixRow(fM, jrow) = v;
1183#else
1184 v = TMatrixFRow_const(fM, irow);
1185 TMatrixFRow(fM, irow) = TMatrixFRow_const(fM, jrow);
1186 TMatrixFRow(fM, jrow) = v;
1187#endif
1188 }
1189
1190 *fLog << warn << GetDescriptor() << ": Attention! Matrix rows have been shuffled." << endl;
1191}
1192
1193// --------------------------------------------------------------------------
1194//
1195// Reduces the number of rows to the given number num by cutting out the
1196// last rows.
1197//
1198void MHMatrix::ReduceRows(UInt_t num)
1199{
1200 if ((Int_t)num>=fM.GetNrows())
1201 {
1202 *fLog << warn << GetDescriptor() << ": Warning - " << num << " >= rows=" << fM.GetNrows() << endl;
1203 return;
1204 }
1205
1206#if ROOT_VERSION_CODE < ROOT_VERSION(3,05,07)
1207 const TMatrix m(fM);
1208#endif
1209 fM.ResizeTo(num, fM.GetNcols());
1210
1211#if ROOT_VERSION_CODE < ROOT_VERSION(3,05,07)
1212 TVector tmp(fM.GetNcols());
1213 for (UInt_t irow=0; irow<num; irow++)
1214 TMatrixRow(fM, irow) = tmp = TMatrixRow(m, irow);
1215#endif
1216}
1217
1218// --------------------------------------------------------------------------
1219//
1220// Remove rows which contains numbers not fullfilling TMath::Finite
1221//
1222Bool_t MHMatrix::RemoveInvalidRows()
1223{
1224 TMatrix m(fM);
1225
1226 const Int_t ncol=fM.GetNcols();
1227
1228#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
1229 TVector vold(ncol);
1230#endif
1231
1232 int irow=0;
1233
1234 for (int i=0; i<m.GetNrows(); i++)
1235 {
1236 const TMatrixRow &row = TMatrixRow(m, i);
1237
1238 // finite (-> math.h) checks for NaN as well as inf
1239 int jcol;
1240 for (jcol=0; jcol<ncol; jcol++)
1241 if (!TMath::Finite(row(jcol)))
1242 break;
1243
1244 if (jcol==ncol)
1245#if ROOT_VERSION_CODE < ROOT_VERSION(4,00,8)
1246 TMatrixRow(fM, irow++) = vold = row;
1247#else
1248 TMatrixFRow(fM, irow++) = row;
1249#endif
1250 else
1251 *fLog << warn << "Warning - MHMatrix::RemoveInvalidRows: row #" << i<< " removed." << endl;
1252 }
1253
1254 // Do not use ResizeTo (in older root versions this doesn't save the contents
1255 ReduceRows(irow);
1256
1257 return kTRUE;
1258}
Note: See TracBrowser for help on using the repository browser.