Library

Types

CALiPPSO.PeriodicNumberType
PeriodicNumber{T}(value::T, L::T)

Define a number whose value is always contained between 0 and L; more precisely ∈ [0,L).

When it is initialized with value>L or value<0 the value is transformed to be its L-modulus.

source
CALiPPSO.MonoParticleType
MonoParticle{d, T}

Struct used to construct a particle in a d-dimensional space.

It's meant to be used as components of monodisperse packings. Hence, while many of its relevant properties are stored, its radius is instead stored in a field in 'MonoPacking' type.

The quantities this struct has are:

  1. X: The particle's center (as a StaticVector{d, PeriodicNumber{T}}})
  2. contact_vecs: The set of contact vectors with the particle's neighbours (as a Vector containing SVector{d,T})
  3. forces: The magnitude of contact forces (Vector{T})
  4. neighbours: The indices of the particle's neighbours (as Vector{Int64})

Naturally, 'contact_vecs', 'forces', and 'neighbours' have the same length.

See also: MonoPacking.

source
CALiPPSO.MonoPackingType
MonoPacking{d,T}

Structure used to store all the relevant info of a monodisperse hard-spheres packing in d-dimensions.

Its fields are:

  1. 'Particles': An vector of MonoParticle that contains all the particles that make the packing.
  2. 'R': The radius of all particles.
  3. 'mechanical_equilibrium': A boolean that specifies whether force balances is satisfied

all across the packing (i.e. for each particle).

  1. 'jammed': A boolean that specifies whether a packing is jammed or not.

Note that a packing might not be necessarily jammed nor in mechanical_equilibrium. That is why specifying such fields is important.

See also: MonoParticle, PolyPacking.

source
CALiPPSO.convergence_infoType

Struct to save convergence information of an CALiPPSO solution. For more info, see docs of its fields: converged, iterations, time, memory, times_LP_optim.

source

Functions and Methods

Main Exported Functions

CALiPPSO.produce_jammed_configuration!Function
produce_jammed_configuration!(Xs::Vector{SVector{d, PeriodicNumber{T}}}, R::T; <keyword arguments>) where {d, T<:Float64, I<:Int64}
produce_jammed_configuration!(Xs::Matrix{T}, R::T, L::T=1.0; <keyword arguments>) where {T<:Float64, I<:Int64}

Use CALiPPSO to generate a jammed packing from a configuration of hard-spheres with positions 'Xs' and radius 'R'.

Because this function is meant to work only with hard-spheres, clearly the initial configuration must be such that no overlaps are present. The dimensionality of the system, 'd', is automatically inferred from the size of the SVectors forming the positions vector. Besides, the periodic boundary conditions are taken into account given that each SVector contains elements of type PeriodicNumber. (If the input is of type Matrix{T}, it is converted to Vector{SVector{d,T}} before CALiPPSO begins.)

The core of this function is the so called main loop, whose essential step consists in using solve_LP_instance! to obtain the maximal inflation factor ($\Gamma$) and set of optimal particles' displacements ($\vec{\mathbf{s}}^\star = \{\mathbf{s}_i^\star\}_{i=1}^N$ –denoted as S⃗ in our scripts.) from a given configuration. The particles' size and position are updated, and a new LP instance is created and solved. (Each of these LP optimizations is considered a single iteration of CALiPPSO.) As explained in our paper, convergence is said to be achieved when both, the packing fraction cannot be further increased (i.e. $\Gamma^\star=1$), and (non-rattler) particles cannot be further displaced ($\vec{\mathbf{s}}^\star = 0$). In practice, a tolerance should be allowed in both quantities, and this is controlled by tol_Γ_convergence and tol_S_conv as described below. Another reason why the main loop might be terminated is that the maximal number of iterations has been exceeded (see below).

In case the packing thus obtained does not satisfy force balance (within a given precision), fine_tune_forces! is called on such packing. In this way, the final packing is guaranteed to be in mechanical equilibrium, within the same precision.

