BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSFieldFactory.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 "BDSArrayReflectionType.hh"
20#include "BDSBeamPipeInfo.hh"
21#include "BDSDebug.hh"
22#include "BDSException.hh"
23#include "BDSPTCOneTurnMap.hh"
24#include "BDSPrimaryGeneratorAction.hh"
25#include "BDSFieldClassType.hh"
26#include "BDSFieldE.hh"
27#include "BDSFieldEGlobal.hh"
28#include "BDSFieldEGlobalPlacement.hh"
29#include "BDSFieldEInterpolated.hh"
30#include "BDSFieldEInterpolated2Layer.hh"
31#include "BDSFieldESinusoid.hh"
32#include "BDSFieldEZero.hh"
33#include "BDSFieldEM.hh"
34#include "BDSFieldEMGlobal.hh"
35#include "BDSFieldEMGlobalPlacement.hh"
36#include "BDSFieldEMInterpolated.hh"
37#include "BDSFieldEMRFCavity.hh"
38#include "BDSFieldEMZero.hh"
39#include "BDSFieldFactory.hh"
40#include "BDSFieldInfo.hh"
41#include "BDSFieldLoader.hh"
42#include "BDSFieldMag.hh"
43#include "BDSFieldMagDecapole.hh"
44#include "BDSFieldMagDipole.hh"
45#include "BDSFieldMagDipoleOuter.hh"
46#include "BDSFieldMagDipoleOuterOld.hh"
47#include "BDSFieldMagDipoleQuadrupole.hh"
48#include "BDSFieldMagGlobal.hh"
49#include "BDSFieldMagGlobalPlacement.hh"
50#include "BDSFieldMagInterpolated.hh"
51#include "BDSFieldMagInterpolated2Layer.hh"
52#include "BDSFieldMagMultipole.hh"
53#include "BDSFieldMagMultipoleOuter.hh"
54#include "BDSFieldMagMultipoleOuterDual.hh"
55#include "BDSFieldMagMultipoleOuterDualOld.hh"
56#include "BDSFieldMagMultipoleOuterOld.hh"
57#include "BDSFieldMagMuonSpoiler.hh"
58#include "BDSFieldMagOctupole.hh"
59#include "BDSFieldMagQuadrupole.hh"
60#include "BDSFieldMagSextupole.hh"
61#include "BDSFieldMagSolenoidSheet.hh"
62#include "BDSFieldMagSkewOwn.hh"
63#include "BDSFieldMagZero.hh"
64#include "BDSFieldObjects.hh"
65#include "BDSFieldType.hh"
66#include "BDSGlobalConstants.hh"
67#include "BDSIntegratorCavityFringe.hh"
68#include "BDSIntegratorDecapole.hh"
69#include "BDSIntegratorDipoleRodrigues.hh"
70#include "BDSIntegratorDipoleRodrigues2.hh"
71#include "BDSIntegratorDipoleFringe.hh"
72#include "BDSIntegratorDipoleFringeScaling.hh"
73#include "BDSIntegratorDipoleQuadrupole.hh"
74#include "BDSIntegratorEuler.hh"
75#include "BDSIntegratorG4RK4MinStep.hh"
76#include "BDSIntegratorKickerThin.hh"
77#include "BDSIntegratorOctupole.hh"
78#include "BDSIntegratorQuadrupole.hh"
79#include "BDSIntegratorMultipoleThin.hh"
80#include "BDSIntegratorParallelTransport.hh"
81#include "BDSIntegratorSextupole.hh"
82#include "BDSIntegratorSolenoid.hh"
83#include "BDSIntegratorTeleporter.hh"
84#include "BDSIntegratorRMatrixThin.hh"
85#include "BDSIntegratorType.hh"
86#include "BDSMagnetOuterFactoryLHC.hh"
87#include "BDSMagnetStrength.hh"
88#include "BDSMagnetType.hh"
89#include "BDSParser.hh"
90#include "BDSParticleDefinition.hh"
91#include "BDSUtilities.hh"
92#include "BDSFieldMagUndulator.hh"
93
94#include "parser/field.h"
95
96#include "globals.hh" // geant4 types / globals
97#include "G4EquationOfMotion.hh"
98#include "G4EqMagElectricField.hh"
99#include "G4MagIntegratorStepper.hh"
100#include "G4Mag_UsualEqRhs.hh"
101#include "G4RotationMatrix.hh"
102#include "G4String.hh"
103#include "G4ThreeVector.hh"
104#include "G4Transform3D.hh"
105#include "G4Version.hh"
106
107// geant4 integrators
108#include "G4CashKarpRKF45.hh"
109#include "G4ClassicalRK4.hh"
110#include "G4ConstRK4.hh"
111#include "G4ExactHelixStepper.hh"
112#include "G4ExplicitEuler.hh"
113#include "G4HelixExplicitEuler.hh"
114#include "G4HelixHeum.hh"
115#include "G4HelixImplicitEuler.hh"
116#include "G4HelixMixedStepper.hh"
117#include "G4HelixSimpleRunge.hh"
118#include "G4ImplicitEuler.hh"
119#include "G4NystromRK4.hh"
120#include "G4RKG3_Stepper.hh"
121#include "G4SimpleHeum.hh"
122#include "G4SimpleRunge.hh"
123#if G4VERSION_NUMBER > 1029
124#include "G4BogackiShampine23.hh"
125#include "G4BogackiShampine45.hh"
126#include "G4DoLoMcPriRK34.hh"
127#include "G4DormandPrince745.hh"
128#include "G4DormandPrinceRK56.hh"
129#include "G4TsitourasRK45.hh"
130#endif
131#if G4VERSION_NUMBER > 1039
132#include "G4DormandPrinceRK78.hh"
133#include "G4RK547FEq1.hh"
134#include "G4RK547FEq2.hh"
135#include "G4RK547FEq3.hh"
136#endif
137
138#include "CLHEP/Units/SystemOfUnits.h"
139#include "CLHEP/Vector/EulerAngles.h"
140
141#include <limits>
142#include <map>
143#include <utility>
144#include <vector>
145
148
150
152{
153 if (!instance)
154 {instance = new BDSFieldFactory();}
155 return instance;
156}
157
159 useOldMultipoleOuterFields(false)
160{
161 G4double defaultRigidity = std::numeric_limits<double>::max();
162 if (designParticle)
163 {defaultRigidity = designParticle->BRho();}
164 PrepareFieldDefinitions(BDSParser::Instance()->GetFields(), defaultRigidity);
165 useOldMultipoleOuterFields = BDSGlobalConstants::Instance()->UseOldMultipoleOuterFields();
166}
167
168BDSFieldFactory::~BDSFieldFactory()
169{
170 for (auto& info : parserDefinitions)
171 {delete info.second;}
172}
173
174void BDSFieldFactory::PrepareFieldDefinitions(const std::vector<GMAD::Field>& definitions,
175 G4double defaultBRho)
176{
177 for (const auto& definition : definitions)
178 {
179 if (definition.type.empty())
180 {
181 G4String msg = "\"type\" not specified in field definition \"";
182 msg += definition.name + "\", but required.";
183 throw BDSException(__METHOD_NAME__, msg);
184 }
185 BDSFieldType fieldType = BDS::DetermineFieldType(definition.type);
186 BDSIntegratorType intType = BDS::DetermineIntegratorType(definition.integrator);
187
188 G4ThreeVector offset = G4ThreeVector(definition.x*CLHEP::m,
189 definition.y*CLHEP::m,
190 definition.z*CLHEP::m);
191
192 G4RotationMatrix rm;
193 if (definition.axisAngle)
194 {
195 G4ThreeVector axis = G4ThreeVector(definition.axisX,
196 definition.axisY,
197 definition.axisZ);
198 rm = G4RotationMatrix(axis, definition.angle*CLHEP::rad);
199 }
200 else
201 {
202 if (BDS::IsFinite(definition.phi) ||
203 BDS::IsFinite(definition.theta) ||
204 BDS::IsFinite(definition.psi))
205 {// only build if finite
206 CLHEP::HepEulerAngles ang = CLHEP::HepEulerAngles(definition.phi*CLHEP::rad,
207 definition.theta*CLHEP::rad,
208 definition.psi*CLHEP::rad);
209 rm = G4RotationMatrix(ang);
210 }
211 // else rm is default rotation matrix
212 }
213
214 G4Transform3D transform = G4Transform3D(rm, offset);
215
216 BDSFieldFormat magFormat = BDSFieldFormat::none;
217 G4String magFile = "";
218 G4bool magFileSpecified = !definition.magneticFile.empty();
219 if (magFileSpecified)
220 {
221 std::pair<G4String, G4String> bf = BDS::SplitOnColon(G4String(definition.magneticFile));
222 magFormat = BDS::DetermineFieldFormat(bf.first);
223 magFile = BDS::GetFullPath(bf.second);
224 }
225
226 BDSFieldFormat eleFormat = BDSFieldFormat::none;
227 G4String eleFile = "";
228 G4bool eleFileSpecified = !definition.electricFile.empty();
229 if (eleFileSpecified)
230 {
231 std::pair<G4String, G4String> ef = BDS::SplitOnColon(G4String(definition.electricFile));
232 eleFormat = BDS::DetermineFieldFormat(ef.first);
233 eleFile = BDS::GetFullPath(ef.second);
234 }
235
236 BDSInterpolatorType magIntType = BDSInterpolatorType::cubic3d;
237 if (magFileSpecified)
238 {// determine and check type of integrator
239 G4int nDimFF = BDS::NDimensionsOfFieldFormat(magFormat);
240 if (!definition.magneticInterpolator.empty())
241 {
242 magIntType = BDS::DetermineInterpolatorType(G4String(definition.magneticInterpolator));
243 // detect if an auto type and match up accordingly, else check it's the right one
244 if (BDS::InterpolatorTypeIsAuto(magIntType))
245 {magIntType = BDS::InterpolatorTypeSpecificFromAuto(nDimFF, magIntType);}
246 else
247 {
248 G4int nDimInt = BDS::NDimensionsOfInterpolatorType(magIntType);
249 if (nDimFF != nDimInt)
250 {
251 G4String message = "mismatch in number of dimensions between magnetic interpolator ";
252 message += "and field map format for field definition \"" + definition.name + "\"";
253 throw BDSException(__METHOD_NAME__, message);
254 }
255 }
256 }
257 else
258 {magIntType = DefaultInterpolatorType(nDimFF);}
259 }
260
261 BDSInterpolatorType eleIntType = BDSInterpolatorType::cubic3d;
262 if (eleFileSpecified)
263 {// determine and check type of integrator
264 G4int nDimFF = BDS::NDimensionsOfFieldFormat(eleFormat);
265 if (!definition.electricInterpolator.empty())
266 {
267 eleIntType = BDS::DetermineInterpolatorType(G4String(definition.electricInterpolator));
268 // detect if an auto type and match up accordingly, else check it's the right one
269 if (BDS::InterpolatorTypeIsAuto(eleIntType))
270 {eleIntType = BDS::InterpolatorTypeSpecificFromAuto(nDimFF, eleIntType);}
271 else
272 {
273 G4int nDimInt = BDS::NDimensionsOfInterpolatorType(eleIntType);
274 if (nDimFF != nDimInt)
275 {
276 G4String message = "mismatch in number of dimensions between electric interpolator ";
277 message += "and field map format for field definition \"" + definition.name + "\"";
278 throw BDSException(__METHOD_NAME__, message);
279 }
280 }
281 }
282 else
283 {eleIntType = DefaultInterpolatorType(nDimFF);}
284 }
285
286 G4UserLimits* fieldLimit = nullptr;
287 if (definition.maximumStepLength > 0)
288 {// only assign if specified
289 auto defaultUL = BDSGlobalConstants::Instance()->DefaultUserLimits();
290 // copy the default and update with the length of the object rather than the default 1m
291 G4double limit = G4double(definition.maximumStepLength) * CLHEP::m;
292 G4UserLimits* ul = BDS::CreateUserLimits(defaultUL, limit, 1.0);
293 // only specify a user limit object if the step length was specified
294 if (ul != defaultUL)
295 {fieldLimit = ul;}
296 }
297
299 G4double poleTipRadius = BDSGlobalConstants::Instance()->DefaultBeamPipeModel()->aper1;
300 if (!definition.fieldParameters.empty())
301 {PrepareFieldStrengthFromParameters(st, definition.fieldParameters, poleTipRadius);}
302
303 BDSFieldInfo* info = new BDSFieldInfo(fieldType,
304 defaultBRho,
305 intType,
306 st,
307 G4bool(definition.globalTransform),
308 transform,
309 magFile,
310 magFormat,
311 magIntType,
312 eleFile,
313 eleFormat,
314 eleIntType,
315 false, /*don't cache transforms*/
316 G4double(definition.eScaling),
317 G4double(definition.bScaling),
318 G4double(definition.t*CLHEP::s),
319 G4bool(definition.autoScale),
320 fieldLimit);
321 info->SetScalingRadius(poleTipRadius);
322
323 if (!definition.magneticSubField.empty())
324 {
325 if (definition.magneticSubField == definition.name)
326 {throw BDSException(__METHOD_NAME__, "error in \"" + definition.name + "\": magneticSubField cannot be the field itself");}
327 info->SetMagneticSubField(G4String(definition.magneticSubField));
328 }
329 if (!definition.electricSubField.empty())
330 {
331 if (definition.electricSubField == definition.name)
332 {throw BDSException(__METHOD_NAME__, "error in \"" + definition.name + "\": electricSubField cannot be the field itself");}
333 info->SetElectricSubField(G4String(definition.electricSubField));
334 }
335 if (!definition.magneticReflection.empty())
336 {
337 G4String magneticReflection = G4String(definition.magneticReflection);
338 BDSArrayReflectionTypeSet mar = BDS::DetermineArrayReflectionTypeSet(magneticReflection);
339 info->SetMagneticArrayReflectionType(mar);
340 }
341 if (!definition.electricReflection.empty())
342 {
343 G4String electricReflection = G4String(definition.electricReflection);
344 BDSArrayReflectionTypeSet ear = BDS::DetermineArrayReflectionTypeSet(electricReflection);
345 info->SetElectricArrayReflectionType(ear);
346 }
347
348 info->SetNameOfParserDefinition(G4String(definition.name));
349 if (BDSGlobalConstants::Instance()->Verbose())
350 {
351 G4cout << "Definition: \"" << definition.name << "\"" << G4endl;
352 G4cout << *info << G4endl;
353 }
354 parserDefinitions[G4String(definition.name)] = info;
355 }
356}
357
359 const G4String& parameterNameForError) const
360{
361 G4double result = 0;
362 try
363 {result = std::stod(value);}
364 catch (std::exception& e)
365 {
366 G4String msg(e.what());
367 G4String baseMsg = "Unable to interpret value (\"" + value + "\" of parameter \"";
368 baseMsg += parameterNameForError + "\" as a number: ";
369 throw BDSException(__METHOD_NAME__, baseMsg + msg);
370 }
371 return result;
372}
373
375 const G4String& fieldParameters,
376 G4double& poleTipRadius) const
377{
378 // use function from BDSUtilities to process user params string into
379 // map of strings and values.
380 std::map<G4String, G4String> map = BDS::GetUserParametersMap(fieldParameters, '=');
381
382 for (const auto& keyValue : map)
383 {
384 if (BDSMagnetStrength::ValidKey(keyValue.first))
385 {
386 G4double value = ConvertToDoubleWithException(keyValue.second, keyValue.first);
387 (*st)[keyValue.first] = value * BDSMagnetStrength::Unit(keyValue.first);
388 }
389 else if (keyValue.first == "poletipradius")
390 {poleTipRadius = ConvertToDoubleWithException(keyValue.second, keyValue.first) * CLHEP::m;}
391 else
392 {
393 G4String msg = "Invalid key \"" + keyValue.first + "\" for field parameters. ";
394 msg += "Acceptable parameters are: \n";
395 const std::vector<G4String>& allKeys = BDSMagnetStrength::AllKeys();
396 for (G4int i = 0; i < (G4int)allKeys.size(); i++)
397 {
398 msg += allKeys[i] + ", ";
399 if ((i % 10 == 0) && (i > 0))
400 {msg += "\n";}
401 }
402 throw BDSException(__METHOD_NAME__, msg);
403 }
404 }
405}
406
408{
409 // Here we test if the string is empty and return nullptr. We do this so
410 // this method can be used without exiting when no key is specified at all.
411 // If a key is given and not found, then that requires the users attention.
412 if (name.empty())
413 {return nullptr;}
414 auto result = parserDefinitions.find(name);
415 if (result == parserDefinitions.end())
416 {// not a valid key
417 G4cerr << __METHOD_NAME__ << "\"" << name << "\" is not a valid field specifier" << G4endl;
418 G4cout << "Defined field specifiers are:" << G4endl;
419 for (const auto& it : parserDefinitions)
420 {G4cout << "\"" << it.first << "\"" << G4endl;}
421 throw BDSException(__METHOD_NAME__, "invalid field name");
422 }
423 return result->second;
424}
425
427 const BDSMagnetStrength* scalingStrength,
428 const G4String& scalingKey)
429{
430#ifdef BDSDEBUG
431 G4cout << __METHOD_NAME__ << info << G4endl;
432#endif
433 // Forward on to delegate functions for the main types of field
434 // such as E, EM and Magnetic
435 BDSFieldObjects* field = nullptr;
436
437 if (info.FieldType() == BDSFieldType::none)
438 {return field;} // as nullptr
439
441 switch (clas.underlying())
442 {
443 case BDSFieldClassType::magnetic:
444 {field = CreateFieldMag(info, scalingStrength, scalingKey); break;}
445 case BDSFieldClassType::electromagnetic:
446 {field = CreateFieldEM(info); break;}
447 case BDSFieldClassType::electric:
448 {field = CreateFieldE(info); break;}
449 case BDSFieldClassType::irregular:
450 {field = CreateFieldIrregular(info); break;}
451 default:
452 {break;} // this will return nullptr
453 }
454 return field;
455}
456
458{
459 BDSInterpolatorType result;
460 switch (numberOfDimensions)
461 {
462 case 1:
463 {result = BDSInterpolatorType::cubic1d; break;}
464 case 2:
465 {result = BDSInterpolatorType::cubic2d; break;}
466 case 3:
467 {result = BDSInterpolatorType::cubic3d; break;}
468 case 4:
469 {result = BDSInterpolatorType::cubic4d; break;}
470 default:
471 {throw BDSException(__METHOD_NAME__, "unsupported number of dimensions " + std::to_string(numberOfDimensions));}
472 }
473 return result;
474}
475
477 const BDSMagnetStrength* scalingStrength,
478 const G4String& scalingKey)
479{
480 const BDSMagnetStrength* strength = info.MagnetStrength();
481 BDSFieldMag* field = CreateFieldMagRaw(info, scalingStrength, scalingKey);
482 if (!field)
483 {return nullptr;} // return nullptr of right type
484
485 BDSFieldMag* resultantField = field;
486 // Optionally provide local to global transform using curvilinear coordinate system.
488 {resultantField = new BDSFieldMagGlobalPlacement(field);}
489 else if (info.ProvideGlobal())
490 {resultantField = new BDSFieldMagGlobal(field);}
491
492 // Always this equation of motion for magnetic (only) fields
493 BDSMagUsualEqRhs* eqOfM = new BDSMagUsualEqRhs(resultantField);
494
495 // Create appropriate integrator
496 G4MagIntegratorStepper* integrator = CreateIntegratorMag(info, eqOfM, strength);
497
498 BDSFieldObjects* completeField = new BDSFieldObjects(&info, resultantField, eqOfM, integrator);
499 return completeField;
500}
501
503 const BDSMagnetStrength* scalingStrength,
504 const G4String& scalingKey)
505{
506 BDSFieldMag* field = nullptr;
507 const BDSMagnetStrength* strength = info.MagnetStrength();
508 G4double brho = info.BRho();
509 G4double poleTipRadius = info.PoleTipRadius();
510 G4double beamPipeRadius = info.BeamPipeRadius();
511 switch (info.FieldType().underlying())
512 {
513 case BDSFieldType::bmap1d:
514 case BDSFieldType::bmap2d:
515 case BDSFieldType::bmap3d:
516 case BDSFieldType::bmap4d:
517 case BDSFieldType::mokka:
518 {
520 scalingStrength,
521 scalingKey);
522 if (ff)
523 {info.UpdateUserLimitsLengthMaximumStepSize(ff->SmallestSpatialStep(), true);}
524 field = ff;
525 break;
526 }
527 case BDSFieldType::bfieldzero:
528 {field = new BDSFieldMagZero(); break;}
529 case BDSFieldType::solenoid:
530 case BDSFieldType::dipole:
531 case BDSFieldType::dipole3d:
532 {field = new BDSFieldMagDipole(strength); break;}
533 case BDSFieldType::solenoidsheet:
534 {field = new BDSFieldMagSolenoidSheet(strength, poleTipRadius); break;}
535 case BDSFieldType::quadrupole:
536 {field = new BDSFieldMagQuadrupole(strength, brho); break;}
537 case BDSFieldType::undulator:
538 {field = new BDSFieldMagUndulator(strength, beamPipeRadius); break;}
539 case BDSFieldType::dipolequadrupole:
540 {field = new BDSFieldMagDipoleQuadrupole(strength, brho); break;}
541 case BDSFieldType::sextupole:
542 {field = new BDSFieldMagSextupole(strength, brho); break;}
543 case BDSFieldType::octupole:
544 {field = new BDSFieldMagOctupole(strength, brho); break;}
545 case BDSFieldType::decapole:
546 {field = new BDSFieldMagDecapole(strength, brho); break;}
547 case BDSFieldType::multipole:
548 {field = new BDSFieldMagMultipole(strength, brho); break;}
549 case BDSFieldType::muonspoiler:
550 {field = new BDSFieldMagMuonSpoiler(strength, brho); break;}
551 case BDSFieldType::skewquadrupole:
552 {field = new BDSFieldMagSkewOwn(new BDSFieldMagQuadrupole(strength, brho), CLHEP::pi/4.); break;}
553 case BDSFieldType::skewsextupole:
554 {field = new BDSFieldMagSkewOwn(new BDSFieldMagSextupole(strength, brho), CLHEP::pi/6.); break;}
555 case BDSFieldType::skewoctupole:
556 {field = new BDSFieldMagSkewOwn(new BDSFieldMagOctupole(strength, brho), CLHEP::pi/8.); break;}
557 case BDSFieldType::skewdecapole:
558 {field = new BDSFieldMagSkewOwn(new BDSFieldMagDecapole(strength, brho), CLHEP::pi/10.); break;}
559 case BDSFieldType::multipoleouterdipole:
560 {// suitable only for querying transversely in x,y - no 3d nature
561 BDSFieldMag* innerField = new BDSFieldMagDipole(strength);
562 G4bool positiveField = (*strength)["field"] < 0; // convention for dipoles - "positive"
563 if (useOldMultipoleOuterFields)
564 {field = new BDSFieldMagMultipoleOuterOld(1, poleTipRadius, innerField, positiveField, GetOuterScaling(strength));}
565 else // positive field for a dipole has no meaning for this class - fix at false
566 {field = new BDSFieldMagMultipoleOuter(1, poleTipRadius, innerField, false, brho, GetOuterScaling(strength));}
567 delete innerField; // no longer required
568 break;
569 }
570 case BDSFieldType::multipoleouterquadrupole:
571 {
572 BDSFieldMag* innerField = new BDSFieldMagQuadrupole(strength, brho);
573 G4bool positiveField = (*strength)["k1"] > 0;
574 if (useOldMultipoleOuterFields)
575 {field = new BDSFieldMagMultipoleOuterOld(2, poleTipRadius, innerField, positiveField, GetOuterScaling(strength));}
576 else
577 {field = new BDSFieldMagMultipoleOuter(2, poleTipRadius, innerField, positiveField, brho, GetOuterScaling(strength));}
578 delete innerField; // no longer required
579 break;
580 }
581 case BDSFieldType::multipoleoutersextupole:
582 {
583 BDSFieldMag* innerField = new BDSFieldMagSextupole(strength, brho);
584 G4bool positiveField = (*strength)["k2"] > 0;
585 if (useOldMultipoleOuterFields)
586 {field = new BDSFieldMagMultipoleOuterOld(3, poleTipRadius, innerField, positiveField, GetOuterScaling(strength));}
587 else
588 {field = new BDSFieldMagMultipoleOuter(3, poleTipRadius, innerField, positiveField, brho, GetOuterScaling(strength));}
589 delete innerField; // no longer required
590 break;
591 }
592 case BDSFieldType::multipoleouteroctupole:
593 {
594 BDSFieldMag* innerField = new BDSFieldMagOctupole(strength, brho);
595 G4bool positiveField = (*strength)["k3"] > 0;
596 if (useOldMultipoleOuterFields)
597 {field = new BDSFieldMagMultipoleOuterOld(4, poleTipRadius, innerField, positiveField, GetOuterScaling(strength));}
598 else
599 {field = new BDSFieldMagMultipoleOuter(4, poleTipRadius, innerField, positiveField, brho, GetOuterScaling(strength));}
600 delete innerField; // no longer required
601 break;
602 }
603 case BDSFieldType::multipoleouterdecapole:
604 {
605 BDSFieldMag* innerField = new BDSFieldMagDecapole(strength, brho);
606 G4bool positiveField = (*strength)["k4"] > 0;
607 if (useOldMultipoleOuterFields)
608 {field = new BDSFieldMagMultipoleOuterOld(5, poleTipRadius, innerField, positiveField, GetOuterScaling(strength));}
609 else
610 {field = new BDSFieldMagMultipoleOuter(5, poleTipRadius, innerField, positiveField, brho, GetOuterScaling(strength));}
611 delete innerField; // no longer required
612 break;
613 }
614 case BDSFieldType::skewmultipoleouterquadrupole:
615 {
616 BDSFieldMag* innerField = new BDSFieldMagQuadrupole(strength, brho);
617 G4bool positiveField = (*strength)["k1"] > 0;
618 BDSFieldMag* normalField;
619 if (useOldMultipoleOuterFields)
620 {normalField = new BDSFieldMagMultipoleOuterOld(2, poleTipRadius, innerField, positiveField, GetOuterScaling(strength));}
621 else
622 {normalField = new BDSFieldMagMultipoleOuter(2, poleTipRadius, innerField, positiveField, brho, GetOuterScaling(strength));}
623 field = new BDSFieldMagSkewOwn(normalField, CLHEP::pi/4.);
624 delete innerField; // no longer required
625 break;
626 }
627 case BDSFieldType::skewmultipoleoutersextupole:
628 {
629 BDSFieldMag* innerField = new BDSFieldMagSextupole(strength, brho);
630 G4bool positiveField = (*strength)["k2"] > 0;
631 BDSFieldMag* normalField;
632 if (useOldMultipoleOuterFields)
633 {normalField = new BDSFieldMagMultipoleOuterOld(3, poleTipRadius, innerField, positiveField, GetOuterScaling(strength));}
634 else
635 {normalField = new BDSFieldMagMultipoleOuter(3, poleTipRadius, innerField, positiveField, brho, GetOuterScaling(strength));}
636 field = new BDSFieldMagSkewOwn(normalField, CLHEP::pi/6.);
637 delete innerField; // no longer required
638 break;
639 }
640 case BDSFieldType::skewmultipoleouteroctupole:
641 {
642 BDSFieldMag* innerField = new BDSFieldMagOctupole(strength, brho);
643 G4bool positiveField = (*strength)["k3"] > 0;
644 BDSFieldMag* normalField;
645 if (useOldMultipoleOuterFields)
646 {normalField = new BDSFieldMagMultipoleOuterOld(4, poleTipRadius, innerField, positiveField, GetOuterScaling(strength));}
647 else
648 {normalField = new BDSFieldMagMultipoleOuter(4, poleTipRadius, innerField, positiveField, brho, GetOuterScaling(strength));}
649 field = new BDSFieldMagSkewOwn(normalField, CLHEP::pi/8.);
650 delete innerField; // no longer required
651 break;
652 }
653 case BDSFieldType::skewmultipoleouterdecapole:
654 {
655 BDSFieldMag* innerField = new BDSFieldMagDecapole(strength, brho);
656 G4bool positiveField = (*strength)["k4"] > 0;
657 BDSFieldMag* normalField;
658 if (useOldMultipoleOuterFields)
659 {normalField = new BDSFieldMagMultipoleOuterOld(5, poleTipRadius, innerField, positiveField, GetOuterScaling(strength));}
660 else
661 {normalField = new BDSFieldMagMultipoleOuter(5, poleTipRadius, innerField, positiveField, brho, GetOuterScaling(strength));}
662 field = new BDSFieldMagSkewOwn(normalField, CLHEP::pi/10.);
663 delete innerField; // no longer required
664 break;
665 }
666 case BDSFieldType::multipoleouterdipole3d:
667 {
668 if (useOldMultipoleOuterFields)
669 {field = new BDSFieldMagDipoleOuterOld(strength, poleTipRadius, GetOuterScaling(strength));}
670 else
671 {field = new BDSFieldMagDipoleOuter(strength, poleTipRadius, GetOuterScaling(strength));}
672 break;
673 }
674 case BDSFieldType::multipoleouterdipolelhc:
675 {
676 BDSFieldMag* innerField = new BDSFieldMagDipole(strength);
677 G4bool positiveField = (*strength)["field"] < 0; // convention for dipoles - "positive"
678 G4bool positiveField2 = (*strength)["angle"] > 0;
680 if (useOldMultipoleOuterFields)
681 {
682 field = new BDSFieldMagMultipoleOuterDualOld(1, poleTipRadius, innerField, positiveField, dx,
683 info.SecondFieldOnLeft(), GetOuterScaling(strength));
684 }
685 else
686 {
687 field = new BDSFieldMagMultipoleOuterDual(1, poleTipRadius, innerField, positiveField2, brho, dx,
688 info.SecondFieldOnLeft(), GetOuterScaling(strength));
689 }
690 delete innerField; // no longer required
691 break;
692 }
693 case BDSFieldType::multipoleouterquadrupolelhc:
694 {
695 BDSFieldMag* innerField = new BDSFieldMagQuadrupole(strength, brho);
696 G4bool positiveField = (*strength)["k1"] > 0;
698 if (useOldMultipoleOuterFields)
699 {
700 field = new BDSFieldMagMultipoleOuterDualOld(2, poleTipRadius, innerField, positiveField, dx,
701 info.SecondFieldOnLeft(), GetOuterScaling(strength));
702 }
703 else
704 {
705 field = new BDSFieldMagMultipoleOuterDual(2, poleTipRadius, innerField, positiveField, brho, dx,
706 info.SecondFieldOnLeft(), GetOuterScaling(strength));
707 }
708 delete innerField; // no longer required
709 break;
710 }
711 case BDSFieldType::multipoleoutersextupolelhc:
712 {
713 BDSFieldMag* innerField = new BDSFieldMagSextupole(strength, brho);
714 G4bool positiveField = (*strength)["k2"] > 0;
716 if (useOldMultipoleOuterFields)
717 {
718 field = new BDSFieldMagMultipoleOuterDualOld(3, poleTipRadius, innerField, positiveField, dx,
719 info.SecondFieldOnLeft(), GetOuterScaling(strength));
720 }
721 else
722 {
723 field = new BDSFieldMagMultipoleOuterDual(3, poleTipRadius, innerField, positiveField, dx, brho,
724 info.SecondFieldOnLeft(), GetOuterScaling(strength));
725 }
726 delete innerField; // no longer required
727 break;
728 }
729 case BDSFieldType::paralleltransporter:
730 default:
731 {// there is no need for case BDSFieldType::none as this won't be used in this function.
732 break;
733 }
734 }
735
736 // Set transform for local geometry offset
737 // Do this before wrapping in global converter BDSFieldMagGlobal so that the sub-field
738 // has it and not the global wrapper.
739 if (field)
740 {field->SetTransform(info.TransformComplete());}
741
742 if (!info.MagneticSubFieldName().empty() && field)
743 {
744 // set the transform of the 'main' field to only the transform defined in that field definition
745 field->SetTransform(info.Transform());
746
747 auto mainField = dynamic_cast<BDSFieldMagInterpolated*>(field);
748 if (!mainField)
749 {throw BDSException(__METHOD_NAME__, "subfield specified for non-field map type field - not supported");}
750
751 BDSFieldInfo* subFieldRecipe = new BDSFieldInfo(*(GetDefinition(info.MagneticSubFieldName())));
752 BDSFieldMag* subFieldRaw = CreateFieldMagRaw(*subFieldRecipe, scalingStrength, scalingKey);
753 auto subField = dynamic_cast<BDSFieldMagInterpolated*>(subFieldRaw);
754 if (!subField)
755 {throw BDSException(__METHOD_NAME__, "subfield type is not a field map type field - not supported");}
756 field = new BDSFieldMagInterpolated2Layer(mainField, subField);
757 // the transform goes beamline transform to the 2Layer class, then inside that the individual field transforms
758 field->SetTransform(info.TransformBeamline());
759 delete subFieldRecipe;
760 }
761
762 return field;
763}
764
766{
767 BDSFieldEM* field = nullptr;
768 switch (info.FieldType().underlying())
769 {
770 case BDSFieldType::rfcavity:
771 {field = new BDSFieldEMRFCavity(info.MagnetStrength(), info.BRho()); break;}
772 case BDSFieldType::ebmap1d:
773 case BDSFieldType::ebmap2d:
774 case BDSFieldType::ebmap3d:
775 case BDSFieldType::ebmap4d:
776 {
778 if (ff)
779 {info.UpdateUserLimitsLengthMaximumStepSize(ff->SmallestSpatialStep(), true);}
780 field = ff;
781 break;
782 }
783 case BDSFieldType::ebfieldzero:
784 {field = new BDSFieldEMZero(); break;}
785 default:
786 return nullptr;
787 break;
788 }
789
790 // Set transform for local geometry offset
791 if (field)
792 {field->SetTransform(info.TransformComplete());}
793
794 if (!field)
795 {return nullptr;}
796
797 // Optionally provide local to global transform using curvilinear coordinate system.
798 BDSFieldEM* resultantField = field;
800 {resultantField = new BDSFieldEMGlobalPlacement(field);}
801 else if (info.ProvideGlobal())
802 {resultantField = new BDSFieldEMGlobal(field);}
803
804 // Equation of motion for em fields
805 G4EqMagElectricField* eqOfM = new G4EqMagElectricField(resultantField);
806
807 // Create appropriate integrator
808 G4MagIntegratorStepper* integrator = CreateIntegratorEM(info, eqOfM);
809
810 BDSFieldObjects* completeField = new BDSFieldObjects(&info, resultantField, eqOfM, integrator);
811 return completeField;
812}
813
815{
816 BDSFieldE* field = CreateFieldERaw(info);
817 if (!field)
818 {return nullptr;}
819
820 // Optionally provide local to global transform using curvilinear coordinate system.
821 BDSFieldE* resultantField = field;
823 {resultantField = new BDSFieldEGlobalPlacement(field);}
824 else if (info.ProvideGlobal())
825 {resultantField = new BDSFieldEGlobal(field);}
826
827 // Equation of motion for em fields
828 G4EqMagElectricField* eqOfM = new G4EqMagElectricField(resultantField);
829
830 // Create appropriate integrator
831 G4MagIntegratorStepper* integrator = CreateIntegratorE(info, eqOfM);
832
833 BDSFieldObjects* completeField = new BDSFieldObjects(&info, resultantField, eqOfM, integrator);
834 return completeField;
835}
836
838{
839 BDSFieldE* field = nullptr;
840 switch (info.FieldType().underlying())
841 {
842 case BDSFieldType::rf:
843 {field = new BDSFieldESinusoid(info.MagnetStrength(), info.BRho()); break;}
844 case BDSFieldType::emap1d:
845 case BDSFieldType::emap2d:
846 case BDSFieldType::emap3d:
847 case BDSFieldType::emap4d:
848 {
850 if (ff)
851 {info.UpdateUserLimitsLengthMaximumStepSize(ff->SmallestSpatialStep(), true);}
852 field = ff;
853 break;
854 }
855 case BDSFieldType::efieldzero:
856 {field = new BDSFieldEZero(); break;}
857 default:
858 return nullptr;
859 break;
860 }
861
862 // Set transform for local geometry offset
863 if (field)
864 {field->SetTransform(info.TransformComplete());}
865
866 if (!info.ElectricSubFieldName().empty() && field)
867 {
868 // set the transform of the 'main' field to only the transform defined in that field definition
869 field->SetTransform(info.Transform());
870
871 auto mainField = dynamic_cast<BDSFieldEInterpolated*>(field);
872 if (!mainField)
873 {throw BDSException(__METHOD_NAME__, "subfield specified for non-field map type field - not supported");}
874
875 BDSFieldInfo* subFieldRecipe = new BDSFieldInfo(*(GetDefinition(info.ElectricSubFieldName())));
876 BDSFieldE* subFieldRaw = CreateFieldERaw(*subFieldRecipe);
877 auto subField = dynamic_cast<BDSFieldEInterpolated*>(subFieldRaw);
878 if (!subField)
879 {throw BDSException(__METHOD_NAME__, "subfield type is not a field map type field - not supported");}
880 field = new BDSFieldEInterpolated2Layer(mainField, subField);
881 // the transform goes beamline transform to the 2Layer class, then inside that the individual field transforms
882 field->SetTransform(info.TransformBeamline());
883 delete subFieldRecipe;
884 }
885
886 return field;
887}
888
890{
891 // special routine for each special / irregular field
892 BDSFieldObjects* result = nullptr;
893 switch (info.FieldType().underlying())
894 {
895 case BDSFieldType::teleporter:
896 {result = CreateTeleporter(info); break;}
897 case BDSFieldType::rmatrix:
898 {result = CreateRMatrix(info); break;}
899 case BDSFieldType::cavityfringe:
900 {result = CreateCavityFringe(info); break;}
901 case BDSFieldType::paralleltransporter:
902 {result = CreateParallelTransport(info); break;}
903 default:
904 {break;}
905 }
906 return result;
907}
908
909G4MagIntegratorStepper* BDSFieldFactory::CreateIntegratorMag(const BDSFieldInfo& info,
910 G4Mag_EqRhs* eqOfM,
911 const BDSMagnetStrength* strength)
912{
913 const G4double minimumRadiusOfCurvature = 10*CLHEP::cm;
914 G4double brho = info.BRho();
915 G4MagIntegratorStepper* integrator = nullptr;
916 // these ones can only be used for magnetic field
917 switch (info.IntegratorType().underlying())
918 {
919 case BDSIntegratorType::solenoid:
920 integrator = new BDSIntegratorSolenoid(strength, brho, eqOfM); break;
921 case BDSIntegratorType::dipolerodrigues:
922 integrator = new BDSIntegratorDipoleRodrigues(strength, brho, eqOfM); break;
923 case BDSIntegratorType::dipolerodrigues2:
924 integrator = new BDSIntegratorDipoleRodrigues2(eqOfM, minimumRadiusOfCurvature); break;
925 case BDSIntegratorType::dipolematrix:
926 integrator = new BDSIntegratorDipoleQuadrupole(strength, brho, eqOfM, minimumRadiusOfCurvature, designParticle, info.Tilt()); break;
927 case BDSIntegratorType::quadrupole:
928 integrator = new BDSIntegratorQuadrupole(strength, brho, eqOfM, minimumRadiusOfCurvature); break;
929 case BDSIntegratorType::sextupole:
930 integrator = new BDSIntegratorSextupole(strength, brho, eqOfM); break;
931 case BDSIntegratorType::octupole:
932 integrator = new BDSIntegratorOctupole(strength, brho, eqOfM); break;
933 case BDSIntegratorType::decapole:
934 integrator = new BDSIntegratorDecapole(strength, brho, eqOfM); break;
935 case BDSIntegratorType::multipolethin:
936 integrator = new BDSIntegratorMultipoleThin(strength, brho, eqOfM); break;
937 case BDSIntegratorType::dipolefringe:
938 integrator = new BDSIntegratorDipoleFringe(strength, brho, eqOfM, minimumRadiusOfCurvature, info.Tilt()); break;
939 case BDSIntegratorType::dipolefringescaling:
940 integrator = new BDSIntegratorDipoleFringeScaling(strength, brho, eqOfM, minimumRadiusOfCurvature, info.Tilt()); break;
941 case BDSIntegratorType::euler:
942 integrator = new BDSIntegratorEuler(eqOfM); break;
943 case BDSIntegratorType::kickerthin:
944 integrator = new BDSIntegratorKickerThin(strength, brho, eqOfM, minimumRadiusOfCurvature); break;
945 case BDSIntegratorType::g4rk4minimumstep:
946 integrator = new BDSIntegratorG4RK4MinStep(eqOfM, BDSGlobalConstants::Instance()->ChordStepMinimumYoke()); break;
947 case BDSIntegratorType::rmatrixthin:
948 integrator = new BDSIntegratorRMatrixThin(strength,eqOfM, info.BeamPipeRadius()); break;
949 case BDSIntegratorType::cavityfringe:
950 integrator = new BDSIntegratorCavityFringe(strength,eqOfM, info.BeamPipeRadius()); break;
951 case BDSIntegratorType::g4constrk4:
952 integrator = new G4ConstRK4(eqOfM); break;
953 case BDSIntegratorType::g4exacthelixstepper:
954 integrator = new G4ExactHelixStepper(eqOfM); break;
955 case BDSIntegratorType::g4helixexpliciteuler:
956 integrator = new G4HelixExplicitEuler(eqOfM); break;
957 case BDSIntegratorType::g4helixheum:
958 integrator = new G4HelixHeum(eqOfM); break;
959 case BDSIntegratorType::g4heliximpliciteuler:
960 integrator = new G4HelixImplicitEuler(eqOfM); break;
961 case BDSIntegratorType::g4helixmixedstepper:
962 integrator = new G4HelixMixedStepper(eqOfM); break;
963 case BDSIntegratorType::g4helixsimplerunge:
964 integrator = new G4HelixSimpleRunge(eqOfM); break;
965 case BDSIntegratorType::g4nystromrk4:
966 integrator = new G4NystromRK4(eqOfM); break;
967 case BDSIntegratorType::g4rkg3stepper:
968 integrator = new G4RKG3_Stepper(eqOfM); break;
969 case BDSIntegratorType::g4cashkarprkf45:
970 case BDSIntegratorType::g4classicalrk4:
971 case BDSIntegratorType::g4expliciteuler:
972 case BDSIntegratorType::g4impliciteuler:
973 case BDSIntegratorType::g4simpleheum:
974 case BDSIntegratorType::g4simplerunge:
975#if G4VERSION_NUMBER > 1029
976 case BDSIntegratorType::g4bogackishampine23:
977 case BDSIntegratorType::g4bogackishampine45:
978 case BDSIntegratorType::g4dolomcprirk34:
979 case BDSIntegratorType::g4dormandprince745:
980 case BDSIntegratorType::g4dormandprincerk56:
981 case BDSIntegratorType::g4tsitourasrk45:
982#endif
983#if G4VERSION_NUMBER > 1039
984 case BDSIntegratorType::g4dormandprincerk78:
985 case BDSIntegratorType::g4rk547feq1:
986 case BDSIntegratorType::g4rk547feq2:
987 case BDSIntegratorType::g4rk547feq3:
988#endif
989 integrator = CreateIntegratorEM(info, (G4EquationOfMotion*)eqOfM); break;
990 default:
991 break; // returns nullptr;
992 }
993
994 return integrator;
995}
996
997G4MagIntegratorStepper* BDSFieldFactory::CreateIntegratorEM(const BDSFieldInfo& info,
998 G4EquationOfMotion* eqOfM)
999{
1000 G4MagIntegratorStepper* integrator = nullptr;
1001 switch (info.IntegratorType().underlying())
1002 {
1003 // do the EM ones first, then complain
1004 case BDSIntegratorType::g4cashkarprkf45:
1005 integrator = new G4CashKarpRKF45(eqOfM, 8); break;
1006 case BDSIntegratorType::g4classicalrk4:
1007 integrator = new G4ClassicalRK4(eqOfM, 8); break;
1008 case BDSIntegratorType::g4expliciteuler:
1009 integrator = new G4ExplicitEuler(eqOfM, 8); break;
1010 case BDSIntegratorType::g4impliciteuler:
1011 integrator = new G4ImplicitEuler(eqOfM, 8); break;
1012 case BDSIntegratorType::g4simpleheum:
1013 integrator = new G4SimpleHeum(eqOfM, 8); break;
1014 case BDSIntegratorType::g4simplerunge:
1015 integrator = new G4SimpleRunge(eqOfM, 8); break;
1016#if G4VERSION_NUMBER > 1029
1017 case BDSIntegratorType::g4bogackishampine23:
1018 {integrator = new G4BogackiShampine45(eqOfM, 8); break;}
1019 case BDSIntegratorType::g4bogackishampine45:
1020 {integrator = new G4BogackiShampine45(eqOfM, 8); break;}
1021 case BDSIntegratorType::g4dolomcprirk34:
1022 {integrator = new G4DoLoMcPriRK34(eqOfM, 8); break;}
1023 case BDSIntegratorType::g4dormandprince745:
1024 {integrator = new G4DormandPrince745(eqOfM, 8); break;}
1025 case BDSIntegratorType::g4dormandprincerk56:
1026 {integrator = new G4DormandPrinceRK56(eqOfM, 8); break;}
1027 case BDSIntegratorType::g4tsitourasrk45:
1028 {integrator = new G4TsitourasRK45(eqOfM, 8); break;}
1029#endif
1030#if G4VERSION_NUMBER > 1039
1031 case BDSIntegratorType::g4dormandprincerk78:
1032 {integrator = new G4DormandPrinceRK78(eqOfM, 8); break;}
1033 case BDSIntegratorType::g4rk547feq1:
1034 {integrator = new G4RK547FEq1(eqOfM, 8); break;}
1035 case BDSIntegratorType::g4rk547feq2:
1036 {integrator = new G4RK547FEq2(eqOfM, 8); break;}
1037 case BDSIntegratorType::g4rk547feq3:
1038 {integrator = new G4RK547FEq3(eqOfM, 8); break;}
1039#endif
1040 case BDSIntegratorType::solenoid:
1041 case BDSIntegratorType::dipolerodrigues:
1042 case BDSIntegratorType::quadrupole:
1043 case BDSIntegratorType::sextupole:
1044 case BDSIntegratorType::octupole:
1045 case BDSIntegratorType::decapole:
1046 case BDSIntegratorType::dipolefringe:
1047 case BDSIntegratorType::g4constrk4:
1048 case BDSIntegratorType::g4exacthelixstepper:
1049 case BDSIntegratorType::g4helixexpliciteuler:
1050 case BDSIntegratorType::g4helixheum:
1051 case BDSIntegratorType::g4heliximpliciteuler:
1052 case BDSIntegratorType::g4helixmixedstepper:
1053 case BDSIntegratorType::g4helixsimplerunge:
1054 case BDSIntegratorType::g4nystromrk4:
1055 case BDSIntegratorType::g4rkg3stepper:
1056 {
1057 G4cerr << "Error: integrator \"" << info.IntegratorType() << "\" is not suitable for an EM field." << G4endl;
1058 G4cout << "Suitable integrators are:" << G4endl;
1059 std::vector<BDSIntegratorType> types = {
1060 BDSIntegratorType::g4cashkarprkf45,
1061 BDSIntegratorType::g4classicalrk4,
1062 BDSIntegratorType::g4expliciteuler,
1063 BDSIntegratorType::g4impliciteuler,
1064 BDSIntegratorType::g4simpleheum,
1065 BDSIntegratorType::g4simplerunge
1066#if G4VERSION_NUMBER > 1029
1067 ,
1068 BDSIntegratorType::g4bogackishampine23,
1069 BDSIntegratorType::g4bogackishampine45,
1070 BDSIntegratorType::g4dolomcprirk34,
1071 BDSIntegratorType::g4dormandprince745,
1072 BDSIntegratorType::g4dormandprincerk56,
1073 BDSIntegratorType::g4tsitourasrk45
1074#endif
1075#if G4VERSION_NUMBER > 1039
1076 ,
1077 BDSIntegratorType::g4dormandprincerk78,
1078 BDSIntegratorType::g4rk547feq1,
1079 BDSIntegratorType::g4rk547feq2,
1080 BDSIntegratorType::g4rk547feq3
1081#endif
1082 };
1083 for (auto type : types)
1084 {G4cout << type << G4endl;}
1085 throw BDSException(__METHOD_NAME__, "invalid integrator type");
1086 }
1087 default:
1088 break; // returns nullptr;
1089 }
1090 return integrator;
1091}
1092
1093G4MagIntegratorStepper* BDSFieldFactory::CreateIntegratorE(const BDSFieldInfo& info,
1094 G4EquationOfMotion* eqOfM)
1095{
1096 return CreateIntegratorEM(info,eqOfM);
1097}
1098
1100{
1101 G4MagneticField* bGlobalField = new BDSFieldMagZero();
1102 G4Mag_EqRhs* bEqOfMotion = new G4Mag_UsualEqRhs(bGlobalField);
1103
1104 G4MagIntegratorStepper* integrator;
1105 auto mapfile = BDSGlobalConstants::Instance()->PTCOneTurnMapFileName();
1106 BDSPTCOneTurnMap* otm = nullptr;
1107
1108 if (!mapfile.empty())
1109 {
1110 otm = new BDSPTCOneTurnMap(mapfile, designParticle);
1112 }
1113
1114 integrator = new BDSIntegratorTeleporter(bEqOfMotion, info.TransformComplete(),
1115 (*info.MagnetStrength())["length"],
1116 otm);
1117
1118 BDSFieldObjects* completeField = new BDSFieldObjects(&info, bGlobalField,
1119 bEqOfMotion, integrator);
1120 return completeField;
1121}
1122
1124{
1125 G4MagneticField* bGlobalField = new BDSFieldMagZero();
1126 G4Mag_EqRhs* bEqOfMotion = new G4Mag_UsualEqRhs(bGlobalField);
1127 G4MagIntegratorStepper* integrator = new BDSIntegratorRMatrixThin(info.MagnetStrength(),bEqOfMotion,0.95*info.BeamPipeRadius());
1128 BDSFieldObjects* completeField = new BDSFieldObjects(&info, bGlobalField,
1129 bEqOfMotion, integrator);
1130 return completeField;
1131}
1132
1134{
1135 BDSFieldMag* bGlobalField = new BDSFieldMagZero();
1136 BDSMagUsualEqRhs* bEqOfMotion = new BDSMagUsualEqRhs(bGlobalField);
1137 G4MagIntegratorStepper* integrator = new BDSIntegratorCavityFringe(info.MagnetStrength(),bEqOfMotion,0.95*info.BeamPipeRadius());
1138 BDSFieldObjects* completeField = new BDSFieldObjects(&info, bGlobalField,
1139 bEqOfMotion, integrator);
1140 return completeField;
1141}
1142
1144{
1145 G4MagneticField* bGlobalField = new BDSFieldMagZero();
1146 G4Mag_EqRhs* bEqOfMotion = new G4Mag_UsualEqRhs(bGlobalField);
1147 G4MagIntegratorStepper* integrator = new BDSIntegratorParallelTransport(bEqOfMotion);
1148 BDSFieldObjects* completeField = new BDSFieldObjects(&info, bGlobalField,
1149 bEqOfMotion, integrator);
1150 return completeField;
1151}
1152
1154{
1155 G4double result = 1.0;
1156 if (st->KeyHasBeenSet("scalingOuter"))
1157 {result = (*st)["scalingOuter"];}
1158 return result;
1159}
G4double aper1
Public member for direct access.
General exception with possible name of object and message.
Definition: BDSException.hh:35
Wrapper class to convert to global coordinates using a navigator for placements.
A base class for electric fields in local to be used in global coordinates.
Two interpolated fields in one. One takes precedence in a subregion.
Class to provide scaling and a base class pointer for interpolator fields.
A base class for electro-magnetic fields in local to be used in global coordinates.
A base class for electro-magnetic fields in local to be used in global coordinates.
Class to provide scaling and a base class pointer for interpolator fields.
Pill box cavity electro-magnetic field.
Null EM field - for special cases where we need a valid object.
Interface for BDSIM electro-magnetic fields that may or may not be local.
Definition: BDSFieldEM.hh:39
virtual void SetTransform(const G4Transform3D &transformIn)
Definition: BDSFieldEM.cc:68
A sinusoidal electric (only) field that doesn't vary with position. Uses cosine.
Null E field - for special cases where we need a valid object.
Interface for BDSIM electric fields that may or may not be local.
Definition: BDSFieldE.hh:37
virtual void SetTransform(const G4Transform3D &transformIn)
Definition: BDSFieldE.cc:66
Factory that produces fields and their associated objects.
G4double GetOuterScaling(const BDSMagnetStrength *st) const
Return the parameter "outerScaling" from strength st, but default to 1.
void PrepareFieldDefinitions(const std::vector< GMAD::Field > &definitions, G4double defaultBRho)
Prepare all required definitions that can be used dynamically.
BDSFieldObjects * CreateFieldE(const BDSFieldInfo &info)
Create an electric field.
BDSFieldObjects * CreateFieldIrregular(const BDSFieldInfo &info)
Create an irregular (special) field.
BDSFieldObjects * CreateRMatrix(const BDSFieldInfo &info)
Create special rmatrix 'field' that applies an rmatrix.
void PrepareFieldStrengthFromParameters(BDSMagnetStrength *st, const G4String &fieldParameters, G4double &poleTipRadius) const
static BDSPrimaryGeneratorAction * primaryGeneratorAction
Cache of primary generator action.
BDSFieldE * CreateFieldERaw(const BDSFieldInfo &info)
Creat just the electric field object.
G4double ConvertToDoubleWithException(const G4String &value, const G4String &parameterNameForError) const
Convert the string 'value' to a double. Throw an exception including the parameterNameForError if it ...
static BDSFieldFactory * Instance()
Public accessor method for singleton pattern.
BDSFieldObjects * CreateField(const BDSFieldInfo &info, const BDSMagnetStrength *scalingStrength=nullptr, const G4String &scalingKey="none")
Main interface to field factory.
BDSFieldObjects * CreateParallelTransport(const BDSFieldInfo &info)
static BDSFieldFactory * instance
Instance - singleton pattern.
BDSFieldInfo * GetDefinition(const G4String &name) const
BDSFieldObjects * CreateCavityFringe(const BDSFieldInfo &info)
Create special rf cavity fringe 'field' that applies an rmatrix.
std::map< G4String, BDSFieldInfo * > parserDefinitions
BDSFieldInfo definitions prepare from parser vector of definitions.
G4MagIntegratorStepper * CreateIntegratorMag(const BDSFieldInfo &info, G4Mag_EqRhs *eqOfM, const BDSMagnetStrength *strength)
G4MagIntegratorStepper * CreateIntegratorE(const BDSFieldInfo &info, G4EquationOfMotion *eqOfM)
BDSFieldObjects * CreateFieldMag(const BDSFieldInfo &info, const BDSMagnetStrength *scalingStrength=nullptr, const G4String &scalingKey="none")
Create a purely magnetic field.
static BDSInterpolatorType DefaultInterpolatorType(G4int numberOfDimensions)
Suggest a default interpolator.
BDSFieldObjects * CreateFieldEM(const BDSFieldInfo &info)
Create a general EM field.
BDSFieldMag * CreateFieldMagRaw(const BDSFieldInfo &info, const BDSMagnetStrength *scalingStrength=nullptr, const G4String &scalingKey="none")
Creat just the magnetic field object.
BDSFieldFactory()
Private default constructor as singleton class.
BDSFieldObjects * CreateTeleporter(const BDSFieldInfo &info)
G4MagIntegratorStepper * CreateIntegratorEM(const BDSFieldInfo &info, G4EquationOfMotion *eqOfM)
static const BDSParticleDefinition * designParticle
Cache of design particle for fields.
All info required to build complete field of any type.
Definition: BDSFieldInfo.hh:65
G4String MagneticSubFieldName() const
Accessor.
void UpdateUserLimitsLengthMaximumStepSize(G4double maximumStepSize, G4bool warn=false) const
G4Transform3D TransformComplete() const
Compound transform of field + beam line transform.
G4bool SecondFieldOnLeft() const
Accessor.
G4double PoleTipRadius() const
Accessor.
G4bool ProvideGlobal() const
Accessor.
BDSIntegratorType IntegratorType() const
Accessor.
G4double BeamPipeRadius() const
Accessor.
G4String ElectricSubFieldName() const
Accessor.
G4Transform3D Transform() const
Transform for the field definition only.
G4bool UsePlacementWorldTransform() const
Accessor.
G4double Tilt() const
Accessor.
BDSFieldType FieldType() const
Accessor.
BDSMagnetStrength * MagnetStrength() const
Accessor.
G4double BRho() const
Accessor.
G4Transform3D TransformBeamline() const
Transform from the curvilinear coordinates to the beam line component.
BDSFieldEMInterpolated * LoadEMField(const BDSFieldInfo &info)
Main interface to load an electro-magnetic field.
BDSFieldMagInterpolated * LoadMagField(const BDSFieldInfo &info, const BDSMagnetStrength *scalingStrength=nullptr, const G4String &scalingKey="none")
Main interface to load a magnetic field.
static BDSFieldLoader * Instance()
Singleton accessor.
BDSFieldEInterpolated * LoadEField(const BDSFieldInfo &info)
Main interface to load an electric field.
Class that provides the magnetic strength in a decapole.
A perfect magetic dipole in 3D, normal field inside 1/2 poleTipRadius.
A perfect magetic dipole in 3D, normal field inside 1/2 poleTipRadius.
Class that provides the magnetic strength in a mixed dipole / quadrupole.
A uniform dipole field.
A base class for magnetic fields in local to be used in global coordinates.
A base class for magnetic fields in local to be used in global coordinates.
Two interpolated fields in one. One takes precedence in a subregion.
Class to provide scaling and a base class pointer for interpolator fields.
Sum of two multipole fields spaced by a distance in x.
Sum of two multipole fields spaced by a distance in x.
A simple paramaterisation of N-Pole outer yoke magnetic field.
A simple parameterisation of N-Pole outer yoke magnetic field.
Class that provides the magnetic strength in a quadrupole.
Field of a Muon Spoiler.
Class that provides the magnetic strength in an octupole.
Class that provides the magnetic strength in a quadrupole.
Class that provides the magnetic strength in a sextupole.
A wrapper class for BDSFieldMagSkew where it owns the field.
Class that provides the magnetic field due to a cylinder of current.
Class that provides the magnetic strength in a quadrupole.
Null B field - for special cases where we need a valid object.
Interface for static magnetic fields that may or may not be local.
Definition: BDSFieldMag.hh:37
virtual void SetTransform(const G4Transform3D &transformIn)
Definition: BDSFieldMag.cc:70
A holder for all the Geant4 field related objects.
static BDSGlobalConstants * Instance()
Access method.
Integrator for RF cavity fringes. Only the transverse momentum kicks are applied, this integrator wil...
Integrator for Decapolar field.
Derived fringe field integrator that does normalise to momentum.
Integrator that ignores the field and uses the analytical solution for a dipole kick.
Integrator for combined dipole and quadrupolar field.
Exact helix through pure dipole field.
Stepper that calculates trajectory through uniform magnetic field.
BDSIM 2nd order Euler integration.
Integrator that wraps a G4ClassicalRK4 and below a minimum step size uses a drift.
Integrator for thin h or v kick.
Integrator that ignores the field and uses the analytical solution to a multipole.
Integrator for octupole field.
Integrator that just moves the particle parallel to the s axis.
Integrator that ignores the field and uses the analytical solution to a quadrupole.
Integrator that just moves the particle parallel to the s axis.
Integrator for sextupole field.
Integrator that ignores the field and uses the analytical solution to a solenoid.
Custom unphysical integrator to advance particle in teleporter.
Override G4Mag_UsualEqRhs, provides BDSIM integrators access to particle attributes.
static const G4double beamSeparation
Used in many places - make it a constant in the code and put here as most relevant.
Efficient storage of magnet strengths.
G4bool KeyHasBeenSet(const G4String &key) const
Whether a key has been set.
static G4bool ValidKey(const G4String &key)
Whether or not the supplied key is a valid magnet strength parameter.
static const std::vector< G4String > & AllKeys()
Accessor to all keys.
static G4double Unit(const G4String &key)
Access a unit factor for a given key.
Class to load and use PTC 1 turn map.
static BDSParser * Instance()
Access method.
Definition: BDSParser.cc:28
Wrapper for particle definition.
G4double BRho() const
Accessor.
Generates primary particle vertices using BDSBunch.
void RegisterPTCOneTurnMap(BDSPTCOneTurnMap *otmIn)
type underlying() const
return underlying value (can be used in switch statement)
BDSArrayReflectionTypeSet DetermineArrayReflectionTypeSet(const G4String &arrayReflectionType)
Return a std::set of reflection types. Split string on white space.
BDSInterpolatorType InterpolatorTypeSpecificFromAuto(G4int nDimension, BDSInterpolatorType autoType)
G4bool InterpolatorTypeIsAuto(BDSInterpolatorType typeIn)
Return true if the type is one containing 'auto'.
BDSFieldFormat DetermineFieldFormat(G4String fieldformat)
Function that gives corresponding enum value for string (case-insensitive)
std::pair< G4String, G4String > SplitOnColon(const G4String &formatAndPath)
BDSFieldType DetermineFieldType(G4String fieldType)
Function that gives corresponding enum value for string (case-insensitive)
Definition: BDSFieldType.cc:88
BDSInterpolatorType DetermineInterpolatorType(G4String interpolatorType)
Function that determines enum from string (case-insensitive).
BDSFieldClassType DetermineFieldClassType(BDSFieldType fieldType)
Function that gives the corresponding enum value for a field type enum.
G4String GetFullPath(G4String filename, bool excludeNameFromPath=false, bool useCWDForPrefix=false)
G4UserLimits * CreateUserLimits(G4UserLimits *defaultUL, G4double length, G4double fraction=1.6)
BDSIntegratorType DetermineIntegratorType(G4String integratorType)
Function that determines enum from string (case-insensitive).
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())
G4int NDimensionsOfInterpolatorType(const BDSInterpolatorType &it)
Report the number of dimensions for that interpolator type.
G4int NDimensionsOfFieldFormat(const BDSFieldFormat &ff)
Report the number of dimensions for that format.
std::map< G4String, G4String > GetUserParametersMap(const G4String &userParameters, char delimiter=':')
Take one long string and split on space and then on colon. "key1:value1 key2:value2" etc.