Geant4 11.4.0
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, const 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 (const 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 61 of file G4VVisCommand.cc.

61{}

◆ ~G4VVisCommand()

G4VVisCommand::~G4VVisCommand ( )
virtual

Definition at line 63 of file G4VVisCommand.cc.

63{}

Member Function Documentation

◆ CheckSceneAndNotifyHandlers()

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

Definition at line 214 of file G4VVisCommand.cc.

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

Referenced by G4VisCommandSceneActivateModel::SetNewValue(), G4VisCommandSceneAddArrow2D::SetNewValue(), G4VisCommandSceneAddArrow::SetNewValue(), G4VisCommandSceneAddAxes::SetNewValue(), G4VisCommandSceneAddDate::SetNewValue(), G4VisCommandSceneAddDigis::SetNewValue(), G4VisCommandSceneAddElectricField::SetNewValue(), G4VisCommandSceneAddEndOfRunMacro::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 244 of file G4VVisCommand.cc.

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

◆ ConvertToColour()

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

Definition at line 126 of file G4VVisCommand.cc.

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

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

◆ ConvertToColourGuidance()

const G4String & G4VVisCommand::ConvertToColourGuidance ( )
protected

Definition at line 116 of file G4VVisCommand.cc.

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

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 92 of file G4VVisCommand.cc.

95{
96 G4double x, y;
97 G4String unit;
98
99 std::istringstream is(paramString);
100 is >> x >> y >> unit;
101
103 xval = x*G4UIcommand::ValueOf(unit);
104 yval = y*G4UIcommand::ValueOf(unit);
105 } else {
106 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
107 if (verbosity >= G4VisManager::errors) {
108 G4warn << "ERROR: Unrecognised unit" << G4endl;
109 }
110 return false;
111 }
112
113 return true;
114}
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 82 of file G4VVisCommand.cc.

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

Referenced by G4VisCommandViewerPan::GetCurrentValue().

◆ CopyCameraParameters()

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

Definition at line 439 of file G4VVisCommand.cc.

441{
442 // Copy view parameters pertaining only to camera
446 target.SetUpVector (from.GetUpVector());
447 target.SetFieldHalfAngle (from.GetFieldHalfAngle());
448 target.SetZoomFactor (from.GetZoomFactor());
449 target.SetScaleFactor (from.GetScaleFactor());
451 target.SetDolly (from.GetDolly());
452}
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 400 of file G4VVisCommand.cc.

402{
403 if (fromCmd && toCmd) {
404 const G4int nGuideEntries = (G4int)fromCmd->GetGuidanceEntries();
405 for (G4int i = startLine; i < nGuideEntries; ++i) {
406 const G4String& guidance = fromCmd->GetGuidanceLine(i);
407 toCmd->SetGuidance(guidance);
408 }
409 }
410}
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 412 of file G4VVisCommand.cc.

414{
415 if (fromCmd && toCmd) {
416 const G4int nParEntries = (G4int)fromCmd->GetParameterEntries();
417 for (G4int i = 0; i < nParEntries; ++i) {
418 G4UIparameter* parameter = new G4UIparameter(*(fromCmd->GetParameter(i)));
419 toCmd->SetParameter(parameter);
420 }
421 }
422}
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 424 of file G4VVisCommand.cc.

425{
426 if (fpVisManager) {
427 const G4double halfX = (extent.GetXmax() - extent.GetXmin()) / 2.;
428 const G4double halfY = (extent.GetYmax() - extent.GetYmin()) / 2.;
429 const G4double halfZ = (extent.GetZmax() - extent.GetZmin()) / 2.;
430 if (halfX > 0. && halfY > 0. && halfZ > 0.) {
431 const G4Box box("vis_extent",halfX,halfY,halfZ);
432 const G4VisAttributes visAtts(G4Color::Red());
433 const G4Point3D& centre = extent.GetExtentCenter();
434 fpVisManager->Draw(box,visAtts,G4Translate3D(centre));
435 }
436 }
437}
HepGeom::Point3D< G4double > G4Point3D
Definition G4Point3D.hh:34
HepGeom::Translate3D G4Translate3D
static G4Colour Red()
Definition G4Colour.hh:180
G4double GetYmin() const
G4double GetXmax() const
const G4Point3D & GetExtentCenter() const
G4double GetYmax() const
G4double GetZmax() const
G4double GetZmin() const
G4double GetXmin() const

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

◆ G4VisCommandsSceneAddUnsuccessful()

void G4VVisCommand::G4VisCommandsSceneAddUnsuccessful ( G4VisManager::Verbosity verbosity)
protected

◆ GetCurrentTextColour()

const G4Colour & G4VVisCommand::GetCurrentTextColour ( )
static

Definition at line 77 of file G4VVisCommand.cc.

78{
79 return fCurrentTextColour;
80}
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 328 of file G4VVisCommand.cc.

335{
336 std::vector<G4ViewParameters> viewVector;
337 viewVector.push_back(oldVP);
338 viewVector.push_back(oldVP);
339 viewVector.push_back(newVP);
340 viewVector.push_back(newVP);
341
343 (currentViewer,
344 viewVector,
345 nInterpolationPoints,
346 waitTimePerPointmilliseconds,
347 exportString);
348}
void InterpolateViews(G4VViewer *currentViewer, const 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,
const std::vector< G4ViewParameters > & viewVector,
const G4int nInterpolationPoints = 50,
const G4int waitTimePerPointmilliseconds = 20,
const G4String & exportString = "" )
protected

Definition at line 295 of file G4VVisCommand.cc.

301{
302 const G4int safety = (G4int)viewVector.size()*nInterpolationPoints;
303 G4int safetyCount = 0;
304 do {
305 // Interpret wait time as a minimum time per interation, i.e., a steady frame rate.
306 auto a = std::chrono::steady_clock::now();
307 G4ViewParameters* vp =
308 G4ViewParameters::CatmullRomCubicSplineInterpolation(viewVector,nInterpolationPoints);
309 if (!vp) break; // Finished.
310 currentViewer->SetViewParameters(*vp);
311 currentViewer->RefreshView();
312 if (exportString == "export" &&
313 G4StrUtil::contains(currentViewer->GetName(), "OpenGL")) {
314 G4UImanager::GetUIpointer()->ApplyCommand("/vis/ogl/export");
315 }
316 currentViewer->ShowView();
317 if (waitTimePerPointmilliseconds > 0) {
318 auto b = std::chrono::steady_clock::now();
319 auto timeTaken = b - a;
320 auto timeLeft = std::chrono::milliseconds(waitTimePerPointmilliseconds) - timeTaken;
321 if (timeLeft > std::chrono::duration<double>::zero()) {
322 std::this_thread::sleep_for(timeLeft);
323 }
324 }
325 } while (safetyCount++ < safety); // Loop checking, 16.02.2016, J.Allison
326}
G4int ApplyCommand(const char *aCommand)
const G4String & GetName() const
void SetViewParameters(const G4ViewParameters &vp)
Definition G4VViewer.cc:216
void RefreshView()
virtual void ShowView()
Definition G4VViewer.cc:110
static G4ViewParameters * CatmullRomCubicSplineInterpolation(const std::vector< G4ViewParameters > &views, G4int nInterpolationPoints=50)
G4bool contains(const G4String &str, std::string_view ss)
Check if a string contains a given substring.

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

◆ ProvideValueOfUnit()

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

Definition at line 182 of file G4VVisCommand.cc.

187{
188 // Return false if there's a problem
189
190 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
191
192 G4bool success = true;
194 if (verbosity >= G4VisManager::warnings) {
195 G4warn << where
196 << "\n Unit \"" << unit << "\" not defined"
197 << G4endl;
198 }
199 success = false;
200 } else if (G4UnitDefinition::GetCategory(unit) != category) {
201 if (verbosity >= G4VisManager::warnings) {
202 G4warn << where
203 << "\n Unit \"" << unit << "\" not a unit of " << category;
204 if (category == "Volumic Mass") G4warn << " (density)";
205 G4warn << G4endl;
206 }
207 success = false;
208 } else {
209 value = G4UnitDefinition::GetValueOf(unit);
210 }
211 return success;
212}
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 278 of file G4VVisCommand.cc.

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

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

◆ SetViewParameters()

◆ SetVisManager()

void G4VVisCommand::SetVisManager ( G4VisManager * pVisManager)
static

Definition at line 72 of file G4VVisCommand.cc.

73{
74 fpVisManager = pVisManager;
75}

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 350 of file G4VVisCommand.cc.

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

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

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(), G4VisCommandsViewerSet::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(), G4VisCommandModelCreate< Factory >::SetNewValue(), G4VisCommandMultithreadingActionOnEventQueueFull::SetNewValue(), G4VisCommandMultithreadingMaxEventQueueSize::SetNewValue(), G4VisCommandOpen::SetNewValue(), G4VisCommandPlot::SetNewValue(), G4VisCommandReviewKeptEvents::SetNewValue(), G4VisCommandReviewPlots::SetNewValue(), G4VisCommandSceneActivateModel::SetNewValue(), G4VisCommandSceneAddArrow2D::SetNewValue(), G4VisCommandSceneAddArrow::SetNewValue(), G4VisCommandSceneAddAxes::SetNewValue(), G4VisCommandSceneAddDate::SetNewValue(), G4VisCommandSceneAddDigis::SetNewValue(), G4VisCommandSceneAddElectricField::SetNewValue(), G4VisCommandSceneAddEndOfRunMacro::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: