PyObject* r2k(PyObject *self, PyObject *args) { Py_complex alpha; PyArrayObject* a; PyArrayObject* b; double beta; PyArrayObject* c; if (!PyArg_ParseTuple(args, "DOOdO", &alpha, &a, &b, &beta, &c)) return NULL; int n = PyArray_DIMS(a)[0]; int k = PyArray_DIMS(a)[1]; for (int d = 2; d < PyArray_NDIM(a); d++) k *= PyArray_DIMS(a)[d]; int ldc = PyArray_STRIDES(c)[0] / PyArray_STRIDES(c)[1]; if (PyArray_DESCR(a)->type_num == NPY_DOUBLE) dsyr2k_("u", "t", &n, &k, (double*)(&alpha), DOUBLEP(a), &k, DOUBLEP(b), &k, &beta, DOUBLEP(c), &ldc); else zher2k_("u", "c", &n, &k, (void*)(&alpha), (void*)COMPLEXP(a), &k, (void*)COMPLEXP(b), &k, &beta, (void*)COMPLEXP(c), &ldc); Py_RETURN_NONE; }
PyObject* dotu(PyObject *self, PyObject *args) { PyArrayObject* a; PyArrayObject* b; if (!PyArg_ParseTuple(args, "OO", &a, &b)) return NULL; int n = PyArray_DIMS(a)[0]; for (int i = 1; i < PyArray_NDIM(a); i++) n *= PyArray_DIMS(a)[i]; int incx = 1; int incy = 1; if (PyArray_DESCR(a)->type_num == NPY_DOUBLE) { double result; result = ddot_(&n, (void*)DOUBLEP(a), &incx, (void*)DOUBLEP(b), &incy); return PyFloat_FromDouble(result); } else { double_complex* ap = COMPLEXP(a); double_complex* bp = COMPLEXP(b); double_complex z = 0.0; for (int i = 0; i < n; i++) z += ap[i] * bp[i]; return PyComplex_FromDoubles(creal(z), cimag(z)); } }
PyObject* scalapack_redist(PyObject *self, PyObject *args) { PyArrayObject* a; // source matrix PyArrayObject* b; // destination matrix PyArrayObject* desca; // source descriptor PyArrayObject* descb; // destination descriptor char uplo; char diag='N'; // copy the diagonal int c_ConTxt; int m; int n; int ia, ja, ib, jb; if (!PyArg_ParseTuple(args, "OOOOiiiiiiic", &desca, &descb, &a, &b, &m, &n, &ia, &ja, &ib, &jb, &c_ConTxt, &uplo)) return NULL; if (uplo == 'G') // General matrix { if (a->descr->type_num == PyArray_DOUBLE) Cpdgemr2d_(m, n, DOUBLEP(a), ia, ja, INTP(desca), DOUBLEP(b), ib, jb, INTP(descb), c_ConTxt); else Cpzgemr2d_(m, n, (void*)COMPLEXP(a), ia, ja, INTP(desca), (void*)COMPLEXP(b), ib, jb, INTP(descb), c_ConTxt); } else // Trapezoidal matrix { if (a->descr->type_num == PyArray_DOUBLE) Cpdtrmr2d_(&uplo, &diag, m, n, DOUBLEP(a), ia, ja, INTP(desca), DOUBLEP(b), ib, jb, INTP(descb), c_ConTxt); else Cpztrmr2d_(&uplo, &diag, m, n, (void*)COMPLEXP(a), ia, ja, INTP(desca), (void*)COMPLEXP(b), ib, jb, INTP(descb), c_ConTxt); } Py_RETURN_NONE; }
PyObject* gemv(PyObject *self, PyObject *args) { Py_complex alpha; PyArrayObject* a; PyArrayObject* x; Py_complex beta; PyArrayObject* y; char trans = 't'; if (!PyArg_ParseTuple(args, "DOODO|c", &alpha, &a, &x, &beta, &y, &trans)) return NULL; int m, n, lda, itemsize, incx, incy; if (trans == 'n') { m = PyArray_DIMS(a)[1]; for (int i = 2; i < PyArray_NDIM(a); i++) m *= PyArray_DIMS(a)[i]; n = PyArray_DIMS(a)[0]; lda = MAX(1, m); } else { n = PyArray_DIMS(a)[0]; for (int i = 1; i < PyArray_NDIM(a)-1; i++) n *= PyArray_DIMS(a)[i]; m = PyArray_DIMS(a)[PyArray_NDIM(a)-1]; lda = MAX(1, m); } if (PyArray_DESCR(a)->type_num == NPY_DOUBLE) itemsize = sizeof(double); else itemsize = sizeof(double_complex); incx = PyArray_STRIDES(x)[0]/itemsize; incy = 1; if (PyArray_DESCR(a)->type_num == NPY_DOUBLE) dgemv_(&trans, &m, &n, &(alpha.real), DOUBLEP(a), &lda, DOUBLEP(x), &incx, &(beta.real), DOUBLEP(y), &incy); else zgemv_(&trans, &m, &n, &alpha, (void*)COMPLEXP(a), &lda, (void*)COMPLEXP(x), &incx, &beta, (void*)COMPLEXP(y), &incy); Py_RETURN_NONE; }
PyObject* gemm(PyObject *self, PyObject *args) { Py_complex alpha; PyArrayObject* a; PyArrayObject* b; Py_complex beta; PyArrayObject* c; char transa = 'n'; if (!PyArg_ParseTuple(args, "DOODO|c", &alpha, &a, &b, &beta, &c, &transa)) return NULL; int m, k, lda, ldb, ldc; if (transa == 'n') { m = PyArray_DIMS(a)[1]; for (int i = 2; i < PyArray_NDIM(a); i++) m *= PyArray_DIMS(a)[i]; k = PyArray_DIMS(a)[0]; lda = MAX(1, PyArray_STRIDES(a)[0] / PyArray_STRIDES(a)[PyArray_NDIM(a) - 1]); ldb = MAX(1, PyArray_STRIDES(b)[0] / PyArray_STRIDES(b)[1]); ldc = MAX(1, PyArray_STRIDES(c)[0] / PyArray_STRIDES(c)[PyArray_NDIM(c) - 1]); } else { k = PyArray_DIMS(a)[1]; for (int i = 2; i < PyArray_NDIM(a); i++) k *= PyArray_DIMS(a)[i]; m = PyArray_DIMS(a)[0]; lda = MAX(1, k); ldb = MAX(1, PyArray_STRIDES(b)[0] / PyArray_STRIDES(b)[PyArray_NDIM(b) - 1]); ldc = MAX(1, PyArray_STRIDES(c)[0] / PyArray_STRIDES(c)[1]); } int n = PyArray_DIMS(b)[0]; if (PyArray_DESCR(a)->type_num == NPY_DOUBLE) dgemm_(&transa, "n", &m, &n, &k, &(alpha.real), DOUBLEP(a), &lda, DOUBLEP(b), &ldb, &(beta.real), DOUBLEP(c), &ldc); else zgemm_(&transa, "n", &m, &n, &k, &alpha, (void*)COMPLEXP(a), &lda, (void*)COMPLEXP(b), &ldb, &beta, (void*)COMPLEXP(c), &ldc); Py_RETURN_NONE; }
PyObject* multi_axpy(PyObject *self, PyObject *args) { PyArrayObject* alpha; PyArrayObject* x; PyArrayObject* y; if (!PyArg_ParseTuple(args, "OOO", &alpha, &x, &y)) return NULL; int n0 = PyArray_DIMS(x)[0]; int n = PyArray_DIMS(x)[1]; for (int d = 2; d < PyArray_NDIM(x); d++) n *= PyArray_DIMS(x)[d]; int incx = 1; int incy = 1; if (PyArray_DESCR(alpha)->type_num == NPY_DOUBLE) { if (PyArray_DESCR(x)->type_num == NPY_CDOUBLE) n *= 2; double *ap = DOUBLEP(alpha); double *xp = DOUBLEP(x); double *yp = DOUBLEP(y); for (int i = 0; i < n0; i++) { daxpy_(&n, &ap[i], (void*)xp, &incx, (void*)yp, &incy); xp += n; yp += n; } } else { double_complex *ap = COMPLEXP(alpha); double_complex *xp = COMPLEXP(x); double_complex *yp = COMPLEXP(y); for (int i = 0; i < n0; i++) { zaxpy_(&n, (void*)(&ap[i]), (void*)xp, &incx, (void*)yp, &incy); xp += n; yp += n; } } Py_RETURN_NONE; }
PyObject *exterior_electron_density_region(PyObject *self, PyObject *args) { PyArrayObject* ai; PyArrayObject* aatom_c; PyArrayObject* beg_c; PyArrayObject* end_c; PyArrayObject* hh_c; PyArrayObject* vdWrad; if (!PyArg_ParseTuple(args, "OOOOOO", &ai, &aatom_c, &beg_c, &end_c, &hh_c, &vdWrad)) return NULL; long *aindex = LONGP(ai); int natoms = PyArray_DIM(aatom_c, 0); double *atom_c = DOUBLEP(aatom_c); long *beg = LONGP(beg_c); long *end = LONGP(end_c); double *h_c = DOUBLEP(hh_c); double *vdWradius = DOUBLEP(vdWrad); int n[3], ij; double pos[3]; for (int c = 0; c < 3; c++) { n[c] = end[c] - beg[c]; } // loop over all points for (int i = 0; i < n[0]; i++) { pos[0] = (beg[0] + i) * h_c[0]; for (int j = 0; j < n[1]; j++) { pos[1] = (beg[1] + j) * h_c[1]; ij = (i*n[1] + j)*n[2]; for (int k = 0; k < n[2]; k++) { pos[2] = (beg[2] + k) * h_c[2]; aindex[ij + k] = (long) 1; /* assume outside the structure */ // loop over all atoms for (int a=0; a < natoms; a++) { double d = distance(atom_c + a*3, pos); if (d < vdWradius[a]) { aindex[ij + k] = (long) 0; /* this is inside */ a = natoms; } } } } } Py_RETURN_NONE; }
PyObject* multi_dotu(PyObject *self, PyObject *args) { PyArrayObject* a; PyArrayObject* b; PyArrayObject* c; if (!PyArg_ParseTuple(args, "OOO", &a, &b, &c)) return NULL; int n0 = PyArray_DIMS(a)[0]; int n = PyArray_DIMS(a)[1]; for (int i = 2; i < PyArray_NDIM(a); i++) n *= PyArray_DIMS(a)[i]; int incx = 1; int incy = 1; if (PyArray_DESCR(a)->type_num == NPY_DOUBLE) { double *ap = DOUBLEP(a); double *bp = DOUBLEP(b); double *cp = DOUBLEP(c); for (int i = 0; i < n0; i++) { cp[i] = ddot_(&n, (void*)ap, &incx, (void*)bp, &incy); ap += n; bp += n; } } else { double_complex* ap = COMPLEXP(a); double_complex* bp = COMPLEXP(b); double_complex* cp = COMPLEXP(c); for (int i = 0; i < n0; i++) { cp[i] = 0.0; for (int j = 0; j < n; j++) cp[i] += ap[j] * bp[j]; ap += n; bp += n; } } Py_RETURN_NONE; }
PyObject* pblas_gemm(PyObject *self, PyObject *args) { char transa; char transb; int m, n, k; Py_complex alpha; Py_complex beta; PyArrayObject *a, *b, *c; PyArrayObject *desca, *descb, *descc; int one = 1; if (!PyArg_ParseTuple(args, "iiiDOODOOOOcc", &m, &n, &k, &alpha, &a, &b, &beta, &c, &desca, &descb, &descc, &transa, &transb)) { return NULL; } // cdesc // int c_ConTxt = INTP(descc)[1]; // If process not on BLACS grid, then return. // if (c_ConTxt == -1) Py_RETURN_NONE; if (c->descr->type_num == PyArray_DOUBLE) pdgemm_(&transa, &transb, &m, &n, &k, &(alpha.real), DOUBLEP(a), &one, &one, INTP(desca), DOUBLEP(b), &one, &one, INTP(descb), &(beta.real), DOUBLEP(c), &one, &one, INTP(descc)); else pzgemm_(&transa, &transb, &m, &n, &k, &alpha, (void*)COMPLEXP(a), &one, &one, INTP(desca), (void*)COMPLEXP(b), &one, &one, INTP(descb), &beta, (void*)COMPLEXP(c), &one, &one, INTP(descc)); Py_RETURN_NONE; }
PyObject* pblas_gemv(PyObject *self, PyObject *args) { char transa; int m, n; Py_complex alpha; Py_complex beta; PyArrayObject *a, *x, *y; int incx = 1, incy = 1; // what should these be? PyArrayObject *desca, *descx, *descy; int one = 1; if (!PyArg_ParseTuple(args, "iiDOODOOOOc", &m, &n, &alpha, &a, &x, &beta, &y, &desca, &descx, &descy, &transa)) { return NULL; } // ydesc // int y_ConTxt = INTP(descy)[1]; // If process not on BLACS grid, then return. // if (y_ConTxt == -1) Py_RETURN_NONE; if (y->descr->type_num == PyArray_DOUBLE) pdgemv_(&transa, &m, &n, &(alpha.real), DOUBLEP(a), &one, &one, INTP(desca), DOUBLEP(x), &one, &one, INTP(descx), &incx, &(beta.real), DOUBLEP(y), &one, &one, INTP(descy), &incy); else pzgemv_(&transa, &m, &n, &alpha, (void*)COMPLEXP(a), &one, &one, INTP(desca), (void*)COMPLEXP(x), &one, &one, INTP(descx), &incx, &beta, (void*)COMPLEXP(y), &one, &one, INTP(descy), &incy); Py_RETURN_NONE; }
/* unpolarised */ PyObject* d2Excdn2(PyObject *self, PyObject *args) { PyArrayObject* den; /* density */ PyArrayObject* res; /* derivative */ if (!PyArg_ParseTuple(args, "OO", &den, &res)) return NULL; int n = den->dimensions[0]; /* printf("<d2Excdn2> nd=%d\n",den->nd); */ for (int d = 1; d < den->nd; d++) n *= den->dimensions[d]; double *denp = DOUBLEP(den); double *resp = DOUBLEP(res); for (int i=0; i<n;i++) { resp[i] = d2exdn2_(denp[i]) + d2ecdrho2_u__(denp+i); } Py_RETURN_NONE; }
PyObject* axpy(PyObject *self, PyObject *args) { Py_complex alpha; PyArrayObject* x; PyArrayObject* y; if (!PyArg_ParseTuple(args, "DOO", &alpha, &x, &y)) return NULL; int n = PyArray_DIMS(x)[0]; for (int d = 1; d < PyArray_NDIM(x); d++) n *= PyArray_DIMS(x)[d]; int incx = 1; int incy = 1; if (PyArray_DESCR(x)->type_num == NPY_DOUBLE) daxpy_(&n, &(alpha.real), DOUBLEP(x), &incx, DOUBLEP(y), &incy); else zaxpy_(&n, &alpha, (void*)COMPLEXP(x), &incx, (void*)COMPLEXP(y), &incy); Py_RETURN_NONE; }
PyObject* pblas_rk(PyObject *self, PyObject *args) { char uplo; int n, k; Py_complex alpha; Py_complex beta; PyArrayObject *a, *c; PyArrayObject *desca, *descc; int one = 1; if (!PyArg_ParseTuple(args, "iiDODOOOc", &n, &k, &alpha, &a, &beta, &c, &desca, &descc, &uplo)) { return NULL; } // cdesc // int c_ConTxt = INTP(descc)[1]; // If process not on BLACS grid, then return. // if (c_ConTxt == -1) Py_RETURN_NONE; if (c->descr->type_num == PyArray_DOUBLE) pdsyrk_(&uplo, "T", &n, &k, &(alpha.real), DOUBLEP(a), &one, &one, INTP(desca), &(beta.real), DOUBLEP(c), &one, &one, INTP(descc)); else pzherk_(&uplo, "C", &n, &k, &alpha, (void*)COMPLEXP(a), &one, &one, INTP(desca), &beta, (void*)COMPLEXP(c), &one, &one, INTP(descc)); Py_RETURN_NONE; }
/* polarised */ PyObject* d2Excdnsdnt(PyObject *self, PyObject *args) { PyArrayObject* dup; /* "up" density */ PyArrayObject* ddn; /* "down" density */ int is; /* i spin (0,1) */ int ks; /* k spin (0,1) */ PyArrayObject* res; /* derivative */ /* printf("<d2Excdnsdnt> is=%p ks=%p\n",is,ks); */ if (!PyArg_ParseTuple(args, "OOiiO", &dup, &ddn, &is, &ks, &res)) return NULL; /* printf("<d2Excdnsdnt> args passed\n"); */ /* printf("<d2Excdnsdnt> is=%d ks=%d\n",is,ks); */ int n = dup->dimensions[0]; for (int d = 1; d < dup->nd; d++) n *= dup->dimensions[d]; /* printf("<d2Excdnsdnt> n=%d\n",n); */ double *dupp = DOUBLEP(dup); double *ddnp = DOUBLEP(ddn); /* printf("<d2Excdnsdnt> dupp=%p ddnp=%p\n",dupp,ddnp); */ double *resp = DOUBLEP(res); /* int iis = *is */ /* int iks = *ks; */ /* printf("<d2Excdnsdnt> is=%i ks=%i\n",iis,iks); */ for (int i=0; i<n;i++) { resp[i] = d2exdnsdnt_(dupp+i,ddnp+i,&is,&ks) + d2ecdnsdnt_(dupp+i,ddnp+i,&is,&ks); /* printf("<d2Excdnsdnt> i=%d dup=%g ddn=%g resp=%g\n",i, */ /* dupp[i],ddnp[i],resp[i]); */ } Py_RETURN_NONE; }
PyObject* rk(PyObject *self, PyObject *args) { double alpha; PyArrayObject* a; double beta; PyArrayObject* c; char trans = 'c'; if (!PyArg_ParseTuple(args, "dOdO|c", &alpha, &a, &beta, &c, &trans)) return NULL; int n = PyArray_DIMS(c)[0]; int k, lda; if (trans == 'c') { k = PyArray_DIMS(a)[1]; for (int d = 2; d < PyArray_NDIM(a); d++) k *= PyArray_DIMS(a)[d]; lda = k; } else { k = PyArray_DIMS(a)[0]; lda = n; } int ldc = PyArray_STRIDES(c)[0] / PyArray_STRIDES(c)[1]; if (PyArray_DESCR(a)->type_num == NPY_DOUBLE) dsyrk_("u", &trans, &n, &k, &alpha, DOUBLEP(a), &lda, &beta, DOUBLEP(c), &ldc); else zherk_("u", &trans, &n, &k, &alpha, (void*)COMPLEXP(a), &lda, &beta, (void*)COMPLEXP(c), &ldc); Py_RETURN_NONE; }
static PyObject* particles_set_cell(particles_t *self, PyObject *args) { PyObject *Abox_obj, *pbc_obj = NULL; PyArrayObject *Abox_arr; PyArrayObject *pbc_arr = NULL; double *Abox; npy_bool *pbc; BOOL pbc_for[3]; int ierror = ERROR_NONE; if (!PyArg_ParseTuple(args, "O|O", &Abox_obj, &pbc_obj)) return NULL; Abox_arr = (PyArrayObject *) PyArray_FROMANY(Abox_obj, NPY_DOUBLE, 2, 2, NPY_C_CONTIGUOUS); if (!Abox_arr) return NULL; Abox = DOUBLEP(Abox_arr); pbc_for[0] = 1; pbc_for[1] = 1; pbc_for[2] = 1; if (pbc_obj) { pbc_arr = (PyArrayObject *) PyArray_FROMANY(pbc_obj, NPY_BOOL, 1, 1, NPY_C_CONTIGUOUS); if (!pbc_arr) return NULL; pbc = (npy_bool *) BOOLP(pbc_arr); pbc_for[0] = pbc[0]; pbc_for[1] = pbc[1]; pbc_for[2] = pbc[2]; } #ifdef DEBUG printf("[particles_set_cell] pbc_for %i, %i, %i\n", pbc_for[0], pbc_for[1], pbc_for[2]); #endif f_particles_set_cell(self->f90obj, Abox, pbc_for, &ierror); if (error_to_py(ierror)) return NULL; Py_RETURN_NONE; }
PyObject* scalapack_set(PyObject *self, PyObject *args) { PyArrayObject* a; // matrix; PyArrayObject* desca; // descriptor Py_complex alpha; Py_complex beta; int m, n; int ia, ja; char uplo; if (!PyArg_ParseTuple(args, "OODDciiii", &a, &desca, &alpha, &beta, &uplo, &m, &n, &ia, &ja)) return NULL; if (a->descr->type_num == PyArray_DOUBLE) pdlaset_(&uplo, &m, &n, &(alpha.real), &(beta.real), DOUBLEP(a), &ia, &ja, INTP(desca)); else pzlaset_(&uplo, &m, &n, &alpha, &beta, (void*)COMPLEXP(a), &ia, &ja, INTP(desca)); Py_RETURN_NONE; }
PyObject* scalapack_diagonalize_dc(PyObject *self, PyObject *args) { // Standard driver for divide and conquer algorithm // Computes all eigenvalues and eigenvectors PyArrayObject* a; // symmetric matrix PyArrayObject* desca; // symmetric matrix description vector PyArrayObject* z; // eigenvector matrix PyArrayObject* w; // eigenvalue array int one = 1; char jobz = 'V'; // eigenvectors also char uplo; if (!PyArg_ParseTuple(args, "OOcOO", &a, &desca, &uplo, &z, &w)) return NULL; // adesc // int a_ConTxt = INTP(desca)[1]; int a_m = INTP(desca)[2]; int a_n = INTP(desca)[3]; // zdesc = adesc; this can be relaxed a bit according to pdsyevd.f // Only square matrices assert (a_m == a_n); int n = a_n; // If process not on BLACS grid, then return. // if (a_ConTxt == -1) Py_RETURN_NONE; // Query part, need to find the optimal size of a number of work arrays int info; int querywork = -1; int* iwork; int liwork; int lwork; int lrwork; int i_work; double d_work; double_complex c_work; if (a->descr->type_num == PyArray_DOUBLE) { pdsyevd_(&jobz, &uplo, &n, DOUBLEP(a), &one, &one, INTP(desca), DOUBLEP(w), DOUBLEP(z), &one, &one, INTP(desca), &d_work, &querywork, &i_work, &querywork, &info); lwork = (int)(d_work); } else { pzheevd_(&jobz, &uplo, &n, (void*)COMPLEXP(a), &one, &one, INTP(desca), DOUBLEP(w), (void*)COMPLEXP(z), &one, &one, INTP(desca), (void*)&c_work, &querywork, &d_work, &querywork, &i_work, &querywork, &info); lwork = (int)(c_work); lrwork = (int)(d_work); } if (info != 0) { PyErr_SetString(PyExc_RuntimeError, "scalapack_diagonalize_dc error in query."); return NULL; } // Computation part liwork = i_work; iwork = GPAW_MALLOC(int, liwork); if (a->descr->type_num == PyArray_DOUBLE) { double* work = GPAW_MALLOC(double, lwork); pdsyevd_(&jobz, &uplo, &n, DOUBLEP(a), &one, &one, INTP(desca), DOUBLEP(w), DOUBLEP(z), &one, &one, INTP(desca), work, &lwork, iwork, &liwork, &info); free(work); }
PyObject * vdw(PyObject* self, PyObject *args) { PyArrayObject* n_obj; PyArrayObject* q0_obj; PyArrayObject* R_obj; PyArrayObject* cell_obj; PyArrayObject* pbc_obj; PyArrayObject* repeat_obj; PyArrayObject* phi_obj; double ddelta; double dD; int iA; int iB; PyArrayObject* rhistogram_obj; double drhist; PyArrayObject* Dhistogram_obj; double dDhist; if (!PyArg_ParseTuple(args, "OOOOOOOddiiOdOd", &n_obj, &q0_obj, &R_obj, &cell_obj, &pbc_obj, &repeat_obj, &phi_obj, &ddelta, &dD, &iA, &iB, &rhistogram_obj, &drhist, &Dhistogram_obj, &dDhist)) return NULL; int ndelta = PyArray_DIMS(phi_obj)[0]; int nD = PyArray_DIMS(phi_obj)[1]; const double* n = (const double*)DOUBLEP(n_obj); const int ni = PyArray_SIZE(n_obj); const double* q0 = (const double*)DOUBLEP(q0_obj); const double (*R)[3] = (const double (*)[3])DOUBLEP(R_obj); const double* cell = (const double*)DOUBLEP(cell_obj); const char* pbc = (const char*)(PyArray_DATA(pbc_obj)); const long* repeat = (const long*)(PyArray_DATA(repeat_obj)); const double (*phi)[nD] = (const double (*)[nD])DOUBLEP(phi_obj); double* rhistogram = (double*)DOUBLEP(rhistogram_obj); double* Dhistogram = (double*)DOUBLEP(Dhistogram_obj); int nbinsr = PyArray_DIMS(rhistogram_obj)[0]; int nbinsD = PyArray_DIMS(Dhistogram_obj)[0]; double energy = 0.0; if (repeat[0] == 0 && repeat[1] == 0 && repeat[2] == 0) for (int i1 = iA; i1 < iB; i1++) { const double* R1 = R[i1]; double q01 = q0[i1]; for (int i2 = 0; i2 <= i1; i2++) { double rr = 0.0; for (int c = 0; c < 3; c++) { double f = R[i2][c] - R1[c]; if (pbc[c]) f = fmod(f + 1.5 * cell[c], cell[c]) - 0.5 * cell[c]; rr += f * f; } double r = sqrt(rr); double d1 = r * q01; double d2 = r * q0[i2]; double D = 0.5 * (d1 + d2); double e12 = (vdwkernel(D, d1, d2, nD, ndelta, dD, ddelta, phi) * n[i1] * n[i2]); if (i1 == i2) e12 /= 2.0; int bin = (int)(r / drhist); if (bin < nbinsr) rhistogram[bin] += e12; bin = (int)(D / dDhist); if (bin < nbinsD) Dhistogram[bin] += e12; energy += e12; } } else for (int i1 = iA; i1 < iB; i1++) { const double* R1 = R[i1]; double q01 = q0[i1]; for (int a1 = -repeat[0]; a1 <= repeat[0]; a1++) for (int a2 = -repeat[1]; a2 <= repeat[1]; a2++) for (int a3 = -repeat[2]; a3 <= repeat[2]; a3++) { double x = 0.5; int i2max = ni-1; if (a1 == 0 && a2 == 0 && a3 == 0) { i2max = i1; x = 1.0; } double R1a[3] = {R1[0] + a1 * cell[0], R1[1] + a2 * cell[1], R1[2] + a3 * cell[2]}; for (int i2 = 0; i2 <= i2max; i2++) { double rr = 0.0; for (int c = 0; c < 3; c++) { double f = R[i2][c] - R1a[c]; rr += f * f; } double r = sqrt(rr); double d1 = r * q01; double d2 = r * q0[i2]; double D = 0.5 * (d1 + d2); double e12 = (vdwkernel(D, d1, d2, nD, ndelta, dD, ddelta, phi) * n[i1] * n[i2] * x); int bin = (int)(r / drhist); if (bin < nbinsr) rhistogram[bin] += e12; bin = (int)(D / dDhist); if (bin < nbinsD) Dhistogram[bin] += e12; energy += e12; } } } return PyFloat_FromDouble(energy); }
// Perform a moving linear least squares interpolation to arrays // Input arguments: // order: order of polynomial used (1 or 2) // cutoff: the cutoff of weight (in grid points) // coords: scaled coords [0,1] for interpolation // N_c: number of grid points // beg_c: first grid point // data: the array used // target: the results are stored in this array PyObject* mlsqr(PyObject *self, PyObject *args) { // The order of interpolation unsigned char order = -1; // The cutoff for moving least squares double cutoff = -1; // The coordinates for interpolation: array of size (3, N) PyArrayObject* coords = 0; // Number of grid points PyArrayObject* N_c = 0; // Beginning of grid PyArrayObject* beg_c = 0; // The 3d-data to be interpolated: array of size (X, Y, Z) PyArrayObject* data; // The interpolation target: array of size (N,) PyArrayObject* target = 0; if (!PyArg_ParseTuple(args, "BdOOOOO", &order, &cutoff, &coords, &N_c, &beg_c, &data, &target)) { return NULL; } int coeffs = -1; if (order == 1) { coeffs = 4; } if (order == 2) { coeffs = 10; // 1 x y z xy yz zx xx yy zz } if (order == 3) { // 1 x y z xy yz zx xx yy zz // xxy xxz yyx yyz zzx zzy // xxx yyy zzz zyz coeffs = 20; } int points = coords->dimensions[0]; double* coord_nc = DOUBLEP(coords); double* grid_points = DOUBLEP(N_c); double* grid_start = DOUBLEP(beg_c); double* target_n = DOUBLEP(target); double* data_g = DOUBLEP(data); // TODO: Calculate fit const int sizex = ceil(cutoff); const int sizey = ceil(cutoff); const int sizez = ceil(cutoff); // Allocate X-matrix and b-vector int source_points = (2*sizex+1)*(2*sizey+1)*(2*sizez+1); double* X = GPAW_MALLOC(double, coeffs*source_points); double* b = GPAW_MALLOC(double, source_points); double* work = GPAW_MALLOC(double, coeffs*source_points); // The multipliers for each dimension int ldx = data->dimensions[1]*data->dimensions[2]; int ldy = data->dimensions[2]; int ldz = 1; // For each point to be interpolated for (int p=0; p< points; p++) { double x = (*coord_nc++)*grid_points[0] - grid_start[0]; double y = (*coord_nc++)*grid_points[1] - grid_start[0]; double z = (*coord_nc++)*grid_points[2] - grid_start[0]; // The grid center point int cx2 = round(x); int cy2 = round(y); int cz2 = round(z); // Scaled to grid int cx = safemod(cx2,data->dimensions[0]); int cy = safemod(cy2,data->dimensions[1]); int cz = safemod(cz2,data->dimensions[2]); double* i_X = X; double* i_b = b; // For each point to take into account for (int dx=-sizex;dx<=sizex;dx++) for (int dy=-sizey;dy<=sizey;dy++) for (int dz=-sizez;dz<=sizez;dz++) { // Coordinates centered on x,y,z double sx = (cx2 + dx) - x; double sy = (cy2 + dy) - y; double sz = (cz2 + dz) - z; // Normalized distance from center double d = sqrt(sx*sx+sy*sy+sz*sz) / cutoff; double w = 0.0; if (d < 1) { w = (1-d)*(1-d); w*=w; w*=(4*d+1); } //double w = exp(-d*d); *i_X++ = w*1.0; *i_X++ = w*sx; *i_X++ = w*sy; *i_X++ = w*sz; if (order > 1) { *i_X++ = w*sx*sy; *i_X++ = w*sy*sz; *i_X++ = w*sz*sx; *i_X++ = w*sx*sx; *i_X++ = w*sy*sy; *i_X++ = w*sz*sz; } if (order > 2) { *i_X++ = w*sx*sy*sz; // xyz *i_X++ = w*sx*sx*sx; // xxx *i_X++ = w*sy*sy*sy; // yyy *i_X++ = w*sz*sz*sz; // zzz *i_X++ = w*sx*sx*sy; // xxy *i_X++ = w*sx*sx*sz; // xxz *i_X++ = w*sy*sy*sx; // yyx *i_X++ = w*sy*sy*sz; // yyz *i_X++ = w*sz*sz*sx; // zzx *i_X++ = w*sz*sz*sy; // zzy } *i_b++ = w*data_g[ safemod(cx+dx, data->dimensions[0]) * ldx + safemod(cy+dy, data->dimensions[1]) * ldy + safemod(cz+dz, data->dimensions[2]) * ldz ]; } int info = 0; int rhs = 1; int worksize = coeffs*source_points; int ldb = source_points; dgels_("T", &coeffs, // ...times 4. &source_points, // lhs is of size sourcepoints... &rhs, // one rhs. X, // provide lhs &coeffs, // Leading dimension of X b, // provide rhs &ldb, // Leading dimension of b work, // work array (and output) &worksize, // the size of work array &info); // info if (info != 0) printf("WARNING: dgels returned %d!", info); // Evaluate the polynomial // Due to centered coordinates, it's just the constant term double value = b[0]; *target_n++ = value; //Nearest neighbour //double value = data_g[ cx*data->dimensions[1]*data->dimensions[2] + cy*data->dimensions[2] + cz ]; //printf("%.5f" , value); } free(work); free(b); free(X); Py_RETURN_NONE; }
static PyObject* lxcXCFunctional_CalculateSpinPaired(lxcXCFunctionalObject *self, PyObject *args) { PyArrayObject* e_array; /* energy per particle*/ PyArrayObject* n_array; /* rho */ PyArrayObject* v_array; /* dE/drho */ PyArrayObject* a2_array = 0; /* |nabla rho|^2*/ PyArrayObject* dEda2_array = 0; /* dE/d|nabla rho|^2 */ PyArrayObject* tau_array = 0; /* tau*/ PyArrayObject* dEdtau_array = 0; /* dE/dtau */ if (!PyArg_ParseTuple(args, "OOO|OOOO", &e_array, &n_array, &v_array, /* object | optional objects*/ &a2_array, &dEda2_array, &tau_array, &dEdtau_array)) return NULL; /* find nspin */ int nspin = self->nspin; assert(nspin == XC_UNPOLARIZED); /* we are spinpaired */ int ng = e_array->dimensions[0]; /* number of grid points */ double* e_g = DOUBLEP(e_array); /* e on the grid */ const double* n_g = DOUBLEP(n_array); /* density on the grid */ double* v_g = DOUBLEP(v_array); /* v on the grid */ const double* a2_g = 0; /* a2 on the grid */ double* tau_g = 0; /* tau on the grid */ double* dEda2_g = 0; /* dEda2 on the grid */ double* dEdtau_g= 0; /* dEdt on the grid */ if (((self->x_functional.family == XC_FAMILY_GGA) || (self->c_functional.family == XC_FAMILY_GGA)) || ((self->x_functional.family == XC_FAMILY_HYB_GGA) || (self->c_functional.family == XC_FAMILY_HYB_GGA))) { a2_g = DOUBLEP(a2_array); dEda2_g = DOUBLEP(dEda2_array); } if ((self->x_functional.family == XC_FAMILY_MGGA) || (self->c_functional.family == XC_FAMILY_MGGA)) { a2_g = DOUBLEP(a2_array); tau_g = DOUBLEP(tau_array); dEda2_g = DOUBLEP(dEda2_array); dEdtau_g = DOUBLEP(dEdtau_array); } assert (self->xc_functional.family == XC_FAMILY_UNKNOWN); /* MDTMP not implemented */ /* find x functional */ switch(self->x_functional.family) { case XC_FAMILY_LDA: self->get_vxc_x = get_vxc_lda; break; case XC_FAMILY_GGA: self->get_vxc_x = get_vxc_gga; break; case XC_FAMILY_HYB_GGA: self->get_vxc_x = get_vxc_gga; break; case XC_FAMILY_MGGA: self->get_vxc_x = get_vxc_mgga; break; /* default: */ /* printf("lxcXCFunctional_CalculateSpinPaired: exchange functional '%d' not found\n", */ /* self->x_functional.family); */ } /* find c functional */ switch(self->c_functional.family) { case XC_FAMILY_LDA: self->get_vxc_c = get_vxc_lda; break; case XC_FAMILY_GGA: self->get_vxc_c = get_vxc_gga; break; case XC_FAMILY_HYB_GGA: self->get_vxc_c = get_vxc_gga; break; case XC_FAMILY_MGGA: self->get_vxc_c = get_vxc_mgga; break; /* default: */ /* printf("lxcXCFunctional_CalculateSpinPaired: correlation functional '%d' not found\n", */ /* self->c_functional.family); */ } /* ################################################################ */ for (int g = 0; g < ng; g++) { double n = n_g[g]; if (n < NMIN) n = NMIN; double a2 = 0.0; /* initialize for lda */ if (((self->x_functional.family == XC_FAMILY_GGA) || (self->c_functional.family == XC_FAMILY_GGA)) || ((self->x_functional.family == XC_FAMILY_HYB_GGA) || (self->c_functional.family == XC_FAMILY_HYB_GGA))) { a2 = a2_g[g]; } double tau = 0.0; if ((self->x_functional.family == XC_FAMILY_MGGA) || (self->c_functional.family == XC_FAMILY_MGGA)) { a2 = a2_g[g]; tau = tau_g[g]; } double dExdn = 0.0; double dExda2 = 0.0; double ex = 0.0; double dEcdn = 0.0; double dEcda2 = 0.0; double ec = 0.0; double dExdtau=0.0; double dEcdtau=0.0; double point[7]; /* generalized point */ // from http://www.tddft.org/programs/octopus/wiki/index.php/Libxc:manual // rhoa rhob sigmaaa sigmaab sigmabb taua taub // \sigma[0] = \nabla n_\uparrow \cdot \nabla n_\uparrow \qquad // \sigma[1] = \nabla n_\uparrow \cdot \nabla n_\downarrow \qquad // \sigma[2] = \nabla n_\downarrow \cdot \nabla n_\downarrow \qquad double derivative_x[7]; /* generalized potential */ double derivative_c[7]; /* generalized potential */ // vrhoa vrhob vsigmaaa vsigmaab vsigmabb dedtaua dedtaub // {\rm vrho}_{\alpha} = \frac{\partial E}{\partial n_\alpha} \qquad // {\rm vsigma}_{\alpha} = \frac{\partial E}{\partial \sigma_\alpha} // {\rm dedtau}_{\alpha} = \frac{\partial E}{\partial \tau_\alpha} for(int j=0; j<7; j++) { point[j] = derivative_x[j] = derivative_c[j] = 0.0; } point[0] = n; /* -> rho */ point[2] = a2; /* -> sigma */ point[5] = tau; /* -> tau */ /* calculate exchange */ if (self->x_functional.family != XC_FAMILY_UNKNOWN) { self->get_vxc_x(&(self->x_functional), point, &ex, derivative_x); dExdn = derivative_x[0]; dExda2 = derivative_x[2]; dExdtau = derivative_x[5]; if (self->c_functional.family == XC_FAMILY_HYB_GGA) { // MDTMP - a hack: HYB_GGA handle h1 internally in c_functional dExdn = 0.0; dExda2 = 0.0; dExdtau = 0.0; } } /* calculate correlation */ if (self->c_functional.family != XC_FAMILY_UNKNOWN) { self->get_vxc_c(&(self->c_functional), point, &ec, derivative_c); dEcdn = derivative_c[0]; dEcda2 = derivative_c[2]; dEcdtau = derivative_c[5]; } if (((self->x_functional.family == XC_FAMILY_GGA) || (self->c_functional.family == XC_FAMILY_GGA)) || ((self->x_functional.family == XC_FAMILY_HYB_GGA) || (self->c_functional.family == XC_FAMILY_HYB_GGA))) { dEda2_g[g] = dExda2 + dEcda2; } if ((self->x_functional.family == XC_FAMILY_MGGA) || (self->c_functional.family == XC_FAMILY_MGGA)) { dEda2_g[g] = dExda2 + dEcda2; dEdtau_g[g] = dExdtau + dEcdtau; } e_g[g] = n * (ex + ec); v_g[g] += dExdn + dEcdn; } Py_RETURN_NONE; }
static PyObject* lxcXCFunctional_CalculateFXCSpinPaired(lxcXCFunctionalObject *self, PyObject *args) { PyArrayObject* n_array; /* rho */ PyArrayObject* v2rho2_array; /* d2E/drho2 */ PyArrayObject* a2_array = 0; /* |nabla rho|^2*/ PyArrayObject* v2rhosigma_array = 0; /* d2E/drhod|nabla rho|^2 */ PyArrayObject* v2sigma2_array = 0; /* d2E/drhod|nabla rho|^2 */ if (!PyArg_ParseTuple(args, "OO|OOO", &n_array, &v2rho2_array, /* object | optional objects*/ &a2_array, &v2rhosigma_array, &v2sigma2_array)) return NULL; /* find nspin */ int nspin = self->nspin; assert(nspin == XC_UNPOLARIZED); /* we are spinpaired */ int ng = n_array->dimensions[0]; /* number of grid points */ const double* n_g = DOUBLEP(n_array); /* density on the grid */ double* v2rho2_g = DOUBLEP(v2rho2_array); /* v on the grid */ const double* a2_g = 0; /* a2 on the grid */ double* v2rhosigma_g = 0; /* d2Ednda2 on the grid */ double* v2sigma2_g = 0; /* d2Eda2da2 on the grid */ if (((self->x_functional.family == XC_FAMILY_GGA) || (self->c_functional.family == XC_FAMILY_GGA)) || ((self->x_functional.family == XC_FAMILY_HYB_GGA) || (self->c_functional.family == XC_FAMILY_HYB_GGA))) { a2_g = DOUBLEP(a2_array); v2rhosigma_g = DOUBLEP(v2rhosigma_array); v2sigma2_g = DOUBLEP(v2sigma2_array); } //assert (self->x_functional.family != XC_FAMILY_MGGA); /* MDTMP - not implemented yet */ //assert (self->c_functional.family != XC_FAMILY_MGGA); /* MDTMP - not implemented yet */ assert (self->xc_functional.family == XC_FAMILY_UNKNOWN); /* MDTMP not implemented */ /* find x functional */ switch(self->x_functional.family) { case XC_FAMILY_LDA: self->get_fxc_x = get_fxc_lda; break; case XC_FAMILY_GGA: self->get_fxc_x = get_fxc_gga; break; /* default: */ /* printf("lxcXCFunctional_CalculateSpinPaired: exchange functional '%d' not found\n", */ /* self->x_functional.family); */ } /* find c functional */ switch(self->c_functional.family) { case XC_FAMILY_LDA: self->get_fxc_c = get_fxc_lda; break; case XC_FAMILY_GGA: self->get_fxc_c = get_fxc_gga; break; /* default: */ /* printf("lxcXCFunctional_CalculateSpinPaired: correlation functional '%d' not found\n", */ /* self->c_functional.family); */ } /* ################################################################ */ for (int g = 0; g < ng; g++) { double n = n_g[g]; if (n < NMIN) n = NMIN; double a2 = 0.0; /* initialize for lda */ if (((self->x_functional.family == XC_FAMILY_GGA) || (self->c_functional.family == XC_FAMILY_GGA)) || ((self->x_functional.family == XC_FAMILY_HYB_GGA) || (self->c_functional.family == XC_FAMILY_HYB_GGA))) { a2 = a2_g[g]; } double point[7]; /* generalized point */ // from http://www.tddft.org/programs/octopus/wiki/index.php/Libxc:manual // rhoa rhob sigmaaa sigmaab sigmabb taua taub // \sigma[0] = \nabla n_\uparrow \cdot \nabla n_\uparrow \qquad // \sigma[1] = \nabla n_\uparrow \cdot \nabla n_\downarrow \qquad // \sigma[2] = \nabla n_\downarrow \cdot \nabla n_\downarrow \qquad double derivative_x[5][5]; /* generalized derivative */ double derivative_c[5][5]; /* generalized derivative */ double v2rho2_x[3]; double v2rhosigma_x[6]; double v2sigma2_x[6]; double v2rho2_c[3]; double v2rhosigma_c[6]; double v2sigma2_c[6]; // one that uses this: please add description of spin derivative order notation // (see c/libxc/src/gga_perdew.c) MDTMP for(int i=0; i<3; i++) v2rho2_x[i] = 0.0; for(int i=0; i<3; i++) v2rho2_c[i] = 0.0; for(int i=0; i<6; i++){ v2rhosigma_x[i] = 0.0; v2sigma2_x[i] = 0.0; v2rhosigma_c[i] = 0.0; v2sigma2_c[i] = 0.0; } for(int j=0; j<7; j++) { point[j] = 0.0; } for(int i=0; i<5; i++) { for(int j=0; j<5; j++) { derivative_x[i][j] = derivative_c[i][j] = 0.0; } } point[0] = n; /* -> rho */ point[2] = a2; /* -> sigma */ /* calculate exchange */ if (self->x_functional.family != XC_FAMILY_UNKNOWN) { self->get_fxc_x(&(self->x_functional), point, derivative_x); //printf("fxc_x '%f'\n", derivative_x[0][0]); // MDTMP v2rho2_x[0] = derivative_x[0][0]; v2rho2_x[1] = derivative_x[0][1]; // XC_POLARIZED v2rho2_x[2] = derivative_x[1][1]; // XC_POLARIZED //printf("fxc_x '%f'\n", derivative_x[0][2]); // MDTMP v2rhosigma_x[0] = derivative_x[0][2]; v2rhosigma_x[1] = derivative_x[0][3]; // XC_POLARIZED v2rhosigma_x[2] = derivative_x[0][4]; // XC_POLARIZED v2rhosigma_x[3] = derivative_x[1][2]; // XC_POLARIZED v2rhosigma_x[4] = derivative_x[1][3]; // XC_POLARIZED v2rhosigma_x[5] = derivative_x[1][4]; // XC_POLARIZED //printf("fxc_x '%f'\n", derivative_x[2][2]); // MDTMP v2sigma2_x[0] = derivative_x[2][2]; /* aa_aa */ v2sigma2_x[1] = derivative_x[2][3]; // XC_POLARIZED /* aa_ab */ v2sigma2_x[2] = derivative_x[2][4]; // XC_POLARIZED /* aa_bb */ v2sigma2_x[3] = derivative_x[3][3]; // XC_POLARIZED /* ab_ab */ v2sigma2_x[4] = derivative_x[3][4]; // XC_POLARIZED /* ab_bb */ v2sigma2_x[5] = derivative_x[4][4]; // XC_POLARIZED /* bb_bb */ } /* calculate correlation */ if (self->c_functional.family != XC_FAMILY_UNKNOWN) { self->get_fxc_c(&(self->c_functional), point, derivative_c); v2rho2_c[0] = derivative_c[0][0]; v2rho2_c[1] = derivative_c[0][1]; // XC_POLARIZED v2rho2_c[2] = derivative_c[1][1]; // XC_POLARIZED v2rhosigma_c[0] = derivative_c[0][2]; v2rhosigma_c[1] = derivative_c[0][3]; // XC_POLARIZED v2rhosigma_c[2] = derivative_c[0][4]; // XC_POLARIZED v2rhosigma_c[3] = derivative_c[1][2]; // XC_POLARIZED v2rhosigma_c[4] = derivative_c[1][3]; // XC_POLARIZED v2rhosigma_c[5] = derivative_c[1][4]; // XC_POLARIZED v2sigma2_c[0] = derivative_c[2][2]; /* aa_aa */ v2sigma2_c[1] = derivative_c[2][3]; // XC_POLARIZED /* aa_ab */ v2sigma2_c[2] = derivative_c[2][4]; // XC_POLARIZED /* aa_bb */ v2sigma2_c[3] = derivative_c[3][3]; // XC_POLARIZED /* ab_ab */ v2sigma2_c[4] = derivative_c[3][4]; // XC_POLARIZED /* ab_bb */ v2sigma2_c[5] = derivative_c[4][4]; // XC_POLARIZED /* bb_bb */ } if (((self->x_functional.family == XC_FAMILY_GGA) || (self->c_functional.family == XC_FAMILY_GGA)) || ((self->x_functional.family == XC_FAMILY_HYB_GGA) || (self->c_functional.family == XC_FAMILY_HYB_GGA))) { v2rhosigma_g[g] = v2rhosigma_x[0] + v2rhosigma_c[0]; v2sigma2_g[g] = v2sigma2_x[0] + v2sigma2_c[0]; } v2rho2_g[g] += v2rho2_x[0] + v2rho2_c[0]; } Py_RETURN_NONE; }
static PyObject* lxcXCFunctional_CalculateSpinPolarized(lxcXCFunctionalObject *self, PyObject *args) { PyArrayObject* e; PyArrayObject* na; PyArrayObject* va; PyArrayObject* nb; PyArrayObject* vb; PyArrayObject* a2 = 0; PyArrayObject* aa2 = 0; PyArrayObject* ab2 = 0; PyArrayObject* dEda2 = 0; PyArrayObject* dEdaa2 = 0; PyArrayObject* dEdab2 = 0; PyArrayObject* taua = 0; /* taua*/ PyArrayObject* taub = 0 ; /* taub*/ PyArrayObject* dEdtaua = 0; /* dE/dtaua */ PyArrayObject* dEdtaub = 0; /* dE/dtaub */ if (!PyArg_ParseTuple(args, "OOOOO|OOOOOOOOOOOOOOO", &e, &na, &va, &nb, &vb, &a2, &aa2, &ab2, &dEda2, &dEdaa2, &dEdab2, &taua, &taub, &dEdtaua, &dEdtaub)) return NULL; /* find nspin */ int nspin = self->nspin; assert(nspin == XC_POLARIZED); /* we are spinpolarized */ int ng = e->dimensions[0]; /* number of grid points */ double* e_g = DOUBLEP(e); /* e on the grid */ const double* na_g = DOUBLEP(na); /* alpha density on the grid */ double* va_g = DOUBLEP(va); /* alpha v on the grid */ const double* nb_g = DOUBLEP(nb); /* beta density on the grid */ double* vb_g = DOUBLEP(vb); /* beta v on the grid */ const double* a2_g = 0; /* sigmaab on the grid */ const double* aa2_g = 0; /* sigmaaa on the grid */ const double* ab2_g = 0; /* sigmabb on the grid */ double* taua_g = 0; /* taua on the grid */ double* taub_g = 0; /* taub on the grid */ double* dEda2_g = 0; /* dEdsigmaab on the grid */ double* dEdaa2_g = 0; /* dEdsigmaaa on the grid */ double* dEdab2_g = 0; /* dEdsigmabb on the grid */ double* dEdtaua_g = 0; /* dEdta on the grid */ double* dEdtaub_g = 0; /* dEdtb on the grid */ if (((self->x_functional.family == XC_FAMILY_GGA) || (self->c_functional.family == XC_FAMILY_GGA)) || ((self->x_functional.family == XC_FAMILY_HYB_GGA) || (self->c_functional.family == XC_FAMILY_HYB_GGA))) { a2_g = DOUBLEP(a2); aa2_g = DOUBLEP(aa2); ab2_g = DOUBLEP(ab2); dEda2_g = DOUBLEP(dEda2); dEdaa2_g = DOUBLEP(dEdaa2); dEdab2_g = DOUBLEP(dEdab2); } if ((self->x_functional.family == XC_FAMILY_MGGA) || (self->c_functional.family == XC_FAMILY_MGGA)) { a2_g = DOUBLEP(a2); aa2_g = DOUBLEP(aa2); ab2_g = DOUBLEP(ab2); taua_g = DOUBLEP(taua); taub_g = DOUBLEP(taub); dEda2_g = DOUBLEP(dEda2); dEdaa2_g = DOUBLEP(dEdaa2); dEdab2_g = DOUBLEP(dEdab2); dEdtaua_g = DOUBLEP(dEdtaua); dEdtaub_g = DOUBLEP(dEdtaub); } assert (self->xc_functional.family == XC_FAMILY_UNKNOWN); /* MDTMP not implemented */ /* find x functional */ switch(self->x_functional.family) { case XC_FAMILY_LDA: self->get_vxc_x = get_vxc_lda; break; case XC_FAMILY_GGA: self->get_vxc_x = get_vxc_gga; break; case XC_FAMILY_HYB_GGA: self->get_vxc_x = get_vxc_gga; break; case XC_FAMILY_MGGA: self->get_vxc_x = get_vxc_mgga; break; /* default: */ /* printf("lxcXCFunctional_CalculateSpinPolarized: exchange functional '%d' not found\n", */ /* self->x_functional.family); */ } /* find c functional */ switch(self->c_functional.family) { case XC_FAMILY_LDA: self->get_vxc_c = get_vxc_lda; break; case XC_FAMILY_GGA: self->get_vxc_c = get_vxc_gga; break; case XC_FAMILY_HYB_GGA: self->get_vxc_c = get_vxc_gga; break; case XC_FAMILY_MGGA: self->get_vxc_c = get_vxc_mgga; break; /* default: */ /* printf("lxcXCFunctional_CalculateSpinPolarized: correlation functional '%d' not found\n", */ /* self->c_functional.family); */ } /* ################################################################ */ for (int g = 0; g < ng; g++) { double na = na_g[g]; if (na < NMIN) na = NMIN; double sigma0 = 0.0; /* initialize for lda */ double sigma1 = 0.0; /* initialize for lda */ double sigma2 = 0.0; /* initialize for lda */ double taua = 0.0, taub =0.0; double dExdtaua = 0.0; double dExdtaub = 0.0; double dEcdtaua = 0.0; double dEcdtaub = 0.0; if (((self->x_functional.family == XC_FAMILY_GGA) || (self->c_functional.family == XC_FAMILY_GGA)) || ((self->x_functional.family == XC_FAMILY_HYB_GGA) || (self->c_functional.family == XC_FAMILY_HYB_GGA))) { sigma0 = a2_g[g]; sigma1 = aa2_g[g]; sigma2 = ab2_g[g]; } double nb = nb_g[g]; if (nb < NMIN) nb = NMIN; if ((self->x_functional.family == XC_FAMILY_MGGA) || (self->c_functional.family == XC_FAMILY_MGGA)) { sigma0 = a2_g[g]; sigma1 = aa2_g[g]; sigma2 = ab2_g[g]; taua = taua_g[g]; taub = taub_g[g]; } double n = na + nb; double dExdna = 0.0; double dExdsigma0 = 0.0; double dExdnb = 0.0; double dExdsigma2 = 0.0; double ex = 0.0; double dExdsigma1 = 0.0; double dEcdna = 0.0; double dEcdsigma0 = 0.0; double dEcdnb = 0.0; double dEcdsigma2 = 0.0; double ec = 0.0; double dEcdsigma1 = 0.0; double point[7]; /* generalized point */ // from http://www.tddft.org/programs/octopus/wiki/index.php/Libxc:manual // rhoa rhob sigmaaa sigmaab sigmabb taua taub // \sigma[0] = \nabla n_\uparrow \cdot \nabla n_\uparrow \qquad // \sigma[1] = \nabla n_\uparrow \cdot \nabla n_\downarrow \qquad // \sigma[2] = \nabla n_\downarrow \cdot \nabla n_\downarrow \qquad double derivative_x[7]; /* generalized potential */ double derivative_c[7]; /* generalized potential */ // vrhoa vrhob vsigmaaa vsigmaab vsigmabb dedtaua dedtaub // {\rm vrho}_{\alpha} = \frac{\partial E}{\partial n_\alpha} \qquad // {\rm vsigma}_{\alpha} = \frac{\partial E}{\partial \sigma_\alpha} // {\rm dedtau}_{\alpha} = \frac{\partial E}{\partial \tau_\alpha} for(int j=0; j<7; j++) { point[j] = derivative_x[j] = derivative_c[j]= 0.0; } point[0] = na; point[1] = nb; point[2] = sigma0; point[3] = sigma1; point[4] = sigma2; point[5] = taua; point[6] = taub; /* calculate exchange */ if (self->x_functional.family != XC_FAMILY_UNKNOWN) { self->get_vxc_x(&(self->x_functional), point, &ex, derivative_x); dExdna = derivative_x[0]; dExdnb = derivative_x[1]; dExdsigma0 = derivative_x[2]; dExdsigma1 = derivative_x[3]; dExdsigma2 = derivative_x[4]; dExdtaua = derivative_x[5]; dExdtaub = derivative_x[6]; if (self->c_functional.family == XC_FAMILY_HYB_GGA) { // MDTMP - a hack: HYB_GGA handle h1 internally in c_functional dExdna = 0.0; dExdnb = 0.0; dExdsigma0 = 0.0; dExdsigma1 = 0.0; dExdsigma2 = 0.0; dExdtaua = 0.0; dExdtaub = 0.0; } } /* calculate correlation */ if (self->c_functional.family != XC_FAMILY_UNKNOWN) { self->get_vxc_c(&(self->c_functional), point, &ec, derivative_c); dEcdna = derivative_c[0]; dEcdnb = derivative_c[1]; dEcdsigma0 = derivative_c[2]; dEcdsigma1 = derivative_c[3]; dEcdsigma2 = derivative_c[4]; dEcdtaua = derivative_c[5]; dEcdtaub = derivative_c[6]; } if (((self->x_functional.family == XC_FAMILY_GGA) || (self->c_functional.family == XC_FAMILY_GGA)) || ((self->x_functional.family == XC_FAMILY_HYB_GGA) || (self->c_functional.family == XC_FAMILY_HYB_GGA))) { dEdaa2_g[g] = dExdsigma1 + dEcdsigma1; dEdab2_g[g] = dExdsigma2 + dEcdsigma2; dEda2_g[g] = dExdsigma0 + dEcdsigma0; } if ((self->x_functional.family == XC_FAMILY_MGGA) || (self->c_functional.family == XC_FAMILY_MGGA)) { dEdaa2_g[g] = dExdsigma1 + dEcdsigma1; dEdab2_g[g] = dExdsigma2 + dEcdsigma2; dEda2_g[g] = dExdsigma0 + dEcdsigma0; dEdtaua_g[g] = dExdtaua + dEcdtaua; dEdtaub_g[g] = dExdtaub + dEcdtaub; } e_g[g] = n* (ex + ec); va_g[g] += dExdna + dEcdna; vb_g[g] += dExdnb + dEcdnb; } Py_RETURN_NONE; }