![]() |
PandaRoot
|
Primary generator with added event filtering capabilities. More...
#include <PndFilteredPrimaryGenerator.h>
Public Member Functions | |
PndFilteredPrimaryGenerator () | |
Default constructor. More... | |
PndFilteredPrimaryGenerator (TString filterset) | |
Constructor with filter-set string. More... | |
virtual | ~PndFilteredPrimaryGenerator () |
Destructor. More... | |
virtual Bool_t | Init () |
Initialize the event generator(s) and the event (veto) filter(s). More... | |
void | AddFilter (TString filterStr) |
Registers a filter as a string to be parsed Each filter set consist of one filter definition or a number of filters connected with logical && (AND) (= a certain event type). All filter sets are combined with a logical || (OR) (all event types to be accepted). Syntax: "Filt_1.1 && Filt_1.2 || Filt_2.1 || Filt_3.1 && Filt_3.2 && Filt_3.2". More... | |
virtual Bool_t | GenerateEvent (FairGenericStack *pStack) |
Calls event generators and the event filters. More... | |
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=nullptr) |
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 | |
PndSmpCandList | CombineList (int pdg, PndSmpCandList *l0, PndSmpCandList *l1, PndSmpCandList *l2=nullptr, PndSmpCandList *l3=nullptr, PndSmpCandList *l4=nullptr) |
Combines upt to five particle lists of PndSmpCand with overlap and double counting prevention. More... | |
void | PrintSmpCandList (PndSmpCandList l, TString name="") |
Prints a candidate lits. More... | |
StrVec | SplitString (TString s, TString delim=" ") |
Splits a TString to substrings. More... | |
bool | CheckKinematic (const PndSmpFilt &f, const TLorentzVector &p4) |
Checks whether P4 kinematics match the criteria of a PndSmpFilt. More... | |
int | AntiPdgCode (int pdg) |
Gets anti-pdg code, if exists. If not returns the code itself (particle is its anti-particle) More... | |
void | GetRangeDouble (TString s, double &a, double &b, TString delim=",", bool forceset=false) |
Turns a string of the form <float><delim><float> (e.g. "124.2,178.3") to a range a,b. More... | |
void | GetRangeInt (TString s, int &a, int &b, TString delim="..") |
Turns a string of the form <int><delim><int> (e.g. "3..8") to a range a,b. More... | |
Protected Attributes | |
std::vector< PndSmpFilterSet > | fFilterSets |
Contains the filter-sets. Each filter set consist of one filter definition or a number of filters connected with logical AND (= a certain event type). All filter sets are combined with a logical OR (all event types to be accepted). Structure is like: (Filt_1.1 AND Filt_1.2) OR (Filt_2.1) OR (Filt_3.1 AND Filt_3.2 AND Filt_3.2) 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... | |
Int_t | fEventNrFiltered |
Event number (Set by the filtered primary generator. More... | |
Int_t | fEventPrintFreq |
Print frequency for filtered events. More... | |
TDatabasePDG * | fdbPdg |
Shortcut to TDatabasePDG. More... | |
std::vector< TString > | fPartNames |
particle names for particles to count (with and w/o charged specification, also simply tracks and neutrals) More... | |
std::vector< int > | fNamePdg |
particle codes for particles to count (with and w/o charged specification, also simply tracks and neutrals) More... | |
std::vector< int > | fCombFsPdg |
particle codes for the lists used for combinatorics; !!! the codes are not selected from MC truth, just the corresponding mass hypos are applied to charged tracks More... | |
std::unordered_set< int > | fSetFsPdg |
set to identify particles from MC truth list which can be combined More... | |
std::map< TString, int > | fNameCodeMap |
mapes names to (pdg) codes More... | |
std::map< int, TString > | fCodeNameMap |
mapes (pdg) codes to names More... | |
Primary generator with added event filtering capabilities.
This class adds event filtering capabilities to FairPrimaryGenerator which is used internally for handling the event generators. The event filtering is performed after the event generation and before the particle transport through the detector model.
The filter is configured with a string of the form "Filt_1.1 && Filt_1.2 || Filt_2.1 || Filt_3.1 && Filt_3.2 && Filt_3.2", where the logical && (AND) has higher priority than the logical || (OR). Each filter can be negated with a leading '!'.
There are two types of filters possible: particle type counters and invariant mass filters. The syntax is
count filter : C(<particle> ; <mult> ; [kinematic requirements]) ==> requires <mult> particles of type <particle> (counted in MC list) to meet certain kinematic conditions (p, pt, pz, theta, phi)
mass filter : M(<particle 1> <particle 2> [<particle 3> ... <particle 5>] ; <mult> ; <mass requirement>=""> ; [kinematic requirements]) ==> takes neutral particles and charged tracks, applies the mass hypotheses from <particle 1..5>, combines them and requires <mult> candidates with invariant mass to be in certain window
<particle> : e+, e-, e+-, mu+, mu-, mu+-, pi+, pi-, pi+-, K+, K-, K+-, p+, p-, p+-, n0, n0b, n00b (=n0 or n0b), gam, nt (neutral = gamma), t+ (positive tracks), t- (neg. tracks), t+- (chg. tracks), any (total multiplicity) <mult> : a..b, a (= a..a), ..a (= 0..a), a.. = (a..10000), not given (>=1 = 1..10000) <mass req.> : m[center,window] : center - window/2 < inv. mass < center + window/2 [kin. req.] : p[min,max] ; pt[min,max] ; pz[min,max] ; tht[min,max] ; phi[min,max]. Can be shortened to [,max] = [0.0,max] and [min,] = [min,1e8]
Example: C(any ; 4..) && C(tr+- ; 2 ; tht[5,20] ; p[1.0,]) && M(pi+ pi- gam gam ; m[0.782,0.1] ; pt[0.5,2.0])
Concerning FairPrimaryGenerator:
The PndFilteredPrimaryGenerator 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 PndFilteredPrimaryGenerator defines position and (optionally) smearing of the primary vertex. This class should be instantiated only once.
Definition at line 72 of file PndFilteredPrimaryGenerator.h.
PndFilteredPrimaryGenerator::PndFilteredPrimaryGenerator | ( | ) |
Default constructor.
PndFilteredPrimaryGenerator::PndFilteredPrimaryGenerator | ( | TString | filterset | ) |
Constructor with filter-set string.
|
inlinevirtual |
Destructor.
Definition at line 82 of file PndFilteredPrimaryGenerator.h.
References AddFilter(), GenerateEvent(), and Init().
void PndFilteredPrimaryGenerator::AddFilter | ( | TString | filterStr | ) |
Registers a filter as a string to be parsed Each filter set consist of one filter definition or a number of filters connected with logical && (AND) (= a certain event type). All filter sets are combined with a logical || (OR) (all event types to be accepted). Syntax: "Filt_1.1 && Filt_1.2 || Filt_2.1 || Filt_3.1 && Filt_3.2 && Filt_3.2".
Referenced by ~PndFilteredPrimaryGenerator().
|
protected |
Gets anti-pdg code, if exists. If not returns the code itself (particle is its anti-particle)
Referenced by SetVerbose().
|
protected |
Checks whether P4 kinematics match the criteria of a PndSmpFilt.
Referenced by SetVerbose().
|
protected |
Combines upt to five particle lists of PndSmpCand with overlap and double counting prevention.
Referenced by SetVerbose().
|
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.
pStack | The particle stack |
Referenced by ~PndFilteredPrimaryGenerator().
|
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 135 of file PndFilteredPrimaryGenerator.h.
References fEvtFilterStat, and FairEvtFilterParams::fFailedFilterEvents.
Referenced by WriteEvtFilterStatsToRootFile().
|
inline |
returns the maximum number of times that this object should try to find an event which suits all event filters.
Definition at line 118 of file PndFilteredPrimaryGenerator.h.
References fEvtFilterStat, and FairEvtFilterParams::fFilterMaxTries.
|
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 122 of file PndFilteredPrimaryGenerator.h.
References fEvtFilterStat, and FairEvtFilterParams::fGeneratedEvents.
Referenced by WriteEvtFilterStatsToRootFile().
|
protected |
Turns a string of the form <float><delim><float> (e.g. "124.2,178.3") to a range a,b.
Referenced by SetVerbose().
|
protected |
Turns a string of the form <int><delim><int> (e.g. "3..8") to a range a,b.
Referenced by SetVerbose().
|
virtual |
Initialize the event generator(s) and the event (veto) filter(s).
Referenced by ~PndFilteredPrimaryGenerator().
|
protected |
Prints a candidate lits.
Referenced by SetVerbose().
|
inline |
Sets the frequency (accepted events) for printout (verbose>0) of accepted and generated events.
Definition at line 125 of file PndFilteredPrimaryGenerator.h.
References fEventPrintFreq.
|
inline |
Define the maximum number of times that this object should try to find an event which suits all event filters.
Definition at line 107 of file PndFilteredPrimaryGenerator.h.
References fEvtFilterStat, and FairEvtFilterParams::fFilterMaxTries.
|
inline |
Set the level of commenting output.
verbose | Level of commenting output, 0 means no output, higher gives more output. |
Definition at line 159 of file PndFilteredPrimaryGenerator.h.
References AntiPdgCode(), CheckKinematic(), CombineList(), f, fVerbose, GetRangeDouble(), GetRangeInt(), PrintSmpCandList(), and SplitString().
|
protected |
Splits a TString to substrings.
Referenced by SetVerbose().
|
inline |
Writes all relevant event filter information to the output root file.
Definition at line 138 of file PndFilteredPrimaryGenerator.h.
References fEvtFilterStat, GetNumberOfFilterFailedEvents(), and GetNumberOfGeneratedEvents().
|
protected |
mapes (pdg) codes to names
Definition at line 229 of file PndFilteredPrimaryGenerator.h.
|
protected |
particle codes for the lists used for combinatorics; !!! the codes are not selected from MC truth, just the corresponding mass hypos are applied to charged tracks
Definition at line 220 of file PndFilteredPrimaryGenerator.h.
|
protected |
Shortcut to TDatabasePDG.
Definition at line 211 of file PndFilteredPrimaryGenerator.h.
|
protected |
Event number (Set by the filtered primary generator.
Definition at line 205 of file PndFilteredPrimaryGenerator.h.
|
protected |
Print frequency for filtered events.
Definition at line 208 of file PndFilteredPrimaryGenerator.h.
Referenced by SetEventPrintFrequency().
|
protected |
Contains the statistics of the event filtering process.
Definition at line 199 of file PndFilteredPrimaryGenerator.h.
Referenced by GetNumberOfFilterFailedEvents(), GetNumberOfFilterMaxTries(), GetNumberOfGeneratedEvents(), SetFilterMaxTries(), and WriteEvtFilterStatsToRootFile().
|
protected |
Contains the filter-sets. Each filter set consist of one filter definition or a number of filters connected with logical AND (= a certain event type). All filter sets are combined with a logical OR (all event types to be accepted). Structure is like: (Filt_1.1 AND Filt_1.2) OR (Filt_2.1) OR (Filt_3.1 AND Filt_3.2 AND Filt_3.2)
Definition at line 196 of file PndFilteredPrimaryGenerator.h.
|
protected |
mapes names to (pdg) codes
Definition at line 226 of file PndFilteredPrimaryGenerator.h.
|
protected |
particle codes for particles to count (with and w/o charged specification, also simply tracks and neutrals)
Definition at line 217 of file PndFilteredPrimaryGenerator.h.
|
protected |
particle names for particles to count (with and w/o charged specification, also simply tracks and neutrals)
Definition at line 214 of file PndFilteredPrimaryGenerator.h.
|
protected |
set to identify particles from MC truth list which can be combined
Definition at line 223 of file PndFilteredPrimaryGenerator.h.
|
protected |
Level of commenting output, 0 means no output, higher gives more output.
Definition at line 202 of file PndFilteredPrimaryGenerator.h.
Referenced by SetVerbose().