Random rough surfaces

The generation of stochatsticly rough surfaces is controlled in Tamaas by two abstract classes: tamaas::SurfaceGenerator and tamaas::Filter. The former provides access lets the user set the surface sizes and random seed, while the latter encodes the information of the spectrum of the surface. Two surface generation methods are provided:

Both of these rely on a tamaas::Filter object to provided the filtering information (usually power spectrum density coefficients). Tamaas provides two such classes and allows for Python subclassing:

Tamaas also provided routines for surface statistics.

Generating a surface

Let us now see how to generate a surface. Frist create a filter object and set the surface sizes:

import tamaas as tm

# Create spectrum object
spectrum = tm.Isopowerlaw2D()

# Set spectrum parameters
spectrum.q0 = 4
spectrum.q1 = 4
spectrum.q2 = 32
spectrum.hurst = 0.8

The spectrum object can be queried for information, such as the root-mean-square of heights, the various statistical moments, the spectrum bandwidth, etc. Then we create a generator object and build the random surface:

generator = tm.SurfaceGeneratorFilter2D([128, 128])
generator.spectrum = spectrum
generator.random_seed = 0

surface = generator.buildSurface()

Important

The surface object is a numpy.ndarray wrapped around internal memory in the generator object, so a subsequent call to buildSurface may change its content. Note that if generator goes out of scope its memory will not be freed if there is still a live reference to the surface data.

Important

If ran in an MPI context, the constructor of SurfaceGeneratorFilter2D (and others) expects the global shape of the surface. The shape can also be changed with generator.shape = [64, 64].

It is common to normalize a surface after it has been generated so that one can scale the surface to a desired stastistic, e.g. to specify the root-mean-square of slopes:

rms_slopes = 0.25
surface /= iso.rmsSlopes()
surface *= rms_slopes

# Compute root-mean-square of slopes in Fourier domain
print(tm.Statistics2D.computeSpectralRMSSlope(surface))
# Compute root-mean-square of slopes with finite differences
# (introduces discretization error)
print(tm.Statistics2D.computeFDRMSSlope(surface))

Note

The spectrum object gives the expected value of surface statistics, but the corresponding value for a surface realization can differ (particularily if the SurfaceGeneratorFilter2D is used). One could normalize a surface with the actual statistic, but this leads to biased quantities when a representative sample of surfaces is used.

Custom filter

Tamaas provides several classes that can be derived directly with Python classes, and tamaas::Filter is one of them. Since it provides a single pure virtual method computeFilter, it is easy to write a sub-class. Here we implement a class that takes a user-defined auto-correlation function and implements the computeFilter virtual function:

import numpy

class AutocorrelationFilter(tm.Filter2D):
    def __init__(self, autocorrelation):
        tm.Filter2D.__init__(self)
        self.autocorrelation = autocorrelation.copy()

    def computeFilter(self, filter_coefficients):
        shifted_ac = numpy.fft.ifftshift(self.autocorrelation)

        # Fill in the PSD coefficients
        filter_coefficients[...] = numpy.sqrt(np.fft.rfft2(shifted_ac))
        # Normalize
        filter_coefficients[...] *= 1 / numpy.sqrt(self.autocorrelation.size)

Here filter_coefficients is also a numpy.ndarray and is therefore easily manipulated. The creation of the surface then follows the same pattern as previously:

# Create spectrum object
autocorrelation = ...  # set your desired autocorrelation
spectrum = AutocorrelationFilter(autocorrelation)

generator = tm.SurfaceGenerator2D()
generator.shape = autocorrelation.shape
generator.spectrum = spectrum

surface = generator.buildSurface()

The lifetime of the spectrum object is associated to the generator when setSpectrum is called.

Surface Satistics

Tamaas provides the C++ class tamaas::Statistics and its wrapper Statistics2D to compute statistics on surfaces, including:

  • power spectrum density

  • autocorrelation

  • spectrum moments

  • root-mean-square of heights \(\sqrt{\langle h^2 \rangle}\)

  • root-mean-square of slopes (computed in Fourier domain) \(\sqrt{\langle |\nabla h|^2\rangle}\)

All these quantities are computed in a discretization-independent manner: increasing the number of points in the surface should not drastically change the computed values (for a given spectrum). This allows to refine the discretization as much as possible to approximate a continuum. Note that the autocorrelation and PSD are fft-shifted. Here is a sample code plotting the PSD and autocorrelation:

psd = tm.Statistics2D.computePowerSpectrum(surface)
psd = numpy.fft.fftshift(psd, axes=0) # Shifting only once axis because of R2C transform

import matplotlib.pyplot as plt
from matplotlib.colors import LogNorm

plt.imshow(psd.real, norm=LogNorm())

acf = tm.Statistics2D.computeAutocorrelation(surface)
acf = numpy.fft.fftshift(acf)

plt.figure()
plt.imshow(acf)

plt.show()

See examples/statistics.py for more usage examples of statistics.