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

Primary generator with added event filtering capabilities. More...

#include <FairFilteredPrimaryGenerator.h>

Inheritance diagram for FairFilteredPrimaryGenerator:

Public Member Functions

 FairFilteredPrimaryGenerator ()
 Default constructor. More...
 
 FairFilteredPrimaryGenerator (const char *name, const char *title="Filtered Generator")
 Constructor with name and title. More...
 
virtual ~FairFilteredPrimaryGenerator ()
 Destructor. More...
 
virtual Bool_t Init ()
 Initialize the event generator(s) and the event (veto) filter(s). More...
 
void AndFilter (FairEvtFilter *filter)
 Register a non-veto event filter using a logical AND to connect with previously defined non-veto event filters. More...
 
void AndNotFilter (FairEvtFilter *filter)
 Register a non-veto event filter using a logical AND NOT to connect with previously defined non-veto event filters. More...
 
void OrFilter (FairEvtFilter *filter)
 Register a non-veto event filter using a logical OR to connect with previously defined non-veto event filters. More...
 
void OrNotFilter (FairEvtFilter *filter)
 Register a non-veto event filter using a logical OR NOT to connect with previously defined non-veto event filters. More...
 
void AddVetoFilter (FairEvtFilter *filter)
 Register an event veto filter. Veto filters have higher priority than regular event filters. If the event matches any veto filter, it will be skipped. More...
 
virtual Bool_t GenerateEvent (FairGenericStack *pStack)
 Calls event generators and the event filters. More...
 
TObjArray * GetListOfFilters ()
 
TObjArray * GetListOfVetoFilters ()
 
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 filters. More...
 
Int_t GetNumberOfFilterMaxTries ()
 returns the maximum number of times that this object should try to find an event which suits all event filters. More...
 
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. More...
 
void SetEventPrintFrequency (int freq)
 Sets the frequency (accepted events) for printout (verbose>0) of accepted and generated events. More...
 
Int_t GetNumberOfFilterFailedEvents ()
 Returns the number of cases in which no matching event was found within the set max. tries. More...
 
void WriteEvtFilterStatsToRootFile (TFile *outputFile=NULL)
 Writes all relevant event filter information to the output root file. More...
 
void SetVerbose (Int_t verbose=12)
 Set the level of commenting output. More...
 

Protected Member Functions

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. See public methods for user interfaces to this method. More...
 

Protected Attributes

TObjArray * fVetoFilterList
 List of registered veto filters. More...
 
TIterator * fVetoFilterIter
 Iterator over veto filter list. More...
 
TObjArray * fFilterList
 List of registered filters. More...
 
TIterator * fFilterIter
 Iterator over filter list. More...
 
FairEvtFilterParams fEvtFilterStat
 Contains the statistics of the event filtering process. More...
 
Int_t fVerbose
 Level of commenting output, 0 means no output, higher gives more output. More...
 
Bool_t fEventVetoFilterActive
 returns kTRUE if any event veto filter is registered. More...
 
Bool_t fEventFilterActive
 returns kTRUE if any non-veto event filter is registerd. More...
 
std::vector< Bool_tfilterAcceptEvent
 Vector containing the results of the EventMatches methods for every registered non-veto event filter in the corresponding order. More...
 
std::vector< UInt_t > fLogicalFilterOperation
 vector containing the logical operations with which the outputs of the non-veto event filters should be combined. More...
 
std::vector< Bool_tfFilterNegation
 vector determining whether the output of a non-veto event filter should be negated or not. More...
 
