PandaRoot
FairEvtFilter.h
Go to the documentation of this file.
1 //****************************************************************************
2 //* This file is part of PandaRoot. *
3 //* *
4 //* PandaRoot is distributed under the terms of the *
5 //* GNU General Public License (GPL) version 3, *
6 //* copied verbatim in the file "LICENSE". *
7 //* *
8 //* Copyright (C) 2006 - 2024 FAIR GmbH and copyright holders of PandaRoot *
9 //* The copyright holders are listed in the file "COPYRIGHTHOLDERS". *
10 //* The authors are listed in the file "AUTHORS". *
11 //****************************************************************************
12 
13 // -------------------------------------------------------------------------
14 // ----- FairEvtFilter header file -----
15 // -------------------------------------------------------------------------
16 
28 #ifndef FAIREVTFILTER_H
29 #define FAIREVTFILTER_H
30 
31 #include "TNamed.h" // for TNamed
32 #include "Rtypes.h" // for Bool_t, etc
33 #include "TClonesArray.h"
34 #include "TParticle.h"
35 #include <iostream>
36 #include <set>
37 #include "TDatabasePDG.h"
38 
39 std::ostream &operator<<(std::ostream &os, const std::set<Int_t> &set);
40 
41 class FairEvtFilter : public TNamed {
42 
43  public:
48 
50  FairEvtFilter();
51 
53  FairEvtFilter(const char *name, const char *title = "FairEvtFilter");
54 
56  virtual ~FairEvtFilter();
57 
58  // contains the list of particles that should be checked in the EventMatches method
59  Bool_t SetListOfParticles(TClonesArray *ParticleList)
60  {
61  fParticleList = ParticleList;
62  return kTRUE;
63  };
64 
65  // prints all TParticle in event
67 
68  // Initialize the filter if needed
69  Bool_t Init() { return kTRUE; }
70 
71  void SetVerbose(Int_t verbose = 12)
72  {
73  if (verbose >= 0) {
74  fVerbose = verbose;
75  std::cout << "FairEvtFilter: fVerbose is now set to " << verbose << "\n";
76  } else {
77  std::cout << "\n\n\n -WARNING from FairEvtFilter: verbose must be a positive number! Check your SetVerbose call!\n\n\n";
78  }
79  }
80 
82  {
83  // for QA
84  // shows fAcceptedEventNumbers that is filled after running a simulation
85  std::cout << "\n fAcceptedEventNumbers"
86  << " = " << fAcceptedEventNumbers;
87  }
88 
90  {
91  // for QA
92  // shows fEvtNrsToAccept that has to be set if you want to run your simulation in test mode
93  std::cout << "\n fEvtNrsToAccept"
94  << " = " << fEvtNrsToAccept;
95  }
96 
97  void SetTestMode(Int_t *arrayPtr, Int_t nEntries)
98  {
99  // turns on the test mode with the declared fEvtNrsToAccept
100  fTestMode = kTRUE;
101  std::set<Int_t> evtNrsToAccept(arrayPtr, arrayPtr + nEntries);
102  fEvtNrsToAccept = evtNrsToAccept;
103  }
104 
105  Bool_t TestPassed()
106  {
107  if (kFALSE == fTestMode) {
108  std::cout << "\n\n\n WARNING from FairEvtFilter: Test mode not set.\n\n\n";
109  return kFALSE;
110  }
112  // std::cout << "\n\n\n FairEvtFilter: Test passed.\n\n\n";
113  return kTRUE;
114  } else {
115  // std::cout << "\n\n\n FairEvtFilter: Test failed. Check your SetTestMode calls. \n\n\n";
116  return kFALSE;
117  }
118  }
119 
120  // returning kTRUE indicates that the event satisfies the filter conditions,
121  // kFALSE otherwise.
122  virtual Bool_t EventMatches(Int_t evtNr) = 0;
123 
124  virtual Bool_t FilterActive() = 0;
125 
126  // returns kTRUE if successful
127  // pdgCodeCharge will contain the charge of the particle with pdg code inPdgCode
128  Bool_t GetCharge(Int_t inPdgCode, Double_t *pdgCodeCharge);
129 
130  protected:
131  TDatabasePDG *fdbPdg;
132  // constant holding a double number which is not a valid charge
133  // this serves to indicate that the value has not been specified by the user
134  static const Double_t kNoChargeSpecified;
135 
136  TClonesArray *fParticleList; // list of particles in the event which was generated
137  Int_t fVerbose; // level of commenting output for your filter, between 0 and 12
138  Bool_t fTestMode; // is kTRUE if the filter should run in QA test mode
139  std::set<Int_t> fAcceptedEventNumbers; // set of event numbers which were accepted by the filter
140  std::set<Int_t> fEvtNrsToAccept; // event numbers that SHOULD be accepted by the filter (QA test mode)
141  Int_t fEventNr; // current event number
142 
143  private:
144  FairEvtFilter(const FairEvtFilter &G);
145  FairEvtFilter &operator=(const FairEvtFilter &) { return *this; }
146 
147  ClassDef(FairEvtFilter, 1);
148 };
149 
150 #endif
void PrintAllTParticleInEvent()
Bool_t GetCharge(Int_t inPdgCode, Double_t *pdgCodeCharge)
void SetVerbose(Int_t verbose=12)
Definition: FairEvtFilter.h:71
Bool_t Init()
Definition: FairEvtFilter.h:69
virtual Bool_t EventMatches(Int_t evtNr)=0
std::set< Int_t > fAcceptedEventNumbers
TClonesArray * fParticleList
void ShowEvtNrsToAccept()
Definition: FairEvtFilter.h:89
static const Double_t kNoChargeSpecified
void SetTestMode(Int_t *arrayPtr, Int_t nEntries)
Definition: FairEvtFilter.h:97
Bool_t TestPassed()
virtual ~FairEvtFilter()
std::set< Int_t > fEvtNrsToAccept
Bool_t SetListOfParticles(TClonesArray *ParticleList)
Definition: FairEvtFilter.h:59
virtual Bool_t FilterActive()=0
basic_ostream< char, char_traits< char > > ostream
void ShowAcceptedEventNumbers()
Definition: FairEvtFilter.h:81
TDatabasePDG * fdbPdg