CMS 3D CMS Logo

RawParticle Class Reference

#include <FastSimulation/Particle/interface/RawParticle.h>

Inheritance diagram for RawParticle:

BaseParticlePropagator ParticlePropagator

List of all members.

Public Types

typedef ROOT::Math::Boost Boost
typedef ROOT::Math::AxisAngle Rotation
typedef ROOT::Math::Rotation3D Rotation3D
typedef ROOT::Math::RotationX RotationX
typedef ROOT::Math::RotationY RotationY
typedef ROOT::Math::RotationZ RotationZ

Public Member Functions

void boost (const Boost &b)
void boost (double bx, double by, double bz)
 Boost the particle.
double charge () const
 get the MEASURED charge
void chargeConjugate ()
 Convert the particle to its charge conjugate state.
double cos2Theta () const
 Cos**2(theta) is faster to determine than eta.
double cos2ThetaV () const
double et () const
 get the transverse energy
double eta () const
 Get the pseudo rapidity of the particle.
int isUsed () const
 Is the particle marked as used.
double mass () const
 get the MEASURED mass
const XYZTLorentzVectormomentum () const
 the momentum fourvector
RawParticleoperator= (const RawParticle &rhs)
 Copy assignment operator.
double PDGcharge () const
 get the THEORETICAL charge
double PDGcTau () const
 get the THEORETICAL lifetime
double PDGmass () const
 get the THEORETICAL mass
std::string PDGname () const
 get the PDG name
int pid () const
 get the HEP particle ID number
void print () const
 Print the formated particle information.
void printName () const
 Print the name of the particle.
double R () const
 vertex radius
double r () const
 vertex radius
double R2 () const
 vertex radius**2
double r2 () const
 vertex radius**2
 RawParticle (const RawParticle &p)
 Copy constructor.
 RawParticle (double px, double py, double pz, double e)
 Construct from fourmomentum components.
 RawParticle (const XYZTLorentzVector &p, const XYZTLorentzVector &xStart)
 Construct from 2 fourvectors.
 RawParticle (const std::string name, const XYZTLorentzVector &p)
 Construct from a fourvector and a name.
 RawParticle (const int id, const XYZTLorentzVector &p)
 Construct from a fourvector and a PID.
 RawParticle (const XYZTLorentzVector &p)
 Construct from a fourvector.
 RawParticle ()
void reUse ()
 Unlock the particle, see isUsed().
void rotate (const RotationZ &r)
void rotate (const RotationY &r)
void rotate (const RotationX &r)
void rotate (const Rotation3D &r)
void rotate (const Rotation &r)
void rotate (double rphi, const XYZVector &raxis)
 Rotate the particle around an axis in space.
void rotateX (double rphi)
 
Warning:
not yet implemented

void rotateY (double rphi)
 Rotate around z axis.
void rotateZ (double rphi)
 Rotate around z axis.
void setCharge (float q)
 set the MEASURED charge
void setID (const std::string name)
 Set identifier for this particle.
void setID (const int id)
 Set identifier for this particle.
void setMass (float m)
 set the RECONSTRUCTED mass
void setStatus (int istat)
 Set the status of this particle.
void setT (const double t)
 set the time of creation
void setVertex (double xv, double yv, double zv, double tv)
void setVertex (const XYZTLorentzVector &vtx)
 set the vertex
int status () const
 get the particle status
double T () const
 vertex time
double t () const
 vertex time
void translate (const XYZVector &t)
 Translate the vertex by a given space amount.
void use ()
 Lock the particle, see isUsed().
const XYZTLorentzVectorvertex () const
 the vertex fourvector
double X () const
 x of vertex
double x () const
 x of vertex
double Y () const
 y of vertex
double y () const
 y of vertex
double Z () const
 z of vertex
double z () const
 z of vertex
virtual ~RawParticle ()

Protected Attributes

double myCharge
 the MEASURED charge
int myId
 the particle id number HEP-PID
const ParticleDatamyInfo
 The pointer to the PDG info.
double myMass
 the RECONSTRUCTED mass
int myStatus
 the status code according to PYTHIA
int myUsed
 status of the locking
XYZTLorentzVector myVertex
 the four vector of the vertex

Private Member Functions

void init ()

Private Attributes

ParticleTabletab


Detailed Description

Definition at line 31 of file RawParticle.h.


Member Typedef Documentation

typedef ROOT::Math::Boost RawParticle::Boost

Definition at line 39 of file RawParticle.h.

typedef ROOT::Math::AxisAngle RawParticle::Rotation

