Pipipi0 Class Reference

#include <Pipipi0.h>

List of all members.

Public Member Functions

 Pipipi0 ()
 ~Pipipi0 ()
bool GetPipipi0md ()
double Gettagmd ()
double Getmass_bc ()
double GetCQtm ()
double GetdelE_tag ()
Vint Gettagtrk1 ()
HepLorentzVector Gettagp1 ()
Vint GettagGam1 ()
void MTotal (double event, SmartDataPtr< EvtRecTrackCol > evtRecTrkCol, Vint iGood, Vint iGam, double Ebeam, int PID_flag, int Charge_candidate_D)

Private Attributes

bool pipipi0md
double tagmd
double mass_bc
double cqtm
double delE_tag
Vint iGoodtag
HepLorentzVector ptag
Vint iGamtag


Detailed Description

Definition at line 20 of file Pipipi0.h.


Constructor & Destructor Documentation

Pipipi0::Pipipi0 (  ) 

Definition at line 28 of file Pipipi0.cxx.

00029 {}

Pipipi0::~Pipipi0 (  ) 

Definition at line 31 of file Pipipi0.cxx.

00032 {}


Member Function Documentation

double Pipipi0::GetCQtm (  )  [inline]

Definition at line 30 of file Pipipi0.h.

References cqtm.

Referenced by Sing::Mdset().

00030 { return cqtm; }

double Pipipi0::GetdelE_tag (  )  [inline]

Definition at line 31 of file Pipipi0.h.

References delE_tag.

Referenced by Sing::Mdset().

00031 { return delE_tag; }

double Pipipi0::Getmass_bc (  )  [inline]

Definition at line 29 of file Pipipi0.h.

References mass_bc.

Referenced by Sing::Mdset().

00029 { return mass_bc; }

bool Pipipi0::GetPipipi0md (  )  [inline]

Definition at line 27 of file Pipipi0.h.

References pipipi0md.

Referenced by Sing::Mdset().

00027 { return pipipi0md; }

Vint Pipipi0::GettagGam1 (  )  [inline]

Definition at line 34 of file Pipipi0.h.

References iGamtag.

Referenced by Sing::Mdset().

00034 { return iGamtag; }

double Pipipi0::Gettagmd (  )  [inline]

Definition at line 28 of file Pipipi0.h.

References tagmd.

Referenced by Sing::Mdset().

00028 { return tagmd; }

HepLorentzVector Pipipi0::Gettagp1 (  )  [inline]

Definition at line 33 of file Pipipi0.h.

References ptag.

Referenced by Sing::Mdset().

00033 { return ptag; }

Vint Pipipi0::Gettagtrk1 (  )  [inline]

Definition at line 32 of file Pipipi0.h.

References iGoodtag.

Referenced by Sing::Mdset().

00032 { return iGoodtag; }

void Pipipi0::MTotal ( double  event,
SmartDataPtr< EvtRecTrackCol evtRecTrkCol,
Vint  iGood,
Vint  iGam,
double  Ebeam,
int  PID_flag,
int  Charge_candidate_D 
)

Definition at line 35 of file Pipipi0.cxx.

References VFHelix::a(), abs, KalmanKinematicFit::AddResonance(), TrackPool::AddTrack(), VertexFit::AddVertex(), KalmanKinematicFit::BuildVirtualParticle(), DstMdcKalTrack::charge(), KalmanKinematicFit::chisq(), cos(), cqtm, delE_tag, ecms, DstEmcShower::energy(), EventModel::EvtRec::EvtRecEvent, VertexFit::Fit(), KalmanKinematicFit::Fit(), RecMdcKalTrack::getZError(), RecMdcKalTrack::getZHelix(), genRecEmupikp::i, iGamtag, iGoodtag, VertexFit::init(), KalmanKinematicFit::init(), VertexFit::instance(), KalmanKinematicFit::instance(), IVertexDbSvc::isVertexValid(), ganga-rec::j, mass_bc, mpi0, KalmanKinematicFit::pfit(), DstEmcShower::phi(), phi1, phi2, DstMdcKalTrack::pion, pipipi0md, VFHelix::pivot(), ISimplePIDSvc::preparePID(), IVertexDbSvc::PrimaryVertex(), ISimplePIDSvc::probKaon(), ISimplePIDSvc::probPion(), ptag, runNo, KalmanKinematicFit::setChisqCut(), VertexParameter::setEvx(), DstMdcKalTrack::setPidType(), VertexParameter::setVx(), IVertexDbSvc::SigmaPrimaryVertex(), sin(), VertexFit::Swim(), tagmd, DstEmcShower::theta(), DstMdcKalTrack::theta(), WTrackParameter::w(), VertexFit::wtrk(), and xmass.

