BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
Rhopi Class Reference

#include <Rhopi.h>

Inheritance diagram for Rhopi:

Public Member Functions

 Rhopi (const std::string &name, ISvcLocator *pSvcLocator)
StatusCode initialize ()
StatusCode execute ()
StatusCode finalize ()

Detailed Description

Definition at line 7 of file Rhopi.h.

Constructor & Destructor Documentation

◆ Rhopi()

Rhopi::Rhopi ( const std::string & name,
ISvcLocator * pSvcLocator )

Definition at line 54 of file Rhopi.cxx.

55 : Algorithm( name, pSvcLocator ) {
56
57 // Declare the properties
58 declareProperty( "Vr0cut", m_vr0cut = 1.0 );
59 declareProperty( "Vz0cut", m_vz0cut = 5.0 );
60 declareProperty( "EnergyThreshold", m_energyThreshold = 0.05 );
61 declareProperty( "GammaPhiCut", m_gammaPhiCut = 20.0 );
62 declareProperty( "GammaThetaCut", m_gammaThetaCut = 20.0 );
63 declareProperty( "GammaAngleCut", m_gammaAngleCut = 10.0 );
64 declareProperty( "Test4C", m_test4C = 1 );
65 declareProperty( "Test5C", m_test5C = 1 );
66 declareProperty( "CheckDedx", m_checkDedx = 1 );
67 declareProperty( "CheckTof", m_checkTof = 1 );
68}

Referenced by Rhopi().

Member Function Documentation

◆ execute()

StatusCode Rhopi::execute ( )

Definition at line 338 of file Rhopi.cxx.

