Molecular Orbitals
Table of Contents
1 Context
The following arrays are stored in the context:
mo_num |
Number of MOs | |
coefficient |
[mo_num][ao_num] |
MO coefficients |
coefficient_t |
[ao_num][mo_num] |
Transposed of the Orbital coefficients |
Computed data:
mo_value |
[point_num][mo_num] |
Value of the MOs at point positions |
mo_value_date |
uint64_t |
Late modification date of the value of the MOs at point positions |
mo_vgl |
[point_num][5][mo_num] |
Value, gradients, Laplacian of the MOs at point positions |
mo_vgl_date |
uint64_t |
Late modification date of Value, gradients, Laplacian of the MOs at point positions |
1.1 Data structure
typedef struct qmckl_mo_basis_struct { int64_t mo_num; double * restrict coefficient; double * restrict coefficient_t; double * restrict mo_vgl; double * restrict mo_value; uint64_t mo_vgl_date; uint64_t mo_value_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*) context; assert (ctx != NULL); ctx->mo_basis.uninitialized = (1 << 2) - 1; return QMCKL_SUCCESS; }
1.2 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);
When the basis set is completely entered, other data structures are computed to accelerate the calculations.
1.3 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);
1.3.1 Fortran interfaces
1.4 Update
Useless MOs can be removed, for instance virtual MOs in a single determinant calculation.
To select a subset of MOs that will be kept, create an array of
integers of size mo_num
. If the integer is zero, the MO is dropped,
otherwise it is kept.
bool qmckl_mo_basis_select_mo (const qmckl_context context, const int32_t* keep, const int64_t size_max);
1.4.1 Fortran interface
2 Computation
2.1 Computation of MOs: values only
2.1.1 Get
qmckl_exit_code qmckl_get_mo_basis_mo_value(qmckl_context context, double* const mo_value, const int64_t size_max);
Uses the given array to compute the values.
qmckl_exit_code qmckl_get_mo_basis_mo_value_inplace (qmckl_context context, double* const mo_value, const int64_t size_max);
2.1.2 Provide
qmckl_exit_code qmckl_provide_mo_basis_mo_value(qmckl_context context);
qmckl_exit_code qmckl_provide_mo_basis_mo_value(qmckl_context context) { qmckl_exit_code rc = QMCKL_SUCCESS; if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return qmckl_failwith( context, QMCKL_INVALID_CONTEXT, "qmckl_provide_mo_basis_mo_value", NULL); } qmckl_context_struct* const ctx = (qmckl_context_struct*) context; assert (ctx != NULL); if (!ctx->mo_basis.provided) { return qmckl_failwith( context, QMCKL_NOT_PROVIDED, "qmckl_provide_mo_basis_mo_value", NULL); } /* Compute if necessary */ if (ctx->point.date > ctx->mo_basis.mo_value_date) { qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = ctx->mo_basis.mo_num * ctx->point.num * sizeof(double); if (ctx->mo_basis.mo_value != NULL) { qmckl_memory_info_struct mem_info_test = qmckl_memory_info_struct_zero; rc = qmckl_get_malloc_info(context, ctx->mo_basis.mo_value, &mem_info_test); /* if rc != QMCKL_SUCCESS, we are maybe in an _inplace function because the memory was not allocated with qmckl_malloc */ if ((rc == QMCKL_SUCCESS) && (mem_info_test.size != mem_info.size)) { rc = qmckl_free(context, ctx->mo_basis.mo_value); assert (rc == QMCKL_SUCCESS); ctx->mo_basis.mo_value = NULL; } } /* Allocate array */ if (ctx->mo_basis.mo_value == NULL) { double* mo_value = (double*) qmckl_malloc(context, mem_info); if (mo_value == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "qmckl_mo_basis_mo_value", NULL); } ctx->mo_basis.mo_value = mo_value; }
if (rc != QMCKL_SUCCESS) { return rc; } ctx->mo_basis.mo_value_date = ctx->date; } return QMCKL_SUCCESS; }
2.1.3 Compute
Variable | Type | In/Out | Description |
---|---|---|---|
context |
qmckl_context |
in | Global state |
ao_num |
int64_t |
in | Number of AOs |
mo_num |
int64_t |
in | Number of MOs |
point_num |
int64_t |
in | Number of points |
coefficient_t |
double[mo_num][ao_num] |
in | Transpose of the AO to MO transformation matrix |
ao_value |
double[point_num][ao_num] |
in | Value of the AOs |
mo_value |
double[point_num][mo_num] |
out | Value of the MOs |
The matrix of AO values is very sparse, so we use a sparse-dense matrix multiplication instead of a dgemm, as exposed in https://dx.doi.org/10.1007/978-3-642-38718-0_14.
integer function qmckl_compute_mo_basis_mo_value_doc_f(context, & ao_num, mo_num, point_num, & coefficient_t, ao_value, mo_value) & result(info) use qmckl implicit none integer(qmckl_context), intent(in) :: context integer*8 , intent(in) :: ao_num, mo_num integer*8 , intent(in) :: point_num double precision , intent(in) :: ao_value(ao_num,point_num) double precision , intent(in) :: coefficient_t(mo_num,ao_num) double precision , intent(out) :: mo_value(mo_num,point_num) integer*8 :: i,j,k double precision :: c1, c2, c3, c4, c5 integer*8 :: LDA, LDB, LDC info = QMCKL_SUCCESS if (.True.) then ! fast algorithm do j=1,point_num mo_value(:,j) = 0.d0 do k=1,ao_num c1 = ao_value(k,j) if (c1 /= 0.d0) then do i=1,mo_num mo_value(i,j) = mo_value(i,j) + coefficient_t(i,k) * c1 end do end if end do end do else ! dgemm for checking LDA = size(coefficient_t,1) LDB = size(ao_value,1) LDC = size(mo_value,1) info = qmckl_dgemm(context,'N', 'N', mo_num, point_num, ao_num, 1.d0, & coefficient_t, LDA, ao_value, LDB, & 0.d0, mo_value, LDC) end if end function qmckl_compute_mo_basis_mo_value_doc_f
qmckl_exit_code qmckl_compute_mo_basis_mo_value ( const qmckl_context context, const int64_t ao_num, const int64_t mo_num, const int64_t point_num, const double* coefficient_t, const double* ao_value, double* const mo_value );
qmckl_exit_code qmckl_compute_mo_basis_mo_value_doc ( const qmckl_context context, const int64_t ao_num, const int64_t mo_num, const int64_t point_num, const double* coefficient_t, const double* ao_value, double* const mo_value );
qmckl_exit_code qmckl_compute_mo_basis_mo_value (const qmckl_context context, const int64_t ao_num, const int64_t mo_num, const int64_t point_num, const double* coefficient_t, const double* ao_value, double* const mo_value ) { #ifdef HAVE_HPC return qmckl_compute_mo_basis_mo_value_hpc (context, ao_num, mo_num, point_num, coefficient_t, ao_value, mo_value); #else return qmckl_compute_mo_basis_mo_value_doc (context, ao_num, mo_num, point_num, coefficient_t, ao_value, mo_value); #endif }
2.1.4 HPC version
#ifdef HAVE_HPC qmckl_exit_code qmckl_compute_mo_basis_mo_value_hpc (const qmckl_context context, const int64_t ao_num, const int64_t mo_num, const int64_t point_num, const double* coefficient_t, const double* ao_value, double* const mo_value ); #endif
#ifdef HAVE_HPC qmckl_exit_code qmckl_compute_mo_basis_mo_value_hpc (const qmckl_context context, const int64_t ao_num, const int64_t mo_num, const int64_t point_num, const double* restrict coefficient_t, const double* restrict ao_value, double* restrict const mo_value ) { assert (context != QMCKL_NULL_CONTEXT); #ifdef HAVE_OPENMP #pragma omp parallel for #endif for (int64_t ipoint=0 ; ipoint < point_num ; ++ipoint) { double* restrict const vgl1 = &(mo_value[ipoint*mo_num]); const double* restrict avgl1 = &(ao_value[ipoint*ao_num]); for (int64_t i=0 ; i<mo_num ; ++i) { vgl1[i] = 0.; } int64_t nidx=0; int64_t idx[ao_num]; double av1[ao_num]; for (int64_t k=0 ; k<ao_num ; ++k) { if (avgl1[k] != 0.) { idx[nidx] = k; av1[nidx] = avgl1[k]; ++nidx; } } int64_t n=0; for (n=0 ; n < nidx-4 ; n+=4) { const double* restrict ck1 = coefficient_t + idx[n ]*mo_num; const double* restrict ck2 = coefficient_t + idx[n+1]*mo_num; const double* restrict ck3 = coefficient_t + idx[n+2]*mo_num; const double* restrict ck4 = coefficient_t + idx[n+3]*mo_num; const double a11 = av1[n ]; const double a21 = av1[n+1]; const double a31 = av1[n+2]; const double a41 = av1[n+3]; #ifdef HAVE_OPENMP #pragma omp simd #endif for (int64_t i=0 ; i<mo_num ; ++i) { vgl1[i] = vgl1[i] + ck1[i] * a11 + ck2[i] * a21 + ck3[i] * a31 + ck4[i] * a41; } } for (int64_t m=n ; m < nidx ; m+=1) { const double* restrict ck = coefficient_t + idx[m]*mo_num; const double a1 = av1[m]; #ifdef HAVE_OPENMP #pragma omp simd #endif for (int64_t i=0 ; i<mo_num ; ++i) { vgl1[i] += ck[i] * a1; } } } return QMCKL_SUCCESS; } #endif
2.2 Computation of MOs: values, gradient, Laplacian
2.2.1 Get
qmckl_exit_code qmckl_get_mo_basis_mo_vgl(qmckl_context context, double* const mo_vgl, const int64_t size_max);
Uses the given array to compute the VGL.
qmckl_exit_code qmckl_get_mo_basis_mo_vgl_inplace (qmckl_context context, double* const mo_vgl, const int64_t size_max);
2.2.2 Provide
qmckl_exit_code qmckl_provide_mo_basis_mo_vgl(qmckl_context context);
qmckl_exit_code qmckl_provide_mo_basis_mo_vgl(qmckl_context context) { qmckl_exit_code rc = QMCKL_SUCCESS; if (qmckl_context_check(context) == QMCKL_NULL_CONTEXT) { return qmckl_failwith( context, QMCKL_INVALID_CONTEXT, "qmckl_provide_mo_basis_mo_vgl", NULL); } qmckl_context_struct* const ctx = (qmckl_context_struct*) context; assert (ctx != NULL); if (!ctx->mo_basis.provided) { return qmckl_failwith( context, QMCKL_NOT_PROVIDED, "qmckl_provide_mo_basis_mo_vgl", NULL); } /* Compute if necessary */ if (ctx->point.date > ctx->mo_basis.mo_vgl_date) { qmckl_memory_info_struct mem_info = qmckl_memory_info_struct_zero; mem_info.size = 5 * ctx->mo_basis.mo_num * ctx->point.num * sizeof(double); if (ctx->mo_basis.mo_vgl != NULL) { qmckl_memory_info_struct mem_info_test = qmckl_memory_info_struct_zero; rc = qmckl_get_malloc_info(context, ctx->mo_basis.mo_vgl, &mem_info_test); /* if rc != QMCKL_SUCCESS, we are maybe in an _inplace function because the memory was not allocated with qmckl_malloc */ if ((rc == QMCKL_SUCCESS) && (mem_info_test.size != mem_info.size)) { rc = qmckl_free(context, ctx->mo_basis.mo_vgl); assert (rc == QMCKL_SUCCESS); ctx->mo_basis.mo_vgl = NULL; } } /* Allocate array */ if (ctx->mo_basis.mo_vgl == NULL) { double* mo_vgl = (double*) qmckl_malloc(context, mem_info); if (mo_vgl == NULL) { return qmckl_failwith( context, QMCKL_ALLOCATION_FAILED, "qmckl_mo_basis_mo_vgl", NULL); } ctx->mo_basis.mo_vgl = mo_vgl; }
if (rc != QMCKL_SUCCESS) { return rc; } ctx->mo_basis.mo_vgl_date = ctx->date; } return QMCKL_SUCCESS; }
2.2.3 Compute
Variable | Type | In/Out | Description |
---|---|---|---|
context |
qmckl_context |
in | Global state |
ao_num |
int64_t |
in | Number of AOs |
mo_num |
int64_t |
in | Number of MOs |
point_num |
int64_t |
in | Number of points |
coefficient_t |
double[mo_num][ao_num] |
in | Transpose of the AO to MO transformation matrix |
ao_vgl |
double[point_num][5][ao_num] |
in | Value, gradients and Laplacian of the AOs |
mo_vgl |
double[point_num][5][mo_num] |
out | Value, gradients and Laplacian of the MOs |
The matrix of AO values is very sparse, so we use a sparse-dense matrix multiplication instead of a dgemm, as exposed in https://dx.doi.org/10.1007/978-3-642-38718-0_14.
integer function qmckl_compute_mo_basis_mo_vgl_doc_f(context, & ao_num, mo_num, point_num, & coefficient_t, 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) :: point_num double precision , intent(in) :: ao_vgl(ao_num,5,point_num) double precision , intent(in) :: coefficient_t(mo_num,ao_num) double precision , intent(out) :: mo_vgl(mo_num,5,point_num) integer*8 :: i,j,k double precision :: c1, c2, c3, c4, c5 do j=1,point_num mo_vgl(:,:,j) = 0.d0 do k=1,ao_num if (ao_vgl(k,1,j) /= 0.d0) then c1 = ao_vgl(k,1,j) c2 = ao_vgl(k,2,j) c3 = ao_vgl(k,3,j) c4 = ao_vgl(k,4,j) c5 = ao_vgl(k,5,j) do i=1,mo_num mo_vgl(i,1,j) = mo_vgl(i,1,j) + coefficient_t(i,k) * c1 mo_vgl(i,2,j) = mo_vgl(i,2,j) + coefficient_t(i,k) * c2 mo_vgl(i,3,j) = mo_vgl(i,3,j) + coefficient_t(i,k) * c3 mo_vgl(i,4,j) = mo_vgl(i,4,j) + coefficient_t(i,k) * c4 mo_vgl(i,5,j) = mo_vgl(i,5,j) + coefficient_t(i,k) * c5 end do end if end do end do info = QMCKL_SUCCESS ! info = qmckl_dgemm(context,'N', 'N', mo_num, point_num, ao_num, 1.d0, & ! coefficient_t, int(size(coefficient_t,1),8), & ! ao_vgl, int(size(ao_vgl,1),8), 0.d0, & ! mo_vgl, int(size(mo_vgl,1),8)) end function qmckl_compute_mo_basis_mo_vgl_doc_f
qmckl_exit_code qmckl_compute_mo_basis_mo_vgl ( const qmckl_context context, const int64_t ao_num, const int64_t mo_num, const int64_t point_num, const double* coefficient_t, const double* ao_vgl, double* const mo_vgl );
qmckl_exit_code qmckl_compute_mo_basis_mo_vgl_doc ( const qmckl_context context, const int64_t ao_num, const int64_t mo_num, const int64_t point_num, const double* coefficient_t, const double* ao_vgl, double* const mo_vgl );
qmckl_exit_code qmckl_compute_mo_basis_mo_vgl (const qmckl_context context, const int64_t ao_num, const int64_t mo_num, const int64_t point_num, const double* coefficient_t, const double* ao_vgl, double* const mo_vgl ) { #ifdef HAVE_HPC return qmckl_compute_mo_basis_mo_vgl_hpc (context, ao_num, mo_num, point_num, coefficient_t, ao_vgl, mo_vgl); #else return qmckl_compute_mo_basis_mo_vgl_doc (context, ao_num, mo_num, point_num, coefficient_t, ao_vgl, mo_vgl); #endif }
2.2.4 HPC version
#ifdef HAVE_HPC qmckl_exit_code qmckl_compute_mo_basis_mo_vgl_hpc (const qmckl_context context, const int64_t ao_num, const int64_t mo_num, const int64_t point_num, const double* coefficient_t, const double* ao_vgl, double* const mo_vgl ); #endif
#ifdef HAVE_HPC qmckl_exit_code qmckl_compute_mo_basis_mo_vgl_hpc (const qmckl_context context, const int64_t ao_num, const int64_t mo_num, const int64_t point_num, const double* restrict coefficient_t, const double* restrict ao_vgl, double* restrict const mo_vgl ) { assert (context != QMCKL_NULL_CONTEXT); #ifdef HAVE_OPENMP #pragma omp parallel for #endif for (int64_t ipoint=0 ; ipoint < point_num ; ++ipoint) { double* restrict const vgl1 = &(mo_vgl[ipoint*5*mo_num]); double* restrict const vgl2 = vgl1 + mo_num; double* restrict const vgl3 = vgl1 + (mo_num << 1); double* restrict const vgl4 = vgl1 + (mo_num << 1) + mo_num; double* restrict const vgl5 = vgl1 + (mo_num << 2); const double* restrict avgl1 = &(ao_vgl[ipoint*5*ao_num]); const double* restrict avgl2 = avgl1 + ao_num; const double* restrict avgl3 = avgl1 + (ao_num << 1); const double* restrict avgl4 = avgl1 + (ao_num << 1) + ao_num; const double* restrict avgl5 = avgl1 + (ao_num << 2); for (int64_t i=0 ; i<mo_num ; ++i) { vgl1[i] = 0.; vgl2[i] = 0.; vgl3[i] = 0.; vgl4[i] = 0.; vgl5[i] = 0.; } int64_t nidx=0; int64_t idx[ao_num]; double av1[ao_num]; double av2[ao_num]; double av3[ao_num]; double av4[ao_num]; double av5[ao_num]; for (int64_t k=0 ; k<ao_num ; ++k) { if (avgl1[k] != 0.) { idx[nidx] = k; av1[nidx] = avgl1[k]; av2[nidx] = avgl2[k]; av3[nidx] = avgl3[k]; av4[nidx] = avgl4[k]; av5[nidx] = avgl5[k]; ++nidx; } } int64_t n=0; for (n=0 ; n < nidx-4 ; n+=4) { const double* restrict ck1 = coefficient_t + idx[n ]*mo_num; const double* restrict ck2 = coefficient_t + idx[n+1]*mo_num; const double* restrict ck3 = coefficient_t + idx[n+2]*mo_num; const double* restrict ck4 = coefficient_t + idx[n+3]*mo_num; const double a11 = av1[n ]; const double a21 = av1[n+1]; const double a31 = av1[n+2]; const double a41 = av1[n+3]; const double a12 = av2[n ]; const double a22 = av2[n+1]; const double a32 = av2[n+2]; const double a42 = av2[n+3]; const double a13 = av3[n ]; const double a23 = av3[n+1]; const double a33 = av3[n+2]; const double a43 = av3[n+3]; const double a14 = av4[n ]; const double a24 = av4[n+1]; const double a34 = av4[n+2]; const double a44 = av4[n+3]; const double a15 = av5[n ]; const double a25 = av5[n+1]; const double a35 = av5[n+2]; const double a45 = av5[n+3]; #ifdef HAVE_OPENMP #pragma omp simd #endif for (int64_t i=0 ; i<mo_num ; ++i) { vgl1[i] = vgl1[i] + ck1[i] * a11 + ck2[i] * a21 + ck3[i] * a31 + ck4[i] * a41; vgl2[i] = vgl2[i] + ck1[i] * a12 + ck2[i] * a22 + ck3[i] * a32 + ck4[i] * a42; vgl3[i] = vgl3[i] + ck1[i] * a13 + ck2[i] * a23 + ck3[i] * a33 + ck4[i] * a43; vgl4[i] = vgl4[i] + ck1[i] * a14 + ck2[i] * a24 + ck3[i] * a34 + ck4[i] * a44; vgl5[i] = vgl5[i] + ck1[i] * a15 + ck2[i] * a25 + ck3[i] * a35 + ck4[i] * a45; } } for (int64_t m=n ; m < nidx ; m+=1) { const double* restrict ck = coefficient_t + idx[m]*mo_num; const double a1 = av1[m]; const double a2 = av2[m]; const double a3 = av3[m]; const double a4 = av4[m]; const double a5 = av5[m]; #ifdef HAVE_OPENMP #pragma omp simd #endif for (int64_t i=0 ; i<mo_num ; ++i) { vgl1[i] += ck[i] * a1; vgl2[i] += ck[i] * a2; vgl3[i] += ck[i] * a3; vgl4[i] += ck[i] * a4; vgl5[i] += ck[i] * a5; } } } return QMCKL_SUCCESS; } #endif