Skip to content

API

impact.Impact

Impact(*args, group=None, always_autophase=False, **kwargs)

Bases: CommandWrapper

Files will be written into a temporary directory within workdir. If workdir=None, a location will be determined by the system.

Source code in impact/impact.py
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
def __init__(self, *args, group=None, always_autophase=False, **kwargs):
    super().__init__(*args, **kwargs)
    # Save init
    self.original_input_file = self.input_file

    self.input = {"header": {}, "lattice": []}
    self.output = {}
    self._units = {}
    self._units.update(EXTRA_UNITS)
    self.group = {}

    # MPI
    self._nnode = 1

    # Convenience lookup of elements in lattice by name
    self.ele = {}

    # Autophase settings to be applied.
    # This will be cleared when actually autophasing
    self._autophase_settings = {}
    self.always_autophase = always_autophase

    # Call configure
    if self.input_file:
        infile = tools.full_path(self.input_file)
        assert os.path.exists(infile), f"Impact input file does not exist: {infile}"
        self.load_input(self.input_file)
        self.configure()

        # Add groups, if any.
        if group:
            for k, v in group.items():
                self.add_group(k, **v)

    else:
        self.vprint("Using default input: 1 m drift lattice")
        self.input = deepcopy(DEFAULT_INPUT)
        self.configure()

Attributes

impact.Impact.cathode_start property writable
cathode_start

Returns a bool if cathode_start is requested. Can also be set as a bool.

impact.Impact.fieldmaps property
fieldmaps

Convenience pointer to .input['fieldmaps']

impact.Impact.header property
header

Convenience pointer to .input['header']

impact.Impact.lattice property
lattice

Convenience pointer to .input['lattice']

impact.Impact.nnode property writable
nnode

Number of MPI nodes to use

impact.Impact.numprocs property writable
numprocs

Number of MPI processors = Npcol*Nprow

impact.Impact.particles property
particles

Convenience pointer to .input['lattice']

impact.Impact.total_charge property writable
total_charge

Returns the total bunch charge in C. Can be set.

Functions

impact.Impact.add_ele
add_ele(ele)

Adds an element to .lattice

Source code in impact/impact.py
109
110
111
112
113
114
115
116
117
def add_ele(self, ele):
    """
    Adds an element to .lattice
    """
    name = ele["name"]
    assert name not in self.lattice, "Element already exists"
    insert_ele_by_s(ele, self.lattice, verbose=self.verbose)
    # Add to the ele dict
    self.ele[name] = ele
impact.Impact.add_group
add_group(name, **kwargs)

Add a control group. See control.py

Source code in impact/impact.py
119
120
121
122
123
124
125
126
127
128
129
130
131
def add_group(self, name, **kwargs):
    """
    Add a control group. See control.py
    """
    assert name not in self.ele
    if name in self.group:
        self.vprint(f"Warning: group {name} already exists, overwriting.")

    g = ControlGroup(**kwargs, name=name)
    g.link(self.ele)
    self.group[name] = g

    return self.group[name]
impact.Impact.archive
archive(h5=None)

Archive all data to an h5 handle or filename.

If no file is given, a file based on the fingerprint will be created.

Source code in impact/impact.py
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
def archive(self, h5=None):
    """
    Archive all data to an h5 handle or filename.

    If no file is given, a file based on the fingerprint will be created.

    """
    if not h5:
        h5 = "impact_" + self.fingerprint() + ".h5"

    if isinstance(h5, str):
        new_h5file = True
        fname = os.path.expandvars(h5)
        g = h5py.File(fname, "w")
        self.vprint(f"Archiving to file {fname}")
    else:
        new_h5file = False
        g = h5

    # Write basic attributes
    archive.impact_init(g)

    # Initial particles
    if self.initial_particles:
        self.initial_particles.write(g, name="initial_particles")

    # All input
    archive.write_input_h5(g, self.input, name="input")

    # All output
    archive.write_output_h5(g, self.output, name="output", units=self._units)

    # Control groups
    if self.group:
        archive.write_control_groups_h5(g, self.group, name="control_groups")

    # Close file if created here.
    if new_h5file:
        g.close()

    return h5
impact.Impact.autophase
autophase(settings=None, full_output=False)

Calculate the relative phases of each rf element by tracking a single particle. This uses a fast method that operates outside of Impact

Parameters:

Name Type Description Default
settings

dict of ele_name:rel_phase_deg

None
full_output

type of output to return (see Returns)

False

Returns:

Type Description
if full_output = True retuns a dict of:

ele_name:info_dict

Otherwise returns a dict of:

ele_name:rel_phase_deg

which is the same format as settings.
Source code in impact/impact.py
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
def autophase(self, settings=None, full_output=False):
    """
    Calculate the relative phases of each rf element
    by tracking a single particle.
    This uses a fast method that operates outside of Impact

    Parameters
    ----------
    settings: dict, optional=None
        dict of ele_name:rel_phase_deg

    full_output: bool, optional = False
        type of output to return (see Returns)


    Returns
    -------
    if full_output = True retuns a dict of:
            ele_name:info_dict

    Otherwise returns a dict of:
        ele_name:rel_phase_deg
    which is the same format as settings.


    """

    if self.initial_particles:
        t0 = self.initial_particles["mean_t"]
        pz0 = self.initial_particles["mean_pz"]
    else:
        t0 = 0
        pz0 = 0

    return fast_autophase_impact(
        self,
        settings=settings,
        t0=t0,
        pz0=pz0,
        full_output=full_output,
        verbose=self.verbose,
    )
impact.Impact.autophase_bookkeeper
autophase_bookkeeper()

Searches for 'autophase_deg' attribute in all eles. If one is found, autophase is called.

If .always_autophase == True, calls autophase is called.

Returns:

Name Type Description
settings dict

Autophase settings found

Source code in impact/impact.py
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
def autophase_bookkeeper(self):
    """
    Searches for `'autophase_deg'` attribute in all eles.
    If one is found, autophase is called.

    If .always_autophase == True, calls autophase is called.

    Returns
    -------
    settings: dict
        Autophase settings found
    """
    if self._autophase_settings or self.always_autophase:
        if self.verbose:
            print("Autophase bookkeeper found settings, applying them")

        # Actual found settings
        settings = self.autophase(settings=self._autophase_settings)

        # Clear
        self._autophase_settings = {}

    else:
        settings = {}

    return settings
impact.Impact.ele_bookkeeper
ele_bookkeeper()

Link .ele = dict to the lattice elements by their 'name' field

Source code in impact/impact.py
231
232
233
234
235
def ele_bookkeeper(self):
    """
    Link .ele = dict to the lattice elements by their 'name' field
    """
    self.ele = ele_dict_from(self.input["lattice"])
impact.Impact.field
field(z=0, t=0, x=0, y=0, component='Ez')

Return the field component at a longitudinal position z at time t.

Warking: This is based on the parsed fieldmaps, and not calculated directly from Impact. Not all elements/parameters are implemented. Currently x, y must be 0.

Source code in impact/impact.py
269
270
271
272
273
274
275
276
277
278
279
280
def field(self, z=0, t=0, x=0, y=0, component="Ez"):
    """
    Return the field component at a longitudinal
    position z at time t.

    Warking: This is based on the parsed fieldmaps,
    and not calculated directly from Impact. Not all elements/parameters
    are implemented. Currently x, y must be 0.
    """
    return lattice_field(
        self.lattice, x=x, y=y, z=z, t=t, component=component, fmaps=self.fieldmaps
    )
impact.Impact.get_executable
get_executable()

Gets the full path of the executable from .command, .command_mpi Will search environmental variables: Impact.command_env='IMPACTT_BIN' Impact.command_mpi_env='IMPACTT_MPI_BIN'

Source code in impact/impact.py
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
def get_executable(self):
    """
    Gets the full path of the executable from .command, .command_mpi
    Will search environmental variables:
            Impact.command_env='IMPACTT_BIN'
            Impact.command_mpi_env='IMPACTT_MPI_BIN'

    """
    if self.use_mpi:
        exe = tools.find_executable(
            exename=self.command_mpi, envname=self.command_mpi_env
        )
    else:
        exe = tools.find_executable(exename=self.command, envname=self.command_env)
    return exe
impact.Impact.get_run_script
get_run_script(write_to_path=False, path=None)

Assembles the run script using self.mpi_run string of the form: 'mpirun -n {n} {command_mpi}'

Optionally writes a file 'run' with this line to path.