Definition at line 34 of file RawParticle.h.

typedef ROOT::Math::Rotation3D RawParticle::Rotation3D

Definition at line 35 of file RawParticle.h.

typedef ROOT::Math::RotationX RawParticle::RotationX

Definition at line 36 of file RawParticle.h.

typedef ROOT::Math::RotationY RawParticle::RotationY

Definition at line 37 of file RawParticle.h.

typedef ROOT::Math::RotationZ RawParticle::RotationZ

Definition at line 38 of file RawParticle.h.


Constructor & Destructor Documentation

RawParticle::RawParticle (  ) 

Definition at line 18 of file RawParticle.cc.

References init().

00018                          {
00019   init();
00020 }

RawParticle::~RawParticle (  )  [virtual]

Definition at line 68 of file RawParticle.cc.

00068                           {
00069   //  nParticles--;
00070 }

RawParticle::RawParticle ( const XYZTLorentzVector p  ) 

Construct from a fourvector.

The fourvector is taken for the particle, the vertex is set to 0.

Definition at line 22 of file RawParticle.cc.

References init().

00023   : XYZTLorentzVector(p) {
00024   init();
00025 }

RawParticle::RawParticle ( const int  id,
const XYZTLorentzVector p 
)

Construct from a fourvector and a PID.

The fourvector and PID are taken for the particle, the vertex is set to 0.

Definition at line 27 of file RawParticle.cc.

References init(), and setID().

00029   : XYZTLorentzVector(p) {
00030   this->init();
00031   this->setID(id);
00032 }

RawParticle::RawParticle ( const std::string  name,
const XYZTLorentzVector p 
)

Construct from a fourvector and a name.

The fourvector and name are taken for the particle, the vertex is set to 0.

Definition at line 34 of file RawParticle.cc.

References init(), and setID().

00036   : XYZTLorentzVector(p) {
00037   this->init();
00038   this->setID(name);
00039 }

RawParticle::RawParticle ( const XYZTLorentzVector p,
const XYZTLorentzVector xStart 
)

Construct from 2 fourvectors.

The first fourvector is taken for the particle, the second for its vertex.

Definition at line 41 of file RawParticle.cc.

References init(), and myVertex.

00042                                                            : 
00043   XYZTLorentzVector(p)
00044 {
00045   init();
00046   myVertex = xStart;
00047 }

RawParticle::RawParticle ( double  px,
double  py,
double  pz,
double  e 
)

Construct from fourmomentum components.

Vertex is set to 0.

Definition at line 49 of file RawParticle.cc.

References init().

00049                                                                   : 
00050   XYZTLorentzVector(px,py,pz,e)
00051 {
00052   init();
00053 }

RawParticle::RawParticle ( const RawParticle p  ) 

Copy constructor.

Definition at line 55 of file RawParticle.cc.

References myCharge, myId, myInfo, myMass, myStatus, myUsed, myVertex, and tab.

00055                                                  : 
00056   XYZTLorentzVector(right.Px(),right.Py(),right.Pz(),right.E()) 
00057 {
00058   myId     = right.myId; 
00059   myStatus = right.myStatus;
00060   myUsed   = right.myUsed;
00061   myCharge = right.myCharge;
00062   myMass   = right.myMass;
00063   myVertex = (right.myVertex);
00064   tab = (right.tab);
00065   myInfo = (right.myInfo);
00066 }


Member Function Documentation

void RawParticle::boost ( const Boost b  )  [inline]

Definition at line 310 of file RawParticle.h.

References b, momentum(), and p.

00310                                                         { 
00311   XYZTLorentzVector p ( b(momentum()) ); SetXYZT(p.Px(),p.Py(),p.Pz(),p.E());
00312 }

void RawParticle::boost ( double  bx,
double  by,
double  bz 
)

Boost the particle.

The arguments are the $\beta$ values of the boost in x, y and z direction.

Warning:
What happens to the vertex?

Definition at line 182 of file RawParticle.cc.

References b, momentum(), and p.

Referenced by NuclearInteractionSimulator::compute().

00182                                                            {
00183   Boost b(betax,betay,betaz);
00184   XYZTLorentzVector p ( b * momentum() );
00185   SetXYZT(p.X(),p.Y(),p.Z(),p.T());
00186 }

double RawParticle::charge ( void   )  const [inline]

get the MEASURED charge

Definition at line 281 of file RawParticle.h.

References myCharge.

