K0pipipi0 Class Reference

#include <K0pipipi0.h>

List of all members.

Public Member Functions

 K0pipipi0 ()
 ~K0pipipi0 ()
bool GetK0pipipi0md ()
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 k0pipipi0md
double tagmd
double mass_bc
double cqtm
double delE_tag
Vint iGoodtag
HepLorentzVector ptag
Vint iGamtag


Detailed Description

Definition at line 20 of file K0pipipi0.h.


Constructor & Destructor Documentation

K0pipipi0::K0pipipi0 (  ) 

Definition at line 27 of file K0pipipi0.cxx.

00028 {}

K0pipipi0::~K0pipipi0 (  ) 

Definition at line 30 of file K0pipipi0.cxx.

00031 {}


Member Function Documentation

double K0pipipi0::GetCQtm (  )  [inline]

Definition at line 30 of file K0pipipi0.h.

References cqtm.

Referenced by Sing::Mdset().

00030 { return cqtm; }

double K0pipipi0::GetdelE_tag (  )  [inline]

Definition at line 31 of file K0pipipi0.h.

References delE_tag.

Referenced by Sing::Mdset().

00031 { return delE_tag; }

bool K0pipipi0::GetK0pipipi0md (  )  [inline]

Definition at line 27 of file K0pipipi0.h.

References k0pipipi0md.

Referenced by Sing::Mdset().

00027 { return k0pipipi0md; }

double K0pipipi0::Getmass_bc (  )  [inline]

Definition at line 29 of file K0pipipi0.h.

References mass_bc.

Referenced by Sing::Mdset().

00029 { return mass_bc; }

Vint K0pipipi0::GettagGam1 (  )  [inline]

Definition at line 34 of file K0pipipi0.h.

References iGamtag.

Referenced by Sing::Mdset().

00034 { return iGamtag; }

double K0pipipi0::Gettagmd (  )  [inline]

Definition at line 28 of file K0pipipi0.h.

References tagmd.

Referenced by Sing::Mdset().

00028 { return tagmd; }

HepLorentzVector K0pipipi0::Gettagp1 (  )  [inline]

Definition at line 33 of file K0pipipi0.h.

References ptag.

Referenced by Sing::Mdset().

00033 { return ptag; }

Vint K0pipipi0::Gettagtrk1 (  )  [inline]

Definition at line 32 of file K0pipipi0.h.

References iGoodtag.

Referenced by Sing::Mdset().

00032 { return iGoodtag; }

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

Definition at line 34 of file K0pipipi0.cxx.

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

Referenced by Sing::Mdset().

