MINT2
TimeDependentGenerator.h
Go to the documentation of this file.
1 #ifndef __TIMEDEPENDENTGENERATOR_H__
2 #define __TIMEDEPENDENTGENERATOR_H__
3 
4 #include <Mint/FitAmpSum.h>
5 #include <Mint/DalitzEvent.h>
6 #include <Mint/SignalGenerator.h>
8 #include <memory>
9 #include <list>
10 #include <map>
11 #include <complex>
12 #include <string>
13 #include <TSpline.h>
14 #include <Mint/Eff3piSymmetric.h>
15 #include <TH1F.h>
16 
17 class TRandom3 ;
18 
19 // Class to generate time dependent phase space events.
21 
22 public :
23  // Class to hold the tag (production flavour), decay time and
24  // Dalitz event that's generated.
25  class GenTimeEvent : public DalitzEvent {
26  public :
27  GenTimeEvent(const IDalitzEvent&, const int, const double,
28  const double smeareddecaytime = -999.) ;
29  GenTimeEvent(const IDalitzEvent&) ;
30  int getTag() const ;
31  void setTag(int) ;
32  double getDecayTime() const ;
33  void setDecayTime(double) ;
34  double getSmearedDecayTime() const ;
35  void setSmearedDecayTime(double) ;
37  static std::map<std::string, unsigned> infoNames ;
38  private :
39  static std::map<std::string, unsigned> makeInfoNames() ;
40  } ;
41 
44 
57  TimeDependentGenerator(const DalitzEventPattern& pattern, double width, double deltam,
58  double deltagamma, double qoverp, double phi,
59  TRandom3* rndm, TH1F* h_efficiency = NULL,
60  float resWidth = 0.05, bool addExpEffects = false,
61  Eff3piSymmetric* sEfficiency = NULL) ;
76  double width, double deltam,
77  double deltagamma, double qoverp, double phi,
78  TRandom3*, TH1F* h_efficiency = NULL,
79  float resWidth = 0.05, bool addExpEffects = false,
80  Eff3piSymmetric* sEfficiency = NULL) ;
81 
83  std::pair<double, double> generate_decay_time() const ;
84 
87 
89  MINT::counted_ptr<IDalitzEvent> generate_event(const double s13, const double s23) ;
90 
92  bool accept_or_reject(const int tag, const double decaytime, IDalitzEvent&) ;
93 
94  double get_scale() const ;
95  float get_gen_efficiency() const ;
96 
98  double pdf_value(int, double, IDalitzEvent&) ;
101  double pdf_value(IDalitzEvent&) ;
104  double envelope_value(const double, IDalitzEvent&) ;
105 
106  typedef std::pair<std::complex<double>, std::complex<double> > AmpPair ;
109  AmpPair amplitude_coefficients(const int tag, const double decaytime) ;
110 
113  TH1F draw_pdf_vs_time(IDalitzEvent&, unsigned, float, float, const std::string& name = "pdf_vs_time") ;
115  TH1F draw_envelope_vs_time(IDalitzEvent&, unsigned, float, float, const std::string& name = "envelope_vs_time") ;
116 
117 private :
118  void init() ;
119 
120  TRandom3* m_rndm ;
127 
128  const double m_width ;
129  const double m_deltam ;
130  const double m_deltagamma ;
131  const std::complex<double> m_qoverp ;
132  double m_scale ;
133  unsigned m_ngen ;
134  unsigned m_naccept ;
135 
137  TSpline3 m_efficiencyFit;
138 
139  float m_resWidth;
141 
143 } ;
144 
145 #endif
FitAmpIncoherentSum m_bothmodel
double envelope_value(const double, IDalitzEvent &)
double pdf_value(int, double, IDalitzEvent &)
static std::map< std::string, unsigned > infoNames
std::pair< double, double > generate_decay_time() const
Generate a decay time, optionally including experimental effects.
static DalitzEventPattern anti(DalitzEventPattern pat)
Take the CP conjugate of the head of the decay pattern.
TimeDependentGenerator(const DalitzEventPattern &pattern, double width, double deltam, double deltagamma, double qoverp, double phi, TRandom3 *rndm, TH1F *h_efficiency=NULL, float resWidth=0.05, bool addExpEffects=false, Eff3piSymmetric *sEfficiency=NULL)
bool accept_or_reject(const int tag, const double decaytime, IDalitzEvent &)
Decide whether the generated event should be accepted or rejected.
static std::map< std::string, unsigned > makeInfoNames()
MINT::counted_ptr< FitAmpSum > m_cpmodel
TH1F draw_envelope_vs_time(IDalitzEvent &, unsigned, float, float, const std::string &name="envelope_vs_time")
TH1F draw_pdf_vs_time(IDalitzEvent &, unsigned, float, float, const std::string &name="pdf_vs_time")
MINT::counted_ptr< IDalitzEvent > generate_event()
Generate a flavour, decay time and Dalitz event.
AmpPair amplitude_coefficients(const int tag, const double decaytime)
const std::complex< double > m_qoverp
const DalitzEventPattern m_pattern
GenTimeEvent(const IDalitzEvent &, const int, const double, const double smeareddecaytime=-999.)
MINT::counted_ptr< FitAmpSum > m_model
std::pair< std::complex< double >, std::complex< double > > AmpPair
const DalitzEventPattern m_cppattern