Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VisCommandSceneAddVolume Class Reference

#include <G4VisCommandsSceneAdd.hh>

+ Inheritance diagram for G4VisCommandSceneAddVolume:

Public Member Functions

 G4VisCommandSceneAddVolume ()
 
virtual ~G4VisCommandSceneAddVolume ()
 
G4String GetCurrentValue (G4UIcommand *command)
 
void SetNewValue (G4UIcommand *command, G4String newValue)
 
- Public Member Functions inherited from G4VVisCommandScene
 G4VVisCommandScene ()
 
virtual ~G4VVisCommandScene ()
 
- Public Member Functions inherited from G4VVisCommand
 G4VVisCommand ()
 
virtual ~G4VVisCommand ()
 
- Public Member Functions inherited from G4UImessenger
 G4UImessenger ()=default
 
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
virtual ~G4UImessenger ()
 
G4bool CommandsShouldBeInMaster () const
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VVisCommand
static G4VisManagerGetVisManager ()
 
static void SetVisManager (G4VisManager *pVisManager)
 
static const G4ColourGetCurrentTextColour ()
 
- Protected Member Functions inherited from G4VVisCommandScene
G4String CurrentSceneName ()
 
- Protected Member Functions inherited from G4VVisCommand
void SetViewParameters (G4VViewer *viewer, const G4ViewParameters &viewParams)
 
void RefreshIfRequired (G4VViewer *viewer)
 
void InterpolateViews (G4VViewer *currentViewer, const std::vector< G4ViewParameters > &viewVector, const G4int nInterpolationPoints=50, const G4int waitTimePerPointmilliseconds=20, const G4String &exportString="")
 
void InterpolateToNewView (G4VViewer *currentViewer, const G4ViewParameters &oldVP, const G4ViewParameters &newVP, const G4int nInterpolationPoints=50, const G4int waitTimePerPointmilliseconds=20, const G4String &exportString="")
 
void Twinkle (G4VViewer *currentViewer, const G4ViewParameters &baseVP, const std::vector< std::vector< G4PhysicalVolumeModel::G4PhysicalVolumeNodeID > > &paths)
 
const G4StringConvertToColourGuidance ()
 
void ConvertToColour (G4Colour &colour, const G4String &redOrString, G4double green, G4double blue, G4double opacity)
 
G4bool ProvideValueOfUnit (const G4String &where, const G4String &unit, const G4String &category, G4double &value)
 
void CopyCameraParameters (G4ViewParameters &target, const G4ViewParameters &from)
 
void CheckSceneAndNotifyHandlers (G4Scene *=nullptr)
 
G4bool CheckView ()
 
void G4VisCommandsSceneAddUnsuccessful (G4VisManager::Verbosity verbosity)
 
void CopyGuidanceFrom (const G4UIcommand *fromCmd, G4UIcommand *toCmd, G4int startLine=0)
 
void CopyParametersFrom (const G4UIcommand *fromCmd, G4UIcommand *toCmd)
 
void DrawExtent (const G4VisExtent &)
 
- Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
 
G4String LtoS (G4long l)
 
G4String DtoS (G4double a)
 
G4String BtoS (G4bool b)
 
G4int StoI (const G4String &s)
 
G4long StoL (const G4String &s)
 
G4double StoD (const G4String &s)
 
G4bool StoB (const G4String &s)
 
void AddUIcommand (G4UIcommand *newCommand)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
template<typename T>
T * CreateCommand (const G4String &cname, const G4String &dsc)
 
- Static Protected Member Functions inherited from G4VVisCommand
static G4String ConvertToString (G4double x, G4double y, const char *unitName)
 
static G4bool ConvertToDoublePair (const G4String &paramString, G4double &xval, G4double &yval)
 
- Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 
- Static Protected Attributes inherited from G4VVisCommand
static G4VisManagerfpVisManager = nullptr
 
static G4int fCurrentArrow3DLineSegmentsPerCircle = 6
 
static G4Colour fCurrentColour = G4Colour::White()
 
static G4double fCurrentLineWidth = 1.
 
static G4Colour fCurrentTextColour = G4Colour::Blue()
 
static G4Text::Layout fCurrentTextLayout = G4Text::left
 
