API Reference

Exported Functions

SPHtoGrid.analytic_synchrotronFunction
analytic_synchrotron(P_cgs::Array{<:Real}, B_cgs::Array{<:Real}, 
                     Mach::Array{<:Real}, θ_B::Union{Nothing,Array{<:Real}}=nothing;
                     dsa_model::Union{Integer,AbstractShockAccelerationEfficiency}=1, 
                     ν0::Real=1.4e9,
                     K_ep::Real=0.01, CR_Emin::Real=1.0,
                     spectrum::Union{Nothing,Function}=nothing,
                     integrate_pitch_angle::Bool=true,
                     polarisation::Bool=false,
                     show_progress::Bool=false)

Computes the analytic synchrotron emission from a spectrum of electrons by explicitly integrating over the distribution function. The integral over the spectrum must be normalized to 1. The total energy density of the relativistic electrons is given by the CR to thermal pressure ratio obtained by employing a DSA model and computing Xcr as in Pfrommer et. al. (2017).

Returns synchrotron emissivity j_nu in units [erg/s/Hz/cm^3].

Arguments

  • P_cgs::Array{<:Real}: Thermal energy density in erg/cm^3.
  • B_cgs::Array{<:Real}: Magnetic field in Gauss.
  • Mach::Array{<:Real}: Mach number.
  • θ_B::Union{Nothing,Array{<:Real}}=nothing: Shock obliquity (optional).

Keyword Arguments

  • ν0::Real=1.4e9: Observation frequency in Hz.
  • dsa_model: Diffusive Shock Acceleration model. Takes values 0...4, or custom model. See next section.
  • K_ep::Real=0.01: Ratio of CR proton to electron energy acceleration.
  • CR_Emin::Real=1: Injection energy of CR electron population in GeV.
  • spectrum::Union{Nothing,Function}=nothing: Spectrum function. Must be normalized so that the integral over it is 1.
  • integrate_pitch_angle::Bool=true: Optional avoid pitch angle integration to reduce computational cost.
  • polarisation::Bool=false: Set to true if you want to compute the polarized emission instead of the total intensity.
  • show_progress::Bool=false: Enables a progress bar if set to true.

DSA Models

Takes either your self-defined AbstractShockAccelerationEfficiency (see DiffusiveShockAccelerationModels.jl for details!) or a numerical value as input. Numerical values correspond to:

Mapping settings

source
SPHtoGrid.analytic_synchrotron_GSFunction
analytic_synchrotron_GS( rho_cgs::Array{<:Real}, B_cgs::Array{<:Real},
                         T_K::Array{<:Real}, Mach::Array{<:Real};
                         xH::Real=0.76, dsa_model::Integer=1, ν0::Real=1.44e9,
                         integrate_pitch_angle::Bool=true )

Computes the analytic synchrotron emission with the simplified approach described in Ginzburg & Syrovatskii 1965, "Cosmic Magnetobremsstrahlung". Uses the implementaion from Donnert et. al. (2016).

Returns synchrotron emissivity j_nu in units [erg/s/Hzcm^3].

Arguments

  • rho_cgs::Array{<:Real}: Density in $g/cm^3$.
  • B_cgs::Array{<:Real}: Magnetic field in Gauss.
  • T_K::Array{<:Real}: Temperature in Kelvin.
  • Mach::Array{<:Real}: Mach number.
  • θ_B::Union{Nothing,Array{<:Real}}=nothing: Shock obliquity (optional).

Keyword Arguments

  • xH::Float64 = 0.76: Hydrogen fraction of the simulation, if run without chemical model.
  • ν0::Real=1.44e9: Observation frequency in $Hz$.
  • dsa_model::Integer=1: Diffusive Shock Acceleration model. Takes values 0...4, see next section.
  • K_ep::Real=0.01: Ratio of CR proton to electron energy acceleration.
  • show_progress::Bool=false: Enables a progress bar if set to true

DSA Models

Takes either your self-defined AbstractShockAccelerationEfficiency (see DiffusiveShockAccelerationModels.jl for details!) or a numerical value as input. Numerical values correspond to:

Mapping settings

source
SPHtoGrid.analytic_synchrotron_HB07Function
analytic_synchrotron_HB07( rho_cgs::Array{<:Real}, m_cgs::Array{<:Real}, hsml_cgs::Array{<:Real},
                           B_cgs::Array{<:Real}, T_keV::Array{<:Real}, Mach::Array{<:Real},
                           θ_B::Union{Nothing, Array{<:Real}}=nothing;
                           xH::Real=0.752, ν0::Real=1.4e9, z::Real=0.0,
                           dsa_model::Union{Nothing,Integer,AbstractShockAccelerationEfficiency}=nothing,
                           ξe::Real=1.e-5,
                           show_progress::Bool=false )

Computes the analytic synchrotron emission with the simplified approach described in Hoeft&Brüggen (2007), following approach by Wittor et. al. (2017).

Returns synchrotron emissivity j_nu in units [erg/s/Hz/cm^3].

Arguments

  • rho_cgs::Array{<:Real}: Density in $g/cm^3$.
  • m_cgs::Array{<:Real}: Particle mass in $g$.
  • hsml_cgs::Array{<:Real}: HSML block in $cm$.
  • B_cgs::Array{<:Real}: Magnetic field in $G$.
  • T_keV::Array{<:Real}: Temperature in $keV$.
  • Mach::Array{<:Real}: Sonic Mach number.
  • θ_B::Union{Nothing,Array{<:Real}}=nothing: Shock obliquity (optional).

Keyword Arguments

  • xH::Float64 = 0.76: Hydrogen fraction of the simulation, if run without chemical model.
  • ν0::Real=1.44e9: Observation frequency in $Hz$.
  • z::Real=0.0: Redshift of the simulation.
  • dsa_model=nothing: Diffusive Shock Acceleration model. If set to a value overwrites the default Hoeft&Brüggen acceleration model. See next section.
  • ξe::Real=1.e-5: Ratio of CR proton to electron energy acceleration. Given as a fraction of thermal gas, essenitally Xcr * Kep. Default value from Nuza+2017. For dsa_model != nothing use something like ξe = 1.e-4.
  • show_progress::Bool=false: Enables a progress bar if set to true.

DSA Models

