Skip to content

Fields

pmd_beamphysics.FieldMesh

FieldMesh(h5=None, data=None)

Class for openPMD External Field Mesh data.

Initialized on on openPMD beamphysics particle group:

  • h5: open h5 handle, or str that is a file
  • data: raw data

The required data is stored in ._data, and consists of dicts:

  • 'attrs'
  • 'components'

Component data is always 3D.

Initialization from openPMD-beamphysics HDF5 file:

  • FieldMesh('file.h5')

Initialization from a data dict:

  • FieldMesh(data=data)

Derived properties:

  • .r, .theta, .z
  • .Br, .Btheta, .Bz
  • .Er, .Etheta, .Ez
  • .E, .B

  • .phase

  • .scale
  • .factor

  • .harmonic

  • .frequency

  • .shape

  • .geometry
  • .mins, .maxs, .deltas
  • .meshgrid
  • .dr, .dtheta, .dz

Booleans:

  • .is_pure_electric
  • .is_pure_magnetic
  • .is_static

Units and labels

  • .units
  • .axis_labels

Plotting:

  • .plot
  • .plot_onaxis

Writers

  • .write
  • .write_astra_1d
  • .write_astra_3d
  • .write_gpt
  • .write_impact_emfield_cartesian
  • .to_cylindrical
  • .to_astra_1d
  • .to_impact_solrf
  • .to_impact_impact_emfield_cartesian
  • .write_gpt
  • .write_superfish

Constructors (class methods):

  • .from_ansys_ascii_3d
  • .from_astra_3d
  • .from_superfish
  • .from_onaxis
  • .expand_onaxis
Source code in pmd_beamphysics/fields/fieldmesh.py
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
def __init__(self, h5=None, data=None):
    if h5:
        # Allow filename
        if isinstance(h5, str):
            fname = os.path.expandvars(os.path.expanduser(h5))
            assert os.path.exists(fname), f"File does not exist: {fname}"

            with File(fname, "r") as hh5:
                fp = field_paths(hh5)
                assert len(fp) == 1, f"Number of field paths in {h5}: {len(fp)}"
                data = load_field_data_h5(hh5[fp[0]])

        else:
            data = load_field_data_h5(h5)
    else:
        data = load_field_data_dict(data)

    # Internal data
    self._data = data

Attributes

pmd_beamphysics.FieldMesh.axis_labels property
axis_labels
pmd_beamphysics.FieldMesh.coord_vecs property
coord_vecs

Uses gridSpacing, gridSize, and gridOriginOffset to return coordinate vectors.

pmd_beamphysics.FieldMesh.factor property
factor

factor to multiply fields by, possibly complex.

factor = scale * exp(i*phase)

pmd_beamphysics.FieldMesh.is_pure_electric property
is_pure_electric

Returns True if there are no non-zero mageneticField components

pmd_beamphysics.FieldMesh.is_pure_magnetic property
is_pure_magnetic

Returns True if there are no non-zero electricField components

pmd_beamphysics.FieldMesh.meshgrid property
meshgrid

Usses coordinate vectors to produce a standard numpy meshgrids.

pmd_beamphysics.FieldMesh.phase property writable
phase

Returns the complex argument phi = -2*pi*RFphase to multiply the oscillating field by.

Can be set.

Functions

pmd_beamphysics.FieldMesh.axis_index
axis_index(key)

Returns axis index for a named axis label key.

Examples:

  • .axis_labels == ('x', 'y', 'z')
  • .axis_index('z') returns 2
Source code in pmd_beamphysics/fields/fieldmesh.py
285
286
287
288
289
290
291
292
293
294
295
296
297
def axis_index(self, key):
    """
    Returns axis index for a named axis label key.

    Examples:

    - `.axis_labels == ('x', 'y', 'z')`
    - `.axis_index('z')` returns `2`
    """
    for i, name in enumerate(self.axis_labels):
        if name == key:
            return i
    raise ValueError(f"Axis not found: {key}")
pmd_beamphysics.FieldMesh.axis_points
axis_points(axis_label)

Returns 3D points for the specified axis to be used by the interpolator.

Parameters:

Name Type Description Default
axis_label str

The label of the coordinate axis. Example: 'r' for cylindrical geometries.

required

Returns:

Type Description
numpy.ndarray of shape (n, 3)

An array of 3D points, where the specified axis is populated, and other axes are zero.

