Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VVisCommand Class Reference

#include <G4VVisCommand.hh>

+ Inheritance diagram for G4VVisCommand:

Public Member Functions

 G4VVisCommand ()
 
virtual ~G4VVisCommand ()
 
- Public Member Functions inherited from G4UImessenger
 G4UImessenger ()=default
 
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
virtual ~G4UImessenger ()
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
virtual void SetNewValue (G4UIcommand *command, G4String newValue)
 
G4bool CommandsShouldBeInMaster () const
 

Static Public Member Functions

static G4VisManagerGetVisManager ()
 
static void SetVisManager (G4VisManager *pVisManager)
 
static const G4ColourGetCurrentTextColour ()
 

Protected Member Functions

void SetViewParameters (G4VViewer *viewer, const G4ViewParameters &viewParams)
 
void RefreshIfRequired (G4VViewer *viewer)
 
void InterpolateViews (G4VViewer *currentViewer, std::vector< G4ViewParameters > viewVector, const G4int nInterpolationPoints=50, const G4int waitTimePerPointmilliseconds=20, const G4String exportString="")
 
void InterpolateToNewView (G4VViewer *currentViewer, const G4ViewParameters &oldVP, const G4ViewParameters &newVP, const G4int nInterpolationPoints=50, const G4int waitTimePerPointmilliseconds=20, const G4String exportString="")
 
void Twinkle (G4VViewer *currentViewer, const G4ViewParameters &baseVP, const std::vector< std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > > &paths)
 
const G4StringConvertToColourGuidance ()
 
void ConvertToColour (G4Colour &colour, const G4String &redOrString, G4double green, G4double blue, G4double opacity)
 
G4bool ProvideValueOfUnit (const G4String &where, const G4String &unit, const G4String &category, G4double &value)
 
void CopyCameraParameters (G4ViewParameters &target, const G4ViewParameters &from)
 
void CheckSceneAndNotifyHandlers (G4Scene *=nullptr)
 
G4bool CheckView ()
 
void G4VisCommandsSceneAddUnsuccessful (G4VisManager::Verbosity verbosity)
 
void CopyGuidanceFrom (const G4UIcommand *fromCmd, G4UIcommand *toCmd, G4int startLine=0)
 
void CopyParametersFrom (const G4UIcommand *fromCmd, G4UIcommand *toCmd)
 
void DrawExtent (const G4VisExtent &)
 
- Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
 
G4String LtoS (G4long l)
 
G4String DtoS (G4double a)
 
G4String BtoS (G4bool b)
 
G4int StoI (const G4String &s)
 
G4long StoL (const G4String &s)
 
G4double StoD (const G4String &s)
 
G4bool StoB (G4String s)
 
void AddUIcommand (G4UIcommand *newCommand)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
template<typename T >
T * CreateCommand (const G4String &cname, const G4String &dsc)
 

Static Protected Member Functions

static G4String ConvertToString (G4double x, G4double y, const char *unitName)
 
static G4bool ConvertToDoublePair (const G4String &paramString, G4double &xval, G4double &yval)
 

Static Protected Attributes

static G4VisManagerfpVisManager = nullptr
 
static G4int fCurrentArrow3DLineSegmentsPerCircle = 6
 
static G4Colour fCurrentColour = G4Colour::White()
 
static G4double fCurrentLineWidth = 1.
 
static G4Colour fCurrentTextColour = G4Colour::Blue()
 
static G4Text::Layout fCurrentTextLayout = G4Text::left
 
static G4double fCurrentTextSize = 12.
 
static G4PhysicalVolumeModel::TouchableProperties fCurrentTouchableProperties
 
static G4VisExtent fCurrentExtentForField
 
static std::vector< G4PhysicalVolumesSearchScene::FindingsfCurrrentPVFindingsForField
 
static G4bool fThereWasAViewer = false
 
static G4ViewParameters fExistingVP
 
static G4SceneTreeItem fExistingSceneTree
 

Additional Inherited Members

- Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Detailed Description

Definition at line 49 of file G4VVisCommand.hh.

Constructor & Destructor Documentation

◆ G4VVisCommand()

G4VVisCommand::G4VVisCommand ( )

Definition at line 60 of file G4VVisCommand.cc.

60{}

◆ ~G4VVisCommand()

G4VVisCommand::~G4VVisCommand ( )
virtual

Definition at line 62 of file G4VVisCommand.cc.

62{}

Member Function Documentation

◆ CheckSceneAndNotifyHandlers()

