FairRoot/PandaRoot
Public Types | Public Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
PndEvtFilterOnInvMassCounts Class Reference

#include <PndEvtFilterOnInvMassCounts.h>

Inheritance diagram for PndEvtFilterOnInvMassCounts:
PndEvtFilter FairEvtFilter

Public Types

enum  ChargeState {
  kNeutral =0, kPlus, kMinus, kCharged,
  kAll, kChargeLastElement
}
 
enum  MomState { kMomTotal =0, kMomTrans, kMomZ, kMomLastElement }
 
enum  GeomState {
  kTheta =0, kPhi, kVertexZ, kVertexRho,
  kVertexRadius, kGeomLastElement
}
 
enum  LogicOp { kAnd =0, kOr, kLogicOpLastElement }
 

Public Member Functions

 PndEvtFilterOnInvMassCounts ()
 
 PndEvtFilterOnInvMassCounts (const char *name, const char *title="PndEvtFilterOnInvMassCounts")
 
virtual ~PndEvtFilterOnInvMassCounts ()
 
Bool_t SetPdgCodesToCombine (Int_t pdgCode1, Int_t pdgCode2, Int_t pdgCode3=kInvalidPdgCode, Int_t pdgCode4=kInvalidPdgCode, Int_t pdgCode5=kInvalidPdgCode)
 
Bool_t SetMinMaxInvMass (Double_t min, Double_t max)
 
Bool_t SetRhoMassParticleSelector (const char *name, Double_t cv, Double_t w, const char *type)
 
Bool_t SetMinMaxCounts (Int_t min, Int_t max)
 
Bool_t SetMinCounts (Int_t min)
 
Bool_t SetMaxCounts (Int_t max)
 
Bool_t EventMatches (Int_t evtNr)
 
Bool_t FilterActive ()
 
Bool_t Init ()
 
Bool_t FillList (RhoCandList &rhoOutList, Int_t inPdgCode, Double_t pdgCodeCharge=kNoChargeSpecified)
 
Bool_t SetListOfParticles (TClonesArray *ParticleList)
 
void PrintAllTParticleInEvent ()
 
void SetVerbose (Int_t verbose=12)
 
void ShowAcceptedEventNumbers ()
 
void ShowEvtNrsToAccept ()
 
void SetTestMode (Int_t *arrayPtr, Int_t nEntries)
 
Bool_t TestPassed ()
 
Bool_t GetCharge (Int_t inPdgCode, Double_t *pdgCodeCharge)
 

Protected Attributes

TDatabasePDG * fdbPdg
 
TClonesArray * fParticleList
 
Int_t fVerbose
 
Bool_t fTestMode
 
std::set< Int_t > fAcceptedEventNumbers
 
std::set< Int_t > fEvtNrsToAccept
 
Int_t fEventNr
 

Static Protected Attributes

static const Double_t kNoChargeSpecified = -999.9
 

Private Member Functions

 ClassDef (PndEvtFilterOnInvMassCounts, 1)
 

Private Attributes

std::pair< Int_t, Int_t > fCountsMinMax
 
RhoMassParticleSelectorfInvMassSel
 
std::vector< std::pair< Int_t,
Double_t > > 
fPdgCodesCharges
 
Bool_t fInvMassRangeSet
 
Bool_t fPgdCodesSet
 
Bool_t fCountRangeSet
 

Static Private Attributes

static const Int_t kInvalidPdgCode = 0
 

Detailed Description

PndEvtFilterOnInvMassCounts.h

Author: Martin Galuska Martin dot J dot Galuska at physik dot uni minus giessen dot de

Purpose: This class is used to filter events right after event generation and before transport through the detector model events will either be accepted (and further simulated) or rejected completely.

This filter class implements such decisions based on multiplicities of invariant mass combinations within a certain mass region

The user has to specify the minimum or maximum multiplicities for the invariant mass combinations and which pdgCodes the particles should correspond to.

only particles which were produced by the event generators are taken into account for the combinations and which have the same charge as the charges corresponding to the input pdgCoded The mass hypotheses of matching particles will be set to the input pdgCodes given by the user

Definition at line 44 of file PndEvtFilterOnInvMassCounts.h.

Member Enumeration Documentation

enum FairEvtFilter::GeomState
inherited
enum FairEvtFilter::LogicOp
inherited
Enumerator
kAnd 
kOr 
kLogicOpLastElement 

Definition at line 42 of file FairEvtFilter.h.

