← Wiki

Building Nanoflakes

Building Nanoflakes

Notes on building 2-D materials with the Atomic Simulation Environment (ASE) and Quantum ESPRESSO (QE).

Background

  • ASE is a Python library which provides a suite of tools for producing atomic simulation files, connecting to various simulation programs, and post-processing the results from these simulations.
  • QE is a free and open source software suite for performing electronic structure calculations and materials modeling at the nanoscale. It is based on density functional theory (DFT) and plane wave basis sets.
  • 2-D materials produce novel quantum phenomena due to topological effects and quantum confinement. These materials are of interest for applications in electronics, photonics, and energy storage.

We will use ASE to build a nanoflake of Silicon Carbide, SiC, a 2-D material with a hexagonal lattice structure, and then set up a simple QE simulation using this nanoflake as the input structure.

Atomistically Flat SiC

Hexagonal lattices have a two-point basis, following the conventions of graphene, the basis vectors for SiC are:

  • α1 = a(3, √3)
  • α2 = a(3, -√3)

where a is the lattice constant distance. This spacing gives the vectors to homogeneous nearest neighbors in the lattice, in this case, Silicon.

Then, from any given Silicon point, the heteroatom nearest neighbors, are given by:

  • δ1 = a(1, √3)/2
  • δ2 = a(1, -√3)/2
  • δ3 = a(-1, 0)

We know from literature that the lattice constant for SiC is approximately 3.08 Å, so we can use this value to build our nanoflake.

We also incorporate hydrogen atoms at the edges of the nanoflake to passivate dangling bonds, which can affect the electronic properties of the material. We can use the vector normal to the plane of the nanoflake to determine the positions of the hydrogen atoms, which will be placed at a distance of 1.48 Å from the edge Si, and 1.09 Å from the edge C atoms, following typical bond lengths for Si-H and C-H bonds.

Finally, we provide 12 Å of vacuum after centering the constructed nanoflake in a simulation cell to prevent interactions between periodic images of the nanoflake in the z-direction.

Building the Nanoflake with ASE

Putting all of the above together, we can use the following Python code to build our SiC nanoflake with ASE:

#!/usr/bin/env python 
from ase import Atoms 
from ase.io import write
from ase.neighborlist import NeighborList 
import ase.data
import numpy as np 

def build_sic_flake(lattice_a=3.08, n_size=4):
    """
    Build a hexagonal SiC nanoflake based on treating it as a modified graphene lattice.

    Parameters:
        lattice_a (float): The Si to Si distance in angstroms. Default is 3.08 Å. 
        n_repeat (int): The 'radius' of the hexagon unit cells. Default is 4. 
    Returns: 
        Atoms: An ASE Atoms object representing the SiC nanoflake.
    """

    bond_length = lattice_a / np.sqrt(3)  # Adjust bond length for SiC 
    bond_length_si_h = 1.48 # Si-H bond length in angstroms 
    bond_length_c_h = 1.09  # C-H bond length in angstroms
    cutoff = 2.1 # Cutoff for neighbor list in angstroms

    a1 = np.array([lattice_a, 0, 0])
    a2 = np.array([lattice_a * 0.5, lattice_a * np.sqrt(3) * 0.5, 0])

    basis_si = np.array([0.,0.,0.])
    basis_c = (a1 + a2) / 3.0 

    atoms_list = []
    grid_range = range(-n_size, n_size + 1) 
    for i in grid_range:
        for j in grid_range:
            disp = i * a1 + j * a2 
            atoms_list.append(('Si', basis_si + disp))
            atoms_list.append(('C', basis_c + disp)) 
    
    symbols = [x[0] for x in atoms_list]
    positions = [x[1] for x in atoms_list]
    flake = Atoms(symbols=symbols, positions=positions, pbc=False) 

    flake.center()
    center_pos = flake.positions.mean(axis=0) 
    flake.translate(-center_pos) 

    nl = NeighborList([cutoff * 0.5] * len(flake), self_interaction=False, bothways=True)
    nl.update(flake) 

    hydrogens = [] 
    for atom in flake:
        indices, offsets = nl.get_neighbors(atom.index)
        coordination = len(indices) 
        if coordination < 3: 
            # the 'surface normal' for the atom; point away from neighbors 
            neighbor_vectors = [flake.positions[n] - atom.position for n in indices]
            sum_vector = np.sum(neighbor_vectors, axis=0)
            if np.linalg.norm(sum_vector) > 1e-6:
                direction = -sum_vector / np.linalg.norm(sum_vector)
            else:
                continue # should not happen 
            d_h = bond_length_si_h if atom.symbol == 'Si' else bond_length_c_h
            h_pos = atom.position + direction * d_h
            hydrogens.append(('H', h_pos))
    if hydrogens:
        h_symbols = [x[0] for x in hydrogens]
        h_positions = [x[1] for x in hydrogens]
        flake += Atoms(symbols=h_symbols, positions=h_positions)
    
    return flake

if __name__ == "__main__":
    flake = build_sic_flake(lattice_a=3.08, n_size=6) 
    flake.center(vacuum=12.0)

    qe_input_data = { 
        'control': {
            'calculation': 'scf',
            'restart_mode': 'from_scratch',
            'pseudo_dir': './',
            'outdir': './outputs',
            'prefix': 'sic_benchmark',
            'disk_io': 'low', 
            'tprnfor': True, # calculate forces 
        },
        'system': {
            'ecutwfc': 50, # standard cutoff 
            'ibrav': 0,    # free cell 
            'nat': len(flake),
            'ntyp': len(set(flake.get_chemical_symbols())), # 3 types: Si, C, H 
            'nbnd': int(len(flake) * 2.5) # generous number of bands
        },
        'electrons': {
            'mixing_beta': 0.3, 
            'diagonalization': 'cg',
            'conv_thr': 1.0e-6, 
        }
    }

    write('sic_flake.pw.in', 
        flake, 
        format='espresso-in',
        input_data=qe_input_data, 
        pseudopotentials={'Si': 'Si.pz-vbc.UPF', 'C': 'C.pz-vbc.UPF', 'H': 'H.pz-vbc.UPF'},
        kpts=None) # Gamma point only 

    # also write an XYZ file for visualization 
    write('sic_flake.xyz', flake) 

    print(f"SiC nanoflake structure and QE input file generated with {len(flake)} atoms.")
    # print the formula name for the produced flake 
    print(f"Formula: {flake.get_chemical_formula()}")
    print("Generated files:")
    print(f"sic_flake.pw.in, sic_flake.xyz")

The generated files can be used to run a (basic) DFT calculation on the SiC nanoflake using QE.