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

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