Source code in pmd_beamphysics/fields/fieldmesh.py
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
def axis_points(self, axis_label):
    """
    Returns 3D points for the specified axis to be used by the interpolator.

    Parameters
    ----------
    axis_label : str
        The label of the coordinate axis. Example: 'r' for cylindrical geometries.

    Returns
    -------
    numpy.ndarray of shape (n, 3)
        An array of 3D points, where the specified axis is populated, and other axes are zero.
    """
    x = self.coord_vec(axis_label)
    points = np.zeros((len(x), 3))
    points[:, self.axis_index(axis_label)] = x
    return points
pmd_beamphysics.FieldMesh.axis_values
axis_values(axis_label, field_key, **kwargs)

Returns the values of the specified field along the given axis, allowing for partial replacement of points.

Parameters:

Name Type Description Default
axis_label str

The label of the coordinate axis.

required
field_key str

The key representing the field data to interpolate.

required
**kwargs dict

Key-value pairs to replace parts of the internal points array. The keys should be axis labels, and the values should be the corresponding values to set. Example: x=0, y=1 will set points along 'x' and 'y' axes.

{}

Returns:

Type Description
tuple(ndarray, ndarray)

A tuple containing: - An array of coordinate values along the specified axis. - An array of interpolated field values at the corresponding points.

Source code in pmd_beamphysics/fields/fieldmesh.py
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
def axis_values(self, axis_label, field_key, **kwargs):
    """
    Returns the values of the specified field along the given axis, allowing for partial replacement of points.

    Parameters
    ----------
    axis_label : str
        The label of the coordinate axis.
    field_key : str
        The key representing the field data to interpolate.
    **kwargs : dict
        Key-value pairs to replace parts of the internal points array.
        The keys should be axis labels, and the values should be the corresponding values to set.
        Example: `x=0, y=1` will set points along 'x' and 'y' axes.

    Returns
    -------
    tuple (numpy.ndarray, numpy.ndarray)
        A tuple containing:
        - An array of coordinate values along the specified axis.
        - An array of interpolated field values at the corresponding points.
    """
    points3d = self.axis_points(axis_label)

    # Replace parts of points3d with the values from kwargs
    for axis, value in kwargs.items():
        points3d[:, self.axis_index(axis)] = value

    vec = points3d[:, self.axis_index(axis_label)]
    values = self.interpolate(field_key, points3d)
    return vec, values
pmd_beamphysics.FieldMesh.component_is_zero
component_is_zero(key)

Returns True if all elements in a component are zero.

Source code in pmd_beamphysics/fields/fieldmesh.py
416
417
418
419
420
421
def component_is_zero(self, key):
    """
    Returns True if all elements in a component are zero.
    """
    a = self[key]
    return not np.any(a)
pmd_beamphysics.FieldMesh.coord_vec
coord_vec(key)

Gets the coordinate vector from a named axis key.

Source code in pmd_beamphysics/fields/fieldmesh.py
360
361
362
363
364
365
def coord_vec(self, key):
    """
    Gets the coordinate vector from a named axis key.
    """
    i = self.axis_index(key)
    return np.linspace(self.mins[i], self.maxs[i], self.shape[i])
pmd_beamphysics.FieldMesh.copy
copy()

Returns a deep copy

Source code in pmd_beamphysics/fields/fieldmesh.py
977
978
979
def copy(self):
    """Returns a deep copy"""
    return deepcopy(self)
pmd_beamphysics.FieldMesh.from_ansys_ascii_3d classmethod
from_ansys_ascii_3d(*, efile=None, hfile=None, frequency=None)

Class method to return a FieldMesh from ANSYS ASCII files.

The format of each file is: header1 (ignored) header2 (ignored) x y z re_fx im_fx re_fy im_fy re_fz im_fz ... in C order, with oscillations as exp(iomegat)

Parameters:

Name Type Description Default
efile

Filename with complex electric field data in V/m

None
hfile

Filename with complex magnetic H field data in A/m

None
frequency

Frequency in Hz

None

Returns:

Type Description
FieldMesh
Source code in pmd_beamphysics/fields/fieldmesh.py
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
@classmethod
def from_ansys_ascii_3d(cls, *, efile=None, hfile=None, frequency=None):
    """
    Class method to return a FieldMesh from ANSYS ASCII files.

    The format of each file is:
    header1 (ignored)
    header2 (ignored)
    x y z re_fx im_fx re_fy im_fy re_fz im_fz
    ...
    in C order, with oscillations as exp(i*omega*t)

    Parameters
    ----------
    efile: str
        Filename with complex electric field data in V/m

    hfile: str
        Filename with complex magnetic H field data in A/m

    frequency: float
        Frequency in Hz

    Returns
    -------
    FieldMesh

    """

    if frequency is None:
        raise ValueError("Please provide a frequency")

    data = read_ansys_ascii_3d_fields(efile, hfile, frequency=frequency)
    return cls(data=data)