Referenced by BaseParticlePropagator::backPropagate(), EnergyLossSimulator::compute(), MultipleScatteringSimulator::compute(), NuclearInteractionSimulator::distanceToPrimary(), ParticlePropagator::fieldMap(), BaseParticlePropagator::helixRadius(), MaterialEffects::interact(), KineParticleFilter::isOKForMe(), TrajectoryManager::makeTrajectoryState(), ConvBremSeedProducer::makeTrajectoryState(), print(), BaseParticlePropagator::propagate(), BaseParticlePropagator::propagateToBeamCylinder(), BaseParticlePropagator::propagateToClosestApproach(), BaseParticlePropagator::propagateToNominalVertex(), TrajectoryManager::updateWithDaughters(), and BaseParticlePropagator::xyImpactParameter().

00281 { return myCharge; }

void RawParticle::chargeConjugate (  ) 

Convert the particle to its charge conjugate state.

This operation resets the particle ID to that of the charge conjugated particle (if one exists). Also the measured charge is multiplied by -1.

Definition at line 143 of file RawParticle.cc.

References myCharge, and myId.

00143                              {
00144   myId = -myId;
00145   myCharge = -1*myCharge;
00146 }

double RawParticle::cos2Theta (  )  const [inline]

Cos**2(theta) is faster to determine than eta.

Definition at line 267 of file RawParticle.h.

Referenced by KineParticleFilter::isOKForMe().

00267 { return Pz()*Pz()/Vect().Mag2(); }

double RawParticle::cos2ThetaV (  )  const [inline]

Definition at line 268 of file RawParticle.h.

References myVertex, and Z().

Referenced by KineParticleFilter::isOKForMe(), BaseParticlePropagator::propagateToEcalEntrance(), BaseParticlePropagator::propagateToHcalEntrance(), and BaseParticlePropagator::propagateToVFcalEntrance().

00268 { return Z()*Z()/myVertex.Vect().Mag2(); }

double RawParticle::et (  )  const

get the transverse energy

Definition at line 306 of file RawParticle.cc.

References momentum(), and funct::sqrt().

00306                       { 
00307   double mypp, tmpEt=-1.;
00308   
00309   mypp = std::sqrt(momentum().mag2());
00310   if ( mypp != 0 ) {
00311     tmpEt = E() * pt() / mypp;
00312   }
00313   return tmpEt;
00314 }

double RawParticle::eta ( void   )  const [inline]

Get the pseudo rapidity of the particle.

$ \eta = -\log ( \tan ( \vartheta/2)) $

Definition at line 266 of file RawParticle.h.

References funct::log(), and funct::tan().

Referenced by EcalHitMaker::preshowerCellLine().

00266 { return -std::log(std::tan(this->theta()/2.)); }

void RawParticle::init ( void   )  [private]

Reimplemented in BaseParticlePropagator.

Definition at line 90 of file RawParticle.cc.

References ParticleTable::instance(), myCharge, myId, myInfo, myMass, myStatus, myUsed, and tab.

Referenced by RawParticle().

00090                   {
00091   myId=0;  
00092   myStatus=99;
00093   myUsed=0;
00094   myCharge=0.;
00095   myMass=0.;
00096   tab = ParticleTable::instance();
00097   myInfo=0;
00098 }

int RawParticle::isUsed (  )  const [inline]

Is the particle marked as used.

The three methods isUsed(), use() and reUse() implement a simple locking mechanism.

Definition at line 233 of file RawParticle.h.

References myUsed.

00233 {return myUsed;}

double RawParticle::mass (  )  const [inline]

get the MEASURED mass

Definition at line 282 of file RawParticle.h.

References myMass.

Referenced by EnergyLossSimulator::compute(), NuclearInteractionSimulator::compute(), MultipleScatteringSimulator::compute(), Pythia6Decays::particleDaughters(), print(), and BaseParticlePropagator::propagate().

00282 { return myMass; }

const XYZTLorentzVector & RawParticle::momentum (  )  const [inline]

the momentum fourvector

Definition at line 285 of file RawParticle.h.

Referenced by PFTrackTransformer::addPoints(), PFTrackTransformer::addPointsAndBrems(), FBaseSimEvent::addSimTrack(), MuonSimHitProducer::applyMaterialEffects(), boost(), BremsstrahlungSimulator::compute(), et(), ConvBremSeedProducer::GoodCluster(), operator<<(), GoodSeedProducer::produce(), TrajectoryManager::propagateToCalorimeters(), and BaseParticlePropagator::propagateToClosestApproach().

00285 { return *this; }

RawParticle & RawParticle::operator= ( const RawParticle rhs  ) 

Copy assignment operator.

