Capped Nanotube Construction

10.3 Capped Nanotube Construction

Single Capped Nanotube Construction

nanocap_single_capped_nanotube.py is an example script to construct and save a capped nanotube. The code is shown below.

'''
-=-=-=-=-=-=-=-=-=-=-=-=-=NanoCap=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
Copyright Marc Robinson  2014
-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-

A script to construct a capped
nanotube.

Input: 
    n,m = Chirality (n,m)
    l = length
    cap_estimate = estimate cap from
                   tube density
    dual_lattice_force_field = force field 
                               for dual lattice
    carbon_force_field = force field 
                        for carbon lattice
    dual_lattice_mintol= energy tolerance for
                         dual lattice optimisation
    dual_lattice_minsteps= steps for dual lattice 
                            optimisation
    carbon_lattice_mintol=as above for carbon lattice
    carbon_lattice_minsteps=as above for carbon lattice
    optimiser=optimsation algorithm
    seed = seed for initial cap generation
                   
Output:
    xyz files containing dual lattice
    and carbon lattice 
    
-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
'''
import sys,os,random,numpy
from nanocap.core import minimisation
from nanocap.structures import cappednanotube
from nanocap.core import output

n,m = 7,3
l = 10.0 
cap_estimate = True
     
dual_force_field = "Thomson"
carbon_force_field = "EDIP"
dual_lattice_mintol=1e-10
dual_lattice_minsteps=100
carbon_lattice_mintol=1e-10
carbon_lattice_minsteps=100
optimiser="LBFGS"
seed = 12345

my_nanotube = cappednanotube.CappedNanotube()

my_nanotube.setup_nanotube(n,m,l=l)

if(cap_estimate):
    NCapDual = my_nanotube.get_cap_dual_lattice_estimate(n,m)

my_nanotube.construct_dual_lattice(N_cap_dual=NCapDual,seed=seed)

my_nanotube.set_Z_cutoff(N_cap_dual=NCapDual)

cap = my_nanotube.cap
outfilename = "n_{}_m_{}_l_{}_cap_{}_dual_lattice_init"
outfilename = outfilename.format(n,m,l,cap.dual_lattice.npoints)
output.write_xyz(outfilename,my_nanotube.dual_lattice)


Dminimiser = minimisation.DualLatticeMinimiser(FFID=dual_force_field,
                                               structure = my_nanotube)
Dminimiser.minimise(my_nanotube.dual_lattice,
                    min_type=optimiser,
                    ftol=dual_lattice_mintol,
                    min_steps=dual_lattice_minsteps)

my_nanotube.update_caps()
outfilename = "n_{}_m_{}_l_{}_cap_{}_dual_lattice"
outfilename = outfilename.format(n,m,l,cap.dual_lattice.npoints)
output.write_xyz(outfilename,my_nanotube.dual_lattice)

my_nanotube.construct_carbon_lattice()

Cminimiser = minimisation.CarbonLatticeMinimiser(FFID=carbon_force_field,
                                                 structure = my_nanotube)

Cminimiser.minimise_scale(my_nanotube.carbon_lattice)
Cminimiser.minimise(my_nanotube.carbon_lattice,
                    min_type=optimiser,
                    ftol=carbon_lattice_mintol,
                    min_steps=carbon_lattice_minsteps)

outfilename = "n_{}_m_{}_l_{}_cap_{}_carbon_atoms"
outfilename = outfilename.format(n,m,l,cap.dual_lattice.npoints)
output.write_xyz(outfilename,my_nanotube.carbon_lattice)
outfilename = "n_{}_m_{}_l_{}_cap_{}_carbon_atoms_constrained"
outfilename = outfilename.format(n,m,l,cap.dual_lattice.npoints)
output.write_xyz(outfilename,my_nanotube.carbon_lattice,constrained=True)

Constructing Multiple Capped Nanotubes

nanocap_bulk_capped_nanotubes.py is an example script to perform a structure search to construct and save multiple capped nanotube. The code is shown below.

'''
-=-=-=-=-=-=-=-=-=-=-=-=-=NanoCap=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
Copyright Marc Robinson  2014
-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-

A script to construct a series of capped
nanotubes of the same chirality

Input: 
    n,m = Chirality (n,m)
    l = length
    cap_estimate = estimate cap from
                   tube density
    dual_lattice_force_field = force field 
                               for dual lattice
    carbon_force_field = force field 
                        for carbon lattice
    dual_lattice_mintol= energy tolerance for
                         dual lattice optimisation
    dual_lattice_minsteps= steps for dual lattice 
                            optimisation
    carbon_lattice_mintol=as above for carbon lattice
    carbon_lattice_minsteps=as above for carbon lattice
    optimiser=optimsation algorithm
    seed = seed for initial cap generation
    N_nanotubes = required number of structures 
    N_max_structures = maximum number of possible 
                        structures to search through
    basin_climb = True/False - climb out of 
                  minima  
    calc_rings = True/False - calculate rings for 
                  each structure
                   
Output:
    -A structure log in myStructures.out
    
    -xyz files containing the carbon lattices 
    
-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-=-
'''
    
import os
from nanocap.core.minimisation import DualLatticeMinimiser
from nanocap.core.minimisation import CarbonLatticeMinimiser  
from nanocap.core.minimasearch import MinimaSearch
from nanocap.structures.cappednanotube import CappedNanotube
from nanocap.core.output import write_points

n,m = 10,10
l = 20.0 
cap_estimate = True

N_nanotubes = 5
N_max_structures = 20
basin_climb = True
calc_rings = True

dual_lattice_minimiser = "Thomson"
carbon_lattice_minimiser = "EDIP"
dual_lattice_mintol=1e-10
dual_lattice_minsteps=100
carbon_lattice_mintol=1e-10
carbon_lattice_minsteps=100
optimiser="LBFGS"
seed = 12345

my_nanotube = CappedNanotube()
my_nanotube.setup_nanotube(n,m,l=l)

if(cap_estimate):
    N_cap_dual = my_nanotube.get_cap_dual_lattice_estimate(n,m)

my_nanotube.construct_dual_lattice(N_cap_dual=N_cap_dual,seed=seed)
my_nanotube.set_Z_cutoff(N_cap_dual=N_cap_dual)

Dminimiser = DualLatticeMinimiser(FFID=dual_lattice_minimiser,
                                  structure = my_nanotube,
                                  min_type= optimiser,
                                  ftol = dual_lattice_mintol,
                                  min_steps = dual_lattice_minsteps)

Cminimiser = CarbonLatticeMinimiser(FFID=carbon_lattice_minimiser,
                                  structure = my_nanotube,
                                  min_type= optimiser,
                                  ftol = carbon_lattice_mintol,
                                  min_steps = carbon_lattice_minsteps)

Searcher = MinimaSearch(Dminimiser,
                        carbon_lattice_minimiser= Cminimiser,
                        basin_climb=basin_climb,
                        calc_rings=calc_rings)

Searcher.start_search(my_nanotube.dual_lattice,
                      N_nanotubes,
                      N_max_structures)

Searcher.structure_log.write_log(os.getcwd(),"myStructures.out")

for i,structure in enumerate(Searcher.structure_log.structures):
    carbon_lattice = structure.carbon_lattice
    filename = "C{}_carbon_atoms_{}".format(carbon_lattice.npoints,i)
    write_points(filename,carbon_lattice,format="xyz")