Phonopy in pyiron

We will use the quasi-harmonic approximation (via PyIron’s implementation of the popular phonopy package) to evaluate look at thermal expansion and self-diffusion in Aluminium

[1]:
# Generic imports
from pyiron.project import Project
import numpy as np
%matplotlib  inline
import matplotlib.pylab as plt
import seaborn as sns
[2]:
pr = Project("PhonopyExample")
pot =  'Al_Mg_Mendelev_eam'
pr.remove_jobs(recursive=True)

Helper functions

Because repeating code is evil.

[3]:
def make_phonopy_job(template_job, name):
    """
    Create a phonopy job from a reference job.

    Args:
        template_job (pyiron job): The job to copy.
        name (str): What to call this new job.

    Returns:
        A new phonopy job.
    """
    project = template_job.project

    # What I want:
    # job_type =  template_job.job_type
    # What I have to do instead:
    job_type = pr.job_type.Lammps

    ref = project.create_job(job_type, name + "_ref")
    ref.structure = template_job.get_final_structure().copy()
    ref.potential = template_job.potential

    phono = project.create_job(pr.job_type.PhonopyJob, name)
    phono.ref_job = ref
    return phono
[4]:
def scale_structure(struct, scale):
    """
    Rescale the atomic positions and cell of a structure simultaneously.
    Accepts rescaling by an arbitrary real-valued 3x3 numpy array, but a float can be given
    for isotropic rescaling.

    Args:
        struct (Structure object): The structure to rescale.
        scale (float or np.array(3,3)): The matrix to rescale by. (float -> isotropic.)

    Returns:
        A rescaled copy of the structure.

    ..TODO: Double check that the scaling matrix still spans 3-space (determinant check?)
    """
    if isinstance(scale, float) or isinstance(scale, int):
        scale_mat = scale * np.eye(3)
    else:
        assert(scale.shape == (3,3))
        scale_mat = scale.T
    new_struct = struct.copy()
    new_struct.cell = np.dot(struct.cell, scale_mat)
    new_struct.positions = np.dot(struct.positions, scale_mat)
    return new_struct
[5]:
def scale_array(arr, scaler=None, new_range=1.):
    """
    Linearly transforms an array so that values equal to the minimum and maximum of the
    `scaler` array are mapped to the range (0, `new_range`). Note that rescaled values can
    still lie outside this range if the orignal values of `arr` are outside the bounds of
    `scaler`.

    Args:
        arr (np.array): Array to rescale.
        scaler (np.array): Array by which to rescale. Default is `arr`.
        new_range (float): New value for data which was the size `np.amax(scaler)`.
          Default is 1.
    """
    if scaler is None:
        scaler = arr
    return new_range * (arr - np.amin(scaler)) / np.ptp(scaler)

Thermal Expansion

What does the QHA say the lattice constant is as a function of temperature?

[6]:
pr_te = pr.create_group("ThermalExpansion")

Relax the unit cell

If we were doing VASP instead it would be important to do the least computation as possible, so here we’ll start by relaxing a simple unit cell to turn into a supercell later.

[7]:
job_unit = pr_te.create_job(pr.job_type.Lammps, "UnitCell")
[8]:
basis = pr_te.create_structure("Al", "fcc", 4.04)
[9]:
job_unit.structure = basis
job_unit.potential = pot
[10]:
job_unit.calc_minimize(pressure=0.0)
job_unit.run()
The job UnitCell was saved and received the ID: 3596380
[11]:
basis_rel = job_unit.get_final_structure()

Relax the bulk supercell

A relaxation which should take zero steps given our starting position!

[12]:
job_bulk_1 = pr_te.create_job(pr.job_type.Lammps, "Bulk_1")
# The _1 here refers to the fact that the volume has been rescaled by a factor of "1.0"
# (i.e. it hasn't been rescaled)
[13]:
n_reps = 3
job_bulk_1.structure = basis_rel.repeat(rep=n_reps)
job_bulk_1.potential = pot
[14]:
job_bulk_1.structure.plot3d();
[15]:
job_bulk_1.calc_minimize(pressure=0.0)
job_bulk_1.run()
The job Bulk_1 was saved and received the ID: 3596381

