Statistiques
| Révision :

root / ase / calculators / vasp.py @ 10

Historique | Voir | Annoter | Télécharger (53,06 ko)

1 1 tkerber
# Copyright (C) 2008 CSC - Scientific Computing Ltd.
2 1 tkerber
"""This module defines an ASE interface to VASP.
3 1 tkerber

4 1 tkerber
Developed on the basis of modules by Jussi Enkovaara and John
5 1 tkerber
Kitchin.  The path of the directory containing the pseudopotential
6 1 tkerber
directories (potpaw,potpaw_GGA, potpaw_PBE, ...) should be set
7 1 tkerber
by the environmental flag $VASP_PP_PATH.
8 1 tkerber

9 1 tkerber
The user should also set the environmental flag $VASP_SCRIPT pointing
10 1 tkerber
to a python script looking something like::
11 1 tkerber

12 1 tkerber
   import os
13 1 tkerber
   exitcode = os.system('vasp')
14 1 tkerber

15 1 tkerber
Alternatively, user can set the environmental flag $VASP_COMMAND pointing
16 1 tkerber
to the command use the launch vasp e.g. 'vasp' or 'mpirun -n 16 vasp'
17 1 tkerber

18 1 tkerber
http://cms.mpi.univie.ac.at/vasp/
19 1 tkerber

20 1 tkerber
-Jonas Bjork j.bjork@liverpool.ac.uk
21 1 tkerber
"""
22 1 tkerber
23 1 tkerber
import os
24 1 tkerber
import sys
25 1 tkerber
from general import Calculator
26 1 tkerber
from os.path import join, isfile, islink
27 1 tkerber
28 1 tkerber
import numpy as np
29 1 tkerber
30 1 tkerber
import ase
31 1 tkerber
32 1 tkerber
# Parameters that can be set in INCAR. The values which are None
33 1 tkerber
# are not written and default parameters of VASP are used for them.
34 1 tkerber
35 1 tkerber
float_keys = [
36 1 tkerber
    'aexx',       # Fraction of exact/DFT exchange
37 1 tkerber
    'aggac',      # Fraction of gradient correction to correlation
38 1 tkerber
    'aggax',      # Fraction of gradient correction to exchange
39 1 tkerber
    'aldac',      # Fraction of LDA correlation energy
40 1 tkerber
    'amix',       #
41 1 tkerber
    'bmix',       # tags for mixing
42 1 tkerber
    'emax',       # energy-range for DOSCAR file
43 1 tkerber
    'emin',       #
44 1 tkerber
    'enaug',      # Density cutoff
45 1 tkerber
    'encut',      # Planewave cutoff
46 1 tkerber
    'encutfock',  # FFT grid in the HF related routines
47 1 tkerber
    'hfscreen',   # attribute to change from PBE0 to HSE
48 1 tkerber
    'potim',      # time-step for ion-motion (fs)
49 1 tkerber
    'nelect',     # total number of electrons
50 1 tkerber
    'pomass',     # mass of ions in am
51 1 tkerber
    'sigma',      # broadening in eV
52 1 tkerber
    'time',       # special control tag
53 1 tkerber
    'zval',       # ionic valence
54 1 tkerber
    ]
55 1 tkerber
56 1 tkerber
exp_keys = [
57 1 tkerber
    'ediff',      # stopping-criterion for electronic upd.
58 1 tkerber
    'ediffg',     # stopping-criterion for ionic upd.
59 1 tkerber
    'symprec',    # precession in symmetry routines
60 1 tkerber
]
61 1 tkerber
62 1 tkerber
string_keys = [
63 1 tkerber
    'algo',       # algorithm: Normal (Davidson) | Fast | Very_Fast (RMM-DIIS)
64 1 tkerber
    'gga',        # xc-type: PW PB LM or 91
65 1 tkerber
    'prec',       # Precission of calculation (Low, Normal, Accurate)
66 1 tkerber
    'system',     # name of System
67 1 tkerber
    'tebeg',      #
68 1 tkerber
    'teend',      # temperature during run
69 1 tkerber
]
70 1 tkerber
71 1 tkerber
int_keys = [
72 1 tkerber
    'ialgo',      # algorithm: use only 8 (CG) or 48 (RMM-DIIS)
73 1 tkerber
    'ibrion',     # ionic relaxation: 0-MD 1-quasi-New 2-CG
74 1 tkerber
    'idipol',     # monopol/dipol and quadropole corrections
75 1 tkerber
    'iniwav',     # initial electr wf. : 0-lowe 1-rand
76 1 tkerber
    'isif',       # calculate stress and what to relax
77 1 tkerber
    'ismear',     # part. occupancies: -5 Blochl -4-tet -1-fermi 0-gaus >0 MP
78 1 tkerber
    'ispin',      # spin-polarized calculation
79 1 tkerber
    'istart',     # startjob: 0-new 1-cont 2-samecut
80 1 tkerber
    'isym',       # symmetry: 0-nonsym 1-usesym
81 1 tkerber
    'iwavpr',     # prediction of wf.: 0-non 1-charg 2-wave 3-comb
82 1 tkerber
    'lmaxmix',    #
83 1 tkerber
    'lorbit',     # create PROOUT
84 1 tkerber
    'ngx',        # FFT mesh for wavefunctions, x
85 1 tkerber
    'ngxf',       # FFT mesh for charges x
86 1 tkerber
    'ngy',        # FFT mesh for wavefunctions, y
87 1 tkerber
    'ngyf',       # FFT mesh for charges y
88 1 tkerber
    'ngz',        # FFT mesh for wavefunctions, z
89 1 tkerber
    'ngzf',       # FFT mesh for charges z
90 1 tkerber
    'nbands',     # Number of bands
91 1 tkerber
    'nblk',       # blocking for some BLAS calls (Sec. 6.5)
92 1 tkerber
    'nbmod',      # specifies mode for partial charge calculation
93 1 tkerber
    'nelm',       #
94 1 tkerber
    'nelmdl',     # nr. of electronic steps
95 1 tkerber
    'nelmin',
96 1 tkerber
    'nfree',      # number of steps per DOF when calculting Hessian using finitite differences
97 1 tkerber
    'nkred',      # define sub grid of q-points for HF with nkredx=nkredy=nkredz
98 1 tkerber
    'nkredx',      # define sub grid of q-points in x direction for HF
99 1 tkerber
    'nkredy',      # define sub grid of q-points in y direction for HF
100 1 tkerber
    'nkredz',      # define sub grid of q-points in z direction for HF
101 1 tkerber
    'npar',       # parallelization over bands
102 1 tkerber
    'nsim',       # evaluate NSIM bands simultaneously if using RMM-DIIS
103 1 tkerber
    'nsw',        # number of steps for ionic upd.
104 1 tkerber
    'nupdown',    # fix spin moment to specified value
105 1 tkerber
    'nwrite',     # verbosity write-flag (how much is written)
106 1 tkerber
    'smass',      # Nose mass-parameter (am)
107 1 tkerber
    'voskown',    # use Vosko, Wilk, Nusair interpolation
108 1 tkerber
]
109 1 tkerber
110 1 tkerber
bool_keys = [
111 1 tkerber
    'addgrid',    # finer grid for augmentation charge density
112 1 tkerber
    'icharg',     # charge: 1-file 2-atom 10-const
113 1 tkerber
    'lasync',     # overlap communcation with calculations
114 1 tkerber
    'lcharg',     #
115 1 tkerber
    'lcorr',      # Harris-correction to forces
116 1 tkerber
    'ldipol',     # potential correction mode
117 1 tkerber
    'lelf',       # create ELFCAR
118 1 tkerber
    'lhfcalc',    # switch to turn on Hartree Fock calculations
119 1 tkerber
    'lpard',      # evaluate partial (band and/or k-point) decomposed charge density
120 1 tkerber
    'lplane',     # parallelisation over the FFT grid
121 1 tkerber
    'lscalapack', # switch off scaLAPACK
122 1 tkerber
    'lscalu',     # switch of LU decomposition
123 1 tkerber
    'lsepb',      # write out partial charge of each band seperately?
124 1 tkerber
    'lsepk',      # write out partial charge of each k-point seperately?
125 1 tkerber
    'lthomas',    #
126 1 tkerber
    'lvtot',      # create WAVECAR/CHGCAR/LOCPOT
127 1 tkerber
    'lwave',      #
128 1 tkerber
]
129 1 tkerber
130 1 tkerber
list_keys = [
131 1 tkerber
    'dipol',      # center of cell for dipol
132 1 tkerber
    'eint',       # energy range to calculate partial charge for
133 1 tkerber
    'ferwe',      # Fixed band occupation
134 1 tkerber
    'iband',      # bands to calculate partial charge for
135 1 tkerber
    'magmom',     # initial magnetic moments
136 1 tkerber
    'kpuse',      # k-point to calculate partial charge for
137 1 tkerber
    'ropt',       # number of grid points for non-local proj in real space
138 1 tkerber
    'rwigs',      # Wigner-Seitz radii
139 1 tkerber
]
140 1 tkerber
141 1 tkerber
special_keys = [
142 1 tkerber
    'lreal',      # non-local projectors in real space
143 1 tkerber
]
144 1 tkerber
145 1 tkerber
keys = [
146 1 tkerber
    # 'NBLOCK' and KBLOCK       inner block; outer block
147 1 tkerber
    # 'NPACO' and APACO         distance and nr. of slots for P.C.
148 1 tkerber
    # 'WEIMIN, EBREAK, DEPER    special control tags
149 1 tkerber
]
150 1 tkerber
151 1 tkerber
class Vasp(Calculator):
152 2 tkerber
    def __init__(self, restart=None, output_template='vasp', track_output=False, write_input = True,
153 1 tkerber
                 **kwargs):
154 10 tkerber
        self.write_input = write_input
155 1 tkerber
        self.name = 'Vasp'
156 1 tkerber
        self.float_params = {}
157 1 tkerber
        self.exp_params = {}
158 1 tkerber
        self.string_params = {}
159 1 tkerber
        self.int_params = {}
160 1 tkerber
        self.bool_params = {}
161 1 tkerber
        self.list_params = {}
162 1 tkerber
        self.special_params = {}
163 1 tkerber
        for key in float_keys:
164 1 tkerber
            self.float_params[key] = None