enum FairEvtFilter::MomState
inherited
Enumerator
kMomTotal 
kMomTrans 
kMomZ 
kMomLastElement 

Definition at line 40 of file FairEvtFilter.h.

Constructor & Destructor Documentation

PndEvtFilterOnInvMassCounts::PndEvtFilterOnInvMassCounts ( )
PndEvtFilterOnInvMassCounts::PndEvtFilterOnInvMassCounts ( const char *  name,
const char *  title = "PndEvtFilterOnInvMassCounts" 
)
PndEvtFilterOnInvMassCounts::~PndEvtFilterOnInvMassCounts ( )
virtual

Definition at line 26 of file PndEvtFilterOnInvMassCounts.cxx.

26 {}

Member Function Documentation

PndEvtFilterOnInvMassCounts::ClassDef ( PndEvtFilterOnInvMassCounts  ,
 
)
private
Bool_t PndEvtFilterOnInvMassCounts::EventMatches ( Int_t  evtNr)
virtual

Implements PndEvtFilter.

Definition at line 169 of file PndEvtFilterOnInvMassCounts.cxx.

References RhoCandList::Combine(), FairEvtFilter::fAcceptedEventNumbers, fCountsMinMax, FairEvtFilter::fdbPdg, PndEvtFilter::FillList(), fInvMassSel, FairEvtFilter::fParticleList, fPdgCodesCharges, FairEvtFilter::fVerbose, RhoCandList::GetLength(), i, p1, p2, FairEvtFilter::PrintAllTParticleInEvent(), and RhoCandList::Select().