void G4VVisCommand::CheckSceneAndNotifyHandlers ( G4Scene * pScene = nullptr)
protected

Definition at line 213 of file G4VVisCommand.cc.

214{
216
217 if (!pScene) {
218 if (verbosity >= G4VisManager::warnings) {
219 G4warn << "WARNING: Scene pointer is null."
220 << G4endl;
221 }
222 return;
223 }
224
225 G4VSceneHandler* pSceneHandler = fpVisManager -> GetCurrentSceneHandler();
226 if (!pSceneHandler) {
227 if (verbosity >= G4VisManager::warnings) {
228 G4warn << "WARNING: Scene handler not found." << G4endl;
229 }
230 return;
231 }
232
233 // Scene has changed. If it is the scene of the currrent scene handler
234 // refresh viewers of all scene handlers using this scene. If not, it may be
235 // a scene that the user is building up before attaching to a scene handler,
236 // so do nothing.
237 if (pScene == pSceneHandler->GetScene()) {
238 G4UImanager::GetUIpointer () -> ApplyCommand ("/vis/scene/notifyHandlers");
239 }
240
241}
#define G4warn
Definition G4Scene.cc:41
#define G4endl
Definition G4ios.hh:67
static G4UImanager * GetUIpointer()
G4Scene * GetScene() const
static G4VisManager * fpVisManager
static Verbosity GetVerbosity()

Referenced by G4VisCommandSceneActivateModel::SetNewValue(), G4VisCommandSceneAddArrow2D::SetNewValue(), G4VisCommandSceneAddArrow::SetNewValue(), G4VisCommandSceneAddAxes::SetNewValue(), G4VisCommandSceneAddDate::SetNewValue(), G4VisCommandSceneAddDigis::SetNewValue(), G4VisCommandSceneAddElectricField::SetNewValue(), G4VisCommandSceneAddEventID::SetNewValue(), G4VisCommandSceneAddExtent::SetNewValue(), G4VisCommandSceneAddFrame::SetNewValue(), G4VisCommandSceneAddGPS::SetNewValue(), G4VisCommandSceneAddHits::SetNewValue(), G4VisCommandSceneAddLine2D::SetNewValue(), G4VisCommandSceneAddLine::SetNewValue(), G4VisCommandSceneAddLocalAxes::SetNewValue(), G4VisCommandSceneAddLogicalVolume::SetNewValue(), G4VisCommandSceneAddLogo2D::SetNewValue(), G4VisCommandSceneAddLogo::SetNewValue(), G4VisCommandSceneAddMagneticField::SetNewValue(), G4VisCommandSceneAddPlotter::SetNewValue(), G4VisCommandSceneAddPSHits::SetNewValue(), G4VisCommandSceneAddScale::SetNewValue(), G4VisCommandSceneAddText2D::SetNewValue(), G4VisCommandSceneAddText::SetNewValue(), G4VisCommandSceneAddTrajectories::SetNewValue(), G4VisCommandSceneAddUserAction::SetNewValue(), G4VisCommandSceneAddVolume::SetNewValue(), G4VisCommandSceneRemoveModel::SetNewValue(), G4VisCommandSceneSelect::SetNewValue(), and G4VisCommandViewerRefresh::SetNewValue().

◆ CheckView()

G4bool G4VVisCommand::CheckView ( )
protected

Definition at line 243 of file G4VVisCommand.cc.

244{
246 G4VViewer* viewer = fpVisManager -> GetCurrentViewer ();
247
248 if (!viewer) {
249 if (verbosity >= G4VisManager::errors) {
250 G4warn <<
251 "ERROR: No current viewer - \"/vis/viewer/list\" to see possibilities."
252 << G4endl;
253 }
254 return false;
255 }
256
257 return true;
258}

◆ ConvertToColour()

void G4VVisCommand::ConvertToColour ( G4Colour & colour,
const G4String & redOrString,
G4double green,
G4double blue,
G4double opacity )
protected

Definition at line 125 of file G4VVisCommand.cc.

