]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TGeant4/TG4StepManager.cxx
updated to AliMC changes - return type of ProdProcess() changed to AliMCProcess,...
[u/mrichter/AliRoot.git] / TGeant4 / TG4StepManager.cxx
CommitLineData
2817d3e2 1// $Id$
2// Category: event
3//
4// See the class description in the header file.
5
6#include "TG4StepManager.h"
7#include "TG4GeometryManager.h"
8#include "TG4PhysicsManager.h"
9#include "TG4VSensitiveDetector.h"
10#include "TG4Globals.h"
11#include "TG3Units.h"
12
2817d3e2 13#include <G4SteppingManager.hh>
14#include <G4UserLimits.hh>
2817d3e2 15#include <G4UImanager.hh>
16#include <G4AffineTransform.hh>
f2510df1 17#include <G4TransportationManager.hh>
18#include <G4Navigator.hh>
2817d3e2 19
20#include <Randomize.hh>
2817d3e2 21#include <TLorentzVector.h>
22
23TG4StepManager* TG4StepManager::fgInstance = 0;
24
25TG4StepManager::TG4StepManager()
f2510df1 26 : fTrack(0),
27 fStep(0),
28 fStepStatus(kNormalStep),
2817d3e2 29 fSteppingManager(0)
30{
31//
32 if (fgInstance) {
33 TG4Globals::Exception(
34 "TG4StepManager: attempt to create two instances of singleton.");
35 }
36
37 fgInstance = this;
38}
39
40TG4StepManager::TG4StepManager(const TG4StepManager& right) {
41//
42 TG4Globals::Exception(
43 "Attempt to copy TG4StepManager singleton.");
44}
45
46TG4StepManager::~TG4StepManager() {
47//
48}
49
50// operators
51
52TG4StepManager& TG4StepManager::operator=(const TG4StepManager& right)
53{
54 // check assignement to self
55 if (this == &right) return *this;
56
57 TG4Globals::Exception(
58 "Attempt to assign TG4StepManager singleton.");
59
60 return *this;
61}
62
63// private methods
64
f2510df1 65void TG4StepManager::CheckTrack() const
66{
67// Gives exception in case the track is not defined.
68// ---
69
70 if (!fTrack)
71 TG4Globals::Exception("TG4StepManager: Track is not defined.");
72}
73
74
12139627 75void TG4StepManager::CheckStep(const G4String& method) const
f2510df1 76{
77// Gives exception in case the step is not defined.
78// ---
79
12139627 80 if (!fStep) {
81 G4String text = "TG4StepManager::";
82 text = text + method + ": Step is not defined.";
83 TG4Globals::Exception(text);
84 }
f2510df1 85}
86
87
88void TG4StepManager::CheckSteppingManager() const
89{
90// Gives exception in case the step is not defined.
91// ---
92
93 if (!fSteppingManager)
94 TG4Globals::Exception("TG4StepManager: Stepping manager is not defined.");
95}
96
97
2817d3e2 98void TG4StepManager::SetTLorentzVector(G4ThreeVector xyz, G4double t,
99 TLorentzVector& lv) const
100{
101// Fills TLorentzVector with G4ThreeVector and G4double.
102// ---
103
104 lv[0] = xyz.x();
105 lv[1] = xyz.y();
106 lv[2] = xyz.z();
107 lv[3] = t;
108}
109
f2510df1 110G4VPhysicalVolume* TG4StepManager::GetCurrentOffPhysicalVolume(G4int off) const
111{
112// Returns the physical volume of the off-th mother's
113// of the current volume.
114// ---
115
116 G4VPhysicalVolume* physVolume = GetCurrentPhysicalVolume();
117
12139627 118 G4VPhysicalVolume* mother = physVolume;
119
f2510df1 120 Int_t level = off;
121 while (level > 0) {
122 if (mother) mother = mother->GetMother();
123 level--;
124 }
125
126 if (!mother) {
127 G4String text = "TG4StepManager::CurrentVolOff: \n";
128 text = text + " Volume ";
129 text = text + physVolume->GetName();
130 text = text + " has not defined mother in the required level.";
131 TG4Globals::Warning(text);
132 }
133
134 return mother;
135}
136
137G4int TG4StepManager::GetVolumeID(G4VPhysicalVolume* physVolume) const
138{
139// Returns the sensitive detector ID of the specified
140// physical volume.
141// ---
142
143 // sensitive detector ID
144 G4VSensitiveDetector* sd
145 = physVolume->GetLogicalVolume()->GetSensitiveDetector();
146 if (sd) {
147 TG4VSensitiveDetector* tsd = dynamic_cast<TG4VSensitiveDetector*>(sd);
148 if (tsd)
149 return tsd->GetID();
150 else {
151 TG4Globals::Exception(
152 "TG4StepManager::GetVolumeID: Unknown sensitive detector type");
153 return 0;
154 }
155 }
156 else {
a6461796 157 G4String text = "TG4StepManager::GetVolumeID: \n";
f2510df1 158 text = text + " Volume " + physVolume->GetName();
159 text = text + " has not a sensitive detector.";
a6461796 160 //TG4Globals::Exception(text);
161 TG4Globals::Warning(text);
f2510df1 162 return 0;
163 }
164}
165
166
2817d3e2 167// public methods
168
169void TG4StepManager::StopTrack()
170{
171// Stops the current track and skips to the next.
172// ?? do we want to invoke rest processes?
173// ?? do we want to stop secondaries too?
174// possible "stop" track status from G4:
175// fStopButAlive, // Invoke active rest physics processes and
176// // and kill the current track afterward
177// fStopAndKill, // Kill the current track
178// fKillTrackAndSecondaries,
179// // Kill the current track and also associated
180// // secondaries.
181// ---
182
12139627 183#ifdef TGEANT4_DEBUG
f2510df1 184 CheckTrack();
12139627 185#endif
f2510df1 186
187 fTrack->SetTrackStatus(fStopAndKill);
188 // fTrack->SetTrackStatus(fStopButAlive);
189 // fTrack->SetTrackStatus(fKillTrackAndSecondaries);
2817d3e2 190}
191
192void TG4StepManager::StopEvent()
193{
194// Aborts the current event processing.
195// ---
196
12139627 197#ifdef TGEANT4_DEBUG
f2510df1 198 CheckTrack();
12139627 199#endif
f2510df1 200
201 fTrack->SetTrackStatus(fKillTrackAndSecondaries);
202 //StopTrack(); // cannot be used as it keeps secondaries
203 G4UImanager::GetUIpointer()->ApplyCommand("/event/abort");
204 G4UImanager::GetUIpointer()->ApplyCommand("/alStacking/clearStack");
2817d3e2 205}
206
207void TG4StepManager::Rndm(Float_t* array, const Int_t size) const
208{
209// Random numbers array of the specified size.
210// ---
211
5025108a 212 G4double* const kpDoubleArray = new G4double[size];
f2510df1 213 RandFlat::shootArray(size, kpDoubleArray);
2817d3e2 214 for (G4int i=0; i<size; i++) {
5025108a 215 array[i] = kpDoubleArray[i];
2817d3e2 216 }
5025108a 217 delete [] kpDoubleArray;
2817d3e2 218}
219
220void TG4StepManager::SetMaxStep(Float_t step)
221{
222// Maximum step allowed in the current logical volume.
223// The maximum value is kept for following tracks - is it ok ??
224// ---
225
f2510df1 226 // check this
227 G4LogicalVolume* curLogVolume
228 = GetCurrentPhysicalVolume()->GetLogicalVolume();
229 G4UserLimits* userLimits
230 = curLogVolume->GetUserLimits();
231 if (userLimits == 0) {
232 userLimits = new G4UserLimits(step);
233 curLogVolume->SetUserLimits(userLimits);
2817d3e2 234 }
f2510df1 235 else
236 userLimits->SetMaxAllowedStep(step);
2817d3e2 237}
238
239void TG4StepManager::SetMaxNStep(Int_t maxNofSteps)
240{
241// Not yet implemented.
242// ---
243
244 TG4Globals::Warning(
245 "TG4StepManager::SetMaxNStep(..) is not yet implemented.");
246}
247
5025108a 248void TG4StepManager::SetUserDecay(Int_t pdg)
2817d3e2 249{
250// Not yet implemented.
251// ---
252
253 TG4Globals::Exception(
254 "TG4StepManager::SetUserDecay(..) is not yet implemented.");
255}
256
f2510df1 257G4VPhysicalVolume* TG4StepManager::GetCurrentPhysicalVolume() const
258{
259// Returns the current physical volume.
260// According to fStepStatus the volume from track vertex,
261// pre step point or post step point is returned.
262// ---
263
264 G4VPhysicalVolume* physVolume;
265 if (fStepStatus == kNormalStep) {
12139627 266
267#ifdef TGEANT4_DEBUG
268 CheckStep("GetCurrentPhysicalVolume");
269#endif
270
f2510df1 271 physVolume = fStep->GetPreStepPoint()->GetPhysicalVolume();
272 }
273 else if (fStepStatus == kBoundary) {
12139627 274
275#ifdef TGEANT4_DEBUG
276 CheckStep("GetCurrentPhysicalVolume");
277#endif
278
f2510df1 279 physVolume = fStep->GetPostStepPoint()->GetPhysicalVolume();
280 }
281 else {
12139627 282
283#ifdef TGEANT4_DEBUG
f2510df1 284 CheckTrack();
12139627 285#endif
286
f2510df1 287 G4ThreeVector position = fTrack->GetPosition();
288 G4Navigator* navigator =
289 G4TransportationManager::GetTransportationManager()->
290 GetNavigatorForTracking();
291 physVolume
292 = navigator->LocateGlobalPointAndSetup(position);
293 }
294
295 return physVolume;
296}
297
2817d3e2 298Int_t TG4StepManager::CurrentVolID(Int_t& copyNo) const
299{
300// Returns the current sensitive detector ID
301// and the copy number of the current physical volume.
302// ---
303
f2510df1 304 G4VPhysicalVolume* physVolume = GetCurrentPhysicalVolume();
feb86ea5 305 copyNo = physVolume->GetCopyNo() + 1;
2817d3e2 306
f2510df1 307 // sensitive detector ID
308 return GetVolumeID(physVolume);
309}
2817d3e2 310
311Int_t TG4StepManager::CurrentVolOffID(Int_t off, Int_t& copyNo) const
312{
313// Returns the off-th mother's of the current volume
314// the sensitive detector ID and the copy number.
315// ---
316
317 if (off == 0) return CurrentVolID(copyNo);
318
f2510df1 319 G4VPhysicalVolume* mother = GetCurrentOffPhysicalVolume(off);
320
321 if (mother) {
feb86ea5 322 copyNo = mother->GetCopyNo() + 1;
f2510df1 323
324 // sensitive detector ID
325 return GetVolumeID(mother);
2817d3e2 326 }
f2510df1 327 else {
328 copyNo = 0;
329 return 0;
330 }
2817d3e2 331}
332
333const char* TG4StepManager::CurrentVolName() const
334{
335// Returns the current physical volume name.
336// ---
337
f2510df1 338 return GetCurrentPhysicalVolume()->GetName();
2817d3e2 339}
340
341const char* TG4StepManager::CurrentVolOffName(Int_t off) const
342{
343// Returns the off-th mother's physical volume name.
344// ---
345
f2510df1 346 if (off == 0) return CurrentVolName();
2817d3e2 347
f2510df1 348 G4VPhysicalVolume* mother = GetCurrentOffPhysicalVolume(off);
349
350 if (mother)
351 return mother->GetName();
352 else
2817d3e2 353 return 0;
2817d3e2 354}
355
356Int_t TG4StepManager::CurrentMaterial(Float_t &a, Float_t &z, Float_t &dens,
357 Float_t &radl, Float_t &absl) const
358{
359// Returns the parameters of the current material during transport
360// the return value is the number of elements in the mixture.
361// ---
362
f2510df1 363 G4VPhysicalVolume* physVolume = GetCurrentPhysicalVolume();
364
365 G4Material* material
366 = physVolume->GetLogicalVolume()->GetMaterial();
367
12139627 368 if (material) {
f2510df1 369 G4int nofElements = material->GetNumberOfElements();
370 TG4GeometryManager* pGeometryManager = TG4GeometryManager::Instance();
371 a = pGeometryManager->GetEffA(material);
372 z = pGeometryManager->GetEffZ(material);
2817d3e2 373
f2510df1 374 // density
375 dens = material->GetDensity();
376 dens /= TG3Units::MassDensity();
2817d3e2 377
f2510df1 378 // radiation length
379 radl = material->GetRadlen();
380 radl /= TG3Units::Length();
2817d3e2 381
f2510df1 382 absl = 0.; // this parameter is not defined in Geant4
383 return nofElements;
384 }
385 else {
386 TG4Globals::Exception(
387 "TG4StepManager::CurrentMaterial(..): material is not defined.");
2817d3e2 388 return 0;
389 }
390}
391
392void TG4StepManager::Gmtod(Float_t* xm, Float_t* xd, Int_t iflag)
393{
394// Transforms a position from the world reference frame
395// to the current volume reference frame.
396//
397// Geant3 desription:
398// ==================
399// Computes coordinates XD (in DRS)
400// from known coordinates XM in MRS
401// The local reference system can be initialized by
402// - the tracking routines and GMTOD used in GUSTEP
403// - a call to GMEDIA(XM,NUMED)
404// - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
405// (inverse routine is GDTOM)
406//
407// If IFLAG=1 convert coordinates
408// IFLAG=2 convert direction cosinus
409//
410// ---
411
12139627 412 G4AffineTransform affineTransform;
f2510df1 413
12139627 414 if (fStepStatus == kVertex) {
415 G4Navigator* navigator =
416 G4TransportationManager::GetTransportationManager()->
417 GetNavigatorForTracking();
418
419 affineTransform = navigator->GetGlobalToLocalTransform();
420 }
421 else {
f2510df1 422
12139627 423#ifdef TGEANT4_DEBUG
424 CheckStep("Gmtod");
425#endif
426
427 affineTransform
428 = fStep->GetPreStepPoint()->GetTouchable()->GetHistory()
429 ->GetTopTransform();
430 }
431
432 G4ThreeVector theGlobalPoint(xm[0],xm[1],xm[2]);
f2510df1 433 G4ThreeVector theLocalPoint;
434 if(iflag == 1)
435 theLocalPoint = affineTransform.TransformPoint(theGlobalPoint);
436 else if ( iflag == 2)
437 theLocalPoint = affineTransform.TransformAxis(theGlobalPoint);
438 else
439 TG4Globals::Exception(
440 "TG4StepManager::Gmtod(..,iflag): iflag is not in 1..2");
2817d3e2 441
f2510df1 442 xd[0] = theLocalPoint.x();
443 xd[1] = theLocalPoint.y();
444 xd[2] = theLocalPoint.z();
2817d3e2 445}
446
447void TG4StepManager::Gdtom(Float_t* xd, Float_t* xm, Int_t iflag)
448{
449// Transforms a position from the current volume reference frame
450// to the world reference frame.
451//
452// Geant3 desription:
453// ==================
454// Computes coordinates XM (Master Reference System
455// knowing the coordinates XD (Detector Ref System)
456// The local reference system can be initialized by
457// - the tracking routines and GDTOM used in GUSTEP
458// - a call to GSCMED(NLEVEL,NAMES,NUMBER)
459// (inverse routine is GMTOD)
460//
461// If IFLAG=1 convert coordinates
462// IFLAG=2 convert direction cosinus
463//
464// ---
465
12139627 466 G4AffineTransform affineTransform;
f2510df1 467
12139627 468 if (fStepStatus == kVertex) {
469 G4Navigator* navigator =
470 G4TransportationManager::GetTransportationManager()->
471 GetNavigatorForTracking();
472
473 affineTransform = navigator->GetLocalToGlobalTransform();
474 }
475 else {
2817d3e2 476
12139627 477#ifdef TGEANT4_DEBUG
478 CheckStep("Gdtom");
479#endif
480
481 // check this
482
483 affineTransform
484 = fStep->GetPreStepPoint()->GetTouchable()->GetHistory()
2817d3e2 485 ->GetTopTransform().Inverse();
12139627 486 }
2817d3e2 487
12139627 488 G4ThreeVector theLocalPoint(xd[0],xd[1],xd[2]);
f2510df1 489 G4ThreeVector theGlobalPoint;
490 if(iflag == 1)
491 theGlobalPoint = affineTransform.TransformPoint(theLocalPoint);
492 else if( iflag == 2)
493 theGlobalPoint = affineTransform.TransformAxis(theLocalPoint);
494 else
495 TG4Globals::Warning(
496 "TG4StepManager::Gdtom(...,iflag): iflag is not in 1..2");
497
498 xm[0] = theGlobalPoint.x();
499 xm[1] = theGlobalPoint.y();
500 xm[2] = theGlobalPoint.z();
2817d3e2 501}
502
503Float_t TG4StepManager::MaxStep() const
504{
505// Returns maximum step allowed in the current logical volume
506// by User Limits.
507// ---
508
12139627 509 G4LogicalVolume* curLogVolume
f2510df1 510 = GetCurrentPhysicalVolume()->GetLogicalVolume();
12139627 511
512 // check this
f2510df1 513 G4UserLimits* userLimits
514 = curLogVolume->GetUserLimits();
515
516 G4double maxStep;
517 if (userLimits == 0) {
518 G4String text = "User Limits are not defined for log volume ";
519 text = text + curLogVolume->GetName();
520 TG4Globals::Warning(text);
2acb5b6d 521 return FLT_MAX;
2817d3e2 522 }
f2510df1 523 else {
524 const G4Track& trackRef = *(fTrack);
525 maxStep = userLimits->GetMaxAllowedStep(trackRef);
526 maxStep /= TG3Units::Length();
527 return maxStep;
528 }
2817d3e2 529}
530
531Int_t TG4StepManager::GetMaxNStep() const
532{
533// Not yet implemented.
534// ---
535
536 TG4Globals::Warning(
537 "Method GetMaxNStep is not yet implemented in TG4StepManager.");
538 return 0;
539}
540
541void TG4StepManager::TrackPosition(TLorentzVector& position) const
542{
543// Current particle position (in the world reference frame)
544// and the local time since the current track is created
545// (position of the PostStepPoint).
546// ---
547
12139627 548#ifdef TGEANT4_DEBUG
f2510df1 549 CheckTrack();
12139627 550#endif
551
f2510df1 552 // get position
553 // check if this is == to PostStepPoint position !!
554 G4ThreeVector positionVector = fTrack->GetPosition();
555 positionVector *= 1./(TG3Units::Length());
2817d3e2 556
f2510df1 557 // local time
558 G4double time = fTrack->GetLocalTime();
559 time /= TG3Units::Time();
560
561 SetTLorentzVector(positionVector, time, position);
2817d3e2 562}
563
564Int_t TG4StepManager::GetMedium() const
565{
566// Returns the second index of the current material (corresponding to
567// G3 tracking medium index).
568// ---
569
f2510df1 570 // current material
571 G4Material* curMaterial
572 = GetCurrentPhysicalVolume()->GetLogicalVolume()->GetMaterial();
2817d3e2 573
f2510df1 574 // medium index
575 TG4GeometryManager* pGeometryManager = TG4GeometryManager::Instance();
576 return pGeometryManager->GetMediumId(curMaterial);
2817d3e2 577}
578
579void TG4StepManager::TrackMomentum(TLorentzVector& momentum) const
580{
581// Current particle "momentum" (px, py, pz, Etot).
582// ---
583
12139627 584#ifdef TGEANT4_DEBUG
f2510df1 585 CheckTrack();
12139627 586#endif
2817d3e2 587
f2510df1 588 G4ThreeVector momentumVector = fTrack->GetMomentum();
589 momentumVector *= 1./(TG3Units::Energy());
2817d3e2 590
f2510df1 591 G4double energy = fTrack->GetDynamicParticle()->GetTotalEnergy();
592 energy /= TG3Units::Energy();
593
594 SetTLorentzVector(momentumVector, energy, momentum);
2817d3e2 595}
596
597void TG4StepManager::TrackVertexPosition(TLorentzVector& position) const
598{
599// The vertex particle position (in the world reference frame)
600// and the local time since the current track is created.
601// ---
602
12139627 603#ifdef TGEANT4_DEBUG
f2510df1 604 CheckTrack();
12139627 605#endif
f2510df1 606
607 // position
608 G4ThreeVector positionVector = fTrack->GetVertexPosition();
609 positionVector *= 1./(TG3Units::Length());
2817d3e2 610
f2510df1 611 // local time
612 // to be checked
613 G4double time = fTrack->GetLocalTime();
614 time /= TG3Units::Time();
2817d3e2 615
f2510df1 616 SetTLorentzVector(positionVector, time, position);
2817d3e2 617}
618
619void TG4StepManager::TrackVertexMomentum(TLorentzVector& momentum) const
620{
621// The vertex particle "momentum" (px, py, pz, Ekin)
622// to do: change Ekin -> Etot
623// ---
624
12139627 625#ifdef TGEANT4_DEBUG
f2510df1 626 CheckTrack();
12139627 627#endif
628
f2510df1 629 G4ThreeVector momentumVector = fTrack->GetVertexMomentumDirection();
630 momentumVector *= 1./(TG3Units::Energy());
2817d3e2 631
f2510df1 632 G4double energy = fTrack->GetVertexKineticEnergy();
633 energy /= TG3Units::Energy();
2817d3e2 634
f2510df1 635 SetTLorentzVector(momentumVector, energy, momentum);
2817d3e2 636}
637
638Float_t TG4StepManager::TrackStep() const
639{
640// Returns the current step length.
641// ---
642
f2510df1 643 G4double length;
644 if (fStepStatus == kNormalStep) {
12139627 645
646#ifdef TGEANT4_DEBUG
647 CheckStep("TrackStep");
648#endif
649
f2510df1 650 length = fStep->GetStepLength();
651 length /= TG3Units::Length();
652 }
653 else
654 length = 0;
655
656 return length;
2817d3e2 657}
658
659Float_t TG4StepManager::TrackLength() const
660{
661// Returns the length of the current track from its origin.
662// ---
663
12139627 664#ifdef TGEANT4_DEBUG
f2510df1 665 CheckTrack();
12139627 666#endif
f2510df1 667
668 G4double length = fTrack->GetTrackLength();
669 length /= TG3Units::Length();
670 return length;
2817d3e2 671}
672
673Float_t TG4StepManager::TrackTime() const
674{
675// Returns the local time since the current track is created.
676// Comment:
677// in Geant4: there is also defined proper time as
678// the proper time of the dynamical particle of the current track.
679// ---
680
12139627 681#ifdef TGEANT4_DEBUG
f2510df1 682 CheckTrack();
12139627 683#endif
f2510df1 684
685 G4double time = fTrack->GetLocalTime();
686 time /= TG3Units::Time();
687 return time;
2817d3e2 688}
689
690Float_t TG4StepManager::Edep() const
691{
692// Returns total energy deposit in this step.
693// ---
694
f2510df1 695 G4double energyDeposit;
696 if (fStepStatus == kNormalStep) {
12139627 697
698#ifdef TGEANT4_DEBUG
699 CheckStep("Edep");
700#endif
701
f2510df1 702 energyDeposit = fStep->GetTotalEnergyDeposit();
703 energyDeposit /= TG3Units::Energy();
2817d3e2 704 }
f2510df1 705 else
706 energyDeposit = 0;
707
708 return energyDeposit;
2817d3e2 709}
710
711Int_t TG4StepManager::TrackPid() const
712{
713// Returns the current particle PDG encoding.
714// ---
715
12139627 716#ifdef TGEANT4_DEBUG
f2510df1 717 CheckTrack();
12139627 718#endif
2817d3e2 719
f2510df1 720 G4ParticleDefinition* particle
721 = fTrack->GetDynamicParticle()->GetDefinition();
722
723 // ask TG4PhysicsManager to get PDG encoding
724 // (in order to get PDG from extended TDatabasePDG
725 // in case the standard PDG code is not defined)
726 TG4PhysicsManager* pPhysicsManager = TG4PhysicsManager::Instance();
727 G4int pdgEncoding = pPhysicsManager->GetPDGEncodingFast(particle);
728
729 return pdgEncoding;
2817d3e2 730}
731
732Float_t TG4StepManager::TrackCharge() const
733{
734// Returns the current particle charge.
735// ---
736
12139627 737#ifdef TGEANT4_DEBUG
f2510df1 738 CheckTrack();
12139627 739#endif
740
f2510df1 741 G4double charge
742 = fTrack->GetDynamicParticle()->GetDefinition()
743 ->GetPDGCharge();
744 charge /= TG3Units::Charge();
745 return charge;
2817d3e2 746}
747
748Float_t TG4StepManager::TrackMass() const
749{
750// Returns current particle rest mass.
751// ---
752
12139627 753#ifdef TGEANT4_DEBUG
f2510df1 754 CheckTrack();
12139627 755#endif
f2510df1 756
757 G4double mass
758 = fTrack->GetDynamicParticle()->GetDefinition()
759 ->GetPDGMass();
760 mass /= TG3Units::Mass();
761 return mass;
2817d3e2 762}
763
764Float_t TG4StepManager::Etot() const
765{
766// Returns total energy of the current particle.
767// ---
768
12139627 769#ifdef TGEANT4_DEBUG
f2510df1 770 CheckTrack();
12139627 771#endif
f2510df1 772
773 G4double energy
774 = fTrack->GetDynamicParticle()->GetTotalEnergy();
775 energy /= TG3Units::Energy();
776 return energy;
2817d3e2 777}
778
779Bool_t TG4StepManager::IsTrackInside() const
780{
781// Returns true if particle does not cross geometrical boundary
f2510df1 782// and is not in vertex.
2817d3e2 783// ---
784
f2510df1 785 if (fStepStatus == kNormalStep && !(IsTrackExiting()) ) {
786 // track is always inside during a normal step
787 return true;
788 }
789
790 return false;
2817d3e2 791}
792
793Bool_t TG4StepManager::IsTrackEntering() const
794{
795// Returns true if particle cross a geometrical boundary
f2510df1 796// or is in vertex.
2817d3e2 797// ---
798
f2510df1 799 if (fStepStatus != kNormalStep) {
800 // track is entering during a vertex or boundary step
801 return true;
802 }
803
804 return false;
2817d3e2 805}
806
807Bool_t TG4StepManager::IsTrackExiting() const
808{
f2510df1 809// Returns true if particle cross a geometrical boundary.
2817d3e2 810// ---
811
f2510df1 812 if (fStepStatus == kNormalStep) {
12139627 813
814#ifdef TGEANT4_DEBUG
815 CheckStep("IsTrackExiting");
816#endif
817
f2510df1 818 if (fStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary)
819 return true;
2817d3e2 820 }
f2510df1 821
822 return false;
2817d3e2 823}
824
825Bool_t TG4StepManager::IsTrackOut() const
826{
68f491b7 827// Returns true if particle cross the world boundary
2817d3e2 828// at post-step point.
829// ---
830
f2510df1 831 if (fStepStatus == kVertex) return false;
832
12139627 833#ifdef TGEANT4_DEBUG
834 CheckStep("IsTrackCut");
835#endif
f2510df1 836
837 // check
838 G4StepStatus status
839 = fStep->GetPostStepPoint()->GetStepStatus();
840 if (status == fWorldBoundary)
841 return true;
842 else
2817d3e2 843 return false;
2817d3e2 844}
845
846Bool_t TG4StepManager::IsTrackStop() const
847{
848// Returns true if particle has stopped
849// or has been killed, suspended or postponed to next event.
850//
851// Possible track status from G4:
852// fAlive, // Continue the tracking
853// fStopButAlive, // Invoke active rest physics processes and
854// // and kill the current track afterward
855// fStopAndKill, // Kill the current track
856// fKillTrackAndSecondaries,
857// // Kill the current track and also associated
858// // secondaries.
859// fSuspend, // Suspend the current track
860// fPostponeToNextEvent
861// // Postpones the tracking of thecurrent track
862// // to the next event.
863// ---
864
12139627 865#ifdef TGEANT4_DEBUG
f2510df1 866 CheckTrack();
12139627 867#endif
f2510df1 868
869 // check
870 G4TrackStatus status
871 = fTrack->GetTrackStatus();
872 if ((status == fStopAndKill) ||
873 (status == fKillTrackAndSecondaries) ||
874 (status == fSuspend) ||
875 (status == fPostponeToNextEvent)) {
876 return true;
2817d3e2 877 }
f2510df1 878 else
879 return false;
2817d3e2 880}
881
882Bool_t TG4StepManager::IsTrackDisappeared() const
883{
884// Returns true if particle has disappeared
885// (due to any physical process)
886// or has been killed, suspended or postponed to next event.
887// ---
888
12139627 889#ifdef TGEANT4_DEBUG
f2510df1 890 CheckTrack();
12139627 891#endif
f2510df1 892
893 // check
894 G4TrackStatus status
895 = fTrack->GetTrackStatus();
896 if ((status == fStopButAlive) ||
897 (status == fKillTrackAndSecondaries) ||
898 (status == fSuspend) ||
899 (status == fPostponeToNextEvent)) {
900 return true;
2817d3e2 901 }
f2510df1 902 else
903 return false;
2817d3e2 904}
905
906Bool_t TG4StepManager::IsTrackAlive() const
907{
908// Returns true if particle continues tracking.
909// ---
910
12139627 911#ifdef TGEANT4_DEBUG
f2510df1 912 CheckTrack();
12139627 913#endif
f2510df1 914
915 G4TrackStatus status
916 = fTrack->GetTrackStatus();
917 if (status == fAlive)
918 return true;
919 else
920 return false;
2817d3e2 921}
922
923Bool_t TG4StepManager::IsNewTrack() const
924{
925// Returns true when track performs the first step.
926// ---
927
f2510df1 928 if (fStepStatus == kVertex)
2817d3e2 929 return true;
930 else
931 return false;
932}
933
934Int_t TG4StepManager::NSecondaries() const
935{
936// Returns the number of secondary particles generated
937// in the current step.
938// ---
939
12139627 940#ifdef TGEANT4_DEBUG
f2510df1 941 CheckSteppingManager();
12139627 942#endif
2817d3e2 943
f2510df1 944 G4int nofSecondaries = 0;
945 nofSecondaries += fSteppingManager->GetfN2ndariesAtRestDoIt();
946 nofSecondaries += fSteppingManager->GetfN2ndariesAlongStepDoIt();
947 nofSecondaries += fSteppingManager->GetfN2ndariesPostStepDoIt();
948
949 return nofSecondaries;
2817d3e2 950}
951
952void TG4StepManager::GetSecondary(Int_t index, Int_t& particleId,
953 TLorentzVector& position, TLorentzVector& momentum)
954{
955// Fills the parameters (particle encoding, position, momentum)
956// of the generated secondary particle which is specified by index.
957// !! Check if indexing of secondaries is same !!
958// ---
959
12139627 960#ifdef TGEANT4_DEBUG
f2510df1 961 CheckSteppingManager();
12139627 962#endif
f2510df1 963
964 G4int nofSecondaries = NSecondaries();
965 G4TrackVector* secondaryTracks = fSteppingManager->GetSecondary();
966
967 if (secondaryTracks){
968 if (index < nofSecondaries) {
969
970 // the index of the first secondary of this step
971 G4int startIndex
972 = secondaryTracks->entries() - nofSecondaries;
973 // (the secondaryTracks vector contains secondaries
974 // produced by the track at previous steps, too)
975 G4Track* track
976 = (*secondaryTracks)[startIndex + index];
2817d3e2 977
f2510df1 978 // particle encoding
979 particleId
980 = track->GetDynamicParticle()->GetDefinition()->GetPDGEncoding();
2817d3e2 981
f2510df1 982 // position & time
983 G4ThreeVector positionVector = track->GetPosition();
984 positionVector *= 1./(TG3Units::Length());
985 G4double time = track->GetLocalTime();
986 time /= TG3Units::Time();
987 SetTLorentzVector(positionVector, time, position);
988
989 // momentum & energy
990 G4ThreeVector momentumVector = track->GetMomentum();
991 G4double energy = track->GetDynamicParticle()->GetTotalEnergy();
992 energy /= TG3Units::Energy();
993 SetTLorentzVector(momentumVector, energy, momentum);
2817d3e2 994 }
995 else {
996 TG4Globals::Exception(
f2510df1 997 "TG4StepManager::GetSecondary(): wrong secondary track index.");
2817d3e2 998 }
999 }
f2510df1 1000 else {
1001 TG4Globals::Exception(
1002 "TG4StepManager::GetSecondary(): secondary tracks vector is empty");
2817d3e2 1003 }
1004}
1005
7bfb4926 1006AliMCProcess TG4StepManager::ProdProcess() const
2817d3e2 1007{
7bfb4926 1008// The process that has produced the secondary particles
1009// in the current step
2817d3e2 1010// ---
1011
12139627 1012 G4int nofSecondaries = NSecondaries();
7bfb4926 1013 if (fStepStatus == kVertex || !nofSecondaries) return kPNoProcess;
1014
12139627 1015#ifdef TGEANT4_DEBUG
1016 CheckStep("ProdProcess");
1017#endif
2817d3e2 1018
7bfb4926 1019/*
f2510df1 1020 const G4VProcess* curProcess
7bfb4926 1021 = fStep->GetPostStepPoint()->GetProcessDefinedStep();
1022*/
1023
1024 G4TrackVector* secondaryTracks = fSteppingManager->GetSecondary();
1025
1026 // should never happen
1027 if (!secondaryTracks) {
1028 TG4Globals::Exception(
1029 "TG4StepManager::ProdProcess(): secondary tracks vector is empty.");
1030 }
f2510df1 1031
7bfb4926 1032 G4Track* firstSecTrack = (*secondaryTracks)[0];
1033 G4Track* lastSecTrack = (*secondaryTracks)[secondaryTracks->entries()-1];
1034
1035 const G4VProcess* firstProcess = firstSecTrack->GetCreatorProcess();
1036 const G4VProcess* lastProcess = lastSecTrack->GetCreatorProcess();
1037
1038 // check if all secondaries comes from the same process
1039 if (firstProcess != lastProcess) {
1040 G4String text = "TG4StepManager: ProdProcess():\n";
1041 text = text + " More than one process has created secondary particles.\n";
1042 text = text + " The creator of the last secondary is returned.";
1043 G4cout << "The first process:" << firstProcess->GetProcessName();
1044 TG4Globals::Warning(text);
1045 }
1046
1047 G4String g4Name = lastProcess->GetProcessName();
1048
1049 TG4PhysicsManager* pPhysicsManager = TG4PhysicsManager::Instance();
1050 AliMCProcess mcProcess = pPhysicsManager->GetMCProcess(g4Name);
1051
1052 // distinguish kPDeltaRay from kPEnergyLoss
1053 if (mcProcess == kPEnergyLoss) mcProcess = kPDeltaRay;
1054
1055 return mcProcess;
2817d3e2 1056}