Skip to content

ImpactZOutput

impact.z.ImpactZOutput

Bases: Mapping, BaseModel

IMPACT-Z command output.

Attributes:

Name Type Description
run RunInfo

Execution information - e.g., how long did it take and what was the output from Impact-Z.

alias dict[str, str]

Dictionary of aliased data keys.

Functions

impact.z.ImpactZOutput.archive
archive(h5)

Dump outputs into the given HDF5 group.

Parameters:

Name Type Description Default
h5 Group

The HDF5 file in which to write the information.

required
Source code in impact/z/output.py
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
def archive(self, h5: h5py.Group) -> None:
    """
    Dump outputs into the given HDF5 group.

    Parameters
    ----------
    h5 : h5py.Group
        The HDF5 file in which to write the information.
    """
    _archive.store_in_hdf5_file(h5, self)
impact.z.ImpactZOutput.from_archive classmethod
from_archive(h5)

Loads output from the given HDF5 group.

Parameters:

Name Type Description Default
h5 str or File

The key to use when restoring the data.

required
Source code in impact/z/output.py
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
@classmethod
def from_archive(cls, h5: h5py.Group) -> ImpactZOutput:
    """
    Loads output from the given HDF5 group.

    Parameters
    ----------
    h5 : str or h5py.File
        The key to use when restoring the data.
    """
    loaded = _archive.restore_from_hdf5_file(h5)
    if not isinstance(loaded, ImpactZOutput):
        raise ValueError(
            f"Loaded {loaded.__class__.__name__} instead of a "
            f"ImpactZOutput instance.  Was the HDF group correct?"
        )
    return loaded
impact.z.ImpactZOutput.from_input_settings classmethod
from_input_settings(input, workdir)

Load ImpactZ output based on the configured input settings.

Returns:

Type Description
ImpactZOutput
Source code in impact/z/output.py
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
@classmethod
def from_input_settings(
    cls,
    input: ImpactZInput,
    workdir: pathlib.Path,
) -> ImpactZOutput:
    """
    Load ImpactZ output based on the configured input settings.

    Returns
    -------
    ImpactZOutput
    """

    species = input.reference_species
    stats = OutputStats.from_stats_files(
        workdir,
        reference_frequency=input.reference_frequency,
        diagnostic_type=input.diagnostic_type,
        reference_particle_mass=input.reference_particle_mass,
    )

    units = stats.units.copy()
    particles_raw = {}
    particles = {}
    slices = {}

    z_start = 0.0
    z_end = 0.0

    for ele in input.lattice:
        z_start = z_end
        z_end += ele.length
        z_end_idx = np.argmin(np.abs(stats.z - z_start))

        if isinstance(ele, WriteSliceInfo):
            key = _get_dict_key(slices, ele.file_id, ele.name)
            slices[ele.file_id] = ImpactZSlices.from_file(
                workdir / f"fort.{ele.file_id}"
            )
        elif ele.class_information().has_output_file and isinstance(
            ele, HasOutputFile
        ):
            raw = ImpactZParticles.from_file(workdir / f"fort.{ele.file_id}")

            key = _get_dict_key(particles_raw, ele.file_id, ele.name)
            particles_raw[key] = raw
            phase_ref = stats.phase_ref[z_end_idx]
            kinetic_energy = stats.kinetic_energy_ref[z_end_idx]
            particles[key] = raw.to_particle_group(
                reference_kinetic_energy=kinetic_energy,
                reference_frequency=input.reference_frequency,
                phase_reference=phase_ref,
            )

    for key, fld in OutputStats.model_fields.items():
        unit = _units_from_metadata(fld.metadata)
        if unit:
            units[key] = unit

    return cls(
        stats=stats,
        key_to_unit=units,
        particles=particles,
        particles_raw=particles_raw,
        reference_frequency=input.reference_frequency,
        slices=slices,
        reference_species=species,
    )