Source code in impact/impact.py
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
def get_run_script(self, write_to_path=False, path=None):
    """
    Assembles the run script using self.mpi_run string of the form:
        'mpirun -n {n} {command_mpi}'

    Optionally writes a file 'run' with this line to path.
    """

    n_procs = self.numprocs

    exe = self.get_executable()

    if self.use_mpi:
        # mpi_exe could be a complicated string like:
        # 'srun -N 1 --cpu_bind=cores {n} {command_mpi}'
        # 'mpirun -n {n} {command_mpi}'

        runscript = self.mpi_run.format(
            nnode=self.nnode, nproc=n_procs, command_mpi=exe
        )

    else:
        if n_procs > 1:
            raise ValueError("Error: n_procs > 1 but use_mpi = False")
        runscript = exe

    if write_to_path:
        if path is None:
            path = self.path
        path = os.path.join(path, "run")
        with open(path, "w") as f:
            f.write(runscript)
        tools.make_executable(path)
    return runscript
impact.Impact.load_archive
load_archive(h5, configure=True)

Loads input and output from archived h5 file.

See: Impact.archive

Source code in impact/impact.py
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
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
def load_archive(self, h5, configure=True):
    """
    Loads input and output from archived h5 file.

    See: Impact.archive
    """
    if isinstance(h5, str):
        fname = os.path.expandvars(h5)
        g = h5py.File(fname, "r")

        glist = archive.find_impact_archives(g)
        n = len(glist)
        if n == 0:
            # legacy: try top level
            message = "legacy"
        elif n == 1:
            gname = glist[0]
            message = f"group {gname} from"
            g = g[gname]
        else:
            raise ValueError(f"Multiple archives found in file {fname}: {glist}")

        self.vprint(f"Reading {message} archive file {h5}")
    else:
        g = h5

    self.input = archive.read_input_h5(g["input"], verbose=self.verbose)
    self.output, self._units = archive.read_output_h5(
        g["output"], verbose=self.verbose
    )
    self._units.update(EXTRA_UNITS)

    if "initial_particles" in g:
        self.initial_particles = ParticleGroup(h5=g["initial_particles"])

    if "control_groups" in g:
        self.group = archive.read_control_groups_h5(
            g["control_groups"], verbose=self.verbose
        )
    self.vprint("Loaded from archive. Note: Must reconfigure to run again.")
    self.configured = False

    if configure:
        self.configure()

        # Re-link groups
        # TODO: cleaner logic
        for _, cg in self.group.items():
            cg.link(self.ele)
impact.Impact.load_output
load_output()

Loads stats, slice_info, and particles.

Source code in impact/impact.py
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
def load_output(self):
    """
    Loads stats, slice_info, and particles.
    """
    self.output["stats"], u = load_stats(
        self.path, species=self.species, types=FORT_STAT_TYPES, verbose=self.verbose
    )
    self._units.update(u)

    # This is not always present
    dipole_stats, u = load_stats(
        self.path,
        species=self.species,
        types=FORT_DIPOLE_STAT_TYPES,
        verbose=self.verbose,
    )
    if dipole_stats:
        self.output["dipole_stats"] = dipole_stats
        self._units.update(u)

    self.output["slice_info"], u = load_slice_info(self.path, self.verbose)
    self._units.update(u)

    self.load_particles()
impact.Impact.old_plot
old_plot(y='sigma_x', x='mean_z', nice=True, include_layout=True)

Simple stat plot

Source code in impact/impact.py
904
905
906
907
908
def old_plot(self, y="sigma_x", x="mean_z", nice=True, include_layout=True):
    """
    Simple stat plot
    """
    return plot_stat(self, y=y, x=x, nice=nice)
