Genesis4 FODO Scan Model¶
This is a variation on the original FODO scan notebook, using more complicated objects to automate some behavior.
In [1]:
Copied!
import warnings
from dataclasses import dataclass
from math import sqrt
import matplotlib.pyplot as plt
import numpy as np
from genesis.version4 import Genesis4, MainInput
import genesis.version4 as g4
%config InlineBackend.figure_format = 'retina'
# Ignore MPI warnings
warnings.filterwarnings("ignore")
import warnings
from dataclasses import dataclass
from math import sqrt
import matplotlib.pyplot as plt
import numpy as np
from genesis.version4 import Genesis4, MainInput
import genesis.version4 as g4
%config InlineBackend.figure_format = 'retina'
# Ignore MPI warnings
warnings.filterwarnings("ignore")
Lattice¶
These are the same parameters as the Genesis4 SASE Benchmark. The lattice is created here in code.
Note that the FODO cell length is 9.5 m, and there are 13 cells, for a total length of 123.5 m.
In [2]:
Copied!
def make_fodo(
    Lcell=9.5,
    kL=None,
    lambdau=15e-3,
    Ltot=150,
    Lpad=0.685 / 2,  # Will be adjusted slightly to make Lcell exact
    Lquad=0.08,
    aw=0.84853,
):
    if kL is None:
        #  Optimal for flat beam
        kL = 1 / (sqrt(2) * Lcell / 2)
        # Optimal for round beam (90 deg phase advance)
        # k1L_optimal = 2*sqrt(2) / Lcell
    k1 = kL / Lquad
    # Length for single wiggler
    Lwig = (Lcell - 4 * Lpad - 2 * Lquad) / 2
    nwig = round(Lwig / lambdau)
    # Set padding exactly
    Lwig = lambdau * nwig
    Lpad = round((Lcell - 2 * Lwig - 2 * Lquad) / 4, 9)
    ncell = round(Ltot // Lcell)
    return g4.Lattice(
        elements={
            "D1": g4.Drift(L=Lpad),
            "D2": g4.Drift(L=Lpad),
            # k1*L = kL
            "QF": g4.Quadrupole(L=Lquad, k1=k1),
            "QD": g4.Quadrupole(L=Lquad, k1=-k1),
            "UND": g4.Undulator(lambdau=lambdau, nwig=nwig, aw=aw),
            "FODO": g4.Line(
                elements=["UND", "D1", "QF", "D2", "UND", "D1", "QD", "D2"]
            ),
            "LAT": g4.Line(elements=["FODO"] * ncell),
        }
    )
print(make_fodo(lambdau=30e-3, Lcell=9.5))
def make_fodo(
    Lcell=9.5,
    kL=None,
    lambdau=15e-3,
    Ltot=150,
    Lpad=0.685 / 2,  # Will be adjusted slightly to make Lcell exact
    Lquad=0.08,
    aw=0.84853,
):
    if kL is None:
        #  Optimal for flat beam
        kL = 1 / (sqrt(2) * Lcell / 2)
        # Optimal for round beam (90 deg phase advance)
        # k1L_optimal = 2*sqrt(2) / Lcell
    k1 = kL / Lquad
    # Length for single wiggler
    Lwig = (Lcell - 4 * Lpad - 2 * Lquad) / 2
    nwig = round(Lwig / lambdau)
    # Set padding exactly
    Lwig = lambdau * nwig
    Lpad = round((Lcell - 2 * Lwig - 2 * Lquad) / 4, 9)
    ncell = round(Ltot // Lcell)
    return g4.Lattice(
        elements={
            "D1": g4.Drift(L=Lpad),
            "D2": g4.Drift(L=Lpad),
            # k1*L = kL
            "QF": g4.Quadrupole(L=Lquad, k1=k1),
            "QD": g4.Quadrupole(L=Lquad, k1=-k1),
            "UND": g4.Undulator(lambdau=lambdau, nwig=nwig, aw=aw),
            "FODO": g4.Line(
                elements=["UND", "D1", "QF", "D2", "UND", "D1", "QD", "D2"]
            ),
            "LAT": g4.Line(elements=["FODO"] * ncell),
        }
    )
print(make_fodo(lambdau=30e-3, Lcell=9.5))
Lattice(
  elements={
    'D1': Drift(L=0.34),
    'D2': Drift(L=0.34),
    'QF': Quadrupole(L=0.08, k1=1.860807318911967),
    'QD': Quadrupole(L=0.08, k1=-1.860807318911967),
    'UND': Undulator(aw=0.84853, lambdau=0.03, nwig=133),
    'FODO': Line(elements=['UND', 'D1', 'QF', 'D2', 'UND', 'D1', 'QD', 'D2']),
    'LAT': Line(
      elements=[
        'FODO',
        'FODO',
        'FODO',
        'FODO',
        'FODO',
        'FODO',
        'FODO',
        'FODO',
        'FODO',
        'FODO',
        'FODO',
        'FODO',
        'FODO',
        'FODO',
        'FODO',
      ],
    ),
  },
)
In [3]:
Copied!
def Krof(*, lambdar, lambdau, gamma):
    """
    K to make lambdar resonant
    """
    Ksq = 2 * (2 * gamma**2 * lambdar / lambdau - 1)
    if Ksq <= 0:
        raise ValueError(
            f"No resonance available, lambdau must be < {2*gamma**2*lambdar*1e3:0.1f}1e-3 m"
        )
    return sqrt(Ksq)
Krof(lambdar=1e-10, lambdau=25e-3, gamma=11357.82) / sqrt(2)
def Krof(*, lambdar, lambdau, gamma):
    """
    K to make lambdar resonant
    """
    Ksq = 2 * (2 * gamma**2 * lambdar / lambdau - 1)
    if Ksq <= 0:
        raise ValueError(
            f"No resonance available, lambdau must be < {2*gamma**2*lambdar*1e3:0.1f}1e-3 m"
        )
    return sqrt(Ksq)
Krof(lambdar=1e-10, lambdau=25e-3, gamma=11357.82) / sqrt(2)
Out[3]:
0.17888711865084084
In [4]:
Copied!
@dataclass
class FODOModel:
    # Lengths
    Lcell: float = 9.5
    Ltot: float = 125
    Lquad: float = 0.08
    Lpad: float = 0.685 / 2
    kL: float = None  # Will be picked automatically
    lambdar: float = 1e-10
    lambdau: float = 15e-3
    gamma: float = 11357.82
    current: float = 3000
    delgam: float = 1.0
    norm_emit_x: float = 0.4e-6
    norm_emit_y: float = 0.4e-6
    nproc = 0  # Auto-select
    seed: int = None  # None will pick a random seed
    def make_lattice(self):
        """Returns the Lattice instance, setting aw for resonance."""
        aw = Krof(lambdar=self.lambdar, lambdau=self.lambdau, gamma=self.gamma) / sqrt(
            2
        )
        return make_fodo(
            Lcell=self.Lcell,
            kL=self.kL,
            lambdau=self.lambdau,
            Ltot=self.Ltot,
            Lpad=self.Lpad,  # Will be adjusted
            Lquad=self.Lquad,
            aw=aw,
        )
    def make_genesis(self):
        if self.seed is None:
            seed = np.random.randint(0, 1e10)
        else:
            seed = self.seed
        input = MainInput(
            namelists=[
                g4.Setup(
                    rootname="Benchmark",
                    beamline="LAT",
                    lambda0=self.lambdar,
                    gamma0=self.gamma,
                    delz=0.045,
                    seed=seed,
                    shotnoise=False,
                    beam_global_stat=True,
                    field_global_stat=True,
                ),
                g4.LatticeNamelist(zmatch=self.Lcell),
                g4.Field(
                    power=5000,
                    dgrid=0.0002,
                    ngrid=255,
                    waist_size=3e-05,
                ),
                g4.Beam(
                    current=self.current,
                    delgam=self.delgam,
                    ex=self.norm_emit_x,
                    ey=self.norm_emit_y,
                ),
                g4.Track(zstop=self.Ltot),
            ]
        )
        lattice = self.make_lattice()
        return Genesis4(input, lattice, verbose=True)
    def run(self):
        G = self.make_genesis()
        G.nproc = self.nproc
        G.run()
        return G
fodo_model = FODOModel(Ltot=20)
@dataclass
class FODOModel:
    # Lengths
    Lcell: float = 9.5
    Ltot: float = 125
    Lquad: float = 0.08
    Lpad: float = 0.685 / 2
    kL: float = None  # Will be picked automatically
    lambdar: float = 1e-10
    lambdau: float = 15e-3
    gamma: float = 11357.82
    current: float = 3000
    delgam: float = 1.0
    norm_emit_x: float = 0.4e-6
    norm_emit_y: float = 0.4e-6
    nproc = 0  # Auto-select
    seed: int = None  # None will pick a random seed
    def make_lattice(self):
        """Returns the Lattice instance, setting aw for resonance."""
        aw = Krof(lambdar=self.lambdar, lambdau=self.lambdau, gamma=self.gamma) / sqrt(
            2
        )
        return make_fodo(
            Lcell=self.Lcell,
            kL=self.kL,
            lambdau=self.lambdau,
            Ltot=self.Ltot,
            Lpad=self.Lpad,  # Will be adjusted
            Lquad=self.Lquad,
            aw=aw,
        )
    def make_genesis(self):
        if self.seed is None:
            seed = np.random.randint(0, 1e10)
        else:
            seed = self.seed
        input = MainInput(
            namelists=[
                g4.Setup(
                    rootname="Benchmark",
                    beamline="LAT",
                    lambda0=self.lambdar,
                    gamma0=self.gamma,
                    delz=0.045,
                    seed=seed,
                    shotnoise=False,
                    beam_global_stat=True,
                    field_global_stat=True,
                ),
                g4.LatticeNamelist(zmatch=self.Lcell),
                g4.Field(
                    power=5000,
                    dgrid=0.0002,
                    ngrid=255,
                    waist_size=3e-05,
                ),
                g4.Beam(
                    current=self.current,
                    delgam=self.delgam,
                    ex=self.norm_emit_x,
                    ey=self.norm_emit_y,
                ),
                g4.Track(zstop=self.Ltot),
            ]
        )
        lattice = self.make_lattice()
        return Genesis4(input, lattice, verbose=True)
    def run(self):
        G = self.make_genesis()
        G.nproc = self.nproc
        G.run()
        return G
fodo_model = FODOModel(Ltot=20)
In [5]:
Copied!
G = fodo_model.run()
G = fodo_model.run()
Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpi4edz4ru Setting use_mpi = True because nproc = 12 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpi4edz4ru /Users/klauer/miniforge3/envs/lume-genesis/bin/mpirun -n 12 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:33:35 2024 MPI-Comm Size: 12 nodes Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 9.8922 alphax : -0.756992 phix (deg): 41.1658 betay (m) : 17.2567 alphay : 1.29774 phiy (deg): 44.1394 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:33:40 2024 Total Wall Clock Time: 4.47621 seconds ------------------------------------- Success - execution took 4.95s.
In [6]:
Copied!
G.plot(
    "field_peak_power", yscale="log", y2=["beam_xsize", "beam_ysize"], ylim2=(0, 50e-6)
)
G.plot(
    "field_peak_power", yscale="log", y2=["beam_xsize", "beam_ysize"], ylim2=(0, 50e-6)
)
Run1 function¶
Make a simple function to run a complete simulation and return a Genesis4 object.
In [7]:
Copied!
def run1(kL):
    G = FODOModel(kL=kL).run()
    return G
G2 = run1(0.136)  # = 1.7 * .08
def run1(kL):
    G = FODOModel(kL=kL).run()
    return G
G2 = run1(0.136)  # = 1.7 * .08
Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp41sp0ts8 Setting use_mpi = True because nproc = 12 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp41sp0ts8 /Users/klauer/miniforge3/envs/lume-genesis/bin/mpirun -n 12 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:33:41 2024 MPI-Comm Size: 12 nodes Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 11.1501 alphax : -0.775418 phix (deg): 37.4688 betay (m) : 18.1031 alphay : 1.24086 phiy (deg): 40.6908 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:34:08 2024 Total Wall Clock Time: 25.9716 seconds ------------------------------------- Success - execution took 27.29s.
In [8]:
Copied!
G2.plot(
    "field_peak_power", yscale="log", y2=["beam_xsize", "beam_ysize"], ylim2=(0, 200e-6)
)
G2.plot(
    "field_peak_power", yscale="log", y2=["beam_xsize", "beam_ysize"], ylim2=(0, 200e-6)
)
Scan kL¶
Now scan use this function to scan various quadrupole k1. Note that Genesis4 is doing matching internally with zmatch in the input.
In [9]:
Copied!
%%time
kLlist = np.linspace(0.1, 0.4, 10)
Glist = [run1(kL) for kL in kLlist]
%%time
kLlist = np.linspace(0.1, 0.4, 10)
Glist = [run1(kL) for kL in kLlist]
Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp6at8y_fz Setting use_mpi = True because nproc = 12 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp6at8y_fz /Users/klauer/miniforge3/envs/lume-genesis/bin/mpirun -n 12 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:34:08 2024 MPI-Comm Size: 12 nodes Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 16.405 alphax : -0.829624 phix (deg): 27.3204 betay (m) : 21.3934 alphay : 1.07148 phiy (deg): 31.5325 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:34:36 2024 Total Wall Clock Time: 26.4122 seconds ------------------------------------- Success - execution took 27.84s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp1y17d5m1 Setting use_mpi = True because nproc = 12 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp1y17d5m1 /Users/klauer/miniforge3/envs/lume-genesis/bin/mpirun -n 12 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:34:36 2024 MPI-Comm Size: 12 nodes Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 11.4414 alphax : -0.779297 phix (deg): 36.7076 betay (m) : 18.2972 alphay : 1.22899 phiy (deg): 39.9866 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:35:01 2024 Total Wall Clock Time: 24.1746 seconds ------------------------------------- Success - execution took 25.25s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpjle4gx3z Setting use_mpi = True because nproc = 12 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpjle4gx3z /Users/klauer/miniforge3/envs/lume-genesis/bin/mpirun -n 12 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:35:01 2024 MPI-Comm Size: 12 nodes Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 8.47309 alphax : -0.73226 phix (deg): 46.3579 betay (m) : 16.3003 alphay : 1.37698 phiy (deg): 49.0476 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:35:32 2024 Total Wall Clock Time: 28.5974 seconds ------------------------------------- Success - execution took 31.22s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp5evwi6k4 Setting use_mpi = True because nproc = 12 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp5evwi6k4 /Users/klauer/miniforge3/envs/lume-genesis/bin/mpirun -n 12 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:35:33 2024 MPI-Comm Size: 12 nodes Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 6.49469 alphax : -0.68826 phix (deg): 56.3707 betay (m) : 15.0411 alphay : 1.53488 phiy (deg): 58.6641 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:36:03 2024 Total Wall Clock Time: 28.8417 seconds ------------------------------------- Success - execution took 30.67s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpg6sadcut Setting use_mpi = True because nproc = 12 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpg6sadcut /Users/klauer/miniforge3/envs/lume-genesis/bin/mpirun -n 12 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:36:03 2024 MPI-Comm Size: 12 nodes Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 5.07466 alphax : -0.647288 phix (deg): 66.8776 betay (m) : 14.3277 alphay : 1.71806 phiy (deg): 68.8959 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:36:28 2024 Total Wall Clock Time: 24.1967 seconds ------------------------------------- Success - execution took 25.24s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpjw0goqyb Setting use_mpi = True because nproc = 12 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpjw0goqyb /Users/klauer/miniforge3/envs/lume-genesis/bin/mpirun -n 12 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:36:29 2024 MPI-Comm Size: 12 nodes Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 3.9959 alphax : -0.609747 phix (deg): 78.0657 betay (m) : 14.0882 alphay : 1.94566 phiy (deg): 79.8946 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:36:58 2024 Total Wall Clock Time: 27.5327 seconds ------------------------------------- Success - execution took 29.14s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmplawn41lw Setting use_mpi = True because nproc = 12 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmplawn41lw /Users/klauer/miniforge3/envs/lume-genesis/bin/mpirun -n 12 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:36:58 2024 MPI-Comm Size: 12 nodes Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 3.13531 alphax : -0.576908 phix (deg): 90.2238 betay (m) : 14.364 alphay : 2.25022 phiy (deg): 91.9346 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:37:23 2024 Total Wall Clock Time: 24.665 seconds ------------------------------------- Success - execution took 25.82s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp2wqcwgcf Setting use_mpi = True because nproc = 12 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp2wqcwgcf /Users/klauer/miniforge3/envs/lume-genesis/bin/mpirun -n 12 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:37:24 2024 MPI-Comm Size: 12 nodes Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 2.4145 alphax : -0.552359 phix (deg): 103.85 betay (m) : 15.3865 alphay : 2.70174 phiy (deg): 105.521 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:37:48 2024 Total Wall Clock Time: 24.0815 seconds ------------------------------------- Success - execution took 25.00s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmptzywm0wt Setting use_mpi = True because nproc = 12 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmptzywm0wt /Users/klauer/miniforge3/envs/lume-genesis/bin/mpirun -n 12 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:37:49 2024 MPI-Comm Size: 12 nodes Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 1.77384 alphax : -0.548072 phix (deg): 119.977 betay (m) : 17.9647 alphay : 3.50197 phiy (deg): 121.738 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:38:17 2024 Total Wall Clock Time: 27.1594 seconds ------------------------------------- Success - execution took 28.74s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpmf_9hwxk Setting use_mpi = True because nproc = 12 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpmf_9hwxk /Users/klauer/miniforge3/envs/lume-genesis/bin/mpirun -n 12 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:38:17 2024 MPI-Comm Size: 12 nodes Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 1.14869 alphax : -0.632829 phix (deg): 141.694 betay (m) : 26.7398 alphay : 5.74149 phiy (deg): 144.007 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:38:43 2024 Total Wall Clock Time: 24.4312 seconds ------------------------------------- Success - execution took 25.50s. CPU times: user 324 ms, sys: 114 ms, total: 438 ms Wall time: 4min 34s
In [10]:
Copied!
fig, ax = plt.subplots()
for k, g in zip(kLlist, Glist):
    x = g.stat("zplot")
    y = g.stat("field_peak_power")
    ax.plot(x, y / 1e6, label=f"{k:0.3f}")
# ax.set_yscale('log')
ax.set_xlabel(r"$z$ (m)")
ax.set_ylabel("power (MW)")
plt.legend(title=r"$k_1L$ (1/m)")
fig, ax = plt.subplots()
for k, g in zip(kLlist, Glist):
    x = g.stat("zplot")
    y = g.stat("field_peak_power")
    ax.plot(x, y / 1e6, label=f"{k:0.3f}")
# ax.set_yscale('log')
ax.set_xlabel(r"$z$ (m)")
ax.set_ylabel("power (MW)")
plt.legend(title=r"$k_1L$ (1/m)")
Out[10]:
<matplotlib.legend.Legend at 0x147e97680>
In [11]:
Copied!
fig, ax = plt.subplots()
y = np.array([g.stat("field_peak_power")[-1] for g in Glist])
ixbest = y.argmax()
Gbest = Glist[ixbest]
kbest = kLlist[ixbest]
ybest = y[ixbest]
ax.plot(kLlist, y / 1e6)
ax.scatter(kbest, ybest / 1e6, marker="*", label=rf"$k_1L$= {kbest:0.1f} 1/m")
ax.set_ylabel("end power (MW)")
ax.set_xlabel(r"$k_1L$ (1/m)")
plt.legend()
fig, ax = plt.subplots()
y = np.array([g.stat("field_peak_power")[-1] for g in Glist])
ixbest = y.argmax()
Gbest = Glist[ixbest]
kbest = kLlist[ixbest]
ybest = y[ixbest]
ax.plot(kLlist, y / 1e6)
ax.scatter(kbest, ybest / 1e6, marker="*", label=rf"$k_1L$= {kbest:0.1f} 1/m")
ax.set_ylabel("end power (MW)")
ax.set_xlabel(r"$k_1L$ (1/m)")
plt.legend()
Out[11]:
<matplotlib.legend.Legend at 0x147de70b0>
In [12]:
Copied!
fig, ax = plt.subplots()
for k, g in zip(kLlist, Glist):
    x = g.stat("zplot")
    y = g.stat("beam_xsize")
    if k == kbest:
        color = "black"
    else:
        color = None
    ax.plot(x, y * 1e6, label=f"{k:0.3f}", color=color)
ax.set_xlabel(r"$z$ (m)")
ax.set_ylabel("Beam xsize (µm)")
ax.set_ylim(0, None)
plt.legend(title=r"$k_1L$ (1/m)");
fig, ax = plt.subplots()
for k, g in zip(kLlist, Glist):
    x = g.stat("zplot")
    y = g.stat("beam_xsize")
    if k == kbest:
        color = "black"
    else:
        color = None
    ax.plot(x, y * 1e6, label=f"{k:0.3f}", color=color)
ax.set_xlabel(r"$z$ (m)")
ax.set_ylabel("Beam xsize (µm)")
ax.set_ylim(0, None)
plt.legend(title=r"$k_1L$ (1/m)");
Scan lambdau¶
In [13]:
Copied!
def run1(lambdau):
    p = {}
    p["lambdau"] = lambdau
    G = FODOModel(**p).make_genesis()
    G.run()
    return G
G2 = run1(10e-3)
G2.plot(
    "field_peak_power", yscale="log", y2=["beam_xsize", "beam_ysize"], ylim2=(0, None)
)
def run1(lambdau):
    p = {}
    p["lambdau"] = lambdau
    G = FODOModel(**p).make_genesis()
    G.run()
    return G
G2 = run1(10e-3)
G2.plot(
    "field_peak_power", yscale="log", y2=["beam_xsize", "beam_ysize"], ylim2=(0, None)
)
Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpwbpvkxx_ Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpwbpvkxx_ /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:38:43 2024 MPI-Comm Size: 1 node Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 9.89977 alphax : -0.757787 phix (deg): 41.1658 betay (m) : 14.1443 alphay : 1.05518 phiy (deg): 54.2739 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:38:53 2024 Total Wall Clock Time: 10.0861 seconds ------------------------------------- Success - execution took 10.27s.
In [14]:
Copied!
%%time
lambdaulist = np.linspace(10e-3, 25e-3, 10)
Glist = [run1(lambdau) for lambdau in lambdaulist]
%%time
lambdaulist = np.linspace(10e-3, 25e-3, 10)
Glist = [run1(lambdau) for lambdau in lambdaulist]
Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpfpw110q3 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpfpw110q3 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:38:54 2024 MPI-Comm Size: 1 node Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 9.89977 alphax : -0.757787 phix (deg): 41.1658 betay (m) : 14.1443 alphay : 1.05518 phiy (deg): 54.2739 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:39:04 2024 Total Wall Clock Time: 10.2502 seconds ------------------------------------- Success - execution took 10.50s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpupyh0fo8 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpupyh0fo8 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:39:04 2024 MPI-Comm Size: 1 node Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 9.8922 alphax : -0.756992 phix (deg): 41.1658 betay (m) : 15.6057 alphay : 1.16926 phiy (deg): 48.9958 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:39:15 2024 Total Wall Clock Time: 10.3603 seconds ------------------------------------- Success - execution took 10.63s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpfeebssjo Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpfeebssjo /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:39:15 2024 MPI-Comm Size: 1 node Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 9.89472 alphax : -0.757257 phix (deg): 41.1658 betay (m) : 16.5956 alphay : 1.24641 phiy (deg): 45.9485 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:39:25 2024 Total Wall Clock Time: 10.2384 seconds ------------------------------------- Success - execution took 10.50s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp655lnff8 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp655lnff8 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:39:26 2024 MPI-Comm Size: 1 node Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 9.8922 alphax : -0.756992 phix (deg): 41.1658 betay (m) : 17.2567 alphay : 1.29774 phiy (deg): 44.1394 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:39:36 2024 Total Wall Clock Time: 10.2317 seconds ------------------------------------- Success - execution took 10.47s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpc1fstlz6 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpc1fstlz6 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:39:36 2024 MPI-Comm Size: 1 node Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 9.89725 alphax : -0.757522 phix (deg): 41.1658 betay (m) : 17.6825 alphay : 1.33097 phiy (deg): 43.0199 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:39:46 2024 Total Wall Clock Time: 10.1847 seconds ------------------------------------- Success - execution took 10.40s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp3qph_uv2 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp3qph_uv2 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:39:46 2024 MPI-Comm Size: 1 node Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 9.90104 alphax : -0.757919 phix (deg): 41.1658 betay (m) : 17.9623 alphay : 1.35283 phiy (deg): 42.3127 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:39:57 2024 Total Wall Clock Time: 10.0746 seconds ------------------------------------- Success - execution took 10.29s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp7kf3k072 Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmp7kf3k072 /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:39:57 2024 MPI-Comm Size: 1 node Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 9.89977 alphax : -0.757787 phix (deg): 41.1658 betay (m) : 18.1549 alphay : 1.36773 phiy (deg): 41.8561 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:40:07 2024 Total Wall Clock Time: 10.0371 seconds ------------------------------------- Success - execution took 10.21s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmph4w6ue7c Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmph4w6ue7c /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:40:07 2024 MPI-Comm Size: 1 node Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 9.89472 alphax : -0.757257 phix (deg): 41.1658 betay (m) : 18.2919 alphay : 1.37818 phiy (deg): 41.5548 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:40:17 2024 Total Wall Clock Time: 10.1753 seconds ------------------------------------- Success - execution took 10.35s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpc44hgnts Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpc44hgnts /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:40:17 2024 MPI-Comm Size: 1 node Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 9.8922 alphax : -0.756992 phix (deg): 41.1658 betay (m) : 18.3836 alphay : 1.3852 phiy (deg): 41.3523 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Calculation: 100% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:40:28 2024 Total Wall Clock Time: 10.1215 seconds ------------------------------------- Success - execution took 10.29s. Configured to run in: /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpes31yz5r Running Genesis4 in /var/folders/vy/s8_hc3m10fddm6n_43cf_m8r0000gn/T/tmpes31yz5r /Users/klauer/miniforge3/envs/lume-genesis/bin/genesis4 -l genesis.lat genesis4.in --------------------------------------------- GENESIS - Version 4.6.6 has started... Compile info: Compiled by runner at 2024-01-11 18:10:26 [UTC] from Git Commit ID: Starting Time: Fri Aug 23 15:40:28 2024 MPI-Comm Size: 1 node Opened input file genesis4.in Parsing lattice file genesis.lat ... Matching for periodic solution between z = 0 and z = 9.5 : betax (m) : 9.90356 alphax : -0.758184 phix (deg): 41.1658 betay (m) : 18.4224 alphay : 1.38862 phiy (deg): 41.2149 Generating input radiation field for HARM = 1 ... Generating input particle distribution... Running Core Simulation... Steady-state run Initial analysis of electron beam and radiation field... Calculation: 0% done Calculation: 10% done Calculation: 20% done Calculation: 30% done Calculation: 40% done Calculation: 50% done Calculation: 60% done Calculation: 70% done Calculation: 80% done Calculation: 90% done Writing output file... Core Simulation done. End of Track Program is terminating... Ending Time: Fri Aug 23 15:40:38 2024 Total Wall Clock Time: 9.9477 seconds ------------------------------------- Success - execution took 10.12s. CPU times: user 211 ms, sys: 75.9 ms, total: 287 ms Wall time: 1min 43s
In [15]:
Copied!
fig, ax = plt.subplots()
for k, g in zip(lambdaulist, Glist):
    x = g.stat("zplot")
    y = g.stat("field_peak_power")
    ax.plot(x, y / 1e6, label=f"{k*1e3:0.1f}")
ax.set_yscale("log")
ax.set_xlabel(r"$z$ (m)")
ax.set_ylabel("power (MW)")
plt.legend(title=r"$\lambda_u$ (mm)")
fig, ax = plt.subplots()
for k, g in zip(lambdaulist, Glist):
    x = g.stat("zplot")
    y = g.stat("field_peak_power")
    ax.plot(x, y / 1e6, label=f"{k*1e3:0.1f}")
ax.set_yscale("log")
ax.set_xlabel(r"$z$ (m)")
ax.set_ylabel("power (MW)")
plt.legend(title=r"$\lambda_u$ (mm)")
Out[15]:
<matplotlib.legend.Legend at 0x1738ef890>