128{
129 // Note: colour is supplied by the caller and some or all of its components
130 // may act as default.
131 //
132 // Note: redOrString is either a number or string. If a string it must be
133 // one of the recognised colours.
134 //
135 // Thus the arguments can be, for example:
136 // (colour,"red",...,...,0.5): will give the colour red with opacity 0.5 (the
137 // third and fourth arguments are ignored), or
138 // (1.,0.,0.,0.5): this also will be red with opacity 0.5.
139
141
142 const std::size_t iPos0 = 0;
143 if (std::isalpha(redOrString[iPos0])) {
144
145 // redOrString is probably alphabetic characters defining the colour
146 if (!G4Colour::GetColour(redOrString, colour)) {
147 // Not a recognised string
148 if (verbosity >= G4VisManager::warnings) {
149 G4warn << "WARNING: Colour \"" << redOrString
150 << "\" not found. Defaulting to " << colour
151 << G4endl;
152 }
153 return;
154 } else {
155 // It was a recognised string. Now add opacity.
156 colour.SetAlpha(opacity);
157 return;
158 }
159
160 } else {
161
162 // redOrString is probably numeric defining the red component
163 std::istringstream iss(redOrString);
164 G4double red;
165 iss >> red;
166 if (iss.fail()) {
167 if (verbosity >= G4VisManager::warnings) {
168 G4warn << "WARNING: String \"" << redOrString
169 << "\" cannot be parsed. Defaulting to " << colour
170 << G4endl;
171 }
172 return;
173 } else {
174 colour = G4Colour(red,green,blue,opacity);
175 return;
176 }
177
178 }
179}
double G4double
Definition G4Types.hh:83
static G4bool GetColour(const G4String &key, G4Colour &result)
Definition G4Colour.cc:155
void SetAlpha(G4double)
Definition G4Colour.cc:70

Referenced by G4VisCommandGeometrySetColour::SetNewValue(), G4VisCommandSceneAddGPS::SetNewValue(), G4VisCommandSetColour::SetNewValue(), G4VisCommandSetTextColour::SetNewValue(), G4VisCommandsTouchableSet::SetNewValue(), and G4VisCommandsViewerSet::SetNewValue().

◆ ConvertToColourGuidance()

const G4String & G4VVisCommand::ConvertToColourGuidance ( )
protected

Definition at line 115 of file G4VVisCommand.cc.

116{
117 static G4String guidance
118 ("Accepts (a) RGB triplet. e.g., \".3 .4 .5\", or"
119 "\n (b) string such as \"white\", \"black\", \"grey\", \"red\"...or"
120 "\n (c) an additional number for opacity, e.g., \".3 .4 .5 .6\""
121 "\n or \"grey ! ! .6\" (note \"!\"'s for unused parameters).");
122 return guidance;
123}

Referenced by G4VisCommandSceneAddGPS::G4VisCommandSceneAddGPS(), G4VisCommandSetColour::G4VisCommandSetColour(), G4VisCommandSetTextColour::G4VisCommandSetTextColour(), G4VisCommandsTouchableSet::G4VisCommandsTouchableSet(), and G4VisCommandsViewerSet::G4VisCommandsViewerSet().

◆ ConvertToDoublePair()

G4bool G4VVisCommand::ConvertToDoublePair ( const G4String & paramString,
G4double & xval,
G4double & yval )
staticprotected

Definition at line 91 of file G4VVisCommand.cc.

94{
95 G4double x, y;
96 G4String unit;
97
98 std::istringstream is(paramString);
99 is >> x >> y >> unit;
100
102 xval = x*G4UIcommand::ValueOf(unit);
103 yval = y*G4UIcommand::ValueOf(unit);
104 } else {
106 if (verbosity >= G4VisManager::errors) {
107 G4warn << "ERROR: Unrecognised unit" << G4endl;
108 }
109 return false;
110 }
111
112 return true;
113}
static G4double ValueOf(const char *unitName)
static G4bool IsUnitDefined(const G4String &)

Referenced by G4VisCommandsViewerSet::SetNewValue(), and G4VisCommandViewerPan::SetNewValue().

◆ ConvertToString()

G4String G4VVisCommand::ConvertToString ( G4double x,
G4double y,
const char * unitName )
staticprotected

Definition at line 81 of file G4VVisCommand.cc.

83{
84 G4double uv = G4UIcommand::ValueOf(unitName);
85
86 std::ostringstream oss;
87 oss << x/uv << " " << y/uv << " " << unitName;
88 return oss.str();
89}

Referenced by G4VisCommandViewerPan::GetCurrentValue().

◆ CopyCameraParameters()

void G4VVisCommand::CopyCameraParameters ( G4ViewParameters & target,
const G4ViewParameters & from )
protected

Definition at line 430 of file G4VVisCommand.cc.

