Calculating the Schottky barrier

From VASP Wiki
Revision as of 10:58, 23 March 2026 by Wolloch (talk | contribs)

The Schottky barrier is the energy barrier that forms at a metal–semiconductor junction. The barrier height is an important characteristic for charge transport across the junction and a critical parameter in the design of semiconductor contacts, transistors, and rectifying diodes. The p-type Schottky barrier height [math]\displaystyle{ \varphi_{\mathrm{p}} }[/math] describes the barrier seen by holes in the semiconductor valence band, and the n-type Schottky barrier height [math]\displaystyle{ \varphi_{\mathrm{n}} }[/math] describes the barrier seen by electrons in the conduction band. In VASP, the potential alignment methodcitation needed can be used to estimate the Schottky barrier.

Three quantities from bulk calculations, the Fermi energy [math]\displaystyle{ E_{\mathrm{F}} }[/math] of the metal, and the valence-band maximum [math]\displaystyle{ E_{\mathrm{VBM}} }[/math] and band gap [math]\displaystyle{ E_{\mathrm{g}} }[/math] of the semiconductor, are needed for the calculation. As well as the macroscopic electrostatic potential difference [math]\displaystyle{ \Delta\bar{V} }[/math] across the interface, which is extracted from a separate interface slab calculation. [1]

Step-by-step instructions

Important: For the calculation of the p- and n-type Schottky barriers, energies, and potentials from different VASP calculations are compared. It is crucial to set a single value for ENCUT for all calculations.

Step 1: Compute the Fermi energy of the metal.

After ensuring that the bulk geometry is optimized, a separate static calculation should be performed with a dense, Γ-centered k-point mesh and ISMEAR = -5 (thedrahedron method with Blöchl corrections), to get an accurate Fermi energy. It can be taken from the OUTCAR and is also easily accessible via py4vasp:
import py4vasp as pv
calc = pv.Calculation.from_path("./Step1/")
E_F = calc.dos.read()["fermi_energy"]
print(f"E_F = {E_F:.4f} eV")

Step 2: Compute the band gap and valence-band maximum of the semiconductor.

Band gaps are notoriously underestimated by standard DFT, so for accurate values, it is usually necessary to use DFT+U, hybrid functionals, or GW calculations. The level of theory to choose depends on the semiconductor of interest. Generally, it is advisable to select BANDGAP = KPOINT for more verbose output, and ideally to make a full bandstructure calculation, since the extrema of the bands will be located at high symmetry lines. The fundamental bandgap [math]\displaystyle{ E_{\mathrm{g}} }[/math] and the valence band maximum [math]\displaystyle{ E_{\mathrm{VBM}} }[/math], can be read from the OUTCAR file, or with py4vasp:
import py4vasp as pv

calc  = pv.Calculation.from_path("./Step2/")

E_g   = calc.bandgap.fundamental()
E_VBM = calc.bandgap.valence_band_maximum()

print(f"E_g   = {E_g:.4f} eV")
print(f"E_VBM = {E_VBM:.4f} eV")


Step 3: Build and relax an interface structure.

Building a lattice-matched metal–semiconductor interface slab is outside the scope of this page. Tools such as the CoherentInterfaceBuilder in pymatgen[2] and the interface construction utilities in ASE[3] can be used to construct the initial geometry. An interface structure for a Schottky barrier calculation should satisfy the following requirements:
  • Low lattice mismatch. The supercell should be chosen such that the in-plane periodicities of the two surfaces match to within ~2%, minimising artificial biaxial strain.
  • Sufficient slab thickness. Each material should contain enough layers that a bulk-like region exists in the centre.
  • No vacuum region. A cell with a vacuum will result in lateral strain, compressing the cell to reduce the surface energy. It is better to have two equivalent interfaces in the cell.
For the relaxation, it is usually prudent to fix the layers in the bulk like regions using selective dynamics in the POSCAR file. Only a couple of layers at the interfaces should be relaxed. Depending on the material pairing and the strain in the interface, it might be beneficial to allow the cell to relax only in the direction normal to the interface plane using the LATTICE_CONSTRAINTS INCAR tag.

Step 4: Compute the electrostatic potential of the interface and calculate the Schottky barrier heights.