Takes either your self-defined AbstractShockAccelerationEfficiency (see DiffusiveShockAccelerationModels.jl for details!) or a numerical value as input. Numerical values correspond to:

Mapping settings

source
SPHtoGrid.analytic_synchrotron_LongairFunction
analytic_synchrotron_emission( rho_cgs::Array{<:Real}, B_cgs::Array{<:Real},
                               T_K::Array{<:Real}, Mach::Array{<:Real};
                               xH::Real=0.76, dsa_model::Integer=1, ν0::Real=1.44e9,
                               integrate_pitch_angle::Bool=true )

Computes the analytic synchrotron emission with the simplified approach described in Longair Eq. 8.128. Returns J_ν in units [erg/cm^3/Hz/s].

Arguments

  • rho_cgs::Array{<:Real}: Density in $g/cm^3$.
  • B_cgs::Array{<:Real}: Magnetic field in Gauss.
  • T_K::Array{<:Real}: Temperature in Kelvin.
  • Mach::Array{<:Real}: Mach number.

Keyword Arguments

  • xH::Float64 = 0.76: Hydrogen fraction of the simulation, if run without chemical model.
  • dsa_model::Integer=1: Diffuse-Shock-Acceleration model. Takes values 0...4, see next section.
  • ν0::Real=1.44e9: Observation frequency in $Hz$.
  • K_ep::Real=0.01: Ratio of CR proton to electron energy density.
  • integrate_pitch_angle::Bool=true: Integrates over the pitch angle as in Longair Eq. 8.87.
  • convert_to_mJy::Bool=false: Convert the result from $[erg/cm^3/Hz/s]$ to $mJy/cm$.

DSA Models

Takes either your self-defined AbstractShockAccelerationEfficiency (see DiffusiveShockAccelerationModels.jl for details!) or a numerical value as input. Numerical values correspond to:

source
SPHtoGrid.beam_in_kpcMethod
beam_in_kpc(θ_beam::Vector{Union{Real, Unitful.AbstractQuantity}}, 
            c::Cosmology.AbstractCosmology, z::Real)

Converts the beam from arcmin to kpc.

source
SPHtoGrid.comptonYMethod
comptonY(n_cm3::Real, T_K::Real, z::Real)

Computes the Compton-Y parameter from electron density n_cm3 and temperature T in Kelvin at redshift z.

Arguments:

  • n_cm3: SPH particle density in [1/cm^3]
  • T_K: SPH particle temperature [K]
  • z: Redshift.

Mapping settings

source
SPHtoGrid.comptonYMethod
comptonY(n_cm3::Vector{<:Real}, T_K::Vector{<:Real}, z::Real)

Computes the Compton-Y parameter from electron density n_cm3 and temperature T in Kelvin at redshift z.

Arguments:

  • n_cm3: SPH particle density in [1/cm^3]
  • T_K: SPH particle temperature [K]
  • z: Redshift

Mapping settings

source
SPHtoGrid.convert_Pnu_map_to_mJy_beamMethod
convert_Pnu_map_to_mJy_beam(map::Matrix{<:Real}, 
                            d_pixel::Real,
                            beam::Union{Real, Unitful.AbstractQuantity}, 
                            c::Cosmology.AbstractCosmology, 
                            z::Real)

Converts a map from units [W / Hz] to [mJy / beam] for a circular beam.

Parameters:

  • map: original map in [W / Hz].
  • d_pixel: size of a pixel in $kpc$.
  • beam: radius of the beam in $arcmin$.
  • c: Cosmology used for conversion.
  • z: Redshift of the image.
source
SPHtoGrid.convert_Pnu_map_to_mJy_beamMethod
convert_Pnu_map_to_mJy_beam(map::Matrix{<:Real}, 
                            d_pixel::Real,
                            beam::Vector{Union{Real, Unitful.AbstractQuantity}}, 
                            c::Cosmology.AbstractCosmology, 
                            z::Real)

Converts a map from units [W / Hz] to [mJy / beam].

Parameters:

  • map: original map in [W / Hz].
  • d_pixel: size of a pixel in $kpc$.
  • beam: dimensions of the beam in $arcmin$.
  • c: Cosmology used for conversion.
  • z: Redshift of the image.
source
SPHtoGrid.convert_Pnu_map_to_mJy_beamMethod
convert_Pnu_map_to_mJy_beam(map::Matrix{<:Real}, 
                            d_pixel::Real,
                            beam::Union{T, Vector{T}}, 
                            h::AbstractGadgetHeader) where T::Union{Real, Unitful.AbstractQuantity}

Converts a map from units [W / Hz] to [mJy / beam] by using a Gadget header.

Parameters:

  • map: original map in [W / Hz].
  • d_pixel: size of a pixel in $kpc$.
  • beam: radius/dimensions of the beam in $arcmin$.
  • h: Gadget header of simulation
source
SPHtoGrid.density_2DFunction
density_2D( rho::Real, pixelSideLength::Real, 
            Mass::Real=1.989e43, Length::Real=3.085678e21)

Computes the density in units of [g/cm^2].

Arguments:

  • rho: SPH particle density in physical code units.
  • pixelSideLength: length of pixel in physical units.
  • Mass: Mass unit in [g].
  • Length: Length unit in [cm]

Mapping settings

source
SPHtoGrid.distributed_allsky_mapMethod
distributed_allsky_map( allsky_filename::String, 
                        Nside::Integer, Nsubfiles::Integer, 
                        mapping_function::Function;
                        reduce_image::Bool=true)

Dynamically dispatches workers to compute one allsky map per subfile, sum up the results and save to a fits file.

Arguments

  • allsky_filename::String: Name of the file under which the image should be saved
  • Nside::Integer: Nside for healpix map, must be a power of 2! Nside = 2^N.
  • Nsubfiles::Integer: Number of subfiles the snapshot is distributed over.
  • mapping_function::Function: The function to be executed per subfile. Must have a call to allsky_map as return value.
  • reduce_image: If the final image should be divided by the weight image set to true
source
SPHtoGrid.distributed_cic_mapFunction
distributed_cic_map(cic_filename::String, Nsubfiles::Integer,
                    mapping_function::Function,
                    Nimages::Integer, param::mappingParameters;
                    reduce_image::Bool=true,
                    Ndim::Integer=2,
                    snap::Integer=0,
                    units::String="",
                    vtk::Bool=true)