Definition at line 73 of file RawParticle.cc.

References myCharge, myId, myInfo, myMass, myStatus, myUsed, myVertex, and tab.

00073                                                    {
00074   //  cout << "Copy assignment " << endl;
00075   if (this != &right) { // don't copy into yourself
00076     this->SetXYZT(right.Px(),right.Py(),right.Pz(),right.E());
00077     myId     = right.myId; 
00078     myStatus = right.myStatus;
00079     myUsed   = right.myUsed;
00080     myCharge = right.myCharge;
00081     myMass   = right.myMass;
00082     myVertex = right.myVertex;
00083     tab      = right.tab;
00084     myInfo   = right.myInfo;
00085   }
00086   return *this;
00087 }

double RawParticle::PDGcharge (  )  const

get the THEORETICAL charge

Definition at line 245 of file RawParticle.cc.

References myInfo.

00245                              { 
00246   double q=-99999;
00247   if ( myInfo ) {
00248     q=myInfo->charge();
00249   }
00250   return q;
00251 }

double RawParticle::PDGcTau (  )  const

get the THEORETICAL lifetime

Definition at line 263 of file RawParticle.cc.

References funct::abs(), ct, myId, myInfo, and w.

Referenced by ParticlePropagator::initProperDecayTime().

00263                            {
00264   double ct=1E99;
00265   if ( myInfo ) {
00266 
00267     // The lifetime is 0. in the Pythia Particle Data Table !
00268     //    ct=tab->theTable()->particle(ParticleID(myId))->lifetime().value();
00269 
00270     // Get it from the width (apparently Gamma/c!)
00271     double w = myInfo->totalWidth().value();
00272     if ( w != 0. && myId != 1000022 ) { 
00273       ct = 6.582119e-25 / w / 10.;   // ctau in cm 
00274     } else {
00275     // Temporary fix of a bug in the particle data table
00276       unsigned amyId = abs(myId);
00277       if ( amyId != 22 &&    // photon 
00278            amyId != 11 &&    // e+/-
00279            amyId != 10 &&    // nu_e
00280            amyId != 12 &&    // nu_mu
00281            amyId != 14 &&    // nu_tau
00282            amyId != 1000022 && // Neutralino
00283            amyId != 1000039 && // Gravitino
00284            amyId != 2112 &&  // neutron/anti-neutron
00285            amyId != 2212 &&  // proton/anti-proton
00286            amyId != 101 &&   // Deutreron etc..
00287            amyId != 102 &&   // Deutreron etc..
00288            amyId != 103 &&   // Deutreron etc..
00289            amyId != 104 ) {  // Deutreron etc.. 
00290         ct = 0.;
00291         /* */
00292       }
00293     }
00294   }
00295 
00296   /*
00297   std::cout << setw(20) << setprecision(18) 
00298        << "myId/ctau/width = " << myId << " " 
00299        << ct << " " << w << endl;  
00300   */
00301 
00302   return ct;
00303 }

double RawParticle::PDGmass (  )  const

get the THEORETICAL mass

Definition at line 254 of file RawParticle.cc.

References m, and myInfo.

Referenced by FBaseSimEvent::addSimTrack().

00254                             { 
00255   double m=-99999;
00256   if ( myInfo ) {
00257     m = myInfo->mass().value();
00258   }
00259   return m;
00260 }

std::string RawParticle::PDGname (  )  const

get the PDG name

Definition at line 188 of file RawParticle.cc.

References myInfo, and tab.

Referenced by printName().

00188                                      {
00189   std::string MyParticleName;
00190   if ( tab && myInfo ) {
00191     MyParticleName = myInfo->name();
00192   } else {
00193     MyParticleName = "unknown  ";
00194   }
00195   return (std::string) MyParticleName;}

int RawParticle::pid (  )  const [inline]

get the HEP particle ID number

Definition at line 264 of file RawParticle.h.

References myId.

Referenced by RawParticleTypeFilter::addAccept(), RawParticleTypeFilter::addReject(), NuclearInteractionSimulator::compute(), ParticlePropagator::initProperDecayTime(), MaterialEffects::interact(), KineParticleFilter::isOKForMe(), RawParticleTypeFilter::isOKForMe(), operator<<(), Pythia6Decays::particleDaughters(), and RawParticleTypeFilter::RawParticleTypeFilter().

00264 { return myId; }

void RawParticle::print ( void   )  const

Print the formated particle information.

The format is: NAME______PX______PY______PZ______E_______Mtheo___Mrec____Qrec____X_______Y_______Z_______T_______

