Tinker9 70bd052 (Thu Nov 9 12:11:35 2023 -0800)
Loading...
Searching...
No Matches
cuinduce.h
1#pragma once
2#include "ff/precision.h"
3
4namespace tinker {
7
36
38__global__
39void pcgUdirV1(int n,
40 const real* polarity,
41 real (*udir)[3],
42 const real (*field)[3]
43);
44
46__global__
47void pcgUdirV2(int n,
48 const real* polarity,
49 real (*udir)[3],
50 real (*udirp)[3],
51 const real (*field)[3],
52 const real (*fieldp)[3]
53);
54
56__global__
57void pcgRsd0V1(int n,
58 const real* polarity_inv,
59 real (*rsd)[3],
60 const real (*udir)[3],
61 const real (*uind)[3],
62 const real (*field)[3]
63);
64
66__global__
67void pcgRsd0V2(int n,
68 const real* polarity_inv,
69 real (*rsd)[3],
70 real (*rsp)[3],
71 const real (*udir)[3],
72 const real (*udip)[3],
73 const real (*uind)[3],
74 const real (*uinp)[3],
75 const real (*field)[3],
76 const real (*fielp)[3]
77);
78
81__global__
82void pcgRsd0V3(int n,
83 const real* polarity_inv,
84 real (*rsd)[3],
85 const real (*udir)[3],
86 const real (*uind)[3],
87 const real (*field)[3],
88 const real (*polscale)[3][3]
89);
90
92__global__
93void pcgRsd0(int n,
94 const real* polarity,
95 real (*rsd)[3],
96 real (*rsdp)[3]
97);
98
100__global__
101void pcgRsd1(int n,
102 const real* polarity,
103 real (*rsd)[3]
104);
105
107__global__
108void pcgP1(int n,
109 const real* polarity_inv,
110 real (*vec)[3],
111 real (*vecp)[3],
112 const real (*conj)[3],
113 const real (*conjp)[3],
114 const real (*field)[3],
115 const real (*fieldp)[3]
116);
117
119__global__
120void pcgP2(int n,
121 const real* polarity,
122 const real* ka,
123 const real* kap,
124 const real* ksum,
125 const real* ksump,
126 real (*uind)[3],
127 real (*uinp)[3],
128 const real (*conj)[3],
129 const real (*conjp)[3],
130 real (*rsd)[3],
131 real (*rsdp)[3],
132 const real (*vec)[3],
133 const real (*vecp)[3]
134);
135
137__global__
138void pcgP3(int n,
139 const real* ksum,
140 const real* ksump,
141 const real* ksum1,
142 const real* ksump1,
143 real (*conj)[3],
144 real (*conjp)[3],
145 const real (*zrsd)[3],
146 const real (*zrsdp)[3]
147);
148
150__global__
151void pcgP4(int n,
152 const real* polarity_inv,
153 real (*vec)[3],
154 const real (*conj)[3],
155 const real (*field)[3]
156);
157
159__global__
160void pcgP5(int n,
161 const real* polarity,
162 const real* ka,
163 const real* ksum,
164 real (*uind)[3],
165 const real (*conj)[3],
166 real (*rsd)[3],
167 const real (*vec)[3]);
168
170__global__
171void pcgP6(int n,
172 const real* ksum,
173 const real* ksum1,
174 real (*conj)[3],
175 const real (*zrsd)[3]
176);
177
179__global__
180void pcgPeek(int n,
181 float pcgpeek,
182 const real* polarity,
183 real (*uind)[3],
184 real (*uinp)[3],
185 const real (*rsd)[3],
186 const real (*rsdp)[3]
187);
188
190__global__
191void pcgPeek1(int n,
192 float pcgpeek,
193 const real* polarity,
194 real (*uind)[3],
195 const real (*rsd)[3]
196);
197
199}
int n
Number of atoms padded by WARP_SIZE.
__global__ void 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.
__global__ void pcgPeek1(int n, float pcgpeek, const real *polarity, real(*uind)[3], const real(*rsd)[3])
Applies the peek step.
__global__ void 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.
__global__ void 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.
__global__ void pcgP4(int n, const real *polarity_inv, real(*vec)[3], const real(*conj)[3], const real(*field)[3])
Calculates the vector Ap – step 5a.
__global__ void pcgRsd1(int n, const real *polarity, real(*rsd)[3])
Sets the residual to zero for atoms with zero polarizability.
__global__ void 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.
__global__ void 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.
__global__ void 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.
__global__ void pcgUdirV1(int n, const real *polarity, real(*udir)[3], const real(*field)[3])
Calculates the direct induced dipoles 2b.
__global__ void 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.
__global__ void 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.
__global__ void 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.
__global__ void 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.
__global__ void pcgRsd0(int n, const real *polarity, real(*rsd)[3], real(*rsdp)[3])
Sets the residual to zero for atoms with zero polarizability.
float real
Definition: precision.h:80
Definition: testrt.h:9
real * polarity_inv
real(* udir)[3]
real(* udirp)[3]
real(* uind)[3]
real * polarity
real(* polscale)[3][3]
real(* uinp)[3]