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"""LAMMPS has the options to use several internal units (which can be different 

2from the ones used in ase). Mapping is therefore necessary. 

3 

4See: https://lammps.sandia.gov/doc/units.html 

5 """ 

6 

7from ase import units 

8from . import unitconvert_constants as u 

9 

10# !TODO add reduced Lennard-Jones units? 

11 

12# NOTE: We assume a three-dimensional simulation here! 

13DIM = 3.0 

14 

15UNITSETS = {} 

16 

17UNITSETS["ASE"] = dict( 

18 mass=1.0 / units.kg, 

19 distance=1.0 / units.m, 

20 time=1.0 / units.second, 

21 energy=1.0 / units.J, 

22 velocity=units.second / units.m, 

23 force=units.m / units.J, 

24 pressure=1.0 / units.Pascal, 

25 charge=1.0 / units.C, 

26) 

27 

28UNITSETS["real"] = dict( 

29 mass=u.gram_per_mole_si, 

30 distance=u.angstrom_si, 

31 time=u.femtosecond_si, 

32 energy=u.kcal_per_mole_si, 

33 velocity=u.angstrom_per_femtosecond_si, 

34 force=u.kcal_per_mole_angstrom_si, 

35 torque=u.kcal_per_mole_si, 

36 temperature=u.kelvin_si, 

37 pressure=u.atmosphere_si, 

38 dynamic_viscosity=u.poise_si, 

39 charge=u.e_si, 

40 dipole=u.electron_angstrom_si, 

41 electric_field=u.volt_per_angstrom_si, 

42 density=u.gram_si / u.centimeter_si ** DIM, 

43) 

44 

45UNITSETS["metal"] = dict( 

46 mass=u.gram_per_mole_si, 

47 distance=u.angstrom_si, 

48 time=u.picosecond_si, 

49 energy=u.ev_si, 

50 velocity=u.angstrom_per_picosecond_si, 

51 force=u.ev_per_angstrom_si, 

52 torque=u.ev_si, 

53 temperature=u.kelvin_si, 

54 pressure=u.bar_si, 

55 dynamic_viscosity=u.poise_si, 

56 charge=u.e_si, 

57 dipole=u.electron_angstrom_si, 

58 electric_field=u.volt_per_angstrom_si, 

59 density=u.gram_si / u.centimeter_si ** DIM, 

60) 

61 

62UNITSETS["si"] = dict( 

63 mass=u.kilogram_si, 

64 distance=u.meter_si, 

65 time=u.second_si, 

66 energy=u.joule_si, 

67 velocity=u.meter_per_second_si, 

68 force=u.newton_si, 

69 torque=u.joule_si, 

70 temperature=u.kelvin_si, 

71 pressure=u.pascal_si, 

72 dynamic_viscosity=u.pascal_si * u.second_si, 

73 charge=u.coulomb_si, 

74 dipole=u.coulomb_meter_si, 

75 electric_field=u.volt_per_meter_si, 

76 density=u.kilogram_si / u.meter_si ** DIM, 

77) 

78 

79UNITSETS["cgs"] = dict( 

80 mass=u.gram_si, 

81 distance=u.centimeter_si, 

82 time=u.second_si, 

83 energy=u.erg_si, 

84 velocity=u.centimeter_per_second_si, 

85 force=u.dyne_si, 

86 torque=u.dyne_centimeter_si, 

87 temperature=u.kelvin_si, 

88 pressure=u.dyne_per_centimetersq_si, # or barye =u. 1.0e-6 bars 

89 dynamic_viscosity=u.poise_si, 

90 charge=u.statcoulomb_si, # or esu (4.8032044e-10 is a proton) 

91 dipole=u.statcoulomb_centimeter_si, # =u. 10^18 debye, 

92 electric_field=u.statvolt_per_centimeter_si, # or dyne / esu 

93 density=u.gram_si / (u.centimeter_si ** DIM), 

94) 

95 

96UNITSETS["electron"] = dict( 

97 mass=u.amu_si, 

98 distance=u.bohr_si, 

99 time=u.femtosecond_si, 

100 energy=u.hartree_si, 

101 velocity=u.bohr_per_atu_si, 

102 force=u.hartree_per_bohr_si, 

103 temperature=u.kelvin_si, 

104 pressure=u.pascal_si, 

105 charge=u.e_si, # multiple of electron charge (1.0 is a proton) 

106 dipole=u.debye_si, 

107 electric_field=u.volt_per_centimeter_si, 

108) 

109 

110UNITSETS["micro"] = dict( 

111 mass=u.picogram_si, 

112 distance=u.micrometer_si, 

113 time=u.microsecond_si, 

114 energy=u.picogram_micrometersq_per_microsecondsq_si, 

115 velocity=u.micrometer_per_microsecond_si, 

116 force=u.picogram_micrometer_per_microsecondsq_si, 

117 torque=u.picogram_micrometersq_per_microsecondsq_si, 

118 temperature=u.kelvin_si, 

119 pressure=u.picogram_per_micrometer_microsecondsq_si, 

120 dynamic_viscosity=u.picogram_per_micrometer_microsecond_si, 

121 charge=u.picocoulomb_si, # (1.6021765e-7 is a proton), 

122 dipole=u.picocoulomb_micrometer_si, 

123 electric_field=u.volt_per_micrometer_si, 

124 density=u.picogram_si / (u.micrometer_si) ** DIM, 

125) 

126 

127UNITSETS["nano"] = dict( 

128 mass=u.attogram_si, 

129 distance=u.nanometer_si, 

130 time=u.nanosecond_si, 

131 energy=u.attogram_nanometersq_per_nanosecondsq_si, 

132 velocity=u.nanometer_per_nanosecond_si, 

133 force=u.attogram_nanometer_per_nanosecondsq_si, 

134 torque=u.attogram_nanometersq_per_nanosecondsq_si, 

135 temperature=u.kelvin_si, 

136 pressure=u.attogram_per_nanometer_nanosecondsq_si, 

137 dynamic_viscosity=u.attogram_per_nanometer_nanosecond_si, 

138 charge=u.e_si, # multiple of electron charge (1.0 is a proton) 

139 dipole=u.electron_nanometer_si, 

140 electric_field=u.volt_si / u.nanometer_si, 

141 density=u.attogram_si / u.nanometer_si ** DIM, 

142) 

143 

144 

145def convert(value, quantity, fromunits, tounits): 

146 """Convert units between LAMMPS and ASE. 

147 

148 :param value: converted value 

149 :param quantity: mass, distance, time, energy, velocity, force, torque, 

150 temperature, pressure, dynamic_viscosity, charge, dipole, 

151 electric_field or density 

152 :param fromunits: ASE, metal, real or other (see lammps docs). 

153 :param tounits: ASE, metal, real or other 

154 :returns: converted value 

155 :rtype: 

156 """ 

157 return UNITSETS[fromunits][quantity] / UNITSETS[tounits][quantity] * value