Module Index
ModeCouplingTheory.ActiveMCTKernelModeCouplingTheory.ActiveMultiComponentKernelModeCouplingTheory.EulerSolverModeCouplingTheory.ExponentiallyDecayingKernelModeCouplingTheory.InterpolatingKernelModeCouplingTheory.MemoryEquationModeCouplingTheory.MemoryEquationCoefficientsModeCouplingTheory.MemoryEquationSolutionModeCouplingTheory.SCGLEKernelModeCouplingTheory.SchematicDiagonalKernelModeCouplingTheory.SchematicF123KernelModeCouplingTheory.SchematicF1KernelModeCouplingTheory.SchematicF2KernelModeCouplingTheory.SchematicMatrixKernelModeCouplingTheory.SjogrenKernelModeCouplingTheory.TaggedActiveMCTKernelModeCouplingTheory.TaggedSchematicF2KernelModeCouplingTheory.TimeDoublingSolverModeCouplingTheory.MSDModeCouplingKernelModeCouplingTheory.MSDMultiComponentModeCouplingKernelModeCouplingTheory.MSDSCGLEKernel_constructorModeCouplingTheory.ModeCouplingKernelModeCouplingTheory.MultiComponentModeCouplingKernelModeCouplingTheory.TaggedModeCouplingKernelModeCouplingTheory.TaggedMultiComponentModeCouplingKernelModeCouplingTheory.allocate_results!ModeCouplingTheory.allocate_temporary_arraysModeCouplingTheory.compute_C3_immutableModeCouplingTheory.compute_C3_mutable!ModeCouplingTheory.convert_multicomponent_structure_factorModeCouplingTheory.do_time_steps!ModeCouplingTheory.evaluate_kernelModeCouplingTheory.evaluate_kernel!ModeCouplingTheory.find_errorModeCouplingTheory.find_relaxation_timeModeCouplingTheory.get_FModeCouplingTheory.get_KModeCouplingTheory.get_tModeCouplingTheory.get_ΔGModeCouplingTheory.initialize_F_temp!ModeCouplingTheory.initialize_K_temp!ModeCouplingTheory.initialize_integrals!ModeCouplingTheory.initialize_output_arraysModeCouplingTheory.mygetindexModeCouplingTheory.mymul!ModeCouplingTheory.new_time_mapping!ModeCouplingTheory.solveModeCouplingTheory.solve_steady_stateModeCouplingTheory.surface_d_dim_unit_sphereModeCouplingTheory.update_F!ModeCouplingTheory.update_Fuchs_parameters!ModeCouplingTheory.update_K!ModeCouplingTheory.update_integrals!ModeCouplingTheory.volume_d_dim_sphere
Detailed API
ModeCouplingTheory.ActiveMCTKernel — TypeActiveMCTKernel(ρ, k_array, wk, w0, Sk, dim=3)Implements the following active MCT kernel:
M(k,t) = ρ / (2(2π)^dim ) ∫ dq w(k) V(k,q)^2 F(q,t) F(k-q,t)
with the vertices:
V(k,q) = c(q) * (k * q)/k + c(p) (k * p)/k
Arguments:
- ρ: number density
- k_array: array of wavevectors at which to evaluate the kernel
- wk: steady-state velocity correlations ( w(k) ) - this depends on k
- w0: local velocity correlations ( w(∞) ) - this is a constant
- Sk: steady-state structure factor
- dim: dimensionality of the problem (the default
dim=3)
ModeCouplingTheory.ActiveMultiComponentKernel — TypeActiveMultiComponentKernel(ρₐ, k_array, wk, w0, Sk, dim)Implements the following multi-component active MCT kernel:
M{αβ}(k,t) = ρ / (2*(2π)^dim) ∑{μνμ'ν'λ} ∫ dq F{μμ'}(q,t) F{νν'}(k-q,t) V{μνα}(k,q) V{μ'ν'λ}(k,q) [w{λβ}(k)]^(-1)
where Greek indices {...} denote species labels and the expression for the vertices V{μνα}(k,q) is given in the documentation. Note: the input data (wk, Sk) should be given in the format Vector{ Matrix }, with the Vector having length Nk and the Matrix having size Ns x Ns. Nk is the number of k-points and Ns is the number of species in the mixture.
Arguments:
- ρₐ: number densities of each species (a vector of length Ns)
- k_array: k-values at which to evaluate the kernel (a vector of length Nk)
- wk: velocity correlations of each species (a vector of Ns x Ns matrices)
- w0: local velocity correlations (a Ns x Ns matrix)
- Sk: partial structure factors (a vector of Ns x Ns matrices)
- dim: dimensionality of the kernel (the default
dim=3)
ModeCouplingTheory.EulerSolver — MethodEulerSolver(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:
equationan instance of MemoryEquationt_maxwhen this time value is reached, the integration returnsΔtfixed time stepverboseiftrue, information will be printed to STDOUT
ModeCouplingTheory.ExponentiallyDecayingKernel — TypeExponentiallyDecayingKernel{T<:Number} <: MemoryKernelScalar kernel with fields ν and τ which when called returns ν exp(-t/τ).
ModeCouplingTheory.InterpolatingKernel — MethodInterpolatingKernel(t, M; k=1)Uses the package Dierckx to provide a kernel that interpolates the data M defined on grid t using spline interpolation of degree k.
ModeCouplingTheory.MemoryEquation — MethodMemoryEquation(α, β, γ, F₀::T, ∂ₜF₀::T, kernel::MemoryKernel) where TArguments:
α: coefficient in front of the second derivative term. IfαandF₀are both vectors,αwill automatically be converted to a diagonal matrix, to make them compatible.β: coefficient in front of the first derivative term. IfβandF₀are both vectors,βwill automatically be converted to a diagonal matrix, to make them compatible.γ: coefficient in front of the second derivative term. IfγandF₀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)kernelinstance of aMemoryKernelthat when called on F₀ and t=0, evaluates to the initial condition of the memory kernel.
ModeCouplingTheory.MemoryEquationCoefficients — MethodMemoryEquationCoefficients(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.
ModeCouplingTheory.MemoryEquationSolution — TypeMemoryEquationSolutionsolution 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.
ModeCouplingTheory.SCGLEKernel — MethodSCGLEKernel(ϕ, k_array, S_array)Constructor of a SCGLEKernel. It implements the kernel K(k,t) = λ(k)Δζ(t) where Δζ(t) = kBT / (36 π ϕ) ∫dq q⁴ M²(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 fractionk_array: vector of wavenumbers at which the structure factor is knownS_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)
ModeCouplingTheory.SchematicDiagonalKernel — TypeSchematicDiagonalKernel{T<:Union{SVector, Vector}} <: MemoryKernelMatrix kernel with field ν which when called returns Diagonal(ν .* F .^ 2), i.e., it implements a non-coupled system of SchematicF2Kernels.
ModeCouplingTheory.SchematicF123Kernel — TypeSchematicF1Kernel{T<:Number} <: MemoryKernelScalar kernel with fields ν1, ν2, and ν3 which when called returns ν1 * F^1 + ν2 * F^2 + ν3 * F^3.
ModeCouplingTheory.SchematicF1Kernel — TypeSchematicF1Kernel{T<:Number} <: MemoryKernelScalar kernel with field ν which when called returns ν F.
ModeCouplingTheory.SchematicF2Kernel — TypeSchematicF2Kernel{T<:Number} <: MemoryKernelScalar kernel with field ν which when called returns ν F^2.
ModeCouplingTheory.SchematicMatrixKernel — TypeSchematicMatrixKernel{T<:Union{SVector, Vector}} <: MemoryKernelMatrix kernel with field ν which when called returns ν * F * Fᵀ, i.e., it implements Kαβ = νFαFβ.
ModeCouplingTheory.SjogrenKernel — TypeSjogrenKernelMemory kernel that implements the kernel K[1] = ν1 F[1]^2, K[2] = ν2 F[1] F[2]. Consider using Static Vectors for performance.
ModeCouplingTheory.TaggedActiveMCTKernel — TypeTaggedActiveMCTkernel(ρ, k_array, wk, w0, Sk, Fsol, dim)This function implements the following kernel for the self-intermediate scattering function:
M(k,t) = ρ w0 / ((2π)^dim ) ∫ dq ( (k*p)/(2k) * c(p) )^2 F(q,t) Fs(k-q,t)
Arguments
- ρ: number density
- k_array: k-values at which to evaluate the kernel
- wk: steady-state velocity correlations ( w(k) ) - this depends on k
- w0: local velocity correlations ( w(∞) ) - this is a constant
- Sk: steady-state structure factor
- Fsol: solution of the collective mode-coupling equation
- dim: dimensionality of the problem (the default
dim=3)
ModeCouplingTheory.TaggedSchematicF2Kernel — TypeSchematicF2Kernel{T<:Number} <: MemoryKernelScalar tagged particle kernel with field ν which when called returns ν F*Fs.
ModeCouplingTheory.TimeDoublingSolver — MethodTimeDoublingSolver(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 MemoryEquationN: The number of time points in the interval is equal to4Nt_max: when this time value is reached, the integration returnsΔt: starting time step, this will be doubled repeatedlymax_iterations: the maximal number of iterations before convergence is reached for each time doubling steptolerance: while the error is bigger than this value, convergence is not reached. The error by default is computed as the absolute sum of squaresverbose: iftrue, information will be printed to STDOUTismutable: iftrueand if the type of F is mutable, the solver will try to avoid allocating many temporariesinit_with_memoryiftrue, the solver will include the memory kernels for the short time expansion to bootstrap the algorithm. Sacrifices performance and memory for accuracy.
ModeCouplingTheory.MSDModeCouplingKernel — MethodMSDModeCouplingKernel(ρ, 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 densitykBT: Thermal energym: particle massk_array: vector of wavenumbers at which the structure factor is knownSₖ: vector with the elements of the structure factorsol: 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)
ModeCouplingTheory.MSDMultiComponentModeCouplingKernel — MethodMSDMultiComponentModeCouplingKernel(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 densitykBT: Thermal energym: particle mass of species sk_array: vector of wavenumbers at which the structure factor is knownSₖ: vector with the elements of the structure factorsol: 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)
ModeCouplingTheory.MSDSCGLEKernel_constructor — MethodMSDSCGLEKernel(ϕ, D⁰, k_array, Sₖ, sol)
Constructor of a MSDSCGLEKernel. It implements the kernel
K(k,t) = D⁰ / (36πϕ) ∫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:
ϕ: volume fractionD⁰: Short times diffusion coefficientk_array: vector of wavenumbers at which the structure factor is knownSₖ: vector with the elements of the structure factorsol: a solution object of an equation with a SCGLEKernel.
Returns:
an instance k of MSDSCGLEKernel <: MemoryKernel, which can be evaluated like: k = evaluate_kernel(kernel, F, t)
ModeCouplingTheory.ModeCouplingKernel — MethodModeCouplingKernel(ρ, 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 densitykBT: Thermal energym: particle massk_array: vector of wavenumbers at which the structure factor is knownSₖ: 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)
ModeCouplingTheory.MultiComponentModeCouplingKernel — MethodMultiComponentModeCouplingKernel(ρ, 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 specieskBT: Thermal energym: vector of particle masses for each speciesk_array: vector of wavenumbers at which the structure factor is knownSₖ: aVectorof NkSMatrixs 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)`
ModeCouplingTheory.TaggedModeCouplingKernel — MethodTaggedModeCouplingKernel(ρ, 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 densitykBT: Thermal energym: particle massk_array: vector of wavenumbers at which the structure factor is knownSₖ: vector with the elements of the structure factorsol: 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)
ModeCouplingTheory.TaggedMultiComponentModeCouplingKernel — MethodTaggedMultiComponentModeCouplingKernel(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 densitykBT: Thermal energym: particle massk_array: vector of wavenumbers at which the structure factor is knownSₖ: vector with the elements of the structure factorsol: 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)
ModeCouplingTheory.allocate_results! — Methodallocate_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.
ModeCouplingTheory.allocate_temporary_arrays — Methodallocate_temporary_arrays(equation::AbstractMemoryEquation, solver::TimeDoublingSolver)Returns a SolverCache containing several arrays that are used for intermediate calculations.
ModeCouplingTheory.compute_C3_immutable — Methodcomputec3immutable(equation::MemoryEquation, solver::TimeDoublingSolver, temp_arrays::SolverCache, it::Int)
computes one of the matrices necessary to propagate F in time.
ModeCouplingTheory.compute_C3_mutable! — Methodcomputec3mutable(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
ModeCouplingTheory.convert_multicomponent_structure_factor — Method`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]ModeCouplingTheory.do_time_steps! — Methoddo_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.
ModeCouplingTheory.evaluate_kernel — Functionevaluate_kernel(kernel::MemoryKernel, F, t)Evaluates the memory kernel out-place. It may mutate the content of kernel.
Returns
outthe kernel evaluated at (F, t)
ModeCouplingTheory.evaluate_kernel! — Functionevaluate_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.
ModeCouplingTheory.find_error — Methodfind_error(F_new::T, F_old::T) where TFinds the error between a new and old iteration of F. The returned scalar will be compared to the tolerance to establish convergence.
ModeCouplingTheory.find_relaxation_time — Methodfind_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 signalFis knownF: the signalthreshold: see abovemode: uses interpolation either in a logarithmic or linear space.
ModeCouplingTheory.get_F — Method`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.
ModeCouplingTheory.get_K — Method`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.
ModeCouplingTheory.get_t — Method`get_t(sol::MemoryEquationSolution)`obtains the time grid t from a MemoryEquationSolution object. Equivalent to sol.t.
ModeCouplingTheory.get_ΔG — MethodΔG = (kBT/60π²)∫dk k⁴[(1/S)(∂S/∂k)]²[(F/S)]² just the G. Naegele formula for shear viscosity relaxation. It must work for the MCT kernel.
ModeCouplingTheory.initialize_F_temp! — Methodinitialize_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.
ModeCouplingTheory.initialize_K_temp! — Methodinitialize_K_temp!(solver::TimeDoublingSolver, kernel::MemoryKernel, temp_arrays::SolverCache)Evaluates the memory kernel at the first 2N time points.
ModeCouplingTheory.initialize_integrals! — Methodinitialize_integrals!(equation::AbstractMemoryEquation, solver::TimeDoublingSolver, temp_arrays::SolverCache)Initializes the integrals on the first 2N time points as prescribed in the literature.
ModeCouplingTheory.initialize_output_arrays — Methodinitialize_output_arrays(equation::AbstractMemoryEquation)initializes arrays that the solver will push results into.
ModeCouplingTheory.mygetindex — Methodmygetindex(A, inds...)Recursively index into an arbitrarily nested array (array of arrays of arrays, etc.). Supports scalars, ranges, and Colon().
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,α,β)
ModeCouplingTheory.new_time_mapping! — Methodnewtimemapping!(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.
ModeCouplingTheory.solve — Methodsolve(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:
tan array of time valuesFThe solution in an array of which the last dimension corresponds to the time.KThe memory kernel corresponding to eachF
ModeCouplingTheory.solve_steady_state — Methodsolve_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 FF₀: initial condition of F. This is also the initial condition of the rootfinding method.kernel: callable memory kernelmax_iterations: the maximal number of iterations before convergence is reachedtolerance: 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 iterationsverbose: iftrue, information will be printed toSTDOUTinplace: 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
ModeCouplingTheory.surface_d_dim_unit_sphere — MethodSurface of a d-dimensional sphere
ModeCouplingTheory.update_F! — Methodupdate_F!(solver::TimeDoublingSolver, temp_arrays::SolverCache, it::Int)updates F using the formula c1F = -KC2 + C3.
ModeCouplingTheory.update_Fuchs_parameters! — Methodupdate_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.
ModeCouplingTheory.update_K! — Methodupdate_K!(solver::TimeDoublingSolver, kernel::MemoryKernel, temp_arrays::SolverCache, it::Int)evaluates the memory kernel, updating the value in solver.temparrays.Ktemp
ModeCouplingTheory.update_integrals! — Methodupdate_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.
ModeCouplingTheory.volume_d_dim_sphere — MethodVolume of a d-dimensional sphere