OPAL (Object Oriented Parallel Accelerator Library)  2024.1
OPAL
OpalWake.cpp
Go to the documentation of this file.
1 //
2 // Class OpalWake
3 // The class for the OPAL WAKE command.
4 //
5 // Copyright (c) 2008 - 2022, 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 //
18 #include "Structure/OpalWake.h"
19 
22 #include "Attributes/Attributes.h"
27 #include "Utilities/OpalFilter.h"
28 
29 #include <unordered_map>
30 
31 extern Inform *gmsg;
32 
33 // The attributes of class OpalWake.
34 namespace {
35  enum {
36  TYPE, // The type of the wake
37  NBIN, // Number of bins for the line density
38  CONST_LENGTH, // True if the length of the Bunch is considered as constant
39  CONDUCT, // Conductivity, either AC or DC
40  Z0,
41  RADIUS, // Radius of the tube
42  SIGMA,
43  TAU,
44  FILTERS, // List of filters to apply on line density
45  FNAME,
46  SIZE
47  };
48 }
49 
51  Definition(SIZE, "WAKE",
52  "The \"WAKE\" statement defines data for the wakefuction "
53  "on an element."),
54  wf_m(0) {
56  ("TYPE", "Specifies the wake function.",
57  {"1D-CSR", "1D-CSR-IGF", "LONG-SHORT-RANGE", "TRANSV-SHORT-RANGE"});
58 
60  ("NBIN", "Number of bins for the line density calculation");
61 
62  itsAttr[CONST_LENGTH] = Attributes::makeBool
63  ("CONST_LENGTH", "True if the length of the Bunch is considered as constant");
64 
65  itsAttr[CONDUCT] = Attributes::makePredefinedString
66  ("CONDUCT", "Conductivity.", {"DC", "AC"});
67 
68  itsAttr[Z0] = Attributes::makeReal
69  ("Z0", "Impedance of the beam pipe ");
70 
71  itsAttr[RADIUS] = Attributes::makeReal
72  ("RADIUS", "The radius of the beam pipe [m]");
73 
74  itsAttr[SIGMA] = Attributes::makeReal
75  ("SIGMA", "Material constant dependant on the beam pipe material");
76 
77  itsAttr[TAU] = Attributes::makeReal
78  ("TAU", "Material constant dependant on the beam pipe material");
79 
80  itsAttr[FILTERS] = Attributes::makeStringArray
81  ("FILTERS", "List of filters to apply on line density");
82 
83  itsAttr[FNAME] = Attributes::makeString
84  ("FNAME", "Filename of the wakefield file");
85 
86  OpalWake* defWake = clone("UNNAMED_WAKE");
87  defWake->builtin = true;
88 
89  try {
90  defWake->update();
91  OpalData::getInstance()->define(defWake);
92  } catch(...) {
93  delete defWake;
94  }
95 
97 }
98 
99 
100 OpalWake::OpalWake(const std::string& name, OpalWake* parent):
101  Definition(name, parent),
102  wf_m(parent->wf_m)
103 {}
104 
105 
107  delete wf_m;
108 }
109 
110 
112  // Can replace only by another WAKE.
113  return dynamic_cast<OpalWake*>(object) != 0;
114 }
115 
116 
117 OpalWake* OpalWake::clone(const std::string& name) {
118  return new OpalWake(name, this);
119 }
120 
121 
123  update();
124 }
125 
126 
127 OpalWake* OpalWake::find(const std::string& name) {
128  OpalWake* wake = dynamic_cast<OpalWake*>(OpalData::getInstance()->find(name));
129  if (wake == 0) {
130  throw OpalException("OpalWake::find()", "Wake \"" + name + "\" not found.");
131  }
132  return wake;
133 }
134 
135 
137  return (int)Attributes::getReal(itsAttr[NBIN]);
138 }
139 
140 
142  // Set default name.
143  if (getOpalName().empty()) setOpalName("UNNAMED_WAKE");
144 }
145 
146 
148  *gmsg << "* ************* W A K E ************************************************************\n";
149  *gmsg << "OpalWake::initWakefunction ";
150  *gmsg << "for element " << element.getName() << "\n";
151  *gmsg << "* **********************************************************************************" << endl;
152 
153  std::vector<std::string> filters_str = Attributes::getStringArray(itsAttr[FILTERS]);
154  std::vector<Filter*> filters;
155 
156  for (std::vector<std::string>::const_iterator fit = filters_str.begin(); fit != filters_str.end(); ++fit) {
157  OpalFilter *f = OpalFilter::find(*fit);
158  if (f) {
159  f->initOpalFilter();
160  filters.push_back(f->filter_m);
161  }
162  }
163 
164  static const std::unordered_map<std::string, OpalWakeType> stringOpalWakeType_s = {
165  {"1D-CSR", OpalWakeType::CSR},
166  {"1D-CSR-IGF", OpalWakeType::CSRIGF},
167  {"LONG-SHORT-RANGE", OpalWakeType::LONGSHORTRANGE},
168  {"TRANSV-SHORT-RANGE", OpalWakeType::TRANSVSHORTRANGE}
169  };
170 
171  if (!itsAttr[TYPE]) {
172  throw OpalException("TrackRun::execute",
173  "The attribute \"TYPE\" isn't set for the \"WAKE\" statement");
174  }
175  OpalWakeType type = stringOpalWakeType_s.at(Attributes::getString(itsAttr[TYPE]));
176  switch (type) {
177  case OpalWakeType::CSR: {
178  if (filters.size() == 0 && Attributes::getReal(itsAttr[NBIN]) <= 7) {
179  throw OpalException("OpalWake::initWakeFunction",
180  "At least 8 bins have to be used, ideally far more");
181  }
182 
183  wf_m = new CSRWakeFunction(getOpalName(), filters,
184  (int)(Attributes::getReal(itsAttr[NBIN])));
185  break;
186  }
187  case OpalWakeType::CSRIGF: {
188  if (filters.size() == 0 && Attributes::getReal(itsAttr[NBIN]) <= 7) {
189  throw OpalException("OpalWake::initWakeFunction",
190  "At least 8 bins have to be used, ideally far more");
191  }
192 
193  wf_m = new CSRIGFWakeFunction(getOpalName(), filters,
194  (int)(Attributes::getReal(itsAttr[NBIN])));
195  break;
196  }
198  int acMode = Attributes::getString(itsAttr[CONDUCT]) == "DC"? 2: 1;
199 
200  wf_m = new GreenWakeFunction(getOpalName(), filters,
203  Attributes::getReal(itsAttr[RADIUS]),
205  acMode,
208  Attributes::getBool(itsAttr[CONST_LENGTH]),
210  break;
211  }
213  int acMode = Attributes::getString(itsAttr[CONDUCT]) == "DC" ? 2: 1;
214 
215  wf_m = new GreenWakeFunction(getOpalName(), filters,
218  Attributes::getReal(itsAttr[RADIUS]),
220  acMode,
223  Attributes::getBool(itsAttr[CONST_LENGTH]),
225  break;
226  }
227  default: {
228  throw OpalException("OpalWake::initWakefunction",
229  "Invalid \"TYPE\" of \"WAKE\" statement");
230  }
231  }
232 }
233 
234 void OpalWake::print(std::ostream& os) const {
235  os << "* ************* W A K E ************************************************************ " << std::endl;
236  os << "* WAKE " << getOpalName() << '\n'
237  << "* BINS " << Attributes::getReal(itsAttr[NBIN]) << '\n'
238  << "* CONST_LENGTH " << Attributes::getReal(itsAttr[CONST_LENGTH]) << '\n'
239  << "* CONDUCT " << Attributes::getReal(itsAttr[CONDUCT]) << '\n'
240  << "* Z0 " << Attributes::getReal(itsAttr[Z0]) << '\n'
241  << "* RADIUS " << Attributes::getReal(itsAttr[RADIUS]) << '\n'
242  << "* SIGMA " << Attributes::getReal(itsAttr[SIGMA]) << '\n'
243  << "* TAU " << Attributes::getReal(itsAttr[TAU]) << '\n';
244  os << "* ********************************************************************************** " << std::endl;
245 }
Attribute makeReal(const std::string &name, const std::string &help)
Make real attribute.
Definition: Attributes.cpp:240
static OpalData * getInstance()
Definition: OpalData.cpp:196
bool builtin
Built-in flag.
Definition: Object.h:233
WakeFunction * wf_m
Definition: OpalWake.h:56
void setOpalName(const std::string &name)
Set object name.
Definition: Object.cpp:331
The base class for all OPAL objects.
Definition: Object.h:48
std::string getString(const Attribute &attr)
Get string value.
Definition: Attributes.cpp:343
OpalWakeType
Definition: OpalWake.h:59
static OpalFilter * find(const std::string &name)
Find named FILTER.
Definition: OpalFilter.cpp:124
void define(Object *newObject)
Define a new object.
Definition: OpalData.cpp:489
virtual ~OpalWake()
Definition: OpalWake.cpp:106
virtual const std::string & getName() const
Get element name.
Inform & endl(Inform &inf)
Definition: Inform.cpp:42
Filter * filter_m
Definition: OpalFilter.h:61
bool getBool(const Attribute &attr)
Return logical value.
Definition: Attributes.cpp:100
Attribute makeBool(const std::string &name, const std::string &help)
Make logical attribute.
Definition: Attributes.cpp:90
void registerOwnership(const AttributeHandler::OwnerType &itsClass) const
Definition: Object.cpp:191
The base class for all OPAL exceptions.
Definition: OpalException.h:28
Attribute makeStringArray(const std::string &name, const std::string &help)
Create a string array attribute.
Definition: Attributes.cpp:473
Attribute makePredefinedString(const std::string &name, const std::string &help, const std::initializer_list< std::string > &predefinedStrings)
Make predefined string attribute.
Definition: Attributes.cpp:409
Definition: Inform.h:42
virtual bool canReplaceBy(Object *object)
Test if replacement is allowed.
Definition: OpalWake.cpp:111
OpalWake()
Exemplar constructor.
Definition: OpalWake.cpp:50
Attribute makeString(const std::string &name, const std::string &help)
Make string attribute.
Definition: Attributes.cpp:332
void initWakefunction(const ElementBase &element)
Definition: OpalWake.cpp:147
const std::string & getOpalName() const
Return object name.
Definition: Object.cpp:310
std::vector< Attribute > itsAttr
The object attributes.
Definition: Object.h:216
virtual void update()
Update the WAKE data.
Definition: OpalWake.cpp:141
void initOpalFilter()
Definition: OpalFilter.cpp:141
const std::string name
void print(std::ostream &os) const
Print the object.
Definition: OpalWake.cpp:234
double getReal(const Attribute &attr)
Return real value.
Definition: Attributes.cpp:252
The base class for all OPAL definitions.
Definition: Definition.h:30
std::vector< std::string > getStringArray(const Attribute &attr)
Get string array value.
Definition: Attributes.cpp:478
int getNumberOfBins()
Definition: OpalWake.cpp:136
Object * find(const std::string &name)
Find entry.
Definition: OpalData.cpp:571
virtual OpalWake * clone(const std::string &name)
Make clone.
Definition: OpalWake.cpp:117
virtual void execute()
Check the WAKE data.
Definition: OpalWake.cpp:122
static OpalWake * find(const std::string &name)
Find named WAKE.
Definition: OpalWake.cpp:127
SDDS1 &description type
Definition: test.stat:4
Inform * gmsg
Definition: Main.cpp:70