예제 #1
0
def gen_test_data(datadir: str, params_fd: dict, supercell):
    from gpaw.elph.electronphonon import ElectronPhononCoupling

    params_gs = copy.deepcopy(params_fd)

    atoms = Cluster(ase.build.bulk('C'))

    calc_gs = GPAW(**params_gs)
    atoms.calc = calc_gs
    atoms.get_potential_energy()
    atoms.calc.write("gs.gpw", mode="all")

    # Make sure the real space grid matches the original.
    # (basically we multiply the number of grid points in each dimension by the supercell dimension)
    params_fd['gpts'] = calc_gs.wfs.gd.N_c * list(supercell)
    if 'h' in params_fd:
        del params_fd['h']
    del calc_gs

    if world.rank == 0:
        os.makedirs(datadir, exist_ok=True)
    calc_fd = GPAW(**params_fd)
    elph = ElectronPhononCoupling(atoms,
                                  calc=calc_fd,
                                  supercell=supercell,
                                  calculate_forces=True,
                                  name=f'{datadir}/elph')
    elph.run()
    calc_fd.wfs.gd.comm.barrier()
    elph = ElectronPhononCoupling(atoms, calc=calc_fd, supercell=supercell)
    elph.set_lcao_calculator(calc_fd)
    elph.calculate_supercell_matrix(dump=1)
예제 #2
0
def main__brute_gpw(structure_path, supercell, log):
    from gpaw.elph.electronphonon import ElectronPhononCoupling
    from gpaw import GPAW

    calc = GPAW(structure_path)
    supercell_atoms = make_gpaw_supercell(calc, supercell, txt=log)

    # NOTE: confusingly, Elph wants primitive atoms, but a calc for the supercell
    elph = ElectronPhononCoupling(calc.atoms, calc=supercell_atoms.calc, supercell=supercell, calculate_forces=True)
    elph.run()
    supercell_atoms.calc.wfs.gd.comm.barrier()
    elph = ElectronPhononCoupling(calc.atoms, calc=supercell_atoms.calc, supercell=supercell)
    elph.set_lcao_calculator(supercell_atoms.calc)
    elph.calculate_supercell_matrix(dump=1)
    return
예제 #3
0
def elph_do_symmetry_expansion(supercell, calc, displacement_dist, phonon, disp_carts, disp_sites, supercell_atoms):
    from gpaw.elph.electronphonon import ElectronPhononCoupling

    natoms_prim = len(calc.get_atoms())
    disp_values = [read_elph_input(f'sym-{index}') for index in range(len(disp_sites))]

    # NOTE: phonon.symmetry includes pure translational symmetries of the supercell
    #       so we use an empty quotient group
    quotient_perms = np.array([np.arange(len(supercell_atoms))])
    super_lattice = supercell_atoms.get_cell()[...]
    super_symmetry = phonon.symmetry.get_symmetry_operations()
    oper_sfrac_rots = super_symmetry['rotations']
    oper_sfrac_trans = super_symmetry['translations']
    oper_cart_rots = np.array([super_lattice.T @ Rfrac @ np.linalg.inv(super_lattice).T for Rfrac in oper_sfrac_rots])
    oper_cart_trans = oper_sfrac_trans @ super_lattice
    oper_phonopy_coperms = phonon.symmetry.get_atomic_permutations()
    oper_phonopy_deperms = np.argsort(oper_phonopy_coperms, axis=1)

    # Convert permutations by composing the following three permutations:   into phonopy order, apply oper, back to ase order
    parprint('phonopy deperms:', oper_phonopy_deperms.shape)
    deperm_phonopy_to_ase = interop.get_deperm_from_phonopy_sc_to_ase_sc(natoms_prim, supercell)
    oper_deperms = [np.argsort(deperm_phonopy_to_ase)[deperm][deperm_phonopy_to_ase] for deperm in oper_phonopy_deperms]
    del oper_phonopy_coperms, oper_phonopy_deperms

    elphsym = symmetry.ElphGpawSymmetrySource.from_setups_and_ops(
        setups=supercell_atoms.calc.wfs.setups,
        lattice=super_lattice,
        oper_cart_rots=oper_cart_rots,
        oper_cart_trans=oper_cart_trans,
        oper_deperms=oper_deperms,
        )

    if world.rank == 0:
        full_derivatives = symmetry.expand_derivs_by_symmetry(
            disp_sites,       # disp -> atom
            disp_carts,       # disp -> 3-vec
            disp_values,      # disp -> T  (displaced value, optionally minus equilibrium value)
            elph_callbacks_2(supercell_atoms.calc.wfs, elphsym, supercell=supercell),        # how to work with T
            oper_cart_rots,   # oper -> 3x3
            oper_perms=oper_deperms,       # oper -> atom' -> atom
            quotient_perms=quotient_perms,
        )

        # NOTE: confusingly, Elph wants primitive atoms, but a calc for the supercell
        elph = ElectronPhononCoupling(calc.atoms, calc=supercell_atoms.calc, supercell=supercell)
        displaced_cell_index = elph.offset
        del elph  # that's all we needed it for

        eq_Vt, eq_dH, eq_forces = read_elph_input('eq')
        for a in range(natoms_prim):
            for c in range(3):
                delta_Vt, delta_dH, delta_forces = full_derivatives[natoms_prim * displaced_cell_index + a][c]
                for sign in [-1, +1]:
                    disp = interop.AseDisplacement(atom=a, axis=c, sign=sign)
                    disp_Vt = eq_Vt + sign * displacement_dist * delta_Vt
                    disp_dH = {k: eq_dH[k] + sign * displacement_dist * delta_dH[k] for k in eq_dH}
                    disp_forces = eq_forces + sign * displacement_dist * delta_forces
                    pickle.dump(disp_forces, paropen(f'phonons.{disp}.pckl', 'wb'), protocol=2)
                    pickle.dump((disp_Vt, disp_dH), paropen(f'elph.{disp}.pckl', 'wb'), protocol=2)
    world.barrier()
