The Plans Module

class ciderpress.dft.plans.FracLaplPlan(settings, nspin, nrho_sl=5, drho_start=1)

Plan for the Fractional Laplacian features. These features do not require additional settings for how to compute them, so this class basically just stores the settings class and nspin. Currently, evaluation of the orbital operations must be implemented in each DFT code separately.

Initialize FracLaplPlan instance

Parameters:

settings (FracLaplSettings) – settings for the fractional laplacian descriptors

class ciderpress.dft.plans.HybridPlan(settings, nspin)

Plan for hybrid DFT. Not yet implemented

Parameters:

settings (HybridSettings)

class ciderpress.dft.plans.NLDFAuxiliaryPlan(nldf_settings, nspin, alpha0, lambd, nalpha, coef_order='gq', alpha_formula='etb', proc_inds=None, rhocut=1e-10, expcut=1e-10, raise_large_expnt_error=True, use_smooth_expnt_cutoff=False)

Plan for expanding the NLDFs using an auxiliary basis or spline. Stores the list of interpolation exponents and how they were constructed, as well as instructions for evaluating the auxiliary coefficients. Also has a method to call the relevant C function to get the coefficients.

Initialize NLDFAuxiliaryPlan

Parameters:
  • nldf_settings (NLDFSettings) – code-universal NLDF settings

  • nspin (int) – 1 (non-spin-polarized) or 2 (spin-polarized)

  • alpha0 (float) – Small-exponent setting for constructing list of alphas.

  • lambd (float) – Inverse density of interpolation. Must be > 1. 1=infinitesimal spacing of Gaussian exponents

  • nalpha (int) – Number of interpolating exponents

  • coef_order (str) – Whether interpolation/expansion coefficients have grid as the slow (first) index (‘gq’) or exponent (‘qg’)

  • alpha_formula (str) – How to generate the interpolating exponents. Must be ‘etb’ (alphas[j] = alpha0 * lambd**j)) or ‘zexp’ (alphas[j] = alpha0 * (lambd**j - 1) / (lambd - 1))

  • proc_inds (list or np.ndarray) – If parallelized, which alpha indexes does this process cover. List or array of integers.

  • rhocut (float) – Small-density cutoff for stability

  • expcut (float) – Small-exponent cutoff for stability

  • raise_large_expnt_error (bool, True) – If True, raise an error if a convolution exponent is larger than the largest interpolation value. It is generally best to set this to true to make sure a calculation doesn’t accidentally lose significant precision to going outside the interpolation range.

  • use_smooth_expnt_cutoff (bool, False) – If True, use a damping function for large exponent to smoothly ensure that as the exponent goes to infinity, the damped exponent goes to alpha_max. The smoothing function is design to contribute insignificantly except for near and beyond alpha_max.

eval_feat_exp(rho_tuple, i=-1)

Evaluate the exponents that determine the length-scale of feature i at each grid point. i can take a range from -1 (theta exponent) to the number of version j/k features minus 1.

Parameters:
  • rho (np.ndarray) – density

  • sigma (np.ndarray) – squared gradient

  • tau (np.ndarray) – kinetic energy density

  • i – feature id. -1 corresponds to theta exponent

Returns:

the cider exponent array (np.ndarray)

eval_rho_full(f, rho_data, spin=0, feat=None, dfeat=None, cache_p=True, apply_transformation=False, coeff_multipliers=None)

Evaluate the raw features. IMPORTANT: For spin-polarized calculations, spin must be passed explicitly and must be different for each spin channel, otherwise cached data will get overwritten, leading to incorrect gradients later on. If apply_transformation is True, f is overwritten.