static G4double fCurrentTextSize = 12.
 
static G4PhysicalVolumeModel::TouchableProperties fCurrentTouchableProperties
 
static G4VisExtent fCurrentExtentForField
 
static std::vector< G4PhysicalVolumesSearchScene::FindingsfCurrrentPVFindingsForField
 
static G4bool fThereWasAViewer = false
 
static G4ViewParameters fExistingVP
 
static G4SceneTreeItem fExistingSceneTree
 

Detailed Description

Definition at line 479 of file G4VisCommandsSceneAdd.hh.

Constructor & Destructor Documentation

◆ G4VisCommandSceneAddVolume()

G4VisCommandSceneAddVolume::G4VisCommandSceneAddVolume ( )

Definition at line 3042 of file G4VisCommandsSceneAdd.cc.

3042 {
3043 G4bool omitable;
3044 fpCommand = new G4UIcommand ("/vis/scene/add/volume", this);
3045 fpCommand -> SetGuidance
3046 ("Adds a physical volume to current scene, with optional clipping volume.");
3047 fpCommand -> SetGuidance
3048 ("If physical-volume-name is \"world\" (the default), the top of the"
3049 "\nmain geometry tree (material world) is added. If \"worlds\", the"
3050 "\ntops of all worlds - material world and parallel worlds, if any - are"
3051 "\nadded. Otherwise a search of all worlds is made.");
3052 fpCommand -> SetGuidance
3053 ("In the last case the names of all volumes in all worlds are matched"
3054 "\nagainst physical-volume-name. If this is of the form \"/regexp/\","
3055 "\nwhere regexp is a regular expression (see C++ regex), the match uses"
3056 "\nthe usual rules of regular expression matching. Otherwise an exact"
3057 "\nmatch is required."
3058 "\nFor example, \"/Shap/\" adds \"Shape1\" and \"Shape2\".");
3059 fpCommand -> SetGuidance
3060 ("It may help to see a textual representation of the geometry hierarchy of"
3061 "\nthe worlds. Try \"/vis/drawTree [worlds]\" or one of the driver/browser"
3062 "\ncombinations that have the required functionality, e.g., HepRepFile.");
3063 fpCommand -> SetGuidance
3064 ("If clip-volume-type is specified, the subsequent parameters are used to"
3065 "\nto define a clipping volume. For example,"
3066 "\n\"/vis/scene/add/volume ! ! ! -box km 0 1 0 1 0 1\" will draw the world"
3067 "\nwith the positive octant cut away. (If the Boolean Processor issues"
3068 "\nwarnings try replacing 0 by 0.000000001 or something.)");
3069 fpCommand -> SetGuidance
3070 ("If clip-volume-type is prepended with '-', the clip-volume is subtracted"
3071 "\n(cutaway). (This is the default if there is no prepended character.)"
3072 "\nIf '*' is prepended, the intersection of the physical-volume and the"
3073 "\nclip-volume is made. (You can make a section through the detector with"
3074 "\na thin box, for example).");
3075 fpCommand -> SetGuidance
3076 ("For \"box\", the parameters are xmin,xmax,ymin,ymax,zmin,zmax."
3077 "\nOnly \"box\" is programmed at present.");
3078 G4UIparameter* parameter;
3079 parameter = new G4UIparameter ("physical-volume-name", 's', omitable = true);
3080 parameter -> SetDefaultValue ("world");
3081 fpCommand -> SetParameter (parameter);
3082 parameter = new G4UIparameter ("copy-no", 'i', omitable = true);
3083 parameter -> SetGuidance ("If negative, matches any copy no.");
3084 parameter -> SetDefaultValue (-1);
3085 fpCommand -> SetParameter (parameter);
3086 parameter = new G4UIparameter ("depth-of-descent", 'i', omitable = true);
3087 parameter -> SetGuidance
3088 ("Depth of descent of geometry hierarchy. Default = unlimited depth.");
3089 parameter -> SetDefaultValue (G4PhysicalVolumeModel::UNLIMITED);
3090 fpCommand -> SetParameter (parameter);
3091 parameter = new G4UIparameter ("clip-volume-type", 's', omitable = true);
3092 parameter -> SetParameterCandidates("none box -box *box");
3093 parameter -> SetDefaultValue ("none");
3094 parameter -> SetGuidance("[-|*]type. See general guidance.");
3095 fpCommand -> SetParameter (parameter);
3096 parameter = new G4UIparameter ("parameter-unit", 's', omitable = true);
3097 parameter -> SetDefaultValue ("m");
3098 fpCommand -> SetParameter (parameter);
3099 parameter = new G4UIparameter ("parameter-1", 'd', omitable = true);
3100 parameter -> SetDefaultValue (0.);
3101 fpCommand -> SetParameter (parameter);
3102 parameter = new G4UIparameter ("parameter-2", 'd', omitable = true);
3103 parameter -> SetDefaultValue (0.);
3104 fpCommand -> SetParameter (parameter);
3105 parameter = new G4UIparameter ("parameter-3", 'd', omitable = true);
3106 parameter -> SetDefaultValue (0.);
3107 fpCommand -> SetParameter (parameter);
3108 parameter = new G4UIparameter ("parameter-4", 'd', omitable = true);
3109 parameter -> SetDefaultValue (0.);
3110 fpCommand -> SetParameter (parameter);
3111 parameter = new G4UIparameter ("parameter-5", 'd', omitable = true);
3112 parameter -> SetDefaultValue (0.);
3113 fpCommand -> SetParameter (parameter);
3114 parameter = new G4UIparameter ("parameter-6", 'd', omitable = true);
3115 parameter -> SetDefaultValue (0.);
3116 fpCommand -> SetParameter (parameter);
3117}
bool G4bool
Definition G4Types.hh:86

