Calibration (TG-51/TRS-398)

Overview

The calibration module actually consists of two submodules: tg51 and trs398, each addressing their respective protocol. Both modules contain functions and classes for calculation the protocol dose. The modules have some overlap, especially with basic functions as well as helper functions. The modules have tried to use the vocabulary of the respective protocol, but occasionally there are differences when we felt that using the same name was clearer. See the vocabulary section for full definitions.

Note

Besides the typical calculations one would expect, the modules also include helper functions, such as a PDD to TPR converter so that a TG-51 calculation can determine kQ from TPR and avoid the tedious PDDx. Additionally, pressure unit converters exist to go from the various units of pressure to kPa which is what pylinac uses.

Vocabulary that may be different than the protocol

  • voltage_reference: Used in both TG-51 and TRS-398 for the voltage used when taking a reference reading; commonly -300V.

  • voltage_reduced: Use in both TG-51 and TRS-398 for the lower voltage used to determine ks/Pion; commonly -150V.

  • m_reference: The ion chamber reading at the reference voltage.

  • m_reduced: The ion chamber reading at the reduced voltage.

  • m_opposite: The ion chamber reading at the opposite polarity of the reference voltage: commonly +300V.

Vocabulary not listed here should be the same as the respective protocol.

Changing a bound

Bounds are placed in the module to prevent catastrophic errors from passing in the wrong values; e.g. the wrong units. If you live in a place that has extreme temperatures or pressures or just otherwise want to change the default bounds, you can change the default range of acceptable values. E.g. to change the minimum allowable temperature that can be passed:

from pylinac import tg51

tg51.p_tp(temp=5, press=100)  # will raise bounds error

# override
tg51.MIN_TEMP = 0

tg51.p_tp(temp=5, press=100)  # no bounds error will be raised

You can override the min/max of temp, pressure, p ion, p elec, p tp, p pol. These bounds are the same for TRS-398. I.e. setting these in either module will set them for both modules.

TG-51

Equation Definitions

