1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // Particles stack class //
21 // Implements the TMCVirtualStack of the Virtual Monte Carlo //
22 // Holds the particles transported during simulation //
23 // Is used to compare results of reconstruction with simulation //
26 ///////////////////////////////////////////////////////////////////////////////
29 #include <TClonesArray.h>
30 #include <TObjArray.h>
32 #include <TParticle.h>
33 #include <TParticlePDG.h>
34 #include <TDatabasePDG.h>
36 #include <TDirectory.h>
43 //_______________________________________________________________________
45 fParticles("TParticle", 1000),
60 // Default constructor
64 //_______________________________________________________________________
65 AliStack::AliStack(Int_t size, const char* /*evfoldname*/):
66 fParticles("TParticle",1000),
85 //_______________________________________________________________________
86 AliStack::AliStack(const AliStack& st):
88 fParticles("TParticle",1000),
89 fParticleMap(*(st.Particles())),
90 fParticleFileMap(st.fParticleFileMap),
93 fTreeK((TTree*)(st.fTreeK->Clone())),
94 fNtrack(st.GetNtrack()),
95 fNprimary(st.GetNprimary()),
106 //_______________________________________________________________________
107 void AliStack::Copy(TObject&) const
109 AliFatal("Not implemented!");
112 //_______________________________________________________________________
113 AliStack::~AliStack()
126 //_____________________________________________________________________________
127 void AliStack::PushTrack(Int_t done, Int_t parent, Int_t pdg, const Float_t *pmom,
128 const Float_t *vpos, const Float_t *polar, Float_t tof,
129 TMCProcess mech, Int_t &ntr, Float_t weight, Int_t is)
132 // Load a track on the stack
134 // done 1 if the track has to be transported
136 // parent identifier of the parent track. -1 for a primary
138 // pmom momentum GeV/c
140 // polar polarisation
141 // tof time of flight in seconds
142 // mecha production mechanism
143 // ntr on output the number of the track stored
146 // const Float_t tlife=0;
149 // Here we get the static mass
150 // For MC is ok, but a more sophisticated method could be necessary
151 // if the calculated mass is required
152 // also, this method is potentially dangerous if the mass
153 // used in the MC is not the same of the PDG database
155 TParticlePDG* pmc = TDatabasePDG::Instance()->GetParticle(pdg);
157 Float_t mass = TDatabasePDG::Instance()->GetParticle(pdg)->Mass();
158 Float_t e=TMath::Sqrt(mass*mass+pmom[0]*pmom[0]+
159 pmom[1]*pmom[1]+pmom[2]*pmom[2]);
161 // printf("Loading mass %f ene %f No %d ip %d parent %d done %d pos %f %f %f mom %f %f %f kS %d m \n",
162 // mass,e,fNtrack,pdg,parent,done,vpos[0],vpos[1],vpos[2],pmom[0],pmom[1],pmom[2],kS);
165 PushTrack(done, parent, pdg, pmom[0], pmom[1], pmom[2], e,
166 vpos[0], vpos[1], vpos[2], tof, polar[0], polar[1], polar[2],
167 mech, ntr, weight, is);
169 AliWarning(Form("Particle type %d not defined in PDG Database !", pdg));
170 AliWarning("Particle skipped !");
174 //_____________________________________________________________________________
175 void AliStack::PushTrack(Int_t done, Int_t parent, Int_t pdg,
176 Double_t px, Double_t py, Double_t pz, Double_t e,
177 Double_t vx, Double_t vy, Double_t vz, Double_t tof,
178 Double_t polx, Double_t poly, Double_t polz,
179 TMCProcess mech, Int_t &ntr, Double_t weight, Int_t is)
182 // Load a track on the stack
184 // done 1 if the track has to be transported
186 // parent identifier of the parent track. -1 for a primary
188 // kS generation status code
189 // px, py, pz momentum GeV/c
190 // vx, vy, vz position
191 // polar polarisation
192 // tof time of flight in seconds
193 // mech production mechanism
194 // ntr on output the number of the track stored
196 // New method interface:
197 // arguments were changed to be in correspondence with TParticle
199 // Note: the energy is not calculated from the static mass but
200 // it is passed by argument e.
202 const Int_t kFirstDaughter=-1;
203 const Int_t kLastDaughter=-1;
207 = new(fParticles[fLoadPoint++])
208 TParticle(pdg, is, parent, -1, kFirstDaughter, kLastDaughter,
209 px, py, pz, e, vx, vy, vz, tof);
211 particle->SetPolarisation(polx, poly, polz);
212 particle->SetWeight(weight);
213 particle->SetUniqueID(mech);
218 particle->SetBit(kDoneBit);
220 particle->SetBit(kTransportBit);
225 // Declare that the daughter information is valid
226 particle->SetBit(kDaughtersBit);
227 // Add the particle to the stack
230 fParticleMap.AddAtAndExpand(particle, fNtrack);//CHECK!!
233 particle = dynamic_cast<TParticle*>(fParticleMap.At(parent));
235 particle->SetLastDaughter(fNtrack);
236 if(particle->GetFirstDaughter()<0) particle->SetFirstDaughter(fNtrack);
239 AliError(Form("Parent %d does not exist",parent));
243 // This is a primary track. Set high water mark for this event
246 // Set also number if primary tracks
247 fNprimary = fHgwmk+1;
253 //_____________________________________________________________________________
254 TParticle* AliStack::PopNextTrack(Int_t& itrack)
257 // Returns next track from stack of particles
261 TParticle* track = GetNextParticle();
265 track->SetBit(kDoneBit);
270 fCurrentTrack = track;
274 //_____________________________________________________________________________
275 TParticle* AliStack::PopPrimaryForTracking(Int_t i)
278 // Returns i-th primary particle if it is flagged to be tracked,
282 TParticle* particle = Particle(i);
284 if (!particle->TestBit(kDoneBit))
290 //_____________________________________________________________________________
291 Bool_t AliStack::PurifyKine()
294 // Compress kinematic tree keeping only flagged particles
295 // and renaming the particle id's in all the hits
298 int nkeep = fHgwmk + 1, parent, i;
299 TParticle *part, *father;
300 fTrackLabelMap.Set(fParticleMap.GetLast()+1);
302 // Save in Header total number of tracks before compression
303 // If no tracks generated return now
304 if(fHgwmk+1 == fNtrack) return kFALSE;
306 // First pass, invalid Daughter information
307 for(i=0; i<fNtrack; i++) {
308 // Preset map, to be removed later
309 if(i<=fHgwmk) fTrackLabelMap[i]=i ;
311 fTrackLabelMap[i] = -99;
312 if((part=dynamic_cast<TParticle*>(fParticleMap.At(i)))) {
314 // Check of this track should be kept for physics reasons
315 if (KeepPhysics(part)) KeepTrack(i);
317 part->ResetBit(kDaughtersBit);
318 part->SetFirstDaughter(-1);
319 part->SetLastDaughter(-1);
323 // Invalid daughter information for the parent of the first particle
324 // generated. This may or may not be the current primary according to
325 // whether decays have been recorded among the primaries
326 part = dynamic_cast<TParticle*>(fParticleMap.At(fHgwmk+1));
327 fParticleMap.At(part->GetFirstMother())->ResetBit(kDaughtersBit);
328 // Second pass, build map between old and new numbering
329 for(i=fHgwmk+1; i<fNtrack; i++) {
330 if(fParticleMap.At(i)->TestBit(kKeepBit)) {
331 // This particle has to be kept
332 fTrackLabelMap[i]=nkeep;
333 // If old and new are different, have to move the pointer
334 if(i!=nkeep) fParticleMap[nkeep]=fParticleMap.At(i);
335 part = dynamic_cast<TParticle*>(fParticleMap.At(nkeep));
336 // as the parent is always *before*, it must be already
337 // in place. This is what we are checking anyway!
338 if((parent=part->GetFirstMother())>fHgwmk) {
339 if(fTrackLabelMap[parent]==-99) Fatal("PurifyKine","fTrackLabelMap[%d] = -99!\n",parent);
340 else part->SetFirstMother(fTrackLabelMap[parent]);}
345 // Fix daughters information
346 for (i=fHgwmk+1; i<nkeep; i++) {
347 part = dynamic_cast<TParticle*>(fParticleMap.At(i));
348 parent = part->GetFirstMother();
350 father = dynamic_cast<TParticle*>(fParticleMap.At(parent));
351 if(father->TestBit(kDaughtersBit)) {
353 if(i<father->GetFirstDaughter()) father->SetFirstDaughter(i);
354 if(i>father->GetLastDaughter()) father->SetLastDaughter(i);
356 // Initialise daughters info for first pass
357 father->SetFirstDaughter(i);
358 father->SetLastDaughter(i);
359 father->SetBit(kDaughtersBit);
364 // Now the output bit, from fHgwmk to nkeep we write everything and we erase
365 if(nkeep > fParticleFileMap.GetSize()) fParticleFileMap.Set(Int_t (nkeep*1.5));
366 for (i=fHgwmk+1; i<nkeep; ++i) {
367 fParticleBuffer = dynamic_cast<TParticle*>(fParticleMap.At(i));
368 fParticleFileMap[i]=static_cast<Int_t>(TreeK()->GetEntries());
370 fParticleMap[i]=fParticleBuffer=0;
373 for (i = nkeep; i < fNtrack; ++i) fParticleMap[i]=0;
375 Int_t toshrink = fNtrack-fHgwmk-1;
376 fLoadPoint-=toshrink;
378 for(i=fLoadPoint; i<fLoadPoint+toshrink; ++i) fParticles.RemoveAt(i);
385 Bool_t AliStack::ReorderKine()
388 // In some transport code children might not come in a continuous sequence.
389 // In this case the stack has to be reordered in order to establish the
390 // mother daughter relation using index ranges.
392 if(fHgwmk+1 == fNtrack) return kFALSE;
395 // Howmany secondaries have been produced ?
396 Int_t nNew = fNtrack - fHgwmk - 1;
402 // Copy pointers to temporary array
403 TParticle** tmp = new TParticle*[nNew];
405 for (i = 0; i < nNew; i++) {
406 if (fParticleMap.At(fHgwmk + 1 + i)) {
407 tmp[i] = (TParticle*) (fParticleMap.At(fHgwmk + 1 + i));
418 Int_t loadPoint = fHgwmk + 1;
420 // Re-Push particles into stack
421 // The outer loop is over parents, the inner over children.
422 // -1 refers to the primary particle
424 for (i = -1; i < nNew-1; i++) {
428 ipa = tmp[0]->GetFirstMother();
429 parP =dynamic_cast<TParticle*>(fParticleMap.At(ipa));
431 ipa = (fHgwmk + 1 + i);
432 // Skip deleted particles
433 if (!tmp[i]) continue;
434 // Skip particles without children
435 if (tmp[i]->GetFirstDaughter() == -1) continue;
438 // Reset daughter information
440 Int_t idaumin = parP->GetFirstDaughter() - fHgwmk - 1;
441 Int_t idaumax = parP->GetLastDaughter() - fHgwmk - 1;
442 parP->SetFirstDaughter(-1);
443 parP->SetLastDaughter(-1);
444 for (j = idaumin; j <= idaumax; j++) {
445 // Skip deleted particles
446 if (!tmp[j]) continue;
447 // Skip particles already handled
448 if (map1[j] != -99) continue;
449 Int_t jpa = tmp[j]->GetFirstMother();
450 // Check if daughter of current parent
452 fParticleMap[loadPoint] = tmp[j];
453 // Re-establish daughter information
454 parP->SetLastDaughter(loadPoint);
455 if (parP->GetFirstDaughter() == -1) parP->SetFirstDaughter(loadPoint);
456 // Set Mother information
458 tmp[j]->SetFirstMother(map1[i]);
462 // Increase load point
471 // Build map for remapping of hits
473 fTrackLabelMap.Set(fNtrack);
474 for (i = 0; i < fNtrack; i ++) {
476 fTrackLabelMap[i] = i;
478 fTrackLabelMap[i] = map1[i - fHgwmk -1];
481 } // new particles poduced
486 Bool_t AliStack::KeepPhysics(const TParticle* part)
489 // Some particles have to kept on the stack for reasons motivated
490 // by physics analysis. Decision is put here.
492 Bool_t keep = kFALSE;
494 // Keep first-generation daughter from primaries with heavy flavor
496 Int_t parent = part->GetFirstMother();
497 if (parent >= 0 && parent <= fHgwmk) {
498 TParticle* father = dynamic_cast<TParticle*>(Particles()->At(parent));
499 Int_t kf = father->GetPdgCode();
503 if (kfl > 10) kfl/=100;
505 if (kfl > 10) kfl/=10;
506 if (kfl > 10) kfl/=10;
514 //_____________________________________________________________________________
515 void AliStack::FinishEvent()
518 // Write out the kinematics that was not yet filled
521 // Update event header
524 // Fatal("FinishEvent", "No kinematics tree is defined.");
525 // Don't panic this is a probably a lego run
530 if(TreeK()->GetEntries() ==0) {
531 // set the fParticleFileMap size for the first time
532 fParticleFileMap.Set(fHgwmk+1);
535 Bool_t allFilled = kFALSE;
537 for(Int_t i=0; i<fHgwmk+1; ++i)
538 if((part=fParticleMap.At(i))) {
539 fParticleBuffer = dynamic_cast<TParticle*>(part);
540 fParticleFileMap[i]= static_cast<Int_t>(TreeK()->GetEntries());
543 fParticleMap.AddAt(0,i);
545 // When all primaries were filled no particle!=0
546 // should be left => to be removed later.
547 if (allFilled) AliWarning(Form("Why != 0 part # %d?\n",i));
551 // // printf("Why = 0 part # %d?\n",i); => We know.
553 // we don't break now in order to be sure there is no
554 // particle !=0 left.
555 // To be removed later and replaced with break.
556 if(!allFilled) allFilled = kTRUE;
559 //_____________________________________________________________________________
561 void AliStack::FlagTrack(Int_t track)
564 // Flags a track and all its family tree to be kept
571 particle=dynamic_cast<TParticle*>(fParticleMap.At(curr));
573 // If the particle is flagged the three from here upward is saved already
574 if(particle->TestBit(kKeepBit)) return;
576 // Save this particle
577 particle->SetBit(kKeepBit);
579 // Move to father if any
580 if((curr=particle->GetFirstMother())==-1) return;
584 //_____________________________________________________________________________
585 void AliStack::KeepTrack(Int_t track)
588 // Flags a track to be kept
591 fParticleMap.At(track)->SetBit(kKeepBit);
594 //_____________________________________________________________________________
595 void AliStack::Clean(Int_t size)
598 // Reset stack data except for fTreeK
609 //_____________________________________________________________________________
610 void AliStack::Reset(Int_t size)
613 // Reset stack data including fTreeK
617 delete fParticleBuffer; fParticleBuffer = 0;
621 //_____________________________________________________________________________
622 void AliStack::ResetArrays(Int_t size)
625 // Resets stack arrays
628 fParticleMap.Clear();
629 if (size>0) fParticleMap.Expand(size);
632 //_____________________________________________________________________________
633 void AliStack::SetHighWaterMark(Int_t)
636 // Set high water mark for last track in event
640 fCurrentPrimary=fHgwmk;
641 // Set also number of primary tracks
642 fNprimary = fHgwmk+1;
645 //_____________________________________________________________________________
646 TParticle* AliStack::Particle(Int_t i)
649 // Return particle with specified ID
651 if(!fParticleMap.At(i)) {
652 Int_t nentries = fParticles.GetEntriesFast();
653 // algorithmic way of getting entry index
654 // (primary particles are filled after secondaries)
655 Int_t entry = TreeKEntry(i);
656 // check whether algorithmic way and
657 // and the fParticleFileMap[i] give the same;
658 // give the fatal error if not
659 if (entry != fParticleFileMap[i]) {
661 "!! The algorithmic way and map are different: !!\n entry: %d map: %d",
662 entry, fParticleFileMap[i]));
664 // Load particle at entry into fParticleBuffer
665 TreeK()->GetEntry(entry);
666 // Add to the TClonesarray
667 new (fParticles[nentries]) TParticle(*fParticleBuffer);
668 // Store a pointer in the TObjArray
669 fParticleMap.AddAt(fParticles[nentries],i);
671 return dynamic_cast<TParticle*>(fParticleMap.At(i));
674 //_____________________________________________________________________________
675 TParticle* AliStack::ParticleFromTreeK(Int_t id) const
678 // return pointer to TParticle with label id
681 if ((entry = TreeKEntry(id)) < 0) return 0;
682 if (fTreeK->GetEntry(entry)<=0) return 0;
683 return fParticleBuffer;
686 //_____________________________________________________________________________
687 Int_t AliStack::TreeKEntry(Int_t id) const
690 // Return entry number in the TreeK for particle with label id
691 // Return negative number if label>fNtrack
693 // The order of particles in TreeK reflects the order of the transport of primaries and production of secondaries:
695 // Before transport there are fNprimary particles on the stack.
696 // They are transported one by one and secondaries (fNtrack - fNprimary) are produced.
697 // After the transport of each particles secondaries are written to the TreeK
698 // They occupy the entries 0 ... fNtrack - fNprimary - 1
699 // The primaries are written after they have been transported and occupy
700 // fNtrack - fNprimary .. fNtrack - 1
704 entry = id+fNtrack-fNprimary;
706 entry = id-fNprimary;
710 //_____________________________________________________________________________
711 Int_t AliStack::GetCurrentParentTrackNumber() const
714 // Return number of the parent of the current track
717 TParticle* current = (TParticle*)fParticleMap.At(fCurrent);
720 return current->GetFirstMother();
722 AliWarning("Current track not found in the stack");
727 //_____________________________________________________________________________
728 Int_t AliStack::GetPrimary(Int_t id)
731 // Return number of primary that has generated track
739 parent=Particle(current)->GetFirstMother();
740 if(parent<0) return current;
744 //_____________________________________________________________________________
745 void AliStack::DumpPart (Int_t i) const
748 // Dumps particle i in the stack
750 dynamic_cast<TParticle*>(fParticleMap.At(i))->Print();
753 //_____________________________________________________________________________
754 void AliStack::DumpPStack ()
757 // Dumps the particle stack
762 printf("\n\n=======================================================================\n");
763 for (i=0;i<fNtrack;i++)
765 TParticle* particle = Particle(i);
767 printf("-> %d ",i); particle->Print();
768 printf("--------------------------------------------------------------\n");
771 Warning("DumpPStack", "No particle with id %d.", i);
774 printf("\n=======================================================================\n\n");
776 // print particle file map
777 printf("\nParticle file map: \n");
778 for (i=0; i<fNtrack; i++)
779 printf(" %d th entry: %d \n",i,fParticleFileMap[i]);
783 //_____________________________________________________________________________
784 void AliStack::DumpLoadedStack() const
787 // Dumps the particle in the stack
788 // that are loaded in memory.
792 "\n\n=======================================================================\n");
793 for (Int_t i=0;i<fNtrack;i++)
795 TParticle* particle = dynamic_cast<TParticle*>(fParticleMap[i]);
797 printf("-> %d ",i); particle->Print();
798 printf("--------------------------------------------------------------\n");
801 printf("-> %d Particle not loaded.\n",i);
802 printf("--------------------------------------------------------------\n");
806 "\n=======================================================================\n\n");
813 //_____________________________________________________________________________
814 void AliStack::CleanParents()
817 // Clean particles stack
818 // Set parent/daughter relations
823 for(i=0; i<fHgwmk+1; i++) {
824 part = dynamic_cast<TParticle*>(fParticleMap.At(i));
825 if(part) if(!part->TestBit(kDaughtersBit)) {
826 part->SetFirstDaughter(-1);
827 part->SetLastDaughter(-1);
832 //_____________________________________________________________________________
833 TParticle* AliStack::GetNextParticle()
836 // Return next particle from stack of particles
839 TParticle* particle = 0;
841 // search secondaries
842 //for(Int_t i=fNtrack-1; i>=0; i--) {
843 for(Int_t i=fNtrack-1; i>fHgwmk; i--) {
844 particle = dynamic_cast<TParticle*>(fParticleMap.At(i));
845 if ((particle) && (!particle->TestBit(kDoneBit))) {
851 // take next primary if all secondaries were done
852 while (fCurrentPrimary>=0) {
853 fCurrent = fCurrentPrimary;
854 particle = dynamic_cast<TParticle*>(fParticleMap.At(fCurrentPrimary--));
855 if ((particle) && (!particle->TestBit(kDoneBit))) {
860 // nothing to be tracked
866 //__________________________________________________________________________________________
868 void AliStack::ConnectTree(TTree* tree)
871 // Creates branch for writing particles
876 AliDebug(1, "Connecting TreeK");
881 AliFatal("Parameter is NULL");//we don't like such a jokes
884 return;//in this case TreeK() calls back this method (ConnectTree)
885 //tree after setting fTreeK, the rest was already executed
886 //it is safe to return now
889 // Create a branch for particles
891 AliDebug(2, Form("Tree name is %s",fTreeK->GetName()));
893 if (fTreeK->GetDirectory())
895 AliDebug(2, Form("and dir is %s",fTreeK->GetDirectory()->GetName()));
898 AliWarning("DIR IS NOT SET !!!");
900 TBranch *branch=fTreeK->GetBranch("Particles");
903 branch = fTreeK->Branch("Particles", "TParticle", &fParticleBuffer, 4000);
904 AliDebug(2, "Creating Branch in Tree");
908 AliDebug(2, "Branch Found in Tree");
909 branch->SetAddress(&fParticleBuffer);
911 if (branch->GetDirectory())
913 AliDebug(1, Form("Branch Dir Name is %s",branch->GetDirectory()->GetName()));
916 AliWarning("Branch Dir is NOT SET");
919 //_____________________________________________________________________________
921 Bool_t AliStack::GetEvent()
924 // Get new event from TreeK
926 // Reset/Create the particle stack
927 Int_t size = (Int_t)TreeK()->GetEntries();
931 //_____________________________________________________________________________
933 Bool_t AliStack::IsStable(Int_t pdg) const
936 // Decide whether particle (pdg) is stable
940 // All ions/nucleons are considered as stable
941 // Nuclear code is 10LZZZAAAI
942 if(pdg>1000000000)return kTRUE;
944 const Int_t kNstable = 15;
947 Int_t pdgStable[kNstable] = {
949 kElectron, // Electron
957 kLambda0, // Lambda_0
958 kSigmaMinus, // Sigma Minus
959 kSigmaPlus, // Sigma Plus
965 Bool_t isStable = kFALSE;
966 for (i = 0; i < kNstable; i++) {
967 if (pdg == TMath::Abs(pdgStable[i])) {
976 Bool_t AliStack::IsPhysicalPrimary(Int_t index)
979 // Test if a particle is a physical primary according to the following definition:
980 // Particles produced in the collision including products of strong and
981 // electromagnetic decay and excluding feed-down from weak decays of strange
984 TParticle* p = Particle(index);
985 Int_t ist = p->GetStatusCode();
988 // Initial state particle
989 if (ist > 1) return kFALSE;
991 Int_t pdg = TMath::Abs(p->GetPdgCode());
993 if (!IsStable(pdg)) return kFALSE;
995 if (index < GetNprimary()) {
997 // Particle produced by generator
1001 // Particle produced during transport
1004 Int_t imo = p->GetFirstMother();
1005 TParticle* pm = Particle(imo);
1006 Int_t mpdg = TMath::Abs(pm->GetPdgCode());
1007 // Check if it comes from a pi0 decay
1009 // What about the pi0 Dalitz ??
1010 // if ((mpdg == kPi0) && (imo < GetNprimary())) return kTRUE;
1012 // Check if this is a heavy flavor decay product
1013 Int_t mfl = Int_t (mpdg / TMath::Power(10, Int_t(TMath::Log10(mpdg))));
1016 if (mfl < 4) return kFALSE;
1019 // Heavy flavor hadron produced by generator
1020 if (imo < GetNprimary()) {
1024 // To be sure that heavy flavor has not been produced in a secondary interaction
1025 // Loop back to the generated mother
1026 while (imo >= GetNprimary()) {
1027 imo = pm->GetFirstMother();
1030 mpdg = TMath::Abs(pm->GetPdgCode());
1031 mfl = Int_t (mpdg / TMath::Power(10, Int_t(TMath::Log10(mpdg))));
1038 } // produced by generator ?