BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSCurvilinearBuilder.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
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 "BDSAcceleratorComponent.hh"
20#include "BDSAcceleratorComponentRegistry.hh"
21#include "BDSBeamline.hh"
22#include "BDSBeamlineElement.hh"
23#include "BDSExtent.hh"
24#include "BDSCurvilinearBuilder.hh"
25#include "BDSCurvilinearFactory.hh"
26#include "BDSGlobalConstants.hh"
27#include "BDSSimpleComponent.hh"
28#include "BDSSamplerType.hh"
29#include "BDSTiltOffset.hh"
30#include "BDSTunnelInfo.hh"
31#include "BDSUtilities.hh"
32#include "BDSWarning.hh"
33
34#include "globals.hh" // geant4 types / globals
35#include "G4ThreeVector.hh"
36#include "G4RotationMatrix.hh"
37
38#include <algorithm>
39#include <cmath>
40#include <iterator>
41
42BDSCurvilinearBuilder::BDSCurvilinearBuilder()
43{
44 const BDSGlobalConstants* globals = BDSGlobalConstants::Instance(); // shortcut
45
47 defaultBridgeLength = paddingLength + 4*globals->LengthSafety();
48 curvilinearRadius = globals->CurvilinearDiameter()*0.5;
49 radiusTolerance = 0.8;
50 if (globals->BuildTunnel() || globals->BuildTunnelStraight())
51 {// query the default tunnel model
52 BDSExtent tunnelExtent = globals->TunnelInfo()->IndicativeExtent();
53 tunnelExtent = tunnelExtent.Translate(globals->TunnelOffsetX(), globals->TunnelOffsetY(), 0);
54 G4double maxTunnelR = tunnelExtent.MaximumAbs();
55 if (curvilinearRadius < maxTunnelR)
56 {
57 if (globals->CurvilinearDiameterShrunkForBends())
58 {
59 G4String message = "Strong bends necessitate smaller curvilinear geometry than the tunnel size, but\n";
60 message += "tunnel hits will need transform for curvilinear coordinates.";
61 BDS::Warning(message);
62 }
63 else
64 {curvilinearRadius = std::max(curvilinearRadius, maxTunnelR);}
65 }
66 }
67 bonusChordLength = 1*CLHEP::m;
68
70}
71
72BDSCurvilinearBuilder::~BDSCurvilinearBuilder()
73{
74 delete factory;
75}
76
78 const G4bool circular)
79{
80 BDSBeamline* result = new BDSBeamline();
81
82 //prepend small sections to machine when non-circular
83 G4bool bonusSections = (!circular && !beamline->empty());
84
85 if (bonusSections)
86 {//prepend small section to machine
87 BDSBeamlineElement* bonusBit = CreateBonusSectionStart(beamline);
88 result->AddBeamlineElement(bonusBit);
89 }
90
91 G4int i = 0;
92 for (BDSBeamline::const_iterator element = beamline->begin(); element != beamline->end(); ++element)
93 {
94 G4String name = (*element)->GetName() + "_cl_" + std::to_string(i);
95 i++;
96 const BDSBeamlineElement* pstEl = nullptr;
97 const BDSBeamlineElement* nxtEl = nullptr;
98 PreviousAndNext(element, beamline->begin(), beamline->end(), pstEl, nxtEl);
99 G4double crRadius = std::min({CurvilinearRadius(pstEl),
100 CurvilinearRadius(*element),
101 CurvilinearRadius(nxtEl)});
102 BDSBeamlineElement* temp = CreateCurvilinearElement(name, element, element, i, crRadius);
103 if (temp)
104 {result->AddBeamlineElement(temp);}
105 }
106
107 if (bonusSections)
108 {// append small section to machine
109 BDSBeamlineElement* bonusBit = CreateBonusSectionEnd(beamline);
110 if (bonusBit)
111 {result->AddBeamlineElement(bonusBit);}
112 }
113 return result;
114}
115
119 const BDSBeamlineElement*& previous,
120 const BDSBeamlineElement*& next) const
121{
122 G4bool isFirst = it == startIt;
123 G4bool isLast = std::next(it) == endIt; // end is one beyond the last
124
125 G4double minimumSize = 0.1*CLHEP::mm;
126 auto pst = std::prev(it);
127 if (isFirst)
128 {previous = nullptr;}
129 else
130 {
131 previous = *pst;
132 if (previous->GetArcLength() < minimumSize && pst != startIt)
133 {
134 pst = std::prev(pst);
135 if (pst != startIt)
136 {previous = *pst;}
137 }
138 }
139 auto nxt = std::next(it);
140 if (isLast)
141 {next = nullptr;}
142 else
143 {
144 next = *nxt;
145 if (next->GetArcLength() < minimumSize && nxt != endIt)
146 {
147 nxt = std::next(nxt);
148 if (nxt != endIt)
149 {next = *nxt;}
150 }
151 }
152}
153
155{
156 BDSBeamline* result = new BDSBeamline();
157
159
160 G4int numberOfUniqueComponents = 0;
161 G4int beamlineIndex = 0;
162 for (BDSBeamline::const_iterator element = beamline->begin(); element != beamline->end(); ++element)
163 {
164 auto nxt = std::next(element);
165 const BDSBeamlineElement* pstEl = nullptr;
166 const BDSBeamlineElement* nxtEl = nullptr;
167 PreviousAndNext(element, beamline->begin(), beamline->end(), pstEl, nxtEl);
168 G4double crRadius = std::min({CurvilinearRadius(pstEl),
169 CurvilinearRadius(*element),
170 CurvilinearRadius(nxtEl)});
171 BDSBeamlineElement* bridgeSection = CreateBridgeSection(defaultBridge, element, nxt,
172 beamline->end(), numberOfUniqueComponents,
173 beamlineIndex, crRadius);
174 if (bridgeSection)
175 {
176 result->AddBeamlineElement(bridgeSection);
177 beamlineIndex++;
178 }
179 }
180 return result;
181}
182
184{
185 if (!el)
186 {return curvilinearRadius;} // default in case there's no element
187 if (BDS::IsFinite(el->GetAngle()))
188 {
189 G4double radiusFromAngleLength = std::abs(el->GetArcLength()/el->GetAngle()); // s = r*theta -> r = s/theta
190 radiusFromAngleLength *= radiusTolerance;
191 return std::min(radiusFromAngleLength, curvilinearRadius);
192 }
193 else // no finite bending angle
194 {return curvilinearRadius;}
195}
196
198 BDSBeamline::const_iterator startElement,
199 BDSBeamline::const_iterator finishElement,
200 G4int index,
201 G4double crRadius)
202{
203 BDSSimpleComponent* component = nullptr;
204
205 // we'll take the tilt from the first element - they should only ever be the same when used here
206 G4bool tilted = BDS::IsFinite((*startElement)->GetTilt());
207 // variables to be defined to create component
208 G4double chordLength=0.0, angle=0.0, arcLength=0.0;
209
210 if (startElement == finishElement)
211 {// build 1:1
212 chordLength = (*startElement)->GetChordLength() + paddingLength;
213 angle = (*startElement)->GetAngle();
214 if (Angled(*startElement))
215 {
216 // Not strictly accurate to add on paddingLength to arcLength, but close for now.
217 arcLength = (*startElement)->GetArcLength() + paddingLength;
218 }
219 }
220 else
221 {// cover a few components
222 G4ThreeVector positionStart = (*startElement)->GetReferencePositionStart();
223 G4ThreeVector positionEnd = (*finishElement)->GetReferencePositionEnd();
224 chordLength = (positionEnd - positionStart).mag() + paddingLength;
225
226 G4double accumulatedAngle = 0;
227 for (auto it = startElement; it < finishElement; it++)
228 {accumulatedAngle += (*it)->GetAngle();}
229
230 angle = accumulatedAngle;
231 if (BDS::IsFinite(angle))
232 {
233 G4double meanBendingRadius = 0.5 * chordLength / std::sin(0.5*std::abs(angle));
234 arcLength = meanBendingRadius * std::abs(angle);
235 }
236 }
237
238 if (!BDS::IsFinite(angle))
239 {// straight
240 component = factory->CreateCurvilinearVolume(elementName,
241 chordLength,
242 crRadius);
243 }
244 else
245 {// angled - tilt matters
246 BDSTiltOffset* to = nullptr;
247 if (tilted)
248 {to = (*startElement)->GetTiltOffset();}
249
250 component = factory->CreateCurvilinearVolume(elementName,
251 arcLength,
252 chordLength,
253 crRadius,
254 angle,
255 to);
256 }
257
259
260 return CreateElementFromComponent(component, startElement, finishElement, index);
261}
262
265 BDSBeamline::const_iterator nextElement,
267 G4int& numberOfUniqueComponents,
268 const G4int beamlineIndex,
269 G4double crRadius)
270{
271 // we can safely assume faces match between two beam line elements so if one's angled, so is the other
272 BDSAcceleratorComponent* component = defaultBridge;
273 if ((*element)->AngledOutputFace()) // angled faces - make one to match to cover the angled gap
274 {component = CreateAngledBridgeComponent(element, numberOfUniqueComponents, crRadius);}
275 else if (BDS::IsFinite((*element)->GetAngle()))
276 {// width may be reduced due to bend - check if required
277 G4double width = (*element)->GetAcceleratorComponent()->GetExtent().DX();
278 width = std::min(width, crRadius);
279 component = CreateStraightBridgeComponent(width, numberOfUniqueComponents);
280 }
281
282 return CreateBridgeElementFromComponent(component, element, nextElement, end, beamlineIndex);
283}
284
285
287{
288 // we're ignoring any possible angled face of the curvilinear geometry
289 BDSSimpleComponent* component = factory->CreateCurvilinearVolume("clb_flat_face",
290 defaultBridgeLength,
292
294
295 return component;
296}
297
299 G4int& numberOfUniqueComponents)
300{
301 BDSSimpleComponent* component = factory->CreateCurvilinearVolume("clb_" + std::to_string(numberOfUniqueComponents),
302 defaultBridgeLength,
303 width);
305 numberOfUniqueComponents++;
306 return component;
307}
308
310 G4int& numberOfUniqueComponents,
311 G4double suggestedRadius)
312{
313 G4ThreeVector outputFaceNormal = (*element)->OutputFaceNormal(); // outgoing face normal
314
315 G4ThreeVector iFNormal = outputFaceNormal;
316 iFNormal *= -1;
317 G4ThreeVector oFNormal = outputFaceNormal; // we assume no angle for the bridge component so this is right.
318
319 G4double width = (*element)->GetAcceleratorComponent()->GetExtent().DX();
320 width = std::min(width, suggestedRadius);
321 // we're ignoring any possible angled face of the curvilinear geometry
322 BDSSimpleComponent* component = factory->CreateCurvilinearVolume("clb_" + std::to_string(numberOfUniqueComponents),
323 defaultBridgeLength,
324 defaultBridgeLength,
325 width*0.5,
326 0, /*angle*/
327 iFNormal,
328 oFNormal);
329
331 numberOfUniqueComponents++;
332 return component;
333}
334
337 BDSBeamline::const_iterator nextElement,
339 const G4int beamlineIndex)
340{
341 G4bool last = nextElement == end;
342 BDSBeamlineElement* pel = (*element); // convenience
343
344 BDSTiltOffset* copyTiltOffset = nullptr;
345 BDSTiltOffset* existingTiltOffset = pel->GetTiltOffset();
346 if (existingTiltOffset)
347 {copyTiltOffset = new BDSTiltOffset(*existingTiltOffset);}
348
349 const G4RotationMatrix* refRotEnd = pel->GetReferenceRotationEnd();
350
351 G4ThreeVector unitZ = G4ThreeVector(0,0,1);
352 unitZ = unitZ.transform(*refRotEnd);
353
354 G4ThreeVector startPosition = pel->GetReferencePositionEnd();
355 G4ThreeVector endPosition;
356 if (!last)
357 {endPosition = (*nextElement)->GetReferencePositionStart();}
358 else
359 {endPosition = startPosition + ( (component->GetChordLength()) * unitZ );}
360 G4ThreeVector midPosition = (startPosition + endPosition) * 0.5;
361
362 G4double startS = pel->GetSPositionEnd();
363 G4double endS;
364 if (!last)
365 {endS = pel->GetSPositionStart();}
366 else
367 {endS = startS + component->GetChordLength();}
368 G4double midS = (startS + endS) * 0.5;
369
370 BDSBeamlineElement* result = new BDSBeamlineElement(component,
371 startPosition, midPosition, endPosition,
372 new G4RotationMatrix(*refRotEnd),
373 new G4RotationMatrix(*refRotEnd),
374 new G4RotationMatrix(*refRotEnd),
375 startPosition, midPosition, endPosition,
376 new G4RotationMatrix(*refRotEnd),
377 new G4RotationMatrix(*refRotEnd),
378 new G4RotationMatrix(*refRotEnd),
379 startS, midS, endS,
380 copyTiltOffset,
381 nullptr, // sampler info
382 beamlineIndex);
383 return result;
384}
385
387 BDSBeamline::const_iterator startElement,
388 BDSBeamline::const_iterator finishElement,
389 G4int index)
390{
391 BDSTiltOffset* copyTiltOffset = nullptr;
392 BDSTiltOffset* existingTiltOffset = (*startElement)->GetTiltOffset();
393 if (existingTiltOffset)
394 {copyTiltOffset = new BDSTiltOffset(*existingTiltOffset);}
395
396 BDSBeamlineElement* result = nullptr;
397
398 if (startElement == finishElement)
399 {// 1:1
400 BDSBeamlineElement* element = *startElement; // convenience
401
402 result = new BDSBeamlineElement(component,
403 element->GetReferencePositionStart(),
405 element->GetReferencePositionEnd(),
406 new G4RotationMatrix(*(element->GetRotationStart())),
407 new G4RotationMatrix(*(element->GetRotationMiddle())),
408 new G4RotationMatrix(*(element->GetRotationEnd())),
409 element->GetReferencePositionStart(),
411 element->GetReferencePositionEnd(),
412 new G4RotationMatrix(*(element->GetReferenceRotationStart())),
413 new G4RotationMatrix(*(element->GetReferenceRotationMiddle())),
414 new G4RotationMatrix(*(element->GetReferenceRotationEnd())),
415 element->GetSPositionStart(),
416 element->GetSPositionMiddle(),
417 element->GetSPositionEnd(),
418 copyTiltOffset,
419 nullptr, // sampler info
420 index);
421 }
422 else
423 {//must cover a few components
424 G4double sStart = (*startElement)->GetSPositionStart();
425 G4double sEnd = (*finishElement)->GetSPositionEnd();
426 G4double sMid = 0.5 * (sEnd + sStart);
427 G4ThreeVector posRefStart = (*startElement)->GetReferencePositionStart();
428 G4ThreeVector posRefEnd = (*finishElement)->GetReferencePositionEnd();
429 G4ThreeVector posRefMid = 0.5 * (posRefStart + posRefEnd);
430 G4RotationMatrix* rotRefStart = (*startElement)->GetReferenceRotationStart();
431 G4RotationMatrix* rotRefEnd = (*finishElement)->GetReferenceRotationEnd();
432
433 G4ThreeVector delta = posRefMid - posRefStart;
434 G4ThreeVector newUnitZ = delta.unit();
435 G4ThreeVector unitXPrevious = G4ThreeVector(1,0,0).transform(*rotRefStart);
436 G4ThreeVector newUnitY = newUnitZ.cross(unitXPrevious).unit();
437 G4ThreeVector unitYPrevious = G4ThreeVector(0,1,0).transform(*rotRefStart);
438 G4ThreeVector newUnitX = unitYPrevious.cross(newUnitZ).unit();
439
440 // create mid point rotation matrix from unit vectors at mid point
441 G4RotationMatrix rotRefMid = G4RotationMatrix(newUnitX, newUnitY, newUnitZ);
442
443 result = new BDSBeamlineElement(component,
444 posRefStart,
445 posRefMid,
446 posRefEnd,
447 new G4RotationMatrix(*rotRefStart),
448 new G4RotationMatrix(rotRefMid),
449 new G4RotationMatrix(*rotRefEnd),
450 posRefStart,
451 posRefMid,
452 posRefEnd,
453 new G4RotationMatrix(*rotRefStart),
454 new G4RotationMatrix(rotRefMid),
455 new G4RotationMatrix(*rotRefEnd),
456 sStart,
457 sMid,
458 sEnd,
459 copyTiltOffset,
460 nullptr, // sampler info
461 index);
462
463 }
464
465 return result;
466}
467
468
470{
471 // we're ignoring any possible angled face of the curvilinear geometry
472 BDSSimpleComponent* component = factory->CreateCurvilinearVolume("cl_start",
475
477
478 const BDSBeamlineElement* firstElement = beamline->GetFirstItem();
479 const G4RotationMatrix* rotStart = firstElement->GetReferenceRotationStart();
480
481 G4ThreeVector unitDir = G4ThreeVector(0,0,1).transform(*rotStart);
482 G4ThreeVector endPos = firstElement->GetReferencePositionStart() - unitDir*0.5*paddingLength;
483 G4ThreeVector midPos = endPos - unitDir*0.5*bonusChordLength;
484 G4ThreeVector staPos = endPos - unitDir*bonusChordLength;
485
486 G4double sStart = firstElement->GetSPositionStart();
487
488 BDSBeamlineElement*result = new BDSBeamlineElement(component,
489 staPos,
490 midPos,
491 endPos,
492 new G4RotationMatrix(*rotStart),
493 new G4RotationMatrix(*rotStart),
494 new G4RotationMatrix(*rotStart),
495 staPos,
496 midPos,
497 endPos,
498 new G4RotationMatrix(*rotStart),
499 new G4RotationMatrix(*rotStart),
500 new G4RotationMatrix(*rotStart),
501 sStart - bonusChordLength,
502 sStart - 0.5*bonusChordLength,
503 sStart,
504 nullptr,
505 nullptr, // sampler info
506 -1); // artificial index of -1 for before beam line
507 return result;
508}
509
511{
512 // we're ignoring any possible angled face of the curvilinear geometry
516
518
519 const BDSBeamlineElement* lastElement = beamline->GetLastItem();
520 const G4RotationMatrix* rotEnd = lastElement->GetReferenceRotationEnd();
521 const G4int lastIndex = lastElement->GetIndex();
522
523 G4ThreeVector unitDir = G4ThreeVector(0,0,1).transform(*rotEnd);
524 G4ThreeVector staPos = lastElement->GetReferencePositionEnd() + unitDir*0.5*paddingLength;
525 G4ThreeVector midPos = staPos + unitDir*0.5*bonusChordLength;
526 G4ThreeVector endPos = staPos + unitDir*bonusChordLength;
527
528 G4double sStart = lastElement->GetSPositionEnd();
529
530 BDSBeamlineElement*result = new BDSBeamlineElement(component,
531 staPos,
532 midPos,
533 endPos,
534 new G4RotationMatrix(*rotEnd),
535 new G4RotationMatrix(*rotEnd),
536 new G4RotationMatrix(*rotEnd),
537 staPos,
538 midPos,
539 endPos,
540 new G4RotationMatrix(*rotEnd),
541 new G4RotationMatrix(*rotEnd),
542 new G4RotationMatrix(*rotEnd),
543 sStart,
544 sStart + 0.5*bonusChordLength,
545 sStart + bonusChordLength,
546 nullptr,
547 nullptr, // sampler info
548 lastIndex + 1); // artificial index of -1 for before beam line
549 return result;
550}
void RegisterCurvilinearComponent(BDSAcceleratorComponent *component)
static BDSAcceleratorComponentRegistry * Instance()
Singleton accessor.
Abstract class that represents a component of an accelerator.
virtual G4double GetChordLength() const
G4ThreeVector OutputFaceNormal() const
A class that holds a fully constructed BDSAcceleratorComponent as well as any information relevant to...
G4double GetAngle() const
Accessor.
G4ThreeVector GetReferencePositionEnd() const
Accessor.
G4RotationMatrix * GetReferenceRotationStart() const
Accessor.
G4RotationMatrix * GetRotationEnd() const
Accessor.
BDSTiltOffset * GetTiltOffset() const
Accessor.
G4RotationMatrix * GetReferenceRotationEnd() const
Accessor.
G4int GetIndex() const
Accessor.
G4ThreeVector GetReferencePositionMiddle() const
Accessor.
G4double GetSPositionEnd() const
Accessor.
G4RotationMatrix * GetRotationStart() const
Accessor.
G4ThreeVector GetReferencePositionStart() const
Accessor.
G4RotationMatrix * GetReferenceRotationMiddle() const
Accessor.
G4double GetSPositionMiddle() const
Accessor.
G4double GetArcLength() const
Accessor.
G4RotationMatrix * GetRotationMiddle() const
Accessor.
G4double GetSPositionStart() const
Accessor.
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
const BDSBeamlineElement * GetLastItem() const
Return a reference to the last element.
Definition: BDSBeamline.hh:126
static G4double PaddingLength()
Access the padding length between each element added to the beamline.
Definition: BDSBeamline.hh:233
const BDSBeamlineElement * GetFirstItem() const
Return a reference to the first element.
Definition: BDSBeamline.hh:123
void AddBeamlineElement(BDSBeamlineElement *element)
Definition: BDSBeamline.cc:571
G4bool empty() const
Iterator mechanics.
Definition: BDSBeamline.hh:175
iterator begin()
Iterator mechanics.
Definition: BDSBeamline.hh:167
iterator end()
Iterator mechanics.
Definition: BDSBeamline.hh:168
BeamlineVector::const_iterator const_iterator
Iterator mechanics.
Definition: BDSBeamline.hh:164
BDSBeamlineElement * CreateElementFromComponent(BDSSimpleComponent *component, BDSBeamline::const_iterator startElement, BDSBeamline::const_iterator finishElement, G4int index)
Create the BDSBeamlineElement by wrapping a BDSSimpleComponent.
BDSBeamlineElement * CreateBridgeSection(BDSAcceleratorComponent *defaultBridge, BDSBeamline::const_iterator element, BDSBeamline::const_iterator nextElement, BDSBeamline::const_iterator end, G4int &numberOfUniqueComponents, const G4int beamlineIndex, G4double crRadius)
BDSAcceleratorComponent * CreateStraightBridgeComponent(G4double width, G4int &numberOfUniqueComponents)
BDSBeamlineElement * CreateBridgeElementFromComponent(BDSAcceleratorComponent *component, BDSBeamline::const_iterator element, BDSBeamline::const_iterator nextElement, BDSBeamline::const_iterator end, const G4int beamlineIndex)
G4double curvilinearRadius
Radius for curvilinear geometry.
G4double paddingLength
Length that was used for padding on the beam line we're building with respesct to.
void PreviousAndNext(BDSBeamline::const_iterator it, BDSBeamline::const_iterator startIt, BDSBeamline::const_iterator endIt, const BDSBeamlineElement *&previous, const BDSBeamlineElement *&next) const
G4double CurvilinearRadius(const BDSBeamlineElement *el) const
BDSBeamline * BuildCurvilinearBeamLine1To1(BDSBeamline const *const beamline, const G4bool circular)
Build a beam line of curvilinear geometry based on another beam line.
BDSAcceleratorComponent * CreateAngledBridgeComponent(BDSBeamline::const_iterator element, G4int &numberOfUniqueComponents, G4double suggestedRadius)
Create a unique accelerator component for an element with angled faces.
G4bool Angled(BDSBeamlineElement const *const element) const
Simple interrogation function to determine if an element has a finite angle or not.
BDSBeamline * BuildCurvilinearBridgeBeamLine(BDSBeamline const *const beamline)
Build bridging volumes to join the curvilinear ones.
BDSCurvilinearFactory * factory
Factory to build curvilinear geometry.
BDSBeamlineElement * CreateBonusSectionStart(BDSBeamline const *const beamline)
Create a small section to extend a linear beam line.
BDSAcceleratorComponent * CreateDefaultBridgeComponent()
G4double bonusChordLength
Length of any possible bonus section added to beginning and end.
BDSBeamlineElement * CreateBonusSectionEnd(BDSBeamline const *const beamline)
Create a small section to extend a linear beam line.
BDSBeamlineElement * CreateCurvilinearElement(const G4String &elementName, BDSBeamline::const_iterator startElement, BDSBeamline::const_iterator finishElement, G4int index, G4double crRadius)
Factory to create curvilinear geometry for parallel world.
BDSSimpleComponent * CreateCurvilinearVolume(const G4String &name, G4double chordLength, G4double radius)
Build a straight curvilinear volume.
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
BDSExtent Translate(const G4ThreeVector &offset) const
Provide a new copy of this extent with an offset applied.
Definition: BDSExtent.hh:125
G4double MaximumAbs() const
Return the maximum absolute value considering all dimensions.
Definition: BDSExtent.cc:157
A class that holds global options and constants.
static BDSGlobalConstants * Instance()
Access method.
A BDSAcceleratorComponent wrapper for BDSGeometryComponent.
A holder for any placement offsets and rotations for a BDSAcceleratorComponent.
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())