OPAL (Object Oriented Parallel Accelerator Library) 2024.2
OPAL
Option.cpp
Go to the documentation of this file.
1//
2// Class Option
3// The OPTION command.
4// The user interface allowing setting of OPAL options.
5// The actual option flags are contained in namespace Options.
6//
7// Copyright (c) 200x - 2020, Paul Scherrer Institut, Villigen PSI, Switzerland
8// All rights reserved
9//
10// This file is part of OPAL.
11//
12// OPAL is free software: you can redistribute it and/or modify
13// it under the terms of the GNU General Public License as published by
14// the Free Software Foundation, either version 3 of the License, or
15// (at your option) any later version.
16//
17// You should have received a copy of the GNU General Public License
18// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
19//
20#include "BasicActions/Option.h"
21
24#include "Parser/FileStream.h"
27#include "Utilities/Options.h"
28#include "Utilities/Util.h"
29
30#include "Utility/Inform.h"
31#include "Utility/IpplInfo.h"
33
34#include <array>
35#include <cstddef>
36#include <ctime>
37#include <iostream>
38#include <limits>
39#include <string_view>
40
41extern Inform* gmsg;
42extern Inform* gmsgALL;
43
44using namespace Options;
45
46namespace {
47 // The attributes of class Option.
48 enum {
49 ECHO,
50 INFO,
51 TRACE,
52 WARN,
53 SEED,
54 TELL,
55 PSDUMPFREQ,
56 STATDUMPFREQ,
57 PSDUMPEACHTURN,
58 PSDUMPFRAME,
59 SPTDUMPFREQ,
60 REPARTFREQ,
61 REBINFREQ,
62 SCSOLVEFREQ,
63 MTSSUBSTEPS,
64 REMOTEPARTDEL,
65 RHODUMP,
66 EBDUMP,
67 CSRDUMP,
68 AUTOPHASE,
69 NUMBLOCKS,
70 RECYCLEBLOCKS,
71 NLHS,
72 CZERO,
73 RNGTYPE,
74 ENABLEHDF5,
75 ENABLEVTK,
76 ASCIIDUMP,
77 BOUNDPDESTROYFQ,
78 BEAMHALOBOUNDARY,
79 CLOTUNEONLY,
80 IDEALIZED,
81 LOGBENDTRAJECTORY,
82 VERSION,
83#ifdef ENABLE_AMR
84 AMR,
85 AMR_YT_DUMP_FREQ,
86 AMR_REGRID_FREQ,
87#endif
88 MEMORYDUMP,
89 HALOSHIFT,
90 DELPARTFREQ,
91 MINBINEMITTED,
92 MINSTEPFORREBIN,
93 COMPUTEPERCENTILES,
94 DUMPBEAMMATRIX,
95 SIZE
96 };
97}
98
99
101 Action(SIZE, "OPTION",
102 "The \"OPTION\" statement defines OPAL execution options.")
103 {
104
106 ("ECHO", "If true, give echo of input", echo);
107
109 ("INFO", "If true, print information messages", info);
110
112 ("TRACE", "If true, print execution trace"
113 "Must be the first option in the inputfile in "
114 "order to render correct results", mtrace);
115
117 ("WARN", "If true, print warning messages", warn);
118
120 ("SEED", "The seed for the random generator, -1 will use time(0) as seed ");
121
123 ("TELL", "If true, print the current settings. "
124 "Must be the last option in the inputfile in "
125 "order to render correct results", false);
126
127 itsAttr[PSDUMPFREQ] = Attributes::makeReal
128 ("PSDUMPFREQ", "The frequency to dump the phase space, "
129 "i.e.dump data when step%psDumpFreq==0, its default value is 10.",
130 psDumpFreq);
131
132 itsAttr[STATDUMPFREQ] = Attributes::makeReal
133 ("STATDUMPFREQ", "The frequency to dump statistical data "
134 "(e.g. RMS beam quantities), i.e. dump data when step%statDumpFreq == 0, "
135 "its default value is 10.", statDumpFreq);
136
137 itsAttr[PSDUMPEACHTURN] = Attributes::makeBool
138 ("PSDUMPEACHTURN", "If true, dump phase space after each "
139 "turn ,only aviable for OPAL-cycl, its default value is false",
141
142 itsAttr[SCSOLVEFREQ] = Attributes::makeReal
143 ("SCSOLVEFREQ", "The frequency to solve space charge fields. its default value is 1");
144
145 itsAttr[MTSSUBSTEPS] = Attributes::makeReal
146 ("MTSSUBSTEPS", "How many small timesteps "
147 "are inside the large timestep used in multiple "
148 "time stepping (MTS) integrator");
149
150 itsAttr[REMOTEPARTDEL] = Attributes::makeReal
151 ("REMOTEPARTDEL", "Artifically delete the remote particle "
152 "if its distance to the beam mass is larger than "
153 "REMOTEPARTDEL times of the beam rms size, "
154 "its default values is 0 (no delete) ",remotePartDel);
155
157 ("PSDUMPFRAME", "Controls the frame of phase space dump in "
158 "stat file and h5 file. If 'GLOBAL' OPAL will dump in the "
159 "lab (global) Cartesian frame; if 'BUNCH_MEAN' OPAL will "
160 "dump in the local Cartesian frame of the beam mean; "
161 "if 'REFERENCE' OPAL will dump in the local Cartesian "
162 "frame of the reference particle 0. Only available for "
163 "OPAL-cycl.", {"BUNCH_MEAN", "REFERENCE", "GLOBAL"}, "GLOBAL");
164
165 itsAttr[SPTDUMPFREQ] = Attributes::makeReal
166 ("SPTDUMPFREQ", "The frequency to dump single "
167 "particle trajectory of particles with ID = 0 & 1, "
168 "its default value is 1.", sptDumpFreq);
169
170 itsAttr[REPARTFREQ] = Attributes::makeReal
171 ("REPARTFREQ", "The frequency to do particles repartition "
172 "for better load balance between nodes, its "
173 "default value is " + std::to_string(repartFreq) + ".", repartFreq);
174
175 itsAttr[MINBINEMITTED] = Attributes::makeReal
176 ("MINBINEMITTED", "The number of bins that have to be emitted before the bins are squashed into "
177 "a single bin; the default value is " + std::to_string(minBinEmitted) + ".", minBinEmitted);
178
179 itsAttr[MINSTEPFORREBIN] = Attributes::makeReal
180 ("MINSTEPFORREBIN", "The number of steps into the simulation before the bins are squashed into "
181 "a single bin; the default value is " + std::to_string(minStepForRebin) + ".", minStepForRebin);
182
183 itsAttr[REBINFREQ] = Attributes::makeReal
184 ("REBINFREQ", "The frequency to reset energy bin ID for "
185 "all particles, its default value is 100.", rebinFreq);
186
188 ("RHODUMP", "If true, in addition to the phase "
189 "space the scalar rho field is also dumped (H5Block)", rhoDump);
190
192 ("EBDUMP", "If true, in addition to the phase space the "
193 "E and B field at each particle is also dumped into the H5 file)", ebDump);
194
196 ("CSRDUMP", "If true, the csr E field, line density "
197 "and the line density derivative is dumped into the "
198 "data directory)", csrDump);
199
200 itsAttr[AUTOPHASE] = Attributes::makeReal
201 ("AUTOPHASE", "If greater than zero OPAL is scanning "
202 "the phases of each rf structure in order to get maximum "
203 "acceleration. Defines the number of refinements of the "
204 "search range", autoPhase);
205
207 ("CZERO", "If set to true a symmetric distribution is "
208 "created -> centroid == 0.0", cZero);
209
211 ("RNGTYPE", "Type of pseudo- or quasi-random number generator, "
212 "see also Quasi-Random Sequences, GSL reference manual.",
213 {"RANDOM", "HALTON", "SOBOL", "NIEDERREITER"}, rngtype);
214
215 itsAttr[CLOTUNEONLY] = Attributes::makeBool
216 ("CLOTUNEONLY", "If set to true stop after "
217 "CLO and tune calculation ", cloTuneOnly);
218
219 itsAttr[NUMBLOCKS] = Attributes::makeReal
220 ("NUMBLOCKS", "Maximum number of vectors in the Krylov "
221 "space (for RCGSolMgr). Default value is 0 and BlockCGSolMgr will be used.");
222
223 itsAttr[RECYCLEBLOCKS] = Attributes::makeReal
224 ("RECYCLEBLOCKS", "Number of vectors in the recycle "
225 "space (for RCGSolMgr). Default value is 0 and BlockCGSolMgr will be used.");
226
228 ("NLHS", "Number of stored old solutions for extrapolating "
229 "the new starting vector. Default value is 1 and just the last solution is used.");
230
231 itsAttr[ENABLEHDF5] = Attributes::makeBool
232 ("ENABLEHDF5", "If true, HDF5 actions are enabled", enableHDF5);
233
234 itsAttr[ENABLEVTK] = Attributes::makeBool
235 ("ENABLEVTK", "If true, writing of VTK files are enabled", enableVTK);
236
237 itsAttr[ASCIIDUMP] = Attributes::makeBool
238 ("ASCIIDUMP", "If true, some of the elements dump in ASCII instead of HDF5", asciidump);
239
240 itsAttr[BOUNDPDESTROYFQ] = Attributes::makeReal
241 ("BOUNDPDESTROYFQ", "The frequency to do boundp_destroy to "
242 "delete lost particles. Default 10", boundpDestroyFreq);
243
244 itsAttr[BEAMHALOBOUNDARY] = Attributes::makeReal
245 ("BEAMHALOBOUNDARY", "Defines in terms of sigma where "
246 "the halo starts. Default 0.0", beamHaloBoundary);
247
248 itsAttr[IDEALIZED] = Attributes::makeBool
249 ("IDEALIZED", "Using the hard edge model for the calculation "
250 "of path length. Default: false", idealized);
251
252 itsAttr[LOGBENDTRAJECTORY] = Attributes::makeBool
253 ("LOGBENDTRAJECTORY", "Writing the trajectory of "
254 "every bend to disk. Default: false", writeBendTrajectories);
255
257 ("VERSION", "Version of OPAL for which input file was written", version);
258
259#ifdef ENABLE_AMR
261 ("AMR", "Use adaptive mesh refinement.", amr);
262
263 itsAttr[AMR_YT_DUMP_FREQ] = Attributes::makeReal("AMR_YT_DUMP_FREQ",
264 "The frequency to dump grid "
265 "and particle data "
266 "(default: 10)", amrYtDumpFreq);
267
268 itsAttr[AMR_REGRID_FREQ] = Attributes::makeReal("AMR_REGRID_FREQ",
269 "The frequency to perform a regrid "
270 "in multi-bunch mode (default: 10)",
272#endif
273
274 itsAttr[MEMORYDUMP] = Attributes::makeBool
275 ("MEMORYDUMP", "If true, write memory to SDDS file", memoryDump);
276
277 itsAttr[HALOSHIFT] = Attributes::makeReal
278 ("HALOSHIFT", "Constant parameter to shift halo value (default: 0.0)", haloShift);
279
280 itsAttr[DELPARTFREQ] = Attributes::makeReal
281 ("DELPARTFREQ", "The frequency to delete particles, "
282 "i.e. delete when step%delPartFreq == 0. Default: 1", delPartFreq);
283
284 itsAttr[COMPUTEPERCENTILES] = Attributes::makeBool
285 ("COMPUTEPERCENTILES", "Flag to control whether the 68.27 "
286 "(1 sigma for normal distribution), the 95.45 (2 sigmas), "
287 "the 99.73 (3 sigmas) and the 99.994 (4 sigmas) percentiles "
288 "for the beam size and the normalized emittance should "
289 "be computed. Default: false", computePercentiles);
290
291 itsAttr[DUMPBEAMMATRIX] = Attributes::makeBool
292 ("DUMPBEAMMATRIX", "Flag to control whether to write "
293 "the 6-dimensional beam matrix (upper triangle only) "
294 "to stat file. Default: false", dumpBeamMatrix);
295
297
299}
300
301
302Option::Option(const std::string& name, Option* parent):
303 Action(name, parent) {
327 Attributes::setPredefinedString(itsAttr[RNGTYPE], std::string(rngtype));
335 Attributes::setReal(itsAttr[BEAMHALOBOUNDARY], beamHaloBoundary);
339#ifdef ENABLE_AMR
341 Attributes::setReal(itsAttr[AMR_YT_DUMP_FREQ], amrYtDumpFreq);
342 Attributes::setReal(itsAttr[AMR_REGRID_FREQ], amrRegridFreq);
343#endif
347 Attributes::setBool(itsAttr[COMPUTEPERCENTILES], computePercentiles);
349}
350
351
354
355
356Option* Option::clone(const std::string& name) {
357 return new Option(name, this);
358}
359
360
362 // Store the option flags.
367 psDumpEachTurn = Attributes::getBool(itsAttr[PSDUMPEACHTURN]);
368 remotePartDel = Attributes::getReal(itsAttr[REMOTEPARTDEL]);
379#ifdef ENABLE_AMR
381 amrYtDumpFreq = int(Attributes::getReal(itsAttr[AMR_YT_DUMP_FREQ]));
382
383 if ( amrYtDumpFreq < 1 ) {
384 amrYtDumpFreq = std::numeric_limits<int>::max();
385 }
386
387 amrRegridFreq = int(Attributes::getReal(itsAttr[AMR_REGRID_FREQ]));
389#endif
390
394 computePercentiles = Attributes::getBool(itsAttr[COMPUTEPERCENTILES]);
395 dumpBeamMatrix = Attributes::getBool(itsAttr[DUMPBEAMMATRIX]);
396 if ( memoryDump ) {
399 memory->sample();
400 }
401
404
405 if (Options::seed == -1)
406 rangen.init55(time(0));
407 else
409
410 gmsgALL->on(info);
411 gmsg->on(info);
414
416
419 if (itsAttr[SEED]) {
420 seed = int(Attributes::getReal(itsAttr[SEED]));
421 if (seed == -1)
422 rangen.init55(time(0));
423 else
425 }
426
427 if (itsAttr[PSDUMPFREQ]) {
428 psDumpFreq = int(Attributes::getReal(itsAttr[PSDUMPFREQ]));
429 if (psDumpFreq==0)
430 psDumpFreq = std::numeric_limits<int>::max();
431 }
432
433 if (itsAttr[STATDUMPFREQ]) {
434 statDumpFreq = int(Attributes::getReal(itsAttr[STATDUMPFREQ]));
435 if (statDumpFreq==0)
436 statDumpFreq = std::numeric_limits<int>::max();
437 }
438
439 if (itsAttr[SPTDUMPFREQ]) {
440 sptDumpFreq = int(Attributes::getReal(itsAttr[SPTDUMPFREQ]));
441 if (sptDumpFreq==0)
442 sptDumpFreq = std::numeric_limits<int>::max();
443 }
444
445 if (itsAttr[SCSOLVEFREQ]) {
446 scSolveFreq = int(Attributes::getReal(itsAttr[SCSOLVEFREQ]));
447 scSolveFreq = ( scSolveFreq < 1 ) ? 1 : scSolveFreq;
448 }
449
450 if (itsAttr[MTSSUBSTEPS]) {
451 mtsSubsteps = int(Attributes::getReal(itsAttr[MTSSUBSTEPS]));
452 }
453
454 if (itsAttr[REPARTFREQ]) {
455 repartFreq = int(Attributes::getReal(itsAttr[REPARTFREQ]));
456 }
457
458 if (itsAttr[MINBINEMITTED]) {
459 minBinEmitted = int(Attributes::getReal(itsAttr[MINBINEMITTED]));
460 }
461
462 if (itsAttr[MINSTEPFORREBIN]) {
463 minStepForRebin = int(Attributes::getReal(itsAttr[MINSTEPFORREBIN]));
464 }
465
466 if (itsAttr[REBINFREQ]) {
467 rebinFreq = int(Attributes::getReal(itsAttr[REBINFREQ]));
468 }
469
470 if (itsAttr[AUTOPHASE]) {
471 autoPhase = int(Attributes::getReal(itsAttr[AUTOPHASE]));
472 }
473
474 if (itsAttr[NUMBLOCKS]) {
475 numBlocks = int(Attributes::getReal(itsAttr[NUMBLOCKS]));
476 }
477
478 if (itsAttr[RECYCLEBLOCKS]) {
479 recycleBlocks = int(Attributes::getReal(itsAttr[RECYCLEBLOCKS]));
480 }
481
482 if (itsAttr[NLHS]) {
483 nLHS = int(Attributes::getReal(itsAttr[NLHS]));
484 }
485
486 if (itsAttr[CZERO]) {
487 cZero = bool(Attributes::getBool(itsAttr[CZERO]));
488 }
489
490 if (itsAttr[RNGTYPE]) {
491 rngtype = std::string(Attributes::getString(itsAttr[RNGTYPE]));
492 } else {
493 rngtype = std::string("RANDOM");
494 }
495
496 if (itsAttr[BEAMHALOBOUNDARY]) {
497 beamHaloBoundary = Attributes::getReal(itsAttr[BEAMHALOBOUNDARY]);
498 } else {
500 }
501
502 if (itsAttr[CLOTUNEONLY]) {
503 cloTuneOnly = bool(Attributes::getBool(itsAttr[CLOTUNEONLY]));
504 } else {
505 cloTuneOnly = false;
506 }
507
508 // Set message flags.
510
511 if (Attributes::getBool(itsAttr[TELL])) {
512 *gmsg << "\nCurrent settings of options:\n" << *this << endl;
513 }
514
515 Option* main = dynamic_cast<Option*>(OpalData::getInstance()->find("OPTION"));
516 if (main) {
518 }
519}
520
521constexpr std::array<std::pair<DumpFrame, std::string_view>, 3> dumpFrameMap {{
522 {DumpFrame::GLOBAL, "GLOBAL"},
523 {DumpFrame::BUNCH_MEAN, "BUNCH_MEAN"},
524 {DumpFrame::REFERENCE, "REFERENCE"}
525}};
526
527void Option::handlePsDumpFrame(std::string_view dumpFrameName) noexcept {
529}
530
531std::string Option::getDumpFrameString(const DumpFrame& df) noexcept {
532 return std::string(Util::enumToString(df, dumpFrameMap, "GLOBAL"));
533}
534
535void Option::update(const std::vector<Attribute>& othersAttributes) {
536 for (int i = 0; i < SIZE; ++ i) {
537 itsAttr[i] = othersAttributes[i];
538 }
539}
DumpFrame
Definition Options.h:26
constexpr std::array< std::pair< DumpFrame, std::string_view >, 3 > dumpFrameMap
Definition Option.cpp:521
Inform * gmsgALL
Definition Main.cpp:71
Inform * gmsg
Definition Main.cpp:70
@ SIZE
Definition IndexMap.cpp:174
Inform * gmsgALL
Definition Main.cpp:71
Inform * gmsg
Definition Main.cpp:70
int main(int argc, char *argv[])
Definition Main.cpp:139
Inform & endl(Inform &inf)
Definition Inform.cpp:42
const std::string name
Some AMR types used a lot.
Definition AmrDefs.h:33
Attribute makeBool(const std::string &name, const std::string &help)
Make logical attribute.
double getReal(const Attribute &attr)
Return real value.
void setBool(Attribute &attr, bool val)
Set logical 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.
void setReal(Attribute &attr, double val)
Set real value.
void setPredefinedString(Attribute &attr, const std::string &val)
Set predefined string value.
std::string getString(const Attribute &attr)
Get string value.
bool computePercentiles
Definition Options.cpp:109
int mtsSubsteps
Definition Options.cpp:57
double remotePartDel
Definition Options.cpp:59
int psDumpFreq
The frequency to dump the phase space, i.e.dump data when steppsDumpFreq==0.
Definition Options.cpp:37
double haloShift
The constant parameter C to shift halo, by < w^4 > / < w^2 > ^2 - C (w=x,y,z)
Definition Options.cpp:105
bool writeBendTrajectories
Definition Options.cpp:93
bool mtrace
Trace flag.
Definition Options.cpp:29
int boundpDestroyFreq
Definition Options.cpp:85
bool memoryDump
Definition Options.cpp:103
bool enableVTK
If true VTK files are written.
Definition Options.cpp:81
int version
opal version of input file
Definition Options.cpp:95
int minBinEmitted
The number of bins that have to be emitted before the bin are squashed into a single bin.
Definition Options.cpp:49
bool enableHDF5
If true HDF5 files are written.
Definition Options.cpp:79
bool psDumpEachTurn
phase space dump flag for OPAL-cycl
Definition Options.cpp:41
bool asciidump
Definition Options.cpp:83
int amrRegridFreq
After how many steps the AMR grid hierarchy is updated.
Definition Options.cpp:101
int numBlocks
RCG: cycle length.
Definition Options.cpp:69
int sptDumpFreq
The frequency to dump single particle trajectory of particles with ID = 0 & 1.
Definition Options.cpp:45
int autoPhase
Definition Options.cpp:67
bool warn
Warn flag.
Definition Options.cpp:31
std::string rngtype
random number generator
Definition Options.cpp:77
bool echo
Echo flag.
Definition Options.cpp:24
bool rhoDump
Definition Options.cpp:61
int minStepForRebin
The number of steps into the simulation before the bins are squashed into a single bin.
Definition Options.cpp:51
bool cloTuneOnly
Do closed orbit and tune calculation only.
Definition Options.cpp:89
bool csrDump
Definition Options.cpp:65
bool ebDump
Definition Options.cpp:63
bool idealized
Definition Options.cpp:91
int amrYtDumpFreq
The frequency to dump AMR grid data and particles into file.
Definition Options.cpp:99
unsigned int delPartFreq
The frequency to delete particles (currently: OPAL-cycl only)
Definition Options.cpp:107
bool dumpBeamMatrix
Definition Options.cpp:111
int scSolveFreq
The frequency to solve space charge fields.
Definition Options.cpp:55
int repartFreq
The frequency to do particles repartition for better load balance between nodes.
Definition Options.cpp:47
bool cZero
If true create symmetric distribution.
Definition Options.cpp:75
double beamHaloBoundary
Definition Options.cpp:87
int nLHS
number of old left hand sides used to extrapolate a new start vector
Definition Options.cpp:73
int rebinFreq
The frequency to reset energy bin ID for all particles.
Definition Options.cpp:53
bool info
Info flag.
Definition Options.cpp:26
Random rangen
Definition Options.cpp:34
int seed
The current random seed.
Definition Options.cpp:35
DumpFrame psDumpFrame
flag to decide in which coordinate frame the phase space will be dumped for OPAL-cycl
Definition Options.cpp:43
int statDumpFreq
The frequency to dump statistical values, e.e. dump data when stepstatDumpFreq==0.
Definition Options.cpp:39
int recycleBlocks
RCG: number of recycle blocks.
Definition Options.cpp:71
constexpr Enum stringToEnum(std::string_view str, const std::array< std::pair< Enum, std::string_view >, N > &map, Enum defaultEnum) noexcept
Definition Util.h:237
constexpr std::string_view enumToString(Enum e, const std::array< std::pair< Enum, std::string_view >, N > &map, std::string_view defaultStr) noexcept
Definition Util.h:247
The base class for all OPAL actions.
Definition Action.h:30
void registerOwnership(const AttributeHandler::OwnerType &itsClass) const
Definition Object.cpp:191
std::vector< Attribute > itsAttr
The object attributes.
Definition Object.h:216
Object * find(const std::string &name)
Find entry.
Definition OpalData.cpp:571
static OpalData * getInstance()
Definition OpalData.cpp:196
virtual ~Option()
Definition Option.cpp:352
virtual void execute()
Execute the command.
Definition Option.cpp:361
void update(const std::vector< Attribute > &)
Definition Option.cpp:535
virtual void update()
Update this object.
Definition Object.cpp:263
static std::string getDumpFrameString(const DumpFrame &df) noexcept
Definition Option.cpp:531
void handlePsDumpFrame(std::string_view dumpFrameName) noexcept
Definition Option.cpp:527
virtual Option * clone(const std::string &name)
Make clone.
Definition Option.cpp:356
Option()
Definition Option.cpp:100
static void setEcho(bool flag)
Set echo flag.
void init55(int seed)
Initialise random number generator.
void on(const bool o)
Definition Inform.h:69
static Inform * Warn
Definition IpplInfo.h:79
static Inform * Info
Definition IpplInfo.h:78
static IpplMemory_p getInstance(Unit unit=Unit::GB, bool reset=true)