The GpumdAtoms Object
The GpumdAtoms class is a subclass of the Atoms class from the popular Python package Atomic
Simulation Environment (ASE) but with additional properties specific to GPUMD. For example, it stores the global
cutoff and max number of neighbors, both of which are needed for an xyz.in file in GPUMD.
To use the GpumdAtoms class, we first need to create an ASE Atoms object. Here, we will use the
graphene_nanoribbon() method to create a graphene nanoribbon and directly
convert it to a GpumdAtoms object:
>>> from ase.build import graphene_nanoribbon
>>> from gpyumd.atoms import GpumdAtoms
>>> gnr = GpumdAtoms(graphene_nanoribbon(60, 36, type='armchair', sheet=True, vacuum=3.35/2, C_C=1.44))
>>> type(gnr)
gpyumd.atoms.GpumdAtoms
>>> gnr.center() # Make sure atoms are in the cell
>>> gnr.set_cutoff(2.1)
>>> gnr.set_max_neighbors(3)
>>> gnr
GpumdAtoms(symbols='C8640', pbc=[True, False, True], cell=[149.64918977395098, 3.35, 155.52])
Note the set_cutoff() and set_max_neighbors() methods
used to set those GPUMD-specific properties.
Grouping Atoms
GPUMD and gpyumd also support grouping atoms. If we want to group atoms, for example, to prepare an NEMD simulation,
we can use the group_by_position() method. Continuing with our graphene nanoribbon
example from above, the code would look like:
>>> lx, ly, lz = gnr.cell.lengths()
>>> split = np.array([lz/100] + [lz/5] + [lz/10]*6)-0.4
>>> split = [0] + list(np.cumsum(split)) + [lz] # Setting the boundaries for each group
>>> print("z-direction boundaries:", [round(z,2) for z in split])
z-direction boundaries: [0, 3.96, 90.83, 134.06, 177.29, 220.52, 263.76, 306.99, 350.22, 436.32]
>>> group_method, ncounts = gnr.group_by_position(split, direction='z') # GpumdAtoms-specific function
>>> ncounts
array([ 400, 8000, 4000, 4000, 4000, 4000, 4000, 4000, 8000])
Here, we defined atom nine groups that can be used to apply thermostats, fix atoms, and measure temperatures. We could
use group 0 (0 indexed) to fix atoms, thermostat groups 1 and 8, and use the remaining groups for temperature
calculations. The second return value for group_by_position(), which we showed here,
tells us how many atoms are in each groups.
We can also group by atomic symbols with group_by_symbol(), where each atomic symbol is
paired with a group in a dictionary. This can be useful when, for example, dealing with heterostructures. A
graphene-MoS2 heterostructure can use groups to differentiate each layer with:
group_method, ncounts = g_mos2_hetero.group_by_symbol({"C":0, "Mo": 1, "S": 1})
meaning that all of the atoms for graphene will be in group 0 and all for MoS2 will be in group 1.
Additional grouping methods will be added as needed.
Sorting Atoms
The GpumdAtoms class also supports sorting the atoms by atom type or by group. This simply
changes the order of the atoms when they are output to an xyz.in file. We can order our graphene-MoS2
structure from the previous section by type using the following:
::
g_mos2_hetero.sort_atoms(sort_key=”type”, order=[“Mo”, “S”, “C”])
In this case, the order parameter takes an ordered list of the atomic symbols in the
GpumdAtoms object.
We can order our graphene nanoribbon structure by group using the following: :: gnr.sort_atoms(sort_key=”group”, order=[0,1,8,2,3,4,5,6,7], group_method=0)
In this case, the order parameter takes an ordered list of integers that directly correspond to the group numbers
of the GpumdAtoms object. Note that we also have to specify the group method as a
GpumdAtoms object may have many different groups.
Warning: In the xyz.in file, the atoms that belong to a single potential must be grouped together although not in
any particular order. For example, if you use the REBO potential for “Mo” and “S” interactions and the Tersoff
potential for “C” interactions, then an order of [“Mo”, “C”, “S”] is not valid. If you use the
Simulation class to build your GPUMD simulation, the ordering is automatically corrected to be
valid.
File Output
The GpumdAtoms class can also directly write relevant GPUMD files. The most common file to write
is the xyz.in, which can be output by using the following:
::
gnr.write_gpumd()
The write_gpumd() method has the optional arguments use_velocity, which is dictates
whether or not the velocities stored by the GpumdAtoms object will be written to the xyz.in
file, gpumd_file, which specifies the file name (xyz.in by default), and directory, which notes where the file
should be output. All other GPUMD-specific properties will be automatically handled by the
GpumdAtoms object.
The other files that the GpumdAtoms class can write are basis.in and kpoints.in, which are
required when using the compute_phonon keyword. Their usage can best be shown through an example where we prepare Si
for a phonon dispersion calculation.
>>> from ase.lattice.cubic import Diamond
>>> from ase.build import bulk
>>> from gpyumd.atoms import GpumdAtoms
>>> a=5.434 # Lattice constant
>>> Si_UC = GpumdAtoms(bulk('Si', 'diamond', a=a))
>>> Si_UC.add_basis() # gpyumd method
>>> Si_UC
GpumdAtoms(symbols='Si2', pbc=True, cell=[[0.0, 2.717, 2.717], [2.717, 0.0, 2.717], [2.717, 2.717, 0.0]])
>>> Si = Si_UC.repeat([2,2,1]) # Create 8 atom diamond structure
>>> Si.set_cell([a, a, a])
>>> Si.wrap()
>>> Si = Si.repeat([2,2,2]) # Complete full supercell
>>> Si.set_max_neighbors(4) # gpyumd method
>>> Si.set_cutoff(3) # gpyumd method
>>> Si
GpumdAtoms(symbols='Si64', pbc=True, cell=[10.868, 10.868, 10.868])
>>> Si.write_basis() # create kpoints.in file
>>> linear_path, sym_points, labels = Si_UC.write_kpoints(path='GXKGL',npoints=400) # create kpoints.in file
The code above creats a Si unit cell and uses the gpyumd.atoms.GpumdAtoms.add_basis() method to assign a basis
index for each atom as is required for the basis.in file. The gpyumd.atoms.GpumdAtoms.repeat() method, which
is an overridden method from the Atoms class in ASE, is then used to create a supercell while keeping the relevant
GPUMD data. Finally, we write the basis.in file with gpyumd.atoms.GpumdAtoms.write_basis() and the
kpoints.in file with gpyumd.atoms.GpumdAtoms.write_kpoints().
List of all methods
- class gpyumd.atoms.GpumdAtoms(symbols=None, positions=None, numbers=None, tags=None, momenta=None, masses=None, magmoms=None, charges=None, scaled_positions=None, cell=None, pbc=None, celldisp=None, constraint=None, calculator=None, info=None, velocities=None)[source]
Bases:
AtomsA thin wrapper class around the ASE Atoms object. Stores additional information
—– Atoms Documentation —– Atoms object.
The Atoms object can represent an isolated molecule, or a periodically repeated structure. It has a unit cell and there may be periodic boundary conditions along any of the three unit cell axes. Information about the atoms (atomic numbers and position) is stored in ndarrays. Optionally, there can be information about tags, momenta, masses, magnetic moments and charges.
In order to calculate energies, forces and stresses, a calculator object has to attached to the atoms object.
Parameters:
- symbols: str (formula) or list of str
Can be a string formula, a list of symbols or a list of Atom objects. Examples: ‘H2O’, ‘COPt12’, [‘H’, ‘H’, ‘O’], [Atom(‘Ne’, (x, y, z)), …].
- positions: list of xyz-positions
Atomic positions. Anything that can be converted to an ndarray of shape (n, 3) will do: [(x1,y1,z1), (x2,y2,z2), …].
- scaled_positions: list of scaled-positions
Like positions, but given in units of the unit cell. Can not be set at the same time as positions.
- numbers: list of int
Atomic numbers (use only one of symbols/numbers).
- tags: list of int
Special purpose tags.
- momenta: list of xyz-momenta
Momenta for all atoms.
- masses: list of float
Atomic masses in atomic units.
- magmoms: list of float or list of xyz-values
Magnetic moments. Can be either a single value for each atom for collinear calculations or three numbers for each atom for non-collinear calculations.
- charges: list of float
Initial atomic charges.
- cell: 3x3 matrix or length 3 or 6 vector
Unit cell vectors. Can also be given as just three numbers for orthorhombic cells, or 6 numbers, where first three are lengths of unit cell vectors, and the other three are angles between them (in degrees), in following order: [len(a), len(b), len(c), angle(b,c), angle(a,c), angle(a,b)]. First vector will lie in x-direction, second in xy-plane, and the third one in z-positive subspace. Default value: [0, 0, 0].
- celldisp: Vector
Unit cell displacement vector. To visualize a displaced cell around the center of mass of a Systems of atoms. Default value = (0,0,0)
- pbc: one or three bool
Periodic boundary conditions flags. Examples: True, False, 0, 1, (1, 1, 0), (True, False, False). Default value: False.
- constraint: constraint object(s)
Used for applying one or more constraints during structure optimization.
- calculator: calculator object
Used to attach a calculator for calculating energies and atomic forces.
- info: dict of key-value pairs
Dictionary of key-value pairs with additional information about the system. The following keys may be used by ase:
spacegroup: Spacegroup instance
unit_cell: ‘conventional’ | ‘primitive’ | int | 3 ints
adsorbate_info: Information about special adsorption sites
Items in the info attribute survives copy and slicing and can be stored in and retrieved from trajectory files given that the key is a string, the value is JSON-compatible and, if the value is a user-defined object, its base class is importable. One should not make any assumptions about the existence of keys.
- add_basis(index: Optional[List[int]] = None, mapping: Optional[List[int]] = None) None[source]
Assigns a basis index for each atom in atoms. Updates atoms.
https://github.com/brucefan1983/GPUMD/tree/master/examples/empirical_potentials/phonon_dispersion https://gpumd.zheyongfan.org/index.php/The_basis.in_input_file
- Parameters
index – Atom indices of those in the unit cell.
mapping – Mapping of all atoms to the relevant basis positions
- add_group_method(group: GroupMethod) int[source]
Add a grouping method to the GpumdAtoms object.
- Parameters
group – The group method to add
- Returns
Index of grouping method
- group_by_position(split: List[float], direction: str) Tuple[int, ndarray][source]
Assigns groups to all atoms based on its position. Only works in one direction as it is used for NEMD. Returns a bookkeeping parameter, but atoms will be udated in-place.
- Parameters
split – List of boundaries in ascending order. First element should be lower boundary of simulation box in specified direction and the last the upper.
direction – Which direction the split will work
- Returns
(index of the grouping method, number of atoms in each group)
- group_by_symbol(symbols: dict) Tuple[int, ndarray][source]
Assigns groups to all atoms based on atom symbols. Returns a bookkeeping parameter, but atoms will be udated in-place.
- Parameters
symbols – Dictionary with symbols for keys and group as a value. Only one group allowed per atom. Assumed groups are integers starting at 0 and increasing in steps of 1.
- Returns
(index of the grouping method, number of atoms in each group)
- group_methods: List[GroupMethod]
- repeat(rep: Union[int, List[int]]) GpumdAtoms[source]
A wrapper of ase.Atoms.repeat that is aware of GPUMD’s basis information.
- Parameters
rep – List of three positive integers or a single integer
- Returns
New repeated GpumdAtoms
- set_cutoff(cutoff: float) None[source]
Sets the global cutoff for the GpumdAtoms object.
- Parameters
cutoff – The cutoff to use in the xyz.in file
- set_max_neighbors(max_neighbors: int) None[source]
Set the maximum size of the neighbor list.
- Parameters
max_neighbors – Maximum number of neighbors
- set_type_dict(type_dict: Dict[str, int], overwrite: bool = False) None[source]
Assigns atomic symbols to GPUMD types
- Parameters
type_dict – Atomic symbol keys and type values
overwrite – To completely change existing type_dict
- sort_atoms(sort_key: Optional[str] = None, order: Optional[Union[List[str], List[int]]] = None, group_method: Optional[int] = None) None[source]
Sorts the atoms according to a specified order. If not order is specified, an ascending order based on the type_dict or group method will be used.
- Parameters
sort_key – How to sort atoms (‘group’, ‘type’)
order – For sort_key==’type’, a list of atomic symbol strings in the desired order. Ex: [“Mo”, “S”, “Si”, “O”]. For sort_key==’group’, a list of ints in desired order for groups at group_index. Ex: [1,3,2,4]
group_method – Selects the group to sort in the output.
- write_basis(filename: str = 'basis.in', directory: Optional[str] = None) None[source]
Creates the basis.in file. Atoms passed to this must already have the basis of every atom defined.
Related: atoms.add_basis, atoms.repeat
- Parameters
filename – File to save the structure data to
directory – Directory to store output
- write_gpumd(use_velocity: bool = False, gpumd_file: str = 'xyz.in', directory: Optional[str] = None) None[source]
Creates and xyz.in file.
- Parameters
use_velocity – Whether or not to set the velocities in the xyz.in file.
gpumd_file – File to save the structure data to
directory – Directory to store output
- write_kpoints(path: str = 'G', npoints: int = 1, special_points: Optional[Mapping[str, Sequence[float]]] = None, filename: str = 'kpoints.in', directory: Optional[str] = None) Tuple[ndarray, None, List[str]][source]
Creates the file “kpoints.in”, which specifies the kpoints needed for the ‘phonon’ keyword
- Parameters
path – String of special point names defining the path, e.g. ‘GXL’
npoints – Number of points in total. Note that at least one point is added for each special point in the path
special_points – Map of special points to scaled kpoint coordinates. For example
{'G': [0, 0, 0], 'X': [1, 0, 0]}filename – File to save the structure data to
directory – Directory to store output
- Returns
- kpoints converted to x-coordinates, x-coordinates of the high
symmetry points, labels of those points.
- class gpyumd.atoms.GroupByPosition(split: List[float], direction: str)[source]
Bases:
GroupMethodStores grouping information for a GpumdAtoms object
- get_group(position: List[float]) int[source]
Gets the group that an atom belongs to based on its position. Only works in one direction as it is used for NEMD.
- Parameters
position – list of floats of length 3 Position of the atom
- Returns
Group of atom
- Return type
int
- update(atoms: GpumdAtoms, order: Optional[List[int]] = None)[source]
- class gpyumd.atoms.GroupBySymbol(symbols: Dict[str, int])[source]
Bases:
GroupMethodStores grouping information for a GpumdAtoms object
- update(atoms: GpumdAtoms, order: Optional[List[int]] = None)[source]
- class gpyumd.atoms.GroupGeneric(groups: List[int])[source]
Bases:
GroupMethodGrouping with no specific guidelines. Mostly used for loaded xyz.in files.
- update(atoms: GpumdAtoms, order: Optional[List[int]] = None)[source]
- class gpyumd.atoms.GroupMethod(group_type: Optional[str] = None)[source]
Bases:
ABCStores grouping information for a GpumdAtoms object
- abstract update(atoms: GpumdAtoms, order: Optional[List[int]] = None)[source]