Tinker9 70bd052 (Thu Nov 9 12:11:35 2023 -0800)
Loading...
Searching...
No Matches
Functions
PCG

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...
 

Detailed Description

Preconditioned Conjugate Gradient.

The PCG method in Tinker9 is implemented as follows

  1. Goal: To solve Au = E, where A = 1/a + T.
    • (a) E: external field;
    • (b) a: polarizability;
    • (c) -Tu: mutual induced field; i.e., ufield(u).
  2. Guess initial induced dipole: u.
    • (a) Predict: u = fromPrediction
    • (b) Direct guess: u = aE
    • (c) Zero: u = 0
  3. Compute initial residual: r = E - Au
    • (a) Predict: r = (aE-u)/a + ufield(u)
    • (b) Direct guess: r = ufield(aE)
    • (c) Zero: r = E
  4. p = Mr; s0 = rMr
    • M is the preconditioner which approximates A.
  5. While not converged (r is still big):
    • (a) g = s0 / pAp, where Ap = p/a - ufield(p)
    • (b) u = u + gp
    • (c) r = r - gAp
    • (d) s = rMr
    • (e) b = s / s0
    • (f) p = Mr + bp; s0 = s
  6. Peek: u = u + peek*a*r – optional, empirical.

Function Documentation

◆ pcgP1()

__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.

Parameters
nNumber of atoms.
polarity_invMultiplicative inverse of polarizabilities.
[out]vecd-vector Ap.
[out]vecpp-vector Ap.
conjd-vector p.
conjpp-vector p.
fieldd-ufield(p).
fieldpp-ufield(p).

◆ pcgP2()

__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.

Parameters
nNumber of atoms.
polarityPolarizabilities.
kaDevice pointer to the dot product d-pAp.
kapDevice pointer to the dot product p-pAp.
ksumDevice pointer to the dot product d-s0.
ksumpDevice pointer to the dot product p-s0.
[in,out]uindd-dipoles.
[in,out]uinpp-dipoles.
conjd-vector p.
conjpp-vector p.
[in,out]rsdd-residuals.
[in,out]rsdpp-residuals.
vecd-vector Ap.
vecpp-vector Ap.

◆ pcgP3()

__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.

Parameters
nNumber of atoms.
ksumDevice pointer to the d-scalar s0.
ksumpDevice pointer to the p-scalar s0.
ksum1Device pointer to the d-scalar s.
ksump1Device pointer to the p-scalar s.
[in,out]conjd-vector p.
[in,out]conjpp-vector p.
zrsdd-vector Mr.
zrsdpp-vector Mr.

◆ pcgP4()

__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.

Parameters
nNumber of atoms.
polarity_invMultiplicative inverse of polarizabilities.
[out]vecVector Ap.
conjVector p.
fieldufield(p).

◆ pcgP5()

__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.

Parameters
nNumber of atoms.
polarityPolarizabilities.
kaDevice pointer to the dot product pAp.
ksumDevice pointer to the dot product s0.
[in,out]uindInduced dipoles.
conjVector p.
[in,out]rsdResiduals.

◆ pcgP6()

__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.

Parameters
nNumber of atoms.
ksumDeivce pointer to the scalar s0.
ksum1Device pointer to the scalar s.
[in,out]conjVector p.
zrsdVector Mr.

◆ pcgPeek()

__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.

Parameters
nNumber of atoms.
pcgpeekThe "peek" parameter.
polarityPolarizabilities.
[in,out]uindd-dipoles.
[in,out]uinpp-dipoles.
rsdFinal d-residual.
rsdpFinal p-residual.

◆ pcgPeek1()

__global__ void tinker::pcgPeek1 ( int  n,
float  pcgpeek,
const real polarity,
real(*)  uind[3],
const real(*)  rsd[3] 
)

Applies the peek step.

Parameters
nNumber of atoms.
pcgpeekThe "peek" parameter.
polarityPolarizabilities.
[in,out]uindDipoles.
rsdFinal residual.

◆ pcgRsd0()

__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.

Parameters
nNumber of atoms.
polaritypolarizabilities.
[out]rsdZero d-residuals.
[out]rsdpZero p-residuals.

◆ pcgRsd0V1()

__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.

Parameters
nNumber of atoms.
polarity_invMultiplicative inverse of polarizabilities.
[out]rsdInitial residuals.
udirDirect induced dipoles 2b.
uindInitial induced dipoles 2a.
fieldMutual field.

◆ pcgRsd0V2()

__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.

Parameters
nNumber of atoms.
polarity_invMultiplicative inverse of polarizabilities.
[out]rsdInitial d-residuals.
[out]rspInitial p-residuals.
udirDirect induced d-dipoles 2b.
udipDirect induced p-dipoles 2b.
uindInitial induced d-dipoles 2a.
uinpInitial induced p-dipoles 2a.
fieldMutual d-field.
fielpMutual p-field.

◆ pcgRsd0V3()

__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.

Parameters
nNumber of atoms.
polarity_invMultiplicative inverse of polarizabilities.
[out]rsdInitial residuals.
udirDirect induced dipoles 2b.
uindInitial induced dipoles 2a.
fieldMutual field.
polscaleThe \( I+kS^2 \) matrices for every atom (Eq. (7d)).

◆ pcgRsd1()

__global__ void tinker::pcgRsd1 ( int  n,
const real polarity,
real(*)  rsd[3] 
)

Sets the residual to zero for atoms with zero polarizability.

Parameters
nNumber of atoms.
polaritypolarizabilities.
[out]rsdZero residuals.

◆ pcgUdirV1()

__global__ void tinker::pcgUdirV1 ( int  n,
const real polarity,
real(*)  udir[3],
const real(*)  field[3] 
)

Calculates the direct induced dipoles 2b.

Parameters
nNumber of atoms.
polarityPolarizabilities.
[out]udirDirect induced dipoles 2b.
fieldExternal field.

◆ pcgUdirV2()

__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.

Parameters
nNumber of atoms.
polarityPolarizabilities.
[out]udirDirect induced d-dipoles 2b.
[out]udirpDirect induced p-dipoles 2b.
fieldExternal d-field.
fieldpExternal p-field.