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