Output

  1. final_packing: A MonoPacking{d,T} corresponding to a jammed state (unless max_iters exceeded).
  2. conv_info: A convergence_info](@ref) struct storing the termination status of CALiPPSO and other useful information.
  3. Γs_vs_t: An array containing the optimal values of \sqrt{\Gamma} obtained after each LP optimization.
  4. smax_vs_t: An analogous array of the largest displacement (in absolute value) of stable particles.
  5. iso_vs_t: A boolean vector whose elements indicate whether the preliminary configurations were isostatic or not.

Keyword arguments

Arguments needed for calling bounds_and_cutoff

  • ℓ0::T=4*R: Upper bound for the radius of influence.
  • sqrΓ0::Real=1.01: Initialization value of $\sqrt{\Gamma}$.
  • thresholds_bounds::Tuple{T, T}=(5e-4, 1e-5): Thresholds that determine the different behaviour of bounds_and_cutoff.
  • sbound::T=0.01: Fraction of 'R' used for very small inflation factor; see bounds_and_cutoff.

Arguments that determine produce_jammed_configuration! termination criteria

The list of default values is specified in this part of the documentation.

  • max_iters::I=1000: Maximum number iterations of the main loop; that is, the maximum number of LP optimizations allowed.
  • tol_Γ_convergence::T=default_tol_Γ_convergence: determines the convergence criterion of the packing fraction as $\sqrt{\Gamma^\star}-1 \leq$ tol_Γ_convergence.
  • tol_S_convergence::T=default_tol_displacements_convergence: determines the convergence criterion for the displacements as $\max |\mathbf{s}_{i,\mu}^\star|_{i=1,\dots,N}^{\mu=1,\dots, d} \leq$ <=tol_S_conv.
  • non_iso_break::I=max_iters: Number of consecutive non-isostatic solutions allowed before produce_jammed_configuration! terminates. The reason is that it is very likely that the final configuration will also be non-isostatic (specially if beginning from a highly compressed state). Note however that every time an isostatic configuration is obtained, this counter resets to 0.

Arguments for controlling the behaviour of the solver/optimizer

  • solver::Module=GLPK: The solver (i.e. the package or library) employed used by JuMP to solve each LP instance.
  • solver_attributes::Dict=default_solver_attributes: The attributes used by the solver. It's used to control some of its features, such as precision, iteration limits, etc. But it depend on which solver is used.
  • solver_args=default_args: Arguments passed to the solver.Optimizer function. It is also used to control the parameters of the optimizer.

More detailed information is available in this part of the documentation.

Arguments to control the tolerance of mechanical equilibrium, overlaps identification, etc.

The list of default values is specified in this part of the documentation.

  • tol_mechanical_equilibrium::Float64=default_tol_force_equilibrium: tolerance to test whether a packing satisfies the force balance condition.
  • zero_force::T=default_tol_zero_forces: threshold for identifying a force as non-zero.
  • tol_overlap::T=default_tol_overlap:: tolerance with which overlaps are identified. That is, if particles are overlapping, but are doing so by a quantity smaller than default_tol_overlap (default 1e-8), no error is thrown. This rather loose tolerance is related to the maximal precision available with Gurobi.

Arguments controlling the screen printed output

  • verbose::Bool=true: Control whether some info about the progress of CALiPPSO and the final packing is printed out or not.
  • non_iso_warn::Bool=false: Print some information whenever a non-isostatic preliminary configuration is obtained, after solving a LOP.
  • monitor_step::I=10: How often info about the progress in the main main loop should be printed out; will only take effect if verbose=true. See print_monitor_progress for more information.
  • initial_monitor::I=monitor_step: print info about the main loop progress during this amount of initial LP optimizations; will only take effect if verbose=true.

Arguments for performing overlaps checks and updating list of distances

