source: trunk/MagicSoft/Mars/mhist/MHMatrix.cc@ 1828

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