00036 {
00037 
00038   int nGood=iGood.size();
00039   int nGam=iGam.size();
00040 
00041   iGoodtag.clear();
00042   iGamtag.clear();
00043   
00044   double mass_bcgg, delE_tag_temp;
00045   int m_chargetag,m_chargepi1,m_chargepi2,m_chargepi3,m_chargepi4;
00046   int ik1_temp,ipi1_temp,ipi2_temp,ipi3_temp,ipi4_temp, iGam1_temp, iGam2_temp;
00047   HepLorentzVector pddd;
00048   HepLorentzVector 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   k0pipipi0md=false;
00061   double  tagmode=0;
00062 
00063   if((evtRecEvent->totalCharged() < 4||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   HepSymMatrix xoriginEx(3,0);
00075   IVertexDbSvc*  vtxsvc;
00076   Gaudi::svcLocator()->service("VertexDbSvc", vtxsvc);
00077   if(vtxsvc->isVertexValid())
00078   {
00079     double* dbv = vtxsvc->PrimaryVertex();
00080     double*  vv = vtxsvc->SigmaPrimaryVertex();
00081     xorigin.setX(dbv[0]);
00082     xorigin.setY(dbv[1]);
00083     xorigin.setZ(dbv[2]);
00084 
00085     xoriginEx[0][0] = vv[0] * vv[0];
00086     xoriginEx[1][1] = vv[1] * vv[1];
00087     xoriginEx[2][2] = vv[2] * vv[2];
00088 
00089   }
00090 
00091   double xv=xorigin.x();
00092   double yv=xorigin.y();
00093   double zv=xorigin.z();
00094 
00095   HepPoint3D point0(0.,0.,0.);
00096   HepPoint3D IP(xorigin[0],xorigin[1],xorigin[2]);
00098   
00099   HepLorentzVector p2gfit;
00100   HepLorentzVector p2gg;
00101 
00102   for(int i = 0; i < evtRecEvent->totalCharged(); i++) {
00103     EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + i;
00104 
00105     int ipi1= (*itTrk)->trackId();
00106 
00107     if(!(*itTrk)->isMdcKalTrackValid()) continue;
00108     RecMdcKalTrack*  mdcKalTrk1 = (*itTrk)->mdcKalTrack();
00109     RecMdcKalTrack::setPidType(RecMdcKalTrack::pion);
00110 
00111     m_chargepi1=mdcKalTrk1->charge();
00112     if(m_chargepi1 != 1) continue;
00113 
00115     HepVector a1 = mdcKalTrk1->getZHelix();
00116     HepSymMatrix Ea1 = mdcKalTrk1->getZError();
00117 
00118     VFHelix helixip3_1(point0,a1,Ea1);
00119     helixip3_1.pivot(IP);
00120     HepVector  vecipa1 = helixip3_1.a();
00121 
00122     double dr1 = fabs(vecipa1[0]);
00123     double dz1 = fabs(vecipa1[3]);
00124     double costheta1 = cos(mdcKalTrk1->theta());
00125 
00126     if (  dr1 >= 15.0) continue;
00127     if (  dz1 >= 25.0) continue; 
00128     if ( fabs(costheta1) >= 0.93) continue; 
00130     WTrackParameter pip1(xmass[2],mdcKalTrk1->getZHelix(),mdcKalTrk1->getZError() );
00131 
00132     //  
00133     // select pi2
00134     //  
00135 
00136     for(int j = 0; j < evtRecEvent->totalCharged();j++) {
00137       EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + j;
00138 
00139       int ipi2= (*itTrk)->trackId();
00140       if(ipi1==ipi2)  continue;
00141 
00142       if(!(*itTrk)->isMdcKalTrackValid()) continue;
00143       RecMdcKalTrack*  mdcKalTrk2 = (*itTrk)->mdcKalTrack();
00144       RecMdcKalTrack::setPidType(RecMdcKalTrack::pion);
00145 
00146       m_chargepi2=mdcKalTrk2->charge();
00147       if((m_chargepi2 + m_chargepi1) != 0) continue;
00148 
00150       HepVector a2 = mdcKalTrk2->getZHelix();
00151       HepSymMatrix Ea2 = mdcKalTrk2->getZError();
00152       VFHelix helixip3_2(point0,a2,Ea2);
00153       helixip3_2.pivot(IP);
00154       HepVector  vecipa2 = helixip3_2.a();
00155 
00156       double dr2 = fabs(vecipa2[0]);
00157       double dz2 = fabs(vecipa2[3]);
00158       double costheta2 = cos(mdcKalTrk2->theta());
00159       if (  dr2 >= 15.0) continue;
00160       if (  dz2 >= 25.0) continue; 
00161       if ( fabs(costheta2) >= 0.93) continue; 
00163       WTrackParameter pim1(xmass[2],mdcKalTrk2->getZHelix(),mdcKalTrk2->getZError() );
00164 
00165 
00166       HepVector pip1_val = HepVector(7,0);
00167       HepVector pim1_val = HepVector(7,0);
00168       pip1_val = pip1.w(); 
00169       pim1_val = pim1.w();
00170       HepLorentzVector ptrktagk0(pip1_val[0]+pim1_val[0],pip1_val[1]+pim1_val[1],pip1_val[2]+pim1_val[2],pip1_val[3]+pim1_val[3]);
00171       double m_xmtagk0_tem = ptrktagk0.mag();
00172       if(fabs(ptrktagk0.m()-0.498)>0.1) continue;
00173 
00174       HepPoint3D vx(xorigin.x(), xorigin.y(), xorigin.z());
00175       HepSymMatrix Evx(3, 0);
00176       double bx = 1E+6; Evx[0][0] = bx*bx;
00177       double by = 1E+6; Evx[1][1] = by*by;
00178       double bz = 1E+6; Evx[2][2] = bz*bz;
00179       VertexParameter vxpar; vxpar.setVx(vx); vxpar.setEvx(Evx);
00180 
00181 
00182       VertexFit *vtxfit0 = VertexFit::instance();
00183       vtxfit0->init();
00184       vtxfit0->AddTrack(0, pip1);
00185       vtxfit0->AddTrack(1, pim1);
00186       vtxfit0->AddVertex(0, vxpar, 0, 1);
00187       if(!(vtxfit0->Fit(0))) continue;
00188       vtxfit0->Swim(0); 
00189       vtxfit0->BuildVirtualParticle(0);
00190       WTrackParameter wksp = vtxfit0->wtrk(0);
00191       WTrackParameter wksm = vtxfit0->wtrk(1);
00192       WTrackParameter wks_Trk = vtxfit0->wVirtualTrack(0);
00193       VertexParameter wks_var = vtxfit0->vpar(0);     
00194 
00195       //        
00196       //select pi3
00197       //        
00198       for(int k = 0; k < evtRecEvent->totalCharged(); k++) {
00199         EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + k;
00200 
00201         int ipi3= (*itTrk)->trackId();
00202         if(ipi2==ipi3 || ipi1==ipi3)  continue;
00203 
00204         if(!(*itTrk)->isMdcKalTrackValid()) continue;
00205         RecMdcKalTrack*  mdcKalTrk3 = (*itTrk)->mdcKalTrack();
00206         RecMdcKalTrack::setPidType(RecMdcKalTrack::pion);
00207 
00208         m_chargepi3=mdcKalTrk3->charge();
00209         if(abs(m_chargepi3) != 1) continue;
00210 
00212         HepVector a3 = mdcKalTrk3->getZHelix();
00213         HepSymMatrix Ea3 = mdcKalTrk3->getZError();
00214         VFHelix helixip3_3(point0,a3,Ea3);
00215         helixip3_3.pivot(IP);
00216         HepVector  vecipa3 = helixip3_3.a();
00217 
00218         double dr3 = fabs(vecipa3[0]);
00219         double dz3 = fabs(vecipa3[3]);
00220         double costheta3 = cos(mdcKalTrk3->theta());
00221         if (  dr3 >= 1.0) continue;
00222         if (  dz3 >= 10.0) continue; 
00223         if ( fabs(costheta3) >= 0.93) continue; 
00225         if(PID_flag == 5) {
00226           simple_pid->preparePID(*itTrk);
00227           if(simple_pid->probPion() < 0.0 || simple_pid->probPion() < simple_pid->probKaon()) continue;  
00228         } 
00230         WTrackParameter pip2(xmass[2],mdcKalTrk3->getZHelix(),mdcKalTrk3->getZError() );
00231 
00232         //
00233         //select pi4
00234         //
00235         for(int l = 0; l < evtRecEvent->totalCharged(); l++) {
00236           EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + l;
00237 
00238           int ipi4= (*itTrk)->trackId();
00239           if(ipi4==ipi3 || ipi4==ipi2 || ipi4 ==ipi1)  continue;
00240 
00241           if(!(*itTrk)->isMdcKalTrackValid()) continue;
00242           RecMdcKalTrack*  mdcKalTrk4 = (*itTrk)->mdcKalTrack();
00243           RecMdcKalTrack::setPidType(RecMdcKalTrack::pion);
00244 
00245           m_chargepi4=mdcKalTrk4->charge();
00246           if((m_chargepi4 + m_chargepi3) != 0) continue;
00247 
00249           HepVector a4 = mdcKalTrk4->getZHelix();
00250           HepSymMatrix Ea4 = mdcKalTrk4->getZError();
00251           VFHelix helixip3_4(point0,a4,Ea4);
00252           helixip3_4.pivot(IP);
00253           HepVector  vecipa4 = helixip3_4.a();
00254 
00255           double dr4 = fabs(vecipa4[0]);
00256           double dz4 = fabs(vecipa4[3]);
00257           double costheta4 = cos(mdcKalTrk4->theta());
00258           if (  dr4 >= 1.0) continue;
00259           if (  dz4 >= 10.0) continue; 
00260           if ( fabs(costheta4) >= 0.93) continue; 
00262           if(PID_flag == 5) {
00263             simple_pid->preparePID(*itTrk);
00264             if(simple_pid->probPion() < 0.0 || simple_pid->probPion() < simple_pid->probKaon()) continue;  
00265           } 
00267           WTrackParameter pim2(xmass[2],mdcKalTrk4->getZHelix(),mdcKalTrk4->getZError() );
00268 
00269           for(int m = 0; m < nGam-1; m++) {
00270             if(iGam[m]==-1) continue;
00271             RecEmcShower *g1Trk = (*(evtRecTrkCol->begin()+iGam[m]))->emcShower();
00272             double eraw1 = g1Trk->energy();
00273             double phi1 = g1Trk->phi();
00274             double the1 = g1Trk->theta();
00275             HepLorentzVector ptrkg1,ptrkg10,ptrkg12;
00276             ptrkg1.setPx(eraw1*sin(the1)*cos(phi1));
00277             ptrkg1.setPy(eraw1*sin(the1)*sin(phi1));
00278             ptrkg1.setPz(eraw1*cos(the1));
00279             ptrkg1.setE(eraw1);
00280             ptrkg10 = ptrkg1;
00281             ptrkg12 = ptrkg1.boost(-0.011,0,0);
00282 
00283             for(int n = m+1; n < nGam; n++) {
00284               if(iGam[n]==-1) continue;
00285               RecEmcShower *g2Trk = (*(evtRecTrkCol->begin()+iGam[n]))->emcShower();
00286               double eraw2 = g2Trk->energy();
00287               double phi2 = g2Trk->phi();
00288               double the2 = g2Trk->theta();
00289               HepLorentzVector ptrkg2,ptrkg20,ptrkg22;
00290               ptrkg2.setPx(eraw2*sin(the2)*cos(phi2));
00291               ptrkg2.setPy(eraw2*sin(the2)*sin(phi2));
00292               ptrkg2.setPz(eraw2*cos(the2));
00293               ptrkg2.setE(eraw2);
00294               ptrkg20 = ptrkg2;
00295               ptrkg22 = ptrkg2.boost(-0.011,0,0);
00296 
00298               HepLorentzVector  ptrkpi0;
00299               ptrkpi0 = ptrkg12+ptrkg22;
00300               double m_xmpi0_tem = ptrkpi0.m();
00301               if(m_xmpi0_tem>0.150||m_xmpi0_tem<0.115)  continue;
00303               bool IsEndcap1 = false; bool IsEndcap2 = false;
00304               if(fabs(cos(the1)) > 0.86 && fabs(cos(the1)) < 0.92) IsEndcap1 = true; 
00305               if(fabs(cos(the2)) > 0.86 && fabs(cos(the2)) < 0.92) IsEndcap2 = true; 
00306               if(IsEndcap1 && IsEndcap2)  continue;
00308               KalmanKinematicFit * kmfit = KalmanKinematicFit::instance();
00309               kmfit->init();
00310               kmfit->setChisqCut(2500);
00311               kmfit->AddTrack(0, 0.0, g1Trk);
00312               kmfit->AddTrack(1, 0.0, g2Trk);
00313               kmfit->AddResonance(0, mpi0, 0, 1);
00314 
00315               kmfit->Fit(0);  // Perform fit
00316               kmfit->BuildVirtualParticle(0);
00317 
00318               double pi0_chisq = kmfit->chisq(0);
00319               if ( pi0_chisq >= 2500) continue;
00320               HepLorentzVector p2gfit = kmfit->pfit(0) + kmfit->pfit(1);
00321               p2gfit.boost(-0.011,0,0);
00323 
00324               VertexFit* vtxfit_2 = VertexFit::instance();
00325               vtxfit_2->init();
00326               vtxfit_2->AddTrack(0,  pip2);
00327               vtxfit_2->AddTrack(1,  pim2);
00328               vtxfit_2->AddVertex(0, vxpar, 0, 1);
00329               if(!vtxfit_2->Fit(0))  continue;
00330               vtxfit_2->Swim(0);
00331 
00332               WTrackParameter wpip2 = vtxfit_2->wtrk(0);
00333               WTrackParameter wpim2 = vtxfit_2->wtrk(1);
00334 
00335               SecondVertexFit *vtxfit = SecondVertexFit::instance();
00336               vtxfit->init();
00337               vxpar.setEvx(xoriginEx);
00338               vtxfit->setPrimaryVertex(vxpar);
00339               vtxfit->AddTrack(0, wks_Trk);
00340               vtxfit->setVpar(wks_var);
00341               if(!vtxfit->Fit()) continue;
00342 
00343               if(vtxfit->chisq() >999.) continue;
00344               if(vtxfit->decayLength()<0.0) continue;
00345 
00346               double m_massks1_tem = vtxfit->p4par().m();
00347               if(m_massks1_tem < 0.485 || m_massks1_tem > 0.515)  continue;
00348               HepLorentzVector  p4kstag =  vtxfit->p4par();
00349               WTrackParameter para_ks = vtxfit0->wVirtualTrack(0);
00350 
00351               HepVector pip2_val = HepVector(7,0);
00352               HepVector pim2_val = HepVector(7,0);
00353               HepVector ksp_val = HepVector(7,0);
00354               HepVector ksm_val = HepVector(7,0);
00355 
00356               pip2_val = wpip2.w();
00357               pim2_val = wpim2.w();
00358               ksp_val = wksp.w();
00359               ksm_val = wksm.w();
00360 
00361               HepLorentzVector P_PIP2(pip2_val[0],pip2_val[1],pip2_val[2],pip2_val[3]);
00362               HepLorentzVector P_PIM2(pim2_val[0],pim2_val[1],pim2_val[2],pim2_val[3]);
00363               HepLorentzVector P_KSP(ksp_val[0],ksp_val[1],ksp_val[2],ksp_val[3]);
00364               HepLorentzVector P_KSM(ksm_val[0],ksm_val[1],ksm_val[2],ksm_val[3]);
00365 
00366               p4kstag.boost(-0.011,0,0);
00367               P_PIP2.boost(-0.011,0,0);
00368               P_PIM2.boost(-0.011,0,0);
00369               P_KSP.boost(-0.011,0,0);
00370               P_KSM.boost(-0.011,0,0);
00371 
00372               //pddd = P_PIP2 + P_PIM2 + P_KSP + P_KSM + p2gfit;
00373               pddd = P_PIP2 + P_PIM2 + p4kstag + p2gfit;
00374 
00375               double pk0pipipi0 = pddd.rho();
00376 
00377               double temp1 = (ecms/2)*(ecms/2)-pk0pipipi0*pk0pipipi0 ;
00378               if(temp1<0) temp1 =0;
00379               double mass_bc_tem  = sqrt(temp1);
00380               if(mass_bc_tem <1.82 || mass_bc_tem > 1.89) continue;
00381 
00382               double    delE_tag_tag = ecms/2-pddd.e();
00383 
00384               if(fabs(delE_tag_tag)<deltaE_tem) {
00385                 deltaE_tem = fabs(delE_tag_tag);
00386                 delE_tag_temp = delE_tag_tag;
00387                 mass_bcgg = mass_bc_tem;
00388 
00389                 pddd_temp = pddd;
00390 
00391                 ipi1_temp=ipi1;
00392                 ipi2_temp=ipi2;
00393                 ipi3_temp=ipi3;
00394                 ipi4_temp=ipi4;
00395                 iGam1_temp = iGam[m];
00396                 iGam2_temp = iGam[n];
00397 
00398                 ncount1 = 1;
00399 
00400               }   
00401             }       
00402           }
00403         }
00404       }
00405     }
00406   }
00407 
00408 
00409   if(ncount1 == 1){  
00410     tagmode=17;
00411     if(m_chargetag<0)  tagmode=-17;
00412     tagmd=tagmode;
00413     mass_bc = mass_bcgg;
00414     delE_tag = delE_tag_temp;  
00415     cqtm     = 0;
00416 
00417     iGoodtag.push_back(ipi1_temp);
00418     iGoodtag.push_back(ipi2_temp);
00419     iGoodtag.push_back(ipi3_temp);
00420     iGoodtag.push_back(ipi4_temp);
00421 
00422     iGamtag.push_back(iGam1_temp);
00423     iGamtag.push_back(iGam2_temp);
00424     iGamtag.push_back(9999);
00425     iGamtag.push_back(9999);
00426 
00427     ptag = pddd_temp;
00428 
00429     k0pipipi0md = true;
00430   }
00431 
00432 }


Member Data Documentation

double K0pipipi0::cqtm [private]

Definition at line 43 of file K0pipipi0.h.

Referenced by GetCQtm(), and MTotal().

double K0pipipi0::delE_tag [private]

Definition at line 44 of file K0pipipi0.h.

Referenced by GetdelE_tag(), and MTotal().

Vint K0pipipi0::iGamtag [private]

Definition at line 47 of file K0pipipi0.h.

Referenced by GettagGam1(), and MTotal().

Vint K0pipipi0::iGoodtag [private]

Definition at line 45 of file K0pipipi0.h.

Referenced by Gettagtrk1(), and MTotal().

bool K0pipipi0::k0pipipi0md [private]

Definition at line 40 of file K0pipipi0.h.

Referenced by GetK0pipipi0md(), and MTotal().

double K0pipipi0::mass_bc [private]

Definition at line 42 of file K0pipipi0.h.

Referenced by Getmass_bc(), and MTotal().

HepLorentzVector K0pipipi0::ptag [private]

Definition at line 46 of file K0pipipi0.h.

Referenced by Gettagp1(), and MTotal().

double K0pipipi0::tagmd [private]

Definition at line 41 of file K0pipipi0.h.

Referenced by Gettagmd(), and MTotal().


Generated on Tue Nov 29 23:19:52 2016 for BOSS_7.0.2 by  doxygen 1.4.7