Defining Magnetic Field#
The user can define either an uniform magnetic field or a custom one.
using Geant4
using Geant4.SystemOfUnits
using Geant4.PhysicalConstants
Uniform magnetic field#
The functionG4JLUniformMagField
creates a uniform magnetic field proving a B vector. For example:
bfield = G4JLUniformMagField(G4ThreeVector(0,0, 1.5tesla))
G4JLUniformMagField("UnifiormB", Geant4.G4JLUniformMagFieldData(G4ThreeVector(0.0,0.0,0.0015)), Geant4.var"#getfield!#33"(), G4JLMagField[])
The G4 toolkit will call the function bfield.getfield(field::G4ThreeVector, pos::G4ThreeVector, data::G4JLFieldData)
each time that is needed. The method getfield(bfield::G4JLMagneticField, pos::G4ThreeVector)
does it as well.
# Will be moved to Geant4.jl
function getfield(bfield::G4JLMagneticField, pos::G4ThreeVector)
B = G4ThreeVector()
bfield.getfield(B, pos, bfield.data)
return B
end
getfield (generic function with 1 method)
@show getfield(bfield, G4ThreeVector(0,0,0))/tesla
getfield(bfield, G4ThreeVector(0, 0, 0)) / tesla = G4ThreeVector(0.0,0.0,1.5)
G4ThreeVector(0.0,0.0,1.5)
Custom magnetic Field#
To define a custom field we need:
define first a user structure for the parameters that inherits from the abstract type
G4JLFieldData
then, define a function with the signature
(result::G4ThreeVector, position::G4ThreeVector, params::G4JLFieldData)::Nothing
and finally, with all this, instantiate the magnetic field calling the function
G4JLMagneticField(<name>, <data>; getfield_method=<function>)
using Geant4.SystemOfUnits: ampere
using Parameters
@with_kw mutable struct WireFieldData <: G4JLFieldData
I::Float64 = 1ampere
wiredir::G4ThreeVector = G4ThreeVector(0,0,1)
end
function getfield!(field::G4ThreeVector, pos::G4ThreeVector, data::WireFieldData)::Nothing
r = cross(data.wiredir, pos)
B = (mu0 * data.I)/(2π*mag2(r)) * r
assign(field, B)
return
end
bfield = G4JLMagneticField("WireField", WireFieldData(); getfield_method=getfield!);
@show getfield(bfield, G4ThreeVector(1,0,0))/tesla
@show getfield(bfield, G4ThreeVector(2,0,0))/tesla
@show getfield(bfield, G4ThreeVector(3,0,0))/tesla
getfield(bfield, G4ThreeVector(1, 0, 0)) / tesla = G4ThreeVector(0.0,0.0002,0.0)
getfield(bfield, G4ThreeVector(2, 0, 0)) / tesla = G4ThreeVector(0.0,0.0001,0.0)
getfield(bfield, G4ThreeVector(3, 0, 0)) / tesla = G4ThreeVector(0.0,6.666666666666667e-5,0.0)
G4ThreeVector(0.0,6.666666666666667e-5,0.0)
using Plots
# Range of distances from the wire
r_values = range(0.01m, stop=1m, length=100) # Distances from 0.01m to 1m
# Calculate magnetic field strengths corresponding to each distance
B_values = [ mag(getfield(bfield, G4ThreeVector(r,0,0)))/tesla for r in r_values]
# Plot
plot(r_values, B_values, xlabel="Distance (mm)", ylabel="Magnetic Field (T)", label="", legend=:bottomright)
title!("Magnetic Field vs. Distance from Wire")