432{
433 // Copy view parameters pertaining only to camera
437 target.SetUpVector (from.GetUpVector());
438 target.SetFieldHalfAngle (from.GetFieldHalfAngle());
439 target.SetZoomFactor (from.GetZoomFactor());
440 target.SetScaleFactor (from.GetScaleFactor());
442 target.SetDolly (from.GetDolly());
443}
void SetViewpointDirection(const G4Vector3D &viewpointDirection)
void SetScaleFactor(const G4Vector3D &scaleFactor)
const G4Vector3D & GetScaleFactor() const
void SetCurrentTargetPoint(const G4Point3D &currentTargetPoint)
const G4Vector3D & GetLightpointDirection() const
void SetFieldHalfAngle(G4double fieldHalfAngle)
const G4Vector3D & GetViewpointDirection() const
const G4Point3D & GetCurrentTargetPoint() const
G4double GetFieldHalfAngle() const
G4double GetZoomFactor() const
void SetDolly(G4double dolly)
const G4Vector3D & GetUpVector() const
void SetZoomFactor(G4double zoomFactor)
void SetUpVector(const G4Vector3D &upVector)
void SetLightpointDirection(const G4Vector3D &lightpointDirection)
void SetLightsMoveWithCamera(G4bool moves)
G4bool GetLightsMoveWithCamera() const
G4double GetDolly() const

Referenced by G4VisCommandViewerCopyViewFrom::SetNewValue(), and G4VisCommandViewerResetCameraParameters::SetNewValue().

◆ CopyGuidanceFrom()

void G4VVisCommand::CopyGuidanceFrom ( const G4UIcommand * fromCmd,
G4UIcommand * toCmd,
G4int startLine = 0 )
protected

Definition at line 391 of file G4VVisCommand.cc.

393{
394 if (fromCmd && toCmd) {
395 const G4int nGuideEntries = (G4int)fromCmd->GetGuidanceEntries();
396 for (G4int i = startLine; i < nGuideEntries; ++i) {
397 const G4String& guidance = fromCmd->GetGuidanceLine(i);
398 toCmd->SetGuidance(guidance);
399 }
400 }
401}
int G4int
Definition G4Types.hh:85
const G4String & GetGuidanceLine(G4int i) const
void SetGuidance(const char *aGuidance)
std::size_t GetGuidanceEntries() const

Referenced by G4VisCommandDrawLogicalVolume::G4VisCommandDrawLogicalVolume(), G4VisCommandDrawVolume::G4VisCommandDrawVolume(), G4VisCommandOpen::G4VisCommandOpen(), G4VisCommandSceneAddMagneticField::G4VisCommandSceneAddMagneticField(), G4VisCommandsTouchable::G4VisCommandsTouchable(), and G4VisCommandViewerCentreOn::G4VisCommandViewerCentreOn().

◆ CopyParametersFrom()

void G4VVisCommand::CopyParametersFrom ( const G4UIcommand * fromCmd,
G4UIcommand * toCmd )
protected

Definition at line 403 of file G4VVisCommand.cc.

405{
406 if (fromCmd && toCmd) {
407 const G4int nParEntries = (G4int)fromCmd->GetParameterEntries();
408 for (G4int i = 0; i < nParEntries; ++i) {
409 G4UIparameter* parameter = new G4UIparameter(*(fromCmd->GetParameter(i)));
410 toCmd->SetParameter(parameter);
411 }
412 }
413}
std::size_t GetParameterEntries() const
G4UIparameter * GetParameter(G4int i) const
void SetParameter(G4UIparameter *const newParameter)

Referenced by G4VisCommandDrawLogicalVolume::G4VisCommandDrawLogicalVolume(), G4VisCommandDrawVolume::G4VisCommandDrawVolume(), G4VisCommandSceneAddMagneticField::G4VisCommandSceneAddMagneticField(), and G4VisCommandViewerCentreOn::G4VisCommandViewerCentreOn().

◆ DrawExtent()

void G4VVisCommand::DrawExtent ( const G4VisExtent & extent)
protected

Definition at line 415 of file G4VVisCommand.cc.

