laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
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 "LauResonanceMaker.hh"
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 "LauFlatNR.hh"
39 #include "LauFlatteRes.hh"
40 #include "LauGaussIncohRes.hh"
41 #include "LauGounarisSakuraiRes.hh"
42 #include "LauKappaRes.hh"
43 #include "LauLASSBWRes.hh"
44 #include "LauLASSNRRes.hh"
45 #include "LauLASSRes.hh"
48 #include "LauNRAmplitude.hh"
49 #include "LauPolNR.hh"
50 #include "LauPolarFormFactorNR.hh"
52 #include "LauPoleRes.hh"
53 #include "LauRelBreitWignerRes.hh"
54 #include "LauRescattering2Res.hh"
55 #include "LauRescatteringRes.hh"
56 #include "LauResonanceInfo.hh"
57 #include "LauRhoOmegaMix.hh"
58 #include "LauSigmaRes.hh"
59 
60 #include <iostream>
61 
63 
65  nResDefMax_( 0 ),
66  bwBarrierType_( LauBlattWeisskopfFactor::BWPrimeBarrier ),
67  bwRestFrame_( LauBlattWeisskopfFactor::ResonanceFrame ),
68  spinFormalism_( LauAbsResonance::Zemach_P ),
69  summaryPrinted_( kFALSE )
70 {
71  this->createResonanceVector();
72  this->setDefaultBWRadius( LauBlattWeisskopfFactor::Parent, 4.0 );
73 }
74 
76 {
77  for ( std::vector<LauBlattWeisskopfFactor*>::iterator iter = bwIndepFactors_.begin();
78  iter != bwIndepFactors_.end();
79  ++iter ) {
80  delete *iter;
81  }
82  bwIndepFactors_.clear();
83  for ( BWFactorCategoryMap::iterator iter = bwFactors_.begin(); iter != bwFactors_.end(); ++iter ) {
84  delete iter->second.bwFactor_;
85  }
86  bwFactors_.clear();
87 }
88 
90 {
91  if ( resonanceMaker_ == 0 ) {
93  }
94 
95  return *resonanceMaker_;
96 }
97 
99 {
100  // Function to create all possible resonances that this class supports.
101  // Also add in the sigma and kappa - but a special paramterisation is used
102  // instead of the PDG "pole mass and width" values.
103 
104  std::cout << "INFO in LauResonanceMaker::createResonanceVector : Setting up possible resonance states..."
105  << 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 =
117  new LauResonanceInfo( "rho0(770)", 0.77526, 0.1478, 1, 0, LauBlattWeisskopfFactor::Light, 5.3 );
118  positve =
119  new LauResonanceInfo( "rho+(770)", 0.77511, 0.1491, 1, 1, LauBlattWeisskopfFactor::Light, 5.3 );
120  negatve = positve->createChargeConjugate();
121  resInfo_.push_back( neutral );
122  resInfo_.push_back( positve );
123  resInfo_.push_back( negatve );
124  // The following two lines of code are placed here in order to allow the following, rather niche, scenario:
125  // 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.
126  // This can be acheived by using the rho0(770) record in one case and the rho0(770)_COPY record in the other.
127  neutral = neutral->createSharedParameterRecord( "rho0(770)_COPY" );
128  resInfo_.push_back( neutral );
129  // rho(1450)
130  neutral = new LauResonanceInfo( "rho0(1450)", 1.465, 0.400, 1, 0, LauBlattWeisskopfFactor::Light );
131  positve = new LauResonanceInfo( "rho+(1450)", 1.465, 0.400, 1, 1, LauBlattWeisskopfFactor::Light );
132  negatve = positve->createChargeConjugate();
133  resInfo_.push_back( neutral );
134  resInfo_.push_back( positve );
135  resInfo_.push_back( negatve );
136  // rho_3(1690)
137  neutral =
138  new LauResonanceInfo( "rho0_3(1690)", 1.686, 0.186, 3, 0, LauBlattWeisskopfFactor::Light );
139  positve =
140  new LauResonanceInfo( "rho+_3(1690)", 1.686, 0.186, 3, 1, LauBlattWeisskopfFactor::Light );
141  negatve = positve->createChargeConjugate();
142  resInfo_.push_back( neutral );
143  resInfo_.push_back( positve );
144  resInfo_.push_back( negatve );
145  // rho(1700)
146  neutral = new LauResonanceInfo( "rho0(1700)", 1.720, 0.250, 1, 0, LauBlattWeisskopfFactor::Light );
147  positve = new LauResonanceInfo( "rho+(1700)", 1.720, 0.250, 1, 1, LauBlattWeisskopfFactor::Light );
148  negatve = positve->createChargeConjugate();
149  resInfo_.push_back( neutral );
150  resInfo_.push_back( positve );
151  resInfo_.push_back( negatve );
152  // rho(1900)
153  neutral = new LauResonanceInfo( "rho0(1900)", 1.909, 0.130, 1, 0, LauBlattWeisskopfFactor::Light );
154  positve = new LauResonanceInfo( "rho+(1900)", 1.909, 0.130, 1, 1, LauBlattWeisskopfFactor::Light );
155  negatve = positve->createChargeConjugate();
156  resInfo_.push_back( neutral );
157  resInfo_.push_back( positve );
158  resInfo_.push_back( negatve );
159  // rho_3(1990)
160  neutral =
161  new LauResonanceInfo( "rho0_3(1990)", 1.982, 0.188, 3, 0, LauBlattWeisskopfFactor::Light );
162  positve =
163  new LauResonanceInfo( "rho+_3(1990)", 1.982, 0.188, 3, 1, LauBlattWeisskopfFactor::Light );
164  negatve = positve->createChargeConjugate();
165  resInfo_.push_back( neutral );
166  resInfo_.push_back( positve );
167  resInfo_.push_back( negatve );
168 
169  // K* resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
170  // K*(892)
171  neutral =
172  new LauResonanceInfo( "K*0(892)", 0.89555, 0.0473, 1, 0, LauBlattWeisskopfFactor::Kstar, 3.0 );
173  positve =
174  new LauResonanceInfo( "K*+(892)", 0.89166, 0.0508, 1, 1, LauBlattWeisskopfFactor::Kstar, 3.0 );
175  negatve = positve->createChargeConjugate();
176  resInfo_.push_back( neutral );
177  resInfo_.push_back( positve );
178  resInfo_.push_back( negatve );
179  // K*(1410)
180  neutral = new LauResonanceInfo( "K*0(1410)", 1.414, 0.232, 1, 0, LauBlattWeisskopfFactor::Kstar );
181  positve = new LauResonanceInfo( "K*+(1410)", 1.414, 0.232, 1, 1, LauBlattWeisskopfFactor::Kstar );
182  negatve = positve->createChargeConjugate();
183  resInfo_.push_back( neutral );
184  resInfo_.push_back( positve );
185  resInfo_.push_back( negatve );
186  // K*_0(1430)
187  neutral = new LauResonanceInfo( "K*0_0(1430)", 1.425, 0.270, 0, 0, LauBlattWeisskopfFactor::Kstar );
188  positve = new LauResonanceInfo( "K*+_0(1430)", 1.425, 0.270, 0, 1, LauBlattWeisskopfFactor::Kstar );
189  negatve = positve->createChargeConjugate();
190  resInfo_.push_back( neutral );
191  resInfo_.push_back( positve );
192  resInfo_.push_back( negatve );
193  // LASS nonresonant model
194  neutral = neutral->createSharedParameterRecord( "LASSNR0" );
195  positve = positve->createSharedParameterRecord( "LASSNR+" );
196  negatve = negatve->createSharedParameterRecord( "LASSNR-" );
197  resInfo_.push_back( neutral );
198  resInfo_.push_back( positve );
199  resInfo_.push_back( negatve );
200  // K*_2(1430)
201  neutral =
202  new LauResonanceInfo( "K*0_2(1430)", 1.4324, 0.109, 2, 0, LauBlattWeisskopfFactor::Kstar );
203  positve =
204  new LauResonanceInfo( "K*+_2(1430)", 1.4273, 0.100, 2, 1, LauBlattWeisskopfFactor::Kstar );
205  negatve = positve->createChargeConjugate();
206  resInfo_.push_back( neutral );
207  resInfo_.push_back( positve );
208  resInfo_.push_back( negatve );
209  // K*(1680)
210  neutral = new LauResonanceInfo( "K*0(1680)", 1.718, 0.322, 1, 0, LauBlattWeisskopfFactor::Kstar );
211  positve = new LauResonanceInfo( "K*+(1680)", 1.718, 0.322, 1, 1, LauBlattWeisskopfFactor::Kstar );
212  negatve = positve->createChargeConjugate();
213  resInfo_.push_back( neutral );
214  resInfo_.push_back( positve );
215  resInfo_.push_back( negatve );
216  // K*(1950)
217  neutral = new LauResonanceInfo( "K*0_0(1950)", 1.945, 0.201, 0, 0, LauBlattWeisskopfFactor::Kstar );
218  positve = new LauResonanceInfo( "K*+_0(1950)", 1.945, 0.201, 0, 1, LauBlattWeisskopfFactor::Kstar );
219  negatve = positve->createChargeConjugate();
220  resInfo_.push_back( neutral );
221  resInfo_.push_back( positve );
222  resInfo_.push_back( negatve );
223 
224  // phi resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
225  // phi(1020)
226  neutral =
227  new LauResonanceInfo( "phi(1020)", 1.019461, 0.004249, 1, 0, LauBlattWeisskopfFactor::Light );
228  resInfo_.push_back( neutral );
229  // phi(1680)
230  neutral = new LauResonanceInfo( "phi(1680)", 1.680, 0.150, 1, 0, LauBlattWeisskopfFactor::Light );
231  resInfo_.push_back( neutral );
232 
233  // f resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
234  // f_0(980)
235  neutral = new LauResonanceInfo( "f_0(980)", 0.990, 0.070, 0, 0, LauBlattWeisskopfFactor::Light );
236  resInfo_.push_back( neutral );
237  // f_2(1270)
238  neutral = new LauResonanceInfo( "f_2(1270)", 1.2755, 0.1867, 2, 0, LauBlattWeisskopfFactor::Light );
239  resInfo_.push_back( neutral );
240  // f_0(1370)
241  neutral = new LauResonanceInfo( "f_0(1370)", 1.370, 0.350, 0, 0, LauBlattWeisskopfFactor::Light );
242  resInfo_.push_back( neutral );
243  // f'_0(1300), from Belle's Kspipi paper
244  neutral = new LauResonanceInfo( "f'_0(1300)", 1.449, 0.126, 0, 0, LauBlattWeisskopfFactor::Light );
245  resInfo_.push_back( neutral );
246  // f_2(1430)
247  neutral = new LauResonanceInfo( "f_2(1430)",
248  1.430,
249  0.150,
250  2,
251  0,
252  LauBlattWeisskopfFactor::Light ); // PDG width in the range 13 - 150
253  resInfo_.push_back( neutral );
254  // f_0(1500)
255  neutral = new LauResonanceInfo( "f_0(1500)", 1.506, 0.112, 0, 0, LauBlattWeisskopfFactor::Light );
256  resInfo_.push_back( neutral );
257  // f'_2(1525)
258  neutral = new LauResonanceInfo( "f'_2(1525)", 1.5174, 0.086, 2, 0, LauBlattWeisskopfFactor::Light );
259  resInfo_.push_back( neutral );
260  // f_2(1565)
261  neutral = new LauResonanceInfo( "f_2(1565)", 1.542, 0.122, 2, 0, LauBlattWeisskopfFactor::Light );
262  resInfo_.push_back( neutral );
263  // f_2(1640)
264  neutral = new LauResonanceInfo( "f_2(1640)", 1.639, 0.099, 2, 0, LauBlattWeisskopfFactor::Light );
265  resInfo_.push_back( neutral );
266  // f_0(1710)
267  neutral = new LauResonanceInfo( "f_0(1710)", 1.704, 0.123, 0, 0, LauBlattWeisskopfFactor::Light );
268  resInfo_.push_back( neutral );
269  // f_2(1810)
270  neutral = new LauResonanceInfo( "f_2(1810)", 1.815, 0.197, 2, 0, LauBlattWeisskopfFactor::Light );
271  resInfo_.push_back( neutral );
272  // f_2(1910)
273  neutral = new LauResonanceInfo( "f_2(1910)", 1.900, 0.167, 2, 0, LauBlattWeisskopfFactor::Light );
274  resInfo_.push_back( neutral );
275  // f_2(1950)
276  neutral = new LauResonanceInfo( "f_2(1950)", 1.936, 0.464, 2, 0, LauBlattWeisskopfFactor::Light );
277  resInfo_.push_back( neutral );
278  // f_2(2010)
279  neutral = new LauResonanceInfo( "f_2(2010)", 2.011, 0.202, 2, 0, LauBlattWeisskopfFactor::Light );
280  resInfo_.push_back( neutral );
281  // f_0(2020)
282  neutral = new LauResonanceInfo( "f_0(2020)", 1.992, 0.442, 0, 0, LauBlattWeisskopfFactor::Light );
283  resInfo_.push_back( neutral );
284  // f_4(2050)
285  neutral = new LauResonanceInfo( "f_4(2050)", 2.018, 0.237, 4, 0, LauBlattWeisskopfFactor::Light );
286  resInfo_.push_back( neutral );
287  // f_0(2100)
288  neutral = new LauResonanceInfo( "f_0(2100)", 2.086, 0.284, 0, 0, LauBlattWeisskopfFactor::Light );
289  resInfo_.push_back( neutral );
290 
291  // omega resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
292  // omega(782)
293  neutral =
294  new LauResonanceInfo( "omega(782)", 0.78265, 0.00849, 1, 0, LauBlattWeisskopfFactor::Light );
295  resInfo_.push_back( neutral );
296 
297  // a resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
298  // a_0(980)
299  neutral = new LauResonanceInfo( "a0_0(980)", 0.980, 0.092, 0, 0, LauBlattWeisskopfFactor::Light );
300  positve = new LauResonanceInfo( "a+_0(980)", 0.980, 0.092, 0, 1, LauBlattWeisskopfFactor::Light );
301  negatve = positve->createChargeConjugate();
302  resInfo_.push_back( neutral );
303  resInfo_.push_back( positve );
304  resInfo_.push_back( negatve );
305  // a_0(1450)
306  neutral = new LauResonanceInfo( "a0_0(1450)", 1.474, 0.265, 0, 0, LauBlattWeisskopfFactor::Light );
307  positve = new LauResonanceInfo( "a+_0(1450)", 1.474, 0.265, 0, 1, LauBlattWeisskopfFactor::Light );
308  negatve = positve->createChargeConjugate();
309  resInfo_.push_back( neutral );
310  resInfo_.push_back( positve );
311  resInfo_.push_back( negatve );
312  // a_2(1320)
313  neutral =
314  new LauResonanceInfo( "a0_2(1320)", 1.3169, 0.1050, 2, 0, LauBlattWeisskopfFactor::Light );
315  positve =
316  new LauResonanceInfo( "a+_2(1320)", 1.3169, 0.1050, 2, 1, LauBlattWeisskopfFactor::Light );
317  negatve = positve->createChargeConjugate();
318  resInfo_.push_back( neutral );
319  resInfo_.push_back( positve );
320  resInfo_.push_back( negatve );
321 
322  // charmonium resonances name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
323  // chi_c0
324  neutral =
325  new LauResonanceInfo( "chi_c0", 3.41471, 0.0108, 0, 0, LauBlattWeisskopfFactor::Charmonium );
326  resInfo_.push_back( neutral );
327  // chi_c1
328  neutral =
329  new LauResonanceInfo( "chi_c1", 3.51067, 0.00084, 0, 0, LauBlattWeisskopfFactor::Charmonium );
330  resInfo_.push_back( neutral );
331  // chi_c2
332  neutral =
333  new LauResonanceInfo( "chi_c2", 3.55617, 0.00197, 2, 0, LauBlattWeisskopfFactor::Charmonium );
334  resInfo_.push_back( neutral );
335  // psi(3770)
336  neutral =
337  new LauResonanceInfo( "psi(3770)", 3.7737, 0.0272, 1, 0, LauBlattWeisskopfFactor::Charmonium );
338  resInfo_.push_back( neutral );
339  // X(3872)
340  neutral =
341  new LauResonanceInfo( "X(3872)", 3.87169, 0.0012, 1, 0, LauBlattWeisskopfFactor::Charmonium );
342  resInfo_.push_back( neutral );
343  // chi_c2(2P)
344  neutral =
345  new LauResonanceInfo( "chi_c2(2P)", 3.9222, 0.0353, 2, 0, LauBlattWeisskopfFactor::Charmonium );
346  resInfo_.push_back( neutral );
347 
348  // unknown scalars name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
349  // sigma
350  neutral = new LauResonanceInfo( "sigma0", 0.475, 0.550, 0, 0, LauBlattWeisskopfFactor::Light );
351  positve = new LauResonanceInfo( "sigma+", 0.475, 0.550, 0, 1, LauBlattWeisskopfFactor::Light );
352  negatve = positve->createChargeConjugate();
353  resInfo_.push_back( neutral );
354  resInfo_.push_back( positve );
355  resInfo_.push_back( negatve );
356  // kappa
357  neutral = new LauResonanceInfo( "kappa0", 0.824, 0.478, 0, 0, LauBlattWeisskopfFactor::Kstar );
358  positve = new LauResonanceInfo( "kappa+", 0.824, 0.478, 0, 1, LauBlattWeisskopfFactor::Kstar );
359  negatve = positve->createChargeConjugate();
360  resInfo_.push_back( neutral );
361  resInfo_.push_back( positve );
362  resInfo_.push_back( negatve );
363  // dabba
364  neutral = new LauResonanceInfo( "dabba0", 2.098, 0.520, 0, 0, LauBlattWeisskopfFactor::Charm );
365  positve = new LauResonanceInfo( "dabba+", 2.098, 0.520, 0, 1, LauBlattWeisskopfFactor::Charm );
366  negatve = positve->createChargeConjugate();
367  resInfo_.push_back( neutral );
368  resInfo_.push_back( positve );
369  resInfo_.push_back( negatve );
370 
371  // excited charm states name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
372  // D*
373  neutral = new LauResonanceInfo( "D*0", 2.00685, 0.0021, 1, 0, LauBlattWeisskopfFactor::Charm );
374  positve = new LauResonanceInfo( "D*+", 2.01026, 83.4e-6, 1, 1, LauBlattWeisskopfFactor::Charm );
375  negatve = positve->createChargeConjugate();
376  resInfo_.push_back( neutral );
377  resInfo_.push_back( positve );
378  resInfo_.push_back( negatve );
379  // D*_0
380  neutral = new LauResonanceInfo( "D*0_0", 2.300, 0.274, 0, 0, LauBlattWeisskopfFactor::Charm );
381  positve = new LauResonanceInfo( "D*+_0", 2.349, 0.221, 0, 1, LauBlattWeisskopfFactor::Charm );
382  negatve = positve->createChargeConjugate();
383  resInfo_.push_back( neutral );
384  resInfo_.push_back( positve );
385  resInfo_.push_back( negatve );
386  // D*_2
387  //AVERAGE--neutral = new LauResonanceInfo("D*0_2", 2.4618, 0.049, 2, 0 );
388  neutral = new LauResonanceInfo( "D*0_2", 2.4607, 0.0475, 2, 0, LauBlattWeisskopfFactor::Charm );
389  positve = new LauResonanceInfo( "D*+_2", 2.4654, 0.0467, 2, 1, LauBlattWeisskopfFactor::Charm );
390  negatve = positve->createChargeConjugate();
391  resInfo_.push_back( neutral );
392  resInfo_.push_back( positve );
393  resInfo_.push_back( negatve );
394  // D1(2420)
395  neutral =
396  new LauResonanceInfo( "D0_1(2420)", 2.4208, 0.0317, 1, 0, LauBlattWeisskopfFactor::Charm );
397  positve = new LauResonanceInfo( "D+_1(2420)", 2.4232, 0.025, 1, 1, LauBlattWeisskopfFactor::Charm );
398  negatve = positve->createChargeConjugate();
399  resInfo_.push_back( neutral );
400  resInfo_.push_back( positve );
401  resInfo_.push_back( negatve );
402  // D(2600)
403  //OLD--neutral = new LauResonanceInfo("D0(2600)", 2.6087, 0.093, 0, 0 );
404  //OLD--positve = new LauResonanceInfo("D+(2600)", 2.6213, 0.093, 0, 1 );
405  neutral = new LauResonanceInfo( "D0(2600)", 2.623, 0.139, 0, 0, LauBlattWeisskopfFactor::Charm );
406  positve = new LauResonanceInfo( "D+(2600)", 2.623, 0.139, 0, 1, LauBlattWeisskopfFactor::Charm );
407  negatve = positve->createChargeConjugate();
408  resInfo_.push_back( neutral );
409  resInfo_.push_back( positve );
410  resInfo_.push_back( negatve );
411  // D(2680)
412  neutral =
413  new LauResonanceInfo( "D*0_1(2680)", 2.6811, 0.1867, 1, 0, LauBlattWeisskopfFactor::Charm );
414  resInfo_.push_back( neutral );
415  // D(2760)
416  //OLD-- neutral = new LauResonanceInfo("D0(2760)", 2.7633, 0.061, 1, 0 );
417  //OLD-- positve = new LauResonanceInfo("D+(2760)", 2.7697, 0.061, 1, 1 );
418  neutral = new LauResonanceInfo( "D0(2760)", 2.761, 0.063, 1, 0, LauBlattWeisskopfFactor::Charm );
419  positve = new LauResonanceInfo( "D+(2760)", 2.761, 0.063, 1, 1, LauBlattWeisskopfFactor::Charm );
420  negatve = positve->createChargeConjugate();
421  resInfo_.push_back( neutral );
422  resInfo_.push_back( positve );
423  resInfo_.push_back( negatve );
424  neutral =
425  new LauResonanceInfo( "D*0_3(2760)", 2.7755, 0.0953, 3, 0, LauBlattWeisskopfFactor::Charm );
426  resInfo_.push_back( neutral );
427  // D(2900)
428  neutral = new LauResonanceInfo( "D0(3000)", 3.214, 0.186, 0, 0, LauBlattWeisskopfFactor::Charm );
429  resInfo_.push_back( neutral );
430  // D(3400)
431  neutral = new LauResonanceInfo( "D0(3400)", 3.4, 0.15, 0, 0, LauBlattWeisskopfFactor::Charm );
432  resInfo_.push_back( neutral );
433 
434  // excited strange charm name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
435  // Ds*
436  positve =
437  new LauResonanceInfo( "Ds*+", 2.1121, 0.0019, 1, 1, LauBlattWeisskopfFactor::StrangeCharm );
438  negatve = positve->createChargeConjugate();
439  resInfo_.push_back( positve );
440  resInfo_.push_back( negatve );
441  // Ds0*(2317)
442  positve = new LauResonanceInfo( "Ds*+_0(2317)",
443  2.3178,
444  0.0038,
445  0,
446  1,
447  LauBlattWeisskopfFactor::StrangeCharm );
448  negatve = positve->createChargeConjugate();
449  resInfo_.push_back( positve );
450  resInfo_.push_back( negatve );
451  // Ds2*(2573)
452  positve = new LauResonanceInfo( "Ds*+_2(2573)",
453  2.5691,
454  0.0169,
455  2,
456  1,
457  LauBlattWeisskopfFactor::StrangeCharm );
458  negatve = positve->createChargeConjugate();
459  resInfo_.push_back( positve );
460  resInfo_.push_back( negatve );
461  // Ds1*(2700)
462  positve = new LauResonanceInfo( "Ds*+_1(2700)",
463  2.7083,
464  0.120,
465  1,
466  1,
467  LauBlattWeisskopfFactor::StrangeCharm );
468  negatve = positve->createChargeConjugate();
469  resInfo_.push_back( positve );
470  resInfo_.push_back( negatve );
471  // Ds1*(2860)
472  positve = new LauResonanceInfo( "Ds*+_1(2860)",
473  2.859,
474  0.159,
475  1,
476  1,
477  LauBlattWeisskopfFactor::StrangeCharm );
478  negatve = positve->createChargeConjugate();
479  resInfo_.push_back( positve );
480  resInfo_.push_back( negatve );
481  // Ds3*(2860)
482  positve = new LauResonanceInfo( "Ds*+_3(2860)",
483  2.860,
484  0.053,
485  3,
486  1,
487  LauBlattWeisskopfFactor::StrangeCharm );
488  negatve = positve->createChargeConjugate();
489  resInfo_.push_back( positve );
490  resInfo_.push_back( negatve );
491 
492  // excited bottom states name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
493  // B*
494  neutral = new LauResonanceInfo( "B*0", 5.3247, 0.00, 1, 0, LauBlattWeisskopfFactor::Beauty, 6.0 );
495  positve = new LauResonanceInfo( "B*+", 5.3247, 0.00, 1, 1, LauBlattWeisskopfFactor::Beauty, 6.0 );
496  negatve = positve->createChargeConjugate();
497  resInfo_.push_back( neutral );
498  resInfo_.push_back( positve );
499  resInfo_.push_back( negatve );
500 
501  // excited strange bottom name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
502  // Bs*
503  neutral =
504  new LauResonanceInfo( "Bs*0", 5.4154, 0.00, 1, 0, LauBlattWeisskopfFactor::StrangeBeauty, 6.0 );
505  resInfo_.push_back( neutral );
506 
507  // nonresonant models name, mass, width, spin, charge, BW category, BW radius parameter (defaults to 4.0)
508  // Phase-space nonresonant model
509  neutral = new LauResonanceInfo( "NonReson", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
510  resInfo_.push_back( neutral );
511  // Theory-based nonresonant model
512  neutral = new LauResonanceInfo( "NRModel", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
513  resInfo_.push_back( neutral );
514  // Belle nonresonant models
515  neutral = new LauResonanceInfo( "BelleSymNR", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
516  resInfo_.push_back( neutral );
517  neutral = new LauResonanceInfo( "BelleNR", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
518  resInfo_.push_back( neutral );
519  positve = new LauResonanceInfo( "BelleNR+", 0.0, 0.0, 0, 1, LauBlattWeisskopfFactor::Light );
520  negatve = positve->createChargeConjugate();
521  resInfo_.push_back( positve );
522  resInfo_.push_back( negatve );
523  neutral = new LauResonanceInfo( "BelleNR_Swave", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
524  resInfo_.push_back( neutral );
525  positve = new LauResonanceInfo( "BelleNR_Swave+", 0.0, 0.0, 0, 1, LauBlattWeisskopfFactor::Light );
526  negatve = positve->createChargeConjugate();
527  resInfo_.push_back( positve );
528  resInfo_.push_back( negatve );
529  neutral = new LauResonanceInfo( "BelleNR_Pwave", 0.0, 0.0, 1, 0, LauBlattWeisskopfFactor::Light );
530  resInfo_.push_back( neutral );
531  positve = new LauResonanceInfo( "BelleNR_Pwave+", 0.0, 0.0, 1, 1, LauBlattWeisskopfFactor::Light );
532  negatve = positve->createChargeConjugate();
533  resInfo_.push_back( positve );
534  resInfo_.push_back( negatve );
535  neutral = new LauResonanceInfo( "BelleNR_Dwave", 0.0, 0.0, 2, 0, LauBlattWeisskopfFactor::Light );
536  resInfo_.push_back( neutral );
537  positve = new LauResonanceInfo( "BelleNR_Dwave+", 0.0, 0.0, 2, 1, LauBlattWeisskopfFactor::Light );
538  negatve = positve->createChargeConjugate();
539  resInfo_.push_back( positve );
540  resInfo_.push_back( negatve );
541  neutral = new LauResonanceInfo( "BelleNR_Fwave", 0.0, 0.0, 3, 0, LauBlattWeisskopfFactor::Light );
542  resInfo_.push_back( neutral );
543  positve = new LauResonanceInfo( "BelleNR_Fwave+", 0.0, 0.0, 3, 1, LauBlattWeisskopfFactor::Light );
544  negatve = positve->createChargeConjugate();
545  resInfo_.push_back( positve );
546  resInfo_.push_back( negatve );
547  // Taylor expansion nonresonant model
548  neutral = new LauResonanceInfo( "NRTaylor", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
549  resInfo_.push_back( neutral );
550  // Polynomial nonresonant models
551  neutral = new LauResonanceInfo( "PolNR_S0", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
552  resInfo_.push_back( neutral );
553  neutral = new LauResonanceInfo( "PolNR_S1", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
554  resInfo_.push_back( neutral );
555  neutral = new LauResonanceInfo( "PolNR_S2", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
556  resInfo_.push_back( neutral );
557  neutral = new LauResonanceInfo( "PolNR_P0", 0.0, 0.0, 1, 0, LauBlattWeisskopfFactor::Light );
558  resInfo_.push_back( neutral );
559  neutral = new LauResonanceInfo( "PolNR_P1", 0.0, 0.0, 1, 0, LauBlattWeisskopfFactor::Light );
560  resInfo_.push_back( neutral );
561  neutral = new LauResonanceInfo( "PolNR_P2", 0.0, 0.0, 1, 0, LauBlattWeisskopfFactor::Light );
562  resInfo_.push_back( neutral );
563 
564  // Fake resonances for S-Wave splines
565  neutral = new LauResonanceInfo( "Spline_S0", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
566  resInfo_.push_back( neutral );
567  neutral = new LauResonanceInfo( "Spline_S0_Bar", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
568  resInfo_.push_back( neutral );
569 
570  // Polar Form Factor nonresonant model
571  neutral = new LauResonanceInfo( "PolarFFSymNR", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
572  resInfo_.push_back( neutral );
573  neutral = new LauResonanceInfo( "PolarFFNR", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
574  resInfo_.push_back( neutral );
575 
576  // PiPi-KK Inelastic Scattering
577  neutral = new LauResonanceInfo( "Rescattering", 0.0, 0.0, 0, 0, LauBlattWeisskopfFactor::Light );
578  resInfo_.push_back( neutral );
579 
580  nResDefMax_ = resInfo_.size();
581 }
582 
584 {
585  // Check whether any BW factors have been created and bail out if so
586  if ( ! bwIndepFactors_.empty() ) {
587  std::cerr << "ERROR in LauResonanceMaker::setBWType : some barrier factors have already been created - cannot change the barrier type now!"
588  << std::endl;
589  return;
590  }
591  for ( BWFactorCategoryMap::const_iterator iter = bwFactors_.begin(); iter != bwFactors_.end();
592  ++iter ) {
593  if ( iter->second.bwFactor_ != 0 ) {
594  std::cerr << "ERROR in LauResonanceMaker::setBWType : some barrier factors have already been created - cannot change the barrier type now!"
595  << std::endl;
596  return;
597  }
598  }
599 
600  bwBarrierType_ = bwType;
601 }
602 
604 {
605  // Check whether any BW factors have been created and bail out if so
606  if ( ! bwIndepFactors_.empty() ) {
607  std::cerr << "ERROR in LauResonanceMaker::setBWBachelorRestFrame : some barrier factors have already been created - cannot change the rest frame now!"
608  << std::endl;
609  return;
610  }
611  for ( BWFactorCategoryMap::const_iterator iter = bwFactors_.begin(); iter != bwFactors_.end();
612  ++iter ) {
613  if ( iter->second.bwFactor_ != 0 ) {
614  std::cerr << "ERROR in LauResonanceMaker::setBWBachelorRestFrame : some barrier factors have already been created - cannot change the rest frame now!"
615  << std::endl;
616  return;
617  }
618  }
619 
620  bwRestFrame_ = restFrame;
621 }
622 
624 {
625  if ( summaryPrinted_ ) {
626  std::cerr << "ERROR in LauResonanceMaker::setSpinFormalism : cannot redefine the spin formalism after creating one or more resonances"
627  << std::endl;
628  return;
629  }
630  spinFormalism_ = spinType;
631 }
632 
635  const Double_t bwRadius )
636 {
637  if ( bwCategory == LauBlattWeisskopfFactor::Default ||
638  bwCategory == LauBlattWeisskopfFactor::Indep ) {
639  std::cerr << "WARNING in LauResonanceMaker::setDefaultBWRadius : cannot set radius values for Default or Indep categories"
640  << std::endl;
641  return;
642  }
643 
644  // Check if we have an information object for this category
645  BWFactorCategoryMap::iterator factor_iter = bwFactors_.find( bwCategory );
646  if ( factor_iter != bwFactors_.end() ) {
647  // If so, we can set the value in the information object
648  BlattWeisskopfCategoryInfo& categoryInfo = factor_iter->second;
649  categoryInfo.defaultRadius_ = bwRadius;
650 
651  // Then we can check if a LauBlattWeisskopfFactor object has been created for this category
652  LauBlattWeisskopfFactor* bwFactor = categoryInfo.bwFactor_;
653  if ( bwFactor != 0 ) {
654  // If it has then we can also set its radius value directly
655  LauParameter* radius = bwFactor->getRadiusParameter();
656  radius->value( bwRadius );
657  radius->initValue( bwRadius );
658  radius->genValue( bwRadius );
659  }
660  } else {
661  // If not then we just store the value to be used later
662  BlattWeisskopfCategoryInfo& categoryInfo = bwFactors_[bwCategory];
663  categoryInfo.bwFactor_ = 0;
664  categoryInfo.defaultRadius_ = bwRadius;
665  categoryInfo.radiusFixed_ = kTRUE;
666  }
667 }
668 
670  const Bool_t fixRadius )
671 {
672  if ( bwCategory == LauBlattWeisskopfFactor::Default ||
673  bwCategory == LauBlattWeisskopfFactor::Indep ) {
674  std::cerr << "WARNING in LauResonanceMaker::fixBWRadius : cannot fix/float radius values for Default or Indep categories"
675  << std::endl;
676  return;
677  }
678 
679  // Check if we have an information object for this category
680  BWFactorCategoryMap::iterator factor_iter = bwFactors_.find( bwCategory );
681  if ( factor_iter != bwFactors_.end() ) {
682  // If so, we can set the value in the information object
683  BlattWeisskopfCategoryInfo& categoryInfo = factor_iter->second;
684  categoryInfo.radiusFixed_ = fixRadius;
685 
686  // Then we can check if a LauBlattWeisskopfFactor object has been created for this category
687  LauBlattWeisskopfFactor* bwFactor = categoryInfo.bwFactor_;
688  if ( bwFactor != 0 ) {
689  // If it has then we can also fix/float its radius value directly
690  LauParameter* radius = bwFactor->getRadiusParameter();
691  radius->fixed( fixRadius );
692  }
693  } else {
694  // If not then we just store the value to be used later
695  BlattWeisskopfCategoryInfo& categoryInfo = bwFactors_[bwCategory];
696  categoryInfo.bwFactor_ = 0;
697  categoryInfo.defaultRadius_ = -1.0;
698  categoryInfo.radiusFixed_ = fixRadius;
699  }
700 }
701 
703  Int_t resSpin,
705 {
706  LauBlattWeisskopfFactor* bwFactor( 0 );
707 
708  // Look up the category in the category information map
709  BWFactorCategoryMap::iterator factor_iter = bwFactors_.find( LauBlattWeisskopfFactor::Parent );
710 
711  if ( factor_iter == bwFactors_.end() ) {
712  // If the category is currently undefined we need to create it
713  bwFactor = new LauBlattWeisskopfFactor( resSpin,
714  4.0,
716  bwRestFrame_,
717  LauBlattWeisskopfFactor::Parent );
718  std::cerr << "WARNING in LauResonanceMaker::getParentBWFactor : Default radius 4.0 set for Blatt-Weisskopf factor category: Parent"
719  << std::endl;
720 
721  BlattWeisskopfCategoryInfo& categoryInfo = bwFactors_[LauBlattWeisskopfFactor::Parent];
722  categoryInfo.bwFactor_ = bwFactor;
723  categoryInfo.defaultRadius_ = bwFactor->getRadiusParameter()->value();
724  categoryInfo.radiusFixed_ = kTRUE;
725  } else {
726  // If it exists, we can check if the factor object has been created
727  BlattWeisskopfCategoryInfo& categoryInfo = factor_iter->second;
728 
729  if ( categoryInfo.bwFactor_ != 0 ) {
730  // If so, simply clone it
731  bwFactor = categoryInfo.bwFactor_->createClone( resSpin, barrierType );
732  } else {
733  // Otherwise we need to create it, using the default value if it has been set
734  if ( categoryInfo.defaultRadius_ >= 0.0 ) {
735  bwFactor = new LauBlattWeisskopfFactor( resSpin,
736  categoryInfo.defaultRadius_,
738  bwRestFrame_,
739  LauBlattWeisskopfFactor::Parent );
740  } else {
741  bwFactor = new LauBlattWeisskopfFactor( resSpin,
742  4.0,
744  bwRestFrame_,
745  LauBlattWeisskopfFactor::Parent );
746  std::cerr << "WARNING in LauResonanceMaker::getParentBWFactor : Default radius 4.0 set for Blatt-Weisskopf factor category: Parent"
747  << std::endl;
748  }
749  categoryInfo.bwFactor_ = bwFactor;
750 
751  // Set whether the radius should be fixed/floated
752  LauParameter* radius = bwFactor->getRadiusParameter();
753  radius->fixed( categoryInfo.radiusFixed_ );
754  }
755  }
756 
757  return bwFactor;
758 }
759 
762  const LauResonanceInfo* resInfo )
763 {
764  LauBlattWeisskopfFactor* bwFactor( 0 );
765 
766  // If this is an independent factor, create it and add it to the list of independent factors, then return it
767  if ( bwCategory == LauBlattWeisskopfFactor::Indep ) {
768  bwFactor = new LauBlattWeisskopfFactor( *resInfo, bwBarrierType_, bwRestFrame_, bwCategory );
769  bwIndepFactors_.push_back( bwFactor );
770  return bwFactor;
771  }
772 
773  // Otherwise, look up the category in the category information map
774  BWFactorCategoryMap::iterator factor_iter = bwFactors_.find( bwCategory );
775 
776  if ( factor_iter == bwFactors_.end() ) {
777  // If the category is currently undefined we need to create it
778  bwFactor = new LauBlattWeisskopfFactor( *resInfo, bwBarrierType_, bwRestFrame_, bwCategory );
779 
780  BlattWeisskopfCategoryInfo& categoryInfo = bwFactors_[bwCategory];
781  categoryInfo.bwFactor_ = bwFactor;
782  categoryInfo.defaultRadius_ = bwFactor->getRadiusParameter()->value();
783  categoryInfo.radiusFixed_ = kTRUE;
784  } else {
785  // If it exists, we can check if the factor object has been created
786  BlattWeisskopfCategoryInfo& categoryInfo = factor_iter->second;
787 
788  if ( categoryInfo.bwFactor_ != 0 ) {
789  // If so, simply clone it
790  const UInt_t resSpin = resInfo->getSpin();
791  bwFactor = categoryInfo.bwFactor_->createClone( resSpin,
792  categoryInfo.bwFactor_->getBarrierType() );
793  } else {
794  // Otherwise we need to create it, using the default value if it has been set
795  if ( categoryInfo.defaultRadius_ >= 0.0 ) {
796  bwFactor = new LauBlattWeisskopfFactor( *resInfo,
797  categoryInfo.defaultRadius_,
799  bwRestFrame_,
800  bwCategory );
801  } else {
802  bwFactor =
803  new LauBlattWeisskopfFactor( *resInfo, bwBarrierType_, bwRestFrame_, bwCategory );
804  }
805  categoryInfo.bwFactor_ = bwFactor;
806 
807  // Set whether the radius should be fixed/floated
808  LauParameter* radius = bwFactor->getRadiusParameter();
809  radius->fixed( categoryInfo.radiusFixed_ );
810  }
811  }
812 
813  return bwFactor;
814 }
815 
817  const LauDaughters* daughters,
818  const TString& resName,
819  const Int_t resPairAmpInt,
822 {
823  // Routine to return the appropriate LauAbsResonance object given the resonance
824  // name (resName), which daughter is the bachelor track (resPairAmpInt = 1,2 or 3),
825  // and the resonance type ("BW" = Breit-Wigner, "Flatte" = Flatte distribution).
826 
827  // If this is the first resonance we are making, first print a summary of the formalism
828  if ( ! summaryPrinted_ ) {
829  std::cout << "INFO in LauResonanceMaker::getResonance : Freezing amplitude formalism:"
830  << std::endl;
831  switch ( spinFormalism_ ) {
833  std::cout << " : Spin factors use Zemach spin tensors, with bachelor momentum in resonance rest frame"
834  << std::endl;
835  break;
837  std::cout << " : Spin factors use Zemach spin tensors, with bachelor momentum in parent rest frame"
838  << std::endl;
839  break;
841  std::cout << " : Spin factors use Covariant spin tensors"
842  << std::endl;
843  break;
845  std::cout << " : Spin factors are just Legendre polynomials"
846  << std::endl;
847  break;
848  }
849  switch ( bwBarrierType_ ) {
851  std::cout << " : Blatt-Weisskopf barrier factors are the 'non-primed' form"
852  << std::endl;
853  break;
855  std::cout << " : Blatt-Weisskopf barrier factors are the 'primed' form"
856  << std::endl;
857  break;
859  std::cout << " : Blatt-Weisskopf barrier factors are the exponential form"
860  << std::endl;
861  break;
862  }
863  switch ( bwRestFrame_ ) {
865  std::cout << " : Blatt-Weisskopf barrier factors use bachelor momentum in parent rest frame"
866  << std::endl;
867  break;
869  std::cout << " : Blatt-Weisskopf barrier factors use bachelor momentum in resonance rest frame"
870  << std::endl;
871  break;
873  std::cout << " : Blatt-Weisskopf barrier factors use covariant expression"
874  << std::endl;
875  break;
876  }
877 
878  summaryPrinted_ = kTRUE;
879  }
880 
881  // Loop over all possible resonance states we have defined in
882  // createResonanceVector() until we get a match with the name of the resonance
883 
884  LauResonanceInfo* resInfo( 0 );
885  for ( std::vector<LauResonanceInfo*>::const_iterator iter = resInfo_.begin();
886  iter != resInfo_.end();
887  ++iter ) {
888 
889  if ( resName == ( *iter )->getName() ) {
890  // We have recognised the resonance name.
891  std::cout << "INFO in LauResonanceMaker::getResonance : Creating resonance: " << resName
892  << std::endl;
893 
894  resInfo = ( *iter );
895 
896  // stop looping
897  break;
898  }
899  }
900 
901  // if we couldn't find the right resonance then we should return a null pointer
902  if ( resInfo == 0 ) {
903  std::cout << "ERROR in LauResonanceMaker::getResonance : Unable to locate resonance info for: "
904  << resName << std::endl;
905  return 0;
906  }
907 
908  // Now construct the resonance using the specified type
909  LauAbsResonance* theResonance( 0 );
910  switch ( resType ) {
911 
912  case LauAbsResonance::BW :
913  // Simple non-relativistic Breit-Wigner
914  std::cout << " : Using simple Breit-Wigner lineshape. "
915  << std::endl;
916  theResonance = new LauBreitWignerRes( resInfo, resPairAmpInt, daughters );
917  break;
918 
919  case LauAbsResonance::RelBW : {
920  // Relativistic Breit-Wigner with Blatt-Weisskopf factors.
921  std::cout << " : Using relativistic Breit-Wigner lineshape. "
922  << std::endl;
923  theResonance = new LauRelBreitWignerRes( resInfo, resPairAmpInt, daughters );
925  LauBlattWeisskopfFactor::Parent;
926  LauBlattWeisskopfFactor::BlattWeisskopfCategory resCategory = bwCategory;
927  if ( bwCategory == LauBlattWeisskopfFactor::Default ) {
928  resCategory = resInfo->getBWCategory();
929  }
930  LauBlattWeisskopfFactor* resBWFactor = this->getBWFactor( resCategory, resInfo );
931  LauBlattWeisskopfFactor* parBWFactor = this->getBWFactor( parCategory, resInfo );
932  theResonance->setBarrierRadii( resBWFactor, parBWFactor );
933  break;
934  }
935 
936  case LauAbsResonance::GS : {
937  // Gounaris-Sakurai function to try and model the rho(770) better
938  std::cout << " : Using Gounaris-Sakurai lineshape. "
939  << std::endl;
940  theResonance = new LauGounarisSakuraiRes( resInfo, resPairAmpInt, daughters );
942  LauBlattWeisskopfFactor::Parent;
943  LauBlattWeisskopfFactor::BlattWeisskopfCategory resCategory = bwCategory;
944  if ( bwCategory == LauBlattWeisskopfFactor::Default ) {
945  resCategory = resInfo->getBWCategory();
946  }
947  LauBlattWeisskopfFactor* resBWFactor = this->getBWFactor( resCategory, resInfo );
948  LauBlattWeisskopfFactor* parBWFactor = this->getBWFactor( parCategory, resInfo );
949  theResonance->setBarrierRadii( resBWFactor, parBWFactor );
950  break;
951  }
952 
954  // Flatte distribution - coupled channel Breit-Wigner
955  std::cout << " : Using Flatte lineshape. "
956  << std::endl;
957  theResonance = new LauFlatteRes( resInfo, resPairAmpInt, daughters );
958  break;
959 
961  // Sigma model - should only be used for the pi-pi system
962  std::cout << " : Using Sigma lineshape. "
963  << std::endl;
964  theResonance = new LauSigmaRes( resInfo, resPairAmpInt, daughters );
965  break;
966 
968  // Kappa model - should only be used for the K-pi system
969  std::cout << " : Using Kappa lineshape. "
970  << std::endl;
971  theResonance = new LauKappaRes( resInfo, resPairAmpInt, daughters );
972  break;
973 
975  // Dabba model - should only be used for the D-pi system
976  std::cout << " : Using Dabba lineshape. "
977  << std::endl;
978  theResonance = new LauDabbaRes( resInfo, resPairAmpInt, daughters );
979  break;
980 
981  case LauAbsResonance::LASS :
982  // LASS function to try and model the K-pi S-wave better
983  std::cout << " : Using LASS lineshape. "
984  << std::endl;
985  theResonance = new LauLASSRes( resInfo, resPairAmpInt, daughters );
986  break;
987 
989  // LASS function to try and model the K-pi S-wave better
990  std::cout << " : Using LASS lineshape (resonant part only). "
991  << std::endl;
992  theResonance = new LauLASSBWRes( resInfo, resPairAmpInt, daughters );
993  break;
994 
996  // LASS function to try and model the K-pi S-wave better
997  std::cout << " : Using LASS lineshape (nonresonant part only). "
998  << std::endl;
999  theResonance = new LauLASSNRRes( resInfo, resPairAmpInt, daughters );
1000  break;
1001 
1003  // EFKLLM form-factor description of the K-pi S-wave
1004  std::cout << " : Using EFKLLM lineshape. "
1005  << std::endl;
1006  theResonance = new LauEFKLLMRes( resInfo, resPairAmpInt, daughters );
1007  break;
1008 
1010  // K-matrix description
1011  std::cerr << "ERROR in LauResonanceMaker::getResonance : K-matrix type specified, which should be separately handled."
1012  << std::endl;
1013  break;
1014 
1016  // uniform NR amplitude - arguments are there to preserve the interface
1017  std::cout << " : Using uniform NR lineshape. "
1018  << std::endl;
1019  // we override resPairAmpInt here
1020  theResonance = new LauFlatNR( resInfo, 0, daughters );
1021  break;
1022 
1024  // NR amplitude model - arguments are there to preserve the interface
1025  std::cout << " : Using NR-model lineshape. "
1026  << std::endl;
1027  // we override resPairAmpInt here
1028  theResonance = new LauNRAmplitude( resInfo, 0, daughters );
1029  break;
1030 
1033  // Belle NR amplitude model - arguments are there to preserve the interface
1034  std::cout << " : Using Belle NR lineshape. "
1035  << std::endl;
1036  theResonance = new LauBelleNR( resInfo, resType, resPairAmpInt, daughters );
1037  break;
1038 
1042  // Belle NR amplitude model - arguments are there to preserve the interface
1043  std::cout << " : Using Belle symmetric NR lineshape. "
1044  << std::endl;
1045  theResonance = new LauBelleSymNR( resInfo, resType, resPairAmpInt, daughters );
1046  break;
1047 
1048  case LauAbsResonance::PolNR :
1049  // Polynomial NR amplitude model - arguments are there to preserve the interface
1050  std::cout << " : Using polynomial NR lineshape. "
1051  << std::endl;
1052  theResonance = new LauPolNR( resInfo, resPairAmpInt, daughters );
1053  break;
1054 
1055  case LauAbsResonance::Pole :
1056  // Scalar pole model
1057  std::cout << " : Using the scalar Pole lineshape.. "
1058  << std::endl;
1059  theResonance = new LauPoleRes( resInfo, resPairAmpInt, daughters );
1060  break;
1061 
1063  // Polar Form Factor NR amplitude model - arguments are there to preserve the interface
1064  std::cout << " : Using Polar FormFactor NR lineshape.. "
1065  << std::endl;
1066  theResonance = new LauPolarFormFactorNR( resInfo, resType, resPairAmpInt, daughters );
1067  break;
1068 
1071  // Polar Form Factor NR amplitude model - arguments are there to preserve the interface
1072  std::cout << " : Using Polar FormFactor Symetric NR lineshape. "
1073  << std::endl;
1074  theResonance = new LauPolarFormFactorSymNR( resInfo, resType, resPairAmpInt, daughters );
1075  break;
1076 
1079  // KKPiPi Inelastic Scattering amplitude - arguments are there to preserve the interface
1080  std::cout << " : KKPiPi Inelastic Scattering amplitude lineshape. "
1081  << std::endl;
1082  theResonance = new LauRescatteringRes( resInfo, resType, resPairAmpInt, daughters );
1083  break;
1084 
1086  // KKPiPi Inelastic Scattering amplitude - arguments are there to preserve the interface
1087  std::cout << " : KKPiPi Inelastic Scattering amplitude lineshape. "
1088  << std::endl;
1089  theResonance = new LauRescattering2Res( resInfo, resPairAmpInt, daughters );
1090  break;
1091 
1093  // Model independent partial wave
1094  std::cout << " : Using model independent partial wave lineshape (magnitude and phase). "
1095  << std::endl;
1096  theResonance = new LauModIndPartWaveMagPhase( resInfo, resPairAmpInt, daughters );
1097  break;
1098 
1100  // Model independent partial wave
1101  std::cout << " : Using model independent partial wave lineshape (real and imaginary part). "
1102  << std::endl;
1103  theResonance = new LauModIndPartWaveRealImag( resInfo, resPairAmpInt, daughters );
1104  break;
1105 
1107  // Incoherent Gaussian
1108  std::cout << " : Using incoherent Gaussian lineshape. "
1109  << std::endl;
1110  theResonance = new LauGaussIncohRes( resInfo, resPairAmpInt, daughters );
1111  break;
1112 
1117  // Rho-omega mass mixing model
1118  std::cout << " : Using rho-omega mass mixing lineshape. "
1119  << std::endl;
1120  theResonance = new LauRhoOmegaMix( resInfo, resType, resPairAmpInt, daughters );
1122  LauBlattWeisskopfFactor::Parent;
1123  LauBlattWeisskopfFactor::BlattWeisskopfCategory resCategory = bwCategory;
1124  if ( bwCategory == LauBlattWeisskopfFactor::Default ) {
1125  resCategory = resInfo->getBWCategory();
1126  }
1127  LauBlattWeisskopfFactor* resBWFactor = this->getBWFactor( resCategory, resInfo );
1128  LauBlattWeisskopfFactor* parBWFactor = this->getBWFactor( parCategory, resInfo );
1129  theResonance->setBarrierRadii( resBWFactor, parBWFactor );
1130  break;
1131  }
1132 
1133  // Set the spin formalism choice
1134  theResonance->setSpinType( spinFormalism_ );
1135 
1136  return theResonance;
1137 }
1138 
1139 Int_t LauResonanceMaker::resTypeInt( const TString& name ) const
1140 {
1141  // Internal function that returns the resonance integer, specified by the
1142  // order of the resonance vector defined in createResonanceVector(),
1143  // for a given resonance name.
1144  Int_t resTypInt( -99 );
1145  Int_t i( 0 );
1146 
1147  for ( std::vector<LauResonanceInfo*>::const_iterator iter = resInfo_.begin();
1148  iter != resInfo_.end();
1149  ++iter ) {
1150 
1151  if ( name.BeginsWith( ( *iter )->getName(), TString::kExact ) == kTRUE ) {
1152  // We have recognised the resonance from those that are available
1153  resTypInt = i;
1154  return resTypInt;
1155  }
1156  ++i;
1157  }
1158 
1159  return resTypInt;
1160 }
1161 
1162 void LauResonanceMaker::printAll( std::ostream& stream ) const
1163 {
1164  for ( std::vector<LauResonanceInfo*>::const_iterator iter = resInfo_.begin();
1165  iter != resInfo_.end();
1166  ++iter ) {
1167  stream << ( **iter ) << std::endl;
1168  }
1169 }
1170 
1171 LauResonanceInfo* LauResonanceMaker::getResInfo( const TString& resName ) const
1172 {
1173  LauResonanceInfo* resInfo( 0 );
1174  for ( std::vector<LauResonanceInfo*>::const_iterator iter = resInfo_.begin();
1175  iter != resInfo_.end();
1176  ++iter ) {
1177 
1178  if ( resName == ( *iter )->getName() ) {
1179  // We have recognised the resonance name.
1180  resInfo = ( *iter );
1181  // stop looping
1182  break;
1183  }
1184  }
1185 
1186  return resInfo;
1187 }
Class for defining the terms of Babar nonresonant model.
Definition: LauPolNR.hh:46
File containing declaration of LauSigmaRes class.
static LauResonanceMaker * resonanceMaker_
The singleton instance.
File containing declaration of LauResonanceInfo class.
static LauResonanceMaker & get()
Get the factory instance.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
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 NR amplitude model.
File containing declaration of LauBreitWignerRes class.
LauBlattWeisskopfFactor::BlattWeisskopfCategory getBWCategory() const
Retrieve the BW category of the resonant particle.
Int_t resTypeInt(const TString &name) const
Retrieve the integer index for the specified resonance.
File containing declaration of LauRescatteringRes class.
Class for defining an alternative rescattering model.
Double_t value() const
The value of the parameter.
void setDefaultBWRadius(const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory, const Double_t bwRadius)
Set the BW radius for the given category.
Class for defininf the Gounaris-Sakurai resonance model.
UInt_t nResDefMax_
The number of known resonances.
LauBlattWeisskopfFactor * getParentBWFactor(Int_t newSpin, LauBlattWeisskopfFactor::BarrierType barrierType)
Retrieve parent Blatt-Weisskopf factor (for use by K-matrix pole/SVP which doesn't have a ‘resInfo’)
Class for defining the relativistic Breit-Wigner resonance model.
Class for defining the rho-omega resonance mixing model.
void setBWBachelorRestFrame(const LauBlattWeisskopfFactor::RestFrame restFrame)
Set the rest frame in which the bachelor momentum should be calculated (for all BW categories)
File containing declaration of LauPolNR class.
LauBlattWeisskopfFactor * createClone(const UInt_t newSpin, const BarrierType newBarrierType)
Method to create a new factor with cloned radius parameter.
File containing declaration of LauBelleNR class.
Class for defining a nonresonant form factor model.
File containing declaration of LauRelBreitWignerRes class.
Bool_t summaryPrinted_
Boolean flag to control printing a summary of the formalism to be used when the first resonance is cr...
File containing declaration of LauAbsResonance class.
LauBlattWeisskopfFactor::BarrierType bwBarrierType_
The type of the Blatt-Weisskopf barrier to use for all resonances.
LauResonanceInfo * createSharedParameterRecord(const TString &name)
Create another record that will share parameters with this one.
void setSpinFormalism(const LauAbsResonance::LauSpinType spinType)
Set the spin formalism to be used for all resonances.
void createResonanceVector()
Create the list of known resonances.
Class for defining a model independent partial wave component where the amplitudes are parameterised ...
File containing declaration of LauPolarFormFactorSymNR class.
std::vector< LauResonanceInfo * > resInfo_
The known resonances.
void setBWType(const LauBlattWeisskopfFactor::BarrierType bwType)
Set the type of BW factor (for all categories)
File containing declaration of LauFlatNR class.
File containing declaration of LauModIndPartWaveMagPhase class.
File containing declaration of LauNRAmplitude class.
void fixBWRadius(const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory, const Bool_t fixRadius)
Fix or release the Blatt-Weisskopf barrier radius for the given category.
Class for defining the simple Breit-Wigner resonance model.
Bool_t radiusFixed_
Whether or not the radius value for this category should be fixed in the fit.
File containing declaration of LauModIndPartWaveRealImag class.
Class for defining the non resonant part of the LASS model.
Definition: LauLASSNRRes.hh:44
void setBarrierRadii(LauBlattWeisskopfFactor *resFactor, LauBlattWeisskopfFactor *parFactor)
Set the form factor model and parameters.
LauAbsResonance::LauSpinType spinFormalism_
The spin formalism that should be used for all resonances.
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.
File containing declaration of LauKappaRes class.
File containing declaration of LauLASSRes class.
Data structure to store information on a given Blatt-Weisskopf category.
UInt_t getSpin() const
Retrieve the spin of the resonant particle.
File containing declaration of LauDaughters class.
Class for defining a model independent partial wave component where the amplitudes are parameterised ...
BarrierType
Define the allowed types of barrier factors.
LauResonanceMaker()
Constructor.
Double_t defaultRadius_
The default value for the radius in this category.
Class for defining the resonant part of the LASS model.
Definition: LauLASSBWRes.hh:44
LauResonanceInfo * createChargeConjugate()
Create the charge conjugate particle info record.
File containing declaration of LauRescattering2Res class.
Class for defining an incoherent resonance with a Gaussian mass dependence.
Class for defining the Belle nonresonant model.
Definition: LauBelleNR.hh:47
Bool_t fixed() const
Check whether the parameter is fixed or floated.
Class for defining the symmetric Belle Non Resonant model.
File containing declaration of LauBelleSymNR class.
Class for defining the properties of a resonant particle.
Class for defining a uniform nonresonant amplitude.
Definition: LauFlatNR.hh:43
File containing declaration of LauLASSBWRes class.
File containing declaration of LauFlatteRes class.
Class for defining the Kappa resonance model.
Definition: LauKappaRes.hh:45
const TString & name() const
The parameter name.
File containing declaration of LauGaussIncohRes class.
Class for defining the Sigma resonance model.
Definition: LauSigmaRes.hh:44
const LauParameter * getRadiusParameter() const
Retrieve the radius parameter.
Singleton factory class for creating resonances.
File containing declaration of LauRhoOmegaMix class.
File containing declaration of LauPolarFormFactorNR class.
Class for defining the EFKLLM K-pi S-wave model.
Definition: LauEFKLLMRes.hh:49
BlattWeisskopfCategory
Define resonance categories that will share common barrier factor radii.
BarrierType getBarrierType() const
Retrieve the barrier type.
LauSpinType
Define the allowed spin formalisms.
Class for defining the Flatte resonance model.
Definition: LauFlatteRes.hh:44
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc....
LauBlattWeisskopfFactor * getBWFactor(const LauBlattWeisskopfFactor::BlattWeisskopfCategory bwCategory, const LauResonanceInfo *resInfo)
Retrieve Blatt-Weisskopf factor for the given category.
File containing declaration of LauGounarisSakuraiRes class.
Class for defining the rescattering model.
Class for defining the Dabba resonance model.
Definition: LauDabbaRes.hh:44
virtual ~LauResonanceMaker()
Destructor.
Double_t genValue() const
The value generated for the parameter.
void setSpinType(const LauSpinType spinType)
Set the spin formalism to be used.
Class for defining a pole-type resonance model.
Definition: LauPoleRes.hh:45
File containing declaration of LauPoleRes class.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
LauResonanceInfo * getResInfo(const TString &resName) const
Get the information for the given resonance name.
File containing declaration of LauEFKLLMRes class.
BWFactorCategoryMap bwFactors_
The Blatt-Weisskopf factor objects (and related information) for each category.
std::vector< LauBlattWeisskopfFactor * > bwIndepFactors_
The Blatt-Weisskopf factor objects for resonances in the independent category.
RestFrame
Define the rest frame in which the momentum should be calculated (only relevant for bachelor)
Double_t initValue() const
The initial value of the parameter.
Class for defining a nonresonant form factor model.
LauResonanceModel
Define the allowed resonance types.
File containing declaration of LauLASSNRRes class.
File containing declaration of LauResonanceMaker class.
LauBlattWeisskopfFactor * bwFactor_
The BW factor object.
LauBlattWeisskopfFactor::RestFrame bwRestFrame_
The rest frame in which the bachelor momentum used in the Blatt-Weisskopf factors should be calculate...
Class for defining the LASS resonance model.
Definition: LauLASSRes.hh:44
Class that implements the Blatt-Weisskopf barrier factor.