Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ViewParameters Class Reference

#include <G4ViewParameters.hh>

Public Types

enum  DrawingStyle { wireframe, hlr, hsr, hlhsr }
 
enum  CutawayMode { cutawayUnion, cutawayIntersection }
 
enum  RotationStyle { constrainUpDirection, freeRotation }
 

Public Member Functions

 G4ViewParameters ()
 
 ~G4ViewParameters ()
 
G4bool operator!= (const G4ViewParameters &) const
 
DrawingStyle GetDrawingStyle () const
 
G4bool IsAuxEdgeVisible () const
 
G4bool IsCulling () const
 
G4bool IsCullingInvisible () const
 
G4bool IsDensityCulling () const
 
G4double GetVisibleDensity () const
 
G4bool IsCullingCovered () const
 
G4bool IsSection () const
 
const G4Plane3DGetSectionPlane () const
 
G4bool IsCutaway () const
 
CutawayMode GetCutawayMode () const
 
const G4PlanesGetCutawayPlanes () const
 
G4bool IsExplode () const
 
G4double GetExplodeFactor () const
 
const G4Point3DGetExplodeCentre () const
 
G4int GetNoOfSides () const
 
const G4Vector3DGetViewpointDirection () const
 
const G4Vector3DGetUpVector () const
 
G4double GetFieldHalfAngle () const
 
G4double GetZoomFactor () const
 
const G4Vector3DGetScaleFactor () const
 
const G4Point3DGetCurrentTargetPoint () const
 
G4double GetDolly () const
 
G4bool GetLightsMoveWithCamera () const
 
const G4Vector3DGetLightpointDirection () const
 
G4Vector3DGetActualLightpointDirection ()
 
const G4VisAttributesGetDefaultVisAttributes () const
 
const G4VisAttributesGetDefaultTextVisAttributes () const
 
const G4VMarkerGetDefaultMarker () const
 
G4double GetGlobalMarkerScale () const
 
G4double GetGlobalLineWidthScale () const
 
G4bool IsMarkerNotHidden () const
 
unsigned int GetWindowSizeHintX () const
 
unsigned int GetWindowSizeHintY () const
 
G4int GetWindowAbsoluteLocationHintX (G4int) const
 
G4int GetWindowAbsoluteLocationHintY (G4int) const
 
G4int GetWindowLocationHintX () const
 
G4int GetWindowLocationHintY () const
 
const G4StringGetXGeometryString () const
 
bool IsWindowSizeHintX () const
 
bool IsWindowSizeHintY () const
 
bool IsWindowLocationHintX () const
 
bool IsWindowLocationHintY () const
 
G4bool IsAutoRefresh () const
 
const G4ColourGetBackgroundColour () const
 
G4bool IsPicking () const
 
RotationStyle GetRotationStyle () const
 
const std::vector
< G4ModelingParameters::VisAttributesModifier > & 
GetVisAttributesModifiers () const
 
G4double GetCameraDistance (G4double radius) const
 
G4double GetNearDistance (G4double cameraDistance, G4double radius) const
 
G4double GetFarDistance (G4double cameraDistance, G4double nearDistance, G4double radius) const
 
G4double GetFrontHalfHeight (G4double nearDistance, G4double radius) const
 
void SetDrawingStyle (G4ViewParameters::DrawingStyle style)
 
void SetAuxEdgeVisible (G4bool)
 
void SetCulling (G4bool)
 
void SetCullingInvisible (G4bool)
 
void SetDensityCulling (G4bool)
 
void SetVisibleDensity (G4double visibleDensity)
 
void SetCullingCovered (G4bool)
 
void SetSectionPlane (const G4Plane3D &sectionPlane)
 
void UnsetSectionPlane ()
 
void SetCutawayMode (CutawayMode)
 
void AddCutawayPlane (const G4Plane3D &cutawayPlane)
 
void ChangeCutawayPlane (size_t index, const G4Plane3D &cutawayPlane)
 
void ClearCutawayPlanes ()
 
void SetExplodeFactor (G4double explodeFactor)
 
void UnsetExplodeFactor ()
 
void SetExplodeCentre (const G4Point3D &explodeCentre)
 
G4int SetNoOfSides (G4int nSides)
 
void SetViewpointDirection (const G4Vector3D &viewpointDirection)
 
void SetViewAndLights (const G4Vector3D &viewpointDirection)
 
void SetUpVector (const G4Vector3D &upVector)
 
void SetFieldHalfAngle (G4double fieldHalfAngle)
 
void SetOrthogonalProjection ()
 
void SetPerspectiveProjection (G4double fieldHalfAngle=30.*CLHEP::deg)
 
void SetZoomFactor (G4double zoomFactor)
 
void MultiplyZoomFactor (G4double zoomFactorMultiplier)
 
void SetScaleFactor (const G4Vector3D &scaleFactor)
 
void MultiplyScaleFactor (const G4Vector3D &scaleFactorMultiplier)
 
void SetCurrentTargetPoint (const G4Point3D &currentTargetPoint)
 
void SetDolly (G4double dolly)
 
void IncrementDolly (G4double dollyIncrement)
 
void SetLightpointDirection (const G4Vector3D &lightpointDirection)
 
void SetLightsMoveWithCamera (G4bool moves)
 
void SetPan (G4double right, G4double up)
 
void IncrementPan (G4double right, G4double up)
 
void IncrementPan (G4double right, G4double up, G4double forward)
 
void SetDefaultVisAttributes (const G4VisAttributes &)
 
void SetDefaultColour (const G4Colour &)
 
void SetDefaultTextVisAttributes (const G4VisAttributes &)
 
void SetDefaultTextColour (const G4Colour &)
 
void SetDefaultMarker (const G4VMarker &defaultMarker)
 
void SetGlobalMarkerScale (G4double globalMarkerScale)
 
void SetGlobalLineWidthScale (G4double globalLineWidthScale)
 
void SetMarkerHidden ()
 
void SetMarkerNotHidden ()
 
void SetWindowSizeHint (G4int xHint, G4int yHint)
 
void SetWindowLocationHint (G4int xHint, G4int yHint)
 
void SetXGeometryString (const G4String &)
 
void SetAutoRefresh (G4bool)
 
void SetBackgroundColour (const G4Colour &)
 
void SetPicking (G4bool)
 
void SetRotationStyle (RotationStyle)
 
void ClearVisAttributesModifiers ()
 
void AddVisAttributesModifier (const G4ModelingParameters::VisAttributesModifier &)
 
G4String CameraAndLightingCommands (const G4Point3D standardTargetPoint) const
 
G4String DrawingStyleCommands () const
 
G4String SceneModifyingCommands () const
 
G4String TouchableCommands () const
 
void PrintDifferences (const G4ViewParameters &v) const
 

Static Public Member Functions

static G4ViewParametersCatmullRomCubicSplineInterpolation (const std::vector< G4ViewParameters > &views, G4int nInterpolationPoints=50)
 

Friends

std::ostream & operator<< (std::ostream &, const DrawingStyle &)
 
std::ostream & operator<< (std::ostream &, const G4ViewParameters &)
 

Detailed Description

Definition at line 90 of file G4ViewParameters.hh.

Member Enumeration Documentation

Enumerator
cutawayUnion 
cutawayIntersection 

Definition at line 101 of file G4ViewParameters.hh.

101  {
102  cutawayUnion, // Union (addition) of result of each cutaway plane.
103  cutawayIntersection // Intersection (multiplication) " .
104  };
Enumerator
wireframe 
hlr 
hsr 
hlhsr 

Definition at line 94 of file G4ViewParameters.hh.

94  {
95  wireframe, // Draw edges - no hidden line removal.
96  hlr, // Draw edges - hidden lines removed.
97  hsr, // Draw surfaces - hidden surfaces removed.
98  hlhsr // Draw surfaces and edges - hidden removed.
99  };
Enumerator
constrainUpDirection 
freeRotation 

Definition at line 106 of file G4ViewParameters.hh.

106  {
107  constrainUpDirection, // Standard, HEP convention.
108  freeRotation // Free, Google-like rotation, using mouse-grab.
109  };

Constructor & Destructor Documentation

G4ViewParameters::G4ViewParameters ( )

Definition at line 43 of file G4ViewParameters.cc.

