PandaRoot
FairFilteredPrimaryGenerator.h
Go to the documentation of this file.
1 
21 #ifndef FairFilteredPrimaryGenerator_H
22 #define FairFilteredPrimaryGenerator_H
23 
24 #include "FairPrimaryGenerator.h"
25 
26 #include "FairEvtFilter.h"
27 #include "FairEvtFilterParams.h"
28 
29 #include "TNamed.h" // for TNamed
30 #include "TFile.h" // for TFile
31 #include "FairRunSim.h"
32 
33 #include "FairGenerator.h" // for FairGenerator
34 
35 #include <iosfwd> // for ostream
36 #include "Rtypes.h" // for Double_t, Bool_t, Int_t, etc
37 #include "TObjArray.h" // for TObjArray
38 #include "TVector3.h" // for TVector3
39 
40 #include <iostream> // for operator<<, basic_ostream, etc
41 #include <vector>
42 
43 class FairGenericStack;
44 class FairMCEventHeader;
45 class TF1;
46 class TIterator;
47 
48 class FairFilteredPrimaryGenerator : public FairPrimaryGenerator {
49 
50  public:
53 
55  FairFilteredPrimaryGenerator(const char *name, const char *title = "Filtered Generator");
56 
59 
61  virtual Bool_t Init();
62 
64  void AndFilter(FairEvtFilter *filter) { AddFilter(filter, FairEvtFilter::kAnd, kFALSE); }
65 
67  void AndNotFilter(FairEvtFilter *filter) { AddFilter(filter, FairEvtFilter::kAnd, kTRUE); }
68 
70  void OrFilter(FairEvtFilter *filter) { AddFilter(filter, FairEvtFilter::kOr, kFALSE); }
71 
73  void OrNotFilter(FairEvtFilter *filter) { AddFilter(filter, FairEvtFilter::kOr, kTRUE); }
74 
77  {
78  if (!fVetoFilterList) {
79  std::cout << "Empty fFilterList pointer ! \n";
80  return;
81  }
82  fVetoFilterList->Add(filter);
83  fEventVetoFilterActive = kTRUE;
84  }
85 
96  virtual Bool_t GenerateEvent(FairGenericStack *pStack);
97 
98  TObjArray *GetListOfFilters() { return fFilterList; }
99 
100  TObjArray *GetListOfVetoFilters() { return fVetoFilterList; }
101 
103  void SetFilterMaxTries(Int_t maxTries = 99999)
104  {
105  if (maxTries > 0) {
106  fEvtFilterStat.fFilterMaxTries = maxTries;
107  std::cout << "FairFilteredPrimaryGenerator: maxTries is now set to " << fEvtFilterStat.fFilterMaxTries << "\n";
108  } else {
109  std::cout << "\n\n\n -WARNING from FairFilteredPrimaryGenerator: maxTries must be a positive number! Check your SetFilterMaxTries call!\n\n\n";
110  }
111  }
112 
115 
119 
121  void SetEventPrintFrequency(int freq) { fEventPrintFreq = freq; }
122 
132 
134  void WriteEvtFilterStatsToRootFile(TFile *outputFile = nullptr)
135  {
136  std::cout << "\n\nGenerated Events = " << GetNumberOfGeneratedEvents() << "\n";
137  if (0 < GetNumberOfFilterFailedEvents()) {
138  std::cout << "WARNING: Number of events where the event filter FAILED " << GetNumberOfFilterFailedEvents() << "\n\n\n";
139  std::cout << "Random events were accepted to avoid infinite loops. \n";
140  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";
141  }
142  if (outputFile == nullptr)
143  outputFile = FairRunSim::Instance()->GetOutputFile();
144  outputFile->cd();
145  outputFile->mkdir("FairEvtFilter");
146  outputFile->cd("FairEvtFilter");
147  fEvtFilterStat.Write();
148  outputFile->cd();
149  }
150 
154  void SetVerbose(Int_t verbose = 12)
155  {
156  if (verbose >= 0) {
157  fVerbose = verbose;
158  std::cout << "FairFilteredPrimaryGenerator: fVerbose is now set to " << verbose << "\n";
159  } else {
160  std::cout << "\n\n\n -WARNING from FairFilteredPrimaryGenerator: verbose must be a positive number! Check your SetVerbose call!\n\n\n";
161  }
162  }
163 
164  protected:
166  TObjArray *fVetoFilterList;
168  TIterator *fVetoFilterIter;
169 
170  TObjArray *fFilterList;
172  TIterator *fFilterIter;
173 
177  Int_t fVerbose;
185  std::vector<Bool_t> filterAcceptEvent;
192  std::vector<UInt_t> fLogicalFilterOperation;
198  std::vector<Bool_t> fFilterNegation;
199 
202 
204  Int_t fEventPrintFreq; // KG, added 03/2017
205 
207  void AddFilter(FairEvtFilter *filter, FairEvtFilter::LogicOp op, Bool_t negateFilter)
208  {
209  if (!fFilterList) {
210  std::cout << "Empty fFilterList pointer ! \n";
211  return;
212  }
213  fFilterList->Add(filter);
214  fEventFilterActive = kTRUE;
215  // add the settings for every new filter
216  if (fFilterList->GetEntriesFast() != 1) {
217  fLogicalFilterOperation.push_back(op);
218  }
219  fFilterNegation.push_back(negateFilter);
220  }
221 
222  private:
225 
226  ClassDef(FairFilteredPrimaryGenerator, 2);
227 };
228 
229 #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.
Bool_t fEventFilterActive
returns kTRUE if any non-veto event filter is registerd.
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 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...
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...
virtual ~FairFilteredPrimaryGenerator()
Destructor.
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.
virtual Bool_t GenerateEvent(FairGenericStack *pStack)
Calls event generators and the event filters.
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...
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...
virtual Bool_t Init()
Initialize the event generator(s) and the event (veto) filter(s).
TIterator * fFilterIter
Iterator over filter list.
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 WriteEvtFilterStatsToRootFile(TFile *outputFile=nullptr)
Writes all relevant event filter information to the output root file.
void SetVerbose(Int_t verbose=12)
Set the level of commenting output.