Geant4 10.7.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 "G4ScaleModel.hh"
44#include "G4TextModel.hh"
45#include "G4ArrowModel.hh"
46#include "G4AxesModel.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
82#include <sstream>
83
84////////////// /vis/scene/add/arrow ///////////////////////////////////////
85
87 fpCommand = new G4UIcommand("/vis/scene/add/arrow", this);
88 fpCommand -> SetGuidance ("Adds arrow to current scene.");
89 G4bool omitable;
90 G4UIparameter* parameter;
91 parameter = new G4UIparameter ("x1", 'd', omitable = false);
92 fpCommand -> SetParameter (parameter);
93 parameter = new G4UIparameter ("y1", 'd', omitable = false);
94 fpCommand -> SetParameter (parameter);
95 parameter = new G4UIparameter ("z1", 'd', omitable = false);
96 fpCommand -> SetParameter (parameter);
97 parameter = new G4UIparameter ("x2", 'd', omitable = false);
98 fpCommand -> SetParameter (parameter);
99 parameter = new G4UIparameter ("y2", 'd', omitable = false);
100 fpCommand -> SetParameter (parameter);
101 parameter = new G4UIparameter ("z2", 'd', omitable = false);
102 fpCommand -> SetParameter (parameter);
103 parameter = new G4UIparameter ("unit", 's', omitable = true);
104 parameter->SetDefaultValue ("m");
105 fpCommand->SetParameter (parameter);
106}
107
109 delete fpCommand;
110}
111
113 return "";
114}
115
117{
119 G4bool warn(verbosity >= G4VisManager::warnings);
120
122 if (!pScene) {
123 if (verbosity >= G4VisManager::errors) {
124 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
125 }
126 return;
127 }
128
129 G4String unitString;
130 G4double x1, y1, z1, x2, y2, z2;
131 std::istringstream is(newValue);
132 is >> x1 >> y1 >> z1 >> x2 >> y2 >> z2 >> unitString;
133 G4double unit = G4UIcommand::ValueOf(unitString);
134 x1 *= unit; y1 *= unit; z1 *= unit;
135 x2 *= unit; y2 *= unit; z2 *= unit;
136
137 // Consult scene for arrow width.
138 const G4VisExtent& sceneExtent = pScene->GetExtent();
139 G4double arrowWidth =
140 0.005 * fCurrentLineWidth * sceneExtent.GetExtentRadius();
141
142 G4VModel* model = new G4ArrowModel
143 (x1, y1, z1, x2, y2, z2,
144 arrowWidth, fCurrentColour, newValue,
146
147 const G4String& currentSceneName = pScene -> GetName ();
148 G4bool successful = pScene -> AddRunDurationModel (model, warn);
149 if (successful) {
150 if (verbosity >= G4VisManager::confirmations) {
151 G4cout << "Arrow has been added to scene \""
152 << currentSceneName << "\"."
153 << G4endl;
154 }
155 }
156 else G4VisCommandsSceneAddUnsuccessful(verbosity);
157
159}
160
161////////////// /vis/scene/add/arrow2D ///////////////////////////////////////
162
164 fpCommand = new G4UIcommand("/vis/scene/add/arrow2D", this);
165 fpCommand -> SetGuidance ("Adds 2D arrow to current scene.");
166 G4bool omitable;
167 G4UIparameter* parameter;
168 parameter = new G4UIparameter ("x1", 'd', omitable = false);
169 fpCommand -> SetParameter (parameter);
170 parameter = new G4UIparameter ("y1", 'd', omitable = false);
171 fpCommand -> SetParameter (parameter);
172 parameter = new G4UIparameter ("x2", 'd', omitable = false);
173 fpCommand -> SetParameter (parameter);
174 parameter = new G4UIparameter ("y2", 'd', omitable = false);
175 fpCommand -> SetParameter (parameter);
176}
177
179 delete fpCommand;
180}
181
183 return "";
184}
185
187{
189 G4bool warn(verbosity >= G4VisManager::warnings);
190
192 if (!pScene) {
193 if (verbosity >= G4VisManager::errors) {
194 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
195 }
196 return;
197 }
198
199 G4double x1, y1, x2, y2;
200 std::istringstream is(newValue);
201 is >> x1 >> y1 >> x2 >> y2;
202
203 Arrow2D* arrow2D = new Arrow2D
204 (x1, y1, x2, y2, fCurrentLineWidth, fCurrentColour);
205 G4VModel* model =
207 model->SetType("Arrow2D");
208 model->SetGlobalTag("Arrow2D");
209 model->SetGlobalDescription("Arrow2D: " + newValue);
210 const G4String& currentSceneName = pScene -> GetName ();
211 G4bool successful = pScene -> AddRunDurationModel (model, warn);
212 if (successful) {
213 if (verbosity >= G4VisManager::confirmations) {
214 G4cout << "A 2D arrow has been added to scene \""
215 << currentSceneName << "\"."
216 << G4endl;
217 }
218 }
219 else G4VisCommandsSceneAddUnsuccessful(verbosity);
220
222}
223
224G4VisCommandSceneAddArrow2D::Arrow2D::Arrow2D
225(G4double x1, G4double y1,
226 G4double x2, G4double y2,
227 G4double width, const G4Colour& colour):
228 fWidth(width), fColour(colour)
229{
230 fShaftPolyline.push_back(G4Point3D(x1,y1,0));
231 fShaftPolyline.push_back(G4Point3D(x2,y2,0));
232 G4Vector3D arrowDirection = G4Vector3D(x2-x1,y2-y1,0).unit();
233 G4Vector3D arrowPointLeftDirection(arrowDirection);
234 arrowPointLeftDirection.rotateZ(150.*deg);
235 G4Vector3D arrowPointRightDirection(arrowDirection);
236 arrowPointRightDirection.rotateZ(-150.*deg);
237 fHeadPolyline.push_back(G4Point3D(x2,y2,0)+0.04*arrowPointLeftDirection);
238 fHeadPolyline.push_back(G4Point3D(x2,y2,0));
239 fHeadPolyline.push_back(G4Point3D(x2,y2,0)+0.04*arrowPointRightDirection);
241 va.SetLineWidth(fWidth);
242 va.SetColour(fColour);
243 fShaftPolyline.SetVisAttributes(va);
244 fHeadPolyline.SetVisAttributes(va);
245}
246
247void G4VisCommandSceneAddArrow2D::Arrow2D::operator()
248 (G4VGraphicsScene& sceneHandler, const G4Transform3D&, const G4ModelingParameters*)
249{
250 sceneHandler.BeginPrimitives2D();
251 sceneHandler.AddPrimitive(fShaftPolyline);
252 sceneHandler.AddPrimitive(fHeadPolyline);
253 sceneHandler.EndPrimitives2D();
254}
255
256////////////// /vis/scene/add/axes //////////////////////////////////
257
259 G4bool omitable;
260 fpCommand = new G4UIcommand ("/vis/scene/add/axes", this);
261 fpCommand -> SetGuidance ("Add axes.");
262 fpCommand -> SetGuidance
263 ("Draws axes at (x0, y0, z0) of given length and colour.");
264 fpCommand -> SetGuidance
265 ("If \"colour-string\" is \"auto\", x, y and z will be red, green and blue"
266 "\n respectively. Otherwise it can be one of the pre-defined text-specified"
267 "\n colours - see information printed by the vis manager at start-up or"
268 "\n use \"/vis/list\".");
269 fpCommand -> SetGuidance
270 ("If \"length\" is negative, it is set to about 25% of scene extent.");
271 fpCommand -> SetGuidance
272 ("If \"showtext\" is false, annotations are suppressed.");
273 G4UIparameter* parameter;
274 parameter = new G4UIparameter ("x0", 'd', omitable = true);
275 parameter->SetDefaultValue (0.);
276 fpCommand->SetParameter (parameter);
277 parameter = new G4UIparameter ("y0", 'd', omitable = true);
278 parameter->SetDefaultValue (0.);
279 fpCommand->SetParameter (parameter);
280 parameter = new G4UIparameter ("z0", 'd', omitable = true);
281 parameter->SetDefaultValue (0.);
282 fpCommand->SetParameter (parameter);
283 parameter = new G4UIparameter ("length", 'd', omitable = true);
284 parameter->SetDefaultValue (-1.);
285 fpCommand->SetParameter (parameter);
286 parameter = new G4UIparameter ("unit", 's', omitable = true);
287 parameter->SetDefaultValue ("m");
288 fpCommand->SetParameter (parameter);
289 parameter = new G4UIparameter ("colour-string", 's', omitable = true);
290 parameter->SetDefaultValue ("auto");
291 fpCommand->SetParameter (parameter);
292 parameter = new G4UIparameter ("showtext", 'b', omitable = true);
293 parameter->SetDefaultValue ("true");
294 fpCommand->SetParameter (parameter);
295}
296
298 delete fpCommand;
299}
300
302 return "";
303}
304
306
308 G4bool warn(verbosity >= G4VisManager::warnings);
309
311 if (!pScene) {
312 if (verbosity >= G4VisManager::errors) {
313 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
314 }
315 return;
316 } else {
317 if (pScene->GetExtent().GetExtentRadius() <= 0.) {
318 if (verbosity >= G4VisManager::errors) {
319 G4cerr
320 << "ERROR: Scene has no extent. Add volumes or use \"/vis/scene/add/extent\"."
321 << G4endl;
322 }
323 return;
324 }
325 }
326
327 G4String unitString, colourString, showTextString;
328 G4double x0, y0, z0, length;
329 std::istringstream is (newValue);
330 is >> x0 >> y0 >> z0 >> length >> unitString
331 >> colourString >> showTextString;
332 G4bool showText = G4UIcommand::ConvertToBool(showTextString);
333
334
335 G4double unit = G4UIcommand::ValueOf(unitString);
336 x0 *= unit; y0 *= unit; z0 *= unit;
337 const G4VisExtent& sceneExtent = pScene->GetExtent(); // Existing extent.
338 if (length < 0.) {
339 const G4double lengthMax = 0.5 * sceneExtent.GetExtentRadius();
340 const G4double intLog10Length = std::floor(std::log10(lengthMax));
341 length = std::pow(10,intLog10Length);
342 if (5.*length < lengthMax) length *= 5.;
343 else if (2.*length < lengthMax) length *= 2.;
344 } else {
345 length *= unit;
346 }
347
348 // Consult scene for arrow width...
349 G4double arrowWidth =
350 0.005 * fCurrentLineWidth * sceneExtent.GetExtentRadius();
351 // ...but limit it to length/50.
352 if (arrowWidth > length/50.) arrowWidth = length/50.;
353
354 G4VModel* model = new G4AxesModel
355 (x0, y0, z0, length, arrowWidth, colourString, newValue,
356 showText, fCurrentTextSize);
357
358 G4bool successful = pScene -> AddRunDurationModel (model, warn);
359 const G4String& currentSceneName = pScene -> GetName ();
360 if (successful) {
361 if (verbosity >= G4VisManager::confirmations) {
362 G4cout << "Axes of length " << G4BestUnit(length,"Length")
363 << "have been added to scene \"" << currentSceneName << "\"."
364 << G4endl;
365 }
366 }
367 else G4VisCommandsSceneAddUnsuccessful(verbosity);
368
370}
371
372////////////// /vis/scene/add/date ///////////////////////////////////////
373
375 G4bool omitable;
376 fpCommand = new G4UIcommand ("/vis/scene/add/date", this);
377 fpCommand -> SetGuidance ("Adds date to current scene.");
378 fpCommand -> SetGuidance
379 ("If \"date\"is omitted, the current date and time is drawn."
380 "\nOtherwise, the string, including the rest of the line, is drawn.");
381 G4UIparameter* parameter;
382 parameter = new G4UIparameter ("size", 'i', omitable = true);
383 parameter -> SetGuidance ("Screen size of text in pixels.");
384 parameter -> SetDefaultValue (18);
385 fpCommand -> SetParameter (parameter);
386 parameter = new G4UIparameter ("x-position", 'd', omitable = true);
387 parameter -> SetGuidance ("x screen position in range -1 < x < 1.");
388 parameter -> SetDefaultValue (0.95);
389 fpCommand -> SetParameter (parameter);
390 parameter = new G4UIparameter ("y-position", 'd', omitable = true);
391 parameter -> SetGuidance ("y screen position in range -1 < y < 1.");
392 parameter -> SetDefaultValue (0.9);
393 fpCommand -> SetParameter (parameter);
394 parameter = new G4UIparameter ("layout", 's', omitable = true);
395 parameter -> SetGuidance ("Layout, i.e., adjustment: left|centre|right.");
396 parameter -> SetDefaultValue ("right");
397 fpCommand -> SetParameter (parameter);
398 parameter = new G4UIparameter ("date", 's', omitable = true);
399 parameter -> SetDefaultValue ("-");
400 fpCommand -> SetParameter (parameter);
401}
402
404 delete fpCommand;
405}
406
408 return "";
409}
410
412{
414 G4bool warn(verbosity >= G4VisManager::warnings);
415
417 if (!pScene) {
418 if (verbosity >= G4VisManager::errors) {
419 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
420 }
421 return;
422 }
423
424 G4int size;
425 G4double x, y;
426 G4String layoutString, dateString;
427 std::istringstream is(newValue);
428 is >> size >> x >> y >> layoutString >> dateString;
429 // Read rest of line, if any.
430 const size_t NREMAINDER = 100;
431 char remainder[NREMAINDER];
432 remainder[0]='\0'; // In case there is nothing remaining.
433 is.getline(remainder, NREMAINDER);
434 dateString += remainder;
436 if (layoutString(0) == 'l') layout = G4Text::left;
437 else if (layoutString(0) == 'c') layout = G4Text::centre;
438 else if (layoutString(0) == 'r') layout = G4Text::right;
439
440 Date* date = new Date(fpVisManager, size, x, y, layout, dateString);
441 G4VModel* model =
443 model->SetType("Date");
444 model->SetGlobalTag("Date");
445 model->SetGlobalDescription("Date: " + newValue);
446 const G4String& currentSceneName = pScene -> GetName ();
447 G4bool successful = pScene -> AddRunDurationModel (model, warn);
448 if (successful) {
449 if (verbosity >= G4VisManager::confirmations) {
450 G4cout << "Date has been added to scene \""
451 << currentSceneName << "\"."
452 << G4endl;
453 }
454 }
455 else G4VisCommandsSceneAddUnsuccessful(verbosity);
456
458}
459
460void G4VisCommandSceneAddDate::Date::operator()
461 (G4VGraphicsScene& sceneHandler, const G4Transform3D&, const G4ModelingParameters*)
462{
463 G4String time;
464 if (fDate == "-") {
465 time = fTimer.GetClockTime();
466 } else {
467 time = fDate;
468 }
469 // Check for \n, starting from back, and erase.
470 std::string::size_type i = time.rfind('\n');
471 if (i != std::string::npos) time.erase(i);
472 G4Text text(time, G4Point3D(fX, fY, 0.));
473 text.SetScreenSize(fSize);
474 text.SetLayout(fLayout);
475 G4VisAttributes textAtts(G4Colour(0.,1.,1));
476 text.SetVisAttributes(textAtts);
477 sceneHandler.BeginPrimitives2D();
478 sceneHandler.AddPrimitive(text);
479 sceneHandler.EndPrimitives2D();
480}
481
482////////////// /vis/scene/add/digis ///////////////////////////////////////
483
485 fpCommand = new G4UIcmdWithoutParameter ("/vis/scene/add/digis", this);
486 fpCommand -> SetGuidance ("Adds digis to current scene.");
487 fpCommand -> SetGuidance
488 ("Digis are drawn at end of event when the scene in which"
489 "\nthey are added is current.");
490}
491
493 delete fpCommand;
494}
495
497 return "";
498}
499
501
503 G4bool warn(verbosity >= G4VisManager::warnings);
504
506 if (!pScene) {
507 if (verbosity >= G4VisManager::errors) {
508 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
509 }
510 return;
511 }
512
513 G4VModel* model = new G4DigiModel;
514 const G4String& currentSceneName = pScene -> GetName ();
515 G4bool successful = pScene -> AddEndOfEventModel (model, warn);
516 if (successful) {
517 if (verbosity >= G4VisManager::confirmations) {
518 G4cout << "Digis, if any, will be drawn at end of run in scene \""
519 << currentSceneName << "\"."
520 << G4endl;
521 }
522 }
523 else G4VisCommandsSceneAddUnsuccessful(verbosity);
524
526}
527
528////////////// /vis/scene/add/electricField ///////////////////////////////////////
529
531 G4bool omitable;
532 fpCommand = new G4UIcommand ("/vis/scene/add/electricField", this);
533 fpCommand -> SetGuidance
534 ("Adds electric field representation to current scene.");
535 fpCommand -> SetGuidance
536 ("The first parameter is no. of data points per half scene. So, possibly, at"
537 "\nmaximum, the number of data points sampled is (2*n+1)^3, which can grow"
538 "\nlarge--be warned!"
539 "\nThe default value is 10, i.e., a 21x21x21 array, i.e., 9,261 sampling points."
540 "\nThat may swamp you scene, but usually, a field is limited to a small part of"
541 "\nthe scene, so it's not a problem. But if it is, here are some of the things"
542 "\nyou can do:"
543 "\n- reduce the number of data points per half scene (first parameter);"
544 "\n- specify \"lightArrow\" (second parameter);"
545 "\n- restrict the region sampled with \"/vis/set/extentForField\";"
546 "\n- restrict the drawing to a specific volume with"
547 "\n \"/vis/set/volumeForField\" or \"/vis/touchable/volumeForField\"."
548 "\nNote: you may have to deactivate existing field models with"
549 "\n \"/vis/scene/activateModel Field false\" and re-issue"
550 "\n \"/vis/scene/add/...Field\" command again.");
551 fpCommand -> SetGuidance
552 ("In the arrow representation, the length of the arrow is proportional"
553 "\nto the magnitude of the field and the colour is mapped onto the range"
554 "\nas a fraction of the maximum magnitude: 0->0.5->1 is red->green->blue.");
555 G4UIparameter* parameter;
556 parameter = new G4UIparameter ("nDataPointsPerHalfScene", 'i', omitable = true);
557 parameter -> SetDefaultValue (10);
558 fpCommand -> SetParameter (parameter);
559 parameter = new G4UIparameter ("representation", 's', omitable = true);
560 parameter -> SetParameterCandidates("fullArrow lightArrow");
561 parameter -> SetDefaultValue ("fullArrow");
562 fpCommand -> SetParameter (parameter);
563}
564
566 delete fpCommand;
567}
568
570 return "";
571}
572
574(G4UIcommand*, G4String newValue) {
575
577 G4bool warn(verbosity >= G4VisManager::warnings);
578
580 if (!pScene) {
581 if (verbosity >= G4VisManager::errors) {
582 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
583 }
584 return;
585 }
586
587 G4int nDataPointsPerHalfScene;
588 G4String representation;
589 std::istringstream iss(newValue);
590 iss >> nDataPointsPerHalfScene >> representation;
592 modelRepresentation = G4ElectricFieldModel::fullArrow;
593 if (representation == "lightArrow") {
594 modelRepresentation = G4ElectricFieldModel::lightArrow;
595 }
596 G4VModel* model;
597 model = new G4ElectricFieldModel
598 (nDataPointsPerHalfScene,modelRepresentation,
602 const G4String& currentSceneName = pScene -> GetName ();
603 G4bool successful = pScene -> AddRunDurationModel (model, warn);
604 if (successful) {
605 if (verbosity >= G4VisManager::confirmations) {
606 G4cout
607 << "Electric field, if any, will be drawn in scene \""
608 << currentSceneName
609 << "\"\n with "
610 << nDataPointsPerHalfScene
611 << " data points per half scene and with representation \""
612 << representation
613 << '\"'
614 << G4endl;
615 }
616 }
617 else G4VisCommandsSceneAddUnsuccessful(verbosity);
618
620}
621
622////////////// /vis/scene/add/eventID ///////////////////////////////////////
623
625 G4bool omitable;
626 fpCommand = new G4UIcommand ("/vis/scene/add/eventID", this);
627 fpCommand -> SetGuidance ("Adds eventID to current scene.");
628 fpCommand -> SetGuidance
629 ("Run and event numbers are drawn at end of event or run when"
630 "\n the scene in which they are added is current.");
631 G4UIparameter* parameter;
632 parameter = new G4UIparameter ("size", 'i', omitable = true);
633 parameter -> SetGuidance ("Screen size of text in pixels.");
634 parameter -> SetDefaultValue (18);
635 fpCommand -> SetParameter (parameter);
636 parameter = new G4UIparameter ("x-position", 'd', omitable = true);
637 parameter -> SetGuidance ("x screen position in range -1 < x < 1.");
638 parameter -> SetDefaultValue (-0.95);
639 fpCommand -> SetParameter (parameter);
640 parameter = new G4UIparameter ("y-position", 'd', omitable = true);
641 parameter -> SetGuidance ("y screen position in range -1 < y < 1.");
642 parameter -> SetDefaultValue (0.9);
643 fpCommand -> SetParameter (parameter);
644 parameter = new G4UIparameter ("layout", 's', omitable = true);
645 parameter -> SetGuidance ("Layout, i.e., adjustment: left|centre|right.");
646 parameter -> SetDefaultValue ("left");
647 fpCommand -> SetParameter (parameter);
648}
649
651 delete fpCommand;
652}
653
655 return "";
656}
657
659{
661 G4bool warn(verbosity >= G4VisManager::warnings);
662
664 if (!pScene) {
665 if (verbosity >= G4VisManager::errors) {
666 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
667 }
668 return;
669 }
670
671 G4int size;
672 G4double x, y;
673 G4String layoutString;
674 std::istringstream is(newValue);
675 is >> size >> x >> y >> layoutString;
676
678 if (layoutString(0) == 'l') layout = G4Text::left;
679 else if (layoutString(0) == 'c') layout = G4Text::centre;
680 else if (layoutString(0) == 'r') layout = G4Text::right;
681
682 // For End of Event (only for reviewing kept events one by one)
683 EventID* eoeEventID
684 = new EventID(forEndOfEvent, fpVisManager, size, x, y, layout);
685 G4VModel* eoeModel =
687 eoeModel->SetType("EoEEventID");
688 eoeModel->SetGlobalTag("EoEEventID");
689 eoeModel->SetGlobalDescription("EoEEventID: " + newValue);
690 G4bool successfulEoE = pScene -> AddEndOfEventModel (eoeModel, warn);
691
692 // For End of Run
693 EventID* eorEventID
694 = new EventID(forEndOfRun, fpVisManager, size, x, y, layout);
695 G4VModel* eorModel =
697 eorModel->SetType("EoREventID");
698 eorModel->SetGlobalTag("EoREventID");
699 eorModel->SetGlobalDescription("EoREventID: " + newValue);
700 G4bool successfulEoR = pScene -> AddEndOfRunModel (eorModel, warn);
701
702 if (successfulEoE && successfulEoR) {
703 if (verbosity >= G4VisManager::confirmations) {
704 const G4String& currentSceneName = pScene -> GetName ();
705 G4cout << "EventID has been added to scene \""
706 << currentSceneName << "\"."
707 << G4endl;
708 }
709 }
710 else G4VisCommandsSceneAddUnsuccessful(verbosity);
711
713}
714
715void G4VisCommandSceneAddEventID::EventID::operator()
716(G4VGraphicsScene& sceneHandler, const G4Transform3D&, const G4ModelingParameters* mp)
717{
719 if(!runManager)
720 return;
721
722 const G4Run* currentRun = runManager->GetCurrentRun();
723 if (!currentRun) return;
724
725 const G4int currentRunID = currentRun->GetRunID();
726
727 std::ostringstream oss;
728 switch (fForWhat) {
729 case forEndOfEvent:
730 {
731 // Only use if reviewing kept events
732 if (!fpVisManager->GetReviewingKeptEvents()) return;
733 const G4Event* currentEvent = mp->GetEvent();
734 if (!currentEvent) return;
735 G4int eventID = currentEvent->GetEventID();
736 oss << "Run " << currentRunID << " Event " << eventID;
737 break;
738 }
739 case forEndOfRun:
740 {
741 // Only use if NOT reviewing kept events
742 if (fpVisManager->GetReviewingKeptEvents()) return;
743 const G4int nEvents = currentRun->GetNumberOfEventToBeProcessed();
744 const auto* events = currentRun->GetEventVector();
745 size_t nKeptEvents = events? events->size(): 0;
746 oss << "Run " << currentRunID << " (" << nEvents << " event";
747 if (nEvents != 1) oss << 's';
748 oss << ", " << nKeptEvents << " kept)";
749 break;
750 }
751 default:
752 return;
753 }
754
755 G4Text text(oss.str(), G4Point3D(fX, fY, 0.));
756 text.SetScreenSize(fSize);
757 text.SetLayout(fLayout);
758 G4VisAttributes textAtts(G4Colour(0.,1.,1));
759 text.SetVisAttributes(textAtts);
760 sceneHandler.BeginPrimitives2D();
761 sceneHandler.AddPrimitive(text);
762 sceneHandler.EndPrimitives2D();
763}
764
765////////////// /vis/scene/add/extent ///////////////////////////////////////
766
768 fpCommand = new G4UIcommand("/vis/scene/add/extent", this);
769 fpCommand -> SetGuidance
770 ("Adds a dummy model with given extent to the current scene."
771 "\nRequires the limits: xmin, xmax, ymin, ymax, zmin, zmax unit."
772 "\nThis can be used to provide an extent to the scene even if"
773 "\nno other models with extent are available. For example,"
774 "\neven if there is no geometry. In that case, for example:"
775 "\n /vis/open OGL"
776 "\n /vis/scene/create"
777 "\n /vis/scene/add/extent -300 300 -300 300 -300 300 cm"
778 "\n /vis/sceneHandler/attach");
779 G4bool omitable;
780 G4UIparameter* parameter;
781 parameter = new G4UIparameter ("xmin", 'd', omitable = true);
782 parameter -> SetDefaultValue (0.);
783 fpCommand -> SetParameter (parameter);
784 parameter = new G4UIparameter ("xmax", 'd', omitable = true);
785 parameter -> SetDefaultValue (0.);
786 fpCommand -> SetParameter (parameter);
787 parameter = new G4UIparameter ("ymin", 'd', omitable = true);
788 parameter -> SetDefaultValue (0.);
789 fpCommand -> SetParameter (parameter);
790 parameter = new G4UIparameter ("ymax", 'd', omitable = true);
791 parameter -> SetDefaultValue (0.);
792 fpCommand -> SetParameter (parameter);
793 parameter = new G4UIparameter ("zmin", 'd', omitable = true);
794 parameter -> SetDefaultValue (0.);
795 fpCommand -> SetParameter (parameter);
796 parameter = new G4UIparameter ("zmax", 'd', omitable = true);
797 parameter -> SetDefaultValue (0.);
798 fpCommand -> SetParameter (parameter);
799 parameter = new G4UIparameter ("unit", 's', omitable = true);
800 parameter -> SetDefaultValue ("m");
801 fpCommand -> SetParameter (parameter);
802}
803
805 delete fpCommand;
806}
807
809 return "";
810}
811
813{
815 G4bool warn(verbosity >= G4VisManager::warnings);
816
818 if (!pScene) {
819 if (verbosity >= G4VisManager::errors) {
820 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
821 }
822 return;
823 }
824
825 G4double xmin, xmax, ymin, ymax, zmin, zmax;
826 G4String unitString;
827 std::istringstream is(newValue);
828 is >> xmin >> xmax >> ymin >> ymax >> zmin >> zmax >> unitString;
829 G4double unit = G4UIcommand::ValueOf(unitString);
830 xmin *= unit; xmax *= unit;
831 ymin *= unit; ymax *= unit;
832 zmin *= unit; zmax *= unit;
833
834 G4VisExtent visExtent(xmin, xmax, ymin, ymax, zmin, zmax);
835 Extent* extent = new Extent(xmin, xmax, ymin, ymax, zmin, zmax);
836 G4VModel* model =
838 model->SetType("Extent");
839 model->SetGlobalTag("Extent");
840 model->SetGlobalDescription("Extent: " + newValue);
841 model->SetExtent(visExtent);
842 const G4String& currentSceneName = pScene -> GetName ();
843 G4bool successful = pScene -> AddRunDurationModel (model, warn);
844 if (successful) {
845 if (verbosity >= G4VisManager::confirmations) {
846 G4cout << "A benign model with extent "
847 << visExtent
848 << " has been added to scene \""
849 << currentSceneName << "\"."
850 << G4endl;
851 }
852 }
853 else G4VisCommandsSceneAddUnsuccessful(verbosity);
854
856}
857
858G4VisCommandSceneAddExtent::Extent::Extent
859(G4double xmin, G4double xmax,
860 G4double ymin, G4double ymax,
861 G4double zmin, G4double zmax):
862fExtent(xmin,xmax,ymin,ymax,zmin,zmax)
863{}
864
865void G4VisCommandSceneAddExtent::Extent::operator()
867{}
868
869////////////// /vis/scene/add/frame ///////////////////////////////////////
870
872 fpCommand = new G4UIcommand("/vis/scene/add/frame", this);
873 fpCommand -> SetGuidance ("Add frame to current scene.");
874 G4bool omitable;
875 G4UIparameter* parameter;
876 parameter = new G4UIparameter ("size", 'd', omitable = true);
877 parameter -> SetGuidance ("Size of frame. 1 = full window.");
878 parameter -> SetParameterRange ("size > 0 && size <=1");
879 parameter -> SetDefaultValue (0.97);
880 fpCommand -> SetParameter (parameter);
881}
882
884 delete fpCommand;
885}
886
888 return "";
889}
890
892{
894 G4bool warn(verbosity >= G4VisManager::warnings);
895
897 if (!pScene) {
898 if (verbosity >= G4VisManager::errors) {
899 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
900 }
901 return;
902 }
903
904 G4double size;
905 std::istringstream is(newValue);
906 is >> size;
907
908 Frame* frame = new Frame(size, fCurrentLineWidth, fCurrentColour);
909 G4VModel* model =
911 model->SetType("Frame");
912 model->SetGlobalTag("Frame");
913 model->SetGlobalDescription("Frame: " + newValue);
914 const G4String& currentSceneName = pScene -> GetName ();
915 G4bool successful = pScene -> AddRunDurationModel (model, warn);
916 if (successful) {
917 if (verbosity >= G4VisManager::confirmations) {
918 G4cout << "Frame has been added to scene \""
919 << currentSceneName << "\"."
920 << G4endl;
921 }
922 }
923 else G4VisCommandsSceneAddUnsuccessful(verbosity);
924
926}
927
928void G4VisCommandSceneAddFrame::Frame::operator()
929 (G4VGraphicsScene& sceneHandler, const G4Transform3D&, const G4ModelingParameters*)
930{
931 G4Polyline frame;
932 frame.push_back(G4Point3D( fSize, fSize, 0.));
933 frame.push_back(G4Point3D(-fSize, fSize, 0.));
934 frame.push_back(G4Point3D(-fSize, -fSize, 0.));
935 frame.push_back(G4Point3D( fSize, -fSize, 0.));
936 frame.push_back(G4Point3D( fSize, fSize, 0.));
938 va.SetLineWidth(fWidth);
939 va.SetColour(fColour);
940 frame.SetVisAttributes(va);
941 sceneHandler.BeginPrimitives2D();
942 sceneHandler.AddPrimitive(frame);
943 sceneHandler.EndPrimitives2D();
944}
945
946////////////// /vis/scene/add/gps ///////////////////////////////////////
947
949 G4bool omitable;
950 G4UIparameter* parameter;
951 fpCommand = new G4UIcommand ("/vis/scene/add/gps", this);
952 fpCommand -> SetGuidance
953 ("A representation of the source(s) of the General Particle Source"
954 "\nwill be added to current scene and drawn, if applicable.");
956 fpCommand->SetGuidance("Default: red and transparent.");
957 parameter = new G4UIparameter("red_or_string", 's', omitable = true);
958 parameter -> SetDefaultValue ("1.");
959 fpCommand -> SetParameter (parameter);
960 parameter = new G4UIparameter("green", 'd', omitable = true);
961 parameter -> SetDefaultValue (0.);
962 fpCommand -> SetParameter (parameter);
963 parameter = new G4UIparameter ("blue", 'd', omitable = true);
964 parameter -> SetDefaultValue (0.);
965 fpCommand -> SetParameter (parameter);
966 parameter = new G4UIparameter ("opacity", 'd', omitable = true);
967 parameter -> SetDefaultValue (0.3);
968 fpCommand -> SetParameter (parameter);
969}
970
972 delete fpCommand;
973}
974
976 return "";
977}
978
980
982 G4bool warn(verbosity >= G4VisManager::warnings);
983
985 if (!pScene) {
986 if (verbosity >= G4VisManager::errors) {
987 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
988 }
989 return;
990 }
991
992 G4String redOrString;
993 G4double green, blue, opacity;
994 std::istringstream iss(newValue);
995 iss >> redOrString >> green >> blue >> opacity;
996 G4Colour colour(1.,0.,0.,0.3); // Default red and transparent.
997 ConvertToColour(colour, redOrString, green, blue, opacity);
998
999 G4VModel* model = new G4GPSModel(colour);
1000 const G4String& currentSceneName = pScene -> GetName ();
1001 G4bool successful = pScene -> AddRunDurationModel (model, warn);
1002 if (successful) {
1003 if (verbosity >= G4VisManager::confirmations) {
1004 G4cout <<
1005 "A representation of the source(s) of the General Particle Source will be drawn"
1006 "\n in colour " << colour << " for scene \""
1007 << currentSceneName << "\" if applicable."
1008 << G4endl;
1009 }
1010 }
1011 else G4VisCommandsSceneAddUnsuccessful(verbosity);
1012
1014}
1015
1016////////////// /vis/scene/add/hits ///////////////////////////////////////
1017
1019 fpCommand = new G4UIcmdWithoutParameter ("/vis/scene/add/hits", this);
1020 fpCommand -> SetGuidance ("Adds hits to current scene.");
1021 fpCommand -> SetGuidance
1022 ("Hits are drawn at end of event when the scene in which"
1023 "\nthey are added is current.");
1024}
1025
1027 delete fpCommand;
1028}
1029
1031 return "";
1032}
1033
1035
1037 G4bool warn(verbosity >= G4VisManager::warnings);
1038
1040 if (!pScene) {
1041 if (verbosity >= G4VisManager::errors) {
1042 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
1043 }
1044 return;
1045 }
1046
1047 G4VModel* model = new G4HitsModel;
1048 const G4String& currentSceneName = pScene -> GetName ();
1049 G4bool successful = pScene -> AddEndOfEventModel (model, warn);
1050 if (successful) {
1051 if (verbosity >= G4VisManager::confirmations) {
1052 G4cout << "Hits, if any, will be drawn at end of run in scene \""
1053 << currentSceneName << "\"."
1054 << G4endl;
1055 }
1056 }
1057 else G4VisCommandsSceneAddUnsuccessful(verbosity);
1058
1060}
1061
1062////////////// /vis/scene/add/line ///////////////////////////////////////
1063
1065 fpCommand = new G4UIcommand("/vis/scene/add/line", this);
1066 fpCommand -> SetGuidance ("Adds line to current scene.");
1067 G4bool omitable;
1068 G4UIparameter* parameter;
1069 parameter = new G4UIparameter ("x1", 'd', omitable = false);
1070 fpCommand -> SetParameter (parameter);
1071 parameter = new G4UIparameter ("y1", 'd', omitable = false);
1072 fpCommand -> SetParameter (parameter);
1073 parameter = new G4UIparameter ("z1", 'd', omitable = false);
1074 fpCommand -> SetParameter (parameter);
1075 parameter = new G4UIparameter ("x2", 'd', omitable = false);
1076 fpCommand -> SetParameter (parameter);
1077 parameter = new G4UIparameter ("y2", 'd', omitable = false);
1078 fpCommand -> SetParameter (parameter);
1079 parameter = new G4UIparameter ("z2", 'd', omitable = false);
1080 fpCommand -> SetParameter (parameter);
1081 parameter = new G4UIparameter ("unit", 's', omitable = true);
1082 parameter->SetDefaultValue ("m");
1083 fpCommand->SetParameter (parameter);
1084}
1085
1087 delete fpCommand;
1088}
1089
1091 return "";
1092}
1093
1095{
1097 G4bool warn(verbosity >= G4VisManager::warnings);
1098
1100 if (!pScene) {
1101 if (verbosity >= G4VisManager::errors) {
1102 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
1103 }
1104 return;
1105 }
1106
1107 G4String unitString;
1108 G4double x1, y1, z1, x2, y2, z2;
1109 std::istringstream is(newValue);
1110 is >> x1 >> y1 >> z1 >> x2 >> y2 >> z2 >> unitString;
1111 G4double unit = G4UIcommand::ValueOf(unitString);
1112 x1 *= unit; y1 *= unit; z1 *= unit;
1113 x2 *= unit; y2 *= unit; z2 *= unit;
1114
1115 Line* line = new Line(x1, y1, z1, x2, y2, z2,
1117 G4VModel* model =
1119 model->SetType("Line");
1120 model->SetGlobalTag("Line");
1121 model->SetGlobalDescription("Line: " + newValue);
1122 const G4String& currentSceneName = pScene -> GetName ();
1123 G4bool successful = pScene -> AddRunDurationModel (model, warn);
1124 if (successful) {
1125 if (verbosity >= G4VisManager::confirmations) {
1126 G4cout << "Line has been added to scene \""
1127 << currentSceneName << "\"."
1128 << G4endl;
1129 }
1130 }
1131 else G4VisCommandsSceneAddUnsuccessful(verbosity);
1132
1134}
1135
1136G4VisCommandSceneAddLine::Line::Line
1137(G4double x1, G4double y1, G4double z1,
1138 G4double x2, G4double y2, G4double z2,
1139 G4double width, const G4Colour& colour):
1140 fWidth(width), fColour(colour)
1141{
1142 fPolyline.push_back(G4Point3D(x1,y1,z1));
1143 fPolyline.push_back(G4Point3D(x2,y2,z2));
1144 G4VisAttributes va;
1145 va.SetLineWidth(fWidth);
1146 va.SetColour(fColour);
1147 fPolyline.SetVisAttributes(va);
1148}
1149
1150void G4VisCommandSceneAddLine::Line::operator()
1151 (G4VGraphicsScene& sceneHandler, const G4Transform3D&, const G4ModelingParameters*)
1152{
1153 sceneHandler.BeginPrimitives();
1154 sceneHandler.AddPrimitive(fPolyline);
1155 sceneHandler.EndPrimitives();
1156}
1157
1158////////////// /vis/scene/add/line2D ///////////////////////////////////////
1159
1161 fpCommand = new G4UIcommand("/vis/scene/add/line2D", this);
1162 fpCommand -> SetGuidance ("Adds 2D line to current scene.");
1163 G4bool omitable;
1164 G4UIparameter* parameter;
1165 parameter = new G4UIparameter ("x1", 'd', omitable = false);
1166 fpCommand -> SetParameter (parameter);
1167 parameter = new G4UIparameter ("y1", 'd', omitable = false);
1168 fpCommand -> SetParameter (parameter);
1169 parameter = new G4UIparameter ("x2", 'd', omitable = false);
1170 fpCommand -> SetParameter (parameter);
1171 parameter = new G4UIparameter ("y2", 'd', omitable = false);
1172 fpCommand -> SetParameter (parameter);
1173}
1174
1176 delete fpCommand;
1177}
1178
1180 return "";
1181}
1182
1184{
1186 G4bool warn(verbosity >= G4VisManager::warnings);
1187
1189 if (!pScene) {
1190 if (verbosity >= G4VisManager::errors) {
1191 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
1192 }
1193 return;
1194 }
1195
1196 G4double x1, y1, x2, y2;
1197 std::istringstream is(newValue);
1198 is >> x1 >> y1 >> x2 >> y2;
1199
1200 Line2D* line2D = new Line2D
1201 (x1, y1, x2, y2, fCurrentLineWidth, fCurrentColour);
1202 G4VModel* model =
1204 model->SetType("Line2D");
1205 model->SetGlobalTag("Line2D");
1206 model->SetGlobalDescription("Line2D: " + newValue);
1207 const G4String& currentSceneName = pScene -> GetName ();
1208 G4bool successful = pScene -> AddRunDurationModel (model, warn);
1209 if (successful) {
1210 if (verbosity >= G4VisManager::confirmations) {
1211 G4cout << "A 2D line has been added to scene \""
1212 << currentSceneName << "\"."
1213 << G4endl;
1214 }
1215 }
1216 else G4VisCommandsSceneAddUnsuccessful(verbosity);
1217
1219}
1220
1221G4VisCommandSceneAddLine2D::Line2D::Line2D
1222(G4double x1, G4double y1,
1223 G4double x2, G4double y2,
1224 G4double width, const G4Colour& colour):
1225 fWidth(width), fColour(colour)
1226{
1227 fPolyline.push_back(G4Point3D(x1,y1,0));
1228 fPolyline.push_back(G4Point3D(x2,y2,0));
1229 G4VisAttributes va;
1230 va.SetLineWidth(fWidth);
1231 va.SetColour(fColour);
1232 fPolyline.SetVisAttributes(va);
1233}
1234
1235void G4VisCommandSceneAddLine2D::Line2D::operator()
1236 (G4VGraphicsScene& sceneHandler, const G4Transform3D&, const G4ModelingParameters*)
1237{
1238 sceneHandler.BeginPrimitives2D();
1239 sceneHandler.AddPrimitive(fPolyline);
1240 sceneHandler.EndPrimitives2D();
1241}
1242
1243////////////// /vis/scene/add/logicalVolume //////////////////////////////////
1244
1246 G4bool omitable;
1247 fpCommand = new G4UIcommand ("/vis/scene/add/logicalVolume", this);
1248 fpCommand -> SetGuidance ("Adds a logical volume to the current scene,");
1249 fpCommand -> SetGuidance
1250 ("Shows boolean components (if any), voxels (if any), readout geometry"
1251 "\n (if any), local axes and overlaps (if any), under control of the"
1252 "\n appropriate flag."
1253 "\n Note: voxels are not constructed until start of run -"
1254 "\n \"/run/beamOn\". (For voxels without a run, \"/run/beamOn 0\".)");
1255 G4UIparameter* parameter;
1256 parameter = new G4UIparameter ("logical-volume-name", 's', omitable = false);
1257 fpCommand -> SetParameter (parameter);
1258 parameter = new G4UIparameter ("depth-of-descent", 'i', omitable = true);
1259 parameter -> SetGuidance ("Depth of descent of geometry hierarchy.");
1260 parameter -> SetDefaultValue (1);
1261 fpCommand -> SetParameter (parameter);
1262 parameter = new G4UIparameter ("booleans-flag", 'b', omitable = true);
1263 parameter -> SetDefaultValue (true);
1264 fpCommand -> SetParameter (parameter);
1265 parameter = new G4UIparameter ("voxels-flag", 'b', omitable = true);
1266 parameter -> SetDefaultValue (true);
1267 fpCommand -> SetParameter (parameter);
1268 parameter = new G4UIparameter ("readout-flag", 'b', omitable = true);
1269 parameter -> SetDefaultValue (true);
1270 fpCommand -> SetParameter (parameter);
1271 parameter = new G4UIparameter ("axes-flag", 'b', omitable = true);
1272 parameter -> SetDefaultValue (true);
1273 parameter -> SetGuidance ("Set \"false\" to suppress axes.");
1274 fpCommand -> SetParameter (parameter);
1275 parameter = new G4UIparameter("check-overlap-flag", 'b', omitable = true);
1276 parameter->SetDefaultValue(true);
1277 parameter -> SetGuidance ("Set \"false\" to suppress overlap check.");
1278 fpCommand->SetParameter(parameter);
1279}
1280
1282 delete fpCommand;
1283}
1284
1286 return "";
1287}
1288
1290 G4String newValue) {
1291
1293 G4bool warn(verbosity >= G4VisManager::warnings);
1294
1296 if (!pScene) {
1297 if (verbosity >= G4VisManager::errors) {
1298 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
1299 }
1300 return;
1301 }
1302
1303 G4String name;
1304 G4int requestedDepthOfDescent;
1305 G4String booleansString, voxelsString, readoutString, axesString;
1306 G4String overlapString;
1307 std::istringstream is (newValue);
1308 is >> name >> requestedDepthOfDescent
1309 >> booleansString >> voxelsString >> readoutString >> axesString
1310 >> overlapString;
1311 G4bool booleans = G4UIcommand::ConvertToBool(booleansString);
1312 G4bool voxels = G4UIcommand::ConvertToBool(voxelsString);
1313 G4bool readout = G4UIcommand::ConvertToBool(readoutString);
1314 G4bool axes = G4UIcommand::ConvertToBool(axesString);
1315 G4bool checkOverlaps = G4UIcommand::ConvertToBool(overlapString);
1316
1318 int nLV = pLVStore -> size ();
1319 int iLV;
1320 G4LogicalVolume* pLV = 0;
1321 for (iLV = 0; iLV < nLV; iLV++ ) {
1322 pLV = (*pLVStore) [iLV];
1323 if (pLV -> GetName () == name) break;
1324 }
1325 if (iLV == nLV) {
1326 if (verbosity >= G4VisManager::errors) {
1327 G4cerr << "ERROR: Logical volume " << name
1328 << " not found in logical volume store." << G4endl;
1329 }
1330 return;
1331 }
1332
1333 const std::vector<G4Scene::Model>& rdModelList =
1334 pScene -> GetRunDurationModelList();
1335 std::vector<G4Scene::Model>::const_iterator i;
1336 for (i = rdModelList.begin(); i != rdModelList.end(); ++i) {
1337 if (i->fpModel->GetGlobalDescription().find("Volume") != std::string::npos) break;
1338 }
1339 if (i != rdModelList.end()) {
1340 if (verbosity >= G4VisManager::errors) {
1341 G4cout << "There is already a volume, \""
1342 << i->fpModel->GetGlobalDescription()
1343 << "\",\n in the run-duration model list of scene \""
1344 << pScene -> GetName()
1345 << "\".\n Your logical volume must be the only volume in the scene."
1346 << "\n Create a new scene and try again:"
1347 << "\n /vis/specify " << name
1348 << "\n or"
1349 << "\n /vis/scene/create"
1350 << "\n /vis/scene/add/logicalVolume " << name
1351 << "\n /vis/sceneHandler/attach"
1352 << "\n (and also, if necessary, /vis/viewer/flush)"
1353 << G4endl;
1354 }
1355 return;
1356 }
1357
1359 (pLV, requestedDepthOfDescent, booleans, voxels, readout, checkOverlaps);
1360 const G4String& currentSceneName = pScene -> GetName ();
1361 G4bool successful = pScene -> AddRunDurationModel (model, warn);
1362
1363 if (successful) {
1364
1365 G4bool axesSuccessful = false;
1366 if (axes) {
1367 const G4double radius = model->GetExtent().GetExtentRadius();
1368 const G4double axisLengthMax = radius / 2.;
1369 const G4double intLog10Length = std::floor(std::log10(axisLengthMax));
1370 G4double axisLength = std::pow(10,intLog10Length);
1371 if (5.*axisLength < axisLengthMax) axisLength *= 5.;
1372 else if (2.*axisLength < axisLengthMax) axisLength *= 2.;
1373 const G4double axisWidth = axisLength / 20.;
1374 G4VModel* axesModel = new G4AxesModel(0.,0.,0.,axisLength,axisWidth);
1375 axesSuccessful = pScene -> AddRunDurationModel (axesModel, warn);
1376 }
1377
1378// if (verbosity >= G4VisManager::warnings) {
1379// const std::map<G4String,G4AttDef>* attDefs = model->GetAttDefs();
1380// std::vector<G4AttValue>* attValues = model->CreateCurrentAttValues();
1381// G4cout << G4AttCheck(attValues, attDefs);
1382// delete attValues;
1383// }
1384
1385 if (verbosity >= G4VisManager::confirmations) {
1386 G4cout << "Logical volume \"" << pLV -> GetName ()
1387 << "\" with requested depth of descent "
1388 << requestedDepthOfDescent
1389 << ",\n with";
1390 if (!booleans) G4cout << "out";
1391 G4cout << " boolean components, with";
1392 if (!voxels) G4cout << "out";
1393 G4cout << " voxels,\n with";
1394 if (!readout) G4cout << "out";
1395 G4cout << " readout geometry and with";
1396 if (!checkOverlaps) G4cout << "out";
1397 G4cout << " overlap checking"
1398 << "\n has been added to scene \"" << currentSceneName << "\".";
1399 if (axes) {
1400 if (axesSuccessful) {
1401 G4cout <<
1402 "\n Axes have also been added at the origin of local cooordinates.";
1403 } else {
1404 G4cout <<
1405 "\n Axes have not been added for some reason possibly stated above.";
1406 }
1407 }
1408 G4cout << G4endl;
1409 }
1410 }
1411 else {
1413 return;
1414 }
1415
1417}
1418
1419
1420////////////// /vis/scene/add/logo //////////////////////////////////
1421
1423 G4bool omitable;
1424 fpCommand = new G4UIcommand ("/vis/scene/add/logo", this);
1425 fpCommand -> SetGuidance ("Adds a G4 logo to the current scene.");
1426 fpCommand -> SetGuidance
1427 ("If \"unit\" is \"auto\", height is roughly one tenth of scene extent.");
1428 fpCommand -> SetGuidance
1429 ("\"direction\" is that of outward-facing normal to front face of logo."
1430 "\nIf \"direction\" is \"auto\", logo faces the user in the current viewer.");
1431 fpCommand -> SetGuidance
1432 ("\nIf \"placement\" is \"auto\", logo is placed at bottom right of screen"
1433 "\n when viewed from logo direction.");
1434 G4UIparameter* parameter;
1435 parameter = new G4UIparameter ("height", 'd', omitable = true);
1436 parameter->SetDefaultValue (1.);
1437 fpCommand->SetParameter (parameter);
1438 parameter = new G4UIparameter ("unit", 's', omitable = true);
1439 parameter->SetDefaultValue ("auto");
1440 fpCommand->SetParameter (parameter);
1441 parameter = new G4UIparameter ("direction", 's', omitable = true);
1442 parameter->SetGuidance ("auto|[-]x|[-]y|[-]z");
1443 parameter->SetDefaultValue ("auto");
1444 fpCommand->SetParameter (parameter);
1445 parameter = new G4UIparameter ("red", 'd', omitable = true);
1446 parameter->SetDefaultValue (0.);
1447 fpCommand->SetParameter (parameter);
1448 parameter = new G4UIparameter ("green", 'd', omitable = true);
1449 parameter->SetDefaultValue (1.);
1450 fpCommand->SetParameter (parameter);
1451 parameter = new G4UIparameter ("blue", 'd', omitable = true);
1452 parameter->SetDefaultValue (0.);
1453 fpCommand->SetParameter (parameter);
1454 parameter = new G4UIparameter ("placement", 's', omitable = true);
1455 parameter -> SetParameterCandidates("auto manual");
1456 parameter->SetDefaultValue ("auto");
1457 fpCommand->SetParameter (parameter);
1458 parameter = new G4UIparameter ("xmid", 'd', omitable = true);
1459 parameter->SetDefaultValue (0.);
1460 fpCommand->SetParameter (parameter);
1461 parameter = new G4UIparameter ("ymid", 'd', omitable = true);
1462 parameter->SetDefaultValue (0.);
1463 fpCommand->SetParameter (parameter);
1464 parameter = new G4UIparameter ("zmid", 'd', omitable = true);
1465 parameter->SetDefaultValue (0.);
1466 fpCommand->SetParameter (parameter);
1467 parameter = new G4UIparameter ("unit", 's', omitable = true);
1468 parameter->SetDefaultValue ("m");
1469 fpCommand->SetParameter (parameter);
1470}
1471
1473 delete fpCommand;
1474}
1475
1477 return "";
1478}
1479
1481
1483 G4bool warn = verbosity >= G4VisManager::warnings;
1484
1486 if (!pScene) {
1487 if (verbosity >= G4VisManager::errors) {
1488 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
1489 }
1490 return;
1491 } else {
1492 if (pScene->GetExtent().GetExtentRadius() <= 0.) {
1493 if (verbosity >= G4VisManager::errors) {
1494 G4cerr
1495 << "ERROR: Scene has no extent. Add volumes or use \"/vis/scene/add/extent\"."
1496 << G4endl;
1497 }
1498 return;
1499 }
1500 }
1501
1503 if (!pViewer) {
1504 if (verbosity >= G4VisManager::errors) {
1505 G4cerr <<
1506 "ERROR: G4VisCommandSceneAddLogo::SetNewValue: no viewer."
1507 "\n Auto direction needs a viewer."
1508 << G4endl;
1509 }
1510 return;
1511 }
1512
1513 G4double userHeight, red, green, blue, xmid, ymid, zmid;
1514 G4String userHeightUnit, direction, auto_manual, positionUnit;
1515 std::istringstream is (newValue);
1516 is >> userHeight >> userHeightUnit >> direction
1517 >> red >> green >> blue
1518 >> auto_manual
1519 >> xmid >> ymid >> zmid >> positionUnit;
1520
1521 G4double height = userHeight;
1522 const G4VisExtent& sceneExtent = pScene->GetExtent(); // Existing extent.
1523 if (userHeightUnit == "auto") {
1524 height *= 0.2 * sceneExtent.GetExtentRadius();
1525 } else {
1526 height *= G4UIcommand::ValueOf(userHeightUnit);
1527 }
1528
1529 G4double unit = G4UIcommand::ValueOf(positionUnit);
1530 xmid *= unit; ymid *= unit; zmid *= unit;
1531
1532 Direction logoDirection = X; // Initialise to keep some compilers happy.
1533 if (direction == "auto") {
1534 // Take cue from viewer
1535 const G4Vector3D& vp =
1537 if (vp.x() > vp.y() && vp.x() > vp.z()) logoDirection = X;
1538 else if (vp.x() < vp.y() && vp.x() < vp.z()) logoDirection = minusX;
1539 else if (vp.y() > vp.x() && vp.y() > vp.z()) logoDirection = Y;
1540 else if (vp.y() < vp.x() && vp.y() < vp.z()) logoDirection = minusY;
1541 else if (vp.z() > vp.x() && vp.z() > vp.y()) logoDirection = Z;
1542 else if (vp.z() < vp.x() && vp.z() < vp.y()) logoDirection = minusZ;
1543 }
1544 else if (direction(0) == 'x') logoDirection = X;
1545 else if (direction(0) == 'y') logoDirection = Y;
1546 else if (direction(0) == 'z') logoDirection = Z;
1547 else if (direction(0) == '-') {
1548 if (direction(1) == 'x') logoDirection = minusX;
1549 else if (direction(1) == 'y') logoDirection = minusY;
1550 else if (direction(1) == 'z') logoDirection = minusZ;
1551 } else {
1552 if (verbosity >= G4VisManager::errors) {
1553 G4cerr << "ERROR: Unrecogniseed direction: \""
1554 << direction << "\"." << G4endl;
1555 return;
1556 }
1557 }
1558
1559 G4bool autoPlacing = false; if (auto_manual == "auto") autoPlacing = true;
1560 // Parameters read and interpreted.
1561
1562 // Current scene extent
1563 const G4double xmin = sceneExtent.GetXmin();
1564 const G4double xmax = sceneExtent.GetXmax();
1565 const G4double ymin = sceneExtent.GetYmin();
1566 const G4double ymax = sceneExtent.GetYmax();
1567 const G4double zmin = sceneExtent.GetZmin();
1568 const G4double zmax = sceneExtent.GetZmax();
1569
1570 // Test existing extent and issue warnings...
1571 G4bool worried = false;
1572 if (sceneExtent.GetExtentRadius() == 0) {
1573 worried = true;
1574 if (verbosity >= G4VisManager::warnings) {
1575 G4cout <<
1576 "WARNING: Existing scene does not yet have any extent."
1577 "\n Maybe you have not yet added any geometrical object."
1578 << G4endl;
1579 }
1580 }
1581
1582 // Useful constants, etc...
1583 const G4double halfHeight(height / 2.);
1584 const G4double comfort(0.01); // 0.15 seems too big. 0.05 might be better.
1585 const G4double freeHeightFraction (1. + 2. * comfort);
1586
1587 // Test existing scene for room...
1588 G4bool room = true;
1589 switch (logoDirection) {
1590 case X:
1591 case minusX:
1592 if (freeHeightFraction * (xmax - xmin) < height) room = false;
1593 break;
1594 case Y:
1595 case minusY:
1596 if (freeHeightFraction * (ymax - ymin) < height) room = false;
1597 break;
1598 case Z:
1599 case minusZ:
1600 if (freeHeightFraction * (zmax - zmin) < height) room = false;
1601 break;
1602 }
1603 if (!room) {
1604 worried = true;
1605 if (verbosity >= G4VisManager::warnings) {
1606 G4cout <<
1607 "WARNING: Not enough room in existing scene. Maybe logo is too large."
1608 << G4endl;
1609 }
1610 }
1611 if (worried) {
1612 if (verbosity >= G4VisManager::warnings) {
1613 G4cout <<
1614 "WARNING: The logo you have asked for is bigger than the existing"
1615 "\n scene. Maybe you have added it too soon. It is recommended that"
1616 "\n you add the logo last so that it can be correctly auto-positioned"
1617 "\n so as not to be obscured by any existing object and so that the"
1618 "\n view parameters can be correctly recalculated."
1619 << G4endl;
1620 }
1621 }
1622
1623 G4double sxmid(xmid), symid(ymid), szmid(zmid);
1624 if (autoPlacing) {
1625 // Aim to place at bottom right of screen when viewed from logoDirection.
1626 // Give some comfort zone.
1627 const G4double xComfort = comfort * (xmax - xmin);
1628 const G4double yComfort = comfort * (ymax - ymin);
1629 const G4double zComfort = comfort * (zmax - zmin);
1630 switch (logoDirection) {
1631 case X: // y-axis up, z-axis to left?
1632 sxmid = xmax + halfHeight + xComfort;
1633 symid = ymin - yComfort;
1634 szmid = zmin - zComfort;
1635 break;
1636 case minusX: // y-axis up, z-axis to right?
1637 sxmid = xmin - halfHeight - xComfort;
1638 symid = ymin - yComfort;
1639 szmid = zmax + zComfort;
1640 break;
1641 case Y: // z-axis up, x-axis to left?
1642 sxmid = xmin - xComfort;
1643 symid = ymax + halfHeight + yComfort;
1644 szmid = zmin - zComfort;
1645 break;
1646 case minusY: // z-axis up, x-axis to right?
1647 sxmid = xmax + xComfort;
1648 symid = ymin - halfHeight - yComfort;
1649 szmid = zmin - zComfort;
1650 break;
1651 case Z: // y-axis up, x-axis to right?
1652 sxmid = xmax + xComfort;
1653 symid = ymin - yComfort;
1654 szmid = zmax + halfHeight + zComfort;
1655 break;
1656 case minusZ: // y-axis up, x-axis to left?
1657 sxmid = xmin - xComfort;
1658 symid = ymin - yComfort;
1659 szmid = zmin - halfHeight - zComfort;
1660 break;
1661 }
1662 }
1663
1664 G4Transform3D transform;
1665 switch (logoDirection) {
1666 case X: // y-axis up, z-axis to left?
1667 transform = G4RotateY3D(halfpi);
1668 break;
1669 case minusX: // y-axis up, z-axis to right?
1670 transform = G4RotateY3D(-halfpi);
1671 break;
1672 case Y: // z-axis up, x-axis to left?
1673 transform = G4RotateX3D(-halfpi) * G4RotateZ3D(pi);
1674 break;
1675 case minusY: // z-axis up, x-axis to right?
1676 transform = G4RotateX3D(halfpi);
1677 break;
1678 case Z: // y-axis up, x-axis to right?
1679 // No transformation required.
1680 break;
1681 case minusZ: // y-axis up, x-axis to left?
1682 transform = G4RotateY3D(pi);
1683 break;
1684 }
1685 transform = G4Translate3D(sxmid,symid,szmid) * transform;
1686
1687 G4VisAttributes visAtts(G4Colour(red, green, blue));
1688 visAtts.SetForceSolid(true); // Always solid.
1689
1690 G4Logo* logo = new G4Logo(height,visAtts);
1691 G4VModel* model =
1693 model->SetType("G4Logo");
1694 model->SetGlobalTag("G4Logo");
1695 model->SetGlobalDescription("G4Logo: " + newValue);
1696 model->SetTransformation(transform);
1697 // Note: it is the responsibility of the model to act upon this, but
1698 // the extent is in local coordinates...
1699 G4double& h = height;
1700 G4double h2 = h/2.;
1701 G4VisExtent extent(-h,h,-h2,h2,-h2,h2);
1702 model->SetExtent(extent);
1703 // This extent gets "added" to existing scene extent in
1704 // AddRunDurationModel below.
1705 const G4String& currentSceneName = pScene -> GetName ();
1706 G4bool successful = pScene -> AddRunDurationModel (model, warn);
1707 if (successful) {
1708 if (verbosity >= G4VisManager::confirmations) {
1709 G4cout << "G4 Logo of height " << userHeight << ' ' << userHeightUnit
1710 << ", " << direction << "-direction, added to scene \""
1711 << currentSceneName << "\"";
1712 if (verbosity >= G4VisManager::parameters) {
1713 G4cout << "\n with extent " << extent
1714 << "\n at " << transform.getRotation()
1715 << transform.getTranslation();
1716 }
1717 G4cout << G4endl;
1718 }
1719 }
1720 else G4VisCommandsSceneAddUnsuccessful(verbosity);
1721
1723}
1724
1725G4VisCommandSceneAddLogo::G4Logo::G4Logo
1726(G4double height, const G4VisAttributes& visAtts):
1727 fVisAtts(visAtts)
1728{
1729 const G4double& h = height;
1730 const G4double h2 = 0.5 * h; // Half height.
1731 const G4double ri = 0.25 * h; // Inner radius.
1732 const G4double ro = 0.5 * h; // Outer radius.
1733 const G4double ro2 = 0.5 * ro; // Half outer radius.
1734 const G4double w = ro - ri; // Width.
1735 const G4double w2 = 0.5 * w; // Half width.
1736 const G4double d2 = 0.2 * h; // Half depth.
1737 const G4double f1 = 0.05 * h; // left edge of stem of "4".
1738 const G4double f2 = -0.3 * h; // bottom edge of cross of "4".
1739 const G4double e = 1.e-4 * h; // epsilon.
1740 const G4double xt = f1, yt = h2; // Top of slope.
1741 const G4double xb = -h2, yb = f2 + w; // Bottom of slope.
1742 const G4double dx = xt - xb, dy = yt - yb;
1743 const G4double angle = std::atan2(dy,dx);
1745 rm.rotateZ(angle*rad);
1746 const G4double d = std::sqrt(dx * dx + dy * dy);
1747 const G4double ss = h; // Half height of square subtractor
1748 const G4double y8 = ss; // Choose y of subtractor for outer slope.
1749 const G4double x8 = ((-ss * d - dx * (yt - y8)) / dy) + xt;
1750 G4double y9 = ss; // Choose y of subtractor for inner slope.
1751 G4double x9 = ((-(ss - w) * d - dx * (yt - y8)) / dy) + xt;
1752 // But to get inner, we make a triangle translated by...
1753 const G4double xtr = ss - f1, ytr = -ss - f2 -w;
1754 x9 += xtr; y9 += ytr;
1755
1756 // The idea here is to create a polyhedron for the G and the 4. To do
1757 // this we use Geant4 geometry solids and make boolean operations.
1758 // Note that we do not need to keep the solids. We use local objects,
1759 // which, of course, are deleted on leaving this function. This
1760 // is contrary to the usual requirement for solids that are part of the
1761 // detector for which solids MUST be created on the heap (with "new").
1762 // Finally we invoke CreatePolyhedron, which creates a polyhedron on the heap
1763 // and returns a pointer. It is the user's responsibility to delete,
1764 // which is done in the destructor of this class. Thus the polyhedra,
1765 // created here, remain on the heap for the lifetime of the job.
1766
1767 // G...
1768 G4Tubs tG("tG",ri,ro,d2,0.15*pi,1.85*pi);
1769 G4Box bG("bG",w2,ro2,d2);
1770 G4UnionSolid logoG("logoG",&tG,&bG,G4Translate3D(ri+w2,-ro2,0.));
1771 fpG = logoG.CreatePolyhedron();
1772 fpG->SetVisAttributes(&fVisAtts);
1773 fpG->Transform(G4Translate3D(-0.55*h,0.,0.));
1774
1775 // 4...
1776 G4Box b1("b1",h2,h2,d2);
1777 G4Box bS("bS",ss,ss,d2+e); // Subtractor.
1778 G4Box bS2("bS2",ss,ss,d2+2.*e); // 2nd Subtractor.
1779 G4SubtractionSolid s1("s1",&b1,&bS,G4Translate3D(f1-ss,f2-ss,0.));
1780 G4SubtractionSolid s2("s2",&s1,&bS,G4Translate3D(f1+ss+w,f2-ss,0.));
1781 G4SubtractionSolid s3("s3",&s2,&bS,G4Translate3D(f1+ss+w,f2+ss+w,0.));
1783 ("s4",&s3,&bS,G4Transform3D(rm,G4ThreeVector(x8,y8,0.)));
1784 G4SubtractionSolid s5 // Triangular hole.
1785 ("s5",&bS,&bS2,G4Transform3D(rm,G4ThreeVector(x9,y9,0.)));
1786 G4SubtractionSolid logo4("logo4",&s4,&s5,G4Translate3D(-xtr,-ytr,0.));
1787 fp4 = logo4.CreatePolyhedron();
1788 /* Experiment with creating own polyhedron...
1789 int nNodes = 4;
1790 int nFaces = 4;
1791 double xyz[][3] = {{0,0,0},{1*m,0,0},{0,1*m,0},{0,0,1*m}};
1792 int faces[][4] = {{1,3,2,0},{1,2,4,0},{1,4,3,0},{2,3,4,0}};
1793 fp4 = new G4Polyhedron();
1794 fp4->createPolyhedron(nNodes,nFaces,xyz,faces);
1795 */
1796 fp4->SetVisAttributes(&fVisAtts);
1797 fp4->Transform(G4Translate3D(0.55*h,0.,0.));
1798}
1799
1800G4VisCommandSceneAddLogo::G4Logo::~G4Logo() {
1801 delete fpG;
1802 delete fp4;
1803}
1804
1805void G4VisCommandSceneAddLogo::G4Logo::operator()
1806 (G4VGraphicsScene& sceneHandler, const G4Transform3D& transform, const G4ModelingParameters*) {
1807 sceneHandler.BeginPrimitives(transform);
1808 sceneHandler.AddPrimitive(*fpG);
1809 sceneHandler.AddPrimitive(*fp4);
1810 sceneHandler.EndPrimitives();
1811}
1812
1813////////////// /vis/scene/add/logo2D ///////////////////////////////////////
1814
1816 G4bool omitable;
1817 fpCommand = new G4UIcommand ("/vis/scene/add/logo2D", this);
1818 fpCommand -> SetGuidance ("Adds 2D logo to current scene.");
1819 G4UIparameter* parameter;
1820 parameter = new G4UIparameter ("size", 'i', omitable = true);
1821 parameter -> SetGuidance ("Screen size of text in pixels.");
1822 parameter -> SetDefaultValue (48);
1823 fpCommand -> SetParameter (parameter);
1824 parameter = new G4UIparameter ("x-position", 'd', omitable = true);
1825 parameter -> SetGuidance ("x screen position in range -1 < x < 1.");
1826 parameter -> SetDefaultValue (-0.9);
1827 fpCommand -> SetParameter (parameter);
1828 parameter = new G4UIparameter ("y-position", 'd', omitable = true);
1829 parameter -> SetGuidance ("y screen position in range -1 < y < 1.");
1830 parameter -> SetDefaultValue (-0.9);
1831 fpCommand -> SetParameter (parameter);
1832 parameter = new G4UIparameter ("layout", 's', omitable = true);
1833 parameter -> SetGuidance ("Layout, i.e., adjustment: left|centre|right.");
1834 parameter -> SetDefaultValue ("left");
1835 fpCommand -> SetParameter (parameter);
1836}
1837
1839 delete fpCommand;
1840}
1841
1843 return "";
1844}
1845
1847{
1849 G4bool warn(verbosity >= G4VisManager::warnings);
1850
1852 if (!pScene) {
1853 if (verbosity >= G4VisManager::errors) {
1854 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
1855 }
1856 return;
1857 }
1858
1859 G4int size;
1860 G4double x, y;
1861 G4String layoutString;
1862 std::istringstream is(newValue);
1863 is >> size >> x >> y >> layoutString;
1865 if (layoutString(0) == 'l') layout = G4Text::left;
1866 else if (layoutString(0) == 'c') layout = G4Text::centre;
1867 else if (layoutString(0) == 'r') layout = G4Text::right;
1868
1869 Logo2D* logo2D = new Logo2D(fpVisManager, size, x, y, layout);
1870 G4VModel* model =
1872 model->SetType("G4Logo2D");
1873 model->SetGlobalTag("G4Logo2D");
1874 model->SetGlobalDescription("G4Logo2D: " + newValue);
1875 const G4String& currentSceneName = pScene -> GetName ();
1876 G4bool successful = pScene -> AddRunDurationModel (model, warn);
1877 if (successful) {
1878 if (verbosity >= G4VisManager::confirmations) {
1879 G4cout << "2D logo has been added to scene \""
1880 << currentSceneName << "\"."
1881 << G4endl;
1882 }
1883 }
1884 else G4VisCommandsSceneAddUnsuccessful(verbosity);
1885
1887}
1888
1889void G4VisCommandSceneAddLogo2D::Logo2D::operator()
1890 (G4VGraphicsScene& sceneHandler, const G4Transform3D&, const G4ModelingParameters*)
1891{
1892 G4Text text("Geant4", G4Point3D(fX, fY, 0.));
1893 text.SetScreenSize(fSize);
1894 text.SetLayout(fLayout);
1895 G4VisAttributes textAtts(G4Colour::Brown());
1896 text.SetVisAttributes(textAtts);
1897 sceneHandler.BeginPrimitives2D();
1898 sceneHandler.AddPrimitive(text);
1899 sceneHandler.EndPrimitives2D();
1900}
1901
1902////////////// /vis/scene/add/magneticField ///////////////////////////////////////
1903
1905 fpCommand = new G4UIcommand ("/vis/scene/add/magneticField", this);
1906 fpCommand -> SetGuidance
1907 ("Adds magnetic field representation to current scene.");
1909 const G4UIcommand* addElecFieldCmd = tree->FindPath("/vis/scene/add/electricField");
1910 // Pick up additional guidance from /vis/scene/add/electricField
1911 CopyGuidanceFrom(addElecFieldCmd,fpCommand,1);
1912 // Pick up parameters from /vis/scene/add/electricField
1913 CopyParametersFrom(addElecFieldCmd,fpCommand);
1914}
1915
1917 delete fpCommand;
1918}
1919
1921 return "";
1922}
1923
1925(G4UIcommand*, G4String newValue) {
1926
1928 G4bool warn(verbosity >= G4VisManager::warnings);
1929
1931 if (!pScene) {
1932 if (verbosity >= G4VisManager::errors) {
1933 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
1934 }
1935 return;
1936 }
1937
1938 G4int nDataPointsPerHalfScene;
1939 G4String representation;
1940 std::istringstream iss(newValue);
1941 iss >> nDataPointsPerHalfScene >> representation;
1943 modelRepresentation = G4ElectricFieldModel::fullArrow;
1944 if (representation == "lightArrow") {
1945 modelRepresentation = G4ElectricFieldModel::lightArrow;
1946 }
1947 G4VModel* model;
1948 model = new G4MagneticFieldModel
1949 (nDataPointsPerHalfScene,modelRepresentation,
1953 const G4String& currentSceneName = pScene -> GetName ();
1954 G4bool successful = pScene -> AddRunDurationModel (model, warn);
1955 if (successful) {
1956 if (verbosity >= G4VisManager::confirmations) {
1957 G4cout
1958 << "Magnetic field, if any, will be drawn in scene \""
1959 << currentSceneName
1960 << "\"\n with "
1961 << nDataPointsPerHalfScene
1962 << " data points per half extent and with representation \""
1963 << representation
1964 << '\"'
1965 << G4endl;
1966 }
1967 }
1968 else G4VisCommandsSceneAddUnsuccessful(verbosity);
1969
1971}
1972
1973////////////// /vis/scene/add/psHits ///////////////////////////////////////
1974
1976 G4bool omitable;
1977 fpCommand = new G4UIcmdWithAString ("/vis/scene/add/psHits", this);
1978 fpCommand -> SetGuidance
1979 ("Adds Primitive Scorer Hits (PSHits) to current scene.");
1980 fpCommand -> SetGuidance
1981 ("PSHits are drawn at end of run when the scene in which"
1982 "\nthey are added is current.");
1983 fpCommand -> SetGuidance
1984 ("Optional parameter specifies name of scoring map. By default all"
1985 "\nscoring maps registered with the G4ScoringManager are drawn.");
1986 fpCommand -> SetParameterName ("mapname", omitable = true);
1987 fpCommand -> SetDefaultValue ("all");
1988}
1989
1991 delete fpCommand;
1992}
1993
1995 return "";
1996}
1997
1999(G4UIcommand*, G4String newValue)
2000{
2002 G4bool warn(verbosity >= G4VisManager::warnings);
2003
2005 if (!pScene) {
2006 if (verbosity >= G4VisManager::errors) {
2007 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
2008 }
2009 return;
2010 }
2011
2012 G4VModel* model = new G4PSHitsModel(newValue);
2013 const G4String& currentSceneName = pScene -> GetName ();
2014 G4bool successful = pScene -> AddEndOfRunModel (model, warn);
2015 if (successful) {
2016 if (verbosity >= G4VisManager::confirmations) {
2017 if (newValue == "all") {
2018 G4cout << "All Primitive Scorer hits";
2019 } else {
2020 G4cout << "Hits of Primitive Scorer \"" << newValue << '"';
2021 }
2022 G4cout << " will be drawn at end of run in scene \""
2023 << currentSceneName << "\"."
2024 << G4endl;
2025 }
2026 }
2027 else G4VisCommandsSceneAddUnsuccessful(verbosity);
2028
2030}
2031
2032////////////// /vis/scene/add/scale //////////////////////////////////
2033
2035 G4bool omitable;
2036 fpCommand = new G4UIcommand ("/vis/scene/add/scale", this);
2037 fpCommand -> SetGuidance
2038 ("Adds an annotated scale line to the current scene.");
2039 fpCommand -> SetGuidance
2040 ("If \"unit\" is \"auto\", length is roughly one tenth of the scene extent.");
2041 fpCommand -> SetGuidance
2042 ("If \"direction\" is \"auto\", scale is roughly in the plane of the current view.");
2043 fpCommand -> SetGuidance
2044 ("If \"placement\" is \"auto\", scale is placed at bottom left of current view."
2045 "\n Otherwise placed at (xmid,ymid,zmid).");
2046 fpCommand -> SetGuidance (G4Scale::GetGuidanceString());
2047 G4UIparameter* parameter;
2048 parameter = new G4UIparameter ("length", 'd', omitable = true);
2049 parameter->SetDefaultValue (1.);
2050 fpCommand->SetParameter (parameter);
2051 parameter = new G4UIparameter ("unit", 's', omitable = true);
2052 parameter->SetDefaultValue ("auto");
2053 fpCommand->SetParameter (parameter);
2054 parameter = new G4UIparameter ("direction", 's', omitable = true);
2055 parameter->SetGuidance ("auto|x|y|z");
2056 parameter->SetDefaultValue ("auto");
2057 fpCommand->SetParameter (parameter);
2058 parameter = new G4UIparameter ("red", 'd', omitable = true);
2059 parameter->SetDefaultValue (1.);
2060 fpCommand->SetParameter (parameter);
2061 parameter = new G4UIparameter ("green", 'd', omitable = true);
2062 parameter->SetDefaultValue (0.);
2063 fpCommand->SetParameter (parameter);
2064 parameter = new G4UIparameter ("blue", 'd', omitable = true);
2065 parameter->SetDefaultValue (0.);
2066 fpCommand->SetParameter (parameter);
2067 parameter = new G4UIparameter ("placement", 's', omitable = true);
2068 parameter -> SetParameterCandidates("auto manual");
2069 parameter->SetDefaultValue ("auto");
2070 fpCommand->SetParameter (parameter);
2071 parameter = new G4UIparameter ("xmid", 'd', omitable = true);
2072 parameter->SetDefaultValue (0.);
2073 fpCommand->SetParameter (parameter);
2074 parameter = new G4UIparameter ("ymid", 'd', omitable = true);
2075 parameter->SetDefaultValue (0.);
2076 fpCommand->SetParameter (parameter);
2077 parameter = new G4UIparameter ("zmid", 'd', omitable = true);
2078 parameter->SetDefaultValue (0.);
2079 fpCommand->SetParameter (parameter);
2080 parameter = new G4UIparameter ("unit", 's', omitable = true);
2081 parameter->SetDefaultValue ("m");
2082 fpCommand->SetParameter (parameter);
2083}
2084
2086 delete fpCommand;
2087}
2088
2090 return "";
2091}
2092
2094
2096 G4bool warn = verbosity >= G4VisManager::warnings;
2097
2099 if (!pScene) {
2100 if (verbosity >= G4VisManager::errors) {
2101 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
2102 }
2103 return;
2104 } else {
2105 if (pScene->GetExtent().GetExtentRadius() <= 0.) {
2106 if (verbosity >= G4VisManager::errors) {
2107 G4cerr
2108 << "ERROR: Scene has no extent. Add volumes or use \"/vis/scene/add/extent\"."
2109 << G4endl;
2110 }
2111 return;
2112 }
2113 }
2114
2115 G4double userLength, red, green, blue, xmid, ymid, zmid;
2116 G4String userLengthUnit, direction, auto_manual, positionUnit;
2117 std::istringstream is (newValue);
2118 is >> userLength >> userLengthUnit >> direction
2119 >> red >> green >> blue
2120 >> auto_manual
2121 >> xmid >> ymid >> zmid >> positionUnit;
2122
2123 G4double length = userLength;
2124 const G4VisExtent& sceneExtent = pScene->GetExtent(); // Existing extent.
2125 if (userLengthUnit == "auto") {
2126 const G4double lengthMax = 0.5 * sceneExtent.GetExtentRadius();
2127 const G4double intLog10Length = std::floor(std::log10(lengthMax));
2128 length = std::pow(10,intLog10Length);
2129 if (5.*length < lengthMax) length *= 5.;
2130 else if (2.*length < lengthMax) length *= 2.;
2131 } else {
2132 length *= G4UIcommand::ValueOf(userLengthUnit);
2133 }
2134 G4String annotation = G4BestUnit(length,"Length");
2135
2136 G4double unit = G4UIcommand::ValueOf(positionUnit);
2137 xmid *= unit; ymid *= unit; zmid *= unit;
2138
2139 G4Scale::Direction scaleDirection (G4Scale::x);
2140 if (direction(0) == 'y') scaleDirection = G4Scale::y;
2141 if (direction(0) == 'z') scaleDirection = G4Scale::z;
2142
2144 if (!pViewer) {
2145 if (verbosity >= G4VisManager::errors) {
2146 G4cerr <<
2147 "ERROR: G4VisCommandSceneAddScale::SetNewValue: no viewer."
2148 "\n Auto direction needs a viewer."
2149 << G4endl;
2150 }
2151 return;
2152 }
2153
2154 const G4Vector3D& vp =
2156 const G4Vector3D& up =
2157 pViewer->GetViewParameters().GetUpVector();
2158
2159 if (direction == "auto") { // Takes cue from viewer.
2160 if (std::abs(vp.x()) > std::abs(vp.y()) &&
2161 std::abs(vp.x()) > std::abs(vp.z())) { // x viewpoint
2162 if (std::abs(up.y()) > std::abs(up.z())) scaleDirection = G4Scale::z;
2163 else scaleDirection = G4Scale::y;
2164 }
2165 else if (std::abs(vp.y()) > std::abs(vp.x()) &&
2166 std::abs(vp.y()) > std::abs(vp.z())) { // y viewpoint
2167 if (std::abs(up.x()) > std::abs(up.z())) scaleDirection = G4Scale::z;
2168 else scaleDirection = G4Scale::x;
2169 }
2170 else if (std::abs(vp.z()) > std::abs(vp.x()) &&
2171 std::abs(vp.z()) > std::abs(vp.y())) { // z viewpoint
2172 if (std::abs(up.y()) > std::abs(up.x())) scaleDirection = G4Scale::x;
2173 else scaleDirection = G4Scale::y;
2174 }
2175 }
2176
2177 G4bool autoPlacing = false; if (auto_manual == "auto") autoPlacing = true;
2178 // Parameters read and interpreted.
2179
2180 // Useful constants, etc...
2181 const G4double halfLength(length / 2.);
2182 const G4double comfort(0.01); // 0.15 seems too big. 0.05 might be better.
2183 const G4double freeLengthFraction (1. + 2. * comfort);
2184
2185 const G4double xmin = sceneExtent.GetXmin();
2186 const G4double xmax = sceneExtent.GetXmax();
2187 const G4double ymin = sceneExtent.GetYmin();
2188 const G4double ymax = sceneExtent.GetYmax();
2189 const G4double zmin = sceneExtent.GetZmin();
2190 const G4double zmax = sceneExtent.GetZmax();
2191
2192 // Test existing extent and issue warnings...
2193 G4bool worried = false;
2194 if (sceneExtent.GetExtentRadius() == 0) {
2195 worried = true;
2196 if (verbosity >= G4VisManager::warnings) {
2197 G4cout <<
2198 "WARNING: Existing scene does not yet have any extent."
2199 "\n Maybe you have not yet added any geometrical object."
2200 << G4endl;
2201 }
2202 }
2203 // Test existing scene for room...
2204 G4bool room = true;
2205 switch (scaleDirection) {
2206 case G4Scale::x:
2207 if (freeLengthFraction * (xmax - xmin) < length) room = false;
2208 break;
2209 case G4Scale::y:
2210 if (freeLengthFraction * (ymax - ymin) < length) room = false;
2211 break;
2212 case G4Scale::z:
2213 if (freeLengthFraction * (zmax - zmin) < length) room = false;
2214 break;
2215 }
2216 if (!room) {
2217 worried = true;
2218 if (verbosity >= G4VisManager::warnings) {
2219 G4cout <<
2220 "WARNING: Not enough room in existing scene. Maybe scale is too long."
2221 << G4endl;
2222 }
2223 }
2224 if (worried) {
2225 if (verbosity >= G4VisManager::warnings) {
2226 G4cout <<
2227 "WARNING: The scale you have asked for is bigger than the existing"
2228 "\n scene. Maybe you have added it too soon. It is recommended that"
2229 "\n you add the scale last so that it can be correctly auto-positioned"
2230 "\n so as not to be obscured by any existing object and so that the"
2231 "\n view parameters can be correctly recalculated."
2232 << G4endl;
2233 }
2234 }
2235
2236 // Let's go ahead a construct a scale and a scale model. Since the
2237 // placing is done here, this G4Scale is *not* auto-placed...
2238 G4Scale scale(length, annotation, scaleDirection,
2239 false, xmid, ymid, zmid,
2241 G4VisAttributes visAttr(G4Colour(red, green, blue));
2242 scale.SetVisAttributes(visAttr);
2243 G4VModel* model = new G4ScaleModel(scale);
2244 G4String globalDescription = model->GetGlobalDescription();
2245 globalDescription += " (" + newValue + ")";
2246 model->SetGlobalDescription(globalDescription);
2247
2248 // Now figure out the extent...
2249 //
2250 // From the G4Scale.hh:
2251 //
2252 // This creates a representation of annotated line in the specified
2253 // direction with tick marks at the end. If autoPlacing is true it
2254 // is required to be centred at the front, right, bottom corner of
2255 // the world space, comfortably outside the existing bounding
2256 // box/sphere so that existing objects do not obscure it. Otherwise
2257 // it is required to be drawn with mid-point at (xmid, ymid, zmid).
2258 //
2259 // The auto placing algorithm might be:
2260 // x = xmin + (1 + comfort) * (xmax - xmin)
2261 // y = ymin - comfort * (ymax - ymin)
2262 // z = zmin + (1 + comfort) * (zmax - zmin)
2263 // if direction == x then (x - length,y,z) to (x,y,z)
2264 // if direction == y then (x,y,z) to (x,y + length,z)
2265 // if direction == z then (x,y,z - length) to (x,y,z)
2266 //
2267 // End of clip from G4Scale.hh:
2268 //
2269 // Implement this in two parts. Here, use the scale's extent to
2270 // "expand" the scene's extent. Then rendering - in
2271 // G4VSceneHandler::AddPrimitive(const G4Scale&) - simply has to
2272 // ensure it's within the new extent.
2273 //
2274
2275 G4double sxmid(xmid), symid(ymid), szmid(zmid);
2276 if (autoPlacing) {
2277 // Aim to place at bottom right of screen in current view.
2278 // Give some comfort zone.
2279 const G4double xComfort = comfort * (xmax - xmin);
2280 const G4double yComfort = comfort * (ymax - ymin);
2281 const G4double zComfort = comfort * (zmax - zmin);
2282 switch (scaleDirection) {
2283 case G4Scale::x:
2284 if (vp.z() > 0.) {
2285 sxmid = xmax + xComfort;
2286 symid = ymin - yComfort;
2287 szmid = zmin - zComfort;
2288 } else {
2289 sxmid = xmin - xComfort;
2290 symid = ymin - yComfort;
2291 szmid = zmax + zComfort;
2292 }
2293 break;
2294 case G4Scale::y:
2295 if (vp.x() > 0.) {
2296 sxmid = xmin - xComfort;
2297 symid = ymax + yComfort;
2298 szmid = zmin - zComfort;
2299 } else {
2300 sxmid = xmax + xComfort;
2301 symid = ymin - yComfort;
2302 szmid = zmin - zComfort;
2303 }
2304 break;
2305 case G4Scale::z:
2306 if (vp.x() > 0.) {
2307 sxmid = xmax + xComfort;
2308 symid = ymin - yComfort;
2309 szmid = zmax + zComfort;
2310 } else {
2311 sxmid = xmin - xComfort;
2312 symid = ymin - yComfort;
2313 szmid = zmax + zComfort;
2314 }
2315 break;
2316 }
2317 }
2318
2319 /* Old code - kept for future reference.
2320 G4double sxmid(xmid), symid(ymid), szmid(zmid);
2321 if (autoPlacing) {
2322 sxmid = xmin + onePlusComfort * (xmax - xmin);
2323 symid = ymin - comfort * (ymax - ymin);
2324 szmid = zmin + onePlusComfort * (zmax - zmin);
2325 switch (scaleDirection) {
2326 case G4Scale::x:
2327 sxmid -= halfLength;
2328 break;
2329 case G4Scale::y:
2330 symid += halfLength;
2331 break;
2332 case G4Scale::z:
2333 szmid -= halfLength;
2334 break;
2335 }
2336 }
2337 */
2338
2339 /* sxmin, etc., not actually used. Comment out to prevent compiler
2340 warnings but keep in case need in future. Extract transform and
2341 scaleExtent into reduced code below.
2342 G4double sxmin(sxmid), sxmax(sxmid);
2343 G4double symin(symid), symax(symid);
2344 G4double szmin(szmid), szmax(szmid);
2345 G4Transform3D transform;
2346 G4VisExtent scaleExtent;
2347 switch (scaleDirection) {
2348 case G4Scale::x:
2349 sxmin = sxmid - halfLength;
2350 sxmax = sxmid + halfLength;
2351 scaleExtent = G4VisExtent(-halfLength,halfLength,0,0,0,0);
2352 break;
2353 case G4Scale::y:
2354 symin = symid - halfLength;
2355 symax = symid + halfLength;
2356 transform = G4RotateZ3D(halfpi);
2357 scaleExtent = G4VisExtent(0,0,-halfLength,halfLength,0,0);
2358 break;
2359 case G4Scale::z:
2360 szmin = szmid - halfLength;
2361 szmax = szmid + halfLength;
2362 transform = G4RotateY3D(halfpi);
2363 scaleExtent = G4VisExtent(0,0,0,0,-halfLength,halfLength);
2364 break;
2365 }
2366 */
2367 G4Transform3D transform;
2368 G4VisExtent scaleExtent;
2369 switch (scaleDirection) {
2370 case G4Scale::x:
2371 scaleExtent = G4VisExtent(-halfLength,halfLength,0,0,0,0);
2372 break;
2373 case G4Scale::y:
2374 transform = G4RotateZ3D(halfpi);
2375 scaleExtent = G4VisExtent(0,0,-halfLength,halfLength,0,0);
2376 break;
2377 case G4Scale::z:
2378 transform = G4RotateY3D(halfpi);
2379 scaleExtent = G4VisExtent(0,0,0,0,-halfLength,halfLength);
2380 break;
2381 }
2382 transform = G4Translate3D(sxmid,symid,szmid) * transform;
2383 ///////// G4VisExtent scaleExtent(sxmin, sxmax, symin, symax, szmin, szmax);
2384
2385 model->SetTransformation(transform);
2386 // Note: it is the responsibility of the model to act upon this, but
2387 // the extent is in local coordinates...
2388 model->SetExtent(scaleExtent);
2389 // This extent gets "added" to existing scene extent in
2390 // AddRunDurationModel below.
2391
2392 const G4String& currentSceneName = pScene -> GetName ();
2393 G4bool successful = pScene -> AddRunDurationModel (model, warn);
2394 if (successful) {
2395 if (verbosity >= G4VisManager::confirmations) {
2396 G4cout << "Scale of " << annotation
2397 << " added to scene \"" << currentSceneName << "\".";
2398 if (verbosity >= G4VisManager::parameters) {
2399 G4cout << "\n with extent " << scaleExtent
2400 << "\n at " << transform.getRotation()
2401 << transform.getTranslation();
2402 }
2403 G4cout << G4endl;
2404 }
2405 }
2406 else G4VisCommandsSceneAddUnsuccessful(verbosity);
2407
2409}
2410
2411
2412////////////// /vis/scene/add/text //////////////////////////////////
2413
2415 G4bool omitable;
2416 fpCommand = new G4UIcommand ("/vis/scene/add/text", this);
2417 fpCommand -> SetGuidance ("Adds text to current scene.");
2418 fpCommand -> SetGuidance
2419 ("Use \"/vis/set/textColour\" to set colour.");
2420 fpCommand -> SetGuidance
2421 ("Use \"/vis/set/textLayout\" to set layout:");
2422 G4UIparameter* parameter;
2423 parameter = new G4UIparameter ("x", 'd', omitable = true);
2424 parameter->SetDefaultValue (0);
2425 fpCommand->SetParameter (parameter);
2426 parameter = new G4UIparameter ("y", 'd', omitable = true);
2427 parameter->SetDefaultValue (0);
2428 fpCommand->SetParameter (parameter);
2429 parameter = new G4UIparameter ("z", 'd', omitable = true);
2430 parameter->SetDefaultValue (0);
2431 fpCommand->SetParameter (parameter);
2432 parameter = new G4UIparameter ("unit", 's', omitable = true);
2433 parameter->SetDefaultValue ("m");
2434 fpCommand->SetParameter (parameter);
2435 parameter = new G4UIparameter ("font_size", 'd', omitable = true);
2436 parameter->SetDefaultValue (12);
2437 parameter->SetGuidance ("pixels");
2438 fpCommand->SetParameter (parameter);
2439 parameter = new G4UIparameter ("x_offset", 'd', omitable = true);
2440 parameter->SetDefaultValue (0);
2441 parameter->SetGuidance ("pixels");
2442 fpCommand->SetParameter (parameter);
2443 parameter = new G4UIparameter ("y_offset", 'd', omitable = true);
2444 parameter->SetDefaultValue (0);
2445 parameter->SetGuidance ("pixels");
2446 fpCommand->SetParameter (parameter);
2447 parameter = new G4UIparameter ("text", 's', omitable = true);
2448 parameter->SetGuidance ("The rest of the line is text.");
2449 parameter->SetDefaultValue ("Hello G4");
2450 fpCommand->SetParameter (parameter);
2451}
2452
2454 delete fpCommand;
2455}
2456
2458 return "";
2459}
2460
2462
2464 G4bool warn = verbosity >= G4VisManager::warnings;
2465
2467 if (!pScene) {
2468 if (verbosity >= G4VisManager::errors) {
2469 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
2470 }
2471 return;
2472 }
2473
2474 G4Tokenizer next(newValue);
2475 G4double x = StoD(next());
2476 G4double y = StoD(next());
2477 G4double z = StoD(next());
2478 G4String unitString = next();
2479 G4double font_size = StoD(next());
2480 G4double x_offset = StoD(next());
2481 G4double y_offset = StoD(next());
2482 G4String text = next("\n");
2483
2484 G4double unit = G4UIcommand::ValueOf(unitString);
2485 x *= unit; y *= unit; z *= unit;
2486
2487 G4Text g4text(text, G4Point3D(x,y,z));
2489 g4text.SetVisAttributes(visAtts);
2491 g4text.SetScreenSize(font_size);
2492 g4text.SetOffset(x_offset,y_offset);
2493 G4VModel* model = new G4TextModel(g4text);
2494 const G4String& currentSceneName = pScene -> GetName ();
2495 G4bool successful = pScene -> AddRunDurationModel (model, warn);
2496 if (successful) {
2497 if (verbosity >= G4VisManager::confirmations) {
2498 G4cout << "Text \"" << text
2499 << "\" has been added to scene \"" << currentSceneName << "\"."
2500 << G4endl;
2501 }
2502 }
2503 else G4VisCommandsSceneAddUnsuccessful(verbosity);
2504
2506}
2507
2508
2509////////////// /vis/scene/add/text2D //////////////////////////////////
2510
2512 G4bool omitable;
2513 fpCommand = new G4UIcommand ("/vis/scene/add/text2D", this);
2514 fpCommand -> SetGuidance ("Adds 2D text to current scene.");
2515 fpCommand -> SetGuidance
2516 ("Use \"/vis/set/textColour\" to set colour.");
2517 fpCommand -> SetGuidance
2518 ("Use \"/vis/set/textLayout\" to set layout:");
2519 G4UIparameter* parameter;
2520 parameter = new G4UIparameter ("x", 'd', omitable = true);
2521 parameter->SetDefaultValue (0);
2522 fpCommand->SetParameter (parameter);
2523 parameter = new G4UIparameter ("y", 'd', omitable = true);
2524 parameter->SetDefaultValue (0);
2525 fpCommand->SetParameter (parameter);
2526 parameter = new G4UIparameter ("font_size", 'd', omitable = true);
2527 parameter->SetDefaultValue (12);
2528 parameter->SetGuidance ("pixels");
2529 fpCommand->SetParameter (parameter);
2530 parameter = new G4UIparameter ("x_offset", 'd', omitable = true);
2531 parameter->SetDefaultValue (0);
2532 parameter->SetGuidance ("pixels");
2533 fpCommand->SetParameter (parameter);
2534 parameter = new G4UIparameter ("y_offset", 'd', omitable = true);
2535 parameter->SetDefaultValue (0);
2536 parameter->SetGuidance ("pixels");
2537 fpCommand->SetParameter (parameter);
2538 parameter = new G4UIparameter ("text", 's', omitable = true);
2539 parameter->SetGuidance ("The rest of the line is text.");
2540 parameter->SetDefaultValue ("Hello G4");
2541 fpCommand->SetParameter (parameter);
2542}
2543
2545 delete fpCommand;
2546}
2547
2549 return "";
2550}
2551
2553
2555 G4bool warn = verbosity >= G4VisManager::warnings;
2556
2558 if (!pScene) {
2559 if (verbosity >= G4VisManager::errors) {
2560 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
2561 }
2562 return;
2563 }
2564
2565 G4Tokenizer next(newValue);
2566 G4double x = StoD(next());
2567 G4double y = StoD(next());
2568 G4double font_size = StoD(next());
2569 G4double x_offset = StoD(next());
2570 G4double y_offset = StoD(next());
2571 G4String text = next("\n");
2572
2573 G4Text g4text(text, G4Point3D(x,y,0.));
2575 g4text.SetVisAttributes(visAtts);
2577 g4text.SetScreenSize(font_size);
2578 g4text.SetOffset(x_offset,y_offset);
2579 G4Text2D* g4text2D = new G4Text2D(g4text);
2580 G4VModel* model =
2582 model->SetType("Text2D");
2583 model->SetGlobalTag("Text2D");
2584 model->SetGlobalDescription("Text2D: " + newValue);
2585 const G4String& currentSceneName = pScene -> GetName ();
2586 G4bool successful = pScene -> AddRunDurationModel (model, warn);
2587 if (successful) {
2588 if (verbosity >= G4VisManager::confirmations) {
2589 G4cout << "2D text \"" << text
2590 << "\" has been added to scene \"" << currentSceneName << "\"."
2591 << G4endl;
2592 }
2593 }
2594 else G4VisCommandsSceneAddUnsuccessful(verbosity);
2595
2597}
2598
2599G4VisCommandSceneAddText2D::G4Text2D::G4Text2D(const G4Text& text):
2600 fText(text)
2601{}
2602
2603void G4VisCommandSceneAddText2D::G4Text2D::operator()
2604 (G4VGraphicsScene& sceneHandler, const G4Transform3D& transform, const G4ModelingParameters*) {
2605 sceneHandler.BeginPrimitives2D(transform);
2606 sceneHandler.AddPrimitive(fText);
2607 sceneHandler.EndPrimitives2D();
2608}
2609
2610
2611////////////// /vis/scene/add/trajectories ///////////////////////////////////
2612
2614 G4bool omitable;
2615 fpCommand = new G4UIcmdWithAString
2616 ("/vis/scene/add/trajectories", this);
2617 fpCommand -> SetGuidance
2618 ("Adds trajectories to current scene.");
2619 fpCommand -> SetGuidance
2620 ("Causes trajectories, if any, to be drawn at the end of processing an"
2621 "\nevent. Switches on trajectory storing and sets the"
2622 "\ndefault trajectory type.");
2623 fpCommand -> SetGuidance
2624 ("The command line parameter list determines the default trajectory type."
2625 "\nIf it contains the string \"smooth\", auxiliary inter-step points will"
2626 "\nbe inserted to improve the smoothness of the drawing of a curved"
2627 "\ntrajectory."
2628 "\nIf it contains the string \"rich\", significant extra information will"
2629 "\nbe stored in the trajectory (G4RichTrajectory) amenable to modeling"
2630 "\nand filtering with \"/vis/modeling/trajectories/create/drawByAttribute\""
2631 "\nand \"/vis/filtering/trajectories/create/attributeFilter\" commands."
2632 "\nIt may contain both strings in any order.");
2633 fpCommand -> SetGuidance
2634 ("\nTo switch off trajectory storing: \"/tracking/storeTrajectory 0\"."
2635 "\nSee also \"/vis/scene/endOfEventAction\".");
2636 fpCommand -> SetGuidance
2637 ("Note: This only sets the default. Independently of the result of this"
2638 "\ncommand, a user may instantiate a trajectory that overrides this default"
2639 "\nin PreUserTrackingAction.");
2640 fpCommand -> SetParameterName ("default-trajectory-type", omitable = true);
2641 fpCommand -> SetDefaultValue ("");
2642}
2643
2645 delete fpCommand;
2646}
2647
2649 return "";
2650}
2651
2653 G4String newValue) {
2654
2656 G4bool warn = verbosity >= G4VisManager::warnings;
2657
2659 if (!pScene) {
2660 if (verbosity >= G4VisManager::errors) {
2661 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
2662 }
2663 return;
2664 }
2665 const G4String& currentSceneName = pScene -> GetName ();
2666
2667 G4bool smooth = false;
2668 G4bool rich = false;
2669 if (newValue.find("smooth") != std::string::npos) smooth = true;
2670 if (newValue.find("rich") != std::string::npos) rich = true;
2671 if (newValue.size() && !(rich || smooth)) {
2672 if (verbosity >= G4VisManager::errors) {
2673 G4cerr << "ERROR: Unrecognised parameter \"" << newValue << "\""
2674 "\n No action taken."
2675 << G4endl;
2676 }
2677 return;
2678 }
2679
2681 G4int keepVerbose = UImanager->GetVerboseLevel();
2682 G4int newVerbose = 2;
2683 UImanager->SetVerboseLevel(newVerbose);
2684 G4String defaultTrajectoryType;
2685 if (smooth && rich) {
2686 UImanager->ApplyCommand("/tracking/storeTrajectory 4");
2687 defaultTrajectoryType = "G4RichTrajectory configured for smooth steps";
2688 } else if (smooth) {
2689 UImanager->ApplyCommand("/tracking/storeTrajectory 2");
2690 defaultTrajectoryType = "G4SmoothTrajectory";
2691 } else if (rich) {
2692 UImanager->ApplyCommand("/tracking/storeTrajectory 3");
2693 defaultTrajectoryType = "G4RichTrajectory";
2694 } else {
2695 UImanager->ApplyCommand("/tracking/storeTrajectory 1");
2696 defaultTrajectoryType = "G4Trajectory";
2697 }
2698 UImanager->SetVerboseLevel(keepVerbose);
2699
2700 if (verbosity >= G4VisManager::errors) {
2701 G4cout <<
2702 "Attributes available for modeling and filtering with"
2703 "\n \"/vis/modeling/trajectories/create/drawByAttribute\" and"
2704 "\n \"/vis/filtering/trajectories/create/attributeFilter\" commands:"
2705 << G4endl;
2707 if (rich) {
2710 } else if (smooth) {
2713 } else {
2716 }
2717 }
2718
2719 const auto& eoeList = pScene->GetEndOfEventModelList();
2720 auto eoeModel = eoeList.begin();
2721 for (; eoeModel != eoeList.end(); ++eoeModel) {
2722 const auto* actualModel = eoeModel->fpModel;
2723 if (dynamic_cast<const G4TrajectoriesModel*>(actualModel)) break;
2724 }
2725 if (eoeModel == eoeList.end()) {
2726 // No trajectories model exists in the scene so create a new one...
2727 G4VModel* model = new G4TrajectoriesModel();
2728 pScene -> AddEndOfEventModel (model, warn);
2729 } // ...else it already exists and there is no need to add a new one
2730 // because G4TrajectoriesModel simply describes trajectories in the
2731 // trajectories store whatever the type.
2732
2733 if (verbosity >= G4VisManager::confirmations) {
2734 G4cout << "Default trajectory type " << defaultTrajectoryType
2735 << "\n will be used to store trajectories for scene \""
2736 << currentSceneName << "\"."
2737 << G4endl;
2738 }
2739
2740 if (verbosity >= G4VisManager::warnings) {
2741 G4cout <<
2742 "WARNING: Trajectory storing has been requested. This action may be"
2743 "\n reversed with \"/tracking/storeTrajectory 0\"."
2744 << G4endl;
2745 }
2746
2748}
2749
2750////////////// /vis/scene/add/userAction ///////////////////////////////////
2751
2753 G4bool omitable;
2754 fpCommand = new G4UIcmdWithAString("/vis/scene/add/userAction",this);
2755 fpCommand -> SetGuidance
2756 ("Add named Vis User Action to current scene.");
2757 fpCommand -> SetGuidance
2758 ("Attempts to match search string to name of action - use unique sub-string.");
2759 fpCommand -> SetGuidance
2760 ("(Use /vis/list to see names of registered actions.)");
2761 fpCommand -> SetGuidance
2762 ("If name == \"all\" (default), all actions are added.");
2763 fpCommand -> SetParameterName("action-name", omitable = true);
2764 fpCommand -> SetDefaultValue("all");
2765}
2766
2768 delete fpCommand;
2769}
2770
2772 return "";
2773}
2774
2776(G4UIcommand*, G4String newValue) {
2777
2779
2781 if (!pScene) {
2782 if (verbosity >= G4VisManager::errors) {
2783 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
2784 }
2785 return;
2786 }
2787
2788 G4bool any = false;
2789
2790 const std::vector<G4VisManager::UserVisAction>& runDurationUserVisActions =
2792 for (size_t i = 0; i < runDurationUserVisActions.size(); i++) {
2793 const G4String& name = runDurationUserVisActions[i].fName;
2794 G4VUserVisAction* visAction = runDurationUserVisActions[i].fpUserVisAction;
2795 if (newValue == "all" || name.find(newValue) != std::string::npos) {
2796 any = true;
2797 AddVisAction(name,visAction,pScene,runDuration,verbosity);
2798 }
2799 }
2800
2801 const std::vector<G4VisManager::UserVisAction>& endOfEventUserVisActions =
2803 for (size_t i = 0; i < endOfEventUserVisActions.size(); i++) {
2804 const G4String& name = endOfEventUserVisActions[i].fName;
2805 G4VUserVisAction* visAction = endOfEventUserVisActions[i].fpUserVisAction;
2806 if (newValue == "all" || name.find(newValue) != std::string::npos) {
2807 any = true;
2808 AddVisAction(name,visAction,pScene,endOfEvent,verbosity);
2809 }
2810 }
2811
2812 const std::vector<G4VisManager::UserVisAction>& endOfRunUserVisActions =
2814 for (size_t i = 0; i < endOfRunUserVisActions.size(); i++) {
2815 const G4String& name = endOfRunUserVisActions[i].fName;
2816 G4VUserVisAction* visAction = endOfRunUserVisActions[i].fpUserVisAction;
2817 if (newValue == "all" || name.find(newValue) != std::string::npos) {
2818 any = true;
2819 AddVisAction(name,visAction,pScene,endOfRun,verbosity);
2820 }
2821 }
2822
2823 if (!any) {
2824 if (verbosity >= G4VisManager::warnings) {
2825 G4cout << "WARNING: No User Vis Action registered." << G4endl;
2826 }
2827 return;
2828 }
2829
2831}
2832
2833void G4VisCommandSceneAddUserAction::AddVisAction
2834(const G4String& name,
2835 G4VUserVisAction* visAction,
2836 G4Scene* pScene,
2837 G4VisCommandSceneAddUserAction::ActionType type,
2838 G4VisManager::Verbosity verbosity)
2839{
2840 G4bool warn = verbosity >= G4VisManager::warnings;
2841
2842 const std::map<G4VUserVisAction*,G4VisExtent>& visExtentMap =
2844 G4VisExtent extent;
2845 std::map<G4VUserVisAction*,G4VisExtent>::const_iterator i =
2846 visExtentMap.find(visAction);
2847 if (i != visExtentMap.end()) extent = i->second;
2848 if (warn) {
2849 if (extent.GetExtentRadius() <= 0.) {
2850 G4cout
2851 << "WARNING: User Vis Action \"" << name << "\" extent is null."
2852 << G4endl;
2853 }
2854 }
2855
2856 G4VModel* model = new G4CallbackModel<G4VUserVisAction>(visAction);
2857 model->SetType("User Vis Action");
2858 model->SetGlobalTag(name);
2859 model->SetGlobalDescription(name);
2860 model->SetExtent(extent);
2861 G4bool successful = false;;
2862 switch (type) {
2863 case runDuration:
2864 successful = pScene -> AddRunDurationModel (model, warn);
2865 break;
2866 case endOfEvent:
2867 successful = pScene -> AddEndOfEventModel (model, warn);
2868 break;
2869 case endOfRun:
2870 successful = pScene -> AddEndOfRunModel (model, warn);
2871 break;
2872 }
2873 if (successful) {
2874 if (verbosity >= G4VisManager::confirmations) {
2875 const G4String& currentSceneName = pScene -> GetName ();
2876 G4cout << "User Vis Action added to scene \""
2877 << currentSceneName << "\"";
2878 if (verbosity >= G4VisManager::parameters) {
2879 G4cout << "\n with extent " << extent;
2880 }
2881 G4cout << G4endl;
2882 }
2883 }
2884 else G4VisCommandsSceneAddUnsuccessful(verbosity);
2885}
2886
2887////////////// /vis/scene/add/volume ///////////////////////////////////////
2888
2890 G4bool omitable;
2891 fpCommand = new G4UIcommand ("/vis/scene/add/volume", this);
2892 fpCommand -> SetGuidance
2893 ("Adds a physical volume to current scene, with optional clipping volume.");
2894 fpCommand -> SetGuidance
2895 ("If physical-volume-name is \"world\" (the default), the top of the"
2896 "\nmain geometry tree (material world) is added. If \"worlds\", the"
2897 "\ntops of all worlds - material world and parallel worlds, if any - are"
2898 "\nadded. Otherwise a search of all worlds is made.");
2899 fpCommand -> SetGuidance
2900 ("In the last case the names of all volumes in all worlds are matched"
2901 "\nagainst physical-volume-name. If this is of the form \"/regexp/\","
2902 "\nwhere regexp is a regular expression (see C++ regex), the match uses"
2903 "\nthe usual rules of regular expression matching. Otherwise an exact"
2904 "\nmatch is required."
2905 "\nFor example, \"/Shap/\" adds \"Shape1\" and \"Shape2\".");
2906 fpCommand -> SetGuidance
2907 ("It may help to see a textual representation of the geometry hierarchy of"
2908 "\nthe worlds. Try \"/vis/drawTree [worlds]\" or one of the driver/browser"
2909 "\ncombinations that have the required functionality, e.g., HepRepFile.");
2910 fpCommand -> SetGuidance
2911 ("If clip-volume-type is specified, the subsequent parameters are used to"
2912 "\nto define a clipping volume. For example,"
2913 "\n\"/vis/scene/add/volume ! ! ! -box km 0 1 0 1 0 1\" will draw the world"
2914 "\nwith the positive octant cut away. (If the Boolean Processor issues"
2915 "\nwarnings try replacing 0 by 0.000000001 or something.)");
2916 fpCommand -> SetGuidance
2917 ("If clip-volume-type is prepended with '-', the clip-volume is subtracted"
2918 "\n(cutaway). (This is the default if there is no prepended character.)"
2919 "\nIf '*' is prepended, the intersection of the physical-volume and the"
2920 "\nclip-volume is made. (You can make a section through the detector with"
2921 "\na thin box, for example).");
2922 fpCommand -> SetGuidance
2923 ("For \"box\", the parameters are xmin,xmax,ymin,ymax,zmin,zmax."
2924 "\nOnly \"box\" is programmed at present.");
2925 G4UIparameter* parameter;
2926 parameter = new G4UIparameter ("physical-volume-name", 's', omitable = true);
2927 parameter -> SetDefaultValue ("world");
2928 fpCommand -> SetParameter (parameter);
2929 parameter = new G4UIparameter ("copy-no", 'i', omitable = true);
2930 parameter -> SetGuidance ("If negative, matches any copy no.");
2931 parameter -> SetDefaultValue (-1);
2932 fpCommand -> SetParameter (parameter);
2933 parameter = new G4UIparameter ("depth-of-descent", 'i', omitable = true);
2934 parameter -> SetGuidance
2935 ("Depth of descent of geometry hierarchy. Default = unlimited depth.");
2936 parameter -> SetDefaultValue (G4Scene::UNLIMITED);
2937 fpCommand -> SetParameter (parameter);
2938 parameter = new G4UIparameter ("clip-volume-type", 's', omitable = true);
2939 parameter -> SetParameterCandidates("none box -box *box");
2940 parameter -> SetDefaultValue ("none");
2941 parameter -> SetGuidance("[-|*]type. See general guidance.");
2942 fpCommand -> SetParameter (parameter);
2943 parameter = new G4UIparameter ("parameter-unit", 's', omitable = true);
2944 parameter -> SetDefaultValue ("m");
2945 fpCommand -> SetParameter (parameter);
2946 parameter = new G4UIparameter ("parameter-1", 'd', omitable = true);
2947 parameter -> SetDefaultValue (0.);
2948 fpCommand -> SetParameter (parameter);
2949 parameter = new G4UIparameter ("parameter-2", 'd', omitable = true);
2950 parameter -> SetDefaultValue (0.);
2951 fpCommand -> SetParameter (parameter);
2952 parameter = new G4UIparameter ("parameter-3", 'd', omitable = true);
2953 parameter -> SetDefaultValue (0.);
2954 fpCommand -> SetParameter (parameter);
2955 parameter = new G4UIparameter ("parameter-4", 'd', omitable = true);
2956 parameter -> SetDefaultValue (0.);
2957 fpCommand -> SetParameter (parameter);
2958 parameter = new G4UIparameter ("parameter-5", 'd', omitable = true);
2959 parameter -> SetDefaultValue (0.);
2960 fpCommand -> SetParameter (parameter);
2961 parameter = new G4UIparameter ("parameter-6", 'd', omitable = true);
2962 parameter -> SetDefaultValue (0.);
2963 fpCommand -> SetParameter (parameter);
2964}
2965
2967 delete fpCommand;
2968}
2969
2971 return "world 0 -1";
2972}
2973
2975 G4String newValue) {
2976
2978 G4bool warn = verbosity >= G4VisManager::warnings;
2979
2981 if (!pScene) {
2982 if (verbosity >= G4VisManager::errors) {
2983 G4cerr << "ERROR: No current scene. Please create one." << G4endl;
2984 }
2985 return;
2986 }
2987
2988 G4String name, clipVolumeType, parameterUnit;
2989 G4int copyNo, requestedDepthOfDescent;
2990 G4double param1, param2, param3, param4, param5, param6;
2991 std::istringstream is (newValue);
2992 is >> name >> copyNo >> requestedDepthOfDescent
2993 >> clipVolumeType >> parameterUnit
2994 >> param1 >> param2 >> param3 >> param4 >> param5 >> param6;
2996 G4PhysicalVolumeModel::subtraction; // Default subtraction mode.
2997 if (clipVolumeType[size_t(0)] == '-') {
2998 clipVolumeType = clipVolumeType.substr(1); // Remove first character.
2999 } else if (clipVolumeType[size_t(0)] == '*') {
3001 clipVolumeType = clipVolumeType.substr(1);
3002 }
3003 G4double unit = G4UIcommand::ValueOf(parameterUnit);
3004 param1 *= unit; param2 *= unit; param3 *= unit;
3005 param4 *= unit; param5 *= unit; param6 *= unit;
3006
3007 G4VSolid* clippingSolid = nullptr;
3008 if (clipVolumeType == "box") {
3009 const G4double dX = (param2 - param1) / 2.;
3010 const G4double dY = (param4 - param3) / 2.;
3011 const G4double dZ = (param6 - param5) / 2.;
3012 const G4double x0 = (param2 + param1) / 2.;
3013 const G4double y0 = (param4 + param3) / 2.;
3014 const G4double z0 = (param6 + param5) / 2.;
3015 clippingSolid = new G4DisplacedSolid
3016 ("_displaced_clipping_box",
3017 new G4Box("_clipping_box",dX,dY,dZ),
3018 G4Translate3D(x0,y0,z0));
3019 }
3020
3021 G4TransportationManager* transportationManager =
3023
3024 size_t nWorlds = transportationManager->GetNoWorlds();
3025 if (nWorlds > 1) { // Parallel worlds in operation...
3026 if (verbosity >= G4VisManager::warnings) {
3027 static G4bool warned = false;
3028 if (!warned && name != "worlds") {
3029 G4cout <<
3030 "WARNING: Parallel worlds in operation. To visualise, specify"
3031 "\n \"worlds\" or the parallel world volume or sub-volume name"
3032 "\n and control visibility with /vis/geometry."
3033 << G4endl;
3034 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
3035 transportationManager->GetWorldsIterator();
3036 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
3037 G4cout << " World " << i << ": " << (*iterWorld)->GetName()
3038 << G4endl;
3039 warned = true;
3040 }
3041 }
3042 }
3043 }
3044
3045 // Get the world (the initial value of the iterator points to the mass world).
3046 G4VPhysicalVolume* world = *(transportationManager->GetWorldsIterator());
3047
3048 if (!world) {
3049 if (verbosity >= G4VisManager::errors) {
3050 G4cerr <<
3051 "ERROR: G4VisCommandSceneAddVolume::SetNewValue:"
3052 "\n No world. Maybe the geometry has not yet been defined."
3053 "\n Try \"/run/initialize\""
3054 << G4endl;
3055 }
3056 return;
3057 }
3058
3059 std::vector<G4PhysicalVolumesSearchScene::Findings> findingsVector;
3060
3061 if (name == "world") {
3062
3063 findingsVector.push_back
3065
3066 } else if (name == "worlds") {
3067
3068 if (nWorlds <= 1) {
3069 if (verbosity >= G4VisManager::warnings) {
3070 G4cout <<
3071 "WARNING: G4VisCommandSceneAddVolume::SetNewValue:"
3072 "\n Parallel worlds requested but none exist."
3073 "\n Just adding material world."
3074 << G4endl;
3075 }
3076 }
3077 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
3078 transportationManager->GetWorldsIterator();
3079 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
3080 findingsVector.push_back
3082 (*iterWorld,*iterWorld));
3083 }
3084
3085 } else { // Search all worlds...
3086
3087 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
3088 transportationManager->GetWorldsIterator();
3089 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
3090 G4PhysicalVolumeModel searchModel (*iterWorld); // Unlimited depth.
3091 G4ModelingParameters mp; // Default - no culling.
3092 searchModel.SetModelingParameters (&mp);
3093 G4PhysicalVolumesSearchScene searchScene (&searchModel, name, copyNo);
3094 searchModel.DescribeYourselfTo (searchScene); // Initiate search.
3095 for (const auto& findings: searchScene.GetFindings()) {
3096 findingsVector.push_back(findings);
3097 }
3098 }
3099 }
3100
3101 for (const auto& findings: findingsVector) {
3102 // Set copy number from search findings for replicas and parameterisations.
3103 findings.fpFoundPV->SetCopyNo(findings.fFoundPVCopyNo);
3104 // Create a physical volume model...
3106 (findings.fpFoundPV,
3107 requestedDepthOfDescent,
3108 findings.fFoundObjectTransformation,
3109 0, // No modelling parameters (these are set later by the scene handler).
3110 false, // Do not use full extent - only use non-culled extent.
3111 findings.fFoundBasePVPath);
3112 if (clippingSolid) {
3113 foundPVModel->SetClippingSolid(clippingSolid);
3114 foundPVModel->SetClippingMode(clippingMode);
3115 }
3116 if (!foundPVModel->Validate(warn)) return; // Calculates extent
3117 // ...so add it to the scene.
3118 G4bool successful = pScene->AddRunDurationModel(foundPVModel,warn);
3119 if (successful) {
3120 if (verbosity >= G4VisManager::confirmations) {
3121 G4cout << "\"" << findings.fpFoundPV->GetName()
3122 << "\", copy no. " << findings.fFoundPVCopyNo
3123 << ",\n found in searched volume \""
3124 << findings.fpSearchPV->GetName()
3125 << "\" at depth " << findings.fFoundDepth
3126 << ",\n base path: \"" << findings.fFoundBasePVPath
3127 << "\",\n with a requested depth of further descent of ";
3128 if (requestedDepthOfDescent < 0) {
3129 G4cout << "<0 (unlimited)";
3130 }
3131 else {
3132 G4cout << requestedDepthOfDescent;
3133 }
3134 G4cout << ",\n has been added to scene \"" << pScene->GetName() << "\"."
3135 << G4endl;
3136 }
3137 } else {
3139 }
3140 }
3141
3142 if (findingsVector.empty()) {
3143 if (verbosity >= G4VisManager::errors) {
3144 G4cerr << "ERROR: Volume \"" << name << "\"";
3145 if (copyNo >= 0) {
3146 G4cerr << ", copy no. " << copyNo << ",";
3147 }
3148 G4cerr << " not found." << G4endl;
3149 }
3151 return;
3152 }
3153
3155}
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:34
CLHEP::Hep3Vector G4ThreeVector
HepGeom::RotateZ3D G4RotateZ3D
HepGeom::Translate3D G4Translate3D
HepGeom::Transform3D G4Transform3D
HepGeom::RotateX3D G4RotateX3D
HepGeom::RotateY3D G4RotateY3D
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
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:87
Definition: G4Box.hh:56
static G4Colour Brown()
Definition: G4Colour.hh:158
G4int GetEventID() const
Definition: G4Event.hh:118
static G4LogicalVolumeStore * GetInstance()
void SetClippingSolid(G4VSolid *pClippingSolid)
G4bool Validate(G4bool warn)
void SetClippingMode(ClippingMode mode)
void DescribeYourselfTo(G4VGraphicsScene &)
const std::vector< Findings > & GetFindings() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
static G4RunManager * GetMasterRunManager()
const G4Run * GetCurrentRun() const
Definition: G4Run.hh:49
G4int GetRunID() const
Definition: G4Run.hh:82
G4int GetNumberOfEventToBeProcessed() const
Definition: G4Run.hh:87
const std::vector< const G4Event * > * GetEventVector() const
Definition: G4Run.hh:119
Direction
Definition: G4Scale.hh:41
static const G4String & GetGuidanceString()
Definition: G4Scale.cc:65
G4bool AddRunDurationModel(G4VModel *, G4bool warn=false)
Definition: G4Scene.cc:50
@ UNLIMITED
Definition: G4Scene.hh:54
const G4VisExtent & GetExtent() const
const G4String & GetName() const
const std::vector< Model > & GetEndOfEventModelList() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
Definition: G4Text.hh:72
void SetLayout(Layout)
void SetOffset(double dx, double dy)
Layout
Definition: G4Text.hh:76
@ centre
Definition: G4Text.hh:76
@ right
Definition: G4Text.hh:76
@ left
Definition: G4Text.hh:76
const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
static G4TransportationManager * GetTransportationManager()
std::vector< G4VPhysicalVolume * >::iterator GetWorldsIterator()
size_t GetNoWorlds() const
Definition: G4Tubs.hh:75
G4UIcommand * FindPath(const char *commandPath) const
static G4double ValueOf(const char *unitName)
Definition: G4UIcommand.cc:348
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:146
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
static G4bool ConvertToBool(const char *st)
Definition: G4UIcommand.cc:530
G4UIcommandTree * GetTree() const
Definition: G4UImanager.hh:179
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:485
G4int GetVerboseLevel() const
Definition: G4UImanager.hh:193
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
void SetVerboseLevel(G4int val)
Definition: G4UImanager.hh:192
G4double StoD(G4String s)
void SetDefaultValue(const char *theDefaultValue)
void SetGuidance(const char *theGuidance)
void SetScreenSize(G4double)
void SetType(const G4String &)
void SetGlobalDescription(const G4String &)
void SetModelingParameters(const G4ModelingParameters *)
void SetGlobalTag(const G4String &)
void SetExtent(const G4VisExtent &)
const G4VisExtent & GetExtent() const
const G4String & GetGlobalDescription() const
void SetTransformation(const G4Transform3D &)
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)
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)
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)
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)
G4String GetCurrentValue(G4UIcommand *command)
void SetNewValue(G4UIcommand *command, G4String newValue)
G4double GetYmin() const
Definition: G4VisExtent.hh:101
G4double GetXmax() const
Definition: G4VisExtent.hh:100
G4double GetExtentRadius() const
Definition: G4VisExtent.cc:75
G4double GetYmax() const
Definition: G4VisExtent.hh:102
G4double GetZmax() const
Definition: G4VisExtent.hh:104
G4double GetZmin() const
Definition: G4VisExtent.hh:103
G4double GetXmin() const
Definition: G4VisExtent.hh:99
G4Scene * GetCurrentScene() const
const std::vector< UserVisAction > & GetRunDurationUserVisActions() const
const std::map< G4VUserVisAction *, G4VisExtent > & GetUserVisActionExtents() const
const std::vector< UserVisAction > & GetEndOfRunUserVisActions() const
G4VViewer * GetCurrentViewer() const
const std::vector< UserVisAction > & GetEndOfEventUserVisActions() const
static Verbosity GetVerbosity()
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:79
BasicVector3D< T > unit() const
CLHEP::HepRotation getRotation() const
CLHEP::Hep3Vector getTranslation() const