See check_for_overlaps and update_distances! for more information

  • interval_overlaps_check::I=10: interval of LP optimizations at which it is verified that no overlaps are present in the system.
  • initial_overlaps_check::I=initial_monitor: number of initial LP optimizations at which it is verified that no overlaps are present; given that for low density initial configurations, CALiPPSO might produce rather large displacements, it is always convenient to keep track of the overlaps in such initial stage.
  • ratio_sℓ_update::T=0.1: This quantity determines the threshold (s_update) after which the list of distances of a particles is updated. Such threshold is calculated as s_update=ratio_sℓ_update*ℓ/sqrt(d), with determined by bounds_and_cutoff.
source
CALiPPSO.network_of_contactsFunction
network_of_contacts(Xs::Vector{SVector{d, PeriodicNumber{T}}}, constraints::Vector{Vector{ConstraintRef}}, neighbours_list::Vector{Vector{Int64}}, images::Vector{SVector{d, T}}) where {d, T<:Float64}

Construct the set of contact vectors, forces magnitudes and list of interacting particles of the full configuration.

In a first loop, the pairs of "interacting" particles are identified from the non-zero dual variables associated to the 'constraints' that are saturated, once the LP instance has been optimized. So every time shadow_price returns a non-zero value, an ordered pair of interacting particles is stored, with the corresponding value of the dual variable identified as the contact force. Finally, the corresponding contact vector for the pair is obtained by calling MIC_vector on the positions ('Xs') of the particles involved. This is why 'images' should also be provided so the MIC contact vector is computed more rapidly.

In the second loop the contact vectors, forces, etc. of the complementary pairs are stored. That is, if in the first loop only pairs (i,j) with i<j are considered, in the second one we assume i>j.

Output

  • all_contact_vectors: A vector whose elements are vectors containing SVector{d,T} entries. So the i-th element is the set of contact vectors of the i-th particle.
  • forces_dual: A Vector{Vector{Float64}} containing the forces magnitudes acting on each particle. So its i-th element is the list of the magnitude of the forces acting on particle i.
  • particles_dual_contact: A Vector{Vector{Int64}} containing the indices of particles in contact with each particle. So its i-th element is the list of indices of particles in contact with the i-th particle.

See also add_non_overlapping_constraints!, @constraint, optimize!, solve_LP_instance!.

source
network_of_contacts(packing::AbstractPacking{d, T}, normalized::Bool=true) where {d, T<:Float64}

Obtain the list of contact indices (as ordered pairs, i.e. [i, j] with j>i), the corresponding contact vectors, and magnitudes of contact forces, from a given 'packing'.

This function is not used in the main CALiPPSO function (i.e. produce_jammed_configuration!), but is useful for extracting the system's micro-structural information once the jammed packings have been generated.

source
CALiPPSO.check_for_overlapsFunction
check_for_overlaps(Xs::Vector{SVector{d, PeriodicNumber{T}}}, R::T, tolerance::T)

Check whether there is an overlap between all pairs of particles centred at 'Xs' and of radius 'R'.

Each element of 'Xs' plays the role of the position of a particle's center, and all of them are assumed to be of the same size, i.e. radius 'R'. A tolerance to determine whether there is an overlap or not should be passed as third argument.

Output

  • 'overlap': a boolean that is 'true' only when an overlap is present
  • 'message': a string that contains some information about the overlapping particles.
  • 'particles': a tuple containing the indices of the overlapping particles; '(0, 0)' if no overlap
source
check_for_overlaps(Xs::Vector{SVector{d, PeriodicNumber{T}}}, Rs::Vector{T}, tolerance::T)

Check whether there is an overlap between all pairs of particles centred at 'Xs' and of radii 'Rs'.

Each element of 'Xs' plays the role of the position of a particle's center, while its associated radius is stored in corresponding element of 'Rs' A tolerance to determine whether there is an overlap or not should be passed as third argument.

