Treatment Plan Pencil Beam source#
Description#
In ion beam therapy, a target’s volume is irradiated by scanning a pencil beam across the transverse plane using steering magnets, while adjusting the beam energy to reach different depths within the target. The Treatment Plan source simulates the delivery of a scanned ion pencil beam treatment plan.
The treatment plan specifies, for each beam and each energy layer, the positions (x,y) of the spots to be irradiated and their weight, represented as either the number of particles or monitor units. Here, the (x,y) coordinates indicate where the beam intersects the plane of the room’s isocenter. The plan also includes the gantry angle for each beam.
The Treatment Plan source reads the treatment plan file provided by the user and generates a Pencil Beam source whose parameters are updated in order to cover the spots described in the plan. To determine the position and orientation of each pencil beam, essential beamline geometrical parameters are required, such as the position of the stearing magnets and the virtual source to isocenter distance. Moreover, the beamline energy dependent optics parameters and the energy spread should be characterized.
Note that the Treatment Plan source simulates only a single beam at a time. If the treatment plan includes multiple beams, each must be simulated individually. You can specify the beam to simulate by using the “beam_nr” property of the source, with beam numbering starting at 1.
image#
To set up a Treatment Plan source, the user shall provide:
a treatment plan file path: DICOM format or Gate 9 text file format. Alternatively, the user can provide a dictionary with the spots data (see opengate.contrib.tps.ionbeamtherapy -> spots_info_from_txt).
a beamline model: including geometrical and energy dependent parametrs’ description. (see opengate.contrib.beamlines.ionbeamline -> BeamlineModel).
the gantry rotation axis. Default is ‘z’.
the number of particles to simulate.
the ion type.
Additional functionality#
flat_generationflag: if True, the same number of particles is generated for each spot and the spot weight is applied to the energy deposition instead. Default is False.sorted_spot_generationflag: if True, the spots are irradiated one by one, following the order in the treatment plan. Otherwise the spots are selected by sampling from a probability density function. Default is False.
Here an example of how to set up a Treatment Plan source in the opengate simulation:
# beamline model
beamline = BeamlineModel()
beamline.name = None
beamline.radiation_types = "proton"
# polinomial coefficients
beamline.energy_mean_coeffs = [1, 0]
beamline.energy_spread_coeffs = [0.4417036946562556]
beamline.sigma_x_coeffs = [2.3335754]
beamline.theta_x_coeffs = [2.3335754e-3]
beamline.epsilon_x_coeffs = [0.00078728e-3]
beamline.sigma_y_coeffs = [1.96433431]
beamline.theta_y_coeffs = [0.00079118e-3]
beamline.epsilon_y_coeffs = [0.00249161e-3]
# tps
n_sim = 80000 # particles to simulate per beam
tps = sim.add_source("TreatmentPlanPBSource", "TPSource")
tps.n = n_sim
tps.beam_model = beamline
tps.plan_path = ref_path / "TreatmentPlan2Spots.txt"
tps.beam_nr = 1
tps.gantry_rot_axis = "x"
tps.sorted_spot_gneration = True
tps.particle = "proton"
To see more examples on the Treatment Plan source usage, the user can refer to test_059*.
Reference#
- class TreatmentPlanPBSource(*args, **kwargs)[source]#
Treatment Plan source Pencil Beam
User input parameters and default values:
activity:
Default value: 0
Description: Activity of the source in Bq (exclusive with ‘n’)
attached_to:
Default value: world
Description: Name of the volume to which the source is attached.
beam_data_dict:
Default value: None
Description: If a plan path is not provided, the source can be initialized by providing custom or plan derived spot data. Check opengate.contrib.tps.ionbeamtherapy.spots_info_from_txt() for more details on the structure of this dictionary.
beam_model:
Default value: None
Description: A BeamModel object instance, containing the geometrical and energy-dependent parameters of the beam model
beam_nr:
Default value: 1
Description: Which beam to simulate. Numbering starts from 1.
dynamic_params (set internally, i.e. read-only):
Default value: None
Description: Dictionary of dictionaries, where each dictionary specifies how the parameters of this object should evolve over time during the simulation. You cannot set this parameter directly. Instead, use the ‘add_dynamic_parametrisation()’ method of your object.If None, the object is static (default).
end_time:
Default value: None
Description: End time of the source
energies:
Default value: []
Description: List of the SPS’s energies. Calculated internally.
energy_sigmas:
Default value: []
Description: List of the SPS’s energy sigmas. Calculated internally.
flat_generation:
Default value: False
Description: If True, the same number of primaries is generated for each spot and the spot weight is applied to the energy deposition instead.
gantry_rot_axis:
Default value: z
Description: By default the source is oriented in +y direction. The source will be rotated of the gantry angle specified in the plan, around the specified axis.
half_life:
Default value: -1
Description: Half-life decay (-1 if no decay). Only when used with ‘activity’
ion:
Default value: {‘Z’: 0, ‘A’: 0, ‘E’: 0}
Description: If the particle is an ion, you must set Z: Atomic Number, A: Atomic Mass (nn + np +nlambda), E: Excitation energy (i.e. for metastable)
mother:
Deprecated: The user input parameter ‘mother’ is deprecated. Use ‘attached_to’ instead.
n:
Default value: 0
Description: Number of particle to generate (exclusive with ‘activity’)
name (must be provided):
Default value: None
partPhSp_xV:
Default value: []
Description: List of phase space parameters for each SPS. Calculated internally.
partPhSp_yV:
Default value: []
Description: List of phase space parameters for each SPS. Calculated internally.
particle:
Default value: None
Description: Name of the particle generated by the source (gamma, e+ … or an ion such as ‘ion 9 18’)
pdf:
Default value: []
Description: Probability distribution function. Calculated internally and used to sample the spots to irradiate.
plan_path:
Default value: None
Description: path of the treatment plan file to simulate. It can be in DICOM or Gate 9 .txt format
position:
Default value: {‘translation’: [0, 0, 0], ‘rotation’: array([[1., 0., 0.],, [0., 1., 0.],, [0., 0., 1.]])}
Description: translation and rotation to be applied to the source. Note that the transform will be applied to the source AFTER the gantry rotation and the source to axis distance have been applied. Note: rotation is not used yet.
positions:
Default value: []
Description: list of the positions to which the single particle source (SPS) will be moved during the simulation to irradiate all the spots. Calculated internally.
rotations:
Default value: []
Description: list of the rotations of each SPS.Calculated internally.
sorted_spot_generation:
Default value: False
Description: If True the source will generate primaries for each spot in the plan one by one, following the order of the plan. Otherwise the spot to fire is randomly sampled at each event from a probability distribution function.
start_time:
Default value: None
Description: Starting time of the source
weights:
Default value: []
Description: List of the weight for each spot. Calculated internally.