impact.Impact.plot
plot(y=['sigma_x', 'sigma_y'], x='mean_z', xlim=None, ylim=None, ylim2=None, y2=[], nice=True, include_layout=True, include_labels=False, include_markers=True, include_particles=True, include_field=True, field_t=None, include_legend=True, return_figure=False, tex=True, **kwargs)
Source code in impact/impact.py
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
def plot(
    self,
    y=["sigma_x", "sigma_y"],
    x="mean_z",
    xlim=None,
    ylim=None,
    ylim2=None,
    y2=[],
    nice=True,
    include_layout=True,
    include_labels=False,
    include_markers=True,
    include_particles=True,
    include_field=True,
    field_t=None,
    include_legend=True,
    return_figure=False,
    tex=True,
    **kwargs,
):
    """ """

    # Just plot fieldmaps if there are no stats
    if "stats" not in self.output:
        return plot_layout(
            self,
            xlim=xlim,
            include_markers=include_markers,
            include_labels=include_labels,
            include_field=include_field,
            field_t=field_t,
            return_figure=return_figure,
            **kwargs,
        )

    return plot_stats_with_layout(
        self,
        ykeys=y,
        ykeys2=y2,
        xkey=x,
        xlim=xlim,
        ylim=ylim,
        ylim2=ylim2,
        nice=nice,
        tex=tex,
        include_layout=include_layout,
        include_labels=include_labels,
        include_field=include_field,
        field_t=field_t,
        include_markers=include_markers,
        include_particles=include_particles,
        include_legend=include_legend,
        return_figure=return_figure,
        **kwargs,
    )
impact.Impact.print_lattice
print_lattice()

Pretty printing of the lattice

Source code in impact/impact.py
966
967
968
969
970
971
972
def print_lattice(self):
    """
    Pretty printing of the lattice
    """
    for ele in self.input["lattice"]:
        line = ele_str(ele)
        print(line)
impact.Impact.run_impact
run_impact(verbose=False, timeout=None)

Runs Impact-T

Source code in impact/impact.py
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
def run_impact(self, verbose=False, timeout=None):
    """
    Runs Impact-T

    """

    # Clear output
    self.output = {}

    # Autophase
    autophase_settings = self.autophase_bookkeeper()
    if autophase_settings:
        self.output["autophase_info"] = autophase_settings

    run_info = self.output["run_info"] = {"error": False}

    # Run script, gets executables
    runscript = self.get_run_script()
    run_info["run_script"] = runscript

    t1 = time()
    run_info["start_time"] = t1

    self.vprint("Running Impact-T in " + self.path)
    self.vprint(runscript)
    # Write input
    self.write_input()

    # Remove previous files
    for f in fort_files(self.path):
        os.remove(f)

    if timeout:
        res = tools.execute2(runscript.split(), timeout=timeout, cwd=self.path)
        log = res["log"]
        self.error = res["error"]
        run_info["error"] = self.error
        run_info["why_run_error"] = res["why_error"]
    else:
        # Interactive output, for Jupyter
        log = []
        counter = 0
        for path in tools.execute(runscript.split(), cwd=self.path):
            # Fancy clearing of old lines
            counter += 1
            if verbose:
                if counter < 15:
                    print(path, end="")
                else:
                    print(
                        "\r",
                        path.strip() + ", elapsed: " + str(time() - t1),
                        end="",
                    )
            log.append(path)
        self.vprint("Finished.")
    self.log = log

    # Load output
    self.load_output()

    run_info["run_time"] = time() - t1

    self.finished = True
impact.Impact.stat
stat(key)

Array from .output['stats'][key]

Additional keys are avalable: 'mean_energy': mean energy 'Ez': z component of the electric field at the centroid particle 'Bz' z component of the magnetic field at the centroid particle 'cov_{a}__{b}': any symmetric covariance matrix term

Source code in impact/impact.py
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
def stat(self, key):
    """
    Array from `.output['stats'][key]`

    Additional keys are avalable:
        'mean_energy': mean energy
        'Ez': z component of the electric field at the centroid particle
        'Bz'  z component of the magnetic field at the centroid particle
        'cov_{a}__{b}': any symmetric covariance matrix term

    """

    if key in ("Ez", "Bz"):
        return self.centroid_field(component=key[0:2])

    if key == "mean_energy":
        return self.stat("mean_kinetic_energy") + self.mc2

    # Allow flipping covariance keys
    if key.startswith("cov_") and key not in self.output["stats"]:
        k1, k2 = key[4:].split("__")
        key = f"cov_{k2}__{k1}"

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

    return self.output["stats"][key]
impact.Impact.track
track(particles, s=None)

Track a ParticleGroup. An optional stopping s can be given.

Source code in impact/impact.py
845
846
847
848
849
850
851
def track(self, particles, s=None):
    """
    Track a ParticleGroup. An optional stopping s can be given.
    """
    if not s:
        s = self.stop
    return track_to_s(self, particles, s)
