A module with procedures for calculating various properties which are interesting when doing numerical percolation experiments.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=dp), | public, | parameter | :: | pc | = | 0.592746 | Known critical probability for a two-dimensional site-percolating system. |
Create a random, binary (logical
) matrix,
which can be used in percolation experiments.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp), | intent(in) | :: | p | Probability for each matrix element to be |
||
integer, | intent(in) | :: | L |
Randomly created matrix, where each element is .true.
or .false.
if a randomly generated number is smaller
or greater than p.
Count the number of sites belonging to each cluster in the labelled matrix.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in), | dimension(:,:) | :: | labelled_matrix | Integer matrix with labelled clusters, resulting from the Hoshen-Kopelman algorithm. |
|
integer, | intent(in) | :: | num_labels | The known number of clusters. |
Array of cluster sizes. The i'th element is the size of the cluster with label i.
Find the label of the percolating cluster, i.e. the one spanning from one side of the system to the opposing side.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in), | dimension(:,:) | :: | labelled_matrix | Labelled matrix of clusters from hoshen_kopelman/label. |
|
integer, | intent(in) | :: | num_labels | Known number of clusters. |
Label of the percolating cluster. -1 if no percolating cluster is found.
Density of the spanning/percolating cluster, i.e. the number of sites on the percolating cluster divided by \(L^2\). Averaged over num_samples Monte Carlo samples (with OpenMP).
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp), | intent(in) | :: | p | The probability for a site to allow transport. |
||
integer, | intent(in) | :: | L | The size of the system. |
||
integer, | intent(in) | :: | num_samples | The number of Monte Carlo samples. |
The number of sites on the spanning cluster divided by \(L^2\).
Calculate the probability of having a spanning/percolating cluster, given a system size L and probability for a site to have transport p.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp), | intent(in) | :: | p | Probability for a site to allow transport. |
||
integer, | intent(in) | :: | L | Size of the system. |
||
integer, | intent(in) | :: | num_samples | Number of Monte Carlo samples. |
The probability of having a percolating cluster, calculated as the number of times a percolating cluster is found, divided by the number of attempts (num_samples).
Find the inverse of spanning_probability by use of the bisection method.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp), | intent(in) | :: | x | The value of spanning_probability for which the inverse is calculated. |
||
integer, | intent(in) | :: | L | Size of the system. |
||
integer, | intent(in) | :: | num_samples | Number of Monte Carlo samples to use when evaluating spanning_probability. |
||
real(kind=dp), | intent(in) | :: | tolerance | Tolerance of approximation. The return value is within tolerance/2 of the correct (but numerical) value. |
The inverse of spanning_probability. If spanning_probability is denoted \(\Pi(p,L)\), this function returns \(p\) such that \(\Pi(p,L)=x\).
Alternative interface to the Hoshen-Kopelman algorithm from the hk module, which uses a binary matrix created by create_binary_matrix.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
logical, | intent(in), | dimension(:,:) | :: | matrix | Binary matrix where clusters will be identified. |
|
integer, | intent(out), | dimension(:,:), allocatable | :: | labelled_matrix | Integer matrix which will store the labels of each cluster. Reallocated if necessary. |
|
integer, | intent(out) | :: | num_labels | Overwritten with the total number of disjoint clusters. |
Calculate the number density of clusters with different sizes. The cluster number density is defined as \begin{equation} n(s,p) = \sum_{\text{MC-samples}} \frac{\text{number of clusters with size }s} {L^2\cdot\text{number of MC-samples}}. \end{equation} Direct calculations will usually give bad results, as there will be very few clusters with large sizes compared to the numbers of clusters with small sizes. This is circumvented by doing logarithmic binning and averaging, i.e. \begin{equation} n\left([s+\Delta s),p\right) = \sum_{\text{MC-samples}} \frac{\text{number of clusters with size }s \in[s,s+\Delta s)} {\Delta s\cdot L^2\cdot\text{number of MC-samples}}, \label{eq:nsp} \end{equation} where \(\Delta s\) are logarithmically distributed. After execution, bin_mids will contain the centres of the bins.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=dp), | intent(in) | :: | p | Probability for a given site to allow transport. |
||
integer, | intent(in) | :: | L | Percolating systems will be \(L\times L\). |
||
integer, | intent(in) | :: | num_samples | Results will be averaged over this number of Monte Carlo-samples. Sampling is parallelised. |
||
real(kind=dp), | intent(out), | dimension(:), allocatable | :: | bin_mids | Centres of bins. |
|
real(kind=dp), | intent(out), | dimension(:), allocatable | :: | results | Cluster number density in \eqref{eq:nsp}. |
|
real(kind=dp), | intent(in), | optional | :: | binsize_base | The edges of the logarithmically distributed bins will be integer powers of this number. Default: 1.5. |