Definition at line 211 of file RawParticle.cc.

References charge(), GenMuonPlsPt100GeV_cfg::cout, lat::endl(), mass(), printName(), status(), T(), X(), Y(), and Z().

00211                          {
00212   printName();
00213   std::cout << std::setw(3) << status();
00214   std::cout.setf(std::ios::fixed, std::ios::floatfield);
00215   std::cout.setf(std::ios::right, std::ios::adjustfield);
00216   std::cout << std::setw(8) << std::setprecision(2) << Px();
00217   std::cout << std::setw(8) << std::setprecision(2) << Py();
00218   std::cout << std::setw(8) << std::setprecision(2) << Pz();
00219   std::cout << std::setw(8) << std::setprecision(2) << E();
00220   std::cout << std::setw(8) << std::setprecision(2) << M();
00221   std::cout << std::setw(8) << std::setprecision(2) << mass();
00222   std::cout << std::setw(8) << std::setprecision(2) << charge();
00223   std::cout << std::setw(8) << std::setprecision(2) << X();
00224   std::cout << std::setw(8) << std::setprecision(2) << Y();
00225   std::cout << std::setw(8) << std::setprecision(2) << Z();
00226   std::cout << std::setw(8) << std::setprecision(2) << T();
00227   std::cout << std::setw(0) << std::endl;
00228 }

void RawParticle::printName (  )  const

Print the name of the particle.

The name is deduced from the particle ID using a particle data table. It is printed with a length of 10 characters. If the id number cannot be found in the table "unknown" is printed as name.

Definition at line 199 of file RawParticle.cc.

References GenMuonPlsPt100GeV_cfg::cout, k, and PDGname().

Referenced by print().

00199                              {
00200   std::string MyParticleName = PDGname();
00201   if (MyParticleName.length() != 0) {
00202     std::cout <<  MyParticleName;
00203     for(unsigned int k=0;k<9-MyParticleName.length() && k<10; k++) 
00204       std::cout << " " ;
00205   } else {
00206     std::cout << "unknown  ";
00207   }
00208 }

double RawParticle::R (  )  const [inline]

vertex radius

Definition at line 277 of file RawParticle.h.

References R2(), and funct::sqrt().

Referenced by TrajectorySeedProducer::compatibleWithBeamAxis(), MaterialEffects::normalVector(), and MaterialEffects::radLengths().

00277 { return std::sqrt(R2()); }

double RawParticle::r (  )  const [inline]

vertex radius

Definition at line 279 of file RawParticle.h.

References r2(), and funct::sqrt().

Referenced by BaseParticlePropagator::propagateToBeamCylinder(), ParticlePropagator::propagateToBoundSurface(), rotate(), rotateX(), rotateY(), and rotateZ().

00279 { return std::sqrt(r2()); }

double RawParticle::R2 (  )  const [inline]

vertex radius**2

Definition at line 278 of file RawParticle.h.

References myVertex.

Referenced by BaseParticlePropagator::inside(), KineParticleFilter::isOKForMe(), BaseParticlePropagator::onBarrel(), BaseParticlePropagator::onEndcap(), BaseParticlePropagator::propagate(), BaseParticlePropagator::propagateToPreshowerLayer1(), BaseParticlePropagator::propagateToPreshowerLayer2(), and R().

00278 { return myVertex.Perp2(); }

double RawParticle::r2 (  )  const [inline]

vertex radius**2

Definition at line 280 of file RawParticle.h.

References myVertex.

Referenced by BaseParticlePropagator::propagateToBeamCylinder(), and r().

00280 { return myVertex.Perp2(); }

void RawParticle::reUse (  )  [inline]

Unlock the particle, see isUsed().

Definition at line 241 of file RawParticle.h.

References myUsed.

00241 { myUsed = 0;}  

void RawParticle::rotate ( const RotationZ r  )  [inline]

Definition at line 306 of file RawParticle.h.

References r(), and v.

00306                                                              { 
00307   XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
00308 }

void RawParticle::rotate ( const RotationY r  )  [inline]

Definition at line 302 of file RawParticle.h.

References r(), and v.

00302                                                              { 
00303   XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
00304 }

void RawParticle::rotate ( const RotationX r  )  [inline]

Definition at line 298 of file RawParticle.h.

References r(), and v.

00298                                                              { 
00299   XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
00300 }

void RawParticle::rotate ( const Rotation3D r  )  [inline]

Definition at line 290 of file RawParticle.h.

References r(), and v.

00290                                                               { 
00291   XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
00292 }

