EvtGen 2.2.0
Monte Carlo generator of particle decays, in particular the weak decays of heavy flavour particles such as B mesons.
Loading...
Searching...
No Matches
EvtFourBodyPhsp.cpp
Go to the documentation of this file.
1
2/***********************************************************************
3* Copyright 1998-2020 CERN for the benefit of the EvtGen authors *
4* *
5* This file is part of EvtGen. *
6* *
7* EvtGen is free software: you can redistribute it and/or modify *
8* it under the terms of the GNU General Public License as published by *
9* the Free Software Foundation, either version 3 of the License, or *
10* (at your option) any later version. *
11* *
12* EvtGen is distributed in the hope that it will be useful, *
13* but WITHOUT ANY WARRANTY; without even the implied warranty of *
14* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15* GNU General Public License for more details. *
16* *
17* You should have received a copy of the GNU General Public License *
18* along with EvtGen. If not, see <https://www.gnu.org/licenses/>. *
19***********************************************************************/
20
22
23#include "EvtGenBase/EvtKine.hh"
24#include "EvtGenBase/EvtPDL.hh"
28
29#include <cmath>
30
31std::string EvtFourBodyPhsp::getName() const
32{
33 return "FOURBODYPHSP";
34}
35
37{
38 return new EvtFourBodyPhsp;
39}
40
41std::array<double, 4> EvtFourBodyPhsp::phspFactor(
42 const double mM, const double m12, const double m34,
43 std::array<double, 4>& daughters ) const
44{
45 std::array<double, 4> result;
46 result[1] = twoBodyMomentum( mM, m12, m34 );
47 result[2] = twoBodyMomentum( m12, daughters[0], daughters[1] );
48 result[3] = twoBodyMomentum( m34, daughters[2], daughters[3] );
49 result[0] = result[1] * result[2] * result[3];
50
51 return result;
52}
53
55{
56 // Check that we have right number of daughters
57 checkNDaug( 4 );
58
59 // Check whether mother is quasi-stable
60 auto parent = getParentId();
61 double width = EvtPDL::getWidth( parent );
62 if ( width > 1e-6 ) {
63 m_stableMother = false;
64 }
65
66 // Check whether all daughters are stable
67 for ( int i = 0; i < 4; ++i ) {
68 auto daughter = getDaug( i );
69 width = EvtPDL::getWidth( daughter );
70 if ( width > 1e-6 ) {
71 m_stableDaughters = false;
72 m_daughterMasses[i] = EvtPDL::getMinMass( daughter );
73 } else {
74 m_daughterMasses[i] = EvtPDL::getMass( daughter );
75 }
76 }
77
78 // check correct number of arguments
79 checkNArg( 0, 2, 4 );
80 double mass1 = m_daughterMasses[0];
81 double mass2 = m_daughterMasses[1];
82 double mass3 = m_daughterMasses[2];
83 double mass4 = m_daughterMasses[3];
84 double mMother = EvtPDL::getMaxMass( parent );
85 if ( getNArg() > 2 ) {
86 m_m12Min = getArg( 0 );
87 m_m12Max = getArg( 1 );
88 m_m34Min = getArg( 2 );
89 m_m34Max = getArg( 3 );
90 } else {
91 if ( getNArg() > 0 ) {
92 m_m12Min = getArg( 0 );
93 m_m12Max = getArg( 1 );
94 } else {
95 m_m12Min = mass1 + mass2;
96 m_m12Max = mMother - mass3 - mass4;
97 }
98 m_m34Min = mass3 + mass4;
99 m_m34Max = mMother - mass1 - mass2;
100 if ( m_stableDaughters == false || m_stableMother == false ) {
101 m_fixedBoundary = false;
102 }
103 }
104 // Make sure that we have correct boundaries
105 if ( m_m12Min < mass1 + mass2 ) {
106 m_m12Min = mass1 + mass2;
107 }
108 if ( m_m12Max > mMother - mass3 - mass4 ) {
109 m_m12Max = mMother - mass3 - mass4;
110 }
111 if ( m_m34Min < mass3 + mass4 ) {
112 m_m34Min = mass3 + mass4;
113 }
114 if ( m_m34Max > mMother - mass1 - mass2 ) {
115 m_m34Max = mMother - mass1 - mass2;
116 }
117
120 mMother );
121 } else {
123 }
124 // If we have fixed boundary, we can precalculate some variables for
125 // m12 and m34 generation
126 if ( m_fixedBoundary ) {
128 const double m12Diff = m_m12Max - m_m12Min;
129 const double minSum = m_m12Min + m_m34Min;
130 m_trapNorm = ( mMother - m_m34Min ) * m12Diff -
131 0.5 * ( m12Diff * ( m_m12Max + m_m12Min ) );
132 m_trapCoeff1 = mMother - m_m34Min;
133 m_trapCoeff2 = mMother * mMother - 2 * mMother * minSum +
134 minSum * minSum;
135 }
137 m_pentagonSplit = mMother - m_m34Max;
138 const double area1 = ( m_pentagonSplit - m_m12Min ) *
139 ( m_m34Max - m_m34Min );
140 const double pm12Diff = m_m12Max - m_pentagonSplit;
141 const double area2 = 0.5 * pm12Diff *
142 ( mMother + m_m34Max - m_m12Max ) -
143 pm12Diff * m_m34Min;
144 m_pentagonFraction = area1 / ( area1 + area2 );
145 const double m12Diff = m_m12Max - m_pentagonSplit;
146 const double minSum = m_pentagonSplit + m_m34Min;
147 m_trapNorm = ( mMother - m_m34Min ) * m12Diff -
148 0.5 * ( m12Diff * ( m_m12Max + m_pentagonSplit ) );
149 m_trapCoeff1 = mMother - m_m34Min;
150 m_trapCoeff2 = mMother * mMother - 2 * mMother * minSum +
151 minSum * minSum;
152 }
153 }
154}
155
157{
158 double startM12 = m_m12Min + ( m_m12Max - m_m12Min ) / 20.;
159 double startM34 = m_m34Min + ( m_m34Max - m_m34Min ) / 20.;
160 bool contCond = true;
161 int iteration = 0;
162
163 auto parent = getParentId();
164 double mMother = EvtPDL::getMaxMass( parent );
165
166 double funcValue = 0;
167 while ( contCond ) {
168 ++iteration;
169 double currentM12 = startM12;
170 double currentM34 = startM34;
171 funcValue = phspFactor( mMother, currentM12, currentM34,
172 m_daughterMasses )[0];
173 // Maximum along m12
174 double step = ( m_m12Max - m_m12Min ) / 100.;
175 while ( step > 1e-4 ) {
176 double point1 = currentM12 + step;
177 if ( point1 > m_m12Max ) {
178 point1 = m_m12Max;
179 }
180 if ( currentM34 > mMother - point1 ) {
181 point1 = mMother - currentM34;
182 }
183 double point2 = currentM12 - step;
184 if ( point2 < m_m12Min ) {
185 point2 = m_m12Min;
186 }
187 double value1 = phspFactor( mMother, point1, currentM34,
188 m_daughterMasses )[0];
189 double value2 = phspFactor( mMother, point2, currentM34,
190 m_daughterMasses )[0];
191 if ( value1 > funcValue && value1 > value2 ) {
192 currentM12 = point1;
193 funcValue = value1;
194 } else if ( value2 > funcValue ) {
195 currentM12 = point2;
196 funcValue = value2;
197 }
198 step /= 2.;
199 }
200 // Maximum along m34
201 step = ( mMother - currentM12 - m_m34Min ) / 100.;
202 while ( step > 1e-4 ) {
203 double point1 = currentM34 + step;
204 if ( point1 > m_m34Max ) {
205 point1 = m_m34Max;
206 }
207 if ( point1 > mMother - currentM12 ) {
208 point1 = mMother - currentM12;
209 }
210 double point2 = currentM34 - step;
211 if ( point2 < m_m34Min ) {
212 point2 = m_m34Min;
213 }
214 double value1 = phspFactor( mMother, currentM12, point1,
215 m_daughterMasses )[0];
216 double value2 = phspFactor( mMother, currentM12, point2,
217 m_daughterMasses )[0];
218 if ( value1 > funcValue && value1 > value2 ) {
219 currentM34 = point1;
220 funcValue = value1;
221 } else if ( value2 > funcValue ) {
222 currentM34 = point2;
223 funcValue = value2;
224 }
225 step /= 2.;
226 }
227
228 // Check termination condition
229 double m12Diff = currentM12 - startM12;
230 double m34Diff = currentM34 - startM34;
231 double distSq = m12Diff * m12Diff + m34Diff * m34Diff;
232 if ( distSq < 1e-8 || iteration > 50 ) {
233 contCond = false;
234 }
235 startM12 = currentM12;
236 startM34 = currentM34;
237 }
238
239 setProbMax( funcValue * 1.05 );
240}
241
243{
244 parent->makeDaughters( getNDaug(), getDaugs() );
245 bool massTreeStatus = parent->generateMassTree();
246 if ( !massTreeStatus ) {
247 EvtGenReport( EVTGEN_ERROR, "EvtGen" )
248 << "Failed to generate daughters masses in EvtFourBodyPhsp."
249 << std::endl;
250 ::abort();
251 }
252
253 double mMother = parent->mass();
254
255 // Need to check whether boundaries are OK and whether we need to work
256 // out boundary shape
257 double cM12Min, cM12Max;
258 double cM34Min, cM34Max;
260 if ( m_fixedBoundary ) {
261 cM12Min = m_m12Min;
262 cM12Max = m_m12Max;
263 cM34Min = m_m34Min;
264 cM34Max = m_m34Max;
265 cShape = m_boundaryShape;
266 } else {
267 // In this case at least one particle has non-zero width and thus
268 // boundaries and shape of the region can change
269 for ( int i = 0; i < 4; ++i ) {
270 auto daughter = parent->getDaug( i );
271 m_daughterMasses[i] = daughter->mass();
272 }
273 cM12Min = m_m12Min > ( m_daughterMasses[0] + m_daughterMasses[1] )
274 ? m_m12Min
276 cM12Max = m_m12Max <
277 ( mMother - m_daughterMasses[2] - m_daughterMasses[3] )
278 ? m_m12Max
279 : mMother - m_daughterMasses[2] - m_daughterMasses[3];
280 cM34Min = m_m34Min > ( m_daughterMasses[2] + m_daughterMasses[3] )
281 ? m_m34Min
283 cM34Max = m_m34Max <
284 ( mMother - m_daughterMasses[0] - m_daughterMasses[1] )
285 ? m_m34Max
286 : mMother - m_daughterMasses[0] - m_daughterMasses[1];
287 cShape = determineBoundaryShape( cM12Min, cM12Max, cM34Max, mMother );
288 }
289
290 // Generate m12 and m34
291 auto masses = generatePairMasses( cM12Min, cM12Max, cM34Min, cM34Max,
292 mMother, cShape );
293 const double m12 = masses.first;
294 const double m34 = masses.second;
295
296 // calculate probability, it will return array with 4 elements with
297 // probability, q, p1 and p3
298 auto probEval = phspFactor( mMother, m12, m34, m_daughterMasses );
299 setProb( probEval[0] );
300
301 // initialise kinematics
302 const double cosTheta1 = EvtRandom::Flat( -1.0, 1.0 );
303 const double sinTheta1 = std::sqrt( 1 - cosTheta1 * cosTheta1 );
304 const double cosTheta3 = EvtRandom::Flat( -1.0, 1.0 );
305 const double sinTheta3 = std::sqrt( 1 - cosTheta3 * cosTheta3 );
306 const double phi = EvtRandom::Flat( 0., EvtConst::twoPi );
307 // m12 and m34 are put along z-axis, 1 and 2 go to x-z plane and 3-4
308 // plane is rotated by phi compared to 1-2 plane. All momenta are set
309 // in 12 and 34 rest frames and then boosted to parent rest frame
310 const double p1x = probEval[2] * sinTheta1;
311 const double p1z = probEval[2] * cosTheta1;
312 const double p1Sq = probEval[2] * probEval[2];
313 const double en1 = std::sqrt( m_daughterMasses[0] * m_daughterMasses[0] +
314 p1Sq );
315 const double en2 = std::sqrt( m_daughterMasses[1] * m_daughterMasses[1] +
316 p1Sq );
317 const double p3T = probEval[3] * sinTheta3;
318 const double p3x = p3T * std::cos( phi );
319 const double p3y = p3T * std::sin( phi );
320 const double p3z = probEval[3] * cosTheta3;
321 const double p3Sq = probEval[3] * probEval[3];
322 const double en3 = std::sqrt( m_daughterMasses[2] * m_daughterMasses[2] +
323 p3Sq );
324 const double en4 = std::sqrt( m_daughterMasses[3] * m_daughterMasses[3] +
325 p3Sq );
326
327 EvtVector4R mom1( en1, p1x, 0.0, p1z );
328 EvtVector4R mom2( en2, -p1x, 0.0, -p1z );
329 EvtVector4R mom3( en3, p3x, p3y, p3z );
330 EvtVector4R mom4( en4, -p3x, -p3y, -p3z );
331
332 const double qSq = probEval[1] * probEval[1];
333 const double en12 = std::sqrt( m12 * m12 + qSq );
334 const double en34 = std::sqrt( m34 * m34 + qSq );
335 EvtVector4R q12( en12, 0.0, 0.0, probEval[1] );
336 EvtVector4R q34( en34, 0.0, 0.0, -probEval[1] );
337 mom1.applyBoostTo( q12 );
338 mom2.applyBoostTo( q12 );
339 mom3.applyBoostTo( q34 );
340 mom4.applyBoostTo( q34 );
341
342 // As final step, rotate everything randomly in space
343 const double euler1 = EvtRandom::Flat( 0., EvtConst::twoPi );
344 const double euler2 = std::acos( EvtRandom::Flat( -1.0, 1.0 ) );
345 const double euler3 = EvtRandom::Flat( 0., EvtConst::twoPi );
346 mom1.applyRotateEuler( euler1, euler2, euler3 );
347 mom2.applyRotateEuler( euler1, euler2, euler3 );
348 mom3.applyRotateEuler( euler1, euler2, euler3 );
349 mom4.applyRotateEuler( euler1, euler2, euler3 );
350
351 // Set momenta for daughters
352 auto daug = parent->getDaug( 0 );
353 daug->init( daug->getId(), mom1 );
354 daug = parent->getDaug( 1 );
355 daug->init( daug->getId(), mom2 );
356 daug = parent->getDaug( 2 );
357 daug->init( daug->getId(), mom3 );
358 daug = parent->getDaug( 3 );
359 daug->init( daug->getId(), mom4 );
360}
361
363 const double m12Min, const double m12Max, const double m34Max,
364 const double mMother ) const
365{
366 double maxY = mMother - m12Min;
367 const bool corner1 = m34Max < maxY;
368 maxY = mMother - m12Max;
369 const bool corner2 = m34Max < maxY;
370
371 if ( corner1 && corner2 ) {
372 return Shape::rectangle;
373 } else if ( !corner1 && !corner2 ) {
374 return Shape::trapezoid;
375 }
376 return Shape::pentagon;
377}
378
379std::pair<double, double> EvtFourBodyPhsp::generatePairMasses(
380 const double m12Min, const double m12Max, const double m34Min,
381 const double m34Max, const double mMother,
382 const EvtFourBodyPhsp::Shape shape ) const
383{
384 switch ( shape ) {
386 return generateRectangle( m12Min, m12Max, m34Min, m34Max );
387 break;
389 return generateTrapezoid( m12Min, m12Max, m34Min, mMother );
390 break;
392 double split, fraction;
393 if ( m_fixedBoundary ) {
394 split = m_pentagonSplit;
395 fraction = m_pentagonFraction;
396 } else {
397 split = mMother - m34Max;
398 const double area1 = ( split - m12Min ) * ( m34Max - m34Min );
399 const double pm12Diff = m12Max - split;
400 const double area2 = 0.5 * pm12Diff *
401 ( mMother + m34Max - m12Max ) -
402 pm12Diff * m34Min;
403 fraction = area1 / ( area1 + area2 );
404 }
405 if ( EvtRandom::Flat() < fraction ) {
406 return generateRectangle( m12Min, split, m34Min, m34Max );
407 } else {
408 return generateTrapezoid( split, m12Max, m34Min, mMother );
409 }
410 break;
411 default:
412 return std::make_pair( m12Min, m34Min );
413 break;
414 }
415}
416
417std::pair<double, double> EvtFourBodyPhsp::generateRectangle(
418 const double m12Min, const double m12Max, const double m34Min,
419 const double m34Max ) const
420{
421 return std::make_pair( EvtRandom::Flat( m12Min, m12Max ),
422 EvtRandom::Flat( m34Min, m34Max ) );
423}
424
425std::pair<double, double> EvtFourBodyPhsp::generateTrapezoid(
426 const double m12Min, const double m12Max, const double m34Min,
427 const double mMother ) const
428{
429 double norm, coeff1, coeff2;
430 if ( m_fixedBoundary ) {
431 norm = m_trapNorm;
432 coeff1 = m_trapCoeff1;
433 coeff2 = m_trapCoeff2;
434 } else {
435 const double m12Diff = m12Max - m12Min;
436 const double minSum = m12Min + m34Min;
437 norm = ( mMother - m34Min ) * m12Diff -
438 0.5 * ( m12Diff * ( m12Max + m12Min ) );
439 coeff1 = mMother - m34Min;
440 coeff2 = mMother * mMother - 2 * mMother * minSum + minSum * minSum;
441 }
442 const double rnd = EvtRandom::Flat();
443 const double m12 = coeff1 - std::sqrt( -2.0 * rnd * norm + coeff2 );
444 const double m34 = EvtRandom::Flat( m34Min, mMother - m12 );
445 return std::make_pair( m12, m34 );
446}
double twoBodyMomentum(const double M, const double m1, const double m2)
Definition EvtKine.cpp:134
std::ostream & EvtGenReport(EvtGenSeverity severity, const char *facility=nullptr)
Definition EvtReport.cpp:32
@ EVTGEN_ERROR
Definition EvtReport.hh:49
static const double twoPi
Definition EvtConst.hh:27
EvtDecayBase()=default
int getNDaug() const
int getNArg() const
double getArg(unsigned int j)
void setProbMax(double prbmx)
EvtId getParentId() const
EvtId getDaug(int i) const
void checkNDaug(int d1, int d2=-1)
void checkNArg(int a1, int a2=-1, int a3=-1, int a4=-1)
const EvtId * getDaugs() const
void setProb(double prob)
std::pair< double, double > generateTrapezoid(const double m12Min, const double m12Max, const double m34Min, const double mMother) const
std::pair< double, double > generatePairMasses(const double m12Min, const double m12Max, const double m34Min, const double m34Max, const double mMother, const EvtFourBodyPhsp::Shape shape) const
std::array< double, 4 > m_daughterMasses
void decay(EvtParticle *parent) override
void init() override
std::array< double, 4 > phspFactor(const double mM, const double m12, const double m34, std::array< double, 4 > &daughters) const
std::pair< double, double > generateRectangle(const double m12Min, const double m12Max, const double m34Min, const double m34Max) const
std::string getName() const override
EvtDecayBase * clone() const override
void initProbMax() override
Shape determineBoundaryShape(const double m12Min, const double m12Max, const double m34Max, const double mMother) const
static double getWidth(EvtId i)
Definition EvtPDL.cpp:346
static double getMaxMass(EvtId i)
Definition EvtPDL.cpp:331
static double getMass(EvtId i)
Definition EvtPDL.cpp:311
static double getMinMass(EvtId i)
Definition EvtPDL.cpp:336
virtual void init(EvtId part_n, const EvtVector4R &p4)=0
EvtParticle * getDaug(const int i)
double mass() const
bool generateMassTree()
void makeDaughters(size_t ndaug, const EvtId *id)
static double Flat()
Definition EvtRandom.cpp:95
void applyRotateEuler(double alpha, double beta, double gamma)
void applyBoostTo(const EvtVector4R &p4, bool inverse=false)