Statistiques
| Révision :

root / ase / calculators / siesta.py @ 20

Historique | Voir | Annoter | Télécharger (24,21 ko)

1
"""This module defines an ASE interface to SIESTA.
2

3
http://www.uam.es/departamentos/ciencias/fismateriac/siesta
4
"""
5

    
6
import os
7
from os.path import join, isfile, islink
8
from cmath import exp
9
import array
10

    
11
import numpy as np
12

    
13
from ase.data import chemical_symbols
14
from ase.units import Rydberg, fs
15

    
16
class Siesta:
17
    """Class for doing SIESTA calculations.
18

19
    The default parameters are very close to those that the SIESTA
20
    Fortran code would use.  These are the exceptions::
21
    
22
      calc = Siesta(label='siesta', xc='LDA', pulay=5, mix=0.1)
23

24
    Use the set_fdf method to set extra FDF parameters::
25
    
26
      calc.set_fdf('PAO.EnergyShift', 0.01 * Rydberg)
27

28
    """
29
    def __init__(self, label='siesta', xc='LDA', kpts=None, nbands=None,
30
                 width=None, meshcutoff=None, charge=None,
31
                 pulay=5, mix=0.1, maxiter=120,
32
                 basis=None, ghosts=[],
33
                 write_fdf=True):
34
        """Construct SIESTA-calculator object.
35

36
        Parameters
37
        ==========
38
        label: str
39
            Prefix to use for filenames (label.fdf, label.txt, ...).
40
            Default is 'siesta'.
41
        xc: str
42
            Exchange-correlation functional.  Must be one of LDA, PBE,
43
            revPBE, RPBE.
44
        kpts: list of three int
45
            Monkhost-Pack sampling.
46
        nbands: int
47
            Number of bands.
48
        width: float
49
            Fermi-distribution width in eV.
50
        meshcutoff: float
51
            Cutoff energy in eV for grid.
52
        charge: float
53
            Total charge of the system.
54
        pulay: int
55
            Number of old densities to use for Pulay mixing.
56
        mix: float
57
            Mixing parameter between zero and one for density mixing.
58
        write_fdf: bool
59
            Use write_fdf=False to use your own fdf-file.
60

61
        Examples
62
        ========
63
        Use default values:
64
        
65
        >>> h = Atoms('H', calculator=Siesta())
66
        >>> h.center(vacuum=3.0)
67
        >>> e = h.get_potential_energy()
68
        
69
        """
70
        
71
        self.label = label#################### != out
72
        self.xc = xc
73
        self.kpts = kpts
74
        self.nbands = nbands
75
        self.width = width
76
        self.meshcutoff = meshcutoff
77
        self.charge = charge
78
        self.pulay = pulay
79
        self.mix = mix
80
        self.maxiter = maxiter
81
        self.basis = basis
82
        self.ghosts = ghosts
83
        self.write_fdf_file = write_fdf
84
    
85
        self.converged = False
86
        self.fdf = {}
87
        self.e_fermi = None
88
        
89
    def update(self, atoms):
90
        if (not self.converged or
91
            len(self.numbers) != len(atoms) or
92
            (self.numbers != atoms.get_atomic_numbers()).any()):
93
            self.initialize(atoms)
94
            self.calculate(atoms)
95
        elif ((self.positions != atoms.get_positions()).any() or
96
              (self.pbc != atoms.get_pbc()).any() or
97
              (self.cell != atoms.get_cell()).any()):
98
            self.calculate(atoms)
99

    
100
    def initialize(self, atoms):
101
        self.numbers = atoms.get_atomic_numbers().copy()
102
        self.species = []
103
        for a, Z in enumerate(self.numbers):
104
            if a in self.ghosts:
105
                Z = -Z
106
            if Z not in self.species:
107
                self.species.append(Z)
108

    
109
        if 'SIESTA_PP_PATH' in os.environ:
110
            pppaths = os.environ['SIESTA_PP_PATH'].split(':')
111
        else:
112
            pppaths = []
113

    
114
        for Z in self.species:
115
            symbol = chemical_symbols[abs(Z)]
116
            name = symbol + '.vps'
117
            name1 = symbol + '.psf'
118
            found = False
119
            for path in pppaths:
120
                filename = join(path, name)
121
                filename1 = join(path, name1)
