Statistiques
| Révision :

root / ase / calculators / dacapo.py @ 20

Historique | Voir | Annoter | Télécharger (6,67 ko)

1 1 tkerber
import numpy as np
2 1 tkerber
3 1 tkerber
from ase.old import OldASEListOfAtomsWrapper
4 1 tkerber
5 1 tkerber
try:
6 1 tkerber
    import Numeric as num
7 1 tkerber
except ImportError:
8 1 tkerber
    pass
9 1 tkerber
10 1 tkerber
def np2num(a, typecode=None):
11 1 tkerber
    if num.__version__ > '23.8':
12 1 tkerber
        return num.array(a, typecode)
13 1 tkerber
    if typecode is None:
14 1 tkerber
        typecode = num.Float
15 1 tkerber
    b = num.fromstring(a.tostring(), typecode)
16 1 tkerber
    b.shape = a.shape
17 1 tkerber
    return b
18 1 tkerber
19 1 tkerber
def restart(filename, **kwargs):
20 1 tkerber
    calc = Dacapo(filename, **kwargs)
21 1 tkerber
    atoms = calc.get_atoms()
22 1 tkerber
    return atoms, calc
23 1 tkerber
24 1 tkerber
class Dacapo:
25 1 tkerber
    def __init__(self, filename=None, stay_alive=False, stress=False,
26 1 tkerber
                 **kwargs):
27 1 tkerber
28 1 tkerber
        self.kwargs = kwargs
29 1 tkerber
        self.stay_alive = stay_alive
30 1 tkerber
        self.stress = stress
31 1 tkerber
32 1 tkerber
        if filename is not None:
33 1 tkerber
            from Dacapo import Dacapo
34 1 tkerber
            self.loa = Dacapo.ReadAtoms(filename, **kwargs)
35 1 tkerber
            self.calc = self.loa.GetCalculator()
36 1 tkerber
        else:
37 1 tkerber
            self.loa = None
38 1 tkerber
            self.calc = None
39 1 tkerber
40 1 tkerber
        self.pps = []
41 1 tkerber
42 1 tkerber
    def set_pp(self, Z, path):
43 1 tkerber
        self.pps.append((Z, path))
44 1 tkerber
45 1 tkerber
    def set_txt(self, txt):
46 1 tkerber
        if self.calc is None:
47 1 tkerber
            self.kwargs['txtout'] = txt
48 1 tkerber
        else:
49 1 tkerber
            self.calc.SetTxtFile(txt)
50 1 tkerber
51 1 tkerber
    def set_nc(self, nc):
52 1 tkerber
        if self.calc is None:
53 1 tkerber
            self.kwargs['out'] = nc
54 1 tkerber
        else:
55 1 tkerber
            self.calc.SetNetCDFFile(nc)
56 1 tkerber
57 1 tkerber
    def update(self, atoms):
58 1 tkerber
        from Dacapo import Dacapo
59 1 tkerber
        if self.calc is None:
60 1 tkerber
            if 'nbands' not in self.kwargs:
61 1 tkerber
                n = sum([valence[atom.symbol] for atom in atoms])
62 1 tkerber
                self.kwargs['nbands'] = int(n * 0.65) + 4
63 1 tkerber
64 1 tkerber
            magmoms = atoms.get_initial_magnetic_moments()
65 1 tkerber
            if magmoms.any():
66 1 tkerber
                self.kwargs['spinpol'] = True
67 1 tkerber
68 1 tkerber
            self.calc = Dacapo(**self.kwargs)
69 1 tkerber
70 1 tkerber
            if self.stay_alive:
71 1 tkerber
                self.calc.StayAliveOn()
72 1 tkerber
            else:
73 1 tkerber
                self.calc.StayAliveOff()
74 1 tkerber
75 1 tkerber
            if self.stress:
76 1 tkerber
                self.calc.CalculateStress()
77 1 tkerber
78 1 tkerber
            for Z, path in self.pps:
79 1 tkerber
                self.calc.SetPseudoPotential(Z, path)
80 1 tkerber
81 1 tkerber
        if self.loa is None:
82 1 tkerber
            from ASE import Atom, ListOfAtoms
83 1 tkerber
            numbers = atoms.get_atomic_numbers()
