16 std::cerr <<
m_className <<
"::SetGeometry: Null pointer.\n";
41 const double ,
double& wx,
double& wy,
42 double& wz,
const std::string& ) {
44 std::cerr <<
m_className <<
"::WeightingField: Function not implemented.\n";
50 const double ,
const double ,
51 double& wx,
double& wy,
double& wz,
52 const std::string& ) {
54 std::cerr <<
m_className <<
"::DelayedWeightingField: Not implemented.\n";
61 const std::string& ) {
63 std::cerr <<
m_className <<
"::WeightingPotential: Not implemented.\n";
72 const std::string& ) {
75 <<
"::DelayedWeightingPotential: Not implemented.\n";
82 double& bx,
double& by,
double& bz,
int& status) {
87 std::cout <<
m_className <<
"::MagneticField: Field at (" << x <<
", " << y
88 <<
", " << z <<
") is (" << bx <<
", " << by <<
", " << bz
100 double& xmax,
double& ymax,
double& zmax) {
106 double& xmax,
double& ymax,
double& zmax) {
111 const double ,
const double ,
112 const double ,
const double ,
113 double& ,
double& ,
double& ,
114 const bool ,
double& ) {
119 const double y0,
const double ,
double& xw,
120 double& yw,
double& rw) {
128 const double r,
const unsigned int nI) {
131 std::cerr <<
m_className <<
"::IntegrateFluxCircle:\n"
132 <<
" Number of intervals must be > 0.\n";
136 constexpr size_t nG = 6;
137 constexpr std::array<double, nG> t = {
138 -0.932469514203152028, -0.661209386466264514, -0.238619186083196909,
139 0.238619186083196909, 0.661209386466264514, 0.932469514203152028};
140 constexpr std::array<double, nG> w = {
141 0.171324492379170345, 0.360761573048138608, 0.467913934572691047,
142 0.467913934572691047, 0.360761573048138608, 0.171324492379170345};
145 const double d = TwoPi / nI;
146 const double h = 0.5 * d;
148 double ex = 0., ey = 0., ez = 0.;
153 for (
size_t i = 0; i < nG; ++i) {
154 const double phi0 = h * (1. + t[i]);
155 for (
unsigned int k = 0; k < nI; ++k) {
156 const double phi = phi0 + k * d;
157 const double cp = cos(phi);
158 const double sp = sin(phi);
159 ElectricField(xc + cp * r, yc + sp * r, 0., ex, ey, ez, m, status);
160 s += w[i] * r * (ex * cp + ey * sp);
163 return h * s * VacuumPermittivity;
167 const double zc,
const double r,
168 const unsigned int nI) {
171 std::cerr <<
m_className <<
"::IntegrateFluxSphere:\n"
172 <<
" Number of intervals must be > 0.\n";
176 constexpr size_t nG = 6;
177 constexpr std::array<double, nG> t = {
178 -0.932469514203152028, -0.661209386466264514, -0.238619186083196909,
179 0.238619186083196909, 0.661209386466264514, 0.932469514203152028};
180 constexpr std::array<double, nG> w = {
181 0.171324492379170345, 0.360761573048138608, 0.467913934572691047,
182 0.467913934572691047, 0.360761573048138608, 0.171324492379170345};
184 const double r2 = r * r;
186 const double dt = Pi / nI;
187 const double ht = 0.5 * dt;
189 const double dp = TwoPi / nI;
190 const double hp = 0.5 * dp;
192 double ex = 0., ey = 0., ez = 0.;
198 for (
size_t i = 0; i < nG; ++i) {
199 const double theta0 = ht * (1. + t[i]) - HalfPi;
200 for (
unsigned int k = 0; k < nI; ++k) {
201 const double theta = theta0 + k * dt;
202 const double ct = cos(theta);
203 const double st = sin(theta);
204 const double z = zc + st * r;
207 for (
size_t ii = 0; ii < nG; ++ii) {
208 const double phi0 = hp * (1. + t[ii]);
209 for (
unsigned int kk = 0; kk < nI; ++kk) {
210 const double phi = phi0 + kk * dp;
211 const double cp = cos(phi);
212 const double sp = sin(phi);
213 const double x = xc + cp * ct * r;
214 const double y = yc + sp * ct * r;
216 s1 += w[ii] * ((ex * cp + ey * sp) * ct + ez * st);
219 s2 += w[i] * r2 * ct * hp * s1;
222 return ht * s2 * VacuumPermittivity;
226 const double x0,
const double y0,
const double z0,
const double dx1,
227 const double dy1,
const double dz1,
const double dx2,
const double dy2,
228 const double dz2,
const unsigned int nU,
const unsigned int nV) {
231 x0, y0, z0, dx1, dy1, dz1, dx2, dy2, dz2, nU, nV,
false,
"");
235 const double x0,
const double y0,
const double z0,
const double dx1,
236 const double dy1,
const double dz1,
const double dx2,
const double dy2,
237 const double dz2,
const unsigned int nU,
const unsigned int nV) {
240 x0, y0, z0, dx1, dy1, dz1, dx2, dy2, dz2, nU, nV,
true,
id);
244 const double x0,
const double y0,
const double z0,
const double dx1,
245 const double dy1,
const double dz1,
const double dx2,
const double dy2,
246 const double dz2,
const unsigned int nU,
const unsigned int nV,
247 const bool wfield,
const std::string& label) {
250 if (nU <= 1 || nV <= 1) {
251 std::cerr <<
m_className <<
"::IntegrateFluxParallelogram:\n"
252 <<
" Number of points to integrate over must be > 1.\n";
256 constexpr size_t nG = 6;
257 constexpr std::array<double, nG> t = {
258 -0.932469514203152028, -0.661209386466264514, -0.238619186083196909,
259 0.238619186083196909, 0.661209386466264514, 0.932469514203152028};
260 constexpr std::array<double, nG> w = {
261 0.171324492379170345, 0.360761573048138608, 0.467913934572691047,
262 0.467913934572691047, 0.360761573048138608, 0.171324492379170345};
265 const double xn = dy1 * dz2 - dz1 * dy2;
266 const double yn = dz1 * dx2 - dx1 * dz2;
267 const double zn = dx1 * dy2 - dy1 * dx2;
269 std::cout <<
m_className <<
"::IntegrateFluxParallelogram:\n"
270 <<
" Normal vector = " << xn <<
", " << yn <<
", " << zn
274 const double d1 = dx1 * dx1 + dy1 * dy1 + dz1 * dz1;
275 const double d2 = dx2 * dx2 + dy2 * dy2 + dz2 * dz2;
276 if (xn * xn + yn * yn + zn * zn < 1.e-10 *
sqrt(d1 * d2) ||
277 d1 < 1.e-10 * d2 || d2 < 1.e-10 * d1) {
278 std::cerr <<
m_className <<
"::IntegrateFluxParallelogram:\n"
279 <<
" Parallelogram does not have non-zero area.\n";
284 const double du = 1. / nU;
285 const double hu = 0.5 * du;
286 const double dv = 1. / nV;
287 const double hv = 0.5 * dv;
289 double fx = 0., fy = 0., fz = 0.;
294 for (
size_t i = 0; i < nG; ++i) {
295 const double v0 = hv * (1. + t[i]);
296 for (
unsigned int k = 0; k < nV; ++k) {
297 const double v = v0 + k * dv;
299 for (
size_t ii = 0; ii < nG; ++ii) {
300 const double u0 = hu * (1. + t[ii]);
301 for (
unsigned int kk = 0; kk < nU; ++kk) {
302 const double u = u0 + kk * du;
303 const double x = x0 + u * dx1 + v * dx2;
304 const double y = y0 + u * dy1 + v * dy2;
305 const double z = z0 + u * dz1 + v * dz2;
311 s1 += w[ii] * (fx * xn + fy * yn + fz * zn);
314 s2 += w[i] * hu * s1;
321 const double z0,
const double x1,
322 const double y1,
const double z1,
323 const double xp,
const double yp,
324 const double zp,
const unsigned int nI,
328 const double pmag2 = xp * xp + yp * yp + zp * zp;
330 std::cerr <<
m_className <<
"::IntegrateFluxLine:\n"
331 <<
" Normal vector has zero length; flux set to 0.\n";
334 const double pmag = sqrt(pmag2);
335 const double xn = xp / pmag;
336 const double yn = yp / pmag;
337 const double zn = zp / pmag;
341 std::cerr <<
m_className <<
"::IntegrateFluxLine:\n"
342 <<
" Number of points to integrate over must be > 1.\n";
346 const double vx = x1 - x0;
347 const double vy = y1 - y0;
348 const double vz = z1 - z0;
349 const double vmag2 = vx * vx + vy * vy + vz * vz;
351 std::cerr <<
m_className <<
"::IntegrateFluxLine:\n"
352 <<
" Segment has zero length; flux set to 0.\n";
355 const double vmag = sqrt(vmag2);
357 if (fabs(vx * xn + vy * yn + vz * zn) > 1.e-4 * vmag) {
358 std::cerr <<
m_className <<
"::IntegrateFluxLine:\n"
359 <<
" Segment is not perpendicular to norm vector.\n";
364 constexpr size_t nG = 6;
365 constexpr std::array<double, nG> t = {
366 -0.932469514203152028, -0.661209386466264514, -0.238619186083196909,
367 0.238619186083196909, 0.661209386466264514, 0.932469514203152028};
368 constexpr std::array<double, nG> w = {
369 0.171324492379170345, 0.360761573048138608, 0.467913934572691047,
370 0.467913934572691047, 0.360761573048138608, 0.171324492379170345};
372 const double d = 1. / nI;
373 const double h = 0.5 * d;
375 double ex = 0., ey = 0., ez = 0.;
379 for (
size_t i = 0; i < nG; ++i) {
380 const double u0 = h * (1. + t[i]);
381 for (
unsigned int k = 0; k < nI; ++k) {
382 const double u = u0 + k * d;
383 const double x = x0 + u * vx;
384 const double y = y0 + u * vy;
385 const double z = z0 + u * vz;
387 double fn = ex * xn + ey * yn + ez * zn;
390 fn = isign * fn > 0 ? fabs(fn) : -1.;
virtual double WeightingPotential(const double x, const double y, const double z, const std::string &label)
virtual bool GetBoundingBox(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
Get the bounding box coordinates.
std::array< bool, 3 > m_rotationSymmetric
Rotation symmetry around x-axis, y-axis, z-axis.
std::array< bool, 3 > m_mirrorPeriodic
Mirror periodicity in x, y, z.
virtual void Reset()=0
Reset the component.
virtual void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status)=0
virtual Medium * GetMedium(const double x, const double y, const double z)
Get the medium at a given location (x, y, z).
virtual void WeightingField(const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label)
bool m_debug
Switch on/off debugging messages.
virtual void DelayedWeightingField(const double x, const double y, const double z, const double t, double &wx, double &wy, double &wz, const std::string &label)
double IntegrateWeightingFluxParallelogram(const std::string &label, const double x0, const double y0, const double z0, const double dx1, const double dy1, const double dz1, const double dx2, const double dy2, const double dz2, const unsigned int nU=20, const unsigned int nV=20)
virtual bool IsWireCrossed(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, double &xc, double &yc, double &zc, const bool centre, double &rc)
Component()=delete
Default constructor.
std::array< bool, 3 > m_periodic
Simple periodicity in x, y, z.
virtual void Clear()
Reset.
std::string m_className
Class name.
double IntegrateFluxCircle(const double xc, const double yc, const double r, const unsigned int nI=50)
Geometry * m_geometry
Pointer to the geometry.
double IntegrateFluxSphere(const double xc, const double yc, const double zc, const double r, const unsigned int nI=20)
virtual void MagneticField(const double x, const double y, const double z, double &bx, double &by, double &bz, int &status)
bool m_ready
Ready for use?
void SetMagneticField(const double bx, const double by, const double bz)
Set a constant magnetic field.
virtual bool IsInTrapRadius(const double q0, const double x0, const double y0, const double z0, double &xw, double &yw, double &rw)
std::array< double, 3 > m_b0
Constant magnetic field.
std::array< bool, 3 > m_axiallyPeriodic
Axial periodicity in x, y, z.
virtual bool GetElementaryCell(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
Get the coordinates of the elementary cell.
virtual double DelayedWeightingPotential(const double x, const double y, const double z, const double t, const std::string &label)
double IntegrateFluxParallelogram(const double x0, const double y0, const double z0, const double dx1, const double dy1, const double dz1, const double dx2, const double dy2, const double dz2, const unsigned int nU=20, const unsigned int nV=20)
virtual void SetGeometry(Geometry *geo)
Define the geometry.
double IntegrateFluxLine(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, const double xp, const double yp, const double zp, const unsigned int nI, const int isign=0)
Abstract base class for geometry classes.
virtual bool GetBoundingBox(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)=0
Get the bounding box (envelope of the geometry).
virtual Medium * GetMedium(const double x, const double y, const double z, const bool tesselated=false) const =0
Retrieve the medium at a given point.
Abstract base class for media.
DoubleAc sqrt(const DoubleAc &f)