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