Calculate phonons

Run phonopy on the bulk supercell

[16]:
phono_bulk_1 = make_phonopy_job(job_bulk_1, "PhonoBulk_1")
[17]:
phono_bulk_1.run()
# Run performs a whole bunch of child calculations
# Each one has the positions slightly deformed in the symmetry-appropriate ways needed
# to get the phonon properties
The job PhonoBulk_1 was saved and received the ID: 3596382
The job supercell_phonon_0 was saved and received the ID: 3596383
[18]:
# Let's see what we got...
T_min = 0
T_max = 800 # a bit below melting
T_step = 25
temperatures = np.linspace(T_min, T_max, int((T_max - T_min) / T_step))
tp_bulk_1 = phono_bulk_1.get_thermal_properties(temperatures=temperatures)
# `get_thermal_properties` uses the displacements and forces to generate phonon information
[19]:
U_bulk_1 = job_bulk_1.output.energy_pot[-1]
Fvib_bulk_1 = tp_bulk_1.free_energies
plt.plot(temperatures, U_bulk_1 + Fvib_bulk_1)
plt.xlabel("Temperature [K]")
plt.ylabel("Free energy  ($U+F_{vib}$)  [eV]")
[19]:
Text(0, 0.5, 'Free energy  ($U+F_{vib}$)  [eV]')
../../_images/source_notebooks_phonopy_example_24_1.png

Calculate thermal expansivity

Above we have the (QHA approximation to the) free energy as a function of temperature at a fixed volume. To evaluate the thermal expansivity, we need to create the entire F(V,T) surface. To get this, we just loop over jobs like the above, but scaled to have different lattice constants.

[20]:
# According to Wikipedia, the thermal expansivity is about 0.0023% / Kelvin
# So at our maximum temperature, we expect around 1.8% expansion
scale_min = 0.995
scale_max = 1.02
scale_step = 0.002
scales = np.linspace(scale_min, scale_max, int((scale_max - scale_min) / scale_step))
[21]:
# Let's keep things clean by making another sub-directory
pr_scales = pr_te.create_group("ScanScales")
[22]:
# Loop the phonon calculation over all the volumes
sc_bulk_rel = job_bulk_1.get_final_structure()
bulk_free_energies = np.zeros((len(scales), len(temperatures)))

for i, scale in enumerate(scales):
    name_tail = "_{}".format(str(scale).replace(".", "_"))

    # Make a bulk job with the rescaled structure
    # (already relaxed, by symmetry won't change, calc static will be enough)
    job_bulk = pr_scales.create_job(pr.job_type.Lammps, "Bulk" + name_tail)
    job_bulk.potential = pot
    job_bulk.structure = scale_structure(sc_bulk_rel, scale)
    job_bulk.calc_static()
    job_bulk.run()
    U = job_bulk.output.energy_tot[-1]

    # Use that job as a reference for a phonopy job
    phono_bulk = make_phonopy_job(job_bulk, "PhonoBulk" + name_tail)
    phono_bulk.run()
    tp_bulk = phono_bulk.get_thermal_properties(temperatures=temperatures)
    Fvib = tp_bulk.free_energies

    # Fill in the row of free energies for this volume
    bulk_free_energies[i] = Fvib + U
The job Bulk_0_995 was saved and received the ID: 3596385
The job PhonoBulk_0_995 was saved and received the ID: 3596386
The job supercell_phonon_0 was saved and received the ID: 3596387
The job Bulk_0_9972727272727273 was saved and received the ID: 3596388
The job PhonoBulk_0_9972727272727273 was saved and received the ID: 3596389
The job supercell_phonon_0 was saved and received the ID: 3596390
The job Bulk_0_9995454545454545 was saved and received the ID: 3596394
The job PhonoBulk_0_9995454545454545 was saved and received the ID: 3596400
The job supercell_phonon_0 was saved and received the ID: 3596405
The job Bulk_1_0018181818181817 was saved and received the ID: 3596418
The job PhonoBulk_1_0018181818181817 was saved and received the ID: 3596419
The job supercell_phonon_0 was saved and received the ID: 3596420
The job Bulk_1_0040909090909091 was saved and received the ID: 3596434
The job PhonoBulk_1_0040909090909091 was saved and received the ID: 3596436
The job supercell_phonon_0 was saved and received the ID: 3596439
The job Bulk_1_0063636363636363 was saved and received the ID: 3596449
The job PhonoBulk_1_0063636363636363 was saved and received the ID: 3596450
The job supercell_phonon_0 was saved and received the ID: 3596451
The job Bulk_1_0086363636363636 was saved and received the ID: 3596456
The job PhonoBulk_1_0086363636363636 was saved and received the ID: 3596458
The job supercell_phonon_0 was saved and received the ID: 3596461
The job Bulk_1_010909090909091 was saved and received the ID: 3596473
The job PhonoBulk_1_010909090909091 was saved and received the ID: 3596475
The job supercell_phonon_0 was saved and received the ID: 3596478
The job Bulk_1_0131818181818182 was saved and received the ID: 3596488
The job PhonoBulk_1_0131818181818182 was saved and received the ID: 3596490
The job supercell_phonon_0 was saved and received the ID: 3596494
The job Bulk_1_0154545454545454 was saved and received the ID: 3596505
The job PhonoBulk_1_0154545454545454 was saved and received the ID: 3596507
The job supercell_phonon_0 was saved and received the ID: 3596510
The job Bulk_1_0177272727272728 was saved and received the ID: 3596521
The job PhonoBulk_1_0177272727272728 was saved and received the ID: 3596523
The job supercell_phonon_0 was saved and received the ID: 3596526
The job Bulk_1_02 was saved and received the ID: 3596537
The job PhonoBulk_1_02 was saved and received the ID: 3596539
The job supercell_phonon_0 was saved and received the ID: 3596542
[23]:
# The lattice constant is probably a more informative value than the 0K-relative strain
latts = basis_rel.cell[0][0] * scales
[24]:
# At each temperature, find the optimal volume by a simple quadratic fit
# ...Wait, which order fit will be good enough? Let's just spot-check
free_en = bulk_free_energies[:, -1]
plt.plot(latts, free_en, color='b', label='data')

# We'll plot the fit on a much denser mesh
fit_deg = 4
p = np.polyfit(latts, free_en, deg=fit_deg)
dense_latts = np.linspace(np.amin(latts), np.amax(latts), 1000)
#dense_latts = np.linspace(0, 10, 1000)
plt.plot(dense_latts, np.polyval(p=p, x=dense_latts), color='r', label='fit')
plt.xlabel('Lattice constant [$\mathrm{\AA}$]')
plt.ylabel('Bulk free energy [eV/supercell]')
plt.legend()
# Ok, a fourth-order fit seems perfectly reasonable
[24]:
<matplotlib.legend.Legend at 0x2b61a42bf278>
../../_images/source_notebooks_phonopy_example_30_1.png
[25]:
# Now find optimal temperatures
best_latts = np.zeros(len(temperatures))
best_latt_guess = basis_rel.cell[0][0]
for i, T in enumerate(temperatures):
    free_en = bulk_free_energies[:, i]
    p = np.polyfit(latts, free_en, deg=fit_deg)
    extrema = np.roots(np.polyder(p, m=1)).real  # Find where first-derivative is zero
    best_latts[i] = extrema[np.argmin(np.abs(extrema - best_latt_guess))]
[26]:
# Check that they're resonable
print(best_latt_guess, '\n', best_latts)
4.045270475668763
 [4.05946291 4.05949371 4.05987201 4.06083718 4.06226186 4.06396024
 4.06579637 4.06768361 4.06956868 4.0714193  4.07321635 4.07494901
 4.07661182 4.07820271 4.07972185 4.08117076 4.08255184 4.08386799
 4.08512237 4.08631821 4.08745877 4.0885472  4.08958656 4.09057975
 4.0915295  4.09243842 4.09330892 4.09414326 4.09494357 4.09571182
 4.09644984 4.09715935]
[27]:
# Let's look at the landscape
fig, ax = plt.subplots()
sns.heatmap(bulk_free_energies, ax=ax, cmap="coolwarm",
            xticklabels=['{:,.0f}'.format(T) for T in temperatures],
            yticklabels=['{:,.2f}'.format(a) for a in latts])
ax.set_xlabel("Temperature [K]")
ax.set_ylabel("Lattice constant [$\mathrm{\AA}$]")

# Overlaying the optimal path takes a couple changes of variables
# since the heatmap is plotting integer cells


ax.plot(scale_array(temperatures, new_range=len(temperatures)),
        scale_array(best_latts, scaler=latts, new_range=len(latts)),
        color='k')

[27]:
[<matplotlib.lines.Line2D at 0x2b61d9e8ee80>]
../../_images/source_notebooks_phonopy_example_33_1.png

Vacancies and diffusion

Another common use of QHA is to calculate the pre-factor for migration in a diffusion event.

In particular, the diffusion jump barrier looks like \(\omega_0 = \nu_0^\star \exp(-H_\mathrm{m} / k_\mathrm{B} T)\), where \(\nu_0^\star = \prod_{i=1}^{3N-3} \nu_i^\mathrm{IS} / \prod_{i=1}^{3N-4} \nu_i^\mathrm{TS}\), with IS and TS indicating the initial and transition states, respectively. Note that the transition state is missing a single frequency, which is from the instability of the transition state. It’s either an imaginary mode, which I think means a negative frequency. Meanwhile, \(H_\mathrm{m}\) is the enthalpic barrier (difference between the initial and transition states) and \(k_\mathrm{B} T\) is the usual thermal energy term.

Typically, these sorts of investigations use the nudged elastic band (NEB) to find the 0K transition state. You can do that with our new flexible jobs, but we’ll save that for later. For now we’ll just “approximate” the transition state with a simple linear interpolation.

Stable vacancy structures

Let’s start by generating and relaxing the initial and final states

[28]:
pr_vac = pr.create_group("Vacancies")
[29]:
# Find two adjacent sites
print(job_bulk_1.structure.positions[0])
print(job_bulk_1.structure.positions[1])
# Yep, 1 and 2 will do
[0. 0. 0.]
[2.02263524 2.02263524 0.        ]
[30]:
job_vac_i = pr_vac.create_job(pr.job_type.Lammps, "VacancyInitial")
job_vac_f = pr_vac.create_job(pr.job_type.Lammps, "VacancyFinal")

job_vac_i.potential = pot
job_vac_f.potential = pot
[31]:
sc_vac_i = sc_bulk_rel.copy()
sc_vac_i.pop(0)
job_vac_i.structure = sc_vac_i

sc_vac_f = sc_bulk_rel.copy()
sc_vac_f.pop(1)
job_vac_f.structure = sc_vac_f
[32]:
# Relax the new vacancy structures
job_vac_i.calc_minimize(pressure=0.0)
job_vac_i.run()

job_vac_f.calc_minimize(pressure=0.0)
job_vac_f.run()
The job VacancyInitial was saved and received the ID: 3596547
The job VacancyFinal was saved and received the ID: 3596549

DOS

The PyIron implementation of phonopy makes it very easy to look at the DOS. Let’s see what the effect is of introducing a vacancy, and confirm that our two vacancies are equivalent.

[33]:
phon_vac_i = make_phonopy_job(job_vac_i, "PhonoVacInitial")
phon_vac_f = make_phonopy_job(job_vac_f, "PhonoVacFinal")
[34]:
phon_vac_i.run()
tp_vac_i = phon_vac_i.get_thermal_properties(temperatures=temperatures)

phon_vac_f.run()
tp_vac_f = phon_vac_i.get_thermal_properties(temperatures=temperatures)

# Note that the vacancy structures spawn many more child processes
# This is because the vacancy structure has lower symmetry
The job PhonoVacInitial was saved and received the ID: 3596552
The job supercell_phonon_0 was saved and received the ID: 3596554
The job supercell_phonon_1 was saved and received the ID: 3596556
The job supercell_phonon_2 was saved and received the ID: 3596559
The job supercell_phonon_3 was saved and received the ID: 3596561
The job supercell_phonon_4 was saved and received the ID: 3596564
The job supercell_phonon_5 was saved and received the ID: 3596566
The job supercell_phonon_6 was saved and received the ID: 3596569
The job supercell_phonon_7 was saved and received the ID: 3596571
The job supercell_phonon_8 was saved and received the ID: 3596573
The job supercell_phonon_9 was saved and received the ID: 3596576
The job supercell_phonon_10 was saved and received the ID: 3596578
The job supercell_phonon_11 was saved and received the ID: 3596580
The job supercell_phonon_12 was saved and received the ID: 3596582
The job supercell_phonon_13 was saved and received the ID: 3596585
The job supercell_phonon_14 was saved and received the ID: 3596587
The job supercell_phonon_15 was saved and received the ID: 3596589
The job supercell_phonon_16 was saved and received the ID: 3596592
The job supercell_phonon_17 was saved and received the ID: 3596594
The job supercell_phonon_18 was saved and received the ID: 3596597
The job supercell_phonon_19 was saved and received the ID: 3596599
The job supercell_phonon_20 was saved and received the ID: 3596601
The job PhonoVacFinal was saved and received the ID: 3596613
The job supercell_phonon_0 was saved and received the ID: 3596616
The job supercell_phonon_1 was saved and received the ID: 3596618
The job supercell_phonon_2 was saved and received the ID: 3596620
The job supercell_phonon_3 was saved and received the ID: 3596623
The job supercell_phonon_4 was saved and received the ID: 3596625
The job supercell_phonon_5 was saved and received the ID: 3596628
The job supercell_phonon_6 was saved and received the ID: 3596630
The job supercell_phonon_7 was saved and received the ID: 3596632
The job supercell_phonon_8 was saved and received the ID: 3596635
The job supercell_phonon_9 was saved and received the ID: 3596637
The job supercell_phonon_10 was saved and received the ID: 3596640
The job supercell_phonon_11 was saved and received the ID: 3596642
The job supercell_phonon_12 was saved and received the ID: 3596644
The job supercell_phonon_13 was saved and received the ID: 3596646
The job supercell_phonon_14 was saved and received the ID: 3596649
The job supercell_phonon_15 was saved and received the ID: 3596651
The job supercell_phonon_16 was saved and received the ID: 3596653
The job supercell_phonon_17 was saved and received the ID: 3596655
The job supercell_phonon_18 was saved and received the ID: 3596658
The job supercell_phonon_19 was saved and received the ID: 3596659
The job supercell_phonon_20 was saved and received the ID: 3596660
[35]:
fig, ax = plt.subplots()
phono_bulk_1.plot_dos(ax=ax, color='b', label='bulk')
phon_vac_i.plot_dos(ax=ax, color='r', label='vac_i')
phon_vac_f.plot_dos(ax=ax, color='orange', label='vac_f')
plt.legend()

[35]:
<matplotlib.legend.Legend at 0x2b61da20bcc0>
../../_images/source_notebooks_phonopy_example_44_1.png

Attack frequency

Now we get the attack frequency by comparing the individual phonon spectra of initial and transition states

[36]:
# Interpolate initial and final positions to guesstimate the transition state
sc_vac_ts = sc_vac_i.copy()
sc_vac_ts.positions = 0.5 * (sc_vac_i.positions + sc_vac_f.positions)
[37]:
job_vac_ts = pr_vac.create_job(pr.job_type.Lammps, "VacancyTransition")
job_vac_ts.potential = pot
job_vac_ts.structure = sc_vac_ts

[38]:
# We _don't_ relax this job, or it would fall into the initial or final state!
job_vac_ts.calc_static()
job_vac_ts.run()
The job VacancyTransition was saved and received the ID: 3596670
[39]:
phon_vac_ts = make_phonopy_job(job_vac_ts, "PhonoVacTransition")
[40]:
phon_vac_ts.run()
tp_vac_ts = phon_vac_ts.get_thermal_properties(temperatures=temperatures)
The job PhonoVacTransition was saved and received the ID: 3596673
The job supercell_phonon_0 was saved and received the ID: 3596676
The job supercell_phonon_1 was saved and received the ID: 3596678
The job supercell_phonon_2 was saved and received the ID: 3596681
The job supercell_phonon_3 was saved and received the ID: 3596683
The job supercell_phonon_4 was saved and received the ID: 3596686
The job supercell_phonon_5 was saved and received the ID: 3596688
The job supercell_phonon_6 was saved and received the ID: 3596691
The job supercell_phonon_7 was saved and received the ID: 3596693
The job supercell_phonon_8 was saved and received the ID: 3596695
The job supercell_phonon_9 was saved and received the ID: 3596698
The job supercell_phonon_10 was saved and received the ID: 3596700
The job supercell_phonon_11 was saved and received the ID: 3596703
The job supercell_phonon_12 was saved and received the ID: 3596705
The job supercell_phonon_13 was saved and received the ID: 3596707
The job supercell_phonon_14 was saved and received the ID: 3596709
The job supercell_phonon_15 was saved and received the ID: 3596712
The job supercell_phonon_16 was saved and received the ID: 3596714
The job supercell_phonon_17 was saved and received the ID: 3596716
The job supercell_phonon_18 was saved and received the ID: 3596719
The job supercell_phonon_19 was saved and received the ID: 3596721
The job supercell_phonon_20 was saved and received the ID: 3596723
The job supercell_phonon_21 was saved and received the ID: 3596726
The job supercell_phonon_22 was saved and received the ID: 3596728
The job supercell_phonon_23 was saved and received the ID: 3596731
The job supercell_phonon_24 was saved and received the ID: 3596733
The job supercell_phonon_25 was saved and received the ID: 3596735
The job supercell_phonon_26 was saved and received the ID: 3596738
The job supercell_phonon_27 was saved and received the ID: 3596740
The job supercell_phonon_28 was saved and received the ID: 3596742
The job supercell_phonon_29 was saved and received the ID: 3596744
The job supercell_phonon_30 was saved and received the ID: 3596746
The job supercell_phonon_31 was saved and received the ID: 3596749
The job supercell_phonon_32 was saved and received the ID: 3596751
The job supercell_phonon_33 was saved and received the ID: 3596754
The job supercell_phonon_34 was saved and received the ID: 3596756
The job supercell_phonon_35 was saved and received the ID: 3596758
The job supercell_phonon_36 was saved and received the ID: 3596760
The job supercell_phonon_37 was saved and received the ID: 3596762
The job supercell_phonon_38 was saved and received the ID: 3596765
The job supercell_phonon_39 was saved and received the ID: 3596767
The job supercell_phonon_40 was saved and received the ID: 3596769
The job supercell_phonon_41 was saved and received the ID: 3596772
The job supercell_phonon_42 was saved and received the ID: 3596774
The job supercell_phonon_43 was saved and received the ID: 3596776
The job supercell_phonon_44 was saved and received the ID: 3596778
The job supercell_phonon_45 was saved and received the ID: 3596781
The job supercell_phonon_46 was saved and received the ID: 3596782
The job supercell_phonon_47 was saved and received the ID: 3596783
The job supercell_phonon_48 was saved and received the ID: 3596784
The job supercell_phonon_49 was saved and received the ID: 3596785
The job supercell_phonon_50 was saved and received the ID: 3596786
The job supercell_phonon_51 was saved and received the ID: 3596788
The job supercell_phonon_52 was saved and received the ID: 3596790
The job supercell_phonon_53 was saved and received the ID: 3596792
The job supercell_phonon_54 was saved and received the ID: 3596794
The job supercell_phonon_55 was saved and received the ID: 3596796
The job supercell_phonon_56 was saved and received the ID: 3596798
The job supercell_phonon_57 was saved and received the ID: 3596799
The job supercell_phonon_58 was saved and received the ID: 3596801
The job supercell_phonon_59 was saved and received the ID: 3596803
The job supercell_phonon_60 was saved and received the ID: 3596806
The job supercell_phonon_61 was saved and received the ID: 3596808
The job supercell_phonon_62 was saved and received the ID: 3596810
The job supercell_phonon_63 was saved and received the ID: 3596813
The job supercell_phonon_64 was saved and received the ID: 3596815
The job supercell_phonon_65 was saved and received the ID: 3596817
The job supercell_phonon_66 was saved and received the ID: 3596819
The job supercell_phonon_67 was saved and received the ID: 3596822
The job supercell_phonon_68 was saved and received the ID: 3596824
The job supercell_phonon_69 was saved and received the ID: 3596826
The job supercell_phonon_70 was saved and received the ID: 3596829
The job supercell_phonon_71 was saved and received the ID: 3596831
The job supercell_phonon_72 was saved and received the ID: 3596833
The job supercell_phonon_73 was saved and received the ID: 3596836
The job supercell_phonon_74 was saved and received the ID: 3596838
The job supercell_phonon_75 was saved and received the ID: 3596842
The job supercell_phonon_76 was saved and received the ID: 3596846
The job supercell_phonon_77 was saved and received the ID: 3596850
The job supercell_phonon_78 was saved and received the ID: 3596855
The job supercell_phonon_79 was saved and received the ID: 3596858
The job supercell_phonon_80 was saved and received the ID: 3596863
The job supercell_phonon_81 was saved and received the ID: 3596866
The job supercell_phonon_82 was saved and received the ID: 3596869
The job supercell_phonon_83 was saved and received the ID: 3596874
The job supercell_phonon_84 was saved and received the ID: 3596878
The job supercell_phonon_85 was saved and received the ID: 3596881
The job supercell_phonon_86 was saved and received the ID: 3596885
The job supercell_phonon_87 was saved and received the ID: 3596890
The job supercell_phonon_88 was saved and received the ID: 3596894
The job supercell_phonon_89 was saved and received the ID: 3596898
The job supercell_phonon_90 was saved and received the ID: 3596902
The job supercell_phonon_91 was saved and received the ID: 3596907
The job supercell_phonon_92 was saved and received the ID: 3596911
The job supercell_phonon_93 was saved and received the ID: 3596915
The job supercell_phonon_94 was saved and received the ID: 3596920
[41]:
# The transition state has an imaginary mode (frequency < 0), let's see it
fig, ax = plt.subplots()
phon_vac_i.plot_dos(ax=ax, color='r', label='initial')
phon_vac_ts.plot_dos(ax=ax, color='b', label='transition')
plt.legend()

[41]:
<matplotlib.legend.Legend at 0x2b61db5474e0>
../../_images/source_notebooks_phonopy_example_51_1.png

To calculate the attack frequency, we’ll ignore both the negative mode of the transition state (which we were warned about in the equation), as well as the three frequencies which correspond to rigid translation and are very near zero, and sometimes dip to be negative. Phonopy sorts the frequencies by magnitude, so we can just skip the first three and four for the initial and transition states, respectively. We take them at q=0.

[42]:
freq_i = phon_vac_i.phonopy.get_frequencies(0)[3:]
freq_ts = phon_vac_i.phonopy.get_frequencies(0)[4:]
[43]:
print(np.prod(freq_i))
6.870675816849329e+236

Recall: \(\nu_0^\star = \prod_{i=1}^{3N-3} \nu_i^\mathrm{IS} / \prod_{i=1}^{3N-4} \nu_i^\mathrm{TS}\)

[44]:
# Products are dangerous beasts, so we'll do a little numeric magic
nu = np.prod(freq_i[:-1] / freq_ts) * freq_i[-1]
print("Attack frequency is ", nu, "THz (10^-12 s)")
Attack frequency is  2.6826827779812032 THz (10^-12 s)

Mantina et al. (PRL 2008) report \(\nu = 19.3\) THz using DFT and NEB, so our linearly-interpolated “transition state” with EAM is actually not doing so poorly.

There are many more things you can do with phonopy, including looking directly at the force constants, the Hessian matrix, etc. But hopefully this is a useful starting point.

[ ]: