rivet is hosted by Hepforge, IPPP Durham
ARGUS_1993_S2653028.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include <iostream>
00003 #include "Rivet/Analysis.hh"
00004 #include "Rivet/Projections/Beam.hh"
00005 #include "Rivet/Projections/UnstableFinalState.hh"
00006 #include "Rivet/ParticleName.hh"
00007 
00008 namespace Rivet {
00009 
00010 
00011   /// @brief BELLE pi+/-, K+/- and proton/antiproton spectrum at Upsilon(4S)
00012   /// @author Peter Richardson
00013   class ARGUS_1993_S2653028 : public Analysis {
00014   public:
00015 
00016     ARGUS_1993_S2653028()
00017       : Analysis("ARGUS_1993_S2653028"), _weightSum(0.)
00018     { }
00019 
00020 
00021     void analyze(const Event& e) {
00022       const double weight = e.weight();
00023 
00024       // Find the upsilons
00025       Particles upsilons;
00026       // First in unstable final state
00027       const UnstableFinalState& ufs = applyProjection<UnstableFinalState>(e, "UFS");
00028       foreach (const Particle& p, ufs.particles()) {
00029         if (p.pdgId() == 300553) upsilons.push_back(p);
00030       }
00031       // Then in whole event if that failed
00032       if (upsilons.empty()) {
00033         foreach (GenParticle* p, Rivet::particles(e.genEvent())) {
00034           if (p->pdg_id() != 300553) continue;
00035           const GenVertex* pv = p->production_vertex();
00036           bool passed = true;
00037           if (pv) {
00038             /// @todo Use better looping
00039             for (GenVertex::particles_in_const_iterator pp = pv->particles_in_const_begin(); pp != pv->particles_in_const_end() ; ++pp) {
00040               if ( p->pdg_id() == (*pp)->pdg_id() ) {
00041                 passed = false;
00042                 break;
00043               }
00044             }
00045           }
00046           if (passed) upsilons.push_back(Particle(*p));
00047         }
00048       }
00049 
00050       // Find an upsilon
00051       foreach (const Particle& p, upsilons) {
00052         _weightSum += weight;
00053         vector<GenParticle *> pionsA,pionsB,protonsA,protonsB,kaons;
00054         // Find the decay products we want
00055         findDecayProducts(p.genParticle(), pionsA, pionsB, protonsA, protonsB, kaons);
00056         LorentzTransform cms_boost;
00057         if (p.momentum().vector3().mod() > 1*MeV)
00058           cms_boost = LorentzTransform(-p.momentum().boostVector());
00059         for (size_t ix = 0; ix < pionsA.size(); ++ix) {
00060           FourMomentum ptemp(pionsA[ix]->momentum());
00061           FourMomentum p2 = cms_boost.transform(ptemp);
00062           double pcm = cms_boost.transform(ptemp).vector3().mod();
00063           _histPiA->fill(pcm,weight);
00064         }
00065         _multPiA->fill(10.58,double(pionsA.size())*weight);
00066         for (size_t ix = 0; ix < pionsB.size(); ++ix) {
00067           double pcm = cms_boost.transform(FourMomentum(pionsB[ix]->momentum())).vector3().mod();
00068           _histPiB->fill(pcm,weight);
00069         }
00070         _multPiB->fill(10.58,double(pionsB.size())*weight);
00071         for (size_t ix = 0; ix < protonsA.size(); ++ix) {
00072           double pcm = cms_boost.transform(FourMomentum(protonsA[ix]->momentum())).vector3().mod();
00073           _histpA->fill(pcm,weight);
00074         }
00075         _multpA->fill(10.58,double(protonsA.size())*weight);
00076         for (size_t ix = 0; ix < protonsB.size(); ++ix) {
00077           double pcm = cms_boost.transform(FourMomentum(protonsB[ix]->momentum())).vector3().mod();
00078           _histpB->fill(pcm,weight);
00079         }
00080         _multpB->fill(10.58,double(protonsB.size())*weight);
00081         for (size_t ix = 0 ;ix < kaons.size(); ++ix) {
00082           double pcm = cms_boost.transform(FourMomentum(kaons[ix]->momentum())).vector3().mod();
00083           _histKA->fill(pcm,weight);
00084           _histKB->fill(pcm,weight);
00085         }
00086         _multK->fill(10.58,double(kaons.size())*weight);
00087       }
00088     }
00089 
00090 
00091     void finalize() {
00092       if (_weightSum > 0.) {
00093         scale(_histPiA, 1./_weightSum);
00094         scale(_histPiB, 1./_weightSum);
00095         scale(_histKA , 1./_weightSum);
00096         scale(_histKB , 1./_weightSum);
00097         scale(_histpA , 1./_weightSum);
00098         scale(_histpB , 1./_weightSum);
00099         scale(_multPiA, 1./_weightSum);
00100         scale(_multPiB, 1./_weightSum);
00101         scale(_multK  , 1./_weightSum);
00102         scale(_multpA , 1./_weightSum);
00103         scale(_multpB , 1./_weightSum);
00104       }
00105     }
00106 
00107 
00108     void init() {
00109       addProjection(UnstableFinalState(), "UFS");
00110 
00111       // spectra
00112       _histPiA = bookHisto1D(1, 1, 1);
00113       _histPiB = bookHisto1D(2, 1, 1);
00114       _histKA  = bookHisto1D(3, 1, 1);
00115       _histKB  = bookHisto1D(6, 1, 1);
00116       _histpA  = bookHisto1D(4, 1, 1);
00117       _histpB  = bookHisto1D(5, 1, 1);
00118       // multiplicities
00119       _multPiA = bookHisto1D( 7, 1, 1);
00120       _multPiB = bookHisto1D( 8, 1, 1);
00121       _multK   = bookHisto1D( 9, 1, 1);
00122       _multpA  = bookHisto1D(10, 1, 1);
00123       _multpB  = bookHisto1D(11, 1, 1);
00124     } // init
00125 
00126 
00127   private:
00128 
00129     //@{
00130     /// Count of weights
00131     double _weightSum;
00132     /// Spectra
00133     Histo1DPtr _histPiA, _histPiB, _histKA, _histKB, _histpA, _histpB;
00134     /// Multiplicities
00135     Histo1DPtr _multPiA, _multPiB, _multK, _multpA, _multpB;
00136     //@}
00137 
00138     void findDecayProducts(const GenParticle* p,
00139                            vector<GenParticle*>& pionsA, vector<GenParticle*>& pionsB,
00140                            vector<GenParticle*>& protonsA, vector<GenParticle*>& protonsB,
00141                            vector<GenParticle*>& kaons)
00142     {
00143       int parentId = p->pdg_id();
00144       const GenVertex* dv = p->end_vertex();
00145       /// @todo Use better looping
00146       for (GenVertex::particles_out_const_iterator pp = dv->particles_out_const_begin(); pp != dv->particles_out_const_end(); ++pp) {
00147         int id = abs((*pp)->pdg_id());
00148         if (id == PID::PIPLUS) {
00149           if (parentId != PID::LAMBDA && parentId != PID::K0S) {
00150             pionsA.push_back(*pp);
00151             pionsB.push_back(*pp);
00152           }
00153           else
00154             pionsB.push_back(*pp);
00155         }
00156         else if (id == PID::PROTON) {
00157           if (parentId != PID::LAMBDA && parentId != PID::K0S) {
00158             protonsA.push_back(*pp);
00159             protonsB.push_back(*pp);
00160           }
00161           else
00162             protonsB.push_back(*pp);
00163         }
00164         else if (id == PID::KPLUS) {
00165           kaons.push_back(*pp);
00166         }
00167         else if ((*pp)->end_vertex())
00168           findDecayProducts(*pp, pionsA, pionsB, protonsA, protonsB, kaons);
00169       }
00170     }
00171   };
00172 
00173 
00174   // The hook for the plugin system
00175   DECLARE_RIVET_PLUGIN(ARGUS_1993_S2653028);
00176 
00177 }