Module Index

Detailed API

ModeCouplingTheory.EulerSolverMethod
EulerSolver(equation::MemoryEquation; t_max=10.0^2, Δt=10^-3, verbose=false)

Constructs a solver object that, when called solve upon will solve an MemoryEquation using a forward Euler method. It will discretise the integral using a Trapezoidal rule. Use this solver only for testing purposes. It is a wildy inefficient way to solve MCT-like equations.

Arguments:

  • equation an instance of MemoryEquation
  • t_max when this time value is reached, the integration returns
  • Δt fixed time step
  • verbose if true, information will be printed to STDOUT
source
ModeCouplingTheory.MemoryEquationMethod
MemoryEquation(α, β, γ, F₀::T, ∂ₜF₀::T, kernel::MemoryKernel) where T

Arguments:

  • α: coefficient in front of the second derivative term. If α and F₀ are both vectors, α will automatically be converted to a diagonal matrix, to make them compatible.
  • β: coefficient in front of the first derivative term. If β and F₀ are both vectors, β will automatically be converted to a diagonal matrix, to make them compatible.
  • γ: coefficient in front of the second derivative term. If γ and F₀ are both vectors, γ will automatically be converted to a diagonal matrix, to make them compatible.
  • δ: coefficient in front of the constant term. δ and F0 must by types that can be added together. If δ is a number and F0 is a vector, this conversion will be done automatically.
  • F₀: initial condition of F(t)
  • ∂ₜF₀ initial condition of the derivative of F(t)
  • kernel instance of a MemoryKernel that when called on F₀ and t=0, evaluates to the initial condition of the memory kernel.
source
ModeCouplingTheory.MemoryEquationCoefficientsMethod
MemoryEquationCoefficients(a, b, c, d, F₀)

Constructor of the MemoryEquationCoefficients struct, which holds the values of the coefficients α, β, γ and δ. It will convert common cases automatically to make the types compatible. For example, if α and F0 are both given as vectors, α is converted to a Diagonal matrix.

source
ModeCouplingTheory.MemoryEquationSolutionType
MemoryEquationSolution

solution object that holds the solution of an AbstractMemoryEquation. It has 4 fields t: array of t values F: array of F for all t K: array of K for all T solver: solver object that holds the solver settings

an MCT object can be indexed such that sol=MemoryEquationSolution sol[2] gives the F[2] for all t.

source
ModeCouplingTheory.SCGLEKernelMethod
SCGLEKernel(ϕ, k_array, S_array)

Constructor of a SCGLEKernel. It implements the kernel K(k,t) = λ(k)Δζ(t) where Δζ(t) = kBT / (144 π^2 ϕ) ∫dq⃗ M^2(k,q) Fs(q,t) F(q,t) in which k and q are vectors. For Fs and F being the self- and collective-intermediate scattering function.

Arguments:

  • ϕ: volume fraction
  • k_array: vector of wavenumbers at which the structure factor is known
  • S_array: concatenated vector for ones(length(structure factor)) and structure factor

Returns:

an instance k of SCGLEKernel <: MemoryKernel, which can be called both in-place and out-of-place: evaluate_kernel!(out, kernel, F, t) out = evaluate_kernel(kernel, F, t)

source
ModeCouplingTheory.SchematicDiagonalKernelType
SchematicDiagonalKernel{T<:Union{SVector, Vector}} <: MemoryKernel

Matrix kernel with field ν which when called returns Diagonal(ν .* F .^ 2), i.e., it implements a non-coupled system of SchematicF2Kernels.

source
ModeCouplingTheory.TimeDoublingSolverMethod
TimeDoublingSolver(N=32, Δt=10^-10, t_max=10.0^10, max_iterations=10^4, tolerance=10^-10, verbose=false, ismutable=true, init_with_memory=true)

Uses the algorithm devised by Fuchs et al.

