Skip to content

Commit

Permalink
Added the details on import with Phonopy an ASE in the user guide
Browse files Browse the repository at this point in the history
  • Loading branch information
mesonepigreco committed Jan 17, 2022
1 parent 694bcad commit d84c02b
Show file tree
Hide file tree
Showing 2 changed files with 76 additions and 4 deletions.
63 changes: 63 additions & 0 deletions UserGuide/gettingstarted.rst
Original file line number Diff line number Diff line change
Expand Up @@ -531,6 +531,69 @@ The structure is read from the first dynamical matrix (usually the gamma point).
An experimental interface to phonopy is under developement.


Load phonons from Phonopy and ASE
---------------------------------

Often, you may have to load a dynamical matrix from other sources than quantum-espresso.
The two most common file formats are the Phonopy and the ASE calculation.

For phonopy, we expect you generate a FORCE_CONSTANTS file and a phonopy.yaml file containing the information about the structure and the atoms in the supercell.

NOTE: Up to version 1.0 the Phonopy importer assumes the FORCE_CONSTANTS and the structure are written in Rydberg atomic units. This will probably change in the near future.

.. code::
import cellconstructor as CC
import cellconstructor.Phonons
# Load the dynamical matrix in phonopy format
dyn = CC.Phonons.Phonons()
dyn.load_phonopy("path/to/phonopy.yaml")
# Save as new_dyn in quantum espresso format
dyn.save_qe("new_dyn")
# Save in the phonopy format
dyn.save_phonopy("path/to/directory")
Optionally, it is possible to provide also the path to the FORCE_CONSTANTS file, which by default is assumed to be in the same directory as the phonopy.yaml.
This script can be employed to convert from Phonopy to quantum espresso file format.

In the same way, it is possible to generate the FORCE_CONSTANTS file with the method save_phonopy. Here, you just need to specify the directory.
Instead of the phonopy.yaml, the unitcell.in in quantum espresso format is generated. Phonopy is able to read this file.


The ASE calculator can be used to directly compute the dynamical matrix from finite displacements.
Once you have your ase.phonons.Phonons object, you can convert it into a CellConstructor Phonons object using the method get_dyn_from_ase_phonons.

.. code::
import ase
from ase.build import bulk
from ase.calculators.emt import EMT
from ase.phonons import Phonons
import cellconstructor as CC, cellconstructor.Phonons
# Here the code to get the ase.phonon.Phonons object
# Setup crystal and EMT calculator
atoms = bulk('Al', 'fcc', a=4.05)
# Phonon calculator with ASE
N = 6
ph = Phonons(atoms, EMT(), supercell=(N, N, N), delta=0.05)
ph.run()
# Read forces and assemble the dynamical matrix
ph.read(acoustic=True)
ph.clean()
# Convert into the cellconstructor and save in quantum espresso format
dyn = CC.get_dyn_from_ase_phonons(ph)
dyn.save_qe("my_dyn_computed_with_ase")
Generate a video of phonon vibrations
-------------------------------------
Expand Down
17 changes: 13 additions & 4 deletions cellconstructor/Phonons.py
Original file line number Diff line number Diff line change
Expand Up @@ -1240,7 +1240,12 @@ def save_qe(self, filename, full_name = False):
fp.write("*" * 75 + "\n")
fp.close()

def save_phononpy(self, supercell_size = (1,1,1)):
def save_phononpy(self, *args, **kwargs):
"Mapping to save_phonopy"
warnings.warn("Deprecated, use save_phonopy instead.")
self.save_phonopy(*args, **kwargs)

def save_phonopy(self, path = ".", supercell_size = None):
"""
EXPORT THE DYN IN THE PHONONPY FORMAT
=====================================
Expand All @@ -1254,16 +1259,20 @@ def save_phononpy(self, supercell_size = (1,1,1)):
Parameters
----------
path: str
Path to the directory in which the FORCE_CONSTANTS and unitcell.in files are created.
supercell_size : list of 3
The supercell that defines the dynamical matrix, note phononpy
works in the supercell.
works in the supercell. If none, it is inferred from the q points
"""
if supercell_size is None:
supercell_size = self.GetSupercell()

# Save it into the phononpy in the supercell
superdyn = self.GenerateSupercellDyn(supercell_size)
filename = "FORCE_CONSTANTS"
filename = os.path.join(path, "FORCE_CONSTANTS")

nat_sc = superdyn.structure.N_atoms
nat = self.structure.N_atoms
Expand Down Expand Up @@ -1314,7 +1323,7 @@ def save_phononpy(self, supercell_size = (1,1,1)):
lines.append("%s %16.8f %16.8f %16.8f\n" % (atm, cov_vect[0], cov_vect[1], cov_vect[2]))


f = open("unitcell.in", "w")
f = open(os.path.join(path, "unitcell.in"), "w")
f.writelines(lines)
f.close()

Expand Down

0 comments on commit d84c02b

Please sign in to comment.