impact.z.ImpactZOutput.plot
plot(y=('sigma_x', 'sigma_y'), x='z', *, y2=(), input=None, xlim=None, ylim=None, ylim2=None, nice=True, tex=True, include_layout=True, include_labels=True, include_markers=True, include_particles=True, include_legend=True, return_figure=False, ax=None, **kwargs)
Source code in impact/z/output.py
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
def plot(
    self,
    y: str | Sequence[str] = ("sigma_x", "sigma_y"),
    x: str = "z",
    *,
    y2: str | Sequence[str] = (),
    input: ImpactZInput | None = None,
    xlim: tuple[float, float] | None = None,
    ylim: tuple[float, float] | None = None,
    ylim2: tuple[float, float] | None = None,
    nice: bool = True,
    tex: bool = True,
    include_layout: bool = True,
    include_labels: bool = True,
    include_markers: bool = True,
    include_particles: bool = True,
    include_legend: bool = True,
    return_figure: bool = False,
    ax: matplotlib.axes.Axes | None = None,
    **kwargs,
):
    """ """

    if "ykeys2" in kwargs:
        y2 = kwargs.pop("ykeys2")

    if input is None:
        # Should we warn?
        include_layout = False

    return plot_stats_with_layout(
        self,
        ykeys=y,
        ykeys2=y2,
        xkey=x,
        xlim=xlim,
        ylim=ylim,
        ylim2=ylim2,
        input=input,
        nice=nice,
        tex=tex,
        include_layout=include_layout,
        include_labels=include_labels,
        # include_field=include_field,
        include_markers=include_markers,
        include_particles=include_particles,
        include_legend=include_legend,
        return_figure=return_figure,
        ax=ax,
        **kwargs,
    )
impact.z.ImpactZOutput.stat
stat(key)

Statistics array from .stats.

Source code in impact/z/output.py
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
def stat(self, key: str):
    """
    Statistics array from `.stats`.
    """
    # Allow flipping covariance keys
    if key.startswith("cov_") and key not in self:
        k1, k2 = key[4:].split("__")
        key = f"cov_{k2}__{k1}"

    if key not in self:
        raise ValueError(f"{key} is not available in the output data")

    return self[self.alias.get(key, key)]

impact.z.output.RunInfo

Bases: BaseModel

Impact-Z run information.

Attributes:

Name Type Description
error bool

True if an error occurred during the Impact-Z run.

error_reason str or None

Error explanation, if error is set.

run_script str

The command-line arguments used to run Impact-Z

output_log str

Impact-Z output log

start_time float

Start time of the process

end_time float

End time of the process

run_time float

Wall clock run time of the process

Attributes

impact.z.output.RunInfo.success property
success

True if the run was successful.

impact.z.output.OutputStats

Bases: BaseModel

Output statistics.

Attributes:

Name Type Description
beta_ref ndarray

Beta of the reference particle.

charge_state_n_particle ndarray

The number of particles for each charge state.

gamma_ref ndarray

Reference particle gamma.

kinetic_energy_ref ndarray

Reference particle kinetic energy. (eV)

loadbalance_max_n_particle ndarray

Maximum number of particles on a processing element (PE).

loadbalance_min_n_particle ndarray

Minimum number of particles on a processing element (PE).

max_amplitude_energy_dev ndarray

Maximum energy deviation (eV).

max_amplitude_gammabeta_x ndarray

Maximum Px (in radians).

max_amplitude_gammabeta_y ndarray

Maximum Py (in radians).

max_amplitude_phase ndarray

Maximum phase (in degrees).

max_amplitude_x ndarray

Maximum X (in meters).

max_amplitude_y ndarray

Maximum Y (in meters).

max_r ndarray

Maximum radius (Rmax) in meters, measured from the axis of the pipe.

mean_phase_deg ndarray

Mean phase (degrees)

mean_px_over_p0 ndarray

Mean \(px / p0\) (unitless).

mean_py_over_p0 ndarray

Mean \(py / p0\) (unitless).

mean_x ndarray

Centroid location in the x-direction (meters).

mean_y ndarray

Centroid location in the y-direction (meters).

moment3_energy ndarray

Third-order central moment for energy deviation (eV).

moment3_phase ndarray

Third-order central moment for phase (degree).

moment3_px_over_p0 ndarray

Third-order central moment for Px (rad).

moment3_py_over_p0 ndarray

Third-order central moment for Py (rad).

moment3_x ndarray