122
                if isfile(filename) or islink(filename):
123
                    found = True
124
                    if path != '.':
125
                        if islink(name) or isfile(name):
126
                            os.remove(name)
127
                        os.symlink(filename, name)
128

    
129
                elif isfile(filename1) or islink(filename1):
130
                    found = True
131
                    if path != '.':
132
                        if islink(name1) or isfile(name1):
133
                            os.remove(name1)
134
                        os.symlink(filename1, name1)
135
            if not found:
136
                raise RuntimeError('No pseudopotential for %s!' % symbol)
137

    
138
        self.converged = False
139
        
140
    def get_potential_energy(self, atoms, force_consistent=False):
141
        self.update(atoms)
142

    
143
        if force_consistent:
144
            return self.efree
145
        else:
146
            # Energy extrapolated to zero Kelvin:
147
            return  (self.etotal + self.efree) / 2
148

    
149
    def get_forces(self, atoms):
150
        self.update(atoms)
151
        return self.forces.copy()
152
    
153
    def get_stress(self, atoms):
154
        self.update(atoms)
155
        return self.stress.copy()
156

    
157
    def get_dipole_moment(self, atoms):
158
        """Returns total dipole moment of the system."""
159
        self.update(atoms)
160
        return self.dipole
161

    
162
    def read_dipole(self):
163
        dipolemoment = np.zeros([1, 3])
164
        for line in open(self.label + '.txt', 'r'):
165
            if line.rfind('Electric dipole (Debye)') > -1:
166
                dipolemoment = np.array([float(f) for f in line.split()[5:8]])
167
        #debye to e*Ang (the units of VASP)
168
        dipolemoment = dipolemoment*0.2081943482534
169
        return dipolemoment
170

    
171
    def calculate(self, atoms):
172
        self.positions = atoms.get_positions().copy()
173
        self.cell = atoms.get_cell().copy()
174
        self.pbc = atoms.get_pbc().copy()
175

    
176
        if self.write_fdf_file:
177
            self.write_fdf(atoms)
178

    
179
        siesta = os.environ['SIESTA_SCRIPT']
180
        locals = {'label': self.label}
181
        execfile(siesta, {}, locals)
182
        exitcode = locals['exitcode']
183
        if exitcode != 0:
184
            raise RuntimeError(('Siesta exited with exit code: %d.  ' +
185
                                'Check %s.txt for more information.') %
186
                               (exitcode, self.label))
187
        
188
        self.dipole = self.read_dipole()
189
        self.read()
190

    
191
        self.converged = True
192
        
193
    def set_fdf(self, key, value):
194
        """Set FDF parameter."""
195
        self.fdf[key] = value
196
                
197
    def write_fdf(self, atoms):
198
        """Write input parameters to fdf-file."""
199
        fh = open(self.label + '.fdf', 'w')
200

    
201
        fdf = {
202
            'SystemLabel': self.label,
203
            'AtomicCoordinatesFormat': 'Ang',
204
            'LatticeConstant': 1.0,
205
            'NumberOfAtoms': len(atoms),
206
            'MeshCutoff': self.meshcutoff,
207
            'NetCharge': self.charge,
208
            'ElectronicTemperature': self.width,
209
            'NumberOfEigenStates': self.nbands,
210
            'DM.UseSaveDM': self.converged,
211
            'PAO.BasisSize': self.basis,
212
            'SolutionMethod': 'diagon',
213
            'DM.NumberPulay': self.pulay,
214
            'DM.MixingWeight': self.mix,
215
            'MaxSCFIterations': self.maxiter
216
            }
217
        
218
        if self.xc != 'LDA':
219
            fdf['xc.functional'] = 'GGA'
220
            fdf['xc.authors'] = self.xc
221

    
222
        magmoms = atoms.get_initial_magnetic_moments()
223
        if magmoms.any():
224
            fdf['SpinPolarized'] = True
225
            fh.write('%block InitSpin\n')
226
            for n, M in enumerate(magmoms):
227
                if M != 0:
228
                    fh.write('%d %.14f\n' % (n + 1, M))
229
            fh.write('%endblock InitSpin\n')
230
        
231
        fdf['Number_of_species'] = len(self.species)
232

    
233
        fdf.update(self.fdf)
234

    
235
        for key, value in fdf.items():
236
            if value is None:
237
                continue
238

    
239
            if isinstance(value, list):
240
                fh.write('%%block %s\n' % key)
241
                for line in value:
242
                    fh.write(line + '\n')
243
                fh.write('%%endblock %s\n' % key)
244
            else:
245
                unit = keys_with_units.get(fdfify(key))
246
                if unit is None:
247
                    fh.write('%s %s\n' % (key, value))
248
                else:
249
                    if 'fs**2' in unit:
250
                        value /= fs**2
251
                    elif 'fs' in unit:
252
                        value /= fs
253
                    fh.write('%s %f %s\n' % (key, value, unit))
254

    
255
        fh.write('%block LatticeVectors\n')
256
        for v in self.cell:
257
            fh.write('%.14f %.14f %.14f\n' % tuple(v))
258
        fh.write('%endblock LatticeVectors\n')
259

    
260
        fh.write('%block Chemical_Species_label\n')
261
        for n, Z in enumerate(self.species):
262
            fh.write('%d %s %s\n' % (n + 1, Z, chemical_symbols[abs(Z)]))
263
        fh.write('%endblock Chemical_Species_label\n')
264

    
265
        fh.write('%block AtomicCoordinatesAndAtomicSpecies\n')
266
        a = 0
267
        for pos, Z in zip(self.positions, self.numbers):
268
            if a in self.ghosts:
269
                Z = -Z
270
            a += 1
271
            fh.write('%.14f %.14f %.14f' %  tuple(pos))
272
            fh.write(' %d\n' % (self.species.index(Z) + 1))
273
        fh.write('%endblock AtomicCoordinatesAndAtomicSpecies\n')
274

    
275
        if self.kpts is not None:
276
            fh.write('%block kgrid_Monkhorst_Pack\n')
277
            for i in range(3):
278
                for j in range(3):
279
                    if i == j:
280
                        fh.write('%d ' % self.kpts[i])
281
                    else:
282
                        fh.write('0 ')
283
                fh.write('%.1f\n' % (((self.kpts[i] + 1) % 2) * 0.5))
284
            fh.write('%endblock kgrid_Monkhorst_Pack\n')
285
        
286
        fh.close()
287
        
288
    def read(self):
289
        """Read results from SIESTA's text-output file."""
290
        text = open(self.label + '.txt', 'r').read().lower()
291
        assert 'error' not in text
292
        lines = iter(text.split('\n'))
293

    
294
        # Get the number of grid points used:
295
        for line in lines:
296
            if line.startswith('initmesh: mesh ='):
297
                self.grid = [int(word) for word in line.split()[3:8:2]]
298
                break
299

    
300
        # Stress (fixed so it's compatible with a MD run from siesta):
301
        for line in lines:
302
            if line.startswith('siesta: stress tensor '):
303
                self.stress = np.empty((3, 3))
304
                for i in range(3):
305
                    tmp = lines.next().split()
306
                    if len(tmp) == 4:
307
                        self.stress[i] = [float(word) for word in tmp[1:]]
308
                    else:
309
                        self.stress[i] = [float(word) for word in tmp]
310
                break
311
        else:
312
            raise RuntimeError
313

    
314
        text = open(self.label + '.txt', 'r').read().lower()
315
        lines = iter(text.split('\n'))
316
        # Energy (again a fix to make it compatible with a MD run from siesta):
317
        counter = 0
318
        for line in lines:
319
            if line.startswith('siesta: etot    =') and counter == 0:
320
                counter += 1
321
            elif line.startswith('siesta: etot    ='):
322
                self.etotal = float(line.split()[-1])
323
                self.efree = float(lines.next().split()[-1])
324
                break
325
        else:
326
            raise RuntimeError
327

    
328
        # Forces (changed so forces smaller than -999eV/A can be fetched):
329
        lines = open(self.label + '.FA', 'r').readlines()
330
        assert int(lines[0]) == len(self.numbers)
331
        assert len(lines) == len(self.numbers) + 1
332
        lines = lines[1:]
333
        self.forces = np.zeros((len(lines), 3))
334
        for i in range(len(lines)):
335
            self.forces[i, 0] = float(lines[i][6:18].strip())
336
            self.forces[i, 1] = float(lines[i][18:30].strip())
337
            self.forces[i, 2] = float(lines[i][30:42].strip())
338

    
339
    def read_eig(self):
340
        if self.e_fermi is not None:
341
            return
342

    
343
        assert os.access(self.label + '.EIG', os.F_OK)
344
        assert os.access(self.label + '.KP', os.F_OK)
345

    
346
        # Read k point weights
347
        text = open(self.label + '.KP', 'r').read()
348
        lines = text.split('\n')
349
        n_kpts = int(lines[0].strip())
350
        self.weights = np.zeros((n_kpts,))
351
        for i in range(n_kpts):
352
            l = lines[i + 1].split()
353
            self.weights[i] = float(l[4])
354

    
355
        # Read eigenvalues and fermi-level
356
        text = open(self.label+'.EIG','r').read()
357
        lines = text.split('\n')
358
        self.e_fermi = float(lines[0].split()[0])
359
        tmp = lines[1].split()
360
        self.n_bands = int(tmp[0])
361
        n_spin_bands = int(tmp[1])
362
        self.spin_pol = n_spin_bands == 2
363
        lines = lines[2:-1]
364
        lines_per_kpt = (self.n_bands * n_spin_bands / 10 +
365
                         int((self.n_bands * n_spin_bands) % 10 != 0))
366
        self.eig = dict()
367
        for i in range(len(self.weights)):
368
            tmp = lines[i * lines_per_kpt:(i + 1) * lines_per_kpt]
369
            v = [float(v) for v in tmp[0].split()[1:]]
370
            for l in tmp[1:]:
371
                v.extend([float(t) for t in l.split()])
372
            if self.spin_pol:
373
                self.eig[(i, 1)] = np.array(v[0:self.n_bands])
374
                self.eig[(i, -1)] = np.array(v[self.n_bands:])
375
            else:
376
                self.eig[(i, 1)] = np.array(v)
377
        
378
    def get_k_point_weights(self):
379
        self.read_eig()
380
        return self.weights
381

    
382
    def get_fermi_level(self):
383
        self.read_eig()
384
        return self.e_fermi
385

    
386
    def get_eigenvalues(self, kpt=0, spin=1):
387
        self.read_eig()
388
        return self.eig[(kpt, spin)]
389

    
390
    def get_number_of_spins(self):
391
        self.read_eig()
392
        if self.spin_pol:
393
            return 2
394
        else:
395
            return 1
396

    
397
    def read_hs(self, filename, is_gamma_only=False, magnus=False):
398
        """Read the Hamiltonian and overlap matrix from a Siesta 
399
           calculation in sparse format. 
400

401
        Parameters
402
        ==========
403
        filename: str
404
            The filename should be on the form jobname.HS  
405
        is_gamma_only: {False, True), optional
406
            Is it a gamma point calculation?
407
        magnus: bool
408
            The fileformat was changed by Magnus in Siesta at some
409
            point around version 2.xxx. 
410
            Use mangus=False, to use the old file format.
411

412
        Note
413
        ====
414
        Data read in is put in self._dat.
415

416
        Examples
417
        ========
418
            >>> calc = Siesta()
419
            >>> calc.read_hs('jobname.HS')
420
            >>> print calc._dat.fermi_level
421
            >>> print 'Number of orbitals: %i' % calc._dat.nuotot 
422
        """
423
        assert not magnus, 'Not implemented; changes by Magnus to file io' 
424
        assert not is_gamma_only, 'Not implemented. Only works for k-points.'
425
        class Dummy:
426
            pass
427
        self._dat = dat = Dummy()
428
        # Try to read supercell and atom data from a jobname.XV file
429
        filename_xv = filename[:-2] + 'XV'
430
        #assert isfile(filename_xv), 'Missing jobname.XV file'
431
        if isfile(filename_xv):
432
            print 'Reading supercell and atom data from ' + filename_xv
433
            fd = open(filename_xv, 'r')
434
            dat.cell = np.zeros((3, 3)) # Supercell
435
            for a_vec in dat.cell:
436
                a_vec[:] = np.array(fd.readline().split()[:3], float)
437
            dat.rcell = 2 * np.pi * np.linalg.inv(dat.cell.T)
438
            dat.natoms = int(fd.readline().split()[0])
439
            dat.symbols = []
