Geant4 VMC Version 6.6
Loading...
Searching...
No Matches
TG4StepManager.cxx
Go to the documentation of this file.
1//------------------------------------------------
2// The Geant4 Virtual Monte Carlo package
3// Copyright (C) 2007 - 2015 Ivana Hrivnacova
4// All rights reserved.
5//
6// For the licensing terms see geant4_vmc/LICENSE.
7// Contact: root-vmc@cern.ch
8//-------------------------------------------------
9
14
15#include "TG4StepManager.h"
16#include "TG4G3Units.h"
17#include "TG4GeometryServices.h"
18#include "TG4Globals.h"
19#include "TG4Limits.h"
20#include "TG4ParticlesManager.h"
21#include "TG4PhysicsManager.h"
22#include "TG4SDServices.h"
23#include "TG4SteppingAction.h"
24#include "TG4TrackInformation.h"
25#include "TG4TrackManager.h"
26
27#include <G4AffineTransform.hh>
28#include <G4Navigator.hh>
29#include <G4OpticalPhoton.hh>
30#include <G4ProcessManager.hh>
31#include <G4ProcessVector.hh>
32#include <G4SteppingManager.hh>
33#include <G4TransportationManager.hh>
34#include <G4TransportationProcessType.hh>
35#include <G4UImanager.hh>
36#include <G4UserLimits.hh>
37#include <G4VProcess.hh>
38#include <G4VTouchable.hh>
39
40#include <TLorentzVector.h>
41#include <TMCParticleStatus.h>
42#include <TMath.h>
43#include <TVector3.h>
44
46
47//_____________________________________________________________________________
48TG4StepManager::TG4StepManager(const TString& userGeometry)
49 : fTrack(0),
50 fStep(0),
51 fGflashSpot(0),
52 fStepStatus(kNormalStep),
53 fLimitsModifiedOnFly(0),
54 fSteppingManager(0),
55 fNameBuffer(),
56 fCopyNoOffset(0),
57 fDivisionCopyNoOffset(0),
58 fTrackManager(0),
59 fInitialVMCTrackStatus(0)
60{
63
64 // G4cout << "TG4StepManager::TG4StepManager " << this << G4endl;
65
66 if (fgInstance) {
67 TG4Globals::Exception("TG4StepManager", "TG4StepManager",
68 "Cannot create two instances of singleton.");
69 }
70
71 fgInstance = this;
72
75 if (userGeometry == "VMCtoGeant4") fCopyNoOffset = 1;
76
80 if (userGeometry == "RootToGeant4" || userGeometry == "Geant4")
82}
83
84//_____________________________________________________________________________
91
92//
93// private methods
94//
95
96//_____________________________________________________________________________
98{
100
101 if (!fTrack)
103 "TG4StepManager", "CheckTrack", "Track is not defined.");
104}
105
106//_____________________________________________________________________________
107void TG4StepManager::CheckStep(const G4String& method) const
108{
110
111 if (!fStep) {
112 TG4Globals::Exception("TG4StepManager", method, "Step is not defined.");
113 }
114}
115
116//_____________________________________________________________________________
117void TG4StepManager::CheckGflashSpot(const G4String& method) const
118{
120
121 if (!fGflashSpot) {
123 "TG4StepManager", method, "Gflash spot is not defined.");
124 }
125}
126
127//_____________________________________________________________________________
129{
131
132 if (!fSteppingManager)
133 TG4Globals::Exception("TG4StepManager", "CheckSteppingManager",
134 "Stepping manager is not defined.");
135}
136
137//_____________________________________________________________________________
139 G4ThreeVector xyz, G4double t, TLorentzVector& lv) const
140{
142
143 lv[0] = xyz.x();
144 lv[1] = xyz.y();
145 lv[2] = xyz.z();
146 lv[3] = t;
147}
148
149//_____________________________________________________________________________
150const G4VTouchable* TG4StepManager::GetCurrentTouchable() const
151{
153
154#ifdef MCDEBUG
155 CheckTrack();
156#endif
157
158 if (fStepStatus == kGflashSpot) {
159 G4ReferenceCountedHandle<G4VTouchable> touchableHandle =
160 fGflashSpot->GetTouchableHandle();
161 return touchableHandle();
162 }
163 else if (fStepStatus != kBoundary)
164 return fTrack->GetTouchable();
165 else
166 return fTrack->GetNextTouchable();
167}
168
169//_____________________________________________________________________________
171 G4int off, G4bool warn) const
172{
175
176 // Get current touchable
177 //
178 const G4VTouchable* touchable = GetCurrentTouchable();
179
180 // Check touchable depth
181 //
182 if (touchable->GetHistoryDepth() < off) {
183 if (warn) {
184 TString text = "level=";
185 text += off;
186 TG4Globals::Warning("TG4StepManager", "GetCurrentOffPhysicalVolume",
187 "Volume " + TString(touchable->GetVolume()->GetName()) +
188 " has not defined mother in " + text + ".");
189 }
190 return 0;
191 }
192
193 return touchable->GetVolume(off);
194}
195
196//
197// public methods
198//
199
200//_____________________________________________________________________________
205
206//_____________________________________________________________________________
208{
216
217 if (fTrack) {
218 fTrack->SetTrackStatus(fStopAndKill);
219 // fTrack->SetTrackStatus(fStopButAlive);
220 // fTrack->SetTrackStatus(fKillTrackAndSecondaries);
221 }
222 else {
223 TG4Globals::Warning("TG4StepManager", "StopTrack()",
224 "There is no current track to be stopped.");
225 }
226}
227
228//_____________________________________________________________________________
230{
232
233 if (fTrack) {
234 fTrack->SetTrackStatus(fStopAndKill);
236 }
237 else {
238 TG4Globals::Warning("TG4StepManager", "InterruptTrack()",
239 "There is no current track to be interrupted.");
240 }
241}
242
243//_____________________________________________________________________________
245{
247
248 if (fTrack) {
249 fTrack->SetTrackStatus(fKillTrackAndSecondaries);
250 // StopTrack(); // cannot be used as it keeps secondaries
251 }
252
253 G4UImanager::GetUIpointer()->ApplyCommand("/event/abort");
254}
255
256//_____________________________________________________________________________
258{
260
262
263 StopEvent();
264 G4UImanager::GetUIpointer()->ApplyCommand("/run/abort");
265}
266
267//_____________________________________________________________________________
268void TG4StepManager::SetMaxStep(Double_t step)
269{
273
274 TG4Limits* userLimits = GetCurrentLimits();
275
276 if (!userLimits) return;
277
278 // G4cout << "TG4StepManager::SetMaxStep in "
279 // << GetCurrentPhysicalVolume()->GetLogicalVolume()->GetName() << " "
280 // << userLimits->GetName() << G4endl;
281
282 // set max step
283 userLimits->SetCurrentMaxAllowedStep(step * TG4G3Units::Length());
284 fLimitsModifiedOnFly = userLimits;
285}
286
287//_____________________________________________________________________________
289{
292
295 "TG4StepManager", "SetMaxStepBack", "No limits modified on fly found.");
296 return;
297 }
298
299 // set max step
302}
303
304//_____________________________________________________________________________
305void TG4StepManager::SetMaxNStep(Int_t maxNofSteps)
306{
308
309 TG4SteppingAction::Instance()->SetMaxNofSteps(TMath::Abs(maxNofSteps));
310}
311
312//_____________________________________________________________________________
313void TG4StepManager::SetCollectTracks(Bool_t collectTracks)
314{
316
318}
319
320//_____________________________________________________________________________
322{
324
325#ifdef MCDEBUG
326 CheckTrack();
327#endif
328
329 G4ParticleDefinition* particle =
330 fTrack->GetDynamicParticle()->GetDefinition();
331
332 // Store the original particle lifetime in track information
333 // (as it has to be set back after track is finished)
334 TG4TrackInformation* trackInformation =
336 trackInformation->SetPDGLifetime(particle->GetPDGLifeTime());
337
338 // Set new lifetime value
339 particle->SetPDGLifeTime(time * TG4G3Units::Time());
340}
341
342//_____________________________________________________________________________
343void TG4StepManager::SetInitialVMCTrackStatus(TMCParticleStatus* status)
344{
347
348 fInitialVMCTrackStatus = status;
349}
350
351//_____________________________________________________________________________
358
359//_____________________________________________________________________________
361{
365
366#ifdef MCDEBUG
367 CheckTrack();
368#endif
369
371 return fGflashSpot->GetTouchableHandle()->GetVolume();
372 else if (fStepStatus != kBoundary)
373 return fTrack->GetVolume();
374 else
375 return fTrack->GetNextVolume();
376}
377
378//_____________________________________________________________________________
380{
382
383#ifdef MCDEBUG
385 GetCurrentPhysicalVolume()->GetLogicalVolume()->GetUserLimits());
386#else
387 TG4Limits* userLimits =
388 (TG4Limits*)GetCurrentPhysicalVolume()->GetLogicalVolume()->GetUserLimits();
389#endif
390
391 if (!userLimits) {
393 "TG4StepManager", "Get current limits", "User limits not defined.");
394 return 0;
395 }
396
397 return userLimits;
398}
399
400//_____________________________________________________________________________
401Int_t TG4StepManager::CurrentVolID(Int_t& copyNo) const
402{
405
406 G4VPhysicalVolume* physVolume = GetCurrentPhysicalVolume();
407 if (!physVolume) {
409 "TG4StepManager", "CurrentVolID", "No current physical volume found");
410 return 0;
411 }
412 copyNo = physVolume->GetCopyNo() + fCopyNoOffset;
413
414 if (physVolume->IsParameterised() || physVolume->IsReplicated())
415 copyNo += fDivisionCopyNoOffset;
416
417 // sensitive detector ID
418 return TG4SDServices::Instance()->GetVolumeID(physVolume->GetLogicalVolume());
419}
420
421//_____________________________________________________________________________
422Int_t TG4StepManager::CurrentVolOffID(Int_t off, Int_t& copyNo) const
423{
426
427 if (off == 0) return CurrentVolID(copyNo);
428#ifdef MCDEBUG
429 G4VPhysicalVolume* mother = GetCurrentOffPhysicalVolume(off, true);
430#else
431 G4VPhysicalVolume* mother = GetCurrentOffPhysicalVolume(off);
432#endif
433
434 if (mother) {
435 copyNo = mother->GetCopyNo() + fCopyNoOffset;
436
437 if (mother->IsParameterised() || mother->IsReplicated())
438 copyNo += fDivisionCopyNoOffset;
439
440 // sensitive detector ID
441 return TG4SDServices::Instance()->GetVolumeID(mother->GetLogicalVolume());
442 }
443 else {
444 copyNo = 0;
445 return 0;
446 }
447}
448
449//_____________________________________________________________________________
451{
453
455 GetCurrentPhysicalVolume()->GetLogicalVolume()->GetName());
456
457 return fNameBuffer.data();
458}
459
460//_____________________________________________________________________________
461const char* TG4StepManager::CurrentVolOffName(Int_t off) const
462{
464
465 if (off == 0) return CurrentVolName();
466
467 G4VPhysicalVolume* mother = GetCurrentOffPhysicalVolume(off);
468
469 if (mother) {
471 mother->GetLogicalVolume()->GetName());
472 }
473 else {
474 fNameBuffer = "";
475 }
476 return fNameBuffer.data();
477}
478
479//_____________________________________________________________________________
481{
483
485
486 // Get current touchable
487 const G4VTouchable* touchable = GetCurrentTouchable();
488
489 // Check touchable depth
490 //
491 G4int depth = touchable->GetHistoryDepth();
492
493 // Compose the path
494 //
495 fNameBuffer = "";
496 for (G4int i = 0; i < depth; i++) {
497 G4VPhysicalVolume* physVolume = touchable->GetHistory()->GetVolume(i);
498 fNameBuffer += "/";
499 fNameBuffer += geometryServices->UserVolumeName(physVolume->GetName());
500 fNameBuffer += "_";
501 TG4Globals::AppendNumberToString(fNameBuffer, physVolume->GetCopyNo());
502 }
503
504 // Add current volume to the path
505 G4VPhysicalVolume* curPhysVolume = GetCurrentPhysicalVolume();
506 fNameBuffer += "/";
507 fNameBuffer += geometryServices->UserVolumeName(curPhysVolume->GetName());
508 fNameBuffer += "_";
509 TG4Globals::AppendNumberToString(fNameBuffer, curPhysVolume->GetCopyNo());
510
511 return fNameBuffer.data();
512}
513
514//_____________________________________________________________________________
516 Double_t& x, Double_t& y, Double_t& z) const
517{
519
520 G4Navigator* theNavigator =
521 G4TransportationManager::GetTransportationManager()
522 ->GetNavigatorForTracking();
523
524 G4bool valid;
525 G4ThreeVector theLocalNormal = theNavigator->GetLocalExitNormal(&valid);
526 if (!valid) return false;
527
528 G4ThreeVector theGlobalNormal =
529 theNavigator->GetLocalToGlobalTransform().TransformAxis(theLocalNormal);
530
531 x = theGlobalNormal.x();
532 y = theGlobalNormal.y();
533 z = theGlobalNormal.z();
534
535 return true;
536}
537
538//_____________________________________________________________________________
540 Float_t& a, Float_t& z, Float_t& dens, Float_t& radl, Float_t& absl) const
541{
549
550 G4VPhysicalVolume* physVolume = GetCurrentPhysicalVolume();
551
552 G4Material* material = physVolume->GetLogicalVolume()->GetMaterial();
553
554 G4int nofElements = material->GetNumberOfElements();
556 a = geometryServices->GetEffA(material);
557 z = geometryServices->GetEffZ(material);
558
559 // density
560 dens = material->GetDensity();
561 dens /= TG4G3Units::MassDensity();
562
563 // radiation length
564 radl = material->GetRadlen();
565 radl /= TG4G3Units::Length();
566
567 absl = 0.; // this parameter is not defined in Geant4
568 return nofElements;
569}
570
571//_____________________________________________________________________________
573{
575
577 GetCurrentPhysicalVolume()->GetLogicalVolume());
578}
579
580//_____________________________________________________________________________
581void TG4StepManager::Gmtod(Float_t* xm, Float_t* xd, Int_t iflag)
582{
591
592 G4double* dxm = TG4GeometryServices::Instance()->CreateG4doubleArray(xm, 3);
593 G4double* dxd = TG4GeometryServices::Instance()->CreateG4doubleArray(xd, 3, false);
594
595 Gmtod(dxm, dxd, iflag);
596
597 // Fill computed xd coordinates
598 for (G4int i = 0; i < 3; i++) {
599 xd[i] = dxd[i];
600 }
601
602 delete[] dxm;
603 delete[] dxd;
604}
605
606//_____________________________________________________________________________
607void TG4StepManager::Gmtod(Double_t* xm, Double_t* xd, Int_t iflag)
608{
617
618#ifdef MCDEBUG
619 if (iflag != 1 && iflag != 2) {
620 TString text = "iflag=";
621 text += iflag;
623 "TG4StepManager", "Gmtod", text + " is different from 1..2.");
624 return;
625 }
626#endif
627
628 const G4AffineTransform& affineTransform =
629 GetCurrentTouchable()->GetHistory()->GetTopTransform();
630
631 G4ThreeVector theGlobalPoint(xm[0] * TG4G3Units::Length(),
632 xm[1] * TG4G3Units::Length(), xm[2] * TG4G3Units::Length());
633 G4ThreeVector theLocalPoint;
634 if (iflag == 1)
635 theLocalPoint = affineTransform.TransformPoint(theGlobalPoint);
636 else {
637 // if ( iflag == 2)
638 theLocalPoint = affineTransform.TransformAxis(theGlobalPoint);
639 }
640
641 xd[0] = theLocalPoint.x() * TG4G3Units::InverseLength();
642 xd[1] = theLocalPoint.y() * TG4G3Units::InverseLength();
643 xd[2] = theLocalPoint.z() * TG4G3Units::InverseLength();
644}
645
646//_____________________________________________________________________________
647void TG4StepManager::Gdtom(Float_t* xd, Float_t* xm, Int_t iflag)
648{
656
657 G4double* dxd = TG4GeometryServices::Instance()->CreateG4doubleArray(xd, 3);
658 G4double* dxm = TG4GeometryServices::Instance()->CreateG4doubleArray(xm, 3, false);
659
660 Gdtom(dxd, dxm, iflag);
661
662 // Fill computed xm coordinates
663 for (G4int i = 0; i < 3; i++) {
664 xm[i] = dxm[i];
665 }
666
667 delete[] dxd;
668 delete[] dxm;
669}
670
671//_____________________________________________________________________________
672void TG4StepManager::Gdtom(Double_t* xd, Double_t* xm, Int_t iflag)
673{
681
682#ifdef MCDEBUG
683 if (iflag != 1 && iflag != 2) {
684 TString text = "iflag=";
685 text += iflag;
687 "TG4StepManager", "Gmtod", text + " is different from 1..2.");
688 return;
689 }
690#endif
691
692 const G4AffineTransform& affineTransform =
693 GetCurrentTouchable()->GetHistory()->GetTopTransform().Inverse();
694
695 G4ThreeVector theLocalPoint(xd[0] * TG4G3Units::Length(),
696 xd[1] * TG4G3Units::Length(), xd[2] * TG4G3Units::Length());
697 G4ThreeVector theGlobalPoint;
698 if (iflag == 1)
699 theGlobalPoint = affineTransform.TransformPoint(theLocalPoint);
700 else {
701 // if( iflag == 2)
702 theGlobalPoint = affineTransform.TransformAxis(theLocalPoint);
703 }
704
705 xm[0] = theGlobalPoint.x() * TG4G3Units::InverseLength();
706 xm[1] = theGlobalPoint.y() * TG4G3Units::InverseLength();
707 xm[2] = theGlobalPoint.z() * TG4G3Units::InverseLength();
708}
709
710//_____________________________________________________________________________
712{
715
716 G4LogicalVolume* curLogVolume =
717 GetCurrentPhysicalVolume()->GetLogicalVolume();
718
719 // check this
720 G4UserLimits* userLimits = curLogVolume->GetUserLimits();
721
722 G4double maxStep;
723 if (userLimits == 0) {
724 TG4Globals::Warning("TG4StepManager", "MaxStep",
725 "User Limits are not defined for the current logical volume " +
726 TString(curLogVolume->GetName()) + ".");
727 return FLT_MAX;
728 }
729 else {
730 const G4Track& trackRef = *(fTrack);
731 maxStep = userLimits->GetMaxAllowedStep(trackRef);
732 maxStep /= TG4G3Units::Length();
733 return maxStep;
734 }
735}
736
737//_____________________________________________________________________________
739{
741
743}
744
745//_____________________________________________________________________________
746void TG4StepManager::TrackPosition(TLorentzVector& position) const
747{
751
752#ifdef MCDEBUG
753 CheckTrack();
754#endif
755
756 G4ThreeVector positionVector;
757 if (fStepStatus == kGflashSpot) {
758 positionVector = fGflashSpot->GetEnergySpot()->GetPosition();
759 }
760 else {
761 // get position
762 // check if this is == to PostStepPoint position !!
763 positionVector = fTrack->GetPosition();
764 }
765 positionVector *= 1. / (TG4G3Units::Length());
766
767 // global time
768 G4double time = fTrack->GetGlobalTime();
769 time /= TG4G3Units::Time();
770
771 SetTLorentzVector(positionVector, time, position);
772}
773
774//_____________________________________________________________________________
775void TG4StepManager::TrackPosition(Double_t& x, Double_t& y, Double_t& z) const
776{
779
780#ifdef MCDEBUG
781 CheckTrack();
782#endif
783
784 G4ThreeVector positionVector;
785 if (fStepStatus == kGflashSpot) {
786 positionVector = fGflashSpot->GetEnergySpot()->GetPosition();
787 }
788 else {
789 // get position
790 // check if this is == to PostStepPoint position !!
791 positionVector = fTrack->GetPosition();
792 }
793 positionVector *= 1. / (TG4G3Units::Length());
794
795 x = positionVector.x();
796 y = positionVector.y();
797 z = positionVector.z();
798}
799
800//_____________________________________________________________________________
801void TG4StepManager::TrackPosition(Float_t& x, Float_t& y, Float_t& z) const
802{
805
806 Double_t dx, dy, dz;
807 TrackPosition(dx, dy, dz);
808
809 x = static_cast<float>(dx);
810 y = static_cast<float>(dy);
811 z = static_cast<float>(dz);
812}
813
814//_____________________________________________________________________________
815void TG4StepManager::TrackMomentum(TLorentzVector& momentum) const
816{
819
820#ifdef MCDEBUG
821 CheckTrack();
822#endif
823
824 G4ThreeVector momentumVector = fTrack->GetMomentum();
825 momentumVector *= 1. / (TG4G3Units::Energy());
826
827 G4double energy = fTrack->GetDynamicParticle()->GetTotalEnergy();
828 energy /= TG4G3Units::Energy();
829
830 SetTLorentzVector(momentumVector, energy, momentum);
831}
832
833//_____________________________________________________________________________
835 Double_t& px, Double_t& py, Double_t& pz, Double_t& etot) const
836{
839
840#ifdef MCDEBUG
841 CheckTrack();
842#endif
843
844 G4ThreeVector momentumVector = fTrack->GetMomentum();
845 momentumVector *= 1. / (TG4G3Units::Energy());
846
847 px = momentumVector.x();
848 py = momentumVector.y();
849 pz = momentumVector.z();
850
851 etot = fTrack->GetDynamicParticle()->GetTotalEnergy();
852 etot /= TG4G3Units::Energy();
853}
854
855//_____________________________________________________________________________
857 Float_t& px, Float_t& py, Float_t& pz, Float_t& etot) const
858{
861
862 Double_t dpx, dpy, dpz, detot;
863 TrackMomentum(dpx, dpy, dpz, detot);
864
865 px = static_cast<float>(dpx);
866 py = static_cast<float>(dpy);
867 pz = static_cast<float>(dpz);
868 etot = static_cast<float>(detot);
869}
870
871//_____________________________________________________________________________
873{
876
877 if (fStepStatus == kNormalStep) {
878#ifdef MCDEBUG
879 CheckStep("TrackStep");
880#endif
881 return fStep->GetStepLength() * TG4G3Units::InverseLength();
882 }
883 else
884 return 0;
885}
886
887//_____________________________________________________________________________
889{
892
893#ifdef MCDEBUG
894 CheckTrack();
895#endif
897 return fTrack->GetTrackLength() * TG4G3Units::InverseLength();
898 }
899 return fTrack->GetTrackLength() * TG4G3Units::InverseLength() +
900 fInitialVMCTrackStatus->fTrackLength;
901}
902
903//_____________________________________________________________________________
905{
911
912#ifdef MCDEBUG
913 CheckTrack();
914#endif
915
916 return fTrack->GetGlobalTime() * TG4G3Units::InverseTime();
917}
918
919//_____________________________________________________________________________
920Double_t TG4StepManager::Edep() const
921{
923
924 if (fStepStatus == kNormalStep) {
925
926#ifdef MCDEBUG
927 CheckStep("Edep");
928#endif
929
930 return fStep->GetTotalEnergyDeposit() * TG4G3Units::InverseEnergy();
931 }
932
933 if (fStepStatus == kBoundary && fTrack->GetTrackStatus() == fStopAndKill) {
934 G4VProcess* proc = fSteppingManager->GetfCurrentProcess();
936 if (proc && physicsManager->GetMCProcess(proc) == kPLightScattering &&
937 physicsManager->GetOpBoundaryStatus() == kPLightDetection) {
938 return fTrack->GetTotalEnergy() * TG4G3Units::InverseEnergy();
939 }
940 }
941
942 if (fStepStatus == kGflashSpot) {
943
944#ifdef MCDEBUG
945 CheckGflashSpot("Edep");
946#endif
947
948 return fGflashSpot->GetEnergySpot()->GetEnergy() * TG4G3Units::InverseEnergy();
949 }
950
951 return 0;
952}
953
954//_____________________________________________________________________________
956{
958
959 if (fStepStatus == kNormalStep) {
960
961#ifdef MCDEBUG
962 CheckStep("NIELEdep");
963#endif
964
965 return fStep->GetNonIonizingEnergyDeposit() * TG4G3Units::InverseEnergy();
966 }
967
968 // return 0. in other cases (including kBoundary, kGflashSpot)
969 return 0;
970}
971
972//_____________________________________________________________________________
974{
976
978 return fTrack->GetCurrentStepNumber();
979 }
980 return fTrack->GetCurrentStepNumber() + fInitialVMCTrackStatus->fStepNumber;
981}
982
983//_____________________________________________________________________________
985{
987
988 return fTrack->GetWeight();
989}
990
991//_____________________________________________________________________________
993 Double_t& polX, Double_t& polY, Double_t& polZ) const
994{
996
997 const G4ThreeVector& pol = fTrack->GetPolarization();
998 polX = pol.x();
999 polY = pol.y();
1000 polZ = pol.z();
1001}
1002
1003//_____________________________________________________________________________
1004void TG4StepManager::TrackPolarization(TVector3& pol) const
1005{
1007
1008 const G4ThreeVector& polG4 = fTrack->GetPolarization();
1009 pol[0] = polG4.x();
1010 pol[1] = polG4.y();
1011 pol[2] = polG4.z();
1012}
1013
1014//_____________________________________________________________________________
1016{
1018
1019#ifdef MCDEBUG
1020 CheckTrack();
1021#endif
1022
1023 G4ParticleDefinition* particle =
1024 fTrack->GetDynamicParticle()->GetDefinition();
1025
1026 // Ask TG4ParticlesManager to get PDG encoding
1027 // (in order to get PDG from extended TDatabasePDG
1028 // in case the standard PDG code is not defined)
1029 G4int pdgEncoding = TG4ParticlesManager::Instance()->GetPDGEncoding(particle);
1030
1031 // Make difference between optical photon from Cerenkov and
1032 // feedback photon generated by user
1033 if (pdgEncoding == 50000050) {
1034 TG4TrackInformation* trackInformation =
1036 if (trackInformation && trackInformation->GetPDGEncoding())
1037 pdgEncoding = trackInformation->GetPDGEncoding();
1038 }
1039
1040 return pdgEncoding;
1041}
1042
1043//_____________________________________________________________________________
1045{
1047
1048#ifdef MCDEBUG
1049 CheckTrack();
1050#endif
1051
1052 return fTrack->GetDynamicParticle()->GetDefinition()->GetPDGCharge() /
1054}
1055
1056//_____________________________________________________________________________
1058{
1060
1061#ifdef MCDEBUG
1062 CheckTrack();
1063#endif
1064
1065 return fTrack->GetDynamicParticle()->GetDefinition()->GetPDGMass() /
1067}
1068
1069//_____________________________________________________________________________
1070Double_t TG4StepManager::Etot() const
1071{
1073
1074#ifdef MCDEBUG
1075 CheckTrack();
1076#endif
1077
1078 return fTrack->GetDynamicParticle()->GetTotalEnergy() * TG4G3Units::InverseEnergy();
1079}
1080
1081// TO DO: revise these with added kGflashSpot status
1082
1083//_____________________________________________________________________________
1085{
1088
1089 if (fStepStatus == kNormalStep && !(IsTrackExiting())) {
1090 // track is always inside during a normal step
1091 return true;
1092 }
1093
1094 return false;
1095}
1096
1097//_____________________________________________________________________________
1099{
1102
1103 if (fStepStatus != kNormalStep) {
1104 // track is entering during a vertex or boundary step
1105 return true;
1106 }
1107
1108 return false;
1109}
1110
1111//_____________________________________________________________________________
1113{
1115
1116 if (fStepStatus == kNormalStep) {
1117
1118#ifdef MCDEBUG
1119 CheckStep("IsTrackExiting");
1120#endif
1121
1122 if (fStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary)
1123 return true;
1124 }
1125
1126 return false;
1127}
1128
1129//_____________________________________________________________________________
1131{
1134
1135 if (fStepStatus == kVertex) return false;
1136
1137#ifdef MCDEBUG
1138 CheckStep("IsTrackOut");
1139#endif
1140
1141 if (fStep->GetPostStepPoint()->GetStepStatus() == fWorldBoundary)
1142 return true;
1143 else
1144 return false;
1145}
1146
1147//_____________________________________________________________________________
1149{
1164
1165#ifdef MCDEBUG
1166 CheckTrack();
1167#endif
1168
1169 // check
1170 G4TrackStatus status = fTrack->GetTrackStatus();
1171 if ((status == fStopAndKill) || (status == fKillTrackAndSecondaries) ||
1172 (status == fSuspend) || (status == fPostponeToNextEvent)) {
1173 return true;
1174 }
1175 else
1176 return false;
1177}
1178
1179//_____________________________________________________________________________
1181{
1185
1186#ifdef MCDEBUG
1187 CheckTrack();
1188#endif
1189
1190 // check
1191 G4TrackStatus status = fTrack->GetTrackStatus();
1192 if ((status == fStopAndKill) || (status == fKillTrackAndSecondaries) ||
1193 (status == fPostponeToNextEvent)) {
1194 return true;
1195 }
1196 else
1197 return false;
1198}
1199
1200//_____________________________________________________________________________
1202{
1204
1205#ifdef MCDEBUG
1206 CheckTrack();
1207#endif
1208
1209 G4TrackStatus status = fTrack->GetTrackStatus();
1210 if ((status == fAlive) || (status == fStopButAlive))
1211 return true;
1212 else
1213 return false;
1214}
1215
1216//_____________________________________________________________________________
1218{
1220
1221 if (fStepStatus == kVertex)
1222 return true;
1223 else
1224 return false;
1225}
1226
1227//_____________________________________________________________________________
1229{
1232
1233 if (fStepStatus == kVertex || fStepStatus == kGflashSpot) return 0;
1234
1235#ifdef MCDEBUG
1237#endif
1238
1239 G4int nofSecondaries = 0;
1240 nofSecondaries += fSteppingManager->GetfN2ndariesAtRestDoIt();
1241 nofSecondaries += fSteppingManager->GetfN2ndariesAlongStepDoIt();
1242 nofSecondaries += fSteppingManager->GetfN2ndariesPostStepDoIt();
1243
1244 return nofSecondaries;
1245}
1246
1247//_____________________________________________________________________________
1248void TG4StepManager::GetSecondary(Int_t index, Int_t& particleId,
1249 TLorentzVector& position, TLorentzVector& momentum)
1250{
1257
1258#ifdef MCDEBUG
1260#endif
1261
1262 G4int nofSecondaries = NSecondaries();
1263 if (!nofSecondaries) return;
1264
1265 const G4TrackVector* secondaryTracks = fSteppingManager->GetSecondary();
1266#ifdef MCDEBUG
1267 if (!secondaryTracks) {
1269 "TG4StepManager", "GetSecondary", "Secondary tracks vector is empty");
1270 }
1271
1272 if (index >= nofSecondaries) {
1274 "TG4StepManager", "GetSecondary", "Wrong secondary track index.");
1275 }
1276#endif
1277
1278 // the index of the first secondary of this step
1279 G4int startIndex = secondaryTracks->size() - nofSecondaries;
1280 // (the secondaryTracks vector contains secondaries
1281 // produced by the track at previous steps, too)
1282 G4Track* track = (*secondaryTracks)[startIndex + index];
1283
1284 // particle encoding
1285 particleId = track->GetDynamicParticle()->GetDefinition()->GetPDGEncoding();
1286
1287 // position & time
1288 G4ThreeVector positionVector = track->GetPosition();
1289 positionVector *= 1. / (TG4G3Units::Length());
1290 G4double time = track->GetGlobalTime();
1291 time /= TG4G3Units::Time();
1292 SetTLorentzVector(positionVector, time, position);
1293
1294 // momentum & energy
1295 G4ThreeVector momentumVector = track->GetMomentum();
1296 momentumVector *= 1. / (TG4G3Units::Energy());
1297 G4double energy = track->GetDynamicParticle()->GetTotalEnergy();
1298 energy /= TG4G3Units::Energy();
1299 SetTLorentzVector(momentumVector, energy, momentum);
1300}
1301
1302//_____________________________________________________________________________
1303TMCProcess TG4StepManager::ProdProcess(Int_t isec) const
1304{
1307
1308 G4int nofSecondaries = NSecondaries();
1309 if (fStepStatus == kVertex || !nofSecondaries) return kPNoProcess;
1310
1311#ifdef MCDEBUG
1312 CheckStep("ProdProcess");
1313#endif
1314
1316 // If this funcion is called from SD, it is earlier than TG4SteppingAction
1317 // fixes the creator processes
1318
1319 const G4TrackVector* secondaryTracks = fSteppingManager->GetSecondary();
1320
1321#ifdef MCDEBUG
1322 // should never happen
1323 if (!secondaryTracks) {
1325 "TG4StepManager", "ProdProcess", "Secondary tracks vector is empty.");
1326
1327 return kPNoProcess;
1328 }
1329
1330 if (isec >= nofSecondaries) {
1332 "TG4StepManager", "ProdProcess", "Wrong secondary track index.");
1333
1334 return kPNoProcess;
1335 }
1336#endif
1337
1338 // the index of the first secondary of this step
1339 G4int startIndex = secondaryTracks->size() - nofSecondaries;
1340 // the secondaryTracks vector contains secondaries
1341 // produced by the track at previous steps, too
1342
1343 // the secondary track with specified isec index
1344 G4Track* track = (*secondaryTracks)[startIndex + isec];
1345
1346 const G4VProcess* kpProcess = track->GetCreatorProcess();
1347
1348 TMCProcess mcProcess = TG4PhysicsManager::Instance()->GetMCProcess(kpProcess);
1349
1350 // distinguish kPDeltaRay from kPEnergyLoss
1351 if (mcProcess == kPEnergyLoss) mcProcess = kPDeltaRay;
1352
1353 return mcProcess;
1354}
1355
1356//_____________________________________________________________________________
1357Int_t TG4StepManager::StepProcesses(TArrayI& processes) const
1358{
1363
1366 G4int nofProcesses = 1;
1367 processes.Set(nofProcesses);
1368 processes[0] = kPNull;
1369 return nofProcesses;
1370 }
1371
1372#ifdef MCDEBUG
1374 CheckStep("StepProcesses");
1375#endif
1376
1378 // If this funcion is called from SD, it is earlier than TG4SteppingAction
1379 // fixes the creator processes
1380
1381 // along step processes
1382 G4ProcessVector* processVector = fStep->GetTrack()
1383 ->GetDefinition()
1384 ->GetProcessManager()
1385 ->GetAlongStepProcessVector();
1386 G4int nofAlongStep = processVector->entries();
1387
1388 // process defined step
1389 const G4VProcess* kpLastProcess =
1390 fStep->GetPostStepPoint()->GetProcessDefinedStep();
1391
1392 // set array size
1393 processes.Set(nofAlongStep + 2);
1394 // maximum number of processes:
1395 // nofAlongStep (along step) - 1 (transportations) + 1 (post step process)
1396 // + possibly 2 (additional processes if OpBoundary )
1397 // => nofAlongStep + 2
1398
1399 // fill array with (nofAlongStep-1) along step processes
1401 G4int counter = 0;
1402 for (G4int i = 0; i < nofAlongStep; i++) {
1403 G4VProcess* g4Process = (*processVector)[i];
1404 // do not fill transportation along step process
1405 if (g4Process && g4Process->GetProcessSubType() != TRANSPORTATION)
1406 processes[counter++] = physicsManager->GetMCProcess(g4Process);
1407 }
1408
1409 // fill array with optical photon information
1410 if (fStep->GetTrack()->GetDefinition() == G4OpticalPhoton::Definition() &&
1411 kpLastProcess->GetProcessSubType() == TRANSPORTATION &&
1412 physicsManager->IsOpBoundaryProcess()) {
1413
1414 // add light scattering anbd reflection/absorption as additional processes
1415 processes[counter++] = kPLightScattering;
1416 processes[counter++] = physicsManager->GetOpBoundaryStatus();
1417 }
1418
1419 // fill array with last process
1420 processes[counter++] = physicsManager->GetMCProcess(kpLastProcess);
1421
1422 return counter;
1423}
Definition of the TG4G3Units class.
Definition of the TG4GeometryServices class.
Definition of the TG4Globals class and basic container types.
Definition of the TG4Limits class.
Definition of the TG4ParticlesManager class.
Definition of the TG4PhysicsManager class.
Definition of the TG4SDServices class.
Definition of the TG4StepManager class.
@ kGflashSpot
in post step point with Gflash
@ kBoundary
when crossing geometrical boundary
@ kVertex
in track vertex
@ kNormalStep
in post step point
Definition of the TG4SteppingAction class.
Definition of the TG4TrackInformation class.
Definition of the TG4TrackManager class.
static G4double Mass()
Definition TG4G3Units.h:111
static G4double Charge()
Definition TG4G3Units.h:99
static G4double Energy()
Definition TG4G3Units.h:105
static G4double Length()
Definition TG4G3Units.h:81
static G4double InverseLength()
Definition TG4G3Units.h:135
static G4double InverseTime()
Definition TG4G3Units.h:145
static G4double InverseEnergy()
Definition TG4G3Units.h:157
static G4double Time()
Definition TG4G3Units.h:93
static G4double MassDensity()
Definition TG4G3Units.h:117
Services for accessing to Geant4 geometry.
G4double * CreateG4doubleArray(Float_t *array, G4int size, G4bool copyValues=true) const
const G4String & UserVolumeName(const G4String &name) const
static TG4GeometryServices * Instance()
TG4Limits * GetLimits(G4UserLimits *limits) const
G4double GetEffA(G4Material *material) const
G4double GetEffZ(G4Material *material) const
static void AppendNumberToString(G4String &string, G4int number)
static void Warning(const TString &className, const TString &methodName, const TString &text)
static void Exception(const TString &className, const TString &methodName, const TString &text)
Extended G4UserLimits class.
Definition TG4Limits.h:38
void SetMaxAllowedStepBack()
void SetCurrentMaxAllowedStep(G4double step)
G4int GetPDGEncoding(G4ParticleDefinition *particle)
static TG4ParticlesManager * Instance()
Geant4 implementation of the TVirtualMC interface methods for building Geant4 physics and access to i...
static TG4PhysicsManager * Instance()
G4bool IsOpBoundaryProcess() const
TMCProcess GetMCProcess(const G4VProcess *process)
TMCProcess GetOpBoundaryStatus()
G4int GetMediumID(G4LogicalVolume *volume) const
void SetIsStopRun(G4bool stopRun)
G4int GetVolumeID(const G4String &volumeName) const
static TG4SDServices * Instance()
Geant4 implementation of the TVirtualMC interface methods for access to Geant4 at step level.
const char * CurrentVolPath()
G4Track * fTrack
current track
Int_t CurrentMaterial(Float_t &a, Float_t &z, Float_t &dens, Float_t &radl, Float_t &absl) const
Double_t Edep() const
Int_t StepProcesses(TArrayI &proc) const
TMCProcess ProdProcess(Int_t isec) const
void CheckTrack() const
G4int fDivisionCopyNoOffset
division copy number offset
Int_t CurrentVolOffID(Int_t off, Int_t &copyNo) const
void TrackMomentum(TLorentzVector &momentum) const
const char * CurrentVolOffName(Int_t off) const
Bool_t IsTrackStop() const
Bool_t IsTrackDisappeared() const
Bool_t IsTrackEntering() const
Double_t TrackStep() const
TG4StepManager(const TString &userGeometry)
Double_t TrackTime() const
void SetMaxNStep(Int_t maxNofSteps)
const char * CurrentVolName() const
G4String fNameBuffer
buffer for current volume name or path
void GetSecondary(Int_t index, Int_t &particleId, TLorentzVector &position, TLorentzVector &momentum)
Double_t MaxStep() const
G4SteppingManager * fSteppingManager
G4SteppingManager.
TG4Limits * GetCurrentLimits() const
Double_t TrackMass() const
Int_t NSecondaries() const
Double_t NIELEdep() const
void CheckGflashSpot(const G4String &method) const
Bool_t CurrentBoundaryNormal(Double_t &x, Double_t &y, Double_t &z) const
Double_t TrackWeight() const
G4GFlashSpot * fGflashSpot
current Gflash spot
TG4TrackManager * fTrackManager
Cached pointer to thread-local track manager.
TG4StepStatus fStepStatus
step status
const G4VTouchable * GetCurrentTouchable() const
Int_t CurrentVolID(Int_t &copyNo) const
Bool_t IsTrackAlive() const
void CheckSteppingManager() const
Double_t TrackCharge() const
Int_t CurrentMedium() const
Double_t TrackLength() const
Bool_t IsTrackInside() const
static G4ThreadLocal TG4StepManager * fgInstance
this instance
Bool_t IsCollectTracks() const
void Gdtom(Double_t *xd, Double_t *xm, Int_t iflag)
G4int fCopyNoOffset
volume copy number offset
Bool_t IsNewTrack() const
Bool_t IsTrackOut() const
void CheckStep(const G4String &method) const
Int_t StepNumber() const
G4VPhysicalVolume * GetCurrentPhysicalVolume() const
TG4Limits * fLimitsModifiedOnFly
limits which step limit was modified during tracking
void SetTLorentzVector(G4ThreeVector xyz, G4double t, TLorentzVector &lv) const
Bool_t IsTrackExiting() const
Int_t TrackPid() const
void Gmtod(Double_t *xm, Double_t *xd, Int_t iflag)
TMCParticleStatus * fInitialVMCTrackStatus
The initial status of a VMC track when it was popped from the VMC stack.
void SetCollectTracks(Bool_t collectTracks)
void ForceDecayTime(Float_t pdg)
G4Step * fStep
current step
void TrackPolarization(Double_t &polX, Double_t &polY, Double_t &polZ) const
Int_t GetMaxNStep() const
G4VPhysicalVolume * GetCurrentOffPhysicalVolume(G4int off, G4bool warn=false) const
Double_t Etot() const
void TrackPosition(TLorentzVector &position) const
void SetInitialVMCTrackStatus(TMCParticleStatus *status)
void SetMaxStep(Double_t step)
void SetMaxNofSteps(G4int number)
void SetCollectTracks(G4bool collectTracks)
void ProcessTrackIfGeneralProcess(const G4Step *step)
G4bool GetCollectTracks() const
static TG4SteppingAction * Instance()
G4int GetMaxNofSteps() const
Defines additional track information.
G4int GetPDGEncoding() const
void SetInterrupt(G4bool interrupt)
void SetPDGLifetime(G4double pdgLifetime)
static TG4TrackManager * Instance()
TG4TrackInformation * GetTrackInformation(const G4Track *track) const