Perform a static calculation with thetrahedron smearing (ISMEAR = -5, and PREC = Accurate as the final step. We need to print out the electrostatic part of the potential (disregarding the XC part) by setting WRT_POTENTIAL = hartree inoic. It is also possible to use LVHAR = .TRUE. to get the same data, but be aware that LVTOT would include the (semi-)local exchange-correlation potential [math]\displaystyle{ V_{\text{xc}}(\mathbf{r}) }[/math], which we want to exclude here.

The planar average [math]\displaystyle{ \bar{V}(z) }[/math] of the electrostatic potential is obtained by averaging the 3D potential over the in-plane (xy) grid at each [math]\displaystyle{ z }[/math] point (assuming the interface is located in the xy-plane). Then, a macroscopic average can be computed by applying a uniform box-car convolution of window length [math]\displaystyle{ L }[/math]. Where [math]\displaystyle{ L }[/math] is a bit less than the thickness of the thinner of the two materials. This removes short-range oscillations within each atomic layer while preserving the long-range potential step across the interface. Read the potential from vaspout.h5, compute both averages, and plot:

#!/usr/bin/env python3
"""
Potential alignment script for Schottky barrier height calculations.

The electrostatic potential (hartree + ionic) is read from vaspout.h5, averaged
in the xy-plane to give the planar average V_bar(z), and then smoothed with a
box-car macroscopic average of width L to remove short-range oscillations.

dV_bar = <V_bar>_M - <V_bar>_SC is the potential offset between the metal (M)
and semiconductor (SC) bulk-like regions of the slab.

L is chosen automatically by minimising the worst-case RMS gradient in the
two plateau regions (M and SC separately). This ensures both regions are flat,
even when the two materials have different lattice periodicities along z.

The plateau centres are given as fractional cell coordinates via --mc (metal)
and --scc (semiconductor). A fixed plateau half-width of 0.08 * c is used.
"""

import argparse
import matplotlib
matplotlib.use("Agg")  # non-interactive backend; remove if you want a GUI window
import matplotlib.pyplot as plt
import h5py
import numpy as np
from scipy.ndimage import uniform_filter1d

parser = argparse.ArgumentParser(description="Compute macroscopic potential alignment.")
parser.add_argument(
    "--mc", "--metal-center",
    type=float,
    default=0.125,
    metavar="F",
    dest="metal_center",
    help="Fractional z-coordinate of the metal plateau centre (default: 0.125)",
)
parser.add_argument(
    "--scc", "--sc-center",
    type=float,
    default=0.5,
    metavar="F",
    dest="sc_center",
    help="Fractional z-coordinate of the semiconductor plateau centre (default: 0.5)",
)
parser.add_argument(
    "--vaspout",
    default="./vaspout.h5",
    metavar="FILE",
    help="Path to vaspout.h5 (default: ./vaspout.h5)",
)
args = parser.parse_args()

# -- load potential and cell geometry ------------------------------------------
with h5py.File(args.vaspout, "r") as f:
    hartree = f["/results/potential/hartree"][0]  # (NGZ, NGY, NGX)
    ionic   = f["/results/potential/ionic"][0]    # (NGZ, NGY, NGX)
    lat     = f["/results/positions/lattice_vectors"][:]  # (3, 3) Ang

c   = np.linalg.norm(lat[2])  # cell length along interface normal (Ang)
NGZ = hartree.shape[0]        # number of grid points along z
dz  = c / NGZ                 # grid spacing (Ang)

# -- planar average: mean over the xy-plane at each z --------------------------
V_planar = (hartree + ionic).mean(axis=(1, 2))  # shape (NGZ,)   [eV]

z = np.arange(NGZ) * dz  # z-coordinates in Ang

# -- plateau masks (used both for auto-selection and for dV_bar extraction) ----
# z_frac maps each grid point to a fractional cell coordinate in [0, 1).
# For each material the mask selects a symmetric window of width 2*HALF_WIDTH
# centred on the user-supplied fractional coordinate (--mc / --scc).  The
# distance is computed with periodic boundary conditions (minimum image), so
# --mc=0.0 correctly selects both the bottom edge (z_frac < 0.08) and the
# top edge (z_frac > 0.92) of the cell, which belong to the same metal region.
# HALF_WIDTH = 0.08 means the window covers 16 % of the cell on each side,
# which is large enough to average over several bulk-like repeat units while
# staying well clear of the interface region for any reasonably converged
# slab (total cell >= ~40 Ang, >= 6 layers per material).
HALF_WIDTH = 0.08
z_frac  = z / c

def _periodic_mask(center):
    d = np.abs(z_frac - center)
    d = np.minimum(d, 1.0 - d)   # minimum-image distance in fractional coords
    return d < HALF_WIDTH

M_mask  = _periodic_mask(args.metal_center)
SC_mask = _periodic_mask(args.sc_center)

# -- auto-select window length -------------------------------------------------
# Sweep L from 2*dz to c/3. For each L the score is the worst-case RMS gradient
# across the two plateau regions (M and SC scored independently). Minimising the
# worst case forces both regions to be flat simultaneously.
L_grid = np.linspace(2 * dz, c / 3, 300)
def _plateau_roughness(L):
    Vm = uniform_filter1d(V_planar, size=max(1, int(round(L / dz))), mode="wrap")
    g  = np.gradient(Vm)
    return max(np.sqrt(np.mean(g[M_mask]**2)),
               np.sqrt(np.mean(g[SC_mask]**2)))
scores = [_plateau_roughness(L) for L in L_grid]
WINDOW_LENGTH = float(L_grid[np.argmin(scores)])
print(f"Auto-selected window length: {WINDOW_LENGTH:.2f} Ang")

# -- macroscopic average: box-car convolution with period L --------------------
# uniform_filter1d computes the running mean over `size` consecutive points.
# mode='wrap' enforces periodic boundary conditions, which is correct for
# a slab calculation where the potential is periodic along z.
n_window = int(round(WINDOW_LENGTH / dz))
V_macro  = uniform_filter1d(V_planar, size=n_window, mode="wrap")

# -- plot ----------------------------------------------------------------------
fig, ax = plt.subplots(figsize=(10, 4))
ax.plot(z, V_planar, lw=0.7, color="steelblue", alpha=0.5, label="Planar average")
ax.plot(z, V_macro,  lw=2.0, color="tomato",
        label=f"Macroscopic average (L = {WINDOW_LENGTH:.2f} Ang, {n_window} pts)")
ax.set_xlabel("z (Ang)")
ax.set_ylabel("Electrostatic potential (eV)")
ax.set_xlim(0, c)
ax.legend()
ax.set_title("M/SC interface - electrostatic potential along interface normal")
plt.tight_layout()
plt.savefig("potential_alignment.png", dpi=150)
print("Figure saved: potential_alignment.png")

# -- extract dV_bar from bulk-like plateaus ------------------------------------
# M_mask and SC_mask are defined above; adjust --mc / --scc if the plateaus shift.

V_M     = V_macro[M_mask].mean()
V_SC    = V_macro[SC_mask].mean()
delta_V = V_M - V_SC

print(f"\n<V_bar>_M   = {V_M:+.4f} eV")
print(f"<V_bar>_SC  = {V_SC:+.4f} eV")
print(f"dV_bar      = <V_bar>_M - <V_bar>_SC = {delta_V:+.4f} eV")

Step 5: Compute the Schottky barrier heights.

With [math]\displaystyle{ E_{\mathrm{F}} }[/math], [math]\displaystyle{ E_{\mathrm{VBM}} }[/math], [math]\displaystyle{ E_{\mathrm{g}} }[/math], and [math]\displaystyle{ \Delta\bar{V} }[/math] in hand, the potential alignment difference [math]\displaystyle{ \Delta\bar{V} = \bar{V}_{\mathrm{m}} - \bar{V}_{\mathrm{sc}} }[/math] (positive if the metal side has the higher average potential) connects the bulk reference frames. The p-type and n-type Schottky barrier heights are then

[math]\displaystyle{ \varphi_{\mathrm{p}} = \Delta\bar{V} + E_{\mathrm{F}} - E_{\mathrm{VBM}}, }[/math]
[math]\displaystyle{ \varphi_{\mathrm{n}} = E_{\mathrm{g}} - \varphi_{\mathrm{p}}. }[/math]

Recommendations and advice

(To be added.)

Related tags and articles

(To be added.)

References

  1. In VASP, the [math]\displaystyle{ \mathbf{G}=0 }[/math] Fourier component of the Hartree potential is set to zero by convention, so all single-particle eigenvalues are already referenced to the average electrostatic potential of their respective unit cell. As a consequence, [math]\displaystyle{ E_{\mathrm{F}} }[/math] and [math]\displaystyle{ E_{\mathrm{VBM}} }[/math] can be read directly from the bulk OUTCAR files without any further potential referencing, and only the interface calculation requires the electrostatic potential written to LOCPOT.
  2. https://pymatgen.org/ (2022).
  3. https://wiki.fysik.dtu.dk/ase/ (2025).