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