1
0
mirror of https://github.com/TREX-CoE/qmckl.git synced 2024-06-30 00:44:52 +02:00
qmckl/org/qmckl_mo.org

29 KiB

Molecular Orbitals

The molecular orbitals (MOs) are defined in the basis of AOs along with a AO to MO coefficient matrix \[C\]. Using these coefficients (e.g. from Hartree Fock SCF method) the MOs are defined as follows:

\[ \phi_i(\mathbf{r}) = C_i * \chi_i (\mathbf{r}) \]

In this section we demonstrate how to use the QMCkl specific DGEMM function to calculate the MOs.

Context

The following arrays are stored in the context:

mo_num Number of MOs
coefficient [mo_num][ao_num] Orbital coefficients

Computed data:

mo_vgl [5][elec_num][mo_num] Value, gradients, Laplacian of the MOs at electron positions
mo_vgl_date uint64_t Late modification date of Value, gradients, Laplacian of the MOs at electron positions

Data structure

typedef struct qmckl_mo_basis_struct {
int64_t   mo_num;
double  * coefficient;

double  * mo_vgl;
int64_t   mo_vgl_date;

int32_t   uninitialized;
bool      provided;
} qmckl_mo_basis_struct;

The uninitialized integer contains one bit set to one for each initialization function which has not been called. It becomes equal to zero after all initialization functions have been called. The struct is then initialized and provided == true. Some values are initialized by default, and are not concerned by this mechanism.

qmckl_exit_code qmckl_init_mo_basis(qmckl_context context);
qmckl_exit_code qmckl_init_mo_basis(qmckl_context context) {

if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) {
 return false;
}

qmckl_context_struct* const ctx = (qmckl_context_struct* const) context;
assert (ctx != NULL);

ctx->mo_basis.uninitialized = (1 << 2) - 1;

return QMCKL_SUCCESS;
}

Access functions

When all the data for the AOs have been provided, the following function returns true.

bool qmckl_mo_basis_provided (const qmckl_context context);

Initialization functions

To set the basis set, all the following functions need to be called.

qmckl_exit_code  qmckl_set_mo_basis_mo_num           (qmckl_context context, const int64_t   mo_num);
qmckl_exit_code  qmckl_set_mo_basis_coefficient      (qmckl_context context, const double  * coefficient);

#+NAME:pre

#+NAME:post

When the basis set is completely entered, other data structures are computed to accelerate the calculations.

Computation

Computation of MOs

Get

qmckl_exit_code qmckl_get_mo_basis_vgl(qmckl_context context, double* const mo_vgl);

Provide

Compute

qmckl_context context in Global state
int64_t ao_num in Number of AOs
int64_t mo_num in Number of MOs
int64_t elec_num in Number of electrons
double coef_normalized[mo_num][ao_num] in AO to MO transformation matrix
double ao_vgl[5][elec_num][ao_num] in Value, gradients and Laplacian of the AOs
double mo_vgl[5][elec_num][mo_num] out Value, gradients and Laplacian of the MOs
integer function qmckl_compute_mo_basis_vgl_f(context, &
 ao_num, mo_num, elec_num, &
 coef_normalized, ao_vgl, mo_vgl) &
 result(info)
use qmckl
implicit none
integer(qmckl_context), intent(in)  :: context
integer*8             , intent(in)  :: ao_num, mo_num
integer*8             , intent(in)  :: elec_num
double precision      , intent(in)  :: ao_vgl(ao_num,elec_num,5)
double precision      , intent(in)  :: coef_normalized(ao_num,mo_num)
double precision      , intent(out) :: mo_vgl(mo_num,elec_num,5)
character                     :: TransA, TransB
double precision,dimension(:,:),allocatable    :: mo_vgl_big
double precision,dimension(:,:),allocatable    :: ao_vgl_big
!double precision,dimension(:,:),allocatable    :: coef_trans
!double precision,dimension(:),allocatable    :: coef_all
double precision              :: alpha, beta
integer                       :: info_qmckl_dgemm_value
integer*8                     :: M, N, K, LDA, LDB, LDC, i,j, idx

integer*8 :: inucl, iprim, iwalk, ielec, ishell
double precision :: x, y, z, two_a, ar2, r2, v, cutoff

allocate(mo_vgl_big(mo_num,elec_num*5))
allocate(ao_vgl_big(ao_num,elec_num*5))
!allocate(coef_all(mo_num*ao_num))
!allocate(coef_trans(mo_num,ao_num))

TransA = 'T'
TransB = 'N'
alpha = 1.0d0
beta  = 0.0d0

info = QMCKL_SUCCESS
info_qmckl_dgemm_value = QMCKL_SUCCESS

! Don't compute exponentials when the result will be almost zero.
! TODO : Use numerical precision here
cutoff = -dlog(1.d-15)
M = mo_num
N = elec_num*5
K = ao_num * 1_8
LDA = size(coef_normalized,1)
idx = 0
!do j = 1,ao_num
!do i = 1,mo_num
!  idx = idx + 1
!  coef_all(idx) = coef_normalized(i,j)
!end do
!end do
!idx = 0
!do j = 1,mo_num
!do i = 1,ao_num
!  idx = idx + 1
!  coef_trans(j,i) = coef_all(idx)
!end do
!end do

ao_vgl_big = reshape(ao_vgl(:, :, :),(/ao_num, elec_num*5_8/))
LDB = size(ao_vgl_big,1)
LDC = size(mo_vgl_big,1)

info = qmckl_dgemm(context,TransA, TransB, M, N, K, alpha,     &
                             coef_normalized,size(coef_normalized,1)*1_8,    &
                             ao_vgl_big, size(ao_vgl_big,1)*1_8, &
                             beta,                                       &
                             mo_vgl_big,LDC)
mo_vgl = reshape(mo_vgl_big,(/mo_num,elec_num,5_8/))

deallocate(mo_vgl_big)
deallocate(ao_vgl_big)

end function qmckl_compute_mo_basis_vgl_f
qmckl_exit_code qmckl_compute_mo_basis_vgl (
const qmckl_context context,
const int64_t ao_num,
const int64_t mo_num,
const int64_t elec_num,
const double* coef_normalized,
const double* ao_vgl,
double* const mo_vgl );

Test