170 {
171  if (fVerbose > 3){
172  std::cout << "\n\n";
173  std::cout << "Internal eventNr:" << evtNr << "\n";
174  std::cout << "Nr. of simulated particles:" << fParticleList->GetEntries()<<"\n";
176  std::cout << "PndEvtFilterOnInvMassCounts " << this->GetTitle() << ": " << this->GetName() << ": Beginning of EventMatches\n";
177  }
178 
179  // sanity checks
180  if (0==fParticleList){
181  std::cout << "\n\n\n PndEvtFilterOnInvMassCounts " << this->GetTitle() << ": " << this->GetName() << ": FATAL ERROR: No particle list! Discard this event.\n\n\n";
182  return kFALSE;
183  }
184 
185 
186  // do not accept events that do not have enough particles
187  if ( fParticleList->GetEntriesFast() < (int)fPdgCodesCharges.size() ){
188  if (0<fVerbose) std::cout << "\n\n\n PndEvtFilterOnInvMassCounts: Event contains less than " << fPdgCodesCharges.size() << " particles. " << this->GetTitle() << ": " << this->GetName() << " will not accept this event.\n\n\n";
189  return kFALSE;
190  }
191 
192  // get rho cand lists of particles and combine them
193  RhoCandList p0, p1, p2, p3, p4;
194  // also care about charged conjugate
195  RhoCandList ap0, ap1, ap2, ap3, ap4;
196 
197  std::vector<int> pdg, apdg;
198 
199  // make list with anti-pdg codes
200  for (UInt_t i=0;i<fPdgCodesCharges.size();++i)
201  {
202  int cpdg = fPdgCodesCharges[i].first;
203  pdg.push_back(cpdg);
204 
205  if (fdbPdg->GetParticle(cpdg)!=0x0 && fdbPdg->GetParticle(cpdg)->AntiParticle()!=0x0)
206  apdg.push_back(fdbPdg->GetParticle(cpdg)->AntiParticle()->PdgCode());
207  else
208  apdg.push_back(cpdg);
209  }
210 
211  // is the cc'd list the same as the original one (= a permutation)?
212  bool isperm = is_permutation(pdg.begin(), pdg.end(), apdg.begin());
213 
214 
215  FillList( p0, fPdgCodesCharges[0].first, fPdgCodesCharges[0].second );
216  FillList( p1, fPdgCodesCharges[1].first, fPdgCodesCharges[1].second );
217 
218  // fill anti-particle lists
219  if (!isperm)
220  {
221  FillList( ap0, apdg[0], -fPdgCodesCharges[0].second );
222  FillList( ap1, apdg[1], -fPdgCodesCharges[1].second );
223  }
224 
225  RhoCandList combinedList, combinedAntiList;
226 
227  switch( fPdgCodesCharges.size() ) {
228  case 2:
229  combinedList.Combine(p0, p1);
230  if (!isperm) combinedAntiList.Combine(ap0, ap1);
231  break;
232  case 3:
233  FillList( p2, fPdgCodesCharges[2].first, fPdgCodesCharges[2].second );
234  combinedList.Combine(p0, p1, p2);
235  if (!isperm)
236  {
237  FillList( ap2, apdg[2], -fPdgCodesCharges[2].second );
238  combinedAntiList.Combine(ap0, ap1, ap2);
239  }
240  break;
241  case 4:
242  FillList( p2, fPdgCodesCharges[2].first, fPdgCodesCharges[2].second );
243  FillList( p3, fPdgCodesCharges[3].first, fPdgCodesCharges[3].second );
244  combinedList.Combine(p0, p1, p2, p3);
245  if (!isperm)
246  {
247  FillList( ap2, apdg[2], -fPdgCodesCharges[2].second );
248  FillList( ap3, apdg[3], -fPdgCodesCharges[3].second );
249  combinedAntiList.Combine(ap0, ap1, ap2, ap3);
250  }
251  break;
252  case 5:
253  FillList( p2, fPdgCodesCharges[2].first, fPdgCodesCharges[2].second );
254  FillList( p3, fPdgCodesCharges[3].first, fPdgCodesCharges[3].second );
255  FillList( p4, fPdgCodesCharges[4].first, fPdgCodesCharges[4].second );
256  combinedList.Combine(p0, p1, p2, p3, p4);
257  if (!isperm)
258  {
259  FillList( p2, apdg[2], -fPdgCodesCharges[2].second );
260  FillList( p3, apdg[3], -fPdgCodesCharges[3].second );
261  FillList( p4, apdg[4], -fPdgCodesCharges[4].second );
262  combinedAntiList.Combine(ap0, ap1, ap2, ap3, ap4);
263  }
264  break;
265  default:
266  std::cout << "FATAL ERROR in PndEvtFilterOnInvMassCounts::EventMatches of " << this->GetTitle() << ": " << this->GetName() << ". The number of pdg codes is not supported. \n";
267  return kFALSE;
268  break;
269  }
270 
271 
272  if (fVerbose > 3){
273  std::cout << "Looking for inv. masses of " << fPdgCodesCharges.size() << " particles.\n";
274  std::cout << p0.GetLength() << " / " << p1.GetLength() << " / " << p2.GetLength() << " / " << p3.GetLength() << " / " << p4.GetLength();
275  std::cout << " particles in initial lists\n";
276  std::cout << combinedList.GetLength() << " particles in combined list before mass selection\n";
277  }
278 
279  // select acceptable invariant mass range
280  combinedList.Select(fInvMassSel);
281  if (!isperm) combinedAntiList.Select(fInvMassSel);
282 
283  int ncand = combinedList.GetLength() + combinedAntiList.GetLength();
284 
285  if (fVerbose > 3){
286  std::cout << combinedList.GetLength()<<" + "<< combinedAntiList.GetLength()<< " particles in combined list after mass selection\n";
287  }
288 
289  // check if number of acceptable combinations is within acceptable limits
290  if ( ncand < fCountsMinMax.first){
291  if (fVerbose >9){
292  std::cout << "Event is not accepted by " << this->GetTitle() << ": " << this->GetName() << " because there are too few matching inv. mass combinations in the event. \n";
293  }
294  return kFALSE;
295  }
296  if ( ncand > fCountsMinMax.second ){
297  if (fVerbose >9){
298  std::cout << "Event is not accepted by " << this->GetTitle() << ": " << this->GetName() << " because there are too many matching inv. mass combinations in the event. \n";
299  }
300  return kFALSE;
301  }
302 
303 
304  // for QA
305  fAcceptedEventNumbers.insert(evtNr);
306  if (fVerbose >5){
307  std::cout << "\n Event matches " << this->GetTitle() << ": " << this->GetName() << "\n\n";
308  }
309  return kTRUE;
310 }
void PrintAllTParticleInEvent()
Int_t i
Definition: run_full.C:25
Bool_t FillList(RhoCandList &rhoOutList, Int_t inPdgCode, Double_t pdgCodeCharge=kNoChargeSpecified)
Int_t GetLength() const
Definition: RhoCandList.h:46
void Combine(RhoCandList &l1, RhoCandList &l2)
RhoMassParticleSelector * fInvMassSel
std::set< Int_t > fAcceptedEventNumbers
TClonesArray * fParticleList
void Select(RhoParticleSelectorBase *pidmgr)
TPad * p2
Definition: hist-t7.C:117
std::pair< Int_t, Int_t > fCountsMinMax
TPad * p1
Definition: hist-t7.C:116
std::vector< std::pair< Int_t, Double_t > > fPdgCodesCharges
TDatabasePDG * fdbPdg
Bool_t PndEvtFilter::FillList ( RhoCandList rhoOutList,
Int_t  inPdgCode,
Double_t  pdgCodeCharge = kNoChargeSpecified 
)
inherited