pmd_beamphysics.FieldMesh.from_astra_3d classmethod
from_astra_3d(common_filename, frequency=0)

Class method to parse multiple 3D astra fieldmap files, based on the common filename.

Source code in pmd_beamphysics/fields/fieldmesh.py
664
665
666
667
668
669
670
671
672
@classmethod
def from_astra_3d(cls, common_filename, frequency=0):
    """
    Class method to parse multiple 3D astra fieldmap files,
    based on the common filename.
    """

    data = read_astra_3d_fieldmaps(common_filename, frequency=frequency)
    return cls(data=data)
pmd_beamphysics.FieldMesh.from_impact_emfield_cartesian classmethod
from_impact_emfield_cartesian(filename, frequency=0, eleAnchorPt='beginning')

Class method to read an Impact-T style 1Tv3.T7 file corresponding to the 111: EMfldCart element.

Parameters:

Name Type Description Default
filename str

Path to the file where the field data will be written.

required
frequency float

Fundamental frequency in Hz This simply adds 'fundamentalFrequency' to attrs default=0

0

Returns:

Type Description
FieldMesh
Source code in pmd_beamphysics/fields/fieldmesh.py
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
@classmethod
def from_impact_emfield_cartesian(
    cls, filename, frequency=0, eleAnchorPt="beginning"
):
    """
    Class method to read an Impact-T style 1Tv3.T7 file corresponding to
    the `111: EMfldCart` element.

    Parameters
    ----------
    filename : str
        Path to the file where the field data will be written.


    frequency : float, optional
        Fundamental frequency in Hz
        This simply adds 'fundamentalFrequency' to attrs
        default=0

    Returns
    -------
        FieldMesh
    """

    attrs, components = parse_impact_emfield_cartesian(filename)

    # These aren't in the file, they must be added
    attrs["fundamentalFrequency"] = frequency
    if frequency == 0:
        attrs["harmonic"] = 0

    attrs["eleAnchorPt"] = eleAnchorPt
    return cls(data=dict(attrs=attrs, components=components))
pmd_beamphysics.FieldMesh.from_onaxis classmethod
from_onaxis(*, z=None, Bz=None, Ez=None, frequency=0, harmonic=None, eleAnchorPt='beginning')

Parameters:

Name Type Description Default
z

z-coordinates. Must be regularly spaced.

None
Bz

magnetic field at r=0 in T Default: None

None
Ez

Electric field at r=0 in V/m Default: None

None
frequency

fundamental frequency in Hz. Default: 0

0
harmonic

Harmonic of the fundamental the field actually oscillates at. Default: 1 if frequency !=0, otherwise 0.

None
eleAnchorPt

Element anchor point. Should be one of 'beginning', 'center', 'end' Default: 'beginning'

'beginning'

Returns:

Name Type Description
field FieldMesh

Instantiated fieldmesh

