BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSMessenger.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#include "BDSAcceleratorModel.hh"
20#include "BDSBeamline.hh"
21#include "BDSBeamlineElement.hh"
22#include "BDSMessenger.hh"
23#include "BDSParser.hh"
24#include "BDSSamplerRegistry.hh"
25#include "BDSUIcmdStringInt.hh"
26#include "BDSUtilities.hh"
27
28#include "parser/query.h"
29
30#include "globals.hh"
31#include "G4UImanager.hh"
32#include "G4UIdirectory.hh"
33#include "G4UIcmdWithoutParameter.hh"
34#include "G4UIcmdWithAString.hh"
35
36#include <iomanip>
37#include <iostream>
38#include <ostream>
39#include <string>
40
41
42BDSMessenger::BDSMessenger()
43{
44 bdsDirectory = new G4UIdirectory("/bds/");
45 bdsDirectory->SetGuidance("BDSIM commands");
46
47 bdsBeamlineDirectory = new G4UIdirectory("/bds/beamline/");
48
49 beamlineListCmd = new G4UIcmdWithoutParameter("/bds/beamline/list",this);
50 beamlineListCmd->SetGuidance("List beamline components");
51
52 elementNameSearchCmd = new G4UIcmdWithAString("/bds/beamline/namesearch",this);
53 elementNameSearchCmd->SetGuidance("Search beamline components for element");
54
55 elementGoToCmd = new BDSUIcmdStringInt("/bds/beamline/goto", this);
56 G4String elementGoToCmdGuidance = "Move the viewpoint to a particular element's location. ";
57 elementGoToCmdGuidance += "Arguments are the name of an element in the beam line and optionally,";
58 elementGoToCmdGuidance += " the instance number (zero-counting).";
59 elementGoToCmd->SetGuidance(elementGoToCmdGuidance);
60
61 bdsSamplersDirectory = new G4UIdirectory("/bds/samplers/");
62 samplerListCmd = new G4UIcmdWithoutParameter("/bds/samplers/list",this);
63 samplerListCmd->SetGuidance("List samplers");
64
65 samplerViewCmd = new G4UIcmdWithoutParameter("/bds/samplers/view",this);
66 samplerViewCmd->SetGuidance("View sampler paralle world");
67
68 queryListCmd = new G4UIcmdWithoutParameter("/bds/field/listQueries", this);
69 queryListCmd->SetGuidance("List all queries defined in input");
70
71 // G4UImanager* UIManager = G4UImanager::GetUIpointer();
72 // UIManager->ApplyCommand("/control/execute " + visMacroFilename);
73}
74
75BDSMessenger::~BDSMessenger()
76{
77 delete bdsDirectory;
78 delete bdsBeamlineDirectory;
79 delete beamlineListCmd;
80 delete elementNameSearchCmd;
81 delete elementGoToCmd;
82 delete bdsSamplersDirectory;
83 delete samplerListCmd;
84 delete samplerViewCmd;
85}
86
87void BDSMessenger::SetNewValue(G4UIcommand* command, G4String newValue)
88{
89 if(command == beamlineListCmd)
90 {BeamLineList();}
91 else if(command == elementNameSearchCmd)
92 {ElementNameSearch(newValue);}
93 else if(command == elementGoToCmd)
94 {GoToElement(newValue);}
95 else if(command == samplerListCmd)
96 {SamplerList();}
97 else if (command == samplerViewCmd)
98 {ViewSamplers();}
99 else if (command == queryListCmd)
100 {ListQueries();}
101}
102
103void BDSMessenger::BeamLineList()
104{
105 const BDSBeamline *beamline = BDSAcceleratorModel::Instance()->BeamlineMain();
106
107 int j = 0;
108 auto flagsCache(G4cout.flags());
109 G4cout << std::right
110 << std::setw(4) << "index" << std::setfill(' ')
111 << std::setw(25) << "name"
112 << std::setw(25) << "placement name"
113 << std::setw(20) << "type"
114 << std::setw(20) << "S-middle[m]" << G4endl;
115 for (auto i = beamline->begin(); i != beamline->end(); ++i, ++j)
116 {G4cout << BDSBeamlineElementToString(j) << G4endl;}
117 G4cout.flags(flagsCache);
118}
119
120std::string BDSMessenger::BDSBeamlineElementToString(G4int iElement)
121{
122 std::stringstream ss;
123
124 const BDSBeamline* beamline = BDSAcceleratorModel::Instance()->BeamlineMain();
125 const BDSBeamlineElement* e = beamline->at(iElement);
126
127 ss << std::setfill('0') << std::right << std::setw(4) << iElement << " " << std::setfill(' ')
128 << std::setw(25) << e->GetName()
129 << std::setw(25) << e->GetPlacementName()
130 << std::setw(20) << e->GetType()
131 << std::setw(20) << std::setprecision(4) << std::fixed << e->GetSPositionMiddle()/CLHEP::m;
132
133 return ss.str();
134}
135
136void BDSMessenger::ElementNameSearch(std::string name)
137{
138 const BDSBeamline* beamline = BDSAcceleratorModel::Instance()->BeamlineMain();
139 int j=0;
140 for (auto i = beamline->begin(); i != beamline->end(); ++i, ++j)
141 {
142 if(BDS::StrContains((*i)->GetName(), name))
143 {G4cout << (*i)->GetName() << G4endl;}
144 }
145}
146
147void BDSMessenger::GoToElement(const std::string& value)
148{
149 if (value.empty())
150 {G4cerr << "empty string given - can't search" << G4endl; return;}
151 std::vector<G4String> words = BDS::SplitOnWhiteSpace(value);
152 G4String name = words[0];
153 G4int instance = 0;
154 if (words.size() > 1)
155 {instance = std::stoi(words[1]);}
156
157 const BDSBeamline* beamline = BDSAcceleratorModel::Instance()->BeamlineMain();
158
159 if (!beamline)
160 {G4cout << "No beam line in this model so not possible to search it."; return;}
161
162 // search for the name exactly
163 const BDSBeamlineElement* e = beamline->GetElement(name, instance);
164
165 G4int count = -1;
166 if (!e)
167 {// search the beam line for any element containing the name at all
168 for (const auto& el : *beamline)
169 {
170 if (BDS::StrContains(el->GetName(), name))
171 {
172 count++;
173 if (count == instance)
174 {
175 e = el;
176 break;
177 }
178 }
179 }
180 if (!e)
181 {
182 G4cout << "No component found by that name" << G4endl;
183 if (count > -1)
184 {G4cout << "only " << count << " instances found." << G4endl;}
185 return;
186 }
187 }
188
189 G4ThreeVector pos = e->GetReferencePositionMiddle();
190 G4cout << "goto> " << name << " at global position> " << pos/CLHEP::m << " m" << G4endl;
191 G4UImanager* UIManager = G4UImanager::GetUIpointer();
192 G4String posStr = std::to_string(pos.x()) + " " + std::to_string(pos.y()) + " " + std::to_string(pos.z());
193 UIManager->ApplyCommand("/vis/viewer/set/targetPoint " + posStr + " mm");
194}
195
196void BDSMessenger::ElementTypeSearch(std::string /*type*/)
197{;}
198
199void BDSMessenger::SamplerList()
200{
201 for (const auto& name : BDSSamplerRegistry::Instance()->GetUniqueNames())
202 {G4cout << name << G4endl;}
203}
204
205void BDSMessenger::ViewSamplers()
206{
207 // we can't use drawVolume because that clears the scene and only draws one volume
208 // inside it really applied add/volume as we do here
209 // some how viewing worlds undoes these trajectory actions in our default vis.mac
210 G4UImanager* UIManager = G4UImanager::GetUIpointer();
211 UIManager->ApplyCommand("/vis/scene/add/volume SamplerWorld_main");
212 UIManager->ApplyCommand("/vis/sceneHandler/attach");
213 UIManager->ApplyCommand("/vis/scene/endOfEventAction accumulate");
214 UIManager->ApplyCommand("/vis/scene/add/trajectories");
215 UIManager->ApplyCommand("/vis/ogl/flushAt endOfEvent");
216}
217
218std::string BDSMessenger::BDSSamplerToString(std::string /*samplerName*/)
219{
220 return std::string();
221}
222
223std::string BDSMessenger::BDSSamplerToString(int iSampler)
224{
225 std::stringstream ss;
226
228
229 ss << std::setfill('0') << std::right << std::setw(4) << iSampler << " " << std::setfill(' ')
230 << std::right << std::setw(20) << sInfo.Name() << " "
231 << std::right << std::setw(20) << sInfo.SPosition();
232 return ss.str();
233}
234
235void BDSMessenger::ListQueries()
236{
237 auto queries = BDSParser::Instance()->GetQuery();
238 for (const auto& qu : queries)
239 {G4cout << qu.name << G4endl;}
240}
const BDSBeamline * BeamlineMain() const
Accessor.
A class that holds a fully constructed BDSAcceleratorComponent as well as any information relevant to...
G4String GetPlacementName() const
Accessor.
G4String GetName() const
Accessor.
G4ThreeVector GetReferencePositionMiddle() const
Accessor.
G4double GetSPositionMiddle() const
Accessor.
G4String GetType() const
Accessor.
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
const BDSBeamlineElement * GetElement(G4String acceleratorComponentName, G4int i=0) const
Get the ith placement of an element in the beam line. Returns null pointer if not found.
Definition: BDSBeamline.cc:739
const BDSBeamlineElement * at(int iElement) const
Return a reference to the element at i.
Definition: BDSBeamline.hh:120
iterator begin()
Iterator mechanics.
Definition: BDSBeamline.hh:167
iterator end()
Iterator mechanics.
Definition: BDSBeamline.hh:168
static BDSParser * Instance()
Access method.
Definition: BDSParser.cc:28
std::vector< GMAD::Query > GetQuery() const
Return the parser list of that object.
Definition: BDSParser.hh:104
Information about a registered sampler.
G4String Name() const
Accessor.
G4double SPosition() const
Accessor.
static BDSSamplerRegistry * Instance()
Accessor for registry.
const BDSSamplerPlacementRecord & GetInfo(G4int index) const
Accessor.
UI command with string, then optional integer argument.
G4bool StrContains(const G4String &str, const G4String &test)
Utility function to simplify lots of syntax changes for pedantic g4 changes.
Definition: BDSUtilities.cc:66
std::vector< G4String > SplitOnWhiteSpace(const G4String &input)
Split a string on whitespace and return a vector of these 'words'.