impact.Impact.track1
track1(x0=0, px0=0, y0=0, py0=0, z0=0, pz0=1e-15, t0=0, s=None, species=None)

Tracks a single particle with starting coordinates: x0, y0, z0 in meters px0, py0, pz0 in eV/c t0 in seconds

to a position 's' in meters

Used for phasing and scaling elements.

If successful, returns a ParticleGroup with the final particle.

Otherwise, returns None

Source code in impact/impact.py
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
def track1(
    self,
    x0=0,
    px0=0,
    y0=0,
    py0=0,
    z0=0,
    pz0=1e-15,
    t0=0,
    s=None,  # final s
    species=None,
):
    """
    Tracks a single particle with starting coordinates:
    x0, y0, z0 in meters
    px0, py0, pz0 in eV/c
    t0 in seconds

    to a position 's' in meters

    Used for phasing and scaling elements.

    If successful, returns a ParticleGroup with the final particle.

    Otherwise, returns None

    """
    if not s:
        s = self.stop

    if not species:
        species = self.species

    # Change to serial exe just for this
    n_procs_save = self.numprocs
    self.numprocs = 1
    result = track1_to_s(
        self,
        s=s,
        x0=x0,
        px0=px0,
        y0=y0,
        py0=py0,
        z0=z0,
        pz0=pz0,
        t0=t0,
        species=species,
    )
    self.numprocs = n_procs_save
    return result
impact.Impact.units
units(key)

pmd_unit of a given key

Source code in impact/impact.py
317
318
319
320
321
322
323
324
325
326
327
328
def units(self, key):
    """pmd_unit of a given key"""

    # Allow flipping covariance keys
    if key.startswith("cov_") and key not in self._units:
        k1, k2 = key[4:].split("__")
        key = f"cov_{k2}__{k1}"

    if key not in self._units:
        raise ValueError(f"Unknown unit for {key}")

    return self._units[key]
impact.Impact.vprint
vprint(*args)

verbose print

Source code in impact/impact.py
974
975
976
977
def vprint(self, *args):
    """verbose print"""
    if self.verbose:
        print(*args)
impact.Impact.write_input
write_input(input_filename='ImpactT.in', path=None)

Write all input.

If .initial_particles are given,

Source code in impact/impact.py
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
def write_input(self, input_filename="ImpactT.in", path=None):
    """
    Write all input.

    If .initial_particles are given,
    """

    if path is None:
        path = self.path

    pathlib.Path(path).mkdir(exist_ok=True, parents=True)

    filePath = os.path.join(path, input_filename)

    # Write fieldmaps
    for name, fieldmap in self.input["fieldmaps"].items():
        file = os.path.join(path, name)
        fieldmaps.write_fieldmap(file, fieldmap)

    # Initial particles (ParticleGroup)
    if self.initial_particles:
        self.write_initial_particles(update_header=True, path=path)

        # Check consistency
        if (
            self.header["Flagimg"] == 1
            and self.header["Nemission"] < 1
            and self.total_charge > 0
        ):
            raise ValueError(
                f"Cathode start with space charge must "
                f"set header['Nemission'] > 0. "
                f"The current value is {self.header['Nemission']}."
            )

    # Symlink
    elif self.header["Flagdist"] == 16:
        src = self.input["input_particle_file"]
        dest = os.path.join(path, "partcl.data")

        # Don't worry about overwriting in temporary directories
        if self._tempdir and os.path.exists(dest):
            os.remove(dest)

        if not os.path.exists(dest):
            writers.write_input_particles_from_file(src, dest, self.header["Np"])
        else:
            self.vprint("partcl.data already exits, will not overwrite.")

    # Check for point-to-point spacechage processor criteria
    if self.numprocs > 1:
        for ele in self.lattice:
            if ele["type"] == "point_to_point_spacecharge":
                Np = self.header["Np"]
                numprocs = self.numprocs
                if Np % numprocs != 0:
                    raise ValueError(
                        f"The number of electrons ({Np}) divided by the number of processors ({numprocs}) must be an integer."
                    )

    # Write main input file. This should come last.
    writers.write_impact_input(filePath, self.header, self.lattice)

    # Write run script
    self.get_run_script(write_to_path=True, path=path)