Dynamically dispatches workers to compute one cic map per subfile, sum up the results and save to a fits file.

Arguments

  • cic_filename::String: Name of the file under which the image should be saved
  • Nsubfiles::Integer: Number of subfiles the snapshot is distributed over.
  • mapping_function::Function: The function to be executed per subfile. Must return a Tuple of the cic image(s) and the weight_image.
  • Nimages::Integer=1: Number of cic image that are constructed simultaniously. Only useful if you map multiple quantities at once, so it defaults to 1.
  • param::mappingParameters: Parameters for the mapping, needed to save the image in the usual way.
  • reduce_image: If the final image should be divided by the weight image set to true.
  • Ndim::Integer=2: Dimensions of image, must be either 2 or 3.
  • snap::Integer=0: Number of the input snapshot, used for saving the image.
  • units::String="": Units of the image, used for saving the image.
source
SPHtoGrid.gamma_flux_pions_PE04Method
gamma_flux_pions_PE04(rho_cgs::Real, m_cgs::Real, T_K::Real, α_p::Real, d::Real;
                      Xcr::Real=0.5,
                      Eγ_π0_min::Real=0.1, Eγ_π0_max::Real=200.0,
                      xH::Real=0.752)

Flux of γ-ray photons produced from a proton spectrum given as a fraction Xcr of the energy density defined by rho_cgs [g/cm^3] and T_K [K], with a powerlaw slope in energy α_p. Integrated over SPH particle volume for particle of mass m_cgs [g]. Flux from a distance d [cm]. Integrated between photon energies Eγ_π0_min and Eγ_π0_max [GeV]. Returns total number of photons in energy band in untis of [γ cm^-2 s^-1]. See Pfrommer&Enßlin (2004), Eq. 25.

Arguments:

  • rho_cgs: SPH particle density in [g/cm^3]
  • m_cgs: SPH particle mass in [g]
  • T_K: SPH particle temperature [K]
  • α_p: Slope of proton energy spectrum S ~ 2.0 - 2.5
  • d: Distance to SPH particle or halo [cm].
  • Xcr: CR proton to thermal pressure ratio.
  • Eγ_π0_min: Minimum photon energy for γ-ray spectrum [GeV]
  • Eγ_π0_max: Maximum photon energy for γ-ray spectrum [GeV]
  • xH: Hydrogen mass fraction in the simulation

Mapping settings

source
SPHtoGrid.gamma_luminosity_pions_PE04Method
gamma_luminosity_pions_PE04(rho_cgs::Real, m_cgs::Real, T_K::Real, α_p::Real;
                            Xcr::Real=0.5,
                            Eγ_π0_min::Real=0.1, Eγ_π0_max::Real=200.0,
                            xH::Real=0.752)

γ-ray luminosity produced from a proton spectrum given as a fraction Xcr of the energy density defined by rho_cgs [g/cm^3] and T_K [K], with a powerlaw slope in energy α_p. Integrated over SPH particle volume for particle of mass m_cgs [g]. Integrated between photon energies Eγ_π0_min and Eγ_π0_max [GeV]. Returns total luminosity integrated over energy band in units of [GeV s^-1]. See Pfrommer&Enßlin (2004), Eq. 25.

Arguments:

  • rho_cgs: SPH particle density in [g/cm^3]
  • m_cgs: SPH particle mass in [g]
  • T_K: SPH particle temperature [K]
  • α_p: Slope of proton energy spectrum S ~ 2.0 - 2.5
  • Xcr: CR proton to thermal pressure ratio.
  • Eγ_π0_min: Minimum photon energy for γ-ray spectrum [GeV]
  • Eγ_π0_max: Maximum photon energy for γ-ray spectrum [GeV]
  • xH: Hydrogen mass fraction in the simulation

Mapping settings

source
SPHtoGrid.healpix_mapMethod
healpix_map(Pos, Hsml, M, Rho, Bin_q, Weights;
            center::Vector{<:Real}=[0.0, 0.0, 0.0],
            radius_limits::Vector{<:Real}=[0.0, Inf],
            Nside::Integer=1024,
            kernel::AbstractSPHKernel,
            calc_mean::Bool=true
            show_progress::Bool=true,
            output_from_all_workers::Bool=false)

Calculate an allsky map from SPH particles. Returns two HealpixMaps: (image, weight_image). To reduce the image afterwards divide image by weight_image.

Arguments:

  • Pos: Positions of particles in physical code units
  • HSML: hsml of particles in physical code units
  • M: Mass of particles in physical code units
  • Bin_Q: Quantitiy for binning in arb. units
  • Weights: Weights for map
  • center: Position from which projection is looking outwards in physical code units
  • radius_limits: Inner and outer radius of the shell that should be mapped
  • Nside: Nside for healpix map, has to be a power of 2
  • calc_mean: Calculate the mean along the line of sight. If set to false only particles with Bin_Q > 0 contribute to the map.
  • show_progress: Print a progress bar
  • output_from_all_workers: Allow output from multiple workers. If set to false only the main process prints a progress bar.
source
SPHtoGrid.jγ_PE04Method
jγ_PE04(rho_cgs::Real, T_K::Real, α_p::Real, Eγ::Real; 
        Xcr::Real=0.5, xH::Real=0.752)

Gamma-ray emissivity at photon energy [GeV] for thermal gas with properties rho_cgs [g/cm^3] and T_K [K]. Sets up a CR proton spectrum with energy slope α_p as a fraction Xcr of the thermal energy density. Returns emissivity in units [GeV cm^-3 s^-1 ]. See Pfrommer&Enßlin (2004), Eq. 19.

Function Arguments:

  • rho_cgs: SPH particle density in [g/cm^3]
  • T_K: SPH particle temperature [K]
  • α_p: Slope of proton energy spectrum S ~ 2.0 - 2.5
  • : Photon energy [GeV]
  • Xcr: CR proton to thermal pressure ratio.
  • xH: Hydrogen mass fraction in the simulation

Mapping settings

For mean value along line-of-sight:

  • weights: rho (weight with density)
  • reduce_image: true

For integral along line-of-sight, aka surface brightness:

source
SPHtoGrid.kinetic_SZFunction
kinetic_SZ(n_cm3::Real, vel_y_cgs::Real, 
           z::Real=0.0, ν::Real=1.e9;
           DI_over_I::Bool=false)

Computes the kinetic Sunyaev-Zel'dovich effect from electron density n_cm3 and velocity in y-direction to the projection plane in cgs units vel_y_cgs. If DI_over_I is set to true you also need to provide an observation frequency ν and redshift z.

Arguments:

  • n_cm3: SPH particle density in [1/cm^3]
  • vel_y_cgs: SPH particle velocity in y-direction in [cm/s]
  • z: Redshift
  • ν: Observing frequency

Mapping settings

source
SPHtoGrid.map_itFunction
map_it(pos_in, hsml, mass, rho, bin_q, weights, RM=nothing;
       param::mappingParameters,
       kernel::AbstractSPHKernel=WendlandC6(2), 
       snap::Integer=0, 
       units::AbstractString="", 
       image_prefix::String="dummy",
       reduce_image::Bool=true, 
       parallel=true,
       calc_mean::Bool=true, show_progress::Bool=true,
       sort_z::Bool=false,
       stokes::Bool=false,
       projection="xy")

Small helper function to copy positions, map particles and save the fits file.

Arguments

  • pos_in: Matrix (3xNpart) with particle positions
  • hsml: Array with particle hsml
  • mass: Array with particle masses
  • rho: Array with particle densities
  • bin_q: Array with particle quantity to be mapped
  • weights: Array with weights. Defaults to density-weighted
  • param: mappingParameters for this map
  • kernel: AbstractSPHKernel to be used for mapping
  • units: Unit string will be saved to FITS file
  • image_prefix: Name of the file to save, withou .fits file-ending
  • reduce_image: If weights need to be applied or not. Set to false for part_weight_physical
  • parallel: Run on multiple processors.
  • calc_mean: Calculates the mean value along the line of sight. If set to false the particle only contributes if its bin_q is larger than 0.
  • show_progress: Show progress bar
  • sort_z: Sort the particles according to their line-of-sight direction. Needed for polarisation mapping.
  • stokes: Set to true of you are mapping Stokes parameters to account for Faraday rotation of the polarisation angle.
  • projection: Which plane the position should be rotated in. Can also be an Array of 3 Euler angles (in [°]) (not used yet!)
source
SPHtoGrid.mass_densityMethod
mass_density(pos::Matrix{<:Real}, mass::Vector{<:Real}; 
             kernel::AbstractSPHKernel=Cubic(Float64, 3),
             Nneighbors::Integer=32,
             boxsize::Vector{<:Real}=zeros(3),
             verbose::Bool=false)

Compute the density of an arbitrary particle distribution using the SPH method. Neighbor searches are performed using a (periodic) BallTree. The density is computed in input units, so additional unit conversion to cgs units is required, if input units are not cgs.

Arguments:

  • pos: particle position in physical code units.
  • mass: particle mass in physical code units.

Keyword Arguments

  • kernel::AbstractSPHKernel=Cubic(Float64, 3): SPH kernel to use for the density estimate. Works with any kernel from SPHKernels.jl.
  • Nneighbors::Integer=32: Number of neighboring particles to use for the density estimate.
  • boxsize::Vector{<:Real}=zeros(3): Boxsize in each dimension. Used for periodic boundary conditions. If set to zero, non-periodic boundary conditions are assumed.
  • verbose::Bool=false: If set to true gives progress reports and progress bar.

Returns

  • rho: mass density at particle position in input units.
  • hsml: smoothing length of each particle in input units.

Mapping settings

source
SPHtoGrid.part_weight_XrayBandFunction
part_weight_XrayBand(T_K::Array{<:Real}, Emin::Real, Emax::Real)

Computes Xray weighted emission of a defined energy band. Emin and Emax are energies in eV.

source
SPHtoGrid.part_weight_emissionMethod
part_weight_emission(rho::Array{<:Real}, T_K::Array{<:Real})

Emission weighted mapping. Takes density in internal untis and temperature in K and computes weights.

source
SPHtoGrid.part_weight_spectroscopicMethod
part_weight_spectroscopic(rho::Array{<:Real}, T_K::Array{<:Real})

Spectroscopic weighted mapping from Mazotta+ 04. Takes density and temperature and computes weights.

source
SPHtoGrid.polarisation_angleFunction
polarisation_angle(Q_image, U_image, Iν_image=nothing, Iν_cutoff= 0.0)

Compute the polarisation fraction image from Stokes Q and U images.

source
SPHtoGrid.project_along_axisFunction
project_along_axis!(x::Array{<:Real}, projection_axis::Integer=3)

Projects and array of 3D along one of the principle axes. projection_axis ∈ {1, 2, 3} => x, y, z axis.

source
SPHtoGrid.project_along_axis!Function
project_along_axis!(x::Array{<:Real}, projection_axis::Integer=3)

Projects and array of 3D along one of the principle axes. projection_axis ∈ {1, 2, 3} => x, y, z axis.

source
SPHtoGrid.project_along_axis!Function
project_along_axis!(x::Array{<:Real}, projection_axis::Integer=3)

Projects and array of 3D along one of the principle axes. projection_axis ∈ {1, 2, 3} => x, y, z axis.

source
SPHtoGrid.read_allsky_fits_imageMethod
read_allsky_fits_image(filename::String)

Read a FITS file containing an allsky image and returns the image, snapshot number and units.

Returns

  • image: A 2D array with the image pixels
  • snap: Number of the mapped snapshot
  • units: A unit string of the image

Example

image, snapnr, unitstring = readallskyfits_image(filename)

source
SPHtoGrid.read_fits_imageFunction
read_fits_image(filename::String)

Read a FITS file and return the image, mappingParameters and the snapshot number.

Returns

  • image: A 2D array with the image pixels
  • par: mappingParameters used for the image
  • snap: Number of the mapped snapshot
  • units: A unit string of the image

Example

image, par, snapnr, unitstring = readfitsimage(filename)

source
SPHtoGrid.read_smac1_binary_infoMethod
read_smac1_binary_info(filename::String)

Returns the image info in a Smac1ImageInfo struct.

