CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/Validation/RecoVertex/src/PrimaryVertexAnalyzer.cc

Go to the documentation of this file.
00001 #include "Validation/RecoVertex/interface/PrimaryVertexAnalyzer.h"
00002 
00003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 #include "FWCore/Version/interface/GetReleaseVersion.h"
00006 
00007 // reco track and vertex 
00008 #include "DataFormats/TrackReco/interface/Track.h"
00009 #include "DataFormats/VertexReco/interface/Vertex.h"
00010 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
00011 
00012 // simulated vertices,..., add <use name=SimDataFormats/Vertex> and <../Track>
00013 #include <SimDataFormats/Vertex/interface/SimVertex.h>
00014 #include <SimDataFormats/Vertex/interface/SimVertexContainer.h>
00015 #include <SimDataFormats/Track/interface/SimTrack.h>
00016 #include <SimDataFormats/Track/interface/SimTrackContainer.h>
00017 
00018 //generator level + CLHEP
00019 #include "HepMC/GenEvent.h"
00020 #include "HepMC/GenVertex.h"
00021 //#include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
00022 
00023 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00024 
00025 // Root
00026 #include <TH1.h>
00027 #include <TH2.h>
00028 #include <TFile.h>
00029 #include <TProfile.h>
00030  
00031 #include <cmath>
00032 #include <gsl/gsl_math.h>
00033 #include <gsl/gsl_eigen.h>
00034 
00035 
00036 using namespace edm;
00037 using namespace reco;
00038 //
00039 // constants, enums and typedefs
00040 //
00041 
00042 //
00043 // static data member definitions
00044 //
00045 
00046 //
00047 // constructors and destructor
00048 //
00049 PrimaryVertexAnalyzer::PrimaryVertexAnalyzer(const ParameterSet& iConfig)
00050 {
00051    //now do what ever initialization is needed
00052   simG4_=iConfig.getParameter<edm::InputTag>( "simG4" );
00053   recoTrackProducer_= iConfig.getUntrackedParameter<std::string>("recoTrackProducer");
00054   // open output file to store histograms}
00055   outputFile_  = iConfig.getUntrackedParameter<std::string>("outputFile");
00056   TString tversion(edm::getReleaseVersion());
00057   tversion = tversion.Remove(0,1);
00058   tversion = tversion.Remove(tversion.Length()-1,tversion.Length());
00059   outputFile_  = std::string(tversion)+"_"+outputFile_;
00060 
00061   vtxSample_   = iConfig.getUntrackedParameter<std::vector< std::string > >("vtxSample");
00062   for(std::vector<std::string>::iterator isample = vtxSample_.begin(); isample!=vtxSample_.end(); ++isample) {
00063     if ( *isample == "offlinePrimaryVertices" ) suffixSample_.push_back("AVF");
00064     if ( *isample == "offlinePrimaryVerticesWithBS" ) suffixSample_.push_back("wBS");
00065   }
00066   if ( suffixSample_.size() == 0 ) throw cms::Exception("NoVertexSamples") << " no known vertex samples given";
00067 
00068   rootFile_ = new TFile(outputFile_.c_str(),"RECREATE");
00069   verbose_= iConfig.getUntrackedParameter<bool>("verbose", false);
00070   simUnit_= 1.0;  // starting with CMSSW_1_2_x ??
00071   if ( (edm::getReleaseVersion()).find("CMSSW_1_1_",0)!=std::string::npos){
00072     simUnit_=0.1;  // for use in  CMSSW_1_1_1 tutorial
00073   }
00074 }
00075 
00076 
00077 PrimaryVertexAnalyzer::~PrimaryVertexAnalyzer()
00078 {
00079     // do anything here that needs to be done at desctruction time
00080    // (e.g. close files, deallocate resources etc.)
00081    delete rootFile_;
00082    /*for(std::map<std::string,TH1*>::const_iterator hist=h.begin(); hist!=h.end(); hist++){
00083     delete hist->second;
00084     }*/
00085    
00086 }
00087 
00088 
00089 
00090 //
00091 // member functions
00092 //
00093 void PrimaryVertexAnalyzer::beginJob(){
00094   std::cout << " PrimaryVertexAnalyzer::beginJob  conversion from sim units to rec units is " << simUnit_ << std::endl;
00095 
00096 
00097   rootFile_->cd();
00098   // release validation histograms used in DoCompare.C
00099   for (std::vector<std::string>::iterator isuffix= suffixSample_.begin(); isuffix!=suffixSample_.end(); isuffix++) {
00100 
00101     hdir[*isuffix] = rootFile_->mkdir(TString(*isuffix));
00102     hdir[*isuffix]->cd();
00103     h["nbvtx"+ *isuffix]        = new TH1F(TString("nbvtx"+ *isuffix),"Reconstructed Vertices in Event",20,-0.5,19.5);
00104     h["nbtksinvtx"+ *isuffix]   = new TH1F(TString("nbtksinvtx"+ *isuffix),"Reconstructed Tracks in Vertex",200,-0.5,199.5); 
00105     h["resx"+ *isuffix]         = new TH1F(TString("resx"+ *isuffix),"Residual X",100,-0.04,0.04);
00106     h["resy"+ *isuffix]         = new TH1F(TString("resy"+ *isuffix),"Residual Y",100,-0.04,0.04);
00107     h["resz"+ *isuffix]         = new TH1F(TString("resz"+ *isuffix),"Residual Z",100,-0.1,0.1);
00108     h["pullx"+ *isuffix]        = new TH1F(TString("pullx"+ *isuffix),"Pull X",100,-25.,25.);
00109     h["pully"+ *isuffix]        = new TH1F(TString("pully"+ *isuffix),"Pull Y",100,-25.,25.);
00110     h["pullz"+ *isuffix]        = new TH1F(TString("pullz"+ *isuffix),"Pull Z",100,-25.,25.);
00111     h["vtxchi2"+ *isuffix]      = new TH1F(TString("vtxchi2"+ *isuffix),"#chi^{2}",100,0.,100.);
00112     h["vtxndf"+ *isuffix]       = new TH1F(TString("vtxndf"+ *isuffix),"ndof",100,0.,100.);
00113     h["tklinks"+ *isuffix]      = new TH1F(TString("tklinks"+ *isuffix),"Usable track links",2,-0.5,1.5);
00114     h["nans"+ *isuffix]         = new TH1F(TString("nans"+ *isuffix),"Illegal values for x,y,z,xx,xy,xz,yy,yz,zz",9,0.5,9.5);
00115     // more histograms
00116     h["vtxprob"+ *isuffix]      = new TH1F(TString("vtxprob"+ *isuffix),"#chi^{2} probability",100,0.,1.);
00117     h["eff"+ *isuffix]          = new TH1F(TString("eff"+ *isuffix),"efficiency",2, -0.5, 1.5);
00118     h["efftag"+ *isuffix]       = new TH1F(TString("efftag"+ *isuffix),"efficiency tagged vertex",2, -0.5, 1.5);
00119     h["effvseta"+ *isuffix]     = new TProfile(TString("effvseta"+ *isuffix),"efficiency vs eta",20, -2.5, 2.5, 0, 1.);
00120     h["effvsptsq"+ *isuffix]    = new TProfile(TString("effvsptsq"+ *isuffix),"efficiency vs ptsq",20, 0., 10000., 0, 1.);
00121     h["effvsntrk"+ *isuffix]    = new TProfile(TString("effvsntrk"+ *isuffix),"efficiency vs # simtracks",200, 0., 200., 0, 1.);
00122     h["effvsnrectrk"+ *isuffix] = new TProfile(TString("effvsnrectrk"+ *isuffix),"efficiency vs # rectracks",200, 0., 200., 0, 1.);
00123     h["effvsz"+ *isuffix]       = new TProfile(TString("effvsz"+ *isuffix),"efficiency vs z",40, -20., 20., 0, 1.);
00124     h["nbsimtksinvtx"+ *isuffix]= new TH1F(TString("nbsimtksinvtx"+ *isuffix),"simulated tracks in vertex",100,-0.5,99.5); 
00125     h["xrec"+ *isuffix]         = new TH1F(TString("xrec"+ *isuffix),"reconstructed x",100,-0.1,0.1);
00126     h["yrec"+ *isuffix]         = new TH1F(TString("yrec"+ *isuffix),"reconstructed y",100,-0.1,0.1);
00127     h["zrec"+ *isuffix]         = new TH1F(TString("zrec"+ *isuffix),"reconstructed z",100,-20.,20.);
00128     h["xsim"+ *isuffix]         = new TH1F(TString("xsim"+ *isuffix),"simulated x",100,-0.1,0.1);
00129     h["ysim"+ *isuffix]         = new TH1F(TString("ysim"+ *isuffix),"simulated y",100,-0.1,0.1);
00130     h["zsim"+ *isuffix]         = new TH1F(TString("zsim"+ *isuffix),"simulated z",100,-20.,20.);
00131     h["nrecvtx"+ *isuffix]      = new TH1F(TString("nrecvtx"+ *isuffix),"# of reconstructed vertices", 50, -0.5, 49.5);
00132     h["nsimvtx"+ *isuffix]      = new TH1F(TString("nsimvtx"+ *isuffix),"# of simulated vertices", 50, -0.5, 49.5);
00133     h["nrectrk"+ *isuffix]      = new TH1F(TString("nrectrk"+ *isuffix),"# of reconstructed tracks", 200, -0.5, 199.5);
00134     h["nsimtrk"+ *isuffix]      = new TH1F(TString("nsimtrk"+ *isuffix),"# of simulated tracks", 200, -0.5, 199.5);
00135     h["xrectag"+ *isuffix]      = new TH1F(TString("xrectag"+ *isuffix),"reconstructed x, signal vtx",100,-0.1,0.1);
00136     h["yrectag"+ *isuffix]      = new TH1F(TString("yrectag"+ *isuffix),"reconstructed y, signal vtx",100,-0.1,0.1);
00137     h["zrectag"+ *isuffix]      = new TH1F(TString("zrectag"+ *isuffix),"reconstructed z, signal vtx",100,-20.,20.);
00138     h["rapidity"+ *isuffix] = new TH1F(TString("rapidity"+ *isuffix),"rapidity ",100,-10., 10.);
00139     h["pt"+ *isuffix] = new TH1F(TString("pt"+ *isuffix),"pt ",100,0., 20.);
00140     h["nrectrk0vtx"+ *isuffix] = new TH1F(TString("nrectrk0vtx"+ *isuffix),"# rec tracks no vertex ",200,0., 200.);
00141     h["zdistancetag"+ *isuffix] = new TH1F(TString("zdistancetag"+ *isuffix),"z-distance between tagged and generated",100, -0.1, 0.1);
00142     h["puritytag"+ *isuffix]    = new TH1F(TString("puritytag"+ *isuffix),"purity of primary vertex tags",2, -0.5, 1.5);
00143   }
00144 
00145 }
00146 
00147 
00148 void PrimaryVertexAnalyzer::endJob() {
00149   rootFile_->cd();
00150   // save all histograms created in beginJob()
00151   for (std::map<std::string, TDirectory*>::const_iterator idir=hdir.begin(); idir!=hdir.end(); idir++){
00152     idir->second->cd();
00153     for(std::map<std::string,TH1*>::const_iterator hist=h.begin(); hist!=h.end(); hist++){
00154       if (TString(hist->first).Contains(idir->first)) hist->second->Write();
00155     }
00156   }
00157 }
00158 
00159 
00160 // helper functions
00161 bool PrimaryVertexAnalyzer::matchVertex(const simPrimaryVertex  &vsim, 
00162                                        const reco::Vertex       &vrec){
00163   return (fabs(vsim.z*simUnit_-vrec.z())<0.0500); // =500um
00164 }
00165 
00166 bool PrimaryVertexAnalyzer::isResonance(const HepMC::GenParticle * p){
00167   double ctau=(pdt->particle( abs(p->pdg_id()) ))->lifetime();
00168   if(verbose_) std::cout << "isResonance   " << p->pdg_id() << " " << ctau << std::endl;
00169   return  ctau >0 && ctau <1e-6;
00170 }
00171 
00172 bool PrimaryVertexAnalyzer::isFinalstateParticle(const HepMC::GenParticle * p){
00173   return ( !p->end_vertex() && p->status()==1 );
00174 }
00175 
00176 bool PrimaryVertexAnalyzer::isCharged(const HepMC::GenParticle * p){
00177   const ParticleData * part = pdt->particle( p->pdg_id() );
00178   if (part){
00179     return part->charge()!=0;
00180   }else{
00181     // the new/improved particle table doesn't know anti-particles
00182     return  pdt->particle( -p->pdg_id() )!=0;
00183   }
00184 }
00185 
00186 void PrimaryVertexAnalyzer::printRecVtxs(const Handle<reco::VertexCollection> recVtxs){
00187     int ivtx=0;
00188     for(reco::VertexCollection::const_iterator v=recVtxs->begin(); 
00189         v!=recVtxs->end(); ++v){
00190       std::cout << "Recvtx "<< std::setw(3) << std::setfill(' ')<<ivtx++
00191                 << "#trk " << std::setw(3) << v->tracksSize() 
00192                 << " chi2 " << std::setw(4) << v->chi2() 
00193                 << " ndof " << std::setw(3) << v->ndof() << std::endl 
00194                 << " x "  << std::setw(8) <<std::fixed << std::setprecision(4) << v->x() 
00195                 << " dx " << std::setw(8) << v->xError()<< std::endl
00196                 << " y "  << std::setw(8) << v->y() 
00197                 << " dy " << std::setw(8) << v->yError()<< std::endl
00198                 << " z "  << std::setw(8) << v->z() 
00199                 << " dz " << std::setw(8) << v->zError()
00200                 << std::endl;
00201     }
00202 }
00203 
00204 
00205 void PrimaryVertexAnalyzer::printSimVtxs(const Handle<SimVertexContainer> simVtxs){
00206     int i=0;
00207     for(SimVertexContainer::const_iterator vsim=simVtxs->begin();
00208         vsim!=simVtxs->end(); ++vsim){
00209       std::cout << i++ << ")" << std::scientific
00210                 << " evtid=" << vsim->eventId().event() 
00211                 << " sim x=" << vsim->position().x()*simUnit_
00212                 << " sim y=" << vsim->position().y()*simUnit_
00213                 << " sim z=" << vsim->position().z()*simUnit_
00214                 << " sim t=" << vsim->position().t()
00215                 << " parent=" << vsim->parentIndex() 
00216                 << std::endl;
00217     }
00218 }
00219 
00220 
00221 void PrimaryVertexAnalyzer::printSimTrks(const Handle<SimTrackContainer> simTrks){
00222   std::cout <<  " simTrks   type, (momentum), vertIndex, genpartIndex"  << std::endl;
00223   int i=1;
00224   for(SimTrackContainer::const_iterator t=simTrks->begin();
00225       t!=simTrks->end(); ++t){
00226     //HepMC::GenParticle* gp=evtMC->GetEvent()->particle( (*t).genpartIndex() );
00227     std::cout << i++ << ")" 
00228               << (*t)
00229               << " index="
00230               << (*t).genpartIndex();
00231     //if (gp) {
00232     //  HepMC::GenVertex *gv=gp->production_vertex();
00233     //  std::cout  <<  " genvertex =" << (*gv);
00234     //}
00235     std::cout << std::endl;
00236   }
00237 }
00238 
00239 
00240 std::vector<PrimaryVertexAnalyzer::simPrimaryVertex> PrimaryVertexAnalyzer::getSimPVs(
00241                                       const Handle<HepMCProduct> evtMC, std::string suffix="")
00242 {
00243   std::vector<PrimaryVertexAnalyzer::simPrimaryVertex> simpv;
00244   const HepMC::GenEvent* evt=evtMC->GetEvent();
00245   if (evt) {
00246     if(verbose_) std::cout << "process id " <<evt->signal_process_id()<<std::endl;
00247     if(verbose_) std::cout <<"signal process vertex "<< ( evt->signal_process_vertex() ?
00248                                              evt->signal_process_vertex()->barcode() : 0 )   <<std::endl;
00249     if(verbose_) std::cout <<"number of vertices " << evt->vertices_size() << std::endl;
00250 
00251 
00252     int idx=0;
00253     for(HepMC::GenEvent::vertex_const_iterator vitr= evt->vertices_begin();
00254         vitr != evt->vertices_end(); ++vitr ) 
00255       { // loop for vertex ...
00256         HepMC::FourVector pos = (*vitr)->position();
00257         //HepLorentzVector pos = (*vitr)->position();
00258 
00259         bool hasMotherVertex=false;
00260         if(verbose_) std::cout << "mothers of vertex[" << ++idx << "]: " << std::endl;
00261         for ( HepMC::GenVertex::particle_iterator
00262               mother  = (*vitr)->particles_begin(HepMC::parents);
00263               mother != (*vitr)->particles_end(HepMC::parents);
00264               ++mother ) {
00265           HepMC::GenVertex * mv=(*mother)->production_vertex();
00266           if (mv) {
00267             hasMotherVertex=true;
00268             if(!verbose_) break; //if verbose_, print all particles of gen vertices
00269           }
00270           if(verbose_)std::cout << "\t";
00271           if(verbose_)(*mother)->print();
00272         }
00273 
00274         if(hasMotherVertex) {continue;}
00275 
00276         // could be a new vertex, check  all primaries found so far to avoid multiple entries
00277         const double mm=0.1;
00278         simPrimaryVertex sv(pos.x()*mm,pos.y()*mm,pos.z()*mm);
00279         simPrimaryVertex *vp=NULL;  // will become non-NULL if a vertex is found and then point to it
00280         for(std::vector<simPrimaryVertex>::iterator v0=simpv.begin();
00281             v0!=simpv.end(); v0++){
00282           if( (fabs(sv.x-v0->x)<1e-5) && (fabs(sv.y-v0->y)<1e-5) && (fabs(sv.z-v0->z)<1e-5)){
00283             vp=&(*v0);
00284             break;
00285           }
00286         }
00287 
00288         if(!vp){
00289           // this is a new vertex
00290           if(verbose_)std::cout << "this is a new vertex " << sv.x << " " << sv.y << " " << sv.z << std::endl;
00291           simpv.push_back(sv);
00292           vp=&simpv.back();
00293         }else{
00294           if(verbose_)std::cout << "this is not new vertex " << sv.x << " " << sv.y << " " << sv.z << std::endl;
00295         }
00296         vp->genVertex.push_back((*vitr)->barcode());
00297         // collect final state descendants
00298         for ( HepMC::GenVertex::particle_iterator
00299               daughter  = (*vitr)->particles_begin(HepMC::descendants);
00300               daughter != (*vitr)->particles_end(HepMC::descendants);
00301               ++daughter ) {
00302           if (isFinalstateParticle(*daughter)){ 
00303             if ( find(vp->finalstateParticles.begin(), vp->finalstateParticles.end(),(*daughter)->barcode())
00304                  == vp->finalstateParticles.end()){
00305               vp->finalstateParticles.push_back((*daughter)->barcode());
00306               HepMC::FourVector m=(*daughter)->momentum();
00307               // the next four lines used to be "vp->ptot+=m;" in the days of CLHEP::HepLorentzVector
00308               // but adding FourVectors seems not to be foreseen
00309               vp->ptot.setPx(vp->ptot.px()+m.px());
00310               vp->ptot.setPy(vp->ptot.py()+m.py());
00311               vp->ptot.setPz(vp->ptot.pz()+m.pz());
00312               vp->ptot.setE(vp->ptot.e()+m.e());
00313               vp->ptsq+=(m.perp())*(m.perp());
00314               if ( (m.perp()>0.8) && (fabs(m.pseudoRapidity())<2.5) && isCharged( *daughter ) ){
00315                 vp->nGenTrk++;
00316               }
00317               
00318               h["rapidity"+suffix]->Fill(m.pseudoRapidity());
00319               h["pt"+suffix]->Fill(m.perp());
00320             }
00321           }
00322         }
00323       }
00324   }
00325   return simpv;
00326 }
00327 
00328 
00329 
00330 std::vector<PrimaryVertexAnalyzer::simPrimaryVertex> PrimaryVertexAnalyzer::getSimPVs(
00331                                       const Handle<HepMCProduct> evtMC, 
00332                                       const Handle<SimVertexContainer> simVtxs, 
00333                                       const Handle<SimTrackContainer> simTrks)
00334 {
00335    // simvertices don't have enough information to decide, 
00336    // genVertices don't have the simulated coordinates  ( with VtxSmeared they might)
00337    // go through simtracks to get the link between simulated and generated vertices
00338    std::vector<PrimaryVertexAnalyzer::simPrimaryVertex> simpv;
00339    int idx=0;
00340    for(SimTrackContainer::const_iterator t=simTrks->begin();
00341        t!=simTrks->end(); ++t){
00342      if ( !(t->noVertex()) && !(t->type()==-99) ){
00343        double ptsq=0;
00344        bool primary=false;   // something coming directly from the primary vertex
00345        bool resonance=false; // resonance
00346        bool track=false;     // undecayed, charged particle
00347        HepMC::GenParticle* gp=evtMC->GetEvent()->barcode_to_particle( (*t).genpartIndex() );
00348        if (gp) {
00349          HepMC::GenVertex * gv=gp->production_vertex();
00350          if (gv) {
00351            for ( HepMC::GenVertex::particle_iterator 
00352                  daughter =  gv->particles_begin(HepMC::descendants);
00353                  daughter != gv->particles_end(HepMC::descendants);
00354                  ++daughter ) {
00355              if (isFinalstateParticle(*daughter)){
00356                ptsq+=(*daughter)->momentum().perp()*(*daughter)->momentum().perp();
00357              }
00358            }
00359            primary =  ( gv->position().t()==0);
00360            //resonance= ( gp->mother() && isResonance(gp->mother()));  // in CLHEP/HepMC days
00361            // no more mother pointer in the improved HepMC GenParticle
00362            resonance= ( isResonance(*(gp->production_vertex()->particles_in_const_begin())));
00363            if (gp->status()==1){
00364              //track=((pdt->particle(gp->pdg_id()))->charge() != 0);
00365              track=not isCharged(gp);
00366            }
00367          }
00368        }
00369 
00370        const HepMC::FourVector & v=(*simVtxs)[t->vertIndex()].position();
00371        //const HepLorentzVector & v=(*simVtxs)[t->vertIndex()].position();
00372        if(primary or resonance){
00373          {
00374            // check all primaries found so far to avoid multiple entries
00375            bool newVertex=true;
00376            for(std::vector<simPrimaryVertex>::iterator v0=simpv.begin();
00377                v0!=simpv.end(); v0++){
00378              if( (fabs(v0->x-v.x())<0.001) && (fabs(v0->y-v.y())<0.001) && (fabs(v0->z-v.z())<0.001) ){
00379                if (track) {
00380                  v0->simTrackIndex.push_back(idx);
00381                  if (ptsq>(*v0).ptsq){(*v0).ptsq=ptsq;}
00382                }
00383                newVertex=false;
00384              }
00385            }
00386            if(newVertex && !resonance){
00387              simPrimaryVertex anotherVertex(v.x(),v.y(),v.z());
00388              if (track) anotherVertex.simTrackIndex.push_back(idx);
00389              anotherVertex.ptsq=ptsq;
00390              simpv.push_back(anotherVertex);
00391            }
00392          }// 
00393        }
00394        
00395      }// simtrack has vertex and valid type
00396      idx++;
00397    }//simTrack loop
00398    return simpv;
00399 }
00400 
00401 
00402 
00403 
00404 // ------------ method called to produce the data  ------------
00405 void
00406 PrimaryVertexAnalyzer::analyze(const Event& iEvent, const EventSetup& iSetup)
00407 {
00408   
00409   Handle<reco::TrackCollection> recTrks;
00410   iEvent.getByLabel(recoTrackProducer_, recTrks);
00411 
00412   Handle<SimVertexContainer> simVtxs;
00413   iEvent.getByLabel( simG4_, simVtxs);
00414   
00415   Handle<SimTrackContainer> simTrks;
00416   iEvent.getByLabel( simG4_, simTrks);
00417 
00418   for (int ivtxSample=0; ivtxSample!= (int)vtxSample_.size(); ++ivtxSample) {
00419     std::string isuffix = suffixSample_[ivtxSample];
00420 
00421     Handle<reco::VertexCollection> recVtxs;
00422     iEvent.getByLabel(vtxSample_[ivtxSample], recVtxs);
00423 
00424     try{
00425       iSetup.getData(pdt);
00426     }catch(const Exception&){
00427       std::cout << "Some problem occurred with the particle data table. This may not work !." <<std::endl;
00428     }
00429 
00430     bool MC=false;
00431     Handle<HepMCProduct> evtMC;
00432     iEvent.getByLabel("generator",evtMC);
00433     if (!evtMC.isValid()) {
00434       MC=false;
00435       if(verbose_){
00436         std::cout << "no HepMCProduct found"<< std::endl;
00437       }
00438     } else {
00439       if(verbose_){
00440         std::cout << "generator HepMCProduct found"<< std::endl;
00441       }
00442       MC=true;
00443     }
00444     
00445     
00446     if(verbose_){
00447       //evtMC->GetEvent()->print();
00448       printRecVtxs(recVtxs);
00449       //printSimVtxs(simVtxs);
00450       //printSimTrks(simTrks);
00451     }
00452     
00453     if(MC){
00454       // make a list of primary vertices:
00455       std::vector<simPrimaryVertex> simpv;
00456       //simpv=getSimPVs(evtMC, simVtxs, simTrks);
00457       simpv=getSimPVs(evtMC,isuffix);
00458       
00459       
00460       // vertex matching and efficiency bookkeeping
00461       h["nsimvtx"+isuffix]->Fill(simpv.size());
00462       int nsimtrk=0;
00463       for(std::vector<simPrimaryVertex>::iterator vsim=simpv.begin();
00464           vsim!=simpv.end(); vsim++){
00465         
00466         hdir[isuffix]->cd();
00467         
00468         h["nbsimtksinvtx"+isuffix]->Fill(vsim->nGenTrk);
00469         h["xsim"+isuffix]->Fill(vsim->x*simUnit_);
00470         h["ysim"+isuffix]->Fill(vsim->y*simUnit_);
00471         h["zsim"+isuffix]->Fill(vsim->z*simUnit_);
00472         nsimtrk+=vsim->nGenTrk;
00473         // look for a matching reconstructed vertex
00474         vsim->recVtx=NULL;
00475         for(reco::VertexCollection::const_iterator vrec=recVtxs->begin(); 
00476             vrec!=recVtxs->end(); ++vrec){
00477           if ( matchVertex(*vsim,*vrec) ){
00478             // if the matching critera are fulfilled, accept the rec-vertex that is closest in z
00479             if(    ((vsim->recVtx) && (fabs(vsim->recVtx->position().z()-vsim->z)>fabs(vrec->z()-vsim->z)))
00480                    || (!vsim->recVtx) )
00481               {
00482                 vsim->recVtx=&(*vrec);
00483               }
00484           }
00485         }
00486         h["nsimtrk"+isuffix]->Fill(float(nsimtrk));
00487         
00488         
00489         // histogram properties of matched vertices
00490         if (vsim->recVtx){
00491           
00492           if(verbose_){std::cout <<"primary matched " << vsim->x << " " << vsim->y << " " << vsim->z << std:: endl;}
00493           
00494           h["resx"+isuffix]->Fill( vsim->recVtx->x()-vsim->x*simUnit_ );
00495           h["resy"+isuffix]->Fill( vsim->recVtx->y()-vsim->y*simUnit_ );
00496           h["resz"+isuffix]->Fill( vsim->recVtx->z()-vsim->z*simUnit_ );
00497           h["pullx"+isuffix]->Fill( (vsim->recVtx->x()-vsim->x*simUnit_)/vsim->recVtx->xError() );
00498           h["pully"+isuffix]->Fill( (vsim->recVtx->y()-vsim->y*simUnit_)/vsim->recVtx->yError() );
00499           h["pullz"+isuffix]->Fill( (vsim->recVtx->z()-vsim->z*simUnit_)/vsim->recVtx->zError() );
00500           h["eff"+isuffix]->Fill( 1.);
00501           if(simpv.size()==1){
00502             if (vsim->recVtx==&(*recVtxs->begin())){
00503               h["efftag"+isuffix]->Fill( 1.); 
00504             }else{
00505               h["efftag"+isuffix]->Fill( 0.); 
00506             }
00507           }
00508           h["effvseta"+isuffix]->Fill(vsim->ptot.pseudoRapidity(),1.);
00509           h["effvsptsq"+isuffix]->Fill(vsim->ptsq,1.);
00510           h["effvsntrk"+isuffix]->Fill(vsim->nGenTrk,1.);
00511           h["effvsnrectrk"+isuffix]->Fill(recTrks->size(),1.);
00512           h["effvsz"+isuffix]->Fill(vsim->z*simUnit_,1.);
00513           
00514         }else{  // no rec vertex found for this simvertex
00515           
00516           if(verbose_){std::cout <<"primary not found " << vsim->x << " " << vsim->y << " " << vsim->z << " nGenTrk=" << vsim->nGenTrk << std::endl;}
00517           
00518           h["eff"+isuffix]->Fill( 0.);
00519           if(simpv.size()==1){ h["efftag"+isuffix]->Fill( 0.); }
00520           h["effvseta"+isuffix]->Fill(vsim->ptot.pseudoRapidity(),0.);
00521           h["effvsptsq"+isuffix]->Fill(vsim->ptsq,0.);
00522           h["effvsntrk"+isuffix]->Fill(float(vsim->nGenTrk),0.);
00523           h["effvsnrectrk"+isuffix]->Fill(recTrks->size(),0.);
00524           h["effvsz"+isuffix]->Fill(vsim->z*simUnit_,0.);
00525         } // no recvertex for this simvertex
00526       }//found MC event
00527       // end of sim/rec matching 
00528       
00529       
00530       // purity of event vertex tags
00531       if (recVtxs->size()>0 && simpv.size()>0){
00532         Double_t dz=(*recVtxs->begin()).z() - (*simpv.begin()).z*simUnit_;
00533         h["zdistancetag"+isuffix]->Fill(dz);
00534         if( fabs(dz)<0.0500){
00535           h["puritytag"+isuffix]->Fill(1.);
00536         }else{
00537           // bad tag: the true primary was more than 500 um away from the tagged primary
00538           h["puritytag"+isuffix]->Fill(0.);
00539         }
00540       }
00541       
00542       
00543     }// MC event
00544     
00545     // test track links, use reconstructed vertices
00546     
00547     h["nrecvtx"+isuffix]->Fill(recVtxs->size());
00548     h["nrectrk"+isuffix]->Fill(recTrks->size());
00549     h["nbvtx"+isuffix]->Fill(recVtxs->size()*1.);
00550     if((recVtxs->size()==0)||recVtxs->begin()->isFake()) {h["nrectrk0vtx"+isuffix]->Fill(recTrks->size());}
00551     
00552     for(reco::VertexCollection::const_iterator v=recVtxs->begin(); 
00553         v!=recVtxs->end(); ++v){
00554 
00555       if(v->isFake()) continue;
00556       
00557       try {
00558         for(reco::Vertex::trackRef_iterator t = v->tracks_begin(); 
00559             t!=v->tracks_end(); t++) {
00560           // illegal charge
00561           if ( (**t).charge() < -1 || (**t).charge() > 1 ) {
00562             h["tklinks"+isuffix]->Fill(0.);
00563           }
00564           else {
00565             h["tklinks"+isuffix]->Fill(1.);
00566           }
00567         }
00568       }
00569       catch (...) {
00570         // exception thrown when trying to use linked track
00571         h["tklinks"+isuffix]->Fill(0.);
00572       }
00573       
00574       h["nbtksinvtx"+isuffix]->Fill(v->tracksSize());
00575       h["vtxchi2"+isuffix]->Fill(v->chi2());
00576       h["vtxndf"+isuffix]->Fill(v->ndof());
00577       h["vtxprob"+isuffix]->Fill(ChiSquaredProbability(v->chi2() ,v->ndof()));
00578       h["xrec"+isuffix]->Fill(v->position().x());
00579       h["yrec"+isuffix]->Fill(v->position().y());
00580       h["zrec"+isuffix]->Fill(v->position().z());
00581       if (v==recVtxs->begin()){
00582         h["xrectag"+isuffix]->Fill(v->position().x());
00583         h["yrectag"+isuffix]->Fill(v->position().y());
00584         h["zrectag"+isuffix]->Fill(v->position().z());
00585       }
00586       
00587       bool problem = false;
00588       h["nans"+isuffix]->Fill(1.,isnan(v->position().x())*1.);
00589       h["nans"+isuffix]->Fill(2.,isnan(v->position().y())*1.);
00590       h["nans"+isuffix]->Fill(3.,isnan(v->position().z())*1.);
00591       
00592       int index = 3;
00593       for (int i = 0; i != 3; i++) {
00594         for (int j = i; j != 3; j++) {
00595           index++;
00596           h["nans"+isuffix]->Fill(index*1., isnan(v->covariance(i, j))*1.);
00597           if (isnan(v->covariance(i, j))) problem = true;
00598           // in addition, diagonal element must be positive
00599           if (j == i && v->covariance(i, j) < 0) {
00600             h["nans"+isuffix]->Fill(index*1., 1.);
00601             problem = true;
00602           }
00603         }
00604       }
00605       
00606       if (problem) {
00607         // analyze track parameter covariance definiteness
00608         double data[25];
00609         try {
00610           int itk = 0;
00611           for(reco::Vertex::trackRef_iterator t = v->tracks_begin(); 
00612               t!=v->tracks_end(); t++) {
00613             std::cout << "Track " << itk++ << std::endl;
00614             int i2 = 0;
00615             for (int i = 0; i != 5; i++) {
00616               for (int j = 0; j != 5; j++) {
00617                 data[i2] = (**t).covariance(i, j);
00618                 std::cout << std:: scientific << data[i2] << " ";
00619                 i2++;
00620               }
00621               std::cout << std::endl;
00622             }
00623             gsl_matrix_view m 
00624               = gsl_matrix_view_array (data, 5, 5);
00625             
00626             gsl_vector *eval = gsl_vector_alloc (5);
00627             gsl_matrix *evec = gsl_matrix_alloc (5, 5);
00628             
00629             gsl_eigen_symmv_workspace * w = 
00630               gsl_eigen_symmv_alloc (5);
00631           
00632             gsl_eigen_symmv (&m.matrix, eval, evec, w);
00633             
00634             gsl_eigen_symmv_free (w);
00635             
00636             gsl_eigen_symmv_sort (eval, evec, 
00637                                   GSL_EIGEN_SORT_ABS_ASC);
00638             
00639             // print sorted eigenvalues
00640             {
00641               int i;
00642               for (i = 0; i < 5; i++) {
00643                 double eval_i 
00644                   = gsl_vector_get (eval, i);
00645                 //gsl_vector_view evec_i 
00646                 //  = gsl_matrix_column (evec, i);
00647                 
00648                 printf ("eigenvalue = %g\n", eval_i);
00649                 //            printf ("eigenvector = \n");
00650                 //            gsl_vector_fprintf (stdout, 
00651                 //                                &evec_i.vector, "%g");
00652               }
00653             }
00654           }
00655         }
00656         catch (...) {
00657           // exception thrown when trying to use linked track
00658           break;
00659         }
00660       }// if (problem)
00661     }
00662   } // for vertex loop
00663 }