Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VisCommandsSceneAdd.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// /vis/scene/add commands - John Allison 9th August 1998
28
30
36#include "G4HitsModel.hh"
37#include "G4DigiModel.hh"
38#include "G4GPSModel.hh"
41#include "G4PSHitsModel.hh"
43#include "G4TextModel.hh"
44#include "G4ArrowModel.hh"
45#include "G4AxesModel.hh"
46#include "G4PlotterModel.hh"
48#include "G4ParticleTable.hh"
50#include "G4ApplicationState.hh"
51#include "G4VUserVisAction.hh"
52#include "G4CallbackModel.hh"
53#include "G4UnionSolid.hh"
54#include "G4SubtractionSolid.hh"
55#include "G4Polyhedron.hh"
56#include "G4UImanager.hh"
57#include "G4UIcommandTree.hh"
58#include "G4UIcommand.hh"
59#include "G4UIcmdWithAString.hh"
62#include "G4Tokenizer.hh"
63#include "G4RunManager.hh"
65#include "G4StateManager.hh"
66#include "G4Run.hh"
67#include "G4Event.hh"
68#include "G4Trajectory.hh"
69#include "G4TrajectoryPoint.hh"
70#include "G4RichTrajectory.hh"
72#include "G4SmoothTrajectory.hh"
74#include "G4AttDef.hh"
75#include "G4AttCheck.hh"
76#include "G4Polyline.hh"
77#include "G4UnitsTable.hh"
79#include "G4SystemOfUnits.hh"
81#include "G4PlotterManager.hh"
82
83#include <sstream>
84
85#define G4warn G4cout
86
87////////////// /vis/scene/add/arrow ///////////////////////////////////////
88
90 fpCommand = new G4UIcommand("/vis/scene/add/arrow", this);
91 fpCommand -> SetGuidance ("Adds arrow to current scene.");
92 G4bool omitable;
93 G4UIparameter* parameter;
94 parameter = new G4UIparameter ("x1", 'd', omitable = false);
95 fpCommand -> SetParameter (parameter);
96 parameter = new G4UIparameter ("y1", 'd', omitable = false);
97 fpCommand -> SetParameter (parameter);
98 parameter = new G4UIparameter ("z1", 'd', omitable = false);
99 fpCommand -> SetParameter (parameter);
100 parameter = new G4UIparameter ("x2", 'd', omitable = false);
101 fpCommand -> SetParameter (parameter);
102 parameter = new G4UIparameter ("y2", 'd', omitable = false);
103 fpCommand -> SetParameter (parameter);
104 parameter = new G4UIparameter ("z2", 'd', omitable = false);
105 fpCommand -> SetParameter (parameter);
106 parameter = new G4UIparameter ("unit", 's', omitable = true);
107 parameter->SetDefaultValue ("m");
108 fpCommand->SetParameter (parameter);
109}
110
114
118
120{
121 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
122 G4bool warn(verbosity >= G4VisManager::warnings);
123
124 G4Scene* pScene = fpVisManager->GetCurrentScene();
125 if (!pScene) {
126 if (verbosity >= G4VisManager::errors) {
127 G4warn << "ERROR: No current scene. Please create one." << G4endl;
128 }
129 return;
130 }
131
132 G4String unitString;
133 G4double x1, y1, z1, x2, y2, z2;
134 std::istringstream is(newValue);
135 is >> x1 >> y1 >> z1 >> x2 >> y2 >> z2 >> unitString;
136 G4double unit = G4UIcommand::ValueOf(unitString);
137 x1 *= unit; y1 *= unit; z1 *= unit;
138 x2 *= unit; y2 *= unit; z2 *= unit;
139
140 // Consult scene for arrow width.
141 const G4VisExtent& sceneExtent = pScene->GetExtent();
142 G4double arrowWidth =
143 0.005 * fCurrentLineWidth * sceneExtent.GetExtentRadius();
144
145 G4VModel* model = new G4ArrowModel
146 (x1, y1, z1, x2, y2, z2,
147 arrowWidth, fCurrentColour, newValue,
149
150 const G4String& currentSceneName = pScene -> GetName ();
151 G4bool successful = pScene -> AddRunDurationModel (model, warn);
152 if (successful) {
153 if (verbosity >= G4VisManager::confirmations) {
154 G4cout << "Arrow has been added to scene \""
155 << currentSceneName << "\"."
156 << G4endl;
157 }
158 }
159 else G4VisCommandsSceneAddUnsuccessful(verbosity);
160
162}
163
164////////////// /vis/scene/add/arrow2D ///////////////////////////////////////
165
167 fpCommand = new G4UIcommand("/vis/scene/add/arrow2D", this);
168 fpCommand -> SetGuidance ("Adds 2D arrow to current scene.");
169 fpCommand -> SetGuidance ("x,y in range [-1,1]");
170 G4bool omitable;
171 G4UIparameter* parameter;
172 parameter = new G4UIparameter ("x1", 'd', omitable = false);
173 fpCommand -> SetParameter (parameter);
174 parameter = new G4UIparameter ("y1", 'd', omitable = false);
175 fpCommand -> SetParameter (parameter);
176 parameter = new G4UIparameter ("x2", 'd', omitable = false);
177 fpCommand -> SetParameter (parameter);
178 parameter = new G4UIparameter ("y2", 'd', omitable = false);
179 fpCommand -> SetParameter (parameter);
180}
181
185
189
191{
192 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
193 G4bool warn(verbosity >= G4VisManager::warnings);
194
195 G4Scene* pScene = fpVisManager->GetCurrentScene();
196 if (!pScene) {
197 if (verbosity >= G4VisManager::errors) {
198 G4warn << "ERROR: No current scene. Please create one." << G4endl;
199 }
200 return;
201 }
202
203 G4double x1, y1, x2, y2;
204 std::istringstream is(newValue);
205 is >> x1 >> y1 >> x2 >> y2;
206
207 Arrow2D* arrow2D = new Arrow2D
208 (x1, y1, x2, y2, fCurrentLineWidth, fCurrentColour);
209 G4VModel* model =
211 model->SetType("Arrow2D");
212 model->SetGlobalTag("Arrow2D");
213 model->SetGlobalDescription("Arrow2D: " + newValue);
214 const G4String& currentSceneName = pScene -> GetName ();
215 G4bool successful = pScene -> AddRunDurationModel (model, warn);
216 if (successful) {
217 if (verbosity >= G4VisManager::confirmations) {
218 G4cout << "A 2D arrow has been added to scene \""
219 << currentSceneName << "\"."
220 << G4endl;
221 }
222 }
223 else G4VisCommandsSceneAddUnsuccessful(verbosity);
224
226}
227
228G4VisCommandSceneAddArrow2D::Arrow2D::Arrow2D
229(G4double x1, G4double y1,
230 G4double x2, G4double y2,
231 G4double width, const G4Colour& colour):
232 fWidth(width), fColour(colour)
233{
234 fShaftPolyline.push_back(G4Point3D(x1,y1,0));
235 fShaftPolyline.push_back(G4Point3D(x2,y2,0));
236 G4Vector3D arrowDirection = G4Vector3D(x2-x1,y2-y1,0).unit();
237 G4Vector3D arrowPointLeftDirection(arrowDirection);
238 arrowPointLeftDirection.rotateZ(150.*deg);
239 G4Vector3D arrowPointRightDirection(arrowDirection);
240 arrowPointRightDirection.rotateZ(-150.*deg);
241 fHeadPolyline.push_back(G4Point3D(x2,y2,0)+0.04*arrowPointLeftDirection);
242 fHeadPolyline.push_back(G4Point3D(x2,y2,0));
243 fHeadPolyline.push_back(G4Point3D(x2,y2,0)+0.04*arrowPointRightDirection);
245 va.SetLineWidth(fWidth);
246 va.SetColour(fColour);
247 fShaftPolyline.SetVisAttributes(va);
248 fHeadPolyline.SetVisAttributes(va);
249}
250
251void G4VisCommandSceneAddArrow2D::Arrow2D::operator()
252 (G4VGraphicsScene& sceneHandler, const G4ModelingParameters*)
253{
254 sceneHandler.BeginPrimitives2D();
255 sceneHandler.AddPrimitive(fShaftPolyline);
256 sceneHandler.AddPrimitive(fHeadPolyline);
257 sceneHandler.EndPrimitives2D();
258}
259
260////////////// /vis/scene/add/axes //////////////////////////////////
261
263 G4bool omitable;
264 fpCommand = new G4UIcommand ("/vis/scene/add/axes", this);
265 fpCommand -> SetGuidance ("Add axes.");
266 fpCommand -> SetGuidance
267 ("Draws axes at (x0, y0, z0) of given length and colour.");
268 fpCommand -> SetGuidance
269 ("If \"colour-string\" is \"auto\", x, y and z will be red, green and blue"
270 "\n respectively. Otherwise it can be one of the pre-defined text-specified"
271 "\n colours - see information printed by the vis manager at start-up or"
272 "\n use \"/vis/list\".");
273 fpCommand -> SetGuidance
274 ("If \"length\" is negative, it is set to about 25% of scene extent.");
275 fpCommand -> SetGuidance
276 ("If \"showtext\" is false, annotations are suppressed.");
277 G4UIparameter* parameter;
278 parameter = new G4UIparameter ("x0", 'd', omitable = true);
279 parameter->SetDefaultValue (0.);
280 fpCommand->SetParameter (parameter);
281 parameter = new G4UIparameter ("y0", 'd', omitable = true);
282 parameter->SetDefaultValue (0.);
283 fpCommand->SetParameter (parameter);
284 parameter = new G4UIparameter ("z0", 'd', omitable = true);
285 parameter->SetDefaultValue (0.);
286 fpCommand->SetParameter (parameter);
287 parameter = new G4UIparameter ("length", 'd', omitable = true);
288 parameter->SetDefaultValue (-1.);
289 fpCommand->SetParameter (parameter);
290 parameter = new G4UIparameter ("unit", 's', omitable = true);
291 parameter->SetDefaultValue ("m");
292 fpCommand->SetParameter (parameter);
293 parameter = new G4UIparameter ("colour-string", 's', omitable = true);
294 parameter->SetDefaultValue ("auto");
295 fpCommand->SetParameter (parameter);
296 parameter = new G4UIparameter ("showtext", 'b', omitable = true);
297 parameter->SetDefaultValue ("true");
298 fpCommand->SetParameter (parameter);
299}
300
304
308
310
311 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
312 G4bool warn(verbosity >= G4VisManager::warnings);
313
314 G4Scene* pScene = fpVisManager->GetCurrentScene();
315 if (!pScene) {
316 if (verbosity >= G4VisManager::errors) {
317 G4warn << "ERROR: No current scene. Please create one." << G4endl;
318 }
319 return;
320 } else {
321 if (pScene->GetExtent().GetExtentRadius() <= 0.) {
322 if (verbosity >= G4VisManager::errors) {
323 G4warn
324 << "ERROR: Scene has no extent. Add volumes or use \"/vis/scene/add/extent\"."
325 << G4endl;
326 }
327 return;
328 }
329 }
330
331 G4String unitString, colourString, showTextString;
332 G4double x0, y0, z0, length;
333 std::istringstream is (newValue);
334 is >> x0 >> y0 >> z0 >> length >> unitString
335 >> colourString >> showTextString;
336 G4bool showText = G4UIcommand::ConvertToBool(showTextString);
337
338
339 G4double unit = G4UIcommand::ValueOf(unitString);
340 x0 *= unit; y0 *= unit; z0 *= unit;
341 const G4VisExtent& sceneExtent = pScene->GetExtent(); // Existing extent.
342 if (length < 0.) {
343 const G4double lengthMax = 0.5 * sceneExtent.GetExtentRadius();
344 const G4double intLog10Length = std::floor(std::log10(lengthMax));
345 length = std::pow(10,intLog10Length);
346 if (5.*length < lengthMax) length *= 5.;
347 else if (2.*length < lengthMax) length *= 2.;
348 } else {
349 length *= unit;
350 }
351
352 // Consult scene for arrow width...
353 G4double arrowWidth =
354 0.05 * fCurrentLineWidth * sceneExtent.GetExtentRadius();
355 // ...but limit it to length/30.
356 if (arrowWidth > length/30.) arrowWidth = length/30.;
357
358 G4VModel* model = new G4AxesModel
359 (x0, y0, z0, length, arrowWidth, colourString, newValue,
360 showText, fCurrentTextSize);
361
362 G4bool successful = pScene -> AddRunDurationModel (model, warn);
363 const G4String& currentSceneName = pScene -> GetName ();
364 if (successful) {
365 if (verbosity >= G4VisManager::confirmations) {
366 G4cout << "Axes of length " << G4BestUnit(length,"Length")
367 << "have been added to scene \"" << currentSceneName << "\"."
368 << G4endl;
369 }
370 }
371 else G4VisCommandsSceneAddUnsuccessful(verbosity);
372
374}
375
376////////////// /vis/scene/add/date ///////////////////////////////////////
377
379 G4bool omitable;
380 fpCommand = new G4UIcommand ("/vis/scene/add/date", this);
381 fpCommand -> SetGuidance ("Adds date to current scene.");
382 fpCommand -> SetGuidance
383 ("If \"date\"is omitted, the current date and time is drawn."
384 "\nOtherwise, the string, including the rest of the line, is drawn.");
385 G4UIparameter* parameter;
386 parameter = new G4UIparameter ("size", 'i', omitable = true);
387 parameter -> SetGuidance ("Screen size of text in pixels.");
388 parameter -> SetDefaultValue (18);
389 fpCommand -> SetParameter (parameter);
390 parameter = new G4UIparameter ("x-position", 'd', omitable = true);
391 parameter -> SetGuidance ("x screen position in range -1 < x < 1.");
392 parameter -> SetDefaultValue (0.95);
393 fpCommand -> SetParameter (parameter);
394 parameter = new G4UIparameter ("y-position", 'd', omitable = true);
395 parameter -> SetGuidance ("y screen position in range -1 < y < 1.");
396 parameter -> SetDefaultValue (0.9);
397 fpCommand -> SetParameter (parameter);
398 parameter = new G4UIparameter ("layout", 's', omitable = true);
399 parameter -> SetGuidance ("Layout, i.e., adjustment: left|centre|right.");
400 parameter -> SetDefaultValue ("right");
401 fpCommand -> SetParameter (parameter);
402 parameter = new G4UIparameter ("date", 's', omitable = true);
403 parameter -> SetDefaultValue ("-");
404 fpCommand -> SetParameter (parameter);
405}
406
410
414
416{
417 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
418 G4bool warn(verbosity >= G4VisManager::warnings);
419
420 G4Scene* pScene = fpVisManager->GetCurrentScene();
421 if (!pScene) {
422 if (verbosity >= G4VisManager::errors) {
423 G4warn << "ERROR: No current scene. Please create one." << G4endl;
424 }
425 return;
426 }
427
428 G4int size;
429 G4double x, y;
430 G4String layoutString, dateString;
431 std::istringstream is(newValue);
432 is >> size >> x >> y >> layoutString >> dateString;
433 // Read rest of line, if any.
434 const size_t NREMAINDER = 100;
435 char remainder[NREMAINDER];
436 remainder[0]='\0'; // In case there is nothing remaining.
437 is.getline(remainder, NREMAINDER);
438 dateString += remainder;
440 if (layoutString[0] == 'l') layout = G4Text::left;
441 else if (layoutString[0] == 'c') layout = G4Text::centre;
442 else if (layoutString[0] == 'r') layout = G4Text::right;
443
444 Date* date = new Date(fpVisManager, size, x, y, layout, dateString);
445 G4VModel* model =
447 model->SetType("Date");
448 model->SetGlobalTag("Date");
449 model->SetGlobalDescription("Date: " + newValue);
450 const G4String& currentSceneName = pScene -> GetName ();
451 G4bool successful = pScene -> AddEndOfEventModel(model, warn);
452 if (successful) {
453 if (verbosity >= G4VisManager::confirmations) {
454 G4cout << "Date has been added to scene \""
455 << currentSceneName << "\"."
456 << G4endl;
457 }
458 }
459 else G4VisCommandsSceneAddUnsuccessful(verbosity);
460
462}
463
464void G4VisCommandSceneAddDate::Date::operator()
465 (G4VGraphicsScene& sceneHandler, const G4ModelingParameters*)
466{
467 G4String time;
468 if (fDate == "-") {
469 time = fTimer.GetClockTime();
470 } else {
471 time = fDate;
472 }
473 // Check for \n, starting from back, and erase.
474 std::string::size_type i = time.rfind('\n');
475 if (i != std::string::npos) time.erase(i);
476 G4Text text(time, G4Point3D(fX, fY, 0.));
477 text.SetScreenSize(fSize);
478 text.SetLayout(fLayout);
479 G4VisAttributes textAtts(G4Colour(0.,1.,1));
480 text.SetVisAttributes(textAtts);
481 sceneHandler.BeginPrimitives2D();
482 sceneHandler.AddPrimitive(text);
483 sceneHandler.EndPrimitives2D();
484}
485
486////////////// /vis/scene/add/digis ///////////////////////////////////////
487
489 fpCommand = new G4UIcmdWithoutParameter ("/vis/scene/add/digis", this);
490 fpCommand -> SetGuidance ("Adds digis to current scene.");
491 fpCommand -> SetGuidance
492 ("Digis are drawn at end of event when the scene in which"
493 "\nthey are added is current.");
494}
495
499
503
505
506 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
507 G4bool warn(verbosity >= G4VisManager::warnings);
508
509 G4Scene* pScene = fpVisManager->GetCurrentScene();
510 if (!pScene) {
511 if (verbosity >= G4VisManager::errors) {
512 G4warn << "ERROR: No current scene. Please create one." << G4endl;
513 }
514 return;
515 }
516
517 G4VModel* model = new G4DigiModel;
518 const G4String& currentSceneName = pScene -> GetName ();
519 G4bool successful = pScene -> AddEndOfEventModel (model, warn);
520 if (successful) {
521 if (verbosity >= G4VisManager::confirmations) {
522 G4cout << "Digis, if any, will be drawn at end of run in scene \""
523 << currentSceneName << "\"."
524 << G4endl;
525 }
526 }
527 else G4VisCommandsSceneAddUnsuccessful(verbosity);
528
530}
531
532////////////// /vis/scene/add/electricField ///////////////////////////////////////
533
535 G4bool omitable;
536 fpCommand = new G4UIcommand ("/vis/scene/add/electricField", this);
537 fpCommand -> SetGuidance
538 ("Adds electric field representation to current scene.");
539 fpCommand -> SetGuidance
540 ("The first parameter is no. of data points per half extent. So, possibly, at"
541 "\nmaximum, the number of data points sampled is (2*n+1)^3, which can grow"
542 "\nlarge--be warned!"
543 "\nThe default value is 10, i.e., a 21x21x21 array, i.e., 9,261 sampling points."
544 "\nThat may swamp your view, but usually, a field is limited to a small part of"
545 "\nthe extent, so it's not a problem. But if it is, here are some of the things"
546 "\nyou can do:"
547 "\n- reduce the number of data points per half extent (first parameter);"
548 "\n- specify \"lightArrow\" (second parameter);"
549 "\n- restrict the region sampled with \"/vis/set/extentForField\";"
550 "\n- restrict the drawing to a specific volume with"
551 "\n \"/vis/set/volumeForField\" or \"/vis/touchable/volumeForField\"."
552 "\nNote: you might have to deactivate existing field models with"
553 "\n \"/vis/scene/activateModel Field false\" and re-issue"
554 "\n \"/vis/scene/add/...Field\" command again.");
555 fpCommand -> SetGuidance
556 ("In the arrow representation, the length of the arrow is proportional"
557 "\nto the magnitude of the field and the colour is mapped onto the range"
558 "\nas a fraction of the maximum magnitude: 0->0.5->1 is red->green->blue.");
559 G4UIparameter* parameter;
560 parameter = new G4UIparameter ("nDataPointsPerHalfExtent", 'i', omitable = true);
561 parameter -> SetDefaultValue (10);
562 fpCommand -> SetParameter (parameter);
563 parameter = new G4UIparameter ("representation", 's', omitable = true);
564 parameter -> SetParameterCandidates("fullArrow lightArrow");
565 parameter -> SetDefaultValue ("fullArrow");
566 fpCommand -> SetParameter (parameter);
567}
568
572
576
578(G4UIcommand*, G4String newValue) {
579
580 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
581 G4bool warn(verbosity >= G4VisManager::warnings);
582
583 G4Scene* pScene = fpVisManager->GetCurrentScene();
584 if (!pScene) {
585 if (verbosity >= G4VisManager::errors) {
586 G4warn << "ERROR: No current scene. Please create one." << G4endl;
587 }
588 return;
589 }
590
591 G4int nDataPointsPerHalfExtent;
592 G4String representation;
593 std::istringstream iss(newValue);
594 iss >> nDataPointsPerHalfExtent >> representation;
596 modelRepresentation = G4ElectricFieldModel::fullArrow;
597 if (representation == "lightArrow") {
598 modelRepresentation = G4ElectricFieldModel::lightArrow;
599 }
600 G4VModel* model;
601 model = new G4ElectricFieldModel
602 (nDataPointsPerHalfExtent,modelRepresentation,
606 const G4String& currentSceneName = pScene -> GetName ();
607 G4bool successful = pScene -> AddRunDurationModel (model, warn);
608 if (successful) {
609 if (verbosity >= G4VisManager::confirmations) {
610 G4cout
611 << "Electric field, if any, will be drawn in scene \""
612 << currentSceneName
613 << "\"\n with "
614 << nDataPointsPerHalfExtent
615 << " data points per half extent and with representation \""
616 << representation
617 << '\"'
618 << G4endl;
619 }
620 }
621 else G4VisCommandsSceneAddUnsuccessful(verbosity);
622
624}
625
626////////////// /vis/scene/add/endOfRunMacro ///////////////////////////////////////
627
629 G4bool omitable;
630 fpCommand = new G4UIcmdWithAString ("/vis/scene/add/endOfRunMacro", this);
631 fpCommand -> SetGuidance ("Macro is executed at end of run and when rebuild required.");
632 fpCommand -> SetGuidance
633 ("WARNING: some vis commands in the macro cause recursion."
634 "\n Stick to simple commmands, e.g., which invoke vis manager Draw() methods.");
635 fpCommand -> SetParameterName ("macro", omitable = false);
636}
637
641
645
647{
648 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
649 G4bool warn(verbosity >= G4VisManager::warnings);
650
651 G4Scene* pScene = fpVisManager->GetCurrentScene();
652 if (!pScene) {
653 if (verbosity >= G4VisManager::errors) {
654 G4warn << "ERROR: No current scene. Please create one." << G4endl;
655 }
656 return;
657 }
658
659 auto endOfRunMacro = new EndOfRunMacro(newValue);
660 G4VModel* model =
662 model->SetType("EndOfRunMacro");
663 model->SetGlobalTag("EndOfRunMacro");
664 model->SetGlobalDescription("EndOfRunMacro: " + newValue);
665 const G4String& currentSceneName = pScene -> GetName ();
666 G4bool successful = pScene -> AddEndOfRunModel(model, warn);
667 if (successful) {
668 if (verbosity >= G4VisManager::confirmations) {
669 G4cout << "EndOfRunMacro has been added to scene \""
670 << currentSceneName << "\"."
671 << G4endl;
672 }
673 }
674 else G4VisCommandsSceneAddUnsuccessful(verbosity);
675
677}
678
679void G4VisCommandSceneAddEndOfRunMacro::EndOfRunMacro::operator()
681{
682 G4UImanager::GetUIpointer()->ApplyCommand("/control/execute " + fMacro);
683}
684
685////////////// /vis/scene/add/eventID ///////////////////////////////////////
686
688 G4bool omitable;
689 fpCommand = new G4UIcommand ("/vis/scene/add/eventID", this);
690 fpCommand -> SetGuidance ("Adds eventID to current scene.");
691 fpCommand -> SetGuidance
692 ("Run and event numbers are drawn at end of event or run when"
693 "\n the scene in which they are added is current.");
694 G4UIparameter* parameter;
695 parameter = new G4UIparameter ("size", 'i', omitable = true);
696 parameter -> SetGuidance ("Screen size of text in pixels.");
697 parameter -> SetDefaultValue (18);
698 fpCommand -> SetParameter (parameter);
699 parameter = new G4UIparameter ("x-position", 'd', omitable = true);
700 parameter -> SetGuidance ("x screen position in range -1 < x < 1.");
701 parameter -> SetDefaultValue (-0.95);
702 fpCommand -> SetParameter (parameter);
703 parameter = new G4UIparameter ("y-position", 'd', omitable = true);
704 parameter -> SetGuidance ("y screen position in range -1 < y < 1.");
705 parameter -> SetDefaultValue (0.9);
706 fpCommand -> SetParameter (parameter);
707 parameter = new G4UIparameter ("layout", 's', omitable = true);
708 parameter -> SetGuidance ("Layout, i.e., adjustment: left|centre|right.");
709 parameter -> SetDefaultValue ("left");
710 fpCommand -> SetParameter (parameter);
711}
712
716
720
722{
723 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
724 G4bool warn(verbosity >= G4VisManager::warnings);
725
726 G4Scene* pScene = fpVisManager->GetCurrentScene();
727 if (!pScene) {
728 if (verbosity >= G4VisManager::errors) {
729 G4warn << "ERROR: No current scene. Please create one." << G4endl;
730 }
731 return;
732 }
733
734 G4int size;
735 G4double x, y;
736 G4String layoutString;
737 std::istringstream is(newValue);
738 is >> size >> x >> y >> layoutString;
739
741 if (layoutString[0] == 'l') layout = G4Text::left;
742 else if (layoutString[0] == 'c') layout = G4Text::centre;
743 else if (layoutString[0] == 'r') layout = G4Text::right;
744
745 // For End of Event (only for reviewing kept events one by one)
746 EventID* eoeEventID
747 = new EventID(forEndOfEvent, fpVisManager, size, x, y, layout);
748 G4VModel* eoeModel =
750 eoeModel->SetType("EoEEventID");
751 eoeModel->SetGlobalTag("EoEEventID");
752 eoeModel->SetGlobalDescription("EoEEventID: " + newValue);
753 G4bool successfulEoE = pScene -> AddEndOfEventModel (eoeModel, warn);
754
755 // For End of Run
756 EventID* eorEventID
757 = new EventID(forEndOfRun, fpVisManager, size, x, y, layout);
758 G4VModel* eorModel =
760 eorModel->SetType("EoREventID");
761 eorModel->SetGlobalTag("EoREventID");
762 eorModel->SetGlobalDescription("EoREventID: " + newValue);
763 G4bool successfulEoR = pScene -> AddEndOfRunModel (eorModel, warn);
764
765 if (successfulEoE && successfulEoR) {
766 if (verbosity >= G4VisManager::confirmations) {
767 const G4String& currentSceneName = pScene -> GetName ();
768 G4cout << "EventID has been added to scene \""
769 << currentSceneName << "\"."
770 << G4endl;
771 }
772 }
773 else G4VisCommandsSceneAddUnsuccessful(verbosity);
774
776}
777
778void G4VisCommandSceneAddEventID::EventID::operator()
779(G4VGraphicsScene& sceneHandler, const G4ModelingParameters* mp)
780{
782 if(!runManager)
783 return;
784
785 const G4Run* currentRun = runManager->GetCurrentRun();
786 if (!currentRun) return;
787
788 const G4int currentRunID = currentRun->GetRunID();
789
790 std::ostringstream oss;
791 switch (fForWhat) {
792 case forEndOfEvent:
793 {
794 // Only use if reviewing kept events
795 if (!fpVisManager->GetReviewingKeptEvents()) return;
796 const G4Event* currentEvent = mp->GetEvent();
797 if (!currentEvent) return;
798 G4int eventID = currentEvent->GetEventID();
799 oss << "Run " << currentRunID << " Event " << eventID;
800 break;
801 }
802 case forEndOfRun:
803 {
804 // Only use if NOT reviewing kept events
805 if (fpVisManager->GetReviewingKeptEvents()) return;
806 const G4int nEvents = currentRun->GetNumberOfEventToBeProcessed();
807 size_t nKeptEvents = (size_t)(currentRun->GetNumberOfKeptEvents());
808 oss << "Run " << currentRunID << " (" << nEvents << " event";
809 if (nEvents != 1) oss << 's';
810 oss << ", " << nKeptEvents << " kept)";
811 break;
812 }
813 default:
814 return;
815 }
816
817 G4Text text(oss.str(), G4Point3D(fX, fY, 0.));
818 text.SetScreenSize(fSize);
819 text.SetLayout(fLayout);
820 G4VisAttributes textAtts(G4Colour(0.,1.,1));
821 text.SetVisAttributes(textAtts);
822 sceneHandler.BeginPrimitives2D();
823 sceneHandler.AddPrimitive(text);
824 sceneHandler.EndPrimitives2D();
825}
826
827////////////// /vis/scene/add/extent ///////////////////////////////////////
828
830 fpCommand = new G4UIcommand("/vis/scene/add/extent", this);
831 fpCommand -> SetGuidance
832 ("Adds a dummy model with given extent to the current scene."
833 "\nRequires the limits: xmin, xmax, ymin, ymax, zmin, zmax unit."
834 "\nThis can be used to provide an extent to the scene even if"
835 "\nno other models with extent are available. For example,"
836 "\neven if there is no geometry. In that case, for example:"
837 "\n /vis/open OGL"
838 "\n /vis/scene/create"
839 "\n /vis/scene/add/extent -300 300 -300 300 -300 300 cm"
840 "\n /vis/sceneHandler/attach");
841 G4bool omitable;
842 G4UIparameter* parameter;
843 parameter = new G4UIparameter ("xmin", 'd', omitable = true);
844 parameter -> SetDefaultValue (0.);
845 fpCommand -> SetParameter (parameter);
846 parameter = new G4UIparameter ("xmax", 'd', omitable = true);
847 parameter -> SetDefaultValue (0.);
848 fpCommand -> SetParameter (parameter);
849 parameter = new G4UIparameter ("ymin", 'd', omitable = true);
850 parameter -> SetDefaultValue (0.);
851 fpCommand -> SetParameter (parameter);
852 parameter = new G4UIparameter ("ymax", 'd', omitable = true);
853 parameter -> SetDefaultValue (0.);
854 fpCommand -> SetParameter (parameter);
855 parameter = new G4UIparameter ("zmin", 'd', omitable = true);
856 parameter -> SetDefaultValue (0.);
857 fpCommand -> SetParameter (parameter);
858 parameter = new G4UIparameter ("zmax", 'd', omitable = true);
859 parameter -> SetDefaultValue (0.);
860 fpCommand -> SetParameter (parameter);
861 parameter = new G4UIparameter ("unit", 's', omitable = true);
862 parameter -> SetDefaultValue ("m");
863 fpCommand -> SetParameter (parameter);
864}
865
869
873
875{
876 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
877 G4bool warn(verbosity >= G4VisManager::warnings);
878
879 G4Scene* pScene = fpVisManager->GetCurrentScene();
880 if (!pScene) {
881 if (verbosity >= G4VisManager::errors) {
882 G4warn << "ERROR: No current scene. Please create one." << G4endl;
883 }
884 return;
885 }
886
887 G4double xmin, xmax, ymin, ymax, zmin, zmax;
888 G4String unitString;
889 std::istringstream is(newValue);
890 is >> xmin >> xmax >> ymin >> ymax >> zmin >> zmax >> unitString;
891 G4double unit = G4UIcommand::ValueOf(unitString);
892 xmin *= unit; xmax *= unit;
893 ymin *= unit; ymax *= unit;
894 zmin *= unit; zmax *= unit;
895
896 G4VisExtent visExtent(xmin, xmax, ymin, ymax, zmin, zmax);
897 Extent* extent = new Extent(xmin, xmax, ymin, ymax, zmin, zmax);
898 G4VModel* model =
900 model->SetType("Extent");
901 model->SetGlobalTag("Extent");
902 model->SetGlobalDescription("Extent: " + newValue);
903 model->SetExtent(visExtent);
904 const G4String& currentSceneName = pScene -> GetName ();
905 G4bool successful = pScene -> AddRunDurationModel (model, warn);
906 if (successful) {
907 if (verbosity >= G4VisManager::confirmations) {
908 G4cout << "A benign model with extent "
909 << visExtent
910 << " has been added to scene \""
911 << currentSceneName << "\"."
912 << G4endl;
913 }
914 }
915 else G4VisCommandsSceneAddUnsuccessful(verbosity);
916
918}
919
920G4VisCommandSceneAddExtent::Extent::Extent
921(G4double xmin, G4double xmax,
922 G4double ymin, G4double ymax,
923 G4double zmin, G4double zmax):
924fExtent(xmin,xmax,ymin,ymax,zmin,zmax)
925{}
926
927void G4VisCommandSceneAddExtent::Extent::operator()
929{}
930
931////////////// /vis/scene/add/frame ///////////////////////////////////////
932
934 fpCommand = new G4UIcommand("/vis/scene/add/frame", this);
935 fpCommand -> SetGuidance ("Add frame to current scene.");
936 G4bool omitable;
937 G4UIparameter* parameter;
938 parameter = new G4UIparameter ("size", 'd', omitable = true);
939 parameter -> SetGuidance ("Size of frame. 1 = full window.");
940 parameter -> SetParameterRange ("size > 0 && size <=1");
941 parameter -> SetDefaultValue (0.97);
942 fpCommand -> SetParameter (parameter);
943}
944
948
952
954{
955 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
956 G4bool warn(verbosity >= G4VisManager::warnings);
957
958 G4Scene* pScene = fpVisManager->GetCurrentScene();
959 if (!pScene) {
960 if (verbosity >= G4VisManager::errors) {
961 G4warn << "ERROR: No current scene. Please create one." << G4endl;
962 }
963 return;
964 }
965
966 G4double size;
967 std::istringstream is(newValue);
968 is >> size;
969
970 Frame* frame = new Frame(size, fCurrentLineWidth, fCurrentColour);
971 G4VModel* model =
973 model->SetType("Frame");
974 model->SetGlobalTag("Frame");
975 model->SetGlobalDescription("Frame: " + newValue);
976 const G4String& currentSceneName = pScene -> GetName ();
977 G4bool successful = pScene -> AddRunDurationModel (model, warn);
978 if (successful) {
979 if (verbosity >= G4VisManager::confirmations) {
980 G4cout << "Frame has been added to scene \""
981 << currentSceneName << "\"."
982 << G4endl;
983 }
984 }
985 else G4VisCommandsSceneAddUnsuccessful(verbosity);
986
988}
989
990void G4VisCommandSceneAddFrame::Frame::operator()
991 (G4VGraphicsScene& sceneHandler, const G4ModelingParameters*)
992{
993 G4Polyline frame;
994 frame.push_back(G4Point3D( fSize, fSize, 0.));
995 frame.push_back(G4Point3D(-fSize, fSize, 0.));
996 frame.push_back(G4Point3D(-fSize, -fSize, 0.));
997 frame.push_back(G4Point3D( fSize, -fSize, 0.));
998 frame.push_back(G4Point3D( fSize, fSize, 0.));
1000 va.SetLineWidth(fWidth);
1001 va.SetColour(fColour);
1002 frame.SetVisAttributes(va);
1003 sceneHandler.BeginPrimitives2D();
1004 sceneHandler.AddPrimitive(frame);
1005 sceneHandler.EndPrimitives2D();
1006}
1007
1008////////////// /vis/scene/add/gps ///////////////////////////////////////
1009
1011 G4bool omitable;
1012 G4UIparameter* parameter;
1013 fpCommand = new G4UIcommand ("/vis/scene/add/gps", this);
1014 fpCommand -> SetGuidance
1015 ("A representation of the source(s) of the General Particle Source"
1016 "\nwill be added to current scene and drawn, if applicable.");
1017 fpCommand->SetGuidance(ConvertToColourGuidance());
1018 fpCommand->SetGuidance("Default: red and transparent.");
1019 parameter = new G4UIparameter("red_or_string", 's', omitable = true);
1020 parameter -> SetDefaultValue ("1.");
1021 fpCommand -> SetParameter (parameter);
1022 parameter = new G4UIparameter("green", 'd', omitable = true);
1023 parameter -> SetDefaultValue (0.);
1024 fpCommand -> SetParameter (parameter);
1025 parameter = new G4UIparameter ("blue", 'd', omitable = true);
1026 parameter -> SetDefaultValue (0.);
1027 fpCommand -> SetParameter (parameter);
1028 parameter = new G4UIparameter ("opacity", 'd', omitable = true);
1029 parameter -> SetDefaultValue (0.3);
1030 fpCommand -> SetParameter (parameter);
1031}
1032
1036
1040
1042
1043 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
1044 G4bool warn(verbosity >= G4VisManager::warnings);
1045
1046 G4Scene* pScene = fpVisManager->GetCurrentScene();
1047 if (!pScene) {
1048 if (verbosity >= G4VisManager::errors) {
1049 G4warn << "ERROR: No current scene. Please create one." << G4endl;
1050 }
1051 return;
1052 }
1053
1054 G4String redOrString;
1055 G4double green, blue, opacity;
1056 std::istringstream iss(newValue);
1057 iss >> redOrString >> green >> blue >> opacity;
1058 G4Colour colour(1.,0.,0.,0.3); // Default red and transparent.
1059 ConvertToColour(colour, redOrString, green, blue, opacity);
1060
1061 G4VModel* model = new G4GPSModel(colour);
1062 const G4String& currentSceneName = pScene -> GetName ();
1063 G4bool successful = pScene -> AddRunDurationModel (model, warn);
1064 if (successful) {
1065 if (verbosity >= G4VisManager::confirmations) {
1066 G4cout <<
1067 "A representation of the source(s) of the General Particle Source will be drawn"
1068 "\n in colour " << colour << " for scene \""
1069 << currentSceneName << "\" if applicable."
1070 << G4endl;
1071 }
1072 }
1073 else G4VisCommandsSceneAddUnsuccessful(verbosity);
1074
1076}
1077
1078////////////// /vis/scene/add/hits ///////////////////////////////////////
1079
1081 fpCommand = new G4UIcmdWithoutParameter ("/vis/scene/add/hits", this);
1082 fpCommand -> SetGuidance ("Adds hits to current scene.");
1083 fpCommand -> SetGuidance
1084 ("Hits are drawn at end of event when the scene in which"
1085 "\nthey are added is current.");
1086}
1087
1091
1095
1097
1098 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
1099 G4bool warn(verbosity >= G4VisManager::warnings);
1100
1101 G4Scene* pScene = fpVisManager->GetCurrentScene();
1102 if (!pScene) {
1103 if (verbosity >= G4VisManager::errors) {
1104 G4warn << "ERROR: No current scene. Please create one." << G4endl;
1105 }
1106 return;
1107 }
1108
1109 G4VModel* model = new G4HitsModel;
1110 const G4String& currentSceneName = pScene -> GetName ();
1111 G4bool successful = pScene -> AddEndOfEventModel (model, warn);
1112 if (successful) {
1113 if (verbosity >= G4VisManager::confirmations) {
1114 G4cout << "Hits, if any, will be drawn at end of run in scene \""
1115 << currentSceneName << "\"."
1116 << G4endl;
1117 }
1118 }
1119 else G4VisCommandsSceneAddUnsuccessful(verbosity);
1120
1122}
1123
1124////////////// /vis/scene/add/line ///////////////////////////////////////
1125
1127 fpCommand = new G4UIcommand("/vis/scene/add/line", this);
1128 fpCommand -> SetGuidance ("Adds line to current scene.");
1129 G4bool omitable;
1130 G4UIparameter* parameter;
1131 parameter = new G4UIparameter ("x1", 'd', omitable = false);
1132 fpCommand -> SetParameter (parameter);
1133 parameter = new G4UIparameter ("y1", 'd', omitable = false);
1134 fpCommand -> SetParameter (parameter);
1135 parameter = new G4UIparameter ("z1", 'd', omitable = false);
1136 fpCommand -> SetParameter (parameter);
1137 parameter = new G4UIparameter ("x2", 'd', omitable = false);
1138 fpCommand -> SetParameter (parameter);
1139 parameter = new G4UIparameter ("y2", 'd', omitable = false);
1140 fpCommand -> SetParameter (parameter);
1141 parameter = new G4UIparameter ("z2", 'd', omitable = false);
1142 fpCommand -> SetParameter (parameter);
1143 parameter = new G4UIparameter ("unit", 's', omitable = true);
1144 parameter->SetDefaultValue ("m");
1145 fpCommand->SetParameter (parameter);
1146}
1147
1151
1155
1157{
1158 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
1159 G4bool warn(verbosity >= G4VisManager::warnings);
1160
1161 G4Scene* pScene = fpVisManager->GetCurrentScene();
1162 if (!pScene) {
1163 if (verbosity >= G4VisManager::errors) {
1164 G4warn << "ERROR: No current scene. Please create one." << G4endl;
1165 }
1166 return;
1167 }
1168
1169 G4String unitString;
1170 G4double x1, y1, z1, x2, y2, z2;
1171 std::istringstream is(newValue);
1172 is >> x1 >> y1 >> z1 >> x2 >> y2 >> z2 >> unitString;
1173 G4double unit = G4UIcommand::ValueOf(unitString);
1174 x1 *= unit; y1 *= unit; z1 *= unit;
1175 x2 *= unit; y2 *= unit; z2 *= unit;
1176
1177 Line* line = new Line(x1, y1, z1, x2, y2, z2,
1179 G4VModel* model =
1181 model->SetType("Line");
1182 model->SetGlobalTag("Line");
1183 model->SetGlobalDescription("Line: " + newValue);
1184 const G4String& currentSceneName = pScene -> GetName ();
1185 G4bool successful = pScene -> AddRunDurationModel (model, warn);
1186 if (successful) {
1187 if (verbosity >= G4VisManager::confirmations) {
1188 G4cout << "Line has been added to scene \""
1189 << currentSceneName << "\"."
1190 << G4endl;
1191 }
1192 }
1193 else G4VisCommandsSceneAddUnsuccessful(verbosity);
1194
1196}
1197
1198G4VisCommandSceneAddLine::Line::Line
1199(G4double x1, G4double y1, G4double z1,
1200 G4double x2, G4double y2, G4double z2,
1201 G4double width, const G4Colour& colour):
1202 fWidth(width), fColour(colour)
1203{
1204 fPolyline.push_back(G4Point3D(x1,y1,z1));
1205 fPolyline.push_back(G4Point3D(x2,y2,z2));
1206 G4VisAttributes va;
1207 va.SetLineWidth(fWidth);
1208 va.SetColour(fColour);
1209 fPolyline.SetVisAttributes(va);
1210}
1211
1212void G4VisCommandSceneAddLine::Line::operator()
1213 (G4VGraphicsScene& sceneHandler, const G4ModelingParameters*)
1214{
1215 sceneHandler.BeginPrimitives();
1216 sceneHandler.AddPrimitive(fPolyline);
1217 sceneHandler.EndPrimitives();
1218}
1219
1220////////////// /vis/scene/add/line2D ///////////////////////////////////////
1221
1223 fpCommand = new G4UIcommand("/vis/scene/add/line2D", this);
1224 fpCommand -> SetGuidance ("Adds 2D line to current scene.");
1225 fpCommand -> SetGuidance ("x,y in range [-1,1]");
1226 G4bool omitable;
1227 G4UIparameter* parameter;
1228 parameter = new G4UIparameter ("x1", 'd', omitable = false);
1229 fpCommand -> SetParameter (parameter);
1230 parameter = new G4UIparameter ("y1", 'd', omitable = false);
1231 fpCommand -> SetParameter (parameter);
1232 parameter = new G4UIparameter ("x2", 'd', omitable = false);
1233 fpCommand -> SetParameter (parameter);
1234 parameter = new G4UIparameter ("y2", 'd', omitable = false);
1235 fpCommand -> SetParameter (parameter);
1236}
1237
1241
1245
1247{
1248 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
1249 G4bool warn(verbosity >= G4VisManager::warnings);
1250
1251 G4Scene* pScene = fpVisManager->GetCurrentScene();
1252 if (!pScene) {
1253 if (verbosity >= G4VisManager::errors) {
1254 G4warn << "ERROR: No current scene. Please create one." << G4endl;
1255 }
1256 return;
1257 }
1258
1259 G4double x1, y1, x2, y2;
1260 std::istringstream is(newValue);
1261 is >> x1 >> y1 >> x2 >> y2;
1262
1263 Line2D* line2D = new Line2D
1264 (x1, y1, x2, y2, fCurrentLineWidth, fCurrentColour);
1265 G4VModel* model =
1267 model->SetType("Line2D");
1268 model->SetGlobalTag("Line2D");
1269 model->SetGlobalDescription("Line2D: " + newValue);
1270 const G4String& currentSceneName = pScene -> GetName ();
1271 G4bool successful = pScene -> AddRunDurationModel (model, warn);
1272 if (successful) {
1273 if (verbosity >= G4VisManager::confirmations) {
1274 G4cout << "A 2D line has been added to scene \""
1275 << currentSceneName << "\"."
1276 << G4endl;
1277 }
1278 }
1279 else G4VisCommandsSceneAddUnsuccessful(verbosity);
1280
1282}
1283
1284G4VisCommandSceneAddLine2D::Line2D::Line2D
1285(G4double x1, G4double y1,
1286 G4double x2, G4double y2,
1287 G4double width, const G4Colour& colour):
1288 fWidth(width), fColour(colour)
1289{
1290 fPolyline.push_back(G4Point3D(x1,y1,0));
1291 fPolyline.push_back(G4Point3D(x2,y2,0));
1292 G4VisAttributes va;
1293 va.SetLineWidth(fWidth);
1294 va.SetColour(fColour);
1295 fPolyline.SetVisAttributes(va);
1296}
1297
1298void G4VisCommandSceneAddLine2D::Line2D::operator()
1299 (G4VGraphicsScene& sceneHandler, const G4ModelingParameters*)
1300{
1301 sceneHandler.BeginPrimitives2D();
1302 sceneHandler.AddPrimitive(fPolyline);
1303 sceneHandler.EndPrimitives2D();
1304}
1305
1306////////////// /vis/scene/add/localAxes ///////////////////////////////////////
1307
1309 G4bool omitable;
1310 fpCommand = new G4UIcommand ("/vis/scene/add/localAxes", this);
1311 fpCommand -> SetGuidance
1312 ("Adds local axes to physical volume(s).");
1313 G4UIparameter* parameter;
1314 parameter = new G4UIparameter ("physical-volume-name", 's', omitable = false);
1315 fpCommand -> SetParameter (parameter);
1316 parameter = new G4UIparameter ("copy-no", 'i', omitable = true);
1317 parameter -> SetGuidance ("If negative, matches any copy no.");
1318 parameter -> SetDefaultValue (-1);
1319 fpCommand -> SetParameter (parameter);
1320}
1321
1325
1329
1331 G4String newValue) {
1332
1333 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
1334 G4bool warn = verbosity >= G4VisManager::warnings;
1335
1336 G4Scene* pScene = fpVisManager->GetCurrentScene();
1337 if (!pScene) {
1338 if (verbosity >= G4VisManager::errors) {
1339 G4warn << "ERROR: No current scene. Please create one." << G4endl;
1340 }
1341 return;
1342 }
1343
1344 G4String name;
1345 G4int copyNo;
1346 std::istringstream is (newValue);
1347 is >> name >> copyNo;
1348
1349 std::vector<G4PhysicalVolumesSearchScene::Findings> findingsVector;
1350
1351 // Search all worlds...
1352 G4TransportationManager* transportationManager =
1354 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
1355 transportationManager->GetWorldsIterator();
1356 size_t nWorlds = transportationManager->GetNoWorlds();
1357 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
1358 G4ModelingParameters mp; // Default - no culling.
1359 G4PhysicalVolumeModel searchModel
1360 (*iterWorld,
1362 G4Transform3D(),
1363 &mp,
1364 true); // Use full extent (avoids initial descent of geometry tree)
1366 (&searchModel, name, copyNo);
1367 searchModel.DescribeYourselfTo (searchScene); // Initiate search.
1368 for (const auto& findings: searchScene.GetFindings()) {
1369 findingsVector.push_back(findings);
1370 }
1371 }
1372
1373 G4int id = 0; // To distinguish axes models by their global description
1374 for (const auto& findings: findingsVector) {
1375
1376 // Create axes model based on size and transformation of found volume(s).
1377 const auto& extent = findings.fpFoundPV->GetLogicalVolume()->GetSolid()->GetExtent();
1378 const auto& transform = findings.fFoundObjectTransformation;
1379
1380 const G4double lengthMax = extent.GetExtentRadius()/2.;
1381 const G4double intLog10LengthMax = std::floor(std::log10(lengthMax));
1382 G4double length = std::pow(10,intLog10LengthMax);
1383 if (5.*length < lengthMax) length *= 5.;
1384 else if (2.*length < lengthMax) length *= 2.;
1385
1386 const auto& axesModel = new G4AxesModel(0.,0.,0.,length,transform);
1387 axesModel->SetGlobalTag("LocalAxesModel");
1388 std::ostringstream oss; oss
1389 << "Local Axes for " << findings.fpFoundPV->GetName()
1390 << ':' << findings.fFoundPVCopyNo << ':' << id++;
1391 axesModel->SetGlobalDescription(oss.str());
1392 // ...so add it to the scene.
1393 G4bool successful = pScene->AddRunDurationModel(axesModel,warn);
1394 if (successful) {
1395 if (verbosity >= G4VisManager::confirmations) {
1396 G4cout << "\"" << findings.fpFoundPV->GetName()
1397 << "\", copy no. " << findings.fFoundPVCopyNo
1398 << ",\n found in searched volume \""
1399 << findings.fpSearchPV->GetName()
1400 << "\" at depth " << findings.fFoundDepth
1401 << ",\n base path: \"" << findings.fFoundBasePVPath
1402 << "\".\n Local axes have been added to scene \""
1403 << pScene->GetName() << "\".";
1404 if (verbosity >= G4VisManager::parameters) {
1405 G4cout << " With extent " << extent
1406 << "\n at " << transform.getRotation()
1407 << " " << transform.getTranslation();
1408 }
1409 G4cout << G4endl;
1410 }
1411 } else {
1413 }
1414 }
1415
1416 if (findingsVector.empty()) {
1417 if (verbosity >= G4VisManager::errors) {
1418 G4warn << "ERROR: Volume \"" << name << "\"";
1419 if (copyNo >= 0) {
1420 G4warn << ", copy no. " << copyNo << ",";
1421 }
1422 G4warn << " not found." << G4endl;
1423 }
1425 return;
1426 }
1427
1429}
1430
1431////////////// /vis/scene/add/logicalVolume //////////////////////////////////
1432
1434 G4bool omitable;
1435 fpCommand = new G4UIcommand ("/vis/scene/add/logicalVolume", this);
1436 fpCommand -> SetGuidance ("Adds a logical volume to the current scene,");
1437 fpCommand -> SetGuidance
1438 ("Shows boolean components (if any), voxels (if any), readout geometry"
1439 "\n (if any), local axes and overlaps (if any), under control of the"
1440 "\n appropriate flag."
1441 "\n Note: voxels are not constructed until start of run -"
1442 "\n \"/run/beamOn\". (For voxels without a run, \"/run/beamOn 0\".)");
1443 G4UIparameter* parameter;
1444 parameter = new G4UIparameter ("logical-volume-name", 's', omitable = false);
1445 fpCommand -> SetParameter (parameter);
1446 parameter = new G4UIparameter ("depth-of-descent", 'i', omitable = true);
1447 parameter -> SetGuidance ("Depth of descent of geometry hierarchy.");
1448 parameter -> SetDefaultValue (1);
1449 fpCommand -> SetParameter (parameter);
1450 parameter = new G4UIparameter ("booleans-flag", 'b', omitable = true);
1451 parameter -> SetDefaultValue (true);
1452 fpCommand -> SetParameter (parameter);
1453 parameter = new G4UIparameter ("voxels-flag", 'b', omitable = true);
1454 parameter -> SetDefaultValue (true);
1455 fpCommand -> SetParameter (parameter);
1456 parameter = new G4UIparameter ("readout-flag", 'b', omitable = true);
1457 parameter -> SetDefaultValue (true);
1458 fpCommand -> SetParameter (parameter);
1459 parameter = new G4UIparameter ("axes-flag", 'b', omitable = true);
1460 parameter -> SetDefaultValue (true);
1461 parameter -> SetGuidance ("Set \"false\" to suppress axes.");
1462 fpCommand -> SetParameter (parameter);
1463 parameter = new G4UIparameter("check-overlap-flag", 'b', omitable = true);
1464 parameter->SetDefaultValue(true);
1465 parameter -> SetGuidance ("Set \"false\" to suppress overlap check.");
1466 fpCommand->SetParameter(parameter);
1467}
1468
1472
1476
1478 G4String newValue) {
1479
1480 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
1481 G4bool warn(verbosity >= G4VisManager::warnings);
1482
1483 G4Scene* pScene = fpVisManager->GetCurrentScene();
1484 if (!pScene) {
1485 if (verbosity >= G4VisManager::errors) {
1486 G4warn << "ERROR: No current scene. Please create one." << G4endl;
1487 }
1488 return;
1489 }
1490
1491 G4String name;
1492 G4int requestedDepthOfDescent;
1493 G4String booleansString, voxelsString, readoutString, axesString;
1494 G4String overlapString;
1495 std::istringstream is (newValue);
1496 is >> name >> requestedDepthOfDescent
1497 >> booleansString >> voxelsString >> readoutString >> axesString
1498 >> overlapString;
1499 G4bool booleans = G4UIcommand::ConvertToBool(booleansString);
1500 G4bool voxels = G4UIcommand::ConvertToBool(voxelsString);
1501 G4bool readout = G4UIcommand::ConvertToBool(readoutString);
1502 G4bool axes = G4UIcommand::ConvertToBool(axesString);
1503 G4bool checkOverlaps = G4UIcommand::ConvertToBool(overlapString);
1504
1506 G4LogicalVolume* pLV = nullptr;
1507 pLV = pLVStore->GetVolume(name);
1508 if (pLV == nullptr) return; // Volume not found; warning message thrown
1509
1510 const std::vector<G4Scene::Model>& rdModelList =
1511 pScene -> GetRunDurationModelList();
1512 std::vector<G4Scene::Model>::const_iterator i;
1513 for (i = rdModelList.begin(); i != rdModelList.end(); ++i) {
1514 if (i->fpModel->GetGlobalDescription().find("Volume") != std::string::npos) break;
1515 }
1516 if (i != rdModelList.end()) {
1517 if (verbosity >= G4VisManager::errors) {
1518 G4warn << "There is already a volume, \""
1519 << i->fpModel->GetGlobalDescription()
1520 << "\",\n in the run-duration model list of scene \""
1521 << pScene -> GetName()
1522 << "\".\n Your logical volume must be the only volume in the scene."
1523 << "\n Create a new scene and try again:"
1524 << "\n /vis/specify " << name
1525 << "\n or"
1526 << "\n /vis/scene/create"
1527 << "\n /vis/scene/add/logicalVolume " << name
1528 << "\n /vis/sceneHandler/attach"
1529 << "\n (and also, if necessary, /vis/viewer/flush)"
1530 << G4endl;
1531 }
1532 return;
1533 }
1534
1536 (pLV, requestedDepthOfDescent, booleans, voxels, readout, checkOverlaps);
1537 const G4String& currentSceneName = pScene -> GetName ();
1538 G4bool successful = pScene -> AddRunDurationModel (model, warn);
1539
1540 if (successful) {
1541
1542 G4bool axesSuccessful = false;
1543 if (axes) {
1544 const G4double radius = model->GetExtent().GetExtentRadius();
1545 const G4double axisLengthMax = radius / 2.;
1546 const G4double intLog10Length = std::floor(std::log10(axisLengthMax));
1547 G4double axisLength = std::pow(10,intLog10Length);
1548 if (5.*axisLength < axisLengthMax) axisLength *= 5.;
1549 else if (2.*axisLength < axisLengthMax) axisLength *= 2.;
1550 const G4double axisWidth = axisLength / 20.;
1551 G4VModel* axesModel = new G4AxesModel(0.,0.,0.,axisLength,axisWidth);
1552 axesSuccessful = pScene -> AddRunDurationModel (axesModel, warn);
1553 }
1554
1555// if (verbosity >= G4VisManager::warnings) {
1556// const std::map<G4String,G4AttDef>* attDefs = model->GetAttDefs();
1557// std::vector<G4AttValue>* attValues = model->CreateCurrentAttValues();
1558// G4warn << G4AttCheck(attValues, attDefs);
1559// delete attValues;
1560// }
1561
1562 if (verbosity >= G4VisManager::confirmations) {
1563 G4cout << "Logical volume \"" << pLV -> GetName ()
1564 << "\" with requested depth of descent "
1565 << requestedDepthOfDescent
1566 << ",\n with";
1567 if (!booleans) G4cout << "out";
1568 G4cout << " boolean components, with";
1569 if (!voxels) G4cout << "out";
1570 G4cout << " voxels,\n with";
1571 if (!readout) G4cout << "out";
1572 G4cout << " readout geometry and with";
1573 if (!checkOverlaps) G4cout << "out";
1574 G4cout << " overlap checking"
1575 << "\n has been added to scene \"" << currentSceneName << "\".";
1576 if (axes) {
1577 if (axesSuccessful) {
1578 G4cout <<
1579 "\n Axes have also been added at the origin of local cooordinates.";
1580 } else {
1581 G4cout <<
1582 "\n Axes have not been added for some reason possibly stated above.";
1583 }
1584 }
1585 G4cout << G4endl;
1586 }
1587 }
1588 else {
1590 return;
1591 }
1592
1594}
1595
1596
1597////////////// /vis/scene/add/logo //////////////////////////////////
1598
1600 G4bool omitable;
1601 fpCommand = new G4UIcommand ("/vis/scene/add/logo", this);
1602 fpCommand -> SetGuidance ("Adds a G4 logo to the current scene.");
1603 fpCommand -> SetGuidance
1604 ("If \"unit\" is \"auto\", height is roughly one tenth of scene extent.");
1605 fpCommand -> SetGuidance
1606 ("\"direction\" is that of outward-facing normal to front face of logo."
1607 "\nIf \"direction\" is \"auto\", logo faces the user in the current viewer.");
1608 fpCommand -> SetGuidance
1609 ("\nIf \"placement\" is \"auto\", logo is placed at bottom right of screen"
1610 "\n when viewed from logo direction.");
1611 G4UIparameter* parameter;
1612 parameter = new G4UIparameter ("height", 'd', omitable = true);
1613 parameter->SetDefaultValue (1.);
1614 fpCommand->SetParameter (parameter);
1615 parameter = new G4UIparameter ("unit", 's', omitable = true);
1616 parameter->SetDefaultValue ("auto");
1617 fpCommand->SetParameter (parameter);
1618 parameter = new G4UIparameter ("direction", 's', omitable = true);
1619 parameter->SetGuidance ("auto|[-]x|[-]y|[-]z");
1620 parameter->SetDefaultValue ("auto");
1621 fpCommand->SetParameter (parameter);
1622 parameter = new G4UIparameter ("red", 'd', omitable = true);
1623 parameter->SetDefaultValue (0.);
1624 fpCommand->SetParameter (parameter);
1625 parameter = new G4UIparameter ("green", 'd', omitable = true);
1626 parameter->SetDefaultValue (1.);
1627 fpCommand->SetParameter (parameter);
1628 parameter = new G4UIparameter ("blue", 'd', omitable = true);
1629 parameter->SetDefaultValue (0.);
1630 fpCommand->SetParameter (parameter);
1631 parameter = new G4UIparameter ("placement", 's', omitable = true);
1632 parameter -> SetParameterCandidates("auto manual");
1633 parameter->SetDefaultValue ("auto");
1634 fpCommand->SetParameter (parameter);
1635 parameter = new G4UIparameter ("xmid", 'd', omitable = true);
1636 parameter->SetDefaultValue (0.);
1637 fpCommand->SetParameter (parameter);
1638 parameter = new G4UIparameter ("ymid", 'd', omitable = true);
1639 parameter->SetDefaultValue (0.);
1640 fpCommand->SetParameter (parameter);
1641 parameter = new G4UIparameter ("zmid", 'd', omitable = true);
1642 parameter->SetDefaultValue (0.);
1643 fpCommand->SetParameter (parameter);
1644 parameter = new G4UIparameter ("unit", 's', omitable = true);
1645 parameter->SetDefaultValue ("m");
1646 fpCommand->SetParameter (parameter);
1647}
1648
1652
1656
1658
1659 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
1660 G4bool warn = verbosity >= G4VisManager::warnings;
1661
1662 G4Scene* pScene = fpVisManager->GetCurrentScene();
1663 if (!pScene) {
1664 if (verbosity >= G4VisManager::errors) {
1665 G4warn << "ERROR: No current scene. Please create one." << G4endl;
1666 }
1667 return;
1668 } else {
1669 if (pScene->GetExtent().GetExtentRadius() <= 0.) {
1670 if (verbosity >= G4VisManager::errors) {
1671 G4warn
1672 << "ERROR: Scene has no extent. Add volumes or use \"/vis/scene/add/extent\"."
1673 << G4endl;
1674 }
1675 return;
1676 }
1677 }
1678
1679 G4VViewer* pViewer = fpVisManager->GetCurrentViewer();
1680 if (!pViewer) {
1681 if (verbosity >= G4VisManager::errors) {
1682 G4warn <<
1683 "ERROR: G4VisCommandSceneAddLogo::SetNewValue: no viewer."
1684 "\n Auto direction needs a viewer."
1685 << G4endl;
1686 }
1687 return;
1688 }
1689
1690 G4double userHeight, red, green, blue, xmid, ymid, zmid;
1691 G4String userHeightUnit, direction, placement, positionUnit;
1692 std::istringstream is (newValue);
1693 is >> userHeight >> userHeightUnit >> direction
1694 >> red >> green >> blue
1695 >> placement
1696 >> xmid >> ymid >> zmid >> positionUnit;
1697
1698 G4double height = userHeight;
1699 const G4VisExtent& sceneExtent = pScene->GetExtent(); // Existing extent.
1700 if (userHeightUnit == "auto") {
1701 height *= 0.2 * sceneExtent.GetExtentRadius();
1702 } else {
1703 height *= G4UIcommand::ValueOf(userHeightUnit);
1704 }
1705
1706 G4double unit = G4UIcommand::ValueOf(positionUnit);
1707 xmid *= unit; ymid *= unit; zmid *= unit;
1708
1709 Direction logoDirection = X; // Initialise to keep some compilers happy.
1710 if (direction == "auto") {
1711 // Take cue from viewer
1712 const G4Vector3D& vp =
1714 if (vp.x() > vp.y() && vp.x() > vp.z()) logoDirection = X;
1715 else if (vp.x() < vp.y() && vp.x() < vp.z()) logoDirection = minusX;
1716 else if (vp.y() > vp.x() && vp.y() > vp.z()) logoDirection = Y;
1717 else if (vp.y() < vp.x() && vp.y() < vp.z()) logoDirection = minusY;
1718 else if (vp.z() > vp.x() && vp.z() > vp.y()) logoDirection = Z;
1719 else if (vp.z() < vp.x() && vp.z() < vp.y()) logoDirection = minusZ;
1720 }
1721 else if (direction[0] == 'x') logoDirection = X;
1722 else if (direction[0] == 'y') logoDirection = Y;
1723 else if (direction[0] == 'z') logoDirection = Z;
1724 else if (direction[0] == '-') {
1725 if (direction[1] == 'x') logoDirection = minusX;
1726 else if (direction[1] == 'y') logoDirection = minusY;
1727 else if (direction[1] == 'z') logoDirection = minusZ;
1728 } else {
1729 if (verbosity >= G4VisManager::errors) {
1730 G4warn << "ERROR: Unrecogniseed direction: \""
1731 << direction << "\"." << G4endl;
1732 return;
1733 }
1734 }
1735
1736 G4bool autoPlacing = false; if (placement == "auto") autoPlacing = true;
1737 // Parameters read and interpreted.
1738
1739 // Current scene extent
1740 const G4double xmin = sceneExtent.GetXmin();
1741 const G4double xmax = sceneExtent.GetXmax();
1742 const G4double ymin = sceneExtent.GetYmin();
1743 const G4double ymax = sceneExtent.GetYmax();
1744 const G4double zmin = sceneExtent.GetZmin();
1745 const G4double zmax = sceneExtent.GetZmax();
1746
1747 // Test existing extent and issue warnings...
1748 G4bool worried = false;
1749 if (sceneExtent.GetExtentRadius() == 0) {
1750 worried = true;
1751 if (verbosity >= G4VisManager::warnings) {
1752 G4warn <<
1753 "WARNING: Existing scene does not yet have any extent."
1754 "\n Maybe you have not yet added any geometrical object."
1755 << G4endl;
1756 }
1757 }
1758
1759 // Useful constants, etc...
1760 const G4double halfHeight(height / 2.);
1761 const G4double comfort(0.01); // 0.15 seems too big. 0.05 might be better.
1762 const G4double freeHeightFraction (1. + 2. * comfort);
1763
1764 // Test existing scene for room...
1765 G4bool room = true;
1766 switch (logoDirection) {
1767 case X:
1768 case minusX:
1769 if (freeHeightFraction * (xmax - xmin) < height) room = false;
1770 break;
1771 case Y:
1772 case minusY:
1773 if (freeHeightFraction * (ymax - ymin) < height) room = false;
1774 break;
1775 case Z:
1776 case minusZ:
1777 if (freeHeightFraction * (zmax - zmin) < height) room = false;
1778 break;
1779 }
1780 if (!room) {
1781 worried = true;
1782 if (verbosity >= G4VisManager::warnings) {
1783 G4warn <<
1784 "WARNING: Not enough room in existing scene. Maybe logo is too large."
1785 << G4endl;
1786 }
1787 }
1788 if (worried) {
1789 if (verbosity >= G4VisManager::warnings) {
1790 G4warn <<
1791 "WARNING: The logo you have asked for is bigger than the existing"
1792 "\n scene. Maybe you have added it too soon. It is recommended that"
1793 "\n you add the logo last so that it can be correctly auto-positioned"
1794 "\n so as not to be obscured by any existing object and so that the"
1795 "\n view parameters can be correctly recalculated."
1796 << G4endl;
1797 }
1798 }
1799
1800 G4double sxmid(xmid), symid(ymid), szmid(zmid);
1801 if (autoPlacing) {
1802 // Aim to place at bottom right of screen when viewed from logoDirection.
1803 // Give some comfort zone.
1804 const G4double xComfort = comfort * (xmax - xmin);
1805 const G4double yComfort = comfort * (ymax - ymin);
1806 const G4double zComfort = comfort * (zmax - zmin);
1807 switch (logoDirection) {
1808 case X: // y-axis up, z-axis to left?
1809 sxmid = xmax + halfHeight + xComfort;
1810 symid = ymin - yComfort;
1811 szmid = zmin - zComfort;
1812 break;
1813 case minusX: // y-axis up, z-axis to right?
1814 sxmid = xmin - halfHeight - xComfort;
1815 symid = ymin - yComfort;
1816 szmid = zmax + zComfort;
1817 break;
1818 case Y: // z-axis up, x-axis to left?
1819 sxmid = xmin - xComfort;
1820 symid = ymax + halfHeight + yComfort;
1821 szmid = zmin - zComfort;
1822 break;
1823 case minusY: // z-axis up, x-axis to right?
1824 sxmid = xmax + xComfort;
1825 symid = ymin - halfHeight - yComfort;
1826 szmid = zmin - zComfort;
1827 break;
1828 case Z: // y-axis up, x-axis to right?
1829 sxmid = xmax + xComfort;
1830 symid = ymin - yComfort;
1831 szmid = zmax + halfHeight + zComfort;
1832 break;
1833 case minusZ: // y-axis up, x-axis to left?
1834 sxmid = xmin - xComfort;
1835 symid = ymin - yComfort;
1836 szmid = zmin - halfHeight - zComfort;
1837 break;
1838 }
1839 }
1840
1841 G4Transform3D transform;
1842 switch (logoDirection) {
1843 case X: // y-axis up, z-axis to left?
1844 transform = G4RotateY3D(halfpi);
1845 break;
1846 case minusX: // y-axis up, z-axis to right?
1847 transform = G4RotateY3D(-halfpi);
1848 break;
1849 case Y: // z-axis up, x-axis to left?
1850 transform = G4RotateX3D(-halfpi) * G4RotateZ3D(pi);
1851 break;
1852 case minusY: // z-axis up, x-axis to right?
1853 transform = G4RotateX3D(halfpi);
1854 break;
1855 case Z: // y-axis up, x-axis to right?
1856 // No transformation required.
1857 break;
1858 case minusZ: // y-axis up, x-axis to left?
1859 transform = G4RotateY3D(pi);
1860 break;
1861 }
1862 transform = G4Translate3D(sxmid,symid,szmid) * transform;
1863
1864 G4VisAttributes visAtts(G4Colour(red, green, blue));
1865 visAtts.SetForceSolid(true); // Always solid.
1866
1867 G4Logo* logo = new G4Logo(height,visAtts,transform);
1868 G4VModel* model =
1870 model->SetType("G4Logo");
1871 model->SetGlobalTag("G4Logo");
1872 model->SetGlobalDescription("G4Logo: " + newValue);
1873 G4double& h = height;
1874 G4double h2 = h/2.;
1875 G4VisExtent extent(-h,h,-h2,h2,-h2,h2);
1876 model->SetExtent(extent.Transform(transform));
1877 // This extent gets "added" to existing scene extent in
1878 // AddRunDurationModel below.
1879 const G4String& currentSceneName = pScene -> GetName ();
1880 G4bool successful = pScene -> AddRunDurationModel (model, warn);
1881 if (successful) {
1882 if (verbosity >= G4VisManager::confirmations) {
1883 G4cout << "G4 Logo of height " << userHeight << ' ' << userHeightUnit
1884 << ", " << direction << "-direction, added to scene \""
1885 << currentSceneName << "\"";
1886 if (verbosity >= G4VisManager::parameters) {
1887 G4cout << "\n with extent " << extent
1888 << "\n at " << transform.getRotation()
1889 << " " << transform.getTranslation();
1890 }
1891 G4cout << G4endl;
1892 }
1893 }
1894 else G4VisCommandsSceneAddUnsuccessful(verbosity);
1895
1897}
1898
1899G4VisCommandSceneAddLogo::G4Logo::G4Logo
1900(G4double height, const G4VisAttributes& visAtts, const G4Transform3D& transform)
1901{
1902 const G4double& h = height;
1903 const G4double h2 = 0.5 * h; // Half height.
1904 const G4double ri = 0.25 * h; // Inner radius.
1905 const G4double ro = 0.5 * h; // Outer radius.
1906 const G4double ro2 = 0.5 * ro; // Half outer radius.
1907 const G4double w = ro - ri; // Width.
1908 const G4double w2 = 0.5 * w; // Half width.
1909 const G4double d2 = 0.2 * h; // Half depth.
1910 const G4double f1 = 0.05 * h; // left edge of stem of "4".
1911 const G4double f2 = -0.3 * h; // bottom edge of cross of "4".
1912 const G4double e = 1.e-4 * h; // epsilon.
1913 const G4double xt = f1, yt = h2; // Top of slope.
1914 const G4double xb = -h2, yb = f2 + w; // Bottom of slope.
1915 const G4double dx = xt - xb, dy = yt - yb;
1916 const G4double angle = std::atan2(dy,dx);
1918 rm.rotateZ(angle*rad);
1919 const G4double d = std::sqrt(dx * dx + dy * dy);
1920 const G4double ss = h; // Half height of square subtractor
1921 const G4double y8 = ss; // Choose y of subtractor for outer slope.
1922 const G4double x8 = ((-ss * d - dx * (yt - y8)) / dy) + xt;
1923 G4double y9 = ss; // Choose y of subtractor for inner slope.
1924 G4double x9 = ((-(ss - w) * d - dx * (yt - y8)) / dy) + xt;
1925 // But to get inner, we make a triangle translated by...
1926 const G4double xtr = ss - f1, ytr = -ss - f2 -w;
1927 x9 += xtr; y9 += ytr;
1928
1929 // The idea here is to create a polyhedron for the G and the 4. To do
1930 // this we use Geant4 geometry solids and make boolean operations.
1931 // Note that we do not need to keep the solids. We use local objects,
1932 // which, of course, are deleted on leaving this function. This
1933 // is contrary to the usual requirement for solids that are part of the
1934 // detector for which solids MUST be created on the heap (with "new").
1935 // Finally we invoke CreatePolyhedron, which creates a polyhedron on the heap
1936 // and returns a pointer. It is the user's responsibility to delete,
1937 // which is done in the destructor of this class. Thus the polyhedra,
1938 // created here, remain on the heap for the lifetime of the job.
1939
1940 // G...
1941 G4Tubs tG("tG",ri,ro,d2,0.15*pi,1.85*pi);
1942 G4Box bG("bG",w2,ro2,d2);
1943 G4UnionSolid logoG("logoG",&tG,&bG,G4Translate3D(ri+w2,-ro2,0.));
1944 // Create with these vis atts (force solid) and current no of sides per circle.
1947 fpG = logoG.CreatePolyhedron();
1949 fpG->SetVisAttributes(visAtts);
1950 fpG->Transform(G4Translate3D(-0.55*h,0.,0.));
1951 fpG->Transform(transform);
1952
1953 // 4...
1954 G4Box b1("b1",h2,h2,d2);
1955 G4Box bS("bS",ss,ss,d2+e); // Subtractor.
1956 G4Box bS2("bS2",ss,ss,d2+2.*e); // 2nd Subtractor.
1957 G4SubtractionSolid s1("s1",&b1,&bS,G4Translate3D(f1-ss,f2-ss,0.));
1958 G4SubtractionSolid s2("s2",&s1,&bS,G4Translate3D(f1+ss+w,f2-ss,0.));
1959 G4SubtractionSolid s3("s3",&s2,&bS,G4Translate3D(f1+ss+w,f2+ss+w,0.));
1961 ("s4",&s3,&bS,G4Transform3D(rm,G4ThreeVector(x8,y8,0.)));
1962 G4SubtractionSolid s5 // Triangular hole.
1963 ("s5",&bS,&bS2,G4Transform3D(rm,G4ThreeVector(x9,y9,0.)));
1964 G4SubtractionSolid logo4("logo4",&s4,&s5,G4Translate3D(-xtr,-ytr,0.));
1965 fp4 = logo4.CreatePolyhedron();
1966 /* Experiment with creating own polyhedron...
1967 int nNodes = 4;
1968 int nFaces = 4;
1969 double xyz[][3] = {{0,0,0},{1*m,0,0},{0,1*m,0},{0,0,1*m}};
1970 int faces[][4] = {{1,3,2,0},{1,2,4,0},{1,4,3,0},{2,3,4,0}};
1971 fp4 = new G4Polyhedron();
1972 fp4->createPolyhedron(nNodes,nFaces,xyz,faces);
1973 */
1974 fp4->SetVisAttributes(visAtts);
1975 fp4->Transform(G4Translate3D(0.55*h,0.,0.));
1976 fp4->Transform(transform);
1977}
1978
1979G4VisCommandSceneAddLogo::G4Logo::~G4Logo() {
1980 delete fpG;
1981 delete fp4;
1982}
1983
1984void G4VisCommandSceneAddLogo::G4Logo::operator()
1985 (G4VGraphicsScene& sceneHandler, const G4ModelingParameters*) {
1986 sceneHandler.BeginPrimitives();
1987 sceneHandler.AddPrimitive(*fpG);
1988 sceneHandler.AddPrimitive(*fp4);
1989 sceneHandler.EndPrimitives();
1990}
1991
1992////////////// /vis/scene/add/logo2D ///////////////////////////////////////
1993
1995 G4bool omitable;
1996 fpCommand = new G4UIcommand ("/vis/scene/add/logo2D", this);
1997 fpCommand -> SetGuidance ("Adds 2D logo to current scene.");
1998 G4UIparameter* parameter;
1999 parameter = new G4UIparameter ("size", 'i', omitable = true);
2000 parameter -> SetGuidance ("Screen size of text in pixels.");
2001 parameter -> SetDefaultValue (48);
2002 fpCommand -> SetParameter (parameter);
2003 parameter = new G4UIparameter ("x-position", 'd', omitable = true);
2004 parameter -> SetGuidance ("x screen position in range -1 < x < 1.");
2005 parameter -> SetDefaultValue (-0.9);
2006 fpCommand -> SetParameter (parameter);
2007 parameter = new G4UIparameter ("y-position", 'd', omitable = true);
2008 parameter -> SetGuidance ("y screen position in range -1 < y < 1.");
2009 parameter -> SetDefaultValue (-0.9);
2010 fpCommand -> SetParameter (parameter);
2011 parameter = new G4UIparameter ("layout", 's', omitable = true);
2012 parameter -> SetGuidance ("Layout, i.e., adjustment: left|centre|right.");
2013 parameter -> SetDefaultValue ("left");
2014 fpCommand -> SetParameter (parameter);
2015}
2016
2020
2024
2026{
2027 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
2028 G4bool warn(verbosity >= G4VisManager::warnings);
2029
2030 G4Scene* pScene = fpVisManager->GetCurrentScene();
2031 if (!pScene) {
2032 if (verbosity >= G4VisManager::errors) {
2033 G4warn << "ERROR: No current scene. Please create one." << G4endl;
2034 }
2035 return;
2036 }
2037
2038 G4int size;
2039 G4double x, y;
2040 G4String layoutString;
2041 std::istringstream is(newValue);
2042 is >> size >> x >> y >> layoutString;
2044 if (layoutString[0] == 'l') layout = G4Text::left;
2045 else if (layoutString[0] == 'c') layout = G4Text::centre;
2046 else if (layoutString[0] == 'r') layout = G4Text::right;
2047
2048 Logo2D* logo2D = new Logo2D(fpVisManager, size, x, y, layout);
2049 G4VModel* model =
2051 model->SetType("G4Logo2D");
2052 model->SetGlobalTag("G4Logo2D");
2053 model->SetGlobalDescription("G4Logo2D: " + newValue);
2054 const G4String& currentSceneName = pScene -> GetName ();
2055 G4bool successful = pScene -> AddRunDurationModel (model, warn);
2056 if (successful) {
2057 if (verbosity >= G4VisManager::confirmations) {
2058 G4cout << "2D logo has been added to scene \""
2059 << currentSceneName << "\"."
2060 << G4endl;
2061 }
2062 }
2063 else G4VisCommandsSceneAddUnsuccessful(verbosity);
2064
2066}
2067
2068void G4VisCommandSceneAddLogo2D::Logo2D::operator()
2069 (G4VGraphicsScene& sceneHandler, const G4ModelingParameters*)
2070{
2071 G4Text text("Geant4", G4Point3D(fX, fY, 0.));
2072 text.SetScreenSize(fSize);
2073 text.SetLayout(fLayout);
2074 G4VisAttributes textAtts(G4Colour::Brown());
2075 text.SetVisAttributes(textAtts);
2076 sceneHandler.BeginPrimitives2D();
2077 sceneHandler.AddPrimitive(text);
2078 sceneHandler.EndPrimitives2D();
2079}
2080
2081////////////// /vis/scene/add/magneticField ///////////////////////////////////////
2082
2084 fpCommand = new G4UIcommand ("/vis/scene/add/magneticField", this);
2085 fpCommand -> SetGuidance
2086 ("Adds magnetic field representation to current scene.");
2088 const G4UIcommand* addElecFieldCmd = tree->FindPath("/vis/scene/add/electricField");
2089 // Pick up additional guidance from /vis/scene/add/electricField
2090 CopyGuidanceFrom(addElecFieldCmd,fpCommand,1);
2091 // Pick up parameters from /vis/scene/add/electricField
2092 CopyParametersFrom(addElecFieldCmd,fpCommand);
2093}
2094
2098
2102
2104(G4UIcommand*, G4String newValue) {
2105
2106 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
2107 G4bool warn(verbosity >= G4VisManager::warnings);
2108
2109 G4Scene* pScene = fpVisManager->GetCurrentScene();
2110 if (!pScene) {
2111 if (verbosity >= G4VisManager::errors) {
2112 G4warn << "ERROR: No current scene. Please create one." << G4endl;
2113 }
2114 return;
2115 }
2116
2117 G4int nDataPointsPerHalfScene;
2118 G4String representation;
2119 std::istringstream iss(newValue);
2120 iss >> nDataPointsPerHalfScene >> representation;
2122 modelRepresentation = G4ElectricFieldModel::fullArrow;
2123 if (representation == "lightArrow") {
2124 modelRepresentation = G4ElectricFieldModel::lightArrow;
2125 }
2126 G4VModel* model;
2127 model = new G4MagneticFieldModel
2128 (nDataPointsPerHalfScene,modelRepresentation,
2132 const G4String& currentSceneName = pScene -> GetName ();
2133 G4bool successful = pScene -> AddRunDurationModel (model, warn);
2134 if (successful) {
2135 if (verbosity >= G4VisManager::confirmations) {
2136 G4cout
2137 << "Magnetic field, if any, will be drawn in scene \""
2138 << currentSceneName
2139 << "\"\n with "
2140 << nDataPointsPerHalfScene
2141 << " data points per half extent and with representation \""
2142 << representation
2143 << '\"'
2144 << G4endl;
2145 }
2146 }
2147 else G4VisCommandsSceneAddUnsuccessful(verbosity);
2148
2150}
2151
2152////////////// /vis/scene/add/psHits ///////////////////////////////////////
2153
2155 G4bool omitable;
2156 fpCommand = new G4UIcmdWithAString ("/vis/scene/add/psHits", this);
2157 fpCommand -> SetGuidance
2158 ("Adds Primitive Scorer Hits (PSHits) to current scene.");
2159 fpCommand -> SetGuidance
2160 ("PSHits are drawn at end of run when the scene in which"
2161 "\nthey are added is current.");
2162 fpCommand -> SetGuidance
2163 ("Optional parameter specifies name of scoring map. By default all"
2164 "\nscoring maps registered with the G4ScoringManager are drawn.");
2165 fpCommand -> SetParameterName ("mapname", omitable = true);
2166 fpCommand -> SetDefaultValue ("all");
2167}
2168
2172
2176
2178(G4UIcommand*, G4String newValue)
2179{
2180 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
2181 G4bool warn(verbosity >= G4VisManager::warnings);
2182
2183 G4Scene* pScene = fpVisManager->GetCurrentScene();
2184 if (!pScene) {
2185 if (verbosity >= G4VisManager::errors) {
2186 G4warn << "ERROR: No current scene. Please create one." << G4endl;
2187 }
2188 return;
2189 }
2190
2191 G4VModel* model = new G4PSHitsModel(newValue);
2192 const G4String& currentSceneName = pScene -> GetName ();
2193 G4bool successful = pScene -> AddEndOfRunModel (model, warn);
2194 if (successful) {
2195 if (verbosity >= G4VisManager::confirmations) {
2196 if (newValue == "all") {
2197 G4cout << "All Primitive Scorer hits";
2198 } else {
2199 G4cout << "Hits of Primitive Scorer \"" << newValue << '"';
2200 }
2201 G4cout << " will be drawn at end of run in scene \""
2202 << currentSceneName << "\"."
2203 << G4endl;
2204 }
2205 }
2206 else G4VisCommandsSceneAddUnsuccessful(verbosity);
2207
2209}
2210
2211////////////// /vis/scene/add/scale //////////////////////////////////
2212
2214 G4bool omitable;
2215 fpCommand = new G4UIcommand ("/vis/scene/add/scale", this);
2216 fpCommand -> SetGuidance
2217 ("Adds an annotated scale line to the current scene.");
2218 fpCommand -> SetGuidance
2219 ("If \"unit\" is \"auto\", length is roughly one tenth of the scene extent.");
2220 fpCommand -> SetGuidance
2221 ("If \"direction\" is \"auto\", scale is roughly in the plane of the current view.");
2222 fpCommand -> SetGuidance
2223 ("If \"placement\" is \"auto\", scale is placed at bottom left of current view."
2224 "\n Otherwise placed at (xmid,ymid,zmid).");
2225 fpCommand -> SetGuidance
2226 ("An annotated line in the specified direction with tick marks at the"
2227 "\nend. If autoPlacing is true it is required to be centred at the"
2228 "\nfront, right, bottom corner of the world space, comfortably outside"
2229 "\nthe existing bounding box/sphere so that existing objects do not"
2230 "\nobscure it. Otherwise it is required to be drawn with mid-point at"
2231 "\n(xmid, ymid, zmid)."
2232 "\n"
2233 "\nThe auto placing algorithm is (approx):"
2234 "\n x = xmin + (1 + comfort) * (xmax - xmin);"
2235 "\n y = ymin - comfort * (ymax - ymin);"
2236 "\n z = zmin + (1 + comfort) * (zmax - zmin);"
2237 "\n if direction == x then (x - length,y,z) to (x,y,z);"
2238 "\n if direction == y then (x,y,z) to (x,y + length,z);"
2239 "\n if direction == z then (x,y,z - length) to (x,y,z);"
2240 );
2241 G4UIparameter* parameter;
2242 parameter = new G4UIparameter ("length", 'd', omitable = true);
2243 parameter->SetDefaultValue (1.);
2244 fpCommand->SetParameter (parameter);
2245 parameter = new G4UIparameter ("unit", 's', omitable = true);
2246 parameter->SetDefaultValue ("auto");
2247 fpCommand->SetParameter (parameter);
2248 parameter = new G4UIparameter ("direction", 's', omitable = true);
2249 parameter->SetGuidance ("auto|x|y|z");
2250 parameter->SetDefaultValue ("auto");
2251 fpCommand->SetParameter (parameter);
2252 parameter = new G4UIparameter ("red", 'd', omitable = true);
2253 parameter->SetDefaultValue (1.);
2254 fpCommand->SetParameter (parameter);
2255 parameter = new G4UIparameter ("green", 'd', omitable = true);
2256 parameter->SetDefaultValue (0.);
2257 fpCommand->SetParameter (parameter);
2258 parameter = new G4UIparameter ("blue", 'd', omitable = true);
2259 parameter->SetDefaultValue (0.);
2260 fpCommand->SetParameter (parameter);
2261 parameter = new G4UIparameter ("placement", 's', omitable = true);
2262 parameter -> SetParameterCandidates("auto manual");
2263 parameter->SetDefaultValue ("auto");
2264 fpCommand->SetParameter (parameter);
2265 parameter = new G4UIparameter ("xmid", 'd', omitable = true);
2266 parameter->SetDefaultValue (0.);
2267 fpCommand->SetParameter (parameter);
2268 parameter = new G4UIparameter ("ymid", 'd', omitable = true);
2269 parameter->SetDefaultValue (0.);
2270 fpCommand->SetParameter (parameter);
2271 parameter = new G4UIparameter ("zmid", 'd', omitable = true);
2272 parameter->SetDefaultValue (0.);
2273 fpCommand->SetParameter (parameter);
2274 parameter = new G4UIparameter ("unit", 's', omitable = true);
2275 parameter->SetDefaultValue ("m");
2276 fpCommand->SetParameter (parameter);
2277}
2278
2282
2286
2288
2289 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
2290 G4bool warn = verbosity >= G4VisManager::warnings;
2291
2292 G4Scene* pScene = fpVisManager->GetCurrentScene();
2293 if (!pScene) {
2294 if (verbosity >= G4VisManager::errors) {
2295 G4warn << "ERROR: No current scene. Please create one." << G4endl;
2296 }
2297 return;
2298 } else {
2299 if (pScene->GetExtent().GetExtentRadius() <= 0.) {
2300 if (verbosity >= G4VisManager::errors) {
2301 G4warn
2302 << "ERROR: Scene has no extent. Add volumes or use \"/vis/scene/add/extent\"."
2303 << G4endl;
2304 }
2305 return;
2306 }
2307 }
2308
2309 G4double userLength, red, green, blue, xmid, ymid, zmid;
2310 G4String userLengthUnit, direction, placement, positionUnit;
2311 std::istringstream is (newValue);
2312 is >> userLength >> userLengthUnit >> direction
2313 >> red >> green >> blue
2314 >> placement
2315 >> xmid >> ymid >> zmid >> positionUnit;
2316
2317 G4double length = userLength;
2318 const G4VisExtent& sceneExtent = pScene->GetExtent(); // Existing extent.
2319 if (userLengthUnit == "auto") {
2320 const G4double lengthMax = 0.5 * sceneExtent.GetExtentRadius();
2321 const G4double intLog10Length = std::floor(std::log10(lengthMax));
2322 length = std::pow(10,intLog10Length);
2323 if (5.*length < lengthMax) length *= 5.;
2324 else if (2.*length < lengthMax) length *= 2.;
2325 } else {
2326 length *= G4UIcommand::ValueOf(userLengthUnit);
2327 }
2328 G4String annotation = G4BestUnit(length,"Length");
2329
2330 G4double unit = G4UIcommand::ValueOf(positionUnit);
2331 xmid *= unit; ymid *= unit; zmid *= unit;
2332
2333 Scale::Direction scaleDirection (Scale::x);
2334 if (direction[0] == 'y') scaleDirection = Scale::y;
2335 if (direction[0] == 'z') scaleDirection = Scale::z;
2336
2337 G4VViewer* pViewer = fpVisManager->GetCurrentViewer();
2338 if (!pViewer) {
2339 if (verbosity >= G4VisManager::errors) {
2340 G4warn <<
2341 "ERROR: G4VisCommandSceneAddScale::SetNewValue: no viewer."
2342 "\n Auto direction needs a viewer."
2343 << G4endl;
2344 }
2345 return;
2346 }
2347
2348 const G4Vector3D& vp =
2350 const G4Vector3D& up =
2351 pViewer->GetViewParameters().GetUpVector();
2352
2353 if (direction == "auto") { // Takes cue from viewer.
2354 if (std::abs(vp.x()) > std::abs(vp.y()) &&
2355 std::abs(vp.x()) > std::abs(vp.z())) { // x viewpoint
2356 if (std::abs(up.y()) > std::abs(up.z())) scaleDirection = Scale::z;
2357 else scaleDirection = Scale::y;
2358 }
2359 else if (std::abs(vp.y()) > std::abs(vp.x()) &&
2360 std::abs(vp.y()) > std::abs(vp.z())) { // y viewpoint
2361 if (std::abs(up.x()) > std::abs(up.z())) scaleDirection = Scale::z;
2362 else scaleDirection = Scale::x;
2363 }
2364 else if (std::abs(vp.z()) > std::abs(vp.x()) &&
2365 std::abs(vp.z()) > std::abs(vp.y())) { // z viewpoint
2366 if (std::abs(up.y()) > std::abs(up.x())) scaleDirection = Scale::x;
2367 else scaleDirection = Scale::y;
2368 }
2369 }
2370
2371 G4bool autoPlacing = false; if (placement == "auto") autoPlacing = true;
2372 // Parameters read and interpreted.
2373
2374 // Useful constants, etc...
2375 const G4double halfLength(length / 2.);
2376 const G4double comfort(0.01); // 0.15 seems too big. 0.05 might be better.
2377 const G4double freeLengthFraction (1. + 2. * comfort);
2378
2379 const G4double xmin = sceneExtent.GetXmin();
2380 const G4double xmax = sceneExtent.GetXmax();
2381 const G4double ymin = sceneExtent.GetYmin();
2382 const G4double ymax = sceneExtent.GetYmax();
2383 const G4double zmin = sceneExtent.GetZmin();
2384 const G4double zmax = sceneExtent.GetZmax();
2385
2386 // Test existing extent and issue warnings...
2387 G4bool worried = false;
2388 if (sceneExtent.GetExtentRadius() == 0) {
2389 worried = true;
2390 if (verbosity >= G4VisManager::warnings) {
2391 G4warn <<
2392 "WARNING: Existing scene does not yet have any extent."
2393 "\n Maybe you have not yet added any geometrical object."
2394 << G4endl;
2395 }
2396 }
2397
2398 // Test existing scene for room...
2399 G4bool room = true;
2400 switch (scaleDirection) {
2401 case Scale::x:
2402 if (freeLengthFraction * (xmax - xmin) < length) room = false;
2403 break;
2404 case Scale::y:
2405 if (freeLengthFraction * (ymax - ymin) < length) room = false;
2406 break;
2407 case Scale::z:
2408 if (freeLengthFraction * (zmax - zmin) < length) room = false;
2409 break;
2410 }
2411 if (!room) {
2412 worried = true;
2413 if (verbosity >= G4VisManager::warnings) {
2414 G4warn <<
2415 "WARNING: Not enough room in existing scene. Maybe scale is too long."
2416 << G4endl;
2417 }
2418 }
2419 if (worried) {
2420 if (verbosity >= G4VisManager::warnings) {
2421 G4warn <<
2422 "WARNING: The scale you have asked for is bigger than the existing"
2423 "\n scene. Maybe you have added it too soon. It is recommended that"
2424 "\n you add the scale last so that it can be correctly auto-positioned"
2425 "\n so as not to be obscured by any existing object and so that the"
2426 "\n view parameters can be correctly recalculated."
2427 << G4endl;
2428 }
2429 }
2430
2431 // Now figure out the extent...
2432 //
2433 // This creates a representation of annotated line in the specified
2434 // direction with tick marks at the end. If autoPlacing is true it
2435 // is required to be centred at the front, right, bottom corner of
2436 // the world space, comfortably outside the existing bounding
2437 // box/sphere so that existing objects do not obscure it. Otherwise
2438 // it is required to be drawn with mid-point at (xmid, ymid, zmid).
2439 //
2440 // The auto placing algorithm might be:
2441 // x = xmin + (1 + comfort) * (xmax - xmin)
2442 // y = ymin - comfort * (ymax - ymin)
2443 // z = zmin + (1 + comfort) * (zmax - zmin)
2444 // if direction == x then (x - length,y,z) to (x,y,z)
2445 // if direction == y then (x,y,z) to (x,y + length,z)
2446 // if direction == z then (x,y,z - length) to (x,y,z)
2447 //
2448 // Implement this in two parts. Here, use the scale's extent to
2449 // "expand" the scene's extent. Then rendering - in
2450 // G4VSceneHandler::AddPrimitive(const G4Scale&) - simply has to
2451 // ensure it's within the new extent.
2452 //
2453
2454 G4double sxmid(xmid), symid(ymid), szmid(zmid);
2455 if (autoPlacing) {
2456 // Aim to place at bottom right of screen in current view.
2457 // Give some comfort zone.
2458 const G4double xComfort = comfort * (xmax - xmin);
2459 const G4double yComfort = comfort * (ymax - ymin);
2460 const G4double zComfort = comfort * (zmax - zmin);
2461 switch (scaleDirection) {
2462 case Scale::x:
2463 if (vp.z() > 0.) {
2464 sxmid = xmax + xComfort;
2465 symid = ymin - yComfort;
2466 szmid = zmin - zComfort;
2467 } else {
2468 sxmid = xmin - xComfort;
2469 symid = ymin - yComfort;
2470 szmid = zmax + zComfort;
2471 }
2472 break;
2473 case Scale::y:
2474 if (vp.x() > 0.) {
2475 sxmid = xmin - xComfort;
2476 symid = ymax + yComfort;
2477 szmid = zmin - zComfort;
2478 } else {
2479 sxmid = xmax + xComfort;
2480 symid = ymin - yComfort;
2481 szmid = zmin - zComfort;
2482 }
2483 break;
2484 case Scale::z:
2485 if (vp.x() > 0.) {
2486 sxmid = xmax + xComfort;
2487 symid = ymin - yComfort;
2488 szmid = zmax + zComfort;
2489 } else {
2490 sxmid = xmin - xComfort;
2491 symid = ymin - yComfort;
2492 szmid = zmax + zComfort;
2493 }
2494 break;
2495 }
2496 }
2497
2498 G4Transform3D transform;
2499 const G4double h = halfLength;
2500 const G4double t = h/5.;
2501 G4VisExtent scaleExtent(-h,h,-t,t,-t,t);
2502 switch (scaleDirection) {
2503 case Scale::x:
2504 break;
2505 case Scale::y:
2506 transform = G4RotateZ3D(halfpi);
2507 break;
2508 case Scale::z:
2509 transform = G4RotateY3D(halfpi);
2510 break;
2511 }
2512 transform = G4Translate3D(sxmid,symid,szmid) * transform;
2513 scaleExtent = scaleExtent.Transform(transform);
2514
2515 G4Colour colour(red, green, blue);
2516 if (direction == "auto") {
2517 switch (scaleDirection) {
2518 case Scale::x:
2519 colour = G4Colour::Red();
2520 break;
2521 case Scale::y:
2522 colour = G4Colour::Green();
2523 break;
2524 case Scale::z:
2525 colour = G4Colour::Blue();
2526 break;
2527 }
2528 }
2529 G4VisAttributes visAttr(colour);
2530
2531 Scale* scale = new Scale
2532 (visAttr, length, transform,
2533 annotation, fCurrentTextSize, colour);
2534 G4VModel* model = new G4CallbackModel<Scale>(scale);
2535 model->SetType("Scale");
2536 model->SetGlobalTag("Scale");
2537 model->SetGlobalDescription("Scale: " + newValue);
2538 model->SetExtent(scaleExtent);
2539
2540 const G4String& currentSceneName = pScene -> GetName ();
2541 G4bool successful = pScene -> AddRunDurationModel (model, warn);
2542 if (successful) {
2543 if (verbosity >= G4VisManager::confirmations) {
2544 G4cout << "Scale of " << annotation
2545 << " added to scene \"" << currentSceneName << "\".";
2546 if (verbosity >= G4VisManager::parameters) {
2547 G4cout << "\n with extent " << scaleExtent
2548 << "\n at " << transform.getRotation()
2549 << " " << transform.getTranslation();
2550 }
2551 G4cout << G4endl;
2552 }
2553 }
2554 else G4VisCommandsSceneAddUnsuccessful(verbosity);
2555
2557}
2558
2559G4VisCommandSceneAddScale::Scale::Scale
2560 (const G4VisAttributes& visAtts,
2561 G4double length, const G4Transform3D& transform,
2562 const G4String& annotation, G4double annotationSize,
2563 const G4Colour& annotationColour):
2564fVisAtts(visAtts)
2565{
2566 // Useful constants...
2567 const G4double halfLength(length / 2.);
2568 const G4double tickLength(length / 20.);
2569
2570 // Create (empty) polylines having the same vis attributes...
2571 // (OK to pass address since fVisAtts is long lived.)
2572 fScaleLine.SetVisAttributes(&fVisAtts);
2573 fTick11.SetVisAttributes(&fVisAtts);
2574 fTick12.SetVisAttributes(&fVisAtts);
2575 fTick21.SetVisAttributes(&fVisAtts);
2576 fTick22.SetVisAttributes(&fVisAtts);
2577
2578 // Add points to the polylines to represent a scale parallel to the
2579 // x-axis centred on the origin...
2580 G4Point3D r1(G4Point3D(-halfLength, 0., 0.));
2581 G4Point3D r2(G4Point3D( halfLength, 0., 0.));
2582 fScaleLine.push_back(r1);
2583 fScaleLine.push_back(r2);
2584 G4Point3D ticky(0., tickLength, 0.);
2585 G4Point3D tickz(0., 0., tickLength);
2586 fTick11.push_back(r1 + ticky);
2587 fTick11.push_back(r1 - ticky);
2588 fTick12.push_back(r1 + tickz);
2589 fTick12.push_back(r1 - tickz);
2590 fTick21.push_back(r2 + ticky);
2591 fTick21.push_back(r2 - ticky);
2592 fTick22.push_back(r2 + tickz);
2593 fTick22.push_back(r2 - tickz);
2594 // ...and transform to chosen position and orientation
2595 fScaleLine.transform(transform);
2596 fTick11.transform(transform);
2597 fTick12.transform(transform);
2598 fTick21.transform(transform);
2599 fTick22.transform(transform);
2600 // Similarly for annotation
2601 G4Point3D textPosition(0., tickLength, 0.);
2602 textPosition.transform(transform);
2603 fText = G4Text(annotation,textPosition);
2604 fText.SetVisAttributes(annotationColour);
2605 fText.SetScreenSize(annotationSize);
2606}
2607
2608void G4VisCommandSceneAddScale::Scale::operator()
2609(G4VGraphicsScene& sceneHandler,const G4ModelingParameters*)
2610{
2611 // Draw...
2612 sceneHandler.BeginPrimitives();
2613 sceneHandler.AddPrimitive(fScaleLine);
2614 sceneHandler.AddPrimitive(fTick11);
2615 sceneHandler.AddPrimitive(fTick12);
2616 sceneHandler.AddPrimitive(fTick21);
2617 sceneHandler.AddPrimitive(fTick22);
2618 sceneHandler.AddPrimitive(fText);
2619 sceneHandler.EndPrimitives();
2620}
2621
2622////////////// /vis/scene/add/text //////////////////////////////////
2623
2625 G4bool omitable;
2626 fpCommand = new G4UIcommand ("/vis/scene/add/text", this);
2627 fpCommand -> SetGuidance ("Adds text to current scene.");
2628 fpCommand -> SetGuidance
2629 ("Use \"/vis/set/textColour\" to set colour.");
2630 fpCommand -> SetGuidance
2631 ("Use \"/vis/set/textLayout\" to set layout:");
2632 G4UIparameter* parameter;
2633 parameter = new G4UIparameter ("x", 'd', omitable = true);
2634 parameter->SetDefaultValue (0);
2635 fpCommand->SetParameter (parameter);
2636 parameter = new G4UIparameter ("y", 'd', omitable = true);
2637 parameter->SetDefaultValue (0);
2638 fpCommand->SetParameter (parameter);
2639 parameter = new G4UIparameter ("z", 'd', omitable = true);
2640 parameter->SetDefaultValue (0);
2641 fpCommand->SetParameter (parameter);
2642 parameter = new G4UIparameter ("unit", 's', omitable = true);
2643 parameter->SetDefaultValue ("m");
2644 fpCommand->SetParameter (parameter);
2645 parameter = new G4UIparameter ("font_size", 'd', omitable = true);
2646 parameter->SetDefaultValue (12);
2647 parameter->SetGuidance ("pixels");
2648 fpCommand->SetParameter (parameter);
2649 parameter = new G4UIparameter ("x_offset", 'd', omitable = true);
2650 parameter->SetDefaultValue (0);
2651 parameter->SetGuidance ("pixels");
2652 fpCommand->SetParameter (parameter);
2653 parameter = new G4UIparameter ("y_offset", 'd', omitable = true);
2654 parameter->SetDefaultValue (0);
2655 parameter->SetGuidance ("pixels");
2656 fpCommand->SetParameter (parameter);
2657 parameter = new G4UIparameter ("text", 's', omitable = true);
2658 parameter->SetGuidance ("The rest of the line is text.");
2659 parameter->SetDefaultValue ("Hello G4");
2660 fpCommand->SetParameter (parameter);
2661}
2662
2666
2670
2672
2673 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
2674 G4bool warn = verbosity >= G4VisManager::warnings;
2675
2676 G4Scene* pScene = fpVisManager->GetCurrentScene();
2677 if (!pScene) {
2678 if (verbosity >= G4VisManager::errors) {
2679 G4warn << "ERROR: No current scene. Please create one." << G4endl;
2680 }
2681 return;
2682 }
2683
2684 G4Tokenizer next(newValue);
2685 G4double x = StoD(next());
2686 G4double y = StoD(next());
2687 G4double z = StoD(next());
2688 G4String unitString = next();
2689 G4double font_size = StoD(next());
2690 G4double x_offset = StoD(next());
2691 G4double y_offset = StoD(next());
2692 G4String text = next("\n");
2693
2694 G4double unit = G4UIcommand::ValueOf(unitString);
2695 x *= unit; y *= unit; z *= unit;
2696
2697 G4Text g4text(text, G4Point3D(x,y,z));
2699 g4text.SetVisAttributes(visAtts);
2701 g4text.SetScreenSize(font_size);
2702 g4text.SetOffset(x_offset,y_offset);
2703 G4VModel* model = new G4TextModel(g4text);
2704 const G4String& currentSceneName = pScene -> GetName ();
2705 G4bool successful = pScene -> AddRunDurationModel (model, warn);
2706 if (successful) {
2707 if (verbosity >= G4VisManager::confirmations) {
2708 G4cout << "Text '" << text
2709 << "' has been added to scene \"" << currentSceneName << "\"."
2710 << G4endl;
2711 }
2712 }
2713 else G4VisCommandsSceneAddUnsuccessful(verbosity);
2714
2716}
2717
2718
2719////////////// /vis/scene/add/text2D //////////////////////////////////
2720
2722 G4bool omitable;
2723 fpCommand = new G4UIcommand ("/vis/scene/add/text2D", this);
2724 fpCommand -> SetGuidance ("Adds 2D text to current scene.");
2725 fpCommand -> SetGuidance ("x,y in range [-1,1]");
2726 fpCommand -> SetGuidance
2727 ("Use \"/vis/set/textColour\" to set colour.");
2728 fpCommand -> SetGuidance
2729 ("Use \"/vis/set/textLayout\" to set layout:");
2730 G4UIparameter* parameter;
2731 parameter = new G4UIparameter ("x", 'd', omitable = true);
2732 parameter->SetDefaultValue (0);
2733 fpCommand->SetParameter (parameter);
2734 parameter = new G4UIparameter ("y", 'd', omitable = true);
2735 parameter->SetDefaultValue (0);
2736 fpCommand->SetParameter (parameter);
2737 parameter = new G4UIparameter ("font_size", 'd', omitable = true);
2738 parameter->SetDefaultValue (12);
2739 parameter->SetGuidance ("pixels");
2740 fpCommand->SetParameter (parameter);
2741 parameter = new G4UIparameter ("x_offset", 'd', omitable = true);
2742 parameter->SetDefaultValue (0);
2743 parameter->SetGuidance ("pixels");
2744 fpCommand->SetParameter (parameter);
2745 parameter = new G4UIparameter ("y_offset", 'd', omitable = true);
2746 parameter->SetDefaultValue (0);
2747 parameter->SetGuidance ("pixels");
2748 fpCommand->SetParameter (parameter);
2749 parameter = new G4UIparameter ("text", 's', omitable = true);
2750 parameter->SetGuidance ("The rest of the line is text.");
2751 parameter->SetDefaultValue ("Hello G4");
2752 fpCommand->SetParameter (parameter);
2753}
2754
2758
2762
2764
2765 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
2766 G4bool warn = verbosity >= G4VisManager::warnings;
2767
2768 G4Scene* pScene = fpVisManager->GetCurrentScene();
2769 if (!pScene) {
2770 if (verbosity >= G4VisManager::errors) {
2771 G4warn << "ERROR: No current scene. Please create one." << G4endl;
2772 }
2773 return;
2774 }
2775
2776 G4Tokenizer next(newValue);
2777 G4double x = StoD(next());
2778 G4double y = StoD(next());
2779 G4double font_size = StoD(next());
2780 G4double x_offset = StoD(next());
2781 G4double y_offset = StoD(next());
2782 G4String text = next("\n");
2783
2784 G4Text g4text(text, G4Point3D(x,y,0.));
2786 g4text.SetVisAttributes(visAtts);
2788 g4text.SetScreenSize(font_size);
2789 g4text.SetOffset(x_offset,y_offset);
2790 G4Text2D* g4text2D = new G4Text2D(g4text);
2791 G4VModel* model =
2793 model->SetType("Text2D");
2794 model->SetGlobalTag("Text2D");
2795 std::ostringstream oss;
2796 oss << "Text2D: '" << g4text.GetText()
2797 << "' at " << g4text.GetPosition().x() << ',' << g4text.GetPosition().y()
2798 << " with size " << g4text.GetScreenSize()
2799 << " with offsets " << g4text.GetXOffset() << ',' << g4text.GetYOffset();
2800 model->SetGlobalDescription(oss.str());
2801 const G4String& currentSceneName = pScene -> GetName ();
2802 G4bool successful = pScene -> AddRunDurationModel (model, warn);
2803 if (successful) {
2804 if (verbosity >= G4VisManager::confirmations) {
2805 G4cout << "2D text '" << text
2806 << "' has been added to scene \"" << currentSceneName << "\"."
2807 << G4endl;
2808 }
2809 }
2810 else G4VisCommandsSceneAddUnsuccessful(verbosity);
2811
2813}
2814
2815G4VisCommandSceneAddText2D::G4Text2D::G4Text2D(const G4Text& text):
2816 fText(text)
2817{}
2818
2819void G4VisCommandSceneAddText2D::G4Text2D::operator()
2820 (G4VGraphicsScene& sceneHandler, const G4ModelingParameters*) {
2821 sceneHandler.BeginPrimitives2D();
2822 sceneHandler.AddPrimitive(fText);
2823 sceneHandler.EndPrimitives2D();
2824}
2825
2826
2827////////////// /vis/scene/add/trajectories ///////////////////////////////////
2828
2830 G4bool omitable;
2831 fpCommand = new G4UIcmdWithAString
2832 ("/vis/scene/add/trajectories", this);
2833 fpCommand -> SetGuidance
2834 ("Adds trajectories to current scene.");
2835 fpCommand -> SetGuidance
2836 ("Causes trajectories, if any, to be drawn at the end of processing an"
2837 "\nevent. Switches on trajectory storing and sets the"
2838 "\ndefault trajectory type.");
2839 fpCommand -> SetGuidance
2840 ("The command line parameter list determines the default trajectory type."
2841 "\nIf it contains the string \"smooth\", auxiliary inter-step points will"
2842 "\nbe inserted to improve the smoothness of the drawing of a curved"
2843 "\ntrajectory."
2844 "\nIf it contains the string \"rich\", significant extra information will"
2845 "\nbe stored in the trajectory (G4RichTrajectory) amenable to modeling"
2846 "\nand filtering with \"/vis/modeling/trajectories/create/drawByAttribute\""
2847 "\nand \"/vis/filtering/trajectories/create/attributeFilter\" commands."
2848 "\nIt may contain both strings in any order.");
2849 fpCommand -> SetGuidance
2850 ("\nTo switch off trajectory storing: \"/tracking/storeTrajectory 0\"."
2851 "\nSee also \"/vis/scene/endOfEventAction\".");
2852 fpCommand -> SetGuidance
2853 ("Note: This only sets the default. Independently of the result of this"
2854 "\ncommand, a user may instantiate a trajectory that overrides this default"
2855 "\nin PreUserTrackingAction.");
2856 fpCommand -> SetParameterName ("default-trajectory-type", omitable = true);
2857 fpCommand -> SetDefaultValue ("");
2858}
2859
2863
2867
2869 G4String newValue) {
2870
2871 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
2872 G4bool warn = verbosity >= G4VisManager::warnings;
2873
2874 G4Scene* pScene = fpVisManager->GetCurrentScene();
2875 if (!pScene) {
2876 if (verbosity >= G4VisManager::errors) {
2877 G4warn << "ERROR: No current scene. Please create one." << G4endl;
2878 }
2879 return;
2880 }
2881 const G4String& currentSceneName = pScene -> GetName ();
2882
2883 G4bool smooth = false;
2884 G4bool rich = false;
2885 if (newValue.find("smooth") != std::string::npos) smooth = true;
2886 if (newValue.find("rich") != std::string::npos) rich = true;
2887 if (newValue.size() && !(rich || smooth)) {
2888 if (verbosity >= G4VisManager::errors) {
2889 G4warn << "ERROR: Unrecognised parameter \"" << newValue << "\""
2890 "\n No action taken."
2891 << G4endl;
2892 }
2893 return;
2894 }
2895
2897 G4String defaultTrajectoryType;
2898 if (smooth && rich) {
2899 UImanager->ApplyCommand("/tracking/storeTrajectory 4");
2900 defaultTrajectoryType = "G4RichTrajectory configured for smooth steps";
2901 } else if (smooth) {
2902 UImanager->ApplyCommand("/tracking/storeTrajectory 2");
2903 defaultTrajectoryType = "G4SmoothTrajectory";
2904 } else if (rich) {
2905 UImanager->ApplyCommand("/tracking/storeTrajectory 3");
2906 defaultTrajectoryType = "G4RichTrajectory";
2907 } else {
2908 UImanager->ApplyCommand("/tracking/storeTrajectory 1");
2909 defaultTrajectoryType = "G4Trajectory";
2910 }
2911
2912 if (verbosity >= G4VisManager::errors) {
2913 G4warn <<
2914 "Attributes available for modeling and filtering with"
2915 "\n \"/vis/modeling/trajectories/create/drawByAttribute\" and"
2916 "\n \"/vis/filtering/trajectories/create/attributeFilter\" commands:"
2917 << G4endl;
2919 if (rich) {
2922 } else if (smooth) {
2925 } else {
2928 }
2929 }
2930
2931 const auto& eoeList = pScene->GetEndOfEventModelList();
2932 auto eoeModel = eoeList.begin();
2933 for (; eoeModel != eoeList.end(); ++eoeModel) {
2934 const auto* actualModel = eoeModel->fpModel;
2935 if (dynamic_cast<const G4TrajectoriesModel*>(actualModel)) break;
2936 }
2937 if (eoeModel == eoeList.end()) {
2938 // No trajectories model exists in the scene so create a new one...
2939 G4VModel* model = new G4TrajectoriesModel();
2940 pScene -> AddEndOfEventModel (model, warn);
2941 } // ...else it already exists and there is no need to add a new one
2942 // because G4TrajectoriesModel simply describes trajectories in the
2943 // trajectories store whatever the type.
2944
2945 if (verbosity >= G4VisManager::confirmations) {
2946 G4cout << "Default trajectory type " << defaultTrajectoryType
2947 << "\n will be used to store trajectories for scene \""
2948 << currentSceneName << "\"."
2949 << G4endl;
2950 }
2951
2952 if (verbosity >= G4VisManager::warnings) {
2953 G4warn <<
2954 "WARNING: Trajectory storing has been requested. This action may be"
2955 "\n reversed with \"/tracking/storeTrajectory 0\"."
2956 << G4endl;
2957 }
2958
2960}
2961
2962////////////// /vis/scene/add/userAction ///////////////////////////////////
2963
2965 G4bool omitable;
2966 fpCommand = new G4UIcmdWithAString("/vis/scene/add/userAction",this);
2967 fpCommand -> SetGuidance
2968 ("Add named Vis User Action to current scene.");
2969 fpCommand -> SetGuidance
2970 ("Attempts to match search string to name of action - use unique sub-string.");
2971 fpCommand -> SetGuidance
2972 ("(Use /vis/list to see names of registered actions.)");
2973 fpCommand -> SetGuidance
2974 ("If name == \"all\" (default), all actions are added.");
2975 fpCommand -> SetParameterName("action-name", omitable = true);
2976 fpCommand -> SetDefaultValue("all");
2977}
2978
2982
2986
2988(G4UIcommand*, G4String newValue) {
2989
2990 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
2991
2992 G4Scene* pScene = fpVisManager->GetCurrentScene();
2993 if (!pScene) {
2994 if (verbosity >= G4VisManager::errors) {
2995 G4warn << "ERROR: No current scene. Please create one." << G4endl;
2996 }
2997 return;
2998 }
2999
3000 G4bool any = false;
3001
3002 const std::vector<G4VisManager::UserVisAction>& runDurationUserVisActions =
3003 fpVisManager->GetRunDurationUserVisActions();
3004 for (const auto& runDurationUserVisAction : runDurationUserVisActions) {
3005 const G4String& name = runDurationUserVisAction.fName;
3006 G4VUserVisAction* visAction = runDurationUserVisAction.fpUserVisAction;
3007 if (newValue == "all" || name.find(newValue) != std::string::npos) {
3008 any = true;
3009 AddVisAction(name,visAction,pScene,runDuration,verbosity);
3010 }
3011 }
3012
3013 const std::vector<G4VisManager::UserVisAction>& endOfEventUserVisActions =
3014 fpVisManager->GetEndOfEventUserVisActions();
3015 for (const auto& endOfEventUserVisAction : endOfEventUserVisActions) {
3016 const G4String& name = endOfEventUserVisAction.fName;
3017 G4VUserVisAction* visAction = endOfEventUserVisAction.fpUserVisAction;
3018 if (newValue == "all" || name.find(newValue) != std::string::npos) {
3019 any = true;
3020 AddVisAction(name,visAction,pScene,endOfEvent,verbosity);
3021 }
3022 }
3023
3024 const std::vector<G4VisManager::UserVisAction>& endOfRunUserVisActions =
3025 fpVisManager->GetEndOfRunUserVisActions();
3026 for (const auto& endOfRunUserVisAction : endOfRunUserVisActions) {
3027 const G4String& name = endOfRunUserVisAction.fName;
3028 G4VUserVisAction* visAction = endOfRunUserVisAction.fpUserVisAction;
3029 if (newValue == "all" || name.find(newValue) != std::string::npos) {
3030 any = true;
3031 AddVisAction(name,visAction,pScene,endOfRun,verbosity);
3032 }
3033 }
3034
3035 if (!any) {
3036 if (verbosity >= G4VisManager::warnings) {
3037 G4warn << "WARNING: No User Vis Action registered." << G4endl;
3038 }
3039 return;
3040 }
3041
3043}
3044
3045void G4VisCommandSceneAddUserAction::AddVisAction
3046(const G4String& name,
3047 G4VUserVisAction* visAction,
3048 G4Scene* pScene,
3049 G4VisCommandSceneAddUserAction::ActionType type,
3050 G4VisManager::Verbosity verbosity)
3051{
3052 G4bool warn = verbosity >= G4VisManager::warnings;
3053
3054 const std::map<G4VUserVisAction*,G4VisExtent>& visExtentMap =
3056 G4VisExtent extent;
3057 std::map<G4VUserVisAction*,G4VisExtent>::const_iterator i =
3058 visExtentMap.find(visAction);
3059 if (i != visExtentMap.end()) extent = i->second;
3060 if (warn) {
3061 if (extent.GetExtentRadius() <= 0.) {
3062 G4warn
3063 << "WARNING: User Vis Action \"" << name << "\" extent is null."
3064 << G4endl;
3065 }
3066 }
3067
3068 G4VModel* model = new G4CallbackModel<G4VUserVisAction>(visAction);
3069 model->SetType("User Vis Action");
3070 model->SetGlobalTag(name);
3071 model->SetGlobalDescription(name);
3072 model->SetExtent(extent);
3073 G4bool successful = false;;
3074 switch (type) {
3075 case runDuration:
3076 successful = pScene -> AddRunDurationModel (model, warn);
3077 break;
3078 case endOfEvent:
3079 successful = pScene -> AddEndOfEventModel (model, warn);
3080 break;
3081 case endOfRun:
3082 successful = pScene -> AddEndOfRunModel (model, warn);
3083 break;
3084 }
3085 if (successful) {
3086 if (verbosity >= G4VisManager::confirmations) {
3087 const G4String& currentSceneName = pScene -> GetName ();
3088 G4cout << "User Vis Action added to scene \""
3089 << currentSceneName << "\"";
3090 if (verbosity >= G4VisManager::parameters) {
3091 G4cout << "\n with extent " << extent;
3092 }
3093 G4cout << G4endl;
3094 }
3095 }
3096 else G4VisCommandsSceneAddUnsuccessful(verbosity);
3097}
3098
3099////////////// /vis/scene/add/volume ///////////////////////////////////////
3100
3102 G4bool omitable;
3103 fpCommand = new G4UIcommand ("/vis/scene/add/volume", this);
3104 fpCommand -> SetGuidance
3105 ("Adds a physical volume to current scene, with optional clipping volume.");
3106 fpCommand -> SetGuidance
3107 ("If physical-volume-name is \"world\" (the default), the top of the"
3108 "\nmain geometry tree (material world) is added. If \"worlds\", the"
3109 "\ntops of all worlds - material world and parallel worlds, if any - are"
3110 "\nadded. Otherwise a search of all worlds is made.");
3111 fpCommand -> SetGuidance
3112 ("In the last case the names of all volumes in all worlds are matched"
3113 "\nagainst physical-volume-name. If this is of the form \"/regexp/\","
3114 "\nwhere regexp is a regular expression (see C++ regex), the match uses"
3115 "\nthe usual rules of regular expression matching. Otherwise an exact"
3116 "\nmatch is required."
3117 "\nFor example, \"/Shap/\" adds \"Shape1\" and \"Shape2\".");
3118 fpCommand -> SetGuidance
3119 ("It may help to see a textual representation of the geometry hierarchy of"
3120 "\nthe worlds. Try \"/vis/drawTree [worlds]\".");
3121 fpCommand -> SetGuidance
3122 ("If clip-volume-type is specified, the subsequent parameters are used to"
3123 "\nto define a clipping volume. For example,"
3124 "\n\"/vis/scene/add/volume ! ! ! -box km 0 1 0 1 0 1\" will draw the world"
3125 "\nwith the positive octant cut away. (If the Boolean Processor issues"
3126 "\nwarnings try replacing 0 by 0.000000001 or something.)");
3127 fpCommand -> SetGuidance
3128 ("If clip-volume-type is prepended with '-', the clip-volume is subtracted"
3129 "\n(cutaway). (This is the default if there is no prepended character.)"
3130 "\nIf '*' is prepended, the intersection of the physical-volume and the"
3131 "\nclip-volume is made. (You can make a section through the detector with"
3132 "\na thin box, for example).");
3133 fpCommand -> SetGuidance
3134 ("For \"box\", the parameters are xmin,xmax,ymin,ymax,zmin,zmax."
3135 "\nOnly \"box\" is programmed at present.");
3136 G4UIparameter* parameter;
3137 parameter = new G4UIparameter ("physical-volume-name", 's', omitable = true);
3138 parameter -> SetDefaultValue ("world");
3139 fpCommand -> SetParameter (parameter);
3140 parameter = new G4UIparameter ("copy-no", 'i', omitable = true);
3141 parameter -> SetGuidance ("If negative, matches any copy no.");
3142 parameter -> SetDefaultValue (-1);
3143 fpCommand -> SetParameter (parameter);
3144 parameter = new G4UIparameter ("depth-of-descent", 'i', omitable = true);
3145 parameter -> SetGuidance
3146 ("Depth of descent of geometry hierarchy. Default = unlimited depth.");
3147 parameter -> SetDefaultValue (G4PhysicalVolumeModel::UNLIMITED);
3148 fpCommand -> SetParameter (parameter);
3149 parameter = new G4UIparameter ("clip-volume-type", 's', omitable = true);
3150 parameter -> SetParameterCandidates("none box -box *box");
3151 parameter -> SetDefaultValue ("none");
3152 parameter -> SetGuidance("[-|*]type. See general guidance.");
3153 fpCommand -> SetParameter (parameter);
3154 parameter = new G4UIparameter ("parameter-unit", 's', omitable = true);
3155 parameter -> SetDefaultValue ("m");
3156 fpCommand -> SetParameter (parameter);
3157 parameter = new G4UIparameter ("parameter-1", 'd', omitable = true);
3158 parameter -> SetDefaultValue (0.);
3159 fpCommand -> SetParameter (parameter);
3160 parameter = new G4UIparameter ("parameter-2", 'd', omitable = true);
3161 parameter -> SetDefaultValue (0.);
3162 fpCommand -> SetParameter (parameter);
3163 parameter = new G4UIparameter ("parameter-3", 'd', omitable = true);
3164 parameter -> SetDefaultValue (0.);
3165 fpCommand -> SetParameter (parameter);
3166 parameter = new G4UIparameter ("parameter-4", 'd', omitable = true);
3167 parameter -> SetDefaultValue (0.);
3168 fpCommand -> SetParameter (parameter);
3169 parameter = new G4UIparameter ("parameter-5", 'd', omitable = true);
3170 parameter -> SetDefaultValue (0.);
3171 fpCommand -> SetParameter (parameter);
3172 parameter = new G4UIparameter ("parameter-6", 'd', omitable = true);
3173 parameter -> SetDefaultValue (0.);
3174 fpCommand -> SetParameter (parameter);
3175}
3176
3180
3184
3186 G4String newValue) {
3187
3188 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
3189 G4bool warn = verbosity >= G4VisManager::warnings;
3190
3191 G4Scene* pScene = fpVisManager->GetCurrentScene();
3192 if (!pScene) {
3193 if (verbosity >= G4VisManager::errors) {
3194 G4warn << "ERROR: No current scene. Please create one." << G4endl;
3195 }
3196 return;
3197 }
3198
3199 G4String name, clipVolumeType, parameterUnit;
3200 G4int copyNo, requestedDepthOfDescent;
3201 G4double param1, param2, param3, param4, param5, param6;
3202 std::istringstream is (newValue);
3203 is >> name >> copyNo >> requestedDepthOfDescent
3204 >> clipVolumeType >> parameterUnit
3205 >> param1 >> param2 >> param3 >> param4 >> param5 >> param6;
3207 G4PhysicalVolumeModel::subtraction; // Default subtraction mode.
3208 if (clipVolumeType[size_t(0)] == '-') {
3209 clipVolumeType = clipVolumeType.substr(1); // Remove first character.
3210 } else if (clipVolumeType[size_t(0)] == '*') {
3212 clipVolumeType = clipVolumeType.substr(1);
3213 }
3214 G4double unit = G4UIcommand::ValueOf(parameterUnit);
3215 param1 *= unit; param2 *= unit; param3 *= unit;
3216 param4 *= unit; param5 *= unit; param6 *= unit;
3217
3218 G4VSolid* clippingSolid = nullptr;
3219 if (clipVolumeType == "box") {
3220 const G4double dX = (param2 - param1) / 2.;
3221 const G4double dY = (param4 - param3) / 2.;
3222 const G4double dZ = (param6 - param5) / 2.;
3223 const G4double x0 = (param2 + param1) / 2.;
3224 const G4double y0 = (param4 + param3) / 2.;
3225 const G4double z0 = (param6 + param5) / 2.;
3226 clippingSolid = new G4DisplacedSolid
3227 ("_displaced_clipping_box",
3228 new G4Box("_clipping_box",dX,dY,dZ),
3229 G4Translate3D(x0,y0,z0));
3230 }
3231
3232 G4TransportationManager* transportationManager =
3234
3235 size_t nWorlds = transportationManager->GetNoWorlds();
3236 if (nWorlds > 1) { // Parallel worlds in operation...
3237 if (verbosity >= G4VisManager::warnings) {
3238 static G4bool warned = false;
3239 if (!warned && name != "worlds") {
3240 G4warn <<
3241 "WARNING: Parallel worlds in operation. To visualise, specify"
3242 "\n \"worlds\" or the parallel world volume or sub-volume name"
3243 "\n and control visibility with /vis/geometry."
3244 << G4endl;
3245 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
3246 transportationManager->GetWorldsIterator();
3247 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
3248 G4warn << " World " << i << ": " << (*iterWorld)->GetName()
3249 << G4endl;
3250 warned = true;
3251 }
3252 }
3253 }
3254 }
3255
3256 // Get the world (the initial value of the iterator points to the mass world).
3257 G4VPhysicalVolume* world = *(transportationManager->GetWorldsIterator());
3258
3259 if (!world) {
3260 if (verbosity >= G4VisManager::errors) {
3261 G4warn <<
3262 "ERROR: G4VisCommandSceneAddVolume::SetNewValue:"
3263 "\n No world. Maybe the geometry has not yet been defined."
3264 "\n Try \"/run/initialize\""
3265 << G4endl;
3266 }
3267 return;
3268 }
3269
3270 std::vector<G4PhysicalVolumesSearchScene::Findings> findingsVector;
3271
3272 // When it comes to determining the extent of a physical volume we normally
3273 // assume the user wishes to ignore "invisible" volumes. For example, most
3274 // users make the world volume invisible. So we ask the physical volume
3275 // model to traverse the geometry hierarchy, starting at the named physical
3276 // volume, until it finds non-invisible ones, whose extents are accumulated
3277 // to determine the overall extent. (Once a non-invisible volume is found,
3278 // the search is curtailed - daughters are always contained within the mother
3279 // so they have no subsequent influence on the extent of the mother - but the
3280 // search continues at the same level until all highest level non-invisible
3281 // volumes are found an their extents accumulated.) So the default is
3282 G4bool useFullExtent = false;
3283 // However, the above procedure can be time consuming in some situations, such
3284 // as a nested parameterisation whose ultimate volumes are the first non-
3285 // visible ones, which are typical of a medical "phantom". So we assume here
3286 // below that if a user specifies a name other than "world" or "worlds" he/she
3287 // wished the extent to be determined by the volume, whether it is visible
3288 // or not. So we set useFullExtent true at that point below.
3289
3290 if (name == "world") {
3291
3292 findingsVector.push_back
3294
3295 } else if (name == "worlds") {
3296
3297 if (nWorlds <= 1) {
3298 if (verbosity >= G4VisManager::warnings) {
3299 G4warn <<
3300 "WARNING: G4VisCommandSceneAddVolume::SetNewValue:"
3301 "\n Parallel worlds requested but none exist."
3302 "\n Just adding material world."
3303 << G4endl;
3304 }
3305 }
3306 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
3307 transportationManager->GetWorldsIterator();
3308 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
3309 findingsVector.push_back
3311 (*iterWorld,*iterWorld));
3312 }
3313
3314 } else { // Search all worlds...
3315
3316 // Use the model's full extent. This assumes the user wants these
3317 // volumes in the findings vector (there could be more than one) to
3318 // determine the scene's extent. Otherwise G4PhysicalVolumeModel would
3319 // re-calculate each volume's extent based on visibility, etc., which
3320 // could be time consuming.
3321 useFullExtent = true;
3322
3323 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
3324 transportationManager->GetWorldsIterator();
3325 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
3326 G4ModelingParameters mp; // Default - no culling.
3327 G4PhysicalVolumeModel searchModel
3328 (*iterWorld,
3329 requestedDepthOfDescent,
3330 G4Transform3D(),
3331 &mp,
3332 useFullExtent);
3333 G4PhysicalVolumesSearchScene searchScene(&searchModel, name, copyNo);
3334 searchModel.DescribeYourselfTo (searchScene); // Initiate search.
3335 for (const auto& findings: searchScene.GetFindings()) {
3336 findingsVector.push_back(findings);
3337 }
3338 }
3339 }
3340
3341 for (const auto& findings: findingsVector) {
3342 // Set copy number from search findings for replicas and parameterisations.
3343 findings.fpFoundPV->SetCopyNo(findings.fFoundPVCopyNo);
3345 (findings.fpFoundPV,
3346 requestedDepthOfDescent,
3347 findings.fFoundObjectTransformation,
3348 0, // No modelling parameters (these are set later by the scene handler).
3349 useFullExtent,
3350 findings.fFoundBasePVPath);
3351 if (clippingSolid) {
3352 foundPVModel->SetClippingSolid(clippingSolid);
3353 foundPVModel->SetClippingMode(clippingMode);
3354 }
3355 if (!foundPVModel->Validate(warn)) return;
3356 // ...so add it to the scene.
3357 G4bool successful = pScene->AddRunDurationModel(foundPVModel,warn);
3358 if (successful) {
3359 if (verbosity >= G4VisManager::confirmations) {
3360 G4cout << "\"" << findings.fpFoundPV->GetName()
3361 << "\", copy no. " << findings.fFoundPVCopyNo
3362 << ",\n found in searched volume \""
3363 << findings.fpSearchPV->GetName()
3364 << "\" at depth " << findings.fFoundDepth
3365 << ",\n base path: \"" << findings.fFoundBasePVPath
3366 << "\",\n with a requested depth of further descent of ";
3367 if (requestedDepthOfDescent < 0) {
3368 G4cout << "<0 (unlimited)";
3369 }
3370 else {
3371 G4cout << requestedDepthOfDescent;
3372 }
3373 G4cout << ",\n has been added to scene \"" << pScene->GetName() << "\"."
3374 << G4endl;
3375 }
3376 } else {
3378 }
3379 }
3380
3381 if (findingsVector.empty()) {
3382 if (verbosity >= G4VisManager::errors) {
3383 G4warn << "ERROR: Volume \"" << name << "\"";
3384 if (copyNo >= 0) {
3385 G4warn << ", copy no. " << copyNo << ",";
3386 }
3387 G4warn << " not found." << G4endl;
3388 }
3390 return;
3391 }
3392
3394}
3395
3396/////////////////////////////////////////////////////////////////////////////
3397////////////// /vis/scene/add/plotter ///////////////////////////////////////
3398/////////////////////////////////////////////////////////////////////////////
3400 fpCommand = new G4UIcommand("/vis/scene/add/plotter", this);
3401 fpCommand -> SetGuidance ("Add a plotter to current scene.");
3402
3403 G4UIparameter* parameter;
3404 parameter = new G4UIparameter ("plotter", 's',false);
3405 fpCommand->SetParameter(parameter);
3406}
3407
3409
3411
3413{
3414 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
3415 G4bool warn(verbosity >= G4VisManager::warnings);
3416
3417 G4Scene* pScene = fpVisManager->GetCurrentScene();
3418 if (!pScene) {
3419 if (verbosity >= G4VisManager::errors) {
3420 G4warn << "ERROR: No current scene. Please create one." << G4endl;
3421 }
3422 return;
3423 }
3424
3425 G4Plotter& _plotter = G4PlotterManager::GetInstance().GetPlotter(newValue);
3426 G4VModel* model = new G4PlotterModel(_plotter,newValue);
3427
3428 const G4String& currentSceneName = pScene -> GetName ();
3429 G4bool successful = pScene -> AddEndOfRunModel(model, warn);
3430 if (successful) {
3431 if (verbosity >= G4VisManager::confirmations) {
3432 G4cout
3433 << "Plotter \"" << model->GetCurrentDescription()
3434 << "\" has been added to scene \"" << currentSceneName << "\"."
3435 << G4endl;
3436 }
3437 }
3438 else G4VisCommandsSceneAddUnsuccessful(verbosity);
3439
3441}
3442
HepGeom::Point3D< G4double > G4Point3D
Definition G4Point3D.hh:34
CLHEP::HepRotation G4RotationMatrix
#define G4warn
Definition G4Scene.cc:41
#define G4BestUnit(a, b)
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Transform3D G4Transform3D
HepGeom::RotateZ3D G4RotateZ3D
HepGeom::RotateX3D G4RotateX3D
HepGeom::RotateY3D G4RotateY3D
HepGeom::Translate3D G4Translate3D
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 G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
HepRotation & rotateZ(double delta)
Definition Rotation.cc:87
G4Box is a cuboid of given half lengths dx,dy,dz. The Box is centred on the origin with sides paralle...
Definition G4Box.hh:58
static G4Colour Green()
Definition G4Colour.hh:181
static G4Colour Red()
Definition G4Colour.hh:180
static G4Colour Brown()
Definition G4Colour.hh:179
static G4Colour Blue()
Definition G4Colour.hh:182
G4DisplacedSolid is a solid that has been shifted from its original frame of reference to a new one....
G4int GetEventID() const
Definition G4Event.hh:126
G4LogicalVolumeStore is a singleton class, acting as container for all logical volumes,...
G4LogicalVolume * GetVolume(const G4String &name, G4bool verbose=true, G4bool reverseSearch=false) const
static G4LogicalVolumeStore * GetInstance()
G4LogicalVolume represents a leaf node or unpositioned subtree in the geometry hierarchy....
void SetClippingSolid(G4VSolid *pClippingSolid)
void SetClippingMode(ClippingMode mode)
void DescribeYourselfTo(G4VGraphicsScene &)
const std::vector< Findings > & GetFindings() const
G4Plotter & GetPlotter(const G4String &a_name)
static G4PlotterManager & GetInstance()
const std::map< G4String, G4AttDef > * GetAttDefs() const override
const std::map< G4String, G4AttDef > * GetAttDefs() const override
static G4RunManager * GetMasterRunManager()
const G4Run * GetCurrentRun() const
Definition G4Run.hh:48
G4int GetRunID() const
Definition G4Run.hh:82
G4int GetNumberOfEventToBeProcessed() const
Definition G4Run.hh:88
G4int GetNumberOfKeptEvents() const
Definition G4Run.cc:79
G4bool AddRunDurationModel(G4VModel *, G4bool warn=false)
Definition G4Scene.cc:160
const G4VisExtent & GetExtent() const
const G4String & GetName() const
const std::vector< Model > & GetEndOfEventModelList() const
const std::map< G4String, G4AttDef > * GetAttDefs() const override
const std::map< G4String, G4AttDef > * GetAttDefs() const override
G4SubtractionSolid is a solid describing the Boolean subtraction of two solids.
void SetLayout(Layout)
G4double GetYOffset() const
G4double GetXOffset() const
void SetOffset(double dx, double dy)
G4String GetText() const
@ centre
Definition G4Text.hh:76
@ right
Definition G4Text.hh:76
@ left
Definition G4Text.hh:76
const std::map< G4String, G4AttDef > * GetAttDefs() const
const std::map< G4String, G4AttDef > * GetAttDefs() const override
const std::map< G4String, G4AttDef > * GetAttDefs() const override
G4TransportationManager is a singleton class which stores the navigator used by the transportation pr...
static G4TransportationManager * GetTransportationManager()
std::vector< G4VPhysicalVolume * >::iterator GetWorldsIterator()
std::size_t GetNoWorlds() const
G4Tubs is a tube or tube segment with curved sides parallel to the Z-axis. The tube has a specified h...
Definition G4Tubs.hh:85
G4UIcommand * FindPath(const char *commandPath) const
static G4double ValueOf(const char *unitName)
static G4bool ConvertToBool(const char *st)
G4UIcommandTree * GetTree() const
G4int ApplyCommand(const char *aCommand)
static G4UImanager * GetUIpointer()
G4double StoD(const G4String &s)
void SetDefaultValue(const char *theDefaultValue)
void SetGuidance(const char *theGuidance)
G4UnionSolid is a solid describing the Boolean union of two solids.
G4double GetScreenSize() const
void SetScreenSize(G4double)
G4Point3D GetPosition() const
void SetType(const G4String &)
void SetGlobalDescription(const G4String &)
virtual G4String GetCurrentDescription() const
Definition G4VModel.cc:58
void SetGlobalTag(const G4String &)
void SetExtent(const G4VisExtent &)
const G4VisExtent & GetExtent() const
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4VSolid is an abstract base class for solids, physical shapes that can be tracked through....
Definition G4VSolid.hh:80
const G4ViewParameters & GetViewParameters() const
static G4double fCurrentTextSize
void G4VisCommandsSceneAddUnsuccessful(G4VisManager::Verbosity verbosity)
static G4Colour fCurrentTextColour
void CheckSceneAndNotifyHandlers(G4Scene *=nullptr)
static std::vector< G4PhysicalVolumesSearchScene::Findings > fCurrrentPVFindingsForField
void ConvertToColour(G4Colour &colour, const G4String &redOrString, G4double green, G4double blue, G4double opacity)
static G4VisManager * fpVisManager
static G4VisExtent fCurrentExtentForField
const G4String & ConvertToColourGuidance()
void CopyParametersFrom(const G4UIcommand *fromCmd, G4UIcommand *toCmd)
static G4int fCurrentArrow3DLineSegmentsPerCircle
static G4Text::Layout fCurrentTextLayout
static G4double fCurrentLineWidth
static G4Colour fCurrentColour
void CopyGuidanceFrom(const G4UIcommand *fromCmd, G4UIcommand *toCmd, G4int startLine=0)
G4int GetNoOfSides() const
const G4Vector3D & GetViewpointDirection() const
const G4Vector3D & GetUpVector() const
void SetColour(const G4Colour &)
void SetLineWidth(G4double)
void SetForceSolid(G4bool=true)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4String GetCurrentValue(G4UIcommand *command)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4double GetYmin() const
G4double GetXmax() const
G4double GetExtentRadius() const
G4VisExtent & Transform(const G4Transform3D &)
G4double GetYmax() const
G4double GetZmax() const
G4double GetZmin() const
G4double GetXmin() const
const std::map< G4VUserVisAction *, G4VisExtent > & GetUserVisActionExtents() const
G4VViewer * GetCurrentViewer() const
void SetVisAttributes(const G4VisAttributes *)
Definition G4Visible.cc:108
BasicVector3D< T > unit() const
CLHEP::HepRotation getRotation() const
CLHEP::Hep3Vector getTranslation() const
static void SetNumberOfRotationSteps(G4int n)
HepPolyhedron & Transform(const G4Transform3D &t)
static void ResetNumberOfRotationSteps()