Tinker9 70bd052 (Thu Nov 9 12:11:35 2023 -0800)
|
Preconditioned Conjugate Gradient. More...
Functions | |
__global__ void | tinker::pcgUdirV1 (int n, const real *polarity, real(*udir)[3], const real(*field)[3]) |
Calculates the direct induced dipoles 2b . More... | |
__global__ void | tinker::pcgUdirV2 (int n, const real *polarity, real(*udir)[3], real(*udirp)[3], const real(*field)[3], const real(*fieldp)[3]) |
Calculates the direct induced dipoles 2b . More... | |
__global__ void | tinker::pcgRsd0V1 (int n, const real *polarity_inv, real(*rsd)[3], const real(*udir)[3], const real(*uind)[3], const real(*field)[3]) |
Calculates the initial residuals for the predictor 3a . More... | |
__global__ void | tinker::pcgRsd0V2 (int n, const real *polarity_inv, real(*rsd)[3], real(*rsp)[3], const real(*udir)[3], const real(*udip)[3], const real(*uind)[3], const real(*uinp)[3], const real(*field)[3], const real(*fielp)[3]) |
Calculates the initial residuals for the predictor 3a . More... | |
__global__ void | tinker::pcgRsd0V3 (int n, const real *polarity_inv, real(*rsd)[3], const real(*udir)[3], const real(*uind)[3], const real(*field)[3], const real(*polscale)[3][3]) |
Calculates the initial residuals for the predictor 3a using the Exchange-Polarization Model. More... | |
__global__ void | tinker::pcgRsd0 (int n, const real *polarity, real(*rsd)[3], real(*rsdp)[3]) |
Sets the residual to zero for atoms with zero polarizability. More... | |
__global__ void | tinker::pcgRsd1 (int n, const real *polarity, real(*rsd)[3]) |
Sets the residual to zero for atoms with zero polarizability. More... | |
__global__ void | tinker::pcgP1 (int n, const real *polarity_inv, real(*vec)[3], real(*vecp)[3], const real(*conj)[3], const real(*conjp)[3], const real(*field)[3], const real(*fieldp)[3]) |
Calculates the vector Ap – step 5a . More... | |
__global__ void | tinker::pcgP2 (int n, const real *polarity, const real *ka, const real *kap, const real *ksum, const real *ksump, real(*uind)[3], real(*uinp)[3], const real(*conj)[3], const real(*conjp)[3], real(*rsd)[3], real(*rsdp)[3], const real(*vec)[3], const real(*vecp)[3]) |
Updates the induced dipoles and residuals – step 5bc . More... | |
__global__ void | tinker::pcgP3 (int n, const real *ksum, const real *ksump, const real *ksum1, const real *ksump1, real(*conj)[3], real(*conjp)[3], const real(*zrsd)[3], const real(*zrsdp)[3]) |
Calculates the scalar b from s0 and s, then updates the vector p – 5ef . More... | |
__global__ void | tinker::pcgP4 (int n, const real *polarity_inv, real(*vec)[3], const real(*conj)[3], const real(*field)[3]) |
Calculates the vector Ap – step 5a . More... | |
__global__ void | tinker::pcgP5 (int n, const real *polarity, const real *ka, const real *ksum, real(*uind)[3], const real(*conj)[3], real(*rsd)[3], const real(*vec)[3]) |
Updates the induced dipoles and residuals – step 5bc . More... | |
__global__ void | tinker::pcgP6 (int n, const real *ksum, const real *ksum1, real(*conj)[3], const real(*zrsd)[3]) |
Calculates the scalar b from s0 and s, then updates the vector p – 5ef . More... | |
__global__ void | tinker::pcgPeek (int n, float pcgpeek, const real *polarity, real(*uind)[3], real(*uinp)[3], const real(*rsd)[3], const real(*rsdp)[3]) |
Applies the peek step. More... | |
__global__ void | tinker::pcgPeek1 (int n, float pcgpeek, const real *polarity, real(*uind)[3], const real(*rsd)[3]) |
Applies the peek step. More... | |
Preconditioned Conjugate Gradient.
ufield(u)
.fromPrediction
ufield(u)
ufield(aE)
ufield(p)
__global__ void tinker::pcgP1 | ( | int | n, |
const real * | polarity_inv, | ||
real(*) | vec[3], | ||
real(*) | vecp[3], | ||
const real(*) | conj[3], | ||
const real(*) | conjp[3], | ||
const real(*) | field[3], | ||
const real(*) | fieldp[3] | ||
) |
Calculates the vector Ap – step 5a
.
n | Number of atoms. | |
polarity_inv | Multiplicative inverse of polarizabilities. | |
[out] | vec | d-vector Ap. |
[out] | vecp | p-vector Ap. |
conj | d-vector p. | |
conjp | p-vector p. | |
field | d-ufield(p) . | |
fieldp | p-ufield(p) . |
__global__ void tinker::pcgP2 | ( | int | n, |
const real * | polarity, | ||
const real * | ka, | ||
const real * | kap, | ||
const real * | ksum, | ||
const real * | ksump, | ||
real(*) | uind[3], | ||
real(*) | uinp[3], | ||
const real(*) | conj[3], | ||
const real(*) | conjp[3], | ||
real(*) | rsd[3], | ||
real(*) | rsdp[3], | ||
const real(*) | vec[3], | ||
const real(*) | vecp[3] | ||
) |
Updates the induced dipoles and residuals – step 5bc
.
n | Number of atoms. | |
polarity | Polarizabilities. | |
ka | Device pointer to the dot product d-pAp. | |
kap | Device pointer to the dot product p-pAp. | |
ksum | Device pointer to the dot product d-s0. | |
ksump | Device pointer to the dot product p-s0. | |
[in,out] | uind | d-dipoles. |
[in,out] | uinp | p-dipoles. |
conj | d-vector p. | |
conjp | p-vector p. | |
[in,out] | rsd | d-residuals. |
[in,out] | rsdp | p-residuals. |
vec | d-vector Ap. | |
vecp | p-vector Ap. |
__global__ void tinker::pcgP3 | ( | int | n, |
const real * | ksum, | ||
const real * | ksump, | ||
const real * | ksum1, | ||
const real * | ksump1, | ||
real(*) | conj[3], | ||
real(*) | conjp[3], | ||
const real(*) | zrsd[3], | ||
const real(*) | zrsdp[3] | ||
) |
Calculates the scalar b from s0 and s, then updates the vector p – 5ef
.
n | Number of atoms. | |
ksum | Device pointer to the d-scalar s0. | |
ksump | Device pointer to the p-scalar s0. | |
ksum1 | Device pointer to the d-scalar s. | |
ksump1 | Device pointer to the p-scalar s. | |
[in,out] | conj | d-vector p. |
[in,out] | conjp | p-vector p. |
zrsd | d-vector Mr. | |
zrsdp | p-vector Mr. |
__global__ void tinker::pcgP4 | ( | int | n, |
const real * | polarity_inv, | ||
real(*) | vec[3], | ||
const real(*) | conj[3], | ||
const real(*) | field[3] | ||
) |
Calculates the vector Ap – step 5a
.
n | Number of atoms. | |
polarity_inv | Multiplicative inverse of polarizabilities. | |
[out] | vec | Vector Ap. |
conj | Vector p. | |
field | ufield(p) . |
__global__ void tinker::pcgP5 | ( | int | n, |
const real * | polarity, | ||
const real * | ka, | ||
const real * | ksum, | ||
real(*) | uind[3], | ||
const real(*) | conj[3], | ||
real(*) | rsd[3], | ||
const real(*) | vec[3] | ||
) |
Updates the induced dipoles and residuals – step 5bc
.
Vector Ap.
n | Number of atoms. | |
polarity | Polarizabilities. | |
ka | Device pointer to the dot product pAp. | |
ksum | Device pointer to the dot product s0. | |
[in,out] | uind | Induced dipoles. |
conj | Vector p. | |
[in,out] | rsd | Residuals. |
__global__ void tinker::pcgP6 | ( | int | n, |
const real * | ksum, | ||
const real * | ksum1, | ||
real(*) | conj[3], | ||
const real(*) | zrsd[3] | ||
) |
Calculates the scalar b from s0 and s, then updates the vector p – 5ef
.
n | Number of atoms. | |
ksum | Deivce pointer to the scalar s0. | |
ksum1 | Device pointer to the scalar s. | |
[in,out] | conj | Vector p. |
zrsd | Vector Mr. |
__global__ void tinker::pcgPeek | ( | int | n, |
float | pcgpeek, | ||
const real * | polarity, | ||
real(*) | uind[3], | ||
real(*) | uinp[3], | ||
const real(*) | rsd[3], | ||
const real(*) | rsdp[3] | ||
) |
Applies the peek step.
n | Number of atoms. | |
pcgpeek | The "peek" parameter. | |
polarity | Polarizabilities. | |
[in,out] | uind | d-dipoles. |
[in,out] | uinp | p-dipoles. |
rsd | Final d-residual. | |
rsdp | Final p-residual. |
__global__ void tinker::pcgPeek1 | ( | int | n, |
float | pcgpeek, | ||
const real * | polarity, | ||
real(*) | uind[3], | ||
const real(*) | rsd[3] | ||
) |
Applies the peek step.
n | Number of atoms. | |
pcgpeek | The "peek" parameter. | |
polarity | Polarizabilities. | |
[in,out] | uind | Dipoles. |
rsd | Final residual. |
Sets the residual to zero for atoms with zero polarizability.
n | Number of atoms. | |
polarity | polarizabilities. | |
[out] | rsd | Zero d-residuals. |
[out] | rsdp | Zero p-residuals. |
__global__ void tinker::pcgRsd0V1 | ( | int | n, |
const real * | polarity_inv, | ||
real(*) | rsd[3], | ||
const real(*) | udir[3], | ||
const real(*) | uind[3], | ||
const real(*) | field[3] | ||
) |
Calculates the initial residuals for the predictor 3a
.
n | Number of atoms. | |
polarity_inv | Multiplicative inverse of polarizabilities. | |
[out] | rsd | Initial residuals. |
udir | Direct induced dipoles 2b . | |
uind | Initial induced dipoles 2a . | |
field | Mutual field. |
__global__ void tinker::pcgRsd0V2 | ( | int | n, |
const real * | polarity_inv, | ||
real(*) | rsd[3], | ||
real(*) | rsp[3], | ||
const real(*) | udir[3], | ||
const real(*) | udip[3], | ||
const real(*) | uind[3], | ||
const real(*) | uinp[3], | ||
const real(*) | field[3], | ||
const real(*) | fielp[3] | ||
) |
Calculates the initial residuals for the predictor 3a
.
n | Number of atoms. | |
polarity_inv | Multiplicative inverse of polarizabilities. | |
[out] | rsd | Initial d-residuals. |
[out] | rsp | Initial p-residuals. |
udir | Direct induced d-dipoles 2b . | |
udip | Direct induced p-dipoles 2b . | |
uind | Initial induced d-dipoles 2a . | |
uinp | Initial induced p-dipoles 2a . | |
field | Mutual d-field. | |
fielp | Mutual p-field. |
__global__ void tinker::pcgRsd0V3 | ( | int | n, |
const real * | polarity_inv, | ||
real(*) | rsd[3], | ||
const real(*) | udir[3], | ||
const real(*) | uind[3], | ||
const real(*) | field[3], | ||
const real(*) | polscale[3][3] | ||
) |
Calculates the initial residuals for the predictor 3a
using the Exchange-Polarization Model.
n | Number of atoms. | |
polarity_inv | Multiplicative inverse of polarizabilities. | |
[out] | rsd | Initial residuals. |
udir | Direct induced dipoles 2b . | |
uind | Initial induced dipoles 2a . | |
field | Mutual field. | |
polscale | The \( I+kS^2 \) matrices for every atom (Eq. (7d)). |
Sets the residual to zero for atoms with zero polarizability.
n | Number of atoms. | |
polarity | polarizabilities. | |
[out] | rsd | Zero residuals. |
__global__ void tinker::pcgUdirV1 | ( | int | n, |
const real * | polarity, | ||
real(*) | udir[3], | ||
const real(*) | field[3] | ||
) |
Calculates the direct induced dipoles 2b
.
n | Number of atoms. | |
polarity | Polarizabilities. | |
[out] | udir | Direct induced dipoles 2b . |
field | External field. |
__global__ void tinker::pcgUdirV2 | ( | int | n, |
const real * | polarity, | ||
real(*) | udir[3], | ||
real(*) | udirp[3], | ||
const real(*) | field[3], | ||
const real(*) | fieldp[3] | ||
) |
Calculates the direct induced dipoles 2b
.
n | Number of atoms. | |
polarity | Polarizabilities. | |
[out] | udir | Direct induced d-dipoles 2b . |
[out] | udirp | Direct induced p-dipoles 2b . |
field | External d-field. | |
fieldp | External p-field. |