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

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