void RawParticle::rotate ( const Rotation r  )  [inline]

Definition at line 294 of file RawParticle.h.

References r(), and v.

00294                                                             { 
00295   XYZVector v ( r(Vect()) ); SetXYZT(v.X(),v.Y(),v.Z(),E());
00296 }

void RawParticle::rotate ( double  rphi,
const XYZVector raxis 
)

Rotate the particle around an axis in space.

The arguments give the amount to rotate rphi in radian and a vector raxis in 3D space around which the rotation is done. The vertex is rotated using the same transformation.

Definition at line 154 of file RawParticle.cc.

References r(), and v.

Referenced by NuclearInteractionSimulator::compute(), BremsstrahlungSimulator::compute(), and MultipleScatteringSimulator::compute().

00154                                                         {
00155   Rotation r(raxis,angle);
00156   XYZVector v(r * Vect());
00157   SetXYZT(v.X(),v.Y(),v.Z(),E());
00158 }

void RawParticle::rotateX ( double  rphi  ) 

Warning:
not yet implemented

Rotate around x axis. Rotate rphi radian around the x axis. The Vertex is rotated as well.

Definition at line 161 of file RawParticle.cc.

References r(), and v.

00161                                 {
00162   RotationX r(rphi);
00163   XYZVector v(r * Vect());
00164   SetXYZT(v.X(),v.Y(),v.Z(),E());
00165 }

void RawParticle::rotateY ( double  rphi  ) 

Rotate around z axis.

Rotate rphi radian around the z axis. The Vertex is rotated as well.

Definition at line 168 of file RawParticle.cc.

References r(), and v.

00168                                 {
00169   RotationY r(rphi);
00170   XYZVector v(r * Vect());
00171   SetXYZT(v.X(),v.Y(),v.Z(),E());
00172 }

void RawParticle::rotateZ ( double  rphi  ) 

Rotate around z axis.

Rotate rphi radian around the z axis. The Vertex is rotated as well.

Definition at line 175 of file RawParticle.cc.

References r(), and v.

00175                                 {
00176   RotationZ r(rphi);
00177   XYZVector v(r * Vect());
00178   SetXYZT(v.X(),v.Y(),v.Z(),E());
00179 }

void RawParticle::setCharge ( float  q  ) 

set the MEASURED charge

Definition at line 138 of file RawParticle.cc.

References myCharge.

Referenced by PFTrackTransformer::addPoints(), PFTrackTransformer::addPointsAndBrems(), BaseParticlePropagator::backPropagate(), FBaseSimEvent::fill(), ParticlePropagator::ParticlePropagator(), GoodSeedProducer::produce(), TrajectorySeedProducer::produce(), ConvBremSeedProducer::produce(), and BaseParticlePropagator::propagateToBeamCylinder().

00138                               {
00139   myCharge = q;
00140 }

void RawParticle::setID ( const std::string  name  ) 

Set identifier for this particle.

This should be a standard HEP-PID name. It will be used to deduce the particle properties from a particle data table.

Definition at line 114 of file RawParticle.cc.

References myCharge, myId, myInfo, myMass, tab, and ParticleTable::theTable().

00114                                        {
00115   if ( tab ) { 
00116     if ( !myInfo ) myInfo = tab->theTable()->particle(name);
00117     if ( myInfo ) { 
00118       myId = myInfo->pid();
00119       myCharge = myInfo->charge();
00120       myMass   = myInfo->mass().value();
00121     } else {
00122       myId = 0;
00123     }
00124   }
00125 }

void RawParticle::setID ( const int  id  ) 

Set identifier for this particle.

This should be a standard HEP-PID number. It will be used to deduce the name and the properties of the particle from a particle data table.

Definition at line 101 of file RawParticle.cc.

References myCharge, myId, myInfo, myMass, gen_jpsi2muons_cfg::ParticleID, tab, and ParticleTable::theTable().

Referenced by MuonSimHitProducer::applyMaterialEffects(), NuclearInteractionSimulator::compute(), PileUpSimulator::produce(), and RawParticle().

00101                                {
00102   myId = id;
00103   if ( tab ) {
00104     if ( !myInfo ) 
00105       myInfo = tab->theTable()->particle(HepPDT::ParticleID(myId));
00106     if ( myInfo ) { 
00107       myCharge = myInfo->charge();
00108       myMass   = myInfo->mass().value();
00109     }
00110   }
00111 }

void RawParticle::setMass ( float  m  ) 

set the RECONSTRUCTED mass

Definition at line 133 of file RawParticle.cc.

References myMass.

