API

Coordinates

FewBodyPhysics.ParticleSystemType
ParticleSystem(masses::Vector{Float64})

A data structure representing a system of particles, storing their masses and associated Jacobi transformation matrices for coordinate transformations.

Fields

  • masses::Vector{Float64}: A vector containing the masses of the particles.
  • J::Matrix{Float64}: The Jacobi transformation matrix, used to convert particle coordinates into Jacobi coordinates.
  • U::Matrix{Float64}: The pseudoinverse of the Jacobi transformation matrix J, used to transform Jacobi coordinates back to the particle coordinate system.

Constructor

  • ParticleSystem(masses::Vector{Float64}): Constructs a new ParticleSystem instance.
    • Arguments:
      • masses: A vector of particle masses. At least two masses are required.
source
FewBodyPhysics.jacobi_transformFunction
jacobi_transform(masses::Vector{Float64})::Tuple{Matrix{Float64}, Matrix{Float64}}

Compute the Jacobi transformation matrices J and U for a system of particles with specified masses.

Arguments

  • masses::Vector{Float64}: A vector of masses for the particles.

Returns

  • (J::Matrix{Float64}, U::Matrix{Float64}): The Jacobi transformation matrix and its pseudoinverse.

Notes

  • The matrices J and U are used to transform between particle coordinates and Jacobi coordinates.
  • The pseudoinverse U is used instead of the inverse to handle cases where J is not square.
source
FewBodyPhysics.generate_A_matrixFunction
generate_A_matrix(bij::Vector{Float64}, w_list::Vector{Vector{Float64}})::Matrix{Float64}

Generate the matrix A for Gaussian basis functions given width parameters bij and weight vectors w_list.

Arguments

  • bij::Vector{Float64}: A vector of width parameters for the Gaussian basis functions.
  • w_list::Vector{Vector{Float64}}: A list of weight vectors.

Returns

  • A::Matrix{Float64}: The sum of weighted outer products of w_list, scaled by bij.

Notes

  • This function constructs the A matrix used in the correlated Gaussian method.
source
FewBodyPhysics.transform_listFunction
transform_list(α::Vector{Float64})::Vector{Matrix{Float64}}

Transform a list of scalar values α into a list of 1x1 matrices.

Arguments

  • α::Vector{Float64}: A list of scalar values.

Returns

  • Array{Matrix{Float64}}: A list of 1x1 matrices where each matrix contains one of the scalar values from α.
source
FewBodyPhysics.shift_vectorsFunction
shift_vectors(a::Matrix{Float64}, b::Matrix{Float64}, mat::Union{Nothing, Matrix{Float64}}=nothing)::Float64

Calculate the weighted sum of the element-wise product of vectors a and b using matrix mat.

Arguments

  • a::Matrix{Float64}: A matrix where each column is a vector a_i.
  • b::Matrix{Float64}: A matrix where each column is a vector b_j.
  • mat::Union{Nothing, Matrix{Float64}}: An optional matrix to weight the product (default is the identity matrix).

Returns

  • sum_val::Float64: The weighted sum of products.

Notes

  • The matrices a and b should have the same dimensions.
source
FewBodyPhysics.generate_weight_vectorFunction
generate_weight_vector(dim::Int, i::Int, j::Int)::Vector{Int}

Generate a weight vector for the i-th and j-th coordinates in a space of dimension dim.

Arguments

  • dim::Int: The dimension of the space.
  • i::Int: The index for the positive element in the weight vector.
  • j::Int: The index for the negative element in the weight vector.

Returns

  • w::Vector{Int}: A vector with 1 at the i-th position, -1 at the j-th position, and 0 elsewhere.
source
FewBodyPhysics.transform_coordinatesFunction
transform_coordinates(J::Matrix{Float64}, r::Vector{Float64})::Vector{Float64}

Transform the coordinates r of a system using the Jacobi matrix J.

Arguments

  • J::Matrix{Float64}: The Jacobi transformation matrix.
  • r::Vector{Float64}: The coordinates to be transformed.