Source code in pmd_beamphysics/fields/fieldmesh.py
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
@classmethod
def from_onaxis(
    cls,
    *,
    z=None,
    Bz=None,
    Ez=None,
    frequency=0,
    harmonic=None,
    eleAnchorPt="beginning",
):
    """


    Parameters
    ----------
    z: array
        z-coordinates. Must be regularly spaced.

    Bz: array, optional
        magnetic field at r=0 in T
        Default: None

    Ez: array, optional
        Electric field at r=0 in V/m
        Default: None

    frequency: float, optional
        fundamental frequency in Hz.
        Default: 0

    harmonic: int, optional
        Harmonic of the fundamental the field actually oscillates at.
        Default: 1 if frequency !=0, otherwise 0.

    eleAnchorPt: str, optional
        Element anchor point.
        Should be one of 'beginning', 'center', 'end'
        Default: 'beginning'


    Returns
    -------
    field: FieldMesh
        Instantiated fieldmesh

    """

    # Get spacing
    nz = len(z)
    dz = np.diff(z)
    if not np.allclose(dz, dz[0]):
        raise NotImplementedError("Irregular spacing not implemented")
    dz = dz[0]

    components = {}
    if Ez is not None:
        Ez = np.squeeze(np.array(Ez))
        if Ez.ndim != 1:
            raise ValueError(f"Ez ndim = {Ez.ndim} must be 1")
        components["electricField/z"] = Ez.reshape(1, 1, len(Ez))

    if Bz is not None:
        Bz = np.squeeze(np.array(Bz))
        if Bz.ndim != 1:
            raise ValueError(f"Bz ndim = {Bz.ndim} must be 1")
        components["magneticField/z"] = Bz.reshape(1, 1, len(Bz))

    if Bz is None and Ez is None:
        raise ValueError("Please enter Ez or Bz")

    # Handle harmonic options
    if frequency == 0:
        harmonic = 0
    elif harmonic is None:
        harmonic = 1

    attrs = {
        "eleAnchorPt": eleAnchorPt,
        "gridGeometry": "cylindrical",
        "axisLabels": np.array(["r", "theta", "z"], dtype="<U5"),
        "gridLowerBound": np.array([0, 0, 0]),
        "gridOriginOffset": np.array([0.0, 0.0, z.min()]),
        "gridSpacing": np.array([0.0, 0.0, dz]),
        "gridSize": np.array([1, 1, nz]),
        "harmonic": harmonic,
        "fundamentalFrequency": frequency,
        "RFphase": 0,
        "fieldScale": 1.0,
    }

    data = dict(attrs=attrs, components=components)
    return cls(data=data)
pmd_beamphysics.FieldMesh.from_superfish classmethod
from_superfish(filename, type=None, geometry='cylindrical')

Class method to parse a superfish T7 style file.

Source code in pmd_beamphysics/fields/fieldmesh.py
708
709
710
711
712
713
714
715
716
@classmethod
@functools.wraps(read_superfish_t7)
def from_superfish(cls, filename, type=None, geometry="cylindrical"):
    """
    Class method to parse a superfish T7 style file.
    """
    data = read_superfish_t7(filename, type=type, geometry=geometry)
    c = cls(data=data)
    return c
pmd_beamphysics.FieldMesh.interpolate
interpolate(key, points)

Interpolates the field data for the given key at specified points.

Parameters:

Name Type Description Default
key str

The key representing the field data to interpolate.

required
points numpy.ndarray of shape (3,) or (n, 3)

An array of n 3d points at which to interpolate the field data. The points should be ordered according to .axis_labels.

required

Returns:

Type Description
ndarray

The interpolated field values at the specified points.

Source code in pmd_beamphysics/fields/fieldmesh.py
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
def interpolate(self, key, points):
    """
    Interpolates the field data for the given key at specified points.

    Parameters
    ----------
    key : str
        The key representing the field data to interpolate.
    points : numpy.ndarray of shape (3,) or (n, 3)
        An array of n 3d points at which to interpolate the field data. The points should
        be ordered according to `.axis_labels`.

    Returns
    -------
    numpy.ndarray
        The interpolated field values at the specified points.
    """
    points = np.array(points)

    # Convenience for a single point
    if len(points.shape) == 1:
        return self.interpolator(key)([points])[0]

    return self.interpolator(key)(points)
pmd_beamphysics.FieldMesh.interpolator
interpolator(key)

Returns an interpolator for a given field key.

Parameters:

Name Type Description Default
key str

The key representing the field data to interpolate. Examples include: - 'Ez' for scaled/phased data - 'magneticField/y' for raw component data

required

Returns:

Type Description
RegularGridInterpolator

An interpolator object that can be used to interpolate points. The points to interpolate should be ordered according to .axis_labels.

Source code in pmd_beamphysics/fields/fieldmesh.py
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
def interpolator(self, key):
    """
    Returns an interpolator for a given field key.

    Parameters
    ----------
    key : str
        The key representing the field data to interpolate. Examples include:
        - 'Ez' for scaled/phased data
        - 'magneticField/y' for raw component data

    Returns
    -------
    RegularGridInterpolator
        An interpolator object that can be used to interpolate points. The points
        to interpolate should be ordered according to `.axis_labels`.
    """
    field = self[key]
    return RegularGridInterpolator(
        tuple(map(self.coord_vec, self.axis_labels)), field
    )
pmd_beamphysics.FieldMesh.scaled_component
scaled_component(key)

Retruns a component scaled by the complex factor factor = scaleexp(iphase)

Source code in pmd_beamphysics/fields/fieldmesh.py
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
def scaled_component(self, key):
    """

    Retruns a component scaled by the complex factor
        factor = scale*exp(i*phase)


    """

    if key in self.components:
        dat = self.components[key]
    # Aliases
    elif key in component_from_alias:
        comp = component_from_alias[key]
        if comp in self.components:
            dat = self.components[comp]
        else:
            # Component not present, make zeros
            return np.zeros(self.shape)
    else:
        raise ValueError(f"Component not available: {key}")

    # Multiply by scale factor
    factor = self.factor

    if factor != 1:
        return factor * dat
    else:
        return dat
pmd_beamphysics.FieldMesh.to_cylindrical
to_cylindrical()

Returns a new FieldMesh in cylindrical geometry.

If the current geometry is rectangular, this will use the y=0 slice.

Source code in pmd_beamphysics/fields/fieldmesh.py
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
def to_cylindrical(self):
    """
    Returns a new FieldMesh in cylindrical geometry.

    If the current geometry is rectangular, this
    will use the y=0 slice.

    """
    if self.geometry == "rectangular":
        return FieldMesh(
            data=fieldmesh_rectangular_to_cylindrically_symmetric_data(self)
        )
    elif self.geometry == "cylindrical":
        return self
    else:
        raise NotImplementedError(f"geometry not implemented: {self.geometry}")
pmd_beamphysics.FieldMesh.units
units(key)

Returns the units of any key

Source code in pmd_beamphysics/fields/fieldmesh.py
516
517
518
519
520
521
522
523
524
525
526
def units(self, key):
    """Returns the units of any key"""

    # Strip any operators
    _, key = get_operator(key)

    # Fill out aliases
    if key in component_from_alias:
        key = component_from_alias[key]

    return pg_units(key)
pmd_beamphysics.FieldMesh.write
write(h5, name=None)

Writes openPMD-beamphysics format to an open h5 handle, or new file if h5 is a str.

Source code in pmd_beamphysics/fields/fieldmesh.py
529
530
531
532
533
534
535
536
537
538
539
540
541
542
def write(self, h5, name=None):
    """
    Writes openPMD-beamphysics format to an open h5 handle, or new file if h5 is a str.

    """
    if isinstance(h5, str):
        fname = os.path.expandvars(os.path.expanduser(h5))
        h5 = File(fname, "w")
        pmd_field_init(h5, externalFieldPath="/ExternalFieldPath/%T/")
        g = h5.create_group("/ExternalFieldPath/1/")
    else:
        g = h5

    write_pmd_field(g, self.data, name=name)
pmd_beamphysics.FieldMesh.write_gpt
write_gpt(filePath, asci2gdf_bin=None, verbose=True)

Writes a GPT field file.

Source code in pmd_beamphysics/fields/fieldmesh.py
581
582
583
584
585
586
587
588
def write_gpt(self, filePath, asci2gdf_bin=None, verbose=True):
    """
    Writes a GPT field file.
    """

    return write_gpt_fieldmap(
        self, filePath, asci2gdf_bin=asci2gdf_bin, verbose=verbose
    )
pmd_beamphysics.FieldMesh.write_impact_emfield_cartesian
write_impact_emfield_cartesian(filename)

Writes Impact-T style 1Tv3.T7 file corresponding to the 111: EMfldCart element.

Parameters:

Name Type Description Default
filename str

Path to the file where the field data will be written.

required
Source code in pmd_beamphysics/fields/fieldmesh.py
590
591
592
593
594
595
596
597
598
599
600
601
602
603
@functools.wraps(write_impact_emfield_cartesian)
def write_impact_emfield_cartesian(self, filename):
    """
    Writes Impact-T style 1Tv3.T7 file corresponding to
    the `111: EMfldCart` element.

    Parameters
    ----------
    filename : str
        Path to the file where the field data will be written.

    """

    return write_impact_emfield_cartesian(self, filename)
pmd_beamphysics.FieldMesh.write_superfish
write_superfish(filePath, verbose=False)

Write a Superfish T7 file.

For static fields, a Poisson T7 file is written.

For dynamic (harmonic /= 0) fields, a Fish T7 file is written

Source code in pmd_beamphysics/fields/fieldmesh.py
606
607
608
609
610
611
612
613
614
615
@functools.wraps(write_superfish_t7)
def write_superfish(self, filePath, verbose=False):
    """
    Write a Superfish T7 file.

    For static fields, a Poisson T7 file is written.

    For dynamic (`harmonic /= 0`) fields, a Fish T7 file is written
    """
    return write_superfish_t7(self, filePath, verbose=verbose)