AmpGen 2.1
Loading...
Searching...
No Matches
CoherentSum.h
Go to the documentation of this file.
1#ifndef AMPGEN_COHERENTSUM_H
2#define AMPGEN_COHERENTSUM_H
3
4#include <memory.h>
5#include <stddef.h>
6#include <complex>
7#include <map>
8#include <memory>
9#include <string>
10#include <utility>
11#include <vector>
12
16#include "AmpGen/EventList.h"
18#include "AmpGen/EventType.h"
19#include "AmpGen/Integrator.h"
20#include "AmpGen/Types.h"
21#include "AmpGen/Event.h"
22#include "AmpGen/Projection.h"
24#include "AmpGen/Store.h"
26namespace AmpGen
27{
30 class FitFraction;
31 class Particle;
32
45 {
46 public:
47 #if ENABLE_AVX
49 #else
51 #endif
53 CoherentSum( const EventType& type, const AmpGen::MinuitParameterSet& mps, const std::string& prefix = "" );
54 virtual ~CoherentSum();
55
56 std::string prefix() const { return m_prefix; }
57
58 auto& operator[]( const size_t& index ) { return m_matrixElements[index]; }
59 auto& operator[]( const size_t& index ) const { return m_matrixElements[index]; }
60 size_t size() const { return m_matrixElements.size(); }
61 real_t getWeight() const { return m_weight; }
62 real_t norm( const Bilinears& norms ) const;
63 real_t norm() const;
64 real_t getNorm( const Bilinears& normalisations );
65
66 complex_t norm( const size_t& x, const size_t& y ) const;
67 complex_t getVal( const Event& evt ) const;
68 complex_t getValNoCache( const Event& evt ) const;
69
71 void prepare();
72 void printVal( const Event& evt );
74 void setWeight( MinuitProxy param ) { m_weight = param; }
76 void reset( bool resetEvents = false );
77 void setEvents( const EventList_type& list );
78 void setMC( const EventList_type& sim );
79 #if ENABLE_AVX
80 void setEvents( const EventList& list) {
81 WARNING("Setting events from a AoS container, will need to make a copy");
82 m_ownEvents = true; setEvents( *(new EventListSIMD(list)) ) ; }
83 void setMC ( const EventList& list) {
84 WARNING("Setting integration events from a AoS container, will need to make a copy");
85 setMC( *(new EventListSIMD(list)) ) ; }
86 //double operator()(const double*, const unsigned) const;
87 #endif
88
89 real_v operator()(const real_v*, const unsigned) const;
90 real_t operator()(const Event& evt ) const { return m_weight*std::norm(getVal(evt))/m_norm; }
91
92 void debug( const Event& evt, const std::string& nameMustContain="");
93 void generateSourceCode( const std::string& fname, const double& normalisation = 1, bool add_mt = false );
94
95 std::vector<FitFraction> fitFractions( const LinearErrorPropagator& linProp );
96 auto matrixElements() const { return m_matrixElements; }
97
98 std::map<std::string, std::vector<unsigned int>> getGroupedAmplitudes();
99 Bilinears norms() const { return m_normalisations ; }
100
101 std::function<real_t(const Event&)> evaluator(const EventList_type* = nullptr) const;
102 std::function<complex_t(const Event&)> amplitudeEvaluator(const EventList_type* = nullptr) const;
103 KeyedFunctors<double(Event)> componentEvaluator(const EventList_type* = nullptr) const;
104 EventType eventType() const { return m_eventType; }
105 const auto& cache() const { return m_cache; }
106 protected:
107 std::vector<MatrixElement> m_matrixElements;
109
111 const EventList_type* m_events = {nullptr};
113
114 bool m_ownEvents = {false};
116 size_t m_prepareCalls = {0};
117 size_t m_lastPrint = {0};
118 size_t m_printFreq = {0};
119 MinuitProxy m_weight = {nullptr, 1};
120 double m_norm = {1};
121 bool m_isConstant = {false};
122 bool m_dbThis = {false};
123 bool m_verbosity = {false};
124 std::string m_objCache = {""};
125 std::string m_prefix = {""};
126 const MinuitParameterSet* m_mps = {nullptr};
127 void addMatrixElement( std::pair<Particle, TotalCoupling>& particleWithCoupling, const MinuitParameterSet& mps );
128 };
129} // namespace AmpGen
130
131#endif
real_t norm(const Bilinears &norms) const
std::string prefix() const
Definition CoherentSum.h:56
void debug(const Event &evt, const std::string &nameMustContain="")
size_t size() const
Definition CoherentSum.h:60
real_t norm() const
auto matrixElements() const
Definition CoherentSum.h:96
std::string m_objCache
Directory that contains (cached) amplitude objects.
Bilinears norms() const
Definition CoherentSum.h:99
bool m_dbThis
Flag to generate amplitude level debugging.
void printVal(const Event &evt)
CoherentSum(const EventType &type, const AmpGen::MinuitParameterSet &mps, const std::string &prefix="")
bool m_ownEvents
Flag as to whether events are owned by this PDF or not.
Integrator m_integrator
Tool to calculate integrals.
std::string m_prefix
Prefix for matrix elements.
complex_t getValNoCache(const Event &evt) const
void setMC(const EventList_type &sim)
std::function< real_t(const Event &)> evaluator(const EventList_type *=nullptr) const
void setWeight(MinuitProxy param)
Definition CoherentSum.h:74
MinuitProxy m_weight
Weight (i.e. the normalised yield)
void setEvents(const EventList_type &list)
std::vector< FitFraction > fitFractions(const LinearErrorPropagator &linProp)
void addMatrixElement(std::pair< Particle, TotalCoupling > &particleWithCoupling, const MinuitParameterSet &mps)
auto & operator[](const size_t &index) const
Definition CoherentSum.h:59
EventList EventList_type
Definition CoherentSum.h:50
bool m_verbosity
Flag for verbose printing.
size_t m_printFreq
Frequency to print verbose PDF info.
Bilinears m_normalisations
Normalisation integrals.
complex_t getVal(const Event &evt) const
const auto & cache() const
EventType eventType() const
EventType m_eventType
Final state for this amplitude.
real_t getNorm(const Bilinears &normalisations)
const EventList_type * m_events
Data events to evaluate PDF on.
size_t m_lastPrint
Last time verbose PDF info was printed.
real_t operator()(const Event &evt) const
Definition CoherentSum.h:90
const MinuitParameterSet * m_mps
virtual ~CoherentSum()
std::vector< MatrixElement > m_matrixElements
Vector of matrix elements.
real_t getWeight() const
Definition CoherentSum.h:61
void reset(bool resetEvents=false)
FunctionCache< EventList_type, complex_v, Alignment::AoS > m_cache
Store of intermediate values for the PDF calculation.
real_v operator()(const real_v *, const unsigned) const
std::map< std::string, std::vector< unsigned int > > getGroupedAmplitudes()
KeyedFunctors< double(Event)> componentEvaluator(const EventList_type *=nullptr) const
double m_norm
Normalisation integral.
bool m_isConstant
Flag for a constant PDF.
auto & operator[](const size_t &index)
Definition CoherentSum.h:58
complex_t norm(const size_t &x, const size_t &y) const
std::function< complex_t(const Event &)> amplitudeEvaluator(const EventList_type *=nullptr) const
size_t m_prepareCalls
Number of times prepare has been called.
void generateSourceCode(const std::string &fname, const double &normalisation=1, bool add_mt=false)
Encapsulates the final state particles of a single event.
Definition Event.h:18
Deals with final state configuration of events, specifically dealing with the ordering of particles i...
Definition EventType.h:22
Propagates uncertainties on functors using either a MinuitParameterSet (thus assuming a diagonal cova...
Describes a particle, its decay process and subsequent decay products, which are also Particles.
Definition Particle.h:103
#define WARNING(X)
Used for printing warning messages, can be switched off using WARNINGLEVEL.
Definition MsgService.h:97
double real_t
Definition Types.h:6
std::complex< real_t > complex_t
Definition Types.h:7
AVX::real_v real_v
Definition utils.h:46