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] |
Orbital coefficients |
Computed data:
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 * 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; }
1.2 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 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.
2 Computation
2.1 Computation of MOs
2.1.1 Get
qmckl_exit_code qmckl_get_mo_basis_vgl(qmckl_context context, double* const mo_vgl);
2.1.2 Provide
2.1.3 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 |
point_num |
in | Number of points |
double |
coef_normalized[mo_num][ao_num] |
in | AO to MO transformation matrix |
double |
ao_vgl[point_num][5][ao_num] |
in | Value, gradients and Laplacian of the AOs |
double |
mo_vgl[point_num][5][mo_num] |
out | Value, gradients and Laplacian of the MOs |
integer function qmckl_compute_mo_basis_vgl_f(context, & ao_num, mo_num, point_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) :: point_num double precision , intent(in) :: ao_vgl(ao_num,5,point_num) double precision , intent(in) :: coef_normalized(ao_num,mo_num) double precision , intent(out) :: mo_vgl(mo_num,5,point_num) character :: TransA, TransB double precision :: alpha, beta integer*8 :: M, N, K, LDA, LDB, LDC, i,j TransA = 'T' TransB = 'N' M = mo_num N = point_num*5_8 K = int(ao_num,8) alpha = 1.d0 beta = 0.d0 LDA = size(coef_normalized,1) LDB = size(ao_vgl,1) LDC = size(mo_vgl,1) info = qmckl_dgemm(context,TransA, TransB, M, N, K, alpha, & coef_normalized, int(size(coef_normalized,1),8), & ao_vgl, int(size(ao_vgl,1),8), beta, & mo_vgl,LDC) info = QMCKL_SUCCESS end function qmckl_compute_mo_basis_vgl_f
qmckl_exit_code qmckl_compute_mo_basis_vgl ( const context qmckl_context, const ao_num int64_t, const mo_num int64_t, const point_num int64_t, const coef_normalized* double, const ao_vgl* double, mo_vgl* const double );