1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38 package ffx.potential.openmm;
39
40 import com.sun.jna.ptr.DoubleByReference;
41 import com.sun.jna.ptr.IntByReference;
42 import edu.uiowa.jopenmm.OpenMMLibrary;
43 import ffx.crystal.Crystal;
44 import ffx.openmm.BondArray;
45 import ffx.openmm.Force;
46 import ffx.openmm.NonbondedForce;
47 import ffx.potential.bonded.Atom;
48 import ffx.potential.bonded.Bond;
49 import ffx.potential.nonbonded.NonbondedCutoff;
50 import ffx.potential.nonbonded.ParticleMeshEwald;
51 import ffx.potential.nonbonded.VanDerWaals;
52 import ffx.potential.nonbonded.VanDerWaalsForm;
53 import ffx.potential.parameters.ForceField;
54 import ffx.potential.parameters.MultipoleType;
55 import ffx.potential.parameters.VDWType;
56
57 import java.util.logging.Level;
58 import java.util.logging.Logger;
59
60 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AngstromsPerNm;
61 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_KJPerKcal;
62 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_NmPerAngstrom;
63 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_Boolean.OpenMM_False;
64 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_Boolean.OpenMM_True;
65 import static ffx.potential.parameters.VDWType.EPSILON_RULE.GEOMETRIC;
66 import static ffx.potential.parameters.VDWType.RADIUS_RULE.ARITHMETIC;
67 import static ffx.potential.parameters.VDWType.RADIUS_SIZE.RADIUS;
68 import static ffx.potential.parameters.VDWType.RADIUS_TYPE.R_MIN;
69 import static ffx.potential.parameters.VDWType.VDW_TYPE.LENNARD_JONES;
70 import static java.lang.String.format;
71 import static org.apache.commons.math3.util.FastMath.abs;
72
73
74
75
76
77
78 public class FixedChargeNonbondedForce extends NonbondedForce {
79
80 private static final Logger logger = Logger.getLogger(FixedChargeNonbondedForce.class.getName());
81
82
83
84
85 private boolean[] chargeExclusion;
86
87
88
89 private boolean[] vdWExclusion;
90
91
92
93 private double[] exceptionChargeProd;
94
95
96
97 private double[] exceptionEps;
98
99 public FixedChargeNonbondedForce(OpenMMEnergy openMMEnergy) {
100 VanDerWaals vdW = openMMEnergy.getVdwNode();
101 if (vdW == null) {
102
103 destroy();
104 return;
105 }
106
107
108
109
110
111 VanDerWaalsForm vdwForm = vdW.getVDWForm();
112 if (vdwForm.vdwType != LENNARD_JONES || vdwForm.radiusRule != ARITHMETIC
113 || vdwForm.epsilonRule != GEOMETRIC) {
114 logger.info(format(" VDW Type: %s", vdwForm.vdwType));
115 logger.info(format(" VDW Radius Rule: %s", vdwForm.radiusRule));
116 logger.info(format(" VDW Epsilon Rule: %s", vdwForm.epsilonRule));
117 logger.log(Level.SEVERE, " Unsupported van der Waals functional form.");
118 return;
119 }
120
121
122 double radScale = 1.0;
123 if (vdwForm.radiusSize == RADIUS) {
124 radScale = 2.0;
125 }
126
127
128 if (vdwForm.radiusType == R_MIN) {
129 radScale /= 1.122462048309372981;
130 }
131
132
133 Atom[] atoms = openMMEnergy.getMolecularAssembly().getAtomArray();
134 for (Atom atom : atoms) {
135 VDWType vdwType = atom.getVDWType();
136 double sigma = OpenMM_NmPerAngstrom * vdwType.radius * radScale;
137 double eps = OpenMM_KJPerKcal * vdwType.wellDepth;
138 double charge = 0.0;
139 MultipoleType multipoleType = atom.getMultipoleType();
140 if (multipoleType != null && atom.getElectrostatics()) {
141 charge = multipoleType.charge;
142 }
143 addParticle(charge, sigma, eps);
144 }
145
146
147 double lj14Scale = vdwForm.getScale14();
148 double coulomb14Scale = 1.0 / 1.2;
149 ParticleMeshEwald pme = openMMEnergy.getPmeNode();
150 if (pme != null) {
151 coulomb14Scale = pme.getScale14();
152 }
153 Bond[] bonds = openMMEnergy.getBonds();
154 if (bonds != null && bonds.length > 0) {
155 BondArray bondArray = new BondArray(0);
156 for (Bond bond : bonds) {
157 int i1 = bond.getAtom(0).getXyzIndex() - 1;
158 int i2 = bond.getAtom(1).getXyzIndex() - 1;
159 bondArray.append(i1, i2);
160 }
161 createExceptionsFromBonds(bondArray, coulomb14Scale, lj14Scale);
162 bondArray.destroy();
163
164 int num = getNumExceptions();
165 chargeExclusion = new boolean[num];
166 vdWExclusion = new boolean[num];
167 exceptionChargeProd = new double[num];
168 exceptionEps = new double[num];
169 IntByReference particle1 = new IntByReference();
170 IntByReference particle2 = new IntByReference();
171 DoubleByReference chargeProd = new DoubleByReference();
172 DoubleByReference sigma = new DoubleByReference();
173 DoubleByReference eps = new DoubleByReference();
174 for (int i = 0; i < num; i++) {
175 getExceptionParameters(i, particle1, particle2, chargeProd, sigma, eps);
176 if (abs(chargeProd.getValue()) > 0.0) {
177 chargeExclusion[i] = false;
178 exceptionChargeProd[i] = chargeProd.getValue();
179 } else {
180 exceptionChargeProd[i] = 0.0;
181 chargeExclusion[i] = true;
182 }
183 if (abs(eps.getValue()) > 0.0) {
184 vdWExclusion[i] = false;
185 exceptionEps[i] = eps.getValue();
186 } else {
187 vdWExclusion[i] = true;
188 exceptionEps[i] = 0.0;
189 }
190 }
191 }
192
193 Crystal crystal = openMMEnergy.getCrystal();
194 if (crystal.aperiodic()) {
195 setNonbondedMethod(OpenMMLibrary.OpenMM_NonbondedForce_NonbondedMethod.OpenMM_NonbondedForce_NoCutoff);
196 } else {
197 setNonbondedMethod(OpenMMLibrary.OpenMM_NonbondedForce_NonbondedMethod.OpenMM_NonbondedForce_PME);
198 if (pme != null) {
199
200
201 double aEwald = OpenMM_AngstromsPerNm * pme.getEwaldCoefficient();
202 int nx = pme.getReciprocalSpace().getXDim();
203 int ny = pme.getReciprocalSpace().getYDim();
204 int nz = pme.getReciprocalSpace().getZDim();
205 setPMEParameters(aEwald, nx, ny, nz);
206 }
207
208 NonbondedCutoff nonbondedCutoff = vdW.getNonbondedCutoff();
209 double off = nonbondedCutoff.off;
210 double cut = nonbondedCutoff.cut;
211 setCutoffDistance(OpenMM_NmPerAngstrom * off);
212 setUseSwitchingFunction(OpenMM_True);
213 if (cut == off) {
214 logger.warning(" OpenMM does not properly handle cutoffs where cut == off!");
215 if (cut == Double.MAX_VALUE || cut == Double.POSITIVE_INFINITY) {
216 logger.info(" Detected infinite or max-value cutoff; setting cut to 1E+40 for OpenMM.");
217 cut = 1E40;
218 } else {
219 logger.info(format(" Detected cut %8.4g == off %8.4g; scaling cut to 0.99 of off for OpenMM.", cut, off));
220 cut *= 0.99;
221 }
222 }
223 setSwitchingDistance(OpenMM_NmPerAngstrom * cut);
224 }
225
226 setUseDispersionCorrection(OpenMM_False);
227
228 ForceField forceField = openMMEnergy.getMolecularAssembly().getForceField();
229 int forceGroup = forceField.getInteger("VDW_FORCE_GROUP", 1);
230 int pmeGroup = forceField.getInteger("PME_FORCE_GROUP", 1);
231 if (forceGroup != pmeGroup) {
232 logger.severe(format(" ERROR: VDW-FORCE-GROUP is %d while PME-FORCE-GROUP is %d. "
233 + "This is invalid for fixed-charge force fields with combined non-bonded forces.",
234 forceGroup, pmeGroup));
235 }
236
237 setForceGroup(forceGroup);
238 logger.log(Level.INFO, format(" Fixed charge non-bonded force \t%1d", forceGroup));
239 }
240
241
242
243
244
245
246
247 public static Force constructForce(OpenMMEnergy openMMEnergy) {
248 VanDerWaals vdW = openMMEnergy.getVdwNode();
249 if (vdW == null) {
250
251 return null;
252 }
253
254
255
256
257
258 VanDerWaalsForm vdwForm = vdW.getVDWForm();
259 if (vdwForm.vdwType != LENNARD_JONES || vdwForm.radiusRule != ARITHMETIC
260 || vdwForm.epsilonRule != GEOMETRIC) {
261 logger.info(format(" VDW Type: %s", vdwForm.vdwType));
262 logger.info(format(" VDW Radius Rule: %s", vdwForm.radiusRule));
263 logger.info(format(" VDW Epsilon Rule: %s", vdwForm.epsilonRule));
264 logger.log(Level.SEVERE, " Unsupported van der Waals functional form.");
265 return null;
266 }
267
268 return new FixedChargeNonbondedForce(openMMEnergy);
269 }
270
271
272
273
274
275
276
277 public void updateForce(Atom[] atoms, OpenMMEnergy openMMEnergy) {
278 VanDerWaals vdW = openMMEnergy.getVdwNode();
279
280
281 VanDerWaalsForm vdwForm = vdW.getVDWForm();
282 if (vdwForm.vdwType != LENNARD_JONES || vdwForm.radiusRule != ARITHMETIC
283 || vdwForm.epsilonRule != GEOMETRIC) {
284 logger.log(Level.SEVERE, " Unsupported van der Waals functional form.");
285 return;
286 }
287
288
289 double radScale = 1.0;
290 if (vdwForm.radiusSize == RADIUS) {
291 radScale = 2.0;
292 }
293
294
295 if (vdwForm.radiusType == R_MIN) {
296 radScale /= 1.122462048309372981;
297 }
298
299
300 for (Atom atom : atoms) {
301 int index = atom.getXyzIndex() - 1;
302 boolean applyLambda = atom.applyLambda();
303
304 double charge = Double.MIN_VALUE;
305 MultipoleType multipoleType = atom.getMultipoleType();
306 if (multipoleType != null && atom.getElectrostatics()) {
307 charge = multipoleType.charge;
308 }
309
310 VDWType vdwType = atom.getVDWType();
311 double sigma = OpenMM_NmPerAngstrom * vdwType.radius * radScale;
312 double eps = OpenMM_KJPerKcal * vdwType.wellDepth;
313
314 if (applyLambda) {
315 OpenMMSystem system = openMMEnergy.getSystem();
316 boolean vdwLambdaTerm = vdW.getLambdaTerm();
317
318 if (vdwLambdaTerm) {
319 eps = 0.0;
320 }
321
322 double permLambda = openMMEnergy.getPmeNode().getAlchemicalParameters().permLambda;
323 charge *= permLambda;
324 }
325
326 if (!atom.getUse()) {
327 eps = 0.0;
328 charge = 0.0;
329 }
330
331 setParticleParameters(index, charge, sigma, eps);
332 }
333
334
335 IntByReference particle1 = new IntByReference();
336 IntByReference particle2 = new IntByReference();
337 DoubleByReference chargeProd = new DoubleByReference();
338 DoubleByReference sigma = new DoubleByReference();
339 DoubleByReference eps = new DoubleByReference();
340
341 int numExceptions = getNumExceptions();
342 for (int i = 0; i < numExceptions; i++) {
343
344
345 if (chargeExclusion[i] && vdWExclusion[i]) {
346 continue;
347 }
348
349 getExceptionParameters(i, particle1, particle2, chargeProd, sigma, eps);
350
351 int i1 = particle1.getValue();
352 int i2 = particle2.getValue();
353
354 double qq = exceptionChargeProd[i];
355 double epsilon = exceptionEps[i];
356
357 Atom atom1 = atoms[i1];
358 Atom atom2 = atoms[i2];
359
360
361
362
363
364 double minEpsilon = 1.0e-12;
365
366 OpenMMSystem system = openMMEnergy.getSystem();
367
368 ParticleMeshEwald pme = openMMEnergy.getPmeNode();
369 double lambdaElec = pme.getAlchemicalParameters().permLambda;
370 boolean vdwLambdaTerm = vdW.getLambdaTerm();
371
372 if (lambdaElec < minEpsilon) {
373 lambdaElec = minEpsilon;
374 }
375
376 if (atom1.applyLambda()) {
377 qq *= lambdaElec;
378 if (vdwLambdaTerm) {
379 epsilon = minEpsilon;
380 }
381 }
382 if (atom2.applyLambda()) {
383 qq *= lambdaElec;
384 if (vdwLambdaTerm) {
385 epsilon = minEpsilon;
386 }
387 }
388 if (!atom1.getUse() || !atom2.getUse()) {
389 qq = minEpsilon;
390 epsilon = minEpsilon;
391 }
392
393 setExceptionParameters(i, i1, i2, qq, sigma.getValue(), epsilon);
394 }
395
396
397 updateParametersInContext(openMMEnergy.getContext());
398 }
399
400 }