00133                             {
00134   myMass = m;
00135 }

void RawParticle::setStatus ( int  istat  ) 

Set the status of this particle.

The coding follows PYTHIAs convention: 1 = stable

Definition at line 128 of file RawParticle.cc.

References myStatus.

00128                                 {
00129   myStatus = istat;
00130 }

void RawParticle::setT ( const double  t  ) 

set the time of creation

Definition at line 149 of file RawParticle.cc.

References myVertex.

00149                                 {
00150   myVertex.SetE(t);
00151 }

void RawParticle::setVertex ( double  xv,
double  yv,
double  zv,
double  tv 
) [inline]

Definition at line 288 of file RawParticle.h.

References myVertex.

00288 { myVertex.SetXYZT(a,b,c,d); }

void RawParticle::setVertex ( const XYZTLorentzVector vtx  )  [inline]

set the vertex

Definition at line 287 of file RawParticle.h.

References myVertex.

Referenced by MuonSimHitProducer::applyMaterialEffects(), CalorimetryManager::EMShowerSimulation(), ParticlePropagator::ParticlePropagator(), BaseParticlePropagator::propagate(), BaseParticlePropagator::propagateToClosestApproach(), and GSPixelHitMatcher::propagateToLayer().

00287 { myVertex = vtx; }

int RawParticle::status ( void   )  const [inline]

get the particle status

Definition at line 265 of file RawParticle.h.

References myStatus.

Referenced by RawStableParticleFilter::isOKForMe(), operator<<(), and print().

00265 { return myStatus; }

double RawParticle::T (  )  const [inline]

vertex time

Definition at line 276 of file RawParticle.h.

References myVertex.

Referenced by Pythia6Decays::particleDaughters(), print(), BaseParticlePropagator::propagate(), and translate().

00276 { return myVertex.E(); }

double RawParticle::t (  )  const [inline]

vertex time

Definition at line 272 of file RawParticle.h.

References myVertex.

Referenced by FBaseSimEvent::addSimTrack().

00272 { return myVertex.E(); }

void RawParticle::translate ( const XYZVector t  )  [inline]

Translate the vertex by a given space amount.

Definition at line 314 of file RawParticle.h.

References myVertex, T(), X(), Y(), and Z().

Referenced by MultipleScatteringSimulator::compute().

00314                                                       { 
00315   myVertex.SetXYZT(X()+tr.X(),Y()+tr.Y(),Z()+tr.Z(),T());
00316 }

void RawParticle::use (  )  [inline]

Lock the particle, see isUsed().

Definition at line 237 of file RawParticle.h.

References myUsed.

00237 { myUsed = 1;}    

const XYZTLorentzVector & RawParticle::vertex (  )  const [inline]

the vertex fourvector

Definition at line 284 of file RawParticle.h.

References myVertex.

Referenced by PFTrackTransformer::addPoints(), PFTrackTransformer::addPointsAndBrems(), MuonSimHitProducer::applyMaterialEffects(), EcalHitMaker::cellLine(), CalorimetryManager::EMShowerSimulation(), FBaseSimEvent::fill(), ConvBremSeedProducer::GoodCluster(), EcalHitMaker::hcalCellLine(), CalorimetryManager::HDShowerSimulation(), MaterialEffects::interact(), KineParticleFilter::isOKForMe(), operator<<(), EcalHitMaker::preshowerCellLine(), GoodSeedProducer::produce(), TrajectoryManager::propagateToCalorimeters(), BaseParticlePropagator::propagateToClosestApproach(), GSPixelHitMatcher::propagateToLayer(), CalorimetryManager::reconstructECAL(), and CalorimetryManager::reconstructHCAL().

00284 { return myVertex ; }

double RawParticle::X (  )  const [inline]

x of vertex

Definition at line 273 of file RawParticle.h.

References myVertex.

Referenced by MuonSimHitProducer::applyMaterialEffects(), BaseParticlePropagator::helixCentreX(), TrajectoryManager::makeTrajectoryState(), ConvBremSeedProducer::makeTrajectoryState(), MaterialEffects::normalVector(), Pythia6Decays::particleDaughters(), ParticlePropagator::ParticlePropagator(), print(), BaseParticlePropagator::propagate(), BaseParticlePropagator::propagateToBeamCylinder(), BaseParticlePropagator::propagateToClosestApproach(), BaseParticlePropagator::propagateToNominalVertex(), translate(), BaseParticlePropagator::xyImpactParameter(), and BaseParticlePropagator::zImpactParameter().

00273 { return myVertex.Px(); }

