source: tags/Mars-V1.0/mhbase/MHMatrix.cc

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