-
Notifications
You must be signed in to change notification settings - Fork 1
/
espresso.py
3682 lines (3180 loc) · 137 KB
/
espresso.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
"""Reads Quantum ESPRESSO files.
Read multiple structures and results from pw.x output files. Read
structures from pw.x input files.
Built for PWSCF v.5.3.0 but should work with earlier and later versions.
Can deal with most major functionality, but might fail with ibrav =/= 0
or crystal_sg positions.
Units are converted using CODATA 2006, as used internally by Quantum
ESPRESSO.
"""
import os
import operator as op
import re
import warnings
from collections import OrderedDict
from os import path
import numpy as np
from ase.atoms import Atoms
from ase.calculators.singlepoint import (SinglePointDFTCalculator,
SinglePointKPoint)
from ase.calculators.calculator import kpts2ndarray, kpts2sizeandoffsets
from ase.dft.kpoints import kpoint_convert
from ase.constraints import FixAtoms, FixCartesian
from ase.data import chemical_symbols, atomic_numbers
from ase.units import create_units
from ase.utils import iofunction
# Quantum ESPRESSO uses CODATA 2006 internally
units = create_units('2006')
# Section identifiers
_PW_START = 'Program PWSCF'
_PW_END = 'End of self-consistent calculation'
_PW_CELL = 'CELL_PARAMETERS'
_PW_POS = 'ATOMIC_POSITIONS'
_PW_MAGMOM = 'Magnetic moment per site'
_PW_FORCE = 'Forces acting on atoms'
_PW_TOTEN = '! total energy'
_PW_STRESS = 'total stress'
_PW_FERMI = 'the Fermi energy is'
_PW_HIGHEST_OCCUPIED = 'highest occupied level'
_PW_HIGHEST_OCCUPIED_LOWEST_FREE = 'highest occupied, lowest unoccupied level'
_PW_KPTS = 'number of k points='
_PW_BANDS = _PW_END
_PW_BANDSTRUCTURE = 'End of band structure calculation'
class Namelist(OrderedDict):
"""Case insensitive dict that emulates Fortran Namelists."""
def __contains__(self, key):
return super(Namelist, self).__contains__(key.lower())
def __delitem__(self, key):
return super(Namelist, self).__delitem__(key.lower())
def __getitem__(self, key):
return super(Namelist, self).__getitem__(key.lower())
def __setitem__(self, key, value):
super(Namelist, self).__setitem__(key.lower(), value)
def get(self, key, default=None):
return super(Namelist, self).get(key.lower(), default)
@iofunction('rU')
def read_espresso_out(fileobj, index=-1, results_required=True):
"""Reads Quantum ESPRESSO output files.
The atomistic configurations as well as results (energy, force, stress,
magnetic moments) of the calculation are read for all configurations
within the output file.
Will probably raise errors for broken or incomplete files.
Parameters
----------
fileobj : file|str
A file like object or filename
index : slice
The index of configurations to extract.
results_required : bool
If True, atomistic configurations that do not have any
associated results will not be included. This prevents double
printed configurations and incomplete calculations from being
returned as the final configuration with no results data.
Yields
------
structure : Atoms
The next structure from the index slice. The Atoms has a
SinglePointCalculator attached with any results parsed from
the file.
"""
# work with a copy in memory for faster random access
pwo_lines = fileobj.readlines()
# TODO: index -1 special case?
# Index all the interesting points
indexes = {
_PW_START: [],
_PW_END: [],
_PW_CELL: [],
_PW_POS: [],
_PW_MAGMOM: [],
_PW_FORCE: [],
_PW_TOTEN: [],
_PW_STRESS: [],
_PW_FERMI: [],
_PW_HIGHEST_OCCUPIED: [],
_PW_HIGHEST_OCCUPIED_LOWEST_FREE: [],
_PW_KPTS: [],
_PW_BANDS: [],
_PW_BANDSTRUCTURE: [],
}
for idx, line in enumerate(pwo_lines):
for identifier in indexes:
if identifier in line:
indexes[identifier].append(idx)
# Configurations are either at the start, or defined in ATOMIC_POSITIONS
# in a subsequent step. Can deal with concatenated output files.
all_config_indexes = sorted(indexes[_PW_START] +
indexes[_PW_POS])
# Slice only requested indexes
# setting results_required argument stops configuration-only
# structures from being returned. This ensures the [-1] structure
# is one that has results. Two cases:
# - SCF of last configuration is not converged, job terminated
# abnormally.
# - 'relax' and 'vc-relax' re-prints the final configuration but
# only 'vc-relax' recalculates.
if results_required:
results_indexes = sorted(indexes[_PW_TOTEN] + indexes[_PW_FORCE] +
indexes[_PW_STRESS] + indexes[_PW_MAGMOM] +
indexes[_PW_BANDS] +
indexes[_PW_BANDSTRUCTURE])
# Prune to only configurations with results data before the next
# configuration
results_config_indexes = []
for config_index, config_index_next in zip(
all_config_indexes,
all_config_indexes[1:] + [len(pwo_lines)]):
if any([config_index < results_index < config_index_next
for results_index in results_indexes]):
results_config_indexes.append(config_index)
# slice from the subset
image_indexes = results_config_indexes[index]
else:
image_indexes = all_config_indexes[index]
# Extract initialisation information each time PWSCF starts
# to add to subsequent configurations. Use None so slices know
# when to fill in the blanks.
pwscf_start_info = dict((idx, None) for idx in indexes[_PW_START])
for image_index in image_indexes:
# Find the nearest calculation start to parse info. Needed in,
# for example, relaxation where cell is only printed at the
# start.
if image_index in indexes[_PW_START]:
prev_start_index = image_index
else:
# The greatest start index before this structure
prev_start_index = [idx for idx in indexes[_PW_START]
if idx < image_index][-1]
# add structure to reference if not there
if pwscf_start_info[prev_start_index] is None:
pwscf_start_info[prev_start_index] = parse_pwo_start(
pwo_lines, prev_start_index)
# Get the bounds for information for this structure. Any associated
# values will be between the image_index and the following one,
# EXCEPT for cell, which will be 4 lines before if it exists.
for next_index in all_config_indexes:
if next_index > image_index:
break
else:
# right to the end of the file
next_index = len(pwo_lines)
# Get the structure
# Use this for any missing data
prev_structure = pwscf_start_info[prev_start_index]['atoms']
if image_index in indexes[_PW_START]:
structure = prev_structure.copy() # parsed from start info
else:
if _PW_CELL in pwo_lines[image_index - 5]:
# CELL_PARAMETERS would be just before positions if present
cell, cell_alat = get_cell_parameters(
pwo_lines[image_index - 5:image_index])
else:
cell = prev_structure.cell
cell_alat = pwscf_start_info[prev_start_index]['alat']
# give at least enough lines to parse the positions
# should be same format as input card
n_atoms = len(prev_structure)
positions_card = get_atomic_positions(
pwo_lines[image_index:image_index + n_atoms + 1],
n_atoms=n_atoms, cell=cell, alat=cell_alat)
# convert to Atoms object
symbols = [label_to_symbol(position[0]) for position in
positions_card]
positions = [position[1] for position in positions_card]
structure = Atoms(symbols=symbols, positions=positions, cell=cell,
pbc=True)
# Extract calculation results
# Energy
energy = None
for energy_index in indexes[_PW_TOTEN]:
if image_index < energy_index < next_index:
energy = float(
pwo_lines[energy_index].split()[-2]) * units['Ry']
# Forces
forces = None
for force_index in indexes[_PW_FORCE]:
if image_index < force_index < next_index:
# Before QE 5.3 'negative rho' added 2 lines before forces
# Use exact lines to stop before 'non-local' forces
# in high verbosity
if not pwo_lines[force_index + 2].strip():
force_index += 4
else:
force_index += 2
# assume contiguous
forces = [
[float(x) for x in force_line.split()[-3:]] for force_line
in pwo_lines[force_index:force_index + len(structure)]]
forces = np.array(forces) * units['Ry'] / units['Bohr']
# Stress
stress = None
for stress_index in indexes[_PW_STRESS]:
if image_index < stress_index < next_index:
sxx, sxy, sxz = pwo_lines[stress_index + 1].split()[:3]
_, syy, syz = pwo_lines[stress_index + 2].split()[:3]
_, _, szz = pwo_lines[stress_index + 3].split()[:3]
stress = np.array([sxx, syy, szz, syz, sxz, sxy], dtype=float)
# sign convention is opposite of ase
stress *= -1 * units['Ry'] / (units['Bohr'] ** 3)
# Magmoms
magmoms = None
for magmoms_index in indexes[_PW_MAGMOM]:
if image_index < magmoms_index < next_index:
magmoms = [
float(mag_line.split()[5]) for mag_line
in pwo_lines[magmoms_index + 1:
magmoms_index + 1 + len(structure)]]
# Fermi level / highest occupied level
efermi = None
for fermi_index in indexes[_PW_FERMI]:
if image_index < fermi_index < next_index:
efermi = float(pwo_lines[fermi_index].split()[-2])
if efermi is None:
for ho_index in indexes[_PW_HIGHEST_OCCUPIED]:
if image_index < ho_index < next_index:
efermi = float(pwo_lines[ho_index].split()[-1])
if efermi is None:
for holf_index in indexes[_PW_HIGHEST_OCCUPIED_LOWEST_FREE]:
if image_index < holf_index < next_index:
efermi = float(pwo_lines[holf_index].split()[-2])
# K-points
ibzkpts = None
weights = None
kpoints_warning = "Number of k-points >= 100: " + \
"set verbosity='high' to print them."
for kpts_index in indexes[_PW_KPTS]:
nkpts = int(pwo_lines[kpts_index].split()[4])
kpts_index += 2
if pwo_lines[kpts_index].strip() == kpoints_warning:
continue
# QE prints the k-points in units of 2*pi/alat
# with alat defined as the length of the first
# cell vector
cell = structure.get_cell()
alat = np.linalg.norm(cell[0])
ibzkpts = []
weights = []
for i in range(nkpts):
L = pwo_lines[kpts_index + i].split()
weights.append(float(L[-1]))
coord = np.array([L[-6], L[-5], L[-4].strip('),')],
dtype=float)
coord *= 2 * np.pi / alat
coord = kpoint_convert(cell, ckpts_kv=coord)
ibzkpts.append(coord)
ibzkpts = np.array(ibzkpts)
weights = np.array(weights)
# Bands
kpts = None
kpoints_warning = "Number of k-points >= 100: " + \
"set verbosity='high' to print the bands."
for bands_index in indexes[_PW_BANDS] + indexes[_PW_BANDSTRUCTURE]:
if image_index < bands_index < next_index:
bands_index += 2
if pwo_lines[bands_index].strip() == kpoints_warning:
continue
assert ibzkpts is not None
spin, bands, eigenvalues = 0, [], [[], []]
while True:
L = pwo_lines[bands_index].replace('-', ' -').split()
if len(L) == 0:
if len(bands) > 0:
eigenvalues[spin].append(bands)
bands = []
elif L == ['occupation', 'numbers']:
# Skip the lines with the occupation numbers
bands_index += len(eigenvalues[spin][0]) // 8 + 1
elif L[0] == 'k' and L[1].startswith('='):
pass
elif 'SPIN' in L:
if 'DOWN' in L:
spin += 1
else:
try:
bands.extend(map(float, L))
except ValueError:
break
bands_index += 1
if spin == 1:
assert len(eigenvalues[0]) == len(eigenvalues[1])
assert len(eigenvalues[0]) == len(ibzkpts), \
(np.shape(eigenvalues), len(ibzkpts))
kpts = []
for s in range(spin + 1):
for w, k, e in zip(weights, ibzkpts, eigenvalues[s]):
kpt = SinglePointKPoint(w, s, k, eps_n=e)
kpts.append(kpt)
# Put everything together
#
# I have added free_energy. Can and should we distinguish
# energy and free_energy? --askhl
calc = SinglePointDFTCalculator(structure, energy=energy,
free_energy=energy,
forces=forces, stress=stress,
magmoms=magmoms, efermi=efermi,
ibzkpts=ibzkpts)
calc.kpts = kpts
structure.calc = calc
yield structure
def parse_pwo_start(lines, index=0):
"""Parse Quantum ESPRESSO calculation info from lines,
starting from index. Return a dictionary containing extracted
information.
- `celldm(1)`: lattice parameters (alat)
- `cell`: unit cell in Angstrom
- `symbols`: element symbols for the structure
- `positions`: cartesian coordinates of atoms in Angstrom
- `atoms`: an `ase.Atoms` object constructed from the extracted data
Parameters
----------
lines : list[str]
Contents of PWSCF output file.
index : int
Line number to begin parsing. Only first calculation will
be read.
Returns
-------
info : dict
Dictionary of calculation parameters, including `celldm(1)`, `cell`,
`symbols`, `positions`, `atoms`.
Raises
------
KeyError
If interdependent values cannot be found (especially celldm(1))
an error will be raised as other quantities cannot then be
calculated (e.g. cell and positions).
"""
# TODO: extend with extra DFT info?
info = {}
for idx, line in enumerate(lines[index:], start=index):
if 'celldm(1)' in line:
# celldm(1) has more digits than alat!!
info['celldm(1)'] = float(line.split()[1]) * units['Bohr']
info['alat'] = info['celldm(1)']
elif 'number of atoms/cell' in line:
info['nat'] = int(line.split()[-1])
elif 'number of atomic types' in line:
info['ntyp'] = int(line.split()[-1])
elif 'crystal axes:' in line:
info['cell'] = info['celldm(1)'] * np.array([
[float(x) for x in lines[idx + 1].split()[3:6]],
[float(x) for x in lines[idx + 2].split()[3:6]],
[float(x) for x in lines[idx + 3].split()[3:6]]])
elif 'positions (alat units)' in line:
info['symbols'], info['positions'] = [], []
for at_line in lines[idx + 1:idx + 1 + info['nat']]:
sym, x, y, z = parse_position_line(at_line)
info['symbols'].append(label_to_symbol(sym))
info['positions'].append([x * info['celldm(1)'],
y * info['celldm(1)'],
z * info['celldm(1)']])
# This should be the end of interesting info.
# Break here to avoid dealing with large lists of kpoints.
# Will need to be extended for DFTCalculator info.
break
# Make atoms for convenience
info['atoms'] = Atoms(symbols=info['symbols'],
positions=info['positions'],
cell=info['cell'], pbc=True)
return info
def parse_position_line(line):
"""Parse a single line from a pw.x output file.
The line must contain information about the atomic symbol and the position,
e.g.
995 Sb tau( 995) = ( 1.4212023 0.7037863 0.1242640 )
Parameters
----------
line : str
Line to be parsed.
Returns
-------
sym : str
Atomic symbol.
x : float
x-position.
y : float
y-position.
z : float
z-position.
"""
pat = re.compile(r'\s*\d+\s*(\S+)\s*tau\(\s*\d+\)\s*='
r'\s*\(\s*(\S+)\s+(\S+)\s+(\S+)\s*\)')
match = pat.match(line)
assert match is not None
sym, x, y, z = match.group(1, 2, 3, 4)
return sym, float(x), float(y), float(z)
@iofunction('rU')
def read_espresso_in(fileobj):
"""Parse a Quantum ESPRESSO input files, '.in', '.pwi'.
ESPRESSO inputs are generally a fortran-namelist format with custom
blocks of data. The namelist is parsed as a dict and an atoms object
is constructed from the included information.
Parameters
----------
fileobj : file | str
A file-like object that supports line iteration with the contents
of the input file, or a filename.
Returns
-------
atoms : Atoms
Structure defined in the input file.
Raises
------
KeyError
Raised for missing keys that are required to process the file
"""
# parse namelist section and extract remaining lines
data, card_lines = read_fortran_namelist(fileobj)
# get the cell if ibrav=0
if 'system' not in data:
raise KeyError('Required section &SYSTEM not found.')
elif 'ibrav' not in data['system']:
raise KeyError('ibrav is required in &SYSTEM')
elif data['system']['ibrav'] == 0:
# celldm(1) is in Bohr, A is in angstrom. celldm(1) will be
# used even if A is also specified.
if 'celldm(1)' in data['system']:
alat = data['system']['celldm(1)'] * units['Bohr']
elif 'A' in data['system']:
alat = data['system']['A']
else:
alat = None
cell, cell_alat = get_cell_parameters(card_lines, alat=alat)
else:
alat, cell = ibrav_to_cell(data['system'])
# species_info holds some info for each element
species_card = get_atomic_species(card_lines, n_species=data['system']['ntyp'])
species_info = {}
for ispec, (label, weight, pseudo) in enumerate(species_card):
symbol = label_to_symbol(label)
valence = get_valence_electrons(symbol, data, pseudo)
# starting_magnetization is in fractions of valence electrons
magnet_key = "starting_magnetization({0})".format(ispec + 1)
magmom = valence * data["system"].get(magnet_key, 0.0)
species_info[symbol] = {"weight": weight, "pseudo": pseudo,
"valence": valence, "magmom": magmom}
positions_card = get_atomic_positions(
card_lines, n_atoms=data['system']['nat'], cell=cell, alat=alat)
symbols = [label_to_symbol(position[0]) for position in positions_card]
positions = [position[1] for position in positions_card]
magmoms = [species_info[symbol]["magmom"] for symbol in symbols]
# TODO: put more info into the atoms object
# e.g magmom, forces.
atoms = Atoms(symbols=symbols, positions=positions, cell=cell, pbc=True,
magmoms=magmoms)
return atoms
def ibrav_to_cell(system):
"""
Convert a value of ibrav to a cell. Any unspecified lattice dimension
is set to 0.0, but will not necessarily raise an error. Also return the
lattice parameter.
Parameters
----------
system : dict
The &SYSTEM section of the input file, containing the 'ibrav' setting,
and either celldm(1)..(6) or a, b, c, cosAB, cosAC, cosBC.
Returns
-------
alat, cell : float, np.array
Cell parameter in Angstrom, and
The 3x3 array representation of the cell.
Raises
------
KeyError
Raise an error if any required keys are missing.
NotImplementedError
Only a limited number of ibrav settings can be parsed. An error
is raised if the ibrav interpretation is not implemented.
"""
if 'celldm(1)' in system and 'a' in system:
raise KeyError('do not specify both celldm and a,b,c!')
elif 'celldm(1)' in system:
# celldm(x) in bohr
alat = system['celldm(1)'] * units['Bohr']
b_over_a = system.get('celldm(2)', 0.0)
c_over_a = system.get('celldm(3)', 0.0)
cosab = system.get('celldm(4)', 0.0)
cosac = system.get('celldm(5)', 0.0)
cosbc = 0.0
if system['ibrav'] == 14:
cosbc = system.get('celldm(4)', 0.0)
cosac = system.get('celldm(5)', 0.0)
cosab = system.get('celldm(6)', 0.0)
elif 'a' in system:
# a, b, c, cosAB, cosAC, cosBC in Angstrom
alat = system['a']
b_over_a = system.get('b', 0.0) / alat
c_over_a = system.get('c', 0.0) / alat
cosab = system.get('cosab', 0.0)
cosac = system.get('cosac', 0.0)
cosbc = system.get('cosbc', 0.0)
else:
raise KeyError("Missing celldm(1) or a cell parameter.")
if system['ibrav'] == 1:
cell = np.identity(3) * alat
elif system['ibrav'] == 2:
cell = np.array([[-1.0, 0.0, 1.0],
[0.0, 1.0, 1.0],
[-1.0, 1.0, 0.0]]) * (alat / 2)
elif system['ibrav'] == 3:
cell = np.array([[1.0, 1.0, 1.0],
[-1.0, 1.0, 1.0],
[-1.0, -1.0, 1.0]]) * (alat / 2)
elif system['ibrav'] == -3:
cell = np.array([[-1.0, 1.0, 1.0],
[1.0, -1.0, 1.0],
[1.0, 1.0, -1.0]]) * (alat / 2)
elif system['ibrav'] == 4:
cell = np.array([[1.0, 0.0, 0.0],
[-0.5, 0.5 * 3**0.5, 0.0],
[0.0, 0.0, c_over_a]]) * alat
elif system['ibrav'] == 5:
tx = ((1.0 - cosab) / 2.0)**0.5
ty = ((1.0 - cosab) / 6.0)**0.5
tz = ((1 + 2 * cosab) / 3.0)**0.5
cell = np.array([[tx, -ty, tz],
[0, 2 * ty, tz],
[-tx, -ty, tz]]) * alat
elif system['ibrav'] == -5:
ty = ((1.0 - cosab) / 6.0)**0.5
tz = ((1 + 2 * cosab) / 3.0)**0.5
a_prime = alat / 3**0.5
u = tz - 2 * 2**0.5 * ty
v = tz + 2**0.5 * ty
cell = np.array([[u, v, v],
[v, u, v],
[v, v, u]]) * a_prime
elif system['ibrav'] == 6:
cell = np.array([[1.0, 0.0, 0.0],
[0.0, 1.0, 0.0],
[0.0, 0.0, c_over_a]]) * alat
elif system['ibrav'] == 7:
cell = np.array([[1.0, -1.0, c_over_a],
[1.0, 1.0, c_over_a],
[-1.0, -1.0, c_over_a]]) * (alat / 2)
elif system['ibrav'] == 8:
cell = np.array([[1.0, 0.0, 0.0],
[0.0, b_over_a, 0.0],
[0.0, 0.0, c_over_a]]) * alat
elif system['ibrav'] == 9:
cell = np.array([[1.0 / 2.0, b_over_a / 2.0, 0.0],
[-1.0 / 2.0, b_over_a / 2.0, 0.0],
[0.0, 0.0, c_over_a]]) * alat
elif system['ibrav'] == -9:
cell = np.array([[1.0 / 2.0, -b_over_a / 2.0, 0.0],
[1.0 / 2.0, b_over_a / 2.0, 0.0],
[0.0, 0.0, c_over_a]]) * alat
elif system['ibrav'] == 10:
cell = np.array([[1.0 / 2.0, 0.0, c_over_a / 2.0],
[1.0 / 2.0, b_over_a / 2.0, 0.0],
[0.0, b_over_a / 2.0, c_over_a / 2.0]]) * alat
elif system['ibrav'] == 11:
cell = np.array([[1.0 / 2.0, b_over_a / 2.0, c_over_a / 2.0],
[-1.0 / 2.0, b_over_a / 2.0, c_over_a / 2.0],
[-1.0, 2.0, -b_over_a / 2.0, c_over_a / 2.0]]) * alat
elif system['ibrav'] == 12:
sinab = (1.0 - cosab**2)**0.5
cell = np.array([[1.0, 0.0, 0.0],
[b_over_a * cosab, b_over_a * sinab, 0.0],
[0.0, 0.0, c_over_a]]) * alat
elif system['ibrav'] == -12:
sinac = (1.0 - cosac**2)**0.5
cell = np.array([[1.0, 0.0, 0.0],
[0.0, b_over_a, 0.0],
[c_over_a * cosac, 0.0, c_over_a * sinac]]) * alat
elif system['ibrav'] == 13:
sinab = (1.0 - cosab**2)**0.5
cell = np.array([[1.0 / 2.0, 0.0, -c_over_a / 2.0],
[b_over_a * cosab, b_over_a * sinab, 0.0],
[1.0 / 2.0, 0.0, c_over_a / 2.0]]) * alat
elif system['ibrav'] == 14:
sinab = (1.0 - cosab**2)**0.5
v3 = [c_over_a * cosac,
c_over_a * (cosbc - cosac * cosab) / sinab,
c_over_a * ((1 + 2 * cosbc * cosac * cosab
- cosbc**2 - cosac**2 - cosab**2)**0.5) / sinab]
cell = np.array([[1.0, 0.0, 0.0],
[b_over_a * cosab, b_over_a * sinab, 0.0],
v3]) * alat
else:
raise NotImplementedError('ibrav = {0} is not implemented'
''.format(system['ibrav']))
return alat, cell
def get_pseudo_dirs(data):
"""Guess a list of possible locations for pseudopotential files.
Parameters
----------
data : Namelist
Namelist representing the quantum espresso input parameters
Returns
-------
pseudo_dirs : list[str]
A list of directories where pseudopotential files could be located.
"""
pseudo_dirs = []
if 'pseudo_dir' in data['control']:
pseudo_dirs.append(data['control']['pseudo_dir'])
if 'ESPRESSO_PSEUDO' in os.environ:
pseudo_dirs.append(os.environ['ESPRESSO_PSEUDO'])
pseudo_dirs.append(path.expanduser('~/espresso/pseudo/'))
return pseudo_dirs
def get_valence_electrons(symbol, data, pseudo=None):
"""The number of valence electrons for a atomic symbol.
Parameters
----------
symbol : str
Chemical symbol
data : Namelist
Namelist representing the quantum espresso input parameters
pseudo : str, optional
File defining the pseudopotential to be used. If missing a fallback
to the number of valence electrons recommended at
http://materialscloud.org/sssp/ is employed.
"""
if pseudo is None:
pseudo = '{}_dummy.UPF'.format(symbol)
for pseudo_dir in get_pseudo_dirs(data):
if path.exists(path.join(pseudo_dir, pseudo)):
valence = grep_valence(path.join(pseudo_dir, pseudo))
break
else: # not found in a file
valence = SSSP_VALENCE[atomic_numbers[symbol]]
return valence
def get_atomic_positions(lines, n_atoms, cell=None, alat=None):
"""Parse atom positions from ATOMIC_POSITIONS card.
Parameters
----------
lines : list[str]
A list of lines containing the ATOMIC_POSITIONS card.
n_atoms : int
Expected number of atoms. Only this many lines will be parsed.
cell : np.array
Unit cell of the crystal. Only used with crystal coordinates.
alat : float
Lattice parameter for atomic coordinates. Only used for alat case.
Returns
-------
positions : list[(str, (float, float, float), (float, float, float))]
A list of the ordered atomic positions in the format:
label, (x, y, z), (if_x, if_y, if_z)
Force multipliers are set to None if not present.
Raises
------
ValueError
Any problems parsing the data result in ValueError
"""
positions = None
# no blanks or comment lines, can the consume n_atoms lines for positions
trimmed_lines = (line for line in lines
if line.strip() and not line[0] == '#')
for line in trimmed_lines:
if line.strip().startswith('ATOMIC_POSITIONS'):
if positions is not None:
raise ValueError('Multiple ATOMIC_POSITIONS specified')
# Priority and behaviour tested with QE 5.3
if 'crystal_sg' in line.lower():
raise NotImplementedError('CRYSTAL_SG not implemented')
elif 'crystal' in line.lower():
cell = cell
elif 'bohr' in line.lower():
cell = np.identity(3) * units['Bohr']
elif 'angstrom' in line.lower():
cell = np.identity(3)
# elif 'alat' in line.lower():
# cell = np.identity(3) * alat
else:
if alat is None:
raise ValueError('Set lattice parameter in &SYSTEM for '
'alat coordinates')
# Always the default, will be DEPRECATED as mandatory
# in future
cell = np.identity(3) * alat
positions = []
for _dummy in range(n_atoms):
split_line = next(trimmed_lines).split()
# These can be fractions and other expressions
position = np.dot((infix_float(split_line[1]),
infix_float(split_line[2]),
infix_float(split_line[3])), cell)
if len(split_line) > 4:
force_mult = (float(split_line[4]),
float(split_line[5]),
float(split_line[6]))
else:
force_mult = None
positions.append((split_line[0], position, force_mult))
return positions
def get_atomic_species(lines, n_species):
"""Parse atomic species from ATOMIC_SPECIES card.
Parameters
----------
lines : list[str]
A list of lines containing the ATOMIC_POSITIONS card.
n_species : int
Expected number of atom types. Only this many lines will be parsed.
Returns
-------
species : list[(str, float, str)]
Raises
------
ValueError
Any problems parsing the data result in ValueError
"""
species = None
# no blanks or comment lines, can the consume n_atoms lines for positions
trimmed_lines = (line.strip() for line in lines
if line.strip() and not line.startswith('#'))
for line in trimmed_lines:
if line.startswith('ATOMIC_SPECIES'):
if species is not None:
raise ValueError('Multiple ATOMIC_SPECIES specified')
species = []
for _dummy in range(n_species):
label_weight_pseudo = next(trimmed_lines).split()
species.append((label_weight_pseudo[0],
float(label_weight_pseudo[1]),
label_weight_pseudo[2]))
return species
def get_cell_parameters(lines, alat=None):
"""Parse unit cell from CELL_PARAMETERS card.
Parameters
----------
lines : list[str]
A list with lines containing the CELL_PARAMETERS card.
alat : float | None
Unit of lattice vectors in Angstrom. Only used if the card is
given in units of alat. alat must be None if CELL_PARAMETERS card
is in Bohr or Angstrom. For output files, alat will be parsed from
the card header and used in preference to this value.
Returns
-------
cell : np.array | None
Cell parameters as a 3x3 array in Angstrom. If no cell is found
None will be returned instead.
cell_alat : float | None
If a value for alat is given in the card header, this is also
returned, otherwise this will be None.
Raises
------
ValueError
If CELL_PARAMETERS are given in units of bohr or angstrom
and alat is not
"""
cell = None
cell_alat = None
# no blanks or comment lines, can take three lines for cell
trimmed_lines = (line for line in lines
if line.strip() and not line[0] == '#')
for line in trimmed_lines:
if line.strip().startswith('CELL_PARAMETERS'):
if cell is not None:
# multiple definitions
raise ValueError('CELL_PARAMETERS specified multiple times')
# Priority and behaviour tested with QE 5.3
if 'bohr' in line.lower():
if alat is not None:
raise ValueError('Lattice parameters given in '
'&SYSTEM celldm/A and CELL_PARAMETERS '
'bohr')
cell_units = units['Bohr']
elif 'angstrom' in line.lower():
if alat is not None:
raise ValueError('Lattice parameters given in '
'&SYSTEM celldm/A and CELL_PARAMETERS '
'angstrom')
cell_units = 1.0
elif 'alat' in line.lower():
# Output file has (alat = value) (in Bohrs)
if '=' in line:
alat = float(line.strip(') \n').split()[-1]) * units['Bohr']
cell_alat = alat
elif alat is None:
raise ValueError('Lattice parameters must be set in '
'&SYSTEM for alat units')
cell_units = alat
elif alat is None:
# may be DEPRECATED in future
cell_units = units['Bohr']
else:
# may be DEPRECATED in future
cell_units = alat
# Grab the parameters; blank lines have been removed
cell = [[ffloat(x) for x in next(trimmed_lines).split()[:3]],
[ffloat(x) for x in next(trimmed_lines).split()[:3]],
[ffloat(x) for x in next(trimmed_lines).split()[:3]]]
cell = np.array(cell) * cell_units
return cell, cell_alat
def str_to_value(string):
"""Attempt to convert string into int, float (including fortran double),
or bool, in that order, otherwise return the string.
Valid (case-insensitive) bool values are: '.true.', '.t.', 'true'
and 't' (or false equivalents).
Parameters
----------
string : str
Test to parse for a datatype
Returns
-------
value : any
Parsed string as the most appropriate datatype of int, float,
bool or string.
"""
# Just an integer
try:
return int(string)
except ValueError:
pass
# Standard float
try:
return float(string)
except ValueError:
pass
# Fortran double
try:
return ffloat(string)
except ValueError:
pass
# possible bool, else just the raw string
if string.lower() in ('.true.', '.t.', 'true', 't'):
return True
elif string.lower() in ('.false.', '.f.', 'false', 'f'):
return False
else:
return string.strip("'")
def read_fortran_namelist(fileobj):
"""Takes a fortran-namelist formatted file and returns nested
dictionaries of sections and key-value data, followed by a list
of lines of text that do not fit the specifications.
Behaviour is taken from Quantum ESPRESSO 5.3. Parses fairly
convoluted files the same way that QE should, but may not get
all the MANDATORY rules and edge cases for very non-standard files:
Ignores anything after '!' in a namelist, split pairs on ','
to include multiple key=values on a line, read values on section
start and end lines, section terminating character, '/', can appear
anywhere on a line.
All of these are ignored if the value is in 'quotes'.
Parameters
----------
fileobj : file
An open file-like object.