Third-order central moment for x (meters).

moment3_y ndarray

Third-order central moment for y (meters).

moment4_energy ndarray

Fourth-order central moment for energy deviation (eV).

moment4_phase ndarray

Fourth-order central moment for phase (degree).

moment4_px_over_p0 ndarray

Fourth-order central moment for Px over p0 (unitless).

moment4_py_over_p0 ndarray

Fourth-order central moment for Py over p0 (unitless).

moment4_x ndarray

Fourth-order central moment for x (meters).

moment4_y ndarray

Fourth-order central moment for y (meters).

n_particle ndarray

Total number of particles in the bunch.

neg_mean_rel_energy ndarray

Negative delta mean energy (eV).

norm_emit_x ndarray

Normalized RMS emittance in x-direction (m-rad).

norm_emit_y ndarray

Normalized RMS emittance in y-direction (m-rad).

norm_emit_z ndarray

Normalized RMS emittance in z-direction (degree-eV).

phase_ref ndarray

Absolute phase in radians.

sigma_energy ndarray

Sigma energy (eV).

sigma_phase_deg ndarray

Sigma phase (degrees).

sigma_px_over_p0 ndarray

Sigma \(px / p0\) (unitless).

sigma_py_over_p0 ndarray

Sigma \(py / p0\) (unitless).

sigma_x ndarray

RMS size in the x-direction (meters).

sigma_y ndarray

RMS size in the y-direction (meters).

twiss_alpha_x ndarray

Twiss parameter alpha for x-direction (unitless).

twiss_alpha_y ndarray

Twiss parameter alpha for y-direction (unitless).

twiss_alpha_z ndarray

Twiss parameter alpha for z-direction (unitless).

z ndarray

Z position (meters)

max_abs_x ndarray

Maximum horizontal displacement from the beam axis: \(\max(|x|)\) (meters)

max_abs_px_over_p0 ndarray

Maximum \(x\)-plane transverse momentum \(\max(|p_x/p_0|)\) (unitless)

max_abs_y ndarray

Maximum vertical displacement from the beam axis \(\max(|y|)\) (meters)

max_abs_py_over_p0 ndarray

Maximum \(y\)-plane transverse momentum \(\max(|p_y/p_0|)\) (unitless)

max_phase ndarray

Maximum deviation in phase (degrees)

max_energy_dev ndarray

Maximum deviation in particle energy (eV)

mean_r ndarray

Mean radius (meters)

sigma_r ndarray

RMS radius (meters)

mean_r_90percent ndarray

90 percent mean radius (meters)

mean_r_95percent ndarray

95 percent mean radius (meters)

mean_r_99percent ndarray

99 percent mean radius (meters)

max_r_rel ndarray

Maximum radius (meters)

max_abs_gammabeta_x ndarray

Maximum \(x\)-plane transverse momentum \(\max(|\gamma\beta_x|)\) (dimensionless)

max_abs_gammabeta_y ndarray

Maximum \(x\)-plane transverse momentum \(\max(|\gamma\beta_y|)\) (dimensionless)

max_gamma_rel ndarray

Maximum deviation in relativistic gamma \(\max(|\gamma - \gamma_0)|)\) (dimensionless)

norm_emit_x_90percent ndarray

90% normalied RMS emittance (meter-rad)

norm_emit_x_95percent ndarray

90% normalied RMS emittance (meter-rad)

norm_emit_x_99percent ndarray

90% normalied RMS emittance (meter-rad)

norm_emit_y_90percent ndarray

90% normalied RMS emittance (meter-rad)

norm_emit_y_95percent ndarray

90% normalied RMS emittance (meter-rad)

norm_emit_y_99percent ndarray

90% normalied RMS emittance (meter-rad)

norm_emit_z_90percent ndarray

90% normalied RMS emittance (meter-rad)

norm_emit_z_95percent ndarray

90% normalied RMS emittance (meter-rad)

norm_emit_z_99percent ndarray

90% normalied RMS emittance (meter-rad)

energy_ref ndarray

Energy reference (eV) (computed)

mean_energy ndarray

Mean energy (eV) (computed)

mean_gamma ndarray

Mean gamma (computed)

mean_px ndarray

