gfields
Loading...
Searching...
No Matches
gfield_multipoles.cc
Go to the documentation of this file.
1// CLHEP units
2#include "CLHEP/Units/PhysicalConstants.h"
3
4// geant4 headers
5#include "G4ThreeVector.hh"
6
7// gfields
8#include "gfield_multipoles.h"
9#include "gfieldConventions.h"
10
11// gemv
12// #include "gutilities.h"
13
14
15// tells the DLL how to create a GFieldFactory in each plugin .so/.dylib
16extern "C" GField* GFieldFactory(const std::shared_ptr<GOptions>& g) { return static_cast<GField*>(new GField_MultipolesFactory(g)); }
17
18
19// Implementation notes (non-Doxygen):
20// - This follows common accelerator multipole conventions for transverse fields.
21// - Strength is defined at a reference radius of 1 m (see r0 below).
22// - References (background reading):
23// - https://cds.cern.ch/record/1333874/files/1.pdf
24// - https://uspas.fnal.gov/materials/12MSU/magnet_elements.pdf
25// - https://cas.web.cern.ch/sites/default/files/lectures/bruges-2009/wolski-1.pdf
26void GField_MultipolesFactory::GetFieldValue(const double pos[3], G4double* bfield) const {
27 // Evaluate an ideal multipole field at the query point in the lab frame.
28 // The computation is performed in a magnet-local frame and then rotated back.
29
30 // ======= Configuration / conventions =======
31 // strength: Tesla at reference radius r0 for all multipole orders
32 const double r0 = CLHEP::m; // reference radius; make this a member if you want
33
34 // ======= Basic checks =======
35 if (pole_number < 2 || (pole_number % 2) != 0) {
37 "Pole number must be an even integer >= 2 (2=dipole,4=quadrupole,...)");
38 }
39
40 // ======= Positions and local frame =======
41 const G4ThreeVector x0(pos[0], pos[1], pos[2]); // query point (lab)
42 const G4ThreeVector x1(origin[0], origin[1], origin[2]); // magnet origin (lab)
43
44 // shift to magnet-centered coordinates and "unroll" the magnet by -rotation_angle
45 G4ThreeVector p = x0 - x1;
46 if (rotaxis == 0) p.rotateX(-rotation_angle);
47 else if (rotaxis == 1) p.rotateY(-rotation_angle);
48 else if (rotaxis == 2) p.rotateZ(-rotation_angle);
49
50 // Identify transverse plane (u,v) ⟂ to the axis, and the axial coordinate w (unused here)
51 double u = 0.0, v = 0.0;
52 switch (rotaxis) {
53 case 0: u = p.y();
54 v = p.z();
55 break; // axis = X ⇒ transverse plane = (Y,Z)
56 case 1: u = p.z();
57 v = p.x();
58 break; // axis = Y ⇒ transverse plane = (Z,X)
59 case 2: u = p.x();
60 v = p.y();
61 break; // axis = Z ⇒ transverse plane = (X,Y)
62 default: break;
63 }
64
65 const double r = std::hypot(u, v); // transverse radius
66 const double phi = std::atan2(v, u); // azimuth in transverse plane
67
68 // ======= Axial (solenoid-like) mode if explicitly requested =======
69 if (longitudinal) {
70 // Uniform axial field aligned with rotaxis; not a multipole.
71 G4ThreeVector B_local(0, 0, 0);
72 switch (rotaxis) {
73 case 0: B_local.setX(strength);
74 break;
75 case 1: B_local.setY(strength);
76 break;
77 case 2: B_local.setZ(strength);
78 break;
79 default: break;
80 }
81
82 // Roll back to lab
83 G4ThreeVector B_lab = B_local;
84 if (rotaxis == 0) B_lab.rotateX(+rotation_angle);
85 else if (rotaxis == 1) B_lab.rotateY(+rotation_angle);
86 else if (rotaxis == 2) B_lab.rotateZ(+rotation_angle);
87
88 bfield[0] = B_lab.x();
89 bfield[1] = B_lab.y();
90 bfield[2] = B_lab.z();
91
92 log->info(2, "Axial field mode (solenoid-like). Strength: ", strength,
93 " T, Field: (", bfield[0], ", ", bfield[1], ", ", bfield[2], ")");
94 return;
95 }
96
97 // ======= Transverse multipole (standard accelerator definition) =======
98 const int n = pole_number / 2; // 1=dipole, 2=quadrupole, 3=sextupole, ...
99 const int a = n - 1; // power of r
100
101 // r^a scaling with reference radius r0; for dipole (a=0) this is 1
102 double ra = 1.0;
103 if (a > 0) {
104 // If r=0 and a>0, |B|=0 (for ideal multipoles).
105 if (r == 0.0) {
106 bfield[0] = bfield[1] = bfield[2] = 0.0;
107 return;
108 }
109 ra = std::pow(r / r0, static_cast<double>(a));
110 }
111
112 // "Normal" multipole angular dependence:
113 // Bu = strength * r^(n-1) * cos((n-1)*phi)
114 // Bv = strength * r^(n-1) * sin((n-1)*phi)
115 // This yields a constant transverse dipole for n=1.
116 const double Bu = strength * ra * std::cos(a * phi);
117 const double Bv = strength * ra * std::sin(a * phi);
118
119 // Place (Bu,Bv) into the correct transverse components of B_local; axial component = 0
120 G4ThreeVector B_local(0, 0, 0);
121 switch (rotaxis) {
122 case 0: B_local.setY(Bu);
123 B_local.setZ(Bv);
124 break; // axis X → (Y,Z)
125 case 1: B_local.setZ(Bu);
126 B_local.setX(Bv);
127 break; // axis Y → (Z,X)
128 case 2: B_local.setX(Bu);
129 B_local.setY(Bv);
130 break; // axis Z → (X,Y)
131 default: break;
132 }
133
134 // Rotate (roll) back to lab
135 G4ThreeVector B_lab = B_local;
136 if (rotaxis == 0) B_lab.rotateX(+rotation_angle);
137 else if (rotaxis == 1) B_lab.rotateY(+rotation_angle);
138 else if (rotaxis == 2) B_lab.rotateZ(+rotation_angle);
139
140 // Output
141 bfield[0] = B_lab.x();
142 bfield[1] = B_lab.y();
143 bfield[2] = B_lab.z();
144
145 log->info(2, "Pole Number: ", pole_number,
146 ", n: ", n,
147 ", Strength: ", strength,
148 ", Requested at: (", pos[0], ", ", pos[1], ", ", pos[2], ")",
149 ", Rotation angle: ", rotation_angle,
150 ", Rotation axis: ", rotaxis,
151 ", longitudinal: ", longitudinal,
152 ", Field: (", bfield[0], ", ", bfield[1], ", ", bfield[2], ")");
153}
154
155
156// Load and cache multipole parameters from the provided field definition (non-Doxygen summary).
158 gfield_definitions = gfd;
159
160 pole_number = get_field_parameter_int("pole_number");
161 origin[0] = get_field_parameter_double("vx");
162 origin[1] = get_field_parameter_double("vy");
163 origin[2] = get_field_parameter_double("vz");
164 rotation_angle = get_field_parameter_double("rotation_angle");
165 auto rot_axis_option = gfield_definitions.field_parameters["rotaxis"];
166 longitudinal = false;
167
168 // Longitudinal mode: return a uniform axial field aligned with rotaxis.
169 if (gfield_definitions.field_parameters["longitudinal"] == "true") {
170 longitudinal = true;
171 log->info(1, "Longitudinal field");
172 }
173 else { log->info(1, "Transverse field"); }
174
175 // Parse rotaxis as a single axis character.
176 if (rot_axis_option == "X" || rot_axis_option == "x") { rotaxis = 0; }
177 else if (rot_axis_option == "Y" || rot_axis_option == "y") { rotaxis = 1; }
178 else if (rot_axis_option == "Z" || rot_axis_option == "z") { rotaxis = 2; }
179 else {
180 log->error(ERR_WRONG_FIELD_ROTATION, "GField_MultipolesFactory::load_field_definitions: Rotation axis " + gfield_definitions.field_parameters["rotaxis"] +
181 " not supported. Exiting.");
182 }
183
184 log->info(1, "Rotation axis: ", rotaxis);
185 strength = get_field_parameter_double("strength");
186}
std::shared_ptr< GLogger > log
Factory class implementing an ideal multipole magnetic field.
void GetFieldValue(const double pos[3], G4double *bfield) const override
Evaluate the magnetic field B at a given position.
void load_field_definitions(GFieldDefinition gfd) override
Load and cache field-definition parameters for fast field evaluation.
Abstract base class representing a magnetic field.
Definition gfield.h:99
GFieldDefinition gfield_definitions
Stored field definition used for configuration and logging.
Definition gfield.h:185
double get_field_parameter_double(const std::string &key)
Convenience accessor for floating-point parameters stored in field_parameters.
Definition gfield.h:151
int get_field_parameter_int(const std::string &key)
Convenience accessor for integer-valued parameters stored in field_parameters.
Definition gfield.h:144
void info(int level, Args &&... args) const
void error(int exit_code, Args &&... args) const
GField * GFieldFactory(const std::shared_ptr< GOptions > &g)
#define ERR_WRONG_FIELD_ROTATION
#define ERR_WRONG_POLE_NUMBER
Lightweight configuration carrier used to load and configure a GField plugin.
Definition gfield.h:26
std::map< std::string, std::string > field_parameters
Field parameters stored as string expressions (often including unit expressions).
Definition gfield.h:45