OPAL (Object Oriented Parallel Accelerator Library)  2024.1
OPAL
StatWriter.cpp
Go to the documentation of this file.
1 //
2 // Class StatWriter
3 // This class writes bunch statistics (*.stat).
4 //
5 // Copyright (c) 2019, Matthias Frey, Paul Scherrer Institut, Villigen PSI, Switzerland
6 // Christof Metzger-Kraus, Open Sourcerer
7 // All rights reserved
8 //
9 // This file is part of OPAL.
10 //
11 // OPAL is free software: you can redistribute it and/or modify
12 // it under the terms of the GNU General Public License as published by
13 // the Free Software Foundation, either version 3 of the License, or
14 // (at your option) any later version.
15 //
16 // You should have received a copy of the GNU General Public License
17 // along with OPAL. If not, see <https://www.gnu.org/licenses/>.
18 //
19 #include "StatWriter.h"
20 
23 #include "Utilities/Timer.h"
24 #include "Physics/Units.h"
25 
26 #include <sstream>
27 
28 StatWriter::StatWriter(const std::string& fname, bool restart)
29  : StatBaseWriter(fname, restart)
30 { }
31 
32 
33 void StatWriter::fillHeader(const losses_t &losses) {
34 
35  if (this->hasColumns()) {
36  return;
37  }
38 
39  columns_m.addColumn("t", "double", "ns", "Time");
40  columns_m.addColumn("s", "double", "m", "Path length");
41  columns_m.addColumn("numParticles", "long", "1", "Number of Macro Particles");
42  columns_m.addColumn("charge", "double", "1", "Bunch Charge");
43  columns_m.addColumn("energy", "double", "MeV", "Mean Bunch Energy");
44 
45  columns_m.addColumn("rms_x", "double", "m", "RMS Beamsize in x");
46  columns_m.addColumn("rms_y", "double", "m", "RMS Beamsize in y");
47  columns_m.addColumn("rms_s", "double", "m", "RMS Beamsize in s");
48 
49  columns_m.addColumn("rms_px", "double", "1", "RMS Normalized Momenta in x");
50  columns_m.addColumn("rms_py", "double", "1", "RMS Normalized Momenta in y");
51  columns_m.addColumn("rms_ps", "double", "1", "RMS Normalized Momenta in s");
52 
53  columns_m.addColumn("emit_x", "double", "m", "Normalized Emittance x");
54  columns_m.addColumn("emit_y", "double", "m", "Normalized Emittance y");
55  columns_m.addColumn("emit_s", "double", "m", "Normalized Emittance s");
56 
57  columns_m.addColumn("mean_x", "double", "m", "Mean Beam Position in x");
58  columns_m.addColumn("mean_y", "double", "m", "Mean Beam Position in y");
59  columns_m.addColumn("mean_s", "double", "m", "Mean Beam Position in s");
60 
61  columns_m.addColumn("ref_x", "double", "m", "x coordinate of reference particle in lab cs");
62  columns_m.addColumn("ref_y", "double", "m", "y coordinate of reference particle in lab cs");
63  columns_m.addColumn("ref_z", "double", "m", "z coordinate of reference particle in lab cs");
64 
65  columns_m.addColumn("ref_px", "double", "1", "x momentum of reference particle in lab cs");
66  columns_m.addColumn("ref_py", "double", "1", "y momentum of reference particle in lab cs");
67  columns_m.addColumn("ref_pz", "double", "1", "z momentum of reference particle in lab cs");
68 
69  columns_m.addColumn("max_x", "double", "m", "Max Beamsize in x");
70  columns_m.addColumn("max_y", "double", "m", "Max Beamsize in y");
71  columns_m.addColumn("max_s", "double", "m", "Max Beamsize in s");
72 
73  columns_m.addColumn("xpx", "double", "1", "Correlation xpx");
74  columns_m.addColumn("ypy", "double", "1", "Correlation ypy");
75  columns_m.addColumn("zpz", "double", "1", "Correlation zpz");
76 
77  columns_m.addColumn("Dx", "double", "m", "Dispersion in x");
78  columns_m.addColumn("DDx", "double", "1", "Derivative of dispersion in x");
79  columns_m.addColumn("Dy", "double", "m", "Dispersion in y");
80  columns_m.addColumn("DDy", "double", "1", "Derivative of dispersion in y");
81 
82  columns_m.addColumn("Bx_ref", "double", "T", "Bx-Field component of ref particle");
83  columns_m.addColumn("By_ref", "double", "T", "By-Field component of ref particle");
84  columns_m.addColumn("Bz_ref", "double", "T", "Bz-Field component of ref particle");
85 
86  columns_m.addColumn("Ex_ref", "double", "MV/m", "Ex-Field component of ref particle");
87  columns_m.addColumn("Ey_ref", "double", "MV/m", "Ey-Field component of ref particle");
88  columns_m.addColumn("Ez_ref", "double", "MV/m", "Ez-Field component of ref particle");
89 
90  columns_m.addColumn("dE", "double", "MeV", "energy spread of the beam");
91  columns_m.addColumn("dt", "double", "ns", "time step size");
92  columns_m.addColumn("partsOutside", "double", "1", "outside n*sigma of the beam");
93  columns_m.addColumn("DebyeLength", "double", "m", "Debye length in the boosted frame");
94  columns_m.addColumn("plasmaParameter", "double", "1", "Plasma parameter that gives no. of particles in a Debye sphere");
95  columns_m.addColumn("temperature", "double", "K", "Temperature of the beam");
96  columns_m.addColumn("rmsDensity", "double", "1", "RMS number density of the beam");
97 
99  columns_m.addColumn("68_Percentile_x", "double", "m",
100  "68.27 percentile (1 sigma of normal distribution) of x-component of position");
101  columns_m.addColumn("68_Percentile_y", "double", "m",
102  "68.27 percentile (1 sigma of normal distribution) of y-component of position");
103  columns_m.addColumn("68_Percentile_z", "double", "m",
104  "68.27 percentile (1 sigma of normal distribution) of z-component of position");
105 
106  columns_m.addColumn("95_Percentile_x", "double", "m",
107  "95.45 percentile (2 sigma of normal distribution) of x-component of position");
108  columns_m.addColumn("95_Percentile_y", "double", "m",
109  "95.45 percentile (2 sigma of normal distribution) of y-component of position");
110  columns_m.addColumn("95_Percentile_z", "double", "m",
111  "95.45 percentile (2 sigma of normal distribution) of z-component of position");
112 
113  columns_m.addColumn("99_Percentile_x", "double", "m",
114  "99.73 percentile (3 sigma of normal distribution) of x-component of position");
115  columns_m.addColumn("99_Percentile_y", "double", "m",
116  "99.73 percentile (3 sigma of normal distribution) of y-component of position");
117  columns_m.addColumn("99_Percentile_z", "double", "m",
118  "99.73 percentile (3 sigma of normal distribution) of z-component of position");
119 
120  columns_m.addColumn("99_99_Percentile_x", "double", "m",
121  "99.994 percentile (4 sigma of normal distribution) of x-component of position");
122  columns_m.addColumn("99_99_Percentile_y", "double", "m",
123  "99.994 percentile (4 sigma of normal distribution) of y-component of position");
124  columns_m.addColumn("99_99_Percentile_z", "double", "m",
125  "99.994 percentile (4 sigma of normal distribution) of z-component of position");
126 
127  columns_m.addColumn("normalizedEps68Percentile_x", "double", "m",
128  "x-component of normalized emittance at 68 percentile (1 sigma of normal distribution)");
129  columns_m.addColumn("normalizedEps68Percentile_y", "double", "m",
130  "y-component of normalized emittance at 68 percentile (1 sigma of normal distribution)");
131  columns_m.addColumn("normalizedEps68Percentile_z", "double", "m",
132  "z-component of normalized emittance at 68 percentile (1 sigma of normal distribution)");
133 
134  columns_m.addColumn("normalizedEps95Percentile_x", "double", "m",
135  "x-component of normalized emittance at 95 percentile (2 sigma of normal distribution)");
136  columns_m.addColumn("normalizedEps95Percentile_y", "double", "m",
137  "y-component of normalized emittance at 95 percentile (2 sigma of normal distribution)");
138  columns_m.addColumn("normalizedEps95Percentile_z", "double", "m",
139  "z-component of normalized emittance at 95 percentile (2 sigma of normal distribution)");
140 
141  columns_m.addColumn("normalizedEps99Percentile_x", "double", "m",
142  "x-component of normalized emittance at 99 percentile (3 sigma of normal distribution)");
143  columns_m.addColumn("normalizedEps99Percentile_y", "double", "m",
144  "y-component of normalized emittance at 99 percentile (3 sigma of normal distribution)");
145  columns_m.addColumn("normalizedEps99Percentile_z", "double", "m",
146  "z-component of normalized emittance at 99 percentile (3 sigma of normal distribution)");
147 
148  columns_m.addColumn("normalizedEps99_99Percentile_x", "double", "m",
149  "x-component of normalized emittance at 99.99 percentile (4 sigma of normal distribution)");
150  columns_m.addColumn("normalizedEps99_99Percentile_y", "double", "m",
151  "y-component of normalized emittance at 99.99 percentile (4 sigma of normal distribution)");
152  columns_m.addColumn("normalizedEps99_99Percentile_z", "double", "m",
153  "z-component of normalized emittance at 99.99 percentile (4 sigma of normal distribution)");
154  }
156  columns_m.addColumn("S11","double","m**2","element 1,1 of 6D beam matrix");
157  columns_m.addColumn("S12","double","m","element 1,2 of 6D beam matrix");
158  columns_m.addColumn("S13","double","m**2","element 1,3 of 6D beam matrix");
159  columns_m.addColumn("S14","double","m","element 1,4 of 6D beam matrix");
160  columns_m.addColumn("S15","double","m**2","element 1,5 of 6D beam matrix");
161  columns_m.addColumn("S16","double","m","element 1,6 of 6D beam matrix");
162  columns_m.addColumn("S22","double","1","element 2,2 of 6D beam matrix");
163  columns_m.addColumn("S23","double","m","element 2,3 of 6D beam matrix");
164  columns_m.addColumn("S24","double","1","element 2,4 of 6D beam matrix");
165  columns_m.addColumn("S25","double","m","element 2,5 of 6D beam matrix");
166  columns_m.addColumn("S26","double","1","element 2,6 of 6D beam matrix");
167  columns_m.addColumn("S33","double","m**2","element 3,3 of 6D beam matrix");
168  columns_m.addColumn("S34","double","m","element 3,4 of 6D beam matrix");
169  columns_m.addColumn("S35","double","m**2","element 3,5 of 6D beam matrix");
170  columns_m.addColumn("S36","double","m","element 3,6 of 6D beam matrix");
171  columns_m.addColumn("S44","double","1","element 4,4 of 6D beam matrix");
172  columns_m.addColumn("S45","double","m","element 4,5 of 6D beam matrix");
173  columns_m.addColumn("S46","double","1","element 4,6 of 6D beam matrix");
174  columns_m.addColumn("S55","double","m**2","element 5,5 of 6D beam matrix");
175  columns_m.addColumn("S56","double","m","element 5,6 of 6D beam matrix");
176  columns_m.addColumn("S66","double","1","element 6,6 of 6D beam matrix");
177  }
178  if (OpalData::getInstance()->isInOPALCyclMode() &&
179  Ippl::getNodes() == 1) {
180  columns_m.addColumn("R0_x", "double", "m", "R0 Particle position in x");
181  columns_m.addColumn("R0_y", "double", "m", "R0 Particle position in y");
182  columns_m.addColumn("R0_s", "double", "m", "R0 Particle position in z");
183 
184  columns_m.addColumn("P0_x", "double", "1", "R0 Particle momentum in x");
185  columns_m.addColumn("P0_y", "double", "1", "R0 Particle momentum in y");
186  columns_m.addColumn("P0_s", "double", "1", "R0 Particle momentum in z");
187  }
188 
189  if (OpalData::getInstance()->isInOPALCyclMode()) {
190  columns_m.addColumn("halo_x", "double", "1", "Halo in x");
191  columns_m.addColumn("halo_y", "double", "1", "Halo in y");
192  columns_m.addColumn("halo_z", "double", "1", "Halo in z");
193 
194  columns_m.addColumn("azimuth", "double", "deg",
195  "Azimuth in global coordinates");
196  }
197 
198  for (size_t i = 0; i < losses.size(); ++ i) {
199  columns_m.addColumn(losses[i].first, "long", "1",
200  "Number of lost particles in element");
201  }
202 
203  if ( mode_m == std::ios::app )
204  return;
205 
206  OPALTimer::Timer simtimer;
207  std::string dateStr(simtimer.date());
208  std::string timeStr(simtimer.time());
209 
210  std::stringstream ss;
211  ss << "Statistics data '"
213  << "' " << dateStr << " " << timeStr;
214 
215  this->addDescription(ss.str(), "stat parameters");
216 
217  this->addDefaultParameters();
218 
219 
220  this->addInfo("ascii", 1);
221 }
222 
223 
225  const losses_t &losses, const double& azimuth,
226  const size_t npOutside)
227 {
228  double Ekin = beam->get_meanKineticEnergy();
229 
230  double pathLength = beam->get_sPos();
231 
234 
235  double Q = beam->getCharge();
236 
237  if (Ippl::myNode() != 0) {
238  return;
239  }
240 
241  fillHeader(losses);
242 
243  this->open();
244 
245  this->writeHeader();
246 
247  columns_m.addColumnValue("t", beam->getT() * Units::s2ns); // 1
248  columns_m.addColumnValue("s", pathLength); // 2
249  columns_m.addColumnValue("numParticles", beam->getTotalNum()); // 3
250  columns_m.addColumnValue("charge", Q); // 4
251  columns_m.addColumnValue("energy", Ekin); // 5
252 
253  columns_m.addColumnValue("rms_x", beam->get_rrms()(0)); // 6
254  columns_m.addColumnValue("rms_y", beam->get_rrms()(1)); // 7
255  columns_m.addColumnValue("rms_s", beam->get_rrms()(2)); // 8
256 
257  columns_m.addColumnValue("rms_px", beam->get_prms()(0)); // 9
258  columns_m.addColumnValue("rms_py", beam->get_prms()(1)); // 10
259  columns_m.addColumnValue("rms_ps", beam->get_prms()(2)); // 11
260 
261  columns_m.addColumnValue("emit_x", beam->get_norm_emit()(0)); // 12
262  columns_m.addColumnValue("emit_y", beam->get_norm_emit()(1)); // 13
263  columns_m.addColumnValue("emit_s", beam->get_norm_emit()(2)); // 14
264 
265  columns_m.addColumnValue("mean_x", beam->get_rmean()(0) ); // 15
266  columns_m.addColumnValue("mean_y", beam->get_rmean()(1) ); // 16
267  columns_m.addColumnValue("mean_s", beam->get_rmean()(2) ); // 17
268 
269  columns_m.addColumnValue("ref_x", beam->RefPartR_m(0)); // 18
270  columns_m.addColumnValue("ref_y", beam->RefPartR_m(1)); // 19
271  columns_m.addColumnValue("ref_z", beam->RefPartR_m(2)); // 20
272 
273  columns_m.addColumnValue("ref_px", beam->RefPartP_m(0)); // 21
274  columns_m.addColumnValue("ref_py", beam->RefPartP_m(1)); // 22
275  columns_m.addColumnValue("ref_pz", beam->RefPartP_m(2)); // 23
276 
277  columns_m.addColumnValue("max_x", beam->get_maxExtent()(0)); // 24
278  columns_m.addColumnValue("max_y", beam->get_maxExtent()(1)); // 25
279  columns_m.addColumnValue("max_s", beam->get_maxExtent()(2)); // 26
280 
281  // Write out Courant Snyder parameters.
282  columns_m.addColumnValue("xpx", beam->get_rprms()(0)); // 27
283  columns_m.addColumnValue("ypy", beam->get_rprms()(1)); // 28
284  columns_m.addColumnValue("zpz", beam->get_rprms()(2)); // 29
285 
286  // Write out dispersion.
287  columns_m.addColumnValue("Dx", beam->get_Dx()); // 30
288  columns_m.addColumnValue("DDx", beam->get_DDx()); // 31
289  columns_m.addColumnValue("Dy", beam->get_Dy()); // 32
290  columns_m.addColumnValue("DDy", beam->get_DDy()); // 33
291 
292  // Write head/reference particle/tail field information.
293  columns_m.addColumnValue("Bx_ref", FDext[0](0)); // 34 B-ref x
294  columns_m.addColumnValue("By_ref", FDext[0](1)); // 35 B-ref y
295  columns_m.addColumnValue("Bz_ref", FDext[0](2)); // 36 B-ref z
296 
297  columns_m.addColumnValue("Ex_ref", FDext[1](0)); // 37 E-ref x
298  columns_m.addColumnValue("Ey_ref", FDext[1](1)); // 38 E-ref y
299  columns_m.addColumnValue("Ez_ref", FDext[1](2)); // 39 E-ref z
300 
301  columns_m.addColumnValue("dE", beam->getdE()); // 40 dE energy spread
302  columns_m.addColumnValue("dt", beam->getdT() * Units::s2ns); // 41 dt time step size
303  columns_m.addColumnValue("partsOutside", npOutside); // 42 number of particles outside n*sigma
304 
305  columns_m.addColumnValue("DebyeLength", beam->get_debyeLength()); // 43 Debye length in the boosted frame
306  columns_m.addColumnValue("plasmaParameter", beam->get_plasmaParameter()); // 43 plasma parameter
307  columns_m.addColumnValue("temperature", beam->get_temperature()); // 44 Temperature
308  columns_m.addColumnValue("rmsDensity", beam->get_rmsDensity()); // 45 RMS number density
309 
311  columns_m.addColumnValue("68_Percentile_x", beam->get_68Percentile()[0]);
312  columns_m.addColumnValue("68_Percentile_y", beam->get_68Percentile()[1]);
313  columns_m.addColumnValue("68_Percentile_z", beam->get_68Percentile()[2]);
314 
315  columns_m.addColumnValue("95_Percentile_x", beam->get_95Percentile()[0]);
316  columns_m.addColumnValue("95_Percentile_y", beam->get_95Percentile()[1]);
317  columns_m.addColumnValue("95_Percentile_z", beam->get_95Percentile()[2]);
318 
319  columns_m.addColumnValue("99_Percentile_x", beam->get_99Percentile()[0]);
320  columns_m.addColumnValue("99_Percentile_y", beam->get_99Percentile()[1]);
321  columns_m.addColumnValue("99_Percentile_z", beam->get_99Percentile()[2]);
322 
323  columns_m.addColumnValue("99_99_Percentile_x", beam->get_99_99Percentile()[0]);
324  columns_m.addColumnValue("99_99_Percentile_y", beam->get_99_99Percentile()[1]);
325  columns_m.addColumnValue("99_99_Percentile_z", beam->get_99_99Percentile()[2]);
326 
327  columns_m.addColumnValue("normalizedEps68Percentile_x", beam->get_normalizedEps_68Percentile()[0]);
328  columns_m.addColumnValue("normalizedEps68Percentile_y", beam->get_normalizedEps_68Percentile()[1]);
329  columns_m.addColumnValue("normalizedEps68Percentile_z", beam->get_normalizedEps_68Percentile()[2]);
330 
331  columns_m.addColumnValue("normalizedEps95Percentile_x", beam->get_normalizedEps_95Percentile()[0]);
332  columns_m.addColumnValue("normalizedEps95Percentile_y", beam->get_normalizedEps_95Percentile()[1]);
333  columns_m.addColumnValue("normalizedEps95Percentile_z", beam->get_normalizedEps_95Percentile()[2]);
334 
335  columns_m.addColumnValue("normalizedEps99Percentile_x", beam->get_normalizedEps_99Percentile()[0]);
336  columns_m.addColumnValue("normalizedEps99Percentile_y", beam->get_normalizedEps_99Percentile()[1]);
337  columns_m.addColumnValue("normalizedEps99Percentile_z", beam->get_normalizedEps_99Percentile()[2]);
338 
339  columns_m.addColumnValue("normalizedEps99_99Percentile_x", beam->get_normalizedEps_99_99Percentile()[0]);
340  columns_m.addColumnValue("normalizedEps99_99Percentile_y", beam->get_normalizedEps_99_99Percentile()[1]);
341  columns_m.addColumnValue("normalizedEps99_99Percentile_z", beam->get_normalizedEps_99_99Percentile()[2]);
342  }
344  columns_m.addColumnValue("S11", beam->getSigmaMatrix()[0][0]);
345  columns_m.addColumnValue("S12", beam->getSigmaMatrix()[0][1]);
346  columns_m.addColumnValue("S13", beam->getSigmaMatrix()[0][2]);
347  columns_m.addColumnValue("S14", beam->getSigmaMatrix()[0][3]);
348  columns_m.addColumnValue("S15", beam->getSigmaMatrix()[0][4]);
349  columns_m.addColumnValue("S16", beam->getSigmaMatrix()[0][5]);
350  columns_m.addColumnValue("S22", beam->getSigmaMatrix()[1][1]);
351  columns_m.addColumnValue("S23", beam->getSigmaMatrix()[1][2]);
352  columns_m.addColumnValue("S24", beam->getSigmaMatrix()[1][3]);
353  columns_m.addColumnValue("S25", beam->getSigmaMatrix()[1][4]);
354  columns_m.addColumnValue("S26", beam->getSigmaMatrix()[1][5]);
355  columns_m.addColumnValue("S33", beam->getSigmaMatrix()[2][2]);
356  columns_m.addColumnValue("S34", beam->getSigmaMatrix()[2][3]);
357  columns_m.addColumnValue("S35", beam->getSigmaMatrix()[2][4]);
358  columns_m.addColumnValue("S36", beam->getSigmaMatrix()[2][5]);
359  columns_m.addColumnValue("S44", beam->getSigmaMatrix()[3][3]);
360  columns_m.addColumnValue("S45", beam->getSigmaMatrix()[3][4]);
361  columns_m.addColumnValue("S46", beam->getSigmaMatrix()[3][5]);
362  columns_m.addColumnValue("S55", beam->getSigmaMatrix()[4][4]);
363  columns_m.addColumnValue("S56", beam->getSigmaMatrix()[4][5]);
364  columns_m.addColumnValue("S66", beam->getSigmaMatrix()[5][5]);
365  }
367  if (Ippl::getNodes() == 1) {
368  if (beam->getLocalNum() > 0) {
369  columns_m.addColumnValue("R0_x", beam->R[0](0));
370  columns_m.addColumnValue("R0_y", beam->R[0](1));
371  columns_m.addColumnValue("R0_s", beam->R[0](2));
372  columns_m.addColumnValue("P0_x", beam->P[0](0));
373  columns_m.addColumnValue("P0_y", beam->P[0](1));
374  columns_m.addColumnValue("P0_s", beam->P[0](2));
375  } else {
376  columns_m.addColumnValue("R0_x", 0.0);
377  columns_m.addColumnValue("R0_y", 0.0);
378  columns_m.addColumnValue("R0_s", 0.0);
379  columns_m.addColumnValue("P0_x", 0.0);
380  columns_m.addColumnValue("P0_y", 0.0);
381  columns_m.addColumnValue("P0_s", 0.0);
382  }
383  }
384  Vector_t halo = beam->get_halo();
385  columns_m.addColumnValue("halo_x", halo(0));
386  columns_m.addColumnValue("halo_y", halo(1));
387  columns_m.addColumnValue("halo_z", halo(2));
388 
389  columns_m.addColumnValue("azimuth", azimuth);
390  }
391 
392  for(size_t i = 0; i < losses.size(); ++ i) {
393  long unsigned int loss = losses[i].second;
394  columns_m.addColumnValue(losses[i].first, loss);
395  }
396 
397  this->writeRow();
398 
399  this->close();
400 }
bool dumpBeamMatrix
Definition: Options.cpp:113
static OpalData * getInstance()
Definition: OpalData.cpp:196
void fillHeader(const losses_t &losses=losses_t())
Definition: StatWriter.cpp:33
double get_rmsDensity() const
void open()
Definition: SDDSWriter.cpp:132
Vector_t get_normalizedEps_99Percentile() const
bool hasColumns() const
Definition: SDDSWriter.h:200
double get_Dx() const
void addDescription(const std::string &text, const std::string &content)
Definition: SDDSWriter.h:159
Vector_t get_rprms() const
ParticleAttrib< Vector_t > P
std::vector< std::pair< std::string, unsigned int > > losses_t
Definition: StatWriter.h:27
static int myNode()
Definition: IpplInfo.cpp:691
SDDSColumnSet columns_m
Definition: SDDSWriter.h:122
double get_DDy() const
std::string getInputFn()
get opals input filename
Definition: OpalData.cpp:670
StatWriter(const std::string &fname, bool restart)
Definition: StatWriter.cpp:28
double get_sPos() const
void write(const PartBunchBase< double, 3 > *beam, Vector_t FDext[], const losses_t &losses=losses_t(), const double &azimuth=-1, const size_t npOutside=0)
Write statistical data.
Definition: StatWriter.cpp:224
double getdE() const
Vector_t get_rrms() const
double get_plasmaParameter() const
double getdT() const
Vector_t get_prms() const
void addDefaultParameters()
Definition: SDDSWriter.cpp:209
static int getNodes()
Definition: IpplInfo.cpp:670
double getT() const
void addInfo(const std::string &mode, const size_t &no_row_counts)
Definition: SDDSWriter.h:178
void writeRow()
Definition: SDDSWriter.h:185
size_t getTotalNum() const
size_t getLocalNum() const
void writeHeader()
Write SDDS header.
Definition: SDDSWriter.cpp:149
constexpr double s2ns
Definition: Units.h:44
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)
FMatrix< double, 2 *Dim, 2 *Dim > getSigmaMatrix() const
double get_DDx() const
double get_Dy() const
Vector_t get_rmean() const
ParticlePos_t & R
Vector_t get_norm_emit() const
std::string time() const
Return time.
Definition: Timer.cpp:42
bool computePercentiles
Definition: Options.cpp:111
std::string date() const
Return date.
Definition: Timer.cpp:35
bool isInOPALCyclMode()
Definition: OpalData.cpp:272
Vector_t get_99Percentile() const
Vector_t RefPartR_m
double getCharge() const
get the total charge per simulation particle
Vector_t get_99_99Percentile() const
double get_meanKineticEnergy() const
Vector_t get_normalizedEps_68Percentile() const
void close()
Definition: SDDSWriter.cpp:142
Vector_t get_maxExtent() const
Vector_t RefPartP_m
void addColumnValue(const std::string &name, const T &val)
Definition: SDDSColumnSet.h:63
Vector_t get_halo() const
double get_temperature() const
Vector_t get_95Percentile() const
Vector_t get_normalizedEps_95Percentile() const
double get_debyeLength() const
Vector_t get_normalizedEps_99_99Percentile() const
std::ios_base::openmode mode_m
First write to the statistics output file.
Definition: SDDSWriter.h:120
Vector_t get_68Percentile() const