double RawParticle::x (  )  const [inline]

x of vertex

Definition at line 269 of file RawParticle.h.

References myVertex.

Referenced by PFSimParticleProducer::produce().

00269 { return myVertex.Px(); }

double RawParticle::Y (  )  const [inline]

y of vertex

Definition at line 274 of file RawParticle.h.

References myVertex.

Referenced by MuonSimHitProducer::applyMaterialEffects(), BaseParticlePropagator::helixCentreY(), ConvBremSeedProducer::makeTrajectoryState(), TrajectoryManager::makeTrajectoryState(), MaterialEffects::normalVector(), Pythia6Decays::particleDaughters(), ParticlePropagator::ParticlePropagator(), print(), BaseParticlePropagator::propagate(), BaseParticlePropagator::propagateToBeamCylinder(), BaseParticlePropagator::propagateToClosestApproach(), BaseParticlePropagator::propagateToNominalVertex(), translate(), BaseParticlePropagator::xyImpactParameter(), and BaseParticlePropagator::zImpactParameter().

00274 { return myVertex.Py(); }

double RawParticle::y (  )  const [inline]

y of vertex

Definition at line 270 of file RawParticle.h.

References myVertex.

Referenced by PFSimParticleProducer::produce().

00270 { return myVertex.Py(); }

double RawParticle::Z (  )  const [inline]

z of vertex

Definition at line 275 of file RawParticle.h.

References myVertex.

Referenced by MuonSimHitProducer::applyMaterialEffects(), TrajectorySeedProducer::compatibleWithBeamAxis(), cos2ThetaV(), TrajectoryManager::createPSimHits(), BaseParticlePropagator::inside(), KineParticleFilter::isOKForMe(), TrajectoryManager::makeTrajectoryState(), ConvBremSeedProducer::makeTrajectoryState(), BaseParticlePropagator::onBarrel(), BaseParticlePropagator::onEndcap(), Pythia6Decays::particleDaughters(), ParticlePropagator::ParticlePropagator(), print(), ConvBremSeedProducer::produce(), BaseParticlePropagator::propagate(), BaseParticlePropagator::propagateToBeamCylinder(), MaterialEffects::radLengths(), translate(), and BaseParticlePropagator::zImpactParameter().

00275 { return myVertex.Pz(); }

double RawParticle::z (  )  const [inline]

z of vertex

Definition at line 271 of file RawParticle.h.

References myVertex.

Referenced by PFSimParticleProducer::produce(), ParticlePropagator::propagateToBoundSurface(), and BaseParticlePropagator::propagateToClosestApproach().

00271 { return myVertex.Pz(); }


Member Data Documentation

double RawParticle::myCharge [protected]

the MEASURED charge

Definition at line 253 of file RawParticle.h.

Referenced by charge(), chargeConjugate(), init(), operator=(), RawParticle(), setCharge(), and setID().

int RawParticle::myId [protected]

the particle id number HEP-PID

Definition at line 250 of file RawParticle.h.

Referenced by chargeConjugate(), init(), operator=(), PDGcTau(), pid(), RawParticle(), and setID().

const ParticleData* RawParticle::myInfo [protected]

The pointer to the PDG info.

Definition at line 255 of file RawParticle.h.

Referenced by init(), operator=(), PDGcharge(), PDGcTau(), PDGmass(), PDGname(), RawParticle(), and setID().

double RawParticle::myMass [protected]

the RECONSTRUCTED mass

Definition at line 254 of file RawParticle.h.

Referenced by init(), mass(), operator=(), RawParticle(), setID(), and setMass().

int RawParticle::myStatus [protected]

the status code according to PYTHIA

Definition at line 251 of file RawParticle.h.

Referenced by init(), operator=(), RawParticle(), setStatus(), and status().

int RawParticle::myUsed [protected]

status of the locking

Definition at line 252 of file RawParticle.h.

Referenced by init(), isUsed(), operator=(), RawParticle(), reUse(), and use().

XYZTLorentzVector RawParticle::myVertex [protected]

the four vector of the vertex

Definition at line 249 of file RawParticle.h.

Referenced by cos2ThetaV(), operator=(), R2(), r2(), RawParticle(), setT(), setVertex(), T(), t(), translate(), vertex(), x(), X(), y(), Y(), z(), and Z().

ParticleTable* RawParticle::tab [private]

Definition at line 258 of file RawParticle.h.

Referenced by init(), operator=(), PDGname(), RawParticle(), and setID().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:30:39 2009 for CMSSW by  doxygen 1.5.4