338 {
339
340 // std::cout << "execute()" << std::endl;
341
342 MsgStream log( msgSvc(), name() );
343 log << MSG::INFO << "in execute()" << endmsg;
344
345 SmartDataPtr<Event::EventHeader> eventHeader( eventSvc(), "/Event/EventHeader" );
346 int runNo = eventHeader->runNumber();
347 int event = eventHeader->eventNumber();
348 log << MSG::DEBUG << "run, evtnum = " << runNo << " , " << event << endmsg;
349 cout << "event " << event << endl;
350 Ncut0++;
351
352 SmartDataPtr<EvtRecEvent> evtRecEvent( eventSvc(), EventModel::EvtRec::EvtRecEvent );
353 // log << MSG::INFO << "get event tag OK" << endmsg;
354 log << MSG::DEBUG << "ncharg, nneu, tottks = " << evtRecEvent->totalCharged() << " , "
355 << evtRecEvent->totalNeutral() << " , " << evtRecEvent->totalTracks() << endmsg;
356
357 SmartDataPtr<EvtRecTrackCol> evtRecTrkCol( eventSvc(), EventModel::EvtRec::EvtRecTrackCol );
358 //
359 // check x0, y0, z0, r0
360 // suggest cut: |z0|<5 && r0<1
361 //
362 Vint iGood, ipip, ipim;
363 iGood.clear();
364 ipip.clear();
365 ipim.clear();
366 Vp4 ppip, ppim;
367 ppip.clear();
368 ppim.clear();
369
370 int nCharge = 0;
371
372 Hep3Vector xorigin( 0, 0, 0 );
373
374 /*VertexDbSvc* vtxsvc;
375 Gaudi::svcLocator()->service("VertexDbSvc", vtxsvc);
376 if(vtxsvc->isVertexValid()){
377 double* dbv = vtxsvc->PrimaryVertex();
378 double* vv = vtxsvc->SigmaPrimaryVertex();
379 xorigin.setX(dbv[0]);
380 xorigin.setY(dbv[1]);
381 xorigin.setZ(dbv[2]);
382 }*/
383 SmartIF<IVertexDbSvc> m_vtxSvc;
384 m_vtxSvc = serviceLocator()->service( "VertexDbSvc" );
385 if ( m_vtxSvc->isVertexValid() )
386 {
387 // if(m_vtxSvc){
388 double* dbv = m_vtxSvc->PrimaryVertex();
389 double* vv = m_vtxSvc->SigmaPrimaryVertex();
390 xorigin.setX( dbv[0] );
391 xorigin.setY( dbv[1] );
392 xorigin.setZ( dbv[2] );
393 }
394
395 for ( int i = 0; i < evtRecEvent->totalCharged(); i++ )
396 {
397 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + i;
398 if ( !( *itTrk )->isMdcTrackValid() ) continue;
399 RecMdcTrack* mdcTrk = ( *itTrk )->mdcTrack();
400 double pch = mdcTrk->p();
401 double x0 = mdcTrk->x();
402 double y0 = mdcTrk->y();
403 double z0 = mdcTrk->z();
404 double theta0 = mdcTrk->theta();
405 double phi0 = mdcTrk->helix( 1 );
406 double xv = xorigin.x();
407 double yv = xorigin.y();
408 double Rxy = ( x0 - xv ) * cos( phi0 ) + ( y0 - yv ) * sin( phi0 );
409 m_vx0 = x0;
410 m_vy0 = y0;
411 m_vz0 = z0;
412 m_vr0 = Rxy;
413
414 HepVector a = mdcTrk->helix();
415 HepSymMatrix Ea = mdcTrk->err();
416 HepPoint3D point0( 0., 0., 0. ); // the initial point for MDC recosntruction
417 HepPoint3D IP( xorigin[0], xorigin[1], xorigin[2] );
418 VFHelix helixip( point0, a, Ea );
419 helixip.pivot( IP );
420 HepVector vecipa = helixip.a();
421 double Rvxy0 = fabs( vecipa[0] ); // the nearest distance to IP in xy plane
422 double Rvz0 = vecipa[3]; // the nearest distance to IP in z direction
423 double Rvphi0 = vecipa[1];
424 m_rvxy0 = Rvxy0;
425 m_rvz0 = Rvz0;
426 m_rvphi0 = Rvphi0;
427
428 m_tuple1->write();
429 // if(fabs(z0) >= m_vz0cut) continue;
430 // if(fabs(Rxy) >= m_vr0cut) continue;
431
432 if ( fabs( Rvz0 ) >= 10.0 ) continue;
433 if ( fabs( Rvxy0 ) >= 1.0 ) continue;
434 if ( fabs( cos( theta0 ) ) >= 0.93 ) continue;
435
436 iGood.push_back( i );
437 nCharge += mdcTrk->charge();
438 }
439
440 //
441 // Finish Good Charged Track Selection
442 //
443 int nGood = iGood.size();
444 log << MSG::DEBUG << "ngood, totcharge = " << nGood << " , " << nCharge << endmsg;
445 if ( ( nGood != 2 ) || ( nCharge != 0 ) ) { return StatusCode::SUCCESS; }
446 Ncut1++;
447
448 Vint iGam;
449 iGam.clear();
450 for ( int i = evtRecEvent->totalCharged(); i < evtRecEvent->totalTracks(); i++ )
451 {
452 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + i;
453 if ( !( *itTrk )->isEmcShowerValid() ) continue;
454 RecEmcShower* emcTrk = ( *itTrk )->emcShower();
455 Hep3Vector emcpos( emcTrk->x(), emcTrk->y(), emcTrk->z() );
456 // find the nearest charged track
457 double dthe = 200.;
458 double dphi = 200.;
459 double dang = 200.;
460 for ( int j = 0; j < evtRecEvent->totalCharged(); j++ )
461 {
462 EvtRecTrackIterator jtTrk = evtRecTrkCol->begin() + j;
463 if ( !( *jtTrk )->isExtTrackValid() ) continue;
464 RecExtTrack* extTrk = ( *jtTrk )->extTrack();
465 if ( extTrk->emcVolumeNumber() == -1 ) continue;
466 Hep3Vector extpos = extTrk->emcPosition();
467 // double ctht = extpos.cosTheta(emcpos);
468 double angd = extpos.angle( emcpos );
469 double thed = extpos.theta() - emcpos.theta();
470 double phid = extpos.deltaPhi( emcpos );
471 thed = fmod( thed + CLHEP::twopi + CLHEP::twopi + pi, CLHEP::twopi ) - CLHEP::pi;
472 phid = fmod( phid + CLHEP::twopi + CLHEP::twopi + pi, CLHEP::twopi ) - CLHEP::pi;
473 if ( angd < dang )
474 {
475 dang = angd;
476 dthe = thed;
477 dphi = phid;
478 }
479 }
480 if ( dang >= 200 ) continue;
481 double eraw = emcTrk->energy();
482 double theta1 = emcTrk->theta();
483 double time = emcTrk->time();
484 dthe = dthe * 180 / ( CLHEP::pi );
485 dphi = dphi * 180 / ( CLHEP::pi );
486 dang = dang * 180 / ( CLHEP::pi );
487 m_dthe = dthe;
488 m_dphi = dphi;
489 m_dang = dang;
490 m_eraw = eraw;
491 m_tuple2->write();
492
493 if ( fabs( cos( theta1 ) ) < 0.80 )
494 {
495 if ( eraw <= ( m_energyThreshold / 2 ) ) continue;
496 }
497 else if ( fabs( cos( theta1 ) ) > 0.86 && fabs( cos( theta1 ) ) < 0.92 )
498 {
499 if ( eraw <= m_energyThreshold ) continue;
500 }
501 else continue;
502
503 // if((fabs(dthe) < m_gammaThetaCut) && (fabs(dphi)<m_gammaPhiCut) ) continue;
504 if ( fabs( dang ) < m_gammaAngleCut ) continue;
505
506 if ( time < 0 || time > 14 ) continue;
507
508 //
509 // good photon cut will be set here
510 //
511 iGam.push_back( i );
512 }
513
514 //
515 // Finish Good Photon Selection
516 //
517 int nGam = iGam.size();
518
519 log << MSG::DEBUG << "num Good Photon " << nGam << " , " << evtRecEvent->totalNeutral()
520 << endmsg;
521 if ( nGam < 2 ) { return StatusCode::SUCCESS; }
522 Ncut2++;
523
524 //
525 //
526 // check dedx infomation
527 //
528 //
529
530 if ( m_checkDedx == 1 )
531 {
532 for ( int i = 0; i < nGood; i++ )
533 {
534 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
535 if ( !( *itTrk )->isMdcTrackValid() ) continue;
536 if ( !( *itTrk )->isMdcDedxValid() ) continue;
537 RecMdcTrack* mdcTrk = ( *itTrk )->mdcTrack();
538 RecMdcDedx* dedxTrk = ( *itTrk )->mdcDedx();
539 m_ptrk = mdcTrk->p();
540
541 m_chie = dedxTrk->chiE();
542 m_chimu = dedxTrk->chiMu();
543 m_chipi = dedxTrk->chiPi();
544 m_chik = dedxTrk->chiK();
545 m_chip = dedxTrk->chiP();
546 m_ghit = dedxTrk->numGoodHits();
547 m_thit = dedxTrk->numTotalHits();
548 m_probPH = dedxTrk->probPH();
549 m_normPH = dedxTrk->normPH();
550 m_tuple7->write();
551 }
552 }
553
554 //
555 // check TOF infomation
556 //
557
558 if ( m_checkTof == 1 )
559 {
560 for ( int i = 0; i < nGood; i++ )
561 {
562 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
563 if ( !( *itTrk )->isMdcTrackValid() ) continue;
564 if ( !( *itTrk )->isTofTrackValid() ) continue;
565
566 RecMdcTrack* mdcTrk = ( *itTrk )->mdcTrack();
567 SmartRefVector<RecTofTrack> tofTrkCol = ( *itTrk )->tofTrack();
568
569 double ptrk = mdcTrk->p();
570
571 SmartRefVector<RecTofTrack>::iterator iter_tof = tofTrkCol.begin();
572 for ( ; iter_tof != tofTrkCol.end(); iter_tof++ )
573 {
574 TofHitStatus* status = new TofHitStatus;
575 status->setStatus( ( *iter_tof )->status() );
576 if ( !( status->is_barrel() ) )
577 { // endcap
578 if ( !( status->is_counter() ) ) continue; // ?
579 if ( status->layer() != 0 ) continue; // layer1
580 double path = ( *iter_tof )->path(); // ?
581 double tof = ( *iter_tof )->tof();
582 double ph = ( *iter_tof )->ph();
583 double rhit = ( *iter_tof )->zrhit();
584 double qual = 0.0 + ( *iter_tof )->quality();
585 double cntr = 0.0 + ( *iter_tof )->tofID();
586 double texp[5];
587 for ( int j = 0; j < 5; j++ )
588 {
589 double gb = ptrk / xmass[j];
590 double beta = gb / sqrt( 1 + gb * gb );
591 texp[j] = 10 * path / beta / velc;
592 }
593 m_cntr_etof = cntr;
594 m_ptot_etof = ptrk;
595 m_ph_etof = ph;
596 m_rhit_etof = rhit;
597 m_qual_etof = qual;
598 m_te_etof = tof - texp[0];
599 m_tmu_etof = tof - texp[1];
600 m_tpi_etof = tof - texp[2];
601 m_tk_etof = tof - texp[3];
602 m_tp_etof = tof - texp[4];
603 m_tuple8->write();
604 }
605 else
606 { // barrel
607 if ( !( status->is_counter() ) ) continue; // ?
608 if ( status->layer() == 1 )
609 { // layer1
610 double path = ( *iter_tof )->path(); // ?
611 double tof = ( *iter_tof )->tof();
612 double ph = ( *iter_tof )->ph();
613 double rhit = ( *iter_tof )->zrhit();
614 double qual = 0.0 + ( *iter_tof )->quality();
615 double cntr = 0.0 + ( *iter_tof )->tofID();
616 double texp[5];
617 for ( int j = 0; j < 5; j++ )
618 {
619 double gb = ptrk / xmass[j];
620 double beta = gb / sqrt( 1 + gb * gb );
621 texp[j] = 10 * path / beta / velc;
622 }
623
624 m_cntr_btof1 = cntr;
625 m_ptot_btof1 = ptrk;
626 m_ph_btof1 = ph;
627 m_zhit_btof1 = rhit;
628 m_qual_btof1 = qual;
629 m_te_btof1 = tof - texp[0];
630 m_tmu_btof1 = tof - texp[1];
631 m_tpi_btof1 = tof - texp[2];
632 m_tk_btof1 = tof - texp[3];
633 m_tp_btof1 = tof - texp[4];
634 m_tuple9->write();
635 }
636
637 if ( status->layer() == 2 )
638 { // layer2
639 double path = ( *iter_tof )->path(); // ?
640 double tof = ( *iter_tof )->tof();
641 double ph = ( *iter_tof )->ph();
642 double rhit = ( *iter_tof )->zrhit();
643 double qual = 0.0 + ( *iter_tof )->quality();
644 double cntr = 0.0 + ( *iter_tof )->tofID();
645 double texp[5];
646 for ( int j = 0; j < 5; j++ )
647 {
648 double gb = ptrk / xmass[j];
649 double beta = gb / sqrt( 1 + gb * gb );
650 texp[j] = 10 * path / beta / velc;
651 }
652
653 m_cntr_btof2 = cntr;
654 m_ptot_btof2 = ptrk;
655 m_ph_btof2 = ph;
656 m_zhit_btof2 = rhit;
657 m_qual_btof2 = qual;
658 m_te_btof2 = tof - texp[0];
659 m_tmu_btof2 = tof - texp[1];
660 m_tpi_btof2 = tof - texp[2];
661 m_tk_btof2 = tof - texp[3];
662 m_tp_btof2 = tof - texp[4];
663 m_tuple10->write();
664 }
665 }
666
667 delete status;
668 }
669 } // loop all charged track
670 } // check tof
671
672 //
673 // Assign 4-momentum to each photon
674 //
675
676 Vp4 pGam;
677 pGam.clear();
678 for ( int i = 0; i < nGam; i++ )
679 {
680 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGam[i];
681 RecEmcShower* emcTrk = ( *itTrk )->emcShower();
682 double eraw = emcTrk->energy();
683 double phi = emcTrk->phi();
684 double the = emcTrk->theta();
685 HepLorentzVector ptrk;
686 ptrk.setPx( eraw * sin( the ) * cos( phi ) );
687 ptrk.setPy( eraw * sin( the ) * sin( phi ) );
688 ptrk.setPz( eraw * cos( the ) );
689 ptrk.setE( eraw );
690
691 // ptrk = ptrk.boost(-0.011,0,0);// boost to cms
692
693 pGam.push_back( ptrk );
694 }
695 // cout<<"before pid"<<endl;
696 //
697 // Assign 4-momentum to each charged track
698 //
699 ParticleID* pid = ParticleID::instance();
700 for ( int i = 0; i < nGood; i++ )
701 {
702 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
703 // if(pid) delete pid;
704 pid->init();
705 pid->setMethod( pid->methodProbability() );
706 // pid->setMethod(pid->methodLikelihood()); //for Likelihood Method
707
708 pid->setChiMinCut( 4 );
709 pid->setRecTrack( *itTrk );
710 // pid->usePidSys(pid->useDedx() | pid->useTof1() | pid->useTof2() | pid->useTofE()); //
711 // for BOSS version < 7.0.0
712 pid->usePidSys( pid->useDedx() | pid->useTofCorr() ); // use PID sub-system
713 pid->identify( pid->onlyPionKaonProton() ); // seperater Pion/Kaon/Proton
714 // pid->identify(pid->onlyPion());
715 // pid->identify(pid->onlyKaon());
716 pid->calculate();
717 if ( !( pid->IsPidInfoValid() ) ) continue;
718 RecMdcTrack* mdcTrk = ( *itTrk )->mdcTrack();
719 m_ptrk_pid = mdcTrk->p();
720 m_cost_pid = cos( mdcTrk->theta() );
721 m_dedx_pid = pid->chiDedx( 2 );
722 m_tof1_pid = pid->chiTof1( 2 );
723 m_tof2_pid = pid->chiTof2( 2 );
724 m_prob_pid = pid->probPion();
725 m_tuple11->write();
726
727 if ( ( pid->probPion() < pid->probProton() ) || ( pid->probPion() < pid->probKaon() ) )
728 continue;
729 // if(pid->probPion() < 0.001) continue;
730 // if(pid->pdf(2)<pid->pdf(3)) continue; // for Likelihood Method(0=electron 1=muon
731 // 2=pion 3=kaon 4=proton)
732
733 RecMdcKalTrack* mdcKalTrk =
734 ( *itTrk )->mdcKalTrack(); // After ParticleID, use RecMdcKalTrack substitute
735 // RecMdcTrack
736 RecMdcKalTrack::setPidType( RecMdcKalTrack::pion ); // PID can set to electron, muon, pion,
737 // kaon and proton;The default setting
738 // is pion
739
740 if ( mdcKalTrk->charge() > 0 )
741 {
742 ipip.push_back( iGood[i] );
743 HepLorentzVector ptrk;
744 ptrk.setPx( mdcKalTrk->px() );
745 ptrk.setPy( mdcKalTrk->py() );
746 ptrk.setPz( mdcKalTrk->pz() );
747 double p3 = ptrk.mag();
748 ptrk.setE( sqrt( p3 * p3 + mpi * mpi ) );
749
750 // ptrk = ptrk.boost(-0.011,0,0);//boost to cms
751
752 ppip.push_back( ptrk );
753 }
754 else
755 {
756 ipim.push_back( iGood[i] );
757 HepLorentzVector ptrk;
758 ptrk.setPx( mdcKalTrk->px() );
759 ptrk.setPy( mdcKalTrk->py() );
760 ptrk.setPz( mdcKalTrk->pz() );
761 double p3 = ptrk.mag();
762 ptrk.setE( sqrt( p3 * p3 + mpi * mpi ) );
763
764 // ptrk = ptrk.boost(-0.011,0,0);//boost to cms
765
766 ppim.push_back( ptrk );
767 }
768 }
769
770 /*
771 for(int i = 0; i < nGood; i++) {//for rhopi without PID
772 EvtRecTrackIterator itTrk = evtRecTrkCol->begin() + iGood[i];
773 RecMdcTrack* mdcTrk = (*itTrk)->mdcTrack();
774 if(mdcTrk->charge() >0) {
775 ipip.push_back(iGood[i]);
776 HepLorentzVector ptrk;
777 ptrk.setPx(mdcTrk->px());
778 ptrk.setPy(mdcTrk->py());
779 ptrk.setPz(mdcTrk->pz());
780 double p3 = ptrk.mag();
781 ptrk.setE(sqrt(p3*p3+mpi*mpi));
782 ppip.push_back(ptrk);
783 } else {
784 ipim.push_back(iGood[i]);
785 HepLorentzVector ptrk;
786 ptrk.setPx(mdcTrk->px());
787 ptrk.setPy(mdcTrk->py());
788 ptrk.setPz(mdcTrk->pz());
789 double p3 = ptrk.mag();
790 ptrk.setE(sqrt(p3*p3+mpi*mpi));
791 ppim.push_back(ptrk);
792 }
793 }// without PID
794 */
795
796 int npip = ipip.size();
797 int npim = ipim.size();
798 if ( npip * npim != 1 ) return StatusCode::SUCCESS;
799
800 Ncut3++;
801
802 //
803 // Loop each gamma pair, check ppi0 and pTot
804 //
805
806 HepLorentzVector pTot;
807 for ( int i = 0; i < nGam - 1; i++ )
808 {
809 for ( int j = i + 1; j < nGam; j++ )
810 {
811 HepLorentzVector p2g = pGam[i] + pGam[j];
812 pTot = ppip[0] + ppim[0];
813 pTot += p2g;
814 m_m2gg = p2g.m();
815 m_etot = pTot.e();
816 m_tuple3->write();
817 }
818 }
819
820 RecMdcKalTrack* pipTrk = ( *( evtRecTrkCol->begin() + ipip[0] ) )->mdcKalTrack();
821 RecMdcKalTrack* pimTrk = ( *( evtRecTrkCol->begin() + ipim[0] ) )->mdcKalTrack();
822
823 WTrackParameter wvpipTrk, wvpimTrk;
824 wvpipTrk = WTrackParameter( mpi, pipTrk->getZHelix(), pipTrk->getZError() );
825 wvpimTrk = WTrackParameter( mpi, pimTrk->getZHelix(), pimTrk->getZError() );
826
827 /* Default is pion, for other particles:
828 wvppTrk = WTrackParameter(mp, pipTrk->getZHelixP(), pipTrk->getZErrorP());//proton
829 wvmupTrk = WTrackParameter(mmu, pipTrk->getZHelixMu(), pipTrk->getZErrorMu());//muon
830 wvepTrk = WTrackParameter(me, pipTrk->getZHelixE(), pipTrk->getZErrorE());//electron
831 wvkpTrk = WTrackParameter(mk, pipTrk->getZHelixK(), pipTrk->getZErrorK());//kaon
832 */
833 //
834 // Test vertex fit
835 //
836
837 HepPoint3D vx( 0., 0., 0. );
838 HepSymMatrix Evx( 3, 0 );
839 double bx = 1E+6;
840 double by = 1E+6;
841 double bz = 1E+6;
842 Evx[0][0] = bx * bx;
843 Evx[1][1] = by * by;
844 Evx[2][2] = bz * bz;
845
846 VertexParameter vxpar;
847 vxpar.setVx( vx );
848 vxpar.setEvx( Evx );
849
850 VertexFit* vtxfit = VertexFit::instance();
851 vtxfit->init();
852 vtxfit->AddTrack( 0, wvpipTrk );
853 vtxfit->AddTrack( 1, wvpimTrk );
854 vtxfit->AddVertex( 0, vxpar, 0, 1 );
855 if ( !vtxfit->Fit( 0 ) ) return StatusCode::SUCCESS;
856 vtxfit->Swim( 0 );
857
858 WTrackParameter wpip = vtxfit->wtrk( 0 );
859 WTrackParameter wpim = vtxfit->wtrk( 1 );
860
861 // KinematicFit * kmfit = KinematicFit::instance();
862 KalmanKinematicFit* kmfit = KalmanKinematicFit::instance();
863
864 //
865 // Apply Kinematic 4C fit
866 //
867 // cout<<"before 4c"<<endl;
868 if ( m_test4C == 1 )
869 {
870 // double ecms = 3.097;
871 HepLorentzVector ecms( 0.034, 0, 0, 3.097 );
872
873 double chisq = 9999.;
874 int ig1 = -1;
875 int ig2 = -1;
876 for ( int i = 0; i < nGam - 1; i++ )
877 {
878 RecEmcShower* g1Trk = ( *( evtRecTrkCol->begin() + iGam[i] ) )->emcShower();
879 for ( int j = i + 1; j < nGam; j++ )
880 {
881 RecEmcShower* g2Trk = ( *( evtRecTrkCol->begin() + iGam[j] ) )->emcShower();
882 kmfit->init();
883 kmfit->AddTrack( 0, wpip );
884 kmfit->AddTrack( 1, wpim );
885 kmfit->AddTrack( 2, 0.0, g1Trk );
886 kmfit->AddTrack( 3, 0.0, g2Trk );
887 kmfit->AddFourMomentum( 0, ecms );
888 bool oksq = kmfit->Fit();
889 if ( oksq )
890 {
891 double chi2 = kmfit->chisq();
892 if ( chi2 < chisq )
893 {
894 chisq = chi2;
895 ig1 = iGam[i];
896 ig2 = iGam[j];
897 }
898 }
899 }
900 }
901
902 if ( chisq < 200 )
903 {
904
905 RecEmcShower* g1Trk = ( *( evtRecTrkCol->begin() + ig1 ) )->emcShower();
906 RecEmcShower* g2Trk = ( *( evtRecTrkCol->begin() + ig2 ) )->emcShower();
907 kmfit->init();
908 kmfit->AddTrack( 0, wpip );
909 kmfit->AddTrack( 1, wpim );
910 kmfit->AddTrack( 2, 0.0, g1Trk );
911 kmfit->AddTrack( 3, 0.0, g2Trk );
912 kmfit->AddFourMomentum( 0, ecms );
913 bool oksq = kmfit->Fit();
914 if ( oksq )
915 {
916 HepLorentzVector ppi0 = kmfit->pfit( 2 ) + kmfit->pfit( 3 );
917 m_mpi0 = ppi0.m();
918 m_chi1 = kmfit->chisq();
919 m_tuple4->write();
920 Ncut4++;
921 }
922 }
923 }
924
925 //
926 // Apply Kinematic 5C Fit
927 //
928
929 // find the best combination over all possible pi+ pi- gamma gamma pair
930 if ( m_test5C == 1 )
931 {
932 // double ecms = 3.097;
933 HepLorentzVector ecms( 0.034, 0, 0, 3.097 );
934 double chisq = 9999.;
935 int ig1 = -1;
936 int ig2 = -1;
937 for ( int i = 0; i < nGam - 1; i++ )
938 {
939 RecEmcShower* g1Trk = ( *( evtRecTrkCol->begin() + iGam[i] ) )->emcShower();
940 for ( int j = i + 1; j < nGam; j++ )
941 {
942 RecEmcShower* g2Trk = ( *( evtRecTrkCol->begin() + iGam[j] ) )->emcShower();
943 kmfit->init();
944 kmfit->AddTrack( 0, wpip );
945 kmfit->AddTrack( 1, wpim );
946 kmfit->AddTrack( 2, 0.0, g1Trk );
947 kmfit->AddTrack( 3, 0.0, g2Trk );
948 kmfit->AddResonance( 0, 0.135, 2, 3 );
949 kmfit->AddFourMomentum( 1, ecms );
950 if ( !kmfit->Fit( 0 ) ) continue;
951 if ( !kmfit->Fit( 1 ) ) continue;
952 bool oksq = kmfit->Fit();
953 if ( oksq )
954 {
955 double chi2 = kmfit->chisq();
956 if ( chi2 < chisq )
957 {
958 chisq = chi2;
959 ig1 = iGam[i];
960 ig2 = iGam[j];
961 }
962 }
963 }
964 }
965
966 log << MSG::INFO << " chisq = " << chisq << endmsg;
967
968 if ( chisq < 200 )
969 {
970 RecEmcShower* g1Trk = ( *( evtRecTrkCol->begin() + ig1 ) )->emcShower();
971 RecEmcShower* g2Trk = ( *( evtRecTrkCol->begin() + ig2 ) )->emcShower();
972
973 kmfit->init();
974 kmfit->AddTrack( 0, wpip );
975 kmfit->AddTrack( 1, wpim );
976 kmfit->AddTrack( 2, 0.0, g1Trk );
977 kmfit->AddTrack( 3, 0.0, g2Trk );
978 kmfit->AddResonance( 0, 0.135, 2, 3 );
979 kmfit->AddFourMomentum( 1, ecms );
980 bool oksq = kmfit->Fit();
981 if ( oksq )
982 {
983 HepLorentzVector ppi0 = kmfit->pfit( 2 ) + kmfit->pfit( 3 );
984 HepLorentzVector prho0 = kmfit->pfit( 0 ) + kmfit->pfit( 1 );
985 HepLorentzVector prhop = ppi0 + kmfit->pfit( 0 );
986 HepLorentzVector prhom = ppi0 + kmfit->pfit( 1 );
987
988 m_chi2 = kmfit->chisq();
989 m_mrh0 = prho0.m();
990 m_mrhp = prhop.m();
991 m_mrhm = prhom.m();
992 double eg1 = ( kmfit->pfit( 2 ) ).e();
993 double eg2 = ( kmfit->pfit( 3 ) ).e();
994 double fcos = abs( eg1 - eg2 ) / ppi0.rho();
995 m_tuple5->write();
996 Ncut5++;
997 //
998 // Measure the photon detection efficiences via
999 // J/psi -> rho0 pi0
1000 //
1001 if ( fabs( prho0.m() - 0.770 ) < 0.150 )
1002 {
1003 if ( fabs( fcos ) < 0.99 )
1004 {
1005 m_fcos = ( eg1 - eg2 ) / ppi0.rho();
1006 m_elow = ( eg1 < eg2 ) ? eg1 : eg2;
1007 m_tuple6->write();
1008 Ncut6++;
1009 }
1010 } // rho0 cut
1011 } // oksq
1012 }
1013 }
1014 return StatusCode::SUCCESS;
1015}
HepGeom::Point3D< double > HepPoint3D
int runNo
Definition DQA_TO_DB.cxx:13
Double_t time
EvtRecTrackCol::iterator EvtRecTrackIterator
double mpi
double pi
std::vector< HepLorentzVector > Vp4
Definition Gam4pikp.cxx:38
const double xmass[5]
Definition Gam4pikp.cxx:35
const double velc
Definition Gam4pikp.cxx:36
std::vector< int > Vint
Definition Gam4pikp.cxx:37
int Ncut2
Definition Ppjrhopi.cxx:35
int Ncut1
Definition Ppjrhopi.cxx:35
int Ncut0
Definition Ppjrhopi.cxx:35
int Ncut5
Definition Ppjrhopi.cxx:35
int Ncut3
Definition Ppjrhopi.cxx:35
int Ncut6
Definition Ppjrhopi.cxx:35
int Ncut4
Definition Ppjrhopi.cxx:35
IMessageSvc * msgSvc()
const HepSymMatrix err() const
const HepVector helix() const
......
void AddFourMomentum(int number, HepLorentzVector p4)
void AddResonance(int number, double mres, std::vector< int > tlis)
static KalmanKinematicFit * instance()
double chiTof2(int n) const
static ParticleID * instance()
bool IsPidInfoValid() const
double chiTof1(int n) const
void calculate()
void init()
double chiDedx(int n) const
void setStatus(unsigned int status)
void AddTrack(const int number, const double mass, const RecMdcTrack *trk)
Definition TrackPool.cxx:21
WTrackParameter wtrk(int n) const
void init()
Definition VertexFit.cxx:27
void AddVertex(int number, VertexParameter vpar, std::vector< int > lis)
Definition VertexFit.cxx:85
static VertexFit * instance()
Definition VertexFit.cxx:15
bool Fit()
const double ecms
Definition inclkstar.cxx:26

◆ finalize()

StatusCode Rhopi::finalize ( )

Definition at line 1018 of file Rhopi.cxx.

1018 {
1019 cout << "total number: " << Ncut0 << endl;
1020 cout << "nGood==2, nCharge==0: " << Ncut1 << endl;
1021 cout << "nGam>=2: " << Ncut2 << endl;
1022 cout << "Pass Pid: " << Ncut3 << endl;
1023 cout << "Pass 4C: " << Ncut4 << endl;
1024 cout << "Pass 5C: " << Ncut5 << endl;
1025 cout << "J/psi->rho0 pi0: " << Ncut6 << endl;
1026 MsgStream log( msgSvc(), name() );
1027 log << MSG::INFO << "in finalize()" << endmsg;
1028 return StatusCode::SUCCESS;
1029}

◆ initialize()

StatusCode Rhopi::initialize ( )

Definition at line 71 of file Rhopi.cxx.

71 {
72 MsgStream log( msgSvc(), name() );
73
74 log << MSG::INFO << "in initialize()" << endmsg;
75
76 StatusCode status;
77 NTuplePtr nt1( ntupleSvc(), "FILE1/vxyz" );
78 if ( nt1 ) m_tuple1 = nt1;
79 else
80 {
81 m_tuple1 = ntupleSvc()->book( "FILE1/vxyz", CLID_ColumnWiseTuple, "ks N-Tuple example" );
82 if ( m_tuple1 )
83 {
84 status = m_tuple1->addItem( "vx0", m_vx0 );
85 status = m_tuple1->addItem( "vy0", m_vy0 );
86 status = m_tuple1->addItem( "vz0", m_vz0 );
87 status = m_tuple1->addItem( "vr0", m_vr0 );
88 status = m_tuple1->addItem( "rvxy0", m_rvxy0 );
89 status = m_tuple1->addItem( "rvz0", m_rvz0 );
90 status = m_tuple1->addItem( "rvphi0", m_rvphi0 );
91 }
92 else
93 {
94 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple1 ) << endmsg;
95 return StatusCode::FAILURE;
96 }
97 }
98
99 NTuplePtr nt2( ntupleSvc(), "FILE1/photon" );
100 if ( nt2 ) m_tuple2 = nt2;
101 else
102 {
103 m_tuple2 = ntupleSvc()->book( "FILE1/photon", CLID_ColumnWiseTuple, "ks N-Tuple example" );
104 if ( m_tuple2 )
105 {
106 status = m_tuple2->addItem( "dthe", m_dthe );
107 status = m_tuple2->addItem( "dphi", m_dphi );
108 status = m_tuple2->addItem( "dang", m_dang );
109 status = m_tuple2->addItem( "eraw", m_eraw );
110 }
111 else
112 {
113 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple2 ) << endmsg;
114 return StatusCode::FAILURE;
115 }
116 }
117
118 NTuplePtr nt3( ntupleSvc(), "FILE1/etot" );
119 if ( nt3 ) m_tuple3 = nt3;
120 else
121 {
122 m_tuple3 = ntupleSvc()->book( "FILE1/etot", CLID_ColumnWiseTuple, "ks N-Tuple example" );
123 if ( m_tuple3 )
124 {
125 status = m_tuple3->addItem( "m2gg", m_m2gg );
126 status = m_tuple3->addItem( "etot", m_etot );
127 }
128 else
129 {
130 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple3 ) << endmsg;
131 return StatusCode::FAILURE;
132 }
133 }
134 if ( m_test4C == 1 )
135 {
136 NTuplePtr nt4( ntupleSvc(), "FILE1/fit4c" );
137 if ( nt4 ) m_tuple4 = nt4;
138 else
139 {
140 m_tuple4 =
141 ntupleSvc()->book( "FILE1/fit4c", CLID_ColumnWiseTuple, "ks N-Tuple example" );
142 if ( m_tuple4 )
143 {
144 status = m_tuple4->addItem( "chi2", m_chi1 );
145 status = m_tuple4->addItem( "mpi0", m_mpi0 );
146 }
147 else
148 {
149 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple4 ) << endmsg;
150 return StatusCode::FAILURE;
151 }
152 }
153 } // test 4C
154
155 if ( m_test5C == 1 )
156 {
157 NTuplePtr nt5( ntupleSvc(), "FILE1/fit5c" );
158 if ( nt5 ) m_tuple5 = nt5;
159 else
160 {
161 m_tuple5 =
162 ntupleSvc()->book( "FILE1/fit5c", CLID_ColumnWiseTuple, "ks N-Tuple example" );
163 if ( m_tuple5 )
164 {
165 status = m_tuple5->addItem( "chi2", m_chi2 );
166 status = m_tuple5->addItem( "mrh0", m_mrh0 );
167 status = m_tuple5->addItem( "mrhp", m_mrhp );
168 status = m_tuple5->addItem( "mrhm", m_mrhm );
169 }
170 else
171 {
172 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple5 ) << endmsg;
173 return StatusCode::FAILURE;
174 }
175 }
176
177 NTuplePtr nt6( ntupleSvc(), "FILE1/geff" );
178 if ( nt6 ) m_tuple6 = nt6;
179 else
180 {
181 m_tuple6 = ntupleSvc()->book( "FILE1/geff", CLID_ColumnWiseTuple, "ks N-Tuple example" );
182 if ( m_tuple6 )
183 {
184 status = m_tuple6->addItem( "fcos", m_fcos );
185 status = m_tuple6->addItem( "elow", m_elow );
186 }
187 else
188 {
189 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple6 ) << endmsg;
190 return StatusCode::FAILURE;
191 }
192 }
193 } // test 5c
194
195 if ( m_checkDedx == 1 )
196 {
197 NTuplePtr nt7( ntupleSvc(), "FILE1/dedx" );
198 if ( nt7 ) m_tuple7 = nt7;
199 else
200 {
201 m_tuple7 = ntupleSvc()->book( "FILE1/dedx", CLID_ColumnWiseTuple, "ks N-Tuple example" );
202 if ( m_tuple7 )
203 {
204 status = m_tuple7->addItem( "ptrk", m_ptrk );
205 status = m_tuple7->addItem( "chie", m_chie );
206 status = m_tuple7->addItem( "chimu", m_chimu );
207 status = m_tuple7->addItem( "chipi", m_chipi );
208 status = m_tuple7->addItem( "chik", m_chik );
209 status = m_tuple7->addItem( "chip", m_chip );
210 status = m_tuple7->addItem( "probPH", m_probPH );
211 status = m_tuple7->addItem( "normPH", m_normPH );
212 status = m_tuple7->addItem( "ghit", m_ghit );
213 status = m_tuple7->addItem( "thit", m_thit );
214 }
215 else
216 {
217 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple7 ) << endmsg;
218 return StatusCode::FAILURE;
219 }
220 }
221 } // check dE/dx
222
223 if ( m_checkTof == 1 )
224 {
225 NTuplePtr nt8( ntupleSvc(), "FILE1/tofe" );
226 if ( nt8 ) m_tuple8 = nt8;
227 else
228 {
229 m_tuple8 = ntupleSvc()->book( "FILE1/tofe", CLID_ColumnWiseTuple, "ks N-Tuple example" );
230 if ( m_tuple8 )
231 {
232 status = m_tuple8->addItem( "ptrk", m_ptot_etof );
233 status = m_tuple8->addItem( "cntr", m_cntr_etof );
234 status = m_tuple8->addItem( "ph", m_ph_etof );
235 status = m_tuple8->addItem( "rhit", m_rhit_etof );
236 status = m_tuple8->addItem( "qual", m_qual_etof );
237 status = m_tuple8->addItem( "te", m_te_etof );
238 status = m_tuple8->addItem( "tmu", m_tmu_etof );
239 status = m_tuple8->addItem( "tpi", m_tpi_etof );
240 status = m_tuple8->addItem( "tk", m_tk_etof );
241 status = m_tuple8->addItem( "tp", m_tp_etof );
242 }
243 else
244 {
245 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple8 ) << endmsg;
246 return StatusCode::FAILURE;
247 }
248 }
249 } // check Tof:endcap
250
251 if ( m_checkTof == 1 )
252 {
253 NTuplePtr nt9( ntupleSvc(), "FILE1/tof1" );
254 if ( nt9 ) m_tuple9 = nt9;
255 else
256 {
257 m_tuple9 = ntupleSvc()->book( "FILE1/tof1", CLID_ColumnWiseTuple, "ks N-Tuple example" );
258 if ( m_tuple9 )
259 {
260 status = m_tuple9->addItem( "ptrk", m_ptot_btof1 );
261 status = m_tuple9->addItem( "cntr", m_cntr_btof1 );
262 status = m_tuple9->addItem( "ph", m_ph_btof1 );
263 status = m_tuple9->addItem( "zhit", m_zhit_btof1 );
264 status = m_tuple9->addItem( "qual", m_qual_btof1 );
265 status = m_tuple9->addItem( "te", m_te_btof1 );
266 status = m_tuple9->addItem( "tmu", m_tmu_btof1 );
267 status = m_tuple9->addItem( "tpi", m_tpi_btof1 );
268 status = m_tuple9->addItem( "tk", m_tk_btof1 );
269 status = m_tuple9->addItem( "tp", m_tp_btof1 );
270 }
271 else
272 {
273 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple9 ) << endmsg;
274 return StatusCode::FAILURE;
275 }
276 }
277 } // check Tof:barrel inner Tof
278
279 if ( m_checkTof == 1 )
280 {
281 NTuplePtr nt10( ntupleSvc(), "FILE1/tof2" );
282 if ( nt10 ) m_tuple10 = nt10;
283 else
284 {
285 m_tuple10 =
286 ntupleSvc()->book( "FILE1/tof2", CLID_ColumnWiseTuple, "ks N-Tuple example" );
287 if ( m_tuple10 )
288 {
289 status = m_tuple10->addItem( "ptrk", m_ptot_btof2 );
290 status = m_tuple10->addItem( "cntr", m_cntr_btof2 );
291 status = m_tuple10->addItem( "ph", m_ph_btof2 );
292 status = m_tuple10->addItem( "zhit", m_zhit_btof2 );
293 status = m_tuple10->addItem( "qual", m_qual_btof2 );
294 status = m_tuple10->addItem( "te", m_te_btof2 );
295 status = m_tuple10->addItem( "tmu", m_tmu_btof2 );
296 status = m_tuple10->addItem( "tpi", m_tpi_btof2 );
297 status = m_tuple10->addItem( "tk", m_tk_btof2 );
298 status = m_tuple10->addItem( "tp", m_tp_btof2 );
299 }
300 else
301 {
302 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple10 ) << endmsg;
303 return StatusCode::FAILURE;
304 }
305 }
306 } // check Tof:barrel outter Tof
307
308 NTuplePtr nt11( ntupleSvc(), "FILE1/pid" );
309 if ( nt11 ) m_tuple11 = nt11;
310 else
311 {
312 m_tuple11 = ntupleSvc()->book( "FILE1/pid", CLID_ColumnWiseTuple, "ks N-Tuple example" );
313 if ( m_tuple11 )
314 {
315 status = m_tuple11->addItem( "ptrk", m_ptrk_pid );
316 status = m_tuple11->addItem( "cost", m_cost_pid );
317 status = m_tuple11->addItem( "dedx", m_dedx_pid );
318 status = m_tuple11->addItem( "tof1", m_tof1_pid );
319 status = m_tuple11->addItem( "tof2", m_tof2_pid );
320 status = m_tuple11->addItem( "prob", m_prob_pid );
321 }
322 else
323 {
324 log << MSG::ERROR << " Cannot book N-tuple:" << long( m_tuple11 ) << endmsg;
325 return StatusCode::FAILURE;
326 }
327 }
328
329 //
330 //--------end of book--------
331 //
332
333 log << MSG::INFO << "successfully return from initialize()" << endmsg;
334 return StatusCode::SUCCESS;
335}
INTupleSvc * ntupleSvc()

The documentation for this class was generated from the following files: