Skip to content

Commit 7a08986

Browse files
committed
Add a doc-page for multigrid
1 parent d5dac46 commit 7a08986

File tree

1 file changed

+242
-0
lines changed

1 file changed

+242
-0
lines changed

source/modules/pbc/multigrid.rst

Lines changed: 242 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,242 @@
1+
.. _pbc_dft_multigrid:
2+
3+
pbc.dft.multigrid --- Multigrid Integration
4+
*******************************************
5+
6+
Multigrid is a numerical integration algorithm optimized for the computation of
7+
the Coulomb matrix and exchange-correlation (XC) functional. It can take
8+
advantage of the locality of density and orbitals when computing electron
9+
density and matrix elements.
10+
Compared to the `get_j` function provided by the FFTDF module and the XC matrix
11+
evaluation functions offered by the `numint` module, the Multigrid algorithm
12+
can achieve an order of magnitude improvement in the computation of the Coulomb
13+
matrix and DFT XC matrix.
14+
15+
Supported Applications
16+
======================
17+
18+
The Multigrid algorithm is designed to accelerate DFT calculations.
19+
It can also be utilized to speed up derived properties based on DFT, such as
20+
TDDFT excited state calculations, stability analysis, analytical nuclear
21+
gradients, and the computation of the orbital Hessian in second-order SCF
22+
(SOSCF) methods. However, it does not support the computation of the HF exchange
23+
matrix or the calculation of MO two-electron integrals in post-HF methods.
24+
25+
In terms of supported systems, the Multigrid algorithm can be applied to both
26+
periodic systems and molecular calculations. For periodic systems, it supports
27+
single k-point calculations as well as k-point sampling in DFT calculations. For
28+
molecular systems, additional configurations are required to enable the multigrid
29+
algorithm, see :ref:`mole_multigrid`.
30+
31+
Generally, the Multigrid algorithm should be used with pseudopotentials and the
32+
corresponding basis sets to reduce computational costs and accuracy.
33+
Using the Multigrid algorithm for all-electron calculations typically results in significant overhead.
34+
35+
36+
How to enable MultiGrid algorithm
37+
=================================
38+
39+
In periodic system calculations, enabling the Multigrid algorithm is
40+
straightforward. The PBC SCF class provides a method `multigrid_numint()`.
41+
This method creates a new SCF instance which is configured with the Multigrid
42+
functionality. This new instance can then be used in the same way as a standard
43+
SCF instance. For example::
44+
45+
cell = pyscf.M(
46+
a = np.eye(3)*3.5668,
47+
atom = '''C 0. 0. 0.
48+
C 0.8917 0.8917 0.8917
49+
C 1.7834 1.7834 0.
50+
C 2.6751 2.6751 0.8917
51+
C 1.7834 0. 1.7834
52+
C 2.6751 0.8917 2.6751
53+
C 0. 1.7834 1.7834
54+
C 0.8917 2.6751 2.6751''',
55+
basis = 'gth-dzv',
56+
pseudo = 'gth-pade',
57+
)
58+
mf = cell.KRKS(xc='pbe', kpts=cell.make_kpts([3,3,3]))
59+
mf = mf.multigrid_numint()
60+
mf.run()
61+
62+
Once the Multigrid algorithm is enabled, it is automatically applied to other
63+
methods based on that SCF instance, such as TDDFT and SOSCF. For example::
64+
65+
mf = cell.KRKS(xc='pbe', kpts=cell.make_kpts([3,3,3]))
66+
mf = mf.multigrid_numint()
67+
68+
mf = mf.soscf().run()
69+
70+
mf.TDA().run()
71+
72+
mf.Gradients().run()
73+
74+
The Multigrid algorithm is enabled through the `._numint` attribute of the SCF
75+
instance You can directly modify the `._numint` attribute of an SCF instance to
76+
apply the Multigrid algorithm.::
77+
78+
from pyscf.pbc.dft.numint import MultiGridNumInt
79+
mf = cell.KRKS(xc='pbe', kpts=cell.make_kpts([3, 3, 3]))
80+
mf._numint = MultiGridNumInt(cell)
81+
mf.run()
82+
83+
Common Options
84+
==============
85+
In most scenarios, the Multigrid module can automatically configure parameters such
86+
as energy cutoff, radius cutoff, sub tasks to balance computational load and
87+
accuracy. If you need more precise control over the Multigrid algorithm,
88+
you can modify the attributes of the `mf._numint`. Here are some configurable options.
89+
90+
* `mesh`: It's a 3-element tuple, which controls the upper limit of the mesh
91+
points utilized by the algorithm. By default, this value is estimated based on
92+
the basis set, the shape, and the size of the lattice. By default, the error in
93+
the electron density at each mesh point is controlled to be below `cell.precision`.
94+
95+
* `ntasks`: The maximum number of subgroups. The Multigrid algorithm divides
96+
orbital pairs into subgroups and processes the integration for each subgroup
97+
separately. Increasing the number of tasks can better utilize the locality of
98+
the orbital functions, reducing the number of floating-point operations
99+
(FLOPS). However, it can also increase the overhead of other operations, such
100+
as FFT. The recommended value for `ntasks` is typically set between 3 and 6.
101+
The default value is 4.
102+
103+
* `xc_with_j`: Determines whether the calculation of the XC matrix should also
104+
include the calculation of the Coulomb energy and the Coulomb matrix.
105+
Combining the calculation of the XC and Coulomb can reduce the computation cost.
106+
For GGA functionals, it can reduce the computational load by 30% - 40%.
107+
The default setting is `True`, and it is generally recommended to enable this option.
108+
109+
* `ke_ratio`: Controls the ratio of the energy cutoffs among the different
110+
subtasks when partitioning the orbital pairs. `ke_ratio` is a number greater
111+
than 1. When approaching 1, the task partitioning is fine, and the orbital
112+
locality is better utilized. However, this also increases the number of tasks.
113+
Regardless of the configuration of this parameter, the maximum number of tasks
114+
generated is still controlled by the `ntasks` attribute.
115+
116+
117+
Compatibility with NumInt
118+
=========================
119+
120+
The `NumInt` class offers numerical integration functions for DFT calculations.
121+
It is assigned to the `_numint` attribute of SCF classes. This module offers
122+
APIs such as `get_rho` (for computing electron density in real space) and
123+
`get_vxc` (for computing the value of the XC energy functional and XC matrices),
124+
for DFT calculations :ref:`dft`. Both molecular and PBC DFT implementations
125+
follow this API design.
126+
127+
The `pyscf.pbc.dft.multigrid` module offers the `MultiGridNumInt` class, which
128+
is compatible with the `NumInt` class. Its methods, such as `get_vxc`, `get_fxc`,
129+
`get_rho`, `cache_xc_kernel`, `nr_rks`, and `nr_uks`, are mostly compatible
130+
with the corresponding methods in NumInt (with only a few additional keyword
131+
arguments for controlling multigrid instances). These methods can be
132+
individually invoked, like those in the `NumInt` class, to compute densities and
133+
XC matrix elements.
134+
135+
The `pyscf.pbc.dft.multigrid` module also provides the `MultiGridNumInt2` class,
136+
which further optimizes the implementations of the `MultiGridNumInt` class.
137+
However, due to differences in the algorithm implementations, the support for
138+
optimized k-points and non-orthogonal lattices is not as comprehensive as that
139+
in the `MultiGridNumInt` class. Currently, the `SCF.multigrid_numint()` method
140+
invokes the `MultiGridNumInt` class. To maximize the multigrid performance, you
141+
can manually assign the `MultiGridNumInt2` instance to the `mf._numint`
142+
attribute.
143+
144+
The two classes will be merged into one in the future release.
145+
146+
147+
.. _mole_multgird:
148+
149+
How to Apply Multigrid in Molecular DFT
150+
=======================================
151+
152+
The Multigrid algorithm is currently designed and implemented for data
153+
structures with periodic boundary conditions. Nevertheless, it can be adapted
154+
for molecular calculations.
155+
156+
First, we need to initialize molecule within a relatively large periodic lattice.
157+
A vacuum space needs to be placed between the box and the molecule to
158+
simulate free boundary conditions::
159+
160+
cell = pyscf.M(atom='''
161+
O 0.000 0.118 0.
162+
H 0.758 -0.470 0.
163+
H -0.758 -0.470 0.''',
164+
a=np.eye(3)*10,
165+
dimension=0,
166+
basis='gth-dzvp', pseudo='gth-pade',)
167+
168+
Here, we apply a 10 x 10 x 10 (unit Angstrom) box. The box does not need to be
169+
excessively large. It only needs to ensure that the electron density of the molecule
170+
does not leak outside the box. If there are no diffused functions, typically, a
171+
5 Angstrom margin around the molecule is sufficient. The lattice is a virtual
172+
box which does not need to be centered on the origin. There is no need to adjust
173+
the molecule's coordinates to the center of the box.
174+
175+
Alternatively, the lattice can be automatically determined. We can just create a
176+
`Mole` instance and then call the `Mole.to_cell()` method to convert a `Mole`
177+
object into a Cell object::
178+
179+
mol = pyscf.M(atom='''
180+
O 0.000 0.118 0.
181+
H 0.758 -0.470 0.
182+
H -0.758 -0.470 0.''',
183+
basis='gth-dzvp', pseudo='gth-pade')
184+
cell = mol.to_cell(dimension=0)
185+
186+
When initializing the Cell, we set `dimension=0`. This setting informs the
187+
program that the system is a 0-dimensional system (molecule), which allows the
188+
program to more effectively utilize this property to truncate the long-range
189+
Coulomb potential and accelerate the computation of certain integrals.
190+
The system can also be treated as a 3-dimensional crystal, with slightly
191+
increased computational load.
192+
193+
In this system, we use the GTH pseudo potential and the GTH basis set. This
194+
basis set does not have very steep GTO functions, making it suitable for the
195+
multigrid method.
196+
197+
With this `Cell` object, we can initialize the DFT instance as we would for
198+
typical PBC DFT calculations. The following example demonstrates another
199+
way to utlize the Multigrid algorithm, whihc integrates
200+
the Multigrid XC matrix functionality into the the molecular DFT instances.::
201+
202+
from pyscf.pbc.dft.multigrid import MultiGridNumInt2
203+
mol = pyscf.M(atom='''
204+
O 0.000 0.118 0.
205+
H 0.758 -0.470 0.
206+
H -0.758 -0.470 0.''',
207+
basis='gth-dzvp', pseudo='gth-pade')
208+
mf = mol.RKS(xc='pbe')
209+
210+
cell = mol.to_cell(dimension=0)
211+
mf._numint = MultiGridNumInt2(cell)
212+
mf._numint.xc_with_j = False
213+
214+
mf.run()
215+
216+
In this example, we use the `MultiGridNumInt` class for `mf._numint` of the
217+
molecular DFT instance. This setup invokes the `MultiGridNumInt` algorithm to
218+
calculate the DFT XC matrix and XC energy, while calls the standard molecular
219+
`get_jk` functions for Coulomb energy.
220+
Note the setting `xc_with_j=False`, which disables the computation of Coulomb
221+
energy by the `MultiGridNumInt` class. This is necessary because the molecular
222+
DFT program employs a different method to compute nuclear repulsion energy
223+
compared to the method used in PBC DFT.
224+
When using the molecular DFT program in conjunction with `MultiGridNumInt`, we
225+
should not use `MultiGridNumInt` to compute the Coulomb energy.
226+
227+
228+
Limitations of Multigrid algorithm
229+
==================================
230+
231+
* The Multigrid algorithm only supports uniform grids. Currently, it cannot be used with Becke grids.
232+
233+
* The `MultiGridNumInt` class does not support the calculation of analytical nuclear gradients.
234+
235+
* The `MultiGridNumInt2` class does not support non-orthogonal lattices, k-point
236+
calculations, meta-GGA functionals, or the computation of the TDDFT fxc kernel.
237+
238+
239+
Examples
240+
========
241+
* :source:`examples/pbc/27-multigrid.py`
242+
* :source:`examples/pbc/27-multigrid2.py`

0 commit comments

Comments
 (0)