root / ase / io / vasp.py @ 14
Historique | Voir | Annoter | Télécharger (11,62 ko)
| 1 |
"""
|
|---|---|
| 2 |
This module contains functionality for reading and writing an ASE
|
| 3 |
Atoms object in VASP POSCAR format.
|
| 4 |
|
| 5 |
"""
|
| 6 |
|
| 7 |
import os |
| 8 |
|
| 9 |
def get_atomtypes(fname): |
| 10 |
"""Given a file name, get the atomic symbols.
|
| 11 |
|
| 12 |
The function can get this information from OUTCAR and POTCAR
|
| 13 |
format files. The files can also be compressed with gzip or
|
| 14 |
bzip2.
|
| 15 |
|
| 16 |
"""
|
| 17 |
atomtypes=[] |
| 18 |
if fname.find('.gz') != -1: |
| 19 |
import gzip |
| 20 |
f = gzip.open(fname) |
| 21 |
elif fname.find('.bz2') != -1: |
| 22 |
import bz2 |
| 23 |
f = bz2.BZ2File(fname) |
| 24 |
else:
|
| 25 |
f = open(fname)
|
| 26 |
for line in f: |
| 27 |
if line.find('TITEL') != -1: |
| 28 |
atomtypes.append(line.split()[3].split('_')[0].split('.')[0]) |
| 29 |
return atomtypes
|
| 30 |
|
| 31 |
def atomtypes_outpot(posfname, numsyms): |
| 32 |
"""Try to retreive chemical symbols from OUTCAR or POTCAR
|
| 33 |
|
| 34 |
If getting atomtypes from the first line in POSCAR/CONTCAR fails, it might
|
| 35 |
be possible to find the data in OUTCAR or POTCAR, if these files exist.
|
| 36 |
|
| 37 |
posfname -- The filename of the POSCAR/CONTCAR file we're trying to read
|
| 38 |
|
| 39 |
numsyms -- The number of symbols we must find
|
| 40 |
|
| 41 |
"""
|
| 42 |
import os.path as op |
| 43 |
import glob |
| 44 |
|
| 45 |
# First check files with exactly same name except POTCAR/OUTCAR instead
|
| 46 |
# of POSCAR/CONTCAR.
|
| 47 |
fnames = [posfname.replace('POSCAR', 'POTCAR').replace('CONTCAR', |
| 48 |
'POTCAR')]
|
| 49 |
fnames.append(posfname.replace('POSCAR', 'OUTCAR').replace('CONTCAR', |
| 50 |
'OUTCAR'))
|
| 51 |
# Try the same but with compressed files
|
| 52 |
fsc = [] |
| 53 |
for fn in fnames: |
| 54 |
fsc.append(fn + '.gz')
|
| 55 |
fsc.append(fn + '.bz2')
|
| 56 |
for f in fsc: |
| 57 |
fnames.append(f) |
| 58 |
# Finally try anything with POTCAR or OUTCAR in the name
|
| 59 |
vaspdir = op.dirname(posfname) |
| 60 |
fs = glob.glob(vaspdir + '*POTCAR*')
|
| 61 |
for f in fs: |
| 62 |
fnames.append(f) |
| 63 |
fs = glob.glob(vaspdir + '*OUTCAR*')
|
| 64 |
for f in fs: |
| 65 |
fnames.append(f) |
| 66 |
|
| 67 |
tried = [] |
| 68 |
files_in_dir = os.listdir('.')
|
| 69 |
for fn in fnames: |
| 70 |
tried.append(fn) |
| 71 |
if fn in files_in_dir: |
| 72 |
at = get_atomtypes(fn) |
| 73 |
if len(at) == numsyms: |
| 74 |
return at
|
| 75 |
|
| 76 |
raise IOError('Could not determine chemical symbols. Tried files ' |
| 77 |
+ str(tried))
|
| 78 |
|
| 79 |
|
| 80 |
def read_vasp(filename='CONTCAR'): |
| 81 |
"""Import POSCAR/CONTCAR type file.
|
| 82 |
|
| 83 |
Reads unitcell, atom positions and constraints from the POSCAR/CONTCAR
|
| 84 |
file and tries to read atom types from POSCAR/CONTCAR header, if this fails
|
| 85 |
the atom types are read from OUTCAR or POTCAR file.
|
| 86 |
"""
|
| 87 |
|
| 88 |
from ase import Atoms, Atom |
| 89 |
from ase.constraints import FixAtoms, FixScaled |
| 90 |
from ase.data import chemical_symbols |
| 91 |
import numpy as np |
| 92 |
|
| 93 |
if isinstance(filename, str): |
| 94 |
f = open(filename)
|
| 95 |
else: # Assume it's a file-like object |
| 96 |
f = filename |
| 97 |
|
| 98 |
# First line should contain the atom symbols , eg. "Ag Ge" in
|
| 99 |
# the same order
|
| 100 |
# as later in the file (and POTCAR for the full vasp run)
|
| 101 |
atomtypes = f.readline().split() |
| 102 |
|
| 103 |
lattice_constant = float(f.readline())
|
| 104 |
|
| 105 |
# Now the lattice vectors
|
| 106 |
a = [] |
| 107 |
for ii in range(3): |
| 108 |
s = f.readline().split() |
| 109 |
floatvect = float(s[0]), float(s[1]), float(s[2]) |
| 110 |
a.append(floatvect) |
| 111 |
|
| 112 |
basis_vectors = np.array(a) * lattice_constant |
| 113 |
|
| 114 |
# Number of atoms. Again this must be in the same order as
|
| 115 |
# in the first line
|
| 116 |
# or in the POTCAR or OUTCAR file
|
| 117 |
atom_symbols = [] |
| 118 |
numofatoms = f.readline().split() |
| 119 |
#vasp5.1 has an additional line which gives the atom types
|
| 120 |
#the following try statement skips this line
|
| 121 |
try:
|
| 122 |
int(numofatoms[0]) |
| 123 |
except ValueError: |
| 124 |
numofatoms = f.readline().split() |
| 125 |
|
| 126 |
numsyms = len(numofatoms)
|
| 127 |
if len(atomtypes) < numsyms: |
| 128 |
# First line in POSCAR/CONTCAR didn't contain enough symbols.
|
| 129 |
atomtypes = atomtypes_outpot(f.name, numsyms) |
| 130 |
else:
|
| 131 |
try:
|
| 132 |
for atype in atomtypes[:numsyms]: |
| 133 |
if not atype in chemical_symbols: |
| 134 |
raise KeyError |
| 135 |
except KeyError: |
| 136 |
atomtypes = atomtypes_outpot(f.name, numsyms) |
| 137 |
|
| 138 |
for i, num in enumerate(numofatoms): |
| 139 |
numofatoms[i] = int(num)
|
| 140 |
[atom_symbols.append(atomtypes[i]) for na in xrange(numofatoms[i])] |
| 141 |
|
| 142 |
# Check if Selective dynamics is switched on
|
| 143 |
sdyn = f.readline() |
| 144 |
selective_dynamics = sdyn[0].lower() == "s" |
| 145 |
|
| 146 |
# Check if atom coordinates are cartesian or direct
|
| 147 |
if selective_dynamics:
|
| 148 |
ac_type = f.readline() |
| 149 |
else:
|
| 150 |
ac_type = sdyn |
| 151 |
cartesian = ac_type[0].lower() == "c" or ac_type[0].lower() == "k" |
| 152 |
tot_natoms = sum(numofatoms)
|
| 153 |
atoms_pos = np.empty((tot_natoms, 3))
|
| 154 |
if selective_dynamics:
|
| 155 |
selective_flags = np.empty((tot_natoms, 3), dtype=bool) |
| 156 |
for atom in xrange(tot_natoms): |
| 157 |
ac = f.readline().split() |
| 158 |
atoms_pos[atom] = (float(ac[0]), float(ac[1]), float(ac[2])) |
| 159 |
if selective_dynamics:
|
| 160 |
curflag = [] |
| 161 |
for flag in ac[3:6]: |
| 162 |
curflag.append(flag == 'F')
|
| 163 |
selective_flags[atom] = curflag |
| 164 |
# Done with all reading
|
| 165 |
if type(filename) == str: |
| 166 |
f.close() |
| 167 |
if cartesian:
|
| 168 |
atoms_pos *= lattice_constant |
| 169 |
atoms = Atoms(symbols = atom_symbols, cell = basis_vectors, pbc = True)
|
| 170 |
if cartesian:
|
| 171 |
atoms.set_positions(atoms_pos) |
| 172 |
else:
|
| 173 |
atoms.set_scaled_positions(atoms_pos) |
| 174 |
if selective_dynamics:
|
| 175 |
constraints = [] |
| 176 |
indices = [] |
| 177 |
for ind, sflags in enumerate(selective_flags): |
| 178 |
if sflags.any() and not sflags.all(): |
| 179 |
constraints.append(FixScaled(atoms.get_cell(), ind, sflags)) |
| 180 |
elif sflags.all():
|
| 181 |
indices.append(ind) |
| 182 |
if indices:
|
| 183 |
constraints.append(FixAtoms(indices)) |
| 184 |
if constraints:
|
| 185 |
atoms.set_constraint(constraints) |
| 186 |
return atoms
|
| 187 |
|
| 188 |
def read_vasp_out(filename='OUTCAR',index = -1): |
| 189 |
"""Import OUTCAR type file.
|
| 190 |
|
| 191 |
Reads unitcell, atom positions, energies, and forces from the OUTCAR file.
|
| 192 |
|
| 193 |
- does not (yet) read any constraints
|
| 194 |
"""
|
| 195 |
import os |
| 196 |
import numpy as np |
| 197 |
from ase.calculators.singlepoint import SinglePointCalculator |
| 198 |
from ase import Atoms, Atom |
| 199 |
|
| 200 |
if isinstance(filename, str): |
| 201 |
f = open(filename)
|
| 202 |
else: # Assume it's a file-like object |
| 203 |
f = filename |
| 204 |
data = f.readlines() |
| 205 |
natoms = 0
|
| 206 |
images = [] |
| 207 |
atoms = Atoms(pbc = True)
|
| 208 |
energy = 0
|
| 209 |
species = [] |
| 210 |
species_num = [] |
| 211 |
symbols = [] |
| 212 |
for n,line in enumerate(data): |
| 213 |
if 'POSCAR:' in line: |
| 214 |
temp = line.split() |
| 215 |
species = temp[1:]
|
| 216 |
if 'ions per type' in line: |
| 217 |
temp = line.split() |
| 218 |
for ispecies in range(len(species)): |
| 219 |
species_num += [int(temp[ispecies+4])] |
| 220 |
natoms += species_num[-1]
|
| 221 |
for iatom in range(species_num[-1]): symbols += [species[ispecies]] |
| 222 |
if 'direct lattice vectors' in line: |
| 223 |
cell = [] |
| 224 |
for i in range(3): |
| 225 |
temp = data[n+1+i].split()
|
| 226 |
cell += [[float(temp[0]), float(temp[1]), float(temp[2])]] |
| 227 |
atoms.set_cell(cell) |
| 228 |
if 'FREE ENERGIE OF THE ION-ELECTRON SYSTEM' in line: |
| 229 |
energy = float(data[n+2].split()[4]) |
| 230 |
if 'POSITION ' in line: |
| 231 |
forces = [] |
| 232 |
for iatom in range(natoms): |
| 233 |
temp = data[n+2+iatom].split()
|
| 234 |
atoms += Atom(symbols[iatom],[float(temp[0]),float(temp[1]),float(temp[2])]) |
| 235 |
forces += [[float(temp[3]),float(temp[4]),float(temp[5])]] |
| 236 |
atoms.set_calculator(SinglePointCalculator(energy,forces,None,None,atoms)) |
| 237 |
images += [atoms] |
| 238 |
atoms = Atoms(pbc = True)
|
| 239 |
|
| 240 |
# return requested images, code borrowed from ase/io/trajectory.py
|
| 241 |
if isinstance(index, int): |
| 242 |
return images[index]
|
| 243 |
else:
|
| 244 |
step = index.step or 1 |
| 245 |
if step > 0: |
| 246 |
start = index.start or 0 |
| 247 |
if start < 0: |
| 248 |
start += len(images)
|
| 249 |
stop = index.stop or len(images) |
| 250 |
if stop < 0: |
| 251 |
stop += len(images)
|
| 252 |
else:
|
| 253 |
if index.start is None: |
| 254 |
start = len(images) - 1 |
| 255 |
else:
|
| 256 |
start = index.start |
| 257 |
if start < 0: |
| 258 |
start += len(images)
|
| 259 |
if index.stop is None: |
| 260 |
stop = -1
|
| 261 |
else:
|
| 262 |
stop = index.stop |
| 263 |
if stop < 0: |
| 264 |
stop += len(images)
|
| 265 |
return [images[i] for i in range(start, stop, step)] |
| 266 |
|
| 267 |
def write_vasp(filename, atoms, label='', direct=False, sort=None, symbol_count = None ): |
| 268 |
"""Method to write VASP position (POSCAR/CONTCAR) files.
|
| 269 |
|
| 270 |
Writes label, scalefactor, unitcell, # of various kinds of atoms,
|
| 271 |
positions in cartesian or scaled coordinates (Direct), and constraints
|
| 272 |
to file. Cartesian coordiantes is default and default label is the
|
| 273 |
atomic species, e.g. 'C N H Cu'.
|
| 274 |
"""
|
| 275 |
|
| 276 |
import numpy as np |
| 277 |
from ase.constraints import FixAtoms, FixScaled |
| 278 |
|
| 279 |
if isinstance(filename, str): |
| 280 |
f = open(filename, 'w') |
| 281 |
else: # Assume it's a 'file-like object' |
| 282 |
f = filename |
| 283 |
|
| 284 |
if isinstance(atoms, (list, tuple)): |
| 285 |
if len(atoms) > 1: |
| 286 |
raise RuntimeError("Don't know how to save more than "+ |
| 287 |
"one image to VASP input")
|
| 288 |
else:
|
| 289 |
atoms = atoms[0]
|
| 290 |
|
| 291 |
# Write atom positions in scaled or cartesian coordinates
|
| 292 |
if direct:
|
| 293 |
coord = atoms.get_scaled_positions() |
| 294 |
else:
|
| 295 |
coord = atoms.get_positions() |
| 296 |
|
| 297 |
if atoms.constraints:
|
| 298 |
sflags = np.zeros((len(atoms), 3), dtype=bool) |
| 299 |
for constr in atoms.constraints: |
| 300 |
if isinstance(constr, FixScaled): |
| 301 |
sflags[constr.a] = constr.mask |
| 302 |
elif isinstance(constr, FixAtoms): |
| 303 |
sflags[constr.index] = [True, True, True] |
| 304 |
|
| 305 |
if sort:
|
| 306 |
ind = np.argsort(atoms.get_chemical_symbols()) |
| 307 |
symbols = np.array(atoms.get_chemical_symbols())[ind] |
| 308 |
coord = coord[ind] |
| 309 |
if atoms.constraints:
|
| 310 |
sflags = sflags[ind] |
| 311 |
else:
|
| 312 |
symbols = atoms.get_chemical_symbols() |
| 313 |
|
| 314 |
# Create a list sc of (symbol, count) pairs
|
| 315 |
if symbol_count:
|
| 316 |
sc = symbol_count |
| 317 |
else:
|
| 318 |
sc = [] |
| 319 |
psym = symbols[0]
|
| 320 |
count = 0
|
| 321 |
for sym in symbols: |
| 322 |
if sym != psym:
|
| 323 |
sc.append((psym, count)) |
| 324 |
psym = sym |
| 325 |
count = 1
|
| 326 |
else:
|
| 327 |
count += 1
|
| 328 |
sc.append((psym, count)) |
| 329 |
|
| 330 |
# Create the label
|
| 331 |
if label == '': |
| 332 |
for sym, c in sc: |
| 333 |
label += '%2s ' % sym
|
| 334 |
f.write(label + '\n')
|
| 335 |
|
| 336 |
# Write unitcell in real coordinates and adapt to VASP convention
|
| 337 |
# for unit cell
|
| 338 |
# ase Atoms doesn't store the lattice constant separately, so always
|
| 339 |
# write 1.0.
|
| 340 |
f.write('%19.16f\n' % 1.0) |
| 341 |
for vec in atoms.get_cell(): |
| 342 |
f.write(' ')
|
| 343 |
for el in vec: |
| 344 |
f.write(' %21.16f' % el)
|
| 345 |
f.write('\n')
|
| 346 |
|
| 347 |
# Numbers of each atom
|
| 348 |
for sym, count in sc: |
| 349 |
f.write(' %3i' % count)
|
| 350 |
f.write('\n')
|
| 351 |
|
| 352 |
if atoms.constraints:
|
| 353 |
f.write('Selective dynamics\n')
|
| 354 |
|
| 355 |
if direct:
|
| 356 |
f.write('Direct\n')
|
| 357 |
else:
|
| 358 |
f.write('Cartesian\n')
|
| 359 |
|
| 360 |
for iatom, atom in enumerate(coord): |
| 361 |
for dcoord in atom: |
| 362 |
f.write(' %19.16f' % dcoord)
|
| 363 |
if atoms.constraints:
|
| 364 |
for flag in sflags[iatom]: |
| 365 |
if flag:
|
| 366 |
s = 'F'
|
| 367 |
else:
|
| 368 |
s = 'T'
|
| 369 |
f.write('%4s' % s)
|
| 370 |
f.write('\n')
|
| 371 |
|
| 372 |
if type(filename) == str: |
| 373 |
f.close() |