Dice interface
Dice is a specialized WF program developed by Prof. Sandeep Sharma. See: Dice Github page and Dice documentation
Dice contains various methods such as SHCI, VMC, GFMC, DMC, FCIQMC, stochastic MRCI and SC-NEVPT2, and AFQMC calculations. The program is best used together with PySCF (see PySCF interface). Primarily it is useful for conveniently performing semi-stochastic heat-bath CI calculations (SHCI) on large systems, allowing CAS-type calculations up to an active space of approximately 100 orbitals.
The program needs to be compiled and set-up together with PySCF.
DiceTheory class:
class DiceTheory:
def __init__(self, dicedir=None, pyscftheoryobject=None, filename='input.dat', printlevel=2, numcores=1,
moreadfile=None, initial_orbitals='MP2', CAS_AO_labels=None,memory=20000, frozencore=True,
SHCI=False, NEVPT2=False, AFQMC=False, label="Dice",
SHCI_stochastic=True, SHCI_PTiter=200, SHCI_sweep_iter= [0,3],
SHCI_DoRDM=False, SHCI_sweep_epsilon = [ 5e-3, 1e-3 ], SHCI_macroiter=0,
SHCI_davidsonTol=5e-05, SHCI_dE=1e-08, SHCI_maxiter=9, SHCI_epsilon2=1e-07,
SHCI_epsilon2Large=1000, SHCI_targetError=0.0001, SHCI_sampleN=200,
SHCI_nroots=1, SHCI_cas_nmin=1.999, SHCI_cas_nmax=0.0, SHCI_active_space=None,
SHCI_active_space_range=None, SHCI_active_all_but_core=None,
Dice_SHCI_direct=None, fcidumpfile=None, Dice_refdeterminant=None,
QMC_trialWF=None, QMC_SHCI_numdets=1000, QMC_dt=0.005, QMC_nsteps=50,
QMC_nblocks=1000, QMC_nwalkers_per_proc=5):
Keyword |
Type |
Default value |
Details |
---|---|---|---|
|
string |
None |
Path to MRCC directory. |
|
integer |
2 |
Printlevel |
|
string |
'Dice' |
String-label used for some output |
|
integer |
1 |
Number of CPU cores that Dice will use (MPI parallelization). |
|
string |
'input.dat' |
Name used for Dice inputfile |
|
integer |
20000 |
Memory in MB for Dice. |
|
PySCFTheory object |
None |
A PySCFTheory object defining a mean-field calculation. |
|
string |
None |
Name of a PySCF checkpoint file (.chk) to be used as input orbitals for Dice. |
|
string |
'MP2' |
|
|
list |
None |
For input-orbitals options (AVAS-CASSCF, DMET-CASSCF) this list selects the active space (see PySCFTheory docs). |
|
Boolean |
True |
For AFQMC and NEVPT2 calculations (not SHCI!), this defines a frozen-core. |
|
Boolean |
False |
Whether to do an SHCI calculation or not |
|
Boolean |
False |
Whether to do a NEVPT2 calculation or not on top of the SHCI CAS WF. |
|
Boolean |
False |
Whether to do an AFQMC calculation or not. |
|
Boolean |
True |
For SHCI: Whether to do the stochastic PT contribution or not |
|
integer |
200 |
For SHCI: How many PT iterations to do. Set to 0 to turn PT off. |
|
list |
[0,3] |
Control the SHCI sweep iterations. See Dice documentation. |
|
list |
[ 5e-3, 1e-3 ] |
The epsilon values to use in the beginning and end of SHCI sweep iterations. |
|
Boolean |
False |
Whether to calculate the density matrices or not. |
|
integer |
0 |
SHCI macro iterations. If > 0 then SHCI-CASSCF is performed (orbital optimization). |
|
float |
5e-05 |
SHCI tolerance for the final Davidson step |
|
float |
1e-08 |
SHCI energy convergence tolerance for the variational part. |
|
float |
1e-08 |
SHCI max iterations for the variational part. |
|
float |
1e-07 |
SHCI: Lower limit for accepting determinants in the selection. Also controls stochastic PT part. |
|
float |
1000 |
SHCI: Activates stochastic PT. Specifies the limit for what PT part will be done deterministically (rest stochastically). |
|
float |
0.0001 |
SHCI: Target standard deviation of the (semi-)stochastic-corrected energy. |
|
integer |
200 |
SHCI: Number of times the determinants outside variational space are sampled in each stochastic iteration. |
|
integer |
1 |
SHCI: Number of states to solve for. |
|
list |
None |
SHCI: Defining active space as n electrons in m orbitals: e.g. [2,4] for 2 electrons in 4 orbitals. |
|
list |
None |
SHCI: Defining active space as a range of orbital indices: e.g. [2,40]. |
|
list |
None |
SHCI: Experimental: |
|
float |
1.999 |
SHCI: Upper limit for selecting natural orbitals for the active space. Requires initial_orbitals to have natural occupations. |
|
float |
0.0 |
SHCI: Lower limit for selecting natural orbitals for the active space. Requires initial_orbitals to have natural occupations. |
|
Boolean |
False |
SHCI: Run Dice directly without pySCF or SHCI plugin. Requires FCIDUMP file and ref determinant (see below) |
|
string |
None |
SHCI: Name of FCIDUMP file containing orbitals and integrals. |
|
string |
None |
SHCI: String specifying reference determinant |
|
string |
None |
QMC: Trial WF for QMC. Set to 'SHCI' to use a SHCI trial WF. |
|
integer |
1000 |
QMC: Number of determinants to use in the SHCI trial WF. |
|
integer |
0.005 |
QMC: Time step for QMC. |
|
integer |
50 |
QMC: Number of steps per block in QMC. |
|
integer |
1000 |
QMC: Number of blocks used in QMC. |
|
integer |
5 |
QMC: Number of walkers per MPI process for QMC. |
Installing Dice
You need to download the Dice source code from Dice Github page and compile it according to the Github instructions. You also need to have installed pyscf (see PySCF interface) and install the SHCI plugin: https://github.com/pyscf/shciscf After some additional settings modification (ASH will prompt you) you should be ready to go.
Using the interface
Typically you first create a pySCFTheory object and then a DiceTheory object pointing to the pySCFTheory object. The default settings for SHCI are mostly sensible with epsilon being the most important parameter.
See the Dice documentation for more details on various settings: https://sanshar.github.io/Dice/gettingstarted.html https://sanshar.github.io/Dice/keywords.html
Parallelization
Parallelization of Dice is possible via MPI if you compiled it with MPI. Just provide the numcores keyword and the MPI environment will need to have been set (PATH, LD_LIBRARY_PATH to the MPI program).
Examples
Example 1: Dice semi-stochastic heat-bath CI calculation
from ash import *
numcores=10
#Fragment
fragment = Fragment(xyzfile="al2h2_mp2.xyz", charge=0, mult=1)
#PySCF object: RHF/cc-pVTZ mean-field calculation
pyscfobject = PySCFTheory(basis="cc-pVTZ", numcores=numcores, scf_type='RHF', conv_tol=1e-9,memory=50000)
#DiceTheory object for an SHCI calculation.
#The pySCFTheory mean-field object has to be provided. The active space is generated using CCSD natural orbitals as input and
#selecting natural orbitals with occupations between 1.999 and 0.0.
eps=1e-4
dicecalc = DiceTheory(pyscftheoryobject=pyscfobject, numcores=numcores, SHCI=True, memory=50000,
initial_orbitals='CCSD', SHCI_cas_nmin=1.999, SHCI_cas_nmax=0.0, SHCI_stochastic=True,
SHCI_PTiter=400, SHCI_sweep_iter= [0,3,6],SHCI_sweep_epsilon = [ 4*eps,2*eps,eps ],
SHCI_davidsonTol=1e-8, SHCI_epsilon2=1e-8, SHCI_epsilon2Large=1e-5, SHCI_macroiter=0)
#Now running Singlepoint job
result = Singlepoint(fragment=fragment, theory=dicecalc)
print(f"Dice eps={eps}: Energy: {result.energy}”)
Example 2: Dice SHCI calculation with decreasing epsilon parameter
from ash import *
numcores=10
#Fragment
fragment = Fragment(xyzfile="al2h2_mp2.xyz", charge=0, mult=1)
#PySCF object: RHF/cc-pVTZ mean-field calculation
pyscfobject = PySCFTheory(basis="cc-pVTZ", numcores=numcores, scf_type='RHF', conv_tol=1e-9)
#Looping over epsilon values
for eps in [1e-2,5e-3,1e-3,5e-4,1e-4,5e-5,1e-5,5e-6,1e-6]:
dicecalc = DiceTheory(pyscftheoryobject=pyscfobject, numcores=numcores, SHCI=True, memory=50000,
initial_orbitals='CCSD', SHCI_cas_nmin=1.999, SHCI_cas_nmax=0.0, SHCI_stochastic=True,
SHCI_PTiter=400, SHCI_sweep_iter= [0,3,6], SHCI_sweep_epsilon = [ 4*eps,2*eps,eps ],
SHCI_davidsonTol=1e-8, SHCI_epsilon2=1e-8, SHCI_epsilon2Large=1e-5, SHCI_macroiter=0)
#Now running Singlepoint job for each epsilon
result = Singlepoint(fragment=fragment, theory=dicecalc)
print(f"Dice eps={eps}: Energy: {result.energy}”)