Output

  • 'overlap': a boolean that is 'true' only when an overlap is present
  • 'message': a string that contains some information about the overlapping particles.
  • 'particles': a tuple containing the indices of the overlapping particles; '(0, 0)' if no overlap
source
check_for_overlaps(packing::MonoPacking, tolerance::AbstractFloat)

Apply check_for_overlaps to all the particles in 'packing'.

source
check_for_overlaps(packing::PolyPacking{d,T}, tolerance::T)

Apply check_for_overlaps to all the particles in 'packing'.

source

Test whether there is an overlap between particles. Also info about displacements (in the possible overlapping pair) is given.

source
check_for_overlaps(packing::AbstractPacking, t::Int64, possible_neighbours::Vector{Vector{Int64}}, jammed::Bool; tolerance=default_tol_overlap)

Check for overlaps in all the particles of a given packing, obtained after CALiPPSO converged.

source
CALiPPSO.generate_random_configurationFunction
generate_random_configuration(d::Int64, N::Int64, ϕ::T, L::T=1.0; max_tries::Int64=5000 )

Generate 'N' random centers of monodisperse particles of radius 'r' and in 'd' dimensions, without any overlaps.

The output is the radius (that corresponds to the packing fraction 'ϕ' used as input), and the vector containing the centers (each as a SVector{d, PeriodicNumber} type). Each center is placed at uniformly random in all space, and when an overlap is detected a new random position is drawn. 'max_tries' (default 5000) attempts are tried for each particle and when this bound surpassed an error is thrown, with the index of the center that was not created.

source
CALiPPSO.total_forceFunction

Compute the sum of forces acting on a given particle. The output is a StaticVector

source

Output an array of StaticVectors corresponding to the sum of forces acting on each particle of the packing.

source

Compute the sum of forces on each particle from the full set of contact vectors and magnitudes of contact forces.

source
CALiPPSO.packing_fractionFunction
packing_fraction(d::Int64, R::AbstractFloat, N::Int64, L=1.0)

Compute the packing fraction of N d-dimensional hyperspheres of the same radius, 'R', inside a box of size 'L'.

See also volume, volume_d_ball

source
packing_fraction(packing::MonoPacking{d, T})

Compute the packing fraction of a monodisperse packing.

See also volume, volume_d_ball

source
packing_fraction(d::Int64, R1::T, N1::Int64, R2::T, N2::Int64, L::T=1.0) where T<:AbstractFloat

Compute the packing fraction of a configuration of d-dimensional hyperspheres, of which 'N1' have radius 'R1', and 'N2' have radius 'R2'. The configuration is assumed to be inside a box of size 'L'.

See also volume, volume_d_ball

source
packing_fraction(d::Int64, Rs::Vector{T}, L::T=1.0) where T<:AbstractFloat

Compute the packing fraction of N d-dimensional hyperspheres with radii 'Rs', inside a box of size 'L'.

See also volume, volume_d_ball

source
packing_fraction(packing::PolyPacking{d, T})

Compute the packing fraction of a monodisperse packing.

See also volume, volume_d_ball

source

Other exported functions

CALiPPSO.PeriodicVectorsFunction
PeriodicVectors(mat::Matrix{T}, L::T=1.0) where {T<:AbstractFloat}

Convert 'mat' to a Vector of elements of type SVector{d, PeriodicNumber{T}}.

mat should be of size d x N. The output is a 1-dimensional array of N elements, each of which consists of StaticVector's of size d, and PeriodicNumber{T} as data. The periodicity of the numbers, 'L', is given as second argument and defaults to 1.0.

source
CALiPPSO.volume_d_ballFunction
volume_d_ball(d::Int64, R::Real)

Compute the volume of a d-dimensional sphere of radius R (that defaults to 1)

source

Constructors

CALiPPSO.PeriodicVectorFunction
PeriodicVector(vec::Vector{T}, L::T=1.0)

