OPAL (Object Oriented Parallel Accelerator Library)  2021.1.99
OPAL
Bend2D.cpp
Go to the documentation of this file.
1 // ------------------------------------------------------------------------
2 // $RCSfile: Bend2D.cpp,v $
3 // ------------------------------------------------------------------------
4 // $Revision: 1.1.1.1 $
5 // ------------------------------------------------------------------------
6 // Copyright: see Copyright.readme
7 // ------------------------------------------------------------------------
8 //
9 // Definitions for class: Bend2D
10 // Defines the abstract interface for a sector bend magnet.
11 //
12 // ------------------------------------------------------------------------
13 // Class category: AbsBeamline
14 // ------------------------------------------------------------------------
15 //
16 // $Date: 2000/03/27 09:32:31 $
17 // $Author: fci $
18 //
19 // ------------------------------------------------------------------------
20 
21 #include "AbsBeamline/Bend2D.h"
24 #include "Utilities/Options.h"
25 #include "Utilities/Util.h"
26 #include "Fields/Fieldmap.h"
29 
30 #include "gsl/gsl_poly.h"
31 
32 #include <iostream>
33 #include <fstream>
34 
35 extern Inform *gmsg;
36 
37 // Class Bend2D
38 // ------------------------------------------------------------------------
39 
41  Bend2D("")
42 {}
43 
44 Bend2D::Bend2D(const Bend2D &right):
45  BendBase(right),
46  messageHeader_m(" * "),
47  pusher_m(right.pusher_m),
48  designRadius_m(right.designRadius_m),
49  exitAngle_m(right.exitAngle_m),
50  fieldIndex_m(right.fieldIndex_m),
51  startField_m(right.startField_m),
52  endField_m(right.endField_m),
53  widthEntranceFringe_m(right.widthEntranceFringe_m),
54  widthExitFringe_m(right.widthExitFringe_m),
55  reinitialize_m(right.reinitialize_m),
56  entranceParameter1_m(right.entranceParameter1_m),
57  entranceParameter2_m(right.entranceParameter2_m),
58  entranceParameter3_m(right.entranceParameter3_m),
59  exitParameter1_m(right.exitParameter1_m),
60  exitParameter2_m(right.exitParameter2_m),
61  exitParameter3_m(right.exitParameter3_m),
62  engeCoeffsEntry_m(right.engeCoeffsEntry_m),
63  engeCoeffsExit_m(right.engeCoeffsExit_m),
64  entryFieldValues_m(NULL),
65  exitFieldValues_m(NULL),
66  entryFieldAccel_m(NULL),
67  exitFieldAccel_m(NULL),
68  deltaBeginEntry_m(right.deltaBeginEntry_m),
69  deltaEndEntry_m(right.deltaEndEntry_m),
70  polyOrderEntry_m(right.polyOrderEntry_m),
71  deltaBeginExit_m(right.deltaBeginExit_m),
72  deltaEndExit_m(right.deltaEndExit_m),
73  polyOrderExit_m(right.polyOrderExit_m),
74  cosEntranceAngle_m(right.cosEntranceAngle_m),
75  sinEntranceAngle_m(right.sinEntranceAngle_m),
76  tanEntranceAngle_m(right.tanEntranceAngle_m),
77  tanExitAngle_m(right.tanExitAngle_m),
78  beginToEnd_m(right.beginToEnd_m),
79  beginToEnd_lcs_m(right.beginToEnd_lcs_m),
80  toEntranceRegion_m(right.toEntranceRegion_m),
81  toExitRegion_m(right.toExitRegion_m),
82  computeAngleTrafo_m(right.computeAngleTrafo_m),
83  maxAngle_m(right.maxAngle_m),
84  nSlices_m(right.nSlices_m){
85 }
86 
87 Bend2D::Bend2D(const std::string &name):
88  BendBase(name),
89  messageHeader_m(" * "),
90  pusher_m(),
91  designRadius_m(0.0),
92  exitAngle_m(0.0),
93  fieldIndex_m(0.0),
94  startField_m(0.0),
95  endField_m(0.0),
96  widthEntranceFringe_m(0.0),
97  widthExitFringe_m(0.0),
98  reinitialize_m(false),
99  entranceParameter1_m(0.0),
100  entranceParameter2_m(0.0),
101  entranceParameter3_m(0.0),
102  exitParameter1_m(0.0),
103  exitParameter2_m(0.0),
104  exitParameter3_m(0.0),
105  entryFieldValues_m(NULL),
106  exitFieldValues_m(NULL),
107  entryFieldAccel_m(NULL),
108  exitFieldAccel_m(NULL),
109  deltaBeginEntry_m(0.0),
110  deltaEndEntry_m(0.0),
111  polyOrderEntry_m(0),
112  deltaBeginExit_m(0.0),
113  deltaEndExit_m(0.0),
114  polyOrderExit_m(0),
115  cosEntranceAngle_m(1.0),
116  sinEntranceAngle_m(0.0),
117  tanEntranceAngle_m(0.0),
118  tanExitAngle_m(0.0),
119  maxAngle_m(0.0),
120  nSlices_m(1){
121 }
122 
124  if (entryFieldAccel_m != NULL) {
125  for (unsigned int i = 0; i < 3u; ++ i) {
126  gsl_spline_free(entryFieldValues_m[i]);
127  gsl_spline_free(exitFieldValues_m[i]);
128  entryFieldValues_m[i] = NULL;
129  exitFieldValues_m[i] = NULL;
130  }
131  delete[] entryFieldValues_m;
132  delete[] exitFieldValues_m;
133 
134  gsl_interp_accel_free(entryFieldAccel_m);
135  gsl_interp_accel_free(exitFieldAccel_m);
136  }
137 }
138 /*
139  * OPAL-T Methods.
140  * ===============
141  */
142 
143 /*
144  * This function merely repackages the field arrays as type Vector_t and calls
145  * the equivalent method but with the Vector_t data types.
146  */
147 bool Bend2D::apply(const size_t &i,
148  const double &t,
149  Vector_t &E,
150  Vector_t &B) {
151  Vector_t P;
152  return apply(RefPartBunch_m->R[i], P, t, E, B);
153 }
154 
156  const Vector_t &/*P*/,
157  const double &/*t*/,
158  Vector_t &/*E*/,
159  Vector_t &B) {
160 
161  if (designRadius_m > 0.0) {
162 
163  Vector_t X = R;
164  Vector_t bField(0.0);
165  if (!calculateMapField(X, bField)) {
166 
167  B += fieldAmplitude_m * bField;
168 
169  return false;
170  }
171 
172  return true;
173  }
174  return false;
175 
176 }
177 
179  const Vector_t &P,
180  const double &t,
181  Vector_t &E,
182  Vector_t &B) {
183  return apply(R, P, t, E, B);
184 }
185 
186 void Bend2D::goOnline(const double &) {
187  online_m = true;
188 }
189 
191  double &startField,
192  double &endField) {
193 
194  Inform msg(messageHeader_m.c_str(), *gmsg);
195 
196  if (initializeFieldMap()) {
197  std::string name = Util::toUpper(getName()) + " ";
198  msg << level2
199  << "======================================================================\n"
200  << "===== " << std::left << std::setw(64) << std::setfill('=') << name << "\n"
201  << "======================================================================\n";
202 
203  setupPusher(bunch);
204  readFieldMap(msg);
205  setupBendGeometry(startField, endField);
206 
207  double bendAngleX = 0.0;
208  double bendAngleY = 0.0;
209  calculateRefTrajectory(bendAngleX, bendAngleY);
210  print(msg, bendAngleX, bendAngleY);
212  // Pass start and end of field to calling function.
213  startField = startField_m;
214  endField = endField_m;
215 
218  elementEdge_m = 0.0;
219 
221  msg << level2
222  << "======================================================================\n"
223  << "======================================================================\n"
224  << endl;
225  } else {
226  ERRORMSG("There is something wrong with your field map \""
227  << fileName_m
228  << "\"");
229  endField = startField - 1.0e-3;
230  }
231 }
232 
233 void Bend2D::adjustFringeFields(double ratio) {
235 
237  entranceParameter1_m = entranceParameter2_m - delta * ratio;
238 
240  entranceParameter3_m = entranceParameter2_m + delta * ratio;
241 
243  exitParameter1_m = exitParameter2_m - delta * ratio;
244 
246  exitParameter3_m = exitParameter2_m + delta * ratio;
247 
249 }
250 
252 
253  const double betaGamma = calcBetaGamma();
254  // const double beta = betaGamma / gamma;
255  const double deltaT = RefPartBunch_m->getdT();
256  const double cdt = Physics::c * deltaT;
257  // Integrate through field for initial angle.
258  Vector_t oldX;
261  double deltaS = 0.0;
262  double bendLength = endField_m - startField_m;
263  const Vector_t eField(0.0);
264 
265  while (deltaS < bendLength) {
266  oldX = X;
267  X /= cdt;
268  pusher_m.push(X, P, deltaT);
269  X *= cdt;
270 
271  Vector_t bField(0.0, 0.0, 0.0);
272  calculateMapField(X, bField);
273  bField = fieldAmplitude_m * bField;
274 
275  X /= cdt;
276  pusher_m.kick(X, P, eField, bField, deltaT);
277 
278  pusher_m.push(X, P, deltaT);
279  X *= cdt;
280 
281  deltaS += euclidean_norm(X - oldX);
282 
283  }
284 
285  double angle = - std::atan2(P(0), P(2)) + entranceAngle_m;
286 
287  return angle;
288 }
289 
290 void Bend2D::calcEngeFunction(double zNormalized,
291  const std::vector<double> &engeCoeff,
292  int polyOrder,
293  double &engeFunc,
294  double &engeFuncDeriv,
295  double &engeFuncSecDerivNorm) {
296 
297  double expSum = 0.0;
298  double expSumDeriv = 0.0;
299  double expSumSecDeriv = 0.0;
300 
301  if (polyOrder >= 2) {
302 
303  expSum = engeCoeff.at(0)
304  + engeCoeff.at(1) * zNormalized;
305  expSumDeriv = engeCoeff.at(1);
306 
307  for (int index = 2; index <= polyOrder; index++) {
308  expSum += engeCoeff.at(index) * std::pow(zNormalized, index);
309  expSumDeriv += index * engeCoeff.at(index)
310  * std::pow(zNormalized, index - 1);
311  expSumSecDeriv += index * (index - 1) * engeCoeff.at(index)
312  * std::pow(zNormalized, index - 2);
313  }
314 
315  } else if (polyOrder == 1) {
316 
317  expSum = engeCoeff.at(0)
318  + engeCoeff.at(1) * zNormalized;
319  expSumDeriv = engeCoeff.at(1);
320 
321  } else
322  expSum = engeCoeff.at(0);
323 
324  double engeExp = std::exp(expSum);
325  engeFunc = 1.0 / (1.0 + engeExp);
326 
327  if (!std::isnan(engeFunc)) {
328 
329  expSumDeriv /= gap_m;
330  expSumSecDeriv /= std::pow(gap_m, 2.0);
331  double engeExpDeriv = expSumDeriv * engeExp;
332  double engeExpSecDeriv = (expSumSecDeriv + std::pow(expSumDeriv, 2.0))
333  * engeExp;
334  double engeFuncSq = std::pow(engeFunc, 2.0);
335 
336  engeFuncDeriv = -engeExpDeriv * engeFuncSq;
337  if (std::isnan(engeFuncDeriv) || std::isinf(engeFuncDeriv))
338  engeFuncDeriv = 0.0;
339 
340  engeFuncSecDerivNorm = -engeExpSecDeriv * engeFunc
341  + 2.0 * std::pow(engeExpDeriv, 2.0)
342  * engeFuncSq;
343  if (std::isnan(engeFuncSecDerivNorm) || std::isinf(engeFuncSecDerivNorm))
344  engeFuncSecDerivNorm = 0.0;
345 
346  } else {
347  engeFunc = 0.0;
348  engeFuncDeriv = 0.0;
349  engeFuncSecDerivNorm = 0.0;
350 
351  }
352 }
353 
354 // :FIXME: is this correct?
355 Vector_t Bend2D::calcCentralField(const Vector_t &/*R*/, double /*deltaX*/) {
356 
357  Vector_t B(0, 0, 0);
358  //double nOverRho = fieldIndex_m / designRadius_m;
359  //double expFactor = exp(-nOverRho * deltaX);
360  //double bxBzFactor = expFactor * nOverRho * R(1);
361  //Vector_t rotationCenter(-designRadius_m, R(1), 0.0);
362  //double cosangle = dot(R - rotationCenter, Vector_t(1, 0, 0)) / euclidean_norm(R - rotationCenter);
363 
364  //B(0) = -bxBzFactor * cosangle;
365  //B(1) = expFactor * (1.0 - std::pow(nOverRho * R(1), 2.0) / 2.0);
366  //B(2) = -bxBzFactor * std::sqrt(1 - std::pow(cosangle, 2));
367 
368  B(1) = 1.0;
369 
370  return B;
371 }
372 
374  double /*deltaX*/) {
375 
376  const CoordinateSystemTrafo toEntranceRegion(Vector_t(0, 0, entranceParameter2_m),
377  Quaternion(0, 0, 1, 0));
378  const Vector_t Rprime = toEntranceRegion.transformTo(R);
379 
380  Vector_t B(0.0);
381 
382  if (Rprime(2) <= -entranceParameter3_m) {
383  B(1) = 1;
384  } else if (Rprime(2) < -entranceParameter1_m) {
385  double engeFunc = gsl_spline_eval(entryFieldValues_m[0], Rprime(2), entryFieldAccel_m);
386  double engeFuncDeriv = gsl_spline_eval(entryFieldValues_m[1], Rprime(2), entryFieldAccel_m);
387  double engeFuncSecDeriv = gsl_spline_eval(entryFieldValues_m[2], Rprime(2), entryFieldAccel_m);
388 
389  //double nOverRho = fieldIndex_m / designRadius_m;
390  //double expFactor = exp(-nOverRho * deltaX);
391  //double trigFactor = pow(nOverRho, 2.0) + engeFuncSecDerivNorm;
392 
393  //double bXEntrance = -nOverRho * expFactor * Rprime(1) * engeFunc;
394  //double bYEntrance = (expFactor * engeFunc *
395  // (1.0 - 0.5 * trigFactor * pow(Rprime(1), 2.0)));
396  //double bZEntrance = expFactor * Rprime(1) * engeFuncDeriv;
397 
398 
399  // B(1) = (engeFunc *
400  // (1.0 - 0.5 * engeFuncSecDerivNorm * pow(Rprime(1), 2.0)));
401 
402  B(1) = (engeFunc - 0.5 * engeFuncSecDeriv * std::pow(Rprime(1), 2.0));
403 
404  B(2) = engeFuncDeriv * Rprime(1);
405  }
406 
407  return toEntranceRegion.rotateFrom(B);
408 }
409 
411  double /*deltaX*/) {
412 
413  const CoordinateSystemTrafo fromEndToExitRegion(Vector_t(0, 0, exitParameter2_m),
414  Quaternion(1, 0, 0, 0));
415  const CoordinateSystemTrafo toExitRegion = (fromEndToExitRegion *
417  const Vector_t Rprime = toExitRegion.transformTo(R);
418 
419  Vector_t B(0.0);
420 
421  if (Rprime(2) <= exitParameter1_m) {
422  B(1) = 1;
423  } else if (Rprime(2) < exitParameter3_m) {
424  double engeFunc = gsl_spline_eval(exitFieldValues_m[0], Rprime(2), exitFieldAccel_m);
425  double engeFuncDeriv = gsl_spline_eval(exitFieldValues_m[1], Rprime(2), exitFieldAccel_m);
426  double engeFuncSecDeriv = gsl_spline_eval(exitFieldValues_m[2], Rprime(2), exitFieldAccel_m);
427 
428  //double nOverRho = fieldIndex_m / designRadius_m;
429  //double expFactor = exp(-nOverRho * deltaX);
430  //double trigFactor = pow(nOverRho, 2.0) + engeFuncSecDerivNorm;
431 
432  //double bXExit = -nOverRho * expFactor * Rprime(1) * engeFunc;
433  //double bYExit = (expFactor * engeFunc *
434  // (1.0 - 0.5 * trigFactor * pow(Rprime(1), 2.0)));
435  //double bZExit = expFactor * Rprime(1) * engeFuncDeriv;
436 
437  //B(1) = (engeFunc *
438  // (1.0 - 0.5 * engeFuncSecDerivNorm * pow(Rprime(1), 2.0)));
439  B(1) = (engeFunc - 0.5 * engeFuncSecDeriv * std::pow(Rprime(1), 2.0));
440 
441  B(2) = engeFuncDeriv * Rprime(1);
442  }
443  return toExitRegion.rotateFrom(B);
444 }
445 
447 
448  B = Vector_t(0.0);
449  bool verticallyInside = (std::abs(R(1)) < 0.5 * gap_m);
450  bool horizontallyInside = false;
452  if (inMagnetCentralRegion(R)) {
453  if (verticallyInside) {
454  double deltaX = 0.0;//euclidean_norm(R - rotationCenter) - designRadius_m;
455  bool inEntranceRegion = isPositionInEntranceField(R);
456  bool inExitRegion = isPositionInExitField(R);
457 
458  if (inEntranceRegion && inExitRegion) {
461 
462  if (std::abs(Rp[2]) < std::abs(Rpp[2])) {
463  inExitRegion = false;
464  } else {
465  inEntranceRegion = false;
466  }
467  }
468  if (inEntranceRegion) {
469  B = calcEntranceFringeField(R, deltaX);
470  } else if (inExitRegion) {
471  B = calcExitFringeField(R, deltaX);
472  } else {
473  B = calcCentralField(R, deltaX);
474  }
475 
476  return false;
477  }
478  return true;
479  }
480 
481  Vector_t BEntrance(0.0), BExit(0.0);
482  verticallyInside = (std::abs(R(1)) < gap_m);
483 
484  bool inEntranceRegion = inMagnetEntranceRegion(R);
485  bool inExitRegion = inMagnetExitRegion(R);
486 
487  if (inEntranceRegion) {
488  horizontallyInside = true;
489  if (verticallyInside) {
490  BEntrance = calcEntranceFringeField(R, R(0));
491  }
492  }
493 
494  if (inExitRegion) {
495  horizontallyInside = true;
496  if (verticallyInside) {
498 
499  BExit = calcExitFringeField(R, Rprime(0));
500  }
501  }
502 
503  B = BEntrance + BExit;
504 
505  bool hitMaterial = (horizontallyInside && (!verticallyInside));
506  return hitMaterial;
507 }
508 
509 void Bend2D::calculateRefTrajectory(double &angleX, double &/*angleY*/) {
510 
511  std::ofstream trajectoryOutput;
513  std::string fname = Util::combineFilePath({
515  OpalData::getInstance()->getInputBasename() + "_" + getName() + "_traj.dat"
516  });
517  trajectoryOutput.open(fname);
518  trajectoryOutput.precision(12);
519  trajectoryOutput << "# " << std::setw(18) << "s"
520  << std::setw(20) << "x"
521  << std::setw(20) << "z"
522  << std::setw(20) << "By"
523  << std::endl;
524  }
525 
526  const double gamma = calcGamma();
527  const double betaGamma = calcBetaGamma();
530 
531  if (!refTrajMap_m.empty())
532  refTrajMap_m.clear();
533 
534  refTrajMap_m.push_back(X);
535 
536  const Vector_t eField(0.0);
537  const double dt = RefPartBunch_m->getdT();
538  const Vector_t scaleFactor(Physics::c * dt);
539  const double stepSize = betaGamma / gamma * Physics::c * dt;
540  const double bendLength = endField_m - startField_m;
541  double deltaS = 0.0;
542  while (deltaS < bendLength) {
543 
544  X /= scaleFactor;
545  pusher_m.push(X, P, dt);
546  X *= scaleFactor;
547 
548  Vector_t bField(0.0, 0.0, 0.0);
549  Vector_t XInBendFrame = X;
550 
551  calculateMapField(XInBendFrame, bField);
552  bField = fieldAmplitude_m * bField;
553 
555  trajectoryOutput << std::setw(20) << deltaS + 0.5 * stepSize
556  << std::setw(20) << X(0)
557  << std::setw(20) << X(2)
558  << std::setw(20) << bField(1)
559  << std::endl;
560  }
561 
562  X /= scaleFactor;
563  pusher_m.kick(X, P, eField, bField, dt);
564 
565  pusher_m.push(X, P, dt);
566  X *= scaleFactor;
567 
568  refTrajMap_m.push_back(X);
569 
570  deltaS += stepSize;
571  }
572 
573  angleX = -std::atan2(P(0), P(2)) + entranceAngle_m;
574 }
575 
576 double Bend2D::estimateFieldAdjustmentStep(double actualBendAngle) {
577 
578  // Estimate field adjustment step.
579  double effectiveLength = angle_m * designRadius_m;
580  double effectiveFieldAmplitude = calcFieldAmplitude(2.0 * effectiveLength);
581  double ratioSquared = std::pow(fieldAmplitude_m / effectiveFieldAmplitude, 2.0);
582  double fieldStep = (angle_m - actualBendAngle) * effectiveFieldAmplitude;
583  if (ratioSquared < 1.0) {
584  fieldStep = (angle_m - actualBendAngle) * effectiveFieldAmplitude * std::sqrt(1.0 - ratioSquared);
585  }
587 
588  return fieldStep;
589 
590 }
591 
592 void Bend2D::findBendEffectiveLength(double startField, double endField) {
593 
594  /*
595  * Use an iterative procedure to set the width of the
596  * default field map for the defined field amplitude
597  * and bend angle.
598  */
599  setEngeOriginDelta(0.0);
601  setFieldBoundaries(startField, endField);
602 
603  double actualBendAngle = calculateBendAngle();
604  double error = std::abs(actualBendAngle - angle_m);
605 
606  if (error > 1.0e-6) {
607 
608  double deltaStep = 0.0;
609  if (std::abs(actualBendAngle) < std::abs(angle_m))
610  deltaStep = -gap_m / 2.0;
611  else
612  deltaStep = gap_m / 2.0;
613 
614  double delta1 = 0.0;
615  double bendAngle1 = actualBendAngle;
616 
617  double delta2 = deltaStep;
618  setEngeOriginDelta(delta2);
620  setFieldBoundaries(startField, endField);
621  double bendAngle2 = calculateBendAngle();
622 
623  if (std::abs(bendAngle1) > std::abs(angle_m)) {
624  while (std::abs(bendAngle2) > std::abs(angle_m)) {
625  delta2 += deltaStep;
626  setEngeOriginDelta(delta2);
628  setFieldBoundaries(startField, endField);
629  bendAngle2 = calculateBendAngle();
630  }
631  } else {
632  while (std::abs(bendAngle2) < std::abs(angle_m)) {
633  delta2 += deltaStep;
634  setEngeOriginDelta(delta2);
636  setFieldBoundaries(startField, endField);
637  bendAngle2 = calculateBendAngle();
638  }
639  }
640 
641  // Now we should have the proper field map width bracketed.
642  unsigned int iterations = 1;
643  error = std::abs(actualBendAngle - angle_m);
644  while (error > 1.0e-6 && iterations < 100) {
645 
646  double delta = (delta1 + delta2) / 2.0;
647  setEngeOriginDelta(delta);
649  setFieldBoundaries(startField, endField);
650  double newBendAngle = calculateBendAngle();
651 
652  error = std::abs(newBendAngle - angle_m);
653 
654  if (error > 1.0e-6) {
655 
656  if (bendAngle1 - angle_m < 0.0) {
657 
658  if (newBendAngle - angle_m < 0.0) {
659  bendAngle1 = newBendAngle;
660  delta1 = delta;
661  } else {
662  // bendAngle2 = newBendAngle;
663  delta2 = delta;
664  }
665 
666  } else {
667 
668  if (newBendAngle - angle_m < 0.0) {
669  // bendAngle2 = newBendAngle;
670  delta2 = delta;
671  } else {
672  bendAngle1 = newBendAngle;
673  delta1 = delta;
674  }
675  }
676  }
677  iterations++;
678  }
679  }
680 }
681 
683 
684  /*
685  * Use an iterative procedure to set the magnet field amplitude
686  * for the defined bend angle.
687  */
688 
689  const double tolerance = 1e-7; // [deg]
690  double actualBendAngle = calculateBendAngle();
691  double error = std::abs(actualBendAngle - angle_m) * Physics::rad2deg;
692  if (error < tolerance)
693  return;
694 
695  double fieldStep = std::copysign(1.0, fieldAmplitude_m) * std::abs(estimateFieldAdjustmentStep(actualBendAngle));
696  double amplitude1 = fieldAmplitude_m;
697  double amplitude2 = amplitude1;
698  double bendAngle1 = actualBendAngle;
699  double bendAngle2 = bendAngle1;
700 
701  int stepSign = std::abs(bendAngle1) > std::abs(angle_m) ? -1: 1;
702  while (true) {
703  amplitude1 = amplitude2;
704  bendAngle1 = bendAngle2;
705 
706  amplitude2 += stepSign * fieldStep;
707  fieldAmplitude_m = amplitude2;
708  bendAngle2 = calculateBendAngle();
709  if (stepSign * (std::abs(bendAngle2) - std::abs(angle_m)) > 0) {
710  break;
711  }
712  }
713 
714  // Now we should have the proper field amplitude bracketed.
715  unsigned int iterations = 1;
716  while (error > tolerance && iterations < 100) {
717 
718  fieldAmplitude_m = (amplitude1 + amplitude2) / 2.0;
719  double newBendAngle = calculateBendAngle();
720 
721  error = std::abs(newBendAngle - angle_m) * Physics::rad2deg;
722 
723  if (error > tolerance) {
724 
725  if (bendAngle1 - angle_m < 0.0) {
726 
727  if (newBendAngle - angle_m < 0.0) {
728  bendAngle1 = newBendAngle;
729  amplitude1 = fieldAmplitude_m;
730  } else {
731  // bendAngle2 = newBendAngle;
732  amplitude2 = fieldAmplitude_m;
733  }
734 
735  } else {
736 
737  if (newBendAngle - angle_m < 0.0) {
738  // bendAngle2 = newBendAngle;
739  amplitude2 = fieldAmplitude_m;
740  } else {
741  bendAngle1 = newBendAngle;
742  amplitude1 = fieldAmplitude_m;
743  }
744  }
745  }
746  iterations++;
747  }
748 }
749 
750 bool Bend2D::findIdealBendParameters(double chordLength) {
751 
752  bool reinitialize = false;
753 
754  if (angle_m != 0.0) {
755 
756  if (angle_m < 0.0) {
757  // Negative angle is a positive bend rotated 180 degrees.
762  }
763  designRadius_m = calcDesignRadius(chordLength, angle_m);
765  reinitialize = true;
766  } else {
767 
769  double refCharge = RefPartBunch_m->getQ();
770  if (refCharge < 0.0) {
772  }
773 
774  fieldAmplitude_m = std::copysign(1.0, refCharge) * std::hypot(fieldAmplitudeX_m, fieldAmplitudeY_m);
776  angle_m = calcBendAngle(chordLength, designRadius_m);
777  reinitialize = false;
778  }
779  return reinitialize;
780 }
781 
783 
785 
786  if (fieldmap_m != NULL) {
787  if (fileName_m != "1DPROFILE1-DEFAULT")
788  return true;
789  else
790  return setupDefaultFieldMap();
791 
792  } else
793  return false;
794 
795 }
796 
798 
800  double distFromRotCenter = euclidean_norm(R - rotationCenter);
803 
804  double effectiveAngle = std::fmod(Physics::two_pi - std::atan2(Rpprime(0), Rpprime(2)), Physics::two_pi);
805 
806  if (std::abs(distFromRotCenter - designRadius_m) < 0.5 * aperture_m.second[0] &&
807  effectiveAngle >= 0.0 && effectiveAngle < maxAngle_m) {
808  if (effectiveAngle < 0.5 * maxAngle_m) return R(2) >= 0.0;
809  return Rprime(2) < 0.0;
810  }
811 
812  return false;
813 }
814 
816 
817  return (R(2) >= entranceParameter1_m &&
818  R(2) < 0.0 &&
819  std::abs(R(0)) < aperture_m.second[0]);
820 }
821 
823 
825 
826  return (Rprime(2) >= 0 &&
827  Rprime(2) < exitParameter3_m &&
828  std::abs(Rprime(0)) < aperture_m.second[0]);
829 }
830 
832 
833  return (polyOrderEntry_m >= 0 &&
834  R(2) >= entranceParameter1_m &&
835  R(2) < entranceParameter3_m);
836 }
837 
840 
841  return (polyOrderExit_m >= 0 &&
842  Rprime(2) >= exitParameter1_m &&
843  Rprime(2) < exitParameter3_m);
844 }
845 
846 void Bend2D::print(Inform &msg, double bendAngleX, double bendAngleY) {
847  msg << level2 << "\n"
848  << "Start of field map: "
849  << startField_m
850  << " m (in s coordinates)"
851  << "\n";
852  msg << "End of field map: "
853  << endField_m
854  << " m (in s coordinates)"
855  << "\n";
856  msg << "Entrance edge of magnet: "
857  << elementEdge_m
858  << " m (in s coordinates)"
859  << "\n";
860  msg << "\n"
861  << "Reference Trajectory Properties"
862  << "\n"
863  << "======================================================================"
864  << "\n\n";
865  msg << "Bend angle magnitude: "
866  << angle_m
867  << " rad ("
869  << " degrees)"
870  << "\n";
871  msg << "Entrance edge angle: "
872  << entranceAngle_m
873  << " rad ("
875  << " degrees)"
876  << "\n";
877  msg << "Exit edge angle: "
878  << exitAngle_m
879  << " rad ("
881  << " degrees)"
882  << "\n";
883  msg << "Bend design radius: "
884  << designRadius_m
885  << " m"
886  << "\n";
887  msg << "Bend design energy: "
888  << designEnergy_m
889  << " eV"
890  << "\n";
891  msg << "\n"
892  << "Bend Field and Rotation Properties"
893  << "\n"
894  << "======================================================================"
895  << "\n\n";
896  msg << "Field amplitude: "
898  << " T"
899  << "\n";
900  msg << "Field index: "
901  << fieldIndex_m
902  << "\n";
903  msg << "Rotation about z axis: "
904  << rotationZAxis_m
905  << " rad ("
907  << " degrees)"
908  << "\n";
909  msg << "\n"
910  << "Reference Trajectory Properties Through Bend Magnet with Fringe Fields"
911  << "\n"
912  << "======================================================================"
913  << "\n\n";
914  msg << "Reference particle is bent: "
915  << bendAngleX
916  << " rad ("
917  << bendAngleX * Physics::rad2deg
918  << " degrees) in x plane"
919  << "\n";
920  msg << "Reference particle is bent: "
921  << bendAngleY
922  << " rad ("
923  << bendAngleY * Physics::rad2deg
924  << " degrees) in y plane"
925  << "\n";
926 
927  if (fileName_m == "1DPROFILE1-DEFAULT") {
928  msg << "\n"
929  << "Effective Field Map\n"
930  << "======================================================================\n"
931  << "\n"
932  << "1DProfile1 " << polyOrderEntry_m << " " << polyOrderExit_m << " " << gap_m * 100 << "\n"
933  << entranceParameter1_m * 100 << " " << entranceParameter2_m * 100 << " " << entranceParameter3_m * 100 << " 0\n"
934  << exitParameter1_m * 100 << " " << exitParameter2_m * 100 << " " << exitParameter3_m * 100 << " 0\n";
935  for (auto coef: engeCoeffsEntry_m) {
936  msg << coef << "\n";
937  }
938  for (auto coef: engeCoeffsExit_m) {
939  msg << coef << "\n";
940  }
941  }
942  msg << endl;
943 }
944 
945 
947  msg << level2 << getName() << " using file ";
948  fieldmap_m->getInfo(&msg);
949 
957 
961  polyOrderEntry_m = engeCoeffsEntry_m.size() - 1;
962  polyOrderExit_m = engeCoeffsExit_m.size() - 1;
963 
964  double stepSize = Physics::c * 1e-12;
965  double entryLength = std::abs(entranceParameter3_m - entranceParameter1_m);
966  unsigned int numStepsEntry = std::ceil(entryLength / stepSize) + 3;
967  double stepSizeEntry = entryLength / (numStepsEntry - 3);
968  std::vector<double> zvalsEntry(numStepsEntry);
969  std::vector<std::vector<double> > fieldValuesEntry(3);
970 
971  double exitLength = std::abs(exitParameter3_m - exitParameter1_m);
972  unsigned int numStepsExit = std::ceil(exitLength / stepSize) + 3;
973  double stepSizeExit = exitLength / (numStepsExit - 3);
974  std::vector<double> zvalsExit(numStepsExit);
975  std::vector<std::vector<double> > fieldValuesExit(3);
976 
977  entryFieldValues_m = new gsl_spline*[3];
978  exitFieldValues_m = new gsl_spline*[3];
979  entryFieldAccel_m = gsl_interp_accel_alloc();
980  exitFieldAccel_m = gsl_interp_accel_alloc();
981 
982  for (unsigned int i = 0; i < 3u; ++ i) {
983  fieldValuesEntry[i].resize(numStepsEntry);
984  fieldValuesExit[i].resize(numStepsExit);
985 
986  entryFieldValues_m[i] = gsl_spline_alloc(gsl_interp_cspline, numStepsEntry);
987  exitFieldValues_m[i] = gsl_spline_alloc(gsl_interp_cspline, numStepsExit);
988  }
989 
990  for (unsigned int j = 0; j < numStepsEntry; ++ j) {
991  if (j == 0) {
992  zvalsEntry[j] = -entranceParameter3_m - stepSizeEntry;
993  } else {
994  zvalsEntry[j] = zvalsEntry[j - 1] + stepSizeEntry;
995  }
996  calcEngeFunction(zvalsEntry[j] / gap_m,
999  fieldValuesEntry[0][j],
1000  fieldValuesEntry[1][j],
1001  fieldValuesEntry[2][j]);
1002  fieldValuesEntry[2][j] *= fieldValuesEntry[0][j];
1003  }
1004 
1005  for (unsigned int j = 0; j < numStepsExit; ++ j) {
1006  if (j == 0) {
1007  zvalsExit[j] = exitParameter1_m - stepSizeExit;
1008  } else {
1009  zvalsExit[j] = zvalsExit[j - 1] + stepSizeExit;
1010  }
1011  calcEngeFunction(zvalsExit[j] / gap_m,
1014  fieldValuesExit[0][j],
1015  fieldValuesExit[1][j],
1016  fieldValuesExit[2][j]);
1017  fieldValuesExit[2][j] *= fieldValuesExit[0][j];
1018  }
1019 
1020  for (unsigned int i = 0; i < 3u; ++ i) {
1021  gsl_spline_init(entryFieldValues_m[i], &(zvalsEntry[0]), &(fieldValuesEntry[i][0]), numStepsEntry);
1022  gsl_spline_init(exitFieldValues_m[i], &(zvalsExit[0]), &(fieldValuesExit[i][0]), numStepsExit);
1023  }
1024 }
1025 
1026 void Bend2D::setBendEffectiveLength(double startField, double endField) {
1027 
1028  // Find initial angle.
1029  double actualBendAngle = calculateBendAngle();
1030 
1031  // Adjust field map to match bend angle.
1032  double error = std::abs(actualBendAngle - angle_m);
1033  if (error > 1.0e-6)
1034  findBendEffectiveLength(startField, endField);
1035 
1036 }
1037 
1039 
1040  // Estimate bend field magnitude.
1042 
1043  // Search for angle if initial guess is not good enough.
1044  findBendStrength();
1045 }
1046 
1047 void Bend2D::setEngeOriginDelta(double delta) {
1048  /*
1049  * This function is used to shift the perpendicular distance of the
1050  * entrance and exit Enge function origins with respect to the entrance
1051  * and exit points in the magnet. A positive delta shifts them towards
1052  * the center of the magnet.
1053  */
1055 
1060  entranceParameter2_m = delta;
1061 
1064  exitParameter2_m = -delta;
1065 
1067 }
1068 
1070 
1073 
1076 
1077  double bendAngle = getBendAngle();
1078  double entranceAngle = getEntranceAngle();
1079 
1080  double rotationAngleAboutZ = getRotationAboutZ();
1081  Quaternion_t rotationAboutZ(std::cos(0.5 * rotationAngleAboutZ),
1082  std::sin(0.5 * rotationAngleAboutZ) * Vector_t(0, 0, 1));
1083 
1084  Vector_t rotationAxis(0, -1, 0);
1085  Quaternion_t halfRotationAboutAxis(std::cos(0.5 * (0.5 * bendAngle - entranceAngle)),
1086  std::sin(0.5 * (0.5 * bendAngle - entranceAngle)) * rotationAxis);
1087  Quaternion_t exitFaceRotation(std::cos(0.5 * (bendAngle - entranceAngle - exitAngle_m)),
1088  std::sin(0.5 * (bendAngle - entranceAngle - exitAngle_m)) * rotationAxis);
1089  Vector_t chord = getChordLength() * halfRotationAboutAxis.rotate(Vector_t(0, 0, 1));
1090  beginToEnd_lcs_m = CoordinateSystemTrafo(chord, exitFaceRotation.conjugate());
1093  Quaternion(0, 0, 1, 0));
1094  const CoordinateSystemTrafo fromEndToExitRegion(Vector_t(0, 0, exitParameter2_m),
1095  Quaternion(1, 0, 0, 0));
1096  toExitRegion_m = CoordinateSystemTrafo(fromEndToExitRegion *
1098 
1100 
1101  Vector_t maxAngleEntranceAperture;
1102  if (rotationCenter(2) < 0.0) {
1103  double tau, tmp;
1104  Vector_t P(0.5 * aperture_m.second[0], 0.0, 0.0), &R = rotationCenter;
1105  gsl_poly_solve_quadratic(dot(P,P),
1106  -2 * dot(P,R),
1107  dot(R,R) - std::pow(designRadius_m + 0.5 * aperture_m.second[0], 2),
1108  &tau,
1109  &tmp);
1110  tau = (std::abs(1.0 - tmp) < std::abs(1.0 - tau)? tmp: tau);
1111  maxAngleEntranceAperture = tau * P;
1112  } else {
1113  double tau, tmp;
1114  Vector_t P(-0.5 * aperture_m.second[0], 0.0, 0.0), &R = rotationCenter;
1115  gsl_poly_solve_quadratic(dot(P,P),
1116  -2 * dot(P,R),
1117  dot(R,R) - std::pow(designRadius_m - 0.5 * aperture_m.second[0], 2),
1118  &tau,
1119  &tmp);
1120  tau = (std::abs(1.0 - tmp) < std::abs(1.0 - tau)? tmp: tau);
1121  maxAngleEntranceAperture = tau * P;
1122  }
1123  Vector_t maxAngleExitAperture;
1124  if (getBeginToEnd_local().transformTo(rotationCenter)(2) > 0.0) {
1125  double tau, tmp;
1126  Vector_t P(0.5 * aperture_m.second[0], 0.0, 0.0), R = getBeginToEnd_local().transformTo(rotationCenter);
1127  gsl_poly_solve_quadratic(dot(P,P),
1128  -2 * dot(P,R),
1129  dot(R,R) - std::pow(designRadius_m + 0.5 * aperture_m.second[0], 2),
1130  &tau,
1131  &tmp);
1132  tau = (std::abs(1.0 - tmp) < std::abs(1.0 - tau)? tmp: tau);
1133  maxAngleExitAperture = getBeginToEnd_local().transformFrom(tau * P);
1134  } else {
1135  double tau, tmp;
1136  Vector_t P(-0.5 * aperture_m.second[0], 0.0, 0.0), R = getBeginToEnd_local().transformTo(rotationCenter);
1137  gsl_poly_solve_quadratic(dot(P,P),
1138  -2 * dot(P,R),
1139  dot(R,R) - std::pow(designRadius_m - 0.5 * aperture_m.second[0], 2),
1140  &tau,
1141  &tmp);
1142  tau = (std::abs(1.0 - tmp) < std::abs(1.0 - tau)? tmp: tau);
1143  maxAngleExitAperture = getBeginToEnd_local().transformFrom(tau * P);
1144  }
1145 
1146  maxAngleEntranceAperture -= rotationCenter;
1147 
1148  Quaternion rotation = getQuaternion(maxAngleEntranceAperture, Vector_t(0, 0, 1));
1149  computeAngleTrafo_m = CoordinateSystemTrafo(rotationCenter,
1150  rotation);
1151  Vector_t tmp = computeAngleTrafo_m.transformTo(maxAngleExitAperture);
1153 
1154  maxAngleExitAperture -= rotationCenter;
1155 }
1156 
1158 
1159  if (gap_m <= 0.0)
1161  else if (gap_m != fieldmap_m->getFieldGap())
1163 
1164 }
1165 
1166 bool Bend2D::setupBendGeometry(double &startField, double &endField) {
1167 
1168  chordLength_m = 0.0;
1170  return false;
1171 
1172  if (isFieldZero()) {
1173  // treat as drift
1174  startField_m = startField;
1175  endField_m = startField + chordLength_m;
1176  return true;
1177  }
1178 
1180  if (getType() == RBEND) {
1182  }
1183 
1184  /*
1185  * Set field map geometry.
1186  */
1187  if (aperture_m.second[0] <= 0.0)
1188  aperture_m.second[0] = designRadius_m / 2.0;
1190 
1191  /*
1192  * If we are using the default field map, then the origins for the
1193  * Enge functions will be shifted so that we get the desired bend
1194  * angle for the given bend strength. (We match the effective length
1195  * of our field map to the ideal bend length.)
1196  *
1197  * If we are not using the default field map, we assume it cannot
1198  * change so we either leave everything alone (if the user defines
1199  * the bend strength) or we adjust the bend field to get the right
1200  * angle.
1201  */
1202  elementEdge_m = startField;
1203  setFieldBoundaries(startField, endField);
1204 
1205  if (fileName_m != "1DPROFILE1-DEFAULT") {
1206  if (reinitialize_m)
1207  setBendStrength();
1208  } else {
1209  setBendEffectiveLength(startField, endField);
1210  }
1211 
1212  startField = startField_m;
1213  endField = endField_m;
1214  return true;
1215 
1216 }
1217 
1219 
1220  if (getElementLength() <= 0.0) {
1221  ERRORMSG("If using \"1DPROFILE1-DEFAULT\" field map you must set the "
1222  "chord length of the bend using the L attribute in the OPAL "
1223  "input file."
1224  << endl);
1225  return false;
1226  } else {
1227  // msg << level2;
1228  // fieldmap_m->getInfo(&msg);
1229  return true;
1230  }
1231 
1232 }
1233 
1234 void Bend2D::setFieldBoundaries(double startField, double /*endField*/) {
1235 
1237  endField_m = (startField + angle_m * designRadius_m +
1239 }
1240 
1242 
1243  RefPartBunch_m = bunch;
1244  pusher_m.initialise(bunch->getReference());
1245 
1246 }
1247 
1249  if (designEnergy_m <= 0.0) {
1250  WARNMSG("Warning: bend design energy is zero. Treating as drift."
1251  << endl);
1252  designRadius_m = 0.0;
1253  return true;
1254  } else if (angle_m == 0.0 &&
1255  std::pow(fieldAmplitudeX_m, 2.0) + std::pow(fieldAmplitudeY_m, 2.0) > 0.0) {
1256 
1257  double radius = calcDesignRadius(fieldAmplitude_m);
1258  double sinArgument = chordLength_m / (2.0 * radius);
1259 
1260  if (std::abs(sinArgument) > 1.0) {
1261  WARNMSG("Warning: bend strength and defined reference trajectory "
1262  << "chord length are not consistent. Treating bend as drift."
1263  << endl);
1264  designRadius_m = 0.0;
1265  return true;
1266  } else
1267  return false;
1268 
1269  } else if (angle_m == 0.0) {
1270 
1271  WARNMSG("Warning bend angle/strength is zero. Treating as drift."
1272  << endl);
1273  designRadius_m = 0.0;
1274  return true;
1275 
1276  } else
1277  return false;
1278 }
1279 
1280 std::vector<Vector_t> Bend2D::getOutline() const {
1281  std::vector<Vector_t> outline;
1283  unsigned int numSteps = 2;
1284 
1285  outline.push_back(Vector_t(-0.5 * widthEntranceFringe_m, 0.0, 0.0));
1289  outline.push_back(Vector_t(0.5 * widthEntranceFringe_m, 0.0, 0.0));
1290 
1291  {
1292  double tau1, tau2;
1293  Vector_t P(0.5 * aperture_m.second[0], 0.0, 0.0), R = rotationCenter;
1294  gsl_poly_solve_quadratic(dot(P,P),
1295  -2 * dot(P,R),
1296  dot(R,R) - std::pow(designRadius_m + 0.5 * aperture_m.second[0], 2),
1297  &tau1,
1298  &tau2);
1299  tau1 = (std::abs(1.0 - tau2) < std::abs(1.0 - tau1)? tau2: tau1);
1300  Vector_t upperCornerAtEntry = tau1 * P;
1301 
1302  R = getBeginToEnd_local().transformTo(rotationCenter);
1303  gsl_poly_solve_quadratic(dot(P,P),
1304  -2 * dot(P,R),
1305  dot(R,R) - std::pow(designRadius_m + 0.5 * aperture_m.second[0], 2),
1306  &tau1,
1307  &tau2);
1308  tau1 = (std::abs(1.0 - tau2) < std::abs(1.0 - tau1)? tau2: tau1);
1309  Vector_t upperCornerAtExit = getBeginToEnd_local().transformFrom(tau1 * P);
1310 
1311  Quaternion rotation = getQuaternion(upperCornerAtEntry - rotationCenter, Vector_t(0,0,1));
1312  Vector_t tmp = CoordinateSystemTrafo(rotationCenter, rotation).transformTo(upperCornerAtExit);
1313  double totalAngle = -std::fmod(Physics::two_pi - std::atan2(tmp(0), tmp(2)), Physics::two_pi);
1314  numSteps = std::max(2.0, std::ceil(-totalAngle / (5.0 * Physics::deg2rad)));
1315  double dAngle = 0.5 * totalAngle / (1.0 * numSteps - 1.0);
1316 
1317  outline.push_back(upperCornerAtEntry);
1318  double angle = 0.0;
1319  for (unsigned int i = 0; i < numSteps; ++ i, angle += dAngle) {
1320  Quaternion rot(std::cos(angle), 0.0, std::sin(angle), 0.0);
1321  outline.push_back(rot.rotate(upperCornerAtEntry - rotationCenter) + rotationCenter);
1322  }
1323  outline.push_back(upperCornerAtExit);
1324  }
1325 
1326  outline.push_back(getBeginToEnd_local().transformFrom(Vector_t(0.5 * widthExitFringe_m, 0.0, 0.0)));
1327  outline.push_back(getBeginToEnd_local().transformFrom(Vector_t(0.5 * widthExitFringe_m - exitParameter3_m * tanExitAngle_m, 0.0, exitParameter3_m)));
1328  outline.push_back(getBeginToEnd_local().transformFrom(Vector_t(-exitParameter3_m * tanExitAngle_m, 0.0, exitParameter3_m)));
1329  outline.push_back(getBeginToEnd_local().transformFrom(Vector_t(-0.5 * widthExitFringe_m - exitParameter3_m * tanExitAngle_m, 0.0, exitParameter3_m)));
1330  outline.push_back(getBeginToEnd_local().transformFrom(Vector_t(-0.5 * widthExitFringe_m, 0.0, 0.0)));
1331 
1332  {
1333  double tau1, tau2;
1334  Vector_t P(-0.5 * aperture_m.second[0], 0.0, 0.0), R = rotationCenter;
1335  gsl_poly_solve_quadratic(dot(P,P),
1336  -2 * dot(P,R),
1337  dot(R,R) - std::pow(designRadius_m - 0.5 * aperture_m.second[0], 2),
1338  &tau1,
1339  &tau2);
1340  tau1 = (std::abs(1.0 - tau2) < std::abs(1.0 - tau1)? tau2: tau1);
1341  Vector_t lowerCornerAtEntry = tau1 * P;
1342 
1343  R = getBeginToEnd_local().transformTo(rotationCenter);
1344  gsl_poly_solve_quadratic(dot(P,P),
1345  -2 * dot(P,R),
1346  dot(R,R) - std::pow(designRadius_m - 0.5 * aperture_m.second[0], 2),
1347  &tau1,
1348  &tau2);
1349  tau1 = (std::abs(1.0 - tau2) < std::abs(1.0 - tau1)? tau2: tau1);
1350  Vector_t lowerCornerAtExit = getBeginToEnd_local().transformFrom(tau1 * P);
1351 
1352  Quaternion rotation = getQuaternion(lowerCornerAtEntry - rotationCenter, Vector_t(0,0,1));
1353  Vector_t tmp = CoordinateSystemTrafo(rotationCenter, rotation).transformTo(lowerCornerAtExit);
1354  double totalAngle = -std::fmod(Physics::two_pi - std::atan2(tmp(0), tmp(2)), Physics::two_pi);
1355  double dAngle = 0.5 * totalAngle / (1.0 * numSteps - 1.0);
1356 
1357  outline.push_back(lowerCornerAtExit);
1358  double angle = 0.5 * totalAngle;
1359  for (unsigned int i = 0; i < numSteps; ++ i, angle -= dAngle) {
1360  Quaternion rot(std::cos(angle), 0.0, std::sin(angle), 0.0);
1361  outline.push_back(rot.rotate(lowerCornerAtEntry - rotationCenter) + rotationCenter);
1362  }
1363  outline.push_back(lowerCornerAtEntry);
1364  }
1365 
1366  std::ofstream outlineOutput;
1368  std::string fname = Util::combineFilePath({
1370  OpalData::getInstance()->getInputBasename() + "_" + getName() + "_outline.dat"
1371  });
1372  outlineOutput.open(fname);
1373  outlineOutput.precision(8);
1374 
1375  for (auto a: outline) {
1376  outlineOutput << std::setw(18) << a(2)
1377  << std::setw(18) << a(0)
1378  << "\n";
1379  }
1380  outlineOutput << std::setw(18) << outline.front()(2)
1381  << std::setw(18) << outline.front()(0)
1382  << std::endl;
1383  }
1384  return outline;
1385 }
1386 
1388  MeshData mesh;
1389  const Vector_t hgap(0, 0.5 * getFullGap(), 0);
1390  std::vector<Vector_t> outline = getOutline();
1391 
1392  unsigned int size = outline.size();
1393  unsigned int last = size - 1;
1394  unsigned int numSteps = (size - 14) / 2;
1395  unsigned int midIdx = numSteps + 7;
1396 
1397  for (unsigned int i = 0; i < 6; ++ i) {
1398  mesh.vertices_m.push_back(outline[i] + 2 * hgap);
1399  }
1400  for (unsigned int i = 6; i < 6 + numSteps; ++ i) {
1401  mesh.vertices_m.push_back(outline[i] + hgap);
1402  }
1403  for (unsigned int i = 6 + numSteps; i < 13 + numSteps; ++ i) {
1404  mesh.vertices_m.push_back(outline[i] + 2 * hgap);
1405  }
1406  for (unsigned int i = 13 + numSteps; i < 13 + 2 * numSteps; ++ i) {
1407  mesh.vertices_m.push_back(outline[i] + hgap);
1408  }
1409  mesh.vertices_m.push_back(outline[last] + 2 * hgap);
1410 
1411  for (unsigned int i = 0; i < 6; ++ i) {
1412  mesh.vertices_m.push_back(outline[i] - 2 * hgap);
1413  }
1414  for (unsigned int i = 6; i < 6 + numSteps; ++ i) {
1415  mesh.vertices_m.push_back(outline[i] - hgap);
1416  }
1417  for (unsigned int i = 6 + numSteps; i < 13 + numSteps; ++ i) {
1418  mesh.vertices_m.push_back(outline[i] - 2 * hgap);
1419  }
1420  for (unsigned int i = 13 + numSteps; i < 13 + 2 * numSteps; ++ i) {
1421  mesh.vertices_m.push_back(outline[i] - hgap);
1422  }
1423  mesh.vertices_m.push_back(outline[last] - 2 * hgap);
1424 
1425  mesh.vertices_m.push_back(outline[0]);
1426  mesh.vertices_m.push_back(outline[4]);
1427  mesh.vertices_m.push_back(outline[midIdx]);
1428  mesh.vertices_m.push_back(outline[midIdx + 4]);
1429 
1430  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(2, 1, 0));
1431  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(2, 4, 3));
1432  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(2, 5, 4));
1433  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(2, 0, last));
1434  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(2, last, 5));
1435  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(last, last - 1, 5));
1436  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(5, last - 1, 6));
1437 
1438  // exit region top
1439  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(midIdx + 2, midIdx + 1, midIdx));
1440  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(midIdx + 2, midIdx + 4, midIdx + 3));
1441  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(midIdx + 2, midIdx + 5, midIdx + 4));
1442  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(midIdx + 2, midIdx, midIdx - 1));
1443  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(midIdx + 2, midIdx - 1, midIdx + 5));
1444  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(midIdx - 1, midIdx - 2, midIdx + 5));
1445  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(midIdx + 5, midIdx - 2, midIdx + 6));
1446 
1447  // entry region bottom
1448  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + 2, size + 0, size + 1));
1449  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + 2, size + 3, size + 4));
1450  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + 2, size + 4, size + 5));
1451  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + 2, size + last, size + 0));
1452  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + 2, size + 5, size + last));
1453  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + last, size + 5, size + last - 1));
1454  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + 5, size + 6, size + last - 1));
1455 
1456  // exit region bottom
1457  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + midIdx + 2, size + midIdx + 0, size + midIdx + 1));
1458  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + midIdx + 2, size + midIdx + 3, size + midIdx + 4));
1459  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + midIdx + 2, size + midIdx + 4, size + midIdx + 5));
1460  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + midIdx + 2, size + midIdx - 1, size + midIdx + 0));
1461  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + midIdx + 2, size + midIdx + 5, size + midIdx - 1));
1462  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + midIdx - 1, size + midIdx + 5, size + midIdx - 2));
1463  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + midIdx + 5, size + midIdx + 6, size + midIdx - 2));
1464 
1465  // central region
1466  for (unsigned int i = 6; i < 5 + numSteps; ++ i) {
1467  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(i, last + 5 - i, i + 1));
1468  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(last + 5 - i, last + 4 - i, i + 1));
1469 
1470  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + i, size + i + 1, size + last + 5 - i));
1471  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + last + 5 - i, size + i + 1, size + last + 4 - i));
1472  }
1473 
1474  // side
1475  for (unsigned int i = 0; i < size - 2; ++ i) {
1476  if (i == 4u || i == 5u ||
1477  i == 5 + numSteps || i == 6 + numSteps ||
1478  i == 11 + numSteps || i == 12 + numSteps) continue;
1479 
1480  unsigned int next = (i + 1) % size;
1481  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(i, next, next + size));
1482  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(i, next + size, i + size));
1483  }
1484 
1485  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(last, 0, last-1));
1486  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(2*size, last - 1, 0));
1487  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(2*size, size + last - 1, last - 1));
1488  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(2*size, size, size + last - 1));
1489  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + last, size + last - 1, size));
1490 
1491  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(4, 5, 6));
1492  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(4, 6, 2*size + 1));
1493  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(2*size + 1, 6, size + 6));
1494  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + 4, 2*size + 1, size + 6));
1495  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + 4, size + 6, size + 5));
1496 
1497  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(6 + numSteps, midIdx, 5 + numSteps));
1498  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(5 + numSteps, midIdx, 2*size + 2));
1499  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(5 + numSteps, 2*size + 2, size + 5 + numSteps));
1500  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + 5 + numSteps, 2*size + 2, size + midIdx));
1501  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + 6 + numSteps, size + 5 + numSteps, size + midIdx));
1502 
1503  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(6 + midIdx, 4 + midIdx, 5 + midIdx));
1504  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(6 + midIdx, 2*size + 3, 4 + midIdx));
1505  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(2*size + 3, 6 + midIdx, size + 6 + midIdx));
1506  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + 4 + midIdx, 2*size + 3, size + 6 + midIdx));
1507  mesh.triangles_m.push_back(Vektor<unsigned int, 3>(size + 4 + midIdx, size + 6 + midIdx, size + 5 + midIdx));
1508 
1510 
1513 
1514  Vector_t T = cross(P1 - rotationCenter, P2 - rotationCenter);
1515  if (T[1] > 0) { // fringe fields are overlapping
1517  if (this->getType() == ElementBase::RBEND ||
1519  mesh.decorations_m.push_back(std::make_pair(0.5 * (P1 + P2) - 0.25 * dir1,
1520  0.5 * (P1 + P2) + 0.25 * dir1));
1521  } else {
1522  Vector_t dir2 = toExitRegion_m.rotateFrom(Vector_t(-1, 0, 0));
1523  Tenzor<double, 3> inv;
1524  double det = -dir1[0] * dir2[2] + dir1[2] * dir2[0];
1525  inv(0, 0) = -dir2[2] / det;
1526  inv(0, 2) = dir2[0] / det;
1527  inv(1,1) = 1.0;
1528  inv(2, 0) = -dir1[2] / det;
1529  inv(2, 2) = dir1[0] / det;
1530  Vector_t Tau = dot(inv, P2 - P1);
1531  Vector_t crossPoint = P1 + Tau[0] * dir1;
1532  double angle = std::asin(cross(dir1, dir2)[1]);
1533  Quaternion halfRot(std::cos(0.25 * angle), std::sin(0.25 * angle) * Vector_t(0, 1, 0));
1534  Vector_t P = halfRot.rotate(dir1);
1535  Vector_t R = crossPoint - rotationCenter;
1536 
1537  double radius = designRadius_m - 0.5 * aperture_m.second[0];
1538  double tau1, tau2;
1539  gsl_poly_solve_quadratic(dot(P,P),
1540  2 * dot(P,R),
1541  dot(R,R) - std::pow(radius, 2),
1542  &tau1,
1543  &tau2);
1544  if (euclidean_norm(crossPoint + tau1 * P - P1) > euclidean_norm(crossPoint + tau2 * P - P1)) {
1545  std::swap(tau1, tau2);
1546  }
1547  Vector_t lowerCornerFringeLimit = crossPoint + tau1 * P;
1548 
1549  radius = designRadius_m + 0.5 * aperture_m.second[0];
1550  gsl_poly_solve_quadratic(dot(P,P),
1551  2 * dot(P,R),
1552  dot(R,R) - std::pow(radius, 2),
1553  &tau1,
1554  &tau2);
1555  if (euclidean_norm(crossPoint + tau1 * P - P1) > euclidean_norm(crossPoint + tau2 * P - P1)) {
1556  std::swap(tau1, tau2);
1557  }
1558  Vector_t upperCornerFringeLimit = crossPoint + tau1 * P;
1559 
1560  mesh.decorations_m.push_back(std::make_pair(lowerCornerFringeLimit,
1561  upperCornerFringeLimit));
1562  }
1563  } else {
1564 
1565  double tau1, tau2;
1566  Vector_t P(-0.5 * aperture_m.second[0], 0.0, 0.0);
1567  Vector_t R = Vector_t(0, 0, entranceParameter3_m) - rotationCenter;
1568  gsl_poly_solve_quadratic(dot(P,P),
1569  2 * dot(P,R),
1570  dot(R,R) - std::pow(designRadius_m - 0.5 * aperture_m.second[0],2),
1571  &tau1,
1572  &tau2);
1573  tau1 = (std::abs(1.0 - tau2) < std::abs(1.0 - tau1)? tau2: tau1);
1574  Vector_t lowerCornerFringeLimitEntrance = Vector_t(0, 0, entranceParameter3_m) + tau1 * P;
1575 
1576  gsl_poly_solve_quadratic(dot(P,P),
1577  -2 * dot(P,R),
1578  dot(R,R) - std::pow(designRadius_m + 0.5 * aperture_m.second[0],2),
1579  &tau1,
1580  &tau2);
1581  tau1 = (std::abs(1.0 - tau2) < std::abs(1.0 - tau1)? tau2: tau1);
1582  Vector_t upperCornerFringeLimitEntrance = Vector_t(0, 0, entranceParameter3_m) - tau1 * P;
1583 
1584  P = Vector_t(0.5 * aperture_m.second[0], 0.0, 0.0);
1585  R = Vector_t(0, 0, exitParameter1_m) - getBeginToEnd_local().transformTo(rotationCenter);
1586  gsl_poly_solve_quadratic(dot(P,P),
1587  2 * dot(P,R),
1588  dot(R,R) - std::pow(designRadius_m + 0.5 * aperture_m.second[0],2),
1589  &tau1,
1590  &tau2);
1591  tau1 = (std::abs(1.0 - tau2) < std::abs(1.0 - tau1)? tau2: tau1);
1592  Vector_t upperCornerFringeLimitExit = getBeginToEnd_local().transformFrom(Vector_t(0, 0, exitParameter1_m) + tau1 * P);
1593 
1594  gsl_poly_solve_quadratic(dot(P,P),
1595  -2 * dot(P,R),
1596  dot(R,R) - std::pow(designRadius_m - 0.5 * aperture_m.second[0],2),
1597  &tau1,
1598  &tau2);
1599  tau1 = (std::abs(1.0 - tau2) < std::abs(1.0 - tau1)? tau2: tau1);
1600  Vector_t lowerCornerFringeLimitExit = getBeginToEnd_local().transformFrom(Vector_t(0, 0, exitParameter1_m) - tau1 * P);
1601 
1602  mesh.decorations_m.push_back(std::make_pair(lowerCornerFringeLimitEntrance, upperCornerFringeLimitEntrance));
1603  mesh.decorations_m.push_back(std::make_pair(lowerCornerFringeLimitExit, upperCornerFringeLimitExit));
1604 
1605  }
1606 
1608  Vector_t(0.0)));
1609  mesh.decorations_m.push_back(std::make_pair(getBeginToEnd_local().transformFrom(Vector_t(0.0)),
1611 
1612  return mesh;
1613 }
1614 
1615 bool Bend2D::isInside(const Vector_t &R) const {
1616  if (std::abs(R(1)) > gap_m) return false;
1617 
1618  if (inMagnetEntranceRegion(R)) {
1619  return true;
1620  }
1621 
1622  if (inMagnetExitRegion(R)) {
1623  return true;
1624  }
1625 
1626  return (std::abs(R(1)) < 0.5 * gap_m && inMagnetCentralRegion(R));
1627 }
1628 
1630 {
1633 }
1634 
1635 //set the number of slices for map tracking
1636 void Bend2D::setNSlices(const std::size_t& nSlices) {
1637  nSlices_m = nSlices;
1638 }
1639 
1640 //get the number of slices for map tracking
1641 std::size_t Bend2D::getNSlices() const {
1642  return nSlices_m;
1643 }
1644 
1646 // Used to create fringe fields in ThickTracker, (before edge[m], after edge[m])
1647 std::array<double,2> Bend2D::getEntranceFringeFieldLength() const {
1648  double entrParam1, entrParam2, entrParam3;
1649 
1651  entrParam2,
1652  entrParam3);
1653 
1654  std::array<double,2> entFFL;
1655  entFFL[0]=entrParam2-entrParam1; //before edge
1656  entFFL[1]=entrParam3-entrParam2; //after edge
1657 
1658 
1659  entFFL[0] = ( entranceParameter2_m - entranceParameter1_m ); //before edge
1660  entFFL[1] = ( entranceParameter3_m - entranceParameter2_m ); //after edge
1661 
1662  return entFFL;
1663 }
1664 
1666 // Used to create fringe fields in ThickTracker, (after edge[m], before edge[m])
1667 std::array<double,2> Bend2D::getExitFringeFieldLength() const {
1668  std::array<double,2> extFFL;
1669 
1670  double exitParam1, exitParam2, exitParam3;
1671 
1672  fieldmap_m->get1DProfile1ExitParam(exitParam1,
1673  exitParam2,
1674  exitParam3);
1675 
1676  extFFL[0]=exitParam3-exitParam2; //after edge
1677  extFFL[1]=exitParam2-exitParam1; //before edge
1678 
1679  extFFL[0] = ( exitParameter3_m-exitParameter2_m ); //after edge
1680  extFFL[1] = ( exitParameter2_m-exitParameter1_m ); //before edge
1681 
1682  return extFFL;
1683 }
1684 
1688 
1689  std::vector<Vector_t> outline = getOutline();
1690 
1691  BoundingBox bb;
1693  bb.upperRightCorner = std::numeric_limits<double>::lowest();
1694 
1695  Vector_t dY(0, 0.5 * getFullGap(), 0);
1696  for (int i : {-1, 1}) {
1697  for (const Vector_t & vec: outline) {
1698  Vector_t vecInLabCoords = csTrafoGlobal2Local_m.transformFrom(vec + i * dY);
1699  for (unsigned int d = 0; d < 3; ++ d) {
1700  if (vecInLabCoords(d) < bb.lowerLeftCorner(d)) {
1701  bb.lowerLeftCorner(d) = vecInLabCoords(d);
1702  }
1703  if (vecInLabCoords(d) > bb.upperRightCorner(d)) {
1704  bb.upperRightCorner(d) = vecInLabCoords(d);
1705  }
1706  }
1707  }
1708  }
1709 
1710  return bb;
1711 }
Inform * gmsg
Definition: Main.cpp:62
Tps< T > cos(const Tps< T > &x)
Cosine.
Definition: TpsMath.h:129
Tps< T > pow(const Tps< T > &x, int y)
Integer power.
Definition: TpsMath.h:76
Tps< T > exp(const Tps< T > &x)
Exponential.
Definition: TpsMath.h:165
Tps< T > sin(const Tps< T > &x)
Sine.
Definition: TpsMath.h:111
Tps< T > sqrt(const Tps< T > &x)
Square root.
Definition: TpsMath.h:91
T euclidean_norm(const Vector< T > &)
Euclidean norm.
Definition: Vector.h:243
Quaternion getQuaternion(Vector_t u, Vector_t ref)
Definition: Quaternion.cpp:34
Vector3D cross(const Vector3D &lhs, const Vector3D &rhs)
Vector cross product.
Definition: Vector3D.cpp:111
double dot(const Vector3D &lhs, const Vector3D &rhs)
Vector dot product.
Definition: Vector3D.cpp:118
T det(const AntiSymTenzor< T, 3 > &)
#define X(arg)
Definition: fftpack.cpp:112
T::PETE_Expr_t::PETE_Return_t max(const PETE_Expr< T > &expr, NDIndex< D > &loc)
Definition: ReductionLoc.h:84
T::PETE_Expr_t::PETE_Return_t min(const PETE_Expr< T > &expr, NDIndex< D > &loc)
Definition: ReductionLoc.h:76
std::complex< double > a
PETE_TUTree< FnAbs, typename T::PETE_Expr_t > abs(const PETE_Expr< T > &l)
PETE_TBTree< FnFmod, PETE_Scalar< Vektor< T1, Dim > >, typename T2::PETE_Expr_t > fmod(const Vektor< T1, Dim > &l, const PETE_Expr< T2 > &r)
PETE_TBTree< FnArcTan2, PETE_Scalar< Vektor< T1, Dim > >, typename T2::PETE_Expr_t > atan2(const Vektor< T1, Dim > &l, const PETE_Expr< T2 > &r)
PETE_TBTree< FnCopysign, PETE_Scalar< Vektor< T1, Dim > >, typename T2::PETE_Expr_t > copysign(const Vektor< T1, Dim > &l, const PETE_Expr< T2 > &r)
PETE_TUTree< FnCeil, typename T::PETE_Expr_t > ceil(const PETE_Expr< T > &l)
Definition: PETE.h:728
PETE_TUTree< FnArcSin, typename T::PETE_Expr_t > asin(const PETE_Expr< T > &l)
Definition: PETE.h:726
Inform & endl(Inform &inf)
Definition: Inform.cpp:42
Inform & level2(Inform &inf)
Definition: Inform.cpp:46
#define ERRORMSG(msg)
Definition: IpplInfo.h:350
#define WARNMSG(msg)
Definition: IpplInfo.h:349
const std::string name
constexpr double two_pi
The value of.
Definition: Physics.h:33
constexpr double deg2rad
The conversion factor from degrees to radians.
Definition: Physics.h:48
constexpr double e
The value of.
Definition: Physics.h:39
constexpr double c
The velocity of light in m/s.
Definition: Physics.h:51
constexpr double pi
The value of.
Definition: Physics.h:30
constexpr double rad2deg
The conversion factor from radians to degrees.
Definition: Physics.h:45
T isinf(T x)
isinf function with adjusted return type
Definition: matheval.hpp:70
T isnan(T x)
isnan function with adjusted return type
Definition: matheval.hpp:66
bool writeBendTrajectories
Definition: Options.cpp:95
std::string combineFilePath(std::initializer_list< std::string > ilist)
Definition: Util.cpp:139
std::string toUpper(const std::string &str)
Definition: Util.cpp:132
ParticlePos_t & R
const PartData * getReference() const
double getQ() const
Access to reference data.
double getdT() const
std::string getInputBasename()
get input file name without extension
Definition: OpalData.cpp:658
static OpalData * getInstance()
Definition: OpalData.cpp:195
std::string getAuxiliaryOutputDirectory() const
get the name of the the additional data directory
Definition: OpalData.cpp:650
Definition: Bend2D.h:51
void calculateRefTrajectory(double &angleX, double &angleY)
Definition: Bend2D.cpp:509
double exitParameter3_m
Definition: Bend2D.h:234
double widthExitFringe_m
Definition: Bend2D.h:208
virtual bool apply(const size_t &i, const double &t, Vector_t &E, Vector_t &B) override
Apply field to particles with coordinates in magnet frame.
Definition: Bend2D.cpp:147
void setupFringeWidths()
Definition: Bend2D.cpp:1629
Vector_t transformToExitRegion(const Vector_t &R) const
Definition: Bend2D.h:375
void readFieldMap(Inform &msg)
Definition: Bend2D.cpp:946
bool reinitialize_m
Definition: Bend2D.h:216
int polyOrderExit_m
function origin that Enge function ends.
Definition: Bend2D.h:259
double exitParameter1_m
Definition: Bend2D.h:232
Vector_t calcExitFringeField(const Vector_t &R, double deltaX)
Definition: Bend2D.cpp:410
virtual bool findChordLength(double &chordLength)=0
gsl_interp_accel * entryFieldAccel_m
Definition: Bend2D.h:242
bool inMagnetCentralRegion(const Vector_t &R) const
Definition: Bend2D.cpp:797
std::vector< double > engeCoeffsEntry_m
Enge coefficients for map entry and exit regions.
Definition: Bend2D.h:237
std::array< double, 2 > getEntranceFringeFieldLength() const
Get entrance fringe field length.
Definition: Bend2D.cpp:1647
virtual void initialise(PartBunchBase< double, 3 > *bunch, double &startField, double &endField) override
Definition: Bend2D.cpp:190
double deltaEndEntry_m
function origin where Enge function starts.
Definition: Bend2D.h:251
double deltaBeginEntry_m
Definition: Bend2D.h:249
void setExitAngle(double exitAngle)
Definition: Bend2D.h:341
double maxAngle_m
Definition: Bend2D.h:272
bool inMagnetExitRegion(const Vector_t &R) const
Definition: Bend2D.cpp:822
CoordinateSystemTrafo beginToEnd_m
Definition: Bend2D.h:266
CoordinateSystemTrafo getBeginToEnd_local() const
Definition: Bend2D.h:354
double deltaBeginExit_m
Enge function order for entry region.
Definition: Bend2D.h:255
Vector_t calcCentralField(const Vector_t &R, double deltaX)
Definition: Bend2D.cpp:355
void setEngeOriginDelta(double delta)
Definition: Bend2D.cpp:1047
double entranceParameter3_m
Definition: Bend2D.h:231
bool calculateMapField(const Vector_t &R, Vector_t &B)
Definition: Bend2D.cpp:446
void setBendStrength()
Definition: Bend2D.cpp:1038
void findBendStrength()
Definition: Bend2D.cpp:682
double widthEntranceFringe_m
End of magnet field map in s coordinates (m).
Definition: Bend2D.h:207
gsl_spline ** exitFieldValues_m
Definition: Bend2D.h:241
std::vector< double > engeCoeffsExit_m
Definition: Bend2D.h:238
double entranceParameter1_m
Definition: Bend2D.h:229
double designRadius_m
through the bend.
Definition: Bend2D.h:199
void setFieldCalcParam()
Definition: Bend2D.cpp:1069
double sinEntranceAngle_m
Definition: Bend2D.h:262
BorisPusher pusher_m
Definition: Bend2D.h:196
void setBendEffectiveLength(double startField, double endField)
Definition: Bend2D.cpp:1026
void setNSlices(const std::size_t &nSlices)
Definition: Bend2D.cpp:1636
Vector_t calcEntranceFringeField(const Vector_t &R, double deltaX)
Definition: Bend2D.cpp:373
bool setupBendGeometry(double &startField, double &endField)
Definition: Bend2D.cpp:1166
bool findIdealBendParameters(double chordLength)
Definition: Bend2D.cpp:750
BoundingBox getBoundingBoxInLabCoords() const override
Definition: Bend2D.cpp:1685
Bend2D()
Definition: Bend2D.cpp:40
std::string messageHeader_m
Definition: Bend2D.h:194
double tanEntranceAngle_m
Definition: Bend2D.h:263
virtual void setEntranceAngle(double entranceAngle) override
Definition: Bend2D.h:332
double startField_m
Dipole field index.
Definition: Bend2D.h:204
virtual void goOnline(const double &kineticEnergy) override
Definition: Bend2D.cpp:186
std::size_t nSlices_m
Definition: Bend2D.h:274
double exitParameter2_m
Definition: Bend2D.h:233
double endField_m
Start of magnet field map in s coordinates (m).
Definition: Bend2D.h:205
CoordinateSystemTrafo toEntranceRegion_m
Definition: Bend2D.h:268
bool isPositionInExitField(const Vector_t &R) const
Definition: Bend2D.cpp:838
void setFieldBoundaries(double startField, double endField)
Definition: Bend2D.cpp:1234
gsl_interp_accel * exitFieldAccel_m
Definition: Bend2D.h:243
double deltaEndExit_m
function origin that Enge function starts.
Definition: Bend2D.h:257
Vector_t transformToEntranceRegion(const Vector_t &R) const
Definition: Bend2D.h:370
bool inMagnetEntranceRegion(const Vector_t &R) const
Definition: Bend2D.cpp:815
void adjustFringeFields(double ratio)
Definition: Bend2D.cpp:233
double exitAngle_m
Bend design radius (m).
Definition: Bend2D.h:201
double fieldIndex_m
and the exit face of the magnet (radians).
Definition: Bend2D.h:203
virtual ~Bend2D()
Definition: Bend2D.cpp:123
MeshData getSurfaceMesh() const
Definition: Bend2D.cpp:1387
bool isFieldZero()
Definition: Bend2D.cpp:1248
std::array< double, 2 > getExitFringeFieldLength() const
Get exit fringe field length.
Definition: Bend2D.cpp:1667
CoordinateSystemTrafo beginToEnd_lcs_m
Definition: Bend2D.h:267
double entranceParameter2_m
Definition: Bend2D.h:230
int polyOrderEntry_m
function origin that Enge function ends.
Definition: Bend2D.h:253
double estimateFieldAdjustmentStep(double actualBendAngle)
Definition: Bend2D.cpp:576
virtual bool isInside(const Vector_t &r) const override
Definition: Bend2D.cpp:1615
virtual CoordinateSystemTrafo getEdgeToEnd() const override
Definition: Bend2D.h:348
bool setupDefaultFieldMap()
Definition: Bend2D.cpp:1218
std::vector< Vector_t > getOutline() const
Definition: Bend2D.cpp:1280
CoordinateSystemTrafo toExitRegion_m
Definition: Bend2D.h:269
gsl_spline ** entryFieldValues_m
Definition: Bend2D.h:240
bool initializeFieldMap()
Definition: Bend2D.cpp:782
CoordinateSystemTrafo computeAngleTrafo_m
Definition: Bend2D.h:271
void calcEngeFunction(double zNormalized, const std::vector< double > &engeCoeff, int polyOrder, double &engeFunc, double &engeFuncDeriv, double &engeFuncSecDerivNorm)
Definition: Bend2D.cpp:290
std::size_t getNSlices() const
Definition: Bend2D.cpp:1641
virtual ElementBase::ElementType getType() const override=0
Get element type std::string.
void setupPusher(PartBunchBase< double, 3 > *bunch)
Definition: Bend2D.cpp:1241
double tanExitAngle_m
Definition: Bend2D.h:264
void print(Inform &msg, double bendAngleX, double bendAngle)
Definition: Bend2D.cpp:846
virtual bool applyToReferenceParticle(const Vector_t &R, const Vector_t &P, const double &t, Vector_t &E, Vector_t &B) override
Definition: Bend2D.cpp:178
void findBendEffectiveLength(double startField, double endField)
Definition: Bend2D.cpp:592
void setGapFromFieldMap()
Definition: Bend2D.cpp:1157
bool isPositionInEntranceField(const Vector_t &R) const
Definition: Bend2D.cpp:831
double cosEntranceAngle_m
Enge function order for entry region.
Definition: Bend2D.h:261
double calculateBendAngle()
Definition: Bend2D.cpp:251
double calcGamma() const
Calculate gamma from design energy.
Definition: BendBase.cpp:89
std::vector< Vector_t > refTrajMap_m
Map of reference particle trajectory.
Definition: BendBase.h:64
const bool fast_m
Flag to turn on fast field calculation.
Definition: BendBase.h:57
double calcBendAngle(double chordLength, double radius) const
Calculate bend angle from chord length and design radius.
Definition: BendBase.cpp:79
double fieldAmplitude_m
Field amplitude.
Definition: BendBase.h:71
double calcBetaGamma() const
Calculate beta*gamma from design energy.
Definition: BendBase.cpp:95
double getChordLength() const
Definition: BendBase.h:82
double designEnergy_m
Bend design energy (eV).
Definition: BendBase.h:61
double entranceAngle_m
Definition: BendBase.h:54
double getBendAngle() const
Definition: BendBase.h:92
double gap_m
Full vertical gap of the magnets.
Definition: BendBase.h:59
double fieldAmplitudeY_m
Definition: BendBase.h:68
Fieldmap * fieldmap_m
Magnet field map.
Definition: BendBase.h:56
std::string fileName_m
Definition: BendBase.h:73
double chordLength_m
Definition: BendBase.h:52
double getEntranceAngle() const
Definition: BendBase.h:103
double calcFieldAmplitude(double radius) const
Calculate field amplitude from design energy and radius.
Definition: BendBase.cpp:70
double calcDesignRadius(double fieldAmplitude) const
Calculate design radius from design energy and field amplitude.
Definition: BendBase.cpp:61
double fieldAmplitudeX_m
Definition: BendBase.h:66
double getFullGap() const
Definition: BendBase.h:113
double angle_m
Bend angle.
Definition: BendBase.h:53
bool online_m
Definition: Component.h:195
PartBunchBase< double, 3 > * RefPartBunch_m
Definition: Component.h:194
virtual const std::string & getName() const
Get element name.
virtual double getElementLength() const
Get design length.
Definition: ElementBase.h:432
virtual void setElementLength(double length)
Set design length.
Definition: ElementBase.h:436
std::pair< ApertureType, std::vector< double > > aperture_m
Definition: ElementBase.h:390
virtual CoordinateSystemTrafo getEdgeToBegin() const
Definition: ElementBase.h:520
double getRotationAboutZ() const
Definition: ElementBase.h:591
double rotationZAxis_m
Definition: ElementBase.h:394
double elementEdge_m
Definition: ElementBase.h:392
CoordinateSystemTrafo csTrafoGlobal2Local_m
Definition: ElementBase.h:387
Vector_t rotateFrom(const Vector_t &r) const
Vector_t transformFrom(const Vector_t &r) const
Vector_t transformTo(const Vector_t &r) const
Vector_t rotate(const Vector_t &) const
Definition: Quaternion.cpp:122
Quaternion conjugate() const
Definition: Quaternion.h:105
virtual void getInfo(Inform *msg)=0
static Fieldmap * getFieldmap(std::string Filename, bool fast=false)
Definition: Fieldmap.cpp:50
virtual void get1DProfile1EntranceParam(double &entranceParameter1, double &entranceParameter2, double &entranceParameter3)
Definition: Fieldmap.cpp:717
virtual double getFieldGap()
Definition: Fieldmap.cpp:729
virtual void get1DProfile1EngeCoeffs(std::vector< double > &engeCoeffsEntry, std::vector< double > &engeCoeffsExit)
Definition: Fieldmap.cpp:712
virtual void get1DProfile1ExitParam(double &exitParameter1, double &exitParameter2, double &exitParameter3)
Definition: Fieldmap.cpp:723
virtual void readMap()=0
std::vector< std::pair< Vector_t, Vector_t > > decorations_m
Definition: MeshGenerator.h:26
std::vector< Vektor< unsigned int, 3 > > triangles_m
Definition: MeshGenerator.h:25
std::vector< Vector_t > vertices_m
Definition: MeshGenerator.h:24
void kick(const Vector_t &R, Vector_t &P, const Vector_t &Ef, const Vector_t &Bf, const double &dt) const
Definition: BorisPusher.h:65
void initialise(const PartData *ref)
Definition: BorisPusher.h:60
void push(Vector_t &R, const Vector_t &P, const double &dt) const
Definition: BorisPusher.h:131
Definition: Inform.h:42
static int myNode()
Definition: IpplInfo.cpp:691
Definition: Vec.h:22
Vektor< double, 3 > Vector_t
Definition: Vektor.h:6