PandaRoot
FairFilteredPrimaryGenerator.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 
33 #ifndef FairFilteredPrimaryGenerator_H
34 #define FairFilteredPrimaryGenerator_H
35 
36 #include "FairPrimaryGenerator.h"
37 
38 #include "FairEvtFilter.h"
39 #include "FairEvtFilterParams.h"
40 
41 #include "TNamed.h" // for TNamed
42 #include "TFile.h" // for TFile
43 #include "FairRunSim.h"
44 
45 #include "FairGenerator.h" // for FairGenerator
46 
47 #include <iosfwd> // for ostream
48 #include "Rtypes.h" // for Double_t, Bool_t, Int_t, etc
49 #include "TObjArray.h" // for TObjArray
50 #include "TVector3.h" // for TVector3
51 
52 #include <iostream> // for operator<<, basic_ostream, etc
53 #include <vector>
54 
55 class FairGenericStack;
56 class FairMCEventHeader;
57 class TF1;
58 class TIterator;
59 
60 class FairFilteredPrimaryGenerator : public FairPrimaryGenerator {
61 
62  public:
65 
67  FairFilteredPrimaryGenerator(const char *name, const char *title = "Filtered Generator");
68 
71 
73  virtual Bool_t Init();
74 
76  void AndFilter(FairEvtFilter *filter) { AddFilter(filter, FairEvtFilter::kAnd, kFALSE); }
77 
79  void AndNotFilter(FairEvtFilter *filter) { AddFilter(filter, FairEvtFilter::kAnd, kTRUE); }
80 
82  void OrFilter(FairEvtFilter *filter) { AddFilter(filter, FairEvtFilter::kOr, kFALSE); }
83 
85  void OrNotFilter(FairEvtFilter *filter) { AddFilter(filter, FairEvtFilter::kOr, kTRUE); }
86 
89  {
90  if (!fVetoFilterList) {
91  std::cout << "Empty fFilterList pointer ! \n";
92  return;
93  }
94  fVetoFilterList->Add(filter);
95  fEventVetoFilterActive = kTRUE;
96  }
97 
108  virtual Bool_t GenerateEvent(FairGenericStack *pStack);
109 
110  TObjArray *GetListOfFilters() { return fFilterList; }
111 
112  TObjArray *GetListOfVetoFilters() { return fVetoFilterList; }
113 
115  void SetFilterMaxTries(Int_t maxTries = 99999)
116  {
117  if (maxTries > 0) {
118  fEvtFilterStat.fFilterMaxTries = maxTries;
119  std::cout << "FairFilteredPrimaryGenerator: maxTries is now set to " << fEvtFilterStat.fFilterMaxTries << "\n";
120  } else {
121  std::cout << "\n\n\n -WARNING from FairFilteredPrimaryGenerator: maxTries must be a positive number! Check your SetFilterMaxTries call!\n\n\n";
122  }
123  }
124 
127 
131 
133  void SetEventPrintFrequency(int freq) { fEventPrintFreq = freq; }
134 
144 
146  void WriteEvtFilterStatsToRootFile(TFile *outputFile = nullptr)
147  {
148  std::cout << "\n\nGenerated Events = " << GetNumberOfGeneratedEvents() << "\n";
149  if (0 < GetNumberOfFilterFailedEvents()) {
150  std::cout << "WARNING: Number of events where the event filter FAILED " << GetNumberOfFilterFailedEvents() << "\n\n\n";
151  std::cout << "Random events were accepted to avoid infinite loops. \n";
152  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";
153  }
154  if (outputFile == nullptr)
155  outputFile = FairRunSim::Instance()->GetOutputFile();
156  outputFile->cd();
157  outputFile->mkdir("FairEvtFilter");
158  outputFile->cd("FairEvtFilter");
159  fEvtFilterStat.Write();
160  outputFile->cd();
161  }
162 
166  void SetVerbose(Int_t verbose = 12)
167  {
168  if (verbose >= 0) {
169  fVerbose = verbose;
170  std::cout << "FairFilteredPrimaryGenerator: fVerbose is now set to " << verbose << "\n";
171  } else {
172  std::cout << "\n\n\n -WARNING from FairFilteredPrimaryGenerator: verbose must be a positive number! Check your SetVerbose call!\n\n\n";
173  }
174  }
175 
176  protected:
178  TObjArray *fVetoFilterList;
180  TIterator *fVetoFilterIter;
181 
182  TObjArray *fFilterList;
184  TIterator *fFilterIter;
185 
189  Int_t fVerbose;
197  std::vector<Bool_t> filterAcceptEvent;
204  std::vector<UInt_t> fLogicalFilterOperation;
210  std::vector<Bool_t> fFilterNegation;
211 
214 
216  Int_t fEventPrintFreq; // KG, added 03/2017
217 
219  void AddFilter(FairEvtFilter *filter, FairEvtFilter::LogicOp op, Bool_t negateFilter)
220  {
221  if (!fFilterList) {
222  std::cout << "Empty fFilterList pointer ! \n";
223  return;
224  }
225  fFilterList->Add(filter);
226  fEventFilterActive = kTRUE;
227  // add the settings for every new filter
228  if (fFilterList->GetEntriesFast() != 1) {
229  fLogicalFilterOperation.push_back(op);
230  }
231  fFilterNegation.push_back(negateFilter);
232  }
233 
234  private:
237 
238  ClassDef(FairFilteredPrimaryGenerator, 2);
239 };
240 
241 #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.