FairRoot/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 
50 #include "PndSmpFilt.h"
51 #include "PndSmpCand.h"
52 
53 #include "TFile.h" // for TFile
54 #include "TLorentzVector.h"
55 #include "TString.h"
56 #include "Rtypes.h" // for Double_t, Bool_t, Int_t, etc
57 
58 //#include <iosfwd> // for ostream
59 #include <iostream> // for operator<<, basic_ostream, etc
60 #include <vector>
61 #include <unordered_set>
62 #include <map>
63 
64 
65 class FairGenericStack;
66 class TDatabasePDG;
67 
68 
69 typedef vector<PndSmpFilt> PndSmpFilterSet;
70 typedef std::vector<PndSmpCand> PndSmpCandList;
71 typedef std::vector<TString> StrVec;
72 
73 
74 class PndFilteredPrimaryGenerator : public FairPrimaryGenerator
75 {
76 
77 public:
78 
81 
84 
87 
89  virtual Bool_t Init();
90 
96  void AddFilter(TString filterStr);
97 
108  virtual Bool_t GenerateEvent(FairGenericStack* pStack);
109 
111  void SetFilterMaxTries(Int_t maxTries=99999)
112  {
113  if ( maxTries > 0 ){
114  fEvtFilterStat.fFilterMaxTries = maxTries;
115  std::cout << "PndFilteredPrimaryGenerator: maxTries is now set to " << fEvtFilterStat.fFilterMaxTries << "\n";
116  } else {
117  std::cout << "\n\n\n -WARNING from PndFilteredPrimaryGenerator: maxTries must be a positive number! Check your SetFilterMaxTries call!\n\n\n";
118  }
119  }
120 
124  }
125 
129  }
130 
132  void SetEventPrintFrequency(int freq)
133  {
134  fEventPrintFreq = freq;
135  }
136 
147  }
148 
150  void WriteEvtFilterStatsToRootFile(TFile* outputFile=NULL){
151  std::cout << "\n\nGenerated Events = " << GetNumberOfGeneratedEvents() << "\n";
152  if (0 < GetNumberOfFilterFailedEvents() ) {
153  std::cout << "WARNING: Number of events where the event filter FAILED " << GetNumberOfFilterFailedEvents() << "\n\n\n";
154  std::cout << "Random events were accepted to avoid infinite loops. \n";
155  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";
156  }
157  if(outputFile==NULL) outputFile = FairRunSim::Instance()->GetOutputFile();
158  outputFile->cd();
159  outputFile->mkdir("FairEvtFilter");
160  outputFile->cd("FairEvtFilter");
161  fEvtFilterStat.Write();
162  outputFile->cd();
163  }
164 
168  void SetVerbose(Int_t verbose=12)
169  {
170  if ( verbose >= 0 )
171  {
172  fVerbose = verbose;
173  std::cout << "PndFilteredPrimaryGenerator: fVerbose is now set to " << verbose << "\n";
174  }
175  else
176  {
177  std::cout << "\n\n\n -WARNING from PndFilteredPrimaryGenerator: verbose must be a positive number! Check your SetVerbose call!\n\n\n";
178  }
179  }
180 
181 
182 
183 protected:
184 
187 
190 
192  StrVec SplitString(TString s, TString delim=" ");
193 
195  bool CheckKinematic(const PndSmpFilt &f, const TLorentzVector &p4);
196 
198  int AntiPdgCode(int pdg);
199 
201  void GetRangeDouble(TString s, double &a, double &b, TString delim=",", bool forceset=false);
202 
204  void GetRangeInt(TString s, int &a, int &b, TString delim="..");
205 
211  std::vector<PndSmpFilterSet> fFilterSets;
212 
215 
217  Int_t fVerbose;
218 
221 
224 
226  TDatabasePDG *fdbPdg;
227 
229  std::vector<TString> fPartNames;
230 
232  std::vector<int> fNamePdg;
233 
235  std::vector<int> fCombFsPdg;
236 
238  std::unordered_set<int> fSetFsPdg;
239 
241  std::map<TString, int> fNameCodeMap;
242 
244  std::map<int, TString> fCodeNameMap;
245 
246 
247 private:
250 
251 
253 };
254 
255 #endif
std::vector< PndSmpCand > PndSmpCandList
TTree * b
void WriteEvtFilterStatsToRootFile(TFile *outputFile=NULL)
Writes all relevant event filter information to the output root file.
#define verbose
StrVec SplitString(TString s, TString delim=" ")
Splits a TString to substrings.
TLorentzVector s
Definition: Pnd2DStar.C:50
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.
void SetEventPrintFrequency(int freq)
Sets the frequency (accepted events) for printout (verbose&gt;0) of accepted and generated events...
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.
virtual Bool_t Init()
Initialize the event generator(s) and the event (veto) filter(s).
std::map< int, TString > fCodeNameMap
mapes (pdg) codes to names
Int_t a
Definition: anaLmdDigi.C:126
TDatabasePDG * fdbPdg
Shortcut to TDatabasePDG.
ClassDef(PndFilteredPrimaryGenerator, 1)
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.
TFile * f
Definition: bump_analys.C:12
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.
PndFilteredPrimaryGenerator & operator=(const PndFilteredPrimaryGenerator &)
vector< PndSmpFilt > PndSmpFilterSet
virtual ~PndFilteredPrimaryGenerator()
Destructor.
std::vector< TString > StrVec
TString name
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.
virtual Bool_t GenerateEvent(FairGenericStack *pStack)
Calls event generators and the event filters.
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...
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 &lt;int&gt;&lt;delim&gt;&lt;int&gt; (e.g. &quot;3..8&quot;) to a range a,b.
PndSmpCandList CombineList(int pdg, PndSmpCandList *l0, PndSmpCandList *l1, PndSmpCandList *l2=0, PndSmpCandList *l3=0, PndSmpCandList *l4=0)
Combines upt to five particle lists of PndSmpCand with overlap and double counting prevention...
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 &lt;float&gt;&lt;delim&gt;&lt;float&gt; (e.g. &quot;124.2,178.3&quot;) to a range a,b.