Int_t fEventNrFiltered
 Event number (Set by the filtered primary generator. More...
 
Int_t fEventPrintFreq
 Print frequency for filtered events. More...
 

Private Member Functions

 FairFilteredPrimaryGenerator (const FairFilteredPrimaryGenerator &)
 
FairFilteredPrimaryGeneratoroperator= (const FairFilteredPrimaryGenerator &)
 
 ClassDef (FairFilteredPrimaryGenerator, 2)
 

Detailed Description

Primary generator with added event filtering capabilities.

Author
Martin J. Galuska <martin [dot] j [dot] galuska (at) physik [dot] uni (minus) giessen [dot] de>
Katja Kleeberg

This class adds event filtering capabilities to FairPrimaryGenerator which is used internally for handling the event generators and so on. The event filtering is performed after the event generation and before the particle transport through the detector model.

From the description of FairPrimaryGenerator:

The FairFilteredPrimaryGenerator is responsible for the handling of the MC input. Several input generators can be registered to it; these have to be derived from the FairGenerator class. The FairFilteredPrimaryGenerator defines position and (optionally) smearing of the primary vertex. This class should be instantiated only once.

Definition at line 52 of file FairFilteredPrimaryGenerator.h.

Constructor & Destructor Documentation

FairFilteredPrimaryGenerator::FairFilteredPrimaryGenerator ( )

Default constructor.

Definition at line 28 of file FairFilteredPrimaryGenerator.cxx.

29 : FairPrimaryGenerator(),
30 
31  fVetoFilterList(new TObjArray()),
32  fVetoFilterIter(fVetoFilterList->MakeIterator()),
33  fFilterList(new TObjArray()),
34  fFilterIter(fFilterList->MakeIterator()),
36  fVerbose(3),
37  fEventVetoFilterActive(kFALSE),
38  fEventFilterActive(kFALSE),
41  fEventPrintFreq(100)
42 {
43 }
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.
Int_t fEventNrFiltered
Event number (Set by the filtered primary generator.
TObjArray * fFilterList
List of registered filters.
std::vector< Bool_t > fFilterNegation
vector determining whether the output of a non-veto event filter should be negated or not...
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 fEventPrintFreq
Print frequency for filtered events.
TIterator * fFilterIter
Iterator over filter list.
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 ...
FairFilteredPrimaryGenerator::FairFilteredPrimaryGenerator ( const char *  name,
const char *  title = "Filtered Generator" 
)

Constructor with name and title.

Definition at line 49 of file FairFilteredPrimaryGenerator.cxx.

50 : FairPrimaryGenerator(name,title),
51 
52  fVetoFilterList(new TObjArray()),
53  fVetoFilterIter(fVetoFilterList->MakeIterator()),
54  fFilterList(new TObjArray()),
55  fFilterIter(fFilterList->MakeIterator()),
57  fVerbose(3),
58  fEventVetoFilterActive(kFALSE),
59  fEventFilterActive(kFALSE),
62  fEventPrintFreq(100)
63 {
64 }
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.
Int_t fEventNrFiltered
Event number (Set by the filtered primary generator.
TObjArray * fFilterList
List of registered filters.
std::vector< Bool_t > fFilterNegation
vector determining whether the output of a non-veto event filter should be negated or not...
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.
TString name
Int_t fEventPrintFreq
Print frequency for filtered events.
TIterator * fFilterIter
Iterator over filter list.
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 ...
FairFilteredPrimaryGenerator::~FairFilteredPrimaryGenerator ( )
virtual

Destructor.

Definition at line 85 of file FairFilteredPrimaryGenerator.cxx.

References fFilterIter, fFilterList, fVetoFilterIter, and fVetoFilterList.

86 {
87  fVetoFilterList->Delete();
88  delete fVetoFilterList;
89  delete fVetoFilterIter;
90 
91  fFilterList->Delete();
92  delete fFilterList;
93  delete fFilterIter;
94 
95  //cout<<" Leave Destructor of FairFilteredPrimaryGenerator"<<endl;
96 }
TObjArray * fVetoFilterList
List of registered veto filters.
TObjArray * fFilterList
List of registered filters.
TIterator * fVetoFilterIter
Iterator over veto filter list.
TIterator * fFilterIter
Iterator over filter list.
FairFilteredPrimaryGenerator::FairFilteredPrimaryGenerator ( const FairFilteredPrimaryGenerator )
private

Member Function Documentation

void FairFilteredPrimaryGenerator::AddFilter ( FairEvtFilter filter,
FairEvtFilter::LogicOp  op,
Bool_t  negateFilter 
)
inlineprotected

Registers a regular (non-veto) filter. This method is not supposed to be directly used by the user. See public methods for user interfaces to this method.

Definition at line 233 of file FairFilteredPrimaryGenerator.h.

References fEventFilterActive, fFilterList, fFilterNegation, and fLogicalFilterOperation.

Referenced by AndFilter(), AndNotFilter(), OrFilter(), and OrNotFilter().

233  {
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  }
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 * fFilterList
List of registered filters.
std::vector< Bool_t > fFilterNegation
vector determining whether the output of a non-veto event filter should be negated or not...
void FairFilteredPrimaryGenerator::AddVetoFilter ( FairEvtFilter filter)
inline

Register an event veto filter. Veto filters have higher priority than regular event filters. If the event matches any veto filter, it will be skipped.

Definition at line 92 of file FairFilteredPrimaryGenerator.h.

References fEventVetoFilterActive, and fVetoFilterList.

Referenced by sim_filter_ex1().

92  {
93  if ( ! fVetoFilterList ) {
94  std::cout << "Empty fFilterList pointer ! \n";
95  return;
96  }
97  fVetoFilterList->Add(filter);
98  fEventVetoFilterActive = kTRUE;
99  }
TObjArray * fVetoFilterList
List of registered veto filters.
Bool_t fEventVetoFilterActive
returns kTRUE if any event veto filter is registered.
void FairFilteredPrimaryGenerator::AndFilter ( FairEvtFilter filter)
inline

Register a non-veto event filter using a logical AND to connect with previously defined non-veto event filters.

Definition at line 72 of file FairFilteredPrimaryGenerator.h.

References AddFilter(), and FairEvtFilter::kAnd.

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

72  {
73  AddFilter(filter, FairEvtFilter::kAnd, kFALSE);
74  }
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 FairFilteredPrimaryGenerator::AndNotFilter ( FairEvtFilter filter)
inline

Register a non-veto event filter using a logical AND NOT to connect with previously defined non-veto event filters.

Definition at line 77 of file FairFilteredPrimaryGenerator.h.

References AddFilter(), and FairEvtFilter::kAnd.

Referenced by sim_filter_ex1().

77  {
78  AddFilter(filter, FairEvtFilter::kAnd, kTRUE);
79  }
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...
FairFilteredPrimaryGenerator::ClassDef ( FairFilteredPrimaryGenerator  ,
 
)
private
Bool_t FairFilteredPrimaryGenerator::GenerateEvent ( FairGenericStack *  pStack)
virtual

Calls event generators and the event filters.

To be called at the beginning of each event from FairMCApplication. Generates an event vertex and calls the ReadEvent methods from the registered generators. Calls defined event (veto) filters to decide whether to process the event or to call the event generators again.

Parameters
pStackThe particle stack
Returns
kTRUE if successful, kFALSE if not

Definition at line 102 of file FairFilteredPrimaryGenerator.cxx.

References Bool_t, FairEvtFilter::EventMatches(), fEventFilterActive, fEventNrFiltered, fEventPrintFreq, fEventVetoFilterActive, fEvtFilterStat, FairEvtFilterParams::fFailedFilterEvents, fFilterIter, fFilterList, FairEvtFilterParams::fFilterMaxTries, fFilterNegation, FairEvtFilterParams::fGeneratedEvents, filterAcceptEvent, FairEvtFilter::FilterActive(), fLogicalFilterOperation, fVerbose, fVetoFilterIter, RhoFactory::Instance(), FairEvtFilter::kAnd, FairEvtFilter::kOr, RhoFactory::Reset(), and FairEvtFilter::SetListOfParticles().

103 {
105 
106  Int_t iTry=0; // number of attempts to find the next event that suits your filter
107  Bool_t acceptEvent = kFALSE; // is kTRUE if the event is finally accepted
108 
109  if(kTRUE == fEventFilterActive){
110  // sanity check settings for logical combinations of event filters
111  if((int)fFilterNegation.size()!=fFilterList->GetEntriesFast()){
112  std::cout << "\n\n\n -WARNING from FairFilteredPrimaryGenerator: size of the negation vector has to be equal to the number of added filters. Check your FilterNegation call!\n\n\n";
113  return kFALSE;
114  }
115  if((int)fLogicalFilterOperation.size()!=fFilterList->GetEntriesFast()-1){
116  std::cout << "\n\n\n -WARNING from FairFilteredPrimaryGenerator: size of the operation vector has to be equal to the number of added filters minus one. Check your LogicalFilterOperation call!\n\n\n";
117  return kFALSE;
118  }
119  }
120 
121  do{
122  ++iTry; // count how often we run the generators before we accept an event
123  ++fEvtFilterStat.fGeneratedEvents; // total number of generated events
124  Bool_t vetoed = kFALSE; // is kTRUE if event matches any of the veto filters (we will skip such events)
125 
126 
127 
128  pStack->Reset(); // Clean the stack
129  FairPrimaryGenerator::GenerateEvent(pStack); // fill the stack
130  fEventNr = fEventNrFiltered; // Fix event numbering in FairPrimaryGenerator (otherwise fRun will stop too early...)
131 
132 
133  if(kTRUE == fEventVetoFilterActive){// skip veto filter checking if no veto filters are set
134  // check veto filters. skip event if the event satisfies any of the veto filters' conditions
135  // Call the FilterAccept methods for all veto filters
136  fVetoFilterIter->Reset();
137  TObject* vetoFilterObject = 0;
138  FairEvtFilter* vetoFilter = 0;
139  while( (vetoFilterObject = fVetoFilterIter->Next()) ) {
140  vetoFilter = dynamic_cast<FairEvtFilter*> (vetoFilterObject);
141  if ( ! vetoFilter ) { // this should never happen
142  cout << " \n\n\n -FATAL ERROR from FairFilteredPrimaryGenerator: Veto filter cast did not work!\n\n\n" << endl;
143  return kFALSE;
144  }
145  if( !( vetoFilter->FilterActive() ) ){
146  cout << " \n\n\n -WARNING from FairFilteredPrimaryGenerator: There are no filter settings for filter " << vetoFilter->GetName() << " : " << vetoFilter->GetTitle() << ". Execution is halted!\n\n\n" << endl;
147  return kFALSE;
148  }
149  Bool_t setParticlesOk = vetoFilter->SetListOfParticles(pStack->GetListOfParticles());//committing the list of particles to be checked
150  if( ! setParticlesOk ) {
151  cout << " \n\n\n -WARNING from FairFilteredPrimaryGenerator: Particles were not pushed correctly to veto filter.\n\n\n" << endl;
152  return kFALSE;
153  }
154  if ( vetoFilter->EventMatches(fEvtFilterStat.fGeneratedEvents) ) {
155  // if event matches veto filter, skip event
156  if( 3 < fVerbose ){
157  cout << "\n Event is NOT accepted because it matches a veto filter \n";
158  }
159  vetoed = kTRUE;
160  break; // no need to check other veto filters as event will be skipped anyway
161  }
162  } // end while
163  } // if veto filtering
164 
165  if ( vetoed ) { continue; }; // skip event if it matches at least one veto filter
166 
167  if ( kFALSE == fEventFilterActive ){ // skip event filtering (i.e. accept all events) if filtering is not requested
168  acceptEvent = kTRUE;
169  if( 3 < fVerbose ){
170  cout << "\n Event is accepted because event filtering is not requested.\n\n";
171  }
172  } else { // do the event filtering (i.e. decide whether to accept or not accept the generated events) if filtering is requested
173  Int_t iCheckFilter=0;
174  // Call the FilterAccept methods for all registered (non-veto) filters
175  fFilterIter->Reset();
176  TObject* eventFilterObject = 0;
177  FairEvtFilter* eventFilter = 0;
178  filterAcceptEvent.clear();
179  while( (eventFilterObject = fFilterIter->Next()) ) {
180  eventFilter = dynamic_cast<FairEvtFilter*> (eventFilterObject);
181  if ( ! eventFilter ) {
182  return kFALSE; } // this should never happen
183  if( !( eventFilter->FilterActive() ) ){
184  cout << " \n\n\n -WARNING from FairFilteredPrimaryGenerator: There are no filter settings for filter " << eventFilter->GetName() << " : " << eventFilter->GetTitle() << ". Execution is halted!\n\n\n" << endl;
185  return kFALSE;
186  }
187  Bool_t setParticlesOk = eventFilter->SetListOfParticles(pStack->GetListOfParticles());//committing the list of particles to be checked
188  if( ! setParticlesOk ) {
189  return kFALSE;
190  }
192  // check whether the simulated event matches the individual event filter
193  iCheckFilter++;
194  }
195 
196  //----------FilterNegation
197  // negate output from the event filters which should be negated
198  for(UInt_t iFil=0; iFil<fFilterNegation.size(); iFil++){
199  if(kTRUE == fFilterNegation[iFil]){
200  filterAcceptEvent[iFil]=(!filterAcceptEvent[iFil]);
201  }
202  }
203 
204 
205  //----------LogicalFilterOperation
206  // combine (possibly negated) outputs from individual event filters using the logical connections stored in fLogicalFilterOperation
207  // in order to check whether the event is accepted
208  acceptEvent=filterAcceptEvent[0];
209  for(UInt_t iOper=0; iOper<fLogicalFilterOperation.size();iOper++){
211  acceptEvent=acceptEvent && filterAcceptEvent[iOper+1];
212  }else{
214  acceptEvent=acceptEvent || filterAcceptEvent[iOper+1];
215  }else{
216  std::cout << "\n\n\n -WARNING from FairFilteredPrimaryGenerator: the entries of the operation vector have to be kAnd or kOr. Check your LogicalFilterOperation call!\n\n\n";
217  return kFALSE;
218  }
219  }
220  }
221  if( !acceptEvent && 3 < fVerbose ){
222  cout << "\n Event is NOT accepted after combining the event filters.\n";
223  }
224  } // if (non-veto) event filtering
225  }while(! acceptEvent && iTry < fEvtFilterStat.fFilterMaxTries);
226 
227 
228  // Set the event number ALWAYS when filtering
230  fEvent->SetEventID(fEventNrFiltered);
231 
232  if ( !acceptEvent && ( fEventFilterActive||fEventVetoFilterActive ) ){
234  cout << "\n -E FairFilteredPrimaryGenerator: No event was found within " << iTry << " tries which satisfies your event filter.\n ";
235  cout << "I accept a random event as evtNr " << fEventNrFiltered << " to avoid infinite loops. \n";
236  cout << "Try increasing the max. number of tries or change your filter\n\n";
237  if (fVerbose > 3 ){
238  if(fFilterList->GetEntriesFast() > 0){cout << "\n Event is accepted after combining the appointed filters.\n";}
239  cout << iTry << " events simulated until I found a good one.\n";
240  cout << fEvtFilterStat.fGeneratedEvents << " events generated for finding " << fEventNr << " accepted events.\n";
241  cout << fEvtFilterStat.fFailedFilterEvents << " unsuccessful attempts in total to find an event that suits your filters\n\n";
242  }
243  }
244 
245  if (0 < fVerbose && (fEventNrFiltered%fEventPrintFreq)==0) cout << fEventNrFiltered << " of " << fEvtFilterStat.fGeneratedEvents << " generated events accepted.\n";
246 
247  return kTRUE;
248 }
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 ...
static void Reset()
Definition: RhoFactory.cxx:28
Int_t fEventNrFiltered
Event number (Set by the filtered primary generator.
TObjArray * fFilterList
List of registered filters.
std::vector< Bool_t > fFilterNegation
vector determining whether the output of a non-veto event filter should be negated or not...
virtual Bool_t EventMatches(Int_t evtNr)=0
FairEvtFilterParams fEvtFilterStat
Contains the statistics of the event filtering process.
TIterator * fVetoFilterIter
Iterator over veto filter list.
static RhoFactory * Instance()
Definition: RhoFactory.cxx:34
Int_t fVerbose
Level of commenting output, 0 means no output, higher gives more output.
Int_t fEventPrintFreq
Print frequency for filtered events.
Bool_t SetListOfParticles(TClonesArray *ParticleList)
Definition: FairEvtFilter.h:57
TIterator * fFilterIter
Iterator over filter list.
Bool_t fEventVetoFilterActive
returns kTRUE if any event veto filter is registered.
virtual Bool_t FilterActive()=0
std::vector< Bool_t > filterAcceptEvent
Vector containing the results of the EventMatches methods for every registered non-veto event filter ...
TObjArray* FairFilteredPrimaryGenerator::GetListOfFilters ( )
inline

Definition at line 115 of file FairFilteredPrimaryGenerator.h.

References fFilterList.

115 { return fFilterList;}
TObjArray * fFilterList
List of registered filters.
TObjArray* FairFilteredPrimaryGenerator::GetListOfVetoFilters ( )
inline

Definition at line 117 of file FairFilteredPrimaryGenerator.h.

References fVetoFilterList.

117 { return fVetoFilterList;}
TObjArray * fVetoFilterList
List of registered veto filters.
Int_t FairFilteredPrimaryGenerator::GetNumberOfFilterFailedEvents ( )
inline

Returns the number of cases in which no matching event was found within the set max. tries.

This method returns 0 if everything works fine. If it returns a value >0 it means that you should set a higher limit in SetFilterMaxTries. If it returns a value which is equal to the number of events that you requested, it means that either the max. number of tries is set way too low or that the generator does not create such events that you are interested in or that your event filters cannot be satisfied at all (logical error).

Definition at line 154 of file FairFilteredPrimaryGenerator.h.

References fEvtFilterStat, and FairEvtFilterParams::fFailedFilterEvents.

Referenced by WriteEvtFilterStatsToRootFile().

154  {
156  }
FairEvtFilterParams fEvtFilterStat
Contains the statistics of the event filtering process.
Int_t FairFilteredPrimaryGenerator::GetNumberOfFilterMaxTries ( )
inline

returns the maximum number of times that this object should try to find an event which suits all event filters.

Definition at line 131 of file FairFilteredPrimaryGenerator.h.

References fEvtFilterStat, and FairEvtFilterParams::fFilterMaxTries.

131  {
133  }
FairEvtFilterParams fEvtFilterStat
Contains the statistics of the event filtering process.
Int_t FairFilteredPrimaryGenerator::GetNumberOfGeneratedEvents ( )
inline

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.

Definition at line 136 of file FairFilteredPrimaryGenerator.h.

References fEvtFilterStat, and FairEvtFilterParams::fGeneratedEvents.

Referenced by prod_fsim(), and WriteEvtFilterStatsToRootFile().

136  {
138  }
FairEvtFilterParams fEvtFilterStat
Contains the statistics of the event filtering process.
Bool_t FairFilteredPrimaryGenerator::Init ( )
virtual

Initialize the event generator(s) and the event (veto) filter(s).

Definition at line 67 of file FairFilteredPrimaryGenerator.cxx.

References fFilterList, fVetoFilterList, FairEvtFilter::Init(), and Init().

68 {
70  // Initialize list of veto filters
71  for(Int_t k=0; k<fVetoFilterList->GetEntries(); ++k ) {
72  FairEvtFilter* vetoFilter= (FairEvtFilter*) fVetoFilterList->At(k);
73  if(vetoFilter) { vetoFilter->Init(); }
74  }
75  // Initialize list of event filters
76  for(Int_t k=0; k<fFilterList->GetEntries(); ++k ) {
77  FairEvtFilter* filter= (FairEvtFilter*) fFilterList->At(k);
78  if(filter) { filter->Init(); }
79  }
80  return kTRUE;
81 }
TObjArray * fVetoFilterList
List of registered veto filters.
TObjArray * fFilterList
List of registered filters.
Bool_t Init()
Definition: FairEvtFilter.h:64
fRun Init()
Definition: NHitsPerEvent.C:20
FairFilteredPrimaryGenerator& FairFilteredPrimaryGenerator::operator= ( const FairFilteredPrimaryGenerator )
private
void FairFilteredPrimaryGenerator::OrFilter ( FairEvtFilter filter)
inline

Register a non-veto event filter using a logical OR to connect with previously defined non-veto event filters.

Definition at line 82 of file FairFilteredPrimaryGenerator.h.

References AddFilter(), and FairEvtFilter::kOr.

82  {
83  AddFilter(filter, FairEvtFilter::kOr, kFALSE);
84  }
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 FairFilteredPrimaryGenerator::OrNotFilter ( FairEvtFilter filter)
inline

Register a non-veto event filter using a logical OR NOT to connect with previously defined non-veto event filters.

Definition at line 87 of file FairFilteredPrimaryGenerator.h.

References AddFilter(), and FairEvtFilter::kOr.

87  {
88  AddFilter(filter, FairEvtFilter::kOr, kTRUE);
89  }
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 FairFilteredPrimaryGenerator::SetEventPrintFrequency ( int  freq)
inline

Sets the frequency (accepted events) for printout (verbose>0) of accepted and generated events.

Definition at line 141 of file FairFilteredPrimaryGenerator.h.

References fEventPrintFreq.

142  {
143  fEventPrintFreq = freq;
144  }
Int_t fEventPrintFreq
Print frequency for filtered events.
void FairFilteredPrimaryGenerator::SetFilterMaxTries ( Int_t  maxTries = 99999)
inline

Define the maximum number of times that this object should try to find an event which suits all event filters.

Definition at line 120 of file FairFilteredPrimaryGenerator.h.

References fEvtFilterStat, and FairEvtFilterParams::fFilterMaxTries.

Referenced by prod_fsim(), quickfsimana(), sim_filter_ex1(), and sim_filter_ex2().

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  }
FairEvtFilterParams fEvtFilterStat
Contains the statistics of the event filtering process.
void FairFilteredPrimaryGenerator::SetVerbose ( Int_t  verbose = 12)
inline

Set the level of commenting output.

Parameters
verboseLevel of commenting output, 0 means no output, higher gives more output.

Definition at line 177 of file FairFilteredPrimaryGenerator.h.

References fVerbose, and verbose.

Referenced by prod_fsim(), quickfsimana(), sim_complete(), sim_complete_newSTT(), sim_day1(), sim_filter_ex1(), sim_filter_ex2(), sim_radlength_complete(), and tut_sim().

177  {
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  }
#define verbose
Int_t fVerbose
Level of commenting output, 0 means no output, higher gives more output.
void FairFilteredPrimaryGenerator::WriteEvtFilterStatsToRootFile ( TFile *  outputFile = NULL)
inline

Writes all relevant event filter information to the output root file.

Definition at line 159 of file FairFilteredPrimaryGenerator.h.

References fEvtFilterStat, GetNumberOfFilterFailedEvents(), and GetNumberOfGeneratedEvents().

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

159  {
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  }
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.
Int_t GetNumberOfFilterFailedEvents()
Returns the number of cases in which no matching event was found within the set max. tries.
FairEvtFilterParams fEvtFilterStat
Contains the statistics of the event filtering process.

Member Data Documentation

Bool_t FairFilteredPrimaryGenerator::fEventFilterActive
protected

returns kTRUE if any non-veto event filter is registerd.

Definition at line 205 of file FairFilteredPrimaryGenerator.h.

Referenced by AddFilter(), and GenerateEvent().

Int_t FairFilteredPrimaryGenerator::fEventNrFiltered
protected

Event number (Set by the filtered primary generator.

Definition at line 226 of file FairFilteredPrimaryGenerator.h.

Referenced by GenerateEvent().

Int_t FairFilteredPrimaryGenerator::fEventPrintFreq
protected

Print frequency for filtered events.

Definition at line 229 of file FairFilteredPrimaryGenerator.h.

Referenced by GenerateEvent(), and SetEventPrintFrequency().

Bool_t FairFilteredPrimaryGenerator::fEventVetoFilterActive
protected

returns kTRUE if any event veto filter is registered.

Definition at line 203 of file FairFilteredPrimaryGenerator.h.

Referenced by AddVetoFilter(), and GenerateEvent().

FairEvtFilterParams FairFilteredPrimaryGenerator::fEvtFilterStat
protected
TIterator* FairFilteredPrimaryGenerator::fFilterIter
protected

Iterator over filter list.

Definition at line 196 of file FairFilteredPrimaryGenerator.h.

Referenced by GenerateEvent(), and ~FairFilteredPrimaryGenerator().

TObjArray* FairFilteredPrimaryGenerator::fFilterList
protected

List of registered filters.

Definition at line 194 of file FairFilteredPrimaryGenerator.h.

Referenced by AddFilter(), GenerateEvent(), GetListOfFilters(), Init(), and ~FairFilteredPrimaryGenerator().

std::vector<Bool_t> FairFilteredPrimaryGenerator::fFilterNegation
protected

vector determining whether the output of a non-veto event filter should be negated or not.

The vector grows automatically with every added non-veto event filter. A kTRUE entry at position i in the vector negates the i. filter's output, kFALSE entries do not negate.

Definition at line 222 of file FairFilteredPrimaryGenerator.h.

Referenced by AddFilter(), and GenerateEvent().

std::vector<Bool_t> FairFilteredPrimaryGenerator::filterAcceptEvent
protected

Vector containing the results of the EventMatches methods for every registered non-veto event filter in the corresponding order.

The content of the vector is overwritten for each generated event.

Definition at line 209 of file FairFilteredPrimaryGenerator.h.

Referenced by GenerateEvent().

std::vector<UInt_t> FairFilteredPrimaryGenerator::fLogicalFilterOperation
protected

vector containing the logical operations with which the outputs of the non-veto event filters should be combined.

The vector grows automatically with every added non-veto event filter. It is used to combine multiple filters via && or ||. The expression is evaluated sequentally from the first registered filter to the last one disregarding operator priorities.

Definition at line 216 of file FairFilteredPrimaryGenerator.h.

Referenced by AddFilter(), and GenerateEvent().

Int_t FairFilteredPrimaryGenerator::fVerbose
protected

Level of commenting output, 0 means no output, higher gives more output.

Definition at line 201 of file FairFilteredPrimaryGenerator.h.

Referenced by GenerateEvent(), and SetVerbose().

TIterator* FairFilteredPrimaryGenerator::fVetoFilterIter
protected

Iterator over veto filter list.

Definition at line 192 of file FairFilteredPrimaryGenerator.h.

Referenced by GenerateEvent(), and ~FairFilteredPrimaryGenerator().

TObjArray* FairFilteredPrimaryGenerator::fVetoFilterList
protected

List of registered veto filters.

Definition at line 190 of file FairFilteredPrimaryGenerator.h.

Referenced by AddVetoFilter(), GetListOfVetoFilters(), Init(), and ~FairFilteredPrimaryGenerator().


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