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

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