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