/home/bes3soft/bes3soft/Boss/7.0.2/dist/7.0.2/Reconstruction/SD0TagAlg/SD0TagAlg-00-00-03/src/Kpipi0.cxx

Go to the documentation of this file.
00001 //
00002 //  Kpipi0.cxx is the single D0 tag code to reconstruct D0 or anti-D0 through the final states of
00003 //  kpipi0 from D0 decays. Kpipi0.cxx was transfered from the Fortran routine "Kpipi0.f"
00004 //  which was orignally used for study of the D0D0-bar production and D0 decays at the BES-II
00005 //  experiment during the time period from 2002 to 2008.
00006 //
00007 //  The orignal Fortran routine "Kpipi0.f" used at the BES-II experiment was coded by G. Rong in 2003.
00008 //
00009 //  Kpipi0.cxx was transfered by G. Rong and J. Liu in December, 2005.
00010 //
00011 //  Since 2008, G. Rong and L.L. Jiang have been working on developing this code to analyze of
00012 //  the data taken at 3.773 GeV with the BES-III detector at the BEPC-II collider.
00013 //
00014 //  During developing this code, many People made significant contributions to this code. These are
00015 //          G. Rong, L.L. Jiang, J. Liu, H.L. Ma, J.C. Chen, D.H. Zhang,
00016 //          M.G. Zhao, B. Zheng, L. Li, Y. Fang, Z.Y. Yi, H.H. Liu, Z.Q. Liu et al.
00017 //
00018 //                                       By G. Rong and L.L. Jiang
00019 //                                       March, 2009
00020 //
00021 //  ==========================================================================================
00022 //
00023 #include "SD0TagAlg/Kpipi0.h"
00024 #include "SD0TagAlg/SingleBase.h"
00025 
00026 Kpipi0::Kpipi0()
00027 {}
00028 
00029 Kpipi0::~Kpipi0()
00030 {}
00031 
00032 
00033 void Kpipi0::MTotal(double event,SmartDataPtr<EvtRecTrackCol> evtRecTrkCol, Vint iGood,Vint
00034     iGam, double Ebeam, int PID_flag, int Charge_candidate_D)
00035 {
00036 
00037   int nGood=iGood.size();
00038   int nGam=iGam.size();
00039 
00040   iGoodtag.clear();
00041   iGamtag.clear();
00042   
00043   double mass_bcgg,delE_tag_temp;
00044   int m_chargetag, m_chargek,m_chargepi;
00045   int ika_temp,ipi_temp, iGam1_temp, iGam2_temp;
00046   HepLorentzVector pddd, pddd_temp;
00047 
00048   int cqtm_temp;
00049   IDataProviderSvc* eventSvc = NULL;
00050   Gaudi::svcLocator()->service("EventDataSvc", eventSvc);
00051   SmartDataPtr<EvtRecEvent> evtRecEvent(eventSvc,EventModel::EvtRec::EvtRecEvent);
00052   SmartDataPtr<Event::EventHeader> eventHeader(eventSvc,"/Event/EventHeader");
00053 
00054   int runNo=eventHeader->runNumber();
00055   int rec=eventHeader->eventNumber();
00056 
00057   double xecm=2*Ebeam;
00058 
00059   kpipi0md=false;
00060   double  tagmode=0;
00061 
00062   if((evtRecEvent->totalCharged() < 2||nGam<2)){     return;  }
00063 
00064   double ecms = xecm;
00065 
00066   ISimplePIDSvc* simple_pid;
00067   Gaudi::svcLocator()->service("SimplePIDSvc", simple_pid);
00068 
00069   double deltaE_tem = 0.20;
00070   int ncount1 = 0; 
00071 
00072   Hep3Vector xorigin(0,0,0);
00073   IVertexDbSvc*  vtxsvc;
00074   Gaudi::svcLocator()->service("VertexDbSvc", vtxsvc);
00075   if(vtxsvc->isVertexValid())
00076   {
00077     double* dbv = vtxsvc->PrimaryVertex();
00078     double*  vv = vtxsvc->SigmaPrimaryVertex();
00079     xorigin.setX(dbv[0]);
00080     xorigin.setY(dbv[1]);
00081     xorigin.setZ(dbv[2]);
00082   }
00083 
00084 
00085   double xv=xorigin.x();
00086   double yv=xorigin.y();
00087   double zv=xorigin.z();
00088 
00089   HepPoint3D point0(0.,0.,0.);
00090   HepPoint3D IP(xorigin[0],xorigin[1],xorigin[2]);
00091 
00092   HepLorentzVector ptrk1_temp, ptrk2_temp, ptrk3_temp, ptrk4_temp, ptrk5_temp;
00094   HepLorentzVector p2gfit;
00095   HepLorentzVector p2gg;
00096   for(int i = 0; i < evtRecEvent->totalCharged(); i++) {
00097     EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + i;
00098 
00099     int ika= (*itTrk)->trackId();
00100 
00101     if(!(*itTrk)->isMdcKalTrackValid()) continue;
00102     RecMdcKalTrack*  mdcKalTrk1 = (*itTrk)->mdcKalTrack();
00103     RecMdcKalTrack::setPidType(RecMdcKalTrack::kaon);
00105     m_chargek=mdcKalTrk1->charge();
00106     if(Charge_candidate_D != 0) {
00107       if(m_chargek != -Charge_candidate_D) continue;
00108     }
00109     if(Charge_candidate_D == 0) {
00110       if(abs(m_chargek) != 1) continue;
00111     }
00113     HepVector a1 = mdcKalTrk1->getZHelixK();
00114     HepSymMatrix Ea1 = mdcKalTrk1->getZErrorK();
00115     VFHelix helixip3_1(point0,a1,Ea1);
00116     helixip3_1.pivot(IP);
00117     HepVector  vecipa1 = helixip3_1.a();
00118 
00119     double dr1 = fabs(vecipa1[0]);
00120     double dz1 = fabs(vecipa1[3]);
00121     double costheta1 = cos(mdcKalTrk1->theta());
00122     if (  dr1 >= 1.0) continue;
00123     if (  dz1 >= 10.0) continue; 
00124     if ( fabs(costheta1) >= 0.93) continue; 
00125 
00127     if(PID_flag == 5) {
00128       simple_pid->preparePID(*itTrk);
00129       if(simple_pid->probKaon() < 0.0 || simple_pid->probKaon() < simple_pid->probPion()) continue;  
00130     } 
00132 
00133     WTrackParameter kam(xmass[3],mdcKalTrk1->getZHelixK(),mdcKalTrk1->getZErrorK() );
00134 
00135     //  
00136     // select pi
00137     //  
00138     for(int j = 0; j< evtRecEvent->totalCharged();j++) {
00139       EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + j;
00140 
00141       int ipi= (*itTrk)->trackId();
00142       if(ipi==ika)  continue;
00143 
00144       if(!(*itTrk)->isMdcKalTrackValid()) continue;
00145       RecMdcKalTrack*  mdcKalTrk2 = (*itTrk)->mdcKalTrack();
00146       RecMdcKalTrack::setPidType(RecMdcKalTrack::pion);
00148       m_chargepi=mdcKalTrk2->charge();
00149       if((m_chargek + m_chargepi) != 0) continue;
00151       HepVector a2 = mdcKalTrk2->getZHelix();
00152       HepSymMatrix Ea2 = mdcKalTrk2->getZError();
00153       VFHelix helixip3_2(point0,a2,Ea2);
00154       helixip3_2.pivot(IP);
00155       HepVector  vecipa2 = helixip3_2.a();
00156 
00157       double dr2 = fabs(vecipa2[0]);
00158       double dz2 = fabs(vecipa2[3]);
00159       double costheta2 = cos(mdcKalTrk2->theta());
00160       if (  dr2 >= 1.0) continue;
00161       if (  dz2 >= 10.0) continue; 
00162       if ( fabs(costheta2) >= 0.93) continue; 
00164       if(PID_flag == 5) {
00165         simple_pid->preparePID(*itTrk);
00166         if(simple_pid->probPion() < 0.0 || simple_pid->probPion() < simple_pid->probKaon()) continue;  
00167       } 
00169 
00170       WTrackParameter pip(xmass[2],mdcKalTrk2->getZHelix(),mdcKalTrk2->getZError() );
00171 
00172       for(int m = 0; m < nGam-1; m++) {
00173         if(iGam[m]==-1) continue;
00174         RecEmcShower *g1Trk = (*(evtRecTrkCol->begin()+iGam[m]))->emcShower();
00175         double eraw1 = g1Trk->energy();
00176         double phi1 = g1Trk->phi();
00177         double the1 = g1Trk->theta();
00178         HepLorentzVector ptrkg1,ptrkg10,ptrkg12;
00179         ptrkg1.setPx(eraw1*sin(the1)*cos(phi1));
00180         ptrkg1.setPy(eraw1*sin(the1)*sin(phi1));
00181         ptrkg1.setPz(eraw1*cos(the1));
00182         ptrkg1.setE(eraw1);
00183         ptrkg10 = ptrkg1;
00184         ptrkg12 = ptrkg1.boost(-0.011,0,0);
00185 
00186         for(int n = m+1; n < nGam; n++) {
00187           if(iGam[n]==-1) continue;
00188           RecEmcShower *g2Trk = (*(evtRecTrkCol->begin()+iGam[n]))->emcShower();
00189           double eraw2 = g2Trk->energy();
00190           double phi2 = g2Trk->phi();
00191           double the2 = g2Trk->theta();
00192           HepLorentzVector ptrkg2,ptrkg20,ptrkg22;
00193           ptrkg2.setPx(eraw2*sin(the2)*cos(phi2));
00194           ptrkg2.setPy(eraw2*sin(the2)*sin(phi2));
00195           ptrkg2.setPz(eraw2*cos(the2));
00196           ptrkg2.setE(eraw2);
00197           ptrkg20 = ptrkg2;
00198           ptrkg22 = ptrkg2.boost(-0.011,0,0);
00199 
00201           HepLorentzVector  ptrkpi0;
00202           ptrkpi0 = ptrkg12+ptrkg22;
00203           double m_xmpi0_tem = ptrkpi0.m();
00204           if(m_xmpi0_tem>0.150||m_xmpi0_tem<0.115)  continue;
00206           bool IsEndcap1 = false; bool IsEndcap2 = false;
00207           if(fabs(cos(the1)) > 0.86 && fabs(cos(the1)) < 0.92) IsEndcap1 = true; 
00208           if(fabs(cos(the2)) > 0.86 && fabs(cos(the2)) < 0.92) IsEndcap2 = true; 
00209           if(IsEndcap1 && IsEndcap2)  continue;
00211           KalmanKinematicFit * kmfit = KalmanKinematicFit::instance();
00212           kmfit->init();
00213           kmfit->setChisqCut(2500);
00214           kmfit->AddTrack(0, 0.0, g1Trk);
00215           kmfit->AddTrack(1, 0.0, g2Trk);
00216           kmfit->AddResonance(0, mpi0, 0, 1);
00217 
00218           kmfit->Fit(0);  // Perform fit
00219           kmfit->BuildVirtualParticle(0);
00220 
00221           double pi0_chisq = kmfit->chisq(0);
00222           if ( pi0_chisq >= 2500) continue;
00223           HepLorentzVector p2gfit = kmfit->pfit(0) + kmfit->pfit(1);
00224           p2gfit.boost(-0.011,0,0);
00225 
00227           HepPoint3D vx(xorigin.x(), xorigin.y(), xorigin.z());
00228           HepSymMatrix Evx(3, 0);
00229           double bx = 1E+6; Evx[0][0] = bx*bx;
00230           double by = 1E+6; Evx[1][1] = by*by;
00231           double bz = 1E+6; Evx[2][2] = bz*bz;
00232           VertexParameter vxpar; vxpar.setVx(vx); vxpar.setEvx(Evx);
00234 
00235           VertexFit* vtxfit = VertexFit::instance();
00236           vtxfit->init();
00237           vtxfit->AddTrack(0,  kam);
00238           vtxfit->AddTrack(1,  pip);
00239           vtxfit->AddVertex(0, vxpar, 0, 1);
00240           if(!vtxfit->Fit(0))  continue;
00241           vtxfit->Swim(0);
00242 
00243           WTrackParameter wkam = vtxfit->wtrk(0);
00244           WTrackParameter wpip = vtxfit->wtrk(1);
00245 
00246           HepVector kam_val = HepVector(7,0);
00247           kam_val = wkam.w();
00248           HepVector pip_val = HepVector(7,0);
00249           pip_val = wpip.w();
00250 
00251           HepLorentzVector P_KAM(kam_val[0],kam_val[1],kam_val[2],kam_val[3]);
00252           HepLorentzVector P_PIP(pip_val[0],pip_val[1],pip_val[2],pip_val[3]);
00253 
00254           P_KAM.boost(-0.011,0,0);
00255           P_PIP.boost(-0.011,0,0);
00256           pddd = P_KAM + P_PIP + p2gfit;
00257 
00258           double pkpipi0=pddd.rho();
00259 
00260           double temp1 = (ecms/2)*(ecms/2)-pkpipi0*pkpipi0 ;
00261           if(temp1<0) temp1 =0;
00262           double mass_bc_tem  = sqrt(temp1);
00263           if(mass_bc_tem < 1.82 || mass_bc_tem > 1.89)   continue;
00264 
00265           double  delE_tag_tag = ecms/2-pddd.e();
00266 
00267 
00268           if(fabs(delE_tag_tag)<deltaE_tem) {
00269             deltaE_tem = fabs(delE_tag_tag);
00270             delE_tag_temp = delE_tag_tag;
00271             mass_bcgg = mass_bc_tem;
00272 
00273             pddd_temp = pddd; 
00274             cqtm_temp = m_chargek; 
00275             
00276             ika_temp=ika;
00277             ipi_temp=ipi;
00278 
00279             iGam1_temp = iGam[m];
00280             iGam2_temp = iGam[n];
00281             ncount1 = 1;
00282 
00283           }
00284         }
00285       }
00286     }
00287   }
00288 
00289   if(ncount1 == 1){
00290     tagmode=12;
00291     if(cqtm_temp <0) tagmode=-12;
00292     tagmd=tagmode;
00293     mass_bc = mass_bcgg;
00294     delE_tag = delE_tag_temp;
00295     cqtm     = -1.0*cqtm_temp;
00296 
00297     iGoodtag.push_back(ipi_temp);
00298     iGoodtag.push_back(ika_temp);
00299     
00300     iGamtag.push_back(iGam1_temp);
00301     iGamtag.push_back(iGam2_temp);
00302     iGamtag.push_back(9999);
00303     iGamtag.push_back(9999);
00304       
00305     ptag = pddd_temp;
00306     
00307     kpipi0md = true;
00308 
00309   }
00310 }
00311 
00312 
00313 

Generated on Tue Nov 29 23:14:03 2016 for BOSS_7.0.2 by  doxygen 1.4.7