OPALX (Object Oriented Parallel Accelerator Library for Exascal) master (dc2a29eed580)
OPALX
Loading...
Searching...
No Matches
StatWriter.cpp
Go to the documentation of this file.
1//
2// Class StatWriter
3// This class writes bunch statistics (*.stat).
4//
5// One StatWriter instance corresponds to one output file. DataSink::init creates
6// statWriters_m.size() == numParticleContainers writers: for a single container the
7// file stem is the input basename (e.g. myjob -> myjob.stat); for multiple containers
8// stems are basename + "_c" + index (run_c0.stat, run_c1.stat, ...), see
9// DataSink::diagnosticStemForContainer. Each write() call must use the same
10// particleContainerIndex as that writer's slot so row data comes from
11// beam->getParticleContainer(particleContainerIndex). Shared beam-level quantities
12// (e.g. time t, dt, rmsDensity, nBins) still come from PartBunch_t regardless of index.
13//
14// Copyright (c) 2019, Matthias Frey, Paul Scherrer Institut, Villigen PSI, Switzerland
15// Christof Metzger-Kraus, Open Sourcerer
16// All rights reserved
17//
18// This file is part of OPAL.
19//
20// OPAL is free software: you can redistribute it and/or modify
21// it under the terms of the GNU General Public License as published by
22// the Free Software Foundation, either version 3 of the License, or
23// (at your option) any later version.
24//
25// You should have received a copy of the GNU General Public License
26// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
27//
28#include "StatWriter.h"
29
30#include <Kokkos_Core.hpp>
31
33#include "BuildInfo.h"
34#include "PartBunch/PartBunch.h"
35#include "Physics/Units.h"
36#include "Utilities/Timer.h"
37#include "Utilities/Util.h"
38
39#include <sstream>
40
41StatWriter::StatWriter(const std::string& fname, bool restart) : StatBaseWriter(fname, restart) {}
42
43void StatWriter::fillHeader(const losses_t& losses, const std::string& species) {
44 if (this->hasColumns()) {
45 return;
46 }
47
48 columns_m.addColumn("t", "double", "ns", "Time");
49 columns_m.addColumn("s", "double", "m", "Path length");
50 columns_m.addColumn("numParticles", "long", "1", "Number of Macro Particles");
51 columns_m.addColumn("charge", "double", "1", "Bunch Charge");
52 columns_m.addColumn("energy", "double", "MeV", "Mean Bunch Energy");
53
54 columns_m.addColumn("rms_x", "double", "m", "RMS Beamsize in x");
55 columns_m.addColumn("rms_y", "double", "m", "RMS Beamsize in y");
56 columns_m.addColumn("rms_s", "double", "m", "RMS Beamsize in s");
57
58 columns_m.addColumn("rms_px", "double", "1", "RMS Normalized Momenta in x");
59 columns_m.addColumn("rms_py", "double", "1", "RMS Normalized Momenta in y");
60 columns_m.addColumn("rms_ps", "double", "1", "RMS Normalized Momenta in s");
61
62 columns_m.addColumn("emit_x", "double", "m", "Normalized Emittance x");
63 columns_m.addColumn("emit_y", "double", "m", "Normalized Emittance y");
64 columns_m.addColumn("emit_s", "double", "m", "Normalized Emittance s");
65
66 columns_m.addColumn("mean_x", "double", "m", "Mean Beam Position in x");
67 columns_m.addColumn("mean_y", "double", "m", "Mean Beam Position in y");
68 columns_m.addColumn("mean_s", "double", "m", "Mean Beam Position in s");
69
70 columns_m.addColumn("ref_x", "double", "m", "x coordinate of reference particle in lab cs");
71 columns_m.addColumn("ref_y", "double", "m", "y coordinate of reference particle in lab cs");
72 columns_m.addColumn("ref_z", "double", "m", "z coordinate of reference particle in lab cs");
73
74 columns_m.addColumn("ref_px", "double", "1", "x momentum of reference particle in lab cs");
75 columns_m.addColumn("ref_py", "double", "1", "y momentum of reference particle in lab cs");
76 columns_m.addColumn("ref_pz", "double", "1", "z momentum of reference particle in lab cs");
77
78 columns_m.addColumn("max_x", "double", "m", "Max Beamsize in x");
79 columns_m.addColumn("max_y", "double", "m", "Max Beamsize in y");
80 columns_m.addColumn("max_s", "double", "m", "Max Beamsize in s");
81
82 columns_m.addColumn("xpx", "double", "1", "Correlation xpx");
83 columns_m.addColumn("ypy", "double", "1", "Correlation ypy");
84 columns_m.addColumn("zpz", "double", "1", "Correlation zpz");
85
86 columns_m.addColumn("Dx", "double", "m", "Dispersion in x");
87 columns_m.addColumn("DDx", "double", "1", "Derivative of dispersion in x");
88 columns_m.addColumn("Dy", "double", "m", "Dispersion in y");
89 columns_m.addColumn("DDy", "double", "1", "Derivative of dispersion in y");
90
91 columns_m.addColumn("Bx_ref", "double", "T", "Bx-Field component of ref particle");
92 columns_m.addColumn("By_ref", "double", "T", "By-Field component of ref particle");
93 columns_m.addColumn("Bz_ref", "double", "T", "Bz-Field component of ref particle");
94
95 columns_m.addColumn("Ex_ref", "double", "MV/m", "Ex-Field component of ref particle");
96 columns_m.addColumn("Ey_ref", "double", "MV/m", "Ey-Field component of ref particle");
97 columns_m.addColumn("Ez_ref", "double", "MV/m", "Ez-Field component of ref particle");
98
99 columns_m.addColumn("dE", "double", "MeV", "energy spread of the beam");
100 columns_m.addColumn("dt", "double", "ns", "time step size");
101 columns_m.addColumn("partsOutside", "double", "1", "outside n*sigma of the beam");
102
103 columns_m.addColumn("DebyeLength", "double", "m", "Debye length in the boosted frame");
105 "plasmaParameter", "double", "1",
106 "Plasma parameter that gives no. of particles in a Debye sphere");
107 columns_m.addColumn("temperature", "double", "K", "Temperature of the beam");
108 columns_m.addColumn("rmsDensity", "double", "1", "RMS number density of the beam");
109
111 "nBins", "int", "1",
112 "Number of field solver bins, potentially after adaptive binning.");
113
115 /*
116 if (Options::computePercentiles) {
117 columns_m.addColumn("68_Percentile_x", "double", "m",
118 "68.27 percentile (1 sigma of normal distribution) of x-component of
119 position"); columns_m.addColumn("68_Percentile_y", "double", "m", "68.27 percentile (1 sigma of
120 normal distribution) of y-component of position"); columns_m.addColumn("68_Percentile_z",
121 "double", "m", "68.27 percentile (1 sigma of normal distribution) of z-component of position");
122
123 columns_m.addColumn("95_Percentile_x", "double", "m",
124 "95.45 percentile (2 sigma of normal distribution) of x-component of
125 position"); columns_m.addColumn("95_Percentile_y", "double", "m", "95.45 percentile (2 sigma of
126 normal distribution) of y-component of position"); columns_m.addColumn("95_Percentile_z",
127 "double", "m", "95.45 percentile (2 sigma of normal distribution) of z-component of position");
128
129 columns_m.addColumn("99_Percentile_x", "double", "m",
130 "99.73 percentile (3 sigma of normal distribution) of x-component of
131 position"); columns_m.addColumn("99_Percentile_y", "double", "m", "99.73 percentile (3 sigma of
132 normal distribution) of y-component of position"); columns_m.addColumn("99_Percentile_z",
133 "double", "m", "99.73 percentile (3 sigma of normal distribution) of z-component of position");
134
135 columns_m.addColumn("99_99_Percentile_x", "double", "m",
136 "99.994 percentile (4 sigma of normal distribution) of x-component of
137 position"); columns_m.addColumn("99_99_Percentile_y", "double", "m", "99.994 percentile (4 sigma
138 of normal distribution) of y-component of position"); columns_m.addColumn("99_99_Percentile_z",
139 "double", "m", "99.994 percentile (4 sigma of normal distribution) of z-component of position");
140
141 columns_m.addColumn("normalizedEps68Percentile_x", "double", "m",
142 "x-component of normalized emittance at 68 percentile (1 sigma of normal
143 distribution)"); columns_m.addColumn("normalizedEps68Percentile_y", "double", "m", "y-component
144 of normalized emittance at 68 percentile (1 sigma of normal distribution)");
145 columns_m.addColumn("normalizedEps68Percentile_z", "double", "m",
146 "z-component of normalized emittance at 68 percentile (1 sigma of normal
147 distribution)");
148
149 columns_m.addColumn("normalizedEps95Percentile_x", "double", "m",
150 "x-component of normalized emittance at 95 percentile (2 sigma of normal
151 distribution)"); columns_m.addColumn("normalizedEps95Percentile_y", "double", "m", "y-component
152 of normalized emittance at 95 percentile (2 sigma of normal distribution)");
153 columns_m.addColumn("normalizedEps95Percentile_z", "double", "m",
154 "z-component of normalized emittance at 95 percentile (2 sigma of normal
155 distribution)");
156
157 columns_m.addColumn("normalizedEps99Percentile_x", "double", "m",
158 "x-component of normalized emittance at 99 percentile (3 sigma of normal
159 distribution)"); columns_m.addColumn("normalizedEps99Percentile_y", "double", "m", "y-component
160 of normalized emittance at 99 percentile (3 sigma of normal distribution)");
161 columns_m.addColumn("normalizedEps99Percentile_z", "double", "m",
162 "z-component of normalized emittance at 99 percentile (3 sigma of normal
163 distribution)");
164
165 columns_m.addColumn("normalizedEps99_99Percentile_x", "double", "m",
166 "x-component of normalized emittance at 99.99 percentile (4 sigma of
167 normal distribution)"); columns_m.addColumn("normalizedEps99_99Percentile_y", "double", "m",
168 "y-component of normalized emittance at 99.99 percentile (4 sigma of
169 normal distribution)"); columns_m.addColumn("normalizedEps99_99Percentile_z", "double", "m",
170 "z-component of normalized emittance at 99.99 percentile (4 sigma of
171 normal distribution)");
172 }
173 */
174 if (OpalData::getInstance()->isInOPALCyclMode() && ippl::Comm->size() == 1) {
175 columns_m.addColumn("R0_x", "double", "m", "R0 Particle position in x");
176 columns_m.addColumn("R0_y", "double", "m", "R0 Particle position in y");
177 columns_m.addColumn("R0_s", "double", "m", "R0 Particle position in z");
178
179 columns_m.addColumn("P0_x", "double", "1", "R0 Particle momentum in x");
180 columns_m.addColumn("P0_y", "double", "1", "R0 Particle momentum in y");
181 columns_m.addColumn("P0_s", "double", "1", "R0 Particle momentum in z");
182 }
183
184 if (OpalData::getInstance()->isInOPALCyclMode()) {
185 columns_m.addColumn("halo_x", "double", "1", "Halo in x");
186 columns_m.addColumn("halo_y", "double", "1", "Halo in y");
187 columns_m.addColumn("halo_z", "double", "1", "Halo in z");
188
189 columns_m.addColumn("azimuth", "double", "deg", "Azimuth in global coordinates");
190 }
191
192 for (size_t i = 0; i < losses.size(); ++i) {
193 columns_m.addColumn(losses[i].first, "long", "1", "Number of lost particles in element");
194 }
195
196 if (mode_m == std::ios::app) return;
197
198 OPALTimer::Timer simtimer;
199 std::string dateStr(simtimer.date());
200 std::string timeStr(simtimer.time());
201
202 std::stringstream ss;
203 ss << "Statistics data '" << OpalData::getInstance()->getInputFn() << "' " << dateStr << " "
204 << timeStr;
205
206 this->addDescription(ss.str(), "stat parameters");
207
208 std::stringstream revision;
209 revision << buildinfo::project_name << " " << buildinfo::project_version << " "
210 << "git rev. #" << Util::getGitRevision();
211
212 addParameter("processors", "long", "Number of Cores used", ippl::Comm->size());
213 addParameter("revision", "string", "git revision of opal", revision.str());
214 addParameter("species", "string", "Particle species of container", species);
215
216 this->addInfo("ascii", 1);
217}
218
220 PartBunch_t& beam, Vector_t<double, 3> FDext[], const losses_t& losses,
221 const double& azimuth, const size_t npOutside, size_t particleContainerIndex) {
223 std::shared_ptr<ParticleContainer_t> pc = beam.getParticleContainer(particleContainerIndex);
224 if (!pc) {
225 return;
226 }
227
228 double pathLength = pc->get_sPos();
229 const std::string species = beam.getParticleName(particleContainerIndex);
230
231 // First write to this writer's .stat file emits SDDS header via fillHeader/writeHeader.
232 // File vs. container: this object was constructed with the stem for particleContainerIndex;
233 // pc must be beam->getParticleContainer(particleContainerIndex) (caller responsibility).
234
235 const size_t numParticles = pc->getTotalNum();
236 double Q = numParticles == 0 ? 0.0 : pc->getTotalCharge();
237
238 if (ippl::Comm->rank() != 0) {
239 return;
240 }
241
242 fillHeader(losses, species);
243
244 this->open();
245
246 this->writeHeader();
247
248 Vector_t<double, 3> rmsR(0.0);
249 Vector_t<double, 3> rmsP(0.0);
250 Vector_t<double, 3> normEmit(0.0);
251 Vector_t<double, 3> meanR(0.0);
252 Vector_t<double, 3> maxR(0.0);
253 Vector_t<double, 3> rmsRP(0.0);
254 double energy = 0.0;
255 double dE = 0.0;
256 double dx = 0.0;
257 double ddx = 0.0;
258 double dy = 0.0;
259 double ddy = 0.0;
260 double debyeLength = 0.0;
261 double plasmaParameter = 0.0;
262 double temperature = 0.0;
263
264 if (numParticles > 0) {
265 energy = pc->getMeanKineticEnergy();
266 rmsR = pc->getRmsR();
267 rmsP = pc->getRmsP();
268 normEmit = pc->getNormEmit();
269 meanR = pc->getMeanR();
270 maxR = pc->getMaxR();
271 rmsRP = pc->getRmsRP();
272 dE = pc->getStdKineticEnergy();
273 dx = pc->getDx();
274 ddx = pc->getDDx();
275 dy = pc->getDy();
276 ddy = pc->getDDy();
277 debyeLength = pc->getDebyeLength();
278 plasmaParameter = pc->getPlasmaParameter();
279 temperature = pc->getTemperature();
280 }
281
282 columns_m.addColumnValue("t", beam.getT() * Units::s2ns); // 1
283 columns_m.addColumnValue("s", pathLength); // 2
284 columns_m.addColumnValue("numParticles", numParticles); // 3
285 columns_m.addColumnValue("charge", Q); // 4
286 columns_m.addColumnValue("energy", energy); // 5
287
288 columns_m.addColumnValue("rms_x", rmsR(0)); // 6
289 columns_m.addColumnValue("rms_y", rmsR(1)); // 7
290 columns_m.addColumnValue("rms_s", rmsR(2)); // 8
291
292 columns_m.addColumnValue("rms_px", rmsP(0)); // 9
293 columns_m.addColumnValue("rms_py", rmsP(1)); // 10
294 columns_m.addColumnValue("rms_ps", rmsP(2)); // 11
295
296 columns_m.addColumnValue("emit_x", normEmit(0)); // 12
297 columns_m.addColumnValue("emit_y", normEmit(1)); // 13
298 columns_m.addColumnValue("emit_s", normEmit(2)); // 14
299
300 columns_m.addColumnValue("mean_x", meanR(0)); // 15
301 columns_m.addColumnValue("mean_y", meanR(1)); // 16
302 columns_m.addColumnValue("mean_s", meanR(2)); // 17
303
304 columns_m.addColumnValue("ref_x", pc->getRefPartR()(0)); // 18
305 columns_m.addColumnValue("ref_y", pc->getRefPartR()(1)); // 19
306 columns_m.addColumnValue("ref_z", pc->getRefPartR()(2)); // 20
307
308 columns_m.addColumnValue("ref_px", pc->getRefPartP()(0)); // 21
309 columns_m.addColumnValue("ref_py", pc->getRefPartP()(1)); // 22
310 columns_m.addColumnValue("ref_pz", pc->getRefPartP()(2)); // 23
311
312 columns_m.addColumnValue("max_x", maxR(0)); // 24
313 columns_m.addColumnValue("max_y", maxR(1)); // 25
314 columns_m.addColumnValue("max_s", maxR(2)); // 26
315
316 // Write out Courant Snyder parameters.
317 columns_m.addColumnValue("xpx", rmsRP(0)); // 27
318 columns_m.addColumnValue("ypy", rmsRP(1)); // 28
319 columns_m.addColumnValue("zpz", rmsRP(2)); // 29
320
321 // Write out dispersion.
322 columns_m.addColumnValue("Dx", dx); // 30
323 columns_m.addColumnValue("DDx", ddx); // 31
324 columns_m.addColumnValue("Dy", dy); // 32
325 columns_m.addColumnValue("DDy", ddy); // 33
326
327 // Write head/reference particle/tail field information.
328 columns_m.addColumnValue("Bx_ref", FDext[0](0)); // 34 B-ref x
329 columns_m.addColumnValue("By_ref", FDext[0](1)); // 35 B-ref y
330 columns_m.addColumnValue("Bz_ref", FDext[0](2)); // 36 B-ref z
331
332 columns_m.addColumnValue("Ex_ref", FDext[1](0)); // 37 E-ref x
333 columns_m.addColumnValue("Ey_ref", FDext[1](1)); // 38 E-ref y
334 columns_m.addColumnValue("Ez_ref", FDext[1](2)); // 39 E-ref z
335
336 columns_m.addColumnValue("dE", dE); // 40 dE energy spread
337 columns_m.addColumnValue("dt", beam.getdT() * Units::s2ns); // 41 dt time step size
338 columns_m.addColumnValue("partsOutside", npOutside); // 42 number of particles outside n*sigma
339
340 columns_m.addColumnValue("DebyeLength", debyeLength); // 43 Debye length
341 columns_m.addColumnValue("plasmaParameter", plasmaParameter); // 43 plasma parameter
342 columns_m.addColumnValue("temperature", temperature); // 44 Temperature
343 columns_m.addColumnValue("rmsDensity", beam.get_rmsDensity()); // 45 RMS number density
345
346 /*
347 if (Options::computePercentiles) {
348 columns_m.addColumnValue("68_Percentile_x", beam->get_68Percentile()[0]);
349 columns_m.addColumnValue("68_Percentile_y", beam->get_68Percentile()[1]);
350 columns_m.addColumnValue("68_Percentile_z", beam->get_68Percentile()[2]);
351
352 columns_m.addColumnValue("95_Percentile_x", beam->get_95Percentile()[0]);
353 columns_m.addColumnValue("95_Percentile_y", beam->get_95Percentile()[1]);
354 columns_m.addColumnValue("95_Percentile_z", beam->get_95Percentile()[2]);
355
356 columns_m.addColumnValue("99_Percentile_x", beam->get_99Percentile()[0]);
357 columns_m.addColumnValue("99_Percentile_y", beam->get_99Percentile()[1]);
358 columns_m.addColumnValue("99_Percentile_z", beam->get_99Percentile()[2]);
359
360 columns_m.addColumnValue("99_99_Percentile_x", beam->get_99_99Percentile()[0]);
361 columns_m.addColumnValue("99_99_Percentile_y", beam->get_99_99Percentile()[1]);
362 columns_m.addColumnValue("99_99_Percentile_z", beam->get_99_99Percentile()[2]);
363
364 columns_m.addColumnValue(
365 "normalizedEps68Percentile_x", beam->get_normalizedEps_68Percentile()[0]);
366 columns_m.addColumnValue(
367 "normalizedEps68Percentile_y", beam->get_normalizedEps_68Percentile()[1]);
368 columns_m.addColumnValue(
369 "normalizedEps68Percentile_z", beam->get_normalizedEps_68Percentile()[2]);
370
371 columns_m.addColumnValue(
372 "normalizedEps95Percentile_x", beam->get_normalizedEps_95Percentile()[0]);
373 columns_m.addColumnValue(
374 "normalizedEps95Percentile_y", beam->get_normalizedEps_95Percentile()[1]);
375 columns_m.addColumnValue(
376 "normalizedEps95Percentile_z", beam->get_normalizedEps_95Percentile()[2]);
377
378 columns_m.addColumnValue(
379 "normalizedEps99Percentile_x", beam->get_normalizedEps_99Percentile()[0]);
380 columns_m.addColumnValue(
381 "normalizedEps99Percentile_y", beam->get_normalizedEps_99Percentile()[1]);
382 columns_m.addColumnValue(
383 "normalizedEps99Percentile_z", beam->get_normalizedEps_99Percentile()[2]);
384
385 columns_m.addColumnValue(
386 "normalizedEps99_99Percentile_x", beam->get_normalizedEps_99_99Percentile()[0]);
387 columns_m.addColumnValue(
388 "normalizedEps99_99Percentile_y", beam->get_normalizedEps_99_99Percentile()[1]);
389 columns_m.addColumnValue(
390 "normalizedEps99_99Percentile_z", beam->get_normalizedEps_99_99Percentile()[2]);
391 }
392 */
393 if (OpalData::getInstance()->isInOPALCyclMode()) {
394 if (ippl::Comm->size() == 1) {
395 if (pc->getLocalNum() > 0) {
396 auto rDev = pc->R.getView();
397 auto pDev = pc->P.getView();
398 auto rHost = Kokkos::create_mirror_view(rDev);
399 auto pHost = Kokkos::create_mirror_view(pDev);
400 Kokkos::deep_copy(rHost, rDev);
401 Kokkos::deep_copy(pHost, pDev);
402 columns_m.addColumnValue("R0_x", rHost(0)[0]);
403 columns_m.addColumnValue("R0_y", rHost(0)[1]);
404 columns_m.addColumnValue("R0_s", rHost(0)[2]);
405 columns_m.addColumnValue("P0_x", pHost(0)[0]);
406 columns_m.addColumnValue("P0_y", pHost(0)[1]);
407 columns_m.addColumnValue("P0_s", pHost(0)[2]);
408 } else {
409 columns_m.addColumnValue("R0_x", 0.0);
410 columns_m.addColumnValue("R0_y", 0.0);
411 columns_m.addColumnValue("R0_s", 0.0);
412 columns_m.addColumnValue("P0_x", 0.0);
413 columns_m.addColumnValue("P0_y", 0.0);
414 columns_m.addColumnValue("P0_s", 0.0);
415 }
416 }
417 Vector_t<double, 3> halo = beam.get_halo();
418 columns_m.addColumnValue("halo_x", halo(0));
419 columns_m.addColumnValue("halo_y", halo(1));
420 columns_m.addColumnValue("halo_z", halo(2));
421
422 columns_m.addColumnValue("azimuth", azimuth);
423 }
424
425 for (size_t i = 0; i < losses.size(); ++i) {
426 long unsigned int loss = losses[i].second;
427 columns_m.addColumnValue(losses[i].first, loss);
428 }
429
430 this->writeRow();
431
432 this->close();
433}
ippl::Vector< T, Dim > Vector_t
Template PIC bunch: IPPL PicManager, shared field mesh/solver, and multiple particle containers.
std::string time() const
Return time.
Definition Timer.cpp:36
std::string date() const
Return date.
Definition Timer.cpp:30
std::string getInputFn()
get opals input filename
Definition OpalData.cpp:569
static OpalData * getInstance()
Definition OpalData.cpp:193
std::string getParticleName(size_t i) const
Particle species name for container i (from BEAM PARTICLE input).
Definition PartBunch.h:603
double getdT() const
Get the global time step.
Definition PartBunch.h:642
double getT() const
Get the current simulation time.
Definition PartBunch.h:651
int getCurrentNBins() const
Effective bin count for diagnostics (1 if binning inactive or still at max bins).
Definition PartBunch.h:574
double get_rmsDensity() const
Legacy RMS density field (may be unused).
Definition PartBunch.h:682
Vector_t< double, Dim > get_halo() const
Stub; logs and returns zero vector.
Definition PartBunch.h:660
Container for all per-particle (and per-simulation) fields tracked during OPALX tracking.
void addColumn(const std::string &name, const std::string &type, const std::string &unit, const std::string &desc, std::ios_base::fmtflags flags=std::ios_base::scientific, unsigned short precision=15)
void addColumnValue(const std::string &name, const T &val)
SDDSColumnSet columns_m
Definition SDDSWriter.h:111
bool hasColumns() const
Definition SDDSWriter.h:168
void addDescription(const std::string &text, const std::string &content)
Definition SDDSWriter.h:141
void writeHeader()
Write SDDS header.
std::ios_base::openmode mode_m
First write to the statistics output file.
Definition SDDSWriter.h:109
void writeRow()
Definition SDDSWriter.h:158
void addInfo(const std::string &mode, const size_t &no_row_counts)
Definition SDDSWriter.h:154
void addParameter(const std::string &name, const std::string &type, const std::string &desc, const T &value)
Definition SDDSWriter.h:146
std::vector< std::pair< std::string, unsigned int > > losses_t
Definition StatWriter.h:26
void fillHeader(const losses_t &losses=losses_t(), const std::string &species="")
void write(PartBunch_t &beam, Vector_t< double, 3 > FDext[], const losses_t &losses=losses_t(), const double &azimuth=-1, const size_t npOutside=0, size_t particleContainerIndex=0)
Write statistical data.
StatWriter(const std::string &fname, bool restart)
constexpr double s2ns
Definition Units.h:44
std::string getGitRevision()
Definition Util.cpp:32