Fork me on GitHub

Ignore:
File:
1 edited

Legend:

Unmodified
Added
Removed
  • examples/ExamplePVtxFind.C

    rb750b0a r46d3442  
    11/*
    22Example of using vertex fitter class to fit primary vertex
    3 assumed to be generated with Pythia8/ee_zh_smr-shf.cmn
     3assumed to be generated in (0,0,0)
    44*/
    55
     
    4444        TH1D* hTrFoundG = new TH1D("hTrFoundG", "Found GOOD primary tracks", 41, -0.5, 40.5);
    4545        TH1D* hTrFoundB = new TH1D("hTrFoundB", "Found BAD  primary tracks", 41, -0.5, 40.5);
    46         TH1D* hTrDiff = new TH1D("hTrDiff", "Found - available tracks", 21, -10.5, 10.5);
    4746        //
    4847        // Loop over all events
     
    8887                                pr[it] = new TVectorD(obsPar);
    8988                                cv[it] = new TMatrixDSym(trk->CovarianceMatrix());
    90                                 // Check covariance
    91                                 /*
    92                                 TMatrixDSymEigen Eign(*cv[it]);
    93                                 TVectorD lambda = Eign.GetEigenValues();
    94                                 Bool_t good = kTRUE;
    95                                 for(Int_t i=0; i<5; i++) if(lambda(i)<0.0) good = kFALSE;
    96                                 if(!good){
    97                                         std::cout<<"Cov["<<it<<"]"<<" eigenvalues: "
    98                                         <<lambda(0)<<", "<<lambda(1)<<", "<<lambda(2)<<", "
    99                                         <<lambda(3)<<", "<<lambda(4)<<std::endl;
    100                                         TMatrixDSym Ctmp = *cv[it];
    101                                         std::cout<<"Error D (dir - cov) "<<trk->ErrorD0
    102                                         <<" - "<<TMath::Sqrt(Ctmp(0,0))<<std::endl;
    103                                 }
    104                                 */
    10589                                //std::cout << "Track loaded it= " << it << std::endl;
    10690                        //
     
    136120                        }               // End loop on tracks
    137121                }
    138                 std::cout<<"PVtxFind true vertex: Nprim= "<<Nprim<<", x,y,z= "<<xpv<<", "<<ypv<<", "<<zpv<<std::endl;
    139122                //
    140123                // Find primary vertex
     
    149132                covpvc(1, 1) = 2.55e-05 * 2.55e-05;
    150133                covpvc(2, 2) = 0.64 * 0.64;
    151                 if(Nprim == 0){
    152                         xpv = xpvc(0);
    153                         ypv = xpvc(1);
    154                         zpv = xpvc(2);
    155                 }
    156134                //
    157135                //
     
    159137                Int_t nSkim = 0;
    160138                Int_t* nSkimmed = new Int_t[NtrG];
    161                 TVectorD** PrSk = new TVectorD * [1];
    162                 TMatrixDSym** CvSk = new TMatrixDSym * [1];
    163                 Double_t MaxChi2 = 2.0+3*2.;
     139                VertexFit* Vskim = new VertexFit();
     140                Vskim->AddVtxConstraint(xpvc, covpvc);
     141                Double_t MaxChi2 = 2.0+3*sqrt(2.);
    164142                for (Int_t n = 0; n < NtrG; n++) {
    165                         PrSk[0] = new TVectorD(*pr[n]);
    166                         CvSk[0] = new TMatrixDSym(*cv[n]);
    167                         VertexFit* Vskim = new VertexFit(1,PrSk, CvSk);
    168                         Vskim->AddVtxConstraint(xpvc, covpvc);
     143                        Vskim->AddTrk(pr[n], cv[n]);
    169144                        Double_t Chi2One = Vskim->GetVtxChi2();
    170145                        //std::cout<<"Track "<<n<<", Chi2 = "<<Chi2One<<std::endl;
     
    174149                                nSkim++;
    175150                        }
    176                         // Cleanup
    177                         delete Vskim;
    178                 }
    179                 delete PrSk[0];
    180                 delete CvSk[0];
    181                 delete[] PrSk;
    182                 delete[] CvSk;
    183                 //
     151                        Vskim->RemoveTrk(0);
     152                }
    184153                // Load tracks for primary fit
    185154                Int_t MinTrk = 1;       // Minumum # tracks for vertex fit
     155                TVectorD** PrFit = new TVectorD * [nSkim];
     156                TMatrixDSym** CvFit = new TMatrixDSym * [nSkim];
    186157                if (nSkim >= MinTrk) {
    187                         TVectorD** PrFit = new TVectorD * [nSkim];
    188                         TMatrixDSym** CvFit = new TMatrixDSym * [nSkim];
    189158                        for (Int_t n = 0; n < nSkim; n++) {
    190159                                PrFit[n] = new TVectorD(*pr[nSkimmed[n]]);
    191160                                CvFit[n] = new TMatrixDSym(*cv[nSkimmed[n]]);
    192                         }
    193                         delete[] nSkimmed;
    194                         Int_t Nfound = nSkim;
    195                         if(entry%500 ==0)std::cout << "Found tracks "<<Nfound << std::endl;
    196                         const Int_t MaxFound = 100; Double_t Chi2LL[MaxFound]; Double_t *Chi2L = Chi2LL;
     161                                TVectorD test = *PrFit[n];
     162                                //std::cout << "Test *PrFit[" << n << "](0) = " << test(0) << std::endl;
     163                        }
     164                }
     165                delete[] nSkimmed;
     166                delete Vskim;
     167                Int_t Nfound = nSkim;
     168                if(entry%500 ==0)std::cout << "Found tracks "<<Nfound << std::endl;
     169                if (nSkim >= MinTrk) {
    197170                        VertexFit* Vtx = new VertexFit(nSkim, PrFit, CvFit);
    198171                        //std::cout << "Vertex fit created " << std::endl;
     
    200173                        //
    201174                        // Remove tracks with large chi2
    202                         Double_t MaxChi2Fit = 4.5;
     175                        Double_t MaxChi2Fit = 9.;
    203176                        Bool_t Done = kFALSE;
    204177                        while (!Done) {
     
    207180                                TVectorD Chi2List = Vtx->GetVtxChi2List();      // Get contributions to Chi2
    208181                                //std::cout << "After Chi2List.  " << std::endl; Chi2List.Print();
    209                                 //Double_t* Chi2L = new Double_t[Nfound];
     182                                Double_t* Chi2L = new Double_t[Nfound];
    210183                                Chi2L = Chi2List.GetMatrixArray();
    211184                                Int_t iMax = TMath::LocMax(Nfound, Chi2L);
     
    213186                                Double_t Chi2Mx = Chi2L[iMax];
    214187                                //std::cout << "Chi2Mx "<<Chi2Mx << std::endl;
    215                                 if (Chi2Mx > MaxChi2Fit && Nfound > 1) {
     188                                if (Chi2Mx > MaxChi2Fit && Nfound > 2) {
    216189                                        //std::cout << "Before remove.  Nfound = "<<Nfound << std::endl;
    217190                                        Vtx->RemoveTrk(iMax);
     
    222195                                        Done = kTRUE;
    223196                                }
     197                                //std::cout << "Done =  " << Done << std::endl;
     198                                //delete[] Chi2L;
     199                                //std::cout << "Array Chi2L removed " << std::endl;
    224200                        }
    225201                        //
     
    227203                        //
    228204                        // Require minimum number of tracks in vertex
    229                         Int_t Nmin = 1;
     205                        Int_t Nmin = 2;
    230206                        if (Nfound >= Nmin) {
    231207                                TVectorD xvtx = Vtx->GetVtx();
    232                                 std::cout << "Found vertex " << xvtx(0)<<", "<<xvtx(1)<<", "<<xvtx(2) << std::endl;
     208                                //std::cout << "Found vertex " << xvtx(0)<<", "<<xvtx(1)<<", "<<xvtx(2) << std::endl;
    233209                                TMatrixDSym covX = Vtx->GetVtxCov();
    234210                                Double_t Chi2 = Vtx->GetVtxChi2();
     
    246222                                hTrPrim->Fill(Nprim);
    247223                                hTrFound->Fill((Double_t)Nfound);
    248                                 hTrDiff->Fill((Double_t)Nfound-Nprim);
    249224                                //std::cout << "Histograms filled " << std::endl;
    250225                        }
    251226                        //
    252                         // Clean
    253227                        delete Vtx;
    254                         for (Int_t i = 0; i < nSkim; i++) delete PrFit[i];
    255                         for (Int_t i = 0; i < nSkim; i++) delete CvFit[i];
    256                         delete[] PrFit;
    257                         delete[] CvFit;
    258228                }
    259229
     
    263233                for (Int_t i = 0; i < NtrG; i++) delete pr[i];
    264234                for (Int_t i = 0; i < NtrG; i++) delete cv[i];
     235                for (Int_t i = 0; i < nSkim; i++) delete PrFit[i];
     236                for (Int_t i = 0; i < nSkim; i++) delete CvFit[i];
    265237                delete[] pr;
    266238                delete[] cv;
     239                delete[] PrFit;
     240                delete[] CvFit;
    267241        }
    268242        //
     
    280254        //
    281255        TCanvas* CnvN = new TCanvas("CnvN", "Primary tracks found", 100, 100, 900, 500);
    282         CnvN->Divide(2, 2);
     256        CnvN->Divide(2, 1);
    283257        CnvN->cd(1);
    284258        hTrPrim->Draw();
     
    286260        hTrFound->SetLineColor(kRed);
    287261        hTrFound->Draw();
    288         CnvN->cd(3);
    289         hTrDiff->Draw();
    290262}
Note: See TracChangeset for help on using the changeset viewer.