Geant4-11
G4ViewParameters.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//
28//
29// John Allison 19th July 1996
30// View parameters and options.
31
32#include "G4ViewParameters.hh"
33
34#include "G4VisManager.hh"
35#include "G4VPhysicalVolume.hh"
36#include "G4UnitsTable.hh"
37#include "G4SystemOfUnits.hh"
38#include "G4Polyhedron.hh"
39
40#include <sstream>
41#include <cmath>
42
44 fDrawingStyle (wireframe),
45 fNumberOfCloudPoints(10000),
46 fAuxEdgeVisible (false),
47 fCulling (true),
48 fCullInvisible (true),
49 fDensityCulling (false),
50 fVisibleDensity (0.01 * g / cm3),
51 fCullCovered (false),
52 fCBDAlgorithmNumber (0),
53 fSection (false),
54 fSectionPlane (),
55 fCutawayMode (cutawayUnion),
56 fCutawayPlanes (),
57 fExplodeFactor (1.),
58 fNoOfSides (),
59 fViewpointDirection (G4Vector3D (0., 0., 1.)), // On z-axis.
60 fUpVector (G4Vector3D (0., 1., 0.)), // y-axis up.
61 fFieldHalfAngle (0.), // Orthogonal projection.
62 fZoomFactor (1.),
63 fScaleFactor (G4Vector3D (1., 1., 1.)),
64 fCurrentTargetPoint (),
65 fDolly (0.),
66 fLightsMoveWithCamera (false),
67 fRelativeLightpointDirection (G4Vector3D (1., 1., 1.)),
68 fActualLightpointDirection (G4Vector3D (1., 1., 1.)),
69 fDefaultVisAttributes (),
70 fDefaultTextVisAttributes (G4Colour (0., 0., 1.)),
71 fDefaultMarker (),
72 fGlobalMarkerScale (1.),
73 fGlobalLineWidthScale (1.),
74 fMarkerNotHidden (true),
75 fWindowSizeHintX (600),
76 fWindowSizeHintY (600),
77 fWindowLocationHintX(0),
78 fWindowLocationHintY(0),
79 fWindowLocationHintXNegative(true),
80 fWindowLocationHintYNegative(false),
81 fGeometryMask(0),
82 fAutoRefresh (false),
83 fBackgroundColour (G4Colour(0.,0.,0.)), // Black
84 fPicking (false),
85 fRotationStyle (constrainUpDirection),
86 fStartTime(-G4VisAttributes::fVeryLongTime),
87 fEndTime(G4VisAttributes::fVeryLongTime),
88 fFadeFactor(0.),
89 fDisplayHeadTime(false),
90 fDisplayHeadTimeX(-0.9),
91 fDisplayHeadTimeY(-0.9),
92 fDisplayHeadTimeSize(24.),
93 fDisplayHeadTimeRed(0.),
94 fDisplayHeadTimeGreen(1.),
95 fDisplayHeadTimeBlue(1.),
96 fDisplayLightFront(false),
97 fDisplayLightFrontX(0.),
98 fDisplayLightFrontY(0.),
99 fDisplayLightFrontZ(0.),
100 fDisplayLightFrontT(0.),
101 fDisplayLightFrontRed(0.),
102 fDisplayLightFrontGreen(1.),
103 fDisplayLightFrontBlue(0.),
104 fSpecialMeshRendering(false)
105{
106 // Pick up default no of sides from G4Polyhedron.
107 // Note that this parameter is variously called:
108 // No of sides
109 // NumberOfRotationSteps
110 // Line segments per circle
111 // It refers to the approximation of a circle by a polygon of
112 // stated number of sides.
114
116 // Markers are 5 pixels "overall" size, i.e., diameter.
117}
118
120
122(const G4Vector3D& scaleFactorMultiplier) {
123 fScaleFactor.setX(fScaleFactor.x() * scaleFactorMultiplier.x());
124 fScaleFactor.setY(fScaleFactor.y() * scaleFactorMultiplier.y());
125 fScaleFactor.setZ(fScaleFactor.z() * scaleFactorMultiplier.z());
126}
127
131}
132
133// Useful quantities - begin snippet.
134// Here Follow functions to evaluate the above algorithms as a
135// function of the radius of the Bounding Sphere of the object being
136// viewed. Call them in the order given - for efficiency, later
137// functions depend on the results of earlier ones (Store the
138// results of earlier functions in your own temporary variables -
139// see, for example, G4OpenGLView::SetView ().)
140
142 G4double cameraDistance;
143 if (fFieldHalfAngle == 0.) {
144 cameraDistance = radius;
145 }
146 else {
147 cameraDistance = radius / std::sin (fFieldHalfAngle) - fDolly;
148 }
149 return cameraDistance;
150}
151
153 G4double radius) const {
154 const G4double small = 1.e-6 * radius;
155 G4double nearDistance = cameraDistance - radius;
156 if (nearDistance < small) nearDistance = small;
157 return nearDistance;
158}
159
161 G4double nearDistance,
162 G4double radius) const {
163 G4double farDistance = cameraDistance + radius;
164 if (farDistance < nearDistance) farDistance = nearDistance;
165 return farDistance;
166}
167
169 G4double radius) const {
170 G4double frontHalfHeight;
171 if (fFieldHalfAngle == 0.) {
172 frontHalfHeight = radius / fZoomFactor;
173 }
174 else {
175 frontHalfHeight = nearDistance * std::tan (fFieldHalfAngle) / fZoomFactor;
176 }
177 return frontHalfHeight;
178}
179// Useful quantities - end snippet.
180
182 if (fCutawayPlanes.size () < 3 ) {
183 fCutawayPlanes.push_back (cutawayPlane);
184 }
185 else {
186 G4cerr <<
187 "ERROR: G4ViewParameters::AddCutawayPlane:"
188 "\n A maximum of 3 cutaway planes supported." << G4endl;
189 }
190}
191
193(size_t index, const G4Plane3D& cutawayPlane) {
194 if (index >= fCutawayPlanes.size()) {
195 G4cerr <<
196 "ERROR: G4ViewParameters::ChangeCutawayPlane:"
197 "\n Plane " << index << " does not exist." << G4endl;
198 } else {
199 fCutawayPlanes[index] = cutawayPlane;
200 }
201}
202
204 const G4double reasonableMaximum = 10.0 * g / cm3;
205 if (visibleDensity < 0) {
206 G4cout << "G4ViewParameters::SetVisibleDensity: attempt to set negative "
207 "density - ignored." << G4endl;
208 }
209 else {
210 if (visibleDensity > reasonableMaximum) {
211 G4cout << "G4ViewParameters::SetVisibleDensity: density > "
212 << G4BestUnit (reasonableMaximum, "Volumic Mass")
213 << " - did you mean this?"
214 << G4endl;
215 }
216 fVisibleDensity = visibleDensity;
217 }
218}
219
222 if (nSides < nSidesMin) {
223 nSides = nSidesMin;
224 G4cout << "G4ViewParameters::SetNoOfSides: attempt to set the"
225 "\nnumber of sides per circle < " << nSidesMin
226 << "; forced to " << nSides << G4endl;
227 }
228 fNoOfSides = nSides;
229 return fNoOfSides;
230}
231
233 const G4int nPointsMin = 100;
234 if (nPoints < nPointsMin) {
235 nPoints = nPointsMin;
236 G4cout << "G4ViewParameters::SetNumberOfCloudPoints:"
237 "\nnumber of points per cloud set to minimum " << nPoints
238 << G4endl;
239 }
240 fNumberOfCloudPoints = nPoints;
242}
243
245(const G4Vector3D& viewpointDirection) {
246
247 fViewpointDirection = viewpointDirection;
248
249 // If the requested viewpoint direction is parallel to the up
250 // vector, the orientation of the view is undefined...
251 if (fViewpointDirection.unit() * fUpVector.unit() > .9999) {
252 static G4bool firstTime = true;
253 if (firstTime) {
254 firstTime = false;
255 G4cout <<
256 "WARNING: Viewpoint direction is very close to the up vector direction."
257 "\n Change the up vector or \"/vis/viewer/set/rotationStyle freeRotation\"."
258 << G4endl;
259 }
260 }
261
262 // Move the lights too if requested...
265 G4Vector3D xprime = (fUpVector.cross (zprime)).unit ();
266 G4Vector3D yprime = zprime.cross (xprime);
268 fRelativeLightpointDirection.x () * xprime +
269 fRelativeLightpointDirection.y () * yprime +
270 fRelativeLightpointDirection.x () * zprime;
271 } else {
273 }
274}
275
277(const G4Vector3D& lightpointDirection) {
278 fRelativeLightpointDirection = lightpointDirection;
280}
281
283 G4Vector3D unitRight = (fUpVector.cross (fViewpointDirection)).unit();
284 G4Vector3D unitUp = (fViewpointDirection.cross (unitRight)).unit();
285 fCurrentTargetPoint = right * unitRight + up * unitUp;
286}
287
289 IncrementPan (right,up, 0);
290}
291
293 G4Vector3D unitRight = (fUpVector.cross (fViewpointDirection)).unit();
294 G4Vector3D unitUp = (fViewpointDirection.cross (unitRight)).unit();
295 fCurrentTargetPoint += right * unitRight + up * unitUp + distance * fViewpointDirection;
296}
297
300 // If target exists with same signifier just change vis attributes.
301 G4bool duplicateTarget = false;
302 auto i = fVisAttributesModifiers.begin();
303 for (; i < fVisAttributesModifiers.end(); ++i) {
304 if (vam.GetPVNameCopyNoPath() == (*i).GetPVNameCopyNoPath() &&
305 vam.GetVisAttributesSignifier() == (*i).GetVisAttributesSignifier()) {
306 duplicateTarget = true;
307 break;
308 }
309 }
310 if (duplicateTarget) (*i).SetVisAttributes(vam.GetVisAttributes());
311 else fVisAttributesModifiers.push_back(vam);
312}
313
315(const G4Point3D standardTargetPoint) const
316{
317 std::ostringstream oss;
318
319 oss << "#\n# Camera and lights commands";
320
321 oss << "\n/vis/viewer/set/viewpointVector "
323 << ' ' << fViewpointDirection.y()
324 << ' ' << fViewpointDirection.z();
325
326 oss << "\n/vis/viewer/set/upVector "
327 << fUpVector.x()
328 << ' ' << fUpVector.y()
329 << ' ' << fUpVector.z();
330
331 oss << "\n/vis/viewer/set/projection ";
332 if (fFieldHalfAngle == 0.) {
333 oss
334 << "orthogonal";
335 } else {
336 oss
337 << "perspective "
339 << " deg";
340 }
341
342 oss << "\n/vis/viewer/zoomTo "
343 << fZoomFactor;
344
345 oss << "\n/vis/viewer/scaleTo "
346 << fScaleFactor.x()
347 << ' ' << fScaleFactor.y()
348 << ' ' << fScaleFactor.z();
349
350 oss << "\n/vis/viewer/set/targetPoint "
351 << G4BestUnit(standardTargetPoint+fCurrentTargetPoint,"Length")
352 << "\n# Note that if you have not set a target point, the vis system sets"
353 << "\n# a target point based on the scene - plus any panning and dollying -"
354 << "\n# so don't be alarmed by strange coordinates here.";
355
356 oss << "\n/vis/viewer/dollyTo "
357 << G4BestUnit(fDolly,"Length");
358
359 oss << "\n/vis/viewer/set/lightsMove ";
361 oss << "camera";
362 } else {
363 oss << "object";
364 }
365
366 oss << "\n/vis/viewer/set/lightsVector "
370
371 oss << "\n/vis/viewer/set/rotationStyle ";
373 oss << "constrainUpDirection";
374 } else {
375 oss << "freeRotation";
376 }
377
379 oss << "\n/vis/viewer/set/background "
380 << c.GetRed()
381 << ' ' << c.GetGreen()
382 << ' ' << c.GetBlue()
383 << ' ' << c.GetAlpha();
384
386 oss << "\n/vis/viewer/set/defaultColour "
387 << c.GetRed()
388 << ' ' << c.GetGreen()
389 << ' ' << c.GetBlue()
390 << ' ' << c.GetAlpha();
391
393 oss << "\n/vis/viewer/set/defaultTextColour "
394 << c.GetRed()
395 << ' ' << c.GetGreen()
396 << ' ' << c.GetBlue()
397 << ' ' << c.GetAlpha();
398
399 oss << std::endl;
400
401 return oss.str();
402}
403
405{
406 std::ostringstream oss;
407
408 oss << "#\n# Drawing style commands";
409
410 oss << "\n/vis/viewer/set/style ";
411 switch (fDrawingStyle) {
412 case wireframe:
413 case hlr:
414 oss << "wireframe";
415 break;
416 case hsr:
417 case hlhsr:
418 oss << "surface";
419 break;
420 case cloud:
421 oss << "cloud";
422 break;
423 }
424
425 oss << "\n/vis/viewer/set/hiddenEdge ";
426 if (fDrawingStyle == hlr || fDrawingStyle == hlhsr) {
427 oss << "true";
428 } else {
429 oss << "false";
430 }
431
432 oss << "\n/vis/viewer/set/auxiliaryEdge ";
433 if (fAuxEdgeVisible) {
434 oss << "true";
435 } else {
436 oss << "false";
437 }
438
439 oss << "\n/vis/viewer/set/hiddenMarker ";
440 if (fMarkerNotHidden) {
441 oss << "false";
442 } else {
443 oss << "true";
444 }
445
446 oss << "\n/vis/viewer/set/globalLineWidthScale "
448
449 oss << "\n/vis/viewer/set/globalMarkerScale "
451
452 oss << "\n/vis/viewer/set/numberOfCloudPoints "
454
455 oss << "\n/vis/viewer/set/specialMeshRendering ";
457 oss << "true";
458 } else {
459 oss << "false";
460 }
461
462 oss << "\n/vis/viewer/set/specialMeshVolumes";
463 for (const auto& volume : fSpecialMeshVolumes) {
464 oss << ' ' << volume.GetName() << ' ' << volume.GetCopyNo();
465 }
466
467 oss << std::endl;
468
469 return oss.str();
470}
471
473{
474 std::ostringstream oss;
475
476 oss << "#\n# Scene-modifying commands";
477
478 oss << "\n/vis/viewer/set/culling global ";
479 if (fCulling) {
480 oss << "true";
481 } else {
482 oss << "false";
483 }
484
485 oss << "\n/vis/viewer/set/culling invisible ";
486 if (fCullInvisible) {
487 oss << "true";
488 } else {
489 oss << "false";
490 }
491
492 oss << "\n/vis/viewer/set/culling density ";
493 if (fDensityCulling) {
494 oss << "true " << fVisibleDensity/(g/cm3) << " g/cm3";
495 } else {
496 oss << "false";
497 }
498
499 oss << "\n/vis/viewer/set/culling coveredDaughters ";
500 if (fCullCovered) {
501 oss << "true";
502 } else {
503 oss << "false";
504 }
505
506 oss << "\n/vis/viewer/colourByDensity "
507 << fCBDAlgorithmNumber << " g/cm3";
508 for (auto p: fCBDParameters) {
509 oss << ' ' << p/(g/cm3);
510 }
511
512 oss << "\n/vis/viewer/set/sectionPlane ";
513 if (fSection) {
514 oss << "on "
515 << G4BestUnit(fSectionPlane.point(),"Length")
516 << fSectionPlane.normal().x()
517 << ' ' << fSectionPlane.normal().y()
518 << ' ' << fSectionPlane.normal().z();
519 } else {
520 oss << "off";
521 }
522
523 oss << "\n/vis/viewer/set/cutawayMode ";
524 if (fCutawayMode == cutawayUnion) {
525 oss << "union";
526 } else {
527 oss << "intersection";
528 }
529
530 oss << "\n/vis/viewer/clearCutawayPlanes";
531 if (fCutawayPlanes.size()) {
532 for (size_t i = 0; i < fCutawayPlanes.size(); i++) {
533 oss << "\n/vis/viewer/addCutawayPlane "
534 << G4BestUnit(fCutawayPlanes[i].point(),"Length")
535 << fCutawayPlanes[i].normal().x()
536 << ' ' << fCutawayPlanes[i].normal().y()
537 << ' ' << fCutawayPlanes[i].normal().z();
538 }
539 } else {
540 oss << "\n# No cutaway planes defined.";
541 }
542
543 oss << "\n/vis/viewer/set/explodeFactor "
545 << ' ' << G4BestUnit(fExplodeCentre,"Length");
546
547 oss << "\n/vis/viewer/set/lineSegmentsPerCircle "
548 << fNoOfSides;
549
550 oss << std::endl;
551
552 return oss.str();
553}
554
556{
557 std::ostringstream oss;
558
559 oss << "#\n# Touchable commands";
560
561 const std::vector<G4ModelingParameters::VisAttributesModifier>& vams =
563
564 if (vams.empty()) {
565 oss
566 << "\n# None"
567 << "\n/vis/viewer/clearVisAttributesModifiers";
568 oss << std::endl;
569 return oss.str();
570 }
571
572 oss
573 << "\n/vis/viewer/clearVisAttributesModifiers";
574
576 std::vector<G4ModelingParameters::VisAttributesModifier>::const_iterator
577 iModifier;
578 for (iModifier = vams.begin();
579 iModifier != vams.end();
580 ++iModifier) {
582 iModifier->GetPVNameCopyNoPath();
583 if (vamPath != lastPath) {
584 lastPath = vamPath;
585 oss << "\n/vis/set/touchable";
587 for (iVAM = vamPath.begin();
588 iVAM != vamPath.end();
589 ++iVAM) {
590 oss << ' ' << iVAM->GetName() << ' ' << iVAM->GetCopyNo();
591 }
592 }
593 const G4VisAttributes& vamVisAtts = iModifier->GetVisAttributes();
594 const G4Colour& c = vamVisAtts.GetColour();
595 switch (iModifier->GetVisAttributesSignifier()) {
597 oss << "\n/vis/touchable/set/visibility ";
598 if (vamVisAtts.IsVisible()) {
599 oss << "true";
600 } else {
601 oss << "false";
602 }
603 break;
605 oss << "\n/vis/touchable/set/daughtersInvisible ";
606 if (vamVisAtts.IsDaughtersInvisible()) {
607 oss << "true";
608 } else {
609 oss << "false";
610 }
611 break;
613 oss << "\n/vis/touchable/set/colour "
614 << c.GetRed()
615 << ' ' << c.GetGreen()
616 << ' ' << c.GetBlue()
617 << ' ' << c.GetAlpha();
618 break;
620 oss << "\n/vis/touchable/set/lineStyle ";
621 switch (vamVisAtts.GetLineStyle()) {
623 oss << "unbroken";
624 break;
626 oss << "dashed";
627 break;
629 oss << "dotted";
630 }
631 break;
633 oss << "\n/vis/touchable/set/lineWidth "
634 << vamVisAtts.GetLineWidth();
635 break;
637 if (vamVisAtts.IsForceDrawingStyle()) {
639 oss << "\n/vis/touchable/set/forceWireframe ";
640 if (vamVisAtts.IsForceDrawingStyle()) {
641 oss << "true";
642 } else {
643 oss << "false";
644 }
645 }
646 }
647 break;
649 if (vamVisAtts.IsForceDrawingStyle()) {
650 if (vamVisAtts.GetForcedDrawingStyle() == G4VisAttributes::solid) {
651 oss << "\n/vis/touchable/set/forceSolid ";
652 if (vamVisAtts.IsForceDrawingStyle()) {
653 oss << "true";
654 } else {
655 oss << "false";
656 }
657 }
658 }
659 break;
661 if (vamVisAtts.IsForceDrawingStyle()) {
662 if (vamVisAtts.GetForcedDrawingStyle() == G4VisAttributes::cloud) {
663 oss << "\n/vis/touchable/set/forceCloud ";
664 if (vamVisAtts.IsForceDrawingStyle()) {
665 oss << "true";
666 } else {
667 oss << "false";
668 }
669 }
670 }
671 break;
673 if (vamVisAtts.IsForceAuxEdgeVisible()) {
674 oss << "\n/vis/touchable/set/forceAuxEdgeVisible ";
675 if (vamVisAtts.IsForcedAuxEdgeVisible()) {
676 oss << "true";
677 } else {
678 oss << "false";
679 }
680 }
681 break;
683 oss << "\n/vis/touchable/set/lineSegmentsPerCircle "
684 << vamVisAtts.GetForcedLineSegmentsPerCircle();
685 break;
687 oss << "\n/vis/touchable/set/numberOfCloudPoints "
688 << vamVisAtts.GetForcedNumberOfCloudPoints();
689 break;
690 }
691 }
692
693 oss << std::endl;
694
695 return oss.str();
696}
697
699{
700 std::ostringstream oss;
701
702 oss << "#\n# Time window commands";
703
704 oss
705 << "\n/vis/viewer/set/timeWindow/startTime "
706 << fStartTime/ns << " ns ";
707
708 oss
709 << "\n/vis/viewer/set/timeWindow/endTime "
710 << fEndTime/ns << " ns ";
711
712 oss << "\n/vis/viewer/set/timeWindow/fadeFactor "
713 << fFadeFactor;
714
715 oss
716 << "\n/vis/viewer/set/timeWindow/displayHeadTime ";
717 if (!fDisplayHeadTime) {
718 oss << "false";
719 } else {
720 oss
721 << "true"
722 << ' ' << fDisplayHeadTimeX
723 << ' ' << fDisplayHeadTimeY
724 << ' ' << fDisplayHeadTimeSize
725 << ' ' << fDisplayHeadTimeRed
726 << ' ' << fDisplayHeadTimeGreen
727 << ' ' << fDisplayHeadTimeBlue;
728 }
729
730 oss
731 << "\n/vis/viewer/set/timeWindow/displayLightFront ";
732 if (!fDisplayLightFront) {
733 oss << "false";
734 } else {
735 oss
736 << "true"
737 << ' ' << fDisplayLightFrontX/mm
738 << ' ' << fDisplayLightFrontY/mm
739 << ' ' << fDisplayLightFrontZ/mm
740 << " mm"
741 << ' ' << fDisplayLightFrontT/ns
742 << " ns"
743 << ' ' << fDisplayLightFrontRed
745 << ' ' << fDisplayLightFrontBlue;
746 }
747
748 oss << std::endl;
749
750 return oss.str();
751}
752
754
755 // Put performance-sensitive parameters first.
756 if (
757 // This first to optimise spin, etc.
759
760 // No particular order from here on.
764 (fCulling != v.fCulling) ||
768 (fCullCovered != v.fCullCovered) ||
770 (fSection != v.fSection) ||
771 (fNoOfSides != v.fNoOfSides) ||
772 (fUpVector != v.fUpVector) ||
774 (fZoomFactor != v.fZoomFactor) ||
775 (fScaleFactor != v.fScaleFactor) ||
777 (fDolly != v.fDolly) ||
790 (fAutoRefresh != v.fAutoRefresh) ||
792 (fPicking != v.fPicking) ||
794 )
795 G4cout << "Difference in 1st batch." << G4endl;
796
797 if (fCBDAlgorithmNumber > 0) {
798 if (fCBDParameters.size() != v.fCBDParameters.size()) {
799 G4cout << "Difference in number of colour by density parameters." << G4endl;
800 } else if (fCBDParameters != v.fCBDParameters) {
801 G4cout << "Difference in values of colour by density parameters." << G4endl;
802 }
803 }
804
805 if (fSection) {
806 if (!(fSectionPlane == v.fSectionPlane))
807 G4cout << "Difference in section planes batch." << G4endl;
808 }
809
810 if (IsCutaway()) {
811 if (fCutawayPlanes.size () != v.fCutawayPlanes.size ()) {
812 G4cout << "Difference in no of cutaway planes." << G4endl;
813 }
814 else {
815 for (size_t i = 0; i < fCutawayPlanes.size (); i++) {
816 if (!(fCutawayPlanes[i] == v.fCutawayPlanes[i]))
817 G4cout << "Difference in cutaway plane no. " << i << G4endl;
818 }
819 }
820 }
821
822 if (IsExplode()) {
824 G4cout << "Difference in explode factor." << G4endl;
826 G4cout << "Difference in explode centre." << G4endl;
827 }
828
830 G4cout << "Difference in vis attributes modifiers." << G4endl;
831 }
832
833 if (fStartTime != v.fStartTime ||
834 fEndTime != v.fEndTime) {
835 G4cout << "Difference in time window." << G4endl;
836 }
837
838 if (fFadeFactor != v.fFadeFactor) {
839 G4cout << "Difference in time window fade factor." << G4endl;
840 }
841
843 G4cout << "Difference in display head time flag." << G4endl;
844 } else {
851 G4cout << "Difference in display head time parameters." << G4endl;
852 }
853 }
854
856 G4cout << "Difference in display light front flag." << G4endl;
857 } else {
865 G4cout << "Difference in display light front parameters." << G4endl;
866 }
867 }
868}
869
870std::ostream& operator <<
871(std::ostream& os, const G4ViewParameters::DrawingStyle& style)
872{
873 switch (style) {
875 os << "wireframe"; break;
877 os << "hlr - hidden lines removed"; break;
879 os << "hsr - hidden surfaces removed"; break;
881 os << "hlhsr - hidden line, hidden surface removed"; break;
883 os << "cloud - draw volume as a cloud of dots"; break;
884 default: os << "unrecognised"; break;
885 }
886 return os;
887}
888
889std::ostream& operator << (std::ostream& os, const G4ViewParameters& v) {
890 os << "View parameters and options:";
891
892 os << "\n Drawing style: " << v.fDrawingStyle;
893
894 os << "\n Number of cloud points: " << v.fNumberOfCloudPoints;
895
896 os << "\n Auxiliary edges: ";
897 if (!v.fAuxEdgeVisible) os << "in";
898 os << "visible";
899
900 os << "\n Culling: ";
901 if (v.fCulling) os << "on";
902 else os << "off";
903
904 os << "\n Culling invisible objects: ";
905 if (v.fCullInvisible) os << "on";
906 else os << "off";
907
908 os << "\n Density culling: ";
909 if (v.fDensityCulling) {
910 os << "on - invisible if density less than "
911 << v.fVisibleDensity / (1. * g / cm3) << " g cm^-3";
912 }
913 else os << "off";
914
915 os << "\n Culling daughters covered by opaque mothers: ";
916 if (v.fCullCovered) os << "on";
917 else os << "off";
918
919 os << "\n Colour by density: ";
920 if (v.fCBDAlgorithmNumber <= 0) {
921 os << "inactive";
922 } else {
923 os << "Algorithm " << v.fCBDAlgorithmNumber << ", Parameters:";
924 for (auto p: v.fCBDParameters) {
925 os << ' ' << G4BestUnit(p,"Volumic Mass");
926 }
927 }
928
929 os << "\n Section flag: ";
930 if (v.fSection) os << "true, section/cut plane: " << v.fSectionPlane;
931 else os << "false";
932
933 if (v.IsCutaway()) {
934 os << "\n Cutaway planes: ";
935 for (size_t i = 0; i < v.fCutawayPlanes.size (); i++) {
936 os << ' ' << v.fCutawayPlanes[i];
937 }
938 }
939 else {
940 os << "\n No cutaway planes";
941 }
942
943 os << "\n Explode factor: " << v.fExplodeFactor
944 << " about centre: " << v.fExplodeCentre;
945
946 os << "\n No. of sides used in circle polygon approximation: "
947 << v.fNoOfSides;
948
949 os << "\n Viewpoint direction: " << v.fViewpointDirection;
950
951 os << "\n Up vector: " << v.fUpVector;
952
953 os << "\n Field half angle: " << v.fFieldHalfAngle;
954
955 os << "\n Zoom factor: " << v.fZoomFactor;
956
957 os << "\n Scale factor: " << v.fScaleFactor;
958
959 os << "\n Current target point: " << v.fCurrentTargetPoint;
960
961 os << "\n Dolly distance: " << v.fDolly;
962
963 os << "\n Light ";
964 if (v.fLightsMoveWithCamera) os << "moves";
965 else os << "does not move";
966 os << " with camera";
967
968 os << "\n Relative lightpoint direction: "
970
971 os << "\n Actual lightpoint direction: "
973
974 os << "\n Derived parameters for standard view of object of unit radius:";
975 G4ViewParameters tempVP = v;
976 tempVP.fDolly = 0.;
977 tempVP.fZoomFactor = 1.;
978 const G4double radius = 1.;
979 const G4double cameraDistance = tempVP.GetCameraDistance (radius);
980 const G4double nearDistance =
981 tempVP.GetNearDistance (cameraDistance, radius);
982 const G4double farDistance =
983 tempVP.GetFarDistance (cameraDistance, nearDistance, radius);
984 const G4double right = tempVP.GetFrontHalfHeight (nearDistance, radius);
985 os << "\n Camera distance: " << cameraDistance;
986 os << "\n Near distance: " << nearDistance;
987 os << "\n Far distance: " << farDistance;
988 os << "\n Front half height: " << right;
989
990 os << "\n Default VisAttributes:\n " << v.fDefaultVisAttributes;
991
992 os << "\n Default TextVisAttributes:\n " << v.fDefaultTextVisAttributes;
993
994 os << "\n Default marker: " << v.fDefaultMarker;
995
996 os << "\n Global marker scale: " << v.fGlobalMarkerScale;
997
998 os << "\n Global lineWidth scale: " << v.fGlobalLineWidthScale;
999
1000 os << "\n Marker ";
1001 if (v.fMarkerNotHidden) os << "not ";
1002 os << "hidden by surfaces.";
1003
1004 os << "\n Window size hint: "
1005 << v.fWindowSizeHintX << 'x'<< v.fWindowSizeHintX;
1006
1007 os << "\n X geometry string: " << v.fXGeometryString;
1008 os << "\n X geometry mask: "
1009 << std::showbase << std::hex << v.fGeometryMask
1010 << std::noshowbase << std::dec;
1011
1012 os << "\n Auto refresh: ";
1013 if (v.fAutoRefresh) os << "true";
1014 else os << "false";
1015
1016 os << "\n Background colour: " << v.fBackgroundColour;
1017
1018 os << "\n Picking requested: ";
1019 if (v.fPicking) os << "true";
1020 else os << "false";
1021
1022 os << "\n Rotation style: ";
1023 switch (v.fRotationStyle) {
1025 os << "constrainUpDirection (conventional HEP view)"; break;
1027 os << "freeRotation (Google-like rotation, using mouse-grab)"; break;
1028 default: os << "unrecognised"; break;
1029 }
1030
1031 os << "\n Vis attributes modifiers: ";
1032 const std::vector<G4ModelingParameters::VisAttributesModifier>& vams =
1034 if (vams.empty()) {
1035 os << "None";
1036 } else {
1037 os << vams;
1038 }
1039
1040 os << "\n Time window parameters:"
1041 << "\n Start time: " << v.fStartTime/ns << " ns"
1042 << "\n End time: " << v.fEndTime/ns << " ns"
1043 << "\n Fade factor: " << v.fFadeFactor;
1044 if (!v.fDisplayHeadTime) {
1045 os << "\n Head time display not requested.";
1046 } else {
1047 os
1048 << "\n Head time position: "
1049 << v.fDisplayHeadTimeX << ' ' << v.fDisplayHeadTimeY
1050 << "\n Head time size: " << v.fDisplayHeadTimeSize
1051 << "\n Head time colour: " << v.fDisplayHeadTimeRed
1052 << ' ' << v.fDisplayHeadTimeGreen << ' ' << v.fDisplayHeadTimeBlue;
1053 }
1054 if (!v.fDisplayLightFront) {
1055 os << "\n Light front display not requested.";
1056 } else {
1057 os
1058 << "\n Light front position: "
1060 << ' ' << v.fDisplayLightFrontZ/mm << " mm"
1061 << "\n Light front time: " << v.fDisplayLightFrontT/ns << " ns"
1062 << "\n Light front colour: " << v.fDisplayLightFrontRed
1063 << ' ' << v.fDisplayLightFrontGreen << ' ' << v.fDisplayLightFrontBlue;
1064 }
1065
1066 os << "\n Special Mesh Rendering: ";
1067 if (v.fSpecialMeshRendering) {
1068 os << "on: ";
1069 if (v.fSpecialMeshVolumes.empty()) {
1070 os << "all meshes";
1071 } else {
1072 os << "selected meshes";
1073 for (const auto& vol: v.fSpecialMeshVolumes) {
1074 os << "\n " << vol.GetName() << ':' << vol.GetCopyNo();
1075 }
1076 }
1077 } else os << "off";
1078
1079 return os;
1080}
1081
1083
1084 // Put performance-sensitive parameters first.
1085 if (
1086 // This first to optimise spin, etc.
1088
1089 // No particular order from here on.
1093 (fCulling != v.fCulling) ||
1096 (fCullCovered != v.fCullCovered) ||
1098 (fSection != v.fSection) ||
1099 (IsCutaway() != v.IsCutaway()) ||
1100 (IsExplode() != v.IsExplode()) ||
1101 (fNoOfSides != v.fNoOfSides) ||
1102 (fUpVector != v.fUpVector) ||
1104 (fZoomFactor != v.fZoomFactor) ||
1105 (fScaleFactor != v.fScaleFactor) ||
1107 (fDolly != v.fDolly) ||
1120 (fAutoRefresh != v.fAutoRefresh) ||
1122 (fPicking != v.fPicking) ||
1125 )
1126 return true;
1127
1128 if (fDensityCulling &&
1129 (fVisibleDensity != v.fVisibleDensity)) return true;
1130
1131 if (fCBDAlgorithmNumber > 0) {
1132 if (fCBDParameters.size() != v.fCBDParameters.size()) return true;
1133 else if (fCBDParameters != v.fCBDParameters) return true;
1134 }
1135
1136 if (fSection &&
1137 (!(fSectionPlane == v.fSectionPlane))) return true;
1138
1139 if (IsCutaway()) {
1140 if (fCutawayPlanes.size () != v.fCutawayPlanes.size ())
1141 return true;
1142 else {
1143 for (size_t i = 0; i < fCutawayPlanes.size (); i++) {
1144 if (!(fCutawayPlanes[i] == v.fCutawayPlanes[i])) return true;
1145 }
1146 }
1147 }
1148
1149 if (IsExplode() &&
1151 (fExplodeCentre != v.fExplodeCentre))) return true;
1152
1154
1155 if (fStartTime != v.fStartTime ||
1156 fEndTime != v.fEndTime ||
1157 fFadeFactor != v.fFadeFactor) return true;
1158
1159 if (fDisplayHeadTime != v.fDisplayHeadTime) return true;
1160 if (fDisplayHeadTime) {
1167 return true;
1168 }
1169 }
1170
1171 if (fDisplayLightFront != v.fDisplayLightFront) return true;
1172 if (fDisplayLightFront) {
1180 return true;
1181 }
1182 }
1183
1186 return true;;
1187 }
1188
1189 return false;
1190}
1191
1193{
1194 G4int x = 0, y = 0;
1195 unsigned int w = 0, h = 0;
1196 G4String geomString = geomStringArg;
1197 // Parse windowSizeHintString for backwards compatibility...
1198 const G4String delimiters("xX+-");
1199 G4String::size_type i = geomString.find_first_of(delimiters);
1200 if (i == G4String::npos) { // Does not contain "xX+-". Assume single number
1201 std::istringstream iss(geomString);
1202 G4int size;
1203 iss >> size;
1204 if (!iss) {
1205 size = 600;
1206 G4cout << "Unrecognised windowSizeHint string: \""
1207 << geomString
1208 << "\". Asuuming " << size << G4endl;
1209 }
1210 std::ostringstream oss;
1211 oss << size << 'x' << size;
1212 geomString = oss.str();
1213 }
1214
1215 fGeometryMask = ParseGeometry( geomString, &x, &y, &w, &h );
1216
1217 // Handle special case :
1218 if ((fGeometryMask & fYValue) == 0)
1219 { // Using default
1221 }
1222 if ((fGeometryMask & fXValue) == 0)
1223 { // Using default
1225 }
1226
1227 // Check errors
1228 // if there is no Width and Height
1229 if ( ((fGeometryMask & fHeightValue) == 0 ) &&
1230 ((fGeometryMask & fWidthValue) == 0 )) {
1231 h = fWindowSizeHintY;
1232 w = fWindowSizeHintX;
1233 } else if ((fGeometryMask & fHeightValue) == 0 ) {
1234
1235 // if there is only Width. Special case to be backward compatible
1236 // We set Width and Height the same to obtain a square windows.
1237
1238 G4cout << "Unrecognised geometry string \""
1239 << geomString
1240 << "\". No Height found. Using Width value instead"
1241 << G4endl;
1242 h = w;
1243 }
1244 if ( ((fGeometryMask & fXValue) == 0 ) ||
1245 ((fGeometryMask & fYValue) == 0 )) {
1246 //Using defaults
1249 }
1250 // Set the string
1251 fXGeometryString = geomString;
1252
1253 // Set values
1254 fWindowSizeHintX = w;
1255 fWindowSizeHintY = h;
1258
1259 if ( ((fGeometryMask & fXValue)) &&
1260 ((fGeometryMask & fYValue))) {
1261
1262 if ( (fGeometryMask & fXNegative) ) {
1264 } else {
1266 }
1267 if ( (fGeometryMask & fYNegative) ) {
1269 } else {
1271 }
1272 }
1273}
1274
1277 return sizeX + fWindowLocationHintX - fWindowSizeHintX;
1278 }
1279 return fWindowLocationHintX;
1280}
1281
1284 return sizeY + fWindowLocationHintY - fWindowSizeHintY;
1285 }
1286 return fWindowLocationHintY;
1287}
1288
1289/* Keep from :
1290 * ftp://ftp.trolltech.com/qt/source/qt-embedded-free-3.0.6.tar.gz/qt-embedded-free-3.0.6/src/kernel/qapplication_qws.cpp
1291 *
1292 * ParseGeometry parses strings of the form
1293 * "=<width>x<height>{+-}<xoffset>{+-}<yoffset>", where
1294 * width, height, xoffset, and yoffset are unsigned integers.
1295 * Example: "=80x24+300-49"
1296 * The equal sign is optional.
1297 * It returns a bitmask that indicates which of the four values
1298 * were actually found in the string. For each value found,
1299 * the corresponding argument is updated; for each value
1300 * not found, the corresponding argument is left unchanged.
1301 */
1302
1304 const char *string,
1305 G4int *x,
1306 G4int *y,
1307 unsigned int *width,
1308 unsigned int *height)
1309{
1310
1311 G4int mask = fNoValue;
1312 char *strind;
1313 unsigned int tempWidth = 0;
1314 unsigned int tempHeight = 0;
1315 G4int tempX = 0;
1316 G4int tempY = 0;
1317 char *nextCharacter;
1318 if ( (string == NULL) || (*string == '\0')) {
1319 return(mask);
1320 }
1321 if (*string == '=')
1322 string++; /* ignore possible '=' at beg of geometry spec */
1323 strind = (char *)string;
1324 if (*strind != '+' && *strind != '-' && *strind != 'x') {
1325 tempWidth = ReadInteger(strind, &nextCharacter);
1326 if (strind == nextCharacter)
1327 return (0);
1328 strind = nextCharacter;
1329 mask |= fWidthValue;
1330 }
1331 if (*strind == 'x' || *strind == 'X') {
1332 strind++;
1333 tempHeight = ReadInteger(strind, &nextCharacter);
1334 if (strind == nextCharacter)
1335 return (0);
1336 strind = nextCharacter;
1337 mask |= fHeightValue;
1338 }
1339
1340 if ((*strind == '+') || (*strind == '-')) {
1341 if (*strind == '-') {
1342 strind++;
1343 tempX = -ReadInteger(strind, &nextCharacter);
1344 if (strind == nextCharacter)
1345 return (0);
1346 strind = nextCharacter;
1347 mask |= fXNegative;
1348
1349 }
1350 else
1351 { strind++;
1352 tempX = ReadInteger(strind, &nextCharacter);
1353 if (strind == nextCharacter)
1354 return(0);
1355 strind = nextCharacter;
1356 }
1357 mask |= fXValue;
1358 if ((*strind == '+') || (*strind == '-')) {
1359 if (*strind == '-') {
1360 strind++;
1361 tempY = -ReadInteger(strind, &nextCharacter);
1362 if (strind == nextCharacter)
1363 return(0);
1364 strind = nextCharacter;
1365 mask |= fYNegative;
1366 }
1367 else
1368 {
1369 strind++;
1370 tempY = ReadInteger(strind, &nextCharacter);
1371 if (strind == nextCharacter)
1372 return(0);
1373 strind = nextCharacter;
1374 }
1375 mask |= fYValue;
1376 }
1377 }
1378 /* If strind isn't at the end of the string the it's an invalid
1379 geometry specification. */
1380 if (*strind != '\0') return (0);
1381 if (mask & fXValue)
1382 *x = tempX;
1383 if (mask & fYValue)
1384 *y = tempY;
1385 if (mask & fWidthValue)
1386 *width = tempWidth;
1387 if (mask & fHeightValue)
1388 *height = tempHeight;
1389 return (mask);
1390}
1391
1392/* Keep from :
1393 * ftp://ftp.trolltech.com/qt/source/qt-embedded-free-3.0.6.tar.gz/qt-embedded-free-3.0.6/src/kernel/qapplication_qws.cpp
1394 *
1395 */
1396G4int G4ViewParameters::ReadInteger(char *string, char **NextString)
1397{
1398 G4int Result = 0;
1399 G4int Sign = 1;
1400
1401 if (*string == '+')
1402 string++;
1403 else if (*string == '-')
1404 {
1405 string++;
1406 Sign = -1;
1407 }
1408 for (; (*string >= '0') && (*string <= '9'); string++)
1409 {
1410 Result = (Result * 10) + (*string - '0');
1411 }
1412 *NextString = string;
1413 if (Sign >= 0)
1414 return (Result);
1415 else
1416 return (-Result);
1417}
1418
1420(const std::vector<G4ViewParameters>& views,
1421 G4int nInterpolationPoints) // No of interpolations points per interval
1422{
1423 // Returns a null pointer when no more to be done. For example:
1424 // do {
1425 // G4ViewParameters* vp =
1426 // G4ViewParameters::CatmullRomCubicSplineInterpolation(viewVector,nInterpolationPoints);
1427 // if (!vp) break;
1428 // ...
1429 // } while (true);
1430
1431 // See https://en.wikipedia.org/wiki/Cubic_Hermite_spline
1432
1433 // Assumes equal intervals
1434
1435 if (views.size() < 2) {
1437 ("G4ViewParameters::CatmullRomCubicSplineInterpolation",
1438 "visman0301", JustWarning,
1439 "There must be at least two views.");
1440 return 0;
1441 }
1442
1443 if (nInterpolationPoints < 1) {
1445 ("G4ViewParameters::CatmullRomCubicSplineInterpolation",
1446 "visman0302", JustWarning,
1447 "Number of interpolation points cannot be zero or negative.");
1448 return 0;
1449 }
1450
1451 const size_t nIntervals = views.size() - 1;
1452 const G4double dt = 1./nInterpolationPoints;
1453
1454 static G4ViewParameters holdingValues;
1455 static G4double t = 0.; // 0. <= t <= 1.
1456 static G4int iInterpolationPoint = 0;
1457 static size_t iInterval = 0;
1458
1459// G4cout << "Interval " << iInterval << ", t = " << t << G4endl;
1460
1461 // Hermite polynomials.
1462 const G4double h00 = 2.*t*t*t - 3.*t*t +1;
1463 const G4double h10 = t*t*t -2.*t*t + t;
1464 const G4double h01 = -2.*t*t*t + 3.*t*t;
1465 const G4double h11 = t*t*t - t*t;
1466
1467 // Aliases (to simplify code)
1468 const size_t& n = nIntervals;
1469 size_t& i = iInterval;
1470 const std::vector<G4ViewParameters>& v = views;
1471
1472 // The Catmull-Rom cubic spline prescription is as follows:
1473 // Slope at first way point is v[1] - v[0].
1474 // Slope at last way point is v[n] - v[n-1].
1475 // Otherwise slope at way point i is 0.5*(v[i+1] - v[i-1]).
1476 // Result = h00*v[i] + h10*m[i] + h01*v[i+1] + h11*m[i+1],
1477 // where m[i] amd m[i+1] are the slopes at the start and end
1478 // of the interval for the particular value.
1479 // If (n == 1), linear interpolation results.
1480 // If (n == 2), quadratic interpolation results.
1481
1482 // Working variables
1483 G4double mi, mi1, real, x, y, z;
1484
1485 // First, a crude interpolation of all parameters. Then, below, a
1486 // smooth interpolation of those for which it makes sense.
1487 holdingValues = t < 0.5? v[i]: v[i+1];
1488
1489 // Catmull-Rom cubic spline interpolation
1490#define INTERPOLATE(param) \
1491/* This works out the interpolated param in i'th interval */ \
1492/* Assumes n >= 1 */ \
1493if (i == 0) { \
1494/* First interval */ \
1495mi = v[1].param - v[0].param; \
1496/* If there is only one interval, make start and end slopes equal */ \
1497/* (This results in a linear interpolation) */ \
1498if (n == 1) mi1 = mi; \
1499/* else the end slope of the interval takes account of the next waypoint along */ \
1500else mi1 = 0.5 * (v[2].param - v[0].param); \
1501} else if (i >= n - 1) { \
1502/* Similarly for last interval */ \
1503mi1 = v[i+1].param - v[i].param; \
1504/* If there is only one interval, make start and end slopes equal */ \
1505if (n == 1) mi = mi1; \
1506/* else the start slope of the interval takes account of the previous waypoint */ \
1507else mi = 0.5 * (v[i+1].param - v[i-1].param); \
1508} else { \
1509/* Full Catmull-Rom slopes use previous AND next waypoints */ \
1510mi = 0.5 * (v[i+1].param - v[i-1].param); \
1511mi1 = 0.5 * (v[i+2].param - v[i ].param); \
1512} \
1513real = h00 * v[i].param + h10 * mi + h01 * v[i+1].param + h11 * mi1;
1514
1515#define INTERPOLATELOG(param) \
1516if (i == 0) { \
1517mi = std::log(v[1].param) - std::log(v[0].param); \
1518if (n == 1) mi1 = mi; \
1519else mi1 = 0.5 * (std::log(v[2].param) - std::log(v[0].param)); \
1520} else if (i >= n - 1) { \
1521mi1 = std::log(v[i+1].param) - std::log(v[i].param); \
1522if (n == 1) mi = mi1; \
1523else mi = 0.5 * (std::log(v[i+1].param) - std::log(v[i-1].param)); \
1524} else { \
1525mi = 0.5 * (std::log(v[i+1].param) - std::log(v[i-1].param)); \
1526mi1 = 0.5 * (std::log(v[i+2].param) - std::log(v[i ].param)); \
1527} \
1528real = std::exp(h00 * std::log(v[i].param) + h10 * mi + h01 * std::log(v[i+1].param) + h11 * mi1);
1529
1530 // Real parameters
1532 if (real < 0.) real = 0.;
1533 holdingValues.fVisibleDensity = real;
1535 holdingValues.fExplodeFactor = real;
1537 if (real < 0.) real = 0.;
1538 holdingValues.fFieldHalfAngle = real;
1540 holdingValues.fZoomFactor = real;
1542 holdingValues.fDolly = real;
1544 if (real < 0.) real = 0.;
1545 holdingValues.fGlobalMarkerScale = real;
1547 if (real < 0.) real = 0.;
1548 holdingValues.fGlobalLineWidthScale = real;
1549
1550 // Unit vectors
1551#define INTERPOLATEUNITVECTOR(vector) \
1552INTERPOLATE(vector.x()); x = real; \
1553INTERPOLATE(vector.y()); y = real; \
1554INTERPOLATE(vector.z()); z = real;
1556 holdingValues.fViewpointDirection = G4Vector3D(x,y,z).unit();
1558 holdingValues.fUpVector = G4Vector3D(x,y,z).unit();
1560 holdingValues.fRelativeLightpointDirection = G4Vector3D(x,y,z).unit();
1562 holdingValues.fActualLightpointDirection = G4Vector3D(x,y,z).unit();
1563
1564 // Un-normalised vectors
1565#define INTERPOLATEVECTOR(vector) \
1566INTERPOLATE(vector.x()); x = real; \
1567INTERPOLATE(vector.y()); y = real; \
1568INTERPOLATE(vector.z()); z = real;
1570 holdingValues.fScaleFactor = G4Vector3D(x,y,z);
1571
1572 // Points
1573#define INTERPOLATEPOINT(point) \
1574INTERPOLATE(point.x()); x = real; \
1575INTERPOLATE(point.y()); y = real; \
1576INTERPOLATE(point.z()); z = real;
1578 holdingValues.fExplodeCentre = G4Point3D(x,y,z);
1580 holdingValues.fCurrentTargetPoint = G4Point3D(x,y,z);
1581
1582 // Colour
1583 G4double red, green, blue, alpha;
1584#define INTERPOLATECOLOUR(colour) \
1585INTERPOLATE(colour.GetRed()); red = real; \
1586INTERPOLATE(colour.GetGreen()); green = real; \
1587INTERPOLATE(colour.GetBlue()); blue = real; \
1588INTERPOLATE(colour.GetAlpha()); alpha = real;
1590 // Components are clamped to 0. <= component <= 1.
1591 holdingValues.fBackgroundColour = G4Colour(red,green,blue,alpha);
1592
1593 // For some parameters we need to check some continuity
1594 G4bool continuous;
1595#define CONTINUITY(quantity) \
1596 continuous = false; \
1597 /* This follows the logic of the INTERPOLATE macro above; see comments therein */ \
1598 if (i == 0) { \
1599 if (v[1].quantity == v[0].quantity) { \
1600 if (n == 1) continuous = true; \
1601 else if (v[2].quantity == v[0].quantity) \
1602 continuous = true; \
1603 } \
1604 } else if (i >= n - 1) { \
1605 if (v[i+1].quantity == v[i].quantity) { \
1606 if (n == 1) continuous = true; \
1607 else if (v[i+1].quantity == v[i-1].quantity) \
1608 continuous = true; \
1609 } \
1610 } else { \
1611 if (v[i-1].quantity == v[i].quantity && \
1612 v[i+1].quantity == v[i].quantity && \
1613 v[i+2].quantity == v[i].quantity) \
1614 continuous = true; \
1615 }
1616
1617 G4double a, b, c, d;
1618#define INTERPOLATEPLANE(plane) \
1619INTERPOLATE(plane.a()); a = real; \
1620INTERPOLATE(plane.b()); b = real; \
1621INTERPOLATE(plane.c()); c = real; \
1622INTERPOLATE(plane.d()); d = real;
1623
1624 // Section plane
1626 if (continuous) {
1628 holdingValues.fSectionPlane = G4Plane3D(a,b,c,d);
1629 }
1630
1631 // Cutaway planes
1632 if (v[i].fCutawayPlanes.size()) {
1633 CONTINUITY(fCutawayPlanes.size());
1634 if (continuous) {
1635 for (size_t j = 0; j < v[i].fCutawayPlanes.size(); ++j) {
1637 holdingValues.fCutawayPlanes[j] = G4Plane3D(a,b,c,d);
1638 }
1639 }
1640 }
1641
1642 // Vis attributes modifiers
1643 // Really, we are only interested in colour - other attributes can follow
1644 // the "crude" interpolation that is guaranteed above.
1645 static G4VisAttributes workingVA;
1646 if (v[i].fVisAttributesModifiers.size()) {
1648 if (continuous) {
1649 for (size_t j = 0; j < v[i].fVisAttributesModifiers.size(); ++j) {
1650 CONTINUITY(fVisAttributesModifiers[j].GetPVNameCopyNoPath());
1651 if (continuous) {
1652 CONTINUITY(fVisAttributesModifiers[j].GetVisAttributesSignifier());
1653 if (continuous) {
1654 if (v[i].fVisAttributesModifiers[j].GetVisAttributesSignifier() ==
1656 INTERPOLATECOLOUR(fVisAttributesModifiers[j].GetVisAttributes().GetColour());
1657 workingVA = v[i].fVisAttributesModifiers[j].GetVisAttributes();
1658 workingVA.SetColour(G4Colour(red,green,blue,alpha));
1659 holdingValues.fVisAttributesModifiers[j].SetVisAttributes(workingVA);
1660 }
1661 }
1662 }
1663 }
1664 }
1665 }
1666
1667 // Time window parameters (for showing particles in flight)
1668 // Only two parameters are interpolated. The others are usually chosen
1669 // once and for all by the user for a given series of views - or at least,
1670 // if not, they will be interpolated by the default "crude" method above.
1672 holdingValues.fStartTime = real;
1674 holdingValues.fEndTime = real;
1675
1676 // Increment counters
1677 iInterpolationPoint++;
1678 t += dt;
1679 if (iInterpolationPoint > nInterpolationPoints) {
1680 iInterpolationPoint = 1; // Ready for next interval.
1681 t = dt;
1682 iInterval++;
1683 }
1684 if (iInterval >= nIntervals) {
1685 iInterpolationPoint = 0; // Ready for a complete restart.
1686 t = 0.;
1687 iInterval = 0;
1688 return 0;
1689 }
1690
1691 return &holdingValues;
1692}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static const G4double alpha
HepGeom::Plane3D< G4double > G4Plane3D
Definition: G4Plane3D.hh:36
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:34
static constexpr double cm3
Definition: G4SIunits.hh:101
static constexpr double mm
Definition: G4SIunits.hh:95
static constexpr double g
Definition: G4SIunits.hh:168
static constexpr double deg
Definition: G4SIunits.hh:132
#define G4BestUnit(a, b)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:34
#define INTERPOLATECOLOUR(colour)
#define INTERPOLATEVECTOR(vector)
#define INTERPOLATE(param)
#define INTERPOLATEPLANE(plane)
#define CONTINUITY(quantity)
#define INTERPOLATEPOINT(point)
#define INTERPOLATEUNITVECTOR(vector)
#define INTERPOLATELOG(param)
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetBlue() const
Definition: G4Colour.hh:154
G4double GetAlpha() const
Definition: G4Colour.hh:155
G4double GetRed() const
Definition: G4Colour.hh:152
G4double GetGreen() const
Definition: G4Colour.hh:153
const PVNameCopyNoPath & GetPVNameCopyNoPath() const
VisAttributesSignifier GetVisAttributesSignifier() const
const G4VisAttributes & GetVisAttributes() const
PVNameCopyNoPath::const_iterator PVNameCopyNoPathConstIterator
std::vector< PVNameCopyNo > PVNameCopyNoPath
void SetScreenSize(G4double)
G4int SetNumberOfCloudPoints(G4int)
G4Point3D fCurrentTargetPoint
static G4ViewParameters * CatmullRomCubicSplineInterpolation(const std::vector< G4ViewParameters > &views, G4int nInterpolationPoints=50)
G4int SetNoOfSides(G4int nSides)
void SetViewAndLights(const G4Vector3D &viewpointDirection)
G4double fDisplayLightFrontT
std::vector< G4ModelingParameters::VisAttributesModifier > fVisAttributesModifiers
G4double fDisplayLightFrontBlue
G4int GetWindowAbsoluteLocationHintY(G4int) const
G4String CameraAndLightingCommands(const G4Point3D standardTargetPoint) const
G4double GetCameraDistance(G4double radius) const
void PrintDifferences(const G4ViewParameters &v) const
void SetVisibleDensity(G4double visibleDensity)
G4bool fWindowLocationHintYNegative
G4double fDisplayHeadTimeBlue
G4bool IsCutaway() const
void AddVisAttributesModifier(const G4ModelingParameters::VisAttributesModifier &)
std::vector< G4ModelingParameters::PVNameCopyNo > fSpecialMeshVolumes
G4Vector3D & GetActualLightpointDirection()
G4double fDisplayLightFrontY
CutawayMode fCutawayMode
void SetXGeometryString(const G4String &)
G4double fDisplayLightFrontRed
G4double fDisplayLightFrontGreen
DrawingStyle fDrawingStyle
G4int ReadInteger(char *string, char **NextString)
G4double GetFarDistance(G4double cameraDistance, G4double nearDistance, G4double radius) const
void MultiplyScaleFactor(const G4Vector3D &scaleFactorMultiplier)
G4double GetFrontHalfHeight(G4double nearDistance, G4double radius) const
G4Vector3D fViewpointDirection
G4int GetWindowAbsoluteLocationHintX(G4int) const
std::vector< G4double > fCBDParameters
G4bool IsExplode() const
G4String SceneModifyingCommands() const
void IncrementPan(G4double right, G4double up)
G4double fDisplayHeadTimeGreen
G4double fDisplayLightFrontX
G4String TimeWindowCommands() const
G4String TouchableCommands() const
G4double fGlobalLineWidthScale
G4VisAttributes fDefaultTextVisAttributes
void ChangeCutawayPlane(size_t index, const G4Plane3D &cutawayPlane)
G4double fDisplayHeadTimeSize
void SetPan(G4double right, G4double up)
G4Vector3D fRelativeLightpointDirection
G4VisAttributes fDefaultVisAttributes
void SetLightpointDirection(const G4Vector3D &lightpointDirection)
G4double fDisplayHeadTimeRed
G4String DrawingStyleCommands() const
G4double fDisplayLightFrontZ
G4bool operator!=(const G4ViewParameters &) const
G4bool fWindowLocationHintXNegative
G4int ParseGeometry(const char *string, G4int *x, G4int *y, unsigned int *width, unsigned int *height)
G4double GetNearDistance(G4double cameraDistance, G4double radius) const
G4Vector3D fActualLightpointDirection
void AddCutawayPlane(const G4Plane3D &cutawayPlane)
RotationStyle fRotationStyle
G4int GetForcedNumberOfCloudPoints() const
G4double GetLineWidth() const
G4bool IsDaughtersInvisible() const
void SetColour(const G4Colour &)
G4int GetForcedLineSegmentsPerCircle() const
LineStyle GetLineStyle() const
const G4Colour & GetColour() const
G4bool IsVisible() const
G4bool IsForceAuxEdgeVisible() const
G4bool IsForcedAuxEdgeVisible() const
ForcedDrawingStyle GetForcedDrawingStyle() const
static G4int GetMinLineSegmentsPerCircle()
G4bool IsForceDrawingStyle() const
BasicVector3D< T > cross(const BasicVector3D< T > &v) const
BasicVector3D< T > unit() const
std::ostream & operator<<(std::ostream &, const BasicVector3D< float > &)
Point3D< T > point(const Point3D< T > &p) const
Definition: Plane3D.h:115
Normal3D< T > normal() const
Definition: Plane3D.h:97
static G4int GetNumberOfRotationSteps()
#define ns
Definition: xmlparse.cc:614