Skip to content

QMol_DFT_Vx_XX_conv

Cristel Chandre edited this page Jul 10, 2024 · 4 revisions

QMol_DFT_Vx_XX_conv

Exact exchange (XX) potential and functional.

Description

Use QMol_DFT_Vx_XX_conv to describe the exact exchange potential and functional. For given Kohn-Sham orbitals $\lbrace \phi_k \rbrace_k$ and population parameters, the exact-exchange energy functional is

$$ H _{{\mathrm{XX}}} = H _{{\mathrm{XX}}} ^{\uparrow } + H _{{\mathrm{XX}}} ^{\downarrow } ~~ ~~ {\mathrm{with}} ~~ ~~ H _{{\mathrm{XX}}} ^{\uparrow } = - \frac{1}{2} \sum _{k,l} n _k^{\uparrow } n _l^{\uparrow } \int \int \phi _k^{\uparrow *} (x)\phi _l^{\uparrow *} (x^{\prime } ){\mathcal{V}} _{{\mathrm{e}\mathrm{e}}} (x-x^{\prime } ) \phi _k^{\uparrow } (x^{\prime } )\phi _l^{\uparrow } (x)~dxdx^{\prime } , ~~~~~~ (1) $$

where ${\mathcal{V}}_{{\mathrm{e}\mathrm{e}}}$ is the (effective) electron-electron interaction potential. The corresponding exact-exchange potential acting on a one-electron wave functional is defined as

$$ {\hat{\mathcal{V}} } _{{\mathrm{XX}}} ~\phi (x) = - \sum _l n _l^{\uparrow } \phi_l^{\uparrow } (x) \int \phi _l^{\uparrow *} (x^{\prime } ) {\mathcal{V}} _{{\mathrm{ee}}} (x-x^{\prime } ) \phi (x^{\prime } )~dx^{\prime } .~~~~~~(2) $$

Self-interaction correction has no effect on QMol_DFT_Vx_XX_conv.

Class properties

Electron-interaction potential model

The QMol_DFT_Vx_XX_conv class defines the following public get-access properties; each can be changed using the set method:

interactionPotential (Vee)

Electron-electron interaction potential [ function handle (default @(x) 1./sqrt(x.^2+2) ) ]

  • This is the potential ${\mathcal{V}}_{{\mathrm{e}\mathrm{e}}}$ used in Eqs. (1) and (2).
  • The signature for the (effective) electron-electron interaction potential should be V = funV(x), where the output V has the same shape as the input vector x and contains element-wise values of the potential at the query points x.
  • Because the computation of the exact-exchange potential and energy is performed over an (internally handled) extended domain, a user-defined discretization of the electron-electron interaction potential is not supported (only function handle).

Other properties

These properties cannot be edited with the set method.

isInitialized (isInit)

Whether the potential object is properly initialized. This is used throughout the QMol-grid package to check that the potential object holds meaningful information and is ready for use. Changing its isSpinPol may cause simulations to fail or produce erroneous results.

type

Flavor of DFT functional [ 'XX' ]

  • This is a constant property, which can be used by other components of the QMol-grid package to determine the flavor/type of functional a given object belongs to.

Class methods

Creation

constructor

Create an exact-exchange object with empty class properties.

obj = QMol_DFT_Vx_XX_conv;

Create an exact-exchange object with the name properties set to the specified value. Several name-value pairs can be specified consecutively. Suitable name is any of the electron-interaction potential model and is case insensitive.

obj = QMol_DFT_Vx_XX_conv(name1,value1);
obj = QMol_DFT_Vx_XX_conv(name1,value1,name2,value2,___);

Changing class properties

set

Update the name properties of an exact-exchange object to the specified value. Several name-value pairs can be specified consecutively. Suitable name is any of the electron-interaction potential model and is case insensitive.

obj.set(name1,value1);
obj.set(name1,value1,name2,value2,___);

This is the common name-value pair assignment method used throughout the QMol-grid package. The set method also reset the class. After running, the set property updates the isInitialized flag to a false value.

reset

Reset the object by deleting/re-initializing all run-time properties of the class and updating the isInitialized flag to false.

obj.reset;
  • This is the common reset method available to all classes throughout the QMol-grid package.

clear

Clear all class properties.

obj.clear;

Clear a specific set of the class properties. Suitable name is any of the electron-interaction potential model and is case insensitive.

obj.clear(name1,name2,___);

This is the common clear method available to all classes throughout the QMol-grid package. The clear method also reset the class. The clear method can be used to delete specific properties before saving an instance of the QMol_DFT_Vx_XX_conv class.

Initializing the object

initialize

Initialize a QMol_DFT_Vx_XX_conv object and set the isInitialized flag to true.

obj.initialize(DFT);
  • DFT is the DFT-model handle object, i.e., QMol_DFT_spinPol or QMol_DFT_spinRes, to which the LDA Slater exchange functional is attached.
  • To avoid any mismatch in internal properties, initialize first reset the object before performing the initialization.

Run-time documentation

getMemoryProfile

Get an estimate of the memory held by a QMol_DFT_Vx_XX_conv object with either

mem = obj.getMemoryProfile;
mem = obj.getMemoryProfile(false);
  • The object must be properly initialized with a domain discretization.
  • The estimate only includes the discretization of member components on the domain grid and ignores other (small) properties.
  • The output mem is the estimated size in bytes.

Additionally display the detail of the memory footprint with

mem = obj.getMemoryProfile(true);

showDocumentation

Display the run-time documentation for the specific configuration of a QMol_DFT_Vx_XX_conv object.

ref = obj.showDocumentation;
  • The output ref is a cell vector containing the list of references to be included in the bibliography.

Exact-exchange functional

Before using any of its exact-exchange functional methods, the QMol_DFT_Vx_XX_conv object must be properly initialized.

applyPotential

For spin-restricted models, apply the exact-exchange potential operator to a wave function as in Eq. (2).

Hp = obj.applyPotential(p);
  • Note: The exact-exchange object must have its potential kernel properly set beforehand with setPotentialKernel. applyPotential does not perform or check for this initialization.
  • p and Hp are both numel(disc.xspan)-by-1 vectors respectively containing the discretization of the orbital to which the potential should be applied and its result.

For spin-polarized models, apply the up- and down-spin exact-exchange potential operators to a wave function, as in Eq. (2), respectively with

Hp = obj.applyPotential(p,true);
Hp = obj.applyPotential(p,false);
  • Note: The exact-exchange object must have its potential kernel properly set beforehand with setPotentialKernel. applyPotential does not perform or check for this initialization.
  • p and Hp are both numel(disc.xspan)-by-1 vectors respectively containing the discretization of the orbital to which the potential should be applied and its result.

applyPotentialDerivative

Warning: applyPotentialDerivative is still an experimental feature and untested. Its use is discouraged.

The derivative of the exact-exchange potential is seldom used as-is, but most often appears in the definition of various observables. applyPotentialDerivative provides a common interface for these calculations with

DVp = obj.applyPotentialDerivative(opt,p);          % Spin restricted
DVp = obj.applyPotentialDerivative(opt,p,true);     % Spin polarized, up-spin component
DVp = obj.applyPotentialDerivative(opt,p,false);    % Spin polarized, up-spin component
  • Note: The exact-exchange object must have its potential kernel properly set beforehand with setPotentialKernel. applyPotentialDerivative does not perform or check for this initialization.
  • opt is a character array specifying the type of calculation to be performed (see next). The type and shape of the output DVp is dictated by this option.
  • p is a numel(disc.xspan)-by-1 vector containing the discretization of the orbital to which the potential should be applied.
  • For spin-polarized models only, the third argument specify whether the up- (true) or down-spin (false) exact-exchange potential derivative should be used.

opt = 'dipacc' computes the unweighted dipole-acceleration contribution from the input orbital, defined as

$$ a_{{\mathrm{X}\mathrm{X}}} [\phi ]=2\Re \left(\int \nabla \phi^* (x)~~{\hat{\mathcal{V}} }_{{\mathrm{X}\mathrm{X}}} ~\phi (x)~dx\right). $$

