File: RoughSurface.py

package info (click to toggle)
bornagain 23.0-5
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 103,948 kB
  • sloc: cpp: 423,131; python: 40,997; javascript: 11,167; awk: 630; sh: 318; ruby: 173; xml: 130; makefile: 51; ansic: 24
file content (108 lines) | stat: -rwxr-xr-x 3,098 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
#!/usr/bin/env python3
"""
Generates random rough surface of the sample with the given
heights statistics and lateral roughness spectrum.
"""

import numpy as np
import bornagain as ba
from bornagain import ba_plot as bp, nm, nm2, nm3
from matplotlib import gridspec
from bornagain.numpyutil import Arrayf64Converter as dac

def plot(h):
    rms = np.std(surface)

    fig = bp.plt.figure(figsize=(10, 6))
    gs = gridspec.GridSpec(3, 2)
    gs.update(left=0.1, right=0.48, wspace=0.15, hspace=0.5)

    # map
    fig.add_subplot(gs[:-1, :])
    bp.plt.imshow(h, extent = [0, Lx, 0, Ly], cmap='inferno', aspect='auto')
    bp.plt.colorbar()

    # psd
    fig.add_subplot(gs[-1, 0])
    factor = X_points*Y_points/np.sqrt(Lx*Ly)
    N = np.floor(X_points/2 + 1).astype(int)
    psd = abs(np.fft.rfft2(h)/factor)**2
    psd_points_x = np.array([psd[0][i] for i in range(1, N)])
    f_points_x = np.array([i/Lx for i in range(1, N)])
    bp.plt.plot(f_points_x, psd_points_x)
    bp.plt.yscale('log')
    bp.plt.xscale('log')
    bp.plt.title('PSD')
    bp.plt.xlabel('$\\nu, nm^{-1}$')

    # hist
    fig.add_subplot(gs[-1, 1])
    bp.plt.hist(h.flatten(), range=[-4*rms, 4*rms],
             bins=300, color='black')
    bp.plt.title('Height histogram')
    bp.plt.xlabel('h, nm')
    bp.plt.tick_params(left = False, labelleft = False)

    bp.plt.show()

def get_sample():
    """
    A sample with one layer on a substrate, with partially
    correlated roughnesses.
    """
    # defining materials
    vacuum = ba.RefractiveMaterial("ambience", 0, 0)
    material_layer = ba.RefractiveMaterial("layer", 5e-6, 0)
    material_substrate = ba.RefractiveMaterial("substrate", 15e-6, 0)

    # defining roughness
    height_distribution = ba.ErfTransient()
    max_spat_freq = 0.5*nm

    # for layer
    omega = 1*nm3
    a1 = 1
    a2 = 10*nm
    a3 = 100*nm2
    a4 = 1000*nm3
    autocorr_layer = ba.LinearGrowthModel(omega, a1, a2, a3, a4, max_spat_freq)
    roughness_layer = ba.Roughness(autocorr_layer, height_distribution)

    # for substrate
    sigma = 0.9*nm
    alpha = 0.5
    xi = 35*nm
    autocorr_sub = ba.SelfAffineFractalModel(sigma, alpha, xi, max_spat_freq)
    roughness_sub = ba.Roughness(autocorr_sub, height_distribution)

    # defining layers
    l_ambience = ba.Layer(vacuum)
    l_layer = ba.Layer(material_layer, 25*nm, roughness_layer)
    l_substrate = ba.Layer(material_substrate, roughness_sub)

    # adding layers
    my_sample = ba.Sample()
    my_sample.addLayer(l_ambience)
    my_sample.addLayer(l_layer)
    my_sample.addLayer(l_substrate)

    return my_sample

if __name__ == '__main__':

    # sample size
    Lx = 1000*nm
    Ly = 1000*nm
    X_points = 10
    Y_points = 10

    sample = get_sample()
    layer_index = 1 # top interface
    seed = -1 # seed < 0 ==> random every time

    # generate roughness map
    roughness_map = ba.RoughnessMap(X_points, Y_points, Lx, Ly,
                                    sample, layer_index, seed)
    surface = dac.npArray(roughness_map.generate())
    print("rms = {:.3}".format(np.std(surface)),"nm")