Parameters:
  • f (np.ndarray) – Feature interpolation coefficients, of shape (ngrids, nalpha) for coef_order==”qg” else (nalpha, ngrids)

  • rho_data (np.ndarray) – Density vector [n, dn/dx, dn/dy, dn/dz, (tau)]

  • spin (int, 0) – Spin index of this density, for caching data.

  • feat (np.ndarray, None) – Shape (nfeat, ngrids) optional buffer for features

  • dfeat (np.ndarray, None) – Shape (nfeat, ngrids) optional buffer for feature derivatives

  • cache_p (bool, True) – Cache nalpha x ngrids-size arrays for interpolation coefficients. Slightly faster, but more memory intensive.

  • apply_transformation (bool, False) – For Gaussian interpolation, apply the linear transformation from the projection values to the interpolation coefficients.

  • coeff_multipliers (np.ndarray, None) – If not None, multiply the interpolation coefficients by this (nalpha,)-shapred array

Returns:

Features (nfeat, ngrids) and derivative of features (nfeat, ngrids) with respect to the feat_params exponent.

Return type:

(np.ndarray, np.ndarray)

eval_vxc_full(vfeat, vrho_data, dfeat, rho_data, spin=0, p_i_qg=None, vf=None)
Parameters:
  • vfeat (np.ndarray) – (nfeat, ngrids) derivative of energy density with respect to features

  • vrho_data (np.ndarray) – (nrho, ngrids) derivative of energy density with respect to density vector [n, dn/dx, dn/dy, dn/dz, (tau)]

  • dfeat (np.ndarray) – (nfeat, ngrids) derivative of features with respect to NLDF exponents.

  • rho_data (np.ndarray) – density vector

  • p_i_qg (np.ndarray, None) – Interpolation coefficients for features. If None, the features must be cached using cache_p=True during eval_rho_full

  • vf (np.ndarray, None) – Optional buffer for vf

Returns:

Shape (nalpha + num_vi_ints + vf_buffer_nslot, ngrids) or transpose if coef_order == ‘gq’. The last vf_buffer_nslot rows (qg) or columns (gq) are zeros. The remaining rows are filled with the functional derivatives with respect to the nonlocal density integrals.

Return type:

(np.ndarray)

get_function_to_convolve(rho_tuple)

Returns the function that gets convolved to make the nonlocal density features. Usually the density rho (settings.rho_mult=’one’), but can also be other quantities like the rho times the theta exponent (rho_mult=’expnt’).

Parameters:
  • rho (np.ndarray) – Density

  • sigma (np.ndarray) – Squared gradient of density

  • tau (np.ndarray) – Kinetic energy density

Returns:

function to convolve, as well as its derivatives with respect to rho, sigma, and tau. Combined, this makes for a tuple of four arrays.

Return type:

tuple(np.ndarray)

get_q2a(q)

Convert (float) q index to exponent alpha.

property num_vj

Number of vj or vk features

class ciderpress.dft.plans.NLDFGaussianPlan(nldf_settings, nspin, alpha0, lambd, nalpha, coef_order='gq', alpha_formula='etb', proc_inds=None, rhocut=1e-10, expcut=1e-10, raise_large_expnt_error=True, use_smooth_expnt_cutoff=False)

Initialize NLDFAuxiliaryPlan

Parameters:
  • nldf_settings (NLDFSettings) – code-universal NLDF settings

  • nspin (int) – 1 (non-spin-polarized) or 2 (spin-polarized)

  • alpha0 (float) – Small-exponent setting for constructing list of alphas.

  • lambd (float) – Inverse density of interpolation. Must be > 1. 1=infinitesimal spacing of Gaussian exponents

  • nalpha (int) – Number of interpolating exponents

  • coef_order (str) – Whether interpolation/expansion coefficients have grid as the slow (first) index (‘gq’) or exponent (‘qg’)

  • alpha_formula (str) – How to generate the interpolating exponents. Must be ‘etb’ (alphas[j] = alpha0 * lambd**j)) or ‘zexp’ (alphas[j] = alpha0 * (lambd**j - 1) / (lambd - 1))

  • proc_inds (list or np.ndarray) – If parallelized, which alpha indexes does this process cover. List or array of integers.

  • rhocut (float) – Small-density cutoff for stability

  • expcut (float) – Small-exponent cutoff for stability

  • raise_large_expnt_error (bool, True) – If True, raise an error if a convolution exponent is larger than the largest interpolation value. It is generally best to set this to true to make sure a calculation doesn’t accidentally lose significant precision to going outside the interpolation range.

  • use_smooth_expnt_cutoff (bool, False) – If True, use a damping function for large exponent to smoothly ensure that as the exponent goes to infinity, the damped exponent goes to alpha_max. The smoothing function is design to contribute insignificantly except for near and beyond alpha_max.

class ciderpress.dft.plans.NLDFSplinePlan(nldf_settings, nspin, alpha0, lambd, nalpha, coef_order='gq', alpha_formula='etb', proc_inds=None, spline_size=None, rhocut=1e-10, expcut=1e-10, raise_large_expnt_error=True, use_smooth_expnt_cutoff=False)

Initialize NLDFAuxiliaryPlan

Parameters:
  • nldf_settings (NLDFSettings) – code-universal NLDF settings

  • nspin (int) – 1 (non-spin-polarized) or 2 (spin-polarized)

  • alpha0 (float) – Small-exponent setting for constructing list of alphas.

  • lambd (float) – Inverse density of interpolation. Must be > 1. 1=infinitesimal spacing of Gaussian exponents

  • nalpha (int) – Number of interpolating exponents

  • coef_order (str) – Whether interpolation/expansion coefficients have grid as the slow (first) index (‘gq’) or exponent (‘qg’)

  • alpha_formula (str) – How to generate the interpolating exponents. Must be ‘etb’ (alphas[j] = alpha0 * lambd**j)) or ‘zexp’ (alphas[j] = alpha0 * (lambd**j - 1) / (lambd - 1))

  • proc_inds (list or np.ndarray) – If parallelized, which alpha indexes does this process cover. List or array of integers.

  • rhocut (float) – Small-density cutoff for stability

  • expcut (float) – Small-exponent cutoff for stability

  • raise_large_expnt_error (bool, True) – If True, raise an error if a convolution exponent is larger than the largest interpolation value. It is generally best to set this to true to make sure a calculation doesn’t accidentally lose significant precision to going outside the interpolation range.

  • use_smooth_expnt_cutoff (bool, False) – If True, use a damping function for large exponent to smoothly ensure that as the exponent goes to infinity, the damped exponent goes to alpha_max. The smoothing function is design to contribute insignificantly except for near and beyond alpha_max.

get_a2q_fast(exp_g)

A fast (parallel C-backend) version of get_a2q that also has a local option to get the q values just for the local alphas.

class ciderpress.dft.plans.SADMPlan(settings, nspin, alpha0, lambd, nalpha, fit_metric='ovlp')

Plan for evaluating exchange energy of spherically averaged density matrix. Stores the settings (basically an empty class since the feature is not tunable, though local range separation might be added later), the nspin, and paramters for an even-tempered Gaussian basis for expanding the spherically averaged density matrix. Currently, evaluation of the integrals must be implemented in each DFT code separately.

Parameters:
  • settings (SADMSettings) – Settings for the descriptor.

  • nspin (int) – Number of spins, 2 for spin-polarized and 1 for non-spin-polarized

  • alpha0 (float) – Minimum exponent. sqrt(1/alpha0) should be at least the width of the system, or altneratively the maximum distance at which the exchange hole is expected to have significant density.

  • lambd (float) – ETB lambda parameter

  • nalpha (int) – Number of ETB basis functions

  • fit_metric (str) – ‘ovlp’ or ‘coul’, the metric for expanding the density matrix in the auxiliary basis.

get_features(p_vag, out=None, l0tmp=None, l1tmp=None)

Should pass out[idm], l0tmp[idm], etc. to this function