snowfake package#

Submodules#

snowfake.snowfake module#

Make Gravner-Griffeath snowfakes.

Implementing: Janko Gravner, David Griffeath, Modeling snow crystal growth II: A mesoscopic lattice map with plausible dynamics, Physica D: Nonlinear Phenomena, Volume 237, Issue 3, 2008, Pages 385-404, https://doi.org/10.1016/j.physd.2007.09.008.

class snowfake.snowfake.Snowfake(size, random=False, **kwargs)[source]#

Bases: object

Simulates a snowfake (fake snowflake), with a slow implementation of Gravner & Griffeath’s algorithm. This class holds all of the state, so you can grow it in stages, eg changing some of the parameters at various times.

Example to reproduce Figure 15b in Gravner & Griffeath (2008):
>>> from snowfake import Snowfake
>>> params =   {
    'ρ': 0.35,
    'β': 1.4,
    'α': 0.001,
    'θ': 0.015,
    'κ': 0.05,
    'μ': 0.015,
    'γ': 0.01,
    'σ': 0.00005,
    'random': False,
}
>>> s = Snowfake(801, **params)
>>> s.grow()
>>> s.plot()
BDY = array([[ 0,  1,  1],        [ 1, -7,  1],        [ 1,  1,  0]])#
NBR = array([[0, 1, 1],        [1, 1, 1],        [1, 1, 0]])#
attachment()[source]#

“This key step in the algorithm decides when a boundary site joins the snowflake.”

What happens depends on how many neighbours a site has.

property boundary#

Updates the boundary sites, if they have changed.

diffusion()[source]#

“Diffusive mass evolves on \(\mathbf{A^c_t}\) (non-crystal) by discrete diffusion with uniform weight 1/7 on the center site and each of its neighbors […] and for \(\mathbf{x ∈ ∂A_t}\) (boundary) any term in the sum corresponding to \(\mathbf{y ∈ ∂A_t}\) is replaced by \(\mathbf{d°(x)}\) (i.e. the value before step).”

freezing()[source]#

“Proportion κ of the diffusive mass at each boundary site crystallizes. The remainder (proportion \(\mathbf{1 - κ}\)) becomes boundary mass.” We control the ‘boundaryness’ with the boundary array.

grow(max_epochs=20000, early_stopping=0.5, snapshot=1000)[source]#

Iterate over the diffusion-freezing-attachment-melting-perturbation cycle, up to max_epochs times. If the flake has previously been grown, calling this method again will continue to grow it. If early_stopping is False, then every epoch will run.

If early_stopping is non-zero (or not False), the snowflake will automatically stop growing soon after it is the given proportion of the width of the space. The default value of 0.5 has proven to be a reasonable heuristic, allowing the flake to grow to a reasonable size, without boundary effects. If you pass True for this parameter, the default value of 0.5 will be used. The growth loop checks for early stopping every 100 epochs.

If snapshot is non-zero (or not False) then a plot of the current state will be produced every N epochs, where N = snapshot. This is on by default and produces a plot every 1000 epochs.

Parameters
  • max_epochs (int) – The maximum number of epochs to run.

  • early_stopping (float) – If False, every epoch will run, even if the flake reaches the boundary of the space. If a float between 0 and 1, the flake will stop growing if it is this proportion of the width of the space. If True, the default of 0.5 is used.

  • snapshot (int) – If False, no plots will be produced. If a positive integer, a plot will be produced every N epochs. Default: 1000.

Returns

None.

melting()[source]#

“Proportion µ of the boundary mass and proportion γ of the crystal mass at each boundary site become diffusive mass. […] Typically µ is small and γ extremely small.”

property neighbours#

Return the number of neighbours each cell has.

property params#

Return a dictionary of the snowfake’s parameters.

perturb()[source]#

“The diffusive mass at each site undergoes a […] perturbation of proportion σ”

In the paper, both random and biased perturbations are used.

plot(ax=None, cmap='Blues')[source]#

Not sure what the prettiest way to plot these things is. For now we’ll do something very simple, which looks a little bit like what the authors do in their paper: Plot the sum of the d and c arrays.

Parameters
  • ax (Axes) – The axes to plot on. If None, a new figure is created.

  • cmap (str) – The colormap to use.

Returns

The axes that were plotted on.

Return type

Axes

rectify(prop='c')[source]#

Transform from hexagonal to Cartesian coordinates with scipy.

Parameters

prop (str) – The array to transform.

Returns

The transformed array.

Return type

ndarray

rectify_skimage(prop='c')[source]#

Transform from hexagonal to Cartesian coordinates using skimage.

Parameters

prop (str) – The array to transform.

Returns

The transformed array.

Return type

ndarray

snapshot()[source]#

Quick look at the arrays during training.

status()[source]#

Return a string representation of the snowfake’s current state.

snowfake.snowfake.main()[source]#

The function that will run if this module is called from the command line.

Not implemented yet.

snowfake.snowfake.random(size=801, seed=None, **kwargs)[source]#

Returns an intialized snowfake with random parameters. Any parameters you pass will override the random ones.

Parameters
  • size (int) – The size in hex pixels of the 2D space.

  • seed (int) – The random seed.

  • kwargs – Any other parameters you want to pass to the snowfake.

Returns

A snowfake with random parameters.

Return type

Snowfake

Module contents#