OPAL (Object Oriented Parallel Accelerator Library) 2024.2
OPAL
OpalElement.cpp
Go to the documentation of this file.
1//
2// Class OpalElement
3// Base class for all beam line elements.
4//
5// Copyright (c) 200x - 2020, Paul Scherrer Institut, Villigen PSI, Switzerland
6// All rights reserved
7//
8// This file is part of OPAL.
9//
10// OPAL is free software: you can redistribute it and/or modify
11// it under the terms of the GNU General Public License as published by
12// the Free Software Foundation, either version 3 of the License, or
13// (at your option) any later version.
14//
15// You should have received a copy of the GNU General Public License
16// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
17//
19
20#include "AbsBeamline/Bend2D.h"
25#include "Parser/Statement.h"
26#include "Physics/Physics.h"
28#include "Utilities/Options.h"
30#include "Utilities/Util.h"
31
32#include <cmath>
33#include <exception>
34#include <ostream>
35#include <regex>
36#include <sstream>
37#include <string>
38#include <utility>
39#include <vector>
40
41OpalElement::OpalElement(int size, const char* name, const char* help):
42 Element(size, name, help), itsSize(size) {
44 ("TYPE", "The element design type.",
45 {"RING",
46 "CARBONCYCL",
47 "CYCIAE",
48 "AVFEQ",
49 "FFA",
50 "BANDRF",
51 "SYNCHROCYCLOTRON",
52 "SINGLEGAP",
53 "STANDING",
54 "TEMPORAL",
55 "SPATIAL"});
56
58 ("L", "The element length [m]");
59
61 ("ELEMEDGE", "The position of the element in path length [m]");
62
64 ("APERTURE", "The element aperture");
65
67 ("WAKEF", "Defines the wake function");
68
70 ("PARTICLEMATTERINTERACTION", "Defines the particle mater interaction handler");
71
73 ("ORIGIN", "The location of the element");
74
76 ("ORIENTATION", "The Tait-Bryan angles for the orientation of the element");
77
79 ("X", "The x-coordinate of the location of the element", 0);
80
82 ("Y", "The y-coordinate of the location of the element", 0);
83
85 ("Z", "The z-coordinate of the location of the element", 0);
86
88 ("THETA", "The rotation about the y-axis of the element", 0);
89
91 ("PHI", "The rotation about the x-axis of the element", 0);
92
94 ("PSI", "The rotation about the z-axis of the element", 0);
95
97 ("DX", "Misalignment in x direction", 0.0);
98
100 ("DY", "Misalignment in y direction", 0.0);
101
103 ("DZ", "Misalignment in z direction", 0.0);
104
106 ("DTHETA", "Misalignment in theta (Tait-Bryan angles)", 0.0);
107
109 ("DPHI", "Misalignment in theta (Tait-Bryan angles)", 0.0);
110
112 ("DPSI", "Misalignment in theta (Tait-Bryan angles)", 0.0);
113
115 ("OUTFN", "Output filename");
116
118 ("DELETEONTRANSVERSEEXIT", "Flag controlling if particles should be deleted if they exit "
119 "the element transversally. Default=TRUE", true);
120
121 const unsigned int end = COMMON;
122 for (unsigned int i = 0; i < end; ++ i) {
124 }
125}
126
127
128OpalElement::OpalElement(const std::string& name, OpalElement* parent):
129 Element(name, parent), itsSize(parent->itsSize)
130{}
131
132
135
136
137std::pair<ApertureType, std::vector<double> > OpalElement::getApert() const {
138
139 std::pair<ApertureType, std::vector<double> > retvalue(ApertureType::ELLIPTICAL,
140 std::vector<double>({0.5, 0.5, 1.0}));
141 if (!itsAttr[APERT]) return retvalue;
142
143 std::string aperture = Attributes::getString(itsAttr[APERT]);
144
145 std::regex square("square *\\((.*)\\)", std::regex::icase);
146 std::regex rectangle("rectangle *\\((.*)\\)", std::regex::icase);
147 std::regex circle("circle *\\((.*)\\)", std::regex::icase);
148 std::regex ellipse("ellipse *\\((.*)\\)", std::regex::icase);
149
150 std::regex twoArguments("([^,]*),([^,]*)");
151 std::regex threeArguments("([^,]*),([^,]*),([^,]*)");
152
153 std::smatch match;
154
155 const double width2HalfWidth = 0.5;
156
157 if (std::regex_search(aperture, match, square)) {
158 std::string arguments = match[1];
159 if (!std::regex_search(arguments, match, twoArguments)) {
160 retvalue.first = ApertureType::RECTANGULAR;
161
162 try {
163 retvalue.second[0] = width2HalfWidth * std::stod(arguments);
164 retvalue.second[1] = retvalue.second[0];
165 } catch (const std::exception &ex) {
166 throw OpalException("OpalElement::getApert()",
167 "could not convert '" + arguments + "' to double");
168 }
169
170 } else {
171 retvalue.first = ApertureType::CONIC_RECTANGULAR;
172
173 try {
174 retvalue.second[0] = width2HalfWidth * std::stod(match[1]);
175 retvalue.second[1] = retvalue.second[0];
176 retvalue.second[2] = std::stod(match[2]);
177 } catch (const std::exception &ex) {
178 throw OpalException("OpalElement::getApert()",
179 "could not convert '" + arguments + "' to doubles");
180 }
181 }
182
183 return retvalue;
184 }
185
186 if (std::regex_search(aperture, match, rectangle)) {
187 std::string arguments = match[1];
188
189 if (!std::regex_search(arguments, match, threeArguments)) {
190 retvalue.first = ApertureType::RECTANGULAR;
191
192 try {
193 size_t sz = 0;
194
195 retvalue.second[0] = width2HalfWidth * std::stod(arguments, &sz);
196 sz = arguments.find_first_of(",", sz) + 1;
197 retvalue.second[1] = width2HalfWidth * std::stod(arguments.substr(sz));
198
199 } catch (const std::exception &ex) {
200 throw OpalException("OpalElement::getApert()",
201 "could not convert '" + arguments + "' to doubles");
202 }
203
204 } else {
205 retvalue.first = ApertureType::CONIC_RECTANGULAR;
206
207 try {
208 retvalue.second[0] = width2HalfWidth * std::stod(match[1]);
209 retvalue.second[1] = width2HalfWidth * std::stod(match[2]);
210 retvalue.second[2] = std::stod(match[3]);
211 } catch (const std::exception &ex) {
212 throw OpalException("OpalElement::getApert()",
213 "could not convert '" + arguments + "' to doubles");
214 }
215 }
216
217 return retvalue;
218 }
219
220 if (std::regex_search(aperture, match, circle)) {
221 std::string arguments = match[1];
222 if (!std::regex_search(arguments, match, twoArguments)) {
223 retvalue.first = ApertureType::ELLIPTICAL;
224
225 try {
226 retvalue.second[0] = width2HalfWidth * std::stod(arguments);
227 retvalue.second[1] = retvalue.second[0];
228 } catch (const std::exception &ex) {
229 throw OpalException("OpalElement::getApert()",
230 "could not convert '" + arguments + "' to double");
231 }
232
233 } else {
234 retvalue.first = ApertureType::CONIC_ELLIPTICAL;
235
236 try {
237 retvalue.second[0] = width2HalfWidth * std::stod(match[1]);
238 retvalue.second[1] = retvalue.second[0];
239 retvalue.second[2] = std::stod(match[2]);
240 } catch (const std::exception &ex) {
241 throw OpalException("OpalElement::getApert()",
242 "could not convert '" + arguments + "' to doubles");
243 }
244 }
245
246 return retvalue;
247 }
248
249 if (std::regex_search(aperture, match, ellipse)) {
250 std::string arguments = match[1];
251
252 if (!std::regex_search(arguments, match, threeArguments)) {
253 retvalue.first = ApertureType::ELLIPTICAL;
254
255 try {
256 size_t sz = 0;
257
258 retvalue.second[0] = width2HalfWidth * std::stod(arguments, &sz);
259 sz = arguments.find_first_of(",", sz) + 1;
260 retvalue.second[1] = width2HalfWidth * std::stod(arguments.substr(sz));
261
262 } catch (const std::exception &ex) {
263 throw OpalException("OpalElement::getApert()",
264 "could not convert '" + arguments + "' to doubles");
265 }
266
267 } else {
268 retvalue.first = ApertureType::CONIC_ELLIPTICAL;
269
270 try {
271 retvalue.second[0] = width2HalfWidth * std::stod(match[1]);
272 retvalue.second[1] = width2HalfWidth * std::stod(match[2]);
273 retvalue.second[2] = std::stod(match[3]);
274 } catch (const std::exception &ex) {
275 throw OpalException("OpalElement::getApert()",
276 "could not convert '" + arguments + "' to doubles");
277 }
278 }
279
280 return retvalue;
281 }
282
283 if (!aperture.empty()) {
284 throw OpalException("OpalElement::getApert()",
285 "Unknown aperture type '" + aperture + "'.");
286 }
287
288 return retvalue;
289}
290
293}
294
295
296const std::string OpalElement::getTypeName() const {
297 const Attribute* attr = findAttribute("TYPE");
298 return attr ? Attributes::getString(*attr) : std::string();
299}
300
304const std::string OpalElement::getWakeF() const {
305 const Attribute* attr = findAttribute("WAKEF");
306 return attr ? Attributes::getString(*attr) : std::string();
307}
308
310 const Attribute* attr = findAttribute("PARTICLEMATTERINTERACTION");
311 return attr ? Attributes::getString(*attr) : std::string();
312}
313
315 while (stat.delimiter(',')) {
316 std::string name = Expressions::parseString(stat, "Attribute name expected.");
318
319 if (attr == 0) {
320 throw OpalException("OpalElement::parse",
321 "unknown attribute \"" + name + "\"");
322 }
323
324 if (stat.delimiter('[')) {
325 int index = int(std::round(Expressions::parseRealConst(stat)));
327
328 if (stat.delimiter('=')) {
329 attr->parseComponent(stat, true, index);
330 } else if (stat.delimiter(":=")) {
331 attr->parseComponent(stat, false, index);
332 } else {
333 throw ParseError("OpalElement::parse()",
334 "Delimiter \"=\" or \":=\" expected.");
335 }
336 } else {
337 if (stat.delimiter('=')) {
338 attr->parse(stat, true);
339 } else if (stat.delimiter(":=")) {
340 attr->parse(stat, false);
341 } else {
342 attr->setDefault();
343 }
344 }
345 }
346}
347
348
349void OpalElement::print(std::ostream& os) const {
350 std::string head = getOpalName();
351
352 Object* parent = getParent();
353 if (parent != 0 && ! parent->getOpalName().empty()) {
354 if (! getOpalName().empty()) head += ':';
355 head += parent->getOpalName();
356 }
357
358 os << head;
359 os << ';'; // << "JMJdebug OPALElement.cc" ;
360 os << std::endl;
361}
362
363
365 int order, int& len,
366 const std::string& sName,
367 const std::string& tName,
368 const Attribute& length,
369 const Attribute& sNorm,
370 const Attribute& sSkew) {
371 // Find out which type of output is required.
372 int flag = 0;
373 if (sNorm) {
374 if (sNorm.getBase().isExpression()) {
375 flag += 2;
376 } else if (Attributes::getReal(sNorm) != 0.0) {
377 flag += 1;
378 }
379 }
380
381 if (sSkew) {
382 if (sSkew.getBase().isExpression()) {
383 flag += 6;
384 } else if (Attributes::getReal(sSkew) != 0.0) {
385 flag += 3;
386 }
387 }
388 // cout << "JMJdebug, OpalElement.cc: flag=" << flag << endl ;
389 // Now do the output.
390 int div = 2 * (order + 1);
391
392 switch (flag) {
393
394 case 0:
395 // No component at all.
396 break;
397
398 case 1:
399 case 2:
400 // Pure normal component.
401 {
402 std::string normImage = sNorm.getImage();
403 if (length) {
404 normImage = "(" + normImage + ")*(" + length.getImage() + ")";
405 }
406 printAttribute(os, sName, normImage, len);
407 }
408 break;
409
410 case 3:
411 case 6:
412 // Pure skew component.
413 {
414 std::string skewImage = sSkew.getImage();
415 if (length) {
416 skewImage = "(" + skewImage + ")*(" + length.getImage() + ")";
417 }
418 printAttribute(os, sName, skewImage, len);
419 double tilt = Physics::pi / double(div);
420 printAttribute(os, tName, tilt, len);
421 }
422 break;
423
424 case 4:
425 // Both components are non-zero constants.
426 {
427 double sn = Attributes::getReal(sNorm);
428 double ss = Attributes::getReal(sSkew);
429 double strength = std::sqrt(sn * sn + ss * ss);
430 if (strength) {
431 std::ostringstream ts;
432 ts << strength;
433 std::string image = ts.str();
434 if (length) {
435 image = "(" + image + ")*(" + length.getImage() + ")";
436 }
437 printAttribute(os, sName, image, len);
438 double tilt = - std::atan2(ss, sn) / double(div);
439 if (tilt) printAttribute(os, tName, tilt, len);
440 }
441 }
442 break;
443
444 case 5:
445 case 7:
446 case 8:
447 // One or both components is/are expressions.
448 {
449 std::string normImage = sNorm.getImage();
450 std::string skewImage = sSkew.getImage();
451 std::string image =
452 "SQRT((" + normImage + ")^2+(" + skewImage + ")^2)";
453 printAttribute(os, sName, image, len);
454 if (length) {
455 image = "(" + image + ")*(" + length.getImage() + ")";
456 }
457 std::string divisor;
458 if (div < 9) {
459 divisor = "0";
460 divisor[0] += div;
461 } else {
462 divisor = "00";
463 divisor[0] += div / 10;
464 divisor[1] += div % 10;
465 }
466 image = "-ATAN2(" + skewImage + ',' + normImage + ")/" + divisor;
467 printAttribute(os, tName, image, len);
468 break;
469 }
470 }
471}
472
474 ElementBase* base = getElement();
475
476 auto apert = getApert();
477 base->setAperture(apert.first, apert.second);
478
480 std::vector<double> ori = Attributes::getRealArray(itsAttr[ORIGIN]);
481 std::vector<double> dir = Attributes::getRealArray(itsAttr[ORIENTATION]);
482 Vector_t origin(0.0);
483 Quaternion rotation;
484
485 if (dir.size() == 3) {
486 Quaternion rotTheta(std::cos(0.5 * dir[0]), 0, std::sin(0.5 * dir[0]), 0);
487 Quaternion rotPhi(std::cos(0.5 * dir[1]), std::sin(0.5 * dir[1]), 0, 0);
488 Quaternion rotPsi(std::cos(0.5 * dir[2]), 0, 0, std::sin(0.5 * dir[2]));
489 rotation = rotTheta * (rotPhi * rotPsi);
490 } else {
491 if (itsAttr[ORIENTATION]) {
492 throw OpalException("OpalElement::update",
493 "Parameter orientation is array of 3 values (theta, phi, psi);\n" +
494 std::to_string(dir.size()) + " values provided");
495 }
496 }
497
498 if (ori.size() == 3) {
499 origin = Vector_t({ori[0], ori[1], ori[2]});
500 } else {
501 if (itsAttr[ORIGIN]) {
502 throw OpalException("OpalElement::update",
503 "Parameter origin is array of 3 values (x, y, z);\n" +
504 std::to_string(ori.size()) + " values provided");
505 }
506 }
507
508 CoordinateSystemTrafo global2local(origin,
509 rotation.conjugate());
510 base->setCSTrafoGlobal2Local(global2local);
511 base->fixPosition();
512
513 } else if (!itsAttr[PSI].defaultUsed() &&
514 itsAttr[X].defaultUsed() &&
515 itsAttr[Y].defaultUsed() &&
516 itsAttr[Z].defaultUsed() &&
517 itsAttr[THETA].defaultUsed() &&
518 itsAttr[PHI].defaultUsed()) {
520 } else if (!itsAttr[X].defaultUsed() ||
521 !itsAttr[Y].defaultUsed() ||
522 !itsAttr[Z].defaultUsed() ||
523 !itsAttr[THETA].defaultUsed() ||
524 !itsAttr[PHI].defaultUsed() ||
525 !itsAttr[PSI].defaultUsed()) {
526 const Vector_t origin({Attributes::getReal(itsAttr[X]),
529
530 const double theta = Attributes::getReal(itsAttr[THETA]);
531 const double phi = Attributes::getReal(itsAttr[PHI]);
532 const double psi = Attributes::getReal(itsAttr[PSI]);
533
534 Quaternion rotTheta(std::cos(0.5 * theta), 0, std::sin(0.5 * theta), 0);
535 Quaternion rotPhi(std::cos(0.5 * phi), std::sin(0.5 * phi), 0, 0);
536 Quaternion rotPsi(std::cos(0.5 * psi), 0, 0, std::sin(0.5 * psi));
537 Quaternion rotation = rotTheta * (rotPhi * rotPsi);
538
539 CoordinateSystemTrafo global2local(origin,
540 rotation.conjugate());
541 base->setCSTrafoGlobal2Local(global2local);
542 base->fixPosition();
544 }
545
546 Vector_t misalignmentShift({Attributes::getReal(itsAttr[DX]),
549 double dtheta = Attributes::getReal(itsAttr[DTHETA]);
550 double dphi = Attributes::getReal(itsAttr[DPHI]);
551 double dpsi = Attributes::getReal(itsAttr[DPSI]);
552 Quaternion rotationY(std::cos(0.5 * dtheta), 0, std::sin(0.5 * dtheta), 0);
553 Quaternion rotationX(std::cos(0.5 * dphi), std::sin(0.5 * dphi), 0, 0);
554 Quaternion rotationZ(std::cos(0.5 * dpsi), 0, 0, std::sin(0.5 * dpsi));
555 Quaternion misalignmentRotation = rotationY * rotationX * rotationZ;
556 CoordinateSystemTrafo misalignment(misalignmentShift,
557 misalignmentRotation.conjugate());
558
559 base->setMisalignment(misalignment);
560
561 if (itsAttr[ELEMEDGE])
563
565}
566
568 for (std::vector<Attribute>::size_type i = itsSize;
569 i < itsAttr.size(); ++i) {
570 Attribute &attr = itsAttr[i];
571 base->setAttribute(attr.getName(), Attributes::getReal(attr));
572
573 }
574}
575
576void OpalElement::printAttribute(std::ostream& os, const std::string& name,
577 const std::string& image, int& len) {
578 len += name.length() + image.length() + 2;
579 if (len > 74) {
580 os << ",&\n ";
581 len = name.length() + image.length() + 3;
582 } else {
583 os << ',';
584 }
585 os << name << '=' << image;
586}
587
589(std::ostream &os, const std::string &name, double value, int &len) {
590 std::ostringstream ss;
591 ss << value << std::ends;
592 printAttribute(os, name, ss.str(), len);
593}
594
595
597 if (getParent() != 0) return;
598
599 const unsigned int end = itsSize;
600 const std::string name = getOpalName();
601 for (unsigned int i = COMMON; i < end; ++ i) {
603 }
604}
PartBunchBase< T, Dim >::ConstIterator end(PartBunchBase< T, Dim > const &bunch)
const std::string name
std::string parseString(Statement &, const char msg[])
Parse string value.
void parseDelimiter(Statement &stat, char delim)
Test for one-character delimiter.
double parseRealConst(Statement &)
Parse real constant.
Attribute makeBool(const std::string &name, const std::string &help)
Make logical attribute.
double getReal(const Attribute &attr)
Return real value.
Attribute makePredefinedString(const std::string &name, const std::string &help, const std::initializer_list< std::string > &predefinedStrings)
Make predefined string attribute.
Attribute makeReal(const std::string &name, const std::string &help)
Make real attribute.
bool getBool(const Attribute &attr)
Return logical value.
Attribute makeRealArray(const std::string &name, const std::string &help)
Create real array attribute.
std::vector< double > getRealArray(const Attribute &attr)
Get array value.
std::string getString(const Attribute &attr)
Get string value.
Attribute makeString(const std::string &name, const std::string &help)
Make string attribute.
constexpr double pi
The value of.
Definition Physics.h:30
A representation of an Object attribute.
Definition Attribute.h:52
AttributeBase & getBase() const
Return reference to polymorphic value.
Definition Attribute.cpp:69
const std::string & getName() const
Return the attribute name.
Definition Attribute.cpp:92
void setDefault()
Assign default value.
void parse(Statement &stat, bool eval)
Parse attribute.
void parseComponent(Statement &stat, bool eval, int index)
Parse array component.
std::string getImage() const
Return printable representation.
Definition Attribute.cpp:87
virtual bool isExpression() const
Test for expression.
static void addAttributeOwner(const std::string &owner, const OwnerType &type, const std::string &name)
ElementBase * getElement() const
Return the embedded CLASSIC element.
Definition Element.h:120
The base class for all OPAL objects.
Definition Object.h:48
Object * getParent() const
Return parent pointer.
Definition Object.cpp:315
const std::string & getOpalName() const
Return object name.
Definition Object.cpp:310
virtual Attribute * findAttribute(const std::string &name)
Find an attribute by name.
Definition Object.cpp:64
std::vector< Attribute > itsAttr
The object attributes.
Definition Object.h:216
void setElementPosition(double elemedge)
Access to ELEMEDGE attribute.
void fixPosition()
void setAperture(const ApertureType &type, const std::vector< double > &args)
void setMisalignment(const CoordinateSystemTrafo &cst)
virtual void setAttribute(const std::string &aKey, double val)
Set value of an attribute.
void setFlagDeleteOnTransverseExit(bool=true)
void setRotationAboutZ(double rotation)
Set rotation about z axis in bend frame.
void setCSTrafoGlobal2Local(const CoordinateSystemTrafo &ori)
Quaternion conjugate() const
Definition Quaternion.h:103
Interface for statements.
Definition Statement.h:38
bool delimiter(char c)
Test for delimiter.
Parse exception.
Definition ParseError.h:32
std::pair< ApertureType, std::vector< double > > getApert() const
static void printMultipoleStrength(std::ostream &os, int order, int &len, const std::string &sName, const std::string &tName, const Attribute &length, const Attribute &vNorm, const Attribute &vSkew)
Print multipole components in OPAL-8 format.
virtual double getLength() const
Return element length.
static void printAttribute(std::ostream &os, const std::string &name, const std::string &image, int &len)
Print an attribute with a OPAL-8 name (as an expression).
@ DELETEONTRANSVERSEEXIT
Definition OpalElement.h:58
@ PARTICLEMATTERINTERACTION
Definition OpalElement.h:42
virtual void parse(Statement &)
Parse the element.
const std::string getParticleMatterInteraction() const
const std::string getWakeF() const
Return the element's type name.
virtual void updateUnknown(ElementBase *)
Transmit the `‘unknown’' (not known to OPAL) attributes to CLASSIC.
const std::string getTypeName() const
Return the element's type name.
virtual void print(std::ostream &) const
Print the object.
virtual ~OpalElement()
virtual void update()
Update the embedded CLASSIC element.
void registerOwnership() const
The base class for all OPAL exceptions.
Vektor< double, 3 > Vector_t
Definition Vektor.h:6