Equation definitions are as follows:

  • Ptp (Temp/Pressure correction) - TG-51 Eqn. 10:

    \[\frac{273.2+T}{273.2+22} * \frac{101.33}{P}\]

    Warning

    Temperature is in Celsius and pressure is in kPa. Use the helper functions fahrenheit2celsius(), mmHg2kPa(), and mbar2kPa() as needed.

  • Ppol (Polarity correction) - Rather than using TG-51 Eqn. 9, we opt instead for TRS-398 Eqn xx, which uses the absolute values of the positive and negative voltages. This is the same results as Eqn. 9 but without worrying about signs.:

    \[\frac{|M^+_{raw}|+|M^-_{raw}|}{2*M_{raw}}\]
  • Pion (Ion collection correction; only for pulsed beams) - TG-51 Eqn. 12:

    \[\frac{1-\frac{V_{H}}{V_{L}}}{\frac{M^H_{raw}}{M^L_{raw}} - \frac{V_{H}}{V_{L}}}\]
  • Dref (Reference electron depth; cm) - TG-51 Eqn. 18:

    \[0.6*R_{50} - 0.1\]
  • R50 (Beam quality specifier; 50% dose depth; cm) - TG-51 Eqn. 16 & 17:

    \[\begin{split}\begin{cases} 1.029*I_{50}-0.06 (cm) & 2\leq I_{50}\leq 10 \\ 1.059*I_{50}-0.37 (cm) & I_{50}\gt 10 \\ \end{cases}\end{split}\]
  • k’R50 (k’R50 for cylindrical chambers) - TG-51 Eqn. 19:

    \[0.9905+0.0710e^\frac{-R_{50}}{3.67}\]
  • PQ_gr (PQ gradient correction for cylindrical chambers) - TG-51 Eqn. 21:

    \[\frac{M_{raw}(d_{ref}+0.5*r_{cav})}{M_{raw}*d_{ref}}\]
  • PDDx (PDD from photons only) - TG-51 Eqns. 13, 14 & 15:

    \[\begin{split}\begin{cases} PDD(10) & energy < 10 \\ PDD(10) & no lead, energy >= 10, PDD(10) <75 \\ 1.267*PDD(10)-20.0 & no lead, 75\leq PDD(10)\leq 89 \\ PDD(10)_{Pb} & lead@50cm, PDD(10)_{Pb} < 73 \\ (0.8905+0.00150*PDD(10)_{Pb})*PDD(10)_{Pb} & lead@50cm, PDD(10)_{Pb} \geq 73 \\ PDD(10)_{Pb} & lead@30cm, PDD(10)_{Pb} < 71 \\ (0.8116+0.00264*PDD(10)_{Pb})*PDD(10)_{Pb} & lead@30cm, PDD(10)_{Pb} \geq 71 \end{cases}\end{split}\]
  • M-corrected (corrected chamber reading) - TG-51 Eqn. 8:

    \[P_{ion}*P_{TP}*P_{elec}*P_{pol}*M_{raw}\]
  • kQ for Photons (cylindrical chamber-specific quality conversion factor) - TG-51 Addendum Eqn 1 & Table I:

    \[\begin{split}\begin{cases} A + B*10^-3*PDD(10)x+C*10^-5*(PDD(10)x)^2 & 63 < PDD(10)x < 86 \\ \end{cases}\end{split}\]

    Where A, B, and C are chamber-specific fitting values as given in Table I. Pylinac automatically retrieves values based on the chamber model passed to the function.

  • kQ for Electrons (cylindrical chamber-specific quality conversion factor) - Muir & Rogers 2014

    The study of Muir & Rogers was to find kecal values that could be determined soley from R50. Through Monte Carlo experiments, the optimal Pgradient was determined as well as fitting parameters for numerous common ion chambers. That study eliminates the need for Pgradient measurements. These kecal values will very likely be incorporated into the next TG-51 addendum (as has their kQ values for photons in the first addendum). From the paper, we can start with the known relationship given in Eqn. 9:

    \[k_{Q} = k_{Q,ecal} * k'_{Q}\]

    where Eqn. 11 states:

    \[k'_{Q} = a + b * R_{50}^{-c}\]

    Where a, b, and c are chamber-specific fitting values as given in Table VII and where \(k_{Q,ecal}\) is given in Table VI.

  • \(D^Q_{w}\) photon (Dose to water at 10cm from a photon beam of quality Q - TG-51 Eqn. 3:

    \[M*k_{Q}*N^{60Co}_{D,w} (Gy)\]
  • \(D^Q_{w}\) electron (Dose to water at 10cm from an electron beam of quality Q - TG-51 Eqn. 6:

    \[M*P^Q_{gr}*k'_{R_{50}}*k_{ecal}*N^{60Co}_{D,w} (Gy)\]

Function-based Use

Using the TG-51 module can be complementary to your existing workflow, or completely replace it. For example, you could use the kQ function to calculate kQ and then calculate the other corrections and values yourself. If you want something a little more complete, you can use the TG51Photon, TG51ElectronLegacy and TG51ElectronModern classes which will calculate all necessary corrections and values.

Note

The Photon class uses kQ values from the TG-51 addendum. The Legacy Electron class will make the user specify a kecal value and measure Pgradient. The Modern Electron class will calculate kQ completely from R50 and the chamber from Muir & Rogers 2014 paper, no kecal or Pgradient needed.

"""A script to calculate TG-51 dose using pylinac functions and following the TG-51 photon form"""
from pylinac.calibration import tg51

ENERGY = 6
TEMP = 22.1
PRESS = tg51.mmHg2kPa(755.0)
CHAMBER = "30013"  # PTW
P_ELEC = 1.000
ND_w = 5.443  # Gy/nC
MU = 200
CLINICAL_PDD = 66.5

# Section 4 (beam quality)
# since energy is 6MV, PDDx == PDD, but we'll run it through anyway just for show
pdd10x = tg51.pddx(pdd=66.4, energy=ENERGY)

# Section 5 (kQ)
kq = tg51.kq_photon_pddx(chamber=CHAMBER, pddx=pdd10x)
# Alternatively, get kQ from TPR (way quicker to measure, without needing to measure TPR!)
tpr = tg51.tpr2010_from_pdd2010(pdd2010=(38.0 / 66.4))
kq = tg51.kq_photon_tpr(chamber=CHAMBER, tpr=tpr)

# Section 6 (Temp/Press)
p_tp = tg51.p_tp(temp=TEMP, press=PRESS)

# Section 7 (polarity)
m_reference = (25.66, 25.67, 25.66)
m_opposite = (25.67, 25.67, 25.68)
p_pol = tg51.p_pol(m_reference=m_reference, m_opposite=m_opposite)

# Section 8 (ionization)
m_reduced = (25.61, 25.62)
p_ion = tg51.p_ion(
    voltage_reference=300,
    voltage_reduced=150,
    m_reference=m_reference,
    m_reduced=m_reduced,
)

# Section 9 (M corrected)
m_corr = tg51.m_corrected(
    p_ion=p_ion, p_tp=p_tp, p_elec=P_ELEC, p_pol=p_pol, m_reference=m_reference
)

# Section 10 (dose to water @ 10cm)
dose_10 = m_corr * kq * ND_w
dose_10_per_mu = dose_10 / MU

# Section 11 (dose/MU to water @ dmax)
dose_ddmax = dose_10_per_mu / CLINICAL_PDD

# Done!
print(dose_ddmax)

Class-based Use

"""A script to calculate TG-51 dose using pylinac classes and following the TG-51 photon form"""
from pylinac.calibration import tg51

ENERGY = 6
TEMP = 22.1
PRESS = tg51.mmHg2kPa(755.0)
CHAMBER = "30013"  # PTW
P_ELEC = 1.000
ND_w = 5.443  # Gy/nC
MU = 200
CLINICAL_PDD = 66.5

tg51_6x = tg51.TG51Photon(
    unit="TrueBeam1",
    chamber=CHAMBER,
    temp=TEMP,
    press=PRESS,
    n_dw=ND_w,
    p_elec=P_ELEC,
    measured_pdd10=66.4,
    lead_foil=None,
    clinical_pdd10=66.5,
    energy=ENERGY,
    voltage_reference=-300,
    voltage_reduced=-150,
    m_reference=(25.65, 25.66, 25.65),
    m_opposite=(25.64, 25.65, 25.65),
    m_reduced=(25.64, 25.63, 25.63),
    mu=MU,
    tissue_correction=1.0,
)

# Done!
print(tg51_6x.dose_mu_dmax)

# examine other parameters
print(tg51_6x.pddx)
print(tg51_6x.kq)
print(tg51_6x.p_ion)

# change readings if you adjust output
tg51_6x.m_reference_adjusted = (25.44, 25.44, 25.43)
# print new dose value
print(tg51_6x.dose_mu_dmax_adjusted)

# generate a PDF for record-keeping
tg51_6x.publish_pdf(
    "TB1 6MV TG-51.pdf",
    notes=["My notes", "I used Pylinac to do this; so easy!"],
    open_file=False,
)

TRS-398

Warning

Pylinac does not calculate electron dose in any other conditions than water; i.e. no solid water.

Equation Definitions

  • Ktp (Temp/Pressure correction):

    \[\frac{273.2+T}{273.2+22} * \frac{101.33}{P}\]

    Warning

    Temperature is in Celsius and pressure is in kPa. Use the helper functions fahrenheit2celsius, mmHg2kPa, and mbar2kPa as needed.

  • Kpol (Polarity correction):

    \[\frac{|M^+_{raw}|+|M^-_{raw}|}{2*M_{raw}}\]
  • Ks (Ion collection correction; only for pulsed beams):

    \[a_{0} + a_{1}*(\frac{M_{1}}{M_{2}}) + a_{2}*(\frac{M_{1}}{M_{2}})^2\]
  • Zref (Reference electron depth; cm) - TRS-398 7.2:

    \[0.6*R_{50} - 0.1\]
  • R50 (Beam quality specifier; 50% dose depth; cm) - TRS-398 7.1:

    \[\begin{split}\begin{cases} 1.029*I_{50}-0.06 (cm) & 2\leq I_{50}\leq 10 \\ 1.059*I_{50}-0.37 (cm) & I_{50}\gt 10 \\ \end{cases}\end{split}\]
  • \(D^Q_{w}\) photon (Dose to water at Zref from a photon or electron beam of quality Q - TRS-398 7.3:

    \[D_{w,Q} = M_{Q}*N_{D,w,Qo}*k_{Q,Qo} (Gy)\]
  • M-corrected (corrected chamber reading):

    \[M_{Q} = k_{s}*k_{TP}*K_{elec}*K_{pol}*M_{1}\]
  • kQ,Qo for Photons (cylindrical chamber-specific quality conversion factor): TRS-398 Table 6.III

  • kQ for Electrons (cylindrical chamber-specific quality conversion factor; calibrated in Co-60): TRS-398 Table 7.III

Function-based Use

"""A script to calculate TRS-398 dose using pylinac functions and following the TRS-398 photon form"""
from pylinac.calibration import trs398

TEMP = 22.1
PRESS = trs398.mmHg2kPa(755.0)
CHAMBER = "30013"  # PTW
K_ELEC = 1.000
ND_w = 5.443  # Gy/nC
MU = 200


# Section 3 (dosimeter corrections)
k_tp = trs398.k_tp(temp=TEMP, press=PRESS)
k_pol = trs398.k_pol(
    m_reference=(25.66, 25.67, 25.66), m_opposite=(25.65, 25.66, 25.66)
)
k_s = trs398.k_s(
    voltage_reference=300,
    voltage_reduced=150,
    m_reference=(25.66, 25.67, 25.66),
    m_reduced=(25.63, 25.65, 25.64),
)
m_corrected = (
    trs398.m_corrected(
        m_reference=(25.66, 25.67, 25.66),
        k_tp=k_tp,
        k_elec=K_ELEC,
        k_pol=k_pol,
        k_s=k_s,
    )
    / MU
)

# Section 4 (kQ + dose at zref)
kq = trs398.kq_photon(chamber=CHAMBER, tpr=(39.2 / 68.1))
dose_mu_zref = m_corrected * ND_w * kq

# Section 5 (Dose at zmax)
# SSD setup
CLINICAL_PDD = 66.5
dose_mu_zmax = dose_mu_zref * 100 / CLINICAL_PDD

# SAD setup
CLINICAL_TMR = 0.666
dose_mu_zmax = dose_mu_zref / CLINICAL_TMR

# Done!
print(dose_mu_zmax)

Class-based Use

"""A script to calculate TRS-398 dose using pylinac classes and following the TRS-398 photon form"""
from pylinac.calibration import trs398

ENERGY = 6
TEMP = 22.1
PRESS = trs398.mmHg2kPa(755.0)
CHAMBER = "30013"  # PTW
K_ELEC = 1.000
ND_w = 5.443  # Gy/nC
MU = 200
CLINICAL_PDD = 66.5

trs398_6x = trs398.TRS398Photon(
    unit="TrueBeam1",
    setup="SSD",
    chamber=CHAMBER,
    temp=TEMP,
    press=PRESS,
    n_dw=ND_w,
    clinical_pdd_zref=CLINICAL_PDD,
    tpr2010=(38.2 / 66.6),
    energy=ENERGY,
    fff=False,
    k_elec=K_ELEC,
    voltage_reference=-300,
    voltage_reduced=-150,
    m_reference=(25.65, 25.66, 25.65),
    m_opposite=(25.64, 25.65, 25.65),
    m_reduced=(25.64, 25.63, 25.63),
    mu=MU,
    tissue_correction=1.0,
)

# Done!
print(trs398_6x.dose_mu_zmax)

# examine other parameters
print(trs398_6x.kq)
print(trs398_6x.k_s)
print(trs398_6x.k_tp)

# change readings if you adjust output
trs398_6x.m_reference_adjusted = (25.44, 25.44, 25.43)
# print new dose value
print(trs398_6x.dose_mu_zmax_adjusted)

# generate a PDF for record-keeping
trs398_6x.publish_pdf(
    "TB1 6MV TRS-398.pdf",
    notes=["My notes", "I used Pylinac to do this; so easy!"],
    open_file=False,
)

TG-51 API Documentation

pylinac.calibration.tg51.mmHg2kPa(mmHg: float) float[source]

Utility function to convert from mmHg to kPa.

pylinac.calibration.tg51.mbar2kPa(mbar: float) float[source]

Utility function to convert from millibars to kPa.

pylinac.calibration.tg51.fahrenheit2celsius(f: float) float[source]

Utility function to convert from Fahrenheit to Celsius.

pylinac.calibration.tg51.tpr2010_from_pdd2010(*, pdd2010: float) float[source]

Calculate TPR20,10 from PDD20,10. From TRS-398 footnote 25, section 6.3.1, p.68 (https://www-pub.iaea.org/MTCD/Publications/PDF/TRS398_scr.pdf), and Followill et al 1998 eqn 1.

pylinac.calibration.tg51.p_tp(*, temp: float, press: float) float[source]

Calculate the temperature & pressure correction.

Parameters

tempfloat (17-27)

The temperature in degrees Celsius.

pressfloat (91-111)

The value of pressure in kPa. Can be converted from mmHg and mbar; see mmHg2kPa() and mbar2kPa().

pylinac.calibration.tg51.p_pol(*, m_reference: float | list | tuple | ndarray, m_opposite: float | list | tuple | ndarray) float[source]

Calculate the polarity correction.

Parameters

m_referencenumber, array

The readings of the ion chamber at the reference polarity and voltage.

m_oppositenumber, array

The readings of the ion chamber at the polarity opposite the reference. The sign does not make a difference.

Raises

BoundsError if calculated Ppol is >1% from 1.0.

pylinac.calibration.tg51.p_ion(*, voltage_reference: int, voltage_reduced: int, m_reference: float | list | tuple | ndarray, m_reduced: float | list | tuple | ndarray) float[source]

Calculate the ion chamber collection correction.

Parameters

voltage_referenceint

The “high” voltage; same as the TG51 measurement voltage.

voltage_reducedint

The “low” voltage; usually half of the high voltage.

m_referencefloat, iterable

The readings of the ion chamber at the “high” voltage.

m_reducedfloat, iterable

The readings of the ion chamber at the “low” voltage.

Raises

BoundsError if calculated Pion is outside the range 1.00-1.05.

pylinac.calibration.tg51.d_ref(*, i_50: float) float[source]

Calculate the dref of an electron beam based on the I50 depth.

Parameters

i_50float

The value of I50 in cm.

pylinac.calibration.tg51.r_50(*, i_50: float) float[source]

Calculate the R50 depth of an electron beam based on the I50 depth.

Parameters

i_50float

The value of I50 in cm.

pylinac.calibration.tg51.kp_r50(*, r_50: float) float[source]

Calculate k’R50 for Farmer-like chambers.

Parameters

r_50float (2-9)

The R50 value in cm.

pylinac.calibration.tg51.pq_gr(*, m_dref_plus: float | list | tuple | ndarray, m_dref: float | list | tuple | ndarray) float[source]

Calculate PQ_gradient for a cylindrical chamber.

Parameters

m_dref_plusfloat, iterable

The readings of the ion chamber at dref + 0.5rcav.

m_dreffloat, iterable

The readings of the ion chamber at dref.

pylinac.calibration.tg51.m_corrected(*, p_ion: float, p_tp: float, p_elec: float, p_pol: float, m_reference: float | list | tuple | ndarray) float[source]

Calculate M_corrected, the ion chamber reading with all corrections applied.

Parameters

p_ionfloat (1.00-1.05)

The ion collection correction.

p_tpfloat (0.92-1.08)

The temperature & pressure correction.

p_elecfloat (0.98-1.02)

The electrometer correction.

p_polfloat (0.98-1.02)

The polarity correction.

m_referencefloat, iterable

The raw ion chamber reading(s).

Returns

float

pylinac.calibration.tg51.pddx(*, pdd: float, energy: int, lead_foil: str | None = None) float[source]

Calculate PDDx based on the PDD.

Parameters

pdd{>62.7, <89.0}

The measured PDD. If lead foil was used, this assumes the pdd as measured with the lead in place.

energyint

The nominal energy in MV.

lead_foil{None, ‘30cm’, ‘50cm’}

Applicable only for energies >10MV. Whether a lead foil was used to acquire the pdd. Use None if no lead foil was used and the interim equation should be used. This is the default Use 50cm if the lead foil was set to 50cm from the phantom surface. Use 30cm if the lead foil was set to 30cm from the phantom surface.

pylinac.calibration.tg51.kq_photon_pddx(*, chamber: str, pddx: float) float[source]

Calculate kQ based on the chamber and clinical measurements of PDD(10)x. This will calculate kQ for photons for CYLINDRICAL chambers only.

Parameters

chamberstr

The chamber of the chamber. Valid values are those listed in Table III of Muir and Rogers and Table I of the TG-51 Addendum.

pddx{>63.0, <86.0}

The PHOTON-ONLY PDD measurement at 10cm depth for a 10x10cm2 field.

Note

Use the pddx() function to convert PDD to PDDx as needed.

Note

Muir and Rogers state limits of 0.627 - 0.861. The TG-51 addendum states them as 0.63 and 0.86. The TG-51 addendum limits are used here.

pylinac.calibration.tg51.kq_photon_tpr(*, chamber: str, tpr: float) float[source]

Calculate kQ based on the chamber and clinical measurements of TPR20,10. This will calculate kQ for photons for CYLINDRICAL chambers only.

Parameters

chamberstr

The chamber of the chamber. Valid values are those listed in Table III of Muir and Rogers and Table I of the TG-51 Addendum.

tpr{>0.630, <0.860}

The TPR(20,10) value.

Note

Use the tpr2010_from_pdd2010() function to convert from PDD without needing to take TPR measurements.

pylinac.calibration.tg51.kq_electron(*, chamber: str, r_50: float) float[source]

Calculate kQ based on the chamber and clinical measurements. This will calculate kQ for electrons for CYLINDRICAL chambers only according to Muir & Rogers.

Parameters

chamberstr

The chamber of the chamber. Valid values are those listed in Tables VI and VII of Muir and Rogers 2014.

r_50float

The R50 value in cm of an electron beam.

class pylinac.calibration.tg51.TG51Photon(*, institution: str = '', physicist: str = '', unit: str, measurement_date: str = '', temp: float, press: float, chamber: str, n_dw: float, p_elec: float, electrometer: str = '', measured_pdd10: float | None = None, lead_foil: str | None = None, clinical_pdd10: float, energy: int, fff: bool = False, voltage_reference: int, voltage_reduced: int, m_reference: float | list | tuple | ndarray, m_opposite: float | list | tuple | ndarray, m_reduced: float | list | tuple | ndarray, mu: int, tissue_correction: float = 1.0, m_reference_adjusted: float | list | tuple | ndarray | None = None)[source]

Bases: TG51Base

Class for calculating absolute dose to water using a cylindrical chamber in a photon beam.

Parameters

institutionstr

Institution name.

physiciststr

Physicist performing calibration.

unitstr

Unit name; e.g. TrueBeam1.

measurement_datestr

Date of measurement. E.g. 10/22/2018.

tempfloat

The temperature in Celsius. Use fahrenheit2celsius() to convert if necessary.

pressfloat

The value of pressure in kPa. Can be converted from mmHg and mbar; see mmHg2kPa() and mbar2kPa().

energyfloat

Nominal energy of the beam in MV.

chamberstr

Chamber model. Must be one of the listed chambers in TG-51 Addendum.

n_dwfloat

NDW value in Gy/nC.

p_elecfloat

Electrometer correction factor; given by the calibration laboratory.

measured_pdd10float

The measured value of PDD(10); will be converted to PDDx(10) and used for calculating kq.

lead_foil{None, ‘50cm’, ‘30cm’}

Whether a lead foil was used to acquire PDD(10)x and where its position was. Used to calculate kq.

clinical_pdd10float

The PDD used to correct the dose at 10cm back to dmax. Usually the TPS PDD(10) value.

voltage_referenceint

Reference voltage; i.e. voltage when taking the calibration measurement.

voltage_reducedint

Reduced voltage; usually half of the reference voltage.

m_referencefloat, tuple

Ion chamber reading(s) at the reference voltage.

m_oppositefloat, tuple

Ion chamber reading(s) at the opposite voltage of reference.

m_reducedfloat, tuple

Ion chamber reading(s) at the reduced voltage.

muint

The MU delivered to measure the reference reading. E.g. 200.

fffbool

Whether the beam is FFF or flat.

tissue_correctionfloat

Correction value to calibration to, e.g., muscle. A value of 1.0 means no correction (i.e. water).

property pddx: float

The photon-only PDD(10) value.

property kq: float

The chamber-specific beam quality correction factor.

property dose_mu_10: float

cGy/MU at a depth of 10cm.

property dose_mu_dmax: float

cGy/MU at a depth of dmax.

property dose_mu_10_adjusted: float

The dose/mu at 10cm depth after adjustment.

property dose_mu_dmax_adjusted: float

The dose/mu at dmax depth after adjustment.

publish_pdf(filename: str, notes: list | None = None, open_file: bool = False, metadata: dict | None = None)[source]

Publish (print) a PDF containing the analysis and quantitative results.

Parameters

filenamestr, file-like object

The file to write the results to.

notesstr, list

Any notes to be added to the report. If a string, adds everything as one line. If a list, must be a list of strings; each string item will be a new line.

open_filebool

Whether to open the file after creation. Will use the default PDF program.

metadatadict

Any data that should be appended to every page of the report. This differs from notes in that metadata is at the top of every page while notes is at the bottom of the report.

class pylinac.calibration.tg51.TG51ElectronLegacy(*, institution: str = '', physicist: str = '', unit: str = '', measurement_date: str = '', energy: int, temp: float, press: float, chamber: str, k_ecal: float, n_dw: float, electrometer: str = '', p_elec: float, clinical_pdd: float, voltage_reference: int, voltage_reduced: int, m_reference: float | list | tuple | ndarray, m_opposite: float | list | tuple | ndarray, m_reduced: float | list | tuple | ndarray, m_gradient: float | list | tuple | ndarray, cone: str, mu: int, i_50: float, tissue_correction: float = 1.0, m_reference_adjusted=None)[source]

Bases: TG51Base

Class for calculating absolute dose to water using a cylindrical chamber in an electron beam.

Parameters

institutionstr

Institution name.

physiciststr

Physicist performing calibration.

unitstr

Unit name; e.g. TrueBeam1.

measurement_datestr

Date of measurement. E.g. 10/22/2018.

tempfloat (17-27)

The temperature in degrees Celsius.

pressfloat (91-111)

The value of pressure in kPa. Can be converted from mmHg and mbar; see mmHg2kPa() and mbar2kPa().

chamberstr

Chamber model; only for bookkeeping.

n_dwfloat

NDW value in Gy/nC. Given by the calibration laboratory.

k_ecalfloat

Kecal value which is chamber specific. This value is the major difference between the legacy class and modern class where no kecal is needed.

p_elecfloat

Electrometer correction factor; given by the calibration laboratory.

clinical_pddfloat

The PDD used to correct the dose back to dref.

voltage_referencefloat

Reference voltage; i.e. voltage when taking the calibration measurement.

voltage_reducedfloat

Reduced voltage; usually half of the reference voltage.

m_referencefloat, tuple

Ion chamber reading(s) at the reference voltage.

m_oppositefloat, tuple

Ion chamber reading(s) at the opposite voltage of reference.

m_reducedfloat, tuple

Ion chamber reading(s) at the reduced voltage.

muint

The MU delivered to measure the reference reading. E.g. 200.

i_50float

Depth of 50% ionization.

tissue_correctionfloat

Correction value to calibration to, e.g., muscle. A value of 1.0 means no correction (i.e. water).

property r_50: float

Depth of the 50% dose value.

property dref: float

Depth of the reference point.

property pq_gr

Gradient factor

property kq: float

The kQ value using classic TG-51

property dose_mu_dref: float

cGy/MU at the depth of Dref.

property dose_mu_dmax: float

cGy/MU at the depth of dmax.

property dose_mu_dref_adjusted: float

cGy/MU at the depth of Dref.

property dose_mu_dmax_adjusted: float

cGy/MU at the depth of dmax.

publish_pdf(filename: str, notes: list | None = None, open_file: bool = False, metadata: dict | None = None)[source]

Publish (print) a PDF containing the analysis and quantitative results.

Parameters

filenamestr, file-like object

The file to write the results to.

notesstr, list

Any notes to be added to the report. If a string, adds everything as one line. If a list, must be a list of strings; each string item will be a new line.

open_filebool

Whether to open the file after creation. Will use the default PDF program.

metadatadict

Any data that should be appended to every page of the report. This differs from notes in that metadata is at the top of every page while notes is at the bottom of the report.

class pylinac.calibration.tg51.TG51ElectronModern(*, institution: str = '', physicist: str = '', unit: str = '', measurement_date: str = '', energy: int, temp: float, press: float, chamber: str, n_dw: float, electrometer: str = '', p_elec: float, clinical_pdd: float, voltage_reference: int, voltage_reduced: int, m_reference: float | list | tuple | ndarray, m_opposite: float | list | tuple | ndarray, m_reduced: float | list | tuple | ndarray, cone: str, mu: int, i_50: float, tissue_correction: float, m_reference_adjusted=None)[source]

Bases: TG51Base

Class for calculating absolute dose to water using a cylindrical chamber in an electron beam.

Warning

This class uses the values of Muir & Rogers. These values are likely to be included in the new TG-51 addendum, but are not official. The results can be up to 1% different. Physicists should use their own judgement when deciding which class to use. To use a manual kecal value, Pgradient and the classic TG-51 equations use the TG51ElectronLegacy class.

Parameters

institutionstr

Institution name.

physiciststr

Physicist performing calibration.

unitstr

Unit name; e.g. TrueBeam1.

measurement_datestr

Date of measurement. E.g. 10/22/2018.

pressfloat

The value of pressure in kPa. Can be converted from mmHg and mbar; see mmHg2kPa() and mbar2kPa().

tempfloat

The temperature in Celsius.

voltage_referenceint

The reference voltage; i.e. the voltage for the calibration reading (e.g. 300V).

voltage_reducedint

The reduced voltage, usually a fraction of the reference voltage (e.g. 150V).

m_referencearray, float

The reading(s) of the chamber at reference voltage.

m_reducedarray, float

The reading(s) of the chamber at the reduced voltage.

m_oppositearray, float

The reading(s) of the chamber at the opposite voltage from reference. Sign of the reading does not matter.

chamberstr

Ion chamber model.

n_dwfloat

NDW value in Gy/nC

p_elecfloat

Electrometer correction given by the calibration laboratory.

clinical_pddfloat

The PDD used to correct the dose back to dref.

muint

MU delivered.

i_50float

Depth of 50% ionization

tissue_correctionfloat

Correction value to calibration to, e.g., muscle. A value of 1.0 means no correction (i.e. water).

property r_50: float

Depth of the 50% dose value.

property dref: float

Depth of the reference point.

property kq: float

The kQ value using the updated Muir & Rogers values from their 2014 paper, equation 11, or classically if kecal is passed.

property dose_mu_dref: float

cGy/MU at the depth of Dref.

property dose_mu_dmax: float

cGy/MU at the depth of dmax.

property dose_mu_dref_adjusted: float

cGy/MU at the depth of Dref.

property dose_mu_dmax_adjusted: float

cGy/MU at the depth of dmax.

publish_pdf(filename: str, notes: list | None = None, open_file: bool = False, metadata: dict | None = None)[source]

Publish (print) a PDF containing the analysis and quantitative results.

Parameters

filenamestr, file-like object

The file to write the results to.

notesstr, list

Any notes to be added to the report. If a string, adds everything as one line. If a list, must be a list of strings; each string item will be a new line.

open_filebool

Whether to open the file after creation. Will use the default PDF program.

metadatadict

Any data that should be appended to every page of the report. This differs from notes in that metadata is at the top of every page while notes is at the bottom of the report.

TRS-398 API Documentation

pylinac.calibration.trs398.k_s(*, voltage_reference: int, voltage_reduced: int, m_reference: float | list | tuple | ndarray, m_reduced: float | list | tuple | ndarray) float[source]

Calculate the ion recombination effect using readings at two voltages. The voltages should have a ratio of 2, 2.5, 3, 3.5, 4, or 5.

Parameters

voltage_referenceint

The voltage at which calibration will be performed (e.g. 300V)

voltage_reducedint

The voltage which is lower than reference (e.g. 150V)

m_referencearray, float

The reading(s) at the reference voltage.

m_reducedarray, float

The reading(s) at the reduced voltage.

Returns

k_sfloat

The ion recombination factor.

Raises

ValueError

If the voltage ratio is not valid.

ValueError

If the calculated ks value is outside the range (1.0, 1.05).

pylinac.calibration.trs398.kq_photon(*, chamber: str, tpr: float) float[source]

Calculate the kQ factor for a photon beam given the chamber model and TPR20/10 using Table 6.III. Linear interpolation is used between given TPR ratios.

Parameters

chamberstr

Allowable chambers are those listed in Table 6.III that are also Farmer-type (e.g. Exradin A14 Farmer).

tprfloat

The ratio of measured TPR(20cm) / TPR(10cm). Note that this can also be calculated from PDD. See tpr2010_from_pdd2010().

Returns

kQfloat

The calculated kQ given table Table 6.III

Raises

KeyError

If the passed chamber is not within the acceptable list.

ValueError

If the TPR is not within the range defined by Table 6.III

pylinac.calibration.trs398.kq_electron(*, chamber: str, r_50: float) float[source]

Calculate the kQ factor for an electron beam given the chamber model and R50 using Table 7.III. Linear interpolation is used between given R50 values.

Parameters

chamberstr

The Farmer-type chambers listed in Table 7.III (e.g. PTW 30004/30012).

r_50float

The depth of R50 in cm in water.

Returns

kQfloat

The calculated kQ from Table 7.III

Raises

KeyError

If the passed chamber is not within the acceptable list.

ValueError

If the R50 is not within the range defined by Table 7.III

pylinac.calibration.trs398.m_corrected(*, m_reference, k_tp, k_elec, k_pol, k_s) float[source]

The fully corrected chamber reading.

Parameters

m_referencearray, float

The chamber reading(s) at the calibration position.

k_tpfloat

Temperature/Pressure correction. See p_tp().

k_elecfloat

Electrometer correction; given by the calibration laboratory.

k_polfloat

Polarity correction. See p_pol().

k_sfloat

Ion recombination correction. See k_s().

Returns

mfloat

The fully corrected chamber reading.

class pylinac.calibration.trs398.TRS398Photon(*, institution: str = '', physicist: str = '', unit: str = '', measurement_date: str = '', electrometer: str = '', setup: str, chamber: str, n_dw: float, mu: int, tpr2010: float, energy: int, fff: bool, press: float, temp: float, voltage_reference: int, voltage_reduced: int, m_reference: tuple | float, m_reduced: tuple | float, m_opposite: tuple | float, k_elec: float, clinical_pdd_zref: float | None = None, clinical_tmr_zref: float | None = None, tissue_correction: float = 1.0)[source]

Bases: TRS398Base

Calculation of dose to water at zmax and zref from a high energy photon beam. Setup can be SSD or SAD.

Parameters

setup{‘SSD’, ‘SAD’}

The physical setup of the calibration.

institutionstr

Institution name.

physiciststr

Physicist performing calibration.

unitstr

Unit name; e.g. TrueBeam1.

measurement_datestr

Date of measurement. E.g. 10/22/2018.

chamberstr

Farmer-type chamber model from Table 6.III.

n_dwfloat

NDw of the chamber given by the calibration laboratory.

mufloat, int

The number of MU given per reading

energyint

Nominal energy of the linac in MV; e.g. 6. Bookkeeping only.

fffbool

Whether the beam is FFF or flat. Bookkeeping only.

tpr2010float

The value of TPR(20)/TPR(10). Can be derived from PDD; see tpr2010_from_pdd2010().

pressfloat

The value of pressure in kPa. Can be converted from mmHg and mbar; see mmHg2kPa() and mbar2kPa().

tempfloat

The temperature in Celsius.

voltage_referenceint

The reference voltage; i.e. the voltage for the calibration reading (e.g. 300V).

voltage_reducedint

The reduced voltage, usually a fraction of the reference voltage (e.g. 150V).

m_referencearray, float

The reading(s) of the chamber at reference voltage.

m_reducedarray, float

The reading(s) of the chamber at the reduced voltage.

m_oppositearray, float

The reading(s) of the chamber at the opposite voltage from reference. Sign of the reading does not matter.

k_elecfloat

The electrometer correction value given by the calibration laboratory.

clinical_pdd_zrefoptional, float

The PDD at the depth of calibration. Use the actual percentage (e.g. 66.7 not 0.667). If not supplied the clinical_tmr_zref value must be supplied.

clinical_tmr_zrefoptional, float

The TMR at the depth of calibration. If not supplied the clinical_pdd_zref value must be supplied.

tissue_correctionfloat

The correction of calibration to a medium other than water. Default value is 1 which is water. E.g. use 0.99 if calibrating to muscle.

property kq

kQ of the chamber and TPR.

property dose_mu_zmax

cGy/MU at a depth of zmax.

property dose_mu_zmax_adjusted

The dose/mu at dmax depth after adjustment.

publish_pdf(filename: str, notes: list | None = None, open_file: bool = False, metadata: dict | None = None)[source]

Publish (print) a PDF containing the analysis and quantitative results.

Parameters

filenamestr, file-like object

The file to write the results to.

notesstr, list

Any notes to be added to the report. If a string, adds everything as one line. If a list, must be a list of strings; each string item will be a new line.

open_filebool

Whether to open the file after creation. Will use the default PDF program.

metadatadict

Any data that should be appended to every page of the report. This differs from notes in that metadata is at the top of every page while notes is at the bottom of the report.

class pylinac.calibration.trs398.TRS398Electron(*, institution: str = '', physicist: str = '', unit: str = '', measurement_date: str = '', electrometer: str = '', energy: str, cone: str, chamber: str, n_dw: float, mu: int, i_50: float, press: float, temp: float, voltage_reference: int, voltage_reduced: int, m_reference: tuple, m_reduced: tuple, m_opposite: tuple, k_elec: float, clinical_pdd_zref: float, tissue_correction: float = 1.0)[source]

Bases: TRS398Base

Calculation of dose to water at zmax and zref from a high energy electron beam.

Parameters

institutionstr

Institution name.

physiciststr

Physicist performing calibration.

unitstr

Unit name; e.g. TrueBeam1.

measurement_datestr

Date of measurement. E.g. 10/22/2018.

chamberstr

Farmer-type chamber model from Table 6.III.

n_dwfloat

NDw of the chamber given by the calibration laboratory.

mufloat, int

The number of MU given per reading.

i_50float

The depth of ionization 50% in cm.

pressfloat

The value of pressure in kPa. Can be converted from mmHg and mbar; see mmHg2kPa() and mbar2kPa().

tempfloat

The temperature in Celsius.

voltage_referenceint

The reference voltage; i.e. the voltage for the calibration reading (e.g. 300V).

voltage_reducedint

The reduced voltage, usually a fraction of the reference voltage (e.g. 150V).

m_referencearray, float

The reading(s) of the chamber at reference voltage.

m_reducedarray, float

The reading(s) of the chamber at the reduced voltage.

m_oppositearray, float

The reading(s) of the chamber at the opposite voltage from reference. Sign of the reading does not matter.

k_elecfloat

The electrometer correction value given by the calibration laboratory.

pdd_zrefoptional, float

The PDD at the depth of calibration. Use the actual percentage (e.g. 66.7 not 0.667). If not supplied the tmr_zref value should be supplied.

tissue_correctionfloat

The correction of calibration to a medium other than water. Default value is 1 which is water. E.g. use 0.99 if calibrating to muscle.

property r_50: float

The depth of R50 in cm, derived from I50.

property zref: float

Depth of the reference point.

property kq

kQ given the chamber and R50.

property dose_mu_zmax

cGy/MU at a depth of zmax.

property dose_mu_zmax_adjusted

The dose/mu at dmax depth after adjustment.

publish_pdf(filename: str, notes: list | None = None, open_file: bool = False, metadata: dict | None = None)[source]

Publish (print) a PDF containing the analysis and quantitative results.

Parameters

filenamestr, file-like object

The file to write the results to.

notesstr, list

Any notes to be added to the report. If a string, adds everything as one line. If a list, must be a list of strings; each string item will be a new line.

open_filebool

Whether to open the file after creation. Will use the default PDF program.

metadatadict

Any data that should be appended to every page of the report. This differs from notes in that metadata is at the top of every page while notes is at the bottom of the report.