416{
417 if (fpVisManager) {
418 const G4double halfX = (extent.GetXmax() - extent.GetXmin()) / 2.;
419 const G4double halfY = (extent.GetYmax() - extent.GetYmin()) / 2.;
420 const G4double halfZ = (extent.GetZmax() - extent.GetZmin()) / 2.;
421 if (halfX > 0. && halfY > 0. && halfZ > 0.) {
422 const G4Box box("vis_extent",halfX,halfY,halfZ);
423 const G4VisAttributes visAtts(G4Color::Red());
424 const G4Point3D& centre = extent.GetExtentCenter();
425 fpVisManager->Draw(box,visAtts,G4Translate3D(centre));
426 }
427 }
428}
HepGeom::Translate3D G4Translate3D
Definition G4Box.hh:56
static G4Colour Red()
Definition G4Colour.hh:161
G4double GetYmin() const
G4double GetXmax() const
const G4Point3D & GetExtentCenter() const
G4double GetYmax() const
G4double GetZmax() const
G4double GetZmin() const
G4double GetXmin() const
void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())

Referenced by G4VisCommandSceneShowExtents::SetNewValue(), G4VisCommandSetVolumeForField::SetNewValue(), and G4VisCommandsTouchable::SetNewValue().

◆ G4VisCommandsSceneAddUnsuccessful()

◆ GetCurrentTextColour()

const G4Colour & G4VVisCommand::GetCurrentTextColour ( )
static

Definition at line 76 of file G4VVisCommand.cc.

77{
78 return fCurrentTextColour;
79}
static G4Colour fCurrentTextColour

◆ GetVisManager()

◆ InterpolateToNewView()

void G4VVisCommand::InterpolateToNewView ( G4VViewer * currentViewer,
const G4ViewParameters & oldVP,
const G4ViewParameters & newVP,
const G4int nInterpolationPoints = 50,
const G4int waitTimePerPointmilliseconds = 20,
const G4String exportString = "" )
protected

Definition at line 319 of file G4VVisCommand.cc.

326{
327 std::vector<G4ViewParameters> viewVector;
328 viewVector.push_back(oldVP);
329 viewVector.push_back(oldVP);
330 viewVector.push_back(newVP);
331 viewVector.push_back(newVP);
332
334 (currentViewer,
335 viewVector,
336 nInterpolationPoints,
337 waitTimePerPointmilliseconds,
338 exportString);
339}
void InterpolateViews(G4VViewer *currentViewer, std::vector< G4ViewParameters > viewVector, const G4int nInterpolationPoints=50, const G4int waitTimePerPointmilliseconds=20, const G4String exportString="")

Referenced by G4VisCommandsTouchable::SetNewValue(), and G4VisCommandViewerCentreOn::SetNewValue().

◆ InterpolateViews()

void G4VVisCommand::InterpolateViews ( G4VViewer * currentViewer,
std::vector< G4ViewParameters > viewVector,
const G4int nInterpolationPoints = 50,
const G4int waitTimePerPointmilliseconds = 20,
const G4String exportString = "" )
protected

Definition at line 294 of file G4VVisCommand.cc.

300{
301 const G4int safety = (G4int)viewVector.size()*nInterpolationPoints;
302 G4int safetyCount = 0;
303 do {
304 G4ViewParameters* vp =
305 G4ViewParameters::CatmullRomCubicSplineInterpolation(viewVector,nInterpolationPoints);
306 if (!vp) break; // Finished.
307 currentViewer->SetViewParameters(*vp);
308 currentViewer->RefreshView();
309 if (exportString == "export" &&
310 G4StrUtil::contains(currentViewer->GetName(), "OpenGL")) {
311 G4UImanager::GetUIpointer()->ApplyCommand("/vis/ogl/export");
312 }
313 currentViewer->ShowView();
314 if (waitTimePerPointmilliseconds > 0)
315 std::this_thread::sleep_for(std::chrono::milliseconds(waitTimePerPointmilliseconds));
316 } while (safetyCount++ < safety); // Loop checking, 16.02.2016, J.Allison
317}
G4int ApplyCommand(const char *aCommand)
const G4String & GetName() const
void SetViewParameters(const G4ViewParameters &vp)
Definition G4VViewer.cc:128
void RefreshView()
virtual void ShowView()
Definition G4VViewer.cc:106
static G4ViewParameters * CatmullRomCubicSplineInterpolation(const std::vector< G4ViewParameters > &views, G4int nInterpolationPoints=50)

Referenced by InterpolateToNewView(), and G4VisCommandViewerInterpolate::SetNewValue().

◆ ProvideValueOfUnit()

G4bool G4VVisCommand::ProvideValueOfUnit ( const G4String & where,
const G4String & unit,
const G4String & category,
G4double & value )
protected

Definition at line 181 of file G4VVisCommand.cc.

