PandaRoot
PndFilteredPrimaryGenerator.h
Go to the documentation of this file.
1 
43 #ifndef PndFilteredPrimaryGenerator_H
44 #define PndFilteredPrimaryGenerator_H
45 
46 #include "FairPrimaryGenerator.h"
47 #include "FairEvtFilterParams.h"
48 #include "FairRunSim.h"
49 #include "FairRootFileSink.h"
50 
51 #include "PndSmpFilt.h"
52 #include "PndSmpCand.h"
53 
54 #include "TFile.h" // for TFile
55 #include "TLorentzVector.h"
56 #include "TString.h"
57 #include "Rtypes.h" // for Double_t, Bool_t, Int_t, etc
58 
59 //#include <iosfwd> // for ostream
60 #include <iostream> // for operator<<, basic_ostream, etc
61 #include <vector>
62 #include <unordered_set>
63 #include <map>
64 
65 class FairGenericStack;
66 class TDatabasePDG;
67 
68 typedef vector<PndSmpFilt> PndSmpFilterSet;
69 typedef std::vector<PndSmpCand> PndSmpCandList;
70 typedef std::vector<TString> StrVec;
71 
72 class PndFilteredPrimaryGenerator : public FairPrimaryGenerator {
73 
74  public:
77 
79  PndFilteredPrimaryGenerator(TString filterset);
80 
83 
85  virtual Bool_t Init();
86 
92  void AddFilter(TString filterStr);
93 
104  virtual Bool_t GenerateEvent(FairGenericStack *pStack);
105 
107  void SetFilterMaxTries(Int_t maxTries = 99999)
108  {
109  if (maxTries > 0) {
110  fEvtFilterStat.fFilterMaxTries = maxTries;
111  std::cout << "PndFilteredPrimaryGenerator: maxTries is now set to " << fEvtFilterStat.fFilterMaxTries << "\n";
112  } else {
113  std::cout << "\n\n\n -WARNING from PndFilteredPrimaryGenerator: maxTries must be a positive number! Check your SetFilterMaxTries call!\n\n\n";
114  }
115  }
116 
119 
123 
125  void SetEventPrintFrequency(int freq) { fEventPrintFreq = freq; }
126 
136 
138  void WriteEvtFilterStatsToRootFile(TFile *outputFile = nullptr)
139  {
140  std::cout << "\n\nGenerated Events = " << GetNumberOfGeneratedEvents() << "\n";
141  if (0 < GetNumberOfFilterFailedEvents()) {
142  std::cout << "WARNING: Number of events where the event filter FAILED " << GetNumberOfFilterFailedEvents() << "\n\n\n";
143  std::cout << "Random events were accepted to avoid infinite loops. \n";
144  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";
145  }
146  if (outputFile == nullptr)
147  outputFile = ((FairRootFileSink*)FairRunSim::Instance()->GetSink())->GetRootFile();
148  //outputFile = FairRunSim::Instance()->GetOutputFile();
149  outputFile->cd();
150  outputFile->mkdir("FairEvtFilter");
151  outputFile->cd("FairEvtFilter");
152  fEvtFilterStat.Write();
153  outputFile->cd();
154  }
155 
159  void SetVerbose(Int_t verbose = 12)
160  {
161  if (verbose >= 0) {
162  fVerbose = verbose;
163  std::cout << "PndFilteredPrimaryGenerator: fVerbose is now set to " << verbose << "\n";
164  } else {
165  std::cout << "\n\n\n -WARNING from PndFilteredPrimaryGenerator: verbose must be a positive number! Check your SetVerbose call!\n\n\n";
166  }
167  }
168 
169  protected:
171  PndSmpCandList CombineList(int pdg, PndSmpCandList *l0, PndSmpCandList *l1, PndSmpCandList *l2 = nullptr, PndSmpCandList *l3 = nullptr, PndSmpCandList *l4 = nullptr);
172 
174  void PrintSmpCandList(PndSmpCandList l, TString name = "");
175 
177  StrVec SplitString(TString s, TString delim = " ");
178 
180  bool CheckKinematic(const PndSmpFilt &f, const TLorentzVector &p4);
181 
183  int AntiPdgCode(int pdg);
184 
186  void GetRangeDouble(TString s, double &a, double &b, TString delim = ",", bool forceset = false);
187 
189  void GetRangeInt(TString s, int &a, int &b, TString delim = "..");
190 
196  std::vector<PndSmpFilterSet> fFilterSets;
197 
200 
202  Int_t fVerbose;
203 
206 
209 
211  TDatabasePDG *fdbPdg;
212 
214  std::vector<TString> fPartNames;
215 
217  std::vector<int> fNamePdg;
218 
220  std::vector<int> fCombFsPdg;
221 
223  std::unordered_set<int> fSetFsPdg;
224 
226  std::map<TString, int> fNameCodeMap;
227 
229  std::map<int, TString> fCodeNameMap;
230 
231  private:
234 
235  ClassDef(PndFilteredPrimaryGenerator, 1);
236 };
237 
238 #endif
std::vector< PndSmpCand > PndSmpCandList
std::vector< TString > StrVec
StrVec SplitString(TString s, TString delim=" ")
Splits a TString to substrings.
std::unordered_set< int > fSetFsPdg
set to identify particles from MC truth list which can be combined
Int_t GetNumberOfFilterMaxTries()
returns the maximum number of times that this object should try to find an event which suits all even...
Simple container for filter definition (criteria) for PndFilteredPrimaryGenerator.
Definition: PndSmpFilt.h:19
PndFilteredPrimaryGenerator()
Default constructor.
virtual Bool_t GenerateEvent(FairGenericStack *pStack)
Calls event generators and the event filters.
void SetEventPrintFrequency(int freq)
Sets the frequency (accepted events) for printout (verbose>0) of accepted and generated events...
void WriteEvtFilterStatsToRootFile(TFile *outputFile=nullptr)
Writes all relevant event filter information to the output root file.
std::map< TString, int > fNameCodeMap
mapes names to (pdg) codes
bool CheckKinematic(const PndSmpFilt &f, const TLorentzVector &p4)
Checks whether P4 kinematics match the criteria of a PndSmpFilt.
Int_t fEventNrFiltered
Event number (Set by the filtered primary generator.
std::map< int, TString > fCodeNameMap
mapes (pdg) codes to names
TDatabasePDG * fdbPdg
Shortcut to TDatabasePDG.
std::vector< TString > fPartNames
particle names for particles to count (with and w/o charged specification, also simply tracks and neu...
Primary generator with added event filtering capabilities.
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.
vector< PndSmpFilt > PndSmpFilterSet
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...
virtual ~PndFilteredPrimaryGenerator()
Destructor.
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.
void AddFilter(TString filterStr)
Registers a filter as a string to be parsed Each filter set consist of one filter definition or a num...
std::vector< PndSmpFilterSet > fFilterSets
Contains the filter-sets. Each filter set consist of one filter definition or a number of filters con...
float f
Definition: P4_F32vec4.h:20
Int_t fEventPrintFreq
Print frequency for filtered events.
void PrintSmpCandList(PndSmpCandList l, TString name="")
Prints a candidate lits.
void SetVerbose(Int_t verbose=12)
Set the level of commenting output.
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.
virtual Bool_t Init()
Initialize the event generator(s) and the event (veto) filter(s).
Int_t fVerbose
Level of commenting output, 0 means no output, higher gives more output.
std::vector< int > fNamePdg
particle codes for particles to count (with and w/o charged specification, also simply tracks and neu...
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...
int AntiPdgCode(int pdg)
Gets anti-pdg code, if exists. If not returns the code itself (particle is its anti-particle) ...
std::vector< int > fCombFsPdg
particle codes for the lists used for combinatorics; !!! the codes are not selected from MC truth...
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.