Statistiques
| Révision :

root / ase / io / aims.py @ 18

Historique | Voir | Annoter | Télécharger (10,59 ko)

1 1 tkerber
2 1 tkerber
def read_aims(filename):
3 1 tkerber
    """Import FHI-aims geometry type files.
4 1 tkerber

5 1 tkerber
    Reads unitcell, atom positions and constraints from
6 1 tkerber
    a geometry.in file.
7 1 tkerber
    """
8 1 tkerber
9 1 tkerber
    from ase import Atoms
10 1 tkerber
    from ase.constraints import FixAtoms, FixCartesian
11 1 tkerber
    import numpy as np
12 1 tkerber
13 1 tkerber
    atoms = Atoms()
14 1 tkerber
    fd = open(filename, 'r')
15 1 tkerber
    lines = fd.readlines()
16 1 tkerber
    fd.close()
17 1 tkerber
    positions = []
18 1 tkerber
    cell = []
19 1 tkerber
    symbols = []
20 1 tkerber
    fix = []
21 1 tkerber
    fix_cart = []
22 1 tkerber
    xyz = np.array([0, 0, 0])
23 1 tkerber
    i = -1
24 1 tkerber
    n_periodic = -1
25 1 tkerber
    periodic = np.array([False, False, False])
26 1 tkerber
    for n, line in enumerate(lines):
27 1 tkerber
        inp = line.split()
28 1 tkerber
        if inp == []:
29 1 tkerber
            continue
30 1 tkerber
        if inp[0] == 'atom':
31 1 tkerber
            if xyz.all():
32 1 tkerber
                fix.append(i)
33 1 tkerber
            elif xyz.any():
34 1 tkerber
                print 1
35 1 tkerber
                fix_cart.append(FixCartesian(i, xyz))
36 1 tkerber
            floatvect = float(inp[1]), float(inp[2]), float(inp[3])
37 1 tkerber
            positions.append(floatvect)
38 1 tkerber
            symbols.append(inp[-1])
39 1 tkerber
            i += 1
40 1 tkerber
            xyz = np.array([0, 0, 0])
41 1 tkerber
        elif inp[0] == 'lattice_vector':
42 1 tkerber
            floatvect = float(inp[1]), float(inp[2]), float(inp[3])
43 1 tkerber
            cell.append(floatvect)
44 1 tkerber
            n_periodic = n_periodic + 1
45 1 tkerber
            periodic[n_periodic] = True
46 1 tkerber
        if inp[0] == 'constrain_relaxation':
47 1 tkerber
            if inp[1] == '.true.':
48 1 tkerber
                fix.append(i)
49 1 tkerber
            elif inp[1] == 'x':
50 1 tkerber
                xyz[0] = 1
51 1 tkerber
            elif inp[1] == 'y':
52 1 tkerber
                xyz[1] = 1
53 1 tkerber
            elif inp[1] == 'z':
54 1 tkerber
                xyz[2] = 1
55 1 tkerber
    if xyz.all():
56 1 tkerber
        fix.append(i)
57 1 tkerber
    elif xyz.any():
58 1 tkerber
        fix_cart.append(FixCartesian(i, xyz))
59 1 tkerber
    atoms = Atoms(symbols, positions)
60 1 tkerber
    if periodic.all():
61 1 tkerber
        atoms.set_cell(cell)
62 1 tkerber
        atoms.set_pbc(periodic)
63 1 tkerber
    if len(fix):
64 1 tkerber
        atoms.set_constraint([FixAtoms(indices=fix)]+fix_cart)
65 1 tkerber
    else:
66 1 tkerber
        atoms.set_constraint(fix_cart)
67 1 tkerber
    return atoms
68 1 tkerber
69 1 tkerber
def write_aims(filename, atoms):
70 1 tkerber
    """Method to write FHI-aims geometry files.
71 1 tkerber

72 1 tkerber
    Writes the atoms positions and constraints (only FixAtoms is
73 1 tkerber
    supported at the moment).
74 1 tkerber
    """
75 1 tkerber
76 1 tkerber
    from ase.constraints import FixAtoms, FixCartesian
77 1 tkerber
    import numpy as np
78 1 tkerber
79 1 tkerber
    if isinstance(atoms, (list, tuple)):
80 1 tkerber
        if len(atoms) > 1:
81 1 tkerber
            raise RuntimeError("Don't know how to save more than "+
82 1 tkerber
                               "one image to FHI-aims input")
83 1 tkerber
        else:
84 1 tkerber
            atoms = atoms[0]
85 1 tkerber
86 1 tkerber
    fd = open(filename, 'w')
87 1 tkerber
    fd.write('#=======================================================\n')
88 1 tkerber
    fd.write('#FHI-aims file: '+filename+'\n')
89 1 tkerber
    fd.write('#Created using the Atomic Simulation Environment (ASE)\n')
90 1 tkerber
    fd.write('#=======================================================\n')
91 1 tkerber
    i = 0
92 1 tkerber
    if atoms.get_pbc().any():
93 1 tkerber
        for n, vector in enumerate(atoms.get_cell()):
94 1 tkerber
            fd.write('lattice_vector ')
95 1 tkerber
            for i in range(3):
96 1 tkerber
                fd.write('%16.16f ' % vector[i])
97 1 tkerber
            fd.write('\n')
98 1 tkerber
    fix_cart = np.zeros([len(atoms),3])
99 1 tkerber
100 1 tkerber
    if atoms.constraints:
101 1 tkerber
        for constr in atoms.constraints:
102 1 tkerber
            if isinstance(constr, FixAtoms):
103 1 tkerber
                fix_cart[constr.index] = [1,1,1]
104 1 tkerber
            elif isinstance(constr, FixCartesian):
105 1 tkerber
                fix_cart[constr.a] = -constr.mask+1
106 1 tkerber
107 1 tkerber
    for i, atom in enumerate(atoms):
108 1 tkerber
        fd.write('atom ')
109 1 tkerber
        for pos in atom.get_position():
110 1 tkerber
            fd.write('%16.16f ' % pos)
111 1 tkerber
        fd.write(atom.symbol)
112 1 tkerber
        fd.write('\n')
113 1 tkerber
# (1) all coords are constrained:
114 1 tkerber
        if fix_cart[i].all():
115 1 tkerber
            fd.write('constrain_relaxation .true.\n')
116 1 tkerber
# (2) some coords are constrained:
117 1 tkerber
        elif fix_cart[i].any():
118 1 tkerber
            xyz = fix_cart[i]
119 1 tkerber
            for n in range(3):
120 1 tkerber
                if xyz[n]:
121 1 tkerber
                    fd.write('constrain_relaxation %s\n' % 'xyz'[n])
122 1 tkerber
        if atom.charge:
123 1 tkerber
            fd.write('initial_charge %16.6f\n' % atom.charge)
124 1 tkerber
        if atom.magmom:
125 1 tkerber
            fd.write('initial_moment %16.6f\n' % atom.magmom)
126 1 tkerber
# except KeyError:
127 1 tkerber
#     continue
128 1 tkerber
129 1 tkerber
def read_energy(filename):
130 1 tkerber
    for line in open(filename, 'r'):
131 1 tkerber
        if line.startswith('  | Total energy corrected'):
132 1 tkerber
            E = float(line.split()[-2])
133 1 tkerber
    return E
134 1 tkerber
135 1 tkerber
def read_aims_output(filename, index = -1):
136 1 tkerber
    """  Import FHI-aims output files with all data available, i.e. relaxations,
137 1 tkerber
    MD information, force information etc etc etc. """
138 1 tkerber
    from ase import Atoms, Atom
139 1 tkerber
    from ase.calculators.singlepoint import SinglePointCalculator
140 1 tkerber
    from ase.units import Ang, fs
141 1 tkerber
    molecular_dynamics = False
142 1 tkerber
    fd = open(filename, 'r')
143 1 tkerber
    cell = []
144 1 tkerber
    images = []
145 1 tkerber
    n_periodic = -1
146 1 tkerber
    f = None
147 1 tkerber
    pbc = False
148 1 tkerber
    found_aims_calculator = False
149 1 tkerber
    v_unit = Ang/(1000.0*fs)
150 1 tkerber
    while True:
151 1 tkerber
        line = fd.readline()
152 1 tkerber
        if not line:
153 1 tkerber
            break
154 1 tkerber
        if "List of parameters used to initialize the calculator:" in line:
155 1 tkerber
            fd.readline()
156 1 tkerber
            calc = read_aims_calculator(fd)
157 1 tkerber
            calc.out = filename
158 1 tkerber
            found_aims_calculator = True
159 1 tkerber
        if "Number of atoms" in line:
160 1 tkerber
            inp = line.split()
161 1 tkerber
            n_atoms = int(inp[5])
162 1 tkerber
        if "| Unit cell:" in line:
163 1 tkerber
            if not pbc:
164 1 tkerber
                pbc = True
165 1 tkerber
                for i in range(3):
166 1 tkerber
                    inp = fd.readline().split()
167 1 tkerber
                    cell.append([inp[1],inp[2],inp[3]])
168 1 tkerber
        if "Atomic structure:" in line and not molecular_dynamics:
169 1 tkerber
            fd.readline()
170 1 tkerber
            atoms = Atoms()
171 1 tkerber
            for i in range(n_atoms):
172 1 tkerber
                inp = fd.readline().split()
173 1 tkerber
                atoms.append(Atom(inp[3],(inp[4],inp[5],inp[6])))
174 1 tkerber
        if "Complete information for previous time-step:" in line:
175 1 tkerber
            molecular_dynamics = True
176 1 tkerber
        if "Updated atomic structure:" in line and not molecular_dynamics:
177 1 tkerber
            fd.readline()
178 1 tkerber
            atoms = Atoms()
179 1 tkerber
            velocities = []
180 1 tkerber
            for i in range(n_atoms):
181 1 tkerber
                inp = fd.readline().split()
182 1 tkerber
                atoms.append(Atom(inp[4],(inp[1],inp[2],inp[3])))
183 1 tkerber
                if molecular_dynamics:
184 1 tkerber
                    inp = fd.readline().split()
185 1 tkerber
        if "Atomic structure (and velocities)" in line:
186 1 tkerber
            fd.readline()
187 1 tkerber
            atoms = Atoms()
188 1 tkerber
            velocities = []
189 1 tkerber
            for i in range(n_atoms):
190 1 tkerber
                inp = fd.readline().split()
191 1 tkerber
                atoms.append(Atom(inp[4],(inp[1],inp[2],inp[3])))
192 1 tkerber
                inp = fd.readline().split()
193 1 tkerber
                velocities += [[float(inp[1])*v_unit,float(inp[2])*v_unit,float(inp[3])*v_unit]]
194 1 tkerber
            atoms.set_velocities(velocities)
195 1 tkerber
            images.append(atoms)
196 1 tkerber
        if "Total atomic forces" in line:
197 1 tkerber
            f = []
198 1 tkerber
            for i in range(n_atoms):
199 1 tkerber
                inp = fd.readline().split()
200 1 tkerber
                f.append([float(inp[2]),float(inp[3]),float(inp[4])])
201 1 tkerber
            if not found_aims_calculator:
202 1 tkerber
                e = images[-1].get_potential_energy()
203 1 tkerber
                images[-1].set_calculator(SinglePointCalculator(e,f,None,None,atoms))
204 1 tkerber
            e = None
205 1 tkerber
            f = None
206 1 tkerber
        if "Total energy corrected" in line:
207 1 tkerber
            e = float(line.split()[5])
208 1 tkerber
            if pbc:
209 1 tkerber
                atoms.set_cell(cell)
210 1 tkerber
                atoms.pbc = True
211 1 tkerber
            if not found_aims_calculator:
212 1 tkerber
                atoms.set_calculator(SinglePointCalculator(e,None,None,None,atoms))
213 1 tkerber
            if not molecular_dynamics:
214 1 tkerber
                images.append(atoms)
215 1 tkerber
            e = None
216 1 tkerber
            if found_aims_calculator:
217 1 tkerber
                calc.set_results(images[-1])
218 1 tkerber
                images[-1].set_calculator(calc)
219 1 tkerber
    fd.close()
220 1 tkerber
    if molecular_dynamics:
221 1 tkerber
        images = images[1:]
222 1 tkerber
223 1 tkerber
    # return requested images, code borrowed from ase/io/trajectory.py
224 1 tkerber
    if isinstance(index, int):
225 1 tkerber
        return images[index]
226 1 tkerber
    else:
227 1 tkerber
        step = index.step or 1
228 1 tkerber
        if step > 0:
229 1 tkerber
            start = index.start or 0
230 1 tkerber
            if start < 0:
231 1 tkerber
                start += len(images)
232 1 tkerber
            stop = index.stop or len(images)
233 1 tkerber
            if stop < 0:
234 1 tkerber
                stop += len(images)
235 1 tkerber
        else:
236 1 tkerber
            if index.start is None:
237 1 tkerber
                start = len(images) - 1
238 1 tkerber
            else:
239 1 tkerber
                start = index.start
240 1 tkerber
                if start < 0:
241 1 tkerber
                    start += len(images)
242 1 tkerber
            if index.stop is None:
243 1 tkerber
                stop = -1
244 1 tkerber
            else:
245 1 tkerber
                stop = index.stop
246 1 tkerber
                if stop < 0:
247 1 tkerber
                    stop += len(images)
248 1 tkerber
        return [images[i] for i in range(start, stop, step)]
249 1 tkerber
250 1 tkerber
def read_aims_calculator(file):
251 1 tkerber
    """  found instructions for building an FHI-aims calculator in the output file,
252 1 tkerber
    read its specifications and return it. """
253 1 tkerber
    from ase.calculators.aims import Aims
254 1 tkerber
    calc = Aims()
255 1 tkerber
    while True:
256 1 tkerber
        line = file.readline()
257 1 tkerber
        if "=======================================================" in line:
258 1 tkerber
            break
259 1 tkerber
        else:
260 1 tkerber
            args = line.split()
261 1 tkerber
            key = args[0]
262 1 tkerber
            if key == '#':
263 1 tkerber
                comment = True
264 1 tkerber
            elif calc.float_params.has_key(key):
265 1 tkerber
                calc.float_params[key] = float(args[1])
266 1 tkerber
            elif calc.exp_params.has_key(key):
267 1 tkerber
                calc.exp_params[key] = float(args[1])
268 1 tkerber
            elif calc.string_params.has_key(key):
269 1 tkerber
                calc.string_params[key] = args[1]
270 1 tkerber
                if len(args) > 2:
271 1 tkerber
                    for s in args[2:]:
272 1 tkerber
                        calc.string_params[key] += " "+s
273 1 tkerber
            elif calc.int_params.has_key(key):
274 1 tkerber
                calc.int_params[key] = int(args[1])
275 1 tkerber
            elif calc.bool_params.has_key(key):
276 1 tkerber
                try:
277 1 tkerber
                    calc.bool_params[key] = bool(args[1])
278 1 tkerber
                except:
279 1 tkerber
                    if key == 'vdw_correction_hirshfeld':
280 1 tkerber
                        calc.bool_params[key] = True
281 1 tkerber
            elif calc.list_params.has_key(key):
282 1 tkerber
                if key == 'output':
283 1 tkerber
                    # build output string from args:
284 1 tkerber
                    out_option = ''
285 1 tkerber
                    for arg in args[1:]:
286 1 tkerber
                        out_option +=str(arg)+' '
287 1 tkerber
                    if calc.list_params['output'] is not None:
288 1 tkerber
                        calc.list_params['output'] += [out_option]
289 1 tkerber
                    else:
290 1 tkerber
                        calc.list_params['output'] = [out_option]
291 1 tkerber
                else:
292 1 tkerber
                    calc.list_params[key] = list(args[1:])
293 1 tkerber
            elif '#' in key:
294 1 tkerber
                key = key[1:]
295 1 tkerber
                if calc.input_parameters.has_key(key):
296 1 tkerber
                    calc.input_parameters[key] = args[1]
297 1 tkerber
                    if len(args) > 2:
298 1 tkerber
                        for s in args[2:]:
299 1 tkerber
                            calc.input_parameters[key] += " "+s
300 1 tkerber
            else:
301 1 tkerber
                raise TypeError('FHI-aims keyword not defined in ASE: ' + key + '. Please check.')
302 1 tkerber
    return calc