Note that while the total dipole acceleration from the exact-exchange vanishes $\sum_k n_n ~a_{{\mathrm{X}\mathrm{X}}} [\phi_k ]=0$ , the individual contributions from each Kohn-Sham orbital generally do not. As can be seen in the equation above, the output DVp is a real-valued scalar.

getEnergy

Get the exact-energy of Eq. (1) for the parent DFT object.

E = obj.getEnergy;
  • This computes the exchange energy associated with the Kohn-Sham orbitals in the parent DFT model.
  • The output scalar E contains the numerical evaluation of the exchange energy of Eq. (2).
  • This is equivalent to, but more efficient than, obj.getEnergy(DFT.getDensity) with DFT being the same DFT-model handle object used to initialize the external-potential object.

getPotential

Get the exact-exchange potential for the parent DFT object.

V = obj.getPotential;
  • This creates a new Kohn-Sham potential object V in which the exact-exchange potential handle is stored. Note that only the handle is returned and the action of the exact-exchange potential on a given one-electron wave function, following Eq. (2), is computed with applyPotential.
  • Note that getPotential does not initialize the output potential object V.

Overwrite the exact-exchange potential in an existing Kohn-Sham potential object with any of

obj.getPotential([],V);
obj.getPotential([],V,false);
  • Note that the first argument [] is required (to provide a common interface with other DFT functionals).
  • This is similar to above without creating a new Kohn-Sham potential object.
  • Any content in the input object V is erased before assigning the exact-exchange potential to it.

Add the exact-exchange potential to an existing Kohn-Sham potential object.

obj.getPotential([],V,true);    % use parent DFT density
  • This is formally equivalent to the in-place addition $\mathcal{V}\gets \mathcal{V}+{\mathcal{V}}_{{\mathrm{X}\mathrm{X}}}$ .
  • Likewise, the first argument [] is required (to provide a common interface with other DFT functionals).

getPotentialDerivative

Get the exact-exchange potential gradient for the parent DFT object with either

DV = obj.getPotentialDerivative(1);
  • This creates a new Kohn-Sham potential gradient object DV in which the exact-exchange potential gradient handle is stored. Note that only the handle is returned and the action of the exact-exchange potential gradient on a given one-electron wave function is computed with applyPotentialDerivative.
  • Note that getPotentialDerivative does not initialize the output potential gradient object DV.
  • Note that the first input 1 is required. This is to provide a uniform signature with higher dimension where the dimension along which the gradient component is considered must be specified.

Overwrite the exact-exchange potential gradient in an existing Kohn-Sham potential gradient object with any of

obj.getPotentialDerivative(1,[],DV);
obj.getPotentialDerivative(1,[],DV,false);
  • Note that the second argument [] is required (to provide a common interface with other DFT functionals).
  • This is similar to above without creating a new Kohn-Sham potential gradient object.
  • Any content in the input object DV is erased before assigning the exact-exchange potential gradient to it.

Add the exact-exchange potential gradient to an existing Kohn-Sham potential gradient object.

obj.getPotentialDerivative(1,[],DV,true);
  • This is formally equivalent to the in-place addition $\nabla \mathcal{V}\gets \nabla \mathcal{V}+\nabla {\mathcal{V}}_{{\mathrm{X}\mathrm{X}}}$ .
  • Likewise, the first argument [] is required (to provide a common interface with other DFT functionals).

setPotentialKernel

Set the kernel for the applyPotential, and applyPotentialDerivative methods.

obj.setPotentialKernel;
  • This creates a local copy of the parent DFT model (used to initialize the object), which is temporarily stored in the object. This copy is then used to when computing the action of the exact-exchange potential (applyPotential) or the action its derivative (applyPotentialDerivative) onto one-electron wave functions.
  • This local copy is required because, e.g., after applying the action of the exact-exchange potential onto a Kohn-Sham orbital with Eq. (2) the exact-exchange potential operator is changed, which is often not the intended feature inside SCF or propagation intermediate steps.

Examples

Create a discretization domain.

disc = QMol_disc('xspan',-20:.1:25);

Create an exact-exchange functional object with default parameters.