◆ ~G4VisCommandSceneAddVolume()

G4VisCommandSceneAddVolume::~G4VisCommandSceneAddVolume ( )
virtual

Definition at line 3119 of file G4VisCommandsSceneAdd.cc.

3119 {
3120 delete fpCommand;
3121}

Member Function Documentation

◆ GetCurrentValue()

G4String G4VisCommandSceneAddVolume::GetCurrentValue ( G4UIcommand * command)
virtual

Reimplemented from G4UImessenger.

Definition at line 3123 of file G4VisCommandsSceneAdd.cc.

3123 {
3124 return "world 0 -1";
3125}

◆ SetNewValue()

void G4VisCommandSceneAddVolume::SetNewValue ( G4UIcommand * command,
G4String newValue )
virtual

Reimplemented from G4UImessenger.

Definition at line 3127 of file G4VisCommandsSceneAdd.cc.

3128 {
3129
3130 G4VisManager::Verbosity verbosity = fpVisManager->GetVerbosity();
3131 G4bool warn = verbosity >= G4VisManager::warnings;
3132
3133 G4Scene* pScene = fpVisManager->GetCurrentScene();
3134 if (!pScene) {
3135 if (verbosity >= G4VisManager::errors) {
3136 G4warn << "ERROR: No current scene. Please create one." << G4endl;
3137 }
3138 return;
3139 }
3140
3141 G4String name, clipVolumeType, parameterUnit;
3142 G4int copyNo, requestedDepthOfDescent;
3143 G4double param1, param2, param3, param4, param5, param6;
3144 std::istringstream is (newValue);
3145 is >> name >> copyNo >> requestedDepthOfDescent
3146 >> clipVolumeType >> parameterUnit
3147 >> param1 >> param2 >> param3 >> param4 >> param5 >> param6;
3149 G4PhysicalVolumeModel::subtraction; // Default subtraction mode.
3150 if (clipVolumeType[size_t(0)] == '-') {
3151 clipVolumeType = clipVolumeType.substr(1); // Remove first character.
3152 } else if (clipVolumeType[size_t(0)] == '*') {
3154 clipVolumeType = clipVolumeType.substr(1);
3155 }
3156 G4double unit = G4UIcommand::ValueOf(parameterUnit);
3157 param1 *= unit; param2 *= unit; param3 *= unit;
3158 param4 *= unit; param5 *= unit; param6 *= unit;
3159
3160 G4VSolid* clippingSolid = nullptr;
3161 if (clipVolumeType == "box") {
3162 const G4double dX = (param2 - param1) / 2.;
3163 const G4double dY = (param4 - param3) / 2.;
3164 const G4double dZ = (param6 - param5) / 2.;
3165 const G4double x0 = (param2 + param1) / 2.;
3166 const G4double y0 = (param4 + param3) / 2.;
3167 const G4double z0 = (param6 + param5) / 2.;
3168 clippingSolid = new G4DisplacedSolid
3169 ("_displaced_clipping_box",
3170 new G4Box("_clipping_box",dX,dY,dZ),
3171 G4Translate3D(x0,y0,z0));
3172 }
3173
3174 G4TransportationManager* transportationManager =
3176
3177 size_t nWorlds = transportationManager->GetNoWorlds();
3178 if (nWorlds > 1) { // Parallel worlds in operation...
3179 if (verbosity >= G4VisManager::warnings) {
3180 static G4bool warned = false;
3181 if (!warned && name != "worlds") {
3182 G4warn <<
3183 "WARNING: Parallel worlds in operation. To visualise, specify"
3184 "\n \"worlds\" or the parallel world volume or sub-volume name"
3185 "\n and control visibility with /vis/geometry."
3186 << G4endl;
3187 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
3188 transportationManager->GetWorldsIterator();
3189 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
3190 G4warn << " World " << i << ": " << (*iterWorld)->GetName()
3191 << G4endl;
3192 warned = true;
3193 }
3194 }
3195 }
3196 }
3197
3198 // Get the world (the initial value of the iterator points to the mass world).
3199 G4VPhysicalVolume* world = *(transportationManager->GetWorldsIterator());
3200
3201 if (!world) {
3202 if (verbosity >= G4VisManager::errors) {
3203 G4warn <<
3204 "ERROR: G4VisCommandSceneAddVolume::SetNewValue:"
3205 "\n No world. Maybe the geometry has not yet been defined."
3206 "\n Try \"/run/initialize\""
3207 << G4endl;
3208 }
3209 return;
3210 }
3211
3212 std::vector<G4PhysicalVolumesSearchScene::Findings> findingsVector;
3213
3214 // When it comes to determining the extent of a physical volume we normally
3215 // assume the user wishes to ignore "invisible" volumes. For example, most
3216 // users make the world volume invisible. So we ask the physical volume
3217 // model to traverse the geometry hierarchy, starting at the named physical
3218 // volume, until it finds non-invisible ones, whose extents are accumulated
3219 // to determine the overall extent. (Once a non-invisible volume is found,
3220 // the search is curtailed - daughters are always contained within the mother
3221 // so they have no subsequent influence on the extent of the mother - but the
3222 // search continues at the same level until all highest level non-invisible
3223 // volumes are found an their extents accumulated.) So the default is
3224 G4bool useFullExtent = false;
3225 // However, the above procedure can be time consuming in some situations, such
3226 // as a nested parameterisation whose ultimate volumes are the first non-
3227 // visible ones, which are typical of a medical "phantom". So we assume here
3228 // below that if a user specifies a name other than "world" or "worlds" he/she
3229 // wished the extent to be determined by the volume, whether it is visible
3230 // or not. So we set useFullExtent true at that point below.
3231
3232 if (name == "world") {
3233
3234 findingsVector.push_back
3235 (G4PhysicalVolumesSearchScene::Findings(world,world));
3236
3237 } else if (name == "worlds") {
3238
3239 if (nWorlds <= 1) {
3240 if (verbosity >= G4VisManager::warnings) {
3241 G4warn <<
3242 "WARNING: G4VisCommandSceneAddVolume::SetNewValue:"
3243 "\n Parallel worlds requested but none exist."
3244 "\n Just adding material world."
3245 << G4endl;
3246 }
3247 }
3248 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
3249 transportationManager->GetWorldsIterator();
3250 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
3251 findingsVector.push_back
3252 (G4PhysicalVolumesSearchScene::Findings
3253 (*iterWorld,*iterWorld));
3254 }
3255
3256 } else { // Search all worlds...
3257
3258 // Use the model's full extent. This assumes the user wants these
3259 // volumes in the findings vector (there could be more than one) to
3260 // determine the scene's extent. Otherwise G4PhysicalVolumeModel would
3261 // re-calculate each volume's extent based on visibility, etc., which
3262 // could be time consuming.
3263 useFullExtent = true;
3264
3265 std::vector<G4VPhysicalVolume*>::iterator iterWorld =
3266 transportationManager->GetWorldsIterator();
3267 for (size_t i = 0; i < nWorlds; ++i, ++iterWorld) {
3268 G4ModelingParameters mp; // Default - no culling.
3269 G4PhysicalVolumeModel searchModel
3270 (*iterWorld,
3271 requestedDepthOfDescent,
3272 G4Transform3D(),
3273 &mp,
3274 useFullExtent);
3275 G4PhysicalVolumesSearchScene searchScene(&searchModel, name, copyNo);
3276 searchModel.DescribeYourselfTo (searchScene); // Initiate search.
3277 for (const auto& findings: searchScene.GetFindings()) {
3278 findingsVector.push_back(findings);
3279 }
3280 }
3281 }
3282
3283 for (const auto& findings: findingsVector) {
3284 // Set copy number from search findings for replicas and parameterisations.
3285 findings.fpFoundPV->SetCopyNo(findings.fFoundPVCopyNo);
3286 G4PhysicalVolumeModel* foundPVModel = new G4PhysicalVolumeModel
3287 (findings.fpFoundPV,
3288 requestedDepthOfDescent,
3289 findings.fFoundObjectTransformation,
3290 0, // No modelling parameters (these are set later by the scene handler).
3291 useFullExtent,
3292 findings.fFoundBasePVPath);
3293 if (clippingSolid) {
3294 foundPVModel->SetClippingSolid(clippingSolid);
3295 foundPVModel->SetClippingMode(clippingMode);
3296 }
3297 if (!foundPVModel->Validate(warn)) return;
3298 // ...so add it to the scene.
3299 G4bool successful = pScene->AddRunDurationModel(foundPVModel,warn);
3300 if (successful) {
3301 if (verbosity >= G4VisManager::confirmations) {
3302 G4cout << "\"" << findings.fpFoundPV->GetName()
3303 << "\", copy no. " << findings.fFoundPVCopyNo
3304 << ",\n found in searched volume \""
3305 << findings.fpSearchPV->GetName()
3306 << "\" at depth " << findings.fFoundDepth
3307 << ",\n base path: \"" << findings.fFoundBasePVPath
3308 << "\",\n with a requested depth of further descent of ";
3309 if (requestedDepthOfDescent < 0) {
3310 G4cout << "<0 (unlimited)";
3311 }
3312 else {
3313 G4cout << requestedDepthOfDescent;
3314 }
3315 G4cout << ",\n has been added to scene \"" << pScene->GetName() << "\"."
3316 << G4endl;
3317 }
3318 } else {
3320 }
3321 }
3322
3323 if (findingsVector.empty()) {
3324 if (verbosity >= G4VisManager::errors) {
3325 G4warn << "ERROR: Volume \"" << name << "\"";
3326 if (copyNo >= 0) {
3327 G4warn << ", copy no. " << copyNo << ",";
3328 }
3329 G4warn << " not found." << G4endl;
3330 }
3332 return;
3333 }
3334
3336}
#define G4warn
Definition G4Scene.cc:41
HepGeom::Transform3D G4Transform3D
HepGeom::Translate3D G4Translate3D
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetClippingSolid(G4VSolid *pClippingSolid)
void SetClippingMode(ClippingMode mode)
G4bool AddRunDurationModel(G4VModel *, G4bool warn=false)
Definition G4Scene.cc:160
const G4String & GetName() const
static G4TransportationManager * GetTransportationManager()
std::vector< G4VPhysicalVolume * >::iterator GetWorldsIterator()
std::size_t GetNoWorlds() const
static G4double ValueOf(const char *unitName)
void G4VisCommandsSceneAddUnsuccessful(G4VisManager::Verbosity verbosity)
void CheckSceneAndNotifyHandlers(G4Scene *=nullptr)
static G4VisManager * fpVisManager
const char * name(G4int ptype)

The documentation for this class was generated from the following files: