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