186{
187 // Return false if there's a problem
188
190
191 G4bool success = true;
193 if (verbosity >= G4VisManager::warnings) {
194 G4warn << where
195 << "\n Unit \"" << unit << "\" not defined"
196 << G4endl;
197 }
198 success = false;
199 } else if (G4UnitDefinition::GetCategory(unit) != category) {
200 if (verbosity >= G4VisManager::warnings) {
201 G4warn << where
202 << "\n Unit \"" << unit << "\" not a unit of " << category;
203 if (category == "Volumic Mass") G4warn << " (density)";
204 G4warn << G4endl;
205 }
206 success = false;
207 } else {
208 value = G4UnitDefinition::GetValueOf(unit);
209 }
210 return success;
211}
bool G4bool
Definition G4Types.hh:86
static G4double GetValueOf(const G4String &)
static G4String GetCategory(const G4String &)

Referenced by G4VisCommandsViewerSet::SetNewValue(), and G4VisCommandViewerColourByDensity::SetNewValue().

◆ RefreshIfRequired()

void G4VVisCommand::RefreshIfRequired ( G4VViewer * viewer)
protected

Definition at line 277 of file G4VVisCommand.cc.

277 {
279 G4VSceneHandler* sceneHandler = viewer->GetSceneHandler();
280 const G4ViewParameters& viewParams = viewer->GetViewParameters();
281 if (sceneHandler && sceneHandler->GetScene()) {
282 if (viewParams.IsAutoRefresh()) {
283 G4UImanager::GetUIpointer()->ApplyCommand("/vis/viewer/refresh");
284 }
285 else {
286 if (verbosity >= G4VisManager::warnings) {
287 G4warn << "Issue /vis/viewer/refresh or flush to see effect."
288 << G4endl;
289 }
290 }
291 }
292}
const G4ViewParameters & GetViewParameters() const
G4VSceneHandler * GetSceneHandler() const
G4bool IsAutoRefresh() const

Referenced by G4VisCommandViewerRebuild::SetNewValue(), G4VisCommandViewerReset::SetNewValue(), G4VisCommandViewerResetCameraParameters::SetNewValue(), G4VisCommandViewerSelect::SetNewValue(), and SetViewParameters().

◆ SetViewParameters()

◆ SetVisManager()

void G4VVisCommand::SetVisManager ( G4VisManager * pVisManager)
static

Definition at line 71 of file G4VVisCommand.cc.

72{
73 fpVisManager = pVisManager;
74}

Referenced by G4VisManager::G4VisManager().

◆ Twinkle()

void G4VVisCommand::Twinkle ( G4VViewer * currentViewer,
const G4ViewParameters & baseVP,
const std::vector< std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > > & paths )
protected

Definition at line 341 of file G4VVisCommand.cc.

346{
347 // Copy view parameters to temporary variables ready for adding VisAttributes Modifiers (VAMs)
348 auto loVP = baseVP; // For black and solid VAMs
349 auto hiVP = baseVP; // For white and solid VAMs
350
351 // Modify them with vis attribute modifiers (VAMs)
352 for (const auto& path: paths) {
353 const auto& touchable = path.back().GetPhysicalVolume();
354 auto loVisAtts
355 = *(currentViewer->GetApplicableVisAttributes
356 (touchable->GetLogicalVolume()->GetVisAttributes()));
357 auto hiVisAtts = loVisAtts;
358 loVisAtts.SetColour(G4Colour::Black());
359 loVisAtts.SetForceSolid();
360 hiVisAtts.SetColour(G4Colour::White());
361 hiVisAtts.SetForceSolid();
362 auto pvNameCopyNoPath
364
366 (loVisAtts, G4ModelingParameters::VASColour, pvNameCopyNoPath);
367 loVP.AddVisAttributesModifier(loVAMColour);
369 (loVisAtts, G4ModelingParameters::VASForceSolid, pvNameCopyNoPath);
370 loVP.AddVisAttributesModifier(loVAMStyle);
371
373 (hiVisAtts, G4ModelingParameters::VASColour, pvNameCopyNoPath);
374 hiVP.AddVisAttributesModifier(hiVAMColour);
376 (hiVisAtts, G4ModelingParameters::VASForceSolid, pvNameCopyNoPath);
377 hiVP.AddVisAttributesModifier(hiVAMStyle);
378 }
379
380 // Twinkle
381 std::vector<G4ViewParameters> viewVector;
382 // Just 5 twinkles is reasonable to get a human's attention
383 for (G4int i = 0; i < 5; i++) {
384 viewVector.push_back(loVP);
385 viewVector.push_back(hiVP);
386 }
387 // Just 5 interpolation points for a reasonable twinkle rate
388 InterpolateViews(currentViewer,viewVector,5);
389}
static G4Colour White()
Definition G4Colour.hh:156
static G4Colour Black()
Definition G4Colour.hh:159
static G4ModelingParameters::PVNameCopyNoPath GetPVNameCopyNoPath(const std::vector< G4PhysicalVolumeNodeID > &)
const G4VisAttributes * GetApplicableVisAttributes(const G4VisAttributes *) const