Arguments:

  • equation: an instance of MemoryEquation
  • N: The number of time points in the interval is equal to 4N
  • t_max: when this time value is reached, the integration returns
  • Δt: starting time step, this will be doubled repeatedly
  • max_iterations: the maximal number of iterations before convergence is reached for each time doubling step
  • tolerance: while the error is bigger than this value, convergence is not reached. The error by default is computed as the absolute sum of squares
  • verbose: if true, information will be printed to STDOUT
  • ismutable: if true and if the type of F is mutable, the solver will try to avoid allocating many temporaries
  • init_with_memory if true, the solver will include the memory kernels for the short time expansion to bootstrap the algorithm. Sacrifices performance and memory for accuracy.
source
ModeCouplingTheory.MSDModeCouplingKernelMethod

MSDModeCouplingKernel(ρ, kBT, m, k_array, Sₖ, sol, taggedsol; dims=3)

Constructor of a MSDModeCouplingKernel. It implements the kernel K(k,t) = ρ kBT / (6π^2m) ∫dq q^4 c(q)^2 F(q,t) Fs(q,t) where the integration runs from 0 to infinity. F and Fs are the coherent and incoherent intermediate scattering functions, and must be passed in as solutions of the corresponding equations.

Arguments:

  • ρ: number density
  • kBT: Thermal energy
  • m : particle mass
  • k_array: vector of wavenumbers at which the structure factor is known
  • Sₖ: vector with the elements of the structure factor
  • sol: a solution object of an equation with a ModeCouplingKernel.
  • taggedsol: a solution object of an equation with a TaggedModeCouplingKernel.

Returns:

an instance k of MSDModeCouplingKernel <: MemoryKernel, which can be evaluated like: k = evaluate_kernel(kernel, F, t)

source
ModeCouplingTheory.MSDMultiComponentModeCouplingKernelMethod

MSDMultiComponentModeCouplingKernel(s, ρ, kBT, m, k_array, Sₖ, sol, taggedsol; dims=3)

Constructor of a MSDModeCouplingKernel. It implements the kernel Kₛ(k,t) = ρ kBT / (6π^2mₛ) Σαβ ∫dq q^4 csα(q)csβ(q) Fαβ(q,t) Fₛ(q,t) where the integration runs from 0 to infinity. F and Fs are the coherent and incoherent intermediate scattering functions, and must be passed in as solutions of the corresponding equations.

Arguments:

  • s: the species for which to solve this equation
  • ρ: number density
  • kBT: Thermal energy
  • m : particle mass of species s
  • k_array: vector of wavenumbers at which the structure factor is known
  • Sₖ: vector with the elements of the structure factor
  • sol: a solution object of an equation with a MultiComponentModeCouplingKernel.
  • taggedsol: a solution object of an equation with a TaggedMultiComponentModeCouplingKernel.

Returns:

an instance k of MSDMultiComponentModeCouplingKernel <: MemoryKernel, which can be evaluated like: k = evaluate_kernel(kernel, F, t)

source
ModeCouplingTheory.ModeCouplingKernelMethod
ModeCouplingKernel(ρ, kBT, m, k_array, Sₖ; dims=3)

Constructor of a ModeCouplingKernel. It implements the kernel K(k,t) = ρ kBT / (16π^3m) ∫dq V^2(k,q) F(q,t) F(k-q,t) in which k and q are vectors.

Arguments:

  • ρ: number density
  • kBT: Thermal energy
  • m : particle mass
  • k_array: vector of wavenumbers at which the structure factor is known
  • Sₖ: structure factor

Returns:

an instance k of ModeCouplingKernel <: MemoryKernel, which can be called both in-place and out-of-place: evaluate_kernel!(out, kernel, F, t) out = evaluate_kernel(kernel, F, t)

source
ModeCouplingTheory.MultiComponentModeCouplingKernelMethod
MultiComponentModeCouplingKernel(ρ, kBT, m, k_array, Sₖ; dims=3)