Referenced by Sing::Mdset().

00037 {
00038 
00039   int nGood=iGood.size();
00040   int nGam=iGam.size();
00041 
00042   iGoodtag.clear();
00043   iGamtag.clear();
00044   
00045   double mass_bcgg,delE_tag_temp;
00046   int m_chargetag,m_chargepi1,m_chargepi2;
00047   int ipi1_temp, ipi2_temp, iGam1_temp, iGam2_temp;
00048   HepLorentzVector pddd, pddd_temp;
00049 
00050   IDataProviderSvc* eventSvc = NULL;
00051   Gaudi::svcLocator()->service("EventDataSvc", eventSvc);
00052   SmartDataPtr<EvtRecEvent> evtRecEvent(eventSvc,EventModel::EvtRec::EvtRecEvent);
00053   SmartDataPtr<Event::EventHeader> eventHeader(eventSvc,"/Event/EventHeader");
00054 
00055   int runNo=eventHeader->runNumber();
00056   int rec=eventHeader->eventNumber();
00057 
00058   double xecm=2*Ebeam;
00059 
00060   pipipi0md=false;
00061   double  tagmode=0;
00062 
00063   if((evtRecEvent->totalCharged() < 2||nGam<2)){    return;  }
00064 
00065   double ecms = xecm; 
00066 
00067   ISimplePIDSvc* simple_pid;
00068   Gaudi::svcLocator()->service("SimplePIDSvc", simple_pid);
00069 
00070   double deltaE_tem = 0.20;
00071   int ncount1 = 0; 
00072 
00073   Hep3Vector xorigin(0,0,0);
00074   IVertexDbSvc*  vtxsvc;
00075   Gaudi::svcLocator()->service("VertexDbSvc", vtxsvc);
00076   if(vtxsvc->isVertexValid())
00077   {
00078     double* dbv = vtxsvc->PrimaryVertex();
00079     double*  vv = vtxsvc->SigmaPrimaryVertex();
00080     xorigin.setX(dbv[0]);
00081     xorigin.setY(dbv[1]);
00082     xorigin.setZ(dbv[2]);
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]);
00092 
00093   HepLorentzVector p2gfit;
00094   HepLorentzVector p2gg;
00095 
00096   for(int i = 0; i < evtRecEvent->totalCharged(); i++) {
00097     EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + i;
00098 
00099     int ipi1= (*itTrk)->trackId();
00100 
00101     if(!(*itTrk)->isMdcKalTrackValid()) continue;
00102     RecMdcKalTrack*  mdcKalTrk1 = (*itTrk)->mdcKalTrack();
00103     RecMdcKalTrack::setPidType(RecMdcKalTrack::pion);
00105     m_chargepi1=mdcKalTrk1->charge();
00106     if(abs(m_chargepi1) != 1) continue;
00108     HepVector a1 = mdcKalTrk1->getZHelix();
00109     HepSymMatrix Ea1 = mdcKalTrk1->getZError();
00110     VFHelix helixip3_1(point0,a1,Ea1);
00111     helixip3_1.pivot(IP);
00112     HepVector  vecipa1 = helixip3_1.a();
00113 
00114     double dr1 = fabs(vecipa1[0]);
00115     double dz1 = fabs(vecipa1[3]);
00116     double costheta1 = cos(mdcKalTrk1->theta());
00117 
00118     if (  dr1 >= 1.0) continue;
00119     if (  dz1 >= 10.0) continue; 
00120     if ( fabs(costheta1) >= 0.93) continue; 
00122     if(PID_flag == 5) {
00123       simple_pid->preparePID(*itTrk);
00124       if(simple_pid->probPion() < 0.0 || simple_pid->probPion() < simple_pid->probKaon()) continue;  
00125     } 
00127     WTrackParameter pip(xmass[2],mdcKalTrk1->getZHelix(),mdcKalTrk1->getZError() );
00128 
00129     //  
00130     // select pi
00131     //  
00132     for(int j = 0; j< evtRecEvent->totalCharged();j++) {
00133       EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + j;
00134 
00135       int ipi2= (*itTrk)->trackId();
00136       if(ipi1==ipi2)  continue;
00137 
00138       if(!(*itTrk)->isMdcKalTrackValid()) continue;
00139       RecMdcKalTrack*  mdcKalTrk2 = (*itTrk)->mdcKalTrack();
00140       RecMdcKalTrack::setPidType(RecMdcKalTrack::pion);
00141 
00142       m_chargepi2=mdcKalTrk2->charge();
00143       if((m_chargepi1 + m_chargepi2) != 0) continue;
00144 
00146       HepVector a2 = mdcKalTrk2->getZHelix();
00147       HepSymMatrix Ea2 = mdcKalTrk2->getZError();
00148       VFHelix helixip3_2(point0,a2,Ea2);
00149       helixip3_2.pivot(IP);
00150       HepVector  vecipa2 = helixip3_2.a();
00151 
00152       double dr2 = fabs(vecipa2[0]);
00153       double dz2 = fabs(vecipa2[3]);
00154       double costheta2 = cos(mdcKalTrk2->theta());
00155       if (  dr2 >= 1.0) continue;
00156       if (  dz2 >= 10.0) continue; 
00157       if ( fabs(costheta2) >= 0.93) continue; 
00159  
00160       if(PID_flag == 5) {
00161         simple_pid->preparePID(*itTrk);
00162         if(simple_pid->probPion() < 0.0 || simple_pid->probPion() < simple_pid->probKaon()) continue;  
00163       } 
00165       WTrackParameter pim(xmass[2],mdcKalTrk2->getZHelix(),mdcKalTrk2->getZError() );
00166 
00167       for(int m = 0; m < nGam-1; m++) {
00168         if(iGam[m]==-1) continue;
00169         int iGam1 = iGam[m];
00170         RecEmcShower *g1Trk = (*(evtRecTrkCol->begin()+iGam[m]))->emcShower();
00171         double eraw1 = g1Trk->energy();
00172         double phi1 = g1Trk->phi();
00173         double the1 = g1Trk->theta();
00174         HepLorentzVector ptrkg1,ptrkg10,ptrkg12;
00175         ptrkg1.setPx(eraw1*sin(the1)*cos(phi1));
00176         ptrkg1.setPy(eraw1*sin(the1)*sin(phi1));
00177         ptrkg1.setPz(eraw1*cos(the1));
00178         ptrkg1.setE(eraw1);
00179         ptrkg10 = ptrkg1;
00180         ptrkg12 = ptrkg1.boost(-0.011,0,0);
00181 
00182         for(int n = m+1; n < nGam; n++) {
00183           if(iGam[n]==-1) continue;
00184           RecEmcShower *g2Trk = (*(evtRecTrkCol->begin()+iGam[n]))->emcShower();
00185           int iGam2 = iGam[n]; 
00186           double eraw2 = g2Trk->energy();
00187           double phi2 = g2Trk->phi();
00188           double the2 = g2Trk->theta();
00189           HepLorentzVector ptrkg2,ptrkg20,ptrkg22;
00190           ptrkg2.setPx(eraw2*sin(the2)*cos(phi2));
00191           ptrkg2.setPy(eraw2*sin(the2)*sin(phi2));
00192           ptrkg2.setPz(eraw2*cos(the2));
00193           ptrkg2.setE(eraw2);
00194           ptrkg20 = ptrkg2;
00195           ptrkg22 = ptrkg2.boost(-0.011,0,0);
00196 
00198           HepLorentzVector  ptrkpi0;
00199           ptrkpi0 = ptrkg12+ptrkg22;
00200           double m_xmpi0_tem = ptrkpi0.m();
00201           if(m_xmpi0_tem>0.150||m_xmpi0_tem<0.115)  continue;
00203           bool IsEndcap1 = false; bool IsEndcap2 = false;
00204           if(fabs(cos(the1)) > 0.86 && fabs(cos(the1)) < 0.92) IsEndcap1 = true; 
00205           if(fabs(cos(the2)) > 0.86 && fabs(cos(the2)) < 0.92) IsEndcap2 = true; 
00206           if(IsEndcap1 && IsEndcap2)  continue;
00208           KalmanKinematicFit * kmfit = KalmanKinematicFit::instance();
00209           kmfit->init();
00210           kmfit->setChisqCut(2500);
00211           kmfit->AddTrack(0, 0.0, g1Trk);
00212           kmfit->AddTrack(1, 0.0, g2Trk);
00213           kmfit->AddResonance(0, mpi0, 0, 1);
00214 
00215           kmfit->Fit(0);  // Perform fit
00216           kmfit->BuildVirtualParticle(0);
00217 
00218           double pi0_chisq = kmfit->chisq(0);
00219           if ( pi0_chisq >= 2500) continue;
00220           HepLorentzVector p2gfit = kmfit->pfit(0) + kmfit->pfit(1);
00221           p2gfit.boost(-0.011,0,0);
00222 
00224           HepPoint3D vx(xorigin.x(), xorigin.y(), xorigin.z());
00225           HepSymMatrix Evx(3, 0);
00226           double bx = 1E+6; Evx[0][0] = bx*bx;
00227           double by = 1E+6; Evx[1][1] = by*by;
00228           double bz = 1E+6; Evx[2][2] = bz*bz;
00229           VertexParameter vxpar; vxpar.setVx(vx); vxpar.setEvx(Evx);
00231 
00232           VertexFit* vtxfit = VertexFit::instance();
00233           vtxfit->init();
00234           vtxfit->AddTrack(0,  pip);
00235           vtxfit->AddTrack(1,  pim);
00236           vtxfit->AddVertex(0, vxpar, 0, 1);
00237           if(!vtxfit->Fit(0))  continue;
00238           vtxfit->Swim(0);
00239 
00240           WTrackParameter wpip = vtxfit->wtrk(0);
00241           WTrackParameter wpim = vtxfit->wtrk(1);
00242 
00243           HepVector pip_val = HepVector(7,0);
00244           HepVector pim_val = HepVector(7,0);
00245           pip_val = wpip.w();
00246           pim_val = wpim.w();
00247 
00248           HepLorentzVector P_PIP(pip_val[0],pip_val[1],pip_val[2],pip_val[3]);
00249           HepLorentzVector P_PIM(pim_val[0],pim_val[1],pim_val[2],pim_val[3]);
00250 
00251           P_PIM.boost(-0.011,0,0);
00252           P_PIP.boost(-0.011,0,0);
00253           HepLorentzVector ppipi = P_PIM + P_PIP;
00254           pddd = P_PIM + P_PIP + p2gfit;
00255 
00256           double   mpipi = ppipi.m();
00257           //if((mpipi-0.497)<0.020)  continue;
00258 
00259           double ppipipi0=pddd.rho();
00260 
00261           double temp1 = (ecms/2)*(ecms/2)-ppipipi0*ppipipi0 ;
00262           if(temp1<0) temp1 =0;
00263           double mass_bc_tem = sqrt(temp1);
00264           if(mass_bc_tem < 1.82 || mass_bc_tem > 1.89) continue;
00265 
00266           double  delE_tag_tag = ecms/2-pddd.e();
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 
00275             ipi1_temp=ipi1;
00276             ipi2_temp=ipi2;
00277             iGam1_temp = iGam1;
00278             iGam2_temp = iGam2;
00279 
00280             ncount1 = 1;
00281           }
00282         }
00283       }
00284     }
00285   }
00286 
00287   if(ncount1 == 1){
00288     tagmode=16;
00289     if(m_chargetag<0)  tagmode=-16;
00290     tagmd=tagmode;
00291     mass_bc  = mass_bcgg;
00292     delE_tag = delE_tag_temp;
00293     cqtm     = 0;  
00294 
00295     iGoodtag.push_back(ipi1_temp);
00296     iGoodtag.push_back(ipi2_temp);
00297     iGamtag.push_back(iGam1_temp);
00298     iGamtag.push_back(iGam2_temp);
00299     iGamtag.push_back(9999);
00300     iGamtag.push_back(9999);
00301 
00302     ptag = pddd_temp;
00303 
00304     pipipi0md=true;
00305   }
00306 }


