//*CMZ : 2.22/08 07/07/99 18.51.24 by Rene Brun
//*CMZ : 2.22/07 01/07/99 19.21.28 by Rene Brun
//*CMZ : 2.22/06 22/06/99 21.02.43 by Rene Brun
//*CMZ : 2.22/02 28/05/99 14.21.32 by Rene Brun
//*CMZ : 2.22/01 03/05/99 09.39.33 by Rene Brun
//*-- Author : Rene Brun , Federico Carminati 26/04/99
//*KEEP,TView.
#include "TView.h"
//*KEEP,TVirtualPad.
#include "TVirtualPad.h"
//*KEEP,TPolyLine3D.
#include "TPolyLine3D.h"
//*KEEP,TDatabasePDG,T=C++.
#include "TDatabasePDG.h"
//*KEEP,TParticle,T=C++.
#include "TParticle.h"
//*KEND.
ClassImp(TParticle)
//______________________________________________________________________________
TParticle::TParticle()
{
fParticlePDG = 0;
}
//______________________________________________________________________________
TParticle::TParticle(Int_t pdg, Int_t status,
Int_t mother1, Int_t mother2,
Int_t daughter1, Int_t daughter2,
Double_t px, Double_t py, Double_t pz, Double_t etot,
Double_t vx, Double_t vy, Double_t vz, Double_t time):
fPdgCode(pdg), fStatusCode(status), fWeight(1.),fPx(px), fPy(py),
fPz(pz), fE(etot), fVx(vx), fVy(vy), fVz(vz), fVt(time)
{
fMother[0] = mother1;
fMother[1] = mother2;
fDaughter[0] = daughter1;
fDaughter[1] = daughter2;
SetPolarisation(0,0,0);
fParticlePDG = TDatabasePDG::Instance()->GetParticle(pdg);
fCalcMass = fParticlePDG->Mass();
}
//______________________________________________________________________________
TParticle::TParticle(Int_t pdg, Int_t status,
Int_t mother1, Int_t mother2,
Int_t daughter1, Int_t daughter2,
const TLorentzVector &p,
const TLorentzVector &v) :
fPdgCode(pdg), fStatusCode(status), fWeight(1.),fPx(p.Px()), fPy(p.Py()),
fPz(p.Pz()), fE(p.E()), fVx(v.X()), fVy(v.Y()), fVz(v.Z()), fVt(v.T())
{
fMother[0] = mother1;
fMother[1] = mother2;
fDaughter[0] = daughter1;
fDaughter[1] = daughter2;
SetPolarisation(0,0,0);
fParticlePDG = TDatabasePDG::Instance()->GetParticle(pdg);
fCalcMass = fParticlePDG->Mass();
}
//______________________________________________________________________________
TParticle::~TParticle() {
}
//______________________________________________________________________________
Int_t TParticle::DistancetoPrimitive(Int_t px, Int_t py)
{
//*-*-*-*-*-*-*-*Compute distance from point px,py to a primary track*-*-*-*
//*-* ====================================================
//*-*
//*-* Compute the closest distance of approach from point px,py to each segment
//*-* of a track.
//*-* The distance is computed in pixels units.
//*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
const Int_t big = 9999;
Float_t xv[3], xe[3], xndc[3];
Float_t rmin[3], rmax[3];
TView *view = gPad->GetView();
if(!view) return big;
// compute first and last point in pad coordinates
Float_t pmom = this->P();
if (pmom == 0) return big;
view->GetRange(rmin,rmax);
Float_t rbox = rmax[2];
xv[0] = fVx;
xv[1] = fVy;
xv[2] = fVz;
xe[0] = xv[0]+rbox*fPx/pmom;
xe[1] = xv[1]+rbox*fPy/pmom;
xe[2] = xv[2]+rbox*fPz/pmom;
view->WCtoNDC(xv, xndc);
Float_t x1 = xndc[0];
Float_t y1 = xndc[1];
view->WCtoNDC(xe, xndc);
Float_t x2 = xndc[0];
Float_t y2 = xndc[1];
return DistancetoLine(px,py,x1,y1,x2,y2);
}
//______________________________________________________________________________
void TParticle::ExecuteEvent(Int_t, Int_t, Int_t)
{
//*-*-*-*-*-*-*-*-*-*-*Execute action corresponding to one event*-*-*-*
//*-* =========================================
gPad->SetCursor(kPointer);
}
//______________________________________________________________________________
const Text_t* TParticle::GetName() const {
static char def[4] = "XXX";
const TParticlePDG *ap = TDatabasePDG::Instance()->GetParticle(fPdgCode);
if (ap) return ap->GetName();
else return def;
}
//______________________________________________________________________________
void TParticle::GetPolarisation(TVector3 &v)
{
if(fPolarTheta == -99 && fPolarPhi == -99)
//No polarisation to return
v.SetXYZ(0.,0.,0.);
else
v.SetXYZ(TMath::Cos(fPolarPhi)*TMath::Sin(fPolarTheta),
TMath::Sin(fPolarPhi)*TMath::Sin(fPolarTheta),
TMath::Cos(fPolarTheta));
}
//______________________________________________________________________________
const Text_t *TParticle::GetTitle() const
{
static char def[4] = "XXX";
const TParticlePDG *ap = TDatabasePDG::Instance()->GetParticle(fPdgCode);
if (ap) return ap->GetTitle();
else return def;
}
//______________________________________________________________________________
void TParticle::Paint(Option_t *option)
{
//
// Paint a primary track
//
Float_t rmin[3], rmax[3];
static TPolyLine3D *pline = 0;
if (!pline) {
pline = new TPolyLine3D(2);
}
Float_t pmom = this->P();
if (pmom == 0) return;
TView *view = gPad->GetView();
if (!view) return;
view->GetRange(rmin,rmax);
Float_t rbox = rmax[2];
pline->SetPoint(0,Vx(), Vy(), Vz());
Float_t xend = Vx()+rbox*Px()/pmom;
Float_t yend = Vy()+rbox*Py()/pmom;
Float_t zend = Vz()+rbox*Pz()/pmom;
pline->SetPoint(1, xend, yend, zend);
pline->SetLineColor(GetLineColor());
pline->SetLineStyle(GetLineStyle());
pline->SetLineWidth(GetLineWidth());
pline->Paint(option);
}
//______________________________________________________________________________
void TParticle::Print(Option_t *)
{
//
// Print the internals of the primary vertex particle
//
TParticlePDG* pdg = this->GetPDG();
Printf("TParticle: %-13s p: %8f %8f %8f Vertex: %8e %8e %8e %5d %5d %s",
GetName(),Px(),Py(),Pz(),Vx(),Vy(),Vz(),
fMother[0],fMother[1],pdg->Type());
}
//______________________________________________________________________________
void TParticle::SetPolarisation(Double_t polx, Double_t poly, Double_t polz)
{
if(polx || poly || polz) {
fPolarTheta = TMath::ACos(polz/TMath::Sqrt(polx*polx+poly*poly+polz*polz));
fPolarPhi = kPI+TMath::ATan2(-poly,-polx);
} else {
fPolarTheta = -99;
fPolarPhi = -99;
}
}
//______________________________________________________________________________
void TParticle::Sizeof3D() const
{
//*-*-*-*-*-*Return total X3D size of this primary*-*-*-*-*-*-*
//*-* =====================================
Float_t pmom = this->P();
if (pmom == 0) return;
Int_t npoints = 2;
gSize3D.numPoints += npoints;
gSize3D.numSegs += (npoints-1);
gSize3D.numPolys += 0;
}
//______________________________________________________________________________
void TParticle::Streamer(TBuffer &R__b)
{
// Stream an object of class TParticle.
if (R__b.IsReading()) {
Version_t R__v = R__b.ReadVersion(); if (R__v) { }
TObject::Streamer(R__b);
TAttLine::Streamer(R__b);
R__b >> fPdgCode;
R__b >> fStatusCode;
R__b.ReadStaticArray(fMother);
R__b.ReadStaticArray(fDaughter);
R__b >> fWeight;
R__b >> fCalcMass;
R__b >> fPx;
R__b >> fPy;
R__b >> fPz;
R__b >> fE;
R__b >> fVx;
R__b >> fVy;
R__b >> fVz;
R__b >> fVt;
R__b >> fPolarTheta;
R__b >> fPolarPhi;
fParticlePDG = TDatabasePDG::Instance()->GetParticle(fPdgCode);
} else {
R__b.WriteVersion(TParticle::IsA());
TObject::Streamer(R__b);
TAttLine::Streamer(R__b);
R__b << fPdgCode;
R__b << fStatusCode;
R__b.WriteArray(fMother, 2);
R__b.WriteArray(fDaughter, 2);
R__b << fWeight;
R__b << fCalcMass;
R__b << fPx;
R__b << fPy;
R__b << fPz;
R__b << fE;
R__b << fVx;
R__b << fVy;
R__b << fVz;
R__b << fVt;
R__b << fPolarTheta;
R__b << fPolarPhi;
}
}
ROOT page - Class index - Top of the page
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.