laura is hosted by Hepforge, IPPP Durham
Laura++  v3r5
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauResonanceMaker.cc
Go to the documentation of this file.
1 
2 /*
3 Copyright 2004 University of Warwick
4 
5 Licensed under the Apache License, Version 2.0 (the "License");
6 you may not use this file except in compliance with the License.
7 You may obtain a copy of the License at
8 
9  http://www.apache.org/licenses/LICENSE-2.0
10 
11 Unless required by applicable law or agreed to in writing, software
12 distributed under the License is distributed on an "AS IS" BASIS,
13 WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 See the License for the specific language governing permissions and
15 limitations under the License.
16 */
17 
18 /*
19 Laura++ package authors:
20 John Back
21 Paul Harrison
22 Thomas Latham
23 */
24 
29 #include <iostream>
30 
31 #include "LauAbsResonance.hh"
32 #include "LauBelleNR.hh"
33 #include "LauBelleSymNR.hh"
34 #include "LauBreitWignerRes.hh"
35 #include "LauDabbaRes.hh"
36 #include "LauDaughters.hh"
37 #include "LauEFKLLMRes.hh"
38 #include "LauFlatteRes.hh"
39 #include "LauFlatNR.hh"
40 #include "LauGaussIncohRes.hh"
41 #include "LauGounarisSakuraiRes.hh"
42 #include "LauKappaRes.hh"
43 #include "LauLASSRes.hh"
44 #include "LauLASSBWRes.hh"
45 #include "LauLASSNRRes.hh"
48 #include "LauNRAmplitude.hh"
49 #include "LauRescatteringRes.hh"
50 #include "LauRescattering2Res.hh"
51 #include "LauPolNR.hh"
52 #include "LauPoleRes.hh"
53 #include "LauPolarFormFactorNR.hh"
55 #include "LauRelBreitWignerRes.hh"
56 #include "LauResonanceInfo.hh"
57 #include "LauResonanceMaker.hh"
58 #include "LauRhoOmegaMix.hh"
59 #include "LauSigmaRes.hh"
60 
62 
63 
65 
66 
68  nResDefMax_(0),
69  bwBarrierType_(LauBlattWeisskopfFactor::BWPrimeBarrier),
70  bwRestFrame_(LauBlattWeisskopfFactor::ResonanceFrame),
71  spinFormalism_(LauAbsResonance::Zemach_P),
72  summaryPrinted_(kFALSE)
73 {
74  this->createResonanceVector();
76 }
77 
79 {
80  for ( std::vector<LauBlattWeisskopfFactor*>::iterator iter = bwIndepFactors_.begin(); iter != bwIndepFactors_.end(); ++iter ) {
81  delete *iter;
82  }
83  bwIndepFactors_.clear();
84  for ( BWFactorCategoryMap::iterator iter = bwFactors_.begin(); iter != bwFactors_.end(); ++iter ) {
85  delete iter->second.bwFactor_;
86  }
87  bwFactors_.clear();
88 }
89 
91 {
92  if ( resonanceMaker_ == 0 ) {
94  }
95 
96  return *resonanceMaker_;
97 }
98 
100 {
101  // Function to create all possible resonances that this class supports.
102  // Also add in the sigma and kappa - but a special paramterisation is used
103  // instead of the PDG "pole mass and width" values.
104 
105  std::cout << "INFO in LauResonanceMaker::createResonanceVector : Setting up possible resonance states..." << std::endl;
106 
107  LauResonanceInfo* neutral(0);
108  LauResonanceInfo* positve(0);
109  LauResonanceInfo* negatve(0);
110 
111  // Define the resonance names and store them in the array
112  resInfo_.clear();
113  resInfo_.reserve(100);
114  // rho resonances name, mass, width, spin, charge, default BW category, BW radius parameter (defaults to 4.0)
115  // rho(770)
116  neutral = new LauResonanceInfo("rho0(770)", 0.77526, 0.1478, 1, 0, LauBlattWeisskopfFactor::Light, 5.3);
117  positve = new LauResonanceInfo("rho+(770)", 0.77511, 0.1491, 1, 1, LauBlattWeisskopfFactor::Light, 5.3);
118  negatve = positve->createChargeConjugate();
119  resInfo_.push_back( neutral );
120  resInfo_.push_back( positve );
121  resInfo_.push_back( negatve );
122  // The following two lines of code are placed here in order to allow the following, rather niche, scenario:
123  // 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.
124  // This can be acheived by using the rho0(770) record in one case and the rho0(770)_COPY record in the other.
125  neutral = neutral->createSharedParameterRecord("rho0(770)_COPY");
126  resInfo_.push_back( neutral );
127  // rho(1450)
128  neutral = new LauResonanceInfo("rho0(1450)", 1.465, 0.400, 1, 0, LauBlattWeisskopfFactor::Light );
129  positve = new LauResonanceInfo("rho+(1450)", 1.465, 0.400, 1, 1, LauBlattWeisskopfFactor::Light );
130  negatve = positve->createChargeConjugate();
131  resInfo_.push_back( neutral );
132  resInfo_.push_back( positve );
133  resInfo_.push_back( negatve );
134  // rho_3(1690)
135  neutral = new LauResonanceInfo("rho0_3(1690)", 1.686, 0.186, 3, 0, LauBlattWeisskopfFactor::Light );
136  positve = new LauResonanceInfo("rho+_3(1690)", 1.686, 0.186, 3, 1, LauBlattWeisskopfFactor::Light );
137  negatve = positve->createChargeConjugate();
138  resInfo_.push_back( neutral );
139  resInfo_.push_back( positve );
140  resInfo_.push_back( negatve );
141  // rho(1700)
142  neutral = new LauResonanceInfo("rho0(1700)", 1.720, 0.250, 1, 0, LauBlattWeisskopfFactor::Light );
143  positve = new LauResonanceInfo("rho+(1700)", 1.720, 0.250, 1, 1, LauBlattWeisskopfFactor::Light );
144  negatve = positve->createChargeConjugate();
145  resInfo_.push_back( neutral );
146  resInfo_.push_back( positve );
147  resInfo_.push_back( negatve );
148  // rho(1900)
149  neutral = new LauResonanceInfo("rho0(1900)", 1.909, 0.130, 1, 0, LauBlattWeisskopfFactor::Light );
150  positve = new LauResonanceInfo("rho+(1900)", 1.909, 0.130, 1, 1, LauBlattWeisskopfFactor::Light );
151  negatve = positve->createChargeConjugate();
152  resInfo_.push_back( neutral );
153  resInfo_.push_back( positve );
154  resInfo_.push_back( negatve );
155  // rho_3(1990)
156  neutral = new LauResonanceInfo("rho0_3(1990)", 1.982, 0.188, 3, 0, LauBlattWeisskopfFactor::Light );
157  positve = new LauResonanceInfo("rho+_3(1990)", 1.982, 0.188, 3, 1, LauBlattWeisskopfFactor::Light );
158  negatve = positve->createChargeConjugate();
159  resInfo_.push_back( neutral );
160  resInfo_.push_back( positve );
161  resInfo_.push_back( negatve );
162 
163  // K* resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
164  // K*(892)
165  neutral = new LauResonanceInfo("K*0(892)", 0.89581, 0.0474, 1, 0, LauBlattWeisskopfFactor::Kstar, 3.0);
166  positve = new LauResonanceInfo("K*+(892)", 0.89166, 0.0508, 1, 1, LauBlattWeisskopfFactor::Kstar, 3.0);
167  negatve = positve->createChargeConjugate();
168  resInfo_.push_back( neutral );
169  resInfo_.push_back( positve );
170  resInfo_.push_back( negatve );
171  // K*(1410)
172  neutral = new LauResonanceInfo("K*0(1410)", 1.414, 0.232, 1, 0, LauBlattWeisskopfFactor::Kstar );
173  positve = new LauResonanceInfo("K*+(1410)", 1.414, 0.232, 1, 1, LauBlattWeisskopfFactor::Kstar );
174  negatve = positve->createChargeConjugate();
175  resInfo_.push_back( neutral );
176  resInfo_.push_back( positve );
177  resInfo_.push_back( negatve );
178  // K*_0(1430)
179  neutral = new LauResonanceInfo("K*0_0(1430)", 1.425, 0.270, 0, 0, LauBlattWeisskopfFactor::Kstar );
180  positve = new LauResonanceInfo("K*+_0(1430)", 1.425, 0.270, 0, 1, LauBlattWeisskopfFactor::Kstar );
181  negatve = positve->createChargeConjugate();
182  resInfo_.push_back( neutral );
183  resInfo_.push_back( positve );
184  resInfo_.push_back( negatve );
185  // LASS nonresonant model
186  neutral = neutral->createSharedParameterRecord("LASSNR0");
187  positve = positve->createSharedParameterRecord("LASSNR+");
188  negatve = negatve->createSharedParameterRecord("LASSNR-");
189  resInfo_.push_back( neutral );
190  resInfo_.push_back( positve );
191  resInfo_.push_back( negatve );
192  // K*_2(1430)
193  neutral = new LauResonanceInfo("K*0_2(1430)", 1.4324, 0.109, 2, 0, LauBlattWeisskopfFactor::Kstar );
194  positve = new LauResonanceInfo("K*+_2(1430)", 1.4256, 0.0985, 2, 1, LauBlattWeisskopfFactor::Kstar );
195  negatve = positve->createChargeConjugate();
196  resInfo_.push_back( neutral );
197  resInfo_.push_back( positve );
198  resInfo_.push_back( negatve );
199  // K*(1680)
200  neutral = new LauResonanceInfo("K*0(1680)", 1.717, 0.322, 1, 0, LauBlattWeisskopfFactor::Kstar );
201  positve = new LauResonanceInfo("K*+(1680)", 1.717, 0.322, 1, 1, LauBlattWeisskopfFactor::Kstar );
202  negatve = positve->createChargeConjugate();
203  resInfo_.push_back( neutral );
204  resInfo_.push_back( positve );
205  resInfo_.push_back( negatve );
206  // K*(1950)
207  neutral = new LauResonanceInfo("K*0_0(1950)", 1.945, 0.201, 0, 0, LauBlattWeisskopfFactor::Kstar );
208  positve = new LauResonanceInfo("K*+_0(1950)", 1.945, 0.201, 0, 1, LauBlattWeisskopfFactor::Kstar );
209  negatve = positve->createChargeConjugate();
210  resInfo_.push_back( neutral );
211  resInfo_.push_back( positve );
212  resInfo_.push_back( negatve );
213 
214  // phi resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
215  // phi(1020)
216  neutral = new LauResonanceInfo("phi(1020)", 1.019461, 0.004266, 1, 0, LauBlattWeisskopfFactor::Light );
217  resInfo_.push_back( neutral );
218  // phi(1680)
219  neutral = new LauResonanceInfo("phi(1680)", 1.680, 0.150, 1, 0, LauBlattWeisskopfFactor::Light );
220  resInfo_.push_back( neutral );
221 
222  // f resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
223  // f_0(980)
224  neutral = new LauResonanceInfo("f_0(980)", 0.990, 0.070, 0, 0, LauBlattWeisskopfFactor::Light );
225  resInfo_.push_back( neutral );
226  // f_2(1270)
227  neutral = new LauResonanceInfo("f_2(1270)", 1.2751, 0.1851, 2, 0, LauBlattWeisskopfFactor::Light );
228  resInfo_.push_back( neutral );
229  // f_0(1370)
230  neutral = new LauResonanceInfo("f_0(1370)", 1.370, 0.350, 0, 0, LauBlattWeisskopfFactor::Light );
231  resInfo_.push_back( neutral );
232  // f'_0(1300), from Belle's Kspipi paper
233  neutral = new LauResonanceInfo("f'_0(1300)", 1.449, 0.126, 0, 0, LauBlattWeisskopfFactor::Light );
234  resInfo_.push_back( neutral );
235  // f_2(1430)
236  neutral = new LauResonanceInfo("f_2(1430)", 1.430, 0.150, 2, 0, LauBlattWeisskopfFactor::Light ); // PDG width in the range 13 - 150
237  resInfo_.push_back( neutral );
238  // f_0(1500)
239  neutral = new LauResonanceInfo("f_0(1500)", 1.505, 0.109, 0, 0, LauBlattWeisskopfFactor::Light );
240  resInfo_.push_back( neutral );
241  // f'_2(1525)
242  neutral = new LauResonanceInfo("f'_2(1525)", 1.525, 0.073, 2, 0, LauBlattWeisskopfFactor::Light );
243  resInfo_.push_back( neutral );
244  // f_2(1565)
245  neutral = new LauResonanceInfo("f_2(1565)", 1.562, 0.134, 2, 0, LauBlattWeisskopfFactor::Light );
246  resInfo_.push_back( neutral );
247  // f_2(1640)
248  neutral = new LauResonanceInfo("f_2(1640)", 1.639, 0.099, 2, 0, LauBlattWeisskopfFactor::Light );
249  resInfo_.push_back( neutral );
250  // f_0(1710)
251  neutral = new LauResonanceInfo("f_0(1710)", 1.722, 0.135, 0, 0, LauBlattWeisskopfFactor::Light );
252  resInfo_.push_back( neutral );
253  // f_2(1810)
254  neutral = new LauResonanceInfo("f_2(1810)", 1.816, 0.197, 2, 0, LauBlattWeisskopfFactor::Light );
255  resInfo_.push_back( neutral );
256  // f_2(1910)
257  neutral = new LauResonanceInfo("f_2(1910)", 1.903, 0.196, 2, 0, LauBlattWeisskopfFactor::Light );
258  resInfo_.push_back( neutral );
259  // f_2(1950)
260  neutral = new LauResonanceInfo("f_2(1950)", 1.944, 0.472, 2, 0, LauBlattWeisskopfFactor::Light );
261  resInfo_.push_back( neutral );
262  // f_2(2010)
263  neutral = new LauResonanceInfo("f_2(2010)", 2.011, 0.202, 2, 0, LauBlattWeisskopfFactor::Light );
264  resInfo_.push_back( neutral );
265  // f_0(2020)
266  neutral = new LauResonanceInfo("f_0(2020)", 1.992, 0.442, 0, 0, LauBlattWeisskopfFactor::Light );
267  resInfo_.push_back( neutral );
268  // f_4(2050)
269  neutral = new LauResonanceInfo("f_4(2050)", 2.018, 0.237, 4, 0, LauBlattWeisskopfFactor::Light );
270  resInfo_.push_back( neutral );
271  // f_0(2100)
272  neutral = new LauResonanceInfo("f_0(2100)", 2.101, 0.224, 0, 0, LauBlattWeisskopfFactor::Light );
273  resInfo_.push_back( neutral );
274 
275  // omega resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
276  // omega(782)
277  neutral = new LauResonanceInfo("omega(782)", 0.78265, 0.00849, 1, 0, LauBlattWeisskopfFactor::Light );
278  resInfo_.push_back( neutral );
279 
280  // a resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
281  // a_0(980)
282  neutral = new LauResonanceInfo("a0_0(980)", 0.980, 0.092, 0, 0, LauBlattWeisskopfFactor::Light );
283  positve = new LauResonanceInfo("a+_0(980)", 0.980, 0.092, 0, 1, LauBlattWeisskopfFactor::Light );
284  negatve = positve->createChargeConjugate();
285  resInfo_.push_back( neutral );
286  resInfo_.push_back( positve );
287  resInfo_.push_back( negatve );
288  // a_0(1450)
289  neutral = new LauResonanceInfo("a0_0(1450)", 1.474, 0.265, 0, 0, LauBlattWeisskopfFactor::Light );
290  positve = new LauResonanceInfo("a+_0(1450)", 1.474, 0.265, 0, 1, LauBlattWeisskopfFactor::Light );
291  negatve = positve->createChargeConjugate();
292  resInfo_.push_back( neutral );
293  resInfo_.push_back( positve );
294  resInfo_.push_back( negatve );
295  // a_2(1320)
296  neutral = new LauResonanceInfo("a0_2(1320)", 1.3190, 0.1050, 2, 0, LauBlattWeisskopfFactor::Light );
297  positve = new LauResonanceInfo("a+_2(1320)", 1.3190, 0.1050, 2, 1, LauBlattWeisskopfFactor::Light );
298  negatve = positve->createChargeConjugate();
299  resInfo_.push_back( neutral );
300  resInfo_.push_back( positve );
301  resInfo_.push_back( negatve );
302 
303  // charmonium resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
304  // chi_c0
305  neutral = new LauResonanceInfo("chi_c0", 3.41475, 0.0105, 0, 0, LauBlattWeisskopfFactor::Charmonium );
306  resInfo_.push_back( neutral );
307  // chi_c1
308  neutral = new LauResonanceInfo("chi_c1", 3.51066, 0.00084, 0, 0, LauBlattWeisskopfFactor::Charmonium );
309  resInfo_.push_back( neutral );
310  // chi_c2
311  neutral = new LauResonanceInfo("chi_c2", 3.55620, 0.00193, 2, 0, LauBlattWeisskopfFactor::Charmonium );
312  resInfo_.push_back( neutral );
313  // psi(3770)
314  neutral = new LauResonanceInfo("psi(3770)", 3.77313, 0.0272, 1, 0, LauBlattWeisskopfFactor::Charmonium );
315  resInfo_.push_back( neutral );
316  // X(3872)
317  neutral = new LauResonanceInfo("X(3872)", 3.87169, 0.0012, 1, 0, LauBlattWeisskopfFactor::Charmonium );
318  resInfo_.push_back( neutral );
319  // chi_c2(2P)
320  neutral = new LauResonanceInfo("chi_c2(2P)", 3.9272, 0.024, 2, 0, LauBlattWeisskopfFactor::Charmonium );
321  resInfo_.push_back( neutral );
322 
323  // unknown scalars name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
324  // sigma
325  neutral = new LauResonanceInfo("sigma0", 0.475, 0.550, 0, 0, LauBlattWeisskopfFactor::Light );
326  positve = new LauResonanceInfo("sigma+", 0.475, 0.550, 0, 1, LauBlattWeisskopfFactor::Light );
327  negatve = positve->createChargeConjugate();
328  resInfo_.push_back( neutral );
329  resInfo_.push_back( positve );
330  resInfo_.push_back( negatve );
331  // kappa
332  neutral = new LauResonanceInfo("kappa0", 0.682, 0.547, 0, 0, LauBlattWeisskopfFactor::Kstar );
333  positve = new LauResonanceInfo("kappa+", 0.682, 0.547, 0, 1, LauBlattWeisskopfFactor::Kstar );
334  negatve = positve->createChargeConjugate();
335  resInfo_.push_back( neutral );
336  resInfo_.push_back( positve );
337  resInfo_.push_back( negatve );
338  // dabba
339  neutral = new LauResonanceInfo("dabba0", 2.098, 0.520, 0, 0, LauBlattWeisskopfFactor::Charm );
340  positve = new LauResonanceInfo("dabba+", 2.098, 0.520, 0, 1, LauBlattWeisskopfFactor::Charm );
341  negatve = positve->createChargeConjugate();
342  resInfo_.push_back( neutral );
343  resInfo_.push_back( positve );
344  resInfo_.push_back( negatve );
345 
346  // excited charm states name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
347  // D*
348  neutral = new LauResonanceInfo("D*0", 2.00696, 0.0021, 1, 0, LauBlattWeisskopfFactor::Charm );
349  positve = new LauResonanceInfo("D*+", 2.01026, 83.4e-6, 1, 1, LauBlattWeisskopfFactor::Charm );
350  negatve = positve->createChargeConjugate();
351  resInfo_.push_back( neutral );
352  resInfo_.push_back( positve );
353  resInfo_.push_back( negatve );
354  // D*_0
355  neutral = new LauResonanceInfo("D*0_0", 2.318, 0.267, 0, 0, LauBlattWeisskopfFactor::Charm );
356  positve = new LauResonanceInfo("D*+_0", 2.403, 0.283, 0, 1, LauBlattWeisskopfFactor::Charm );
357  negatve = positve->createChargeConjugate();
358  resInfo_.push_back( neutral );
359  resInfo_.push_back( positve );
360  resInfo_.push_back( negatve );
361  // D*_2
362  //AVERAGE--neutral = new LauResonanceInfo("D*0_2", 2.4618, 0.049, 2, 0 );
363  neutral = new LauResonanceInfo("D*0_2", 2.4626, 0.049, 2, 0, LauBlattWeisskopfFactor::Charm );
364  positve = new LauResonanceInfo("D*+_2", 2.4643, 0.037, 2, 1, LauBlattWeisskopfFactor::Charm );
365  negatve = positve->createChargeConjugate();
366  resInfo_.push_back( neutral );
367  resInfo_.push_back( positve );
368  resInfo_.push_back( negatve );
369  // D1(2420)
370  neutral = new LauResonanceInfo("D0_1(2420)", 2.4214, 0.0274, 1, 0, LauBlattWeisskopfFactor::Charm );
371  positve = new LauResonanceInfo("D+_1(2420)", 2.4232, 0.025, 1, 1, LauBlattWeisskopfFactor::Charm );
372  negatve = positve->createChargeConjugate();
373  resInfo_.push_back( neutral );
374  resInfo_.push_back( positve );
375  resInfo_.push_back( negatve );
376  // D(2600)
377  //OLD--neutral = new LauResonanceInfo("D0(2600)", 2.6087, 0.093, 0, 0 );
378  //OLD--positve = new LauResonanceInfo("D+(2600)", 2.6213, 0.093, 0, 1 );
379  neutral = new LauResonanceInfo("D0(2600)", 2.612, 0.093, 0, 0, LauBlattWeisskopfFactor::Charm );
380  positve = new LauResonanceInfo("D+(2600)", 2.612, 0.093, 0, 1, LauBlattWeisskopfFactor::Charm );
381  negatve = positve->createChargeConjugate();
382  resInfo_.push_back( neutral );
383  resInfo_.push_back( positve );
384  resInfo_.push_back( negatve );
385  // D(2760)
386  //OLD-- neutral = new LauResonanceInfo("D0(2760)", 2.7633, 0.061, 1, 0 );
387  //OLD-- positve = new LauResonanceInfo("D+(2760)", 2.7697, 0.061, 1, 1 );
388  neutral = new LauResonanceInfo("D0(2760)", 2.761, 0.063, 1, 0, LauBlattWeisskopfFactor::Charm );
389  positve = new LauResonanceInfo("D+(2760)", 2.761, 0.063, 1, 1, LauBlattWeisskopfFactor::Charm );
390  negatve = positve->createChargeConjugate();
391  resInfo_.push_back( neutral );
392  resInfo_.push_back( positve );
393  resInfo_.push_back( negatve );
394  // D(2900)
395  neutral = new LauResonanceInfo("D0(3000)", 3.0, 0.15, 0, 0, LauBlattWeisskopfFactor::Charm );
396  resInfo_.push_back( neutral );
397  // D(3400)
398  neutral = new LauResonanceInfo("D0(3400)", 3.4, 0.15, 0, 0, LauBlattWeisskopfFactor::Charm );
399  resInfo_.push_back( neutral );
400 
401  // excited strange charm name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
402  // Ds*
403  positve = new LauResonanceInfo("Ds*+", 2.1121, 0.0019, 1, 1, LauBlattWeisskopfFactor::StrangeCharm );
404  negatve = positve->createChargeConjugate();
405  resInfo_.push_back( positve );
406  resInfo_.push_back( negatve );
407  // Ds0*(2317)
408  positve = new LauResonanceInfo("Ds*+_0(2317)", 2.3177, 0.0038, 0, 1, LauBlattWeisskopfFactor::StrangeCharm );
409  negatve = positve->createChargeConjugate();
410  resInfo_.push_back( positve );
411  resInfo_.push_back( negatve );
412  // Ds2*(2573)
413  positve = new LauResonanceInfo("Ds*+_2(2573)", 2.5719, 0.017, 2, 1, LauBlattWeisskopfFactor::StrangeCharm );
414  negatve = positve->createChargeConjugate();
415  resInfo_.push_back( positve );
416  resInfo_.push_back( negatve );
417  // Ds1*(2700)
418  positve = new LauResonanceInfo("Ds*+_1(2700)", 2.709, 0.117, 1, 1, LauBlattWeisskopfFactor::StrangeCharm );
419  negatve = positve->createChargeConjugate();
420  resInfo_.push_back( positve );
421  resInfo_.push_back( negatve );
422  // Ds1*(2860)
423  positve = new LauResonanceInfo("Ds*+_1(2860)", 2.862, 0.180, 1, 1, LauBlattWeisskopfFactor::StrangeCharm );
424  negatve = positve->createChargeConjugate();
425  resInfo_.push_back( positve );
426  resInfo_.push_back( negatve );
427  // Ds3*(2860)
428  positve = new LauResonanceInfo("Ds*+_3(2860)", 2.862, 0.058, 3, 1, LauBlattWeisskopfFactor::StrangeCharm );
429  negatve = positve->createChargeConjugate();
430  resInfo_.push_back( positve );
431  resInfo_.push_back( negatve );
432 
433  // excited bottom states name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
434  // B*
435  neutral = new LauResonanceInfo("B*0", 5.3252, 0.00, 1, 0, LauBlattWeisskopfFactor::Beauty, 6.0);
436  positve = new LauResonanceInfo("B*+", 5.3252, 0.00, 1, 1, LauBlattWeisskopfFactor::Beauty, 6.0);
437  negatve = positve->createChargeConjugate();
438  resInfo_.push_back( neutral );
439  resInfo_.push_back( positve );
440  resInfo_.push_back( negatve );
441 
442  // excited strange bottom name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
443  // Bs*
444  neutral = new LauResonanceInfo("Bs*0", 5.4154, 0.00, 1, 0, LauBlattWeisskopfFactor::StrangeBeauty, 6.0);
445  resInfo_.push_back( neutral );
446 
447  // nonresonant models name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
448  // Phase-space nonresonant model
449  neutral = new LauResonanceInfo("NonReson", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
450  resInfo_.push_back( neutral );
451  // Theory-based nonresonant model
452  neutral = new LauResonanceInfo("NRModel", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
453  resInfo_.push_back( neutral );
454  // Belle nonresonant models
455  neutral = new LauResonanceInfo("BelleSymNR", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
456  resInfo_.push_back( neutral );
457  neutral = new LauResonanceInfo("BelleNR", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
458  resInfo_.push_back( neutral );
459  positve = new LauResonanceInfo("BelleNR+", 0.0, 0.0, 0, 1, LauBlattWeisskopfFactor::Light );
460  negatve = positve->createChargeConjugate();
461  resInfo_.push_back( positve );
462  resInfo_.push_back( negatve );
463  neutral = new LauResonanceInfo("BelleNR_Swave", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
464  resInfo_.push_back( neutral );
465  positve = new LauResonanceInfo("BelleNR_Swave+",0.0, 0.0, 0, 1, LauBlattWeisskopfFactor::Light );
466  negatve = positve->createChargeConjugate();
467  resInfo_.push_back( positve );
468  resInfo_.push_back( negatve );
469  neutral = new LauResonanceInfo("BelleNR_Pwave", 0.0, 0.0, 1, 0, LauBlattWeisskopfFactor::Light );
470  resInfo_.push_back( neutral );
471  positve = new LauResonanceInfo("BelleNR_Pwave+",0.0, 0.0, 1, 1, LauBlattWeisskopfFactor::Light );
472  negatve = positve->createChargeConjugate();
473  resInfo_.push_back( positve );
474  resInfo_.push_back( negatve );
475  neutral = new LauResonanceInfo("BelleNR_Dwave", 0.0, 0.0, 2, 0, LauBlattWeisskopfFactor::Light );
476  resInfo_.push_back( neutral );
477  positve = new LauResonanceInfo("BelleNR_Dwave+",0.0, 0.0, 2, 1, LauBlattWeisskopfFactor::Light );
478  negatve = positve->createChargeConjugate();
479  resInfo_.push_back( positve );
480  resInfo_.push_back( negatve );
481  neutral = new LauResonanceInfo("BelleNR_Fwave", 0.0, 0.0, 3, 0, LauBlattWeisskopfFactor::Light );
482  resInfo_.push_back( neutral );
483  positve = new LauResonanceInfo("BelleNR_Fwave+",0.0, 0.0, 3, 1, LauBlattWeisskopfFactor::Light );
484  negatve = positve->createChargeConjugate();
485  resInfo_.push_back( positve );
486  resInfo_.push_back( negatve );
487  // Taylor expansion nonresonant model
488  neutral = new LauResonanceInfo("NRTaylor", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
489  resInfo_.push_back( neutral );
490  // Polynomial nonresonant models
491  neutral = new LauResonanceInfo("PolNR_S0", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
492  resInfo_.push_back( neutral );
493  neutral = new LauResonanceInfo("PolNR_S1", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
494  resInfo_.push_back( neutral );
495  neutral = new LauResonanceInfo("PolNR_S2", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
496  resInfo_.push_back( neutral );
497  neutral = new LauResonanceInfo("PolNR_P0", 0.0, 0.0, 1, 0, LauBlattWeisskopfFactor::Light );
498  resInfo_.push_back( neutral );
499  neutral = new LauResonanceInfo("PolNR_P1", 0.0, 0.0, 1, 0, LauBlattWeisskopfFactor::Light );
500  resInfo_.push_back( neutral );
501  neutral = new LauResonanceInfo("PolNR_P2", 0.0, 0.0, 1, 0, LauBlattWeisskopfFactor::Light );
502  resInfo_.push_back( neutral );
503 
504  // Fake resonances for S-Wave splines
505  neutral = new LauResonanceInfo("Spline_S0", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
506  resInfo_.push_back( neutral );
507  neutral = new LauResonanceInfo("Spline_S0_Bar", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
508  resInfo_.push_back( neutral );
509 
510  // Polar Form Factor nonresonant model
511  neutral = new LauResonanceInfo("PolarFFSymNR", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
512  resInfo_.push_back( neutral );
513  neutral = new LauResonanceInfo("PolarFFNR", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
514  resInfo_.push_back( neutral );
515 
516  // PiPi-KK Inelastic Scattering
517  neutral = new LauResonanceInfo("Rescattering", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
518  resInfo_.push_back( neutral );
519 
520  nResDefMax_ = resInfo_.size();
521 }
522 
524 {
525  // Check whether any BW factors have been created and bail out if so
526  if ( ! bwIndepFactors_.empty() ) {
527  std::cerr << "ERROR in LauResonanceMaker::setBWType : some barrier factors have already been created - cannot change the barrier type now!" << std::endl;
528  return;
529  }
530  for ( BWFactorCategoryMap::const_iterator iter = bwFactors_.begin(); iter != bwFactors_.end(); ++iter ) {
531  if ( iter->second.bwFactor_ != 0 ) {
532  std::cerr << "ERROR in LauResonanceMaker::setBWType : some barrier factors have already been created - cannot change the barrier type now!" << std::endl;
533  return;
534  }
535  }
536 
537  bwBarrierType_ = bwType;
538 }
539 
541 {
542  // Check whether any BW factors have been created and bail out if so
543  if ( ! bwIndepFactors_.empty() ) {
544  std::cerr << "ERROR in LauResonanceMaker::setBWBachelorRestFrame : some barrier factors have already been created - cannot change the rest frame now!" << std::endl;
545  return;
546  }
547  for ( BWFactorCategoryMap::const_iterator iter = bwFactors_.begin(); iter != bwFactors_.end(); ++iter ) {
548  if ( iter->second.bwFactor_ != 0 ) {
549  std::cerr << "ERROR in LauResonanceMaker::setBWBachelorRestFrame : some barrier factors have already been created - cannot change the rest frame now!" << std::endl;
550  return;
551  }
552  }
553 
554  bwRestFrame_ = restFrame;
555 }
556 
558 {
559  if ( summaryPrinted_ ) {
560  std::cerr << "ERROR in LauResonanceMaker::setSpinFormalism : cannot redefine the spin formalism after creating one or more resonances" << std::endl;
561  return;
562  }
563  spinFormalism_ = spinType;
564 }
565 
567 {
568  if ( bwCategory == LauBlattWeisskopfFactor::Default || bwCategory == LauBlattWeisskopfFactor::Indep ) {
569  std::cerr << "WARNING in LauResonanceMaker::setDefaultBWRadius : cannot set radius values for Default or Indep categories" << std::endl;
570  return;
571  }
572 
573  // Check if we have an information object for this category
574  BWFactorCategoryMap::iterator factor_iter = bwFactors_.find( bwCategory );
575  if ( factor_iter != bwFactors_.end() ) {
576  // If so, we can set the value in the information object
577  BlattWeisskopfCategoryInfo& categoryInfo = factor_iter->second;
578  categoryInfo.defaultRadius_ = bwRadius;
579 
580  // Then we can check if a LauBlattWeisskopfFactor object has been created for this category
581  LauBlattWeisskopfFactor* bwFactor = categoryInfo.bwFactor_;
582  if ( bwFactor != 0 ) {
583  // If it has then we can also set its radius value directly
584  LauParameter* radius = bwFactor->getRadiusParameter();
585  radius->value(bwRadius);
586  radius->initValue(bwRadius);
587  radius->genValue(bwRadius);
588  }
589  } else {
590  // If not then we just store the value to be used later
591  BlattWeisskopfCategoryInfo& categoryInfo = bwFactors_[bwCategory];
592  categoryInfo.bwFactor_ = 0;
593  categoryInfo.defaultRadius_ = bwRadius;
594  categoryInfo.radiusFixed_ = kTRUE;
595  }
596 }
597 
599 {
600  if ( bwCategory == LauBlattWeisskopfFactor::Default || bwCategory == LauBlattWeisskopfFactor::Indep ) {
601  std::cerr << "WARNING in LauResonanceMaker::fixBWRadius : cannot fix/float radius values for Default or Indep categories" << std::endl;
602  return;
603  }
604 
605  // Check if we have an information object for this category
606  BWFactorCategoryMap::iterator factor_iter = bwFactors_.find( bwCategory );
607  if ( factor_iter != bwFactors_.end() ) {
608  // If so, we can set the value in the information object
609  BlattWeisskopfCategoryInfo& categoryInfo = factor_iter->second;
610  categoryInfo.radiusFixed_ = fixRadius;
611 
612  // Then we can check if a LauBlattWeisskopfFactor object has been created for this category
613  LauBlattWeisskopfFactor* bwFactor = categoryInfo.bwFactor_;
614  if ( bwFactor != 0 ) {
615  // If it has then we can also fix/float its radius value directly
616  LauParameter* radius = bwFactor->getRadiusParameter();
617  radius->fixed(fixRadius);
618  }
619  } else {
620  // If not then we just store the value to be used later
621  BlattWeisskopfCategoryInfo& categoryInfo = bwFactors_[bwCategory];
622  categoryInfo.bwFactor_ = 0;
623  categoryInfo.defaultRadius_ = -1.0;
624  categoryInfo.radiusFixed_ = fixRadius;
625  }
626 }
627 
629 {
630  LauBlattWeisskopfFactor* bwFactor(0);
631 
632  // If this is an independent factor, create it and add it to the list of independent factors, then return it
633  if ( bwCategory == LauBlattWeisskopfFactor::Indep ) {
634  bwFactor = new LauBlattWeisskopfFactor( *resInfo, bwBarrierType_, bwRestFrame_, bwCategory );
635  bwIndepFactors_.push_back(bwFactor);
636  return bwFactor;
637  }
638 
639  // Otherwise, look up the category in the category information map
640  BWFactorCategoryMap::iterator factor_iter = bwFactors_.find( bwCategory );
641 
642  if ( factor_iter == bwFactors_.end() ) {
643  // If the category is currently undefined we need to create it
644  bwFactor = new LauBlattWeisskopfFactor( *resInfo, bwBarrierType_, bwRestFrame_, bwCategory );
645 
646  BlattWeisskopfCategoryInfo& categoryInfo = bwFactors_[bwCategory];
647  categoryInfo.bwFactor_ = bwFactor;
648  categoryInfo.defaultRadius_ = bwFactor->getRadiusParameter()->value();
649  categoryInfo.radiusFixed_ = kTRUE;
650  } else {
651  // If it exists, we can check if the factor object has been created
652  BlattWeisskopfCategoryInfo& categoryInfo = factor_iter->second;
653 
654  if ( categoryInfo.bwFactor_ != 0 ) {
655  // If so, simply clone it
656  const UInt_t resSpin = resInfo->getSpin();
657  bwFactor = categoryInfo.bwFactor_->createClone( resSpin );
658  } else {
659  // Otherwise we need to create it, using the default value if it has been set
660  if ( categoryInfo.defaultRadius_ >= 0.0 ) {
661  bwFactor = new LauBlattWeisskopfFactor( *resInfo, categoryInfo.defaultRadius_, bwBarrierType_, bwRestFrame_, bwCategory );
662  } else {
663  bwFactor = new LauBlattWeisskopfFactor( *resInfo, bwBarrierType_, bwRestFrame_, bwCategory );
664  }
665  categoryInfo.bwFactor_ = bwFactor;
666 
667  // Set whether the radius should be fixed/floated
668  LauParameter* radius = bwFactor->getRadiusParameter();
669  radius->fixed( categoryInfo.radiusFixed_ );
670  }
671  }
672 
673  return bwFactor;
674 }
675 
676 LauAbsResonance* LauResonanceMaker::getResonance(const LauDaughters* daughters, const TString& resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType, const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory)
677 {
678  // Routine to return the appropriate LauAbsResonance object given the resonance
679  // name (resName), which daughter is the bachelor track (resPairAmpInt = 1,2 or 3),
680  // and the resonance type ("BW" = Breit-Wigner, "Flatte" = Flatte distribution).
681 
682  // If this is the first resonance we are making, first print a summary of the formalism
683  if ( ! summaryPrinted_ ) {
684  std::cout << "INFO in LauResonanceMaker::getResonance : Freezing amplitude formalism:" << std::endl;
685  switch ( spinFormalism_ ) {
687  std::cout << " : Spin factors use Zemach spin tensors, with bachelor momentum in resonance rest frame" << std::endl;
688  break;
690  std::cout << " : Spin factors use Zemach spin tensors, with bachelor momentum in parent rest frame" << std::endl;
691  break;
693  std::cout << " : Spin factors use Covariant spin tensors" << std::endl;
694  break;
696  std::cout << " : Spin factors are just Legendre polynomials" << std::endl;
697  break;
698  }
699  switch ( bwBarrierType_ ) {
701  std::cout << " : Blatt-Weisskopf barrier factors are the 'non-primed' form" << std::endl;
702  break;
704  std::cout << " : Blatt-Weisskopf barrier factors are the 'primed' form" << std::endl;
705  break;
707  std::cout << " : Blatt-Weisskopf barrier factors are the exponential form" << std::endl;
708  break;
709  }
710  switch ( bwRestFrame_ ) {
712  std::cout << " : Blatt-Weisskopf barrier factors use bachelor momentum in parent rest frame" << std::endl;
713  break;
715  std::cout << " : Blatt-Weisskopf barrier factors use bachelor momentum in resonance rest frame" << std::endl;
716  break;
718  std::cout << " : Blatt-Weisskopf barrier factors use covariant expression" << std::endl;
719  break;
720  }
721 
722  summaryPrinted_ = kTRUE;
723  }
724 
725  // Loop over all possible resonance states we have defined in
726  // createResonanceVector() until we get a match with the name of the resonance
727 
728  LauResonanceInfo* resInfo(0);
729  for (std::vector<LauResonanceInfo*>::const_iterator iter=resInfo_.begin(); iter!=resInfo_.end(); ++iter) {
730 
731  if (resName == (*iter)->getName()) {
732  // We have recognised the resonance name.
733  std::cout<<"INFO in LauResonanceMaker::getResonance : Creating resonance: "<<resName<<std::endl;
734 
735  resInfo = (*iter);
736 
737  // stop looping
738  break;
739  }
740  }
741 
742  // if we couldn't find the right resonance then we should return a null pointer
743  if ( resInfo == 0 ) {
744  std::cout<<"ERROR in LauResonanceMaker::getResonance : Unable to locate resonance info for: "<<resName<<std::endl;
745  return 0;
746  }
747 
748  // Now construct the resonance using the specified type
749  LauAbsResonance* theResonance(0);
750  switch ( resType ) {
751 
752  case LauAbsResonance::BW :
753  // Simple non-relativistic Breit-Wigner
754  std::cout<<" : Using simple Breit-Wigner lineshape. "<<std::endl;
755  theResonance = new LauBreitWignerRes(resInfo, resPairAmpInt, daughters);
756  break;
757 
759  {
760  // Relativistic Breit-Wigner with Blatt-Weisskopf factors.
761  std::cout<<" : Using relativistic Breit-Wigner lineshape. "<<std::endl;
762  theResonance = new LauRelBreitWignerRes(resInfo, resPairAmpInt, daughters);
764  LauBlattWeisskopfFactor::BlattWeisskopfCategory resCategory = bwCategory;
765  if ( bwCategory == LauBlattWeisskopfFactor::Default ) {
766  resCategory = resInfo->getBWCategory();
767  }
768  LauBlattWeisskopfFactor* resBWFactor = this->getBWFactor( resCategory, resInfo );
769  LauBlattWeisskopfFactor* parBWFactor = this->getBWFactor( parCategory, resInfo );
770  theResonance->setBarrierRadii( resBWFactor, parBWFactor );
771  break;
772  }
773 
774  case LauAbsResonance::GS :
775  {
776  // Gounaris-Sakurai function to try and model the rho(770) better
777  std::cout<<" : Using Gounaris-Sakurai lineshape. "<<std::endl;
778  theResonance = new LauGounarisSakuraiRes(resInfo, resPairAmpInt, daughters);
780  LauBlattWeisskopfFactor::BlattWeisskopfCategory resCategory = bwCategory;
781  if ( bwCategory == LauBlattWeisskopfFactor::Default ) {
782  resCategory = resInfo->getBWCategory();
783  }
784  LauBlattWeisskopfFactor* resBWFactor = this->getBWFactor( resCategory, resInfo );
785  LauBlattWeisskopfFactor* parBWFactor = this->getBWFactor( parCategory, resInfo );
786  theResonance->setBarrierRadii( resBWFactor, parBWFactor );
787  break;
788  }
789 
791  // Flatte distribution - coupled channel Breit-Wigner
792  std::cout<<" : Using Flatte lineshape. "<<std::endl;
793  theResonance = new LauFlatteRes(resInfo, resPairAmpInt, daughters);
794  break;
795 
797  // Sigma model - should only be used for the pi-pi system
798  std::cout<<" : Using Sigma lineshape. "<<std::endl;
799  theResonance = new LauSigmaRes(resInfo, resPairAmpInt, daughters);
800  break;
801 
803  // Kappa model - should only be used for the K-pi system
804  std::cout<<" : Using Kappa lineshape. "<<std::endl;
805  theResonance = new LauKappaRes(resInfo, resPairAmpInt, daughters);
806  break;
807 
809  // Dabba model - should only be used for the D-pi system
810  std::cout<<" : Using Dabba lineshape. "<<std::endl;
811  theResonance = new LauDabbaRes(resInfo, resPairAmpInt, daughters);
812  break;
813 
814  case LauAbsResonance::LASS :
815  // LASS function to try and model the K-pi S-wave better
816  std::cout<<" : Using LASS lineshape. "<<std::endl;
817  theResonance = new LauLASSRes(resInfo, resPairAmpInt, daughters);
818  break;
819 
821  // LASS function to try and model the K-pi S-wave better
822  std::cout<<" : Using LASS lineshape (resonant part only). "<<std::endl;
823  theResonance = new LauLASSBWRes(resInfo, resPairAmpInt, daughters);
824  break;
825 
827  // LASS function to try and model the K-pi S-wave better
828  std::cout<<" : Using LASS lineshape (nonresonant part only). "<<std::endl;
829  theResonance = new LauLASSNRRes(resInfo, resPairAmpInt, daughters);
830  break;
831 
833  // EFKLLM form-factor description of the K-pi S-wave
834  std::cout<<" : Using EFKLLM lineshape. "<<std::endl;
835  theResonance = new LauEFKLLMRes(resInfo, resPairAmpInt, daughters);
836  break;
837 
839  // K-matrix description of S-wave
840  std::cerr<<"ERROR in LauResonanceMaker::getResonance : K-matrix type specified, which should be separately handled."<<std::endl;
841  break;
842 
844  // uniform NR amplitude - arguments are there to preserve the interface
845  std::cout<<" : Using uniform NR lineshape. "<<std::endl;
846  // we override resPairAmpInt here
847  theResonance = new LauFlatNR(resInfo, 0, daughters);
848  break;
849 
851  // NR amplitude model - arguments are there to preserve the interface
852  std::cout<<" : Using NR-model lineshape. "<<std::endl;
853  // we override resPairAmpInt here
854  theResonance = new LauNRAmplitude(resInfo, 0, daughters);
855  break;
856 
859  // Belle NR amplitude model - arguments are there to preserve the interface
860  std::cout<<" : Using Belle NR lineshape. "<<std::endl;
861  theResonance = new LauBelleNR(resInfo, resType, resPairAmpInt, daughters);
862  break;
863 
867  // Belle NR amplitude model - arguments are there to preserve the interface
868  std::cout<<" : Using Belle symmetric NR lineshape. "<<std::endl;
869  theResonance = new LauBelleSymNR(resInfo, resType, resPairAmpInt, daughters);
870  break;
871 
873  // Polynomial NR amplitude model - arguments are there to preserve the interface
874  std::cout<<" : Using polynomial NR lineshape. "<<std::endl;
875  theResonance = new LauPolNR(resInfo, resPairAmpInt, daughters);
876  break;
877 
878  case LauAbsResonance::Pole :
879  // Scalar pole model
880  std::cout<<" : Using the scalar Pole lineshape.. "<<std::endl;
881  theResonance = new LauPoleRes(resInfo, resPairAmpInt, daughters);
882  break;
883 
885  // Polar Form Factor NR amplitude model - arguments are there to preserve the interface
886  std::cout<<" : Using Polar FormFactor NR lineshape.. "<<std::endl;
887  theResonance = new LauPolarFormFactorNR(resInfo, resType, resPairAmpInt, daughters);
888  break;
889 
892  // Polar Form Factor NR amplitude model - arguments are there to preserve the interface
893  std::cout<<" : Using Polar FormFactor Symetric NR lineshape. "<<std::endl;
894  theResonance = new LauPolarFormFactorSymNR(resInfo, resType, resPairAmpInt, daughters);
895  break;
896 
899  // KKPiPi Inelastic Scattering amplitude - arguments are there to preserve the interface
900  std::cout<<" : KKPiPi Inelastic Scattering amplitude lineshape. "<<std::endl;
901  theResonance = new LauRescatteringRes(resInfo, resType, resPairAmpInt, daughters);
902  break;
903 
905  // KKPiPi Inelastic Scattering amplitude - arguments are there to preserve the interface
906  std::cout<<" : KKPiPi Inelastic Scattering amplitude lineshape. "<<std::endl;
907  theResonance = new LauRescattering2Res(resInfo, resPairAmpInt, daughters);
908  break;
909 
911  // Model independent partial wave
912  std::cout<<" : Using model independent partial wave lineshape (magnitude and phase). "<<std::endl;
913  theResonance = new LauModIndPartWaveMagPhase(resInfo, resPairAmpInt, daughters);
914  break;
915 
917  // Model independent partial wave
918  std::cout<<" : Using model independent partial wave lineshape (real and imaginary part). "<<std::endl;
919  theResonance = new LauModIndPartWaveRealImag(resInfo, resPairAmpInt, daughters);
920  break;
921 
923  // Incoherent Gaussian
924  std::cout<<" : Using incoherent Gaussian lineshape. "<<std::endl;
925  theResonance = new LauGaussIncohRes(resInfo, resPairAmpInt, daughters);
926  break;
927 
932  // Rho-omega mass mixing model
933  std::cout<<" : Using rho-omega mass mixing lineshape. "<<std::endl;
934  theResonance = new LauRhoOmegaMix(resInfo, resType, resPairAmpInt, daughters);
936  LauBlattWeisskopfFactor::BlattWeisskopfCategory resCategory = bwCategory;
937  if ( bwCategory == LauBlattWeisskopfFactor::Default ) {
938  resCategory = resInfo->getBWCategory();
939  }
940  LauBlattWeisskopfFactor* resBWFactor = this->getBWFactor( resCategory, resInfo );
941  LauBlattWeisskopfFactor* parBWFactor = this->getBWFactor( parCategory, resInfo );
942  theResonance->setBarrierRadii( resBWFactor, parBWFactor );
943  break;
944 
945  }
946 
947  // Set the spin formalism choice
948  theResonance->setSpinType( spinFormalism_ );
949 
950  return theResonance;
951 }
952 
953 Int_t LauResonanceMaker::resTypeInt(const TString& name) const
954 {
955  // Internal function that returns the resonance integer, specified by the
956  // order of the resonance vector defined in createResonanceVector(),
957  // for a given resonance name.
958  Int_t resTypInt(-99);
959  Int_t i(0);
960 
961  for (std::vector<LauResonanceInfo*>::const_iterator iter=resInfo_.begin(); iter!=resInfo_.end(); ++iter) {
962 
963  if (name.BeginsWith((*iter)->getName(), TString::kExact) == kTRUE) {
964  // We have recognised the resonance from those that are available
965  resTypInt = i;
966  return resTypInt;
967  }
968  ++i;
969  }
970 
971  return resTypInt;
972 }
973 
974 void LauResonanceMaker::printAll( std::ostream& stream ) const
975 {
976  for ( std::vector<LauResonanceInfo*>::const_iterator iter = resInfo_.begin(); iter != resInfo_.end(); ++iter ) {
977  stream << (**iter) << std::endl;
978  }
979 }
980 
981 LauResonanceInfo* LauResonanceMaker::getResInfo(const TString& resName) const
982 {
983  LauResonanceInfo* resInfo(0);
984  for (std::vector<LauResonanceInfo*>::const_iterator iter=resInfo_.begin(); iter!=resInfo_.end(); ++iter) {
985 
986  if (resName == (*iter)->getName()) {
987  // We have recognised the resonance name.
988  resInfo = (*iter);
989  // stop looping
990  break;
991  }
992  }
993 
994  return resInfo;
995 
996 }
File containing declaration of LauNRAmplitude class.
Bool_t fixed() const
Check whether the parameter is fixed or floated.
File containing declaration of LauRescatteringRes class.
LauBlattWeisskopfFactor * getBWFactor(const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory, const LauResonanceInfo *resInfo)
Retrieve Blatt-Weisskopf factor for the given category.
LauBlattWeisskopfFactor::BlattWeisskopfCategory getBWCategory() const
Retrieve the BW category of the resonant particle.
Class for defining the LASS resonance model.
Definition: LauLASSRes.hh:45
void setBarrierRadii(LauBlattWeisskopfFactor *resFactor, LauBlattWeisskopfFactor *parFactor)
Set the form factor model and parameters.
File containing declaration of LauResonanceInfo class.
void setSpinFormalism(const LauAbsResonance::LauSpinType spinType)
Set the spin formalism to be used for all resonances.
ClassImp(LauAbsCoeffSet)
LauBlattWeisskopfFactor::RestFrame bwRestFrame_
The rest frame in which the bachelor momentum used in the Blatt-Weisskopf factors should be calculate...
BWFactorCategoryMap bwFactors_
The Blatt-Weisskopf factor objects (and related information) for each category.
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:47
static LauResonanceMaker * resonanceMaker_
The singleton instance.
std::vector< LauResonanceInfo * > resInfo_
The known resonances.
File containing declaration of LauBelleNR class.
void setBWType(const LauBlattWeisskopfFactor::BarrierType bwType)
Set the type of BW factor (for all categories)
File containing declaration of LauPolarFormFactorNR 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.
const LauParameter * getRadiusParameter() const
Retrieve the radius parameter.
File containing declaration of LauModIndPartWaveRealImag class.
LauBlattWeisskopfFactor * bwFactor_
The BW factor object.
File containing declaration of LauGounarisSakuraiRes class.
Double_t defaultRadius_
The default value for the radius in this category.
File containing declaration of LauBelleSymNR class.
File containing declaration of LauModIndPartWaveMagPhase class.
File containing declaration of LauBreitWignerRes class.
LauSpinType
Define the allowed spin formalisms.
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:50
Class for defining the non resonant part of the LASS model.
Definition: LauLASSNRRes.hh:45
File containing declaration of LauEFKLLMRes class.
Class for defining the Sigma resonance model.
Definition: LauSigmaRes.hh:45
Singleton factory class for creating resonances.
void setBWBachelorRestFrame(const LauBlattWeisskopfFactor::RestFrame restFrame)
Set the rest frame in which the bachelor momentum should be calculated (for all BW categories) ...
void createResonanceVector()
Create the list of known resonances.
UInt_t getSpin() const
Retrieve the spin of the resonant particle.
File containing declaration of LauPolarFormFactorSymNR class.
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:48
Class for defining a nonresonant form factor model.
Class for defining the NR amplitude model.
Class for defining the simple Breit-Wigner resonance model.
Bool_t summaryPrinted_
Boolean flag to control printing a summary of the formalism to be used when the first resonance is cr...
Class for defining a model independent partial wave component where the amplitudes are parameterised ...
LauBlattWeisskopfFactor::BarrierType bwBarrierType_
The type of the Blatt-Weisskopf barrier to use for all resonances.
Class for defining the terms of Babar nonresonant model.
Definition: LauPolNR.hh:47
Class for defining the Kappa resonance model.
Definition: LauKappaRes.hh:46
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:49
File containing declaration of LauRescattering2Res class.
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.
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.
LauAbsResonance::LauSpinType spinFormalism_
The spin formalism that should be used for all resonances.
LauResonanceMaker()
Constructor.
BarrierType
Define the allowed types of barrier factors.
Class for defining the rescattering model.
LauAbsResonance * getResonance(const LauDaughters *daughters, const TString &resName, const Int_t resPairAmpInt, const LauAbsResonance::LauResonanceModel resType, const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory=LauBlattWeisskopfFactor::Default)
Create a resonance.
Class for defining an incoherent resonance with a Gaussian mass dependence.
Bool_t radiusFixed_
Whether or not the radius value for this category should be fixed in the fit.
Class for defining the resonant part of the LASS model.
Definition: LauLASSBWRes.hh:45
Class for defining a pole-type resonance model.
Definition: LauPoleRes.hh:46
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc...
Class for defining the Dabba resonance model.
Definition: LauDabbaRes.hh:45
Double_t initValue() const
The initial value of the parameter.
File containing declaration of LauAbsResonance class.
RestFrame
Define the rest frame in which the momentum should be calculated (only relevant for bachelor) ...
Class that implements the Blatt-Weisskopf barrier factor.
Class for defining an alternative rescattering model.
void setSpinType(const LauSpinType spinType)
Set the spin formalism to be used.
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.
File containing declaration of LauPoleRes class.
Class for defining a uniform nonresonant amplitude.
Definition: LauFlatNR.hh:44
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:45
Double_t genValue() const
The value generated for the parameter.
Class for defining a nonresonant form factor model.
File containing declaration of LauGaussIncohRes class.
LauResonanceInfo * getResInfo(const TString &resName) const
Get the information for the given resonance name.
Data structure to store information on a given Blatt-Weisskopf category.
File containing declaration of LauKappaRes class.
File containing declaration of LauFlatteRes class.
File containing declaration of LauFlatNR class.