| 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 | !
|
|---|
| 20 | ! Copyright: MAGIC Software Development, 2000-2002
|
|---|
| 21 | !
|
|---|
| 22 | !
|
|---|
| 23 | \* ======================================================================== */
|
|---|
| 24 |
|
|---|
| 25 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 26 | //
|
|---|
| 27 | // MHMatrix
|
|---|
| 28 | //
|
|---|
| 29 | // This is a histogram container which holds a matrix with one column per
|
|---|
| 30 | // data variable. The data variable can be a complex rule (MDataChain).
|
|---|
| 31 | // Each event for wich Fill is called (by MFillH) is added as a new
|
|---|
| 32 | // row to the matrix.
|
|---|
| 33 | //
|
|---|
| 34 | // For example:
|
|---|
| 35 | // MHMatrix m;
|
|---|
| 36 | // m.AddColumn("MHillas.fSize");
|
|---|
| 37 | // m.AddColumn("MMcEvt.fImpact/100");
|
|---|
| 38 | // m.AddColumn("HillasSource.fDist*MGeomCam.fConvMm2Deg");
|
|---|
| 39 | // MFillH fillm(&m);
|
|---|
| 40 | // taskliost.AddToList(&fillm);
|
|---|
| 41 | // [...]
|
|---|
| 42 | // m.Print();
|
|---|
| 43 | //
|
|---|
| 44 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 45 | #include "MHMatrix.h"
|
|---|
| 46 |
|
|---|
| 47 | #include <fstream.h>
|
|---|
| 48 |
|
|---|
| 49 | #include <TList.h>
|
|---|
| 50 | #include <TArrayF.h>
|
|---|
| 51 | #include <TArrayD.h>
|
|---|
| 52 | #include <TArrayI.h>
|
|---|
| 53 |
|
|---|
| 54 | #include "MLog.h"
|
|---|
| 55 | #include "MLogManip.h"
|
|---|
| 56 |
|
|---|
| 57 | #include "MFillH.h"
|
|---|
| 58 | #include "MEvtLoop.h"
|
|---|
| 59 | #include "MParList.h"
|
|---|
| 60 | #include "MTaskList.h"
|
|---|
| 61 |
|
|---|
| 62 | #include "MData.h"
|
|---|
| 63 | #include "MDataArray.h"
|
|---|
| 64 |
|
|---|
| 65 | ClassImp(MHMatrix);
|
|---|
| 66 |
|
|---|
| 67 | const TString MHMatrix::gsDefName = "MHMatrix";
|
|---|
| 68 | const TString MHMatrix::gsDefTitle = "Multidimensional Matrix";
|
|---|
| 69 |
|
|---|
| 70 | // --------------------------------------------------------------------------
|
|---|
| 71 | //
|
|---|
| 72 | // Default Constructor
|
|---|
| 73 | //
|
|---|
| 74 | MHMatrix::MHMatrix(const char *name, const char *title)
|
|---|
| 75 | : fNumRow(0), fData(NULL)
|
|---|
| 76 | {
|
|---|
| 77 | fName = name ? name : gsDefName.Data();
|
|---|
| 78 | fTitle = title ? title : gsDefTitle.Data();
|
|---|
| 79 | }
|
|---|
| 80 |
|
|---|
| 81 | // --------------------------------------------------------------------------
|
|---|
| 82 | //
|
|---|
| 83 | // Constructor. Initializes the columns of the matrix with the entries
|
|---|
| 84 | // from a MDataArray
|
|---|
| 85 | //
|
|---|
| 86 | MHMatrix::MHMatrix(MDataArray *mat, const char *name, const char *title)
|
|---|
| 87 | : fNumRow(0), fData(mat)
|
|---|
| 88 | {
|
|---|
| 89 | fName = name ? name : gsDefName.Data();
|
|---|
| 90 | fTitle = title ? title : gsDefTitle.Data();
|
|---|
| 91 | }
|
|---|
| 92 |
|
|---|
| 93 | // --------------------------------------------------------------------------
|
|---|
| 94 | //
|
|---|
| 95 | // Destructor. Does not deleted a user given MDataArray, except IsOwner
|
|---|
| 96 | // was called.
|
|---|
| 97 | //
|
|---|
| 98 | MHMatrix::~MHMatrix()
|
|---|
| 99 | {
|
|---|
| 100 | if (TestBit(kIsOwner) && fData)
|
|---|
| 101 | delete fData;
|
|---|
| 102 | }
|
|---|
| 103 |
|
|---|
| 104 | // --------------------------------------------------------------------------
|
|---|
| 105 | //
|
|---|
| 106 | // Add a new column to the matrix. This can only be done before the first
|
|---|
| 107 | // event (row) was filled into the matrix. For the syntax of the rule
|
|---|
| 108 | // see MDataChain.
|
|---|
| 109 | // Returns the index of the new column, -1 in case of failure.
|
|---|
| 110 | // (0, 1, 2, ... for the 1st, 2nd, 3rd, ...)
|
|---|
| 111 | //
|
|---|
| 112 | Int_t MHMatrix::AddColumn(const char *rule)
|
|---|
| 113 | {
|
|---|
| 114 | if (fM.IsValid())
|
|---|
| 115 | {
|
|---|
| 116 | *fLog << warn << "Warning - matrix is already in use. Can't add a new column... skipped." << endl;
|
|---|
| 117 | return -1;
|
|---|
| 118 | }
|
|---|
| 119 |
|
|---|
| 120 | if (TestBit(kIsLocked))
|
|---|
| 121 | {
|
|---|
| 122 | *fLog << warn << "Warning - matrix is locked. Can't add new column... skipped." << endl;
|
|---|
| 123 | return -1;
|
|---|
| 124 | }
|
|---|
| 125 |
|
|---|
| 126 | if (!fData)
|
|---|
| 127 | {
|
|---|
| 128 | fData = new MDataArray;
|
|---|
| 129 | SetBit(kIsOwner);
|
|---|
| 130 | }
|
|---|
| 131 |
|
|---|
| 132 | fData->AddEntry(rule);
|
|---|
| 133 | return fData->GetNumEntries()-1;
|
|---|
| 134 | }
|
|---|
| 135 |
|
|---|
| 136 | void MHMatrix::AddColumns(MDataArray *matrix)
|
|---|
| 137 | {
|
|---|
| 138 | if (fM.IsValid())
|
|---|
| 139 | {
|
|---|
| 140 | *fLog << warn << "Warning - matrix is already in use. Can't add new columns... skipped." << endl;
|
|---|
| 141 | return;
|
|---|
| 142 | }
|
|---|
| 143 |
|
|---|
| 144 | if (TestBit(kIsLocked))
|
|---|
| 145 | {
|
|---|
| 146 | *fLog << warn << "Warning - matrix is locked. Can't add new columns... skipped." << endl;
|
|---|
| 147 | return;
|
|---|
| 148 | }
|
|---|
| 149 |
|
|---|
| 150 | if (fData)
|
|---|
| 151 | *fLog << warn << "Warning - columns already added... replacing." << endl;
|
|---|
| 152 |
|
|---|
| 153 | if (fData && TestBit(kIsOwner))
|
|---|
| 154 | {
|
|---|
| 155 | delete fData;
|
|---|
| 156 | ResetBit(kIsOwner);
|
|---|
| 157 | }
|
|---|
| 158 |
|
|---|
| 159 | fData = matrix;
|
|---|
| 160 | }
|
|---|
| 161 |
|
|---|
| 162 | // --------------------------------------------------------------------------
|
|---|
| 163 | //
|
|---|
| 164 | // Checks whether at least one column is available and PreProcesses all
|
|---|
| 165 | // data chains.
|
|---|
| 166 | //
|
|---|
| 167 | Bool_t MHMatrix::SetupFill(const MParList *plist)
|
|---|
| 168 | {
|
|---|
| 169 | if (!fData)
|
|---|
| 170 | {
|
|---|
| 171 | *fLog << err << "Error - No Columns initialized... aborting." << endl;
|
|---|
| 172 | return kFALSE;
|
|---|
| 173 | }
|
|---|
| 174 |
|
|---|
| 175 | return fData->PreProcess(plist);
|
|---|
| 176 | }
|
|---|
| 177 |
|
|---|
| 178 | // --------------------------------------------------------------------------
|
|---|
| 179 | //
|
|---|
| 180 | // If the matrix has not enough rows double the number of available rows.
|
|---|
| 181 | //
|
|---|
| 182 | void MHMatrix::AddRow()
|
|---|
| 183 | {
|
|---|
| 184 | fNumRow++;
|
|---|
| 185 |
|
|---|
| 186 | if (fM.GetNrows() > fNumRow)
|
|---|
| 187 | return;
|
|---|
| 188 |
|
|---|
| 189 | if (!fM.IsValid())
|
|---|
| 190 | {
|
|---|
| 191 | fM.ResizeTo(1, fData->GetNumEntries());
|
|---|
| 192 | return;
|
|---|
| 193 | }
|
|---|
| 194 |
|
|---|
| 195 | TMatrix m(fM);
|
|---|
| 196 |
|
|---|
| 197 | fM.ResizeTo(fM.GetNrows()*2, fData->GetNumEntries());
|
|---|
| 198 |
|
|---|
| 199 | for (int x=0; x<m.GetNcols(); x++)
|
|---|
| 200 | for (int y=0; y<m.GetNrows(); y++)
|
|---|
| 201 | fM(y, x) = m(y, x);
|
|---|
| 202 | }
|
|---|
| 203 |
|
|---|
| 204 | // --------------------------------------------------------------------------
|
|---|
| 205 | //
|
|---|
| 206 | // Add the values correspoding to the columns to the new row
|
|---|
| 207 | //
|
|---|
| 208 | Bool_t MHMatrix::Fill(const MParContainer *par)
|
|---|
| 209 | {
|
|---|
| 210 | AddRow();
|
|---|
| 211 |
|
|---|
| 212 | for (int col=0; col<fData->GetNumEntries(); col++)
|
|---|
| 213 | fM(fNumRow-1, col) = (*fData)(col);
|
|---|
| 214 |
|
|---|
| 215 | return kTRUE;
|
|---|
| 216 | }
|
|---|
| 217 |
|
|---|
| 218 | // --------------------------------------------------------------------------
|
|---|
| 219 | //
|
|---|
| 220 | // Resize the matrix to a number of rows which corresponds to the number of
|
|---|
| 221 | // rows which have really been filled with values.
|
|---|
| 222 | //
|
|---|
| 223 | Bool_t MHMatrix::Finalize()
|
|---|
| 224 | {
|
|---|
| 225 | //
|
|---|
| 226 | // It's not a fatal error so we don't need to stop PostProcessing...
|
|---|
| 227 | //
|
|---|
| 228 | if (fData->GetNumEntries()==0 || fNumRow<1)
|
|---|
| 229 | return kTRUE;
|
|---|
| 230 |
|
|---|
| 231 | TMatrix m(fM);
|
|---|
| 232 |
|
|---|
| 233 | fM.ResizeTo(fNumRow, fData->GetNumEntries());
|
|---|
| 234 |
|
|---|
| 235 | for (int x=0; x<fM.GetNcols(); x++)
|
|---|
| 236 | for (int y=0; y<fM.GetNrows(); y++)
|
|---|
| 237 | fM(y, x) = m(y, x);
|
|---|
| 238 |
|
|---|
| 239 | return kTRUE;
|
|---|
| 240 | }
|
|---|
| 241 | /*
|
|---|
| 242 | // --------------------------------------------------------------------------
|
|---|
| 243 | //
|
|---|
| 244 | // Draw clone of histogram. So that the object can be deleted
|
|---|
| 245 | // and the histogram is still visible in the canvas.
|
|---|
| 246 | // The cloned object are deleted together with the canvas if the canvas is
|
|---|
| 247 | // destroyed. If you want to handle destroying the canvas you can get a
|
|---|
| 248 | // pointer to it from this function
|
|---|
| 249 | //
|
|---|
| 250 | TObject *MHMatrix::DrawClone(Option_t *opt) const
|
|---|
| 251 | {
|
|---|
| 252 | TCanvas &c = *MH::MakeDefCanvas(fHist);
|
|---|
| 253 |
|
|---|
| 254 | //
|
|---|
| 255 | // This is necessary to get the expected bahviour of DrawClone
|
|---|
| 256 | //
|
|---|
| 257 | gROOT->SetSelectedPad(NULL);
|
|---|
| 258 |
|
|---|
| 259 | fHist->DrawCopy(opt);
|
|---|
| 260 |
|
|---|
| 261 | TString str(opt);
|
|---|
| 262 | if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2)
|
|---|
| 263 | {
|
|---|
| 264 | TProfile *p = ((TH2*)fHist)->ProfileX();
|
|---|
| 265 | p->Draw("same");
|
|---|
| 266 | p->SetBit(kCanDelete);
|
|---|
| 267 | }
|
|---|
| 268 | if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2)
|
|---|
| 269 | {
|
|---|
| 270 | TProfile *p = ((TH2*)fHist)->ProfileY();
|
|---|
| 271 | p->Draw("same");
|
|---|
| 272 | p->SetBit(kCanDelete);
|
|---|
| 273 | }
|
|---|
| 274 |
|
|---|
| 275 | c.Modified();
|
|---|
| 276 | c.Update();
|
|---|
| 277 |
|
|---|
| 278 | return &c;
|
|---|
| 279 | }
|
|---|
| 280 |
|
|---|
| 281 | // --------------------------------------------------------------------------
|
|---|
| 282 | //
|
|---|
| 283 | // Creates a new canvas and draws the histogram into it.
|
|---|
| 284 | // Be careful: The histogram belongs to this object and won't get deleted
|
|---|
| 285 | // together with the canvas.
|
|---|
| 286 | //
|
|---|
| 287 | void MHMatrix::Draw(Option_t *opt)
|
|---|
| 288 | {
|
|---|
| 289 | if (!gPad)
|
|---|
| 290 | MH::MakeDefCanvas(fHist);
|
|---|
| 291 |
|
|---|
| 292 | fHist->Draw(opt);
|
|---|
| 293 |
|
|---|
| 294 | TString str(opt);
|
|---|
| 295 | if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2)
|
|---|
| 296 | {
|
|---|
| 297 | TProfile *p = ((TH2*)fHist)->ProfileX();
|
|---|
| 298 | p->Draw("same");
|
|---|
| 299 | p->SetBit(kCanDelete);
|
|---|
| 300 | }
|
|---|
| 301 | if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2)
|
|---|
| 302 | {
|
|---|
| 303 | TProfile *p = ((TH2*)fHist)->ProfileY();
|
|---|
| 304 | p->Draw("same");
|
|---|
| 305 | p->SetBit(kCanDelete);
|
|---|
| 306 | }
|
|---|
| 307 |
|
|---|
| 308 | gPad->Modified();
|
|---|
| 309 | gPad->Update();
|
|---|
| 310 | }
|
|---|
| 311 | */
|
|---|
| 312 |
|
|---|
| 313 | // --------------------------------------------------------------------------
|
|---|
| 314 | //
|
|---|
| 315 | // Prints the meaning of the columns and the contents of the matrix.
|
|---|
| 316 | // Becareful, this can take a long time for matrices with many rows.
|
|---|
| 317 | // Use the option 'size' to print the size of the matrix.
|
|---|
| 318 | // Use the option 'cols' to print the culumns
|
|---|
| 319 | // Use the option 'data' to print the contents
|
|---|
| 320 | //
|
|---|
| 321 | void MHMatrix::Print(Option_t *o) const
|
|---|
| 322 | {
|
|---|
| 323 | TString str(o);
|
|---|
| 324 |
|
|---|
| 325 | *fLog << all << flush;
|
|---|
| 326 |
|
|---|
| 327 | if (str.Contains("size", TString::kIgnoreCase))
|
|---|
| 328 | {
|
|---|
| 329 | *fLog << GetDescriptor() << ": NumColumns=" << fM.GetNcols();
|
|---|
| 330 | *fLog << " NumRows=" << fM.GetNrows() << endl;
|
|---|
| 331 | }
|
|---|
| 332 |
|
|---|
| 333 | if (!fData && str.Contains("cols", TString::kIgnoreCase))
|
|---|
| 334 | *fLog << "Sorry, no column information available." << endl;
|
|---|
| 335 |
|
|---|
| 336 | if (fData && str.Contains("cols", TString::kIgnoreCase))
|
|---|
| 337 | fData->Print();
|
|---|
| 338 |
|
|---|
| 339 | if (str.Contains("data", TString::kIgnoreCase))
|
|---|
| 340 | fM.Print();
|
|---|
| 341 | }
|
|---|
| 342 |
|
|---|
| 343 | const TMatrix *MHMatrix::InvertPosDef()
|
|---|
| 344 | {
|
|---|
| 345 | TMatrix m(fM);
|
|---|
| 346 |
|
|---|
| 347 | const Int_t rows = m.GetNrows();
|
|---|
| 348 | const Int_t cols = m.GetNcols();
|
|---|
| 349 |
|
|---|
| 350 | for (int x=0; x<cols; x++)
|
|---|
| 351 | {
|
|---|
| 352 | Double_t avg = 0;
|
|---|
| 353 | for (int y=0; y<rows; y++)
|
|---|
| 354 | avg += fM(y, x);
|
|---|
| 355 |
|
|---|
| 356 | avg /= rows;
|
|---|
| 357 |
|
|---|
| 358 | for (int y=0; y<rows; y++)
|
|---|
| 359 | m(y, x) -= avg;
|
|---|
| 360 | }
|
|---|
| 361 |
|
|---|
| 362 | TMatrix *m2 = new TMatrix(m, TMatrix::kTransposeMult, m);
|
|---|
| 363 |
|
|---|
| 364 | Double_t det;
|
|---|
| 365 | m2->Invert(&det);
|
|---|
| 366 | if (det==0)
|
|---|
| 367 | {
|
|---|
| 368 | *fLog << err << "ERROR - MHMatrix::InvertPosDef failed (Matrix is singular)." << endl;
|
|---|
| 369 | delete m2;
|
|---|
| 370 | return NULL;
|
|---|
| 371 | }
|
|---|
| 372 |
|
|---|
| 373 | // m2->Print();
|
|---|
| 374 |
|
|---|
| 375 | return m2;
|
|---|
| 376 | }
|
|---|
| 377 |
|
|---|
| 378 | // --------------------------------------------------------------------------
|
|---|
| 379 | //
|
|---|
| 380 | // Calculated the distance of vector evt from the reference sample
|
|---|
| 381 | // represented by the covariance metrix m.
|
|---|
| 382 | // - If n<0 the kernel method is applied and
|
|---|
| 383 | // -log(sum(epx(-d/h))/n) is returned.
|
|---|
| 384 | // - For n>0 the n nearest neighbors are summed and
|
|---|
| 385 | // sqrt(sum(d)/n) is returned.
|
|---|
| 386 | // - if n==0 all distances are summed
|
|---|
| 387 | //
|
|---|
| 388 | Double_t MHMatrix::CalcDist(const TMatrix &m, const TVector &evt, Int_t num) const
|
|---|
| 389 | {
|
|---|
| 390 | if (num==0) // may later be used for another method
|
|---|
| 391 | {
|
|---|
| 392 | TVector d = evt;
|
|---|
| 393 | d *= m;
|
|---|
| 394 | return TMath::Sqrt(d*evt);
|
|---|
| 395 | }
|
|---|
| 396 |
|
|---|
| 397 | const Int_t rows = fM.GetNrows();
|
|---|
| 398 | const Int_t cols = fM.GetNcols();
|
|---|
| 399 |
|
|---|
| 400 | TArrayD dists(rows);
|
|---|
| 401 |
|
|---|
| 402 | //
|
|---|
| 403 | // Calculate: v^T * M * v
|
|---|
| 404 | //
|
|---|
| 405 | for (int i=0; i<rows; i++)
|
|---|
| 406 | {
|
|---|
| 407 | TVector col(cols);
|
|---|
| 408 | col = TMatrixRow(fM, i);
|
|---|
| 409 |
|
|---|
| 410 | TVector d = evt;
|
|---|
| 411 | d -= col;
|
|---|
| 412 |
|
|---|
| 413 | TVector d2 = d;
|
|---|
| 414 | d2 *= m;
|
|---|
| 415 |
|
|---|
| 416 | dists[i] = d2*d; // square of distance
|
|---|
| 417 |
|
|---|
| 418 | //
|
|---|
| 419 | // This corrects for numerical uncertanties in cases of very
|
|---|
| 420 | // small distances...
|
|---|
| 421 | //
|
|---|
| 422 | if (dists[i]<0)
|
|---|
| 423 | dists[i]=0;
|
|---|
| 424 | }
|
|---|
| 425 |
|
|---|
| 426 | TArrayI idx(rows);
|
|---|
| 427 | TMath::Sort(dists.GetSize(), dists.GetArray(), idx.GetArray(), kFALSE);
|
|---|
| 428 |
|
|---|
| 429 | Int_t from = 0;
|
|---|
| 430 | Int_t to = TMath::Abs(num)<rows ? TMath::Abs(num) : rows;
|
|---|
| 431 | //
|
|---|
| 432 | // This is a zero-suppression for the case a test- and trainings
|
|---|
| 433 | // sample is identical. This would result in an unwanted leading
|
|---|
| 434 | // zero in the array. To suppress also numerical uncertanties of
|
|---|
| 435 | // zero we cut at 1e-5. Due to Rudy this should be enough. If
|
|---|
| 436 | // you encounter problems we can also use (eg) 1e-25
|
|---|
| 437 | //
|
|---|
| 438 | if (dists[idx[0]]<1e-5)
|
|---|
| 439 | {
|
|---|
| 440 | from++;
|
|---|
| 441 | to ++;
|
|---|
| 442 | if (to>rows)
|
|---|
| 443 | to = rows;
|
|---|
| 444 | }
|
|---|
| 445 |
|
|---|
| 446 | if (num<0)
|
|---|
| 447 | {
|
|---|
| 448 | //
|
|---|
| 449 | // Kernel function sum (window size h set according to literature)
|
|---|
| 450 | //
|
|---|
| 451 | const Double_t h = TMath::Power(rows, -1./(cols+4));
|
|---|
| 452 | const Double_t hwin = h*h*2;
|
|---|
| 453 |
|
|---|
| 454 | Double_t res = 0;
|
|---|
| 455 | for (int i=from; i<to; i++)
|
|---|
| 456 | res += TMath::Exp(-dists[idx[i]]/hwin);
|
|---|
| 457 |
|
|---|
| 458 | return -TMath::Log(res/(to-from));
|
|---|
| 459 | }
|
|---|
| 460 | else
|
|---|
| 461 | {
|
|---|
| 462 | //
|
|---|
| 463 | // Nearest Neighbor sum
|
|---|
| 464 | //
|
|---|
| 465 | Double_t res = 0;
|
|---|
| 466 | for (int i=from; i<to; i++)
|
|---|
| 467 | res += dists[idx[i]];
|
|---|
| 468 |
|
|---|
| 469 | return TMath::Sqrt(res/(to-from));
|
|---|
| 470 | }
|
|---|
| 471 | }
|
|---|
| 472 |
|
|---|
| 473 | // --------------------------------------------------------------------------
|
|---|
| 474 | //
|
|---|
| 475 | // Calls calc dist. In the case of the first call the covariance matrix
|
|---|
| 476 | // fM2 is calculated.
|
|---|
| 477 | // - If n<0 it is divided by (nrows-1)/h while h is the kernel factor.
|
|---|
| 478 | //
|
|---|
| 479 | Double_t MHMatrix::CalcDist(const TVector &evt, Int_t num)
|
|---|
| 480 | {
|
|---|
| 481 | if (!fM2.IsValid())
|
|---|
| 482 | {
|
|---|
| 483 | const TMatrix *m = InvertPosDef();
|
|---|
| 484 | if (!m)
|
|---|
| 485 | return -1;
|
|---|
| 486 |
|
|---|
| 487 | fM2.ResizeTo(*m);
|
|---|
| 488 | fM2 = *m;
|
|---|
| 489 | fM2 *= fM.GetNrows()-1;
|
|---|
| 490 | delete m;
|
|---|
| 491 | }
|
|---|
| 492 |
|
|---|
| 493 | return CalcDist(fM2, evt, num);
|
|---|
| 494 | }
|
|---|
| 495 |
|
|---|
| 496 | void MHMatrix::Reassign()
|
|---|
| 497 | {
|
|---|
| 498 | TMatrix m = fM;
|
|---|
| 499 | fM.ResizeTo(1,1);
|
|---|
| 500 | fM.ResizeTo(m);
|
|---|
| 501 | fM = m;
|
|---|
| 502 | }
|
|---|
| 503 |
|
|---|
| 504 | // --------------------------------------------------------------------------
|
|---|
| 505 | //
|
|---|
| 506 | // Implementation of SavePrimitive. Used to write the call to a constructor
|
|---|
| 507 | // to a macro. In the original root implementation it is used to write
|
|---|
| 508 | // gui elements to a macro-file.
|
|---|
| 509 | //
|
|---|
| 510 | void MHMatrix::StreamPrimitive(ofstream &out) const
|
|---|
| 511 | {
|
|---|
| 512 | Bool_t data = fData && !TestBit(kIsOwner);
|
|---|
| 513 |
|
|---|
| 514 | if (data)
|
|---|
| 515 | {
|
|---|
| 516 | fData->SavePrimitive(out);
|
|---|
| 517 | out << endl;
|
|---|
| 518 | }
|
|---|
| 519 |
|
|---|
| 520 | out << " MHMatrix " << GetUniqueName();
|
|---|
| 521 |
|
|---|
| 522 | if (data || fName!=gsDefName || fTitle!=gsDefTitle)
|
|---|
| 523 | {
|
|---|
| 524 | out << "(";
|
|---|
| 525 | if (data)
|
|---|
| 526 | out << "&" << fData->GetUniqueName();
|
|---|
| 527 | if (fName!=gsDefName || fTitle!=gsDefTitle)
|
|---|
| 528 | {
|
|---|
| 529 | if (data)
|
|---|
| 530 | out << ", ";
|
|---|
| 531 | out << "\"" << fName << "\"";
|
|---|
| 532 | if (fTitle!=gsDefTitle)
|
|---|
| 533 | out << ", \"" << fTitle << "\"";
|
|---|
| 534 | }
|
|---|
| 535 | }
|
|---|
| 536 | out << ");" << endl;
|
|---|
| 537 |
|
|---|
| 538 | if (fData && TestBit(kIsOwner))
|
|---|
| 539 | for (int i=0; i<fData->GetNumEntries(); i++)
|
|---|
| 540 | out << " " << GetUniqueName() << ".AddColumn(\"" << (*fData)[i].GetRule() << "\");" << endl;
|
|---|
| 541 | }
|
|---|
| 542 |
|
|---|
| 543 | const TArrayI MHMatrix::GetIndexOfSortedColumn(Int_t ncol, Bool_t desc) const
|
|---|
| 544 | {
|
|---|
| 545 | TMatrixColumn col(fM, ncol);
|
|---|
| 546 |
|
|---|
| 547 | const Int_t n = fM.GetNrows();
|
|---|
| 548 |
|
|---|
| 549 | TArrayF array(n);
|
|---|
| 550 |
|
|---|
| 551 | for (int i=0; i<n; i++)
|
|---|
| 552 | array[i] = col(i);
|
|---|
| 553 |
|
|---|
| 554 | TArrayI idx(n);
|
|---|
| 555 | TMath::Sort(n, array.GetArray(), idx.GetArray(), desc);
|
|---|
| 556 |
|
|---|
| 557 | return idx;
|
|---|
| 558 | }
|
|---|
| 559 |
|
|---|
| 560 | void MHMatrix::SortMatrixByColumn(Int_t ncol, Bool_t desc)
|
|---|
| 561 | {
|
|---|
| 562 | TArrayI idx = GetIndexOfSortedColumn(ncol, desc);
|
|---|
| 563 |
|
|---|
| 564 | const Int_t n = fM.GetNrows();
|
|---|
| 565 |
|
|---|
| 566 | TMatrix m(n, fM.GetNcols());
|
|---|
| 567 | for (int i=0; i<n; i++)
|
|---|
| 568 | {
|
|---|
| 569 | TVector vold(n);
|
|---|
| 570 | vold = TMatrixRow(fM, idx[i]);
|
|---|
| 571 |
|
|---|
| 572 | TMatrixRow rownew(m, i);
|
|---|
| 573 | rownew = vold;
|
|---|
| 574 | }
|
|---|
| 575 |
|
|---|
| 576 | fM = m;
|
|---|
| 577 | }
|
|---|
| 578 |
|
|---|
| 579 | Bool_t MHMatrix::Fill(MParList *plist, MTask *read)
|
|---|
| 580 | {
|
|---|
| 581 | //
|
|---|
| 582 | // Read data into Matrix
|
|---|
| 583 | //
|
|---|
| 584 | const Bool_t is = plist->IsOwner();
|
|---|
| 585 | plist->SetOwner(kFALSE);
|
|---|
| 586 |
|
|---|
| 587 | MTaskList tlist;
|
|---|
| 588 | plist->Replace(&tlist);
|
|---|
| 589 |
|
|---|
| 590 | MFillH fillh(this);
|
|---|
| 591 |
|
|---|
| 592 | tlist.AddToList(read);
|
|---|
| 593 | tlist.AddToList(&fillh);
|
|---|
| 594 |
|
|---|
| 595 | MEvtLoop evtloop;
|
|---|
| 596 | evtloop.SetParList(plist);
|
|---|
| 597 |
|
|---|
| 598 | if (!evtloop.Eventloop())
|
|---|
| 599 | return kFALSE;
|
|---|
| 600 |
|
|---|
| 601 | plist->Remove(&tlist);
|
|---|
| 602 | plist->SetOwner(is);
|
|---|
| 603 |
|
|---|
| 604 | return kTRUE;
|
|---|
| 605 | }
|
|---|
| 606 |
|
|---|
| 607 | // --------------------------------------------------------------------------
|
|---|
| 608 | //
|
|---|
| 609 | // Return a comma seperated list of all data members used in the matrix.
|
|---|
| 610 | // This is mainly used in MTask::AddToBranchList
|
|---|
| 611 | //
|
|---|
| 612 | TString MHMatrix::GetDataMember() const
|
|---|
| 613 | {
|
|---|
| 614 | return fData ? fData->GetDataMember() : TString("");
|
|---|
| 615 | }
|
|---|
| 616 |
|
|---|
| 617 | // --------------------------------------------------------------------------
|
|---|
| 618 | //
|
|---|
| 619 | // Return a comma seperated list of all data members used in the matrix.
|
|---|
| 620 | // This is mainly used in MTask::AddToBranchList
|
|---|
| 621 | //
|
|---|
| 622 | void MHMatrix::ReduceNumberOfRows(UInt_t numrows, const TString opt)
|
|---|
| 623 | {
|
|---|
| 624 | UInt_t rows = fM.GetNrows();
|
|---|
| 625 |
|
|---|
| 626 | if (rows==numrows)
|
|---|
| 627 | {
|
|---|
| 628 | *fLog << warn << "Matrix has already the correct number of rows..." << endl;
|
|---|
| 629 | return;
|
|---|
| 630 | }
|
|---|
| 631 |
|
|---|
| 632 | Float_t ratio = (Float_t)numrows/fM.GetNrows();
|
|---|
| 633 |
|
|---|
| 634 | if (ratio>=1)
|
|---|
| 635 | {
|
|---|
| 636 | *fLog << warn << "Matrix cannot be enlarged..." << endl;
|
|---|
| 637 | return;
|
|---|
| 638 | }
|
|---|
| 639 |
|
|---|
| 640 | Double_t sum = 0;
|
|---|
| 641 |
|
|---|
| 642 | UInt_t oldrow = 0;
|
|---|
| 643 | UInt_t newrow = 0;
|
|---|
| 644 |
|
|---|
| 645 | while (oldrow<rows)
|
|---|
| 646 | {
|
|---|
| 647 | sum += ratio;
|
|---|
| 648 |
|
|---|
| 649 | if (newrow<=(unsigned int)sum)
|
|---|
| 650 | {
|
|---|
| 651 | TVector vold(fM.GetNcols());
|
|---|
| 652 | vold = TMatrixRow(fM, oldrow);
|
|---|
| 653 |
|
|---|
| 654 | TMatrixRow rownew(fM, newrow);
|
|---|
| 655 | rownew = vold;
|
|---|
| 656 | newrow++;
|
|---|
| 657 | }
|
|---|
| 658 |
|
|---|
| 659 | oldrow++;
|
|---|
| 660 | }
|
|---|
| 661 |
|
|---|
| 662 | /*
|
|---|
| 663 |
|
|---|
| 664 | TMatrix m(n, fM.GetNcols());
|
|---|
| 665 | for (int i=0; i<n; i++)
|
|---|
| 666 | {
|
|---|
| 667 | TVector vold(n);
|
|---|
| 668 | vold = TMatrixRow(fM, idx[i]);
|
|---|
| 669 |
|
|---|
| 670 | TMatrixRow rownew(m, i);
|
|---|
| 671 | rownew = vold;
|
|---|
| 672 | }
|
|---|
| 673 | */
|
|---|
| 674 |
|
|---|
| 675 | }
|
|---|