Constructor of a MultiComponentModeCouplingKernel. It implements the kernel Kαβ(k,t) = ρ / (2 xα xβ (2π)³) Σμνμ'ν' ∫dq Vμ'ν'(k,q) Fμμ'(q,t) Fνν'(k-q,t) Vμν(k,q) in which k and q are vectors and α and β species labels.

Arguments:

  • ρ: vector of number densities for each species
  • kBT: Thermal energy
  • m : vector of particle masses for each species
  • k_array: vector of wavenumbers at which the structure factor is known
  • Sₖ: a Vector of Nk SMatrixs containing the structure factor of each component at each wave number

Returns:

an instance k of ModeCouplingKernel <: MemoryKernel, which can be called both in-place and out-of-place: k = MultiComponentModeCouplingKernel(ρ, kBT, m, karray, Sₖ) `evaluatekernel!(out, kernel, F, t)out = evaluate_kernel(kernel, F, t)`

source
ModeCouplingTheory.TaggedModeCouplingKernelMethod

TaggedModeCouplingKernel(ρ, kBT, m, k_array, Sₖ, sol; dims=3)

Constructor of a Tagged ModeCouplingKernel. It implements the kernel K(k,t) = ρ kBT / (8π^3m) ∫dq V^2(k,q) F(q,t) Fs(k-q,t) in which k and q are vectors. Here V(k,q) = c(q) (k dot q)/k.

Arguments:

  • ρ: number density
  • kBT: Thermal energy
  • m : particle mass
  • k_array: vector of wavenumbers at which the structure factor is known
  • Sₖ: vector with the elements of the structure factor
  • sol: a solution object of an equation with a ModeCouplingKernel.

Returns:

an instance k of TaggedModeCouplingKernel <: MemoryKernel, which can be called both in-place and out-of-place: evaluate_kernel!(out, kernel, Fs, t) out = evaluate_kernel(kernel, Fs, t)

source
ModeCouplingTheory.TaggedMultiComponentModeCouplingKernelMethod

TaggedMultiComponentModeCouplingKernel(s, ρ, kBT, m, k_array, Sₖ, sol; dims=3)

Constructor of a Tagged ModeCouplingKernel. It implements the kernel Kₛ(k,t) = ρ kBT / (8π^3 mₛ) Σαβ ∫dq V^2sαβ(k,q) Fαβ(q,t) Fₛ(k-q,t) in which k and q are vectors.

Arguments:

  • s: the species for which to solve this equation
  • ρ: number density
  • kBT: Thermal energy
  • m : particle mass
  • k_array: vector of wavenumbers at which the structure factor is known
  • Sₖ: vector with the elements of the structure factor
  • sol: a solution object of an equation with a MultiComponentModeCouplingKernel.

Returns:

an instance k of TaggedMultiComponentModeCouplingKernel <: MemoryKernel, which can be called both in-place and out-of-place: evaluate_kernel!(out, kernel, Fs, t) out = evaluate_kernel(kernel, Fs, t)

source
ModeCouplingTheory.allocate_results!Method
allocate_results!(t_array, F_array, K_array, equation::AbstractMemoryEquation, solver::TimeDoublingSolver, temp_arrays::SolverCache; istart=2(solver.N)+1, iend=4(solver.N))

pushes the found solution, stored in temp_arrays with indices istart until istop to the output arrays.

source
ModeCouplingTheory.compute_C3_immutableMethod

computec3immutable(equation::MemoryEquation, solver::TimeDoublingSolver, temp_arrays::SolverCache, it::Int)

computes one of the matrices necessary to propagate F in time.

source
ModeCouplingTheory.compute_C3_mutable!Method

computec3mutable(equation::MemoryEquation, solver::TimeDoublingSolver, temp_arrays::SolverCache, it::Int)

computes one of the matrices necessary to propagate F in time, by mutating it inplace. The commented code is the corresponding scalar equivalent