84 1 tkerber
            positions = atoms.get_positions()
85 1 tkerber
            magmoms = atoms.get_initial_magnetic_moments()
86 1 tkerber
            self.loa = ListOfAtoms([Atom(Z=numbers[a],
87 1 tkerber
                                         position=positions[a],
88 1 tkerber
                                         magmom=magmoms[a])
89 1 tkerber
                                    for a in range(len(atoms))],
90 1 tkerber
                                   cell=np2num(atoms.get_cell()),
91 1 tkerber
                                   periodic=tuple(atoms.get_pbc()))
92 1 tkerber
            self.loa.SetCalculator(self.calc)
93 1 tkerber
        else:
94 1 tkerber
            self.loa.SetCartesianPositions(np2num(atoms.get_positions()))
95 1 tkerber
            self.loa.SetUnitCell(np2num(atoms.get_cell()), fix=True)
96 1 tkerber
97 1 tkerber
    def get_atoms(self):
98 1 tkerber
        atoms = OldASEListOfAtomsWrapper(self.loa).copy()
99 1 tkerber
        atoms.set_calculator(self)
100 1 tkerber
        return atoms
101 1 tkerber
102 1 tkerber
    def get_potential_energy(self, atoms):
103 1 tkerber
        self.update(atoms)
104 1 tkerber
        return self.calc.GetPotentialEnergy()
105 1 tkerber
106 1 tkerber
    def get_forces(self, atoms):
107 1 tkerber
        self.update(atoms)
108 1 tkerber
        return np.array(self.calc.GetCartesianForces())
109 1 tkerber
110 1 tkerber
    def get_stress(self, atoms):
111 1 tkerber
        self.update(atoms)
112 1 tkerber
        stress = np.array(self.calc.GetStress())
113 1 tkerber
        if stress.ndim == 2:
114 1 tkerber
            return stress.ravel()[[0, 4, 8, 5, 2, 1]]
115 1 tkerber
        else:
116 1 tkerber
            return stress
117 1 tkerber
118 1 tkerber
    def calculation_required(self, atoms, quantities):
119 1 tkerber
        if self.calc is None:
120 1 tkerber
            return True
121 1 tkerber
122 1 tkerber
        if atoms != self.get_atoms():
123 1 tkerber
            return True
124 1 tkerber
125 1 tkerber
        return False
126 1 tkerber
127 1 tkerber
    def get_number_of_bands(self):
128 1 tkerber
        return self.calc.GetNumberOfBands()
129 1 tkerber
130 1 tkerber
    def get_k_point_weights(self):
131 1 tkerber
        return np.array(self.calc.GetIBZKPointWeights())
132 1 tkerber
133 1 tkerber
    def get_number_of_spins(self):
134 1 tkerber
        return 1 + int(self.calc.GetSpinPolarized())
135 1 tkerber
136 1 tkerber
    def get_eigenvalues(self, kpt=0, spin=0):
137 1 tkerber
        return np.array(self.calc.GetEigenvalues(kpt, spin))
138 1 tkerber
139 1 tkerber
    def get_fermi_level(self):
140 1 tkerber
        return self.calc.GetFermiLevel()
141 1 tkerber
142 1 tkerber
    def get_number_of_grid_points(self):
143 1 tkerber
        return np.array(self.get_pseudo_wave_function(0, 0, 0).shape)
144 1 tkerber
145 1 tkerber
    def get_pseudo_density(self, spin=0):
146 1 tkerber
        return np.array(self.calc.GetDensityArray(s))
147 1 tkerber
148 1 tkerber
    def get_pseudo_wave_function(self, band=0, kpt=0, spin=0, pad=True):
149 1 tkerber
        kpt_c = self.get_bz_k_points()[kpt]
150 1 tkerber
        state = self.calc.GetElectronicStates().GetState(band=band, spin=spin,
151 1 tkerber
                                                         kptindex=kpt)
152 1 tkerber
153 1 tkerber
        # Get wf, without bloch phase (Phase = True doesn't do anything!)
154 1 tkerber
        wave = state.GetWavefunctionOnGrid(phase=False)
155 1 tkerber
156 1 tkerber
        # Add bloch phase if this is not the Gamma point
157 1 tkerber
        if np.all(kpt_c == 0):
158 1 tkerber
            return wave
159 1 tkerber
        coord = state.GetCoordinates()
160 1 tkerber
        phase = coord[0] * kpt_c[0] + coord[1] * kpt_c[1] + coord[2] * kpt_c[2]
161 1 tkerber
        return np.array(wave) * np.exp(-2.j * np.pi * phase) # sign! XXX
162 1 tkerber
163 1 tkerber
        #return np.array(self.calc.GetWaveFunctionArray(n, k, s)) # No phase!
164 1 tkerber
165 1 tkerber
    def get_bz_k_points(self):
166 1 tkerber
        return np.array(self.calc.GetBZKPoints())
167 1 tkerber
168 1 tkerber
    def get_ibz_k_points(self):
169 1 tkerber
        return np.array(self.calc.GetIBZKPoints())
170 1 tkerber
171 1 tkerber
    def get_wannier_localization_matrix(self, nbands, dirG, kpoint,
172 1 tkerber
                                        nextkpoint, G_I, spin):
173 1 tkerber
        return np.array(self.calc.GetWannierLocalizationMatrix(
174 1 tkerber
            G_I=G_I.tolist(), nbands=nbands, dirG=dirG.tolist(),
175 1 tkerber
            kpoint=kpoint, nextkpoint=nextkpoint, spin=spin))
176 1 tkerber
177 1 tkerber
    def initial_wannier(self, initialwannier, kpointgrid, fixedstates,
178 1 tkerber
                        edf, spin):
179 1 tkerber
        # Use initial guess to determine U and C
180 1 tkerber
        init = self.calc.InitialWannier(initialwannier, self.atoms,
181 1 tkerber
                                        np2num(kpointgrid, num.Int))
182 1 tkerber
183 1 tkerber
        states = self.calc.GetElectronicStates()
184 1 tkerber
        waves = [[state.GetWaveFunction()
185 1 tkerber
                  for state in states.GetStatesKPoint(k, spin)]
186 1 tkerber
                 for k in self.calc.GetIBZKPoints()]
187 1 tkerber
188 1 tkerber
        init.SetupMMatrix(waves, self.calc.GetBZKPoints())
189 1 tkerber
        c, U = init.GetListOfCoefficientsAndRotationMatrices(
190 1 tkerber
            (self.calc.GetNumberOfBands(), fixedstates, edf))
191 1 tkerber
        U = np.array(U)
192 1 tkerber
        for k in range(len(c)):
193 1 tkerber
            c[k] = np.array(c[k])
194 1 tkerber
        return c, U
195 1 tkerber
196 1 tkerber
valence = {
197 1 tkerber
'H':   1,
198 1 tkerber
'B':   3,
199 1 tkerber
'C':   4,
200 1 tkerber
'N':   5,
201 1 tkerber
'O':   6,
202 1 tkerber
'Li':  1,
203 1 tkerber
'Na':  1,
204 1 tkerber
'K':   9,
205 1 tkerber
'Mg':  8,
206 1 tkerber
'Ca': 10,
207 1 tkerber
'Sr': 10,
208 1 tkerber
'Al':  3,
209 1 tkerber
'Ga': 13,
210 1 tkerber
'Sc': 11,
211 1 tkerber
'Ti': 12,
212 1 tkerber
'V':  13,
213 1 tkerber
'Cr': 14,
214 1 tkerber
'Mn':  7,
215 1 tkerber
'Fe':  8,
216 1 tkerber
'Co':  9,
217 1 tkerber
'Ni': 10,
218 1 tkerber
'Cu': 11,
219 1 tkerber
'Zn': 12,
220 1 tkerber
'Y':  11,
221 1 tkerber
'Zr': 12,
222 1 tkerber
'Nb': 13,
223 1 tkerber
'Mo':  6,
224 1 tkerber
'Ru':  8,
225 1 tkerber
'Rh':  9,
226 1 tkerber
'Pd': 10,
227 1 tkerber
'Ag': 11,
228 1 tkerber
'Cd': 12,
229 1 tkerber
'Ir': 9,
230 1 tkerber
'Pt': 10,
231 1 tkerber
'Au': 11,
232 1 tkerber
}