StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StarParticleFilter.cxx
1 #include "StarParticleFilter.h"
2 
3 #include "StarGenerator/EVENT/StarGenEvent.h"
4 #include "StarGenerator/EVENT/StarGenParticle.h"
5 
6 //______________________________________________________________________________________________
7 StarParticleFilter::StarParticleFilter( const char* name ) : StarFilterMaker(name)
8 {
9 
10 }
11 //______________________________________________________________________________________________
12 void StarParticleFilter::AddTrigger( int id, double mnpt, double mxpt, double mneta, double mxeta, int idm )
13 {
14  mTriggers.push_back( Trigger_t( id, mnpt, mxpt, mneta, mxeta, idm ) );
15 }
16 //______________________________________________________________________________________________
18 {
19  // Get a reference to the current event
20  StarGenEvent& event = (_event)? *_event : *mEvent;
21 
22  // Loop over tracks to find particles of interest
23  int npart = event.GetNumberOfParticles();
24  for ( int ipart=1 /*skip header*/;
25  ipart<npart;
26  ipart++ )
27  {
28  StarGenParticle *part = event[ipart];
29 
30  // Make sure we're looking at something with no color
31  if ( TMath::Abs(part->GetId()) < 10 ) continue;
32 
33  if ( part->GetStatus() != StarGenParticle::kFinal ) continue; // Must be a final state particle
34 
35  TLorentzVector momentum = part->momentum();
36  double pT = momentum.Perp();
37  if ( pT < 0.05 ) continue;
38  double eta = momentum.Eta();
39 
40  for ( auto trig : mTriggers )
41  {
42  if ( trig.pdgid != part->GetId() ) continue;
43  if ( trig.pdgid_parent ) {
44  int idx = part->GetFirstMother();
45  int pdgidm = (event[idx])? event[idx]->GetId() : 0;
46  if ( pdgidm != trig.pdgid_parent ) continue;
47  }
48  if ( pT < trig.ptmn ) continue;
49  if ( pT > trig.ptmx && trig.ptmn < trig.ptmx ) continue;
50  if ( eta < trig.etamn ) continue;
51  if ( eta > trig.etamx ) continue;
52 
53  return StarGenEvent::kAccept;
54  }
55 
56  }
57 
58  return StarGenEvent::kReject;
59 
60 }
int Filter(StarGenEvent *event=0)
Int_t GetNumberOfParticles()
Obtain the number of particles in the event record.
Definition: StarGenEvent.h:187
Yet another particle class.
Int_t GetStatus()
Get the status code of the particle according to the HEPEVT standard.
void AddTrigger(int _pdgid, double _ptmn=0, double _ptmx=-1, double _etamn=-1, double _etamx=1, int _pdgidParent=0)
Int_t GetId()
Get the id code of the particle according to the PDG standard.
Main filter class. Goes anywhere in the chain, filters StarGenEvent objects.
TLorentzVector momentum()
Return the 4-momentum of the particle.
Int_t GetFirstMother()
Get the first mother particle.
Base class for event records.
Definition: StarGenEvent.h:81