Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1"""Module to read and write atoms in xtl file format for the prismatic and 

2computem software. 

3 

4See https://prism-em.com/docs-inputs for an example of this format and the 

5documentation of prismatic. 

6 

7See https://sourceforge.net/projects/computem/ for the source code of the 

8computem software. 

9""" 

10 

11import numpy as np 

12 

13from ase.atoms import Atoms, symbols2numbers 

14from ase.utils import reader 

15from .utils import verify_cell_for_export, verify_dictionary 

16 

17 

18@reader 

19def read_prismatic(fd): 

20 r"""Import prismatic and computem xyz input file as an Atoms object. 

21 

22 Reads cell, atom positions, occupancies and Debye Waller factor. 

23 The occupancy values and the Debye Waller factors are obtained using the 

24 `get_array` method and the `occupancies` and `debye_waller_factors` keys, 

25 respectively. The root means square (RMS) values from the 

26 prismatic/computem xyz file are converted to Debye-Waller factors (B) in Ų 

27 by: 

28 

29 .. math:: 

30 

31 B = RMS^2 * 8\pi^2 

32 

33 """ 

34 # Read comment: 

35 fd.readline() 

36 

37 # Read unit cell parameters: 

38 cellpar = [float(i) for i in fd.readline().split()] 

39 

40 # Read all data at once 

41 # Use genfromtxt instead of loadtxt to skip last line 

42 read_data = np.genfromtxt(fname=fd, skip_footer=1) 

43 # Convert from RMS to Debye-Waller factor 

44 RMS = read_data[:, 5]**2 * 8 * np.pi**2 

45 

46 atoms = Atoms(symbols=read_data[:, 0], 

47 positions=read_data[:, 1:4], 

48 cell=cellpar, 

49 ) 

50 atoms.set_array('occupancies', read_data[:, 4]) 

51 atoms.set_array('debye_waller_factors', RMS) 

52 

53 return atoms 

54 

55 

56class XYZPrismaticWriter: 

57 """ See the docstring of the `write_prismatic` function. 

58 """ 

59 

60 def __init__(self, atoms, debye_waller_factors=None, comments=None): 

61 verify_cell_for_export(atoms.get_cell()) 

62 

63 self.atoms = atoms.copy() 

64 self.atom_types = set(atoms.symbols) 

65 self.comments = comments 

66 

67 self.occupancies = self._get_occupancies() 

68 debye_waller_factors = self._get_debye_waller_factors( 

69 debye_waller_factors) 

70 self.RMS = np.sqrt(debye_waller_factors / (8 * np.pi**2)) 

71 

72 def _get_occupancies(self): 

73 if 'occupancies' in self.atoms.arrays: 

74 occupancies = self.atoms.get_array('occupancies', copy=False) 

75 else: 

76 occupancies = np.ones_like(self.atoms.numbers) 

77 

78 return occupancies 

79 

80 def _get_debye_waller_factors(self, DW): 

81 if np.isscalar(DW): 

82 if len(self.atom_types) > 1: 

83 raise ValueError('This cell contains more then one type of ' 

84 'atoms and the Debye-Waller factor needs to ' 

85 'be provided for each atom using a ' 

86 'dictionary.') 

87 DW = np.ones_like(self.atoms.numbers) * DW 

88 elif isinstance(DW, dict): 

89 verify_dictionary(self.atoms, DW, 'DW') 

90 # Get the arrays of DW from mapping the DW defined by symbol 

91 DW = {symbols2numbers(k)[0]: v for k, v in DW.items()} 

92 DW = np.vectorize(DW.get)(self.atoms.numbers) 

93 else: 

94 for name in ['DW', 'debye_waller_factors']: 

95 if name in self.atoms.arrays: 

96 DW = self.atoms.get_array(name) 

97 

98 if DW is None: 

99 raise ValueError('Missing Debye-Waller factors. It can be ' 

100 'provided as a dictionary with symbols as key or ' 

101 'can be set for each atom by using the ' 

102 '`set_array("debye_waller_factors", values)` of ' 

103 'the `Atoms` object.') 

104 

105 return DW 

106 

107 def _get_file_header(self): 

108 # 1st line: comment line 

109 if self.comments is None: 

110 s = "{0} atoms with chemical formula: {1}.".format( 

111 len(self.atoms), 

112 self.atoms.get_chemical_formula()) 

113 else: 

114 s = self.comments 

115 

116 s = s.strip() 

117 s += " generated by the ase library.\n" 

118 # 2nd line: lattice parameter 

119 s += "{} {} {}".format( 

120 *self.atoms.cell.cellpar()[:3]) 

121 

122 return s 

123 

124 def write_to_file(self, f): 

125 data_array = np.vstack((self.atoms.numbers, 

126 self.atoms.positions.T, 

127 self.occupancies, 

128 self.RMS) 

129 ).T 

130 

131 np.savetxt(fname=f, 

132 X=data_array, 

133 fmt='%.6g', 

134 header=self._get_file_header(), 

135 newline='\n', 

136 footer='-1', 

137 comments='' 

138 ) 

139 

140 

141def write_prismatic(fd, *args, **kwargs): 

142 r"""Write xyz input file for the prismatic and computem software. The cell 

143 needs to be orthorhombric. If the cell contains the `occupancies` and 

144 `debye_waller_factors` arrays (see the `set_array` method to set them), 

145 these array will be written to the file. 

146 If the occupancies is not specified, the default value will be set to 0. 

147 

148 Parameters: 

149 

150 atoms: Atoms object 

151 

152 debye_waller_factors: float or dictionary of float or None (optional). 

153 If float, set this value to all 

154 atoms. If dictionary, each atom needs to specified with the symbol as 

155 key and the corresponding Debye-Waller factor as value. 

156 If None, the `debye_waller_factors` array of the Atoms object needs to 

157 be set (see the `set_array` method to set them), otherwise raise a 

158 ValueError. Since the prismatic/computem software use root means square 

159 (RMS) displacements, the Debye-Waller factors (B) needs to be provided 

160 in Ų and these values are converted to RMS displacement by: 

161 

162 .. math:: 

163 

164 RMS = \sqrt{\frac{B}{8\pi^2}} 

165 

166 Default is None. 

167 

168 comment: str (optional) 

169 Comments to be written in the first line of the file. If not 

170 provided, write the total number of atoms and the chemical formula. 

171 

172 """ 

173 

174 writer = XYZPrismaticWriter(*args, **kwargs) 

175 writer.write_to_file(fd)