Geant4 10.7.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
42
44
46 (const G4VTrajectory& traj,
47 const G4VisTrajContext& context,
48 G4Polyline& trajectoryLine,
49 G4Polymarker& auxiliaryPoints,
50 G4Polymarker& stepPoints,
51 std::vector<G4double>& trajectoryLineTimes,
52 std::vector<G4double>& auxiliaryPointTimes,
53 std::vector<G4double>& stepPointTimes)
54 {
55 TimesValidity validity = InvalidTimes;
56 if (context.GetTimeSliceInterval()) validity = ValidTimes;
57
58 // Memory for last trajectory point position for auxiliary point
59 // time interpolation algorithm. There are no auxiliary points
60 // for the first trajectory point, so its initial value is
61 // immaterial.
62 G4ThreeVector lastTrajectoryPointPosition;
63
64 // Keep positions. Don't store unless first or different.
65 std::vector<G4ThreeVector> positions;
66
67 for (G4int iPoint=0; iPoint<traj.GetPointEntries(); iPoint++) {
68
69 G4VTrajectoryPoint* aTrajectoryPoint = traj.GetPoint(iPoint);
70 const G4ThreeVector& trajectoryPointPosition =
71 aTrajectoryPoint->GetPosition();
72
73 // Only store if first or if different
74 if (positions.size() == 0 ||
75 trajectoryPointPosition != positions[positions.size()-1]) {
76
77 // Pre- and Post-Point times from the trajectory point...
78 G4double trajectoryPointPreTime = -std::numeric_limits<double>::max();
79 G4double trajectoryPointPostTime = std::numeric_limits<double>::max();
80
81 if (context.GetTimeSliceInterval() && validity == ValidTimes) {
82
83 std::vector<G4AttValue>* trajectoryPointAttValues =
84 aTrajectoryPoint->CreateAttValues();
85 if (!trajectoryPointAttValues) {
86 static G4bool warnedNoAttValues = false;
87 if (!warnedNoAttValues) {
88 G4cout <<
89 "*************************************************************************"
90 "\n* WARNING: G4TrajectoryDrawerUtils::GetPointsAndTimes: no att values."
91 "\n*************************************************************************"
92 << G4endl;
93 warnedNoAttValues = true;
94 }
95 validity = InvalidTimes;
96 } else {
97 G4bool foundPreTime = false, foundPostTime = false;
98 for (std::vector<G4AttValue>::iterator i =
99 trajectoryPointAttValues->begin();
100 i != trajectoryPointAttValues->end(); ++i) {
101 if (i->GetName() == "PreT") {
102 trajectoryPointPreTime =
104 foundPreTime = true;
105 }
106 if (i->GetName() == "PostT") {
107 trajectoryPointPostTime =
109 foundPostTime = true;
110 }
111 }
112 if (!foundPreTime || !foundPostTime) {
113 static G4bool warnedTimesNotFound = false;
114 if (!warnedTimesNotFound) {
115 G4cout <<
116 "*************************************************************************"
117 "\n* WARNING: G4TrajectoryDrawerUtils::GetPointsAndTimes: times not found."
118 "\n You need to specify \"/vis/scene/add/trajectories rich\""
119 "\n*************************************************************************"
120 << G4endl;
121 warnedTimesNotFound = true;
122 }
123 validity = InvalidTimes;
124 }
125 }
126 delete trajectoryPointAttValues; // (Must be deleted after use.)
127 }
128
129 const std::vector<G4ThreeVector>* auxiliaries
130 = aTrajectoryPoint->GetAuxiliaryPoints();
131 if (0 != auxiliaries) {
132 for (size_t iAux=0; iAux<auxiliaries->size(); ++iAux) {
133 const G4ThreeVector& auxPointPosition = (*auxiliaries)[iAux];
134 if (positions.size() == 0 ||
135 auxPointPosition != positions[positions.size()-1]) {
136 // Only store if first or if different
137 positions.push_back(trajectoryPointPosition);
138 trajectoryLine.push_back(auxPointPosition);
139 auxiliaryPoints.push_back(auxPointPosition);
140 if (validity == ValidTimes) {
141 // Interpolate time for auxiliary points...
142 G4double s1 =
143 (auxPointPosition - lastTrajectoryPointPosition).mag();
144 G4double s2 =
145 (trajectoryPointPosition - auxPointPosition).mag();
146 G4double t = trajectoryPointPreTime +
147 (trajectoryPointPostTime - trajectoryPointPreTime) *
148 (s1 / (s1 + s2));
149 trajectoryLineTimes.push_back(t);
150 auxiliaryPointTimes.push_back(t);
151 }
152 }
153 }
154 }
155
156 positions.push_back(trajectoryPointPosition);
157 trajectoryLine.push_back(trajectoryPointPosition);
158 stepPoints.push_back(trajectoryPointPosition);
159 if (validity == ValidTimes) {
160 trajectoryLineTimes.push_back(trajectoryPointPostTime);
161 stepPointTimes.push_back(trajectoryPointPostTime);
162 }
163 lastTrajectoryPointPosition = trajectoryPointPosition;
164 }
165 }
166 return validity;
167 }
168
169 static void SliceLine(G4double timeIncrement,
170 G4Polyline& trajectoryLine,
171 std::vector<G4double>& trajectoryLineTimes)
172 {
173 // Assumes valid arguments from GetPoints and GetTimes.
174
175 G4Polyline newTrajectoryLine;
176 std::vector<G4double> newTrajectoryLineTimes;
177
178 newTrajectoryLine.push_back(trajectoryLine[0]);
179 newTrajectoryLineTimes.push_back(trajectoryLineTimes[0]);
180 size_t lineSize = trajectoryLine.size();
181 if (lineSize > 1) {
182 for (size_t i = 1; i < trajectoryLine.size(); ++i) {
183 G4double deltaT = trajectoryLineTimes[i] - trajectoryLineTimes[i - 1];
184 if (deltaT > 0.) {
185 G4double practicalTimeIncrement =
186 std::max(timeIncrement, deltaT / 100.);
187 for (G4double t =
188 (int(trajectoryLineTimes[i - 1]/practicalTimeIncrement) + 1) *
189 practicalTimeIncrement;
190 t <= trajectoryLineTimes[i];
191 t += practicalTimeIncrement) {
192 G4ThreeVector pos = trajectoryLine[i - 1] +
193 (trajectoryLine[i] - trajectoryLine[i - 1]) *
194 ((t - trajectoryLineTimes[i - 1]) / deltaT);
195 newTrajectoryLine.push_back(pos);
196 newTrajectoryLineTimes.push_back(t);
197 }
198 }
199 newTrajectoryLine.push_back(trajectoryLine[i]);
200 newTrajectoryLineTimes.push_back(trajectoryLineTimes[i]);
201 }
202 }
203
204 trajectoryLine = newTrajectoryLine;
205 trajectoryLineTimes = newTrajectoryLineTimes;
206 }
207
208 static void DrawWithoutTime(const G4VisTrajContext& myContext,
209 G4Polyline& trajectoryLine,
210 G4Polymarker& auxiliaryPoints,
211 G4Polymarker& stepPoints)
212 {
213 // Draw without time slice information
214
216 if (0 == pVVisManager) return;
217
218 if (myContext.GetDrawLine() && myContext.GetLineVisible()) {
219 G4VisAttributes trajectoryLineAttribs(myContext.GetLineColour());
220 trajectoryLineAttribs.SetLineWidth(myContext.GetLineWidth());
221 trajectoryLine.SetVisAttributes(&trajectoryLineAttribs);
222
223 pVVisManager->Draw(trajectoryLine);
224 }
225
226 if (myContext.GetDrawAuxPts() && myContext.GetAuxPtsVisible()
227 && (auxiliaryPoints.size() > 0)) {
228 auxiliaryPoints.SetMarkerType(myContext.GetAuxPtsType());
229 auxiliaryPoints.SetSize(myContext.GetAuxPtsSizeType(), myContext.GetAuxPtsSize());
230 auxiliaryPoints.SetFillStyle(myContext.GetAuxPtsFillStyle());
231
232 G4VisAttributes auxiliaryPointsAttribs(myContext.GetAuxPtsColour());
233 auxiliaryPoints.SetVisAttributes(&auxiliaryPointsAttribs);
234
235 pVVisManager->Draw(auxiliaryPoints);
236 }
237
238 if (myContext.GetDrawStepPts() && myContext.GetStepPtsVisible()
239 && (stepPoints.size() > 0)) {
240 stepPoints.SetMarkerType(myContext.GetStepPtsType());
241 stepPoints.SetSize(myContext.GetStepPtsSizeType(), myContext.GetStepPtsSize());
242 stepPoints.SetFillStyle(myContext.GetStepPtsFillStyle());
243
244 G4VisAttributes stepPointsAttribs(myContext.GetStepPtsColour());
245 stepPoints.SetVisAttributes(&stepPointsAttribs);
246
247 pVVisManager->Draw(stepPoints);
248 }
249 }
250
251 static void DrawWithTime(const G4VisTrajContext& myContext,
252 G4Polyline& trajectoryLine,
253 G4Polymarker& auxiliaryPoints,
254 G4Polymarker& stepPoints,
255 std::vector<G4double>& trajectoryLineTimes,
256 std::vector<G4double>& auxiliaryPointTimes,
257 std::vector<G4double>& stepPointTimes)
258 {
259 // Draw with time slice information
260
262 if (0 == pVVisManager) return;
263
264 if (myContext.GetDrawLine() && myContext.GetLineVisible()) {
265 G4VisAttributes trajectoryLineAttribs(myContext.GetLineColour());
266 trajectoryLineAttribs.SetLineWidth(myContext.GetLineWidth());
267
268 for (size_t i = 1; i < trajectoryLine.size(); ++i ) {
269 G4Polyline slice;
270 slice.push_back(trajectoryLine[i -1]);
271 slice.push_back(trajectoryLine[i]);
272 trajectoryLineAttribs.SetStartTime(trajectoryLineTimes[i - 1]);
273 trajectoryLineAttribs.SetEndTime(trajectoryLineTimes[i]);
274 slice.SetVisAttributes(&trajectoryLineAttribs);
275 pVVisManager->Draw(slice);
276 }
277 }
278
279 if (myContext.GetDrawAuxPts() && myContext.GetAuxPtsVisible()
280 && (auxiliaryPoints.size() > 0)) {
281 G4VisAttributes auxiliaryPointsAttribs(myContext.GetAuxPtsColour());
282
283 for (size_t i = 0; i < auxiliaryPoints.size(); ++i ) {
284 G4Polymarker point;
285 point.push_back(auxiliaryPoints[i]);
286 point.SetMarkerType(myContext.GetAuxPtsType());
287 point.SetSize(myContext.GetAuxPtsSizeType(), myContext.GetAuxPtsSize());
288 point.SetFillStyle(myContext.GetAuxPtsFillStyle());
289 auxiliaryPointsAttribs.SetStartTime(auxiliaryPointTimes[i]);
290 auxiliaryPointsAttribs.SetEndTime(auxiliaryPointTimes[i]);
291 point.SetVisAttributes(&auxiliaryPointsAttribs);
292 pVVisManager->Draw(point);
293 }
294 }
295
296 if (myContext.GetDrawStepPts() && myContext.GetStepPtsVisible()
297 && (stepPoints.size() > 0)) {
298 G4VisAttributes stepPointsAttribs(myContext.GetStepPtsColour());
299
300 for (size_t i = 0; i < stepPoints.size(); ++i ) {
301 G4Polymarker point;
302 point.push_back(stepPoints[i]);
303 point.SetMarkerType(myContext.GetStepPtsType());
304 point.SetSize(myContext.GetStepPtsSizeType(), myContext.GetStepPtsSize());
305 point.SetFillStyle(myContext.GetStepPtsFillStyle());
306 stepPointsAttribs.SetStartTime(stepPointTimes[i]);
307 stepPointsAttribs.SetEndTime(stepPointTimes[i]);
308 point.SetVisAttributes(&stepPointsAttribs);
309 pVVisManager->Draw(point);
310 }
311 }
312 }
313
314 void DrawLineAndPoints(const G4VTrajectory& traj, const G4VisTrajContext& context)
315 {
316 // Return if don't need to do anything
317 if (!context.GetDrawLine() && !context.GetDrawAuxPts() && !context.GetDrawStepPts()) return;
318
319 // Get points and times (times are returned only if time-slicing
320 // is requested).
321 G4Polyline trajectoryLine;
322 G4Polymarker stepPoints;
323 G4Polymarker auxiliaryPoints;
324 std::vector<G4double> trajectoryLineTimes;
325 std::vector<G4double> stepPointTimes;
326 std::vector<G4double> auxiliaryPointTimes;
327
329 (traj, context,
330 trajectoryLine, auxiliaryPoints, stepPoints,
331 trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
332
333 if (validity == ValidTimes) {
334
335 SliceLine(context.GetTimeSliceInterval(),
336 trajectoryLine, trajectoryLineTimes);
337
338 DrawWithTime(context,
339 trajectoryLine, auxiliaryPoints, stepPoints,
340 trajectoryLineTimes, auxiliaryPointTimes, stepPointTimes);
341
342 } else {
343
344 DrawWithoutTime(context, trajectoryLine, auxiliaryPoints, stepPoints);
345
346 }
347 }
348}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetMarkerType(MarkerType)
static G4double ConvertToDimensionedDouble(const char *st)
Definition: G4UIcommand.cc:570
void SetSize(SizeType, G4double)
Definition: G4VMarker.cc:118
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:79
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 &)