BOSS 7.1.2
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtDecay.cxx
Go to the documentation of this file.
1//*****************************************************************************
2//
3// EvtDecay.cxx
4//
5// EvtDecay algorithm takes generated event from transient store, and decay
6// particles, including weak decays of its secondary particles.
7// EvtDecay can be used as a TopAlg in conjunction with algorithms Pythia,
8// KKMC or a SingleParticleGun.
9//
10// History:
11// Original LHCb code by Witold Pokorski and Patric Robbe.
12// August 2002: Malte Muller, adopted for ATHENA to work inside algorithm PythiaBModule (only private version within 3.0.0)
13// Sept 2003: James Catmore, adopted for 6.0.3 (not reeased in 6.0.3 ony private version) to work inside PythiaBModule.
14// Nov 2003: M.Smizanska, rewritten a) to be standalone EvtDecay algorithm so it can be combined
15// with any type of Pythia generator b) to match to new LHCb(CDF) code dedicated to LHC, c) to work in 7.3.0.
16// EvtDecay first time released in 7.3.0 20.Nov. 2003
17// Dec 2003: M.Smizanska: define user's input - decay table
18// Feb 2004 J Catmore, addition of random seed facility. TEMPORARY FIX
19// Oct 2005 A. Zhemchugov, adapted for BES3
20// May 2006 K.L He, set spin density matrix for e+e- -> V
21// Sept.2007, R.G.Ping, change to the BesEvtGen
22// Jan. 2008, R.G.Ping, to print McTruth table to the file DecayTopology when only doing generation, not simulation
23// Feb. 2008, R.G.Ping, add an option to only run the BesEvtGen
24// Mar. 2008, R.G. Ping, user options "DecayDecDir" and "PdtTableDir" are changed.
25// Mar. 2008-03, R.G. Ping, dump the user decay card to the bosslog file
26//*****************************************************************************
27
28// Header for this module:-
29#include "EvtDecay.h"
30#include "ReadME.h"
35#include "EvtGen.hh"
40#include "EvtGlobalSet.hh"
44
45// Framework Related Headers:-
46#include "HepMC/GenEvent.h"
47#include "HepMC/GenVertex.h"
48#include "HepMC/GenParticle.h"
50#include "GaudiKernel/SmartDataPtr.h"
53#include "GaudiKernel/MsgStream.h"
54#include "GaudiKernel/ISvcLocator.h"
55#include "GaudiKernel/AlgFactory.h"
56#include "GaudiKernel/DataSvc.h"
57#include "GaudiKernel/SmartDataPtr.h"
58#include "GaudiKernel/IPartPropSvc.h"
59
62#include "CLHEP/Random/Ranlux64Engine.h"
63#include "CLHEP/Random/RandFlat.h"
65
66#include <stdlib.h>
67#include <string.h>
68
69static string mydecayrec;
70static double _amps_max;
71int nchr = 0;
72int nchr_e = 0;
73int nchr_mu= 0;
74int nchr_pi= 0;
75int nchr_k = 0;
76int nchr_p = 0;
77int N_gamma=0;
78int N_gammaFSR = 0;
79int fst[50];
80int EvtDecay::m_runNo=0;
81double EvtConExc::SetMthr=0;
82//////these three are allocate memory for EvtParticle //new modified
86DECLARE_COMPONENT(EvtDecay)
87
88EvtDecay::EvtDecay(const string& name, ISvcLocator* pSvcLocator):Algorithm( name, pSvcLocator )
89{
90
91
92
93 // these can be used to specify alternative locations or filenames
94 // for the EvtGen particle and channel definition files.
95
96 declareProperty("DecayDecDir", m_DecayDec = "");
97 declareProperty("PdtTableDir", m_PdtTable = "");
98 declareProperty("userDecayTableName", userDecFileName = "NONE");
99 declareProperty("DecayTopology", m_DecayTop = ""); // output decay topology to a file specified by user
100 declareProperty("DecayRecord", m_DecayRec = ""); // output decay topology of one particle specified by user to a file
101 declareProperty("ParentParticle", m_ParentPart = "NONE");// Mothor particle name in pdt.table for only running BesEvtGen
102 declareProperty("Multiplicity", m_Ncharge = false); // output ncharge number of an event
103 declareProperty("NutupleFile",m_NtupleFile = false); // output Ntuple file
104 declareProperty("mDIY",_mDIY= false); // mDIY model
105 declareProperty("FDPdata",m_SB3File = ""); // Fit Dalitz Plot (FDP) data file name (*.root)
106 declareProperty("FDPHisT",m_SB3HT = ""); // histogram title of Dalitz plot in data file
107 declareProperty("FDPpart",m_FDPparticle =""); //to assign the FDP parent particle name (string)
108 declareProperty("TruthFile",m_truthFile ="");
109 declareProperty("TruthPart",m_truthPart ="");
110 declareProperty("Psi3SopenCharm",m_Psi4040OpenCharm=false);
111 declareProperty("Psi2openCharm", m_Psi2openCharm=false);
112 declareProperty("SetMthrForConExc",m_SetMthr=0.0);
113 declareProperty("statDecays",m_statDecays=false);
114 declareProperty("TagLundCharmModel", m_tagLundModel=false);
115 m_mystring.clear();
116 declareProperty("FileForTrackGen", m_mystring);
117 //for polarized charmonium production
118 m_polarization.clear();//= diag{1+Pe, 0, 1-Pe} with electron polarization Pe
119 declareProperty("polarization", m_polarization);
120 declareProperty("eBeamPolarization", m_eBeamPolarization=-1);
121 m_cluster0.clear();m_wind0.clear();
122 m_cluster1.clear();m_wind1.clear();
123 m_cluster2.clear();m_wind2.clear();
124 declareProperty("cluster0",m_cluster0);
125 declareProperty("cluster1",m_cluster1);
126 declareProperty("cluster2",m_cluster2);
127 declareProperty("masswin0",m_wind0);
128 declareProperty("masswin1",m_wind1);
129 declareProperty("masswin2",m_wind2);
130 declareProperty("OutputP4",m_outputp4="");
131 declareProperty("ReadMeasuredEcms", m_RdMeasuredEcms = false);
132 declareProperty("beamEnergySpread", m_beamEnergySpread = 0);
133 m_ReadTruth.clear();
134 declareProperty("ReadTruthAtCM",_truthAtCM=false);
135 declareProperty("ReadTruth",m_ReadTruth);
136 //ReadTruth={{"ParentName"},{"i0","i1","i2"},{"j0","j1","j2","j3"}}, where the first part. is Parent->getDaug(i0)->getDaug(i1)->getDaug(i2),
137 //and second particle is Parent ->getDaug(j0)->getDaug(j1)->getDaug(j2)->getDaug(j3);
138 declareProperty("RvalueTag",_RvalueTag=false);
139 declareProperty("PrintFSFor",m_printfs="");
140}
141
142
144
145 MsgStream log(messageService(), name());
146 //system("cat $BESEVTGENROOT/share/phokhara_9.1.fferr>phokhara_9.1.fferr");
147 //system("cat $BESEVTGENROOT/share/phokhara_9.1.ffwarn>phokhara_9.1.ffwarn");
148 system("cat $BESEVTGENROOT/share/phokhara_10.0.fferr>phokhara_10.0.fferr"); //modified by luojiashun 2022.2.21
149 system("cat $BESEVTGENROOT/share/phokhara_10.0.ffwarn>phokhara_10.0.ffwarn"); //modified by luojiashun 2022.2.21
150
151 if(m_truthFile!=""){truth.open(m_truthFile.c_str());}
152 log << MSG::INFO << "EvtDecay initialize" << endreq;
153 static const bool CREATEIFNOTTHERE(true);
154 StatusCode RndmStatus = service("BesRndmGenSvc", p_BesRndmGenSvc, CREATEIFNOTTHERE);
155 if (!RndmStatus.isSuccess() || 0 == p_BesRndmGenSvc)
156 {
157 log << MSG::ERROR << " Could not initialize Random Number Service" << endreq;
158 return RndmStatus;
159 }
160
161 EvtGlobalSet::SV.clear();
162 for(int i=0;i<m_mystring.size();i++){EvtGlobalSet::SV.push_back(m_mystring[i]);}
163
164 EvtGlobalSet::iVV.clear();
165 EvtGlobalSet::dVV.clear();
166 std::vector<std::vector<int> >myivv;
167 std::vector<std::vector<double> >mydvv;
168
169 EvtConExc::SetMthr=m_SetMthr;
170
171 myivv.push_back(m_cluster0);
172 myivv.push_back(m_cluster1);
173 myivv.push_back(m_cluster2);
174
175 mydvv.push_back(m_wind0);
176 mydvv.push_back(m_wind1);
177 mydvv.push_back(m_wind2);
178
179 for(int i=0;i<myivv.size();i++){
180 std::vector<int> theivv;
181 for(int j=0;j<myivv[i].size();j++){
182 theivv.push_back(myivv[i][j]);
183 }
184 if(theivv.size()) EvtGlobalSet::iVV.push_back(theivv);
185 }
186
187 for(int i=0;i<mydvv.size();i++){
188 std::vector<double> thedvv;
189 for(int j=0;j<mydvv[i].size();j++){
190 thedvv.push_back(mydvv[i][j]);
191 }
192 if(thedvv.size()) EvtGlobalSet::dVV.push_back(thedvv);
193 }
194
195
196 if(m_eBeamPolarization>1) {std::cout<<"The e-beam polaziation should be in 0 to 1"<<std::endl;abort();}
197 m_numberEvent=0;
198 first = true;
199 m_ampscalflag = true;
200 // Save the random number seeds in the event
201 CLHEP::HepRandomEngine* engine = p_BesRndmGenSvc->GetEngine("EVTGEN");
202 const long s = engine->getSeed();
203 p_BesRndmGenSvc->setGenseed(s+1);
204
205 m_seeds.clear();
206 m_seeds.push_back(s);
207 totalChannels = 0;
208
209 // open an interface to EvtGen
210
211 if ( m_DecayDec == "") { //use default DECAY.DEC and pdt.table if not defined by user
212 m_DecayDec=getenv("BESEVTGENROOT");
213 m_DecayDec +="/share/DECAY.DEC";
214 }
215
216 if ( m_PdtTable == "") {
217 m_PdtTable =getenv("BESEVTGENROOT");
218 m_PdtTable +="/share/pdt.table";
219 }
220
221 char catcmd[300]; //output user decay cards to log file
222 std::cout<<"===================== user decay card ========================"<<std::endl;
223 strcpy(catcmd, "cat ");
224 strcat(catcmd, userDecFileName.c_str());
225 system(catcmd);
226
229
230 // write decay cards to the root file: pingrg@2009-09-09
231 StatusCode status;
232 status = service("DataInfoSvc",tmpInfoSvc);
233 if (status.isSuccess()) {
234 log << MSG::INFO << "got the DataInfoSvc" << endreq;
235 dataInfoSvc=dynamic_cast<DataInfoSvc *>(tmpInfoSvc);
236 dataInfoSvc->setDecayCard(userDecFileName);
237 } else {
238 log << MSG::ERROR << "could not get the DataInfoSvc" << endreq;
239 return StatusCode::FAILURE;
240 }
241
242
243
244 m_RandomEngine = new EvtBesRandom(engine);
245 m_Gen = new EvtGen( m_DecayDec.c_str(), m_PdtTable.c_str(), m_RandomEngine );
246
247 if(userDecFileName =="NONE") log << MSG::INFO << "EvtDecay User did not define his Decay table EvtGen will use standart" << endreq;
248 if(!(userDecFileName =="NONE")) m_Gen->readUDecay(userDecFileName.c_str());
249
250 if(!(m_DecayTop==" ")) {
251 outfile.open(m_DecayTop.c_str());
252 }
253
254 //for output Ntuple file, pingrg-2009-09-07
255
256
257 if(m_NtupleFile) {
258 NTuplePtr nt1(ntupleSvc(),"MYROOTFILE/Trk");
259 if(nt1) {m_tuple=nt1;}
260 else {
261 m_tuple = ntupleSvc()->book ("MYROOTFILE/Trk", CLID_ColumnWiseTuple, "Generator-trk-Ntuple");
262 if(m_tuple){
263 status = m_tuple->addItem ("TotNumTrk", TotNumTrk, 0,100);
264 status = m_tuple->addIndexedItem ("Trk_index", TotNumTrk, m_Trk_index);
265 status = m_tuple->addIndexedItem ("m_px_trk", TotNumTrk, m_px_trk);
266 status = m_tuple->addIndexedItem ("m_py_trk", TotNumTrk, m_py_trk);
267 status = m_tuple->addIndexedItem ("m_pz_trk", TotNumTrk, m_pz_trk);
268 status = m_tuple->addIndexedItem ("m_en_trk", TotNumTrk, m_en_trk);
269 status = m_tuple->addIndexedItem ("FST", TotNumTrk, m_fst);
270
271 status = m_tuple->addItem ("nchr", m_nchr);
272 status = m_tuple->addItem ("nchr_e", m_nchr_e);
273 status = m_tuple->addItem ("nchr_mu", m_nchr_mu);
274 status = m_tuple->addItem ("nchr_pi", m_nchr_pi);
275 status = m_tuple->addItem ("nchr_k", m_nchr_k);
276 status = m_tuple->addItem ("nchr_p", m_nchr_p);
277 status = m_tuple->addItem ("N_gamma", m_gamma);
278 status = m_tuple->addItem ("N_gammaFSR", m_gammaFSR);
279 status = m_tuple->addItem ("Flag1", m_flag1);
280 } else {
281 log << MSG::ERROR << " Cannot book N-tuple:"<< long(m_tuple) << endmsg;
282 return StatusCode::FAILURE;
283 }
284 }
285
286 NTuplePtr nt2(ntupleSvc(),"MYROOTFILE/mass");
287 if(nt2) {mass_tuple=nt2;}
288 else {
289 mass_tuple = ntupleSvc()->book ("MYROOTFILE/mass", CLID_ColumnWiseTuple, "Generator-mass-Ntuple");
290 if(mass_tuple){
291 status = mass_tuple->addItem ("m12", m_m12);
292 status = mass_tuple->addItem ("m13", m_m13);
293 status = mass_tuple->addItem ("m23", m_m23);
294 status = mass_tuple->addItem ("m1", m_m1);
295 status = mass_tuple->addItem ("m2", m_m2);
296 status = mass_tuple->addItem ("m3", m_m3);
297 status = mass_tuple->addItem ("costheta1", m_cos1);
298 status = mass_tuple->addItem ("costheta2", m_cos2);
299 status = mass_tuple->addItem ("costheta3", m_cos3);
300 status = mass_tuple->addItem ("ichannel", m_ich);
301 } else {
302 log << MSG::ERROR << " Cannot book N-tuple:"<< long(m_tuple) << endmsg;
303 return StatusCode::FAILURE;
304 }
305 }
306
307 NTuplePtr nt3(ntupleSvc(),"MYROOTFILE/massGen");
308 if(nt3) {massgen_tuple=nt3;}
309 else {
310 massgen_tuple = ntupleSvc()->book ("MYROOTFILE/massGen", CLID_ColumnWiseTuple, "Generator-mass-Ntuple");
311 if(massgen_tuple){
312 status = massgen_tuple->addItem ("m12", _m12);
313 status = massgen_tuple->addItem ("m13", _m13);
314 status = massgen_tuple->addItem ("m23", _m23);
315 status = massgen_tuple->addItem ("m1", _m1);
316 status = massgen_tuple->addItem ("m2", _m2);
317 status = massgen_tuple->addItem ("m3", _m3);
318 status = massgen_tuple->addItem ("costheta1", _cos1);
319 status = massgen_tuple->addItem ("costheta2", _cos2);
320 status = massgen_tuple->addItem ("costheta3", _cos3);
321 status = massgen_tuple->addItem ("ichannel", _ich);
322 } else {
323 log << MSG::ERROR << " Cannot book N-tuple:"<< long(m_tuple) << endmsg;
324 return StatusCode::FAILURE;
325 }
326 }
327
328
329 }
330 for(int i=0;i<500;i++){br[i]=0;}
331 if(!(m_SB3File=="" && m_SB3HT=="")) {
332 const char *filename, *histitle;
333 filename=(char*)m_SB3File.c_str();
334 histitle=(char*)m_SB3HT.c_str();
335 SuperBody3decay.setFile(filename);
336 SuperBody3decay.setHTitle(histitle);
337 SuperBody3decay.init();
338 }
339
340 return StatusCode::SUCCESS;
341}
342
343StatusCode EvtDecay::execute()
344{
345 //new Added
346 EvtParticle::_NextLevelDauNum=0; //////////add for Model:FromParent
347
348 MsgStream log(messageService(), name());
349 // log << MSG::INFO << "EvtDecay executing" << endreq;
350 //------------------
351 string key = "GEN_EVENT";
352 // retrieve event from Transient Store (Storegate)
353 SmartDataPtr<McGenEventCol> McEvtColl(eventSvc(), "/Event/Gen");
354
355 m_numberEvent += 1;
356 writeFlag = 0;
357 //std::cout<<"EvtNumber= "<<m_numberEvent<<std::endl;
358 if ( McEvtColl == 0 )
359 {
360 //std::cout<<"McEvtColl == 0 INFO: Using the callBesEvtGen!!"<<std::endl;
361 HepMC::GenEvent* evt=new GenEvent();
362 evt->set_event_number(m_numberEvent);
363
364 //read Ecms from the database
365 SmartDataPtr<Event::EventHeader> eventHeader(eventSvc(),"/Event/EventHeader");
366 int runNo=eventHeader->runNumber();
367 int event=eventHeader->eventNumber();
368 bool newRunFlag=0;
369 if(runNo != 0 && runNo != m_runNo){m_runNo=runNo;newRunFlag = true;}else{newRunFlag=false;}
370 if(m_RdMeasuredEcms&& newRunFlag){// using cms energy of beam read from database
371 runNo = std::abs(runNo);
372 ReadME theME(runNo);
373 if(theME.isRunNoValid()){
374 dbEcms=theME.getEcms();
375 std::cout<<"Read Ecms= "<<dbEcms<<std::endl;
376 if(dbEcms!=0){
377 std::cout<<"INFO: Read the MeasuredEcms"<<std::endl;
378 }
379 else{
380 std::cout<<"ERROR: Can not read the MeasuredEcms, try to turn off the ReadEcmsDB"<<std::endl;
381 return StatusCode::FAILURE;
382 }
383 }
384 }
385 //end of read Ecms
386
387 callBesEvtGen( evt );
388 if(!(m_DecayTop=="")) evt->print(outfile);
389
390 //create a Transient store
391 McGenEventCol *mcColl = new McGenEventCol;
392 McGenEvent* mcEvent = new McGenEvent(evt);
393 mcColl->push_back(mcEvent);
394 StatusCode sc = eventSvc()->registerObject("/Event/Gen",mcColl);
395 if(sc != SUCCESS) return StatusCode::FAILURE;
396
397 } else // (McEvtColl != 0)
398 {
399
400 //std::cout<<"McEvtColl == 1 INFO: Using the callEvtGen!!"<<std::endl;
401 McGenEventCol::iterator mcItr;
402 for( mcItr = McEvtColl->begin(); mcItr != McEvtColl->end(); mcItr++ )
403 {
404 HepMC::GenEvent* hEvt = (*mcItr)->getGenEvt();
405 // MeVToGeV( hEvt );
406
407 callEvtGen( hEvt );
408
409 if(!(m_DecayTop=="")) hEvt->print(outfile);
410 // GeVToMeV( hEvt );
411 // if(!(m_DecayRec=="")) outfile2<<std::endl;
412 if(!(m_DecayRec=="")) std::cout<<std::endl;
413 };
414 }
415 if( m_NtupleFile ){ m_tuple->write();}
416
417 if(_mDIY){
418 std::cout<<"helangs ";
419 for(int i=0;i<m_vangs.size();i++) std::cout<<m_vangs[i]<<" ";
420 std::cout<<std::endl;
421 }
422 return StatusCode::SUCCESS;
423}
424
425// main procedure, loops over all particles in given event,
426// converts them to EvtGenParticles,
427// feeds them to EvtGen,
428// converts back to HepMC particles and puts them in the event.
429
430StatusCode EvtDecay::callEvtGen( HepMC::GenEvent* hepMCevt )
431{
432 MsgStream log(messageService(), name());
433 nchr = 0;
434 nchr_e = 0;
435 nchr_mu= 0;
436 nchr_pi= 0;
437 nchr_k = 0;
438 nchr_p = 0;
439 N_gamma=0;
440 N_gammaFSR = 0;
441 SmartDataPtr<Event::EventHeader> eventHeader(eventSvc(),"/Event/EventHeader");
442
443 for (int i=0;i<13;i++) fst[i]=0;
444 HepMC::GenEvent::particle_const_iterator itp;
445 AllTrk_index=0;
446 for( itp=hepMCevt->particles_begin(); itp!=hepMCevt->particles_end(); ++itp )
447 {
448 // This is the main loop.
449 // We iterate over particles and we look for ones that
450 // should be decayed with EvtGen.
451 //
452 // status == 1 - undecayed Pythia particle (also for particles that are not supposed to decay)
453 // status == 999 - particle already decayed by EvtGen
454 // status == 899 - particle was supposed to decay by EvtGen - but found no model
455 // this may be also intentional - such events are then excluded from
456 // being written into persistent output.
457 // ISTHEP(IHEP):
458 // status code for entry IHEP, with the following meanings:
459 //
460 // = 0 : null entry.
461 // = 1 : an existing entry, which has not decayed or fragmented. This is the main class of entries, which represents the `final state' given by the generator.
462 // = 2 : an entry which has decayed or fragmented and is therefore not appearing in the final state, but is retained for event history information.
463 // = 3 : a documentation line, defined separately from the event history. This could include the two incoming reacting particles, etc.
464 // = 4 - 10 : undefined, but reserved for future standards.
465 // = 11 - 200 : at the disposal of each model builder for constructs specific to his program, but equivalent to a null line in the context of any other program.
466 // = 201 - : at the disposal of users, in particular for event tracking in the detector.
467
468 HepMC::GenParticle* hepMCpart=*itp;
469 int stat = hepMCpart->status();
470
471
472 if ( stat != 1 ) continue;
473
474 loop_to_select_eventsB:
475 int id = hepMCpart->pdg_id();
476 parentPDGcode=id;
477 hepMCpart->set_status(899);
479 string parentName=EvtPDL::name(eid);
480
481 double en =(hepMCpart->momentum()).e();
482 double pex=(hepMCpart->momentum()).px();
483 double pey=(hepMCpart->momentum()).py();
484 double pez=(hepMCpart->momentum()).pz();
485
486 EvtVector4R p_init(en,pex,pey,pez);
487 parentMass=p_init.mass();
488 EvtPDL::reSetMass(eid,parentMass);
489
491
492 // set spin density matrix for e+ e- -> V
493 bool parentFlag=false;
494 if( writeFlag==0 && part->getSpinType() == EvtSpinType::VECTOR) {
495 if(m_polarization.size()>0) {
496 part->setPolarizedSpinDensity(m_polarization[0],m_polarization[1],m_polarization[2]);
497 }else if(m_eBeamPolarization>0){
498 part->setPolarizedSpinDensity(1+ m_eBeamPolarization,0,1- m_eBeamPolarization);
499 } else{
500 part->setVectorSpinDensity();
501 }
502 parentFlag=true;
503 writeFlag++;
504 } else {
505 int id = hepMCpart->pdg_id();
506 if( id == 22 ) {N_gamma ++;fst[0]++;} //fst[0]:photon
507 if( id == -22 ){N_gammaFSR ++;}
508 if( id == 11 ) {fst[1]++;} else if(id == -11){fst[2]++;} // fst[1]: electron ; fst[]: positron
509 else if(id == 13){fst[3]++;} // fst[3]: mu-
510 else if(id == -13){fst[4]++;} // fst[4]: mu+
511 else if(id == 211){fst[5]++;} // fst[5]: pi+
512 else if(id ==-211){fst[6]++;} // fst[6]: pi-
513 else if(id == 321){fst[7]++;} // fst[7]: K+
514 else if(id ==-321){fst[8]++;} // fst[8]: K-
515 else if(id ==2212){fst[9]++;} // fst[9]: pronton
516 else if(id ==-2212){fst[10]++;} // fst[10]: anti-proton
517 else if(id == 2112){fst[11]++;} // fst[11]: nucleon
518 else if(id ==-2112){fst[12]++;} // fst[12]: ant-nucleon
519 if( fabs(id) == 11) {nchr_e++;}
520 if( fabs(id) == 13) {nchr_mu++;}
521 if( fabs(id) == 211) {nchr_pi++;}
522 if( fabs(id) == 321) {nchr_k++;}
523 if( fabs(id) == 2212) {nchr_p++;}
524
525 //std::cout<<"id= "<<id<<" AllTrk_index= "<<AllTrk_index<<std::endl;
526 if( m_NtupleFile==true ){
527 m_Trk_index[AllTrk_index]=id;
528 m_px_trk[AllTrk_index]=pex;
529 m_py_trk[AllTrk_index]=pey;
530 m_pz_trk[AllTrk_index]=pez;
531 m_en_trk[AllTrk_index]=en ;
532
533 AllTrk_index++;
534 Trk_index[AllTrk_index]=0;
535 }
536
537
538 }
539
540 // for SuperBody3decay generator
541 // EvtVector4R pd1,pd2,pd3;
542 EvtVector4R fdp_init;
543 EvtId fdp_ppid;
544 if(m_FDPparticle!=""){
545 findPart(part);
546 fdp_ppid = FDP_id;
547 fdp_init = FDP_init;
548 } else{
549 fdp_ppid = eid;
550 fdp_init = p_init;
551 }
552
553 if(!(m_SB3File=="" || m_SB3HT=="") && parentFlag && first ){
554 SuperBody3decay_make(fdp_ppid, fdp_init );
555 }
556
557 loop_to_select_eventsA:
558 m_Gen->generateDecay(part);
559 if(m_truthFile!="" && m_truthPart!=""){
560 if(EvtPDL::name(part->getId())==m_truthPart ){
561 int ndaug = part->getNDaug();
562 for(int id=0;id<ndaug;id++){
563 EvtParticle* sonp=part->getDaug(id);
564 EvtVector4R son=sonp->getP4();
565 truth<<son.get(0)<<" "<<son.get(1)<<" "<<son.get(2)<<" "<<son.get(3)<<std::endl;
566 }
567 }
568 }
569 double ratio,rdm;
570 if(_mDIY && parentFlag && m_ampscalflag ) _amps_max=CalAmpsMax(part);
571 if(_mDIY && parentFlag) {
572 rdm=EvtRandom::Flat(0.0, 1.0);
573 double amps=CalAmpsMDIY( part );
574 ratio=amps/_amps_max;
575 }
576 if(_mDIY && parentFlag && ratio<=rdm){ goto loop_to_select_eventsA;}
577
578 //std::cout<<"m_vangs.size= "<<m_vangs.size()<<std::endl;
579 //-------- For superbody3------------------------------
580 bool accept;
581 if(m_FDPparticle==""){FDP_part=part;}
582 if(!(m_SB3File=="" || m_SB3HT=="") && parentFlag ){ accept=SuperBody3decay_judge( FDP_part); }
583 if(!(m_SB3File=="" || m_SB3HT=="") && parentFlag && !accept){
584 part->deleteTree();
585 writeFlag=0;
586 goto loop_to_select_eventsB;
587 }else if(!(m_SB3File=="" || m_SB3HT=="") && parentFlag && accept){
588 countChannel(FDP_part);
589 }
590 //-------- for psi4040 open charm inclusive decay
591 if((m_Psi2openCharm || m_Psi4040OpenCharm) && parentFlag ){
592 if(getModel(part) == "OPENCHARM"){
593 std::cout<<"OPENCHARM model need to set Psi2openCharm and Psi3SopenCharm to be false!"<<std::endl;
594 abort();
595 }
596 EvtPsi3Sdecay opencharm(en, part);
597 bool theDecay = opencharm.choseDecay();
598 if(!theDecay ) {
599 part->deleteTree();
600 writeFlag=0;
601 goto loop_to_select_eventsB;
602 }
603 }
604 //------------- dump the decay tree to the event model
605 // if(Nchannel>=0) part->dumpTree();
606 // part->printTree();
607
608 if(parentFlag){
609 EvtDecayTag theTag(part);
610 unsigned int modeTag, multiplicityTag;
611 modeTag = theTag.getModeTag();
612 if( getModel(part) == "OPENCHARM"||getModel(part) == "LUNDCHARM" && m_tagLundModel){
613 multiplicityTag = decayType(part);
614 //debugging
615 // std::cout<<"Opencharm Mode: "<<decayType(part)<<std::endl;
616 } else {
617 multiplicityTag = theTag.getMultTag() + decayType(part);
618 }
619
620 //std::cout<<"StdHep= "<<EvtPDL::getStdHep(part->getId())<<" getChannel "<<part->getChannel()<<" flag1,flag2= "<<modeTag<<" "<<multiplicityTag<<std::endl;
621 if(eventHeader == 0) std::cout << "event header unavailable: " << std::endl;
622 if (eventHeader != 0) {
623 eventHeader->setFlag1(modeTag);
624 eventHeader->setFlag2(multiplicityTag);
625 //std::cout<<"flag1,flag2= "<<modeTag<<" "<<multiplicityTag<<std::endl;
626 }
627 }
628
629 if(!(m_DecayRec=="")) mydecayrec=part->writeTreeRec(m_DecayRec);
630 //-- decay statistics
631 if(m_statDecays && parentFlag ) countChannel(part);
632 // ----------- pingrg 2008-03-25 ---------
633
634 if ( log.level() <= MSG::DEBUG ) part->printTree() ;
635 log << MSG::DEBUG << "Converting particles to HepMC" << endreq;
636
637 makeHepMC(part, hepMCevt, hepMCpart);
638 if(part->getNDaug()!=0)
639 hepMCpart->set_status(999);
640
641 //debug
642
643 if(part->getId()==EvtPDL::getId(m_outputp4)){
644 int nds=part->getNDaug();
645 for(int i=0;i<nds;i++){
646 if(EvtPDL::name(part->getDaug(i)->getId())=="gammaFSR") continue;
647 EvtVector4R vp1=part->getDaug(i)->getP4Lab();
648 std::cout<<"vvpp "<<vp1.get(0)<<" "<<vp1.get(1)<<" "<<vp1.get(2)<<" "<<vp1.get(3)<<std::endl;
649 }
650 }
651
652 if(m_ReadTruth.size())ReadTruth(part,m_ReadTruth);
653 //if(part->getId()==EvtPDL::getId("psi(4260)")) std::cout<<"mydebug "<<part->getP4().mass()<<std::endl;
654 part->deleteTree();
655 };
656
657 for( itp=hepMCevt->particles_begin(); itp!=hepMCevt->particles_end(); ++itp )
658 {
659 if (!(*itp)->end_vertex() && ( (*itp)->status()==899 || (*itp)->status()==999 ))
660 (*itp)->set_status(1);
661 }
662 //nchr = nchr_e + nchr_mu + nchr_pi + nchr_k +nchr_p;
664 if(m_Ncharge == true ) {std::cout<<"Multiplicity: TotoalCharge_e_mu_pi_K_p_gamma_gammaFSR: "
665 <<nchr <<" "
666 <<nchr_e <<" "
667 << nchr_mu <<" "
668 << nchr_pi <<" "
669 << nchr_k <<" "
670 <<nchr_p <<" "
671 <<N_gamma <<" "
672 <<N_gammaFSR<<endl;}
673 if(m_Ncharge == true ){std::cout<<"FinalStates: "
674 <<fst[0]<<" "
675 <<fst[1]<<" "
676 <<fst[2]<<" "
677 <<fst[3]<<" "
678 <<fst[4]<<" "
679 <<fst[5]<<" "
680 <<fst[6]<<" "
681 <<fst[7]<<" "
682 <<fst[8]<<" "
683 <<fst[9]<<" "
684 <<fst[10]<<" "
685 <<fst[11]<<" "
686 <<fst[12]<<std::endl;}
687 if( m_NtupleFile ){
688 // TotNumTrk=500;
689 m_nchr = nchr;
690 m_nchr_e = nchr_e;
691 m_nchr_mu = nchr_mu;
692 m_nchr_pi = nchr_pi;
693 m_nchr_k = nchr_k;
694 m_nchr_p = nchr_p;
695 m_gamma=N_gamma;
696 m_gammaFSR=N_gammaFSR;
697 TotNumTrk = AllTrk_index;// nchr + N_gamma + N_gammaFSR;
698 }
699 // std::cout<<"Flag= "<<eventHeader->flag1()<<" "<<eventHeader->flag2()<<std::endl;
700
701 return StatusCode::SUCCESS;
702}
703
704//****** CallBesEvtGen ************
705// This part is responsible for only ruuing BesEvtGen. //pingrg Feb. 3,2008
706//***************************************************
707StatusCode EvtDecay::callBesEvtGen( HepMC::GenEvent* hepMCevt )
708{
709 MsgStream log(messageService(), name());
710
711 nchr = -1;
712 nchr_e = 0;
713 nchr_mu= 0;
714 nchr_pi= 0;
715 nchr_k = 0;
716 nchr_p = 0;
717 N_gamma=0;
718 N_gammaFSR = 0;
719 SmartDataPtr<Event::EventHeader> eventHeader(eventSvc(),"/Event/EventHeader");
720 for (int i=0;i<13;i++) fst[i]=0;
721 HepMC::GenEvent::particle_const_iterator itp;
722 // This is the main loop.
723 // We iterate over particles and we look for ones that
724 // should be decayed with EvtGen.
725 //
726 // status == 1 - undecayed Pythia particle (also for particles that are not supposed to decay)
727 // status == 999 - particle already decayed by EvtGen
728 // status == 899 - particle was supposed to decay by EvtGen - but found no model
729 // this may be also intentional - such events are then excluded from
730 // being written into persistent output.
731 // ISTHEP(IHEP):
732 // status code for entry IHEP, with the following meanings:
733 //
734 // = 0 : null entry.
735 // = 1 : an existing entry, which has not decayed or fragmented. This is the main class of entries, which represents the `final state' given by the generator.
736 // = 2 : an entry which has decayed or fragmented and is therefore not appearing in the final state, but is retained for event history information.
737 // = 3 : a documentation line, defined separately from the event history. This could include the two incoming reacting particles, etc.
738 // = 4 - 10 : undefined, but reserved for future standards.
739 // = 11 - 200 : at the disposal of each model builder for constructs specific to his program, but equivalent to a null line in the context of any other program.
740 // = 201 - : at the disposal of users, in particular for event tracking in the detector.
741
742 loop_to_select_eventsB:
743 EvtId ppid=EvtPDL::getId(m_ParentPart); //parent particle name
744 double ppmass=0;
745 if(m_RdMeasuredEcms && m_beamEnergySpread==0) { //read Ecms, no energy spread
746 EvtPDL::reSetMass(ppid,dbEcms);
747 ppmass=EvtPDL::getMass(ppid);
748 }else if(m_RdMeasuredEcms && m_beamEnergySpread!=0){ //read Ecms, with nergy spread
749 double cmssig=sqrt(2.)*m_beamEnergySpread;
750 double myppmass = energySpread(dbEcms, cmssig);
751 EvtPDL::reSetMass(ppid,myppmass);
752 ppmass=EvtPDL::getMass(ppid);
753 }else if( (!m_RdMeasuredEcms) && m_beamEnergySpread!=0){//not read Ecms, with energy spread
754 if(m_numberEvent==1 )dbEcms = EvtPDL::getMass(ppid);
755 double cmssig=sqrt(2.)*m_beamEnergySpread;
756 double myppmass = energySpread(dbEcms, cmssig);
757 EvtPDL::reSetMass(ppid,myppmass);
758 ppmass=EvtPDL::getMass(ppid);
759 }else if((!m_RdMeasuredEcms) && m_beamEnergySpread==0 ){//not read Ecms, no energy spread
760 ppmass=EvtPDL::getMass(ppid);
761 }else {std::cout<<"EvtDecay::callBesEvtGen:: bad option"<<std::endl; abort();}
762
763 // std::cout<<"testMass "<<ppmass<<std::endl;
764 //using Ecms from data base, for XYZ data sets
765 parentMass=ppmass;
766 int pid=EvtPDL::getStdHep(ppid);
767 parentPDGcode=pid;
768
769 EvtVector4R p_init(ppmass,0.0,0.0,0.0);
770
772 HepMC::GenParticle* hepMCpart=new GenParticle(HepLorentzVector(0,0,0,ppmass),pid,3);
773
774 bool parentFlag = false;
775 // set spin density matrix for e+ e- -> V
776 if(writeFlag ==0 && part->getSpinType() == EvtSpinType::VECTOR) {
777 if(m_polarization.size()>0) {
778 part->setPolarizedSpinDensity(m_polarization[0],m_polarization[1],m_polarization[2]);
779 } else if(m_eBeamPolarization>0){
780 part->setPolarizedSpinDensity(1+ m_eBeamPolarization,0,1- m_eBeamPolarization);
781 } else{
782 part->setVectorSpinDensity();
783 }
784 parentFlag = true;
785 writeFlag++;
786 }
787
788 // for SuperBody3decay generator
789 EvtVector4R fdp_init;
790 EvtId fdp_ppid;
791 if(m_FDPparticle!=""){
792 findPart(part);
793 fdp_ppid = FDP_id;
794 fdp_init = FDP_init;
795 } else{
796 fdp_ppid = ppid;
797 fdp_init = p_init;
798 }
799
800 if(!(m_SB3File=="" || m_SB3HT=="") && first ){ SuperBody3decay_make( fdp_ppid, fdp_init);}
801 // -----------------------------------
802
803 loop_to_select_events:
804 m_Gen->generateDecay(part); if(m_numberEvent<=1){ std::cout<<"after m_Gen->decay "<<std::endl; part->printTree();}
805
806 double ratio,rdm;
807 if(_mDIY && m_ampscalflag ) _amps_max=CalAmpsMax(part);
808
809 if(_mDIY ) {
810 for(int myloop=0;myloop<100;myloop++){
811 m_Gen->generateDecay(part);
812 double amps=CalAmpsMDIY( part);
813 ratio=amps/_amps_max;
814 rdm=EvtRandom::Flat(0.0, 1.0);
815 if( !isNumber(amps) || !isNumber(_amps_max) || ratio<=0 ) {
816 part->deleteTree();
817 part=EvtParticleFactory::particleFactory(ppid,p_init);
818 continue;
819 }else if( rdm<=ratio) {
820 break;
821 }
822 }
823 }
824
825 if(m_truthFile!="" && m_truthPart!=""){
826 if(EvtPDL::name(part->getId())==m_truthPart ){
827 int ndaug = part->getNDaug();
828 for(int id=0;id<ndaug;id++){
829 EvtParticle* sonp=part->getDaug(id);
830 EvtVector4R son=sonp->getP4();
831 truth<<son.get(0)<<" "<<son.get(1)<<" "<<son.get(2)<<" "<<son.get(3)<<std::endl;
832 }
833 }
834 }
835//--------- For superbody3
836 bool accept;
837 if(m_FDPparticle==""){FDP_part=part;}
838 if(!(m_SB3File=="" || m_SB3HT=="")){
839 accept=SuperBody3decay_judge( FDP_part);
840 }
841 if(!(m_SB3File=="" || m_SB3HT=="") && !accept) {
842 delete hepMCpart;
843 part->deleteTree();
844 goto loop_to_select_eventsB;
845 }else if(!(m_SB3File=="" || m_SB3HT=="") && parentFlag && accept){
846 countChannel(FDP_part);
847 }
848
849 int Nchannel=part->getChannel();
850
851 if(m_statDecays && parentFlag) countChannel(part);
852 //------------- dump the decay tree to the event model
853 // int Nchannel=part->getChannel();
854 // if(Nchannel>=0) part->dumpTree();
855
856 if(parentFlag){
857 EvtDecayTag theTag(part);
858 unsigned int modeTag, multiplicityTag;
859 modeTag = theTag.getModeTag();
860 if( getModel(part) == "OPENCHARM"|| getModel(part) == "LUNDCHARM" && m_tagLundModel ){
861 multiplicityTag = decayType(part);
862 //debugging
863 // std::cout<<"Opencharm Mode: "<<decayType(part)<<std::endl;
864 }else if(_RvalueTag){
865 multiplicityTag = decayType(part);
866 } else {
867 multiplicityTag = theTag.getMultTag() + decayType(part);
868 }
869 if(eventHeader == 0) std::cout << "event header unavailable: " << std::endl;
870 if (eventHeader != 0) {
871 eventHeader->setFlag1(modeTag);
872 eventHeader->setFlag2(multiplicityTag);
873 }
874 //debugg
875 // std::cout<<modeTag<<" "<<multiplicityTag<<std::endl;
876 if(m_NtupleFile)m_flag1 = modeTag;
877 //std::cout<<"Flag1 "<<modeTag<<std::endl;
878 }
879
880 if(!(m_DecayRec=="")) {mydecayrec=part->writeTreeRec(m_DecayRec);std::cout<<std::endl;};
881 // ----------- pingrg 2008-03-25 ---------
882
883 if ( log.level() <= MSG::DEBUG ) part->printTree() ;
884 log << MSG::DEBUG << "Converting particles to HepMC" << endreq;
885
886 makeHepMC(part, hepMCevt, hepMCpart);
887
888 if(part->getNDaug()!=0) hepMCpart->set_status(999);
889
890 //-------------
891 if(EvtPDL::getId(m_printfs)==part->getId()) {
892 int myd=part->getNDaug();
893 std::cout<<"_FinalState "<<myd;
894 for(int ii=0;ii<myd;ii++) {
895 std::cout<<" "<<EvtPDL::name(part->getDaug(ii)->getId());
896 }
897 std::cout<<std::endl;
898 }
899
900 AllTrk_index=0;
901 for( itp=hepMCevt->particles_begin(); itp!=hepMCevt->particles_end(); ++itp )
902 {
903 if (!(*itp)->end_vertex() && ( (*itp)->status()==899 || (*itp)->status()==999 ))
904 (*itp)->set_status(1);
905
906 HepMC::GenParticle* hepMCpart=*itp;
907 int stat = hepMCpart->status();
908 int id = hepMCpart->pdg_id();
909 if(abs(id)==411 ||abs(id)==413 ) { (*itp)->set_status(999);stat=999;}
910
911 if ( stat != 1 ) continue;
912 if( id == 22 ) {N_gamma ++;fst[0]++;}
913 if( id == -22 ){N_gammaFSR ++;}
914 if(id == 11 ) {fst[1]++;}
915 else if(id == -11) {fst[2]++;}
916 else if(id == 13 ) {fst[3]++;}
917 else if(id ==-13) {fst[4]++;}
918 else if(id==211) {fst[5]++;}
919 else if(id==-211) {fst[6]++;}
920 else if(id== 321) {fst[7]++;}
921 else if(id ==-321) {fst[8]++;}
922 else if(id ==2212) {fst[9]++;}
923 else if(id==-2212) {fst[10]++;}
924 else if(id == 2112){fst[11]++;}
925 else if(id==-2112) {fst[12]++;}
926 if( fabs(id) == 11) {nchr_e++;}
927 if( fabs(id) == 13) {nchr_mu++;}
928 if( fabs(id) == 211) {nchr_pi++;}
929 if( fabs(id) == 321) {nchr_k++;}
930 if( fabs(id) == 2212) {nchr_p++;}
931
932 //for Nuple
933 double en =(hepMCpart->momentum()).e();
934 double pex=(hepMCpart->momentum()).px();
935 double pey=(hepMCpart->momentum()).py();
936 double pez=(hepMCpart->momentum()).pz();
937
938 if( m_NtupleFile==true && id!=0){
939 m_Trk_index[AllTrk_index]=id;
940 m_px_trk[AllTrk_index]=pex;
941 m_py_trk[AllTrk_index]=pey;
942 m_pz_trk[AllTrk_index]=pez;
943 m_en_trk[AllTrk_index]=en ; /*
944 std::cout<<"trk# " <<AllTrk_index
945 <<" PID:" <<id
946 <<" px: " <<pex
947 <<" py: " <<pey
948 <<" pz: " <<pez
949 <<" E: " <<en<<std::endl; */
950 AllTrk_index++;
951 }
952
953 }
954
955 itp=hepMCevt->particles_begin(); // parent decaying particle status set
956 (*itp)->set_status(2);
957
958 //nchr = nchr_e + nchr_mu + nchr_pi + nchr_k +nchr_p;
960 if(m_Ncharge == true ) {std::cout<<"Multiplicity: TotoalCharge_e_mu_pi_K_p_gamma_gammaFSR: "
961 <<nchr<<" "
962 <<nchr_e<<" "
963 << nchr_mu <<" "
964 << nchr_pi <<" "
965 << nchr_k <<" "
966 <<nchr_p<<" "
967 <<N_gamma<<" "
968 <<N_gammaFSR<<endl;}
969 if(m_Ncharge == true ){std::cout<<"FinalStates: "
970 <<fst[0]<<" "
971 <<fst[1]<<" "
972 <<fst[2]<<" "
973 <<fst[3]<<" "
974 <<fst[4]<<" "
975 <<fst[5]<<" "
976 <<fst[6]<<" "
977 <<fst[7]<<" "
978 <<fst[8]<<" "
979 <<fst[9]<<" "
980 <<fst[10]<<" "
981 <<fst[11]<<" "
982 <<fst[12]<<std::endl;}
983
984 //if(nchr==3) part->printTree();
985 if( m_NtupleFile ){
986 m_nchr = nchr;
987 m_nchr_e = nchr_e;
988 m_nchr_mu = nchr_mu;
989 m_nchr_pi = nchr_pi;
990 m_nchr_k = nchr_k;
991 m_nchr_p = nchr_p;
992 m_gamma=N_gamma;
993 m_gammaFSR=N_gammaFSR;
994 TotNumTrk = AllTrk_index; //nchr + N_gamma + N_gammaFSR;
995 }
996
997 //debug
998 if(m_ReadTruth.size())ReadTruth(part,m_ReadTruth);
999
1000
1001 //if(part->getId()==EvtPDL::getId("psi(4260)")) std::cout<<"mydebug "<<part->getP4().mass()<<std::endl;
1002
1003 part->deleteTree();
1004 return StatusCode::SUCCESS;
1005}
1006
1007//************ end of callBesEvtGen *********************
1009{
1010
1011 if(EvtCalHelAmp::nevt>0){
1012 double H2=EvtCalHelAmp::_H2;
1013 double H1=EvtCalHelAmp::_H1;
1014 double H0=EvtCalHelAmp::_H0;
1015 double H2err=EvtCalHelAmp::_H2err;
1016 double H1err=EvtCalHelAmp::_H1err;
1017 double H0err=EvtCalHelAmp::_H0err;
1018 int nd = EvtCalHelAmp::nevt;
1019 std::cout<<"Calculated helicity amplitude sqaured are:"<<std::endl;
1020 std::cout<<"|H_{2}|^2=|H_{-2}|^2= "<<H2/nd<<" +/- "<<sqrt(H2err/nd)<<endl;
1021 std::cout<<"|H_{1}|^2=|H_{-1}|^2= "<<H1/nd<<" +/- "<<sqrt(H1err/nd)<<endl;
1022 //std::cout<<"|H_{0}|^2= "<<H0/nd <<" +/- "<<sqrt(H0err/nd)<<endl;
1023 }
1024 MsgStream log(messageService(), name());
1025 truth.close();
1026 log << MSG::INFO << "EvtDecay finalized" << endreq;
1027 fstream Forfile;
1028 Forfile.open("fort.16",ios::in);
1029 char delcmd[300];
1030 strcpy(delcmd,"rm -f ");
1031 strcat(delcmd,"fort.16");
1032 // if(Forfile)system(delcmd);
1033
1034 fstream Forfile2;
1035 Forfile2.open("fort.22",ios::in);
1036 strcpy(delcmd,"rm -f ");
1037 strcat(delcmd,"fort.22");
1038 // if(Forfile2)system(delcmd);
1039
1040 // To get the branching fraction of the specified mode in Lundcharm Model
1041 /*
1042 EvtLundCharm lundcharmEvt;
1043 int nevt=lundcharmEvt.getTotalEvt();
1044 std::cout<<"The total number of the specified mode is "<<nevt<<std::endl;
1045 */
1046 int totalEvt=0;
1047 if(!(m_SB3File=="" || m_SB3HT=="")){
1048 for(int i=0;i<500;i++){totalEvt=totalEvt+br[i];}
1049 for(int i=0;i<500;i++){
1050 double bi=1.0*br[i]/totalEvt;
1051 std::cout<<"Branching fraction of channel "<<i<<" is: "<<bi<<std::endl;
1052 if(br[i]==0) break;
1053 }
1054 }
1055
1056 if(m_statDecays){
1057 int totalEvt=0;
1058 for(int i=0;i<=totalChannels;i++){totalEvt=totalEvt+br[i];}
1059 std::cout<<"==================Statistical first chain decay ==============================="<<std::endl;
1060 std::cout<<"i-th channel Events Generated Branching fraction generated "<<std::endl;
1061 for(int i=0;i<=totalChannels;i++){
1062 std::cout<<"| "<<i<<" "<<br[i]<<" "<<br[i]*1.0/totalEvt<<std::endl;
1063
1064 }
1065 std::cout<<"--------------------------------------------------------------------------------"<<std::endl;
1066 std::cout<<" sum: "<<totalEvt<<" "<<"1"<<std::endl;
1067 std::cout<<"================== End of statistical first chain decay ========================"<<std::endl;
1068 }
1069
1070 return StatusCode::SUCCESS;
1071}
1072
1073
1074StatusCode EvtDecay::makeHepMC(EvtParticle* part, HepMC::GenEvent* hEvt,
1075 HepMC::GenParticle* hPart)
1076{
1077 MsgStream log(messageService(), name());
1078
1079 if(part->getNDaug()!=0)
1080 {
1081 double ct=(part->getDaug(0)->get4Pos()).get(0);
1082 double x=(part->getDaug(0)->get4Pos()).get(1);
1083 double y=(part->getDaug(0)->get4Pos()).get(2);
1084 double z=(part->getDaug(0)->get4Pos()).get(3);
1085
1086 HepMC::GenVertex* end_vtx = new HepMC::GenVertex(HepLorentzVector(x,y,z,ct));
1087 hEvt->add_vertex( end_vtx );
1088 end_vtx->add_particle_in(hPart);
1089
1090 int ndaug=part->getNDaug();
1091
1092 for(int it=0;it<ndaug;it++)
1093 {
1094
1095 double e=(part->getDaug(it)->getP4Lab()).get(0);
1096 double px=(part->getDaug(it)->getP4Lab()).get(1);
1097 double py=(part->getDaug(it)->getP4Lab()).get(2);
1098 double pz=(part->getDaug(it)->getP4Lab()).get(3);
1099 int id=EvtPDL::getStdHep(part->getDaug(it)->getId());
1100 int status=1;
1101
1102 if(part->getDaug(it)->getNDaug()!=0) status=999;
1103
1104 HepMC::GenParticle* prod_part = new HepMC::GenParticle(HepLorentzVector(px,py,pz,e),id,status);
1105 end_vtx->add_particle_out(prod_part);
1106 makeHepMC(part->getDaug(it),hEvt,prod_part);
1107
1108
1109
1110 if( log.level()<MSG::INFO )
1111 prod_part->print();
1112 }
1113 }
1114
1115 return StatusCode::SUCCESS;
1116}
1117
1118
1119void EvtDecay::MeVToGeV (HepMC::GenEvent* evt)
1120{
1121 for ( HepMC::GenEvent::particle_iterator p = evt->particles_begin(); p != evt->particles_end(); ++p ) {
1122 // cout << " PDG, BAR " << (*p)->pdg_id() << " " << (*p)->barcode() << endl;
1123 // (*p)->set_momentum( (*p)->momentum() / 1000. );
1124
1125 HepMC::FourVector tmpfv((*p)->momentum().x()/1000.,
1126 (*p)->momentum().y()/1000.,
1127 (*p)->momentum().z()/1000.,
1128 (*p)->momentum().t()/1000.
1129 );
1130
1131 (*p)->set_momentum( tmpfv );
1132 }
1133}
1134
1135void EvtDecay::GeVToMeV (HepMC::GenEvent* evt)
1136{
1137 for ( HepMC::GenEvent::particle_iterator p = evt->particles_begin(); p != evt->particles_end(); ++p ) {
1138 // cout << " PDG, BAR " << (*p)->pdg_id() << " " << (*p)->barcode() << endl;
1139 // (*p)->set_momentum( (*p)->momentum() * 1000. );
1140 HepMC::FourVector tmpfv((*p)->momentum().x()*1000.,
1141 (*p)->momentum().y()*1000.,
1142 (*p)->momentum().z()*1000.,
1143 (*p)->momentum().t()*1000.
1144 );
1145
1146 (*p)->set_momentum( tmpfv );
1147 }
1148}
1149
1150
1151double EvtDecay::CalAmpsMax( EvtParticle* part )
1152{
1153 double amps_max=0;
1154 for(int i=0;i<100000;i++){
1155 m_Gen->generateDecay(part);
1156 double amps=CalAmpsMDIY(part);
1157 if(amps > amps_max) amps_max=amps;
1158 }
1159 amps_max=amps_max*1.05;
1160 m_ampscalflag = false;
1161 return amps_max;
1162}
1163
1164// EvtBesRandom class implementation, basically BesRandomSvc -> EvtGen random engine interface
1165
1166EvtBesRandom::EvtBesRandom(HepRandomEngine* engine)
1167{
1168 m_engine = engine;
1169 m_engine->showStatus();
1170}
1171
1174
1176{
1177 //double rdm=EvtRandom::Flat(0.0, 1.0);;
1178 double rdm= RandFlat::shoot(m_engine);
1179 // std::cout<<"rdm= "<<rdm<<std::endl;
1180 return rdm;
1181 //return RandFlat::shoot(m_engine);
1182}
1183
1184
1185void EvtDecay::SuperBody3decay_make(EvtId ppid, EvtVector4R p_init ){
1186 double xmass2,ymass2;
1187 bool accept;
1188 EvtVector4R pd1,pd2,pd3;
1189
1190 //---- find out the pdg codes of final state and count the identical particles
1191 FinalState_make( ppid, p_init);
1192
1193 //begin to loop with events
1194 for(int i=0;i<100000;i++){
1195 regen:
1197 HepMC::GenParticle* hepMCpart=new GenParticle(HepLorentzVector(0,0,0,parentMass),parentPDGcode,3); //this line make the decay to select different channel
1198
1199 if(part->getSpinType() == EvtSpinType::VECTOR) {
1200 part->setVectorSpinDensity();
1201 }
1202 m_Gen->generateDecay(part);
1203
1204
1205 FinalState_sort(part);
1206 pd1=son0;
1207 pd2=son1;
1208 pd3=son2;
1209
1210 // std::cout<<"multi, pdg, pi= "<<multi<<" "<<pdg0<<" "<<pdg1<<" "<<pdg2<<" "<<son0<<son1<<son2<<std::endl;
1211
1212 xmass2=(pd1+pd2).mass2(); // Dalitz plot m12^2 ~ m13^2
1213 ymass2=(pd1+pd3).mass2();
1214 double xlow=SuperBody3decay.getXlow();
1215 double xup=SuperBody3decay.getXup();
1216 double ylow=SuperBody3decay.getYlow();
1217 double yup=SuperBody3decay.getYup();
1218 if(xmass2<xlow || xmass2>xup || ymass2<ylow || ymass2>yup) { part->deleteTree();goto regen;}
1219 SuperBody3decay.HFill(xmass2,ymass2);
1220
1221 if( m_NtupleFile ){
1222 m_ich=part->getChannel();
1223 m_m1=pd1.mass();
1224 m_m2=pd2.mass();
1225 m_m3=pd3.mass();
1226 m_m12=(pd1+pd2).mass();
1227 m_m23=(pd2+pd3).mass();
1228 m_m13=(pd1+pd3).mass();
1229 m_cos1=pd1.get(3)/pd1.d3mag();
1230 m_cos2=pd2.get(3)/pd2.d3mag();
1231 m_cos3=pd3.get(3)/pd3.d3mag();
1232 mass_tuple->write();
1233 }
1234 double m1=pd1.mass(); double m2=pd2.mass();double m3=pd3.mass();
1235 double mj=(pd2+pd3).mass();
1236 // std::cout<<"mass 1 2 3 chicj <<=="<<m1<<" "<<m2<<" "<<m3<<" "<<mj<<std::endl;
1237 // delete hepMCpart;
1238 }
1239
1240 SuperBody3decay.HReweight(); //reweighting Dalitz plot
1241 SuperBody3decay.setZmax(); // find out the maximum value after reweighting
1242 first=false;
1243}
1244
1245bool EvtDecay::SuperBody3decay_judge(EvtParticle* part){
1246 double xmass2,ymass2;
1247 bool accept;
1248 EvtVector4R pd1,pd2,pd3;
1249
1250
1251 FinalState_sort( part);
1252
1253 pd1=son0;
1254 pd2=son1;
1255 pd3=son2;
1256
1257
1258 xmass2=(pd1+pd2).mass2(); // Dalitz plot m12^2 ~ m13^2
1259 ymass2=(pd1+pd3).mass2();
1260 accept=SuperBody3decay.AR(xmass2,ymass2);
1261
1262 //debugging
1263 if(accept && m_NtupleFile) {
1264 _ich=part->getChannel();
1265 _m1=pd1.mass();
1266 _m2=pd2.mass();
1267 _m3=pd3.mass();
1268 _m12=(pd1+pd2).mass();
1269 _m23=(pd2+pd3).mass();
1270 _m13=(pd1+pd3).mass();
1271 _cos1=pd1.get(3)/pd1.d3mag();
1272 _cos2=pd2.get(3)/pd2.d3mag();
1273 _cos3=pd3.get(3)/pd3.d3mag();
1274 massgen_tuple->write();
1275 }
1276 // std::cout<<"mass 1 2 3 chicj>> "<<_m1<<" "<<_m2<<" "<<_m3<<" "<<_m23<<std::endl;
1277
1278 return accept;
1279}
1280
1281
1282void EvtDecay:: FinalState_make(EvtId ppid, EvtVector4R p_init){
1283 // get the final state pdg cond and count the identicle particle
1284 multi=1;
1285 identical_flag=true;
1286
1287 for(int i=1;i<10000;i++){ //sigle out the final state
1289 HepMC::GenParticle* hepMCpart=new GenParticle(HepLorentzVector(0,0,0,parentMass),parentPDGcode,3);
1290
1291 m_Gen->generateDecay(part);
1292 int nson=part->getNDaug();
1293
1294 if(nson == 2) {continue;} else
1295 if(nson ==3){
1296 EvtId xid0,xid1,xid2;
1297 xid0=part->getDaug(0)->getId();
1298 xid1=part->getDaug(1)->getId();
1299 xid2=part->getDaug(2)->getId();
1300
1301//-- pdg code
1302 pdg0=EvtPDL::getStdHep(xid0);
1303 pdg1=EvtPDL::getStdHep(xid1);
1304 pdg2=EvtPDL::getStdHep(xid2);
1305
1306 if(pdg0==pdg1 && pdg1==pdg2) {multi=3;}
1307 else if(pdg0==pdg1 ){multi=2;} // two identical particle list as 0,1 index
1308 else if(pdg0==pdg2 ){multi=2;pdg=pdg1; pdg1=pdg2; pdg2=pdg;}
1309 else if(pdg1==pdg2) {multi=2;pdg=pdg0; pdg0=pdg1; pdg1=pdg2; pdg2=pdg;}
1310 else {multi=1;}
1311 break;
1312 } else{
1313 std::cout<<"No direct three body decay channel found in decay card, I stop run"<<std::endl;
1314 ::abort();
1315 }
1316 }
1317 // std::cout<<"pdg_i "<<pdg0<<" "<<pdg1<<" "<<pdg2<<std::endl;
1318 // std::cout<<"identical particle "<<multi<<std::endl;
1319}
1320
1321void EvtDecay::FinalState_sort(EvtParticle* part){
1322 // sort the momentum to aon0, son1, son2
1323 EvtVector4R pd0,pd1,pd2;
1324 EvtId xid0,xid1,xid2;
1325 int id0,id1,id2;
1326
1327 int nson=part->getNDaug();
1328 if(nson==2){
1329 pd0=part->getDaug(0)->getP4Lab();
1330 pd1=part->getDaug(1)->getDaug(0)->getP4Lab();
1331 pd2=part->getDaug(1)->getDaug(1)->getP4Lab();
1332
1333 xid0=part->getDaug(0)->getId();
1334 xid1=part->getDaug(1)->getDaug(0)->getId();
1335 xid2=part->getDaug(1)->getDaug(1)->getId();
1336
1337 } else if(nson==3){
1338 pd0=part->getDaug(0)->getP4Lab();
1339 pd1=part->getDaug(1)->getP4Lab();
1340 pd2=part->getDaug(2)->getP4Lab();
1341
1342 xid0=part->getDaug(0)->getId();
1343 xid1=part->getDaug(1)->getId();
1344 xid2=part->getDaug(2)->getId();
1345 }
1346 //-- pdg code
1347 id0=EvtPDL::getStdHep(xid0);
1348 id1=EvtPDL::getStdHep(xid1);
1349 id2=EvtPDL::getStdHep(xid2);
1350
1351 // std::cout<<"pid before match: "<<id0<<" "<<id1<<" "<<id2<<std::endl;
1352 //-- assign sons momentum
1353 if(multi==1){
1354 assign_momentum(id0,pd0);
1355 assign_momentum(id1,pd1);
1356 assign_momentum(id2,pd2);
1357 } else if(multi==2){
1358 assign_momentum2(id0,pd0);
1359 assign_momentum2(id1,pd1);
1360 assign_momentum2(id2,pd2);
1361 if(son0.get(0)>son1.get(0)){son=son0;son0=son1;son1=son;} // change into E_0 < E_1
1362 } else if(multi==3){ // sort sons according to energy E_0 < E_1 < E_2
1363 double En0=pd0.get(0);
1364 double En1=pd1.get(0);
1365 double En2=pd2.get(0);
1366 if(En0 < En1 && En1 < En2) {son0=pd0;son1=pd1;son2=pd2;}
1367 if(En0 < En2 && En2 < En1) {son0=pd0;son1=pd2;son2=pd1;}
1368 if(En1 < En0 && En0 < En2) {son0=pd1;son1=pd0;son2=pd2;}
1369 if(En1 < En2 && En2 < En0) {son0=pd1;son1=pd2;son2=pd0;}
1370 if(En2 < En0 && En0 < En1) {son0=pd2;son1=pd0;son2=pd1;}
1371 if(En2 < En1 && En1 < En0) {son0=pd2;son1=pd1;son2=pd0;}
1372 }
1373
1374}
1375
1376
1377void EvtDecay::assign_momentum(int id0, EvtVector4R pd0){
1378 if(id0==pdg0){son0=pd0;}
1379 else if(id0==pdg1){son1=pd0;}
1380 else if(id0==pdg2){son2=pd0;}
1381 else {std::cout<<"PID= "<<id0<<" not machted, I stop"<<std::endl; ::abort();}
1382}
1383
1384void EvtDecay::assign_momentum2(int id0, EvtVector4R pd0){ // for two identicle particle case
1385 if(id0==pdg0 && identical_flag){son0=pd0;identical_flag=false;}
1386 else if(id0==pdg1){son1=pd0;identical_flag=true;}
1387 else if(id0==pdg2){son2=pd0;}
1388 else {std::cout<<"PID not machted, I stop"<<std::endl; ::abort();}
1389}
1390
1391void EvtDecay::findPart(EvtParticle* part){
1392 FDP_id=EvtPDL::getId(m_FDPparticle);
1393 int FDP_pdg=EvtPDL::getStdHep(FDP_id);
1394
1395 m_Gen->generateDecay(part);
1396 int dn=part->getNDaug();
1397
1398 if(dn!=0){
1399 for(int i=0;i<dn;i++){
1400 EvtParticle* part1=part->getDaug(i);
1401 EvtId sonid=part1->getId();
1402 int son_pdg = EvtPDL::getStdHep(sonid);
1403 if(son_pdg == FDP_pdg){
1404 FDP_init=part1->getP4Lab();
1405 FDP_part=part1;
1406 return;
1407 } else{
1408 findPart(part1);
1409 }
1410 }
1411 } //if (dn block
1412
1413}
1414
1415void EvtDecay::countChannel(EvtParticle* part){
1416 int ich=part->getChannel();
1417
1418 //debuging
1419 //if(getModel(part,ich)=="OPENCHARM" &&EvtPDL::name( part->getId() )=="psi(4260)") ich = part->getGeneratorFlag();
1420 //std::cout<<"the channel is "<<ich<<std::endl;
1421
1422 br[ich]++;
1423 if(ich>totalChannels) totalChannels = ich;
1424 //std::cout<<"EVENT IN br_i "<<br[ich]<<std::endl;
1425}
1426
1427bool EvtDecay::isCharmonium(EvtId xid){
1428EvtId psi4415 = EvtPDL::getId(std::string("psi(4415)"));
1429EvtId psi4160 = EvtPDL::getId(std::string("psi(4160)"));
1430EvtId psi4040 = EvtPDL::getId(std::string("psi(4040)"));
1431EvtId psi3770 = EvtPDL::getId(std::string("psi(3770)"));
1432EvtId psiprim = EvtPDL::getId(std::string("psi(2S)"));
1433EvtId jpsi = EvtPDL::getId(std::string("J/psi"));
1434EvtId etac = EvtPDL::getId(std::string("eta_c"));
1435EvtId etac2s = EvtPDL::getId(std::string("eta_c(2S)"));
1436EvtId hc = EvtPDL::getId(std::string("h_c"));
1437EvtId chic0 = EvtPDL::getId(std::string("chi_c0"));
1438EvtId chic1 = EvtPDL::getId(std::string("chi_c1"));
1439EvtId chic2 = EvtPDL::getId(std::string("chi_c2"));
1440EvtId chic0p = EvtPDL::getId(std::string("chi_c0p"));
1441EvtId chic1p = EvtPDL::getId(std::string("chi_c1p"));
1442EvtId chic2p = EvtPDL::getId(std::string("chi_c1p"));
1443 std::vector<EvtId> Vid; Vid.clear();
1444 Vid.push_back(psi4415);
1445 Vid.push_back(psi4160);
1446 Vid.push_back(psi4040);
1447 Vid.push_back(psi3770);
1448 Vid.push_back(psiprim);
1449 Vid.push_back(jpsi);
1450 Vid.push_back(etac);
1451 Vid.push_back(etac2s);
1452 Vid.push_back(hc);
1453 Vid.push_back(chic0);
1454 Vid.push_back(chic1);
1455 Vid.push_back(chic2);
1456 Vid.push_back(chic0p);
1457 Vid.push_back(chic1p);
1458 Vid.push_back(chic2p);
1459
1460 bool flag=true;
1461 for(int i=0;i<Vid.size();i++){ if(xid == Vid[i]) return flag;}
1462 return false;
1463}
1464
1465
1466bool EvtDecay::isCharm(EvtId xid){
1467EvtId d0 = EvtPDL::getId(std::string("D0"));
1468EvtId d0bar = EvtPDL::getId(std::string("anti-D0"));
1469EvtId dp = EvtPDL::getId(std::string("D+"));
1470EvtId dm = EvtPDL::getId(std::string("D-"));
1471EvtId d0h = EvtPDL::getId(std::string("D0H"));
1472EvtId d0l = EvtPDL::getId(std::string("D0L"));
1473EvtId dstp = EvtPDL::getId(std::string("D*+"));
1474EvtId dstm = EvtPDL::getId(std::string("D*-"));
1475EvtId ds0 = EvtPDL::getId(std::string("D*0"));
1476EvtId ds0bar = EvtPDL::getId(std::string("anti-D*0"));
1477EvtId dsp = EvtPDL::getId(std::string("D_s+"));
1478EvtId dsm = EvtPDL::getId(std::string("D_s-"));
1479EvtId dsstp = EvtPDL::getId(std::string("D_s*+"));
1480EvtId dsstm = EvtPDL::getId(std::string("D_s*-"));
1481EvtId ds0stp = EvtPDL::getId(std::string("D_s0*+"));
1482EvtId ds0stm = EvtPDL::getId(std::string("D_s0*-"));
1483
1484 std::vector<EvtId> Vid; Vid.clear();
1485 Vid.push_back(d0);
1486 Vid.push_back(d0bar);
1487 Vid.push_back(dp);
1488 Vid.push_back(dm);
1489 Vid.push_back(d0h);
1490 Vid.push_back(d0l);
1491 Vid.push_back(dstp);
1492 Vid.push_back(dstm);
1493 Vid.push_back(ds0);
1494 Vid.push_back(ds0bar );
1495 Vid.push_back(dsp );
1496 Vid.push_back(dsm );
1497 Vid.push_back(dsstp );
1498 Vid.push_back(dsstm );
1499 Vid.push_back(ds0stp );
1500 Vid.push_back(ds0stm );
1501
1502 bool flag=true;
1503 for(int i=0;i<Vid.size();i++){ if(xid == Vid[i]) return flag;}
1504 return false;
1505}
1506
1507bool EvtDecay::isRadecay(EvtParticle *par){
1508 int ndg = par->getNDaug();
1509 for(int i=0;i<ndg;i++){
1510 EvtId xid = par->getDaug(i)->getId();
1511 if(EvtPDL::getStdHep(xid) == 22){return true;}
1512 }
1513 return false;
1514}
1515
1516int EvtDecay::decayType(EvtParticle *par){
1517 /*************************************
1518 * 0: gamma light_hadrons
1519 * 1: gamma charmonium
1520 * 2: DDbar (D0 D0bar or D+D-)
1521 * 3: ligh_hadrons
1522 * 4: others
1523 *************************************/
1524 int Nchannel=par->getChannel();
1525 int ndg = par->getNDaug();
1526 int nfsr=0;
1527
1528 if(_RvalueTag){return Nchannel;}
1529 // std::cout<<"decayType= "<<Nchannel<<std::endl;
1530
1531 std::string model = getModel(par,Nchannel);
1532 if( model == "OPENCHARM" || model == "LUNDCHARM" && m_tagLundModel){ // lundcharm model tag
1533 int ldm = par->getGeneratorFlag();
1534 // std::cout<<"LundCharmFlag = "<<ldm<<std::endl;
1535 return ldm;
1536 //return 9;
1537 }
1538
1539 for(int i=0;i<ndg;i++){ //remove the FSR photon
1540 EvtId xid =par->getDaug(i)->getId();
1541 if(EvtPDL::getStdHep(xid) == -22 ){nfsr++;}
1542 }
1543
1544 if( isRadecay(par) ){
1545 for(int i=0;i<ndg;i++){
1546 EvtId xid = par->getDaug(i)->getId();
1547 if( isCharmonium(xid) ) return 1;
1548 }
1549 return 0;
1550 } else{
1551 if(ndg-nfsr ==2 ){
1552 int ndg = par->getNDaug();
1553 EvtId xd1 = par->getDaug(0)->getId();
1554 EvtId xd2 = par->getDaug(1)->getId();
1555 if( isCharm(xd1) && isCharm(xd2) ) {return 2;} else
1556 if( !isCharmonium(xd1) && !isCharmonium(xd2) && !isCharm(xd1) && !isCharm(xd2) ){
1557 return 3;
1558 }
1559 } else{ // ndg>=3
1560 bool flag = false;
1561 for(int i=0;i<ndg;i++){
1562 EvtId xid = par->getDaug(i)->getId();
1563 if( isCharmonium(xid) ) {flag = true; break;}
1564 if( isCharm(xid) ) {flag = true; break;}
1565 } // for_i block
1566 if( !flag ){return 3;} else {return 4;}
1567 }// if_ndg block
1568 }
1569
1570}
1571
1572
1573std::string EvtDecay::getModel(EvtParticle *par, int mode){
1574 EvtDecayBase* thedecay = EvtDecayTable::gettheDecay(par->getId(),mode);
1575 return thedecay->getModelName();
1576}
1577
1578std::string EvtDecay::getModel(EvtParticle *par){
1579 int mode = par->getChannel();
1580 EvtDecayBase* thedecay = EvtDecayTable::gettheDecay(par->getId(),mode);
1581 return thedecay->getModelName();
1582}
1583
1584
1585void EvtDecay::ReadTruth(EvtParticle* part,std::vector<std::vector<string> > mylist){
1586 if(mylist.size()<2) {std::cout<<" Empty list"<<std::endl;abort();}
1587 EvtVector4R vp1;
1588 for(int k=0;k<mylist[0].size();k++){
1589 if(part->getId()==EvtPDL::getId(mylist[0][k])){
1590 if(part->getDaug(1)->getId()==EvtPDL::getId("vhdr")){part=part->getDaug(1);continue;}
1591 if(part->getDaug(0)->getId()==EvtPDL::getId("vgam")){part=part->getDaug(1);continue;}
1592 for(int i=1;i<mylist.size();i++){
1593 EvtParticle* mypart=part;
1594 for(int j=0;j<mylist[i].size();j++){
1595 mypart=mypart->getDaug(atoi(mylist[i][j].c_str()));
1596 //std::cout<<atoi(mylist[i][j].c_str());
1597 }
1598 //std::cout<<std::endl;
1599 std::cout<<EvtPDL::name(mypart->getId())<<std::endl;
1600 if(_truthAtCM){
1601 vp1=mypart->getP4();
1602 }else{
1603 vp1=mypart->getP4Lab();
1604 }
1605 if( !isNumber(vp1.get(0)) ) {part->printTree(); abort();}
1606 //std::cout<<"truth "<<vp1.get(0)<<" "<<vp1.get(1)<<" "<<vp1.get(2)<<" "<<vp1.get(3)<<std::endl;
1607 printf("truth %.12lf %.12lf %.12lf %.12lf \n",vp1.get(0),vp1.get(1),vp1.get(2),vp1.get(3) );
1608 }
1609 }
1610 }
1611}
1612
1613int EvtDecay::isNumber(double d){return (d==d);}//if d=nan, return 0, otherwise, return 1
1614
1615
1616double EvtDecay::energySpread(double mu,double sigma){
1617 //mu: mean value in Gaussian
1618 //sigma: variance in Gaussian
1619 rloop:
1620 int n=12;
1621 double ri=0;
1622 for(int i=0;i<n;i++){
1623 double pm= EvtRandom::Flat(0.,1);
1624 ri += pm;
1625 }
1626 double eta=sqrt(n*12.0)*(ri/12-0.5);
1627 double xsig= sigma*eta+mu;
1628 double mx0=EvtConExc::mlow;
1629 double mx1=EvtConExc::mup;
1630 if(xsig<mx0 || xsig>mx1) goto rloop;
1631 return xsig;
1632}
1633
1634#include "../user/Lenu.inc"
double mass
TTree * sigma
int runNo
Definition DQA_TO_DB.cxx:12
const Int_t n
Double_t x[10]
int fst[50]
Definition EvtDecay.cxx:79
int nchr
Definition EvtDecay.cxx:71
int nchr_mu
Definition EvtDecay.cxx:73
int nchr_pi
Definition EvtDecay.cxx:74
int N_gamma
Definition EvtDecay.cxx:77
int nchr_k
Definition EvtDecay.cxx:75
int N_gammaFSR
Definition EvtDecay.cxx:78
int nchr_p
Definition EvtDecay.cxx:76
int nchr_e
Definition EvtDecay.cxx:72
DOUBLE_PRECISION xup[20]
DOUBLE_PRECISION xlow[20]
XmlRpcServer s
ObjectVector< McGenEvent > McGenEventCol
Definition McGenEvent.h:39
INTupleSvc * ntupleSvc()
*************DOUBLE PRECISION m_pi *DOUBLE PRECISION m_HvecTau2 DOUBLE PRECISION m_HvClone2 DOUBLE PRECISION m_gamma1 DOUBLE PRECISION m_gamma2 DOUBLE PRECISION m_thet1 DOUBLE PRECISION m_thet2 INTEGER m_IFPHOT *COMMON c_Taupair $ !Spin Polarimeter vector first Tau $ !Spin Polarimeter vector second Tau $ !Clone Spin Polarimeter vector first Tau $ !Clone Spin Polarimeter vector second Tau $ !Random Euler angle for cloning st tau $ !Random Euler angle for cloning st tau $ !Random Euler angle for cloning st tau $ !Random Euler angle for cloning nd tau $ !Random Euler angle for cloning nd tau $ !Random Euler angle for cloning nd tau $ !phi of HvecTau1 $ !theta of HvecTau1 $ !phi of HvecTau2 $ !theta of HvecTau2 $ !super key
Definition Taupair.h:42
void setDecayCard(string card)
double random()
EvtBesRandom(HepRandomEngine *engine)
virtual ~EvtBesRandom()
static double _H0err
static double _H1err
static double _H0
static int nevt
static double _H1
static double _H2
static double _H2err
static double SetMthr
Definition EvtConExc.hh:192
static double mlow
Definition EvtConExc.hh:235
static double mup
Definition EvtConExc.hh:235
std::string getModelName()
static EvtDecayBase * gettheDecay(EvtId parent, int imode)
DataInfoSvc * dataInfoSvc
Definition EvtDecay.h:75
IDataInfoSvc * tmpInfoSvc
Definition EvtDecay.h:74
double energySpread(double mu, double sigma)
StatusCode initialize()
Definition EvtDecay.cxx:143
StatusCode finalize()
StatusCode execute()
Definition EvtDecay.cxx:343
void readUDecay(const char *const udecay_name)
Definition EvtGen.cc:126
void generateDecay(int stdhepid, EvtVector4R P, EvtVector4R D, EvtStdHep *evtStdHep, EvtSpinDensity *spinDensity=0)
Definition EvtGen.cc:152
static std::vector< std::string > SV
static std::vector< std::vector< double > > dVV
static std::vector< std::vector< int > > iVV
void HReweight()
Definition EvtHis2F.cc:137
void init()
Definition EvtHis2F.cc:98
double getYup()
Definition EvtHis2F.hh:69
bool AR(double xmass2, double ymass2)
Definition EvtHis2F.cc:205
void HFill(double xmass2, double ymass2)
Definition EvtHis2F.cc:129
void setHTitle(const char *htitle)
Definition EvtHis2F.cc:40
double getYlow()
Definition EvtHis2F.hh:67
void setZmax()
Definition EvtHis2F.cc:177
double getXup()
Definition EvtHis2F.hh:68
double getXlow()
Definition EvtHis2F.hh:66
void setFile(const char *dtfile)
Definition EvtHis2F.cc:35
Definition EvtId.hh:27
static int getStdHep(EvtId id)
Definition EvtPDL.hh:56
static void reSetMass(EvtId i, double mass)
Definition EvtPDL.hh:68
static EvtId evtIdFromStdHep(int stdhep)
Definition EvtPDL.cc:244
static std::string name(EvtId i)
Definition EvtPDL.hh:64
static double getMass(EvtId i)
Definition EvtPDL.hh:46
static EvtId getId(const std::string &name)
Definition EvtPDL.cc:287
static EvtParticle * particleFactory(EvtSpinType::spintype spinType)
EvtVector4R getP4Lab()
EvtId getId() const
static EvtId _NextLevelId[20]
void setVectorSpinDensity()
static int _NextLevelDauNum
void setPolarizedSpinDensity(double r00, double r11, double r22)
EvtSpinType::spintype getSpinType() const
const EvtVector4R & getP4() const
void printTree() const
int getGeneratorFlag()
int getNDaug() const
EvtParticle * getDaug(int i)
static EvtVector4R _NextLevelP4[20]
int getChannel() const
std::string writeTreeRec(std::string) const
void deleteTree()
static double Flat()
Definition EvtRandom.cc:74
double mass() const
double get(int i) const
double d3mag() const
virtual void setGenseed(long)=0
virtual CLHEP::HepRandomEngine * GetEngine(const std::string &StreamName)=0
Interface to the CLHEP engine.
double y[1000]
double complex pa0 double complex pb0ij double complex pc0ijk double complex pd0
char * c_str(Index i)
const double hc
Definition TConstant.h:41
double double double double double * m3
Definition qcdloop1.h:76
double * m1
Definition qcdloop1.h:75
double double * m2
Definition qcdloop1.h:75