Fork me on GitHub

Changeset 7dac4ea in git for examples/ExamplePVtxFind.C


Ignore:
Timestamp:
May 17, 2021, 6:05:45 PM (3 years ago)
Author:
GitHub <noreply@…>
Branches:
master
Children:
1c1c9c2, 33a6b3a, f8e61b2
Parents:
4afb18d (diff), 4acf2fd (diff)
Note: this is a merge changeset, the changes displayed below correspond to the merge itself.
Use the (diff) links above to see all the changes relative to each parent.
git-author:
Michele Selvaggi <michele.selvaggi@…> (05/17/21 18:05:45)
git-committer:
GitHub <noreply@…> (05/17/21 18:05:45)
Message:

Merge pull request #95 from fbedesch/master

Fix for backward tracks and vertex fitting improvements

File:
1 edited

Legend:

Unmodified
Added
Removed
  • examples/ExamplePVtxFind.C

    r4afb18d r7dac4ea  
    5151                // Load selected branches with data from specified event
    5252                treeReader->ReadEntry(entry);
    53                 Int_t Ntr = 0;  // # of starting tracks used in primary vertex
    5453                Int_t NtrG = branchTrack->GetEntries();
    55                 //std::cout << "Event opened containing " << NtrG << " tracks" << std::endl;
     54                if(entry%500 ==0)std::cout << "Event "<<entry<<" opened containing " << NtrG << " tracks" << std::endl;
    5655                TVectorD** pr = new TVectorD * [NtrG];
    5756                TMatrixDSym** cv = new TMatrixDSym * [NtrG];
     
    6362                //
    6463                Double_t Nprim = 0.0;
     64                Double_t xpv = 0.0;             // Init true primary vertex position
     65                Double_t ypv = 0.0;
     66                Double_t zpv = 0.0;
    6567                if (branchTrack->GetEntries() > 0)
    6668                {
     
    8082                                //std::cout << "Got track parameters for track " << it << std::endl;
    8183                                //
    82                                 // Load tracks for vertex fit if impact parameters is not ridiculous
    83                                 Double_t Dmax = 1.0;    // max is 1 mm
    84                                 if (TMath::Abs(obsD0) < Dmax) {
    85                                         Double_t oPar[5] = { obsD0, obsPhi, obsC, obsZ0, obsCtg };
    86                                         TVectorD obsPar(5, oPar);       // Fill observed parameters
    87                                         pr[Ntr] = new TVectorD(obsPar);
    88                                         cv[Ntr] = new TMatrixDSym(trk->CovarianceMatrix());
    89                                         Ntr++;
    90                                         //std::cout << "Track loaded Ntr= " << Ntr << std::endl;
    91                                 }
    92                                 //
    93                                 // Get associated generated particle
     84                                // Load all tracks for vertex fit
     85                                Double_t oPar[5] = { obsD0, obsPhi, obsC, obsZ0, obsCtg };
     86                                TVectorD obsPar(5, oPar);       // Fill observed parameters
     87                                pr[it] = new TVectorD(obsPar);
     88                                cv[it] = new TMatrixDSym(trk->CovarianceMatrix());
     89                                //std::cout << "Track loaded it= " << it << std::endl;
     90                        //
     91                        // Find true primary vertex
    9492                                GenParticle* gp = (GenParticle*)trk->Particle.GetObject();
    9593                                //std::cout << "GenParticle pointer "<<gp << std::endl;
     
    9997                                Double_t y = gp->Y;
    10098                                Double_t z = gp->Z;
    101                                 //std::cout << "Got position of origin " << std::endl;
    102                                 //
    103                                 // Count tracks originating from the primary vertex
    104                                 if (x == 0.0 && y == 0.0) Nprim++;
     99                                Bool_t prim = kTRUE;    // Is primary?
     100                                Int_t mp = gp->M1;      // Mother
     101                                while (mp > 0) {
     102                                        GenParticle* gm =
     103                                                (GenParticle*)branchGenPart->At(mp);
     104                                        Double_t xm = gm->X;
     105                                        Double_t ym = gm->Y;
     106                                        Double_t zm = gm->Z;
     107                                        if (x != xm || y != ym || z != zm) {
     108                                                prim = kFALSE;
     109                                                break;
     110                                        }
     111                                        else mp = gm->M1;
     112                                }
     113                                if (prim) {             // It's a primary track
     114                                        Nprim++;
     115                                        xpv = x;        // Store true primary
     116                                        ypv = y;
     117                                        zpv = z;
     118                                }
    105119
    106120                        }               // End loop on tracks
    107121                }
    108122                //
    109                 // Fit primary vertex
    110                 //
    111                 Int_t Nfound = Ntr;
    112                 //std::cout << "Found tracks "<<Nfound << std::endl;
    113                 Int_t MinTrk = 2;       // Minumum # tracks for vertex fit
    114                 Double_t MaxChi2 = 6.;
    115                 if (Ntr >= MinTrk) {
    116                         VertexFit* Vtx = new VertexFit(Ntr, pr, cv);
     123                // Find primary vertex
     124                //
     125                //Beam constraint
     126                TVectorD xpvc(3);
     127                xpvc(0) = 1.0;
     128                xpvc(1) = -2.0;
     129                xpvc(2) = 10.0;
     130                TMatrixDSym covpvc(3); covpvc.Zero();
     131                covpvc(0, 0) = 0.0097 * 0.0097;
     132                covpvc(1, 1) = 2.55e-05 * 2.55e-05;
     133                covpvc(2, 2) = 0.64 * 0.64;
     134                //
     135                //
     136                // Skim tracks
     137                Int_t nSkim = 0;
     138                Int_t* nSkimmed = new Int_t[NtrG];
     139                VertexFit* Vskim = new VertexFit();
     140                Vskim->AddVtxConstraint(xpvc, covpvc);
     141                Double_t MaxChi2 = 2.0+3*sqrt(2.);
     142                for (Int_t n = 0; n < NtrG; n++) {
     143                        Vskim->AddTrk(pr[n], cv[n]);
     144                        Double_t Chi2One = Vskim->GetVtxChi2();
     145                        //std::cout<<"Track "<<n<<", Chi2 = "<<Chi2One<<std::endl;
     146                        if (Chi2One < MaxChi2) {
     147                                nSkimmed[nSkim] = n;
     148                                //std::cout << "nSkimmed[" << nSkim << "] = " << n << std::endl;
     149                                nSkim++;
     150                        }
     151                        Vskim->RemoveTrk(0);
     152                }
     153                // Load tracks for primary fit
     154                Int_t MinTrk = 1;       // Minumum # tracks for vertex fit
     155                TVectorD** PrFit = new TVectorD * [nSkim];
     156                TMatrixDSym** CvFit = new TMatrixDSym * [nSkim];
     157                if (nSkim >= MinTrk) {
     158                        for (Int_t n = 0; n < nSkim; n++) {
     159                                PrFit[n] = new TVectorD(*pr[nSkimmed[n]]);
     160                                CvFit[n] = new TMatrixDSym(*cv[nSkimmed[n]]);
     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) {
     170                        VertexFit* Vtx = new VertexFit(nSkim, PrFit, CvFit);
    117171                        //std::cout << "Vertex fit created " << std::endl;
     172                        Vtx->AddVtxConstraint(xpvc, covpvc);
    118173                        //
    119174                        // Remove tracks with large chi2
     175                        Double_t MaxChi2Fit = 9.;
    120176                        Bool_t Done = kFALSE;
    121                         while (!Done){
     177                        while (!Done) {
    122178                                //std::cout << "After while " << std::endl;
    123179                                // Find largest Chi2 contribution
     
    130186                                Double_t Chi2Mx = Chi2L[iMax];
    131187                                //std::cout << "Chi2Mx "<<Chi2Mx << std::endl;
    132                                 if (Chi2Mx > MaxChi2 && Nfound > 2) {
    133                                         Vtx->RemoveTrk(iMax);
     188                                if (Chi2Mx > MaxChi2Fit && Nfound > 2) {
     189                                        //std::cout << "Before remove.  Nfound = "<<Nfound << std::endl;
     190                                        Vtx->RemoveTrk(iMax);
     191                                        //std::cout << "After remove." << std::endl;
    134192                                        Nfound--;
    135193                                }
     
    145203                        //
    146204                        // Require minimum number of tracks in vertex
    147                         Int_t Nmin = 4;
     205                        Int_t Nmin = 2;
    148206                        if (Nfound >= Nmin) {
    149207                                TVectorD xvtx = Vtx->GetVtx();
     
    151209                                TMatrixDSym covX = Vtx->GetVtxCov();
    152210                                Double_t Chi2 = Vtx->GetVtxChi2();
    153                                 Double_t Ndof = 2 * (Double_t)Nfound - 3;
    154                                 Double_t PullX = xvtx(0) / TMath::Sqrt(covX(0, 0));
    155                                 Double_t PullY = xvtx(1) / TMath::Sqrt(covX(1, 1));
    156                                 Double_t PullZ = xvtx(2) / TMath::Sqrt(covX(2, 2));
     211                                Double_t Ndof = 2 * (Double_t)Nfound;
     212                                Double_t PullX = (xvtx(0) - xpv) / TMath::Sqrt(covX(0, 0));
     213                                Double_t PullY = (xvtx(1) - ypv) / TMath::Sqrt(covX(1, 1));
     214                                Double_t PullZ = (xvtx(2) - zpv) / TMath::Sqrt(covX(2, 2));
    157215                                //
    158216                                // Fill histograms
     
    173231                //
    174232                // Cleanup
    175                 for (Int_t i = 0; i < Ntr; i++) delete pr[i];
    176                 for (Int_t i = 0; i < Ntr; i++) delete cv[i];
     233                for (Int_t i = 0; i < NtrG; i++) delete pr[i];
     234                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];
    177237                delete[] pr;
    178238                delete[] cv;
     239                delete[] PrFit;
     240                delete[] CvFit;
    179241        }
    180242        //
     
    195257        CnvN->cd(1);
    196258        hTrPrim->Draw();
    197         CnvN-> cd(2);
     259        CnvN->cd(2);
    198260        hTrFound->SetLineColor(kRed);
    199261        hTrFound->Draw();
Note: See TracChangeset for help on using the changeset viewer.