Referenced by G4VisCommandsTouchable::SetNewValue(), and G4VisCommandViewerCentreOn::SetNewValue().

Member Data Documentation

◆ fCurrentArrow3DLineSegmentsPerCircle

◆ fCurrentColour

◆ fCurrentExtentForField

◆ fCurrentLineWidth

◆ fCurrentTextColour

◆ fCurrentTextLayout

◆ fCurrentTextSize

G4double G4VVisCommand::fCurrentTextSize = 12.
staticprotected

◆ fCurrentTouchableProperties

◆ fCurrrentPVFindingsForField

◆ fExistingSceneTree

G4SceneTreeItem G4VVisCommand::fExistingSceneTree
staticprotected

◆ fExistingVP

G4ViewParameters G4VVisCommand::fExistingVP
staticprotected

◆ fpVisManager

G4VisManager * G4VVisCommand::fpVisManager = nullptr
staticprotected

Definition at line 148 of file G4VVisCommand.hh.

Referenced by CheckSceneAndNotifyHandlers(), CheckView(), ConvertToColour(), ConvertToDoublePair(), G4VVisCommandScene::CurrentSceneName(), DrawExtent(), G4VisCommandSceneHandlerCreate::G4VisCommandSceneHandlerCreate(), G4VisCommandOpen::GetCurrentValue(), G4VisCommandSceneHandlerAttach::GetCurrentValue(), G4VisCommandSceneHandlerCreate::GetCurrentValue(), G4VisCommandViewerClear::GetCurrentValue(), G4VisCommandViewerClearTransients::GetCurrentValue(), G4VisCommandViewerClone::GetCurrentValue(), G4VisCommandViewerCreate::GetCurrentValue(), G4VisCommandViewerFlush::GetCurrentValue(), G4VisCommandViewerRebuild::GetCurrentValue(), G4VisCommandViewerRefresh::GetCurrentValue(), G4VisCommandViewerReset::GetCurrentValue(), G4VisCommandViewerResetCameraParameters::GetCurrentValue(), G4VisCommandViewerUpdate::GetCurrentValue(), GetVisManager(), ProvideValueOfUnit(), RefreshIfRequired(), G4VVisCommandGeometrySet::Set(), G4VVisCommandGeometrySet::SetLVVisAtts(), G4VisCommandAbortReviewKeptEvents::SetNewValue(), G4VisCommandAbortReviewPlots::SetNewValue(), G4VisCommandDrawLogicalVolume::SetNewValue(), G4VisCommandDrawOnlyToBeKeptEvents::SetNewValue(), G4VisCommandDrawTree::SetNewValue(), G4VisCommandDrawView::SetNewValue(), G4VisCommandDrawVolume::SetNewValue(), G4VisCommandEnable::SetNewValue(), G4VisCommandGeometryList::SetNewValue(), G4VisCommandGeometryRestore::SetNewValue(), G4VisCommandGeometrySetDaughtersInvisible::SetNewValue(), G4VisCommandGeometrySetVisibility::SetNewValue(), G4VisCommandInitialize::SetNewValue(), G4VisCommandList::SetNewValue(), G4VisCommandOpen::SetNewValue(), G4VisCommandPlot::SetNewValue(), G4VisCommandReviewKeptEvents::SetNewValue(), G4VisCommandReviewPlots::SetNewValue(), G4VisCommandSceneActivateModel::SetNewValue(), G4VisCommandSceneAddArrow2D::SetNewValue(), G4VisCommandSceneAddArrow::SetNewValue(), G4VisCommandSceneAddAxes::SetNewValue(), G4VisCommandSceneAddDate::SetNewValue(), G4VisCommandSceneAddDigis::SetNewValue(), G4VisCommandSceneAddElectricField::SetNewValue(), G4VisCommandSceneAddEventID::SetNewValue(), G4VisCommandSceneAddExtent::SetNewValue(), G4VisCommandSceneAddFrame::SetNewValue(), G4VisCommandSceneAddGPS::SetNewValue(), G4VisCommandSceneAddHits::SetNewValue(), G4VisCommandSceneAddLine2D::SetNewValue(), G4VisCommandSceneAddLine::SetNewValue(), G4VisCommandSceneAddLocalAxes::SetNewValue(), G4VisCommandSceneAddLogicalVolume::SetNewValue(), G4VisCommandSceneAddLogo2D::SetNewValue(), G4VisCommandSceneAddLogo::SetNewValue(), G4VisCommandSceneAddMagneticField::SetNewValue(), G4VisCommandSceneAddPlotter::SetNewValue(), G4VisCommandSceneAddPSHits::SetNewValue(), G4VisCommandSceneAddScale::SetNewValue(), G4VisCommandSceneAddText2D::SetNewValue(), G4VisCommandSceneAddText::SetNewValue(), G4VisCommandSceneAddTrajectories::SetNewValue(), G4VisCommandSceneAddUserAction::SetNewValue(), G4VisCommandSceneAddVolume::SetNewValue(), G4VisCommandSceneCreate::SetNewValue(), G4VisCommandSceneEndOfEventAction::SetNewValue(), G4VisCommandSceneEndOfRunAction::SetNewValue(), G4VisCommandSceneHandlerAttach::SetNewValue(), G4VisCommandSceneHandlerCreate::SetNewValue(), G4VisCommandSceneHandlerList::SetNewValue(), G4VisCommandSceneHandlerSelect::SetNewValue(), G4VisCommandSceneList::SetNewValue(), G4VisCommandSceneNotifyHandlers::SetNewValue(), G4VisCommandSceneRemoveModel::SetNewValue(), G4VisCommandSceneSelect::SetNewValue(), G4VisCommandSceneShowExtents::SetNewValue(), G4VisCommandSetArrow3DLineSegmentsPerCircle::SetNewValue(), G4VisCommandSetColour::SetNewValue(), G4VisCommandSetExtentForField::SetNewValue(), G4VisCommandSetLineWidth::SetNewValue(), G4VisCommandSetTextColour::SetNewValue(), G4VisCommandSetTextLayout::SetNewValue(), G4VisCommandSetTextSize::SetNewValue(), G4VisCommandSetTouchable::SetNewValue(), G4VisCommandSetVolumeForField::SetNewValue(), G4VisCommandSpecify::SetNewValue(), G4VisCommandsTouchable::SetNewValue(), G4VisCommandsTouchableSet::SetNewValue(), G4VisCommandsViewerSet::SetNewValue(), G4VisCommandVerbose::SetNewValue(), G4VisCommandViewerAddCutawayPlane::SetNewValue(), G4VisCommandViewerCentreOn::SetNewValue(), G4VisCommandViewerChangeCutawayPlane::SetNewValue(), G4VisCommandViewerClear::SetNewValue(), G4VisCommandViewerClearCutawayPlanes::SetNewValue(), G4VisCommandViewerClearTransients::SetNewValue(), G4VisCommandViewerClearVisAttributesModifiers::SetNewValue(), G4VisCommandViewerClone::SetNewValue(), G4VisCommandViewerColourByDensity::SetNewValue(), G4VisCommandViewerCopyViewFrom::SetNewValue(), G4VisCommandViewerCreate::SetNewValue(), G4VisCommandViewerDefaultHiddenEdge::SetNewValue(), G4VisCommandViewerDefaultStyle::SetNewValue(), G4VisCommandViewerDolly::SetNewValue(), G4VisCommandViewerFlush::SetNewValue(), G4VisCommandViewerInterpolate::SetNewValue(), G4VisCommandViewerList::SetNewValue(), G4VisCommandViewerPan::SetNewValue(), G4VisCommandViewerRebuild::SetNewValue(), G4VisCommandViewerRefresh::SetNewValue(), G4VisCommandViewerReset::SetNewValue(), G4VisCommandViewerResetCameraParameters::SetNewValue(), G4VisCommandViewerSave::SetNewValue(), G4VisCommandViewerScale::SetNewValue(), G4VisCommandViewerSelect::SetNewValue(), G4VisCommandViewerUpdate::SetNewValue(), G4VisCommandViewerZoom::SetNewValue(), G4VisCommandGeometrySetVisibility::SetNewValueOnLV(), and SetVisManager().

◆ fThereWasAViewer

G4bool G4VVisCommand::fThereWasAViewer = false
staticprotected

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