165 1 tkerber
        for key in exp_keys:
166 1 tkerber
            self.exp_params[key] = None
167 1 tkerber
        for key in string_keys:
168 1 tkerber
            self.string_params[key] = None
169 1 tkerber
        for key in int_keys:
170 1 tkerber
            self.int_params[key] = None
171 1 tkerber
        for key in bool_keys:
172 1 tkerber
            self.bool_params[key] = None
173 1 tkerber
        for key in list_keys:
174 1 tkerber
            self.list_params[key] = None
175 1 tkerber
        for key in special_keys:
176 1 tkerber
            self.special_params[key] = None
177 1 tkerber
        self.string_params['prec'] = 'Normal'
178 1 tkerber
179 1 tkerber
        self.input_params = {
180 1 tkerber
            'xc':     'PW91',  # exchange correlation potential
181 1 tkerber
            'setups': None,    # Special setups (e.g pv, sv, ...)
182 1 tkerber
            'txt':    '-',     # Where to send information
183 1 tkerber
            'kpts':   (1,1,1), # k-points
184 1 tkerber
            'gamma':  False,   # Option to use gamma-sampling instead
185 1 tkerber
                               # of Monkhorst-Pack
186 1 tkerber
            }
187 1 tkerber
188 1 tkerber
        self.restart = restart
189 1 tkerber
        self.track_output = track_output
190 1 tkerber
        self.output_template = output_template
191 1 tkerber
        if restart:
192 1 tkerber
            self.restart_load()
193 1 tkerber
            return
194 1 tkerber
195 1 tkerber
        if self.input_params['xc'] not in ['PW91','LDA','PBE']:
196 1 tkerber
            raise ValueError(
197 1 tkerber
                '%s not supported for xc! use one of: PW91, LDA or PBE.' %
198 1 tkerber
                kwargs['xc'])
199 1 tkerber
        self.nbands = self.int_params['nbands']
200 1 tkerber
        self.atoms = None
201 1 tkerber
        self.positions = None
202 1 tkerber
        self.run_counts = 0
203 1 tkerber
        self.set(**kwargs)
204 1 tkerber
205 1 tkerber
    def set(self, **kwargs):
206 1 tkerber
        for key in kwargs:
207 1 tkerber
            if self.float_params.has_key(key):
208 1 tkerber
                self.float_params[key] = kwargs[key]
209 1 tkerber
            elif self.exp_params.has_key(key):
210 1 tkerber
                self.exp_params[key] = kwargs[key]
211 1 tkerber
            elif self.string_params.has_key(key):
212 1 tkerber
                self.string_params[key] = kwargs[key]
213 1 tkerber
            elif self.int_params.has_key(key):
214 1 tkerber
                self.int_params[key] = kwargs[key]
215 1 tkerber
            elif self.bool_params.has_key(key):
216 1 tkerber
                self.bool_params[key] = kwargs[key]
217 1 tkerber
            elif self.list_params.has_key(key):
218 1 tkerber
                self.list_params[key] = kwargs[key]
219 1 tkerber
            elif self.special_params.has_key(key):
220 1 tkerber
                self.special_params[key] = kwargs[key]
221 1 tkerber
            elif self.input_params.has_key(key):
222 1 tkerber
                self.input_params[key] = kwargs[key]
223 1 tkerber
            else:
224 1 tkerber
                raise TypeError('Parameter not defined: ' + key)
225 1 tkerber
226 1 tkerber
    def update(self, atoms):
227 1 tkerber
        if self.calculation_required(atoms, ['energy']):
228 1 tkerber
            if (self.atoms is None or
229 1 tkerber
                self.atoms.positions.shape != atoms.positions.shape):
230 1 tkerber
                # Completely new calculation just reusing the same
231 1 tkerber
                # calculator, so delete any old VASP files found.
232 1 tkerber
                self.clean()
233 1 tkerber
            self.calculate(atoms)
234 1 tkerber
235 1 tkerber
    def initialize(self, atoms):
236 1 tkerber
        """Initialize a VASP calculation
237 1 tkerber

238 1 tkerber
        Constructs the POTCAR file. User should specify the PATH
239 1 tkerber
        to the pseudopotentials in VASP_PP_PATH environment variable"""
240 1 tkerber
241 1 tkerber
        p = self.input_params
242 1 tkerber
243 1 tkerber
        self.all_symbols = atoms.get_chemical_symbols()
244 1 tkerber
        self.natoms = len(atoms)
245 1 tkerber
        self.spinpol = atoms.get_initial_magnetic_moments().any()
246 1 tkerber
        atomtypes = atoms.get_chemical_symbols()
247 1 tkerber
248 1 tkerber
        # Determine the number of atoms of each atomic species
249 1 tkerber
        # sorted after atomic species
250 1 tkerber
        special_setups = []
251 1 tkerber
        if self.input_params['setups']:
252 1 tkerber
            for m in self.input_params['setups']:
253 1 tkerber
                try :
254 1 tkerber
                    #special_setup[self.input_params['setups'][m]] = int(m)
255 1 tkerber
                    special_setups.append(int(m))
256 1 tkerber
                except:
257 1 tkerber
                    #print 'setup ' + m + ' is a groups setup'
258 1 tkerber
                    continue
259 1 tkerber
            #print 'special_setups' , special_setups
260 1 tkerber
261 10 tkerber
        #symbols = {}
262 10 tkerber
        #for m,atom in enumerate(atoms):
263 10 tkerber
        #    symbol = atom.get_symbol()
264 10 tkerber
        #    if m in special_setups:
265 10 tkerber
        #        pass
266 10 tkerber
        #    else:
267 10 tkerber
        #        if not symbols.has_key(symbol):
268 10 tkerber
        #            symbols[symbol] = 1
269 10 tkerber
        #        else:
270 10 tkerber
        #            symbols[symbol] += 1
271 10 tkerber
272 10 tkerber
        symbols = []
273 1 tkerber
        for m,atom in enumerate(atoms):
274 1 tkerber
            symbol = atom.get_symbol()
275 1 tkerber
            if m in special_setups:
276 1 tkerber
                pass
277 1 tkerber
            else:
278 10 tkerber
                ind = 0
279 10 tkerber
                bfound = False
280 10 tkerber
                for symbolX, iatom in symbols:
281 10 tkerber
                    if symbol == symbolX:
282 10 tkerber
                        bfound = True
283 10 tkerber
                        break
284 10 tkerber
                    ind += 1
285 10 tkerber
                if bfound:
286 10 tkerber
                    symbolX, iatom  = symbols[ind]
287 10 tkerber
                    symbols[ind] = symbolX, iatom+1
288 1 tkerber
                else:
289 10 tkerber
                    symbols.append((symbol, 1))
290 1 tkerber
291 1 tkerber
        # Build the sorting list
292 1 tkerber
        self.sort = []
293 1 tkerber
        self.sort.extend(special_setups)
294 1 tkerber
295 10 tkerber
        for symbol, iatom in symbols:
296 1 tkerber
            for m,atom in enumerate(atoms):
297 1 tkerber
                if m in special_setups:
298 1 tkerber
                    pass
299 1 tkerber
                else:
300 1 tkerber
                    if atom.get_symbol() == symbol:
301 1 tkerber
                        self.sort.append(m)
302 1 tkerber
        self.resort = range(len(self.sort))
303 1 tkerber
        for n in range(len(self.resort)):
304 1 tkerber
            self.resort[self.sort[n]] = n
305 1 tkerber
        self.atoms_sorted = atoms[self.sort]
306 10 tkerber
307 1 tkerber
        # Check if the necessary POTCAR files exists and
308 1 tkerber
        # create a list of their paths.
309 1 tkerber
        self.symbol_count = []
310 1 tkerber
        for m in special_setups:
311 10 tkerber
            for symbol, iatom in symbols:
312 10 tkerber
                if symbol == m:
313 10 tkerber
                    self.symbol_count.append([atomtypes[m],1])
314 10 tkerber
            #self.symbol_count.append([atomtypes[m],1])
315 1 tkerber
        for m in symbols:
316 10 tkerber
            for symbol, iatom in symbols:
317 10 tkerber
                if symbol == m:
318 10 tkerber
                    self.symbol_count.append([m,iatom])
319 10 tkerber
            #self.symbol_count.append([m,symbols[m]])
320 1 tkerber
        sys.stdout.flush()
321 10 tkerber
        if self.write_input:
322 10 tkerber
            xc = '/'
323 10 tkerber
            #print 'p[xc]',p['xc']
324 10 tkerber
            if p['xc'] == 'PW91':
325 10 tkerber
                xc = '_gga/'
326 10 tkerber
            elif p['xc'] == 'PBE':
327 10 tkerber
                xc = '_pbe/'
328 10 tkerber
            if 'VASP_PP_PATH' in os.environ:
329 10 tkerber
                pppaths = os.environ['VASP_PP_PATH'].split(':')
330 10 tkerber
            else:
331 10 tkerber
                pppaths = []
332 10 tkerber
            self.ppp_list = []
333 10 tkerber
            # Setting the pseudopotentials, first special setups and
334 10 tkerber
            # then according to symbols
335 10 tkerber
            for m in special_setups:
336 10 tkerber
                name = 'potpaw'+xc.upper() + p['setups'][str(m)] + '/POTCAR'
337 10 tkerber
                found = False
338 10 tkerber
                for path in pppaths:
339 10 tkerber
                    filename = join(path, name)
340 10 tkerber
                    #print 'filename', filename
341 10 tkerber
                    if isfile(filename) or islink(filename):
342 10 tkerber
                        found = True
343 10 tkerber
                        self.ppp_list.append(filename)
344 10 tkerber
                        break
345 10 tkerber
                    elif isfile(filename + '.Z') or islink(filename + '.Z'):
346 10 tkerber
                        found = True
347 10 tkerber
                        self.ppp_list.append(filename+'.Z')
348 10 tkerber
                        break
349 10 tkerber
                if not found:
350 10 tkerber
                    raise RuntimeError('No pseudopotential for %s!' % symbol)
351 10 tkerber
            #print 'symbols', symbols
352 10 tkerber
            for symbol in symbols:
353 10 tkerber
                try:
354 10 tkerber
                    name = 'potpaw'+xc.upper()+symbol + p['setups'][symbol]
355 10 tkerber
                except (TypeError, KeyError):
356 10 tkerber
                    name = 'potpaw' + xc.upper() + symbol
357 10 tkerber
                name += '/POTCAR'
358 10 tkerber
                found = False
359 10 tkerber
                for path in pppaths:
360 10 tkerber
                    filename = join(path, name)
361 10 tkerber
                    #print 'filename', filename
362 10 tkerber
                    if isfile(filename) or islink(filename):
363 10 tkerber
                        found = True
364 10 tkerber
                        self.ppp_list.append(filename)
365 10 tkerber
                        break
366 10 tkerber
                    elif isfile(filename + '.Z') or islink(filename + '.Z'):
367 10 tkerber
                        found = True
368 10 tkerber
                        self.ppp_list.append(filename+'.Z')
369 10 tkerber
                        break
370 10 tkerber
                if not found:
371 10 tkerber
                    raise RuntimeError('No pseudopotential for %s!' % symbol)
372 1 tkerber
        self.converged = None
373 1 tkerber
        self.setups_changed = None
374 1 tkerber
375 1 tkerber
    def calculate(self, atoms):
376 1 tkerber
        """Generate necessary files in the working directory and run VASP.
377 1 tkerber

378 1 tkerber
        The method first write VASP input files, then calls the method
379 1 tkerber
        which executes VASP. When the VASP run is finished energy, forces,
380 1 tkerber
        etc. are read from the VASP output.
381 1 tkerber
        """
382 1 tkerber
383 1 tkerber
        # Write input
384 1 tkerber
        from ase.io.vasp import write_vasp
385 1 tkerber
        self.initialize(atoms)
386 1 tkerber
        write_vasp('POSCAR', self.atoms_sorted, symbol_count = self.symbol_count)
387 10 tkerber
        if self.write_input:
388 2 tkerber
            self.write_incar(atoms)
389 2 tkerber
            self.write_potcar()
390 2 tkerber
            self.write_kpoints()
391 1 tkerber
        self.write_sort_file()
392 1 tkerber
393 1 tkerber
        # Execute VASP
394 1 tkerber
        self.run()
395 1 tkerber
        # Read output
396 1 tkerber
        atoms_sorted = ase.io.read('CONTCAR', format='vasp')
397 1 tkerber
        if self.int_params['ibrion']>-1 and self.int_params['nsw']>0:
398 1 tkerber
            # Replace entire atoms object with the one read from CONTCAR.
399 1 tkerber
            atoms.__dict__ = atoms_sorted[self.resort].__dict__
400 1 tkerber
        self.converged = self.read_convergence()
401 1 tkerber
        self.set_results(atoms)
402 1 tkerber
403 1 tkerber
    def set_results(self, atoms):
404 1 tkerber
        self.read(atoms)
405 1 tkerber
        if self.spinpol:
406 1 tkerber
            self.magnetic_moment = self.read_magnetic_moment()
407 1 tkerber
            if self.int_params['lorbit']>=10 or (self.int_params['lorbit']!=None and self.list_params['rwigs']):
408 1 tkerber
                self.magnetic_moments = self.read_magnetic_moments(atoms)
409 1 tkerber
        self.old_float_params = self.float_params.copy()
410 1 tkerber
        self.old_exp_params = self.exp_params.copy()
411 1 tkerber
        self.old_string_params = self.string_params.copy()
412 1 tkerber
        self.old_int_params = self.int_params.copy()
413 1 tkerber
        self.old_input_params = self.input_params.copy()
414 1 tkerber
        self.old_bool_params = self.bool_params.copy()
415 1 tkerber
        self.old_list_params = self.list_params.copy()
416 1 tkerber
        self.atoms = atoms.copy()
417 1 tkerber
418 1 tkerber
    def run(self):
419 1 tkerber
        """Method which explicitely runs VASP."""
420 1 tkerber
421 1 tkerber
        if self.track_output:
422 1 tkerber
            self.out = self.output_template+str(self.run_counts)+'.out'
423 1 tkerber
            self.run_counts += 1
424 1 tkerber
        else:
425 1 tkerber
            self.out = self.output_template+'.out'
426 1 tkerber
        stderr = sys.stderr
427 1 tkerber
        p=self.input_params
428 1 tkerber
        if p['txt'] is None:
429 1 tkerber
            sys.stderr = devnull
430 1 tkerber
        elif p['txt'] == '-':
431 1 tkerber
            pass
432 1 tkerber
        elif isinstance(p['txt'], str):
433 1 tkerber
            sys.stderr = open(p['txt'], 'w')
434 1 tkerber
        if os.environ.has_key('VASP_COMMAND'):
435 1 tkerber
            vasp = os.environ['VASP_COMMAND']
436 1 tkerber
            exitcode = os.system('%s > %s' % (vasp, self.out))
437 1 tkerber
        elif os.environ.has_key('VASP_SCRIPT'):
438 1 tkerber
            vasp = os.environ['VASP_SCRIPT']
439 1 tkerber
            locals={}
440 1 tkerber
            execfile(vasp, {}, locals)
441 1 tkerber
            exitcode = locals['exitcode']
442 1 tkerber
        else:
443 1 tkerber
            raise RuntimeError('Please set either VASP_COMMAND or VASP_SCRIPT environment variable')
444 1 tkerber
        sys.stderr = stderr
445 1 tkerber
        if exitcode != 0:
446 1 tkerber
            raise RuntimeError('Vasp exited with exit code: %d.  ' % exitcode)
447 1 tkerber
448 1 tkerber
    def restart_load(self):
449 1 tkerber
        """Method which is called upon restart."""
450 1 tkerber
451 1 tkerber
        # Try to read sorting file
452 1 tkerber
        if os.path.isfile('ase-sort.dat'):
453 1 tkerber
            self.sort = []
454 1 tkerber
            self.resort = []
455 1 tkerber
            file = open('ase-sort.dat', 'r')
456 1 tkerber
            lines = file.readlines()
457 1 tkerber
            file.close()
458 1 tkerber
            for line in lines:
459 1 tkerber
                data = line.split()
460 1 tkerber
                self.sort.append(int(data[0]))
461 1 tkerber
                self.resort.append(int(data[1]))
462 1 tkerber
            atoms = ase.io.read('CONTCAR', format='vasp')[self.resort]
463 1 tkerber
        else:
464 1 tkerber
            atoms = ase.io.read('CONTCAR', format='vasp')
465 1 tkerber
            self.sort = range(len(atoms))
466 1 tkerber
            self.resort = range(len(atoms))
467 1 tkerber
        self.atoms = atoms.copy()
468 1 tkerber
        self.read_incar()
469 1 tkerber
        self.read_outcar()
470 1 tkerber
        self.set_results(atoms)
471 1 tkerber
        self.read_kpoints()
472 1 tkerber
        self.read_potcar()
473 1 tkerber
#        self.old_incar_params = self.incar_params.copy()
474 1 tkerber
        self.old_input_params = self.input_params.copy()
475 1 tkerber
        self.converged = self.read_convergence()
476 1 tkerber
477 1 tkerber
    def clean(self):
478 1 tkerber
        """Method which cleans up after a calculation.
479 1 tkerber

480 1 tkerber
        The default files generated by Vasp will be deleted IF this
481 1 tkerber
        method is called.
482 1 tkerber

483 1 tkerber
        """
484 1 tkerber
        files = ['CHG', 'CHGCAR', 'POSCAR', 'INCAR', 'CONTCAR', 'DOSCAR',
485 1 tkerber
                 'EIGENVAL', 'IBZKPT', 'KPOINTS', 'OSZICAR', 'OUTCAR', 'PCDAT',
486 1 tkerber
                 'POTCAR', 'vasprun.xml', 'WAVECAR', 'XDATCAR',
487 1 tkerber
                 'PROCAR', 'ase-sort.dat']
488 1 tkerber
        for f in files:
489 1 tkerber
            try:
490 1 tkerber
                os.remove(f)
491 1 tkerber
            except OSError:
492 1 tkerber
                pass
493 1 tkerber
494 1 tkerber
    def set_atoms(self, atoms):
495 1 tkerber
        if (atoms != self.atoms):
496 1 tkerber
            self.converged = None
497 1 tkerber
        self.atoms = atoms.copy()
498 1 tkerber
499 1 tkerber
    def get_atoms(self):
500 1 tkerber
        atoms = self.atoms.copy()
501 1 tkerber
        atoms.set_calculator(self)
502 1 tkerber
        return atoms
503 1 tkerber
504 1 tkerber
    def get_potential_energy(self, atoms, force_consistent=False):
505 1 tkerber
        self.update(atoms)
506 1 tkerber
        if force_consistent:
507 1 tkerber
            return self.energy_free
508 1 tkerber
        else:
509 1 tkerber
            return self.energy_zero
510 1 tkerber
511 1 tkerber
    def get_forces(self, atoms):
512 1 tkerber
        self.update(atoms)
513 1 tkerber
        return self.forces
514 1 tkerber
515 1 tkerber
    def get_stress(self, atoms):
516 1 tkerber
        self.update(atoms)
517 1 tkerber
        return self.stress
518 1 tkerber
519 1 tkerber
    def read_stress(self):
520 1 tkerber
        for line in open('OUTCAR'):
521 1 tkerber
            if line.find(' in kB  ') != -1:
522 1 tkerber
                stress = -np.array([float(a) for a in line.split()[2:]]) \
523 1 tkerber
                         [[0, 1, 2, 4, 5, 3]] \
524 1 tkerber
                         * 1e-1 * ase.units.GPa
525 1 tkerber
        return stress
526 1 tkerber
527 1 tkerber
    def calculation_required(self, atoms, quantities):
528 1 tkerber
        if (self.positions is None or
529 1 tkerber
            (self.atoms != atoms) or
530 1 tkerber
            (self.float_params != self.old_float_params) or
531 1 tkerber
            (self.exp_params != self.old_exp_params) or
532 1 tkerber
            (self.string_params != self.old_string_params) or
533 1 tkerber
            (self.int_params != self.old_int_params) or
534 1 tkerber
            (self.bool_params != self.old_bool_params) or
535 1 tkerber
            (self.list_params != self.old_list_params) or
536 1 tkerber
            (self.input_params != self.old_input_params)
537 1 tkerber
            or not self.converged):
538 1 tkerber
            return True
539 1 tkerber
        if 'magmom' in quantities:
540 1 tkerber
            return not hasattr(self, 'magnetic_moment')
541 1 tkerber
        return False
542 1 tkerber
543 1 tkerber
    def get_number_of_bands(self):
544 1 tkerber
        return self.nbands
545 1 tkerber
546 1 tkerber
    def get_k_point_weights(self):
547 1 tkerber
        self.update(self.atoms)
548 1 tkerber
        return self.read_k_point_weights()
549 1 tkerber
550 1 tkerber
    def get_number_of_spins(self):
551 1 tkerber
        return 1 + int(self.spinpol)
552 1 tkerber
553 1 tkerber
    def get_eigenvalues(self, kpt=0, spin=0):
554 1 tkerber
        self.update(self.atoms)
555 1 tkerber
        return self.read_eigenvalues(kpt, spin)
556 1 tkerber
557 1 tkerber
    def get_fermi_level(self):
558 1 tkerber
        return self.fermi
559 1 tkerber
560 1 tkerber
    def get_number_of_grid_points(self):
561 1 tkerber
        raise NotImplementedError
562 1 tkerber
563 1 tkerber
    def get_pseudo_density(self):
564 1 tkerber
        raise NotImplementedError
565 1 tkerber
566 1 tkerber
    def get_pseudo_wavefunction(self, n=0, k=0, s=0, pad=True):
567 1 tkerber
        raise NotImplementedError
568 1 tkerber
569 1 tkerber
    def get_bz_k_points(self):
570 1 tkerber
        raise NotImplementedError
571 1 tkerber
572 1 tkerber
    def get_ibz_kpoints(self):
573 1 tkerber
        self.update(self.atoms)
574 1 tkerber
        return self.read_ibz_kpoints()
575 1 tkerber
576 1 tkerber
    def get_spin_polarized(self):
577 1 tkerber
        if not hasattr(self, 'spinpol'):
578 1 tkerber
            self.spinpol = self.atoms.get_initial_magnetic_moments().any()
579 1 tkerber
        return self.spinpol
580 1 tkerber
581 1 tkerber
    def get_magnetic_moment(self, atoms):
582 1 tkerber
        self.update(atoms)
583 1 tkerber
        return self.magnetic_moment
584 1 tkerber
585 1 tkerber
    def get_magnetic_moments(self, atoms):
586 1 tkerber
        if self.int_params['lorbit']>=10 or self.list_params['rwigs']:
587 1 tkerber
            self.update(atoms)
588 1 tkerber
            return self.magnetic_moments
589 1 tkerber
        else:
590 1 tkerber
            raise RuntimeError(
591 1 tkerber
                "The combination %s for lorbit with %s for rwigs not supported to calculate magnetic moments" % (p['lorbit'], p['rwigs']))
592 1 tkerber
593 1 tkerber
    def get_dipole_moment(self, atoms):
594 1 tkerber
        """Returns total dipole moment of the system."""
595 1 tkerber
        self.update(atoms)
596 1 tkerber
        return self.dipole
597 1 tkerber
598 1 tkerber
    def get_xc_functional(self):
599 1 tkerber
        return self.input_params['xc']
600 1 tkerber
601 1 tkerber
    def write_incar(self, atoms, **kwargs):
602 1 tkerber
        """Writes the INCAR file."""
603 1 tkerber
        incar = open('INCAR', 'w')
604 1 tkerber
        incar.write('INCAR created by Atomic Simulation Environment\n')
605 1 tkerber
        for key, val in self.float_params.items():
606 1 tkerber
            if val is not None:
607 1 tkerber
                incar.write(' %s = %5.6f\n' % (key.upper(), val))
608 1 tkerber
        for key, val in self.exp_params.items():
609 1 tkerber
            if val is not None:
610 1 tkerber
                incar.write(' %s = %5.2e\n' % (key.upper(), val))
611 1 tkerber
        for key, val in self.string_params.items():
612 1 tkerber
            if val is not None:
613 1 tkerber
                incar.write(' %s = %s\n' % (key.upper(), val))
614 1 tkerber
        for key, val in self.int_params.items():
615 1 tkerber
            if val is not None:
616 1 tkerber
                incar.write(' %s = %d\n' % (key.upper(), val))
617 1 tkerber
        for key, val in self.list_params.items():
618 1 tkerber
            if val is not None:
619 1 tkerber
                incar.write(' %s = ' % key.upper())
620 1 tkerber
                if key in ('dipol', 'eint', 'ferwe', 'ropt', 'rwigs'):
621 1 tkerber
                    [incar.write('%.4f ' % x) for x in val]
622 1 tkerber
                elif key in ('iband', 'kpuse'):
623 1 tkerber
                    [incar.write('%i ' % x) for x in val]
624 1 tkerber
                elif key == 'magmom':
625 1 tkerber
                    list = [[1, val[0]]]
626 1 tkerber
                    for n in range(1, len(val)):
627 1 tkerber
                        if val[n] == val[n-1]:
628 1 tkerber
                            list[-1][0] += 1
629 1 tkerber
                        else:
630 1 tkerber
                            list.append([1, val[n]])
631 1 tkerber
                    [incar.write('%i*%.4f ' % (mom[0], mom[1])) for mom in list]
632 1 tkerber
                incar.write('\n')
633 1 tkerber
        for key, val in self.bool_params.items():
634 1 tkerber
            if val is not None:
635 1 tkerber
                incar.write(' %s = ' % key.upper())
636 1 tkerber
                if val:
637 1 tkerber
                    incar.write('.TRUE.\n')
638 1 tkerber
                else:
639 1 tkerber
                    incar.write('.FALSE.\n')
640 1 tkerber
        for key, val in self.special_params.items():
641 1 tkerber
            if val is not None:
642 1 tkerber
                incar.write(' %s = ' % key.upper())
643 1 tkerber
                if key is 'lreal':
644 1 tkerber
                    if type(val)==type('str'):
645 1 tkerber
                        incar.write(val+'\n')
646 1 tkerber
                    elif type(val)==type(True):
647 1 tkerber
                       if val:
648 1 tkerber
                           incar.write('.TRUE.\n')
649 1 tkerber
                       else:
650 1 tkerber
                           incar.write('.FALSE.\n')
651 1 tkerber
        if self.spinpol and not self.int_params['ispin']:
652 1 tkerber
            incar.write(' ispin = 2\n'.upper())
653 1 tkerber
            # Write out initial magnetic moments
654 1 tkerber
            magmom = atoms.get_initial_magnetic_moments()[self.sort]
655 1 tkerber
            list = [[1, magmom[0]]]
656 1 tkerber
            for n in range(1, len(magmom)):
657 1 tkerber
                if magmom[n] == magmom[n-1]:
658 1 tkerber
                    list[-1][0] += 1
659 1 tkerber
                else:
660 1 tkerber
                    list.append([1, magmom[n]])
661 1 tkerber
            incar.write(' magmom = '.upper())
662 1 tkerber
            [incar.write('%i*%.4f ' % (mom[0], mom[1])) for mom in list]
663 1 tkerber
            incar.write('\n')
664 1 tkerber
        incar.close()
665 1 tkerber
666 1 tkerber
    def write_kpoints(self, **kwargs):
667 1 tkerber
        """Writes the KPOINTS file."""
668 1 tkerber
        p = self.input_params
669 1 tkerber
        kpoints = open('KPOINTS', 'w')
670 1 tkerber
        kpoints.write('KPOINTS created by Atomic Simulation Environemnt\n')
671 1 tkerber
        shape=np.array(p['kpts']).shape
672 1 tkerber
        if len(shape)==1:
673 1 tkerber
            kpoints.write('0\n')
674 1 tkerber
            if p['gamma']:
675 1 tkerber
                kpoints.write('Gamma\n')
676 1 tkerber
            else:
677 1 tkerber
                kpoints.write('Monkhorst-Pack\n')
678 1 tkerber
            [kpoints.write('%i ' % kpt) for kpt in p['kpts']]
679 1 tkerber
            kpoints.write('\n0 0 0')
680 1 tkerber
        elif len(shape)==2:
681 1 tkerber
            kpoints.write('%i \n' % (len(p['kpts'])))
682 1 tkerber
            kpoints.write('Cartesian\n')
683 1 tkerber
            for n in range(len(p['kpts'])):
684 1 tkerber
                [kpoints.write('%f ' % kpt) for kpt in p['kpts'][n]]
685 1 tkerber
                if shape[1]==4:
686 1 tkerber
                    kpoints.write('\n')
687 1 tkerber
                elif shape[1]==3:
688 1 tkerber
                    kpoints.write('1.0 \n')
689 1 tkerber
        kpoints.close()
690 1 tkerber
691 1 tkerber
    def write_potcar(self,suffix = ""):
692 1 tkerber
        """Writes the POTCAR file."""
693 1 tkerber
        import tempfile
694 1 tkerber
        potfile = open('POTCAR'+suffix,'w')
695 1 tkerber
        for filename in self.ppp_list:
696 1 tkerber
            if filename.endswith('R'):
697 1 tkerber
                for line in open(filename, 'r'):
698 1 tkerber
                    potfile.write(line)
699 1 tkerber
            elif filename.endswith('.Z'):
700 1 tkerber
                file_tmp = tempfile.NamedTemporaryFile()
701 1 tkerber
                os.system('gunzip -c %s > %s' % (filename, file_tmp.name))
702 1 tkerber
                for line in file_tmp.readlines():
703 1 tkerber
                    potfile.write(line)
704 1 tkerber
                file_tmp.close()
705 1 tkerber
        potfile.close()
706 1 tkerber
707 1 tkerber
    def write_sort_file(self):
708 1 tkerber
        """Writes a sortings file.
709 1 tkerber

710 1 tkerber
        This file contains information about how the atoms are sorted in
711 1 tkerber
        the first column and how they should be resorted in the second
712 1 tkerber
        column. It is used for restart purposes to get sorting right
713 1 tkerber
        when reading in an old calculation to ASE."""
714 1 tkerber
715 1 tkerber
        file = open('ase-sort.dat', 'w')
716 1 tkerber
        for n in range(len(self.sort)):
717 1 tkerber
            file.write('%5i %5i \n' % (self.sort[n], self.resort[n]))
718 1 tkerber
719 1 tkerber
    # Methods for reading information from OUTCAR files:
720 1 tkerber
    def read_energy(self, all=None):
721 1 tkerber
        [energy_free, energy_zero]=[0, 0]
722 1 tkerber
        if all:
723 1 tkerber
            energy_free = []
724 1 tkerber
            energy_zero = []
725 1 tkerber
        for line in open('OUTCAR', 'r'):
726 1 tkerber
            # Free energy
727 1 tkerber
            if line.startswith('  free  energy   toten'):
728 1 tkerber
                if all:
729 1 tkerber
                    energy_free.append(float(line.split()[-2]))
730 1 tkerber
                else:
731 1 tkerber
                    energy_free = float(line.split()[-2])
732 1 tkerber
            # Extrapolated zero point energy
733 1 tkerber
            if line.startswith('  energy  without entropy'):
734 1 tkerber
                if all:
735 1 tkerber
                    energy_zero.append(float(line.split()[-1]))
736 1 tkerber
                else:
737 1 tkerber
                    energy_zero = float(line.split()[-1])
738 1 tkerber
        return [energy_free, energy_zero]
739 1 tkerber
740 1 tkerber
    def read_forces(self, atoms, all=False):
741 1 tkerber
        """Method that reads forces from OUTCAR file.
742 1 tkerber

743 1 tkerber
        If 'all' is switched on, the forces for all ionic steps
744 1 tkerber
        in the OUTCAR file be returned, in other case only the
745 1 tkerber
        forces for the last ionic configuration is returned."""
746 1 tkerber
747 1 tkerber
        file = open('OUTCAR','r')
748 1 tkerber
        lines = file.readlines()
749 1 tkerber
        file.close()
750 1 tkerber
        n=0
751 1 tkerber
        if all:
752 1 tkerber
            all_forces = []
753 1 tkerber
        for line in lines:
754 1 tkerber
            if line.rfind('TOTAL-FORCE') > -1:
755 1 tkerber
                forces=[]
756 1 tkerber
                for i in range(len(atoms)):
757 1 tkerber
                    forces.append(np.array([float(f) for f in lines[n+2+i].split()[3:6]]))
758 1 tkerber
                if all:
759 1 tkerber
                    all_forces.append(np.array(forces)[self.resort])
760 1 tkerber
            n+=1
761 1 tkerber
        if all:
762 1 tkerber
            return np.array(all_forces)
763 1 tkerber
        else:
764 1 tkerber
            return np.array(forces)[self.resort]
765 1 tkerber
766 1 tkerber
    def read_fermi(self):
767 1 tkerber
        """Method that reads Fermi energy from OUTCAR file"""
768 1 tkerber
        E_f=None
769 1 tkerber
        for line in open('OUTCAR', 'r'):
770 1 tkerber
            if line.rfind('E-fermi') > -1:
771 1 tkerber
                E_f=float(line.split()[2])
772 1 tkerber
        return E_f
773 1 tkerber
774 1 tkerber
    def read_dipole(self):
775 1 tkerber
        dipolemoment=np.zeros([1,3])
776 1 tkerber
        for line in open('OUTCAR', 'r'):
777 1 tkerber
            if line.rfind('dipolmoment') > -1:
778 1 tkerber
                dipolemoment=np.array([float(f) for f in line.split()[1:4]])
779 1 tkerber
        return dipolemoment
780 1 tkerber
781 1 tkerber
    def read_magnetic_moments(self, atoms):
782 1 tkerber
        magnetic_moments = np.zeros(len(atoms))
783 1 tkerber
        n = 0
784 1 tkerber
        lines = open('OUTCAR', 'r').readlines()
785 1 tkerber
        for line in lines:
786 1 tkerber
            if line.rfind('magnetization (x)') > -1:
787 1 tkerber
                for m in range(len(atoms)):
788 1 tkerber
                    magnetic_moments[m] = float(lines[n + m + 4].split()[4])
789 1 tkerber
            n += 1
790 1 tkerber
        return np.array(magnetic_moments)[self.resort]
791 1 tkerber
792 1 tkerber
    def read_magnetic_moment(self):
793 1 tkerber
        n=0
794 1 tkerber
        for line in open('OUTCAR','r'):
795 1 tkerber
            if line.rfind('number of electron  ') > -1:
796 1 tkerber
                magnetic_moment=float(line.split()[-1])
797 1 tkerber
            n+=1
798 1 tkerber
        return magnetic_moment
799 1 tkerber
800 1 tkerber
    def read_nbands(self):
801 1 tkerber
        for line in open('OUTCAR', 'r'):
802 1 tkerber
            if line.rfind('NBANDS') > -1:
803 1 tkerber
                return int(line.split()[-1])
804 1 tkerber
805 1 tkerber
    def read_convergence(self):
806 1 tkerber
        """Method that checks whether a calculation has converged."""
807 1 tkerber
        converged = None
808 1 tkerber
        # First check electronic convergence
809 1 tkerber
        for line in open('OUTCAR', 'r'):
810 1 tkerber
            if line.rfind('EDIFF  ') > -1:
811 1 tkerber
                ediff = float(line.split()[2])
812 1 tkerber
            if line.rfind('total energy-change')>-1:
813 1 tkerber
                split = line.split(':')
814 1 tkerber
                a = float(split[1].split('(')[0])
815 1 tkerber
                b = float(split[1].split('(')[1][0:-2])
816 1 tkerber
                if [abs(a), abs(b)] < [ediff, ediff]:
817 1 tkerber
                    converged = True
818 1 tkerber
                else:
819 1 tkerber
                    converged = False
820 1 tkerber
                    continue
821 1 tkerber
        # Then if ibrion > 0, check whether ionic relaxation condition been fulfilled
822 1 tkerber
        if self.int_params['ibrion'] > 0:
823 1 tkerber
            if not self.read_relaxed():
824 1 tkerber
                converged = False
825 1 tkerber
            else:
826 1 tkerber
                converged = True
827 1 tkerber
        return converged
828 1 tkerber
829 1 tkerber
    def read_ibz_kpoints(self):
830 1 tkerber
        lines = open('OUTCAR', 'r').readlines()
831 1 tkerber
        ibz_kpts = []
832 1 tkerber
        n = 0
833 1 tkerber
        i = 0
834 1 tkerber
        for line in lines:
835 1 tkerber
            if line.rfind('Following cartesian coordinates')>-1:
836 1 tkerber
                m = n+2
837 1 tkerber
                while i==0:
838 1 tkerber
                    ibz_kpts.append([float(lines[m].split()[p]) for p in range(3)])
839 1 tkerber
                    m += 1
840 1 tkerber
                    if lines[m]==' \n':
841 1 tkerber
                        i = 1
842 1 tkerber
            if i == 1:
843 1 tkerber
                continue
844 1 tkerber
            n += 1
845 1 tkerber
        ibz_kpts = np.array(ibz_kpts)
846 1 tkerber
        return np.array(ibz_kpts)
847 1 tkerber
848 1 tkerber
    def read_k_point_weights(self):
849 1 tkerber
        file = open('IBZKPT')
850 1 tkerber
        lines = file.readlines()
851 1 tkerber
        file.close()
852 1 tkerber
        kpt_weights = []
853 1 tkerber
        for n in range(3, len(lines)):
854 1 tkerber
            kpt_weights.append(float(lines[n].split()[3]))
855 1 tkerber
        kpt_weights = np.array(kpt_weights)
856 1 tkerber
        kpt_weights /= np.sum(kpt_weights)
857 1 tkerber
        return kpt_weights
858 1 tkerber
859 1 tkerber
    def read_eigenvalues(self, kpt=0, spin=0):
860 1 tkerber
        file = open('EIGENVAL', 'r')
861 1 tkerber
        lines = file.readlines()
862 1 tkerber
        file.close()
863 1 tkerber
        eigs = []
864 1 tkerber
        for n in range(8+kpt*(self.nbands+2), 8+kpt*(self.nbands+2)+self.nbands):
865 1 tkerber
            eigs.append(float(lines[n].split()[spin+1]))
866 1 tkerber
        return np.array(eigs)
867 1 tkerber
868 1 tkerber
    def read_relaxed(self):
869 1 tkerber
        for line in open('OUTCAR', 'r'):
870 1 tkerber
            if line.rfind('reached required accuracy') > -1:
871 1 tkerber
                return True
872 1 tkerber
        return False
873 1 tkerber
874 1 tkerber
# The below functions are used to restart a calculation and are under early constructions
875 1 tkerber
876 1 tkerber
    def read_incar(self, filename='INCAR'):
877 1 tkerber
        """Method that imports settings from INCAR file."""
878 1 tkerber
879 1 tkerber
        self.spinpol = False
880 1 tkerber
        file=open(filename, 'r')
881 1 tkerber
        file.readline()
882 1 tkerber
        lines=file.readlines()
883 1 tkerber
        for line in lines:
884 1 tkerber
            try:
885 1 tkerber
                data = line.split()
886 1 tkerber
                # Skip empty and commented lines.
887 1 tkerber
                if len(data) == 0:
888 1 tkerber
                    continue
889 1 tkerber
                elif data[0][0] in ['#', '!']:
890 1 tkerber
                    continue
891 1 tkerber
                key = data[0].lower()
892 1 tkerber
                if key in float_keys:
893 1 tkerber
                    self.float_params[key] = float(data[2])
894 1 tkerber
                elif key in exp_keys:
895 1 tkerber
                    self.exp_params[key] = float(data[2])
896 1 tkerber
                elif key in string_keys:
897 1 tkerber
                    self.string_params[key] = str(data[2])
898 1 tkerber
                elif key in int_keys:
899 1 tkerber
                    if key == 'ispin':
900 1 tkerber
                        if int(data[2]) == 2:
901 1 tkerber
                            self.spinpol = True
902 1 tkerber
                    else:
903 1 tkerber
                        self.int_params[key] = int(data[2])
904 1 tkerber
                elif key in bool_keys:
905 1 tkerber
                    if 'true' in data[2].lower():
906 1 tkerber
                        self.bool_params[key] = True
907 1 tkerber
                    elif 'false' in data[2].lower():
908 1 tkerber
                        self.bool_params[key] = False
909 1 tkerber
                elif key in list_keys:
910 1 tkerber
                    list = []
911 1 tkerber
                    if key in ('dipol', 'eint', 'ferwe', 'ropt', 'rwigs'):
912 1 tkerber
                        for a in data[2:]:
913 1 tkerber
                            list.append(float(a))
914 1 tkerber
                    elif key in ('iband', 'kpuse'):
915 1 tkerber
                        for a in data[2:]:
916 1 tkerber
                            list.append(int(a))
917 1 tkerber
                    self.list_params[key] = list
918 1 tkerber
                    if key == 'magmom':
919 1 tkerber
                        done = False
920 1 tkerber
                        for a in data[2:]:
921 1 tkerber
                            if '!' in a or done:
922 1 tkerber
                                done = True
923 1 tkerber
                            elif '*' in a:
924 1 tkerber
                                a = a.split('*')
925 1 tkerber
                                for b in range(int(a[0])):
926 1 tkerber
                                    list.append(float(a[1]))
927 1 tkerber
                            else:
928 1 tkerber
                                list.append(float(a))
929 1 tkerber
                        list = np.array(list)
930 1 tkerber
                        self.atoms.set_initial_magnetic_moments(list[self.resort])
931 1 tkerber
            except KeyError:
932 1 tkerber
                raise IOError('Keyword "%s" in INCAR is not known by calculator.' % key)
933 1 tkerber
            except IndexError:
934 1 tkerber
                raise IOError('Value missing for keyword "%s".' % key)
935 1 tkerber
936 1 tkerber
    def read_outcar(self):
937 1 tkerber
        # Spin polarized calculation?
938 1 tkerber
        file = open('OUTCAR', 'r')
939 1 tkerber
        lines = file.readlines()
940 1 tkerber
        file.close()
941 1 tkerber
        for line in lines:
942 1 tkerber
            if line.rfind('ISPIN') > -1:
943 1 tkerber
                if int(line.split()[2])==2:
944 1 tkerber
                    self.spinpol = True
945 1 tkerber
                else:
946 1 tkerber
                    self.spinpol = None
947 1 tkerber
        self.energy_free, self.energy_zero = self.read_energy()
948 1 tkerber
        self.forces = self.read_forces(self.atoms)
949 1 tkerber
        self.dipole = self.read_dipole()
950 1 tkerber
        self.fermi = self.read_fermi()
951 1 tkerber
        self.stress = self.read_stress()
952 1 tkerber
        self.nbands = self.read_nbands()
953 1 tkerber
        p=self.int_params
954 1 tkerber
        q=self.list_params
955 1 tkerber
        if self.spinpol:
956 1 tkerber
            self.magnetic_moment = self.read_magnetic_moment()
957 1 tkerber
            if p['lorbit']>=10 or (p['lorbit']!=None and q['rwigs']):
958 1 tkerber
                self.magnetic_moments = self.read_magnetic_moments(self.atoms)
959 1 tkerber
        self.set(nbands=self.nbands)
960 1 tkerber
961 1 tkerber
    def read_kpoints(self, filename='KPOINTS'):
962 1 tkerber
        file = open(filename, 'r')
963 1 tkerber
        lines = file.readlines()
964 1 tkerber
        file.close()
965 1 tkerber
        type = lines[2].split()[0].lower()[0]
966 1 tkerber
        if type in ['g', 'm']:
967 1 tkerber
            if type=='g':
968 1 tkerber
                self.set(gamma=True)
969 1 tkerber
            kpts = np.array([int(lines[3].split()[i]) for i in range(3)])
970 1 tkerber
            self.set(kpts=kpts)
971 1 tkerber
        elif type in ['c', 'k']:
972 1 tkerber
            raise NotImplementedError('Only Monkhorst-Pack and gamma centered grid supported for restart.')
973 1 tkerber
        else:
974 1 tkerber
            raise NotImplementedError('Only Monkhorst-Pack and gamma centered grid supported for restart.')
975 1 tkerber
976 1 tkerber
    def read_potcar(self):
977 1 tkerber
        """ Method that reads the Exchange Correlation functional from POTCAR file.
978 1 tkerber
        """
979 1 tkerber
        file = open('POTCAR', 'r')
980 1 tkerber
        lines = file.readlines()
981 1 tkerber
        file.close()
982 1 tkerber
983 1 tkerber
        # Search for key 'LEXCH' in POTCAR
984 1 tkerber
        xc_flag = None
985 1 tkerber
        for line in lines:
986 1 tkerber
            key = line.split()[0].upper()
987 1 tkerber
            if key == 'LEXCH':
988 1 tkerber
                xc_flag = line.split()[-1].upper()
989 1 tkerber
                break
990 1 tkerber
991 1 tkerber
        if xc_flag is None:
992 1 tkerber
            raise ValueError('LEXCH flag not found in POTCAR file.')
993 1 tkerber
994 1 tkerber
        # Values of parameter LEXCH and corresponding XC-functional
995 1 tkerber
        xc_dict = {'PE':'PBE', '91':'PW91', 'CA':'LDA'}
996 1 tkerber
997 1 tkerber
        if xc_flag not in xc_dict.keys():
998 1 tkerber
            raise ValueError(
999 1 tkerber
                'Unknown xc-functional flag found in POTCAR, LEXCH=%s' % xc_flag)
1000 1 tkerber
1001 1 tkerber
        self.input_params['xc'] = xc_dict[xc_flag]
1002 1 tkerber
1003 1 tkerber
1004 1 tkerber
class VaspChargeDensity(object):
1005 1 tkerber
    """Class for representing VASP charge density"""
1006 1 tkerber
1007 1 tkerber
    def __init__(self, filename='CHG'):
1008 1 tkerber
        # Instance variables
1009 1 tkerber
        self.atoms = []   # List of Atoms objects
1010 1 tkerber
        self.chg = []     # Charge density
1011 1 tkerber
        self.chgdiff = [] # Charge density difference, if spin polarized
1012 1 tkerber
        self.aug = ''     # Augmentation charges, not parsed just a big string
1013 1 tkerber
        self.augdiff = '' # Augmentation charge differece, is spin polarized
1014 1 tkerber
1015 1 tkerber
        # Note that the augmentation charge is not a list, since they
1016 1 tkerber
        # are needed only for CHGCAR files which store only a single
1017 1 tkerber
        # image.
1018 1 tkerber
        if filename != None:
1019 1 tkerber
            self.read(filename)
1020 1 tkerber
1021 1 tkerber
    def is_spin_polarized(self):
1022 1 tkerber
        if len(self.chgdiff) > 0:
1023 1 tkerber
            return True
1024 1 tkerber
        return False
1025 1 tkerber
1026 1 tkerber
    def _read_chg(self, fobj, chg, volume):
1027 1 tkerber
        """Read charge from file object
1028 1 tkerber

1029 1 tkerber
        Utility method for reading the actual charge density (or
1030 1 tkerber
        charge density difference) from a file object. On input, the
1031 1 tkerber
        file object must be at the beginning of the charge block, on
1032 1 tkerber
        output the file position will be left at the end of the
1033 1 tkerber
        block. The chg array must be of the correct dimensions.
1034 1 tkerber

1035 1 tkerber
        """
1036 1 tkerber
        # VASP writes charge density as
1037 1 tkerber
        # WRITE(IU,FORM) (((C(NX,NY,NZ),NX=1,NGXC),NY=1,NGYZ),NZ=1,NGZC)
1038 1 tkerber
        # Fortran nested implied do loops; innermost index fastest
1039 1 tkerber
        # First, just read it in
1040 1 tkerber
        for zz in range(chg.shape[2]):
1041 1 tkerber
            for yy in range(chg.shape[1]):
1042 1 tkerber
                chg[:, yy, zz] = np.fromfile(fobj, count = chg.shape[0],
1043 1 tkerber
                                             sep=' ')
1044 1 tkerber
        chg /= volume
1045 1 tkerber
1046 1 tkerber
    def read(self, filename='CHG'):
1047 1 tkerber
        """Read CHG or CHGCAR file.
1048 1 tkerber

1049 1 tkerber
        If CHG contains charge density from multiple steps all the
1050 1 tkerber
        steps are read and stored in the object. By default VASP
1051 1 tkerber
        writes out the charge density every 10 steps.
1052 1 tkerber

1053 1 tkerber
        chgdiff is the difference between the spin up charge density
1054 1 tkerber
        and the spin down charge density and is thus only read for a
1055 1 tkerber
        spin-polarized calculation.
1056 1 tkerber

1057 1 tkerber
        aug is the PAW augmentation charges found in CHGCAR. These are
1058 1 tkerber
        not parsed, they are just stored as a string so that they can
1059 1 tkerber
        be written again to a CHGCAR format file.
1060 1 tkerber

1061 1 tkerber
        """
1062 1 tkerber
        import ase.io.vasp as aiv
1063 1 tkerber
        f = open(filename)
1064 1 tkerber
        self.atoms = []
1065 1 tkerber
        self.chg = []
1066 1 tkerber
        self.chgdiff = []
1067 1 tkerber
        self.aug = ''
1068 1 tkerber
        self.augdiff = ''
1069 1 tkerber
        while True:
1070 1 tkerber
            try:
1071 1 tkerber
                atoms = aiv.read_vasp(f)
1072 1 tkerber
            except ValueError, e:
1073 1 tkerber
                # Probably an empty line, or we tried to read the
1074 1 tkerber
                # augmentation occupancies in CHGCAR
1075 1 tkerber
                break
1076 1 tkerber
            f.readline()
1077 1 tkerber
            ngr = f.readline().split()
1078 1 tkerber
            ng = (int(ngr[0]), int(ngr[1]), int(ngr[2]))
1079 1 tkerber
            chg = np.empty(ng)
1080 1 tkerber
            self._read_chg(f, chg, atoms.get_volume())
1081 1 tkerber
            self.chg.append(chg)
1082 1 tkerber
            self.atoms.append(atoms)
1083 1 tkerber
            # Check if the file has a spin-polarized charge density part, and
1084 1 tkerber
            # if so, read it in.
1085 1 tkerber
            fl = f.tell()
1086 1 tkerber
            # First check if the file has an augmentation charge part (CHGCAR file.)
1087 1 tkerber
            line1 = f.readline()
1088 1 tkerber
            if line1=='':
1089 1 tkerber
                break
1090 1 tkerber
            elif line1.find('augmentation') != -1:
1091 1 tkerber
                augs = [line1]
1092 1 tkerber
                while True:
1093 1 tkerber
                    line2 = f.readline()
1094 1 tkerber
                    if line2.split() == ngr:
1095 1 tkerber
                        self.aug = ''.join(augs)
1096 1 tkerber
                        augs = []
1097 1 tkerber
                        chgdiff = np.empty(ng)
1098 1 tkerber
                        self._read_chg(f, chgdiff, atoms.get_volume())
1099 1 tkerber
                        self.chgdiff.append(chgdiff)
1100 1 tkerber
                    elif line2 == '':
1101 1 tkerber
                        break
1102 1 tkerber
                    else:
1103 1 tkerber
                        augs.append(line2)
1104 1 tkerber
                if len(self.aug) == 0:
1105 1 tkerber
                    self.aug = ''.join(augs)
1106 1 tkerber
                    augs = []
1107 1 tkerber
                else:
1108 1 tkerber
                    self.augdiff = ''.join(augs)
1109 1 tkerber
                    augs = []
1110 1 tkerber
            elif line1.split() == ngr:
1111 1 tkerber
                chgdiff = np.empty(ng)
1112 1 tkerber
                self._read_chg(f, chgdiff, atoms.get_volume())
1113 1 tkerber
                self.chgdiff.append(chgdiff)
1114 1 tkerber
            else:
1115 1 tkerber
                f.seek(fl)
1116 1 tkerber
        f.close()
1117 1 tkerber
1118 1 tkerber
    def _write_chg(self, fobj, chg, volume, format='chg'):
1119 1 tkerber
        """Write charge density
1120 1 tkerber

1121 1 tkerber
        Utility function similar to _read_chg but for writing.
1122 1 tkerber

1123 1 tkerber
        """
1124 1 tkerber
        # Make a 1D copy of chg, must take transpose to get ordering right
1125 1 tkerber
        chgtmp=chg.T.ravel()
1126 1 tkerber
        # Multiply by volume
1127 1 tkerber
        chgtmp=chgtmp*volume
1128 1 tkerber
        # Must be a tuple to pass to string conversion
1129 1 tkerber
        chgtmp=tuple(chgtmp)
1130 1 tkerber
        # CHG format - 10 columns
1131 1 tkerber
        if format.lower() == 'chg':
1132 1 tkerber
            # Write all but the last row
1133 1 tkerber
            for ii in range((len(chgtmp)-1)/10):
1134 1 tkerber
                fobj.write(' %#11.5G %#11.5G %#11.5G %#11.5G %#11.5G\
1135 1 tkerber
 %#11.5G %#11.5G %#11.5G %#11.5G %#11.5G\n' % chgtmp[ii*10:(ii+1)*10]
1136 1 tkerber
                           )
1137 1 tkerber
            # If the last row contains 10 values then write them without a newline
1138 1 tkerber
            if len(chgtmp)%10==0:
1139 1 tkerber
                fobj.write(' %#11.5G %#11.5G %#11.5G %#11.5G %#11.5G\
1140 1 tkerber
 %#11.5G %#11.5G %#11.5G %#11.5G %#11.5G' % chgtmp[len(chgtmp)-10:len(chgtmp)])
1141 1 tkerber
            # Otherwise write fewer columns without a newline
1142 1 tkerber
            else:
1143 1 tkerber
                for ii in range(len(chgtmp)%10):
1144 1 tkerber
                    fobj.write((' %#11.5G') % chgtmp[len(chgtmp)-len(chgtmp)%10+ii])
1145 1 tkerber
        # Other formats - 5 columns
1146 1 tkerber
        else:
1147 1 tkerber
            # Write all but the last row
1148 1 tkerber
            for ii in range((len(chgtmp)-1)/5):
1149 1 tkerber
                fobj.write(' %17.10E %17.10E %17.10E %17.10E %17.10E\n' % chgtmp[ii*5:(ii+1)*5])
1150 1 tkerber
            # If the last row contains 5 values then write them without a newline
1151 1 tkerber
            if len(chgtmp)%5==0:
1152 1 tkerber
                fobj.write(' %17.10E %17.10E %17.10E %17.10E %17.10E' % chgtmp[len(chgtmp)-5:len(chgtmp)])
1153 1 tkerber
            # Otherwise write fewer columns without a newline
1154 1 tkerber
            else:
1155 1 tkerber
                for ii in range(len(chgtmp)%5):
1156 1 tkerber
                    fobj.write((' %17.10E') % chgtmp[len(chgtmp)-len(chgtmp)%5+ii])
1157 1 tkerber
        # Write a newline whatever format it is
1158 1 tkerber
        fobj.write('\n')
1159 1 tkerber
        # Clean up
1160 1 tkerber
        del chgtmp
1161 1 tkerber
1162 1 tkerber
    def write(self, filename='CHG', format=None):
1163 1 tkerber
        """Write VASP charge density in CHG format.
1164 1 tkerber

1165 1 tkerber
        filename: str
1166 1 tkerber
            Name of file to write to.
1167 1 tkerber
        format: str
1168 1 tkerber
            String specifying whether to write in CHGCAR or CHG
1169 1 tkerber
            format.
1170 1 tkerber

1171 1 tkerber
        """
1172 1 tkerber
        import ase.io.vasp as aiv
1173 1 tkerber
        if format == None:
1174 1 tkerber
            if filename.lower().find('chgcar') != -1:
1175 1 tkerber
                format = 'chgcar'
1176 1 tkerber
            elif filename.lower().find('chg') != -1:
1177 1 tkerber
                format = 'chg'
1178 1 tkerber
            elif len(self.chg) == 1:
1179 1 tkerber
                format = 'chgcar'
1180 1 tkerber
            else:
1181 1 tkerber
                format = 'chg'
1182 1 tkerber
        f = open(filename, 'w')
1183 1 tkerber
        for ii, chg in enumerate(self.chg):
1184 1 tkerber
            if format == 'chgcar' and ii != len(self.chg) - 1:
1185 1 tkerber
                continue # Write only the last image for CHGCAR
1186 1 tkerber
            aiv.write_vasp(f, self.atoms[ii], direct=True)
1187 1 tkerber
            f.write('\n')
1188 1 tkerber
            for dim in chg.shape:
1189 1 tkerber
                f.write(' %4i' % dim)
1190 1 tkerber
            f.write('\n')
1191 1 tkerber
            vol = self.atoms[ii].get_volume()
1192 1 tkerber
            self._write_chg(f, chg, vol, format)
1193 1 tkerber
            if format == 'chgcar':
1194 1 tkerber
                f.write(self.aug)
1195 1 tkerber
            if self.is_spin_polarized():
1196 1 tkerber
                if format == 'chg':
1197 1 tkerber
                    f.write('\n')
1198 1 tkerber
                for dim in chg.shape:
1199 1 tkerber
                    f.write(' %4i' % dim)
1200 1 tkerber
                self._write_chg(f, self.chgdiff[ii], vol, format)
1201 1 tkerber
                if format == 'chgcar':
1202 1 tkerber
                    f.write('\n')
1203 1 tkerber
                    f.write(self.augdiff)
1204 1 tkerber
            if format == 'chg' and len(self.chg) > 1:
1205 1 tkerber
                f.write('\n')
1206 1 tkerber
        f.close()
1207 1 tkerber
1208 1 tkerber
1209 1 tkerber
class VaspDos(object):
1210 1 tkerber
    """Class for representing density-of-states produced by VASP
1211 1 tkerber

1212 1 tkerber
    The energies are in property self.energy
1213 1 tkerber

1214 1 tkerber
    Site-projected DOS is accesible via the self.site_dos method.
1215 1 tkerber

1216 1 tkerber
    Total and integrated DOS is accessible as numpy.ndarray's in the
1217 1 tkerber
    properties self.dos and self.integrated_dos. If the calculation is
1218 1 tkerber
    spin polarized, the arrays will be of shape (2, NDOS), else (1,
1219 1 tkerber
    NDOS).
1220 1 tkerber

1221 1 tkerber
    The self.efermi property contains the currently set Fermi
1222 1 tkerber
    level. Changing this value shifts the energies.
1223 1 tkerber

1224 1 tkerber
    """
1225 1 tkerber
1226 1 tkerber
    def __init__(self, doscar='DOSCAR', efermi=0.0):
1227 1 tkerber
        """Initialize"""
1228 1 tkerber
        self._efermi = 0.0
1229 1 tkerber
        self.read_doscar(doscar)
1230 1 tkerber
        self.efermi = efermi
1231 1 tkerber
1232 1 tkerber
    def _set_efermi(self, efermi):
1233 1 tkerber
        """Set the Fermi level."""
1234 1 tkerber
        ef = efermi - self._efermi
1235 1 tkerber
        self._efermi = efermi
1236 1 tkerber
        self._total_dos[0, :] = self._total_dos[0, :] - ef
1237 1 tkerber
        try:
1238 1 tkerber
            self._site_dos[:, 0, :] = self._site_dos[:, 0, :] - ef
1239 1 tkerber
        except IndexError:
1240 1 tkerber
            pass
1241 1 tkerber
1242 1 tkerber
    def _get_efermi(self):
1243 1 tkerber
        return self._efermi
1244 1 tkerber
1245 1 tkerber
    efermi = property(_get_efermi, _set_efermi, None, "Fermi energy.")
1246 1 tkerber
1247 1 tkerber
    def _get_energy(self):
1248 1 tkerber
        """Return the array with the energies."""
1249 1 tkerber
        return self._total_dos[0, :]
1250 1 tkerber
    energy = property(_get_energy, None, None, "Array of energies")
1251 1 tkerber
1252 1 tkerber
    def site_dos(self, atom, orbital):
1253 1 tkerber
        """Return an NDOSx1 array with dos for the chosen atom and orbital.
1254 1 tkerber

1255 1 tkerber
        atom: int
1256 1 tkerber
            Atom index
1257 1 tkerber
        orbital: int or str
1258 1 tkerber
            Which orbital to plot
1259 1 tkerber

1260 1 tkerber
        If the orbital is given as an integer:
1261 1 tkerber
        If spin-unpolarized calculation, no phase factors:
1262 1 tkerber
        s = 0, p = 1, d = 2
1263 1 tkerber
        Spin-polarized, no phase factors:
1264 1 tkerber
        s-up = 0, s-down = 1, p-up = 2, p-down = 3, d-up = 4, d-down = 5
1265 1 tkerber
        If phase factors have been calculated, orbitals are
1266 1 tkerber
        s, py, pz, px, dxy, dyz, dz2, dxz, dx2
1267 1 tkerber
        double in the above fashion if spin polarized.
1268 1 tkerber

1269 1 tkerber
        """
1270 1 tkerber
        # Integer indexing for orbitals starts from 1 in the _site_dos array
1271 1 tkerber
        # since the 0th column contains the energies
1272 1 tkerber
        if isinstance(orbital, int):
1273 1 tkerber
            return self._site_dos[atom, orbital + 1, :]
1274 1 tkerber
        n = self._site_dos.shape[1]
1275 1 tkerber
        if n == 4:
1276 1 tkerber
            norb = {'s':1, 'p':2, 'd':3}
1277 1 tkerber
        elif n == 7:
1278 1 tkerber
            norb = {'s+':1, 's-up':1, 's-':2, 's-down':2,
1279 1 tkerber
                    'p+':3, 'p-up':3, 'p-':4, 'p-down':4,
1280 1 tkerber
                    'd+':5, 'd-up':5, 'd-':6, 'd-down':6}
1281 1 tkerber
        elif n == 10:
1282 1 tkerber
            norb = {'s':1, 'py':2, 'pz':3, 'px':4,
1283 1 tkerber
                    'dxy':5, 'dyz':6, 'dz2':7, 'dxz':8,
1284 1 tkerber
                    'dx2':9}
1285 1 tkerber
        elif n == 19:
1286 1 tkerber
            norb = {'s+':1, 's-up':1, 's-':2, 's-down':2,
1287 1 tkerber
                    'py+':3, 'py-up':3, 'py-':4, 'py-down':4,
1288 1 tkerber
                    'pz+':5, 'pz-up':5, 'pz-':6, 'pz-down':6,
1289 1 tkerber
                    'px+':7, 'px-up':7, 'px-':8, 'px-down':8,
1290 1 tkerber
                    'dxy+':9, 'dxy-up':9, 'dxy-':10, 'dxy-down':10,
1291 1 tkerber
                    'dyz+':11, 'dyz-up':11, 'dyz-':12, 'dyz-down':12,
1292 1 tkerber
                    'dz2+':13, 'dz2-up':13, 'dz2-':14, 'dz2-down':14,
1293 1 tkerber
                    'dxz+':15, 'dxz-up':15, 'dxz-':16, 'dxz-down':16,
1294 1 tkerber
                    'dx2+':17, 'dx2-up':17, 'dx2-':18, 'dx2-down':18}
1295 1 tkerber
        return self._site_dos[atom, norb[orbital.lower()], :]
1296 1 tkerber
1297 1 tkerber
    def _get_dos(self):
1298 1 tkerber
        if self._total_dos.shape[0] == 3:
1299 1 tkerber
            return self._total_dos[1, :]
1300 1 tkerber
        elif self._total_dos.shape[0] == 5:
1301 1 tkerber
            return self._total_dos[1:3, :]
1302 1 tkerber
    dos = property(_get_dos, None, None, 'Average DOS in cell')
1303 1 tkerber
1304 1 tkerber
    def _get_integrated_dos(self):
1305 1 tkerber
        if self._total_dos.shape[0] == 3:
1306 1 tkerber
            return self._total_dos[2, :]
1307 1 tkerber
        elif self._total_dos.shape[0] == 5:
1308 1 tkerber
            return self._total_dos[3:5, :]
1309 1 tkerber
    integrated_dos = property(_get_integrated_dos, None, None,
1310 1 tkerber
                              'Integrated average DOS in cell')
1311 1 tkerber
1312 1 tkerber
    def read_doscar(self, fname="DOSCAR"):
1313 1 tkerber
        """Read a VASP DOSCAR file"""
1314 1 tkerber
        f = open(fname)
1315 1 tkerber
        natoms = int(f.readline().split()[0])
1316 1 tkerber
        [f.readline() for nn in range(4)]  # Skip next 4 lines.
1317 1 tkerber
        # First we have a block with total and total integrated DOS
1318 1 tkerber
        ndos = int(f.readline().split()[2])
1319 1 tkerber
        dos = []
1320 1 tkerber
        for nd in xrange(ndos):
1321 1 tkerber
            dos.append(np.array([float(x) for x in f.readline().split()]))
1322 1 tkerber
        self._total_dos = np.array(dos).T
1323 1 tkerber
        # Next we have one block per atom, if INCAR contains the stuff
1324 1 tkerber
        # necessary for generating site-projected DOS
1325 1 tkerber
        dos = []
1326 1 tkerber
        for na in xrange(natoms):
1327 1 tkerber
            line = f.readline()
1328 1 tkerber
            if line == '':
1329 1 tkerber
                # No site-projected DOS
1330 1 tkerber
                break
1331 1 tkerber
            ndos = int(line.split()[2])
1332 1 tkerber
            line = f.readline().split()
1333 1 tkerber
            cdos = np.empty((ndos, len(line)))
1334 1 tkerber
            cdos[0] = np.array(line)
1335 1 tkerber
            for nd in xrange(1, ndos):
1336 1 tkerber
                line = f.readline().split()
1337 1 tkerber
                cdos[nd] = np.array([float(x) for x in line])
1338 1 tkerber
            dos.append(cdos.T)
1339 1 tkerber
        self._site_dos = np.array(dos)
1340 1 tkerber
1341 1 tkerber
1342 1 tkerber
import pickle
1343 1 tkerber
1344 1 tkerber
class xdat2traj:
1345 1 tkerber
    def __init__(self, trajectory=None, atoms=None, poscar=None,
1346 1 tkerber
                 xdatcar=None, sort=None, calc=None):
1347 1 tkerber
        if not poscar:
1348 1 tkerber
            self.poscar = 'POSCAR'
1349 1 tkerber
        else:
1350 1 tkerber
            self.poscar = poscar
1351 1 tkerber
        if not atoms:
1352 1 tkerber
            self.atoms = ase.io.read(self.poscar, format='vasp')
1353 1 tkerber
        else:
1354 1 tkerber
            self.atoms = atoms
1355 1 tkerber
        if not xdatcar:
1356 1 tkerber
            self.xdatcar = 'XDATCAR'
1357 1 tkerber
        else:
1358 1 tkerber
            self.xdatcar = xdatcar
1359 1 tkerber
        if not trajectory:
1360 1 tkerber
            self.trajectory = 'out.traj'
1361 1 tkerber
        else:
1362 1 tkerber
            self.trajectory = trajectory
1363 1 tkerber
        if not calc:
1364 1 tkerber
            self.calc = Vasp()
1365 1 tkerber
        else:
1366 1 tkerber
            self.calc = calc
1367 1 tkerber
        if not sort:
1368 1 tkerber
            if not hasattr(self.calc, 'sort'):
1369 1 tkerber
                self.calc.sort = range(len(self.atoms))
1370 1 tkerber
        else:
1371 1 tkerber
            self.calc.sort = sort
1372 1 tkerber
        self.calc.resort = range(len(self.calc.sort))
1373 1 tkerber
        for n in range(len(self.calc.resort)):
1374 1 tkerber
            self.calc.resort[self.calc.sort[n]] = n
1375 1 tkerber
        self.out = ase.io.trajectory.PickleTrajectory(self.trajectory, mode='w')
1376 1 tkerber
        self.energies = self.calc.read_energy(all=True)[1]
1377 1 tkerber
        self.forces = self.calc.read_forces(self.atoms, all=True)
1378 1 tkerber
1379 1 tkerber
    def convert(self):
1380 1 tkerber
        lines = open(self.xdatcar).readlines()
1381 1 tkerber
        if len(lines[5].split())==0:
1382 1 tkerber
            del(lines[0:6])
1383 1 tkerber
        elif len(lines[4].split())==0:
1384 1 tkerber
            del(lines[0:5])
1385 1 tkerber
        step = 0
1386 1 tkerber
        iatom = 0
1387 1 tkerber
        scaled_pos = []
1388 1 tkerber
        for line in lines:
1389 1 tkerber
            if iatom == len(self.atoms):
1390 1 tkerber
                if step == 0:
1391 1 tkerber
                    self.out.write_header(self.atoms[self.calc.resort])
1392 1 tkerber
                scaled_pos = np.array(scaled_pos)
1393 1 tkerber
                self.atoms.set_scaled_positions(scaled_pos)
1394 1 tkerber
                d = {'positions': self.atoms.get_positions()[self.calc.resort],
1395 1 tkerber
                     'cell': self.atoms.get_cell(),
1396 1 tkerber
                     'momenta': None,
1397 1 tkerber
                     'energy': self.energies[step],
1398 1 tkerber
                     'forces': self.forces[step],
1399 1 tkerber
                     'stress': None}
1400 1 tkerber
                pickle.dump(d, self.out.fd, protocol=-1)
1401 1 tkerber
                scaled_pos = []
1402 1 tkerber
                iatom = 0
1403 1 tkerber
                step += 1
1404 1 tkerber
            else:
1405 1 tkerber
1406 1 tkerber
                iatom += 1
1407 1 tkerber
                scaled_pos.append([float(line.split()[n]) for n in range(3)])
1408 1 tkerber
1409 1 tkerber
        # Write also the last image
1410 1 tkerber
        # I'm sure there is also more clever fix...
1411 1 tkerber
        scaled_pos = np.array(scaled_pos)
1412 1 tkerber
        self.atoms.set_scaled_positions(scaled_pos)
1413 1 tkerber
        d = {'positions': self.atoms.get_positions()[self.calc.resort],
1414 1 tkerber
             'cell': self.atoms.get_cell(),
1415 1 tkerber
             'momenta': None,
1416 1 tkerber
             'energy': self.energies[step],
1417 1 tkerber
             'forces': self.forces[step],
1418 1 tkerber
             'stress': None}
1419 1 tkerber
        pickle.dump(d, self.out.fd, protocol=-1)
1420 1 tkerber
1421 1 tkerber
        self.out.fd.close()