YAP
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Model.h
Go to the documentation of this file.
1 /* YAP - Yet another PWA toolkit
2  Copyright 2015, Technische Universitaet Muenchen,
3  Authors: Daniel Greenwald, Johannes Rauch
4 
5  This program is free software: you can redistribute it and/or modify
6  it under the terms of the GNU General Public License as published by
7  the Free Software Foundation, either version 3 of the License, or
8  (at your option) any later version.
9 
10  This program is distributed in the hope that it will be useful,
11  but WITHOUT ANY WARRANTY; without even the implied warranty of
12  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  GNU General Public License for more details.
14 
15  You should have received a copy of the GNU General Public License
16  along with this program. If not, see <http://www.gnu.org/licenses/>.
17  */
18 
20 
21 #ifndef yap_Model_h
22 #define yap_Model_h
23 
24 #include "fwd/DataAccessor.h"
25 #include "fwd/DataPartition.h"
26 #include "fwd/DataPoint.h"
27 #include "fwd/DataSet.h"
28 #include "fwd/DecayingParticle.h"
29 #include "fwd/Filters.h"
30 #include "fwd/FinalStateParticle.h"
31 #include "fwd/FourMomenta.h"
32 #include "fwd/FourVector.h"
33 #include "fwd/FreeAmplitude.h"
34 #include "fwd/HelicityAngles.h"
35 #include "fwd/MassAxes.h"
36 #include "fwd/MeasuredBreakupMomenta.h"
37 #include "fwd/Parameter.h"
38 #include "fwd/Particle.h"
39 #include "fwd/RecalculableDataAccessor.h"
40 #include "fwd/StaticDataAccessor.h"
41 #include "fwd/StatusManager.h"
42 
43 #include "CoordinateSystem.h"
45 #include "SpinAmplitudeCache.h"
46 
47 #include <complex>
48 #include <memory>
49 #include <vector>
50 
51 namespace yap {
52 
56 class Model
57 {
58 public:
59 
63  Model(std::unique_ptr<SpinAmplitudeCache> SAC, bool include_phsp_factors = true);
64 
67  Model(const Model&) = delete;
68 
71  Model& operator=(const Model&) = delete;
72 
75  Model(Model&&) = delete;
76 
79  Model& operator=(Model&&) = delete;
80 
84  void calculate(DataPartition& D) const;
85 
87  virtual bool consistent() const;
88 
90  bool locked() const
91  { return Locked_; }
92 
96  void lock();
97 
100 
102  const bool includePhaseSpaceFactors() const
103  { return IncludePhaseSpaceFactors_; }
104 
106  const CoordinateSystem<double, 3>& coordinateSystem() const
107  { return CoordinateSystem_; }
108 
110  std::shared_ptr<FourMomenta> fourMomenta()
111  { return FourMomenta_; }
112 
114  const std::shared_ptr<FourMomenta> fourMomenta() const
115  { return FourMomenta_; }
116 
118  std::shared_ptr<MeasuredBreakupMomenta> measuredBreakupMomenta()
119  { return MeasuredBreakupMomenta_; }
120 
122  const std::shared_ptr<MeasuredBreakupMomenta> measuredBreakupMomenta() const
123  { return MeasuredBreakupMomenta_; }
124 
127  std::shared_ptr<HelicityAngles> helicityAngles()
128  { return HelicityAngles_; }
129 
132  const std::shared_ptr<HelicityAngles> helicityAngles() const
133  { return HelicityAngles_; }
134 
137  { return ParticleCombinationCache_; }
138 
141  { return ParticleCombinationCache_; }
142 
145  { return SpinAmplitudeCache_.get(); }
146 
149  { return SpinAmplitudeCache_.get(); }
150 
152  const FinalStateParticleVector& finalStateParticles() const
153  { return FinalStateParticles_; }
154 
156  const InitialStateParticleMap& initialStateParticles() const
157  { return InitialStateParticles_; }
158 
160  const DataAccessorSet& dataAccessors() const
161  { return DataAccessors_; }
162 
164 
167 
173  void setFinalState(const FinalStateParticleVector& FSP);
174 
180  template <typename ... Types>
181  void setFinalState(Types ... FSPs)
182  { FinalStateParticleVector V{FSPs...}; setFinalState(V); }
183 
188  const InitialStateParticleMap::value_type& addInitialStateParticle(std::shared_ptr<DecayingParticle> bg);
189 
193  void setFinalStateMomenta(DataPoint& d, const std::vector<FourVector<double> >& P, StatusManager& sm) const;
194 
196  void setCoordinateSystem(const CoordinateSystem<double, 3>& cs);
197 
201 
203 
206 
213  const MassAxes massAxes(std::vector<std::vector<unsigned> > pcs = {});
214 
216 
218  void requireHelicityAngles();
219 
222 
225  DataSet createDataSet(size_t n = 0);
226 
229 
231  void printDataAccessors(bool printParticleCombinations = true) const;
232 
234  void printFlags(const StatusManager& sm) const;
235 
237  friend class DataAccessor;
238 
241 
243  friend class StaticDataAccessor;
244 
246  friend class DecayingParticle;
247 
248 protected:
249 
253  virtual void addParticleCombination(std::shared_ptr<ParticleCombination> pc);
254 
255  /* /// remove a DataAccessor from this Model */
256  /* virtual void removeDataAccessor(DataAccessorSet::value_type da); */
257 
258 private:
259 
262  bool Locked_;
263 
265  CoordinateSystem<double, 3> CoordinateSystem_;
266 
269 
272 
274  std::unique_ptr<SpinAmplitudeCache> SpinAmplitudeCache_;
275 
277  DataAccessorSet DataAccessors_;
278 
281  StaticDataAccessorVector StaticDataAccessors_;
282 
284  RecalculableDataAccessorSet RecalculableDataAccessors_;
285 
288  InitialStateParticleMap InitialStateParticles_;
289 
291  FinalStateParticleVector FinalStateParticles_;
292 
294  std::shared_ptr<FourMomenta> FourMomenta_;
295 
297  std::shared_ptr<MeasuredBreakupMomenta> MeasuredBreakupMomenta_;
298 
300  std::shared_ptr<HelicityAngles> HelicityAngles_;
301 
302 };
303 
305 std::string to_string(const AdmixtureMap& mix);
306 
308 const bool decays_to_full_final_state(const Particle& p);
309 
313 std::vector<std::shared_ptr<DecayingParticle> > full_final_state_isp(const Model& M);
314 
316 const double intensity(const InitialStateParticleMap::value_type& isp_mix, const DataPoint& d);
317 
319 const double intensity(const InitialStateParticleMap& isp_map, const DataPoint& d);
320 
325 const double sum_of_log_intensity(const Model& M, DataPartition& D, double ped = 0);
326 
330 const double sum_of_log_intensity(const Model& M, DataPartitionVector& DP, double ped = 0);
331 
333 FreeAmplitudeSet free_amplitudes(const Model& M);
334 
336 template <typename Last, typename ... UnaryPredicates>
337 FreeAmplitudeSet free_amplitudes(const Model& M, Last p, UnaryPredicates ... P)
338 { return filter(free_amplitudes(M), p, P...); }
339 
342 template <typename Last, typename ... UnaryPredicates>
343 FreeAmplitudeSet::value_type free_amplitude(const Model& M, Last p, UnaryPredicates ... P)
344 { return lone_elt(filter(free_amplitudes(M), p, P...)); }
345 
347 ParticleSet particles(const Model& M);
348 
350 template <typename Last, typename ... UnaryPredicates>
351 ParticleSet particles(const Model& M, Last p, UnaryPredicates ... P)
352 { return filter(particles(M), p, P...); }
353 
356 template <typename Last, typename ... UnaryPredicates>
357 ParticleSet::value_type particle(const Model& M, Last p, UnaryPredicates ... P)
358 { return lone_elt(filter(particles(M), p, P...)); }
359 
360 }
361 
362 #endif
SpinAmplitudeCache * spinAmplitudeCache()
Definition: Model.h:144
const std::shared_ptr< FourMomenta > fourMomenta() const
Definition: Model.h:114
void setParameterFlagsToUnchanged()
Set VariableStatus'es of all Parameter's to unchanged, or leave as fixed.
Definition: Model.cxx:464
const double sum_of_log_intensity(const Model &M, DataPartition &D, double ped=0)
Definition: Model.cxx:99
ParticleCombinationCache & particleCombinationCache()
Definition: Model.h:136
virtual bool consistent() const
Check consistency of object.
Definition: Model.cxx:139
Model(std::unique_ptr< SpinAmplitudeCache > SAC, bool include_phsp_factors=true)
Definition: Model.cxx:36
const ParticleCombinationCache & particleCombinationCache() const
Definition: Model.h:140
DataSet createDataSet(size_t n=0)
Definition: Model.cxx:425
virtual void addParticleCombination(std::shared_ptr< ParticleCombination > pc)
Definition: Model.cxx:156
container::value_type lone_elt(container &C)
Definition: Filter.h:39
const InitialStateParticleMap::value_type & addInitialStateParticle(std::shared_ptr< DecayingParticle > bg)
Definition: Model.cxx:229
const FinalStateParticleVector & finalStateParticles() const
Definition: Model.h:152
Four-Vector handling.
Definition: FourVector.h:41
FinalStateParticleVector FinalStateParticles_
vector of final state particles
Definition: Model.h:291
ParticleCombinationVector specialized to contain axes for defining a phase-space coordinate.
Definition: MassAxes.h:31
const std::shared_ptr< MeasuredBreakupMomenta > measuredBreakupMomenta() const
Definition: Model.h:122
ParticleCombinationCache ParticleCombinationCache_
ParticleCombination cache.
Definition: Model.h:271
void printFlags(const StatusManager &sm) const
Print all VariableStatus'es and CalculationStatus'es.
Definition: Model.cxx:502
ParticleSet particles(DecayingParticle &dp)
Definition: DecayingParticle.cxx:357
void setIncludePhaseSpaceFactors(bool b)
turn on/off inclusion of phase-space factors
Definition: Model.h:199
std::shared_ptr< MeasuredBreakupMomenta > measuredBreakupMomenta()
Definition: Model.h:118
Class for holding data and cached values per data point for fast calculation.
Definition: DataPoint.h:35
Base class for all data accessors that will only write to DataPoint once at initial data loading...
Definition: StaticDataAccessor.h:36
const bool includePhaseSpaceFactors() const
Definition: Model.h:102
ParticleSet::value_type particle(const Model &M, Last p, UnaryPredicates...P)
Definition: Model.h:357
bool Locked_
Definition: Model.h:262
void setCoordinateSystem(const CoordinateSystem< double, 3 > &cs)
set coordinate system
Definition: Model.cxx:292
std::shared_ptr< HelicityAngles > helicityAngles()
Definition: Model.h:127
RecalculableDataAccessorSet RecalculableDataAccessors_
set of pointers to RecalculableDataAccessors
Definition: Model.h:284
void requireMeasuredBreakupMomenta()
tell model to calculate and store measured breakup momenta
Definition: Model.cxx:308
std::unique_ptr< SpinAmplitudeCache > SpinAmplitudeCache_
SpinAmplitude cache.
Definition: Model.h:274
void requireHelicityAngles()
tell model to calculate and store helicity angles
Definition: Model.cxx:301
std::shared_ptr< HelicityAngles > HelicityAngles_
helicity angles manager
Definition: Model.h:300
const SpinAmplitudeCache * spinAmplitudeCache() const
Definition: Model.h:148
DataAccessorSet DataAccessors_
Set of all DataAccessor's registered to this model.
Definition: Model.h:277
const DataAccessorSet & dataAccessors() const
Definition: Model.h:160
Class defining a partition of the DataSet.
Definition: DataPartition.h:158
Caches SpinAmplitudes.
Definition: SpinAmplitudeCache.h:39
Definition: StatusManager.h:32
CoordinateSystem< double, 3 > CoordinateSystem_
Lab coordinate system to use in calculating helicity angles.
Definition: Model.h:265
Class implementing a PWA model.
Definition: Model.h:56
std::shared_ptr< FourMomenta > fourMomenta()
Definition: Model.h:110
const double intensity(const DecayTreeVector &dtv, const DataPoint &d)
Definition: DecayTree.h:167
Model & operator=(const Model &)=delete
std::shared_ptr< FourMomenta > FourMomenta_
four momenta manager
Definition: Model.h:294
void setFinalState(Types...FSPs)
Definition: Model.h:181
FreeAmplitudeSet free_amplitudes(const DecayingParticle &dp)
Definition: DecayingParticle.cxx:377
const CoordinateSystem< double, 3 > & coordinateSystem() const
Definition: Model.h:106
const bool decays_to_full_final_state(const Particle &p)
Definition: Particle.cxx:45
const std::shared_ptr< HelicityAngles > helicityAngles() const
Definition: Model.h:132
Caches list of ParticleCombination's.
Definition: ParticleCombinationCache.h:36
const InitialStateParticleMap & initialStateParticles() const
Definition: Model.h:156
Abstract base class for all objects accessing DataPoint's.
Definition: DataAccessor.h:38
FreeAmplitudeSet::value_type free_amplitude(const DecayingParticle &dp, Last p, UnaryPredicates...P)
Definition: DecayingParticle.h:208
void calculate(DataPartition &D) const
Definition: Model.cxx:51
void setFinalState(const FinalStateParticleVector &FSP)
Definition: Model.cxx:180
std::set< std::shared_ptr< T > > filter(const std::set< std::shared_ptr< T > > &S, Last p, UnaryPredicates...P)
Definition: Filter.h:60
std::string to_string(const CachedValue::Status &S)
streaming operator for CachedValue::Status
Definition: CachedValue.cxx:27
bool IncludePhaseSpaceFactors_
whether to include phase-space factors
Definition: Model.h:268
void setFinalStateMomenta(DataPoint &d, const std::vector< FourVector< double > > &P, StatusManager &sm) const
Definition: Model.cxx:208
const MassAxes massAxes(std::vector< std::vector< unsigned > > pcs={})
Definition: Model.cxx:357
Definition: DecayingParticle.h:51
StaticDataAccessorVector StaticDataAccessors_
Definition: Model.h:281
Abstract Particle base class.
Definition: Particle.h:47
std::vector< std::shared_ptr< DecayingParticle > > full_final_state_isp(const Model &M)
Definition: Model.cxx:276
bool locked() const
Definition: Model.h:90
void lock()
Definition: Model.cxx:315
InitialStateParticleMap InitialStateParticles_
Definition: Model.h:288
std::shared_ptr< MeasuredBreakupMomenta > MeasuredBreakupMomenta_
Breakup momenta manager.
Definition: Model.h:297
void printDataAccessors(bool printParticleCombinations=true) const
Print the list of DataAccessor's.
Definition: Model.cxx:471
Definition: RecalculableDataAccessor.h:44