Correlated roughness

Scattering from a multilayered sample with correlated roughness.

  • The sample is composed of a substrate on which is sitting a stack of layers. These layers consist in a repetition of 5 times two different superimposed layers (from bottom to top):

- layer A: 2.5 nm thick with a real refractive index n = 1-5e-6.
- layer B: 5 nm thick with a real refractive index n = 1-1e-5.

  • There is no added particle.
  • All layers present the same type of roughness on the top surface, which is characterized by:

- a rms roughness of the interfaces σ =1 nm,
- a Hurst parameter H equal to 0.3,
- a lateral correlation length ξ of 5 nm,
- a cross correlation length ξ equal to 1e-4 nm.

  • The incident beam is characterized by a wavelength of 1 Å.
  • The incident angles are αi = 0.2° and Φi = 0°.
Note:
The roughness profile is described by a normally-distributed random function. The
roughness correlation function at jth interface is expressed as
<Uj (x, y)Uj (x', y')>= σ2exp(-(τ/ξ)2H), τ= [(x-x')2+(y-y')2]1/2,
where Uj(x, y) is the height deviation of the jth interface at position (x, y).

σ gives the rms roughness of the interface.
The Hurst parameter H, comprised between 0 and 1 is connected to the fractal dimension
D=3-H of the interface. The smaller H is, the more serrate the surface profile
looks. If H = 1, the interface has a non fractal nature.

The lateral correlation length ξ acts as a cut-off for the lateral length scale on which an
interface begins to look smooth. If ξ » τ the surface looks smooth.

The cross correlation length ξis the vertical distance over which the correlation between layers is damped by a factor 1/e. It is assumed to be the same for all interfaces.
If ξ = 0 there is no correlations between layers. If ξ is much larger than the layer
thickness, the layers are perfectly correlated.
Real-space model: 
Intensity Image: 
Python Script: 
"""
MultiLayer with correlated roughness
"""
import numpy
import bornagain as ba
from bornagain import deg, angstrom, nm

phi_min, phi_max = -0.5, 0.5
alpha_min, alpha_max = 0.0, 1.0


def get_sample():
    """
    Returns a sample with two layers on a substrate, with correlated roughnesses.
    """
    # defining materials
    m_ambience = ba.HomogeneousMaterial("ambience", 0.0, 0.0)
    m_part_a = ba.HomogeneousMaterial("PartA", 5e-6, 0.0)
    m_part_b = ba.HomogeneousMaterial("PartB", 10e-6, 0.0)
    m_substrate = ba.HomogeneousMaterial("substrate", 15e-6, 0.0)

    # defining layers
    l_ambience = ba.Layer(m_ambience)
    l_part_a = ba.Layer(m_part_a, 2.5*nm)
    l_part_b = ba.Layer(m_part_b, 5.0*nm)
    l_substrate = ba.Layer(m_substrate)

    roughness = ba.LayerRoughness()
    roughness.setSigma(1.0*nm)
    roughness.setHurstParameter(0.3)
    roughness.setLatteralCorrLength(5.0*nm)

    my_sample = ba.MultiLayer()

    # adding layers
    my_sample.addLayer(l_ambience)

    n_repetitions = 5
    for i in range(n_repetitions):
        my_sample.addLayerWithTopRoughness(l_part_a, roughness)
        my_sample.addLayerWithTopRoughness(l_part_b, roughness)

    my_sample.addLayerWithTopRoughness(l_substrate, roughness)
    my_sample.setCrossCorrLength(1e-4)

    print(my_sample.treeToString())

    return my_sample


def get_simulation():
    """
    Characterizing the input beam and output detector
    """
    simulation = ba.GISASSimulation()
    simulation.setDetectorParameters(200, phi_min*deg, phi_max*deg,
                                     200, alpha_min*deg, alpha_max*deg)
    simulation.setBeamParameters(1.0*angstrom, 0.2*deg, 0.0*deg)
    simulation.setBeamIntensity(5e11)
    return simulation


def run_simulation():
    """
    Runs simulation and returns intensity map.
    """
    sample = get_sample()
    simulation = get_simulation()
    simulation.setSample(sample)
    simulation.runSimulation()
    return simulation.getIntensityData()


if __name__ == '__main__':
    result = run_simulation()
    ba.plot_intensity_data(result)