BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBendBuilder.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 "globals.hh" // geant4 globals / types
20
21#include "BDSAcceleratorComponent.hh"
22#include "BDSBeamPipeInfo.hh"
23#include "BDSBendBuilder.hh"
24#include "BDSComponentFactory.hh"
25#include "BDSDebug.hh"
26#include "BDSException.hh"
27#include "BDSFieldType.hh"
28#include "BDSGlobalConstants.hh"
29#include "BDSIntegratorSet.hh"
30#include "BDSIntegratorType.hh"
31#include "BDSIntegratorSetType.hh"
32#include "BDSMagnetStrength.hh"
33#include "BDSLine.hh"
34#include "BDSMagnet.hh"
35#include "BDSMagnetOuterInfo.hh"
36#include "BDSUtilities.hh"
37
38#include "parser/element.h"
39#include "parser/elementtype.h"
40
41#include "G4Transform3D.hh"
42
43using namespace GMAD;
44
46{
47 G4bool finiteAngle = BDS::IsFinite((*st)["angle"]);
48 G4bool finiteField = BDS::IsFinite((*st)["field"]);
49 return !(finiteAngle || finiteField);
50}
51
52BDSAcceleratorComponent* BDS::BuildSBendLine(const G4String& elementName,
53 const Element* element,
55 G4double brho,
56 const BDSIntegratorSet* integratorSet,
57 G4double incomingFaceAngle,
58 G4double outgoingFaceAngle,
59 G4bool buildFringeFields,
60 const GMAD::Element* prevElement,
61 const GMAD::Element* nextElement)
62{
63 const G4String baseName = elementName;
64 const G4double thinElementArcLength = BDSGlobalConstants::Instance()->ThinElementLength();
65 const G4bool yokeOnLeft = BDSComponentFactory::YokeOnLeft(element,st);
66 const G4double arcLength = element->l * CLHEP::m;
67 const G4double angle = (*st)["angle"];
68 G4double fint = element->fint;
69 G4double fintx = element->fintx;
70 G4double hgap = element->hgap * CLHEP::m;
71
72 // Note for tilted dipoles, the geometry is tilted but the curvilinear world isn't,
73 // therefore we tilt the field to match the geometry.
74 G4Transform3D fieldTiltOffset = BDSComponentFactory::CreateFieldTransform(element);
75
76 G4bool buildFringeIncoming = buildFringeFields;
77 G4bool buildFringeOutgoing = buildFringeFields;
78
79 G4bool finiteK1 = BDS::IsFinite((*st)["k1"]);
80 BDSFieldType dipoleFieldType = finiteK1 ? BDSFieldType::dipolequadrupole : BDSFieldType::dipole;
81
82 if (buildFringeFields)
83 {
84 // perform series of checks on fringe field parameters
85 // nominally, don't build the fringe element if there's no face angle
86 if (!BDS::IsFinite(element->e1))
87 {buildFringeIncoming = false;} // by default false
88 if (!BDS::IsFinite(element->e2))
89 {buildFringeOutgoing = false;}
90
91 // fintx = -1 is how madx writes out that fintx should default to fint
92 // it's also our default. If by default fint=0 and fintx=0, no fringes
93 // would be built. If finite fint and -1 for fintx both will be built.
94 // if fint=0 and fintx != -1, only the outgoin will be built
95 if (fintx == -1)
96 {fintx = fint;}
97
98 // however if finite hgap and fint or fintx specified, there is an effect
99 if (BDS::IsFinite(fint) && BDS::IsFinite(hgap))
100 {buildFringeIncoming = true;}
101 if (BDS::IsFinite(fintx) && BDS::IsFinite(hgap))
102 {buildFringeOutgoing = true;}
103
104 // don't build fringe elements if a non-matrix integrator set is used and
105 // only a pole face effect is specified. The angled pole face geometry is
106 // constructed so a thin integrated pole face kick shouldn't be applied as well
107 if (!integratorSet->IsMatrixIntegratorSet() && BDSGlobalConstants::Instance()->BuildPoleFaceGeometry())
108 {
109 // no hgap check - finite hgap but 0 fint/x will not produce a fringe kick so fringes won't be built anyway
110 if (BDS::IsFinite(element->e1) && !BDS::IsFinite(fint))
111 {buildFringeIncoming = false;}
112 if (BDS::IsFinite(element->e2) && !BDS::IsFinite(fintx))
113 {buildFringeOutgoing = false;}
114 }
115
116 // overriding checks - don't build fringe field if we're butt against another
117 // sbend.
118 if (prevElement)
119 {
120 if (prevElement->type == ElementType::_SBEND)
121 {
122 buildFringeIncoming = false;
123 if (BDS::IsFinite(prevElement->e2 - element->e1))
124 {throw BDSException( __METHOD_NAME__, prevElement->name + " e2 clashes with " + elementName + " e1");}
125 }
126 }
127 if (nextElement)
128 {
129 if (nextElement->type == ElementType::_SBEND)
130 {buildFringeOutgoing = false;}
131 // check above on clashing sbends with e1 e2 will be used for next bend.
132 }
133
134 if (!BDS::IsFinite(angle))
135 {// no fringe fields if magnet of zero angle, ie strength
136 buildFringeIncoming = false;
137 buildFringeOutgoing = false;
138 }
139 }// end of checks
140
141 // Calculate number of sbends to split parent into
142 G4int nSBends = BDS::CalculateNSBendSegments(arcLength, angle, incomingFaceAngle,
143 outgoingFaceAngle);
144
145 G4bool zeroStrength = BDS::ZeroStrengthDipole(st);
146 if (zeroStrength)
147 {
148 buildFringeIncoming = false;
149 buildFringeOutgoing = false;
150 }
151
152 // If no angle, means no fringe strength, so even if there are pole face rotations, we don't
153 // need the fringe elements. If there's only 1 segment we may still have a finite
154 // angle and therefore finite strength and the pole face rotations will require fringe
155 // elements
156 if (!BDS::IsFinite(angle))
157 {
158 BDSIntegratorType intType = BDS::GetDipoleIntegratorType(integratorSet, element);
159 BDSFieldInfo* vacuumField = nullptr;
160 BDSFieldInfo* outerField = nullptr;
161 // prepare one sbend segment
162 auto bpInfo = BDSComponentFactory::PrepareBeamPipeInfo(element, -incomingFaceAngle,
163 -outgoingFaceAngle);
164 auto mgInfo = BDSComponentFactory::PrepareMagnetOuterInfo(baseName, element,
165 -incomingFaceAngle,
166 -outgoingFaceAngle,
167 bpInfo,
168 yokeOnLeft);
169 if (!zeroStrength)
170 {
171 vacuumField = new BDSFieldInfo(dipoleFieldType,
172 brho,
173 intType,
174 st,
175 true,
176 fieldTiltOffset);
178 BDSFieldType::dipole,
179 bpInfo,
180 mgInfo,
181 fieldTiltOffset,
182 integratorSet,
183 brho,
184 BDSComponentFactory::ScalingFieldOuter(element));
185 }
186 BDSMagnet* oneBend = new BDSMagnet(BDSMagnetType::sectorbend,
187 baseName,
188 arcLength,
189 bpInfo,
190 mgInfo,
191 vacuumField,
192 -angle,
193 outerField); // angle 0, but should be minus for 3d cart.
194 return oneBend;
195 }
196
197 // We're definitely building a line now.
198 BDSLine* sbendline = new BDSLine(baseName); // line for resultant sbend
199
200 // length of bends minus the fringe segments
201 G4double centralArcLength = arcLength;
202 G4double centralAngle = angle;
203 G4double oneFringeAngle = 0;
204 if (BDS::IsFinite(angle))
205 {oneFringeAngle = (thinElementArcLength / arcLength) * angle;}
206
207 if (buildFringeIncoming)
208 {
209 centralArcLength -= thinElementArcLength;
210 centralAngle -= oneFringeAngle;
211 }
212 if (buildFringeOutgoing)
213 {
214 centralArcLength -= thinElementArcLength;
215 centralAngle -= oneFringeAngle;
216 }
217
218 // calculate their angles and length
219 G4double semiAngle = centralAngle / (G4double) nSBends;
220 G4double semiArcLength = centralArcLength / (G4double) nSBends;
221
222 BDSMagnetStrength* semiStrength = new BDSMagnetStrength(*st); // the copy is crucial to copy the field strength
223 (*semiStrength)["angle"] = semiAngle;
224 (*semiStrength)["length"] = semiArcLength;
225
226 G4double zExtentIn = 0;
227 G4double zExtentOut = 0;
228 G4bool fadeIn = true;
229 G4bool fadeOut = true;
230
231 // unlike an rbend, the sbend will mostly likely be split up into segments.
232 // we must check that the faces of each segment (varying from e1 to e2) will
233 // not overlap given the outer diameter.
234 // calculate extent along z due poleface rotation at half the horizontal width.
235 // note this is only w.r.t. pole face rotation and not the general bending angle
236 G4double horizontalWidth = BDSComponentFactory::PrepareHorizontalWidth(element);
237 if (incomingFaceAngle > 0)
238 {zExtentIn = 0.5*horizontalWidth*std::tan(incomingFaceAngle - 0.5*std::abs(semiAngle));}
239 else if (incomingFaceAngle < 0)
240 {zExtentIn = 0.5*horizontalWidth*std::tan(0.5*std::abs(semiAngle) + incomingFaceAngle);}
241 if (outgoingFaceAngle > 0)
242 {zExtentOut = 0.5*horizontalWidth*std::tan(outgoingFaceAngle - 0.5*std::abs(semiAngle));}
243 else if (outgoingFaceAngle < 0)
244 {zExtentOut = 0.5*horizontalWidth*std::tan(0.5*std::abs(semiAngle) + outgoingFaceAngle);}
245
246 // decide if segment angles fade or not depending on the extents
247 if (std::abs(zExtentIn) < semiArcLength/4.0)
248 {fadeIn = false;}
249 if (std::abs(zExtentOut) < semiArcLength/4.0)
250 {fadeOut = false;}
251
252 // field recipe for one segment of the sbend
253 G4String centralName = baseName + "_even_ang";
254 BDSIntegratorType intType = BDS::GetDipoleIntegratorType(integratorSet, element);
255 BDSFieldInfo* semiVacuumField = nullptr;
256 BDSFieldInfo* semiOuterField = nullptr;
257 auto bpInfo = BDSComponentFactory::PrepareBeamPipeInfo(element, 0.5*semiAngle, 0.5*semiAngle);
258 auto mgInfo = BDSComponentFactory::PrepareMagnetOuterInfo(centralName, element,
259 0.5*semiAngle, 0.5*semiAngle, bpInfo,
260 yokeOnLeft);
261 if (!zeroStrength)
262 {
263 semiVacuumField = new BDSFieldInfo(dipoleFieldType,
264 brho,
265 intType,
266 semiStrength,
267 true,
268 fieldTiltOffset);
269 semiOuterField = BDSComponentFactory::PrepareMagnetOuterFieldInfo(semiStrength,
270 BDSFieldType::dipole,
271 bpInfo,
272 mgInfo,
273 fieldTiltOffset,
274 integratorSet,
275 brho,
276 BDSComponentFactory::ScalingFieldOuter(element));
277 }
278 mgInfo->name = centralName;
279 BDSMagnet* centralWedge = new BDSMagnet(BDSMagnetType::sectorbend,
280 centralName,
281 semiArcLength,
282 bpInfo,
283 mgInfo,
284 semiVacuumField,
285 -semiAngle,
286 semiOuterField); // minus for 3d cartesian conversion
287
288 // check magnet outer info
289 auto magnetOuterInfoCheck = BDSComponentFactory::PrepareMagnetOuterInfo("checking", element,
290 -incomingFaceAngle,
291 -outgoingFaceAngle,
292 bpInfo, yokeOnLeft);
293 // minus for conversion to 3d cartesian
294 G4double minimalRadius = 2*magnetOuterInfoCheck->MinimumIntersectionRadiusRequired();
295 // extra pedantic check for dipoles with certain geometry types
296 if (!magnetOuterInfoCheck->hStyle)
297 {// in the case of C-shaped poled dipoles, the full 'horizontalWidth' is shifted to the yoke side
298 switch (magnetOuterInfoCheck->geometryType.underlying())
299 {
300 case BDSMagnetGeometryType::polescircular:
301 case BDSMagnetGeometryType::polesfacet:
302 case BDSMagnetGeometryType::polesfacetcrop:
303 case BDSMagnetGeometryType::polessquare:
304 {
305 minimalRadius *= element->yokeOnInside ? 2.0 : 0.5;
306 break;
307 }
308 default:
309 {break;}
310 }
311 }
312 BDSComponentFactory::CheckBendLengthAngleWidthCombo(semiArcLength, -semiAngle, minimalRadius, element->name);
313 delete magnetOuterInfoCheck; // clean up
314
315 // build incoming fringe field if required
316 if (buildFringeIncoming)
317 {
318 G4double e1 = element->e1;
319 // if using a non-matrix integrator set, this code should only be reached if a finite fint is specified.
320 // As the angled face geometry is constructed, a thin pole face kick shouldn't also be applied, so set
321 // e1 to 0 in the magnet strength object so only fringe effects are applied
322 if (!integratorSet->IsMatrixIntegratorSet() && BDSGlobalConstants::Instance()->BuildPoleFaceGeometry()
323 && BDS::IsFinite(fint))
324 {e1 = 0;}
325
326 BDSMagnetStrength* fringeStIn = BDS::GetFringeMagnetStrength(element, st, oneFringeAngle,
327 e1, element->e2, fintx, true);
328 G4String segmentName = baseName + "_e1_fringe";
329 G4double fringeAngleIn = 0.5*oneFringeAngle - incomingFaceAngle;
330 G4double fringeAngleOut = 0.5*oneFringeAngle + incomingFaceAngle;
331 BDSMagnet* startfringe = BDS::BuildDipoleFringe(element, fringeAngleIn, fringeAngleOut,
332 segmentName, fringeStIn, brho,
333 integratorSet, dipoleFieldType);
334 sbendline->AddComponent(startfringe);
335 }
336
337 // logic for elements in the sbend line:
338 // reuse central segment for all segment of in/out half if no poleface angle(s)
339 // if small poleface, new first/last segment, reuse central wedge for remainder of in/out half
340 // otherwise fade in/out faces for all segments in app. halves.
341 // 'central' one is definitely used for the central part, but also it's just a segment
342 // with even incoming and outgoing face angles w.r.t. the chord.
343 G4double segmentAngleIn = 0;
344 G4double segmentAngleOut = 0;
345 G4int numberOfUniqueComponents = 1; // used for naming purposes
346 BDSMagnet* oneBend = nullptr;
347 G4bool centralWedgeUsed = false; // keep track to avoid memory leak
348 const G4int numSegmentsEitherSide = (nSBends - 1) / 2;
349 for (G4int i = 0; i < nSBends; ++i)
350 {
351 G4String name = baseName;
352 if (i < numSegmentsEitherSide)
353 {// first half of magnet
354 if (!BDS::IsFinite(incomingFaceAngle)) // no pole face rotation so just repeat central segment
355 {oneBend = centralWedge;}
356 else if (fadeIn) // build incremented angled segment
357 {
358 name += "_"+std::to_string(numberOfUniqueComponents);
359 numberOfUniqueComponents++;
360 BDS::UpdateSegmentAngles(i,nSBends,semiAngle,incomingFaceAngle,outgoingFaceAngle,segmentAngleIn,segmentAngleOut);
361 oneBend = BDS::BuildSingleSBend(element, name, semiArcLength, semiAngle,
362 segmentAngleIn, segmentAngleOut, semiStrength,
363 brho, integratorSet, yokeOnLeft, semiOuterField);
364 }
365 else
366 {// finite pole face, but not strong so build one angled, then repeat the rest to save memory
367 if (i == 0) // the first one is unique
368 {
369 name += "_"+std::to_string(numberOfUniqueComponents);
370 numberOfUniqueComponents++;
371 segmentAngleIn = 0.5*semiAngle - incomingFaceAngle; // whole pole face angle
372 segmentAngleOut = 0.5*semiAngle; // even matching angle
373 oneBend = BDS::BuildSingleSBend(element, name, semiArcLength, semiAngle,
374 segmentAngleIn, segmentAngleOut, semiStrength,
375 brho, integratorSet, yokeOnLeft, semiOuterField);
376 }
377 else // others afterwards are a repeat of the even angled one
378 {oneBend = centralWedge;}
379 }
380 }
381 else if (i > numSegmentsEitherSide) // i is zero counting
382 {// second half of magnet
383 if (!BDS::IsFinite(outgoingFaceAngle)) // no pole face rotation so just repeat central segment
384 {oneBend = centralWedge;}
385 else if (fadeOut) // build incremented angled segment
386 {
387 name += "_"+std::to_string(numberOfUniqueComponents);
388 numberOfUniqueComponents++;
389 BDS::UpdateSegmentAngles(i,nSBends,semiAngle,incomingFaceAngle,outgoingFaceAngle,segmentAngleIn,segmentAngleOut);
390 oneBend = BDS::BuildSingleSBend(element, name, semiArcLength, semiAngle,
391 segmentAngleIn, segmentAngleOut, semiStrength,
392 brho, integratorSet, yokeOnLeft, semiOuterField);
393 }
394 else
395 {// finite pole face, but not strong so build only one unique angled on output face
396 if (i == (nSBends-1)) // last one
397 {
398 name += "_"+std::to_string(numberOfUniqueComponents);
399 numberOfUniqueComponents++;
400 segmentAngleIn = 0.5*semiAngle;
401 segmentAngleOut = 0.5*semiAngle - outgoingFaceAngle;
402 oneBend = BDS::BuildSingleSBend(element, name, semiArcLength, semiAngle,
403 segmentAngleIn, segmentAngleOut, semiStrength,
404 brho, integratorSet, yokeOnLeft, semiOuterField);
405 }
406 else // after central, but before unique end piece - even angled.
407 {oneBend = centralWedge;}
408 }
409 }
410 else // the middle piece
411 {oneBend = centralWedge;}
412
413 // append to the line
414 sbendline->AddComponent(oneBend);
415
416 centralWedgeUsed = centralWedgeUsed || (oneBend == centralWedge);
417 }
418
419 if (!centralWedgeUsed)
420 {delete centralWedge;}
421
422 //Last element should be fringe if poleface specified
423 if (buildFringeOutgoing)
424 {
425 G4double e2 = element->e2;
426 // if using a non-matrix integrator set, this code should only be reached if a finite fintx is specified.
427 // As the angled face geometry is constructed, a thin pole face kick shouldn't also be applied, so set
428 // e2 to 0 in the magnet strength object so only fringe effects are applied
429 if (!integratorSet->IsMatrixIntegratorSet() && BDSGlobalConstants::Instance()->BuildPoleFaceGeometry()
430 && BDS::IsFinite(fintx))
431 {e2 = 0;}
432
433 BDSMagnetStrength* fringeStOut = BDS::GetFringeMagnetStrength(element, st, oneFringeAngle,
434 element->e1, e2, fintx, false);
435 G4double fringeAngleIn = 0.5*oneFringeAngle + outgoingFaceAngle;
436 G4double fringeAngleOut = 0.5*oneFringeAngle - outgoingFaceAngle;
437 G4String segmentName = baseName + "_e2_fringe";
438
439 BDSMagnet* endfringe = BDS::BuildDipoleFringe(element, fringeAngleIn, fringeAngleOut,
440 segmentName, fringeStOut, brho,
441 integratorSet, dipoleFieldType);
442 sbendline->AddComponent(endfringe);
443 }
444 return sbendline;
445}
446
447void BDS::UpdateSegmentAngles(G4int index,
448 G4int nSBends,
449 G4double semiAngle,
450 G4double incomingFaceAngle,
451 G4double outgoingFaceAngle,
452 G4double& segmentAngleIn,
453 G4double& segmentAngleOut)
454{
455 // nSBends is always odd.
456 // index assumed 0 counting.
457 G4bool central = index == (nSBends - 1) / 2;
458 if (central)
459 {// shouldn't happen, but just in case
460 segmentAngleIn = 0.5*semiAngle;
461 segmentAngleOut = 0.5*semiAngle;
462 return;
463 }
464
465 G4bool firstHalf = index + 1 < (nSBends + 1) / 2;
466 G4int numberToFadeOver = (nSBends - 1) / 2;
467 if (firstHalf)
468 {// fade from incomingFaceAngle to 0
469 G4double delta = incomingFaceAngle / (G4double) numberToFadeOver;
470 G4double inputFaceAngle = incomingFaceAngle - ((G4double)index * delta);
471 G4double outputFaceAngle = incomingFaceAngle - ((G4double)(index+1) * delta);
472 // if incomingFaceAngle was 0, each would have 1/2 of semiAngle
473 segmentAngleIn = 0.5*semiAngle - inputFaceAngle;
474 segmentAngleOut = 0.5*semiAngle + outputFaceAngle;
475 }
476 else
477 {// second half - fade from 0 to outgoingFaceAngle
478 G4double delta = outgoingFaceAngle / (G4double) numberToFadeOver;
479 G4int secondHalfIndex = index - ((nSBends + 1) / 2);
480 // here we fade from 0 contribution to maximum (opposite of above if statement)
481 G4double inputFaceAngle = outgoingFaceAngle - ( (G4double)(numberToFadeOver - secondHalfIndex) * delta);
482 G4double outputFaceAngle = outgoingFaceAngle - ( (G4double)(numberToFadeOver - (secondHalfIndex + 1)) * delta);
483 // if incomingFaceAngle was 0, each would have 1/2 of semiAngle
484 segmentAngleIn = 0.5*semiAngle + inputFaceAngle;
485 segmentAngleOut = 0.5*semiAngle - outputFaceAngle;
486 }
487}
488
490 const G4String& name,
491 G4double arcLength,
492 G4double angle,
493 G4double angleIn,
494 G4double angleOut,
495 const BDSMagnetStrength* strength,
496 G4double brho,
497 const BDSIntegratorSet* integratorSet,
498 G4bool yokeOnLeft,
499 const BDSFieldInfo* outerFieldIn)
500{
501 auto bpInfo = BDSComponentFactory::PrepareBeamPipeInfo(element, angleIn, angleOut);
502
503 BDSMagnetStrength* strengthCopy = new BDSMagnetStrength(*strength); // the copy is crucial to copy the field strength
504 auto magnetOuterInfo = BDSComponentFactory::PrepareMagnetOuterInfo(name, element, angleIn, angleOut, bpInfo, yokeOnLeft);
505 // set the name to the desired one rather than the one from the element
506 magnetOuterInfo->name = name;
507
508 G4Transform3D fieldTiltOffset = BDSComponentFactory::CreateFieldTransform(element);
509
510 G4bool finiteK1 = BDS::IsFinite((*strength)["k1"]);
511 BDSFieldType dipoleFieldType = finiteK1 ? BDSFieldType::dipolequadrupole : BDSFieldType::dipole;
512 BDSIntegratorType intType = BDS::GetDipoleIntegratorType(integratorSet, element);
513 BDSFieldInfo* vacuumField = nullptr;
514 BDSFieldInfo* outerField = nullptr;
515 G4bool zeroStrength = BDS::ZeroStrengthDipole(strengthCopy);
516 if (!zeroStrength)
517 {
518 vacuumField = new BDSFieldInfo(dipoleFieldType,
519 brho,
520 intType,
521 strengthCopy,
522 true,
523 fieldTiltOffset);
524 outerField = new BDSFieldInfo(*outerFieldIn);
525 }
526
527 BDSMagnet* magnet = new BDSMagnet(BDSMagnetType::sectorbend,
528 name,
529 arcLength,
530 bpInfo,
531 magnetOuterInfo,
532 vacuumField,
533 -angle,
534 outerField);
535
536 return magnet;
537}
538
539BDSLine* BDS::BuildRBendLine(const G4String& elementName,
540 const Element* element,
541 const Element* prevElement,
542 const Element* nextElement,
543 G4double brho,
545 const BDSIntegratorSet* integratorSet,
546 G4double incomingFaceAngle,
547 G4double outgoingFaceAngle,
548 G4bool buildFringeFields)
549{
550 const G4String name = elementName;
551 BDSLine* rbendline = new BDSLine(name); // line for resultant rbend
552
553 const G4double thinElementArcLength = BDSGlobalConstants::Instance()->ThinElementLength();
554 const G4bool yokeOnLeft = BDSComponentFactory::YokeOnLeft(element, st);
555 G4bool buildFringeIncoming = buildFringeFields;
556 G4bool buildFringeOutgoing = buildFringeFields;
557 G4double fint = element->fint;
558 G4double fintx = element->fintx;
559 G4double hgap = element->hgap * CLHEP::m;
560
561 // Angle here is in the 'strength' convention of +ve angle -> -ve x deflection
562 const G4double angle = (*st)["angle"];
563 const G4double arcLength = (*st)["length"];
564
565 G4bool zeroStrength = BDS::ZeroStrengthDipole(st);
566
567 G4Transform3D fieldTiltOffset = BDSComponentFactory::CreateFieldTransform(element);
568
569 G4bool finiteK1 = BDS::IsFinite((*st)["k1"]);
570 BDSFieldType dipoleFieldType = finiteK1 ? BDSFieldType::dipolequadrupole : BDSFieldType::dipole;
571
572 // Here, 'no face angle' really means that the rbend becomes an sbend.
573 // Calculate how far away we are from an sbend.
574 G4double incomingFaceAngleWRTSBend = incomingFaceAngle + angle*0.5;
575 G4double outgoingFaceangleWRTSBend = outgoingFaceAngle + angle*0.5;
576 if (!BDS::IsFinite(incomingFaceAngleWRTSBend) && (!integratorSet->IsMatrixIntegratorSet()))
577 {buildFringeIncoming = false;}
578 if (!BDS::IsFinite(outgoingFaceangleWRTSBend) && (!integratorSet->IsMatrixIntegratorSet()))
579 {buildFringeOutgoing = false;}
580
581 // fintx = -1 is how madx writes out that fintx should default to fint
582 // it's also our default. If by default fint=0 and fintx=0, no fringes
583 // would be built. If finite fint and -1 for fintx both will be built.
584 // if fint=0 and fintx != -1, only the outgoin will be built
585 if (fintx == -1)
586 {fintx = fint;}
587
588 // however if finite hgap and fint or fintx specified, there is an effect
589 if (BDS::IsFinite(fint) && BDS::IsFinite(hgap))
590 {buildFringeIncoming = true;}
591 if (BDS::IsFinite(fintx) && BDS::IsFinite(hgap))
592 {buildFringeOutgoing = true;}
593
594 // the poleface angles to be used in tracking only.
595 G4double trackingPolefaceAngleIn = element->e1;
596 G4double trackingPolefaceAngleOut = element->e2;
597 if (integratorSet->IsMatrixIntegratorSet())
598 {// for this integrator set we track the rbend as an sbend with extra pole face rotation
599 trackingPolefaceAngleIn += angle*0.5;
600 trackingPolefaceAngleOut += angle*0.5;
601 }
602
603 G4double e1 = -incomingFaceAngle;
604 G4double e2 = -outgoingFaceAngle;
605
606 if (prevElement)
607 {
608 if (BDS::IsFinite(prevElement->e2) && BDS::IsFinite(element->e1))
609 {
610 throw BDSException(__METHOD_NAME__, prevElement->name + " has finite e2!\n Clashes with " + elementName + " with finite e1");
611 }
612 if (prevElement->type == ElementType::_RBEND)
613 {
614 buildFringeIncoming = false;
615 e1 = angle*0.5;
616 }
617 }
618 if (nextElement)
619 {
620 if (BDS::IsFinite(nextElement->e1) && BDS::IsFinite(element->e2))
621 {
622 throw BDSException(__METHOD_NAME__ + nextElement->name + " has finite e1!\n Clashes with " + elementName + " with finite e2");
623 }
624 if (nextElement->type == ElementType::_RBEND)
625 {
626 buildFringeOutgoing = false;
627 e2 = angle*0.5;
628 }
629 }
630
631 if (zeroStrength)
632 {
633 buildFringeIncoming = false;
634 buildFringeOutgoing = false;
635 }
636
637 // used for debugging purposes to forefully try out one and not the other fringe
638 //buildFringeIncoming = false;
639 //buildFringeOutgoing = false;
640
641 // if we're building the fringe elements, we reduce the length of the central section
642 // and proportion the bending by their length w.r.t. the full length of the total rbend.
643 // angles of the faces are in order:
644 // angleIn / fringeInOutputAngle / centralInputFaceAngle / centralOutputFaceAngle
645 // fringeOutInputAngle / angleOut
646 // default face angles for an rbend are 0 - ie parallel faces, plus any pole face rotation
647 // angle in and out of total rbend are nominally the face angles.
648 G4double angleIn = incomingFaceAngle;
649 G4double angleOut = outgoingFaceAngle;
650 G4double fringeInOutputAngle = 0;
651 G4double centralInputFaceAngle = angleIn;
652 G4double centralOutputFaceAngle = angleOut;
653 G4double fringeOutInputAngle = 0;
654 G4double centralArcLength = arcLength; // full length to start with
655 G4double centralAngle = angle;
656 G4double oneFringeAngle = 0;
657 if (BDS::IsFinite(angle))
658 {oneFringeAngle = (thinElementArcLength / arcLength) * angle;}
659
660 if (buildFringeIncoming && buildFringeOutgoing)
661 {// both
662 centralArcLength -= 2*thinElementArcLength;
663 centralAngle -= 2*oneFringeAngle;
664 angleIn = e1 + (0.5*oneFringeAngle - 0.5*angle);
665 fringeInOutputAngle = -angleIn;
666 centralInputFaceAngle = e1;
667 centralOutputFaceAngle = e2;
668 fringeOutInputAngle = - (e2 + (0.5*oneFringeAngle - 0.5*angle));
669 angleOut = -fringeOutInputAngle;
670 }
671 else if (buildFringeIncoming)
672 {// incoming only
673 centralArcLength -= thinElementArcLength;
674 centralAngle -= oneFringeAngle;
675 angleIn = e1 + (0.5*oneFringeAngle - 0.5*angle);
676 fringeInOutputAngle = -angleIn; // fringe built the same
677 centralInputFaceAngle = e1 + 0.5*oneFringeAngle;
678 centralOutputFaceAngle = e2 - 0.5*oneFringeAngle;
679 }
680 else if (buildFringeOutgoing)
681 {// outgoing only
682 centralArcLength -= thinElementArcLength;
683 centralAngle -= oneFringeAngle;
684 centralInputFaceAngle = e1 - 0.5*oneFringeAngle;
685 centralOutputFaceAngle = e2 + 0.5*oneFringeAngle;
686 fringeOutInputAngle = - (e2 + (0.5*oneFringeAngle - 0.5*angle));;
687 angleOut = e2 + (0.5*oneFringeAngle - 0.5*angle);;
688 }
689 else
690 {// neither
691 // centralArcLength the same
692 // centralAngle the same
693 centralInputFaceAngle = e1;
694 centralOutputFaceAngle = e2;
695 }
696
697 if (buildFringeIncoming)
698 {
699 G4double trackingPolefaceAngle = trackingPolefaceAngleIn;
700 // if using a non-matrix integrator set, this code should only be reached if a finite fint is specified.
701 // As the angled face geometry is constructed, a thin pole face kick shouldn't also be applied, so subtract
702 // e1 for in the magnet strength object so only the fringe & the natural rbend angled face effects are applied
703 if (!integratorSet->IsMatrixIntegratorSet() && BDSGlobalConstants::Instance()->BuildPoleFaceGeometry()
704 && BDS::IsFinite(fint))
705 {trackingPolefaceAngle -= element->e1;}
706
707 BDSMagnetStrength* fringeStIn = BDS::GetFringeMagnetStrength(element, st, oneFringeAngle,
708 trackingPolefaceAngle, trackingPolefaceAngleOut,
709 fintx, true);
710 G4String fringeName = name + "_e1_fringe";
711
712 // element used for beam pipe materials etc - not strength, angle or length.
713 BDSMagnet* startfringe = BDS::BuildDipoleFringe(element, angleIn, fringeInOutputAngle,
714 fringeName,
715 fringeStIn, brho,
716 integratorSet, dipoleFieldType);
717 rbendline->AddComponent(startfringe);
718 }
719
720 // override copied length and angle
721 (*st)["length"] = centralArcLength;
722 (*st)["angle"] = centralAngle;
723
724 BDSIntegratorType intType = BDS::GetDipoleIntegratorType(integratorSet, element);
725 BDSFieldInfo* vacuumField = nullptr;
726 BDSFieldInfo* outerField = nullptr;
727 auto bpInfo = BDSComponentFactory::PrepareBeamPipeInfo(element, centralInputFaceAngle, centralOutputFaceAngle);
728 auto mgInfo = BDSComponentFactory::PrepareMagnetOuterInfo(elementName, element, centralInputFaceAngle, centralOutputFaceAngle, bpInfo, yokeOnLeft);
729 mgInfo->name = elementName;
730 if (!zeroStrength)
731 {
732 vacuumField = new BDSFieldInfo(dipoleFieldType,
733 brho,
734 intType,
735 st,
736 true,
737 fieldTiltOffset);
739 BDSFieldType::dipole,
740 bpInfo,
741 mgInfo,
742 fieldTiltOffset,
743 integratorSet,
744 brho,
745 BDSComponentFactory::ScalingFieldOuter(element));
746 }
747
748 // Here we change from the strength angle convention of +ve angle corresponds to
749 // deflection in negative x, to correct 3d +ve angle corresponds to deflection in
750 // positive x. Hence, angle sign flip for construction.
751 BDSMagnet* oneBend = new BDSMagnet(BDSMagnetType::rectangularbend,
752 elementName+"_centre",
753 centralArcLength,
754 bpInfo,
755 mgInfo,
756 vacuumField,
757 -centralAngle,
758 outerField);
759
760 rbendline->AddComponent(oneBend);
761
762 //Last element should be fringe if poleface specified
763 if (buildFringeOutgoing)
764 {
765 G4double trackingPolefaceAngle = trackingPolefaceAngleOut;
766 // if using a non-matrix integrator set, this code should only be reached if a finite fintx is specified.
767 // As the angled face geometry is constructed, a thin pole face kick shouldn't also be applied, so subtract
768 // e2 for in the magnet strength object so only the fringe & the natural rbend angled face effects are applied
769 if (!integratorSet->IsMatrixIntegratorSet() && BDSGlobalConstants::Instance()->BuildPoleFaceGeometry()
770 && BDS::IsFinite(fintx))
771 {trackingPolefaceAngle -= element->e2;}
772 BDSMagnetStrength* fringeStOut = BDS::GetFringeMagnetStrength(element, st, oneFringeAngle,
773 trackingPolefaceAngleIn, trackingPolefaceAngle,
774 fintx, false);
775 G4String fringeName = name + "_e2_fringe";
776
777 BDSMagnet* endfringe = BDS::BuildDipoleFringe(element, fringeOutInputAngle, angleOut,
778 fringeName,
779 fringeStOut, brho,
780 integratorSet, dipoleFieldType);
781 rbendline->AddComponent(endfringe);
782 }
783
784 return rbendline;
785}
786
788 G4double angleIn,
789 G4double angleOut,
790 const G4String& name,
792 G4double brho,
793 const BDSIntegratorSet* integratorSet,
794 BDSFieldType dipoleFieldType)
795{
796 BDSBeamPipeInfo* beamPipeInfo = BDSComponentFactory::PrepareBeamPipeInfo(element, angleIn, angleOut);
797 beamPipeInfo->beamPipeType = BDSBeamPipeType::circularvacuum;
798 auto magnetOuterInfo = BDSComponentFactory::PrepareMagnetOuterInfo(name, element,
799 angleIn, angleOut, beamPipeInfo);
800 magnetOuterInfo->geometryType = BDSMagnetGeometryType::none;
801 magnetOuterInfo->name = name;
802 magnetOuterInfo->buildEndPieces = false;
803
804 G4Transform3D fieldTiltOffset = BDSComponentFactory::CreateFieldTransform(element);
805
806 BDSIntegratorType intType = integratorSet->dipoleFringe;
807 BDSFieldInfo* vacuumField = new BDSFieldInfo(dipoleFieldType,
808 brho,
809 intType,
810 st,
811 true,
812 fieldTiltOffset);
813
814 return new BDSMagnet(BDSMagnetType::dipolefringe,
815 name,
816 (*st)["length"],
817 beamPipeInfo,
818 magnetOuterInfo,
819 vacuumField,
820 -(*st)["angle"],
821 nullptr,
822 true);
823}
824
825G4int BDS::CalculateNSBendSegments(G4double length,
826 G4double angle,
827 G4double incomingFaceAngle,
828 G4double outgoingFaceAngle,
829 G4double aperturePrecision)
830{
831 // Split a bend into equal segments such that the maximum distance between the
832 // chord and arc is 1mm.
833
834 // from formula: L/2 / N tan (angle/N) < precision. (L=physical length)
835 // add poleface rotations onto angle as absolute number (just to be safe)
836 G4double totalAngle = std::abs(angle) + std::abs(incomingFaceAngle) + std::abs(outgoingFaceAngle);
837 G4int nSBends = (G4int) ceil(std::sqrt(length*totalAngle/2/aperturePrecision));
838 if (nSBends==0)
839 {nSBends = 1;} // can happen in case angle = 0
840 if (BDSGlobalConstants::Instance()->DontSplitSBends())
841 {nSBends = 1;} // use for debugging
842 if (nSBends % 2 == 0)
843 {nSBends += 1;} // always have odd number of poles for poleface rotations
844#ifdef BDSDEBUG
845 G4cout << __METHOD_NAME__ << " splitting sbend into " << nSBends << " sbends" << G4endl;
846#endif
847 return nSBends;
848}
849
851 const Element* element)
852{
853 // default is dipole integrator
854 BDSIntegratorType intType = integratorSet->Integrator(BDSFieldType::dipole);
855
856 // check for finite k1 and change integrator type if needed
857 if (BDS::IsFinite(element->k1))
858 {intType = integratorSet->Integrator(BDSFieldType::dipolequadrupole);}
859
860 //if (BDS::IsFinite(element->tilt))
861 // {intType = BDSIntegratorType::dipolerodrigues2;}
862
863 return intType;
864}
865
866BDSMagnetStrength* BDS::GetFringeMagnetStrength(const Element* element,
867 const BDSMagnetStrength* st,
868 G4double fringeAngle,
869 G4double e1,
870 G4double e2,
871 G4double fintx,
872 G4bool isEntrance)
873{
874 BDSMagnetStrength* fringeSt = new BDSMagnetStrength(*st);
875 (*fringeSt)["length"] = BDSGlobalConstants::Instance()->ThinElementLength();
876 (*fringeSt)["angle"] = fringeAngle;
877 (*fringeSt)["e1"] = e1; // supply separately as it may be modified for rbends
878 (*fringeSt)["e2"] = e2; // supply separately as it may be modified for rbends
879 (*fringeSt)["fint"] = element->fint;
880 (*fringeSt)["fintx"] = fintx; // supply separately as it may be modified to match madx behaviour
881 (*fringeSt)["fintk2"] = element->fintK2;
882 (*fringeSt)["fintk2"] = element->fintxK2;
883 (*fringeSt)["hgap"] = element->hgap * CLHEP::m;
884 (*fringeSt)["isentrance"] = isEntrance;
885 (*fringeSt)["h1"] = element->h1 / CLHEP::m;
886 (*fringeSt)["h2"] = element->h2 / CLHEP::m;
887 return fringeSt;
888}
Abstract class that represents a component of an accelerator.
Holder class for all information required to describe a beam pipe model.
BDSBeamPipeType beamPipeType
Public member for direct access.
static G4bool YokeOnLeft(const GMAD::Element *el, const BDSMagnetStrength *st)
static G4double PrepareHorizontalWidth(GMAD::Element const *el, G4double defaultHorizontalWidth=-1)
Prepare the element horizontal width in Geant4 units - if not set, use the global default.
static G4Transform3D CreateFieldTransform(const BDSTiltOffset *tiltOffset)
Create a transform from a tilt offset. If nullptr, returns identity transform.
static BDSMagnetOuterInfo * PrepareMagnetOuterInfo(const G4String &elementNameIn, const GMAD::Element *el, const BDSMagnetStrength *st, const BDSBeamPipeInfo *beamPipe, G4double defaultHorizontalWidth=-1, G4double defaultVHRatio=1.0, G4double defaultCoilWidthFraction=-1, G4double defaultCoilHeightFraction=-1)
static void CheckBendLengthAngleWidthCombo(G4double arcLength, G4double angle, G4double horizontalWidth, const G4String &name="not given")
static BDSFieldInfo * PrepareMagnetOuterFieldInfo(const BDSMagnetStrength *vacuumSt, const BDSFieldType &fieldType, const BDSBeamPipeInfo *bpInfo, const BDSMagnetOuterInfo *outerInfo, const G4Transform3D &fieldTransform, const BDSIntegratorSet *integratorSetIn, G4double brhoIn, G4double outerFieldScaling=1.0)
Prepare the field definition for the yoke of a magnet.
static BDSBeamPipeInfo * PrepareBeamPipeInfo(GMAD::Element const *el, const G4ThreeVector &inputFaceNormal=G4ThreeVector(0, 0,-1), const G4ThreeVector &outputFaceNormal=G4ThreeVector(0, 0, 1))
General exception with possible name of object and message.
Definition: BDSException.hh:35
All info required to build complete field of any type.
Definition: BDSFieldInfo.hh:65
static BDSGlobalConstants * Instance()
Access method.
Which integrator to use for each type of magnet / field object.
G4bool IsMatrixIntegratorSet() const
Accessor for bool of is the integrator set matrix style.
BDSIntegratorType Integrator(const BDSFieldType field) const
Get appropriate integrator based on the field type.
A class that hold multiple accelerator components.
Definition: BDSLine.hh:38
void AddComponent(BDSAcceleratorComponent *component)
Add a component to the line.
Definition: BDSLine.cc:28
Efficient storage of magnet strengths.
Abstract base class that implements features common to all magnets.
Definition: BDSMagnet.hh:45
G4bool ZeroStrengthDipole(const BDSMagnetStrength *st)
Return whether finite angle or field for a dipole.
G4int CalculateNSBendSegments(G4double length, G4double angle, G4double incomingFaceAngle=0, G4double outgoingFaceAngle=0, G4double aperturePrecision=1.0)
BDSMagnet * BuildDipoleFringe(const GMAD::Element *element, G4double angleIn, G4double angleOut, const G4String &name, BDSMagnetStrength *st, G4double brho, const BDSIntegratorSet *integratorSet, BDSFieldType dipoleFieldType)
BDSLine * BuildRBendLine(const G4String &elementName, const GMAD::Element *element, const GMAD::Element *prevElement, const GMAD::Element *nextElement, G4double brho, BDSMagnetStrength *st, const BDSIntegratorSet *integratorSet, G4double incomingFaceAngle, G4double outgoingFaceAngle, G4bool buildFringeFields)
BDSAcceleratorComponent * BuildSBendLine(const G4String &elementName, const GMAD::Element *element, BDSMagnetStrength *st, G4double brho, const BDSIntegratorSet *integratorSet, G4double incomingFaceAngle, G4double outgoingFaceAngle, G4bool buildFringeFields, const GMAD::Element *prevElement, const GMAD::Element *nextElement)
BDSIntegratorType GetDipoleIntegratorType(const BDSIntegratorSet *integratorSet, const GMAD::Element *element)
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())
BDSMagnet * BuildSingleSBend(const GMAD::Element *element, const G4String &name, G4double arcLength, G4double angle, G4double angleIn, G4double angleOut, const BDSMagnetStrength *strength, G4double brho, const BDSIntegratorSet *integratorSet, G4bool yokeOnLeft, const BDSFieldInfo *outerFieldIn)
Function to return a single sector bend section.
Parser namespace for GMAD language. Combination of Geant4 and MAD.
Element class.
Definition: element.h:43
double hgap
half distance of pole separation for purposes of fringe fields - 'half gap'
Definition: element.h:66
double fintx
fringe field integral at the dipole exit
Definition: element.h:63
double h2
output pole face curvature for bends
Definition: element.h:68
double h1
input pole face curvature for bends
Definition: element.h:67
double fintxK2
second fringe field integral at the dipole exit - for TRANSPORT matching
Definition: element.h:65
double fint
fringe field integral at the dipole entrance
Definition: element.h:62
double e2
output pole face rotation for bends
Definition: element.h:61
double l
length in metres
Definition: element.h:49
double e1
input pole face rotation for bends
Definition: element.h:60
double k1
quadrupole
Definition: element.h:54
double fintK2
second fringe field integral at the dipole entrance - for TRANSPORT matching
Definition: element.h:64
ElementType type
element enum
Definition: element.h:44