Definition at line 31 of file PndEvtFilter.cxx.

References RhoCandList::Add(), RhoCandList::Cleanup(), Double_t, FairEvtFilter::fdbPdg, FairEvtFilter::fParticleList, FairEvtFilter::GetCharge(), FairEvtFilter::kNoChargeSpecified, particle, RhoCandidate::SetMcTruth(), RhoCandidate::SetPos(), and RhoCandidate::SetType().

Referenced by EventMatches().

32 {
33  rhoOutList.Cleanup();
34 
35 
36  // Get charge corresponding to pdgCode if it has not been passed as an argument by the user
37  if ( kNoChargeSpecified == pdgCodeCharge ){
38  if ( kFALSE == GetCharge( inPdgCode, &pdgCodeCharge ) ) {
39  return kFALSE;
40  }
41  }
42 
43 
44 
45  // find all TParticle on stack with a charge(particle) == charge(pdgCode)
46  for (Int_t iPart=0; iPart<fParticleList->GetEntries(); ++iPart) {
47  TParticle *particle = (TParticle*)fParticleList->At(iPart);
48 
49  // get charge for particle
50  TParticlePDG* pdt = particle->GetPDG();
51  if (0==pdt) continue; // unknown particle type (KG, 10/2015)
52 
53  Double_t pCharge = pdt->Charge()/3.; // TParticlePDG contains charge in units of |e|/3
54 
55  if ( pdgCodeCharge != pCharge ){ continue; } // skip all particles with different charge
56 
57 
58  // // mark stable particles
59  // Bool_t isStable = kFALSE;
60  // switch(abs(particle->GetPdgCode())) {
61  // case 22: isStable = true; break;
62  // case 11: isStable = true; break;
63  // case 13: isStable = true; break;
64  // case 211: isStable = true; break;
65  // case 321: isStable = true; break;
66  // case 2212: isStable = true; break;
67  // }
68  //
69  // if (kFALSE == isStable) { continue; } // skip all particles which are not considered stable
70 
71 
72  // store TParticle as RhoCandidate in output TCA and modify 4 momentum to match particle hypothesis
73  TVector3 pVertex(particle->Vx(),particle->Vy(),particle->Vz());
74  TLorentzVector p4(particle->Px(),particle->Py(),particle->Pz(),particle->Energy());
75  // assign mass hypothesis according to pdgCode
76  p4.SetVectM(p4.Vect(),fdbPdg->GetParticle(inPdgCode)->Mass());
77  // construct RhoCandidate
78  RhoCandidate rhoCand(p4,pCharge);
79  rhoCand.SetMcTruth(&rhoCand);
80  rhoCand.SetPos(pVertex);
81  rhoCand.SetType(inPdgCode);
82  // save to output
83  rhoOutList.Add(&rhoCand);
84  }//trackloop
85 
86 
87  return kTRUE;
88 }
void Add(const RhoCandidate *c)
Definition: RhoCandList.h:49
void Cleanup()
Definition: RhoCandList.cxx:62
Bool_t GetCharge(Int_t inPdgCode, Double_t *pdgCodeCharge)
const int particle
Double_t
TClonesArray * fParticleList
static const Double_t kNoChargeSpecified
TDatabasePDG * fdbPdg
Bool_t PndEvtFilterOnInvMassCounts::FilterActive ( )
inlinevirtual
Bool_t FairEvtFilter::GetCharge ( Int_t  inPdgCode,
Double_t pdgCodeCharge 
)
inherited

Definition at line 62 of file FairEvtFilter.cxx.

References FairEvtFilter::fdbPdg, FairEvtFilter::fVerbose, and FairEvtFilter::kNoChargeSpecified.

Referenced by FairEvtFilterOnSingleParticleCounts::CountCharge(), PndEvtFilter::FillList(), and SetPdgCodesToCombine().

63 {
64  // Try to find the pdg code
65  TParticlePDG *ptrToPdg = fdbPdg->GetParticle(inPdgCode);
66  if ( 0 == ptrToPdg) {
67  // ignore particles with unknown charges
68  std::cout << "WARNING from FairEvtFilter::GetCharge Charge of pdgCode " << inPdgCode << " is unknown and will be ignored!\n";
69  *pdgCodeCharge = kNoChargeSpecified;
70  return kFALSE;
71  }
72  *pdgCodeCharge = ptrToPdg->Charge()/3.; // TParticlePDG contains charge in units of |e|/3
73  if ( fVerbose > 1 ) std::cout << "Found pdgCodeCharge = " << *pdgCodeCharge << " for inPdgCode " << inPdgCode << '\n';
74  return kTRUE;
75 }
static const Double_t kNoChargeSpecified
TDatabasePDG * fdbPdg
Bool_t PndEvtFilter::Init ( )
inlineinherited

Definition at line 54 of file PndEvtFilter.h.

54 { return kTRUE;}
void FairEvtFilter::PrintAllTParticleInEvent ( )
inherited

Definition at line 52 of file FairEvtFilter.cxx.

References FairEvtFilter::fParticleList, and particle.

Referenced by EventMatches(), and FairEvtFilterOnSingleParticleCounts::EventMatches().

52  {
53  // sanity checks
54  if (0==fParticleList){ return; }
55  for (Int_t iPart=0; iPart<fParticleList->GetEntries(); ++iPart) {
56  TParticle *particle = (TParticle*)fParticleList->At(iPart);
57  particle->Print();
58  }
59 }
const int particle
TClonesArray * fParticleList
Bool_t FairEvtFilter::SetListOfParticles ( TClonesArray *  ParticleList)
inlineinherited

Definition at line 57 of file FairEvtFilter.h.

References FairEvtFilter::fParticleList.

Referenced by FairFilteredPrimaryGenerator::GenerateEvent().

57 {fParticleList=ParticleList; return kTRUE;};
TClonesArray * fParticleList
Bool_t PndEvtFilterOnInvMassCounts::SetMaxCounts ( Int_t  max)
inline

Definition at line 69 of file PndEvtFilterOnInvMassCounts.h.

References SetMinMaxCounts().

69 { return SetMinMaxCounts(0, max); };
Bool_t SetMinMaxCounts(Int_t min, Int_t max)
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:26
Bool_t PndEvtFilterOnInvMassCounts::SetMinCounts ( Int_t  min)
inline

Definition at line 68 of file PndEvtFilterOnInvMassCounts.h.

References SetMinMaxCounts().

68 { return SetMinMaxCounts(min, 9999); };
Bool_t SetMinMaxCounts(Int_t min, Int_t max)
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:25
Bool_t PndEvtFilterOnInvMassCounts::SetMinMaxCounts ( Int_t  min,
Int_t  max 
)

Definition at line 29 of file PndEvtFilterOnInvMassCounts.cxx.

References fCountRangeSet, fCountsMinMax, max(), and min().

Referenced by prod_fsim(), quickfsimana(), SetMaxCounts(), SetMinCounts(), and sim_filter_inv_mass().

30 {
31  if ( kTRUE == fCountRangeSet ){
32  std::cout << "\n\n\n -WARNING from PndEvtFilterOnInvMassCounts of " << this->GetTitle() << ": " << this->GetName() << ": \n";
33  std::cout << "You are trying to set min and max counts for the events more than once! \n";
34  std::cout << "I take the first setting of min = " << fCountsMinMax.first << " max = " << fCountsMinMax.second << " and ignore the later ones! \n\n\n";
35  std::cout << "Your current input of min = " << min << " max = " << max << " is ignored! \n\n\n";
36  return kFALSE;
37  }
38 
39  // check if filter can be applied (make sure that min and max make sense)
40  if ( min < 0 ){
41  std::cout << "\n\n\n -WARNING from PndEvtFilterOnInvMassCounts of " << this->GetTitle() << ": " << this->GetName() << ": Filter could not be added. You are trying to request that your events should have at least a number <= 0 of some particles. That makes no sense. \n\n\n";
42  return kFALSE;
43  }
44  if ( max < min ){
45  std::cout << "\n\n\n -WARNING from PndEvtFilterOnInvMassCounts of " << this->GetTitle() << ": " << this->GetName() << ": Filter could not be added. You are trying to request that your events should have at most a number of some particles which is less than the minimum number that you request. That makes no sense. \n\n\n";
46  return kFALSE;
47  }
48 
49  // set the minimum and maximum number of particle combinations
50  fCountsMinMax = std::pair<Int_t, Int_t> (min,max);
51 
52  // turn filter on
53  fCountRangeSet = kTRUE;
54 
55  return kTRUE;
56 }
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:26
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:25
std::pair< Int_t, Int_t > fCountsMinMax
Bool_t PndEvtFilterOnInvMassCounts::SetMinMaxInvMass ( Double_t  min,
Double_t  max 
)

Definition at line 59 of file PndEvtFilterOnInvMassCounts.cxx.

References Double_t, fInvMassRangeSet, fInvMassSel, max(), and min().

Referenced by prod_fsim(), quickfsimana(), and sim_filter_inv_mass().

60 {
61  if ( kTRUE == fInvMassRangeSet ){
62  std::cout << "\n\n\n -WARNING from PndEvtFilterOnInvMassCounts of " << this->GetTitle() << ": " << this->GetName() << ": \n";
63  std::cout << "You are trying to set the limits for the invariant mass selection set more than once! \n";
64  std::cout << "I take the first setting and ignore the later ones! \n\n\n";
65  return kFALSE;
66  }
67 
68  // check if filter can be applied (make sure that min and max make sense)
69  if ( min < 0. ){
70  std::cout << "\n\n\n -WARNING from PndEvtFilterOnInvMassCounts " << this->GetTitle() << ": " << this->GetName() << ": Filter could not be added. Min inv. mass too low. \n\n\n";
71  return kFALSE;
72  }
73  if ( max < min ){
74  std::cout << "\n\n\n -WARNING from PndEvtFilterOnInvMassCounts " << this->GetTitle() << ": " << this->GetName() << ": Filter could not be added. Max inv. mass lower than min inv. mass. \n\n\n";
75  return kFALSE;
76  }
77 
78  // set up the inv. mass selector
79  Double_t mu = (min + max)/2.0;
80  Double_t fullWidth = max - min;
81 
82  fInvMassSel = new RhoMassParticleSelector("invMassFilter", mu, fullWidth);
83 
84  // turn filter on
85  fInvMassRangeSet = kTRUE;
86 
87  return kTRUE;
88 }
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:26
Double_t
RhoMassParticleSelector * fInvMassSel
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:25
Bool_t PndEvtFilterOnInvMassCounts::SetPdgCodesToCombine ( Int_t  pdgCode1,
Int_t  pdgCode2,
Int_t  pdgCode3 = kInvalidPdgCode,
Int_t  pdgCode4 = kInvalidPdgCode,
Int_t  pdgCode5 = kInvalidPdgCode 
)

Definition at line 119 of file PndEvtFilterOnInvMassCounts.cxx.

References Double_t, fPdgCodesCharges, fPgdCodesSet, FairEvtFilter::GetCharge(), kInvalidPdgCode, and FairEvtFilter::kNoChargeSpecified.

Referenced by prod_fsim(), quickfsimana(), and sim_filter_inv_mass().

119  {
120  if ( kTRUE == fPgdCodesSet ){
121  std::cout << "\n\n\n -WARNING from PndEvtFilterOnInvMassCounts of " << this->GetTitle() << ": " << this->GetName() << ": \n";
122  std::cout << "You are trying to set the pdgCodes to combine into invariant mass combinations more than once! \n";
123  std::cout << "I take the first setting and ignore the later ones! \n\n\n";
124  return kFALSE;
125  }
126 
127 
128  // create temp. vector of all input pdg codes
129  std::vector< Int_t > pdgCodes, apdgCodes;
130  pdgCodes.push_back(pdgCode1);
131  pdgCodes.push_back(pdgCode2);
132  pdgCodes.push_back(pdgCode3);
133  pdgCodes.push_back(pdgCode4);
134  pdgCodes.push_back(pdgCode5);
135 
136 
137 // std::cout <<"[PndEvtFilterOnInvMassCounts]"<<pdgCode1<<" "<<pdgCode2<<" "<<pdgCode3<<" "<<pdgCode4<<" "<<pdgCode5<<std::endl;
138 
139  // Copy valid entries from pdgCodes to fPdgCodesCharges
140  for (UInt_t iPdgCodes = 0; iPdgCodes < pdgCodes.size(); ++iPdgCodes){
141  // skip kInvalidPdgCode entries in pdgCodes
142  if ( kInvalidPdgCode == pdgCodes[iPdgCodes] ) { continue; }
143  // If we have invalid pdg codes warn user
144  Double_t pdgCodeCharge = kNoChargeSpecified;
145  if ( kFALSE == GetCharge( pdgCodes[iPdgCodes], &pdgCodeCharge ) ) {
146  std::cout << "\n\n\nFATAL ERROR FROM PndEvtFilterOnInvMassCounts::SetPdgCodesToCombine of " << this->GetTitle() << ": " << this->GetName() << " Charge for pdg code " << pdgCodes[iPdgCodes] << " was not found. \n";
147  std::cout << "Filter was not set. \n\n\n”";
148  fPgdCodesSet = kFALSE;
149  return kFALSE;
150  }
151  // If we get to this point, we have successfully added a pdg code
152  fPdgCodesCharges.push_back(std::pair<Int_t, Double_t> (pdgCodes[iPdgCodes], pdgCodeCharge));
153  }
154 
155  // Check that we had at least two correctly set pdg codes
156  if ( 1 < fPdgCodesCharges.size() ){
157  fPgdCodesSet = kTRUE;
158  return kTRUE;
159  } else {
160  std::cout << "\n\n\nFATAL ERROR FROM PndEvtFilterOnInvMassCounts::SetPdgCodesToCombine of " << this->GetTitle() << ": " << this->GetName() << " You have to specify at least 2 valid pdg codes. \n";
161  std::cout << "Filter was not set. \n\n\n”";
162  fPgdCodesSet = kFALSE;
163  return kFALSE;
164  }
165 }
Bool_t GetCharge(Int_t inPdgCode, Double_t *pdgCodeCharge)
Double_t
static const Double_t kNoChargeSpecified
std::vector< std::pair< Int_t, Double_t > > fPdgCodesCharges
Bool_t PndEvtFilterOnInvMassCounts::SetRhoMassParticleSelector ( const char *  name,
Double_t  cv,
Double_t  w,
const char *  type 
)

Definition at line 91 of file PndEvtFilterOnInvMassCounts.cxx.

References fInvMassRangeSet, and fInvMassSel.

92 {
93  if ( kTRUE == fInvMassRangeSet ){
94  std::cout << "\n\n\n -WARNING from PndEvtFilterOnInvMassCounts of " << this->GetTitle() << ": " << this->GetName() << ": \n";
95  std::cout << "You are trying to set the limits for the invariant mass selection set more than once! \n";
96  std::cout << "I take the first setting and ignore the later ones! \n\n\n";
97  return kFALSE;
98  }
99 
100  // check if filter can be applied (make sure that min and max make sense)
101  if ( cv <= 0. ){
102  std::cout << "\n\n\n -WARNING from PndEvtFilterOnInvMassCounts " << this->GetTitle() << ": " << this->GetName() << ": Filter could not be added. cv has to be positive. \n\n\n";
103  return kFALSE;
104  }
105  if ( w < 0. ){
106  std::cout << "\n\n\n -WARNING from PndEvtFilterOnInvMassCounts " << this->GetTitle() << ": " << this->GetName() << ": Filter could not be added. w has to be non-negative. \n\n\n";
107  return kFALSE;
108  }
109 
110  fInvMassSel = new RhoMassParticleSelector(name, cv, w, type);
111 
112  // turn filter on
113  fInvMassRangeSet = kTRUE;
114 
115  return kTRUE;
116 }
RhoMassParticleSelector * fInvMassSel
TString name
void FairEvtFilter::SetTestMode ( Int_t *  arrayPtr,
Int_t  nEntries 
)
inlineinherited

Definition at line 91 of file FairEvtFilter.h.

References FairEvtFilter::fEvtNrsToAccept, and FairEvtFilter::fTestMode.

91  {
92  //turns on the test mode with the declared fEvtNrsToAccept
93  fTestMode=kTRUE;
94  std::set<Int_t> evtNrsToAccept (arrayPtr,arrayPtr+nEntries);
95  fEvtNrsToAccept=evtNrsToAccept;
96  }
std::set< Int_t > fEvtNrsToAccept
void FairEvtFilter::SetVerbose ( Int_t  verbose = 12)
inlineinherited

Definition at line 67 of file FairEvtFilter.h.

References FairEvtFilter::fVerbose, and verbose.

Referenced by sim_filter_inv_mass().