Convert 'vec' to a Static Vector of the same size, but with elements of 'PeriodicNumber{T}' type.

The periodicity of the numbers, 'L', is given as second argument and defaults to 1.0.

source
CALiPPSO.MonoPackingMethod
MonoPacking(Xs::Vector{SVector{d, PeriodicNumber{T}}}, contact_vecs::Vector{Vector{SVector{d,T}}}, fs::Vector{Vector{T}}, neighbours::Vector{Vector{Int64}}, R::T, jammed::Bool=false; <keyword arguments>)
MonoPacking(Xs::Vector{SVector{d, PeriodicNumber{T}}}, contact_vecs::Vector{Matrix{T}}, fs::Vector{Vector{T}}, neighbours::Vector{Vector{Int64}}, R::T, jammed::Bool=false; <keyword arguments>)

Create a d-dimensional packing, where the attributes of each particles are inferred from the elements of the input arrays.

The number of particles (N) is inferred from the length of the input arrays. Then a Vector{MonoParticle{d,T}} of size N but undefined elements is constructed. Each of its elements is then defined by calling 'MonoParticle(Xs[i], contact_vecs[i], fs[i], neighbours[i])'.

The constructors also asses whether force balance for each particle is satisfied, within a given precision 'tolmechanicalequilibrium' (that defaults to default_tol_force_equilibrium=1e-12). When this condition is not met, it throws a warning, but the packing is created.

source
CALiPPSO.MonoPackingMethod
MonoPacking(Xs::Vector{SVector{d, PeriodicNumber{T}}}, constraints::Vector{Vector{ConstraintRef}}, neighbours_list::Vector{Vector{Int64}}, R::T, images::Vector{SVector{d, T}}, jammed::Bool=false; <keyword arguments>) where {d, T<:AbstractFloat}

Construct a MonoPacking from the set of particles' position ('Xs'), set of all constraints defined in the LP model ('constraints'), list of possible neighbours ('neighbours_list'), and virtual images ('images') needed for the MIC contact vectors. "

source

Secondary functions

Missing docstring.

Missing docstring for CALiPPSO.solve_LP_instance. Check Documenter's build log for details.

CALiPPSO.fine_tune_forces!Function
fine_tune_forces!(Packing::MonoPacking{d, T}, force_mismatch::T, sqrΓ::T, ℓ0::T, images::Vector{SVector{d, T}};      
<keyword arguments> ) where {d, T<:Float64}

Update the forces of all particles in 'Packing' so the global force balance condition is more closely satisfied.

This function is needed because when CALiPPSO has converged (or terminated due to reaching maximum iterations) it might happen that mechanical equilibrium is rather inaccurate. This is caused by the displacements in the last LP optimization. So the total force on particles whose position was updated is about as large as the convergence tolerance on the displacements. Thus, to produce a packing that satisfies force balance much more precisely an extra LP optimization is done (so the dual variables are recalculated with the final positions), but the particles' position and radius are NOT updated. Also importantly, to avoid the introduction of extra constraints that might ruin force balance (mainly when CALiPPSO terminated without producing a jammed packing), the optimization performed here is done without bounding the displacements.

Output

Besides updating the forces of each particle in 'Packing', this function produces the following output

  1. t_solve: The time elapsed during the extra LP optimization.
  2. isostatic: A boolean that specifies whether the updated packing is isostatic or not.
  3. Nc: The updated number of contacts.
  4. Nnr: The updated number of non-rattlers.

Keyword arguments

  • solver::Symbol=:Gurobi: the solver used to call the optimizer when creating the JuMP model.
  • solver_attributes::Dict=default_solver_attributes: The attributes (i.e. parameters) passed to the solver after creating model, using set_optimizer_attributes.
  • solver_args=default_args: The arguments passed to Optimizer of the chosen solver. It should be either nothing or a NamedTuple. Choose the former if testing a solver other than Gurobi, GLPK, or Hypatia.
  • thresholds::Tuple{T, T}=(5e-4, 1e-5): thresholds that define the different criteria to determine the radius of influence, ℓ, and the displacements' bound when calling bounds_and_cutoff.
  • tol_mechanical_equilibrium=default_tol_force_equilibrium: The tolerance to determine whether force balance is fulfilled in each particle.

See also add_non_overlapping_constraints!, optimize!, produce_jammed_configuration!, solve_LP_instance!.

source
CALiPPSO.update_packing_forces!Function
update_packing_forces!(Packing::MonoPacking{d,T}, constraints::Vector{Vector{ConstraintRef}}, neighbours_list::Vector{Vector{Int64}}, images::Vector{SVector{d, T}}; tol_mechanical_equilibrium::Float64=default_tol_force_equilibrium) where {d, T<:Float64}
update_packing_forces!(Packing::PolyPacking{d,T}, constraints::Vector{Vector{ConstraintRef}}, neighbours_list::Vector{Vector{Int64}}, images::Vector{SVector{d, T}}; tol_mechanical_equilibrium::Float64=default_tol_force_equilibrium) where {d, T<:Float64}

Update the forces field of each MonoParticle in 'Packing', (or the analogous field of each Particle in case of a polydispersity) from a new set of 'constraints'.

The arguments needed are the ones needed to call network_of_contacts, which is the main function used here. The kwarg 'tolmechanicalequilibrium' defines the tolerance to assess whether the force balance condition is satisfied in each particle, and consequently update the mechanical_equilibrium field of 'Packing'.

See also network_of_contacts, MonoPacking, PolyPacking, total_force, fine_tune_forces!.

source
CALiPPSO.add_non_overlapping_constraints!Function
add_non_overlapping_constraints!(model::JuMP.Model, Xs::Vector{SVector{d, PeriodicNumber{T}}}, R::T, ℓ::T, images::Vector{SVector{d, T}})
add_non_overlapping_constraints!(model::JuMP.Model, Xs::Vector{SVector{d, PeriodicNumber{T}}}, Rs::Vector{T}, ℓ::T, images::Vector{SVector{d, T}})

Assign the linearized non-overlapping constraints to 'model', according to particle positions, 'Xs', and radius, 'R' (or radii 'Rs' in polydisperse systems).

A pair of particles is included in the set of constraints only if their distance is smaller than the cutoff 'ℓ'. Because the MIC distance is considered, also the set of virtual images 'images' should also be provided as input.

Output

  1. 'constraints': a Vector whose elements are themselves arrays of 'ConstraintRef' type.
  2. nearby_particles_list: a vector whose i-th entry is the list of particles' indices that induce a constraint on particle i.

Note that the constraint associated to the pair (i,j) is only counted once (is associated to 'i' if j>i; or to 'j' otherwise). Thus the i-th entry of the output arrays only contains constraints and indices associated to particles of index greater than 'i'.

See also MIC_distance, MIC_vector, @constraint, solve_LP_instance!, fine_tune_forces!.

source
CALiPPSO.bounds_and_cutoffFunction
bounds_and_cutoff(sqrΓ::T, R::T, ℓ0::T, d::Int64; thresholds::Vector{T}=[5e-4, 1e-5], sbound::T=0.01)

Obtain displacements' bound and radius of influence (or cutoff) for particles of radius 'R' and after an inflation factor with square root 'sqrΓ'. In case of a polydisperse packing the largest radius is used.

Output

  1. The cutoff for assigning constraints, i.e. the radius of influence, ℓ.
  2. The bounds to be imposed on the each component of the displacements.

The value of the displacements' bound and ℓ are chosen depending on how close the input value of the growth factor, Γ (or better √Γ), is to 1. Thus, as explained in the paper, even when the final φ_J is not known, the displacements' bounds and ℓ are estimated with the value of Γ from the previous LP optimization (or an initial suitable guess if needed). Three different criteria are selected for large, small, and very small values of sqrΓ-1, respectively

The other arguments are the current value of 'R', an upper bound for the radius of influence 'ℓ0', and the dimensionality of the system, 'd'. Keywords arguments are used to control when different criteria are triggered 'thresholds', and the fraction of 'R' used as displacement bound when Γ is already very close to 1, 'sbound' (default 0.01).

source
CALiPPSO.MIC_vectorFunction
MIC_vector(V1::SVector{d, PeriodicNumber{T}}, V2::SVector{d, PeriodicNumber{T}}, images::Vector{SVector{d, T}})

Compute the vector joining V2 to V1 (i.e. V1-V2) considering periodic boundary conditions and using the so called Minimum Image Convention (MIC).

To avoid unnecessary calls in other functions, the system's set of virtual images should be given as third argument, as an array of SVector's. It is assumed that both vectors are contained in a d-dimensional box of size L. The MIC guarantees that all the possible periodic shifts or 'virtual images' of the vectors are considered when calculating their distance.

The output is the MIC vector (as a SVector type) and the index of the image of minimum distance (0 if no such image is needed).

See also: MIC_distance.

source
CALiPPSO.MIC_distanceFunction
MIC_distance(V1::SVector{d, PeriodicNumber{T}}, V2::SVector{d, PeriodicNumber{T}}, images::Vector{SVector{d, T}})

Compute the distance between 'V1' and 'V2' considering periodic boundary conditions and using the so called Minimum Image Convention (MIC).

To avoid unnecessary calls in other functions, the system's set of virtual images should be given as third argument, as an array of SVector's. It is assumed that both vectors are contained in a d-dimensional box of size L. The MIC guarantees that all the possible periodic shifts or 'virtual images' of the relative vectors are considered when calculating their distance.

The output is the index of the image of minimum distance (0 if no shift is needed) and the value of such distance.

See also: MIC_vector.

source
CALiPPSO.volumeFunction
volume(P::Particle)

Compute the volume of 'P' (i.e. a hypersphere of d dimensions and radius R=P.R.

source
CALiPPSO.update_distances!Function
function update_distances!(S⃗_cum::Vector{SVector{d, T}}, s_update::T, Xs::Vector{SVector{d, PeriodicNumber{T}}}, distances::Symmetric{T, Matrix{T}}; verbose=Bool=true) where {d, T<:AbstractFloat}

Update (in place) the list of distances between pairs of particles, depending on the cumulative displacements S⃗_cum and a threshold s_update. If less than 10% of the particles have a displacement, whose norm is smaller than s_update, only the list of distances of such particles is recalculated. Otherwise, the distances between of all pairs is recomputed.

For convenience, if this function is called with s_update=0.0, the distances are always recalculated, without the need to first find the particles which have moved a lot. This feature might be useful if you want to be sure that before each LP instance all the relevant constrictions are considered.

source

Functions for printing CALiPPSO info and progress

CALiPPSO.print_monitor_progressFunction

Print info about CALiPPSO progress, including new values of $|s_i|$, density and $R$, statistics of constraints and contacts, forces, etc.

source

Print info about CALiPPSO progress, including new values of $|s_i|$, density and some statistics of $R$, statistics of constraints and contacts, forces, etc.

source
CALiPPSO.print_info_convergenceFunction
print_info_convergence(final_packing::MonoPacking, isostatic::Bool, time, memory; digs::Int64=2)

Print extra details about the configuration obtained once CALiPPSO converges.

The maximum mismatch in the force balance condition is computed (and printed), and also the isostaticity of the packing is assessed. Besides, some information about performance (execution time and allocated memory) is also printed out.

See also print_converged, print_monitor_progress.

source
CALiPPSO.print_non_isostaticFunction

Print information about CALiPPSO progress and current status of the system whenever a non-isostatic configuration is created.

source