1 // ****************************************************************************** 2 // 3 // Title: Force Field X. 4 // Description: Force Field X - Software for Molecular Biophysics. 5 // Copyright: Copyright (c) Michael J. Schnieders 2001-2025. 6 // 7 // This file is part of Force Field X. 8 // 9 // Force Field X is free software; you can redistribute it and/or modify it 10 // under the terms of the GNU General Public License version 3 as published by 11 // the Free Software Foundation. 12 // 13 // Force Field X is distributed in the hope that it will be useful, but WITHOUT 14 // ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS 15 // FOR A PARTICULAR PURPOSE. See the GNU General Public License for more 16 // details. 17 // 18 // You should have received a copy of the GNU General Public License along with 19 // Force Field X; if not, write to the Free Software Foundation, Inc., 59 Temple 20 // Place, Suite 330, Boston, MA 02111-1307 USA 21 // 22 // Linking this library statically or dynamically with other modules is making a 23 // combined work based on this library. Thus, the terms and conditions of the 24 // GNU General Public License cover the whole combination. 25 // 26 // As a special exception, the copyright holders of this library give you 27 // permission to link this library with independent modules to produce an 28 // executable, regardless of the license terms of these independent modules, and 29 // to copy and distribute the resulting executable under terms of your choice, 30 // provided that you also meet, for each linked independent module, the terms 31 // and conditions of the license of that module. An independent module is a 32 // module which is not derived from or based on this library. If you modify this 33 // library, you may extend this exception to your version of the library, but 34 // you are not obligated to do so. If you do not wish to do so, delete this 35 // exception statement from your version. 36 // 37 // ****************************************************************************** 38 package ffx.openmm.amoeba; 39 40 import com.sun.jna.ptr.DoubleByReference; 41 import com.sun.jna.ptr.IntByReference; 42 import com.sun.jna.ptr.PointerByReference; 43 import ffx.openmm.Context; 44 import ffx.openmm.Force; 45 46 import java.nio.DoubleBuffer; 47 import java.nio.IntBuffer; 48 49 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_addException; 50 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_addParticle; 51 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_create; 52 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_destroy; 53 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getCutoffDistance; 54 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getDPMEParameters; 55 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getDPMEParametersInContext; 56 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getEwaldErrorTolerance; 57 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getExceptionParameters; 58 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getExtrapolationCoefficients; 59 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getInducedDipoles; 60 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getLabFramePermanentDipoles; 61 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getNonbondedMethod; 62 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getNumExceptions; 63 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getNumParticles; 64 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getPMEParameters; 65 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getPMEParametersInContext; 66 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getParticleParameters; 67 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_getSwitchingDistance; 68 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_setCutoffDistance; 69 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_setDPMEParameters; 70 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_setEwaldErrorTolerance; 71 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_setExceptionParameters; 72 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_setExtrapolationCoefficients; 73 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_setNonbondedMethod; 74 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_setPMEParameters; 75 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_setParticleParameters; 76 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_setSwitchingDistance; 77 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_updateParametersInContext; 78 import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_HippoNonbondedForce_usesPeriodicBoundaryConditions; 79 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_Boolean.OpenMM_True; 80 81 /** 82 * This class implements all nonbonded interactions in the HIPPO force field: electrostatics, 83 * induction, charge transfer, dispersion, and repulsion. Although some of these are 84 * conceptually distinct, they share parameters in common and are most efficiently computed 85 * together. For example, the same multipole definitions are used for both electrostatics 86 * and Pauli repulsion. Therefore, all of them are computed by a single Force object. 87 * <p> 88 * To use it, create a HippoNonbondedForce object, then call addParticle() once for each particle. 89 * After an entry has been added, you can modify its force field parameters by calling setParticleParameters(). 90 * This will have no effect on Contexts that already exist unless you call updateParametersInContext(). 91 * <p> 92 * You also can specify "exceptions", particular pairs of particles whose interactions should be 93 * reduced or completely omitted. Call addException() to define exceptions. 94 */ 95 public class HippoNonbondedForce extends Force { 96 97 /** 98 * Create a new HippoNonbondedForce. 99 */ 100 public HippoNonbondedForce() { 101 super(OpenMM_HippoNonbondedForce_create()); 102 } 103 104 /** 105 * Add an interaction to the list of exceptions that should be calculated differently from other interactions. 106 * If all scale factors are set to 0, this will cause the interaction to be completely omitted from 107 * force and energy calculations. 108 * 109 * @param particle1 the index of the first particle involved in the interaction 110 * @param particle2 the index of the second particle involved in the interaction 111 * @param multipoleMultipoleScale the factor by which to scale the Coulomb interaction between fixed multipoles 112 * @param dipoleMultipoleScale the factor by which to scale the Coulomb interaction between an induced dipole and a fixed multipole 113 * @param dipoleDipoleScale the factor by which to scale the Coulomb interaction between induced dipoles 114 * @param dispersionScale the factor by which to scale the dispersion interaction 115 * @param repulsionScale the factor by which to scale the Pauli repulsion 116 * @param chargeTransferScale the factor by which to scale the charge transfer interaction 117 * @param replace determines the behavior if there is already an exception for the same two particles. If true, the existing one is replaced. If false, an exception is thrown. 118 * @return the index of the exception that was added 119 */ 120 public int addException(int particle1, int particle2, double multipoleMultipoleScale, double dipoleMultipoleScale, 121 double dipoleDipoleScale, double dispersionScale, double repulsionScale, double chargeTransferScale, int replace) { 122 return OpenMM_HippoNonbondedForce_addException(pointer, particle1, particle2, multipoleMultipoleScale, 123 dipoleMultipoleScale, dipoleDipoleScale, dispersionScale, repulsionScale, chargeTransferScale, replace); 124 } 125 126 /** 127 * Add the nonbonded force parameters for a particle. This should be called once for each particle 128 * in the System. When it is called for the i'th time, it specifies the parameters for the i'th particle. 129 * 130 * @param charge the particle's charge 131 * @param dipole the particle's molecular dipole (vector of size 3) 132 * @param quadrupole the particle's molecular quadrupole (vector of size 9) 133 * @param coreCharge the charge of the atomic core 134 * @param alpha controls the width of the particle's electron density 135 * @param epsilon sets the magnitude of charge transfer 136 * @param damping sets the length scale for charge transfer 137 * @param c6 the coefficient of the dispersion interaction 138 * @param pauliK the coefficient of the Pauli repulsion interaction 139 * @param pauliQ the charge used in computing the Pauli repulsion interaction 140 * @param pauliAlpha the width of the particle's electron density for computing the Pauli repulsion interaction 141 * @param polarizability atomic polarizability 142 * @param axisType the particle's axis type 143 * @param multipoleAtomZ index of first atom used in defining the local coordinate system for multipoles 144 * @param multipoleAtomX index of second atom used in defining the local coordinate system for multipoles 145 * @param multipoleAtomY index of third atom used in defining the local coordinate system for multipoles 146 * @return the index of the particle that was added 147 */ 148 public int addParticle(double charge, PointerByReference dipole, PointerByReference quadrupole, double coreCharge, 149 double alpha, double epsilon, double damping, double c6, double pauliK, double pauliQ, double pauliAlpha, 150 double polarizability, int axisType, int multipoleAtomZ, int multipoleAtomX, int multipoleAtomY) { 151 return OpenMM_HippoNonbondedForce_addParticle(pointer, charge, dipole, quadrupole, coreCharge, 152 alpha, epsilon, damping, c6, pauliK, pauliQ, pauliAlpha, 153 polarizability, axisType, multipoleAtomZ, multipoleAtomX, multipoleAtomY); 154 } 155 156 /** 157 * Destroy the force. 158 */ 159 @Override 160 public void destroy() { 161 if (pointer != null) { 162 OpenMM_HippoNonbondedForce_destroy(pointer); 163 pointer = null; 164 } 165 } 166 167 /** 168 * Get the cutoff distance. 169 * 170 * @return The cutoff distance, measured in nm. 171 */ 172 public double getCutoffDistance() { 173 return OpenMM_HippoNonbondedForce_getCutoffDistance(pointer); 174 } 175 176 /** 177 * Get the DPME parameters. 178 * 179 * @param alpha The Ewald alpha parameter (output). 180 * @param nx The number of grid points along the x axis (output). 181 * @param ny The number of grid points along the y axis (output). 182 * @param nz The number of grid points along the z axis (output). 183 */ 184 public void getDPMEParameters(DoubleByReference alpha, IntByReference nx, IntByReference ny, IntByReference nz) { 185 OpenMM_HippoNonbondedForce_getDPMEParameters(pointer, alpha, nx, ny, nz); 186 } 187 188 /** 189 * Get the DPME parameters. 190 * 191 * @param alpha The Ewald alpha parameter (output). 192 * @param nx The number of grid points along the x axis (output). 193 * @param ny The number of grid points along the y axis (output). 194 * @param nz The number of grid points along the z axis (output). 195 */ 196 public void getDPMEParameters(DoubleBuffer alpha, IntBuffer nx, IntBuffer ny, IntBuffer nz) { 197 OpenMM_HippoNonbondedForce_getDPMEParameters(pointer, alpha, nx, ny, nz); 198 } 199 200 /** 201 * Get the DPME parameters in context. 202 * 203 * @param context The context. 204 * @param alpha The Ewald alpha parameter (output). 205 * @param nx The number of grid points along the x axis (output). 206 * @param ny The number of grid points along the y axis (output). 207 * @param nz The number of grid points along the z axis (output). 208 */ 209 public void getDPMEParametersInContext(Context context, DoubleByReference alpha, IntByReference nx, 210 IntByReference ny, IntByReference nz) { 211 OpenMM_HippoNonbondedForce_getDPMEParametersInContext(pointer, context.getPointer(), alpha, nx, ny, nz); 212 } 213 214 /** 215 * Get the DPME parameters in context. 216 * 217 * @param context The context. 218 * @param alpha The Ewald alpha parameter (output). 219 * @param nx The number of grid points along the x axis (output). 220 * @param ny The number of grid points along the y axis (output). 221 * @param nz The number of grid points along the z axis (output). 222 */ 223 public void getDPMEParametersInContext(Context context, DoubleBuffer alpha, IntBuffer nx, 224 IntBuffer ny, IntBuffer nz) { 225 OpenMM_HippoNonbondedForce_getDPMEParametersInContext(pointer, context.getPointer(), alpha, nx, ny, nz); 226 } 227 228 /** 229 * Get the Ewald error tolerance. 230 * 231 * @return The Ewald error tolerance. 232 */ 233 public double getEwaldErrorTolerance() { 234 return OpenMM_HippoNonbondedForce_getEwaldErrorTolerance(pointer); 235 } 236 237 /** 238 * Get the scale factors for an interaction that should be calculated differently from others. 239 * 240 * @param index the index of the interaction for which to get parameters 241 * @param particle1 the index of the first particle involved in the interaction 242 * @param particle2 the index of the second particle involved in the interaction 243 * @param multipoleMultipoleScale the factor by which to scale the Coulomb interaction between fixed multipoles 244 * @param dipoleMultipoleScale the factor by which to scale the Coulomb interaction between an induced dipole and a fixed multipole 245 * @param dipoleDipoleScale the factor by which to scale the Coulomb interaction between induced dipoles 246 * @param dispersionScale the factor by which to scale the dispersion interaction 247 * @param repulsionScale the factor by which to scale the Pauli repulsion 248 * @param chargeTransferScale the factor by which to scale the charge transfer interaction 249 */ 250 public void getExceptionParameters(int index, IntByReference particle1, IntByReference particle2, 251 DoubleByReference multipoleMultipoleScale, DoubleByReference dipoleMultipoleScale, 252 DoubleByReference dipoleDipoleScale, DoubleByReference dispersionScale, 253 DoubleByReference repulsionScale, DoubleByReference chargeTransferScale) { 254 OpenMM_HippoNonbondedForce_getExceptionParameters(pointer, index, particle1, particle2, 255 multipoleMultipoleScale, dipoleMultipoleScale, dipoleDipoleScale, dispersionScale, repulsionScale, chargeTransferScale); 256 } 257 258 /** 259 * Get the scale factors for an interaction that should be calculated differently from others. 260 * 261 * @param index the index of the interaction for which to get parameters 262 * @param particle1 the index of the first particle involved in the interaction 263 * @param particle2 the index of the second particle involved in the interaction 264 * @param multipoleMultipoleScale the factor by which to scale the Coulomb interaction between fixed multipoles 265 * @param dipoleMultipoleScale the factor by which to scale the Coulomb interaction between an induced dipole and a fixed multipole 266 * @param dipoleDipoleScale the factor by which to scale the Coulomb interaction between induced dipoles 267 * @param dispersionScale the factor by which to scale the dispersion interaction 268 * @param repulsionScale the factor by which to scale the Pauli repulsion 269 * @param chargeTransferScale the factor by which to scale the charge transfer interaction 270 */ 271 public void getExceptionParameters(int index, IntBuffer particle1, IntBuffer particle2, 272 DoubleBuffer multipoleMultipoleScale, DoubleBuffer dipoleMultipoleScale, 273 DoubleBuffer dipoleDipoleScale, DoubleBuffer dispersionScale, 274 DoubleBuffer repulsionScale, DoubleBuffer chargeTransferScale) { 275 OpenMM_HippoNonbondedForce_getExceptionParameters(pointer, index, particle1, particle2, 276 multipoleMultipoleScale, dipoleMultipoleScale, dipoleDipoleScale, dispersionScale, repulsionScale, chargeTransferScale); 277 } 278 279 /** 280 * Get the extrapolation coefficients. 281 * 282 * @return The extrapolation coefficients. 283 */ 284 public PointerByReference getExtrapolationCoefficients() { 285 return OpenMM_HippoNonbondedForce_getExtrapolationCoefficients(pointer); 286 } 287 288 /** 289 * Get the induced dipoles. 290 * 291 * @param context The context. 292 * @param dipoles The induced dipoles (output). 293 */ 294 public void getInducedDipoles(Context context, PointerByReference dipoles) { 295 OpenMM_HippoNonbondedForce_getInducedDipoles(pointer, context.getPointer(), dipoles); 296 } 297 298 /** 299 * Get the lab frame permanent dipoles. 300 * 301 * @param context The context. 302 * @param dipoles The lab frame permanent dipoles (output). 303 */ 304 public void getLabFramePermanentDipoles(Context context, PointerByReference dipoles) { 305 OpenMM_HippoNonbondedForce_getLabFramePermanentDipoles(pointer, context.getPointer(), dipoles); 306 } 307 308 /** 309 * Get the nonbonded method. 310 * 311 * @return The nonbonded method. 312 */ 313 public int getNonbondedMethod() { 314 return OpenMM_HippoNonbondedForce_getNonbondedMethod(pointer); 315 } 316 317 /** 318 * Get the number of exceptions. 319 * 320 * @return The number of exceptions. 321 */ 322 public int getNumExceptions() { 323 return OpenMM_HippoNonbondedForce_getNumExceptions(pointer); 324 } 325 326 /** 327 * Get the number of particles. 328 * 329 * @return The number of particles. 330 */ 331 public int getNumParticles() { 332 return OpenMM_HippoNonbondedForce_getNumParticles(pointer); 333 } 334 335 /** 336 * Get the PME parameters. 337 * 338 * @param alpha The Ewald alpha parameter (output). 339 * @param nx The number of grid points along the x axis (output). 340 * @param ny The number of grid points along the y axis (output). 341 * @param nz The number of grid points along the z axis (output). 342 */ 343 public void getPMEParameters(DoubleByReference alpha, IntByReference nx, IntByReference ny, IntByReference nz) { 344 OpenMM_HippoNonbondedForce_getPMEParameters(pointer, alpha, nx, ny, nz); 345 } 346 347 /** 348 * Get the PME parameters. 349 * 350 * @param alpha The Ewald alpha parameter (output). 351 * @param nx The number of grid points along the x axis (output). 352 * @param ny The number of grid points along the y axis (output). 353 * @param nz The number of grid points along the z axis (output). 354 */ 355 public void getPMEParameters(DoubleBuffer alpha, IntBuffer nx, IntBuffer ny, IntBuffer nz) { 356 OpenMM_HippoNonbondedForce_getPMEParameters(pointer, alpha, nx, ny, nz); 357 } 358 359 /** 360 * Get the PME parameters in context. 361 * 362 * @param context The context. 363 * @param alpha The Ewald alpha parameter (output). 364 * @param nx The number of grid points along the x axis (output). 365 * @param ny The number of grid points along the y axis (output). 366 * @param nz The number of grid points along the z axis (output). 367 */ 368 public void getPMEParametersInContext(Context context, DoubleByReference alpha, IntByReference nx, 369 IntByReference ny, IntByReference nz) { 370 OpenMM_HippoNonbondedForce_getPMEParametersInContext(pointer, context.getPointer(), alpha, nx, ny, nz); 371 } 372 373 /** 374 * Get the PME parameters in context. 375 * 376 * @param context The context. 377 * @param alpha The Ewald alpha parameter (output). 378 * @param nx The number of grid points along the x axis (output). 379 * @param ny The number of grid points along the y axis (output). 380 * @param nz The number of grid points along the z axis (output). 381 */ 382 public void getPMEParametersInContext(Context context, DoubleBuffer alpha, IntBuffer nx, 383 IntBuffer ny, IntBuffer nz) { 384 OpenMM_HippoNonbondedForce_getPMEParametersInContext(pointer, context.getPointer(), alpha, nx, ny, nz); 385 } 386 387 /** 388 * Get the nonbonded force parameters for a particle. 389 * 390 * @param index the index of the particle for which to get parameters 391 * @param charge the particle's charge 392 * @param dipole the particle's molecular dipole (vector of size 3) 393 * @param quadrupole the particle's molecular quadrupole (vector of size 9) 394 * @param coreCharge the charge of the atomic core 395 * @param alpha controls the width of the particle's electron density 396 * @param epsilon sets the magnitude of charge transfer 397 * @param damping sets the length scale for charge transfer 398 * @param c6 the coefficient of the dispersion interaction 399 * @param pauliK the coefficient of the Pauli repulsion interaction 400 * @param pauliQ the charge used in computing the Pauli repulsion interaction 401 * @param pauliAlpha the width of the particle's electron density for computing the Pauli repulsion interaction 402 * @param polarizability atomic polarizability 403 * @param axisType the particle's axis type 404 * @param multipoleAtomZ index of first atom used in defining the local coordinate system for multipoles 405 * @param multipoleAtomX index of second atom used in defining the local coordinate system for multipoles 406 * @param multipoleAtomY index of third atom used in defining the local coordinate system for multipoles 407 */ 408 public void getParticleParameters(int index, DoubleByReference charge, PointerByReference dipole, 409 PointerByReference quadrupole, DoubleByReference coreCharge, 410 DoubleByReference alpha, DoubleByReference epsilon, 411 DoubleByReference damping, DoubleByReference c6, 412 DoubleByReference pauliK, DoubleByReference pauliQ, DoubleByReference pauliAlpha, 413 DoubleByReference polarizability, IntByReference axisType, 414 IntByReference multipoleAtomZ, IntByReference multipoleAtomX, IntByReference multipoleAtomY) { 415 OpenMM_HippoNonbondedForce_getParticleParameters(pointer, index, charge, dipole, quadrupole, 416 coreCharge, alpha, epsilon, damping, c6, pauliK, pauliQ, pauliAlpha, 417 polarizability, axisType, multipoleAtomZ, multipoleAtomX, multipoleAtomY); 418 } 419 420 /** 421 * Get the nonbonded force parameters for a particle. 422 * 423 * @param index the index of the particle for which to get parameters 424 * @param charge the particle's charge 425 * @param dipole the particle's molecular dipole (vector of size 3) 426 * @param quadrupole the particle's molecular quadrupole (vector of size 9) 427 * @param coreCharge the charge of the atomic core 428 * @param alpha controls the width of the particle's electron density 429 * @param epsilon sets the magnitude of charge transfer 430 * @param damping sets the length scale for charge transfer 431 * @param c6 the coefficient of the dispersion interaction 432 * @param pauliK the coefficient of the Pauli repulsion interaction 433 * @param pauliQ the charge used in computing the Pauli repulsion interaction 434 * @param pauliAlpha the width of the particle's electron density for computing the Pauli repulsion interaction 435 * @param polarizability atomic polarizability 436 * @param axisType the particle's axis type 437 * @param multipoleAtomZ index of first atom used in defining the local coordinate system for multipoles 438 * @param multipoleAtomX index of second atom used in defining the local coordinate system for multipoles 439 * @param multipoleAtomY index of third atom used in defining the local coordinate system for multipoles 440 */ 441 public void getParticleParameters(int index, DoubleBuffer charge, PointerByReference dipole, 442 PointerByReference quadrupole, DoubleBuffer coreCharge, 443 DoubleBuffer alpha, DoubleBuffer epsilon, 444 DoubleBuffer damping, DoubleBuffer c6, 445 DoubleBuffer pauliK, DoubleBuffer pauliQ, DoubleBuffer pauliAlpha, 446 DoubleBuffer polarizability, IntBuffer axisType, 447 IntBuffer multipoleAtomZ, IntBuffer multipoleAtomX, IntBuffer multipoleAtomY) { 448 OpenMM_HippoNonbondedForce_getParticleParameters(pointer, index, charge, dipole, quadrupole, 449 coreCharge, alpha, epsilon, damping, c6, pauliK, pauliQ, pauliAlpha, 450 polarizability, axisType, multipoleAtomZ, multipoleAtomX, multipoleAtomY); 451 } 452 453 /** 454 * Get the switching distance. 455 * 456 * @return The switching distance, measured in nm. 457 */ 458 public double getSwitchingDistance() { 459 return OpenMM_HippoNonbondedForce_getSwitchingDistance(pointer); 460 } 461 462 /** 463 * Set the cutoff distance. 464 * 465 * @param distance The cutoff distance, measured in nm. 466 */ 467 public void setCutoffDistance(double distance) { 468 OpenMM_HippoNonbondedForce_setCutoffDistance(pointer, distance); 469 } 470 471 /** 472 * Set the DPME parameters. 473 * 474 * @param alpha The Ewald alpha parameter. 475 * @param nx The number of grid points along the x axis. 476 * @param ny The number of grid points along the y axis. 477 * @param nz The number of grid points along the z axis. 478 */ 479 public void setDPMEParameters(double alpha, int nx, int ny, int nz) { 480 OpenMM_HippoNonbondedForce_setDPMEParameters(pointer, alpha, nx, ny, nz); 481 } 482 483 /** 484 * Set the Ewald error tolerance. 485 * 486 * @param tolerance The Ewald error tolerance. 487 */ 488 public void setEwaldErrorTolerance(double tolerance) { 489 OpenMM_HippoNonbondedForce_setEwaldErrorTolerance(pointer, tolerance); 490 } 491 492 /** 493 * Set the scale factors for an interaction that should be calculated differently from others. 494 * 495 * @param index the index of the interaction for which to set parameters 496 * @param particle1 the index of the first particle involved in the interaction 497 * @param particle2 the index of the second particle involved in the interaction 498 * @param multipoleMultipoleScale the factor by which to scale the Coulomb interaction between fixed multipoles 499 * @param dipoleMultipoleScale the factor by which to scale the Coulomb interaction between an induced dipole and a fixed multipole 500 * @param dipoleDipoleScale the factor by which to scale the Coulomb interaction between induced dipoles 501 * @param dispersionScale the factor by which to scale the dispersion interaction 502 * @param repulsionScale the factor by which to scale the Pauli repulsion 503 * @param chargeTransferScale the factor by which to scale the charge transfer interaction 504 */ 505 public void setExceptionParameters(int index, int particle1, int particle2, double multipoleMultipoleScale, 506 double dipoleMultipoleScale, double dipoleDipoleScale, double dispersionScale, 507 double repulsionScale, double chargeTransferScale) { 508 OpenMM_HippoNonbondedForce_setExceptionParameters(pointer, index, particle1, particle2, 509 multipoleMultipoleScale, dipoleMultipoleScale, dipoleDipoleScale, dispersionScale, repulsionScale, chargeTransferScale); 510 } 511 512 /** 513 * Set the extrapolation coefficients. 514 * 515 * @param coefficients The extrapolation coefficients. 516 */ 517 public void setExtrapolationCoefficients(PointerByReference coefficients) { 518 OpenMM_HippoNonbondedForce_setExtrapolationCoefficients(pointer, coefficients); 519 } 520 521 /** 522 * Set the nonbonded method. 523 * 524 * @param method The nonbonded method. 525 */ 526 public void setNonbondedMethod(int method) { 527 OpenMM_HippoNonbondedForce_setNonbondedMethod(pointer, method); 528 } 529 530 /** 531 * Set the PME parameters. 532 * 533 * @param alpha The Ewald alpha parameter. 534 * @param nx The number of grid points along the x axis. 535 * @param ny The number of grid points along the y axis. 536 * @param nz The number of grid points along the z axis. 537 */ 538 public void setPMEParameters(double alpha, int nx, int ny, int nz) { 539 OpenMM_HippoNonbondedForce_setPMEParameters(pointer, alpha, nx, ny, nz); 540 } 541 542 /** 543 * Set the nonbonded force parameters for a particle. 544 * 545 * @param index the index of the particle for which to set parameters 546 * @param charge the particle's charge 547 * @param dipole the particle's molecular dipole (vector of size 3) 548 * @param quadrupole the particle's molecular quadrupole (vector of size 9) 549 * @param coreCharge the charge of the atomic core 550 * @param alpha controls the width of the particle's electron density 551 * @param epsilon sets the magnitude of charge transfer 552 * @param damping sets the length scale for charge transfer 553 * @param c6 the coefficient of the dispersion interaction 554 * @param pauliK the coefficient of the Pauli repulsion interaction 555 * @param pauliQ the charge used in computing the Pauli repulsion interaction 556 * @param pauliAlpha the width of the particle's electron density for computing the Pauli repulsion interaction 557 * @param polarizability atomic polarizability 558 * @param axisType the particle's axis type 559 * @param multipoleAtomZ index of first atom used in defining the local coordinate system for multipoles 560 * @param multipoleAtomX index of second atom used in defining the local coordinate system for multipoles 561 * @param multipoleAtomY index of third atom used in defining the local coordinate system for multipoles 562 */ 563 public void setParticleParameters(int index, double charge, PointerByReference dipole, 564 PointerByReference quadrupole, double coreCharge, double alpha, double epsilon, 565 double damping, double c6, double pauliK, double pauliQ, double pauliAlpha, 566 double polarizability, int axisType, int multipoleAtomZ, 567 int multipoleAtomX, int multipoleAtomY) { 568 OpenMM_HippoNonbondedForce_setParticleParameters(pointer, index, charge, dipole, quadrupole, 569 coreCharge, alpha, epsilon, damping, c6, pauliK, pauliQ, pauliAlpha, 570 polarizability, axisType, multipoleAtomZ, multipoleAtomX, multipoleAtomY); 571 } 572 573 /** 574 * Set the switching distance. 575 * 576 * @param distance The switching distance, measured in nm. 577 */ 578 public void setSwitchingDistance(double distance) { 579 OpenMM_HippoNonbondedForce_setSwitchingDistance(pointer, distance); 580 } 581 582 /** 583 * Update the parameters in a Context to match those stored in this Force object. 584 * 585 * @param context The Context in which to update the parameters. 586 */ 587 public void updateParametersInContext(Context context) { 588 if (context.hasContextPointer()) { 589 OpenMM_HippoNonbondedForce_updateParametersInContext(pointer, context.getPointer()); 590 } 591 } 592 593 /** 594 * Check if the force uses periodic boundary conditions. 595 * 596 * @return True if the force uses periodic boundary conditions. 597 */ 598 @Override 599 public boolean usesPeriodicBoundaryConditions() { 600 int pbc = OpenMM_HippoNonbondedForce_usesPeriodicBoundaryConditions(pointer); 601 return pbc == OpenMM_True; 602 } 603 }