Smac1ImageInfo fields

  • snap::Int32: number of input snapshot
  • z::Float32: redshift of snapshot
  • m_vir::Float32: virial mass of halo
  • r_vir::Float32: virial radius of halo
  • xcm::Float32: x coordinate of image center
  • ycm::Float32: y coordinate of image center
  • zcm::Float32: z coordinate of image center
  • z_slice_kpc::Float32: depth of the image in kpc
  • boxsize_kpc::Float32: xy-size of the image in kpc
  • boxsize_pix::Float32: xy-size of the image in pixels
  • pixsize_kpc::Float32: size of one pixel in kpc
  • xlim::Array{Float64,1}: x limits of image
  • ylim::Array{Float64,1}: y limits of image
  • zlim::Array{Float64,1}: z limits of image
  • units::String: unitstring of image
source
SPHtoGrid.rotate_3D!Method
rotate_3D!(x::Array{<:Real}, alpha::Real, beta::Real, gamma::Real)

Rotates and array of 3D positions around the euler angles α, β and γ corresponding to rotations around the x, y, and z-axis respectively. α, β and γ need to be given in degrees.

source
SPHtoGrid.rotate_3DMethod
rotate_3D(x::Array{<:Real}, alpha::Real, beta::Real, gamma::Real)

Rotates and array of 3D positions around the euler angles α, β and γ corresponding to rotations around the x, y, and z-axis respectively. α, β and γ need to be given in degrees.

source
SPHtoGrid.rotation_measureMethod
rotation_measure(n_cm3::Real, B_los::Real, dz::Real, ν_obs::Real)

Computes the rotation measure of the parallel magnetic field along the LOS at a given frequency. To be used for continuous rotation of polarized emission along the LOS.

Arguments

  • n_cm3::Real: Electron number density in [cm^-3].
  • B_los::Real: Magnetic field strength along the LOS in [G].
  • dz::Real: Depth along the LOS in [cm]. See get_dz for a convenient helper function.
  • ν_obs::Real: Observing frequency in [Hz].

Returns

  • RM::Real: Rotation measure in [rad/cm^2].

Mapping settings

  • weight function: part_weight_physical
  • reduce image: false
  • stokes: true
  • sort_z: true
  • RM: use output of this function.
source
SPHtoGrid.rotation_measureMethod
rotation_measure(n_cm3::Real, B_los::Real, dz::Real; ν_obs=nothing)

Computes the rotation measure of the parallel magnetic field along the LOS.

Arguments

  • n_cm3::Real: Electron number density in [cm^-3].
  • B_los::Real: Magnetic field strength along the LOS in [G].

Returns

  • RM::Real: Rotation measure in [rad/m^2].

Mapping settings

source
SPHtoGrid.sphMappingFunction
sphMapping( Pos::Array{<:Real}, HSML::Array{<:Real}, M::Array{<:Real}, 
            Rho::Array{<:Real}, Bin_Quant::Array{<:Real}, 
            Weights::Array{<:Real}=Rho;
            param::mappingParameters,
            kernel::AbstractSPHKernel,
            show_progress::Bool=true,
            parallel::Bool=false,
            reduce_image::Bool=true,
            return_both_maps::Bool=false,
            filter_particles::Bool=true,
            dimensions::Int=2,
            calc_mean::Bool=false)

Maps the data in Bin_Quant to a grid. Parameters of mapping are supplied in param and the kernel to be used in kernel.

Arguments

  • Pos: Matrix (3xNpart) with particle positions
  • HSML: Array with particle hsml
  • M: Array with particle masses
  • Rho: Array with particle densities
  • Bin_Quant: Array with particle quantity to be mapped
  • Weights: Array with weights. Defaults to density-weighted
  • param: mappingParameters for this map
  • kernel: AbstractSPHKernel to be used for mapping
  • show_progress: Show progress bar
  • parallel: Run on multiple processors
  • reduce_image: If weights need to be applied or not. Set to false for part_weight_physical
  • return_both_maps: Returns the full image array. To be used with parallel mapping of subfiles
  • filter_particles: Find the particles that are actually contained in the image
  • dimensions: Number of mapping dimensions (2 = to grid, 3 = to cube)
  • calc_mean: Calculates the mean value along the line of sight. If set to false the particle only contributes if its Bin_Quant is larger than 0
source
SPHtoGrid.sphMappingMethod
sphMapping( Pos, HSML, Bin_Quant;
            param::mappingParameters,
            kernel::AbstractSPHKernel,
            show_progress::Bool=true,
            parallel::Bool=false,
            reduce_image::Bool=true,
            filter_particles::Bool=true,
            dimensions::Int=2)

Maps the data in Bin_Quant to a grid. Parameters of mapping are supplied in param and the kernel to be used in kernel.

Arguments

  • Pos: Matrix (3xNpart) with particle positions.
  • HSML: Array with particle hsml.
  • Bin_Quant: Array with particle quantity to be mapped.
  • kernel::AbstractSPHKernel: Kernel object to be used.
  • show_progress::Bool=true: Show progress bar.
  • parallel::Bool=true: Run on multiple processors.
  • reduce_image::Bool=true: If weights need to be applied or not. Set to false for part_weight_one and part_weight_physical.
  • filter_particles::Bool=true: Find the particles that are actually contained in the image.
  • dimensions::Int=2: Number of mapping dimensions (2 = to grid, 3 = to cube).
source
SPHtoGrid.surface_brightness_to_luminosityMethod
surface_brightness_to_luminosity(map::Matrix{<:Real}, pixelSideLength::Real; unit_factor::Real=1.0)

Converts a map of surface brightness to luminosity per pixel. Uses pixelSideLength as the diameter of a pixel in [kpc]. If unit_factor is provided it is multiplied to every pixel to perform unit conversion.

source
SPHtoGrid.surface_brightness_to_luminosityMethod
surface_brightness_to_luminosity(map::Matrix{<:Real}, par::mappingParameters; unit_factor::Real=1.0)

Converts a map of surface brightness to luminosity per pixel. If unit_factor is provided it is multiplied to every pixel to perform unit conversion.

source
SPHtoGrid.synchrotron_SB_to_luminosityMethod
synchrotron_SB_to_luminosity(map, pixelSideLength::Real)

Converts a map of synchrotron surface brightness [erg/s/Hz/cm^2] to synchrotron luminosity [W/Hz]. Uses pixelSideLength as the diameter of a pixel in [kpc].