예제 #4
0
def elph_do_supercell_matrix(log, calc, supercell):
    from gpaw.elph.electronphonon import ElectronPhononCoupling

    # calculate_supercell_matrix breaks if parallelized over domains so parallelize over kpt instead
    # (note: it prints messages from all processes but it DOES run faster with more processes)
    supercell_atoms = GPAW('supercell.eq.gpw', txt=log, parallel={'domain': (1,1,1), 'band': 1, 'kpt': world.size}).get_atoms()

    elph = ElectronPhononCoupling(calc.atoms, supercell=supercell, calc=supercell_atoms.calc)
    elph.set_lcao_calculator(supercell_atoms.calc)
    # to initialize bfs.M_a
    ensure_gpaw_setups_initialized(supercell_atoms.calc, supercell_atoms)
    elph.calculate_supercell_matrix(dump=1)

    world.barrier()
예제 #5
0
def do_elph_symmetry(
    data_subdir: str,
    params_fd: dict,
    supercell,
    all_displacements: tp.Iterable[AseDisplacement],
    symmetry_type: tp.Optional[str],
):
    atoms = Cluster(ase.build.bulk('C'))

    # a supercell exactly like ElectronPhononCoupling makes
    supercell_atoms = atoms * supercell
    quotient_perms = list(
        interop.ase_repeat_translational_symmetry_perms(len(atoms), supercell))

    # Make sure the grid matches our calculations (we repeated the grid of the groundstate)
    params_fd = copy.deepcopy(params_fd)
    params_fd['gpts'] = GPAW('gs.gpw').wfs.gd.N_c * list(supercell)
    if 'h' in params_fd:
        del params_fd['h']

    wfs_with_sym = get_wfs_with_sym(params_fd=params_fd,
                                    supercell_atoms=supercell_atoms,
                                    symmetry_type=symmetry_type)
    calc_fd = GPAW(**params_fd)

    # GPAW displaces the center cell for some reason instead of the first cell
    elph = ElectronPhononCoupling(atoms,
                                  calc=calc_fd,
                                  supercell=supercell,
                                  calculate_forces=True)
    displaced_cell_index = elph.offset
    del elph  # just showing that we don't use these anymore
    del calc_fd

    get_displaced_index = lambda prim_atom: displaced_cell_index * len(
        atoms) + prim_atom

    all_displacements = list(all_displacements)
    disp_atoms = [get_displaced_index(disp.atom) for disp in all_displacements]
    disp_carts = [
        disp.cart_displacement(DISPLACEMENT_DIST) for disp in all_displacements
    ]
    disp_values = [
        read_elph_input(data_subdir, disp) for disp in all_displacements
    ]

    full_Vt = np.empty((len(supercell_atoms), 3) + disp_values[0][0].shape)
    full_dH = np.empty((len(supercell_atoms), 3), dtype=object)
    full_forces = np.empty((len(supercell_atoms), 3) + disp_values[0][2].shape)

    lattice = supercell_atoms.get_cell()[...]
    oper_cart_rots = interop.gpaw_op_scc_to_cart_rots(
        wfs_with_sym.kd.symmetry.op_scc, lattice)
    if world.rank == 0:
        full_values = symmetry.expand_derivs_by_symmetry(
            disp_atoms,  # disp -> atom
            disp_carts,  # disp -> 3-vec
            disp_values,  # disp -> T  (displaced value, optionally minus equilibrium value)
            elph_callbacks(wfs_with_sym, supercell),  # how to work with T
            oper_cart_rots,  # oper -> 3x3
            oper_perms=wfs_with_sym.kd.symmetry.a_sa,  # oper -> atom' -> atom
            quotient_perms=quotient_perms,
        )
        for a in range(len(full_values)):
            for c in range(3):
                full_Vt[a][c] = full_values[a][c][0]
                full_dH[a][c] = full_values[a][c][1]
                full_forces[a][c] = full_values[a][c][2]
    else:
        # FIXME
        # the symmetry part is meant to be done in serial but we should return back to
        # our original parallel state after it...
        pass

    return full_Vt, full_dH, full_forces
예제 #6
0
def get_elph_elements(atoms,
                      gpw_name,
                      calc_fd,
                      sc=(1, 1, 1),
                      basename=None,
                      phononname='phonon'):
    """
        Evaluates the dipole transition matrix elements

        Input
        ----------
        params_fd : Calculation parameters used for the phonon calculation
        sc (tuple): Supercell, default is (1,1,1) used for gamma phonons
        basename  : If you want give a specific name (gqklnn_{}.pckl)

        Output
        ----------
        gqklnn.pckl, the electron-phonon matrix elements
    """
    from ase.phonons import Phonons
    from gpaw.elph.electronphonon import ElectronPhononCoupling

    calc_gs = GPAW(gpw_name)
    world = calc_gs.wfs.world

    #calc_fd = GPAW(**params_fd)
    calc_gs.initialize_positions(atoms)
    kpts = calc_gs.get_ibz_k_points()
    nk = len(kpts)
    gamma_kpt = [[0, 0, 0]]
    nbands = calc_gs.wfs.bd.nbands
    qpts = gamma_kpt

    # calc_fd.get_potential_energy()  # XXX needed to initialize C_nM ??????

    # Phonon calculation, We'll read the forces from the elph.run function
    # This only looks at gamma point phonons
    ph = Phonons(atoms=atoms, name=phononname, supercell=sc)
    ph.read()
    frequencies, modes = ph.band_structure(qpts, modes=True)

    if world.rank == 0:
        print("Phonon frequencies are loaded.")

    # Find el-ph matrix in the LCAO basis
    elph = ElectronPhononCoupling(atoms, calc=None, supercell=sc)

    elph.set_lcao_calculator(calc_fd)
    elph.load_supercell_matrix()
    if world.rank == 0:
        print("Supercell matrix is loaded")

    # Non-root processes on GD comm seem to be missing kpoint data.
    assert calc_gs.wfs.gd.comm.size == 1, "domain parallelism not supported"  # not sure how to fix this, sorry

    gcomm = calc_gs.wfs.gd.comm
    kcomm = calc_gs.wfs.kd.comm
    if gcomm.rank == 0:
        # Find the bloch expansion coefficients
        c_kn = np.empty((nk, nbands, calc_gs.wfs.setups.nao), dtype=complex)
        for k in range(calc_gs.wfs.kd.nibzkpts):
            c_k = calc_gs.wfs.collect_array("C_nM", k, 0)
            if kcomm.rank == 0:
                c_kn[k] = c_k
        kcomm.broadcast(c_kn, 0)

        # And we finally find the electron-phonon coupling matrix elements!
        g_qklnn = elph.bloch_matrix(c_kn=c_kn,
                                    kpts=kpts,
                                    qpts=qpts,
                                    u_ql=modes)

    if world.rank == 0:
        print("Saving the elctron-phonon coupling matrix")
        np.save("gqklnn{}.npy".format(make_suffix(basename)),
                np.array(g_qklnn))
예제 #7
0
        'gamma': True
    },
    'txt': None,
    'basis': 'dzp',
    'symmetry': {
        'point_group': False
    },
    'xc': 'PBE'
}
elph_calc = GPAW(**parameters)
atoms.set_calculator(elph_calc)
atoms.get_potential_energy()
gamma_bands = elph_calc.wfs.kpt_u[0].C_nM

elph = ElectronPhononCoupling(atoms,
                              elph_calc,
                              supercell=supercell,
                              calculate_forces=True)
elph.run()

parameters['parallel'] = {'domain': 1}
elph_calc = GPAW(**parameters)
elph = ElectronPhononCoupling(atoms, calc=None, supercell=supercell)
elph.set_lcao_calculator(elph_calc)
elph.calculate_supercell_matrix(dump=1)

ph = Phonons(atoms=atoms, name='phonons', supercell=supercell, calc=None)
ph.read()
kpts = [[0, 0, 0]]
frequencies, modes = ph.band_structure(kpts, modes=True)

c_kn = np.array([[gamma_bands[0]]])