Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TrajectoryDrawerUtils.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// Jane Tinslay, John Allison, Joseph Perl November 2005
28//
30#include "G4Colour.hh"
31#include "G4Polyline.hh"
32#include "G4Polymarker.hh"
33#include "G4VTrajectory.hh"
34#include "G4VTrajectoryPoint.hh"
35#include "G4VisAttributes.hh"
36#include "G4VisTrajContext.hh"
37#include "G4VVisManager.hh"
38#include "G4UIcommand.hh"
39#include "G4AttValue.hh"
40
41#include <utility>
42
43#define G4warn G4cout
44
46
48
50 (const G4VTrajectory& traj,
51 const G4VisTrajContext& context,
52 G4Polyline& trajectoryLine,
53 G4Polymarker& auxiliaryPoints,
54 G4Polymarker& stepPoints,
55 std::vector<G4double>& trajectoryLineTimes,
56 std::vector<G4double>& auxiliaryPointTimes,
57 std::vector<G4double>& stepPointTimes)
58 {
59 TimesValidity validity = InvalidTimes;
60 if (context.GetTimeSliceInterval()) validity = ValidTimes;
61
62 // Memory for last trajectory point position for auxiliary point
63 // time interpolation algorithm. There are no auxiliary points
64 // for the first trajectory point, so its initial value is
65 // immaterial.
66 G4ThreeVector lastTrajectoryPointPosition;
67
68 // Keep positions. Don't store unless first or different.
69 std::vector<G4ThreeVector> positions;
70
71 for (G4int iPoint=0; iPoint<traj.GetPointEntries(); iPoint++) {
72
73 G4VTrajectoryPoint* aTrajectoryPoint = traj.GetPoint(iPoint);
74 const G4ThreeVector& trajectoryPointPosition =
75 aTrajectoryPoint->GetPosition();
76
77 // Only store if first or if different
78 if (positions.size() == 0 ||
79 trajectoryPointPosition != positions[positions.size()-1]) {
80
81 // Pre- and Post-Point times from the trajectory point...
82 G4double trajectoryPointPreTime = -std::numeric_limits<double>::max();
83 G4double trajectoryPointPostTime = std::numeric_limits<double>::max();
84
85 if (context.GetTimeSliceInterval() && validity == ValidTimes) {
86
87 std::vector<G4AttValue>* trajectoryPointAttValues =
88 aTrajectoryPoint->CreateAttValues();
89 if (!trajectoryPointAttValues) {
90 static G4bool warnedNoAttValues = false;
91 if (!warnedNoAttValues) {
92 G4warn <<
93 "*************************************************************************"
94 "\n* WARNING: G4TrajectoryDrawerUtils::GetPointsAndTimes: no att values."
95 "\n*************************************************************************"
96 << G4endl;
97 warnedNoAttValues = true;
98 }
99 validity = InvalidTimes;
100 } else {
101 G4bool foundPreTime = false, foundPostTime = false;
102 for (std::vector<G4AttValue>::iterator i =
103 trajectoryPointAttValues->begin();
104 i != trajectoryPointAttValues->end(); ++i) {
105 if (i->GetName() == "PreT") {
106 trajectoryPointPreTime =
108 foundPreTime = true;
109 }
110 if (i->GetName() == "PostT") {
111 trajectoryPointPostTime =
113 foundPostTime = true;
114 }
115 }
116 if (!foundPreTime || !foundPostTime) {
117 static G4bool warnedTimesNotFound = false;
118 if (!warnedTimesNotFound) {
119 G4warn <<
120 "*************************************************************************"
121 "\n* WARNING: G4TrajectoryDrawerUtils::GetPointsAndTimes: times not found."
122 "\n You need to specify \"/vis/scene/add/trajectories rich\""
123 "\n*************************************************************************"
124 << G4endl;
125 warnedTimesNotFound = true;
126 }
127 validity = InvalidTimes;
128 }
129 }
130 delete trajectoryPointAttValues; // (Must be deleted after use.)
131 }
132
133 const std::vector<G4ThreeVector>* auxiliaries
134 = aTrajectoryPoint->GetAuxiliaryPoints();
135 if (0 != auxiliaries) {
136 for (size_t iAux=0; iAux<auxiliaries->size(); ++iAux) {
137 const G4ThreeVector& auxPointPosition = (*auxiliaries)[iAux];
138 if (positions.size() == 0 ||
139 auxPointPosition != positions[positions.size()-1]) {
140 // Only store if first or if different
141 positions.push_back(trajectoryPointPosition);
142 trajectoryLine.push_back(auxPointPosition);
143 auxiliaryPoints.push_back(auxPointPosition);
144 if (validity == ValidTimes) {
145 // Interpolate time for auxiliary points...
146 G4double s1 =
147 (auxPointPosition - lastTrajectoryPointPosition).mag();
148 G4double s2 =
149 (trajectoryPointPosition - auxPointPosition).mag();
150 G4double t = trajectoryPointPreTime +
151 (trajectoryPointPostTime - trajectoryPointPreTime) *
152 (s1 / (s1 + s2));
153 trajectoryLineTimes.push_back(t);
154 auxiliaryPointTimes.push_back(t);
155 }
156 }
157 }
158 }
159
160 positions.push_back(trajectoryPointPosition);
161 trajectoryLine.push_back(trajectoryPointPosition);
162 stepPoints.push_back(trajectoryPointPosition);
163 if (validity == ValidTimes) {
164 trajectoryLineTimes.push_back(trajectoryPointPostTime);
165 stepPointTimes.push_back(trajectoryPointPostTime);
166 }
167 lastTrajectoryPointPosition = trajectoryPointPosition;
168 }
169 }
170 return validity;
171 }
172
173 static void SliceLine(G4double timeIncrement,
174 G4Polyline& trajectoryLine,
175 std::vector<G4double>& trajectoryLineTimes)
176 {
177 // Assumes valid arguments from GetPoints and GetTimes.
178
179 G4Polyline newTrajectoryLine;
180 std::vector<G4double> newTrajectoryLineTimes;
181
182 newTrajectoryLine.push_back(trajectoryLine[0]);
183 newTrajectoryLineTimes.push_back(trajectoryLineTimes[0]);
184 size_t lineSize = trajectoryLine.size();
185 if (lineSize > 1) {
186 for (size_t i = 1; i < trajectoryLine.size(); ++i) {
187 G4double deltaT = trajectoryLineTimes[i] - trajectoryLineTimes[i - 1];
188 if (deltaT > 0.) {
189 G4double practicalTimeIncrement =
190 std::max(timeIncrement, deltaT / 100.);
191 for (G4double t =
192 (int(trajectoryLineTimes[i - 1]/practicalTimeIncrement) + 1) *
193 practicalTimeIncrement;
194 t <= trajectoryLineTimes[i];
195 t += practicalTimeIncrement) {
196 G4ThreeVector pos = trajectoryLine[i - 1] +
197 (trajectoryLine[i] - trajectoryLine[i - 1]) *
198 ((t - trajectoryLineTimes[i - 1]) / deltaT);
199 newTrajectoryLine.push_back(pos);
200 newTrajectoryLineTimes.push_back(t);
201 }
202 }
203 newTrajectoryLine.push_back(trajectoryLine[i]);
204 newTrajectoryLineTimes.push_back(trajectoryLineTimes[i]);
205 }
206 }
207
208 trajectoryLine = std::move(newTrajectoryLine);
209 trajectoryLineTimes = std::move(newTrajectoryLineTimes);
210 }
211
212 static void DrawWithoutTime(const G4VisTrajContext& myContext,
213 G4Polyline& trajectoryLine,
214 G4Polymarker& auxiliaryPoints,
215 G4Polymarker& stepPoints)
216 {
217 // Draw without time slice information
218
219 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
220 if (0 == pVVisManager) return;
221
222 if (myContext.GetDrawLine() && myContext.GetLineVisible()) {
223 G4VisAttributes trajectoryLineAttribs(myContext.GetLineColour());
224 trajectoryLineAttribs.SetLineWidth(myContext.GetLineWidth());
225 trajectoryLine.SetVisAttributes(&trajectoryLineAttribs);
226
227 pVVisManager->Draw(trajectoryLine);
228 }
229
230 if (myContext.GetDrawAuxPts() && myContext.GetAuxPtsVisible()
231 && (auxiliaryPoints.size() > 0)) {
232 auxiliaryPoints.SetMarkerType(myContext.GetAuxPtsType());
233 auxiliaryPoints.SetSize(myContext.GetAuxPtsSizeType(), myContext.GetAuxPtsSize());
234 auxiliaryPoints.SetFillStyle(myContext.GetAuxPtsFillStyle());
235
236 G4VisAttributes auxiliaryPointsAttribs(myContext.GetAuxPtsColour());
237 auxiliaryPoints.SetVisAttributes(&auxiliaryPointsAttribs);
238
239 pVVisManager->Draw(auxiliaryPoints);
240 }
241
242 if (myContext.GetDrawStepPts() && myContext.GetStepPtsVisible()
243 && (stepPoints.size() > 0)) {
244 stepPoints.SetMarkerType(myContext.GetStepPtsType());
245 stepPoints.SetSize(myContext.GetStepPtsSizeType(), myContext.GetStepPtsSize());
246 stepPoints.SetFillStyle(myContext.GetStepPtsFillStyle());
247
248 G4VisAttributes stepPointsAttribs(myContext.GetStepPtsColour());
249 stepPoints.SetVisAttributes(&stepPointsAttribs);
250
251 pVVisManager->Draw(stepPoints);
252 }
253 }
254
255 static void DrawWithTime(const G4VisTrajContext& myContext,
256 G4Polyline& trajectoryLine,
257 G4Polymarker& auxiliaryPoints,
258 G4Polymarker& stepPoints,
259 std::vector<G4double>& trajectoryLineTimes,
260 std::vector<G4double>& auxiliaryPointTimes,
261 std::vector<G4double>& stepPointTimes)
262 {
263 // Draw with time slice information
264
265 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
266 if (0 == pVVisManager) return;
267
268 if (myContext.GetDrawLine() && myContext.GetLineVisible()) {
269 G4VisAttributes trajectoryLineAttribs(myContext.GetLineColour());
270 trajectoryLineAttribs.SetLineWidth(myContext.GetLineWidth());
271
272 for (size_t i = 1; i < trajectoryLine.size(); ++i ) {
273 G4Polyline slice;
274 slice.push_back(trajectoryLine[i -1]);
275 slice.push_back(trajectoryLine[i]);
276 trajectoryLineAttribs.SetStartTime(trajectoryLineTimes[i - 1]);
277 trajectoryLineAttribs.SetEndTime(trajectoryLineTimes[i]);
278 slice.SetVisAttributes(&trajectoryLineAttribs);
279 pVVisManager->Draw(slice);
280 }
281 }
282
283 if (myContext.GetDrawAuxPts() && myContext.GetAuxPtsVisible()
284 && (auxiliaryPoints.size() > 0)) {
285 G4VisAttributes auxiliaryPointsAttribs(myContext.GetAuxPtsColour());
286
287 for (size_t i = 0; i < auxiliaryPoints.size(); ++i ) {
288 G4Polymarker point;
289 point.push_back(auxiliaryPoints[i]);
290 point.SetMarkerType(myContext.GetAuxPtsType());
291 point.SetSize(myContext.GetAuxPtsSizeType(), myContext.GetAuxPtsSize());
292 point.SetFillStyle(myContext.GetAuxPtsFillStyle());
293 auxiliaryPointsAttribs.SetStartTime(auxiliaryPointTimes[i]);
294 auxiliaryPointsAttribs.SetEndTime(auxiliaryPointTimes[i]);
295 point.SetVisAttributes(&auxiliaryPointsAttribs);
296 pVVisManager->Draw(point);
297 }
298 }
299
300 if (myContext.GetDrawStepPts() && myContext.GetStepPtsVisible()
301 && (stepPoints.size() > 0)) {
302 G4VisAttributes stepPointsAttribs(myContext.GetStepPtsColour());
303
304 for (size_t i = 0; i < stepPoints.size(); ++i ) {
305 G4Polymarker point;
306 point.push_back(stepPoints[i]);
307 point.SetMarkerType(myContext.GetStepPtsType());
308 point.SetSize(myContext.GetStepPtsSizeType(), myContext.GetStepPtsSize());
309 point.SetFillStyle(myContext.GetStepPtsFillStyle());
310 stepPointsAttribs.SetStartTime(stepPointTimes[i]);
311 stepPointsAttribs.SetEndTime(stepPointTimes[i]);
312 point.SetVisAttributes(&stepPointsAttribs);
313 pVVisManager->Draw(point);
314 }
315 }
316 }
317
318 void DrawLineAndPoints(const G4VTrajectory& traj, const G4VisTrajContext& context)
319 {
320 // Return if don't need to do anything
321 if (!context.GetDrawLine() && !context.GetDrawAuxPts() && !context.GetDrawStepPts()) return;
322
323 // Get points and times (times are returned only if time-slicing
324 // is requested).
325 G4Polyline trajectoryLine;
326 G4Polymarker stepPoints;
327 G4Polymarker auxiliaryPoints;
328 std::vector<G4double> trajectoryLineTimes;
329 std::vector<G4double> stepPointTimes;
330 std::vector<G4double> auxiliaryPointTimes;
331
333 (traj, context,
334 trajectoryLine, auxiliaryPoints, stepPoints,
335 trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
336
337 if (validity == ValidTimes) {
338
339 SliceLine(context.GetTimeSliceInterval(),
340 trajectoryLine, trajectoryLineTimes);
341
342 DrawWithTime(context,
343 trajectoryLine, auxiliaryPoints, stepPoints,
344 trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
345
346 } else {
347
348 DrawWithoutTime(context, trajectoryLine, auxiliaryPoints, stepPoints);
349
350 }
351 }
352}
#define G4warn
Definition G4Scene.cc:41
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
void SetMarkerType(MarkerType)
static G4double ConvertToDimensionedDouble(const char *st)
void SetSize(SizeType, G4double)
Definition G4VMarker.cc:86
void SetFillStyle(FillStyle)
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual const std::vector< G4ThreeVector > * GetAuxiliaryPoints() const
virtual const G4ThreeVector GetPosition() const =0
virtual G4VTrajectoryPoint * GetPoint(G4int i) const =0
virtual G4int GetPointEntries() const =0
static G4VVisManager * GetConcreteInstance()
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
G4bool GetDrawAuxPts() const
G4Colour GetStepPtsColour() const
G4double GetLineWidth() const
G4double GetStepPtsSize() const
G4bool GetLineVisible() const
G4double GetTimeSliceInterval() const
G4bool GetDrawLine() const
G4VMarker::SizeType GetStepPtsSizeType() const
G4Polymarker::MarkerType GetAuxPtsType() const
G4double GetAuxPtsSize() const
G4VMarker::SizeType GetAuxPtsSizeType() const
G4Colour GetAuxPtsColour() const
G4VMarker::FillStyle GetStepPtsFillStyle() const
G4VMarker::FillStyle GetAuxPtsFillStyle() const
G4bool GetAuxPtsVisible() const
G4bool GetStepPtsVisible() const
G4Colour GetLineColour() const
G4Polymarker::MarkerType GetStepPtsType() const
G4bool GetDrawStepPts() const
void SetVisAttributes(const G4VisAttributes *)
Definition G4Visible.cc:98
TimesValidity GetPointsAndTimes(const G4VTrajectory &traj, const G4VisTrajContext &context, G4Polyline &trajectoryLine, G4Polymarker &auxiliaryPoints, G4Polymarker &stepPoints, std::vector< G4double > &trajectoryLineTimes, std::vector< G4double > &auxiliaryPointTimes, std::vector< G4double > &stepPointTimes)
void DrawLineAndPoints(const G4VTrajectory &traj, const G4VisTrajContext &)