source
SPHtoGrid.synchrotron_SB_to_luminosityMethod
synchrotron_SB_to_luminosity(map, par::mappingParameters)

Converts a map of synchrotron surface brightness [erg/s/Hz/cm^2] to synchrotron luminosity [W/Hz]. Uses par as the mappingParameters of the original map.

source
SPHtoGrid.thermal_SZFunction
thermal_SZ( n_cm3::Vector{<:Real}, T_K::Vector{<:Real},
            z::Real=0.0, ν::Real=1.44e9; 
            DI_over_I::Bool=false )

Computes the thermal Sunyaev-Zel'dovich effect for electron density n_cm3 and temperature T_K in Kelvin at redshift z and observer frequency ν. DI_over_I outputs in units of $dI/I$ if set to true and dT/T otherwise.

Arguments:

  • n_cm3: SPH particle density in [1/cm^3]
  • T_K: SPH particle temperature [K]
  • z: Redshift
  • ν: Observing frequency

Mapping settings

source
SPHtoGrid.total_synch_luminosity_from_SBMethod
total_synch_luminosity_from_SB(map::Matrix{<:Real}, pixelSideLength::Real)

Computes the total synchrotron luminosity in [W/Hz] from a map of synchrotron surface brightness in [erg/s/Hz/cm^2].

source
SPHtoGrid.total_synch_luminosity_from_SBMethod
total_synch_luminosity_from_SB(map::Matrix{<:Real}, par::mappingParameters)

Computes the total synchrotron luminosity in [W/Hz] from a map of synchrotron surface brightness in [erg/s/Hz/cm^2].

source
SPHtoGrid.write_fits_imageMethod
write_fits_image(filename::String, image::Array{<:Real}, 
                        par::mappingParameters; 
                        units::String="[i.u.]",
                        snap::Integer=0)

Writes a mapped image to a FITS file and stored the essential mapping parameters in the header.

source
SPHtoGrid.write_fits_imageMethod
write_fits_image(filename::String, image::Array{<:Real};
                units::String = "[i.u.]",
                snap::Integer = 0)

Writes a mapped allsky image to a FITS file and stored units and snapshot number in the header.

source
SPHtoGrid.write_smac1_parFunction
write_smac1_par( path="./"; kwargs...)

Writes a Smac parameter file. Add keyword arguments according to the smac1 parameter file.

source
SPHtoGrid.write_smac2_parFunction
write_smac2_par(x, y, z,
                euler_angle_0, euler_angle_1, euler_angle_2,
                xy_size, z_depth, xy_pix::Integer,
                input_file, output_file, path,
                effect_module::Integer=0, effect_flag::Integer=0,
                ν_obs::Real=1.44e9, cosmology::Integer=0,
                CR_pmin::Real=10.0, CR_pmax::Real=1.e7)

Writes a P-Smac2 parameter file. Not all relevant parameters implemented yet!

source
SPHtoGrid.write_vtk_imageMethod
write_vtk_image(filename::String, image::Array{<:Real},
                image_name::String, par::mappingParameters;
                units::String = "[i.u.]", snap::Integer = 0)

Writes a mapped image to a vti file.

source
SPHtoGrid.x_ray_emissivityFunction
x_ray_emissivity(T_keV::Vector{<:Real}, 
                 rho_cgs::Vector{<:Real},
                 metalicity::Union{Vector{Float64, Nothing}}=nothing; 
                 E0::Real=0.1, E1::Real=2.4, 
                 xH::Real=0.752,
                 cooling_function::Bool=false,
                 z::Real=0.0)

X-Ray emissivity for particles with temperature T_keV in $keV$, and density rho_cgs in $g/cm^3$. If available you can also add the metalicity in the gas. Emin and Emax give the minimum and maximum energy of the oberservation. xH gives the hydrogen fraction used in the simulation.

Returns

X-Ray emissivity in units of [erg/s/cm^3].

Arguments:

  • T_keV: SPH particle temperature [keV]
  • m_cgs: SPH particle mass in [g]
  • rho_cgs: SPH particle density in [g/cm^3]
  • E0: Minimum photon energy for Xray spectrum [keV]
  • E1: Maximum photon energy for Xray spectrum [keV]
  • xH: Hydrogen mass fraction in the simulation

Mapping settings

source
SPHtoGrid.λγ_PE04Method
λγ_PE04(rho_cgs::Real, T_K::Real, α_p::Real; 
        Xcr::Real=0.5,
        Eγ_π0_min::Real=0.1, Eγ_π0_max::Real=200.0,
        xH::Real=0.752)

Number of γ-ray photons produced per time and volume from a proton spectrum given as a fraction Xcr of the energy density defined by rho_cgs [g/cm^3] and T_K [K], with a powerlaw slope in energy α_p. Integrated between photon energies Eγ_π0_min and Eγ_π0_max [GeV]. Returns number of photons in energy band in units of [γ cm^-3 s^-1]. See Pfrommer&Enßlin (2004), Eq. 25.

Function Arguments:

  • rho_cgs: SPH particle density in [g/cm^3]
  • T_K: SPH particle temperature [K]
  • α_p: Slope of proton energy spectrum S ~ 2.0 - 2.5
  • Xcr: CR proton to thermal pressure ratio.
  • Eγ_π0_min: Minimum photon energy for γ-ray spectrum [GeV]
  • Eγ_π0_max: Maximum photon energy for γ-ray spectrum [GeV]
  • xH: Hydrogen mass fraction in the simulation

Mapping settings

For mean value along line-of-sight:

  • weights: rho (weight with density)
  • reduce_image: true

For integral along line-of-sight, aka surface brightness:

source

Exported Types

SPHtoGrid.mappingParametersType
mappingParameters( T::DataType=Float64;
                   x_lim::Vector{Float64}   = [-1.0, -1.0],
                   y_lim::Vector{Float64}   = [-1.0, -1.0],
                   z_lim::Vector{Float64}   = [-1.0, -1.0],
                   center::Vector{Float64}  = [-1.0, -1.0, -1.0],
                   x_size::Float64          =  -1.0,
                   y_size::Float64          =  -1.0,
                   z_size::Float64          =  -1.0,
                   pixelSideLength::Float64 =  -1.0,
                   Npixels::Int64           =   0)

Parameter object for sph to grid mapping. Define either *_lim, or center and *_size. Resolution is defined by pixelSideLength or Npixels.

source

Private Functions

SPHtoGrid.B_cmbMethod
B_cmb(z::Real)

Magnetic field equivalent of the CMB for a given redshift z in [G].

source
SPHtoGrid.X_crMethod
X_cr(Mach::Real, ξ::Real)

Ratio between CR and thermal pressure for an injected CR spectrum following Pfrommer et. al. (2017). Uses the assumption that the CRs do not significantly alter the shock properties, which is not valid in the case of efficient CR proton injection!

source
SPHtoGrid.add_subtr_boxMethod
add_subtr_box(pos::T, boxsize::T) where T

Helper function to add or subtract the boxsize from the given position

source
SPHtoGrid.calculate_indexMethod
function calculate_index(i::Integer, j::Integer, x_pixels::Integer)

Calculates the index of a flattened 3D image array.

source
SPHtoGrid.calculate_indexMethod
function calculate_index(i::Integer, j::Integer, x_pixels::Integer)

Calculates the index of a flattened 2D image array.

source
SPHtoGrid.calculate_weightsMethod
function calculate_weights_2D(  wk::Array{<:Real,1}, 
                                iMin::Integer, iMax::Integer, 
                                jMin::Integer, jMax::Integer,
                                x::Real, y::Real, hsml::Real, hsml_inv::Real,
                                kernel::AbstractSPHKernel,
                                x_pixels::Integer )

Calculates the kernel- and geometric weights of the pixels a particle contributes to.

source
SPHtoGrid.calculate_weightsMethod
function calculate_weights_3D(  wk::Array{<:Real,1}, 
                                iMin::Integer, iMax::Integer, 
                                jMin::Integer, jMax::Integer,
                                kMin::Integer, kMax::Integer,
                                x::Real, y::Real, z::Real, 
                                hsml::Real, hsml_inv::Real,
                                kernel::AbstractSPHKernel,
                                x_pixels::Integer, y_pixels::Integer )

Calculates the kernel- and geometric weights of the pixels a particle contributes to.

source
SPHtoGrid.calculate_weightsMethod
calculate_weights(wk::Vector{Float64}, A::Vector{Float64}, 
                _pos::Vector{Float64}, _hsml::Float64,
                _Δx::Float64, res::Resolution,
                pixidx::Vector{Int64}, ang_pix::Float64,
                kernel::AbstractSPHKernel)
source
SPHtoGrid.calculate_weights_splashMethod
function calculate_weights_2D(  wk::Array{<:Real,1}, 
                                iMin::Integer, iMax::Integer, 
                                jMin::Integer, jMax::Integer,
                                x::Real, y::Real, hsml::Real, hsml_inv::Real,
                                kernel::AbstractSPHKernel,
                                x_pixels::Integer )

Calculates the kernel- and geometric weights of the pixels a particle contributes to.

source
SPHtoGrid.calculate_weights_splashMethod
function calculate_weights_3D(  wk::Array{<:Real,1}, 
                                iMin::Integer, iMax::Integer, 
                                jMin::Integer, jMax::Integer,
                                kMin::Integer, kMax::Integer,
                                x::Real, y::Real, z::Real, 
                                hsml::Real, hsml_inv::Real,
                                kernel::AbstractSPHKernel,
                                x_pixels::Integer, y_pixels::Integer )

Calculates the kernel- and geometric weights of the pixels a particle contributes to.

source
SPHtoGrid.cic_mapping_2DFunction

cicmapping2D( Pos, HSML, M, Rho, BinQ, Weights; param::mappingParameters, kernel::AbstractSPHKernel, showprogress::Bool=false, calc_mean::Bool=true )

Underlying function to map SPH data to a 2D grid.

source
SPHtoGrid.contributing_pixelsMethod
contributing_pixels(pos::Vector{T}, radius::T, 
                    res::Resolution, allsky_map) where {T}

Computes the pixels the particle contributes to.

source
SPHtoGrid.cre_spec_norm_particleMethod
cre_spec_norm_particle(M::Real, η_model::AbstractShockAccelerationEfficiency)

Computes the CR electron norm of the particles. This depends on the Mach number M and the acceleration efficiency given by η_model.

source
SPHtoGrid.distance_to_pixel_centerMethod
distance_to_pixel_center(r::T, pos::Vector{T}, pixel_center::Vector{T}) where T<:Real

Compute the distance between particle vector and pixel center in radians.

source
SPHtoGrid.faraday_rotate_pixel!Method
faraday_rotate_pixel!(image::Array{Float64}, idx::Integer, 
                      pRM::Float64, pix_weight::Float64)

Faraday rotates the current pixel state based on the contribution from particle p.

source
SPHtoGrid.filter_sort_particlesMethod
filter_sort_particles(Pos, Hsml, Bin_q, Weights, center, radius_limits)

Filters the particles that are in the shell that should be mapped and sorts them by according to their position along the line of sight. Returns arrays with only relevant particles

source
SPHtoGrid.gamma_source_PE04Method
gamma_source_PE04(Eγ::T, f_p::T, q_p::T, nH::T) where {T}

Source function of gamma-ray photons at energy in units of N_photons erg^-1 s^-1 cm^-3 as given in Pfrommer&Enßlin (2004), Eq.19.

source
SPHtoGrid.get_d_hsmlMethod
get_d_hsml_3D( dx::Real, dy::Real, dz::Real, hsml_inv::Real )

Computes the distance in 3D to the pixel center in units of the kernel support.

source
SPHtoGrid.get_d_hsmlMethod
get_d_hsml( dx::Real, dy::Real, hsml_inv::Real )

Computes the distance in 2D to the pixel center in units of the kernel support.

source
SPHtoGrid.get_dxyzMethod
function get_dxyz(x::Real, hsml::Real, i::Integer)

Calculates the extent of the current particle size in units of pixels.

source
SPHtoGrid.get_dzMethod
get_dz(M, rho, hsml)

Helper function to get the particle depth along the LOS in the given units. All additional unit conversion must be performed by hand.

source
SPHtoGrid.get_number_densitiesMethod
get_number_densities(rho_cgs::Real, T_K::Real, α_p::Real, Xcr::Real, xH::Real)

