File indexing completed on 2026-04-03 07:52:59
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029 #include "LXeTrajectory.hh"
0030
0031 #include "G4Circle.hh"
0032 #include "G4Colour.hh"
0033 #include "G4ParticleTable.hh"
0034 #include "G4ParticleTypes.hh"
0035 #include "G4Polyline.hh"
0036 #include "G4Polymarker.hh"
0037 #include "G4ThreeVector.hh"
0038 #include "G4Trajectory.hh"
0039 #include "G4TrajectoryPoint.hh"
0040 #include "G4VVisManager.hh"
0041 #include "G4VisAttributes.hh"
0042
0043 G4ThreadLocal G4Allocator<LXeTrajectory>* LXeTrajectoryAllocator = nullptr;
0044
0045
0046
0047 LXeTrajectory::LXeTrajectory(const G4Track* aTrack) : G4Trajectory(aTrack)
0048 {
0049 fParticleDefinition = aTrack->GetDefinition();
0050 }
0051
0052
0053
0054 LXeTrajectory::LXeTrajectory(LXeTrajectory& right)
0055 : G4Trajectory(right), fWls(right.fWls), fDrawit(right.fDrawit)
0056 {
0057 fParticleDefinition = right.fParticleDefinition;
0058 }
0059
0060
0061
0062 void LXeTrajectory::DrawTrajectory() const
0063 {
0064
0065
0066
0067 if (!fForceDraw && (!fDrawit || fForceNoDraw)) return;
0068
0069 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
0070 if (!pVVisManager) return;
0071
0072 const G4double markerSize = 0.05;
0073 G4bool lineRequired = true;
0074 G4bool markersRequired = true;
0075
0076 G4Polyline trajectoryLine;
0077 G4Polymarker stepPoints;
0078 G4Polymarker auxiliaryPoints;
0079
0080 for (G4int i = 0; i < GetPointEntries(); ++i) {
0081 G4VTrajectoryPoint* aTrajectoryPoint = GetPoint(i);
0082 const std::vector<G4ThreeVector>* auxiliaries = aTrajectoryPoint->GetAuxiliaryPoints();
0083 if (auxiliaries) {
0084 for (size_t iAux = 0; iAux < auxiliaries->size(); ++iAux) {
0085 const G4ThreeVector pos((*auxiliaries)[iAux]);
0086 if (lineRequired) {
0087 trajectoryLine.push_back(pos);
0088 }
0089 if (markersRequired) {
0090 auxiliaryPoints.push_back(pos);
0091 }
0092 }
0093 }
0094 const G4ThreeVector pos(aTrajectoryPoint->GetPosition());
0095 if (lineRequired) {
0096 trajectoryLine.push_back(pos);
0097 }
0098 if (markersRequired) {
0099 stepPoints.push_back(pos);
0100 }
0101 }
0102
0103 if (lineRequired) {
0104 G4Colour colour;
0105
0106 if (fParticleDefinition == G4OpticalPhoton::OpticalPhotonDefinition()) {
0107 if (fWls)
0108 colour = G4Colour(1., 0., 0.);
0109 else {
0110 colour = G4Colour(0., 1., 0.);
0111 }
0112 }
0113 else
0114 colour = G4Colour(0., 0., 1.);
0115
0116 G4VisAttributes trajectoryLineAttribs(colour);
0117 trajectoryLine.SetVisAttributes(&trajectoryLineAttribs);
0118 pVVisManager->Draw(trajectoryLine);
0119 }
0120 if (markersRequired) {
0121 auxiliaryPoints.SetMarkerType(G4Polymarker::squares);
0122 auxiliaryPoints.SetScreenSize(markerSize);
0123 auxiliaryPoints.SetFillStyle(G4VMarker::filled);
0124 G4VisAttributes auxiliaryPointsAttribs(G4Colour(0., 1., 1.));
0125 auxiliaryPoints.SetVisAttributes(&auxiliaryPointsAttribs);
0126 pVVisManager->Draw(auxiliaryPoints);
0127
0128 stepPoints.SetMarkerType(G4Polymarker::circles);
0129 stepPoints.SetScreenSize(markerSize);
0130 stepPoints.SetFillStyle(G4VMarker::filled);
0131 G4VisAttributes stepPointsAttribs(G4Colour(1., 1., 0.));
0132 stepPoints.SetVisAttributes(&stepPointsAttribs);
0133 pVVisManager->Draw(stepPoints);
0134 }
0135 }