43  :
44  fDrawingStyle (wireframe),
45  fAuxEdgeVisible (false),
46  fCulling (true),
47  fCullInvisible (true),
48  fDensityCulling (false),
49  fVisibleDensity (0.01 * g / cm3),
50  fCullCovered (false),
51  fSection (false),
52  fSectionPlane (),
53  fCutawayMode (cutawayUnion),
54  fCutawayPlanes (),
55  fExplodeFactor (1.),
56  fNoOfSides (72),
57  fViewpointDirection (G4Vector3D (0., 0., 1.)), // On z-axis.
58  fUpVector (G4Vector3D (0., 1., 0.)), // y-axis up.
59  fFieldHalfAngle (0.), // Orthogonal projection.
60  fZoomFactor (1.),
61  fScaleFactor (G4Vector3D (1., 1., 1.)),
62  fCurrentTargetPoint (),
63  fDolly (0.),
64  fLightsMoveWithCamera (false),
65  fRelativeLightpointDirection (G4Vector3D (1., 1., 1.)),
66  fActualLightpointDirection (G4Vector3D (1., 1., 1.)),
67  fDefaultVisAttributes (),
68  fDefaultTextVisAttributes (G4Colour (0., 0., 1.)),
69  fDefaultMarker (),
70  fGlobalMarkerScale (1.),
71  fGlobalLineWidthScale (1.),
72  fMarkerNotHidden (true),
73  fWindowSizeHintX (600),
74  fWindowSizeHintY (600),
75  fWindowLocationHintX(0),
76  fWindowLocationHintY(0),
77  fWindowLocationHintXNegative(true),
78  fWindowLocationHintYNegative(false),
79  fGeometryMask(0),
80  fAutoRefresh (false),
81  fBackgroundColour (G4Colour(0.,0.,0.)), // Black
82  fPicking (false),
83  fRotationStyle (constrainUpDirection)
84 {
85  fDefaultMarker.SetScreenSize (5.);
86  // Markers are 5 pixels "overall" size, i.e., diameter.
87 }
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
static constexpr double cm3
Definition: G4SIunits.hh:121
void SetScreenSize(G4double)

Here is the call graph for this function:

G4ViewParameters::~G4ViewParameters ( )

Definition at line 89 of file G4ViewParameters.cc.

89 {}

Member Function Documentation

void G4ViewParameters::AddCutawayPlane ( const G4Plane3D cutawayPlane)

Definition at line 151 of file G4ViewParameters.cc.

151  {
152  if (fCutawayPlanes.size () < 3 ) {
153  fCutawayPlanes.push_back (cutawayPlane);
154  }
155  else {
156  G4cerr <<
157  "ERROR: G4ViewParameters::AddCutawayPlane:"
158  "\n A maximum of 3 cutaway planes supported." << G4endl;
159  }
160 }
#define G4endl
Definition: G4ios.hh:61
G4GLOB_DLL std::ostream G4cerr

Here is the caller graph for this function:

void G4ViewParameters::AddVisAttributesModifier ( const G4ModelingParameters::VisAttributesModifier vam)

Definition at line 257 of file G4ViewParameters.cc.

257  {
258  // If target exists with same signifier just change vis attributes.
259  G4bool duplicateTarget = false;
260  auto i = fVisAttributesModifiers.begin();
261  for (; i < fVisAttributesModifiers.end(); ++i) {
262  if (vam.GetPVNameCopyNoPath() == (*i).GetPVNameCopyNoPath() &&
263  vam.GetVisAttributesSignifier() == (*i).GetVisAttributesSignifier()) {
264  duplicateTarget = true;
265  break;
266  }
267  }
268  if (duplicateTarget) (*i).SetVisAttributes(vam.GetVisAttributes());
269  else fVisAttributesModifiers.push_back(vam);
270 }
const G4VisAttributes & GetVisAttributes() const
bool G4bool
Definition: G4Types.hh:79
const PVNameCopyNoPath & GetPVNameCopyNoPath() const
VisAttributesSignifier GetVisAttributesSignifier() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4String G4ViewParameters::CameraAndLightingCommands ( const G4Point3D  standardTargetPoint) const

Definition at line 273 of file G4ViewParameters.cc.

274 {
275  std::ostringstream oss;
276 
277  oss << "#\n# Camera and lights commands";
278 
279  oss << "\n/vis/viewer/set/viewpointVector "
280  << fViewpointDirection.x()
281  << ' ' << fViewpointDirection.y()
282  << ' ' << fViewpointDirection.z();
283 
284  oss << "\n/vis/viewer/set/upVector "
285  << fUpVector.x()
286  << ' ' << fUpVector.y()
287  << ' ' << fUpVector.z();
288 
289  oss << "\n/vis/viewer/set/projection ";
290  if (fFieldHalfAngle == 0.) {
291  oss
292  << "orthogonal";
293  } else {
294  oss
295  << "perspective "
296  << fFieldHalfAngle/deg
297  << " deg";
298  }
299 
300  oss << "\n/vis/viewer/zoomTo "
301  << fZoomFactor;
302 
303  oss << "\n/vis/viewer/scaleTo "
304  << fScaleFactor.x()
305  << ' ' << fScaleFactor.y()
306  << ' ' << fScaleFactor.z();
307 
308  oss << "\n/vis/viewer/set/targetPoint "
309  << G4BestUnit(standardTargetPoint+fCurrentTargetPoint,"Length")
310  << "\n# Note that if you have not set a target point, the vis system sets"
311  << "\n# a target point based on the scene - plus any panning and dollying -"
312  << "\n# so don't be alarmed by strange coordinates here.";
313 
314  oss << "\n/vis/viewer/dollyTo "
315  << G4BestUnit(fDolly,"Length");
316 
317  oss << "\n/vis/viewer/set/lightsMove ";
318  if (fLightsMoveWithCamera) {
319  oss << "camera";
320  } else {
321  oss << "object";
322  }
323 
324  oss << "\n/vis/viewer/set/lightsVector "
325  << fRelativeLightpointDirection.x()
326  << ' ' << fRelativeLightpointDirection.y()
327  << ' ' << fRelativeLightpointDirection.z();
328 
329  oss << "\n/vis/viewer/set/rotationStyle ";
330  if (fRotationStyle == constrainUpDirection) {
331  oss << "constrainUpDirection";
332  } else {
333  oss << "freeRotation";
334  }
335 
336  G4Colour c = fBackgroundColour;
337  oss << "\n/vis/viewer/set/background "
338  << c.GetRed()
339  << ' ' << c.GetGreen()
340  << ' ' << c.GetBlue()
341  << ' ' << c.GetAlpha();
342 
343  c = fDefaultVisAttributes.GetColour();
344  oss << "\n/vis/viewer/set/defaultColour "
345  << c.GetRed()
346  << ' ' << c.GetGreen()
347  << ' ' << c.GetBlue()
348  << ' ' << c.GetAlpha();
349 
350  c = fDefaultTextVisAttributes.GetColour();
351  oss << "\n/vis/viewer/set/defaultTextColour "
352  << c.GetRed()
353  << ' ' << c.GetGreen()
354  << ' ' << c.GetBlue()
355  << ' ' << c.GetAlpha();
356 
357  oss << std::endl;
358 
359  return oss.str();
360 }
G4double GetAlpha() const
Definition: G4Colour.hh:142
const G4Colour & GetColour() const
G4double GetBlue() const
Definition: G4Colour.hh:141
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4double GetRed() const
Definition: G4Colour.hh:139
G4double GetGreen() const
Definition: G4Colour.hh:140
static constexpr double deg
Definition: G4SIunits.hh:152
tuple c
Definition: test.py:13

Here is the call graph for this function:

G4ViewParameters * G4ViewParameters::CatmullRomCubicSplineInterpolation ( const std::vector< G4ViewParameters > &  views,
G4int  nInterpolationPoints = 50 
)
static

Definition at line 1146 of file G4ViewParameters.cc.

1148 {
1149  // Returns a null pointer when no more to be done. For example:
1150  // do {
1151  // G4ViewParameters* vp =
1152  // G4ViewParameters::CatmullRomCubicSplineInterpolation(viewVector,nInterpolationPoints);
1153  // if (!vp) break;
1154  // ...
1155  // } while (true);
1156 
1157  // See https://en.wikipedia.org/wiki/Cubic_Hermite_spline
1158 
1159  // Assumes equal intervals
1160 
1161  if (views.size() < 2) {
1162  G4Exception
1163  ("G4ViewParameters::CatmullRomCubicSplineInterpolation",
1164  "visman0301", JustWarning,
1165  "There must be at least two views.");
1166  return 0;
1167  }
1168 
1169  if (nInterpolationPoints < 1) {
1170  G4Exception
1171  ("G4ViewParameters::CatmullRomCubicSplineInterpolation",
1172  "visman0302", JustWarning,
1173  "Number of interpolation points cannot be zero or negative.");
1174  return 0;
1175  }
1176 
1177  const size_t nIntervals = views.size() - 1;
1178  const G4double dt = 1./nInterpolationPoints;
1179 
1180  static G4ViewParameters holdingValues;
1181  static G4double t = 0.; // 0. <= t <= 1.
1182  static G4int iInterpolationPoint = 0;
1183  static size_t iInterval = 0;
1184 
1185 // G4cout << "Interval " << iInterval << ", t = " << t << G4endl;
1186 
1187  // Hermite polynomials.
1188  const G4double h00 = 2.*t*t*t - 3.*t*t +1;
1189  const G4double h10 = t*t*t -2.*t*t + t;
1190  const G4double h01 = -2.*t*t*t + 3.*t*t;
1191  const G4double h11 = t*t*t - t*t;
1192 
1193  // Aliases (to simplify code)
1194  const size_t& n = nIntervals;
1195  size_t& i = iInterval;
1196  const std::vector<G4ViewParameters>& v = views;
1197 
1198  // The Catmull-Rom cubic spline prescription is as follows:
1199  // Slope at first way point is v[1] - v[0].
1200  // Slope at last way point is v[n] - v[n-1].
1201  // Otherwise slope at way point i is 0.5*(v[i+1] - v[i-1]).
1202  // Result = h00*v[i] + h10*m[i] + h01*v[i+1] + h11*m[i+1],
1203  // where m[i] amd m[i+1] are the slopes at the start and end
1204  // of the interval for the particular value.
1205  // If (n == 1), linear interpolation results.
1206  // If (n == 2), quadratic interpolation results.
1207 
1208  // Working variables
1209  G4double mi, mi1, real, x, y, z;
1210 
1211  // First, a crude interpolation of all parameters. Then, below, a
1212  // smooth interpolation of those for which it makes sense.
1213  holdingValues = t < 0.5? v[i]: v[i+1];
1214 
1215  // Catmull-Rom cubic spline interpolation
1216 #define INTERPOLATE(param) \
1217  /* This works out the interpolated param in i'th interval */ \
1218  /* Assumes n >= 1 */ \
1219  if (i == 0) { \
1220  /* First interval */ \
1221  mi = v[1].param - v[0].param; \
1222  /* If there is only one interval, make start and end slopes equal */ \
1223  /* (This results in a linear interpolation) */ \
1224  if (n == 1) mi1 = mi; \
1225  /* else the end slope of the interval takes account of the next waypoint along */ \
1226  else mi1 = 0.5 * (v[2].param - v[0].param); \
1227  } else if (i >= n - 1) { \
1228  /* Similarly for last interval */ \
1229  mi1 = v[i+1].param - v[i].param; \
1230  /* If there is only one interval, make start and end slopes equal */ \
1231  if (n == 1) mi = mi1; \
1232  /* else the start slope of the interval takes account of the previous waypoint */ \
1233  else mi = 0.5 * (v[i+1].param - v[i-1].param); \
1234  } else { \
1235  /* Full Catmull-Rom slopes use previous AND next waypoints */ \
1236  mi = 0.5 * (v[i+1].param - v[i-1].param); \
1237  mi1 = 0.5 * (v[i+2].param - v[i ].param); \
1238  } \
1239  real = h00 * v[i].param + h10 * mi + h01 * v[i+1].param + h11 * mi1;
1240 
1241  // Real parameters
1242  INTERPOLATE(fVisibleDensity);
1243  if (real < 0.) real = 0.;
1244  holdingValues.fVisibleDensity = real;
1245  INTERPOLATE(fExplodeFactor);
1246  if (real < 0.) real = 0.;
1247  holdingValues.fExplodeFactor = real;
1248  INTERPOLATE(fFieldHalfAngle);
1249  if (real < 0.) real = 0.;
1250  holdingValues.fFieldHalfAngle = real;
1251  INTERPOLATE(fZoomFactor);
1252  if (real < 0.) real = 0.;
1253  holdingValues.fZoomFactor = real;
1254  INTERPOLATE(fDolly);
1255  holdingValues.fDolly = real;
1256  INTERPOLATE(fGlobalMarkerScale);
1257  if (real < 0.) real = 0.;
1258  holdingValues.fGlobalMarkerScale = real;
1259  INTERPOLATE(fGlobalLineWidthScale);
1260  if (real < 0.) real = 0.;
1261  holdingValues.fGlobalLineWidthScale = real;
1262 
1263  // Unit vectors
1264 #define INTERPOLATEUNITVECTOR(vector) \
1265 INTERPOLATE(vector.x()); x = real; \
1266 INTERPOLATE(vector.y()); y = real; \
1267 INTERPOLATE(vector.z()); z = real;
1268  INTERPOLATEUNITVECTOR(fViewpointDirection);
1269  holdingValues.fViewpointDirection = G4Vector3D(x,y,z).unit();
1270  INTERPOLATEUNITVECTOR(fUpVector);
1271  holdingValues.fUpVector = G4Vector3D(x,y,z).unit();
1272  INTERPOLATEUNITVECTOR(fRelativeLightpointDirection);
1273  holdingValues.fRelativeLightpointDirection = G4Vector3D(x,y,z).unit();
1274  INTERPOLATEUNITVECTOR(fActualLightpointDirection);
1275  holdingValues.fActualLightpointDirection = G4Vector3D(x,y,z).unit();
1276 
1277  // Un-normalised vectors
1278 #define INTERPOLATEVECTOR(vector) \
1279 INTERPOLATE(vector.x()); x = real; \
1280 INTERPOLATE(vector.y()); y = real; \
1281 INTERPOLATE(vector.z()); z = real;
1282  INTERPOLATEVECTOR(fScaleFactor);
1283  holdingValues.fScaleFactor = G4Vector3D(x,y,z);
1284 
1285  // Points
1286 #define INTERPOLATEPOINT(point) \
1287 INTERPOLATE(point.x()); x = real; \
1288 INTERPOLATE(point.y()); y = real; \
1289 INTERPOLATE(point.z()); z = real;
1290  INTERPOLATEPOINT(fExplodeCentre);
1291  holdingValues.fExplodeCentre = G4Point3D(x,y,z);
1292  INTERPOLATEPOINT(fCurrentTargetPoint);
1293  holdingValues.fCurrentTargetPoint = G4Point3D(x,y,z);
1294 
1295  // Colour
1296  G4double red, green, blue, alpha;
1297 #define INTERPOLATECOLOUR(colour) \
1298 INTERPOLATE(colour.GetRed()); red = real; \
1299 INTERPOLATE(colour.GetGreen()); green = real; \
1300 INTERPOLATE(colour.GetBlue()); blue = real; \
1301 INTERPOLATE(colour.GetAlpha()); alpha = real;
1302  INTERPOLATECOLOUR(fBackgroundColour);
1303  // Components are clamped to 0. <= component <= 1.
1304  holdingValues.fBackgroundColour = G4Colour(red,green,blue,alpha);
1305 
1306  // For some parameters we need to check some continuity
1307  G4bool continuous;
1308 #define CONTINUITY(quantity) \
1309  continuous = false; \
1310  /* This follows the logic of the INTERPOLATE macro above; see comments therein */ \
1311  if (i == 0) { \
1312  if (v[1].quantity == v[0].quantity) { \
1313  if (n == 1) continuous = true; \
1314  else if (v[2].quantity == v[0].quantity) \
1315  continuous = true; \
1316  } \
1317  } else if (i >= n - 1) { \
1318  if (v[i+1].quantity == v[i].quantity) { \
1319  if (n == 1) continuous = true; \
1320  else if (v[i+1].quantity == v[i-1].quantity) \
1321  continuous = true; \
1322  } \
1323  } else { \
1324  if (v[i-1].quantity == v[i].quantity && \
1325  v[i+1].quantity == v[i].quantity && \
1326  v[i+2].quantity == v[i].quantity) \
1327  continuous = true; \
1328  }
1329 
1330  G4double a, b, c, d;
1331 #define INTERPOLATEPLANE(plane) \
1332 INTERPOLATE(plane.a()); a = real; \
1333 INTERPOLATE(plane.b()); b = real; \
1334 INTERPOLATE(plane.c()); c = real; \
1335 INTERPOLATE(plane.d()); d = real;
1336 
1337  // Section plane
1338  CONTINUITY(fSection);
1339  if (continuous) {
1340  INTERPOLATEPLANE(fSectionPlane);
1341  holdingValues.fSectionPlane = G4Plane3D(a,b,c,d);
1342  }
1343 
1344  // Cutaway planes
1345  if (v[i].fCutawayPlanes.size()) {
1346  CONTINUITY(fCutawayPlanes.size());
1347  if (continuous) {
1348  for (size_t j = 0; j < v[i].fCutawayPlanes.size(); ++j) {
1349  INTERPOLATEPLANE(fCutawayPlanes[j]);
1350  holdingValues.fCutawayPlanes[j] = G4Plane3D(a,b,c,d);
1351  }
1352  }
1353  }
1354 
1355  // Vis attributes modifiers
1356  // Really, we are only interested in colour - other attributes can follow
1357  // the "crude" interpolation that is guaranteed above.
1358  static G4VisAttributes workingVA;
1359  if (v[i].fVisAttributesModifiers.size()) {
1360  CONTINUITY(fVisAttributesModifiers.size());
1361  if (continuous) { \
1362  for (size_t j = 0; j < v[i].fVisAttributesModifiers.size(); ++j) { \
1363  CONTINUITY(fVisAttributesModifiers[j].GetPVNameCopyNoPath());
1364  if (continuous) {
1365  CONTINUITY(fVisAttributesModifiers[j].GetVisAttributesSignifier());
1366  if (continuous) {
1367  if (v[i].fVisAttributesModifiers[j].GetVisAttributesSignifier() ==
1369  INTERPOLATECOLOUR(fVisAttributesModifiers[j].GetVisAttributes().GetColour());
1370  workingVA = v[i].fVisAttributesModifiers[j].GetVisAttributes();
1371  workingVA.SetColour(G4Colour(red,green,blue,alpha));
1372  holdingValues.fVisAttributesModifiers[j].SetVisAttributes(workingVA);
1373  }
1374  }
1375  }
1376  }
1377  }
1378  }
1379 
1380  // Increment counters
1381  iInterpolationPoint++;
1382  t += dt;
1383  if (iInterpolationPoint > nInterpolationPoints) {
1384  iInterpolationPoint = 1; // Ready for next interval.
1385  t = dt;
1386  iInterval++;
1387  }
1388  if (iInterval >= nIntervals) {
1389  iInterpolationPoint = 0; // Ready for a complete restart.
1390  t = 0.;
1391  iInterval = 0;
1392  return 0;
1393  }
1394 
1395  return &holdingValues;
1396 }
void SetColour(const G4Colour &)
Definition: test07.cc:36
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
#define INTERPOLATEUNITVECTOR(vector)
BasicVector3D< T > unit() const
#define INTERPOLATECOLOUR(colour)
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:35
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:35
tuple x
Definition: test.py:50
Definition: test07.cc:36
int G4int
Definition: G4Types.hh:78
#define INTERPOLATEPLANE(plane)
tuple b
Definition: test.py:12
bool G4bool
Definition: G4Types.hh:79
#define INTERPOLATEVECTOR(vector)
const G4int n
tuple v
Definition: test.py:18
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define INTERPOLATE(param)
tuple z
Definition: test.py:28
HepGeom::Plane3D< G4double > G4Plane3D
Definition: G4Plane3D.hh:37
double G4double
Definition: G4Types.hh:76
tuple c
Definition: test.py:13
static const G4double alpha
for(G4int i1=0;i1< theStableOnes.GetNumberOfIsotopes(static_cast< G4int >(anE->GetZ()));i1++)
#define INTERPOLATEPOINT(point)
#define CONTINUITY(quantity)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ViewParameters::ChangeCutawayPlane ( size_t  index,
const G4Plane3D cutawayPlane 
)

Definition at line 163 of file G4ViewParameters.cc.

163  {
164  if (index >= fCutawayPlanes.size()) {
165  G4cerr <<
166  "ERROR: G4ViewParameters::ChangeCutawayPlane:"
167  "\n Plane " << index << " does not exist." << G4endl;
168  } else {
169  fCutawayPlanes[index] = cutawayPlane;
170  }
171 }
#define G4endl
Definition: G4ios.hh:61
G4GLOB_DLL std::ostream G4cerr

Here is the caller graph for this function:

void G4ViewParameters::ClearCutawayPlanes ( )

Here is the caller graph for this function:

void G4ViewParameters::ClearVisAttributesModifiers ( )

Here is the caller graph for this function:

G4String G4ViewParameters::DrawingStyleCommands ( ) const

Definition at line 362 of file G4ViewParameters.cc.

363 {
364  std::ostringstream oss;
365 
366  oss << "#\n# Drawing style commands";
367 
368  oss << "\n/vis/viewer/set/style ";
369  if (fDrawingStyle == wireframe || fDrawingStyle == hlr) {
370  oss << "wireframe";
371  } else {
372  oss << "surface";
373  }
374 
375  oss << "\n/vis/viewer/set/hiddenEdge ";
376  if (fDrawingStyle == hlr || fDrawingStyle == hlhsr) {
377  oss << "true";
378  } else {
379  oss << "false";
380  }
381 
382  oss << "\n/vis/viewer/set/auxiliaryEdge ";
383  if (fAuxEdgeVisible) {
384  oss << "true";
385  } else {
386  oss << "false";
387  }
388 
389  oss << "\n/vis/viewer/set/hiddenMarker ";
390  if (fMarkerNotHidden) {
391  oss << "false";
392  } else {
393  oss << "true";
394  }
395 
396  oss << "\n/vis/viewer/set/globalLineWidthScale "
397  << fGlobalLineWidthScale;
398 
399  oss << "\n/vis/viewer/set/globalMarkerScale "
400  << fGlobalMarkerScale;
401 
402  oss << std::endl;
403 
404  return oss.str();
405 }
G4Vector3D & G4ViewParameters::GetActualLightpointDirection ( )

Definition at line 98 of file G4ViewParameters.cc.

98  {
99  SetViewAndLights (fViewpointDirection);
100  return fActualLightpointDirection;
101 }
void SetViewAndLights(const G4Vector3D &viewpointDirection)

Here is the call graph for this function:

Here is the caller graph for this function:

const G4Colour& G4ViewParameters::GetBackgroundColour ( ) const

Here is the caller graph for this function:

G4double G4ViewParameters::GetCameraDistance ( G4double  radius) const

Definition at line 111 of file G4ViewParameters.cc.

111  {
112  G4double cameraDistance;
113  if (fFieldHalfAngle == 0.) {
114  cameraDistance = radius;
115  }
116  else {
117  cameraDistance = radius / std::sin (fFieldHalfAngle) - fDolly;
118  }
119  return cameraDistance;
120 }
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

const G4Point3D& G4ViewParameters::GetCurrentTargetPoint ( ) const

Here is the caller graph for this function:

CutawayMode G4ViewParameters::GetCutawayMode ( ) const

Here is the caller graph for this function:

const G4Planes& G4ViewParameters::GetCutawayPlanes ( ) const

Here is the caller graph for this function:

const G4VMarker& G4ViewParameters::GetDefaultMarker ( ) const
const G4VisAttributes* G4ViewParameters::GetDefaultTextVisAttributes ( ) const

Here is the caller graph for this function:

const G4VisAttributes* G4ViewParameters::GetDefaultVisAttributes ( ) const

Here is the caller graph for this function:

G4double G4ViewParameters::GetDolly ( ) const

Here is the caller graph for this function:

DrawingStyle G4ViewParameters::GetDrawingStyle ( ) const

Here is the caller graph for this function:

const G4Point3D& G4ViewParameters::GetExplodeCentre ( ) const

Here is the caller graph for this function:

G4double G4ViewParameters::GetExplodeFactor ( ) const

Here is the caller graph for this function:

G4double G4ViewParameters::GetFarDistance ( G4double  cameraDistance,
G4double  nearDistance,
G4double  radius 
) const

Definition at line 130 of file G4ViewParameters.cc.

132  {
133  G4double farDistance = cameraDistance + radius;
134  if (farDistance < nearDistance) farDistance = nearDistance;
135  return farDistance;
136 }
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

G4double G4ViewParameters::GetFieldHalfAngle ( ) const

Here is the caller graph for this function:

G4double G4ViewParameters::GetFrontHalfHeight ( G4double  nearDistance,
G4double  radius 
) const

Definition at line 138 of file G4ViewParameters.cc.

139  {
140  G4double frontHalfHeight;
141  if (fFieldHalfAngle == 0.) {
142  frontHalfHeight = radius / fZoomFactor;
143  }
144  else {
145  frontHalfHeight = nearDistance * std::tan (fFieldHalfAngle) / fZoomFactor;
146  }
147  return frontHalfHeight;
148 }
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

G4double G4ViewParameters::GetGlobalLineWidthScale ( ) const

Here is the caller graph for this function:

G4double G4ViewParameters::GetGlobalMarkerScale ( ) const

Here is the caller graph for this function:

const G4Vector3D& G4ViewParameters::GetLightpointDirection ( ) const

Here is the caller graph for this function:

G4bool G4ViewParameters::GetLightsMoveWithCamera ( ) const

Here is the caller graph for this function:

G4double G4ViewParameters::GetNearDistance ( G4double  cameraDistance,
G4double  radius 
) const

Definition at line 122 of file G4ViewParameters.cc.

123  {
124  const G4double small = 1.e-6 * radius;
125  G4double nearDistance = cameraDistance - radius;
126  if (nearDistance < small) nearDistance = small;
127  return nearDistance;
128 }
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

G4int G4ViewParameters::GetNoOfSides ( ) const

Here is the caller graph for this function:

RotationStyle G4ViewParameters::GetRotationStyle ( ) const

Here is the caller graph for this function:

const G4Vector3D& G4ViewParameters::GetScaleFactor ( ) const

Here is the caller graph for this function:

const G4Plane3D& G4ViewParameters::GetSectionPlane ( ) const

Here is the caller graph for this function:

const G4Vector3D& G4ViewParameters::GetUpVector ( ) const

Here is the caller graph for this function:

const G4Vector3D& G4ViewParameters::GetViewpointDirection ( ) const

Here is the caller graph for this function:

const std::vector<G4ModelingParameters::VisAttributesModifier>& G4ViewParameters::GetVisAttributesModifiers ( ) const

Here is the caller graph for this function:

G4double G4ViewParameters::GetVisibleDensity ( ) const

Here is the caller graph for this function:

G4int G4ViewParameters::GetWindowAbsoluteLocationHintX ( G4int  sizeX) const

Definition at line 1001 of file G4ViewParameters.cc.

1001  {
1002  if ( fWindowLocationHintXNegative ) {
1003  return sizeX + fWindowLocationHintX - fWindowSizeHintX;
1004  }
1005  return fWindowLocationHintX;
1006 }
G4int G4ViewParameters::GetWindowAbsoluteLocationHintY ( G4int  sizeY) const

Definition at line 1008 of file G4ViewParameters.cc.

1008  {
1009  if ( fWindowLocationHintYNegative ) {
1010  return sizeY + fWindowLocationHintY - fWindowSizeHintY;
1011  }
1012  return fWindowLocationHintY;
1013 }
G4int G4ViewParameters::GetWindowLocationHintX ( ) const
G4int G4ViewParameters::GetWindowLocationHintY ( ) const
unsigned int G4ViewParameters::GetWindowSizeHintX ( ) const
unsigned int G4ViewParameters::GetWindowSizeHintY ( ) const
const G4String& G4ViewParameters::GetXGeometryString ( ) const

Here is the caller graph for this function:

G4double G4ViewParameters::GetZoomFactor ( ) const

Here is the caller graph for this function:

void G4ViewParameters::IncrementDolly ( G4double  dollyIncrement)

Here is the caller graph for this function:

void G4ViewParameters::IncrementPan ( G4double  right,
G4double  up 
)

Definition at line 246 of file G4ViewParameters.cc.

246  {
247  IncrementPan (right,up, 0);
248 }
void IncrementPan(G4double right, G4double up)

Here is the caller graph for this function:

void G4ViewParameters::IncrementPan ( G4double  right,
G4double  up,
G4double  forward 
)

Definition at line 250 of file G4ViewParameters.cc.

250  {
251  G4Vector3D unitRight = (fUpVector.cross (fViewpointDirection)).unit();
252  G4Vector3D unitUp = (fViewpointDirection.cross (unitRight)).unit();
253  fCurrentTargetPoint += right * unitRight + up * unitUp + distance * fViewpointDirection;
254 }
BasicVector3D< T > cross(const BasicVector3D< T > &v) const

Here is the call graph for this function:

G4bool G4ViewParameters::IsAutoRefresh ( ) const

Here is the caller graph for this function:

G4bool G4ViewParameters::IsAuxEdgeVisible ( ) const

Here is the caller graph for this function:

G4bool G4ViewParameters::IsCulling ( ) const

Here is the caller graph for this function:

G4bool G4ViewParameters::IsCullingCovered ( ) const

Here is the caller graph for this function:

G4bool G4ViewParameters::IsCullingInvisible ( ) const

Here is the caller graph for this function:

G4bool G4ViewParameters::IsCutaway ( ) const

Here is the caller graph for this function:

G4bool G4ViewParameters::IsDensityCulling ( ) const

Here is the caller graph for this function:

G4bool G4ViewParameters::IsExplode ( ) const

Here is the caller graph for this function:

G4bool G4ViewParameters::IsMarkerNotHidden ( ) const

Here is the caller graph for this function:

G4bool G4ViewParameters::IsPicking ( ) const

Here is the caller graph for this function:

G4bool G4ViewParameters::IsSection ( ) const

Here is the caller graph for this function:

bool G4ViewParameters::IsWindowLocationHintX ( ) const
bool G4ViewParameters::IsWindowLocationHintY ( ) const
bool G4ViewParameters::IsWindowSizeHintX ( ) const
bool G4ViewParameters::IsWindowSizeHintY ( ) const
void G4ViewParameters::MultiplyScaleFactor ( const G4Vector3D scaleFactorMultiplier)

Definition at line 92 of file G4ViewParameters.cc.

92  {
93  fScaleFactor.setX(fScaleFactor.x() * scaleFactorMultiplier.x());
94  fScaleFactor.setY(fScaleFactor.y() * scaleFactorMultiplier.y());
95  fScaleFactor.setZ(fScaleFactor.z() * scaleFactorMultiplier.z());
96 }

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ViewParameters::MultiplyZoomFactor ( G4double  zoomFactorMultiplier)

Here is the caller graph for this function:

G4bool G4ViewParameters::operator!= ( const G4ViewParameters v) const

Definition at line 850 of file G4ViewParameters.cc.

850  {
851 
852  // Put performance-sensitive parameters first.
853  if (
854  // This first to optimise spin, etc.
855  (fViewpointDirection != v.fViewpointDirection) ||
856 
857  // No particular order from here on.
858  (fDrawingStyle != v.fDrawingStyle) ||
859  (fAuxEdgeVisible != v.fAuxEdgeVisible) ||
860  (fCulling != v.fCulling) ||
861  (fCullInvisible != v.fCullInvisible) ||
862  (fDensityCulling != v.fDensityCulling) ||
863  (fCullCovered != v.fCullCovered) ||
864  (fSection != v.fSection) ||
865  (IsCutaway() != v.IsCutaway()) ||
866  (IsExplode() != v.IsExplode()) ||
867  (fNoOfSides != v.fNoOfSides) ||
868  (fUpVector != v.fUpVector) ||
869  (fFieldHalfAngle != v.fFieldHalfAngle) ||
870  (fZoomFactor != v.fZoomFactor) ||
871  (fScaleFactor != v.fScaleFactor) ||
872  (fCurrentTargetPoint != v.fCurrentTargetPoint) ||
873  (fDolly != v.fDolly) ||
874  (fRelativeLightpointDirection != v.fRelativeLightpointDirection) ||
875  (fLightsMoveWithCamera != v.fLightsMoveWithCamera) ||
876  (fDefaultVisAttributes != v.fDefaultVisAttributes) ||
877  (fDefaultTextVisAttributes != v.fDefaultTextVisAttributes) ||
878  (fDefaultMarker != v.fDefaultMarker) ||
879  (fGlobalMarkerScale != v.fGlobalMarkerScale) ||
880  (fGlobalLineWidthScale != v.fGlobalLineWidthScale) ||
881  (fMarkerNotHidden != v.fMarkerNotHidden) ||
882  (fWindowSizeHintX != v.fWindowSizeHintX) ||
883  (fWindowSizeHintY != v.fWindowSizeHintY) ||
884  (fXGeometryString != v.fXGeometryString) ||
885  (fGeometryMask != v.fGeometryMask) ||
886  (fAutoRefresh != v.fAutoRefresh) ||
887  (fBackgroundColour != v.fBackgroundColour) ||
888  (fPicking != v.fPicking) ||
889  (fRotationStyle != v.fRotationStyle)
890  )
891  return true;
892 
893  if (fDensityCulling &&
894  (fVisibleDensity != v.fVisibleDensity)) return true;
895 
896  if (fSection &&
897  (!(fSectionPlane == v.fSectionPlane))) return true;
898 
899  if (IsCutaway()) {
900  if (fCutawayPlanes.size () != v.fCutawayPlanes.size ())
901  return true;
902  else {
903  for (size_t i = 0; i < fCutawayPlanes.size (); i++) {
904  if (!(fCutawayPlanes[i] == v.fCutawayPlanes[i])) return true;
905  }
906  }
907  }
908 
909  if (IsExplode() &&
910  ((fExplodeFactor != v.fExplodeFactor) ||
911  (fExplodeCentre != v.fExplodeCentre))) return true;
912 
913  if (fVisAttributesModifiers != v.fVisAttributesModifiers) return true;
914 
915  return false;
916 }
G4bool IsExplode() const
G4bool IsCutaway() const

Here is the call graph for this function:

void G4ViewParameters::PrintDifferences ( const G4ViewParameters v) const

Definition at line 610 of file G4ViewParameters.cc.

610  {
611 
612  // Put performance-sensitive parameters first.
613  if (
614  // This first to optimise spin, etc.
615  (fViewpointDirection != v.fViewpointDirection) ||
616 
617  // No particular order from here on.
618  (fDrawingStyle != v.fDrawingStyle) ||
619  (fAuxEdgeVisible != v.fAuxEdgeVisible) ||
620  (fCulling != v.fCulling) ||
621  (fCullInvisible != v.fCullInvisible) ||
622  (fDensityCulling != v.fDensityCulling) ||
623  (fVisibleDensity != v.fVisibleDensity) ||
624  (fCullCovered != v.fCullCovered) ||
625  (fSection != v.fSection) ||
626  (fNoOfSides != v.fNoOfSides) ||
627  (fUpVector != v.fUpVector) ||
628  (fFieldHalfAngle != v.fFieldHalfAngle) ||
629  (fZoomFactor != v.fZoomFactor) ||
630  (fScaleFactor != v.fScaleFactor) ||
631  (fCurrentTargetPoint != v.fCurrentTargetPoint) ||
632  (fDolly != v.fDolly) ||
633  (fRelativeLightpointDirection != v.fRelativeLightpointDirection) ||
634  (fLightsMoveWithCamera != v.fLightsMoveWithCamera) ||
635  (fDefaultVisAttributes != v.fDefaultVisAttributes) ||
636  (fDefaultTextVisAttributes != v.fDefaultTextVisAttributes) ||
637  (fDefaultMarker != v.fDefaultMarker) ||
638  (fGlobalMarkerScale != v.fGlobalMarkerScale) ||
639  (fGlobalLineWidthScale != v.fGlobalLineWidthScale) ||
640  (fMarkerNotHidden != v.fMarkerNotHidden) ||
641  (fWindowSizeHintX != v.fWindowSizeHintX) ||
642  (fWindowSizeHintY != v.fWindowSizeHintY) ||
643  (fXGeometryString != v.fXGeometryString) ||
644  (fGeometryMask != v.fGeometryMask) ||
645  (fAutoRefresh != v.fAutoRefresh) ||
646  (fBackgroundColour != v.fBackgroundColour) ||
647  (fPicking != v.fPicking) ||
648  (fRotationStyle != v.fRotationStyle)
649  )
650  G4cout << "Difference in 1st batch." << G4endl;
651 
652  if (fSection) {
653  if (!(fSectionPlane == v.fSectionPlane))
654  G4cout << "Difference in section planes batch." << G4endl;
655  }
656 
657  if (IsCutaway()) {
658  if (fCutawayPlanes.size () != v.fCutawayPlanes.size ()) {
659  G4cout << "Difference in no of cutaway planes." << G4endl;
660  }
661  else {
662  for (size_t i = 0; i < fCutawayPlanes.size (); i++) {
663  if (!(fCutawayPlanes[i] == v.fCutawayPlanes[i]))
664  G4cout << "Difference in cutaway plane no. " << i << G4endl;
665  }
666  }
667  }
668 
669  if (IsExplode()) {
670  if (fExplodeFactor != v.fExplodeFactor)
671  G4cout << "Difference in explode factor." << G4endl;
672  if (fExplodeCentre != v.fExplodeCentre)
673  G4cout << "Difference in explode centre." << G4endl;
674  }
675 
676  if (fVisAttributesModifiers != v.fVisAttributesModifiers) {
677  G4cout << "Difference in vis attributes modifiers." << G4endl;
678  }
679 }
G4GLOB_DLL std::ostream G4cout
G4bool IsExplode() const
G4bool IsCutaway() const
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4String G4ViewParameters::SceneModifyingCommands ( ) const

Definition at line 407 of file G4ViewParameters.cc.

408 {
409  std::ostringstream oss;
410 
411  oss << "#\n# Scene-modifying commands";
412 
413  oss << "\n/vis/viewer/set/culling global ";
414  if (fCulling) {
415  oss << "true";
416  } else {
417  oss << "false";
418  }
419 
420  oss << "\n/vis/viewer/set/culling invisible ";
421  if (fCullInvisible) {
422  oss << "true";
423  } else {
424  oss << "false";
425  }
426 
427  oss << "\n/vis/viewer/set/culling density ";
428  if (fDensityCulling) {
429  oss << "true " << fVisibleDensity/(g/cm3) << " g/cm3";
430  } else {
431  oss << "false";
432  }
433 
434  oss << "\n/vis/viewer/set/culling coveredDaughters ";
435  if (fCullCovered) {
436  oss << "true";
437  } else {
438  oss << "false";
439  }
440 
441  oss << "\n/vis/viewer/set/sectionPlane ";
442  if (fSection) {
443  oss << "on "
444  << G4BestUnit(fSectionPlane.point(),"Length")
445  << fSectionPlane.normal().x()
446  << ' ' << fSectionPlane.normal().y()
447  << ' ' << fSectionPlane.normal().z();
448  } else {
449  oss << "off";
450  }
451 
452  oss << "\n/vis/viewer/set/cutawayMode ";
453  if (fCutawayMode == cutawayUnion) {
454  oss << "union";
455  } else {
456  oss << "intersection";
457  }
458 
459  oss << "\n/vis/viewer/clearCutawayPlanes";
460  if (fCutawayPlanes.size()) {
461  for (size_t i = 0; i < fCutawayPlanes.size(); i++) {
462  oss << "\n/vis/viewer/addCutawayPlane "
463  << G4BestUnit(fCutawayPlanes[i].point(),"Length")
464  << fCutawayPlanes[i].normal().x()
465  << ' ' << fCutawayPlanes[i].normal().y()
466  << ' ' << fCutawayPlanes[i].normal().z();
467  }
468  } else {
469  oss << "\n# No cutaway planes defined.";
470  }
471 
472  oss << "\n/vis/viewer/set/explodeFactor "
473  << fExplodeFactor
474  << ' ' << G4BestUnit(fExplodeCentre,"Length");
475 
476  oss << "\n/vis/viewer/set/lineSegmentsPerCircle "
477  << fNoOfSides;
478 
479  oss << std::endl;
480 
481  return oss.str();
482 }
Normal3D< T > normal() const
Definition: Plane3D.h:90
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
Point3D< T > point(const Point3D< T > &p) const
Definition: Plane3D.h:108
static constexpr double cm3
Definition: G4SIunits.hh:121

Here is the call graph for this function:

void G4ViewParameters::SetAutoRefresh ( G4bool  )

Here is the caller graph for this function:

void G4ViewParameters::SetAuxEdgeVisible ( G4bool  )

Here is the caller graph for this function:

void G4ViewParameters::SetBackgroundColour ( const G4Colour )

Here is the caller graph for this function:

void G4ViewParameters::SetCulling ( G4bool  )

Here is the caller graph for this function:

void G4ViewParameters::SetCullingCovered ( G4bool  )

Here is the caller graph for this function:

void G4ViewParameters::SetCullingInvisible ( G4bool  )

Here is the caller graph for this function:

void G4ViewParameters::SetCurrentTargetPoint ( const G4Point3D currentTargetPoint)

Here is the caller graph for this function:

void G4ViewParameters::SetCutawayMode ( CutawayMode  )

Here is the caller graph for this function:

void G4ViewParameters::SetDefaultColour ( const G4Colour )
void G4ViewParameters::SetDefaultMarker ( const G4VMarker defaultMarker)
void G4ViewParameters::SetDefaultTextColour ( const G4Colour )
void G4ViewParameters::SetDefaultTextVisAttributes ( const G4VisAttributes )

Here is the caller graph for this function:

void G4ViewParameters::SetDefaultVisAttributes ( const G4VisAttributes )

Here is the caller graph for this function:

void G4ViewParameters::SetDensityCulling ( G4bool  )

Here is the caller graph for this function:

void G4ViewParameters::SetDolly ( G4double  dolly)

Here is the caller graph for this function:

void G4ViewParameters::SetDrawingStyle ( G4ViewParameters::DrawingStyle  style)

Here is the caller graph for this function:

void G4ViewParameters::SetExplodeCentre ( const G4Point3D explodeCentre)

Here is the caller graph for this function:

void G4ViewParameters::SetExplodeFactor ( G4double  explodeFactor)

Here is the caller graph for this function:

void G4ViewParameters::SetFieldHalfAngle ( G4double  fieldHalfAngle)

Here is the caller graph for this function:

void G4ViewParameters::SetGlobalLineWidthScale ( G4double  globalLineWidthScale)

Here is the caller graph for this function:

void G4ViewParameters::SetGlobalMarkerScale ( G4double  globalMarkerScale)

Here is the caller graph for this function:

void G4ViewParameters::SetLightpointDirection ( const G4Vector3D lightpointDirection)

Definition at line 235 of file G4ViewParameters.cc.

235  {
236  fRelativeLightpointDirection = lightpointDirection;
237  SetViewAndLights (fViewpointDirection);
238 }
void SetViewAndLights(const G4Vector3D &viewpointDirection)

Here is the caller graph for this function:

void G4ViewParameters::SetLightsMoveWithCamera ( G4bool  moves)

Here is the caller graph for this function:

void G4ViewParameters::SetMarkerHidden ( )

Here is the caller graph for this function:

void G4ViewParameters::SetMarkerNotHidden ( )

Here is the caller graph for this function:

G4int G4ViewParameters::SetNoOfSides ( G4int  nSides)

Definition at line 190 of file G4ViewParameters.cc.