Returns

  • x::Vector{Float64}: The transformed coordinates.
source
FewBodyPhysics.inverse_transform_coordinatesFunction
inverse_transform_coordinates(U::Matrix{Float64}, x::Vector{Float64})::Vector{Float64}

Transform the coordinates x back to the original system using the inverse of the Jacobi matrix.

Arguments

  • U::Matrix{Float64}: The inverse Jacobi transformation matrix.
  • x::Vector{Float64}: The coordinates in Jacobi space.

Returns

  • r::Vector{Float64}: The coordinates transformed back to the original system.
source

Matrix elements

FewBodyPhysics.S_elementsFunction
S_elements(A::Matrix{Float64}, B::Matrix{Float64}, K::Matrix{Float64}, w::Union{Nothing, Vector{Vector{Float64}}}=nothing)

Calculate matrix elements for overlap, kinetic energy, and optionally the Coulomb term.

Arguments

  • A: Width matrix of Gaussian basis functions for state i.
  • B: Width matrix of Gaussian basis functions for state j.
  • K: Kinetic energy matrix.
  • w (optional): List of weight vectors for the particles involved.

Returns

  • M0: The overlap matrix element between the two states.
  • trace: The trace used in the kinetic energy calculation.
  • Coulomb_term (optional): The Coulomb interaction term, if weight vectors w are provided.
source
FewBodyPhysics.S_waveFunction
S_wave(α::Vector{Matrix{Float64}}, K::Matrix{Float64}, w::Union{Nothing, Vector{Vector{Float64}}}=nothing)

Calculate the overlap, kinetic energy, and optionally Coulomb interaction matrices for a given set of basis functions.

Arguments

  • α: A list of width matrices for the Gaussian basis functions.
  • K: Kinetic energy matrix.
  • w (optional): List of weight vectors for the particles involved.

Returns

  • overlap: The overlap matrix for the basis functions.
  • kinetic: The kinetic energy matrix for the basis functions.
  • Coulomb: The Coulomb interaction matrix, if weight vectors w are specified.
source
FewBodyPhysics.S_energyFunction
S_energy(bij::Vector{Float64}, K::Matrix{Float64}, w::Vector{Vector{Float64}})

Compute the ground state energy of the system using the basis functions specified by the width parameters bij.

Arguments

  • bij: A vector of width parameters for the Gaussian basis functions.
  • K: Kinetic energy matrix.
  • w: List of weight vectors for the particles involved.

Returns

  • E0: The lowest eigenvalue computed from the Hamiltonian.
source
FewBodyPhysics.P_elementsFunction
P_elements(a::Vector{Float64}, b::Vector{Float64}, A::PositiveDefiniteSymmetricMatrix, B::PositiveDefiniteSymmetricMatrix, K::Matrix{Float64}, w::Union{Nothing, Vector{Vector{Float64}}}=nothing)

Calculate the perturbation matrix elements given two basis states represented by vectors a and b, and their respective width matrices A and B.

Arguments

  • a: The coefficient vector for basis state i.
  • b: The coefficient vector for basis state j.
  • A: Width matrix for state i.
  • B: Width matrix for state j.
  • K: Kinetic energy matrix.
  • w (optional): List of weight vectors for the particles involved.

Returns

  • M1: The overlap perturbation term.
  • kinetic: The kinetic energy perturbation term.
  • Coulomb_term (optional): The Coulomb interaction perturbation term, if weight vectors w are provided.
source
FewBodyPhysics.pion_nucleonFunction
pion_nucleon(alphas, masses, params)

Calculate the overlap and kinetic matrices for a pion-nucleon system.

Arguments

  • alphas: A vector of alpha values, which are parameters related to the Gaussian basis functions.
  • masses: A 2-element vector containing the masses of the nucleon and the pion, respectively.
  • params: A 2-element vector containing the parameters b and S.

Returns

  • overlap: A matrix representing the overlap between the basis functions.
  • kinetic: A matrix representing the kinetic energy elements.

Description

The function calculates the overlap and kinetic matrices for a pion-nucleon system using Gaussian basis functions. The overlap matrix elements are calculated as integrals of the product of two basis functions, and the kinetic matrix elements are calculated as integrals of the product of the derivatives of two basis functions.

The function uses the reduced mass of the pion-nucleon system, the parameter b related to the width of the Gaussian functions, and the parameter S related to the amplitude of the Gaussian functions.

The matrices are symmetric, and the diagonal elements of the overlap matrix are 1. The off-diagonal elements are calculated using the alpha parameters and the b and S parameters.

source
FewBodyPhysics.ComputeEigenSystemFunction
ComputeEigenSystem(bs, masses, params)

Calculate the eigenvalues and eigenvectors of a system defined by parameters bs, masses, and params.

Arguments

  • bs: Array of parameter values used in the computation.
  • masses: Array of masses, representing physical properties of the system.
  • params: Additional parameters required for the calculation.

Returns

  • Tuple of eigenvalues (E) and eigenvectors (c).
source
FewBodyPhysics.GetMinimumEnergyFunction
GetMinimumEnergy(bs, masses, params)

Compute the minimum energy of a system characterized by bs, masses, and params.

Arguments

  • bs: Array of parameter values used in the computation.
  • masses: Array of masses, representing physical properties of the system.
  • params: Additional parameters required for the calculation.

Returns

  • Minimum energy value of the system.
source
FewBodyPhysics.OptimizeGlobalParametersFunction
OptimizeGlobalParameters(ngauss, dim, bmax, masses, params)

Perform global optimization over a given parameter space to find optimal parameters for a physical system.

Arguments

  • ngauss: Number of Gaussian functions used in the optimization.
  • dim: Dimensionality of the parameter space.
  • bmax: Maximum value of the parameter b.
  • masses: Array of masses, representing physical properties of the system.
  • params: Additional parameters used in the optimization.

Returns

  • E_list: List of optimized energies.
  • gaussians: List of Gaussian functions used.
  • coords: Optimized coordinates in the parameter space.
  • eigenvectors: Eigenvectors corresponding to the optimized coordinates.
  • masses: Updated masses array.
source

Sampling

FewBodyPhysics.corputFunction
corput(n::Int, b::Int=3) -> Float64

Generates the nth term of the van der Corput sequence in the given base b, which is often used for quasi-random number generation.

Arguments

  • n::Int: The position of the term in the sequence to calculate.
  • b::Int=3: The base for the sequence. Defaults to 3 if not provided.

Returns

  • Float64: The nth term in the van der Corput sequence for the given base b.

Example

```julia corput(1, 2) # Returns 0.5 for base 2

source
FewBodyPhysics.haltonFunction
halton(n::Int, d::Int) -> Vector{Float64}

Generates a point in the Halton sequence with d dimensions, used in quasi-random sampling.

# Arguments
- `n::Int`: The index of the sequence point to generate.  
- `d::Int`: The dimensionality of the Halton sequence (i.e., the number of bases to use).
# Returns
Vector{Float64}: A vector of length d representing the nth point in the Halton sequence.
source
FewBodyPhysics.calculate_energiesFunction
calculate_energies(num_gauss::Int, w_transformed::Vector{Vector{Float64}}, K_transformed::Matrix{Float64}, b1::Float64, method::Symbol) -> Tuple{Vector{Float64}, Vector{Vector{Float64}}, Vector{Int}, Float64}

Calculates and refines energies for a set of Gaussian basis functions using quasi-random or pseudo-random methods.

source
FewBodyPhysics.run_simulationFunction
run_simulation(num_gauss::Int, method::Symbol, w_transformed::Vector{Vector{Float64}}, K_transformed::Matrix{Float64}, plot_result::Bool=true) -> Tuple{Plots.Plot, Float64, Vector{Vector{Float64}}}

Runs a simulation to calculate energy values for Gaussian basis functions and optionally plots the results.

source

Constants

FewBodyPhysics.μConstant

μ = (mbare * mpi0) / (mbare + mpi0)

Reduced mass of the nucleon-pion system in MeV/c^2.

source