/* ======================================================================== *\ ! ! * ! * This file is part of MARS, the MAGIC Analysis and Reconstruction ! * Software. It is distributed to you in the hope that it can be a useful ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes. ! * It is distributed WITHOUT ANY WARRANTY. ! * ! * Permission to use, copy, modify and distribute this software and its ! * documentation for any purpose is hereby granted without fee, ! * provided that the above copyright notice appear in all copies and ! * that both that copyright notice and this permission notice appear ! * in supporting documentation. It is provided "as is" without express ! * or implied warranty. ! * ! ! ! Author(s): Thomas Bretz 2002 ! Rudy Boeck 2003 ! ! Copyright: MAGIC Software Development, 2000-2003 ! ! \* ======================================================================== */ ///////////////////////////////////////////////////////////////////////////// // // MHMatrix // // This is a histogram container which holds a matrix with one column per // data variable. The data variable can be a complex rule (MDataChain). // Each event for wich Fill is called (by MFillH) is added as a new // row to the matrix. // // For example: // MHMatrix m; // m.AddColumn("MHillas.fSize"); // m.AddColumn("MMcEvt.fImpact/100"); // m.AddColumn("HillasSource.fDist*MGeomCam.fConvMm2Deg"); // MFillH fillm(&m); // taskliost.AddToList(&fillm); // [...] // m.Print(); // ///////////////////////////////////////////////////////////////////////////// #include "MHMatrix.h" #include #include #include #include #include #include #include #include #include "MLog.h" #include "MLogManip.h" #include "MFillH.h" #include "MEvtLoop.h" #include "MParList.h" #include "MTaskList.h" #include "MProgressBar.h" #include "MData.h" #include "MDataArray.h" #include "MFilter.h" ClassImp(MHMatrix); using namespace std; const TString MHMatrix::gsDefName = "MHMatrix"; const TString MHMatrix::gsDefTitle = "Multidimensional Matrix"; // -------------------------------------------------------------------------- // // Default Constructor // MHMatrix::MHMatrix(const char *name, const char *title) : fNumRows(0), fData(NULL) { fName = name ? name : gsDefName.Data(); fTitle = title ? title : gsDefTitle.Data(); } // -------------------------------------------------------------------------- // // Default Constructor // MHMatrix::MHMatrix(const TMatrix &m, const char *name, const char *title) : fNumRows(m.GetNrows()), fM(m), fData(NULL) { fName = name ? name : gsDefName.Data(); fTitle = title ? title : gsDefTitle.Data(); } // -------------------------------------------------------------------------- // // Constructor. Initializes the columns of the matrix with the entries // from a MDataArray // MHMatrix::MHMatrix(MDataArray *mat, const char *name, const char *title) : fNumRows(0), fData(mat) { fName = name ? name : gsDefName.Data(); fTitle = title ? title : gsDefTitle.Data(); } // -------------------------------------------------------------------------- // // Destructor. Does not deleted a user given MDataArray, except IsOwner // was called. // MHMatrix::~MHMatrix() { if (TestBit(kIsOwner) && fData) delete fData; } // -------------------------------------------------------------------------- // Bool_t MHMatrix::SetNumRow(Int_t row) { if (row>=fM.GetNrows() || row<0) return kFALSE; fRow = row; return kTRUE; } // -------------------------------------------------------------------------- // // Add a new column to the matrix. This can only be done before the first // event (row) was filled into the matrix. For the syntax of the rule // see MDataChain. // Returns the index of the new column, -1 in case of failure. // (0, 1, 2, ... for the 1st, 2nd, 3rd, ...) // Int_t MHMatrix::AddColumn(const char *rule) { if (fM.IsValid()) { *fLog << warn << "Warning - matrix is already in use. Can't add a new column... skipped." << endl; return -1; } if (TestBit(kIsLocked)) { *fLog << warn << "Warning - matrix is locked. Can't add new column... skipped." << endl; return -1; } if (!fData) { fData = new MDataArray; SetBit(kIsOwner); } fData->AddEntry(rule); return fData->GetNumEntries()-1; } // -------------------------------------------------------------------------- // void MHMatrix::AddColumns(MDataArray *matrix) { if (fM.IsValid()) { *fLog << warn << "Warning - matrix is already in use. Can't add new columns... skipped." << endl; return; } if (TestBit(kIsLocked)) { *fLog << warn << "Warning - matrix is locked. Can't add new columns... skipped." << endl; return; } if (fData) *fLog << warn << "Warning - columns already added... replacing." << endl; if (fData && TestBit(kIsOwner)) { delete fData; ResetBit(kIsOwner); } fData = matrix; } // -------------------------------------------------------------------------- // // Checks whether at least one column is available and PreProcesses all // data chains. // Bool_t MHMatrix::SetupFill(const MParList *plist) { if (!fData) { *fLog << err << "Error - No Columns initialized... aborting." << endl; return kFALSE; } return fData->PreProcess(plist); } // -------------------------------------------------------------------------- // // If the matrix has not enough rows double the number of available rows. // void MHMatrix::AddRow() { fNumRows++; if (fM.GetNrows() > fNumRows) return; if (!fM.IsValid()) { fM.ResizeTo(1, fData->GetNumEntries()); return; } #if ROOT_VERSION_CODE < ROOT_VERSION(3,05,07) TMatrix m(fM); #endif fM.ResizeTo(fM.GetNrows()*2, fData->GetNumEntries()); #if ROOT_VERSION_CODE < ROOT_VERSION(3,05,07) TVector vold(fM.GetNcols()); for (int x=0; xGetNumEntries(); col++) fM(fNumRows-1, col) = (*fData)(col); return kTRUE; } // -------------------------------------------------------------------------- // // Resize the matrix to a number of rows which corresponds to the number of // rows which have really been filled with values. // Bool_t MHMatrix::Finalize() { // // It's not a fatal error so we don't need to stop PostProcessing... // if (fData->GetNumEntries()==0 || fNumRows<1) return kTRUE; if (fNumRows != fM.GetNrows()) { TMatrix m(fM); CopyCrop(fM, m, fNumRows); } return kTRUE; } /* // -------------------------------------------------------------------------- // // Draw clone of histogram. So that the object can be deleted // and the histogram is still visible in the canvas. // The cloned object are deleted together with the canvas if the canvas is // destroyed. If you want to handle destroying the canvas you can get a // pointer to it from this function // TObject *MHMatrix::DrawClone(Option_t *opt) const { TCanvas &c = *MH::MakeDefCanvas(fHist); // // This is necessary to get the expected bahviour of DrawClone // gROOT->SetSelectedPad(NULL); fHist->DrawCopy(opt); TString str(opt); if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2) { TProfile *p = ((TH2*)fHist)->ProfileX(); p->Draw("same"); p->SetBit(kCanDelete); } if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2) { TProfile *p = ((TH2*)fHist)->ProfileY(); p->Draw("same"); p->SetBit(kCanDelete); } c.Modified(); c.Update(); return &c; } // -------------------------------------------------------------------------- // // Creates a new canvas and draws the histogram into it. // Be careful: The histogram belongs to this object and won't get deleted // together with the canvas. // void MHMatrix::Draw(Option_t *opt) { if (!gPad) MH::MakeDefCanvas(fHist); fHist->Draw(opt); TString str(opt); if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2) { TProfile *p = ((TH2*)fHist)->ProfileX(); p->Draw("same"); p->SetBit(kCanDelete); } if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2) { TProfile *p = ((TH2*)fHist)->ProfileY(); p->Draw("same"); p->SetBit(kCanDelete); } gPad->Modified(); gPad->Update(); } */ // -------------------------------------------------------------------------- // // Prints the meaning of the columns and the contents of the matrix. // Becareful, this can take a long time for matrices with many rows. // Use the option 'size' to print the size of the matrix. // Use the option 'cols' to print the culumns // Use the option 'data' to print the contents // void MHMatrix::Print(Option_t *o) const { TString str(o); *fLog << all << flush; if (str.Contains("size", TString::kIgnoreCase)) { *fLog << GetDescriptor() << ": NumColumns=" << fM.GetNcols(); *fLog << " NumRows=" << fM.GetNrows() << endl; } if (!fData && str.Contains("cols", TString::kIgnoreCase)) *fLog << "Sorry, no column information available." << endl; if (fData && str.Contains("cols", TString::kIgnoreCase)) fData->Print(); if (str.Contains("data", TString::kIgnoreCase)) fM.Print(); } // -------------------------------------------------------------------------- // const TMatrix *MHMatrix::InvertPosDef() { TMatrix m(fM); const Int_t rows = m.GetNrows(); const Int_t cols = m.GetNcols(); for (int x=0; xInvert(&det); if (det==0) { *fLog << err << "ERROR - MHMatrix::InvertPosDef failed (Matrix is singular)." << endl; delete m2; return NULL; } // m2->Print(); return m2; } // -------------------------------------------------------------------------- // // Calculated the distance of vector evt from the reference sample // represented by the covariance metrix m. // - If n<0 the kernel method is applied and // -log(sum(epx(-d/h))/n) is returned. // - For n>0 the n nearest neighbors are summed and // sqrt(sum(d)/n) is returned. // - if n==0 all distances are summed // Double_t MHMatrix::CalcDist(const TMatrix &m, const TVector &evt, Int_t num) const { if (num==0) // may later be used for another method { TVector d = evt; d *= m; return TMath::Sqrt(d*evt); } const Int_t rows = fM.GetNrows(); const Int_t cols = fM.GetNcols(); TArrayD dists(rows); // // Calculate: v^T * M * v // for (int i=0; irows) to = rows; } if (num<0) { // // Kernel function sum (window size h set according to literature) // const Double_t h = TMath::Power(rows, -1./(cols+4)); const Double_t hwin = h*h*2; Double_t res = 0; for (int i=from; iSavePrimitive(out); out << endl; } out << " MHMatrix " << GetUniqueName(); if (data || fName!=gsDefName || fTitle!=gsDefTitle) { out << "("; if (data) out << "&" << fData->GetUniqueName(); if (fName!=gsDefName || fTitle!=gsDefTitle) { if (data) out << ", "; out << "\"" << fName << "\""; if (fTitle!=gsDefTitle) out << ", \"" << fTitle << "\""; } } out << ");" << endl; if (fData && TestBit(kIsOwner)) for (int i=0; iGetNumEntries(); i++) out << " " << GetUniqueName() << ".AddColumn(\"" << (*fData)[i].GetRule() << "\");" << endl; } // -------------------------------------------------------------------------- // const TArrayI MHMatrix::GetIndexOfSortedColumn(Int_t ncol, Bool_t desc) const { TMatrixColumn col(fM, ncol); const Int_t n = fM.GetNrows(); TArrayF array(n); for (int i=0; iIsOwner(); plist->SetOwner(kFALSE); MTaskList tlist; plist->Replace(&tlist); MFillH fillh(this); tlist.AddToList(read); if (filter) { tlist.AddToList(filter); fillh.SetFilter(filter); } tlist.AddToList(&fillh); //MProgressBar bar; MEvtLoop evtloop("MHMatrix::Fill-EvtLoop"); evtloop.SetParList(plist); //evtloop.SetProgressBar(&bar); if (!evtloop.Eventloop()) return kFALSE; tlist.PrintStatistics(); plist->Remove(&tlist); plist->SetOwner(is); return kTRUE; } // -------------------------------------------------------------------------- // // Return a comma seperated list of all data members used in the matrix. // This is mainly used in MTask::AddToBranchList // TString MHMatrix::GetDataMember() const { return fData ? fData->GetDataMember() : TString(""); } // -------------------------------------------------------------------------- // // void MHMatrix::ReduceNumberOfRows(UInt_t numrows, const TString opt) { UInt_t rows = fM.GetNrows(); if (rows==numrows) { *fLog << warn << "Matrix has already the correct number of rows..." << endl; return; } Float_t ratio = (Float_t)numrows/fM.GetNrows(); if (ratio>=1) { *fLog << warn << "Matrix cannot be enlarged..." << endl; return; } Double_t sum = 0; UInt_t oldrow = 0; UInt_t newrow = 0; TVector vold(fM.GetNcols()); while (oldrowDivide(2,2); th1->cd(1); ((TH1&)hth).DrawCopy(); // real histogram before th1->cd(2); ((TH1&)hta).DrawCopy(); // histogram after th1->cd(3); ((TH1&)hthd).DrawCopy(); // correction factors th1->cd(4); ((TH1&)thsh).DrawCopy(); // target } // ------------------------------------------------------------------------ // // Resizes th etarget matrix to rows*source.GetNcol() and copies // the data from the first (n)rows or the source into the target matrix. // void MHMatrix::CopyCrop(TMatrix &target, const TMatrix &source, Int_t rows) { TVector v(source.GetNcols()); target.ResizeTo(rows, source.GetNcols()); for (Int_t ir=0; irfM.GetNrows()) { *fLog << warn << dbginf << "No.requested (" << nmaxevts; *fLog << ") > available events (" << fM.GetNrows() << ")... "; *fLog << "setting equal." << endl; nmaxevts = fM.GetNrows(); } if (nmaxevts<0) { *fLog << err << dbginf << "Number of requested events < 0" << endl; return kFALSE; } if (nmaxevts==0) nmaxevts = fM.GetNrows(); // // refcol is the column number starting at 0; it is >= 0 // // number of the column (count from 0) containing // the variable for which the target distribution is given // // // Calculate normalization factors // //const int nbins = thsh.GetNbinsX(); //const double frombin = thsh.GetBinLowEdge(1); //const double tobin = thsh.GetBinLowEdge(nbins+1); //const double dbin = thsh.GetBinWidth(1); const Int_t nrows = fM.GetNrows(); const Int_t ncols = fM.GetNcols(); // // set up the real histogram (distribution before) // //TH1F hth("th", "Distribution before reduction", nbins, frombin, tobin); TH1F hth; hth.SetNameTitle("th", "Distribution before reduction"); SetBinning(&hth, &thsh); hth.SetDirectory(NULL); for (Int_t j=0; j= nmaxevts) break; TMatrixRow(mnewtmp, evtcount1) = vold = TMatrixRow(fM, ind[ir]); } for (/*empty*/; irfM.GetNrows()) { *fLog << dbginf << "No.of requested events (" << nmaxevts << ") exceeds no.of available events (" << fM.GetNrows() << ")" << endl; *fLog << dbginf << " set no.of requested events = no.of available events" << endl; nmaxevts = fM.GetNrows(); } if (nmaxevts<0) { *fLog << err << dbginf << "Number of requested events < 0" << endl; return kFALSE; } if (nmaxevts==0) nmaxevts = fM.GetNrows(); const Int_t nrows = fM.GetNrows(); const Int_t ncols = fM.GetNcols(); // // get random access // TArrayI ind(nrows); GetRandomArrayI(ind); // // define new matrix // Int_t evtcount1 = 0; Int_t evtcount2 = 0; TMatrix mnewtmp(nrows, ncols); TMatrix mrest(nrows, ncols); // // select events (distribution after renormalization) // TVector vold(fM.GetNcols()); for (Int_t ir=0; irDelete(); } else { fData = new MDataArray; SetBit(kIsOwner); } fData->AddEntry(name); i++; } return i!=0; } // -------------------------------------------------------------------------- // // ShuffleEvents. Shuffles the order of the matrix rows. // // void MHMatrix::ShuffleRows(UInt_t seed) { TRandom rnd(seed); TVector v(fM.GetNcols()); TVector tmp(fM.GetNcols()); for (Int_t irow = 0; irow=fM.GetNrows()) { *fLog << warn << GetDescriptor() << ": Warning - " << num << " >= rows=" << fM.GetNrows() << endl; return; } #if ROOT_VERSION_CODE < ROOT_VERSION(3,05,07) const TMatrix m(fM); #endif fM.ResizeTo(num, fM.GetNcols()); #if ROOT_VERSION_CODE < ROOT_VERSION(3,05,07) TVector tmp(fM.GetNcols()); for (UInt_t irow=0; irow math.h) checks for NaN as well as inf int jcol; for (jcol=0; jcol