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