Thermal proton number density and CR proton normalisation.

source
SPHtoGrid.get_quantities_2DMethod
get_quantities_2D( pos, weight, hsml, 
                   rho, m, len2pix::T) where T

Helper function to convert quantities to pixel units and the correct data type.

source
SPHtoGrid.get_quantities_splashMethod
get_quantities_splash( pos, weight, hsml, 
                   rho, m, len2pix::T) where T

Helper function to convert quantities to pixel units and the correct data type.

source
SPHtoGrid.get_weight_per_pixelMethod
get_weight_per_pixel(distr_area, distr_weight,
                            n_tot_pix, n_distr_pix,
                            dxdy, u, hsml_inv,
                            kernel)

Computes the area/volume contribution for each pixel and evaluates the kernel at the pixel center.

source
SPHtoGrid.get_xyzMethod
function get_xyz( pos, hsml, k::Integer,
                  par::mappingParameters)

Calculates x, y, z position in units of pixels and performs periodic mapping, if required.

source
SPHtoGrid.integrate_θFunction
integrate_θ(x_in::Real, synch_F::Function, θ_steps::Integer=50)

Pitch angle integration in Donnert+16, Eq. 17. synch_F can be either the first or second synchrotron function.

source
SPHtoGrid.kSzPrefacMethod
kSzPrefac(ν::Real, z::Real, DI_over_I::Bool)

Prefactor for the kinetic Sunyaev-Zel'dovich effect.

source
SPHtoGrid.particle_area_and_depthMethod
particle_area_and_depth(d, hsml, pix_radian)

Computes the particle area and depth. Caution: This has to be represented as a cylinder instead of a sphere, which introduces an error by design. Also computes the pixel radius at the particle horizon.

source
SPHtoGrid.particle_in_imageMethod
particle_in_image(x::Real, y::Real, z::Real, hsml::Real,
                            halfXsize::Real, halfYsize::Real, halfZsize::Real)

Checks if a periodically mapped particle is still in the image frame.

source
SPHtoGrid.particle_in_imageMethod
particle_in_image( pos::Array{T}, hsml::T,
                   halfsize::Array{Float64}) where T

Checks if a periodically mapped particle is still in the image frame.

source
SPHtoGrid.pix_index_min_maxMethod
function pix_index_min_max(x::T, hsml::T, n_pixels::Int64) where T

Calculates the minimum and maximum pixel to which a particle contributes.

source
SPHtoGrid.qγ_PE04Method
qγ_PE04(rho_cgs, T_K, α_p, Eγ; 
        Xcr::Real=0.5, 
        xH=0.752)

Gamma-ray sorce function at photon energy for thermal gas with properties rho_cgs [g/cm^3] and T_K [K]. Sets up a CR proton spectrum with energy slope α_p as a fraction Xcr of the thermal energy density. Returns source function in units [γ cm^-3 s^-1 GeV^-1]. See Pfrommer&Enßlin (2004), Eq. 19.

source
SPHtoGrid.reduce_image_2DMethod
function reduce_image_2D( image::Array{<:Real},
                          x_pixels::Int64, y_pixels::Int64)

Unflattens an image array to a 2D array of pixels.

source
SPHtoGrid.reduce_image_2DMethod
function reduce_image_2D( image::Array{<:Real},
                          x_pixels::Int64, y_pixels::Int64)

Unflattens an image array to a 2D array of pixels.

source
SPHtoGrid.reduce_image_3DMethod
function reduce_image_3D( image::Array{<:Real}, w_image::Array{<:Real},
                                        x_pixels::Int64, y_pixels::Int64, z_pixels::Int64)

Unflattens an image array to a 3D array of pixels.

source
SPHtoGrid.reduce_image_3DMethod
function reduce_image_3D( image::Array{<:Real}, w_image::Array{<:Real},
                                        x_pixels::Int64, y_pixels::Int64, z_pixels::Int64)

Unflattens an image array to a 3D array of pixels.

source
SPHtoGrid.select_dsa_modelMethod
select_dsa_model(dsa_model::AbstractShockAccelerationEfficiency)

Helper function to assign the self-defined DSA model.

source
SPHtoGrid.splash_mapping_2DMethod

splashmapping2D( Pos, HSML, M, Rho, BinQ; param::mappingParameters, kernel::AbstractSPHKernel, showprogress::Bool=false )

Underlying function to map SPH data to a 2D grid.

source
SPHtoGrid.splash_mapping_3DMethod

splashmapping3D( Pos, HSML, BinQ; param::mappingParameters, kernel::AbstractSPHKernel, showprogress::Bool=false )

Underlying function to map SPH data to a 3D grid.

source
SPHtoGrid.tSzPrefacMethod
tSzPrefac(ν::Real, z::Real, DI_over_I::Bool)

Computes the prefactor for the thermal Sunyaev-Zel'dovich effect.

source
SPHtoGrid.update_image!Method
update_image!(allsky_map, weight_map, wk,
            pixidx,
            area, dz,
            A, N, weight_per_pix,
            quantity_weight, bin_quantity)

Updates the images.

source
SPHtoGrid.update_image!Method
function update_image( image::Real, w_image::Real, 
                       wk::Real, bin_q::Real, 
                       geometry_norm::Real )

Applies the different contributions to the image and the weight image.

source
SPHtoGrid.weight_per_indexMethod
weight_per_index(_Δx, _pos, _hsml, _hsml_inv,
                _pixidx,
                distr_area, distr_weight, 
                n_tot_pix, n_distr_pix, 
                res, ang_pix,
                kernel)

Helper function to compute the weights for one pixel

source
SPHtoGrid.ys_fMethod
ys_f(xs::Real, ξ::Real)

Helper function to get ys (see Pfrommer+2017)

source
SPHtoGrid.ñ_crpFunction
ñ_crp(rho_cgs::Real, T_K::Real, α_p::Real, Xcr::Real=0.01; 
        xH::Real=0.76)

CR proton normalisation in [1/cm^3]

source
SPHtoGrid.α_γMethod
α_γ(slope::T) where T

Slope of γ-ray spectrum as a function of the slope of a proton spectrum in energy space.

source
SPHtoGrid.σ_ppMethod
σ_pp(α::T) where T

Scattering cross-section of proton-proton scatterin in $cm^2$

source

Private Types