rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

BESIII_2019_I1691850

Analysis of $J/\psi$ decays to $\Lambda^0\bar\Lambda^0$
Experiment: BESIII (BEPC)
Inspire ID: 1691850
Status: VALIDATED
Authors:
  • Peter Richardson
References:
  • Nature Phys. 15 (2019) 631-634
Beams: e- e+
Beam energies: (1.6, 1.6) GeV
Run details:
  • e+e- > J/psi

Analysis of the angular distribution of the baryons, and decay products, produced in $e^+e^-\to J/\psi \to \Lambda^0\bar\Lambda^0$. Gives information about the decay and is useful for testing correlations in hadron decays. N.B. the moment data is not corrected and should only be used qualatively.

Source code: BESIII_2019_I1691850.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/Beam.hh"
  4#include "Rivet/Projections/FinalState.hh"
  5#include "Rivet/Projections/UnstableParticles.hh"
  6
  7namespace Rivet {
  8
  9
 10  /// @brief JPsi > Lambda, Lambdabar
 11  class BESIII_2019_I1691850 : public Analysis {
 12  public:
 13
 14    /// Constructor
 15    RIVET_DEFAULT_ANALYSIS_CTOR(BESIII_2019_I1691850);
 16
 17    /// @name Analysis methods
 18    //@{
 19
 20    /// Book histograms and initialise projections before the run
 21    void init() {
 22
 23      // Initialise and register projections
 24      declare(Beam(), "Beams");
 25      declare(UnstableParticles(), "UFS");
 26      declare(FinalState(), "FS");
 27
 28      // Book histograms
 29      book(_h_T1_p, "T1_p",20,-1.,1.);
 30      book(_h_T2_p, "T2_p",20,-1.,1.);
 31      book(_h_T3_p, "T3_p",20,-1.,1.);
 32      book(_h_T4_p, "T4_p",20,-1.,1.);
 33      book(_h_T5_p, "T5_p",20,-1.,1.);
 34           		       
 35      book(_h_T1_n, "T1_n",20,-1.,1.);
 36      book(_h_T2_n, "T2_n",20,-1.,1.);
 37      book(_h_T3_n, "T3_n",20,-1.,1.);
 38      book(_h_T4_n, "T4_n",20,-1.,1.);
 39      book(_h_T5_n, "T5_n",20,-1.,1.);
 40      
 41      book(_h_cThetaL,"cThetaL",20,-1.,1.);
 42      
 43      book(_h_mu_p, 2,1,1);
 44      book(_h_mu_n, 2,1,2);
 45      book(_wsum_p,"TMP/wsum_p");
 46      book(_wsum_n,"TMP/wsum_n");
 47    }
 48
 49    void findChildren(const Particle & p,map<long,int> & nRes, int &ncount) {
 50      for( const Particle &child : p.children()) {
 51	if(child.children().empty()) {
 52	  nRes[child.pid()]-=1;
 53	  --ncount;
 54	}
 55	else
 56	  findChildren(child,nRes,ncount);
 57      }
 58    }
 59
 60    /// Perform the per-event analysis
 61    void analyze(const Event& event) {
 62      // get the axis, direction of incoming electron
 63      const ParticlePair& beams = apply<Beam>(event, "Beams").beams();
 64      Vector3 axis;
 65      if(beams.first.pid()>0)
 66	axis = beams.first .momentum().p3().unit();
 67      else
 68	axis = beams.second.momentum().p3().unit();
 69      // types of final state particles
 70      const FinalState& fs = apply<FinalState>(event, "FS");
 71      map<long,int> nCount;
 72      int ntotal(0);
 73      for (const Particle& p :  fs.particles()) {
 74	nCount[p.pid()] += 1;
 75	++ntotal;
 76      }
 77      // loop over lambda0 baryons
 78      const UnstableParticles & ufs = apply<UnstableParticles>(event, "UFS");
 79      Particle Lambda,LamBar;
 80      bool matched(false);
 81      for (const Particle& p :  ufs.particles(Cuts::abspid==3122)) {
 82       	if(p.children().empty()) continue;
 83       	map<long,int> nRes=nCount;
 84       	int ncount = ntotal;
 85       	findChildren(p,nRes,ncount);
 86       	matched=false;
 87       	// check for antiparticle
 88      	for (const Particle& p2 :  ufs.particles(Cuts::pid==-p.pid())) {
 89      	  if(p2.children().empty()) continue;
 90      	  map<long,int> nRes2=nRes;
 91      	  int ncount2 = ncount;
 92      	  findChildren(p2,nRes2,ncount2);
 93      	  if(ncount2==0) {
 94      	    matched = true;
 95      	    for(auto const & val : nRes2) {
 96      	      if(val.second!=0) {
 97      		matched = false;
 98      		break;
 99      	      }
100      	    }
101      	    // fond baryon and antibaryon
102      	    if(matched) {
103	      if(p.pid()>0) {
104		Lambda = p;
105		LamBar = p2;
106	      }
107	      else {
108		Lambda = p2;
109		LamBar = p;
110	      }	
111       	      break;
112       	    }
113       	  }
114       	}
115      	if(matched) break;
116      }
117      if(!matched) vetoEvent;
118      Particle proton;
119      matched = false;
120      for (const Particle & p : Lambda.children()) {
121	if(p.pid()==2212) {
122	  matched=true;
123	  proton=p;
124	}
125      }
126      if(!matched) vetoEvent;
127      Particle baryon;
128      int mode(-1);
129      for (const Particle & p : LamBar.children()) {
130	if(p.pid()==-2212) {
131	  baryon=p;
132	  mode=0;
133	}
134	else if(p.pid()==-2112) {
135	  baryon=p;
136	  mode=1;
137	}
138      }
139      if(mode<0) vetoEvent;
140      // boost to the Lambda rest frame
141      LorentzTransform boost1 = LorentzTransform::mkFrameTransformFromBeta(Lambda.momentum().betaVec());
142      Vector3 e1z = Lambda.momentum().p3().unit();
143      Vector3 e1y = e1z.cross(axis).unit();
144      Vector3 e1x = e1y.cross(e1z).unit();
145      Vector3 axis1 = boost1.transform(proton.momentum()).p3().unit();
146      double n1x(e1x.dot(axis1)),n1y(e1y.dot(axis1)),n1z(e1z.dot(axis1));
147      // boost to the Lambda bar
148      LorentzTransform boost2 = LorentzTransform::mkFrameTransformFromBeta(LamBar.momentum().betaVec());
149      Vector3 axis2 = boost2.transform(baryon.momentum()).p3().unit();
150      double n2x(e1x.dot(axis2)),n2y(e1y.dot(axis2)),n2z(e1z.dot(axis2));
151      double cosL = axis.dot(Lambda.momentum().p3().unit());
152      double sinL = sqrt(1.-sqr(cosL));
153      double T1 = sqr(sinL)*n1x*n2x+sqr(cosL)*n1z*n2z;
154      double T2 = -sinL*cosL*(n1x*n2z+n1z*n2x);
155      double T3 = -sinL*cosL*n1y;
156      double T4 = -sinL*cosL*n2y;
157      double T5 = n1z*n2z-sqr(sinL)*n1y*n2y;
158      double mu = n1y-n2y;
159      if(mode==0) {
160	_h_T1_p->fill(cosL,T1);
161	_h_T2_p->fill(cosL,T2);
162	_h_T3_p->fill(cosL,T3);
163	_h_T4_p->fill(cosL,T4);
164	_h_T5_p->fill(cosL,T5);
165	_h_mu_p->fill(cosL,mu);
166	_wsum_p->fill();
167      }
168      else {
169	_h_T1_n->fill(cosL,T1);
170	_h_T2_n->fill(cosL,T2);
171	_h_T3_n->fill(cosL,T3);
172	_h_T4_n->fill(cosL,T4);
173	_h_T5_n->fill(cosL,T5);
174	_h_mu_n->fill(cosL,mu);
175	_wsum_n->fill();
176      }
177      _h_cThetaL->fill(cosL);
178    }
179    
180    pair<double,pair<double,double> > calcAlpha0(Histo1DPtr hist) {
181      if(hist->numEntries()==0.) return make_pair(0.,make_pair(0.,0.));
182      double d = 3./(pow(hist->xMax(),3)-pow(hist->xMin(),3));
183      double c = 3.*(hist->xMax()-hist->xMin())/(pow(hist->xMax(),3)-pow(hist->xMin(),3));
184      double sum1(0.),sum2(0.),sum3(0.),sum4(0.),sum5(0.);
185      for (auto bin : hist->bins() ) {
186       	double Oi = bin.area();
187	if(Oi==0.) continue;
188	double a =  d*(bin.xMax() - bin.xMin());
189	double b = d/3.*(pow(bin.xMax(),3) - pow(bin.xMin(),3));
190       	double Ei = bin.areaErr();
191	sum1 +=   a*Oi/sqr(Ei);
192	sum2 +=   b*Oi/sqr(Ei);
193	sum3 += sqr(a)/sqr(Ei);
194	sum4 += sqr(b)/sqr(Ei);
195	sum5 +=    a*b/sqr(Ei);
196      }
197      // calculate alpha
198      double alpha = (-c*sum1 + sqr(c)*sum2 + sum3 - c*sum5)/(sum1 - c*sum2 + c*sum4 - sum5);
199      // and error
200      double cc = -pow((sum3 + sqr(c)*sum4 - 2*c*sum5),3);
201      double bb = -2*sqr(sum3 + sqr(c)*sum4 - 2*c*sum5)*(sum1 - c*sum2 + c*sum4 - sum5);
202      double aa =  sqr(sum1 - c*sum2 + c*sum4 - sum5)*(-sum3 - sqr(c)*sum4 + sqr(sum1 - c*sum2 + c*sum4 - sum5) + 2*c*sum5);      
203      double dis = sqr(bb)-4.*aa*cc;
204      if(dis>0.) {
205	dis = sqrt(dis);
206	return make_pair(alpha,make_pair(0.5*(-bb+dis)/aa,-0.5*(-bb-dis)/aa));
207      }
208      else {
209	return make_pair(alpha,make_pair(0.,0.));
210      }
211    }
212    
213    pair<double,double> calcCoeff(unsigned int imode,Histo1DPtr hist) {
214      if(hist->numEntries()==0.) return make_pair(0.,0.);
215      double sum1(0.),sum2(0.);
216      for (auto bin : hist->bins() ) {
217	double Oi = bin.area();
218	if(Oi==0.) continue;
219	double ai(0.),bi(0.);
220	if(imode==0) {
221	  bi = (pow(1.-sqr(bin.xMin()),1.5) - pow(1.-sqr(bin.xMax()),1.5))/3.;
222	}
223	else if(imode>=2 && imode<=4) {
224	  bi = ( pow(bin.xMin(),3)*( -5. + 3.*sqr(bin.xMin()))  +
225		 pow(bin.xMax(),3)*(  5. - 3.*sqr(bin.xMax())))/15.;
226	}
227	else
228	  assert(false);
229	double Ei = bin.areaErr();
230	sum1 += sqr(bi/Ei);
231	sum2 += bi/sqr(Ei)*(Oi-ai);
232      }
233      return make_pair(sum2/sum1,sqrt(1./sum1));
234    }
235
236    /// Normalise histograms etc., after the run
237    void finalize() {
238      normalize(_h_cThetaL);
239      scale(_h_T1_p,1./ *_wsum_p);
240      scale(_h_T2_p,1./ *_wsum_p);
241      scale(_h_T3_p,1./ *_wsum_p);
242      scale(_h_T4_p,1./ *_wsum_p);
243      scale(_h_T5_p,1./ *_wsum_p);
244      scale(_h_mu_p,0.04/ *_wsum_p);
245      
246      scale(_h_T1_n,1./ *_wsum_n);
247      scale(_h_T2_n,1./ *_wsum_n);
248      scale(_h_T3_n,1./ *_wsum_n);
249      scale(_h_T4_n,1./ *_wsum_n);
250      scale(_h_T5_n,1./ *_wsum_n);
251      scale(_h_mu_n,0.04/ *_wsum_n);
252
253      // calculate alpha0
254      pair<double,pair<double,double> > alpha0 = calcAlpha0(_h_cThetaL);
255      Scatter2DPtr _h_alpha0;
256      book(_h_alpha0,1,1,1);
257      _h_alpha0->addPoint(0.5, alpha0.first, make_pair(0.5,0.5),
258			  make_pair(alpha0.second.first,alpha0.second.second) );
259      double s2 = -1. + sqr(alpha0.first);
260      double s3 = 3 + alpha0.first;
261      double s1 = sqr(s3);
262      // alpha- and alpha+ from proton data
263      pair<double,double> c_T2_p = calcCoeff(2,_h_T2_p);
264      pair<double,double> c_T3_p = calcCoeff(3,_h_T3_p);
265      pair<double,double> c_T4_p = calcCoeff(4,_h_T4_p);
266      double s4 = sqr(c_T2_p.first);
267      double s5 = sqr(c_T3_p.first);
268      double s6 = sqr(c_T4_p.first);
269      double disc = s1*s5*s6*(-9.*s2*s4 + 4.*s1*s5*s6);
270      if(disc>=0.) {
271	disc = sqrt(disc);
272	double aM = sqrt(-1./s2/s6*(2.*s1*s5*s6+disc));
273	double aP = c_T4_p.first/c_T3_p.first*aM;
274	double aM_M = (2*(alpha0.first*c_T4_p.first*alpha0.second.first + c_T4_p.second*s2)*(disc + 2*s1*s5*s6)
275		       - c_T4_p.first*s2*(4*s3*c_T3_p.first*c_T4_p.first*(c_T3_p.first*c_T4_p.first*alpha0.second.first +s3*c_T4_p.first*c_T3_p.second +s3*c_T3_p.first*c_T4_p.second) +
276				(disc*(- 9*s2*s3*c_T2_p.first*c_T3_p.first*c_T4_p.first* c_T2_p.second
277				       + 9*((1 -  alpha0.first*(3 + 2*alpha0.first))* c_T3_p.first*c_T4_p.first*alpha0.second.first -  s2*s3*c_T4_p.first*c_T3_p.second
278					     - s2*s3*c_T3_p.first*c_T4_p.second)* s4
279				       + 8*(c_T3_p.first*c_T4_p.first*alpha0.second.first +  s3*c_T4_p.first*c_T3_p.second +  s3*c_T3_p.first*c_T4_p.second)* s1*s5*s6))
280				/(4*pow(3 + alpha0.first,3)*pow(c_T3_p.first,3)*pow(c_T4_p.first,3) -9*s2*s3*c_T3_p.first*c_T4_p.first*s4)))/
281	  (2.*pow(c_T4_p.first,3)*pow(s2,2)*sqrt(-((disc + 2*s1*s5*s6)/(s2*s6))));
282	double aM_P = (2*(alpha0.first*c_T4_p.first*alpha0.second.second + c_T4_p.second*s2)*(disc + 2*s1*s5*s6)
283		       - c_T4_p.first*s2*(4*s3*c_T3_p.first*c_T4_p.first*(c_T3_p.first*c_T4_p.first*alpha0.second.second +s3*c_T4_p.first*c_T3_p.second +s3*c_T3_p.first*c_T4_p.second) +
284				(disc*(- 9*s2*s3*c_T2_p.first*c_T3_p.first*c_T4_p.first* c_T2_p.second
285				       + 9*((1 -  alpha0.first*(3 + 2*alpha0.first))* c_T3_p.first*c_T4_p.first*alpha0.second.second -  s2*s3*c_T4_p.first*c_T3_p.second
286					     - s2*s3*c_T3_p.first*c_T4_p.second)* s4
287				       + 8*(c_T3_p.first*c_T4_p.first*alpha0.second.second +  s3*c_T4_p.first*c_T3_p.second +  s3*c_T3_p.first*c_T4_p.second)* s1*s5*s6))
288				/(4*pow(3 + alpha0.first,3)*pow(c_T3_p.first,3)*pow(c_T4_p.first,3) -9*s2*s3*c_T3_p.first*c_T4_p.first*s4)))/
289	  (2.*pow(c_T4_p.first,3)*pow(s2,2)*sqrt(-((disc + 2*s1*s5*s6)/(s2*s6))));
290	double aP_M = (c_T4_p.first*sqrt(-((disc + 2*s1*s5*s6)/   (s2*s6)))*
291		       (-2*c_T3_p.second -  (2*alpha0.first*c_T3_p.first*alpha0.second.first)/s2 +  (c_T3_p.first*(4*s3*c_T3_p.first*c_T4_p.first*(c_T3_p.first*c_T4_p.first*alpha0.second.first +  s3*c_T4_p.first*c_T3_p.second +  s3*c_T3_p.first*c_T4_p.second)
292							    + (disc*(-9*s2*s3*c_T2_p.first*c_T3_p.first*c_T4_p.first* c_T2_p.second
293								      +  9*((1 -  alpha0.first*(3 + 2*alpha0.first))* c_T3_p.first*c_T4_p.first*alpha0.second.first -  s2*s3*c_T4_p.first*c_T3_p.second
294									     -  s2*s3*c_T3_p.first*c_T4_p.second)* s4 +
295								      8*(c_T3_p.first*c_T4_p.first*alpha0.second.first +  s3*c_T4_p.first*c_T3_p.second +  s3*c_T3_p.first*c_T4_p.second)* s1*s5*s6))/
296							    (4* pow(3 + alpha0.first,3)* pow(c_T3_p.first,3)* pow(c_T4_p.first,3) -  9*s2*s3*c_T3_p.first*c_T4_p.first*s4)))/
297			(disc + 2*s1*s5*s6)))/(2.*pow(c_T3_p.first,2));
298	double aP_P = (c_T4_p.first*sqrt(-((disc + 2*s1*s5*s6)/   (s2*s6)))*
299		       (-2*c_T3_p.second -  (2*alpha0.first*c_T3_p.first*alpha0.second.second)/s2 +  (c_T3_p.first*(4*s3*c_T3_p.first*c_T4_p.first*(c_T3_p.first*c_T4_p.first*alpha0.second.second +  s3*c_T4_p.first*c_T3_p.second +  s3*c_T3_p.first*c_T4_p.second)
300							    + (disc*(-9*s2*s3*c_T2_p.first*c_T3_p.first*c_T4_p.first* c_T2_p.second
301								      +  9*((1 -  alpha0.first*(3 + 2*alpha0.first))* c_T3_p.first*c_T4_p.first*alpha0.second.second -  s2*s3*c_T4_p.first*c_T3_p.second
302									     -  s2*s3*c_T3_p.first*c_T4_p.second)* s4 +
303								      8*(c_T3_p.first*c_T4_p.first*alpha0.second.second +  s3*c_T4_p.first*c_T3_p.second +  s3*c_T3_p.first*c_T4_p.second)* s1*s5*s6))/
304							    (4* pow(3 + alpha0.first,3)* pow(c_T3_p.first,3)* pow(c_T4_p.first,3) -  9*s2*s3*c_T3_p.first*c_T4_p.first*s4)))/
305			(disc + 2*s1*s5*s6)))/(2.*pow(c_T3_p.first,2));
306	Scatter2DPtr _h_alphaM;
307	book(_h_alphaM,1,3,1);
308	_h_alphaM->addPoint(0.5, aM, make_pair(0.5,0.5),
309			    make_pair(-aM_M , -aM_P ) );
310
311	Scatter2DPtr _h_alphaP;
312	book(_h_alphaP,1,4,1);
313	_h_alphaP->addPoint(0.5, aP, make_pair(0.5,0.5),
314			    make_pair(-aP_M , -aP_P  ) );
315	// now for Delta
316	double sDelta = (-2.*(3. + alpha0.first)*c_T3_p.first)/(aM*sqrt(1 - sqr(alpha0.first)));
317	double cDelta = (-3*(3 + alpha0.first)*c_T2_p.first)/(aM*aP*sqrt(1 - sqr(alpha0.first)));
318
319	double Delta = asin(sDelta);
320	if(cDelta<0.) Delta = M_PI-Delta;
321	double ds_P = (-9*c_T2_p.first*((-1 + alpha0.first)*(1 + alpha0.first)*  (3 + alpha0.first)*c_T3_p.first*c_T4_p.first*c_T2_p.second +  c_T2_p.first*c_T4_p.first*(c_T3_p.first*(alpha0.second.first + 3*alpha0.first*alpha0.second.first) -(-1 + alpha0.first)*(1 + alpha0.first)*(3 + alpha0.first)*c_T3_p.second)
322			      -  (-1 + alpha0.first)*(1 + alpha0.first)*  (3 + alpha0.first)*c_T2_p.first*c_T3_p.first*c_T4_p.second)*disc)/
323	  (pow(1 - pow(alpha0.first,2),1.5)*pow(c_T4_p.first,3)*pow(-((disc + 2*s1*s5*s6)/   (s2*s6)),1.5)*(-9*s2*s4 + 4*s1*s5*s6));
324	double ds_M = (-9*c_T2_p.first*((-1 + alpha0.first)*(1 + alpha0.first)*  (3 + alpha0.first)*c_T3_p.first*c_T4_p.first*c_T2_p.second +  c_T2_p.first*c_T4_p.first*(c_T3_p.first*(alpha0.second.second + 3*alpha0.first*alpha0.second.second) -(-1 + alpha0.first)*(1 + alpha0.first)*(3 + alpha0.first)*c_T3_p.second)
325			      -  (-1 + alpha0.first)*(1 + alpha0.first)*  (3 + alpha0.first)*c_T2_p.first*c_T3_p.first*c_T4_p.second)*disc)/
326	  (pow(1 - pow(alpha0.first,2),1.5)*pow(c_T4_p.first,3)*pow(-((disc + 2*s1*s5*s6)/   (s2*s6)),1.5)*(-9*s2*s4 + 4*s1*s5*s6));
327	ds_P /= sqrt(1.-sqr(sDelta));
328	ds_M /= sqrt(1.-sqr(sDelta));
329	Scatter2DPtr _h_sin;
330	book(_h_sin,1,2,1);
331	_h_sin->addPoint(0.5, Delta/M_PI*180., make_pair(0.5,0.5), make_pair( -ds_P/M_PI*180., -ds_M/M_PI*180. ) );
332      }
333      // alpha 0
334      pair<double,double> c_T2_n = calcCoeff(2,_h_T2_n);
335      pair<double,double> c_T3_n = calcCoeff(3,_h_T3_n);
336      pair<double,double> c_T4_n = calcCoeff(4,_h_T4_n);
337      s4 = sqr(c_T2_n.first);
338      s5 = sqr(c_T3_n.first);
339      s6 = sqr(c_T4_n.first);
340      disc = s1*s5*s6*(-9.*s2*s4 + 4.*s1*s5*s6);
341      if(disc>=0.) {
342	disc = sqrt(disc);
343	double aM = sqrt(-1./s2/s6*(2.*s1*s5*s6+disc));
344	double a0 = c_T4_n.first/c_T3_n.first*aM;
345	double a0_M = (c_T4_n.first*sqrt(-((disc + 2*s1*s5*s6)/   (s2*s6)))*
346		       (-2*c_T3_n.second -  (2*alpha0.first*c_T3_n.first*alpha0.second.first)/s2 +  (c_T3_n.first*(4*s3*c_T3_n.first*c_T4_n.first*(c_T3_n.first*c_T4_n.first*alpha0.second.first +  s3*c_T4_n.first*c_T3_n.second +  s3*c_T3_n.first*c_T4_n.second)
347							    + (disc*(-9*s2*s3*c_T2_n.first*c_T3_n.first*c_T4_n.first* c_T2_n.second
348								      +  9*((1 -  alpha0.first*(3 + 2*alpha0.first))* c_T3_n.first*c_T4_n.first*alpha0.second.first -  s2*s3*c_T4_n.first*c_T3_n.second
349									     -  s2*s3*c_T3_n.first*c_T4_n.second)* s4 +
350								      8*(c_T3_n.first*c_T4_n.first*alpha0.second.first +  s3*c_T4_n.first*c_T3_n.second +  s3*c_T3_n.first*c_T4_n.second)* s1*s5*s6))/
351							    (4* pow(3 + alpha0.first,3)* pow(c_T3_n.first,3)* pow(c_T4_n.first,3) -  9*s2*s3*c_T3_n.first*c_T4_n.first*s4)))/
352			(disc + 2*s1*s5*s6)))/(2.*pow(c_T3_n.first,2));
353	double a0_P = (c_T4_n.first*sqrt(-((disc + 2*s1*s5*s6)/   (s2*s6)))*
354		       (-2*c_T3_n.second -  (2*alpha0.first*c_T3_n.first*alpha0.second.second)/s2 +  (c_T3_n.first*(4*s3*c_T3_n.first*c_T4_n.first*(c_T3_n.first*c_T4_n.first*alpha0.second.second +  s3*c_T4_n.first*c_T3_n.second +  s3*c_T3_n.first*c_T4_n.second)
355							    + (disc*(-9*s2*s3*c_T2_n.first*c_T3_n.first*c_T4_n.first* c_T2_n.second
356								      +  9*((1 -  alpha0.first*(3 + 2*alpha0.first))* c_T3_n.first*c_T4_n.first*alpha0.second.second -  s2*s3*c_T4_n.first*c_T3_n.second
357									     -  s2*s3*c_T3_n.first*c_T4_n.second)* s4 +
358								      8*(c_T3_n.first*c_T4_n.first*alpha0.second.second +  s3*c_T4_n.first*c_T3_n.second +  s3*c_T3_n.first*c_T4_n.second)* s1*s5*s6))/
359							    (4* pow(3 + alpha0.first,3)* pow(c_T3_n.first,3)* pow(c_T4_n.first,3) -  9*s2*s3*c_T3_n.first*c_T4_n.first*s4)))/
360			(disc + 2*s1*s5*s6)))/(2.*pow(c_T3_n.first,2));
361	Scatter2DPtr _h_alpha0;
362	book(_h_alpha0,1,5,1);
363	_h_alpha0->addPoint(0.5, a0, make_pair(0.5,0.5),
364			    make_pair(-a0_M , -a0_P ) );
365      }
366    }
367
368    //@}
369
370    /// @name Histograms
371    //@{
372    Histo1DPtr _h_T1_p,_h_T2_p,_h_T3_p,_h_T4_p,_h_T5_p;
373    Histo1DPtr _h_T1_n,_h_T2_n,_h_T3_n,_h_T4_n,_h_T5_n;
374    Histo1DPtr _h_cThetaL;
375    Histo1DPtr _h_mu_p,_h_mu_n;
376    CounterPtr _wsum_p,_wsum_n;
377    //@}
378
379  };
380
381
382  // The hook for the plugin system
383  RIVET_DECLARE_PLUGIN(BESIII_2019_I1691850);
384
385
386}