laura is hosted by Hepforge, IPPP Durham
Laura++  v3r1
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauResonanceMaker.cc
Go to the documentation of this file.
1 
2 // Copyright University of Warwick 2004 - 2014.
3 // Distributed under the Boost Software License, Version 1.0.
4 // (See accompanying file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
5 
6 // Authors:
7 // Thomas Latham
8 // John Back
9 // Paul Harrison
10 
15 #include <iostream>
16 
17 #include "LauAbsResonance.hh"
18 #include "LauBelleNR.hh"
19 #include "LauBelleSymNR.hh"
20 #include "LauBreitWignerRes.hh"
21 #include "LauDabbaRes.hh"
22 #include "LauDaughters.hh"
23 #include "LauEFKLLMRes.hh"
24 #include "LauFlatteRes.hh"
25 #include "LauFlatNR.hh"
26 #include "LauGaussIncohRes.hh"
27 #include "LauGounarisSakuraiRes.hh"
28 #include "LauKappaRes.hh"
29 #include "LauLASSRes.hh"
30 #include "LauLASSBWRes.hh"
31 #include "LauLASSNRRes.hh"
34 #include "LauNRAmplitude.hh"
35 #include "LauPolNR.hh"
36 #include "LauRelBreitWignerRes.hh"
37 #include "LauResonanceInfo.hh"
38 #include "LauResonanceMaker.hh"
39 #include "LauRhoOmegaMix.hh"
40 #include "LauSigmaRes.hh"
41 
43 
44 
46 
47 
49  nResDefMax_(0)
50 {
51  this->createResonanceVector();
53 }
54 
56 {
57  for ( std::vector<LauBlattWeisskopfFactor*>::iterator iter = bwIndepFactors_.begin(); iter != bwIndepFactors_.end(); ++iter ) {
58  delete *iter;
59  }
60  bwIndepFactors_.clear();
61  for ( BWFactorCategoryMap::iterator iter = bwFactors_.begin(); iter != bwFactors_.end(); ++iter ) {
62  delete iter->second;
63  }
64  bwFactors_.clear();
65 }
66 
68 {
69  if ( resonanceMaker_ == 0 ) {
71  }
72 
73  return *resonanceMaker_;
74 }
75 
77 {
78  // Function to create all possible resonances that this class supports.
79  // Also add in the sigma and kappa - but a special paramterisation is used
80  // instead of the PDG "pole mass and width" values.
81 
82  std::cout << "INFO in LauResonanceMaker::createResonanceVector : Setting up possible resonance states..." << std::endl;
83 
84  LauResonanceInfo* neutral(0);
85  LauResonanceInfo* positve(0);
86  LauResonanceInfo* negatve(0);
87 
88  // Define the resonance names and store them in the array
89  resInfo_.clear();
90  resInfo_.reserve(100);
91  // rho resonances name, mass, width, spin, charge, default BW category, BW radius parameter (defaults to 4.0)
92  // rho(770)
93  neutral = new LauResonanceInfo("rho0(770)", 0.77526, 0.1478, 1, 0, LauBlattWeisskopfFactor::Light, 5.3);
94  positve = new LauResonanceInfo("rho+(770)", 0.77511, 0.1491, 1, 1, LauBlattWeisskopfFactor::Light, 5.3);
95  negatve = positve->createChargeConjugate();
96  resInfo_.push_back( neutral );
97  resInfo_.push_back( positve );
98  resInfo_.push_back( negatve );
99  // The following two lines of code are placed here in order to allow the following, rather niche, scenario:
100  // The LauRhoOmegaMix code permits (through the use of the optional independentPar argument of LauResonanceInfo::addExtraParameter) the magnitude and phase of the rho/omega mixing to potentially differ between the decay of the parent particle to rho0 X and the parent antiparticle to rho0 Xbar.
101  // This can be acheived by using the rho0(770) record in one case and the rho0(770)_COPY record in the other.
102  neutral = neutral->createSharedParameterRecord("rho0(770)_COPY");
103  resInfo_.push_back( neutral );
104  // rho(1450)
105  neutral = new LauResonanceInfo("rho0(1450)", 1.465, 0.400, 1, 0, LauBlattWeisskopfFactor::Light );
106  positve = new LauResonanceInfo("rho+(1450)", 1.465, 0.400, 1, 1, LauBlattWeisskopfFactor::Light );
107  negatve = positve->createChargeConjugate();
108  resInfo_.push_back( neutral );
109  resInfo_.push_back( positve );
110  resInfo_.push_back( negatve );
111  // rho(1700)
112  neutral = new LauResonanceInfo("rho0(1700)", 1.720, 0.250, 1, 0, LauBlattWeisskopfFactor::Light );
113  positve = new LauResonanceInfo("rho+(1700)", 1.720, 0.250, 1, 1, LauBlattWeisskopfFactor::Light );
114  negatve = positve->createChargeConjugate();
115  resInfo_.push_back( neutral );
116  resInfo_.push_back( positve );
117  resInfo_.push_back( negatve );
118 
119  // K* resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
120  // K*(892)
121  neutral = new LauResonanceInfo("K*0(892)", 0.89581, 0.0474, 1, 0, LauBlattWeisskopfFactor::Kstar, 3.0);
122  positve = new LauResonanceInfo("K*+(892)", 0.89166, 0.0508, 1, 1, LauBlattWeisskopfFactor::Kstar, 3.0);
123  negatve = positve->createChargeConjugate();
124  resInfo_.push_back( neutral );
125  resInfo_.push_back( positve );
126  resInfo_.push_back( negatve );
127  // K*(1410)
128  neutral = new LauResonanceInfo("K*0(1410)", 1.414, 0.232, 1, 0, LauBlattWeisskopfFactor::Kstar );
129  positve = new LauResonanceInfo("K*+(1410)", 1.414, 0.232, 1, 1, LauBlattWeisskopfFactor::Kstar );
130  negatve = positve->createChargeConjugate();
131  resInfo_.push_back( neutral );
132  resInfo_.push_back( positve );
133  resInfo_.push_back( negatve );
134  // K*_0(1430)
135  neutral = new LauResonanceInfo("K*0_0(1430)", 1.425, 0.270, 0, 0, LauBlattWeisskopfFactor::Kstar );
136  positve = new LauResonanceInfo("K*+_0(1430)", 1.425, 0.270, 0, 1, LauBlattWeisskopfFactor::Kstar );
137  negatve = positve->createChargeConjugate();
138  resInfo_.push_back( neutral );
139  resInfo_.push_back( positve );
140  resInfo_.push_back( negatve );
141  // LASS nonresonant model
142  neutral = neutral->createSharedParameterRecord("LASSNR0");
143  positve = positve->createSharedParameterRecord("LASSNR+");
144  negatve = negatve->createSharedParameterRecord("LASSNR-");
145  resInfo_.push_back( neutral );
146  resInfo_.push_back( positve );
147  resInfo_.push_back( negatve );
148  // K*_2(1430)
149  neutral = new LauResonanceInfo("K*0_2(1430)", 1.4324, 0.109, 2, 0, LauBlattWeisskopfFactor::Kstar );
150  positve = new LauResonanceInfo("K*+_2(1430)", 1.4256, 0.0985, 2, 1, LauBlattWeisskopfFactor::Kstar );
151  negatve = positve->createChargeConjugate();
152  resInfo_.push_back( neutral );
153  resInfo_.push_back( positve );
154  resInfo_.push_back( negatve );
155  // K*(1680)
156  neutral = new LauResonanceInfo("K*0(1680)", 1.717, 0.322, 1, 0, LauBlattWeisskopfFactor::Kstar );
157  positve = new LauResonanceInfo("K*+(1680)", 1.717, 0.322, 1, 1, LauBlattWeisskopfFactor::Kstar );
158  negatve = positve->createChargeConjugate();
159  resInfo_.push_back( neutral );
160  resInfo_.push_back( positve );
161  resInfo_.push_back( negatve );
162  // K*(1950)
163  neutral = new LauResonanceInfo("K*0_0(1950)", 1.945, 0.201, 0, 0, LauBlattWeisskopfFactor::Kstar );
164  positve = new LauResonanceInfo("K*+_0(1950)", 1.945, 0.201, 0, 1, LauBlattWeisskopfFactor::Kstar );
165  negatve = positve->createChargeConjugate();
166  resInfo_.push_back( neutral );
167  resInfo_.push_back( positve );
168  resInfo_.push_back( negatve );
169 
170  // phi resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
171  // phi(1020)
172  neutral = new LauResonanceInfo("phi(1020)", 1.019461, 0.004266, 1, 0, LauBlattWeisskopfFactor::Light );
173  resInfo_.push_back( neutral );
174  // phi(1680)
175  neutral = new LauResonanceInfo("phi(1680)", 1.680, 0.150, 1, 0, LauBlattWeisskopfFactor::Light );
176  resInfo_.push_back( neutral );
177 
178  // f resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
179  // f_0(980)
180  neutral = new LauResonanceInfo("f_0(980)", 0.990, 0.070, 0, 0, LauBlattWeisskopfFactor::Light );
181  resInfo_.push_back( neutral );
182  // f_2(1270)
183  neutral = new LauResonanceInfo("f_2(1270)", 1.2751, 0.1851, 2, 0, LauBlattWeisskopfFactor::Light );
184  resInfo_.push_back( neutral );
185  // f_0(1370)
186  neutral = new LauResonanceInfo("f_0(1370)", 1.370, 0.350, 0, 0, LauBlattWeisskopfFactor::Light );
187  resInfo_.push_back( neutral );
188  // f'_0(1300), from Belle's Kspipi paper
189  neutral = new LauResonanceInfo("f'_0(1300)", 1.449, 0.126, 0, 0, LauBlattWeisskopfFactor::Light );
190  resInfo_.push_back( neutral );
191  // f_0(1500)
192  neutral = new LauResonanceInfo("f_0(1500)", 1.505, 0.109, 0, 0, LauBlattWeisskopfFactor::Light );
193  resInfo_.push_back( neutral );
194  // f'_2(1525)
195  neutral = new LauResonanceInfo("f'_2(1525)", 1.525, 0.073, 2, 0, LauBlattWeisskopfFactor::Light );
196  resInfo_.push_back( neutral );
197  // f_0(1710)
198  neutral = new LauResonanceInfo("f_0(1710)", 1.722, 0.135, 0, 0, LauBlattWeisskopfFactor::Light );
199  resInfo_.push_back( neutral );
200  // f_2(2010)
201  neutral = new LauResonanceInfo("f_2(2010)", 2.011, 0.202, 2, 0, LauBlattWeisskopfFactor::Light );
202  resInfo_.push_back( neutral );
203 
204  // omega resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
205  // omega(782)
206  neutral = new LauResonanceInfo("omega(782)", 0.78265, 0.00849, 1, 0, LauBlattWeisskopfFactor::Light );
207  resInfo_.push_back( neutral );
208 
209  // a resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
210  // a_0(980)
211  neutral = new LauResonanceInfo("a0_0(980)", 0.980, 0.092, 0, 0, LauBlattWeisskopfFactor::Light );
212  positve = new LauResonanceInfo("a+_0(980)", 0.980, 0.092, 0, 1, LauBlattWeisskopfFactor::Light );
213  negatve = positve->createChargeConjugate();
214  resInfo_.push_back( neutral );
215  resInfo_.push_back( positve );
216  resInfo_.push_back( negatve );
217  // a_0(1450)
218  neutral = new LauResonanceInfo("a0_0(1450)", 1.474, 0.265, 0, 0, LauBlattWeisskopfFactor::Light );
219  positve = new LauResonanceInfo("a+_0(1450)", 1.474, 0.265, 0, 1, LauBlattWeisskopfFactor::Light );
220  negatve = positve->createChargeConjugate();
221  resInfo_.push_back( neutral );
222  resInfo_.push_back( positve );
223  resInfo_.push_back( negatve );
224  // a_2(1320)
225  neutral = new LauResonanceInfo("a0_2(1320)", 1.3190, 0.1050, 2, 0, LauBlattWeisskopfFactor::Light );
226  positve = new LauResonanceInfo("a+_2(1320)", 1.3190, 0.1050, 2, 1, LauBlattWeisskopfFactor::Light );
227  negatve = positve->createChargeConjugate();
228  resInfo_.push_back( neutral );
229  resInfo_.push_back( positve );
230  resInfo_.push_back( negatve );
231 
232  // charmonium resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
233  // chi_c0
234  neutral = new LauResonanceInfo("chi_c0", 3.41475, 0.0105, 0, 0, LauBlattWeisskopfFactor::Charmonium );
235  resInfo_.push_back( neutral );
236  // chi_c1
237  neutral = new LauResonanceInfo("chi_c1", 3.51066, 0.00084, 0, 0, LauBlattWeisskopfFactor::Charmonium );
238  resInfo_.push_back( neutral );
239  // chi_c2
240  neutral = new LauResonanceInfo("chi_c2", 3.55620, 0.00193, 2, 0, LauBlattWeisskopfFactor::Charmonium );
241  resInfo_.push_back( neutral );
242  // X(3872)
243  neutral = new LauResonanceInfo("X(3872)", 3.87169, 0.0012, 1, 0, LauBlattWeisskopfFactor::Charmonium );
244  resInfo_.push_back( neutral );
245 
246  // unknown scalars name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
247  // sigma
248  neutral = new LauResonanceInfo("sigma0", 0.475, 0.550, 0, 0, LauBlattWeisskopfFactor::Light );
249  positve = new LauResonanceInfo("sigma+", 0.475, 0.550, 0, 1, LauBlattWeisskopfFactor::Light );
250  negatve = positve->createChargeConjugate();
251  resInfo_.push_back( neutral );
252  resInfo_.push_back( positve );
253  resInfo_.push_back( negatve );
254  // kappa
255  neutral = new LauResonanceInfo("kappa0", 0.682, 0.547, 0, 0, LauBlattWeisskopfFactor::Kstar );
256  positve = new LauResonanceInfo("kappa+", 0.682, 0.547, 0, 1, LauBlattWeisskopfFactor::Kstar );
257  negatve = positve->createChargeConjugate();
258  resInfo_.push_back( neutral );
259  resInfo_.push_back( positve );
260  resInfo_.push_back( negatve );
261  // dabba
262  neutral = new LauResonanceInfo("dabba0", 2.098, 0.520, 0, 0, LauBlattWeisskopfFactor::Charm );
263  positve = new LauResonanceInfo("dabba+", 2.098, 0.520, 0, 1, LauBlattWeisskopfFactor::Charm );
264  negatve = positve->createChargeConjugate();
265  resInfo_.push_back( neutral );
266  resInfo_.push_back( positve );
267  resInfo_.push_back( negatve );
268 
269  // excited charm states name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
270  // D*
271  neutral = new LauResonanceInfo("D*0", 2.00696, 0.0021, 1, 0, LauBlattWeisskopfFactor::Charm );
272  positve = new LauResonanceInfo("D*+", 2.01026, 83.4e-6, 1, 1, LauBlattWeisskopfFactor::Charm );
273  negatve = positve->createChargeConjugate();
274  resInfo_.push_back( neutral );
275  resInfo_.push_back( positve );
276  resInfo_.push_back( negatve );
277  // D*_0
278  neutral = new LauResonanceInfo("D*0_0", 2.318, 0.267, 0, 0, LauBlattWeisskopfFactor::Charm );
279  positve = new LauResonanceInfo("D*+_0", 2.403, 0.283, 0, 1, LauBlattWeisskopfFactor::Charm );
280  negatve = positve->createChargeConjugate();
281  resInfo_.push_back( neutral );
282  resInfo_.push_back( positve );
283  resInfo_.push_back( negatve );
284  // D*_2
285  //AVERAGE--neutral = new LauResonanceInfo("D*0_2", 2.4618, 0.049, 2, 0 );
286  neutral = new LauResonanceInfo("D*0_2", 2.4626, 0.049, 2, 0, LauBlattWeisskopfFactor::Charm );
287  positve = new LauResonanceInfo("D*+_2", 2.4643, 0.037, 2, 1, LauBlattWeisskopfFactor::Charm );
288  negatve = positve->createChargeConjugate();
289  resInfo_.push_back( neutral );
290  resInfo_.push_back( positve );
291  resInfo_.push_back( negatve );
292  // D1(2420)
293  neutral = new LauResonanceInfo("D0_1(2420)", 2.4214, 0.0274, 1, 0, LauBlattWeisskopfFactor::Charm );
294  positve = new LauResonanceInfo("D+_1(2420)", 2.4232, 0.025, 1, 1, LauBlattWeisskopfFactor::Charm );
295  negatve = positve->createChargeConjugate();
296  resInfo_.push_back( neutral );
297  resInfo_.push_back( positve );
298  resInfo_.push_back( negatve );
299  // D(2600)
300  //OLD--neutral = new LauResonanceInfo("D0(2600)", 2.6087, 0.093, 0, 0 );
301  //OLD--positve = new LauResonanceInfo("D+(2600)", 2.6213, 0.093, 0, 1 );
302  neutral = new LauResonanceInfo("D0(2600)", 2.612, 0.093, 0, 0, LauBlattWeisskopfFactor::Charm );
303  positve = new LauResonanceInfo("D+(2600)", 2.612, 0.093, 0, 1, LauBlattWeisskopfFactor::Charm );
304  negatve = positve->createChargeConjugate();
305  resInfo_.push_back( neutral );
306  resInfo_.push_back( positve );
307  resInfo_.push_back( negatve );
308  // D(2760)
309  //OLD-- neutral = new LauResonanceInfo("D0(2760)", 2.7633, 0.061, 1, 0 );
310  //OLD-- positve = new LauResonanceInfo("D+(2760)", 2.7697, 0.061, 1, 1 );
311  neutral = new LauResonanceInfo("D0(2760)", 2.761, 0.063, 1, 0, LauBlattWeisskopfFactor::Charm );
312  positve = new LauResonanceInfo("D+(2760)", 2.761, 0.063, 1, 1, LauBlattWeisskopfFactor::Charm );
313  negatve = positve->createChargeConjugate();
314  resInfo_.push_back( neutral );
315  resInfo_.push_back( positve );
316  resInfo_.push_back( negatve );
317  // D(2900)
318  neutral = new LauResonanceInfo("D0(3000)", 3.0, 0.15, 0, 0, LauBlattWeisskopfFactor::Charm );
319  resInfo_.push_back( neutral );
320  // D(3400)
321  neutral = new LauResonanceInfo("D0(3400)", 3.4, 0.15, 0, 0, LauBlattWeisskopfFactor::Charm );
322  resInfo_.push_back( neutral );
323 
324  // excited strange charm name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
325  // Ds*
326  positve = new LauResonanceInfo("Ds*+", 2.1121, 0.0019, 1, 1, LauBlattWeisskopfFactor::StrangeCharm );
327  negatve = positve->createChargeConjugate();
328  resInfo_.push_back( positve );
329  resInfo_.push_back( negatve );
330  // Ds0*(2317)
331  positve = new LauResonanceInfo("Ds*+_0(2317)", 2.3177, 0.0038, 0, 1, LauBlattWeisskopfFactor::StrangeCharm );
332  negatve = positve->createChargeConjugate();
333  resInfo_.push_back( positve );
334  resInfo_.push_back( negatve );
335  // Ds2*(2573)
336  positve = new LauResonanceInfo("Ds*+_2(2573)", 2.5719, 0.017, 2, 1, LauBlattWeisskopfFactor::StrangeCharm );
337  negatve = positve->createChargeConjugate();
338  resInfo_.push_back( positve );
339  resInfo_.push_back( negatve );
340  // Ds1*(2700)
341  positve = new LauResonanceInfo("Ds*+_1(2700)", 2.709, 0.117, 1, 1, LauBlattWeisskopfFactor::StrangeCharm );
342  negatve = positve->createChargeConjugate();
343  resInfo_.push_back( positve );
344  resInfo_.push_back( negatve );
345  // Ds1*(2860)
346  positve = new LauResonanceInfo("Ds*+_1(2860)", 2.862, 0.180, 1, 1, LauBlattWeisskopfFactor::StrangeCharm );
347  negatve = positve->createChargeConjugate();
348  resInfo_.push_back( positve );
349  resInfo_.push_back( negatve );
350  // Ds3*(2860)
351  positve = new LauResonanceInfo("Ds*+_3(2860)", 2.862, 0.058, 3, 1, LauBlattWeisskopfFactor::StrangeCharm );
352  negatve = positve->createChargeConjugate();
353  resInfo_.push_back( positve );
354  resInfo_.push_back( negatve );
355 
356  // excited bottom states name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
357  // B*
358  neutral = new LauResonanceInfo("B*0", 5.3252, 0.00, 1, 0, LauBlattWeisskopfFactor::Beauty, 6.0);
359  positve = new LauResonanceInfo("B*+", 5.3252, 0.00, 1, 1, LauBlattWeisskopfFactor::Beauty, 6.0);
360  negatve = positve->createChargeConjugate();
361  resInfo_.push_back( neutral );
362  resInfo_.push_back( positve );
363  resInfo_.push_back( negatve );
364 
365  // excited strange bottom name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
366  // Bs*
367  neutral = new LauResonanceInfo("Bs*0", 5.4154, 0.00, 1, 0, LauBlattWeisskopfFactor::StrangeBeauty, 6.0);
368  resInfo_.push_back( neutral );
369 
370  // nonresonant models name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
371  // Phase-space nonresonant model
372  neutral = new LauResonanceInfo("NonReson", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
373  resInfo_.push_back( neutral );
374  // Theory-based nonresonant model
375  neutral = new LauResonanceInfo("NRModel", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
376  resInfo_.push_back( neutral );
377  // Belle nonresonant models
378  neutral = new LauResonanceInfo("BelleSymNR", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
379  resInfo_.push_back( neutral );
380  neutral = new LauResonanceInfo("BelleNR", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
381  resInfo_.push_back( neutral );
382  positve = new LauResonanceInfo("BelleNR+", 0.0, 0.0, 0, 1, LauBlattWeisskopfFactor::Light );
383  negatve = positve->createChargeConjugate();
384  resInfo_.push_back( positve );
385  resInfo_.push_back( negatve );
386  neutral = new LauResonanceInfo("BelleNR_Swave", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
387  resInfo_.push_back( neutral );
388  positve = new LauResonanceInfo("BelleNR_Swave+",0.0, 0.0, 0, 1, LauBlattWeisskopfFactor::Light );
389  negatve = positve->createChargeConjugate();
390  resInfo_.push_back( positve );
391  resInfo_.push_back( negatve );
392  neutral = new LauResonanceInfo("BelleNR_Pwave", 0.0, 0.0, 1, 0, LauBlattWeisskopfFactor::Light );
393  resInfo_.push_back( neutral );
394  positve = new LauResonanceInfo("BelleNR_Pwave+",0.0, 0.0, 1, 1, LauBlattWeisskopfFactor::Light );
395  negatve = positve->createChargeConjugate();
396  resInfo_.push_back( positve );
397  resInfo_.push_back( negatve );
398  neutral = new LauResonanceInfo("BelleNR_Dwave", 0.0, 0.0, 2, 0, LauBlattWeisskopfFactor::Light );
399  resInfo_.push_back( neutral );
400  positve = new LauResonanceInfo("BelleNR_Dwave+",0.0, 0.0, 2, 1, LauBlattWeisskopfFactor::Light );
401  negatve = positve->createChargeConjugate();
402  resInfo_.push_back( positve );
403  resInfo_.push_back( negatve );
404  // Taylor expansion nonresonant model
405  neutral = new LauResonanceInfo("NRTaylor", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
406  resInfo_.push_back( neutral );
407  // Polynomial nonresonant models
408  neutral = new LauResonanceInfo("PolNR_S0", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
409  resInfo_.push_back( neutral );
410  neutral = new LauResonanceInfo("PolNR_S1", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
411  resInfo_.push_back( neutral );
412  neutral = new LauResonanceInfo("PolNR_S2", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
413  resInfo_.push_back( neutral );
414  neutral = new LauResonanceInfo("PolNR_P0", 0.0, 0.0, 1, 0, LauBlattWeisskopfFactor::Light );
415  resInfo_.push_back( neutral );
416  neutral = new LauResonanceInfo("PolNR_P1", 0.0, 0.0, 1, 0, LauBlattWeisskopfFactor::Light );
417  resInfo_.push_back( neutral );
418  neutral = new LauResonanceInfo("PolNR_P2", 0.0, 0.0, 1, 0, LauBlattWeisskopfFactor::Light );
419  resInfo_.push_back( neutral );
420 
421  nResDefMax_ = resInfo_.size();
422 }
423 
425 {
426  if ( bwCategory == LauBlattWeisskopfFactor::Default || bwCategory == LauBlattWeisskopfFactor::Indep ) {
427  std::cerr << "WARNING in LauResonanceMaker::setDefaultBWRadius : cannot set radius values for Default or Indep categories" << std::endl;
428  return;
429  }
430 
431  // Check if a LauBlattWeisskopfFactor object has been created for this category
432  // If it has then we can set it directly
433  // If not then we just store the value to be used when it is created later
434  BWFactorCategoryMap::iterator factor_iter = bwFactors_.find( bwCategory );
435  if ( factor_iter != bwFactors_.end() ) {
436  LauBlattWeisskopfFactor* bwFactor = factor_iter->second;
437  LauParameter* radius = bwFactor->getRadiusParameter();
438  radius->value(bwRadius);
439  radius->initValue(bwRadius);
440  radius->genValue(bwRadius);
441  }
442 
443  bwDefaultRadii_[bwCategory] = bwRadius;
444 }
445 
447 {
448  if ( bwCategory == LauBlattWeisskopfFactor::Default || bwCategory == LauBlattWeisskopfFactor::Indep ) {
449  std::cerr << "WARNING in LauResonanceMaker::fixBWRadius : cannot fix/float radius values for Default or Indep categories" << std::endl;
450  return;
451  }
452 
453  // Check if a LauBlattWeisskopfFactor object has been created for this category
454  // If it has then we can set it directly
455  // If not then we just store the value to be used when it is created later
456  BWFactorCategoryMap::iterator factor_iter = bwFactors_.find( bwCategory );
457  if ( factor_iter != bwFactors_.end() ) {
458  LauBlattWeisskopfFactor* bwFactor = factor_iter->second;
459  LauParameter* radius = bwFactor->getRadiusParameter();
460  radius->fixed(fixRadius);
461  }
462 
463  bwFixRadii_[bwCategory] = fixRadius;
464 }
465 
467 {
468  LauBlattWeisskopfFactor* bwFactor(0);
469 
470  BWFactorCategoryMap::iterator factor_iter = bwFactors_.find( bwCategory );
471  BWRadiusCategoryMap::const_iterator radius_iter = bwDefaultRadii_.find( bwCategory );
472 
473  if ( bwCategory == LauBlattWeisskopfFactor::Indep ) {
474  bwFactor = new LauBlattWeisskopfFactor( *resInfo, bwType, bwCategory );
475  bwIndepFactors_.push_back(bwFactor);
476  } else if ( factor_iter == bwFactors_.end() ) {
477  if ( radius_iter == bwDefaultRadii_.end() ) {
478  bwFactor = new LauBlattWeisskopfFactor( *resInfo, bwType, bwCategory );
479  } else {
480  bwFactor = new LauBlattWeisskopfFactor( *resInfo, radius_iter->second, bwType, bwCategory );
481  }
482  bwFactors_[bwCategory] = bwFactor;
483  } else {
484  const UInt_t resSpin = resInfo->getSpin();
485  bwFactor = factor_iter->second->createClone( resSpin );
486  if ( bwFactor->getBarrierType() != bwType ) {
487  std::cerr << "WARNING in LauResonanceMaker::getBWFactor : A barrier factor already exists for the specified category but does not have the specified barrier type.\n";
488  std::cerr << " : If you want to use that type you will need to use a different category for this resonance." << std::endl;
489  }
490  }
491 
492  BWRadiusFixedCategoryMap::const_iterator radius_fixed_iter = bwFixRadii_.find( bwCategory );
493 
494  if ( radius_fixed_iter != bwFixRadii_.end() ) {
495  const Bool_t fixed = radius_fixed_iter->second;
496 
497  LauParameter* radius = bwFactor->getRadiusParameter();
498  radius->fixed( fixed );
499  }
500 
501  return bwFactor;
502 }
503 
504 LauAbsResonance* LauResonanceMaker::getResonance(const LauDaughters* daughters, const TString& resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType, const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory, const LauBlattWeisskopfFactor::BarrierType bwType)
505 {
506  // Routine to return the appropriate LauAbsResonance object given the resonance
507  // name (resName), which daughter is the bachelor track (resPairAmpInt = 1,2 or 3),
508  // and the resonance type ("BW" = Breit-Wigner, "Flatte" = Flatte distribution).
509 
510  // Loop over all possible resonance states we have defined in
511  // createResonanceVector() until we get a match with the name of the resonance
512 
513  LauResonanceInfo* resInfo(0);
514  for (std::vector<LauResonanceInfo*>::const_iterator iter=resInfo_.begin(); iter!=resInfo_.end(); ++iter) {
515 
516  if (resName == (*iter)->getName()) {
517  // We have recognised the resonance name.
518  std::cout<<"INFO in LauResonanceMaker::getResonance : Creating resonance: "<<resName<<std::endl;
519 
520  resInfo = (*iter);
521 
522  // stop looping
523  break;
524  }
525  }
526 
527  // if we couldn't find the right resonance then we should return a null pointer
528  if ( resInfo == 0 ) {
529  std::cout<<"ERROR in LauResonanceMaker::getResonance : Unable to locate resonance info for: "<<resName<<std::endl;
530  return 0;
531  }
532 
533  LauAbsResonance* theResonance(0);
534 
535  // Now construct the resonnace using the right type.
536  // If we don't recognise the resonance model name, just use a simple Breit-Wigner.
537 
538  switch ( resType ) {
539 
540  case LauAbsResonance::BW :
541  // Simple non-relativistic Breit-Wigner
542  std::cout<<" : Using simple Breit-Wigner lineshape. "<<std::endl;
543  theResonance = new LauBreitWignerRes(resInfo, resPairAmpInt, daughters);
544  break;
545 
547  {
548  // Relativistic Breit-Wigner with Blatt-Weisskopf factors.
549  std::cout<<" : Using relativistic Breit-Wigner lineshape. "<<std::endl;
550  theResonance = new LauRelBreitWignerRes(resInfo, resPairAmpInt, daughters);
552  LauBlattWeisskopfFactor::BlattWeisskopfCategory resCategory = bwCategory;
553  if ( bwCategory == LauBlattWeisskopfFactor::Default ) {
554  resCategory = resInfo->getBWCategory();
555  }
556  LauBlattWeisskopfFactor* resBWFactor = this->getBWFactor( resCategory, resInfo, bwType );
557  LauBlattWeisskopfFactor* parBWFactor = this->getBWFactor( parCategory, resInfo, bwType );
558  theResonance->setBarrierRadii( resBWFactor, parBWFactor );
559  break;
560  }
561 
562  case LauAbsResonance::GS :
563  {
564  // Gounaris-Sakurai function to try and model the rho(770) better
565  std::cout<<" : Using Gounaris-Sakurai lineshape. "<<std::endl;
566  theResonance = new LauGounarisSakuraiRes(resInfo, resPairAmpInt, daughters);
568  LauBlattWeisskopfFactor::BlattWeisskopfCategory resCategory = bwCategory;
569  if ( bwCategory == LauBlattWeisskopfFactor::Default ) {
570  resCategory = resInfo->getBWCategory();
571  }
572  LauBlattWeisskopfFactor* resBWFactor = this->getBWFactor( resCategory, resInfo, bwType );
573  LauBlattWeisskopfFactor* parBWFactor = this->getBWFactor( parCategory, resInfo, bwType );
574  theResonance->setBarrierRadii( resBWFactor, parBWFactor );
575  break;
576  }
577 
579  // Flatte distribution - coupled channel Breit-Wigner
580  std::cout<<" : Using Flatte lineshape. "<<std::endl;
581  theResonance = new LauFlatteRes(resInfo, resPairAmpInt, daughters);
582  break;
583 
585  // Sigma model - should only be used for the pi-pi system
586  std::cout<<" : Using Sigma lineshape. "<<std::endl;
587  theResonance = new LauSigmaRes(resInfo, resPairAmpInt, daughters);
588  break;
589 
591  // Kappa model - should only be used for the K-pi system
592  std::cout<<" : Using Kappa lineshape. "<<std::endl;
593  theResonance = new LauKappaRes(resInfo, resPairAmpInt, daughters);
594  break;
595 
597  // Dabba model - should only be used for the D-pi system
598  std::cout<<" : Using Dabba lineshape. "<<std::endl;
599  theResonance = new LauDabbaRes(resInfo, resPairAmpInt, daughters);
600  break;
601 
602  case LauAbsResonance::LASS :
603  // LASS function to try and model the K-pi S-wave better
604  std::cout<<" : Using LASS lineshape. "<<std::endl;
605  theResonance = new LauLASSRes(resInfo, resPairAmpInt, daughters);
606  break;
607 
609  // LASS function to try and model the K-pi S-wave better
610  std::cout<<" : Using LASS lineshape (resonant part only). "<<std::endl;
611  theResonance = new LauLASSBWRes(resInfo, resPairAmpInt, daughters);
612  break;
613 
615  // LASS function to try and model the K-pi S-wave better
616  std::cout<<" : Using LASS lineshape (nonresonant part only). "<<std::endl;
617  theResonance = new LauLASSNRRes(resInfo, resPairAmpInt, daughters);
618  break;
619 
621  // EFKLLM form-factor description of the K-pi S-wave
622  std::cout<<" : Using EFKLLM lineshape. "<<std::endl;
623  theResonance = new LauEFKLLMRes(resInfo, resPairAmpInt, daughters);
624  break;
625 
627  // K-matrix description of S-wave
628  std::cerr<<"ERROR in LauResonanceMaker::getResonance : K-matrix type specified, which should be separately handled."<<std::endl;
629  break;
630 
632  // uniform NR amplitude - arguments are there to preserve the interface
633  std::cout<<" : Using uniform NR lineshape. "<<std::endl;
634  // we override resPairAmpInt here
635  theResonance = new LauFlatNR(resInfo, 0, daughters);
636  break;
637 
639  // NR amplitude model - arguments are there to preserve the interface
640  std::cout<<" : Using NR-model lineshape. "<<std::endl;
641  // we override resPairAmpInt here
642  theResonance = new LauNRAmplitude(resInfo, 0, daughters);
643  break;
644 
647  // Belle NR amplitude model - arguments are there to preserve the interface
648  std::cout<<" : Using Belle NR lineshape. "<<std::endl;
649  theResonance = new LauBelleNR(resInfo, resType, resPairAmpInt, daughters);
650  break;
651 
655  // Belle NR amplitude model - arguments are there to preserve the interface
656  std::cout<<" : Using Belle symmetric NR lineshape. "<<std::endl;
657  theResonance = new LauBelleSymNR(resInfo, resType, resPairAmpInt, daughters);
658  break;
659 
661  // Polynomial NR amplitude model - arguments are there to preserve the interface
662  std::cout<<" : Using polynomial NR lineshape. "<<std::endl;
663  theResonance = new LauPolNR(resInfo, resPairAmpInt, daughters);
664  break;
665 
667  // Model independent partial wave
668  std::cout<<" : Using model independent partial wave lineshape (magnitude and phase). "<<std::endl;
669  theResonance = new LauModIndPartWaveMagPhase(resInfo, resPairAmpInt, daughters);
670  break;
671 
673  // Model independent partial wave
674  std::cout<<" : Using model independent partial wave lineshape (real and imaginary part). "<<std::endl;
675  theResonance = new LauModIndPartWaveRealImag(resInfo, resPairAmpInt, daughters);
676  break;
677 
679  // Incoherent Gaussian
680  std::cout<<" : Using incoherent Gaussian lineshape. "<<std::endl;
681  theResonance = new LauGaussIncohRes(resInfo, resPairAmpInt, daughters);
682  break;
683 
688  // Rho-omega mass mixing model
689  std::cout<<" : Using rho-omega mass mixing lineshape. "<<std::endl;
690  theResonance = new LauRhoOmegaMix(resInfo, resType, resPairAmpInt, daughters);
692  LauBlattWeisskopfFactor::BlattWeisskopfCategory resCategory = bwCategory;
693  if ( bwCategory == LauBlattWeisskopfFactor::Default ) {
694  resCategory = resInfo->getBWCategory();
695  }
696  LauBlattWeisskopfFactor* resBWFactor = this->getBWFactor( resCategory, resInfo, bwType );
697  LauBlattWeisskopfFactor* parBWFactor = this->getBWFactor( parCategory, resInfo, bwType );
698  theResonance->setBarrierRadii( resBWFactor, parBWFactor );
699  break;
700 
701  }
702 
703  return theResonance;
704 }
705 
706 Int_t LauResonanceMaker::resTypeInt(const TString& name) const
707 {
708  // Internal function that returns the resonance integer, specified by the
709  // order of the resonance vector defined in createResonanceVector(),
710  // for a given resonance name.
711  Int_t resTypInt(-99);
712  Int_t i(0);
713 
714  for (std::vector<LauResonanceInfo*>::const_iterator iter=resInfo_.begin(); iter!=resInfo_.end(); ++iter) {
715 
716  if (name.BeginsWith((*iter)->getName(), TString::kExact) == kTRUE) {
717  // We have recognised the resonance from those that are available
718  resTypInt = i;
719  return resTypInt;
720  }
721  ++i;
722  }
723 
724  return resTypInt;
725 }
726 
727 void LauResonanceMaker::printAll( std::ostream& stream ) const
728 {
729  for ( std::vector<LauResonanceInfo*>::const_iterator iter = resInfo_.begin(); iter != resInfo_.end(); ++iter ) {
730  stream << (**iter) << std::endl;
731  }
732 }
733 
734 LauResonanceInfo* LauResonanceMaker::getResInfo(const TString& resName) const
735 {
736  LauResonanceInfo* resInfo(0);
737  for (std::vector<LauResonanceInfo*>::const_iterator iter=resInfo_.begin(); iter!=resInfo_.end(); ++iter) {
738 
739  if (resName == (*iter)->getName()) {
740  // We have recognised the resonance name.
741  resInfo = (*iter);
742  // stop looping
743  break;
744  }
745  }
746 
747  return resInfo;
748 
749 }
File containing declaration of LauNRAmplitude class.
Bool_t fixed() const
Check whether the parameter is fixed or floated.
LauBlattWeisskopfFactor::BlattWeisskopfCategory getBWCategory() const
Retrieve the BW category of the resonant particle.
Class for defining the LASS resonance model.
Definition: LauLASSRes.hh:31
BWRadiusFixedCategoryMap bwFixRadii_
void setBarrierRadii(LauBlattWeisskopfFactor *resFactor, LauBlattWeisskopfFactor *parFactor)
Set the form factor model and parameters.
BWRadiusCategoryMap bwDefaultRadii_
File containing declaration of LauResonanceInfo class.
ClassImp(LauAbsCoeffSet)
BWFactorCategoryMap bwFactors_
Class for defining the properties of a resonant particle.
const TString & name() const
The parameter name.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:33
static LauResonanceMaker * resonanceMaker_
The singleton instance.
std::vector< LauResonanceInfo * > resInfo_
The known resonances.
File containing declaration of LauBelleNR class.
Class for defininf the Gounaris-Sakurai resonance model.
void fixBWRadius(const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory, const Bool_t fixRadius)
Fix or release the Blatt-Weisskopf barrier radius for the given category.
File containing declaration of LauDaughters class.
LauAbsResonance * getResonance(const LauDaughters *daughters, const TString &resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType, const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory=LauBlattWeisskopfFactor::Default, const LauBlattWeisskopfFactor::BarrierType bwType=LauBlattWeisskopfFactor::BWPrimeBarrier)
Create a resonance.
const LauParameter * getRadiusParameter() const
Retrieve the radius parameter.
File containing declaration of LauModIndPartWaveRealImag class.
File containing declaration of LauGounarisSakuraiRes class.
File containing declaration of LauBelleSymNR class.
File containing declaration of LauModIndPartWaveMagPhase class.
File containing declaration of LauBreitWignerRes class.
Int_t resTypeInt(const TString &name) const
Retrieve the integer index for the specified resonance.
File containing declaration of LauRelBreitWignerRes class.
File containing declaration of LauLASSBWRes class.
void setDefaultBWRadius(const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory, const Double_t bwRadius)
Set the BW radius for the given category.
Class for defining the EFKLLM K-pi S-wave model.
Definition: LauEFKLLMRes.hh:36
Class for defining the non resonant part of the LASS model.
Definition: LauLASSNRRes.hh:31
File containing declaration of LauEFKLLMRes class.
Class for defining the Sigma resonance model.
Definition: LauSigmaRes.hh:31
Singleton factory class for creating resonances.
void createResonanceVector()
Create the list of known resonances.
UInt_t getSpin() const
Retrieve the spin of the resonant particle.
File containing declaration of LauResonanceMaker class.
File containing declaration of LauRhoOmegaMix class.
File containing declaration of LauSigmaRes class.
Class for defining the Belle nonresonant model.
Definition: LauBelleNR.hh:34
Class for defining the NR amplitude model.
Class for defining the simple Breit-Wigner resonance model.
Class for defining a model independent partial wave component where the amplitudes are parameterised ...
Class for defining the terms of Babar nonresonant model.
Definition: LauPolNR.hh:33
Class for defining the Kappa resonance model.
Definition: LauKappaRes.hh:32
File containing declaration of LauDabbaRes class.
void printAll(std::ostream &stream) const
Print the information records, one per line, to the requested stream.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:35
Class for defining the relativistic Breit-Wigner resonance model.
Class for defining the rho-omega resonance mixing model.
std::vector< LauBlattWeisskopfFactor * > bwIndepFactors_
The Blatt-Weisskopf factor objects for resonances in the independent category.
BarrierType getBarrierType() const
Retrieve the barrier type.
File containing declaration of LauLASSRes class.
File containing declaration of LauLASSNRRes class.
virtual ~LauResonanceMaker()
Destructor.
LauResonanceModel
Define the allowed resonance types.
Class for defining the symmetric Belle Non Resonant model.
LauResonanceMaker()
Constructor.
BarrierType
Define the allowed types of barrier factors.
LauBlattWeisskopfFactor * getBWFactor(const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory, const LauResonanceInfo *resInfo, const LauBlattWeisskopfFactor::BarrierType bwType)
Retrieve Blatt-Weisskopf factor for the given category.
Class for defining an incoherent resonance with a Gaussian mass dependence.
Class for defining the resonant part of the LASS model.
Definition: LauLASSBWRes.hh:31
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc...
Class for defining the Dabba resonance model.
Definition: LauDabbaRes.hh:31
Double_t initValue() const
The initial value of the parameter.
File containing declaration of LauAbsResonance class.
Class that implements the Blatt-Weisskopf barrier factor.
LauBlattWeisskopfFactor * createClone(const UInt_t newSpin)
Method to create a new factor with cloned radius parameter.
Double_t value() const
The value of the parameter.
File containing declaration of LauPolNR class.
BlattWeisskopfCategory
Define resonance categories that will share common barrier factor radii.
UInt_t nResDefMax_
The number of known resonances.
Class for defining a model independent partial wave component where the amplitudes are parameterised ...
static LauResonanceMaker & get()
Get the factory instance.
Class for defining a uniform nonresonant amplitude.
Definition: LauFlatNR.hh:30
LauResonanceInfo * createChargeConjugate()
Create the charge conjugate particle info record.
LauResonanceInfo * createSharedParameterRecord(const TString &name)
Create another record that will share parameters with this one.
Class for defining the Flatte resonance model.
Definition: LauFlatteRes.hh:31
Double_t genValue() const
The value generated for the parameter.
File containing declaration of LauGaussIncohRes class.
LauResonanceInfo * getResInfo(const TString &resName) const
Get the information for the given resonance name.
File containing declaration of LauKappaRes class.
File containing declaration of LauFlatteRes class.
File containing declaration of LauFlatNR class.