Mean px (eV/c) (computed)

mean_py ndarray

Mean py (eV/c) (computed)

mean_t ndarray

Mean time (s) (computed)

mean_t_rel ndarray

Mean time relative (s) (computed)

p0c ndarray

Momentum reference (eV) (computed)

sigma_px ndarray

Sigma px (eV/c) (computed)

sigma_py ndarray

Sigma py (eV/c) (computed)

sigma_t ndarray

RMS size in time (rad) (computed)

t_ref ndarray

Reference time (sec) (computed)

twiss_beta_x ndarray

Twiss beta x (m) (computed)

twiss_beta_y ndarray

Twiss beta y (m) (computed)

units dict[str, pmd_unit]

Mapping of attribute name to pmd_unit.

extra dict[str, ndarray]

Additional Impact-Z output data. This is a future-proofing mechanism in case Impact-Z changes and LUME-Impact is not yet ready for it.

impact.z.output.ImpactZSlices

Bases: BaseModel

A class to represent the impact Z slices.

Attributes:

Name Type Description
bunch_length NDArray

Bunch length coordinate (m).

particles_per_slice NDArray

Number of particles per slice.

current_per_slice NDArray

Current (A) per slice.

normalized_emittance_x NDArray

X normalized emittance (m-rad) per slice.

normalized_emittance_y NDArray

Y normalized emittance (m-rad) per slice.

dE_E NDArray

dE/E.

uncorrelated_energy_spread NDArray

Uncorrelated energy spread (eV) per slice.

mean_x NDArray

(m) of each slice.

mean_y NDArray

(m) of each slice.

mismatch_factor_x NDArray

X mismatch factor.

mismatch_factor_y NDArray

Y mismatch factor.

Functions

impact.z.output.ImpactZSlices.from_contents classmethod
from_contents(contents, filename=None)

Load main input from its file contents.

Parameters:

Name Type Description Default
contents str

The contents of the main input file.

required
filename AnyPath or None

The filename, if known.

None

Returns:

Type Description
ImpactZSlices
Source code in impact/z/output.py
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
@classmethod
def from_contents(
    cls, contents: str, filename: AnyPath | None = None
) -> ImpactZSlices:
    """
    Load main input from its file contents.

    Parameters
    ----------
    contents : str
        The contents of the main input file.
    filename : AnyPath or None, optional
        The filename, if known.

    Returns
    -------
    ImpactZSlices
    """

    contents = parsers.fix_line(contents)
    fields = list(cls.model_fields)[:9]
    dtype = np.dtype(
        {
            "names": fields,
            "formats": [np.float64] * 9,
        }
    )

    lines = contents.splitlines()
    if not lines:
        return ImpactZSlices(
            filename=pathlib.Path(filename) if filename else None,
        )

    arrays = np.loadtxt(
        lines,
        dtype=dtype,
        usecols=range(len(fields)),
        unpack=True,
    )

    data = {field: arr for field, arr in zip(fields, arrays)}
    return ImpactZSlices(
        **data,
        filename=pathlib.Path(filename) if filename else None,
    )
impact.z.output.ImpactZSlices.from_file classmethod
from_file(filename)

Load a main input file from disk.

Parameters:

Name Type Description Default
filename AnyPath

The filename to load.

required

Returns:

Type Description
ImpactZSlices
Source code in impact/z/output.py
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
@classmethod
def from_file(cls, filename: AnyPath) -> ImpactZSlices:
    """
    Load a main input file from disk.

    Parameters
    ----------
    filename : AnyPath
        The filename to load.

    Returns
    -------
    ImpactZSlices
    """
    with open(filename) as fp:
        contents = fp.read()
    return cls.from_contents(contents, filename=filename)

Internal / Raw file data classes

impact.z.output.FortranOutputFileData

Bases: SequenceBaseModel

Base class for representing a single row of file data from fort.{file_id}.

Subclasses of this are used to: 1. Match output file IDs with the type of data they contain 2. Give names to each column 3. Track units (by way of annotated attributes)

Raw file data classes

impact.z.output.ReferenceParticles

Bases: FortranOutputFileData

Reference particle information from an output file.

Attributes:

