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