Member Data Documentation

double Pipipi0::cqtm [private]

Definition at line 44 of file Pipipi0.h.

Referenced by GetCQtm(), and MTotal().

double Pipipi0::delE_tag [private]

Definition at line 45 of file Pipipi0.h.

Referenced by GetdelE_tag(), and MTotal().

Vint Pipipi0::iGamtag [private]

Definition at line 48 of file Pipipi0.h.

Referenced by GettagGam1(), and MTotal().

Vint Pipipi0::iGoodtag [private]

Definition at line 46 of file Pipipi0.h.

Referenced by Gettagtrk1(), and MTotal().

double Pipipi0::mass_bc [private]

Definition at line 43 of file Pipipi0.h.

Referenced by Getmass_bc(), and MTotal().

bool Pipipi0::pipipi0md [private]

Definition at line 41 of file Pipipi0.h.

Referenced by GetPipipi0md(), and MTotal().

HepLorentzVector Pipipi0::ptag [private]

Definition at line 47 of file Pipipi0.h.

Referenced by Gettagp1(), and MTotal().

double Pipipi0::tagmd [private]

Definition at line 42 of file Pipipi0.h.

Referenced by Gettagmd(), and MTotal().


Generated on Tue Nov 29 23:20:42 2016 for BOSS_7.0.2 by  doxygen 1.4.7