Name Type Description
z float

Distance in meters (1st column).

phase_ref float

Absolute phase in radians (2nd column).

gamma_ref float

Reference particle gamma (3rd column).

kinetic_energy_ref float

Reference particle kinetic energy in MeV (4th column). LUME-ImpactZ converts this automatically to eV.

beta_ref float

Beta (5th column).

max_r float

Rmax in meters, measured from the axis of pipe (6th column).

impact.z.output.RmsX

Bases: FortranOutputFileData

RMS size information in X.

Attributes:

Name Type Description
z float

Mean z distance (m)

mean_x float

Centroid location (m)

sigma_x float

RMS size (m)

mean_px_over_p0 float

Mean \(px / p0\) (unitless)

sigma_px_over_p0 float

Sigma \(px / p0\) (unitless)

twiss_alpha_x float

Twiss parameter, alpha

norm_emit_x float

normalized RMS emittance [m-rad]

norm_emit_x_90percent float

90% normalized RMS horizontal emittance (meter-rad) Only available in diagnostic_types=extended mode.

norm_emit_x_95percent float

95% normalized RMS horizontal emittance (meter-rad) Only available in diagnostic_types=extended mode.

norm_emit_x_99percent float

99% normalized RMS horizontal emittance (meter-rad) Only available in diagnostic_types=extended mode.

impact.z.output.RmsY

Bases: FortranOutputFileData

RMS size information in Y.

Attributes:

Name Type Description
z float

z distance (m)

mean_y float

centroid location (m)

sigma_y float

RMS size (m)

mean_py_over_p0 float

Mean \(py / p0\) [unitless]

sigma_py_over_p0 float

\(py / p0\) [unitless]

twiss_alpha_y float

Twiss parameter, alpha

norm_emit_y float

normalized RMS emittance [m-rad]

norm_emit_y_90percent float

90% normalized RMS vertical emittance (meter-rad) Only available in diagnostic_type=extended mode.

norm_emit_y_95percent float

95% normalized RMS vertical emittance (meter-rad) Only available in diagnostic_type=extended mode.

norm_emit_y_99percent float

99% normalized RMS vertical emittance (meter-rad) Only available in diagnostic_type=extended mode.

impact.z.output.RmsZ

Bases: FortranOutputFileData

RMS size information in Z.

Attributes:

Name Type Description
z float

z distance (m)

mean_phase_deg float

Mean phase (degrees)

sigma_phase_deg float

RMS phase in degrees.

neg_mean_rel_energy float

Negative delta mean energy, used to convert to mean energy [eV] where neg_mean_rel_energy = (kinetic_energy_ref - mean_energy) + reference_particle_mass In the file, this is stored as MeV and LUME-Impact converts to eV automatically.

sigma_energy float

RMS momentum [eV] In the file, this is stored as MeV and LUME-Impact converts to eV automatically.

twiss_alpha_z float

Twiss parameter, alpha

norm_emit_z float

normalized RMS emittance [degree-MeV]

norm_emit_z_90percent float

90% normalized RMS longitudinal emittance (degree-MeV) Only available in diagnostic_type=extended mode.

norm_emit_z_95percent float

95% normalized RMS longitudinal emittance (degree-MeV) Only available in diagnostic_type=extended mode.

norm_emit_z_99percent float

99% normalized RMS longitudinal emittance (degree-MeV) Only available in diagnostic_type=extended mode.

impact.z.output.MaxAmplitudeStandard

Bases: FortranOutputFileData

File fort.27: maximum amplitude information (standard)

Attributes:

Name Type Description
z float

Longitudinal position along the beamline (meters)

max_abs_x float

Maximum horizontal displacement from the beam axis: \(\max(|x|)\) (meters)

max_abs_px_over_p0 float

Maximum \(x\)-plane transverse momentum \(\max(|p_x/p_0|)\) (unitless)

max_abs_y float

Maximum vertical displacement from the beam axis: \(\max(|y|)\) (meters)

max_abs_py_over_p0 float

Maximum \(y\)-plane transverse momentum \(\max(|p_y/p_0|)\) (unitless)

max_phase float

Maximum deviation in phase (deg)

max_energy_dev float

Maximum deviation in particle energy (MeV)

impact.z.output.MaxAmplitudeExtended

Bases: FortranOutputFileData

File fort.27: maximum amplitude information (extended)

Attributes:

Name Type Description
z float

longitudinal position along the beamline (meters)

max_abs_x float

Maximum horizontal displacement from the beam axis \(\max(|x|)\) (meters)

max_abs_gammabeta_x float

Maximum \(x\)-plane transverse momentum \(\max(|\gamma\beta_x|)\) (dimensionless)

max_abs_y float

Maximum vertical displacement from the beam axis \(\max(|y|)\) (meters)

max_abs_gammabeta_y float

Maximum \(x\)-plane transverse momentum \(\max(|\gamma\beta_y|)\) (dimensionless)

max_phase float

Maximum deviation in phase (deg???)

max_gamma_rel float

Maximum deviation in relativistic gamma \(\max(|\gamma - \gamma_0)|)\) (dimensionless)

impact.z.output.LoadBalanceLossDiagnostic

Bases: FortranOutputFileData

File fort.28: Load balance and loss diagnostic.

Attributes:

Name Type Description
z float

z distance (m)

loadbalance_min_n_particle float

Minimum number of particles on a PE

loadbalance_max_n_particle float

Maximum number of particles on a PE

n_particle float

Total number of particles in the bunch

impact.z.output.BeamDistribution3rdStandard

Bases: FortranOutputFileData

File fort.29: cubic root of 3rd moments of the beam distribution

Attributes:

Name Type Description
z float

z distance (m)

moment3_x float

Cubic root of the third moment of the horizontal position \(M_3(x)\) (meters)

moment3_px_over_p0 float

Cubic root of the third moment of the horizontal momentum \(M_3(p_x/p_0)\) (unitless)

moment3_y float

Cubic root of the third moment of the vertical position \(M_3(y)\) (meters)

moment3_py_over_p0 float

Cubic root of the third moment of the vertical momentum \(M_3(p_y/p_0)\) (unitless)

moment3_phase float

Cubic root of the third moment of the phase \(M_3(\phi)\) (deg)

moment3_energy float

Cubic root of the energy \(M_3(E)\) (MeV) In the file, this is stored as MeV and LUME-Impact converts to eV automatically.

impact.z.output.BeamDistribution3rdExtended

Bases: FortranOutputFileData

File fort.29: contains radius moments of the beam distribution.

Attributes:

Name Type Description
z float

z distance (m)

mean_r float

Mean radius (meters)

sigma_r float

RMS radius (meters)

mean_r_90percent float

90 percent mean radius (meters)

mean_r_95percent float

95 percent mean radius (meters)

mean_r_99percent float

99 percent mean radius (meters)

max_r_rel float

Maximum radius (meters)

impact.z.output.BeamDistribution4th

Bases: FortranOutputFileData

File fort.30 with diagnostic_type=1 contains the cubic root of the third moments of the beam distribution.

Here $ M_4(x) \equiv\left< (x-\left< x \right>)^4 \right>^{1/4} $, averaging over all particles.

Attributes:

Name Type Description
z float

Longitudinal position along the beamline (meters)

moment4_x float

Fourth root of the third moment of the horizontal position \(M_4(x)\) (meters)

moment4_px_over_p0 float

Fourth root of the third moment of the horizontal momentum \(M_4(p_x/p_0)\) (rad)

moment4_y float

Fourth root of the third moment of the vertical position \(M_4(y)\) (meters)

moment4_py_over_p0 float

Fourth root of the third moment of the vertical momentum \(M_4(p_y/p_0)\) (rad)

moment4_phase float

Fourth root of the third moment of the phase \(M_4(\phi)\) (deg)

moment4_energy float

Fourth root of the energy \(M_4(E)\) (MeV) In the file, this is stored as MeV and LUME-Impact converts to eV automatically.

impact.z.output.ParticlesAtChargedState

Bases: FortranOutputFileData

File fort.32: number of particles for each charge state

This file contains data about the number of particles for each charge state at different z distances.

Attributes:

Name Type Description
z float

The z distance in meters.

charge_state_n_particle int

The number of particles for each charge state.