V_X = QMol_DFT_Vx_XX_conv;

Create a minimal DFT-model object required to initialize the exchange functional class and display the run-time documentation.

DFT = QMol_DFT_spinRes('discretization',disc,'occupation',[1 1 1]);
disc.initialize(DFT);
V_X.initialize(DFT);
V_X.showDocumentation;

yielding

  * Exact-exchange functional                         explicit convolution
    Interaction pot. = @(x)1./sqrt(x.^2+2) (elec.-elec.)
    V-01.21.000 (06/17/2024)                                     F. Mauger

Display the estimated memory footprint for the object.

V_X.getMemoryProfile(true);
  * Exact-exchange functional (conv.)                            
    > interaction potential                                         7.0 KB
    > kernel                                                       21.1 KB

Looking at Eq. (2) we see that, at minimum, the DFT model should define occupation parameters and a set of Kohn-Sham orbitals.

% Create Kohn-Sham orbitals
x   =   disc.xspan(:);
dx  =   x(2)-x(1);
p   =  [exp(-x.^2*.25/4),x.*exp(-x.^2*.25/4),x.^2.*exp(-x.^2*.25/4)];

for k = 1:3
for l = 1:k-1
    p(:,k)  =   p(:,k) - p(:,l)*sum(p(:,k).*p(:,l))*dx; % orthonormalize
end
    p(:,k)  =   DFT.disc.DFT_normalizeOrbital(p(:,k));  % normalize
end

KSO = DFT.discretization.DFT_allocateOrbital(3);        % create orbital object
KSO.set('orbital',p);                                   % assign orbitals
KSO.initialize(disc);

% Attach the hand-made Kohn-Sham orbitals to the DFT object
DFT.set('occupation',[1 1 1],'orbital',KSO);
disc.initialize(DFT);                                   % for good measure, re-initialize everything
V_X.initialize(DFT);

For good measure, the beginning of the first part of the code above orthonormalizes the set of "hand-made" Kohn-Sham orbitals. Finally, we compute the action of the exact-exchange potential on a given one-electron wave function and plot the result.

% Apply the exact-exchange potential on a wave function
v   =   (x(:)-1).*exp(-(x(:)-2).^2*.25/1.5^2);          % define one-electron wave function
v   =   DFT.disc.DFT_normalizeOrbital(v);

V_X.setPotentialKernel;                                 % initialize the potential kernel

% Plot the results
figure; hold on
plot(disc.xspan,v,'-','LineWidth',2,'DisplayName','\phi')
plot(disc.xspan,V_X.applyPotential(v),'-','LineWidth',2','DisplayName','V_{XX}\phi')
xlabel('position (a.u.)'); xlim(disc.xspan([1 end]));
ylabel('wave function')
legend show
exact-exchange potential action

Test suite

Run the test suite for the class in normal or summary mode respectively with

QMol_test.test('DFT_Vx_XX_conv');
QMol_test.test('-summary','DFT_Vx_XX_conv');

For developers

Other hidden class properties

QMol_DFT_Vx_XX_conv defines a handful of additional transient/constant and hidden properties to facilitate and speed up computations. These properties cannot be edited with the set method, nor by any function outside of the object (SetAccess=private attribute).

DFT

DFT-model object [ [] (default) | QMol_DFT_spinPol handle object | QMol_DFT_spinRes handle object ]

  • This is a copy of the DFT-model handle object passed to initialize.
  • Un-initialized QMol_DFT_Vx_XX_conv objects, i.e., isInitialized == false , have empty DFT.
  • For practical reasons, DFT is editable by QMol_DFT classes.

KSO

For spin-restricted models, Kohn-Sham orbital kernel [ [] (default) | numel(DFT.disc.x))-by-numel(DFT.occ) matrix]

  • This is a copy of the discretization of Kohn-Sham orbitals on the domain.
  • For basis-set models, the kernel holds a reconstruction of the Kohn-Sham orbitals on the underlying discretization grid.
  • It is set by setPotentialKernel.

KSOup

For spin-polarized models, up-spin Kohn-Sham orbital kernel [ [] (default) | numel(DFT.disc.x))-by-numel(DFT.occ{1}) matrix]

  • This is a copy of the discretization of Kohn-Sham orbitals on the domain.
  • For basis-set models, the kernel holds a reconstruction of the Kohn-Sham orbitals on the underlying discretization grid.
  • It is set by setPotentialKernel.

KSOdw

For spin-polarized models, up-spin Kohn-Sham orbital kernel [ [] (default) | numel(DFT.disc.x))-by-numel(DFT.occ{1}) matrix]

  • This is a copy of the discretization of Kohn-Sham orbitals on the domain.
  • For basis-set models, the kernel holds a reconstruction of the Kohn-Sham orbitals on the underlying discretization grid.
  • It is set by setPotentialKernel.

tol

Population threshold [ nonnegative scalar (default 1e-10) ]

  • Orbitals that have a population coefficient smaller than tol are ignored in computations of the energy of Eq. (1) with getEnergy and of the exact-exchange potential of Eq. (2) with getPotential (and this feature gets automatically included in getPotentialDerivative).

Other class methods

initialize

While QMol_DFT_Vx_XX_conv does not define any self-interaction correction scheme (SIC) -- XX does not suffer from self-interaction errors -- to provide a common interface with other DFT functional objects the initialize method supports passing a SIC scheme as a second argument

obj.initialize(DFT,SIC);
  • The SIC flag is ignored in the initialization (and all calculations in the class).

getEnergy

The exact-exchange energy of Eq. (1) is defined from Kohn-Sham orbitals, not from the one-body density as is common in DFT functionals. To provide a common interface with other DFT functional objects, the getEnergy method supports passing a one-body density object argument.

E = obj.getEnergy(rho);
  • The one-body density object rho is ignored and the exact-exchange energy is computed like in E = obj.getEnergy(); i.e., using the Kohn-Sham orbitals in the parent DFT model.

getPotential

The exact-exchange potential of Eq. (2) is defined from Kohn-Sham orbitals, not from the one-body density as is common in DFT functionals. To provide a common interface with other DFT functional objects, the getPotential method supports passing a one-body density object argument.

V = obj.getPotential(rho);      % create new potential object
obj.getPotential(rho,V);        % use existing potential object
obj.getPotential(rho,V,false);
  • The one-body density object rho is ignored and the exact-exchange potential is returned like in V = obj.getPotential().

getPotentialDerivative

The exact-exchange potential of Eq. (2), and thus its gradient, is defined from Kohn-Sham orbitals, not from the one-body density as is common in DFT functionals. To provide a common interface with other DFT functional objects, the getPotentialDerivative method supports passing a one-body density object argument.

DV = obj.getPotentialDerivative(1,rho);     % create new potential object
obj.getPotentialDerivative(1,rho,DV);       % use existing potential object
obj.getPotentialDerivative(1,rho,DV,false);
  • The one-body density object rho is ignored and the exact-exchange potential is returned like in DV = obj.getPotentialDerivative().

References

[Baker 2015] T.E. Baker, E.M. Stoudenmire, L.O. Wagner, K. Burke, and S.R. White, "One-dimensional mimicking of electronic structure: The case for exponentials," Physical Review B 91, 235141 (2015).

[Mauger 2024] F. Mauger, C. Chandre, M.B. Gaarde, K. Lopata, and K.J. Schafer, "Hamiltonian formulation and symplectic split-operator schemes for time-dependent density-functional-theory equations of electron dynamics in molecules," Communications in Nonlinear Science and Numerical Simulation 129, 107685 (2024).

Notes

The results displayed in this documentation page were generated using version 01.21 of the QMol-grid package.

  • QMol_DFT_Vx_XX_conv was introduced in version 01.00.
  • getMemoryProfile was introduced in version 01.10.
Density-functional theory (DFT)

$~~$ Hartree-Fock theory (HF)

Schrödinger equation (SE)
Time-dependent density-functional theory (TDDFT)
Time-dependent Schrödinger equation (TDSE)
Discretization
Pseudopotentials
External field
External components
Tutorials

$~~$ Documentation

$~~$ Test suite

For developers
Clone this wiki locally