FairRoot/PandaRoot
FairFilteredPrimaryGenerator.h
Go to the documentation of this file.
1 
22 #ifndef FairFilteredPrimaryGenerator_H
23 #define FairFilteredPrimaryGenerator_H
24 
25 #include "FairPrimaryGenerator.h"
26 
27 #include "FairEvtFilter.h"
28 #include "FairEvtFilterParams.h"
29 
30 #include "TNamed.h" // for TNamed
31 #include "TFile.h" // for TFile
32 #include "FairRunSim.h"
33 
34 #include "FairGenerator.h" // for FairGenerator
35 
36 #include <iosfwd> // for ostream
37 #include "Rtypes.h" // for Double_t, Bool_t, Int_t, etc
38 #include "TObjArray.h" // for TObjArray
39 #include "TVector3.h" // for TVector3
40 
41 #include <iostream> // for operator<<, basic_ostream, etc
42 #include <vector>
43 
44 
45 class FairGenericStack;
46 class FairMCEventHeader;
47 class TF1;
48 class TIterator;
49 
50 
51 
52 class FairFilteredPrimaryGenerator : public FairPrimaryGenerator
53 {
54 
55 public:
56 
59 
60 
62  FairFilteredPrimaryGenerator(const char* name, const char* title="Filtered Generator");
63 
64 
67 
69  virtual Bool_t Init();
70 
72  void AndFilter(FairEvtFilter* filter) {
73  AddFilter(filter, FairEvtFilter::kAnd, kFALSE);
74  }
75 
77  void AndNotFilter(FairEvtFilter* filter) {
78  AddFilter(filter, FairEvtFilter::kAnd, kTRUE);
79  }
80 
82  void OrFilter(FairEvtFilter* filter) {
83  AddFilter(filter, FairEvtFilter::kOr, kFALSE);
84  }
85 
87  void OrNotFilter(FairEvtFilter* filter) {
88  AddFilter(filter, FairEvtFilter::kOr, kTRUE);
89  }
90 
92  void AddVetoFilter(FairEvtFilter* filter) {
93  if ( ! fVetoFilterList ) {
94  std::cout << "Empty fFilterList pointer ! \n";
95  return;
96  }
97  fVetoFilterList->Add(filter);
98  fEventVetoFilterActive = kTRUE;
99  }
100 
101 
112  virtual Bool_t GenerateEvent(FairGenericStack* pStack);
113 
114 
115  TObjArray* GetListOfFilters() { return fFilterList;}
116 
117  TObjArray* GetListOfVetoFilters() { return fVetoFilterList;}
118 
120  void SetFilterMaxTries(Int_t maxTries=99999)
121  {
122  if ( maxTries > 0 ){
123  fEvtFilterStat.fFilterMaxTries = maxTries;
124  std::cout << "FairFilteredPrimaryGenerator: maxTries is now set to " << fEvtFilterStat.fFilterMaxTries << "\n";
125  } else {
126  std::cout << "\n\n\n -WARNING from FairFilteredPrimaryGenerator: maxTries must be a positive number! Check your SetFilterMaxTries call!\n\n\n";
127  }
128  }
129 
133  }
134 
138  }
139 
141  void SetEventPrintFrequency(int freq)
142  {
143  fEventPrintFreq = freq;
144  }
145 
156  }
157 
159  void WriteEvtFilterStatsToRootFile(TFile* outputFile=NULL){
160  std::cout << "\n\nGenerated Events = " << GetNumberOfGeneratedEvents() << "\n";
161  if (0 < GetNumberOfFilterFailedEvents() ) {
162  std::cout << "WARNING: Number of events where the event filter FAILED " << GetNumberOfFilterFailedEvents() << "\n\n\n";
163  std::cout << "Random events were accepted to avoid infinite loops. \n";
164  std::cout << "Try increasing the max. number of tries or change your filter (maybe the generators do not produce such events as you want).\n\n";
165  }
166  if(outputFile==NULL) outputFile = FairRunSim::Instance()->GetOutputFile();
167  outputFile->cd();
168  outputFile->mkdir("FairEvtFilter");
169  outputFile->cd("FairEvtFilter");
170  fEvtFilterStat.Write();
171  outputFile->cd();
172  }
173 
177  void SetVerbose(Int_t verbose=12){
178  if ( verbose >= 0 ){
179  fVerbose = verbose;
180  std::cout << "FairFilteredPrimaryGenerator: fVerbose is now set to " << verbose << "\n";
181  } else {
182  std::cout << "\n\n\n -WARNING from FairFilteredPrimaryGenerator: verbose must be a positive number! Check your SetVerbose call!\n\n\n";
183  }
184  }
185 
186 
187 
188 protected:
190  TObjArray* fVetoFilterList;
192  TIterator* fVetoFilterIter;
193 
194  TObjArray* fFilterList;
196  TIterator* fFilterIter;
197 
201  Int_t fVerbose;
209  std::vector<Bool_t> filterAcceptEvent;
216  std::vector<UInt_t> fLogicalFilterOperation;
222  std::vector<Bool_t> fFilterNegation;
223 
224 
227 
229  Int_t fEventPrintFreq; // KG, added 03/2017
230 
231 
233  void AddFilter(FairEvtFilter* filter, FairEvtFilter::LogicOp op, Bool_t negateFilter) {
234  if ( ! fFilterList ) {
235  std::cout << "Empty fFilterList pointer ! \n";
236  return;
237  }
238  fFilterList->Add(filter);
239  fEventFilterActive = kTRUE;
240  // add the settings for every new filter
241  if(fFilterList->GetEntriesFast()!=1){
242  fLogicalFilterOperation.push_back(op);
243  }
244  fFilterNegation.push_back(negateFilter);
245  }
246 
247 
248 private:
251 
252 
254 };
255 
256 #endif
Int_t GetNumberOfGeneratedEvents()
returns the total (accepted + rejected) number of events generated by the event generators. If no event filters are used this number is equal to the number of simulated events.
virtual Bool_t Init()
Initialize the event generator(s) and the event (veto) filter(s).
Bool_t fEventFilterActive
returns kTRUE if any non-veto event filter is registerd.
FairFilteredPrimaryGenerator & operator=(const FairFilteredPrimaryGenerator &)
std::vector< UInt_t > fLogicalFilterOperation
vector containing the logical operations with which the outputs of the non-veto event filters should ...
TObjArray * fVetoFilterList
List of registered veto filters.
void WriteEvtFilterStatsToRootFile(TFile *outputFile=NULL)
Writes all relevant event filter information to the output root file.
void AndNotFilter(FairEvtFilter *filter)
Register a non-veto event filter using a logical AND NOT to connect with previously defined non-veto ...
void AndFilter(FairEvtFilter *filter)
Register a non-veto event filter using a logical AND to connect with previously defined non-veto even...
#define verbose
Int_t fEventNrFiltered
Event number (Set by the filtered primary generator.
Int_t GetNumberOfFilterFailedEvents()
Returns the number of cases in which no matching event was found within the set max. tries.
void SetFilterMaxTries(Int_t maxTries=99999)
Define the maximum number of times that this object should try to find an event which suits all event...
TObjArray * fFilterList
List of registered filters.
void SetEventPrintFrequency(int freq)
Sets the frequency (accepted events) for printout (verbose>0) of accepted and generated events...
ClassDef(FairFilteredPrimaryGenerator, 2)
Primary generator with added event filtering capabilities.
void OrNotFilter(FairEvtFilter *filter)
Register a non-veto event filter using a logical OR NOT to connect with previously defined non-veto e...
std::vector< Bool_t > fFilterNegation
vector determining whether the output of a non-veto event filter should be negated or not...
FairFilteredPrimaryGenerator()
Default constructor.
FairEvtFilterParams fEvtFilterStat
Contains the statistics of the event filtering process.
TIterator * fVetoFilterIter
Iterator over veto filter list.
Int_t fVerbose
Level of commenting output, 0 means no output, higher gives more output.
Int_t GetNumberOfFilterMaxTries()
returns the maximum number of times that this object should try to find an event which suits all even...
TString name
Int_t fEventPrintFreq
Print frequency for filtered events.
void OrFilter(FairEvtFilter *filter)
Register a non-veto event filter using a logical OR to connect with previously defined non-veto event...
TIterator * fFilterIter
Iterator over filter list.
virtual Bool_t GenerateEvent(FairGenericStack *pStack)
Calls event generators and the event filters.
void AddVetoFilter(FairEvtFilter *filter)
Register an event veto filter. Veto filters have higher priority than regular event filters...
Bool_t fEventVetoFilterActive
returns kTRUE if any event veto filter is registered.
std::vector< Bool_t > filterAcceptEvent
Vector containing the results of the EventMatches methods for every registered non-veto event filter ...
void AddFilter(FairEvtFilter *filter, FairEvtFilter::LogicOp op, Bool_t negateFilter)
Registers a regular (non-veto) filter. This method is not supposed to be directly used by the user...
void SetVerbose(Int_t verbose=12)
Set the level of commenting output.
virtual ~FairFilteredPrimaryGenerator()
Destructor.