440
            dat.pos_ac = np.zeros((dat.natoms, 3))
441
            for a in range(dat.natoms):
442
                line = fd.readline().split()
443
                dat.symbols.append(chemical_symbols[int(line[1])])
444
                dat.pos_ac[a, :] = [float(line[i]) for i in range(2, 2 + 3)]
445
        # Read in the jobname.HS file
446
        fileobj = file(filename, 'rb')
447
        fileobj.seek(0)
448
        dat.fermi_level = float(open(filename[:-3] + '.EIG', 'r').readline())
449
        dat.is_gammay_only = is_gamma_only 
450
        dat.nuotot, dat.ns, dat.mnh = getrecord(fileobj, 'l')
451
        nuotot, ns, mnh = dat.nuotot, dat.ns, dat.mnh
452
        print 'Number of orbitals found: %i' % nuotot
453
        dat.numh = numh = np.array([getrecord(fileobj, 'l')
454
                                    for i in range(nuotot)], 'l')
455
        dat.maxval = max(numh)
456
        dat.listhptr = listhptr = np.zeros(nuotot, 'l')
457
        listhptr[0] = 0
458
        for oi in xrange(1, nuotot):
459
            listhptr[oi] = listhptr[oi - 1] + numh[oi - 1]
460
        dat.listh = listh = np.zeros(mnh, 'l')
461
        
462
        print 'Reading sparse info'
463
        for oi in xrange(nuotot):
464
            for mi in xrange(numh[oi]):
465
                listh[listhptr[oi] + mi] = getrecord(fileobj, 'l')
466

    
467
        dat.nuotot_sc = max(listh)
468
        dat.h_sparse = h_sparse = np.zeros((mnh, ns), float)
469
        dat.s_sparse = s_sparse = np.zeros(mnh, float)
470
        print 'Reading H'
471
        for si in xrange(ns):
472
            for oi in xrange(nuotot):
473
                for mi in xrange(numh[oi]):
474
                    h_sparse[listhptr[oi] + mi, si] = getrecord(fileobj, 'd')
475
        print 'Reading S'
476
        for oi in xrange(nuotot):
477
            for mi in xrange(numh[oi]):
478
                s_sparse[listhptr[oi] + mi] = getrecord(fileobj, 'd')
479

    
480
        dat.qtot, dat.temperature = getrecord(fileobj, 'd')
481
        if not is_gamma_only:
482
            print 'Reading X'
483
            dat.xij_sparse = xij_sparse = np.zeros([3, mnh], float)
484
            for oi in xrange(nuotot):
485
                for mi in xrange(numh[oi]):
486
                    xij_sparse[:, listhptr[oi] + mi] = getrecord(fileobj, 'd')
487
        fileobj.close()
488

    
489
    def get_hs(self, kpt=(0, 0, 0), spin=0, remove_pbc=None, kpt_scaled=True):
490
        """Hamiltonian and overlap matrices for an arbitrary k-point.
491
       
492
        The default values corresponds to the Gamma point for 
493
        spin 0 and periodic boundary conditions.
494

495
        Parameters
496
        ==========
497
        kpt : {(0, 0, 0), (3,) array_like}, optional
498
            k-point in scaled or absolute coordinates.
499
            For the latter the units should be Bohr^-1.
500
        spin : {0, 1}, optional
501
            Spin index 
502
        remove_pbc : {None, ({'x', 'y', 'z'}, basis)}, optional
503
            Use remove_pbc to truncate h and s along a cartesian
504
            axis. 
505
        basis: {str, dict}
506
            The basis specification as either a string or a dictionary.
507
        kpt_scaled : {True, bool}, optional
508
            Use kpt_scaled=False if `kpt` is in absolute units (Bohr^-1).
509

510
        Note
511
        ====
512
        read_hs should be called before get_hs gets called.
513

514
        Examples
515
        ========
516
        >>> calc = Siesta()
517
        >>> calc.read_hs('jobname.HS')
518
        >>> h, s = calc.get_hs((0.0, 0.375, 0.375))
519
        >>> h -= s * calc._dat.fermi_level # fermi level is now at 0.0
520
        >>> basis = 'szp'
521
        >>> h, s = calc.get_hs((0.0, 0.375, 0.375), remove_pbc=('x', basis))
522
        >>> basis = {'Au:'sz}', 'C':'dzp', None:'szp'}
523
        >>> h, s = calc.get_hs((0.0, 0.375, 0.375), remove_pbc=('x', basis))
524

525
        """
526
        if not hasattr(self, '_dat'):# XXX Crude check if data is avail.
527
            print 'Please read in data first by calling the method read_hs.'
528
            return None, None
529
        dot = np.dot
530
        dat = self._dat            
531
        kpt_c = np.array(kpt, float)
532
        if kpt_scaled:
533
            kpt_c = dot(kpt_c, dat.rcell)
534

    
535
        h_MM = np.zeros((dat.nuotot, dat.nuotot), complex)
536
        s_MM = np.zeros((dat.nuotot, dat.nuotot), complex)
537
        h_sparse, s_sparse = dat.h_sparse, dat.s_sparse
538
        x_sparse = dat.xij_sparse
539
        numh, listhptr, listh = dat.numh, dat.listhptr, dat.listh
540
        indxuo = np.mod(np.arange(dat.nuotot_sc), dat.nuotot)
541
        
542
        for iuo in xrange(dat.nuotot):
543
            for j in range(numh[iuo]):
544
                ind =  listhptr[iuo] + j 
545
                jo = listh[ind] - 1
546
                juo = indxuo[jo]
547
                kx = dot(kpt_c, x_sparse[:, ind])
548
                phasef = exp(1.0j * kx)
549
                h_MM[iuo, juo] += phasef * h_sparse[ind, spin] 
550
                s_MM[iuo, juo] += phasef * s_sparse[ind]
551
        
552
        if remove_pbc is not None:
553
            direction, basis = remove_pbc
554
            centers_ic = get_bf_centers(dat.symbols, dat.pos_ac, basis)
555
            d = 'xyz'.index(direction)
556
            cutoff = dat.cell[d, d] * 0.5
557
            truncate_along_axis(h_MM, s_MM, direction, centers_ic, cutoff)
558
        
559
        h_MM *= complex(Rydberg)
560
        return h_MM, s_MM
561

    
562

    
563
def getrecord(fileobj, dtype):
564
    """Used to read in binary files. 
565
    """
566
    typetosize = {'l':4, 'f':4, 'd':8}# XXX np.int, np.float32, np.float64
567
    assert dtype in typetosize # XXX
568
    size = typetosize[dtype]
569
    record = array.array('l')
570
    trunk = array.array(dtype)
571
    record.fromfile(fileobj, 1)
572
    nofelements = int(record[-1]) / size
573
    trunk.fromfile(fileobj, nofelements)
574
    record.fromfile(fileobj, 1)
575
    data = np.array(trunk, dtype=dtype)
576
    if len(data)==1:
577
        data = data[0]
578
    return data
579

    
580
def truncate_along_axis(h, s, direction, centers_ic, cutoff):
581
    """Truncate h and s such along a cartesian axis.
582

583
    Parameters:
584

585
    h: (N, N) ndarray
586
        Hamiltonian matrix.
587
    s: (N, N) ndarray
588
        Overlap matrix.
589
    direction: {'x', 'y', 'z'} 
590
        Truncate allong a cartesian axis.
591
    centers_ic: (N, 3) ndarray
592
        Centers of the basis functions.
593
    cutoff: float
594
        The (direction-axis projected) cutoff distance.
595
    """
596
    dtype = h.dtype
597
    ni = len(centers_ic)
598
    d = 'xyz'.index(direction)
599
    pos_i = centers_ic[:, d]
600
    for i in range(ni):
601
        dpos_i = abs(pos_i - pos_i[i])
602
        mask_i = (dpos_i < cutoff).astype(dtype)
603
        h[i, :] *= mask_i
604
        h[:, i] *= mask_i
605
        s[i, :] *= mask_i
606
        s[:, i] *= mask_i
607

    
608
def get_nao(symbol, basis):
609
    """Number of basis functions. 
610

611
    Parameters
612
    ==========
613
    symbol: str
614
        The chemical symbol.
615
    basis: str
616
        Basis function type.
617
    """
618
    ls = valence_config[symbol]
619
    nao = 0
620
    zeta = {'s':1, 'd':2, 't':3, 'q':4}
621
    nzeta = zeta[basis[0]]
622
    is_pol = 'p' in basis
623
    for l in ls: 
624
        nao += (2 * l + 1) * nzeta
625
    if is_pol:
626
        l_pol = None
627
        l = -1 
628
        while l_pol is None:
629
            l += 1
630
            if not l in ls:
631
                l_pol = l
632
        nao += 2 * l_pol + 1
633
    return nao        
634

    
635
def get_bf_centers(symbols, positions, basis):
636
    """Centers of basis functions.
637

638
    Parameters
639
    ==========
640
    symbols: str, (N, ) array_like 
641
        chemical symbol for each atom. 
642
    positions: float, (N, 3) array_like
643
        Positions of the atoms.
644
    basis: {str,  dict}
645
        Basis set specification as either a string or a dictionary
646

647
    Examples
648
    ========
649
    >>> symbols = ['O', 'H']
650
    >>> positions = [(0, 0, 0), (0, 0, 1)]
651
    >>> basis = 'sz'
652
    >>> print get_bf_centers(symbols, positions, basis)
653
    [[0 0 0]
654
     [0 0 0]
655
     [0 0 0]
656
     [0 0 0]
657
     [0 0 1]]
658
    >>> basis = {'H':'dz', None:'sz'}
659
    >>> print get_bf_centers(symbols, positions, basis)
660
    [[0 0 0]
661
     [0 0 0]
662
     [0 0 0]
663
     [0 0 0]
664
     [0 0 1]
665
     [0 0 1]]
666

667
    """
668
    centers_ic = []
669
    dict_basis = False
670
    if type(basis)==dict:
671
        dict_basis = True
672
    for symbol, pos in zip(symbols, positions):
673
        if dict_basis:
674
            if symbol not in basis:
675
                bas = basis[None]
676
            else:
677
                bas = basis[symbol]
678
        else:
679
            bas = basis
680
        for i in range(get_nao(symbol, bas)):
681
            centers_ic.append(pos)
682
    return np.asarray(centers_ic)
683

    
684
def fdfify(key):
685
    return key.lower().replace('_', '').replace('.', '').replace('-', '')
686

    
687
valence_config = {
688
    'H': (0,),
689
    'C': (0, 1),
690
    'N': (0, 1),
691
    'O': (0, 1),
692
    'S': (0, 1), 
693
    'Li': (0,),
694
    'Na': (0,),
695
    'Ni': (0, 2),
696
    'Cu': (0, 2), 
697
    'Pd': (0, 2), 
698
    'Ag': (0, 2), 
699
    'Pt': (0, 2), 
700
    'Au': (0, 2)}
701

    
702
keys_with_units = {
703
    'paoenergyshift': 'eV',
704
    'zmunitslength': 'Bohr',
705
    'zmunitsangle': 'rad',
706
    'zmforcetollength': 'eV/Ang',
707
    'zmforcetolangle': 'eV/rad',
708
    'zmmaxdispllength': 'Ang',
709
    'zmmaxdisplangle': 'rad',
710
    'meshcutoff': 'eV',
711
    'dmenergytolerance': 'eV',
712
    'electronictemperature': 'eV',
713
    'oneta': 'eV',
714
    'onetaalpha': 'eV',
715
    'onetabeta': 'eV',
716
    'onrclwf': 'Ang',
717
    'onchemicalpotentialrc': 'Ang',
718
    'onchemicalpotentialtemperature': 'eV',
719
    'mdmaxcgdispl': 'Ang',
720
    'mdmaxforcetol': 'eV/Ang',
721
    'mdmaxstresstol': 'eV/Ang**3',
722
    'mdlengthtimestep': 'fs',
723
    'mdinitialtemperature': 'eV',
724
    'mdtargettemperature': 'eV',
725
    'mdtargetpressure': 'eV/Ang**3',
726
    'mdnosemass': 'eV*fs**2',
727
    'mdparrinellorahmanmass': 'eV*fs**2',
728
    'mdtaurelax': 'fs',
729
    'mdbulkmodulus': 'eV/Ang**3',
730
    'mdfcdispl': 'Ang',
731
    'warningminimumatomicdistance': 'Ang',
732
    'rcspatial': 'Ang',
733
    'kgridcutoff': 'Ang',
734
    'latticeconstant': 'Ang'}