67  {
68  if ( verbose >= 0 ){
69  fVerbose = verbose;
70  std::cout << "FairEvtFilter: fVerbose is now set to " << verbose << "\n";
71  } else {
72  std::cout << "\n\n\n -WARNING from FairEvtFilter: verbose must be a positive number! Check your SetVerbose call!\n\n\n";
73  }
74  }
#define verbose
void FairEvtFilter::ShowAcceptedEventNumbers ( )
inlineinherited

Definition at line 77 of file FairEvtFilter.h.

References FairEvtFilter::fAcceptedEventNumbers.

77  {
78  // for QA
79  // shows fAcceptedEventNumbers that is filled after running a simulation
80  std::cout << "\n fAcceptedEventNumbers" << " = " << fAcceptedEventNumbers;
81  }
std::set< Int_t > fAcceptedEventNumbers
void FairEvtFilter::ShowEvtNrsToAccept ( )
inlineinherited

Definition at line 84 of file FairEvtFilter.h.

References FairEvtFilter::fEvtNrsToAccept.

84  {
85  // for QA
86  //shows fEvtNrsToAccept that has to be set if you want to run your simulation in test mode
87  std::cout << "\n fEvtNrsToAccept" << " = " << fEvtNrsToAccept;
88  }
std::set< Int_t > fEvtNrsToAccept
Bool_t FairEvtFilter::TestPassed ( )
inlineinherited

Definition at line 99 of file FairEvtFilter.h.

References FairEvtFilter::fAcceptedEventNumbers, FairEvtFilter::fEvtNrsToAccept, and FairEvtFilter::fTestMode.

99  {
100  if(kFALSE==fTestMode){
101  std::cout << "\n\n\n WARNING from FairEvtFilter: Test mode not set.\n\n\n";
102  return kFALSE;
103  }
105  //std::cout << "\n\n\n FairEvtFilter: Test passed.\n\n\n";
106  return kTRUE;
107  }else{
108  //std::cout << "\n\n\n FairEvtFilter: Test failed. Check your SetTestMode calls. \n\n\n";
109  return kFALSE;
110  }
111  }
std::set< Int_t > fAcceptedEventNumbers
std::set< Int_t > fEvtNrsToAccept

Member Data Documentation

std::set<Int_t> FairEvtFilter::fAcceptedEventNumbers
protectedinherited
Bool_t PndEvtFilterOnInvMassCounts::fCountRangeSet
private

Definition at line 95 of file PndEvtFilterOnInvMassCounts.h.

Referenced by FilterActive(), and SetMinMaxCounts().

std::pair<Int_t,Int_t> PndEvtFilterOnInvMassCounts::fCountsMinMax
private

Definition at line 80 of file PndEvtFilterOnInvMassCounts.h.

Referenced by EventMatches(), and SetMinMaxCounts().

TDatabasePDG* FairEvtFilter::fdbPdg
protectedinherited
Int_t FairEvtFilter::fEventNr
protectedinherited

Definition at line 136 of file FairEvtFilter.h.

std::set<Int_t> FairEvtFilter::fEvtNrsToAccept
protectedinherited
Bool_t PndEvtFilterOnInvMassCounts::fInvMassRangeSet
private
RhoMassParticleSelector* PndEvtFilterOnInvMassCounts::fInvMassSel
private
TClonesArray* FairEvtFilter::fParticleList
protectedinherited
std::vector< std::pair< Int_t, Double_t > > PndEvtFilterOnInvMassCounts::fPdgCodesCharges
private

Definition at line 89 of file PndEvtFilterOnInvMassCounts.h.

Referenced by EventMatches(), and SetPdgCodesToCombine().

Bool_t PndEvtFilterOnInvMassCounts::fPgdCodesSet
private

Definition at line 94 of file PndEvtFilterOnInvMassCounts.h.

Referenced by FilterActive(), and SetPdgCodesToCombine().

Bool_t FairEvtFilter::fTestMode
protectedinherited

Definition at line 133 of file FairEvtFilter.h.

Referenced by FairEvtFilter::SetTestMode(), and FairEvtFilter::TestPassed().

Int_t FairEvtFilter::fVerbose
protectedinherited
const Int_t PndEvtFilterOnInvMassCounts::kInvalidPdgCode = 0
staticprivate

Definition at line 100 of file PndEvtFilterOnInvMassCounts.h.

Referenced by SetPdgCodesToCombine().

const Double_t FairEvtFilter::kNoChargeSpecified = -999.9
staticprotectedinherited

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