Source code for pyscf.pbc.dft.uks

#!/usr/bin/env python
# Copyright 2014-2019 The PySCF Developers. All Rights Reserved.
#
# Licensed under the Apache License, Version 2.0 (the "License");
# you may not use this file except in compliance with the License.
# You may obtain a copy of the License at
#
#     http://www.apache.org/licenses/LICENSE-2.0
#
# Unless required by applicable law or agreed to in writing, software
# distributed under the License is distributed on an "AS IS" BASIS,
# WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
# See the License for the specific language governing permissions and
# limitations under the License.
#
# Author: Qiming Sun <osirpt.sun@gmail.com>
#

'''
Non-relativistic unrestricted Kohn-Sham for periodic systems at a single k-point

See Also:
    pyscf.pbc.dft.krks.py : Non-relativistic Restricted Kohn-Sham for periodic
                            systems with k-point sampling
'''


import numpy
import pyscf.dft
from pyscf import lib
from pyscf.pbc.scf import uhf as pbcuhf
from pyscf.lib import logger
from pyscf.dft import uks as mol_uks
from pyscf.pbc.dft import gen_grid
from pyscf.pbc.dft import rks
from pyscf.pbc.dft import multigrid
from pyscf import __config__

get_rho = rks.get_rho

[docs] def get_veff(ks, cell=None, dm=None, dm_last=0, vhf_last=0, hermi=1, kpt=None, kpts_band=None): '''Coulomb + XC functional for UKS. See pyscf/pbc/dft/uks.py :func:`get_veff` fore more details. ''' if cell is None: cell = ks.cell if dm is None: dm = ks.make_rdm1() if kpt is None: kpt = ks.kpt t0 = (logger.process_clock(), logger.perf_counter()) ni = ks._numint hybrid = ni.libxc.is_hybrid_xc(ks.xc) if not hybrid and isinstance(ks.with_df, multigrid.MultiGridFFTDF): if ks.nlc or ni.libxc.is_nlc(ks.xc): raise NotImplementedError(f'MultiGrid for NLC functional {ks.xc} + {ks.nlc}') n, exc, vxc = multigrid.nr_uks(ks.with_df, ks.xc, dm, hermi, kpt.reshape(1,3), kpts_band, with_j=True, return_j=False) logger.debug(ks, 'nelec by numeric integration = %s', n) t0 = logger.timer(ks, 'vxc', *t0) return vxc if not isinstance(dm, numpy.ndarray): dm = numpy.asarray(dm) if dm.ndim == 2: # RHF DM dm = numpy.asarray((dm*.5,dm*.5)) # ndim = 3 : dm.shape = ([alpha,beta], nao, nao) ground_state = (dm.ndim == 3 and dm.shape[0] == 2 and kpts_band is None) ks.initialize_grids(cell, dm, kpt, ground_state) if hermi == 2: # because rho = 0 n, exc, vxc = (0,0), 0, 0 else: max_memory = ks.max_memory - lib.current_memory()[0] n, exc, vxc = ni.nr_uks(cell, ks.grids, ks.xc, dm, 0, hermi, kpt, kpts_band, max_memory=max_memory) if ks.nlc or ni.libxc.is_nlc(ks.xc): if ni.libxc.is_nlc(ks.xc): xc = ks.xc else: assert ni.libxc.is_nlc(ks.nlc) xc = ks.nlc n, enlc, vnlc = ni.nr_nlc_vxc(cell, ks.nlcgrids, xc, dm[0]+dm[1], 0, hermi, kpt, max_memory=max_memory) exc += enlc vxc += vnlc logger.debug(ks, 'nelec by numeric integration = %s', n) t0 = logger.timer(ks, 'vxc', *t0) if not hybrid: vj = ks.get_j(cell, dm[0]+dm[1], hermi, kpt, kpts_band) vxc += vj else: omega, alpha, hyb = ni.rsh_and_hybrid_coeff(ks.xc, spin=cell.spin) vj, vk = ks.get_jk(cell, dm, hermi, kpt, kpts_band) vj = vj[0] + vj[1] vk *= hyb if omega != 0: vklr = ks.get_k(cell, dm, hermi, kpt, kpts_band, omega=omega) vklr *= (alpha - hyb) vk += vklr vxc += vj - vk if ground_state: exc -=(numpy.einsum('ij,ji', dm[0], vk[0]) + numpy.einsum('ij,ji', dm[1], vk[1])).real * .5 if ground_state: ecoul = numpy.einsum('ij,ji', dm[0]+dm[1], vj).real * .5 else: ecoul = None vxc = lib.tag_array(vxc, ecoul=ecoul, exc=exc, vj=None, vk=None) return vxc
[docs] class UKS(rks.KohnShamDFT, pbcuhf.UHF): '''UKS class adapted for PBCs. This is a literal duplication of the molecular UKS class with some `mol` variables replaced by `cell`. ''' get_rho = get_rho get_vsap = mol_uks.UKS.get_vsap init_guess_by_vsap = mol_uks.UKS.init_guess_by_vsap get_veff = get_veff energy_elec = pyscf.dft.uks.energy_elec def __init__(self, cell, kpt=numpy.zeros(3), xc='LDA,VWN', exxdiv=getattr(__config__, 'pbc_scf_SCF_exxdiv', 'ewald')): pbcuhf.UHF.__init__(self, cell, kpt, exxdiv=exxdiv) rks.KohnShamDFT.__init__(self, xc)
[docs] def dump_flags(self, verbose=None): pbcuhf.UHF.dump_flags(self, verbose) rks.KohnShamDFT.dump_flags(self, verbose) return self
[docs] def to_hf(self): '''Convert to UHF object.''' from pyscf.pbc import scf return self._transfer_attrs_(scf.UHF(self.cell, self.kpt))
if __name__ == '__main__': from pyscf.pbc import gto cell = gto.Cell() cell.unit = 'A' cell.atom = 'C 0., 0., 0.; C 0.8917, 0.8917, 0.8917' cell.a = '''0. 1.7834 1.7834 1.7834 0. 1.7834 1.7834 1.7834 0. ''' cell.basis = 'gth-szv' cell.pseudo = 'gth-pade' cell.verbose = 7 cell.output = '/dev/null' cell.build() mf = UKS(cell) print(mf.kernel())