source
ModeCouplingTheory.convert_multicomponent_structure_factorMethod
`convert_multicomponent_structure_factor(Sk_in::Matrix{Vector{T}})``

This function takes a matrix of vectors Sk_in where each vector contains structure factor information for different species at different k-points and converts it into a vector of SMatrix.

Input:

  • Sk_in::Matrix{Vector{T}} - A square matrix of vectors, where each vector contains structure factor information for a specific set of species at different k-points.

Output:

  • Sk_out - A Vector of SMatrix where each SMatrix contains structure factor information for different species at a single k-point.

Example

Ns = 2;
Nk = 3;
S11 = [1,2,3]; S21 = [4,5,6]; S22 = [8,9,10];
S = [zeros(Nk) for i=1:2, j=1:2];
S[1,1] = S11; S[1,2] = S21; S[2,1] = S21; S[2,2] = S22;
convert_multicomponent_structure_factor(S)

    3-element Vector{StaticArraysCore.SMatrix{2, 2, Float64, 4}}:
    [1.0 4.0; 4.0 8.0]
    [2.0 5.0; 5.0 9.0]
    [3.0 6.0; 6.0 10.0]
source
ModeCouplingTheory.do_time_steps!Method
do_time_steps!(equation::MemoryEquation, solver::TimeDoublingSolver, kernel::MemoryKernel, temp_arrays::SolverCache)

Solves the equation on the time points with index 2N+1 until 4N, for each point doing a recursive iteration to find the solution to the nonlinear equation C1 F = -C2 M(F) + C3.

source
ModeCouplingTheory.evaluate_kernelFunction
evaluate_kernel(kernel::MemoryKernel, F, t)

Evaluates the memory kernel out-place. It may mutate the content of kernel.

Returns

  • out the kernel evaluated at (F, t)
source
ModeCouplingTheory.evaluate_kernel!Function
evaluate_kernel!(out, kernel::MemoryKernel, F, t)

Evaluates the memory kernel in-place, overwriting the elements of the out variable. It may mutate the content of kernel.

source
ModeCouplingTheory.find_errorMethod
find_error(F_new::T, F_old::T) where T

Finds the error between a new and old iteration of F. The returned scalar will be compared to the tolerance to establish convergence.

source
ModeCouplingTheory.find_relaxation_timeMethod
find_relaxation_time(t, F; threshold=exp(-1), mode=:log)

Finds the time at which a signal first decreases below some threshold. Uses interpolation to find an accurate result. if it never decreases below the threshold, returns Inf.

Arguments:

  • t: a vector of time values at which the signal F is known
  • F: the signal
  • threshold: see above
  • mode: uses interpolation either in a logarithmic or linear space.
source
ModeCouplingTheory.get_FMethod
`get_F(sol::MemoryEquationSolution)`

obtains the solution F from a MemoryEquationSolution object. Equivalent to sol.F.

`get_F(sol::MemoryEquationSolution, I...)`

obtains the solution F from a MemoryEquationSolution object and indexes into it. Enables convenient indexing into the multidimensional object. I is interpreted as a set of indices.

Examples:

If sol is the solution to a scalar equation get_F(sol, 2:4) gets the elements 2:4 If sol is the solution to a vector-valued equation get_F(sol, 5, 2:43) gets the solution at the 5th time point for vector indices 2:43. If sol is the solution to a vector-valued multicomponent equation get_F(sol, 5, 2:43, (1,2)) gets the solution at the 5th time point for vector indices 2:43, for species 1 and 2.

source
ModeCouplingTheory.get_KMethod
`get_K(sol::MemoryEquationSolution)`

obtains the kernel K from a MemoryEquationSolution object. Equivalent to sol.K.

`get_K(sol::MemoryEquationSolution, I...)`

obtains the solution K from a MemoryEquationSolution object and indexes into it. Enables convenient indexing into the multidimensional object. See get_F for examples.

source
ModeCouplingTheory.get_tMethod
`get_t(sol::MemoryEquationSolution)`

obtains the time grid t from a MemoryEquationSolution object. Equivalent to sol.t.

source
ModeCouplingTheory.initialize_F_temp!Method
initialize_F_temp!(equation::MemoryEquation, solver::TimeDoublingSolver, temp_arrays::SolverCache)

Fills the first 2N entries of the temporary arrays of F using forward Euler without a memory kernel in order to kickstart the algorithm' scheme.

source
ModeCouplingTheory.initialize_integrals!Method
initialize_integrals!(equation::AbstractMemoryEquation, solver::TimeDoublingSolver, temp_arrays::SolverCache)

Initializes the integrals on the first 2N time points as prescribed in the literature.

source
ModeCouplingTheory.mymul!Method

" mymul!(c,a,b,α,β)

prescribes how types used in this solver should be multiplied in place. In particular, it performs C.= βC .+ αa*b. defaults to mul!(c,a,b,α,β)

source
ModeCouplingTheory.new_time_mapping!Method

newtimemapping!(equation::AbstractMemoryEquation, solver::TimeDoublingSolver, temp_arrays::SolverCache)

Performs the time mapping central to Fuchs' algorithm with the conventions prescribed in "Flenner, Elijah, and Grzegorz Szamel. Physical Review E 72.3 (2005): 031508". Performs them inplace if solver.inplace = true in order to avoid unnecessary allocations.

source
ModeCouplingTheory.solveMethod
solve(equation::AbstractMemoryEquation, solver::Solver)

Solves the AbstractMemoryEquation with the provided kernel using solver. Search for a specific solver or kernel object to find more specific information.

If no solver is provided, it uses the default TimeDoublingSolver.

Returns:

  • t an array of time values
  • F The solution in an array of which the last dimension corresponds to the time.
  • K The memory kernel corresponding to each F
source
ModeCouplingTheory.solve_steady_stateMethod
solve_steady_state(γ, F₀, kernel; tolerance=10^-8, max_iterations=10^6, verbose=false, inplace=true)

Finds the steady-state solution (non-ergodicity parameter) of the generalized Langevin equation by recursive iteration of F = (K + γ)⁻¹ * K(F) * F₀ This function assumes δ = 0, since that would lead to a divergence.

Arguments:

  • γ: parameter in front of the linear term in F
  • F₀: initial condition of F. This is also the initial condition of the rootfinding method.
  • kernel: callable memory kernel
  • max_iterations: the maximal number of iterations before convergence is reached
  • tolerance: while the error is bigger than this value, convergence is not reached. The error by default is computed as the absolute sum of squares between successive iterations
  • verbose: if true, information will be printed to STDOUT
  • inplace: if true and if the type of F is mutable, the solver will try to avoid allocating many temporaries. default =true``

Returns:

  • The steady state solution
source
ModeCouplingTheory.update_Fuchs_parameters!Method
update_Fuchs_parameters!(equation::MemoryEquation, solver::TimeDoublingSolver, temp_arrays::SolverCache, it::Int)

Updates the parameters c1, c2, c3, according to the appendix of "Flenner, Elijah, and Grzegorz Szamel. Physical Review E 72.3 (2005): 031508" using the naming conventions from that paper. If F is mutable (and therefore also c1,c2,c3), it will update the variables in place, otherwise it will create new copies. This is controlled by the solver.inplace setting.

source
ModeCouplingTheory.update_K!Method
update_K!(solver::TimeDoublingSolver, kernel::MemoryKernel, temp_arrays::SolverCache, it::Int)

evaluates the memory kernel, updating the value in solver.temparrays.Ktemp

source
ModeCouplingTheory.update_integrals!Method
update_integrals!(solver::TimeDoublingSolver, equation::AbstractMemoryEquation, temp_arrays::SolverCache, it::Int)

Update the discretisation of the integral of F and K, see the literature for details.

source