一个用于生成高效率通用 Monkhorst-Pack K 点网格的库,可加速电子结构计算,如 DFT。
项目描述
KpLib
KpLib 是一个用于寻找最优通用 Monkhorst-Pack K 点网格的 C++ 库。它可以作为高效通用 k 点网格的生成器导入到电子结构软件包中,或通过 Python 接口 kpLib
集成到用户脚本中。
有关 KpLib 和底层算法的问题,欢迎您查看我们的论文在此或发送电子邮件至 kpoints@jhu.edu。
注意:我们在 v.1.1.0
(2023.04.08)版本中修复了Python接口的一些错误,这些错误导致脚本在 $hR$ 表示法中的单斜晶系和三角晶系中挂起。我们还将外部Python包的依赖降到最低,以便用户可以选择他们喜欢的结构IO包。请查看 RELEASE.md 文件以获取详细信息。
使用方法
路线 I:将 KpLib 作为 C++ 库集成到模拟软件包中
编译 KpLib 为库
我们使用 cmake
来检测本地构建环境和生成本地构建文件。对于类Unix操作系统,用户可以通过以下步骤构建项目:
$ git clone https://gitlab.com/muellergroup/kplib.git
$ cd kplib
$ mkdir build
$ cd build
$ cmake ..
$ make
然后你可以在 build/libkpoints.a
找到静态库,在 build/libkpoints.so
找到动态库。
在 C++ 代码中使用 KpLib
基本上有两个步骤
-
将头文件
src/kPointLatticeGenerator.h
复制到你的include
文件夹,并在你的源代码中添加以下行:#include "kPointLatticeGenerator.h"
-
在链接阶段链接库。
例如,要将静态库链接到对象
myapp.o
并编译成最终的可执行文件myapp
,你可以:$ g++ myapp.o -L /path/to/lib libkpoints.a -o myapp
如果你想使用动态库,你可以:
$ export LD_LIBRARY_PATH=$LD_LIBRARY_PATH:/path/to/lib $ g++ myapp.o -L /path/to/lib -lkpoints -o myapp
第一行告诉加载器
ld
在运行时在哪里找到共享库,因为动态链接只在可执行文件中放置库的引用。
然后你就可以开始了!
下面将详细介绍API,同时还有一个用于演示在C++环境中使用此API的示例C++程序(demo_kplib/src
)以及API假定的约定。
路线 II:将 KpLib 作为 Python 包使用
安装
使用 pip
从 PyPI 安装
这将安装核心 KpLib 模块,第三方依赖最小化。用户将需要使用他们喜欢的材料科学包(例如 ase
和 pymatgen
)解析输入结构并输出生成的 K 点网格。
$ pip install --user pybind11 build wheel
$ pip install --user kplib # Install only the core kplib package.
$ pip install --user kplib[cli] # To install the CLI tool `kpgen`.
# This will also install `pymatgen` and `click`
使用 pip
从源码安装
Python 接口是 C++ 库的轻量级包装。从源码构建应该是直截了当的。它还提供了一种在终端中调用 CLI 命令的方法(需要 pymatgen
进行 IO)。
第一步。从 https://gitlab.com/muellergroup/kplib 下载源代码。
$ git clone https://gitlab.com/muellergroup/kplib.git
第二步。C++-Python 接口使用 pybind11
编写,而 setup.py
在开始时导入它以将 C++ 代码编译成一个 Python 外部模块。因此,在安装 kpLib
之前应该先安装它。
$ pip install --user pybind11
$ pip install --user build wheel # Python build toolchain that you might already have.
第三步。如果你想选择自己的结构解析器(如 pymatgen
、ase
和其他优秀的包),你可以仅安装 kpLib
的核心模块和最小依赖项。
$ cd kplib
$ pip install --user .
如果你想使用 CLI 命令工具 kpgen
,可以使用以下代码安装 kpLib
。这将额外安装 pymatgen
(作为解析器)和 click
(提供友好的 CLI 界面)。
$ cd kplib
$ pip install --user .[cli]
$ which kpgen
如果你想测试你的安装,可以使用启用测试套件的依赖项安装 kpLib
(pymatgen
和 pytest
)。安装后,你可以运行测试以确保一切安装正确。
$ cd kplib
$ pip install --user .[tests]
$ pytest
在 Python 脚本中使用 KpLib
《kpLib》Python库以提供核心KpLib函数的接口并允许最大化的使用灵活性为编写方式。解析输入结构的包和将生成的k点网格写入所需格式的选择留给了用户。以下部分将首先描述《get_kpoints》函数的签名,然后给出两个使用两个流行的包《pymatgen》和《ase》的示例脚本,用于I/O。
《kpLib.interface::get_kpoints》是《kpLib》Python库提供的核心和唯一函数。其签名如下
def get_kpoints(
lattice: List[List[float]] or numpy.ndarray,
fractional_coords: List[List[float]] or numpy.ndarray,
atomic_numbers: List[int] or Tuple[int],
min_distance: float = 0.1,
min_total_kpoints: int = 1,
include_gamma: str = "auto",
symprec: float = 1e-5,
use_scale_factor: bool = False,
) -> Dict
《get_kpoints》使用3个位置参数来定义一个结构。其余的关键字参数定义了输出K点网格的要求。K点网格的大小要求由min_distance
和min_total_kpoints
定义。它们可以单独使用,也可以一起使用。我们推荐使用min_distance
。有关将传统的通过三个整数m1 x m2 x m3
指定的Monkhorst-Pack网格与min_distance
关联的示例,请参阅下面的部分。
所有参数的含义如下
lattice
:一个3x3矩阵,表示输入结构的单位细胞的晶格向量。每一行表示一个向量,即格式为[a,b,c]。它可以是列表的列表或一个numpy.ndarray
。fractional_coords
:一个分数格式的三维坐标列表,相对于晶格向量。它可以是列表的列表或一个numpy.ndarray
。atomic_numbers
:一个原子序数列表,其中第i个数定义了位于fractional_coords[i]
位置的i第原子的元素。它可以是列表或元组。min_distance
:返回的K点网格对应的超晶格中任意两个晶格点在实空间中的最小所需距离。它指定了所需K点网格的大小,具有埃格单位。min_distance
的值越大,返回的网格越密集(即,更多对称不同的K点)。min_total_kpoints
:返回的K点网格中所需的最小总K点数。注意,这不是在计算中实际使用的对称不同的K点的数量。所有的K点集合将被对称地减少到最小的一组不同的K点,这些点将用于计算。include_gamma
:它可以有true
或false
或auto
的值。使用true
将Gamma点(即,倒空间中的[0, 0, 0])包括在返回的网格中。使用false
排除Gamma点。使用auto
让kplib在这两种类型的网格之间选择,并返回具有较少对称不同K点的网格。如果它们具有相同数量的不同K点,则返回具有更大的有效min_distance
的网格。symprec
:在确定结构对称性时实空间的精度。默认值为1E-5
。差异小于symprec
的原子坐标将被视为等效。如果您想包含更多的对称性,请使用更大的值。use_scale_factor
:一种以牺牲返回次优K点网格为代价来加速搜索的方案。默认行为是不使用该方案。当代码返回网格需要非常长的时间时,可以考虑使用此方案。有关详细说明,请参阅此README的末尾。
《get_kpoints》函数返回一个K点网格作为dict
。键是
min_periodic_distance
:返回的K点网格对应的实空间超晶格中任意两个晶格点之间的最小距离。num_distinct_kpts
:对称不同的K-点的数量。这是在计算中实际使用的K-点的数量。distinct_coords
:对称不同的K-点的分数坐标。distinct_weights
:在distinct_coords
中对应于对称不同的K-点的权重。num_total_kpts
:总K-点的数量。coords
:所有K-点的分数坐标。weights
:在coords
中对应于K-点的权重。
通常,对称不同的K-点的坐标和权重足以定义一个计算的输入文件。
理解"min_distance
"的含义
"最小距离"是指真实空间晶格上任意一对晶格点之间的最小周期距离。原始单元胞的倒格是真实空间超单元胞的倒格的超胞。这意味着每个K-点网格对应于真实空间的超晶格。因此,我们可以使用真实空间超晶格的"最小距离"来指定其对应的倒空间中K-点网格的大小。min_distance
越大,超胞就越大,因此K-点网格越密集。
例如,考虑一个具有晶格$[[2, 0, 0], [0, 2, 0], [0, 0, 3]]$的人工四方结构,以及位于原点的单个原子。其空间群为$P4/mmm$。带有Γ点的传统Monkhorst-Pack 4x4x4网格将对应于$[[8, 0, 0], [0, 8, 0], [0, 0, 12]]$的超胞,并且有18个对称不同的K-点。在这种情况下,最小周期距离是8埃。
当使用kpLib
且min_distance
=8时,以Γ为中心的结果将只有9个对称不同的K-点和一个有效最小距离为8.485。对于平移网格,结果将是一个具有6个对称不同的K-点和一个有效最小距离为8的网格。计算精度相似,但不同的K-点数量是原来的三分之一。
Python脚本示例
在这里,我们给出了两个在Python脚本中使用kpLib
的示例,分别使用流行的pymatgen
和ase
包进行输入输出。
from kplib import get_kpoints
from pymatgen.core import Structure
from pymatgen.io.vasp.inputs import Kpoints, Kpoints_supported_modes
struct = Structure.from_file("./POSCAR")
kpts = get_kpoints(struct.lattice.matrix,
struct.frac_coords,
struct.atomic_numbers,
min_distance=25.0,
include_gamma="auto")
kpoints_file = Kpoints(style=Kpoints_supported_modes.Reciprocal,
num_kpts=kpts["num_distinct_kpts"]
kpts=kpts["distinct_coords"],
kpts_weights=["distinct_weights"])
kpoints_file.write_file("./KPOINTS")
from kpLib.interface import get_kpoints
import ase
struct = ase.io.read("./POSCAR", format="vasp")
kpts = get_kpoints(struct.get_cell()[:],
struct.get_scaled_positions(), # Note: the atomic coordinates are in fractional format, not Cartesian.
struct.get_atomic_numbers(),
min_distance=25.0,
include_gamma="gamma")
# ... User can then use the kpts dict to build VASP calculator and perform other tasks
命令行界面kpgen
为了方便用户,Python包还提供了一个可以在终端中调用的命令行界面。它使用pymatgen
作为结构解析器,并将K-点网格输出到VASP格式的KPOINTS
文件。
一个简单的示例是
$ kpgen -g auto -d 25.0 ./POSCAR ./KPOINTS
-g
是--gamma
的简称,它定义了K-点网格的类型,即以Γ为中心、平移或自动确定的网格。-d
是--distance
的简称,指定了min_distance
。这是一个上述讨论的get_kpoints
函数的简单包装器。此CLI接受get_kpoints
接受的全部选项。使用kpgen --help
获取允许的参数和选项的完整列表。
默认的pip install kplib
不会安装此CLI。要启用此功能,在克隆源代码后,运行
$ pip install --user .[cli]
"[cli]
"部分指示setup.py
安装pymatgen
和click
。
如何引用
如果您在生成广义k-点网格时使用了kplib或K-点网格生成器,请引用以下论文
王,杨,维萨,巴拉苏布拉曼南,德瓦卡南特,穆勒。快速生成最优的广义Monkhorst-Pack网格。计算材料科学,110100,doi:https://doi.org/10.1016/j.commatsci.2020.110100 (2020)。
文档
约定
本节说明了在C++中,KPointLatticeGenerator
构造函数的参数中我们假设的约定。Python接口kpLib
是以与此一致的方式编写的。
晶格矢量
晶格矢量用晶格矩阵中的行向量表示
double primtiveVectors[3][3] = {{a_x, a_y, a_z},
{b_x, b_y, b_z},
{c_a, c_z, c_y}}
点算符
每个算符是一个3x3整数矩阵,表示在原始晶格基中的分数坐标在这个操作下如何变换。
int latticePointOperations[][3][3];
由于晶格矢量用行表示,每个对称操作都是通过 ${x'}^T = {x}^T\cdot R$ 来完成的。
传统晶格矢量
我们开发的算法可以有效地迭代保持对称性的超晶格,它对输入的传统晶格矢量假设以下条件,除了通常的约定(例如,参见spglib
文档中的标准单元)。
-
对于所有晶体系统,除了单斜晶系之外,$\mathbf{c}$-矢量应沿着最高阶旋转操作的轴,即对于立方、六方、三方、四方、正交和单斜晶格,分别是4重、6重、3重、4重、任何2重和2重旋转。这个方向通常被称为晶体学教科书中“主要对称方向”。国际晶体学表(见第9.1章)中的约定已满足这一条件,除了单斜晶格。
-
对于单斜晶系,国际公约假设$\mathbf{b}$-独特设置,其中$\beta$角是非锐角,两重旋转轴沿$\mathbf{b}$矢量,并且在$\mathbf{a}$-$\mathbf{b}$平面中的一面中心晶格点(即空间群符号中的“C”)。但kplib库假设两重旋转沿$\mathbf{c}$-矢量进行,换句话说,就是$\mathbf{b}$-独特设置。一个简单的转换就是旋转向量,如下所示:$\mathbf{b} \rightarrow \mathbf{c}$,$\mathbf{c} \rightarrow \mathbf{a}$,$\mathbf{a} \rightarrow \mathbf{b}$。
-
对于三方晶格,无论是$hR$还是$hP$表示,传统晶格应该是原始六方晶格,即三棱柱中心六方晶格。这与国际公约相同。
该算法不对三斜晶系或正交晶格的中心类型施加约束。
(注意:用户可以从点对称操作中获得主要方向。)
C++ API
类:KPointLatticeGenerator
template <typename T>
using Tensor = std::vector<std::vector<T>>;
/**
* Constructor.
*
* To see how the variables are defined, check the "Conventions" section below)
*
* @param primVectorsArray 3x3 matrix with primitive lattice vectors in rows.
* @param conventionalVectorsArray 3x3 matrix with conventional lattice vectors in rows.
* @param latticePointOperatorsArray point operators of the Laue Class
of the input structure, expressed
in the basis of primitive lattice vectors
* @param numOperators number of point operators in above array
* @param isConventionalHexaognal whether the conventional lattice is
hexagonal
*/
KPointLatticeGenerator(const double primVectorsArray[][3],
const double conventionalVectorsArray[][3],
const int latticePointOperatorsArray[][3][3],
const int numOperators,
const bool isConventionalHexagonal);
/*
* Specify whether to generate a gamma-centered grid or a shifted grid.
* The available shifts are:
* {{0.0, 0.0, 0.5}, {0.0, 0.5, 0.0}, {0.5, 0.0, 0.0}, {0.5, 0.5, 0.0},
* {0.5, 0.0, 0.5}, {0.0, 0.5, 0.5}, {0.5, 0.5, 0.5}}
* Basiclly, side centers, face centers and the body center.
*
* @param includeGamma TRUE: gamma-centered grid
* FALSE: grid with one of the above shift
* AUTO: search both shifted and gamma-centered grid
* and return the best one.
*/
enum INCLUDE_GAMMA { TRUE, FALSE, AUTO };
void includeGamma(INCLUDE_GAMMA includeGamma);
/*
* @param minDistance The returned grid should have a corresponding
* real-space superlattice whose "minimum periodic distance"
* is no smaller than this value.
* @param minSize Minimum number of total k-points of grids returned.
*/
KPointLattice getKPointLattice(const double minDistance,
const int minSize);
类:KPointLattice
它用于保存找到的k点网格并提供查询函数。此类的主要查询例程
double getMinPeriodicDistance();
int getNumDistinctKPoints();
int numTotalKPoints();
/*
* @return Tensor<double> 2D arrays of coordinates. It's basically a wrapper
* of "double coords[][3]".
*/
Tensor<double> getKPointCoordinates();
/*
* @return vector<int> 1D array of k-points weights.
*/
std::vector<int> getKPointWeights();
示例:使用C++代码中的KpLib -- demo_kplib
为了展示kplib的使用方法,我们实现了一个简单的C++应用程序来生成优化的通用 Monkhorst-Pack k点网格。它使用spglib
来查找对称性(Togo and Tanaka, 2010)并以VASP KPOINTS文件的格式输出k点网格。
它位于demo_kplib
文件夹中。我们已经包含了预构建的spglib库的二进制文件。用户也可以选择从源代码构建spglib
的最新版本并替换它。可以按照以下命令构建可执行的demo_kplib
:
$ cd demo_kplib
$ mkdir build
$ cd build
$ cmake ..
$ make
二进制文件将放在./build/demo_kplib
。要调用它,请使用
$ ./demo_kplib /path/to/POSCAR /path/to/PRECALC > KPOINTS
POSCAR
是VASP标准输入文件之一。《PRECALC
》文件是demo_kplib
的输入文件,用户可以在我们的网站上找到其规范。由于它是演示目的,只有MINDISTANCE
、MINTOTALKPOINTS
和INCLUDEGAMMA
参数有效。
不同晶体系统使用此应用示例可以在 demo_kplib/examples
中找到。KPOINTS_ref
是由我们的服务器生成的 KPOINTS
文件。KPOINTS_kplib
是由 demo_kplib
应用创建的。
对于具有更多功能(如自动检测片状结构)的独立应用程序,请检查我们的 K-点网格生成器 和 K-点服务器。
demo_kpib 代码片段
以下是来自 demo_kplib
应用程序的摘录,展示了如何使用 KpLib API。
demo_kplib/src/main.cpp
:
#include "kPointLatticeGenerator.h"
#include "utils.h"
#include "precalc.h"
#include "poscar.h"
#include <iostream>
int main(int argc, char **argv) {
if (argc < 3) {
std::cerr << "Usage: ./main /path/to/POSCAR /path/to/PRECALC"
<< std::endl;
return 1;
}
// Parse POSCAR and PRECALC.
Poscar poscar;
poscar.readFromPoscar(std::string(argv[1]));
Precalc precalc(argv[2]);
// Execute the main routines.
KPointLatticeGenerator generator = initializeKPointLatticeGeneratorObject(
poscar.primitiveLattice, poscar.coordinates, poscar.atomTypes);
if (precalc.getIncludeGamma() == "TRUE") {
generator.includeGamma(TRUE);
} else if (precalc.getIncludeGamma() == "FALSE") {
generator.includeGamma(FALSE);
} else if (precalc.getIncludeGamma() == "AUTO") {
generator.includeGamma(AUTO);
}
KPointLattice latticeGamma = generator.getKPointLattice(precalc.getMinDistance(),
precalc.getMinTotalKpoints());
outputLattice(latticeGamma);
}
demo_kplib/src/utils.cpp
:
#include "utils.h"
#include "spglib.h"
... // other includes and functions
// Wrapper of the KPointLatticeGenerator constructor.
KPointLatticeGenerator initializeKPointLatticeGeneratorObject(
Tensor<double> primitiveLattice,
Tensor<double> coordinates,
std::vector<int> atomTypes) {
double primLatticeArray[3][3] = {0};
double conventionalLatticeArray[3][3] = {0};
int rotation[192][3][3] = {0};
int size = 0;
bool isConventionalHexagonal = false;
// use spglib to get necessary parameters for the consturctor
// of KPointLatticeGenerator.
...
KPointLatticeGenerator generator = KPointLatticeGenerator(primLatticeArray,
conventionalLatticeArray, rotation, size, isConventionalHexagonal);
if (useScaleFactor == "TRUE") {
generator.useScaleFactor(spaceGroup); // Otherwise, use the fully dynamic search.
}
return generator;
}
使用比例因子加速搜索
比例因子是一种为了加快优化广义网格搜索而实施的方法,但会以返回次优 K-点网格为代价。它在 C++ 库中的 KPointsLatticeGenerator::useScaleFactor
和 Python 接口中的 use_scale_factor
标志控制。
当在 kpLib.interface::get_kpoints
函数中将 use_scale_factor
设置为 True
或在 CLI 命令 kpgen
中,代码将全面搜索,直到达到预定义的总 K-点数阈值。如果找不到有效网格,其有效最小距离至少为 min_distance
,则用户提供的 min_distance
将除以 n
,然后 kplib 再次搜索满足 min_distance/n
的网格。如果仍然没有有效网格,kplib 将 n
增加 1 并再次搜索,直到找到有效网格。当最终找到有效网格时,网格将按 n
缩放(即总 K-点数增加 n^3
倍),以满足 min_distance
要求。目前 n
的最大值设置为 3。穷举搜索的最大总 K-点数阈值分别为:三斜结构 729(9x9x9),单斜结构 1728(12x12x12),正交、四方、三角和六方结构 5832(18x18x18),立方结构 46656(36x36x36)。用户可以在 KPointLatticeGenerator::useScaleFactor
中增加 n
的最大值和这些阈值,然后重新编译。
C++ 库通过 KPointLatticeGenerator::useScaleFactor
函数开启和关闭此方案。它接受输入结构的空间群编号,以提供选择性地为用户选择的子集 lattice 类型使用此方案的功能。
项目详情
kpLib-1.1.1.tar.gz 的哈希值
算法 | 哈希摘要 | |
---|---|---|
SHA256 | a7b8d51e77e383d69ab5896a6a79a7fc02def93a9ab0f2e32407afb2a7060731 |
|
MD5 | 4700a953a0231cc18205a879c298902d |
|
BLAKE2b-256 | e3ce820745322bf79f17bb258d1affe0ad43d483312418bd3387f1c22911911c |