190  {
191  const G4int nSidesMin = fDefaultVisAttributes.GetMinLineSegmentsPerCircle();
192  if (nSides < nSidesMin) {
193  nSides = nSidesMin;
194  G4cout << "G4ViewParameters::SetNoOfSides: attempt to set the"
195  "\nnumber of sides per circle < " << nSidesMin
196  << "; forced to " << nSides << G4endl;
197  }
198  fNoOfSides = nSides;
199  return fNoOfSides;
200 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
static G4int GetMinLineSegmentsPerCircle()
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ViewParameters::SetOrthogonalProjection ( )
void G4ViewParameters::SetPan ( G4double  right,
G4double  up 
)

Definition at line 240 of file G4ViewParameters.cc.

240  {
241  G4Vector3D unitRight = (fUpVector.cross (fViewpointDirection)).unit();
242  G4Vector3D unitUp = (fViewpointDirection.cross (unitRight)).unit();
243  fCurrentTargetPoint = right * unitRight + up * unitUp;
244 }
BasicVector3D< T > cross(const BasicVector3D< T > &v) const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ViewParameters::SetPerspectiveProjection ( G4double  fieldHalfAngle = 30.*CLHEP::deg)
void G4ViewParameters::SetPicking ( G4bool  )

Here is the caller graph for this function:

void G4ViewParameters::SetRotationStyle ( RotationStyle  )

Here is the caller graph for this function:

void G4ViewParameters::SetScaleFactor ( const G4Vector3D scaleFactor)

Here is the caller graph for this function:

void G4ViewParameters::SetSectionPlane ( const G4Plane3D sectionPlane)

Here is the caller graph for this function:

void G4ViewParameters::SetUpVector ( const G4Vector3D upVector)

Here is the caller graph for this function:

void G4ViewParameters::SetViewAndLights ( const G4Vector3D viewpointDirection)

Definition at line 203 of file G4ViewParameters.cc.

203  {
204 
205  fViewpointDirection = viewpointDirection;
206 
207  // If the requested viewpoint direction is parallel to the up
208  // vector, the orientation of the view is undefined...
209  if (fViewpointDirection.unit() * fUpVector.unit() > .9999) {
210  static G4bool firstTime = true;
211  if (firstTime) {
212  firstTime = false;
213  G4cout <<
214  "WARNING: Viewpoint direction is very close to the up vector direction."
215  "\n Consider setting the up vector to obtain definable behaviour."
216  << G4endl;
217  }
218  }
219 
220  // Move the lights too if requested...
221  if (fLightsMoveWithCamera) {
222  G4Vector3D zprime = fViewpointDirection.unit ();
223  G4Vector3D xprime = (fUpVector.cross (zprime)).unit ();
224  G4Vector3D yprime = zprime.cross (xprime);
225  fActualLightpointDirection =
226  fRelativeLightpointDirection.x () * xprime +
227  fRelativeLightpointDirection.y () * yprime +
228  fRelativeLightpointDirection.x () * zprime;
229  } else {
230  fActualLightpointDirection = fRelativeLightpointDirection;
231  }
232 }
BasicVector3D< T > unit() const
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
#define G4endl
Definition: G4ios.hh:61
BasicVector3D< T > cross(const BasicVector3D< T > &v) const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ViewParameters::SetViewpointDirection ( const G4Vector3D viewpointDirection)

Here is the caller graph for this function:

void G4ViewParameters::SetVisibleDensity ( G4double  visibleDensity)

Definition at line 173 of file G4ViewParameters.cc.

173  {
174  const G4double reasonableMaximum = 10.0 * g / cm3;
175  if (visibleDensity < 0) {
176  G4cout << "G4ViewParameters::SetVisibleDensity: attempt to set negative "
177  "density - ignored." << G4endl;
178  }
179  else {
180  if (visibleDensity > reasonableMaximum) {
181  G4cout << "G4ViewParameters::SetVisibleDensity: density > "
182  << G4BestUnit (reasonableMaximum, "Volumic Mass")
183  << " - did you mean this?"
184  << G4endl;
185  }
186  fVisibleDensity = visibleDensity;
187  }
188 }
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
G4GLOB_DLL std::ostream G4cout
static constexpr double cm3
Definition: G4SIunits.hh:121
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ViewParameters::SetWindowLocationHint ( G4int  xHint,
G4int  yHint 
)
void G4ViewParameters::SetWindowSizeHint ( G4int  xHint,
G4int  yHint 
)
void G4ViewParameters::SetXGeometryString ( const G4String geomStringArg)

Definition at line 918 of file G4ViewParameters.cc.

919 {
920  G4int x = 0, y = 0;
921  unsigned int w = 0, h = 0;
922  G4String geomString = geomStringArg;
923  // Parse windowSizeHintString for backwards compatibility...
924  const G4String delimiters("xX+-");
925  G4String::size_type i = geomString.find_first_of(delimiters);
926  if (i == G4String::npos) { // Does not contain "xX+-". Assume single number
927  std::istringstream iss(geomString);
928  G4int size;
929  iss >> size;
930  if (!iss) {
931  size = 600;
932  G4cout << "Unrecognised windowSizeHint string: \""
933  << geomString
934  << "\". Asuuming " << size << G4endl;
935  }
936  std::ostringstream oss;
937  oss << size << 'x' << size;
938  geomString = oss.str();
939  }
940 
941  fGeometryMask = ParseGeometry( geomString, &x, &y, &w, &h );
942 
943  // Handle special case :
944  if ((fGeometryMask & fYValue) == 0)
945  { // Using default
946  y = fWindowLocationHintY;
947  }
948  if ((fGeometryMask & fXValue) == 0)
949  { // Using default
950  x = fWindowLocationHintX;
951  }
952 
953  // Check errors
954  // if there is no Width and Height
955  if ( ((fGeometryMask & fHeightValue) == 0 ) &&
956  ((fGeometryMask & fWidthValue) == 0 )) {
957  h = fWindowSizeHintY;
958  w = fWindowSizeHintX;
959  } else if ((fGeometryMask & fHeightValue) == 0 ) {
960 
961  // if there is only Width. Special case to be backward compatible
962  // We set Width and Height the same to obtain a square windows.
963 
964  G4cout << "Unrecognised geometry string \""
965  << geomString
966  << "\". No Height found. Using Width value instead"
967  << G4endl;
968  h = w;
969  }
970  if ( ((fGeometryMask & fXValue) == 0 ) ||
971  ((fGeometryMask & fYValue) == 0 )) {
972  //Using defaults
973  x = fWindowLocationHintX;
974  y = fWindowLocationHintY;
975  }
976  // Set the string
977  fXGeometryString = geomString;
978 
979  // Set values
980  fWindowSizeHintX = w;
981  fWindowSizeHintY = h;
982  fWindowLocationHintX = x;
983  fWindowLocationHintY = y;
984 
985  if ( ((fGeometryMask & fXValue)) &&
986  ((fGeometryMask & fYValue))) {
987 
988  if ( (fGeometryMask & fXNegative) ) {
989  fWindowLocationHintXNegative = true;
990  } else {
991  fWindowLocationHintXNegative = false;
992  }
993  if ( (fGeometryMask & fYNegative) ) {
994  fWindowLocationHintYNegative = true;
995  } else {
996  fWindowLocationHintYNegative = false;
997  }
998  }
999 }
tuple x
Definition: test.py:50
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the caller graph for this function:

void G4ViewParameters::SetZoomFactor ( G4double  zoomFactor)

Here is the caller graph for this function:

G4String G4ViewParameters::TouchableCommands ( ) const

Definition at line 484 of file G4ViewParameters.cc.

485 {
486  std::ostringstream oss;
487 
488  oss
489  << "#\n# Touchable commands"
490  << "\n/vis/viewer/clearVisAttributesModifiers";
491 
492  const std::vector<G4ModelingParameters::VisAttributesModifier>& vams =
493  fVisAttributesModifiers;
494 
495  if (vams.empty()) {
496  oss << "\n# None";
497  oss << std::endl;
498  return oss.str();
499  }
500 
502  std::vector<G4ModelingParameters::VisAttributesModifier>::const_iterator
503  iModifier;
504  for (iModifier = vams.begin();
505  iModifier != vams.end();
506  ++iModifier) {
508  iModifier->GetPVNameCopyNoPath();
509  if (vamPath != lastPath) {
510  lastPath = vamPath;
511  oss << "\n/vis/set/touchable";
513  for (iVAM = vamPath.begin();
514  iVAM != vamPath.end();
515  ++iVAM) {
516  oss << ' ' << iVAM->GetName() << ' ' << iVAM->GetCopyNo();
517  }
518  }
519  const G4VisAttributes& vamVisAtts = iModifier->GetVisAttributes();
520  const G4Colour& c = vamVisAtts.GetColour();
521  switch (iModifier->GetVisAttributesSignifier()) {
523  oss << "\n/vis/touchable/set/visibility ";
524  if (vamVisAtts.IsVisible()) {
525  oss << "true";
526  } else {
527  oss << "false";
528  }
529  break;
531  oss << "\n/vis/touchable/set/daughtersInvisible ";
532  if (vamVisAtts.IsDaughtersInvisible()) {
533  oss << "true";
534  } else {
535  oss << "false";
536  }
537  break;
539  oss << "\n/vis/touchable/set/colour "
540  << c.GetRed()
541  << ' ' << c.GetGreen()
542  << ' ' << c.GetBlue()
543  << ' ' << c.GetAlpha();
544  break;
546  oss << "\n/vis/touchable/set/lineStyle ";
547  switch (vamVisAtts.GetLineStyle()) {
549  oss << "unbroken";
550  break;
552  oss << "dashed";
553  break;
555  oss << "dotted";
556  }
557  break;
559  oss << "\n/vis/touchable/set/lineWidth "
560  << vamVisAtts.GetLineWidth();
561  break;
563  if (vamVisAtts.IsForceDrawingStyle()) {
565  oss << "\n/vis/touchable/set/forceWireframe ";
566  if (vamVisAtts.IsForceDrawingStyle()) {
567  oss << "true";
568  } else {
569  oss << "false";
570  }
571  }
572  }
573  break;
575  if (vamVisAtts.IsForceDrawingStyle()) {
576  if (vamVisAtts.GetForcedDrawingStyle() == G4VisAttributes::solid) {
577  oss << "\n/vis/touchable/set/forceSolid ";
578  if (vamVisAtts.IsForceDrawingStyle()) {
579  oss << "true";
580  } else {
581  oss << "false";
582  }
583  }
584  }
585  break;
587  if (vamVisAtts.IsForceAuxEdgeVisible()) {
588  oss << "\n/vis/touchable/set/forceAuxEdgeVisible ";
589  if (vamVisAtts.IsForcedAuxEdgeVisible()) {
590  oss << "true";
591  } else {
592  oss << "false";
593  }
594  }
595  break;
597  if (vamVisAtts.GetForcedLineSegmentsPerCircle() > 0) {
598  oss << "\n/vis/touchable/set/lineSegmentsPerCircle "
599  << vamVisAtts.GetForcedLineSegmentsPerCircle();
600  }
601  break;
602  }
603  }
604 
605  oss << std::endl;
606 
607  return oss.str();
608 }
G4bool IsForceAuxEdgeVisible() const
G4double GetAlpha() const
Definition: G4Colour.hh:142
G4double GetLineWidth() const
G4bool IsVisible() const
const G4Colour & GetColour() const
G4double GetBlue() const
Definition: G4Colour.hh:141
LineStyle GetLineStyle() const
G4double GetRed() const
Definition: G4Colour.hh:139
G4bool IsDaughtersInvisible() const
G4double GetGreen() const
Definition: G4Colour.hh:140
G4bool IsForcedAuxEdgeVisible() const
G4int GetForcedLineSegmentsPerCircle() const
std::vector< PVNameCopyNo > PVNameCopyNoPath
G4bool IsForceDrawingStyle() const
tuple c
Definition: test.py:13
ForcedDrawingStyle GetForcedDrawingStyle() const
PVNameCopyNoPath::const_iterator PVNameCopyNoPathConstIterator

Here is the call graph for this function:

void G4ViewParameters::UnsetExplodeFactor ( )
void G4ViewParameters::UnsetSectionPlane ( )

Here is the caller graph for this function:

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  os,
const DrawingStyle style 
)
friend

Definition at line 681 of file G4ViewParameters.cc.

682  {
683  switch (style) {
685  os << "wireframe"; break;
687  os << "hlr - hidden lines removed"; break;
689  os << "hsr - hidden surfaces removed"; break;
691  os << "hlhsr - hidden line, hidden surface removed"; break;
692  default: os << "unrecognised"; break;
693  }
694  return os;
695 }
std::ostream& operator<< ( std::ostream &  os,
const G4ViewParameters v 
)
friend

Definition at line 697 of file G4ViewParameters.cc.

697  {
698  os << "View parameters and options:";
699 
700  os << "\n Drawing style: ";
701  switch (v.fDrawingStyle) {
703  os << "edges, wireframe"; break;
705  os << "edges, hidden line removal"; break;
707  os << "surfaces, hidden surface removal"; break;
709  os << "surfaces and edges, hidden line and surface removal"; break;
710  default: os << "unrecognised"; break;
711  }
712 
713  os << "\n Auxiliary edges: ";
714  if (!v.fAuxEdgeVisible) os << "in";
715  os << "visible";
716 
717  os << "\n Culling: ";
718  if (v.fCulling) os << "on";
719  else os << "off";
720 
721  os << "\n Culling invisible objects: ";
722  if (v.fCullInvisible) os << "on";
723  else os << "off";
724 
725  os << "\n Density culling: ";
726  if (v.fDensityCulling) {
727  os << "on - invisible if density less than "
728  << v.fVisibleDensity / (1. * g / cm3) << " g cm^-3";
729  }
730  else os << "off";
731 
732  os << "\n Culling daughters covered by opaque mothers: ";
733  if (v.fCullCovered) os << "on";
734  else os << "off";
735 
736  os << "\n Section flag: ";
737  if (v.fSection) os << "true, section/cut plane: " << v.fSectionPlane;
738  else os << "false";
739 
740  if (v.IsCutaway()) {
741  os << "\n Cutaway planes: ";
742  for (size_t i = 0; i < v.fCutawayPlanes.size (); i++) {
743  os << ' ' << v.fCutawayPlanes[i];
744  }
745  }
746  else {
747  os << "\n No cutaway planes";
748  }
749 
750  os << "\n Explode factor: " << v.fExplodeFactor
751  << " about centre: " << v.fExplodeCentre;
752 
753  os << "\n No. of sides used in circle polygon approximation: "
754  << v.fNoOfSides;
755 
756  os << "\n Viewpoint direction: " << v.fViewpointDirection;
757 
758  os << "\n Up vector: " << v.fUpVector;
759 
760  os << "\n Field half angle: " << v.fFieldHalfAngle;
761 
762  os << "\n Zoom factor: " << v.fZoomFactor;
763 
764  os << "\n Scale factor: " << v.fScaleFactor;
765 
766  os << "\n Current target point: " << v.fCurrentTargetPoint;
767 
768  os << "\n Dolly distance: " << v.fDolly;
769 
770  os << "\n Light ";
771  if (v.fLightsMoveWithCamera) os << "moves";
772  else os << "does not move";
773  os << " with camera";
774 
775  os << "\n Relative lightpoint direction: "
776  << v.fRelativeLightpointDirection;
777 
778  os << "\n Actual lightpoint direction: "
779  << v.fActualLightpointDirection;
780 
781  os << "\n Derived parameters for standard view of object of unit radius:";
782  G4ViewParameters tempVP = v;
783  tempVP.fDolly = 0.;
784  tempVP.fZoomFactor = 1.;
785  const G4double radius = 1.;
786  const G4double cameraDistance = tempVP.GetCameraDistance (radius);
787  const G4double nearDistance =
788  tempVP.GetNearDistance (cameraDistance, radius);
789  const G4double farDistance =
790  tempVP.GetFarDistance (cameraDistance, nearDistance, radius);
791  const G4double right = tempVP.GetFrontHalfHeight (nearDistance, radius);
792  os << "\n Camera distance: " << cameraDistance;
793  os << "\n Near distance: " << nearDistance;
794  os << "\n Far distance: " << farDistance;
795  os << "\n Front half height: " << right;
796 
797  os << "\n Default VisAttributes:\n " << v.fDefaultVisAttributes;
798 
799  os << "\n Default TextVisAttributes:\n " << v.fDefaultTextVisAttributes;
800 
801  os << "\n Default marker: " << v.fDefaultMarker;
802 
803  os << "\n Global marker scale: " << v.fGlobalMarkerScale;
804 
805  os << "\n Global lineWidth scale: " << v.fGlobalLineWidthScale;
806 
807  os << "\n Marker ";
808  if (v.fMarkerNotHidden) os << "not ";
809  os << "hidden by surfaces.";
810 
811  os << "\n Window size hint: "
812  << v.fWindowSizeHintX << 'x'<< v.fWindowSizeHintX;
813 
814  os << "\n X geometry string: " << v.fXGeometryString;
815  os << "\n X geometry mask: "
816  << std::showbase << std::hex << v.fGeometryMask
817  << std::noshowbase << std::dec;
818 
819  os << "\n Auto refresh: ";
820  if (v.fAutoRefresh) os << "true";
821  else os << "false";
822 
823  os << "\n Background colour: " << v.fBackgroundColour;
824 
825  os << "\n Picking requested: ";
826  if (v.fPicking) os << "true";
827  else os << "false";
828 
829  os << "\n Rotation style: ";
830  switch (v.fRotationStyle) {
832  os << "constrainUpDirection (conventional HEP view)"; break;
834  os << "freeRotation (Google-like rotation, using mouse-grab)"; break;
835  default: os << "unrecognised"; break;
836  }
837 
838  os << "\n Vis attributes modifiers: ";
839  const std::vector<G4ModelingParameters::VisAttributesModifier>& vams =
840  v.fVisAttributesModifiers;
841  if (vams.empty()) {
842  os << "None";
843  } else {
844  os << vams;
845  }
846 
847  return os;
848 }
G4double GetFarDistance(G4double cameraDistance, G4double nearDistance, G4double radius) const
G4double GetCameraDistance(G4double radius) const
G4double GetNearDistance(G4double cameraDistance, G4double radius) const
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
G4double GetFrontHalfHeight(G4double nearDistance, G4double radius) const
static constexpr double cm3
Definition: G4SIunits.hh:121
tuple v
Definition: test.py:18
G4bool IsCutaway() const
double G4double
Definition: G4Types.hh:76

The documentation for this class was generated from the following files: