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