View Javadoc
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.DoubleArray;
45  import ffx.openmm.Force;
46  import ffx.openmm.IntArray;
47  
48  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_addMultipole;
49  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_create;
50  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_destroy;
51  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getAEwald;
52  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getCovalentMap;
53  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getCovalentMaps;
54  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getCutoffDistance;
55  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getElectrostaticPotential;
56  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getEwaldErrorTolerance;
57  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getExtrapolationCoefficients;
58  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getInducedDipoles;
59  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getLabFramePermanentDipoles;
60  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getMultipoleParameters;
61  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getMutualInducedMaxIterations;
62  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getMutualInducedTargetEpsilon;
63  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getNonbondedMethod;
64  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getNumMultipoles;
65  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getPMEParameters;
66  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getPMEParametersInContext;
67  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getPmeBSplineOrder;
68  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getPmeGridDimensions;
69  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getPolarizationType;
70  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getSystemMultipoleMoments;
71  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_getTotalDipoles;
72  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_setAEwald;
73  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_setCovalentMap;
74  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_setCutoffDistance;
75  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_setEwaldErrorTolerance;
76  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_setExtrapolationCoefficients;
77  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_setMultipoleParameters;
78  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_setMutualInducedMaxIterations;
79  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_setMutualInducedTargetEpsilon;
80  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_setNonbondedMethod;
81  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_setPMEParameters;
82  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_setPmeGridDimensions;
83  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_setPolarizationType;
84  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_updateParametersInContext;
85  import static edu.uiowa.jopenmm.OpenMMAmoebaLibrary.OpenMM_AmoebaMultipoleForce_usesPeriodicBoundaryConditions;
86  import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_Boolean.OpenMM_True;
87  
88  /**
89   * This class implements the Amoeba multipole interaction.
90   * <p>
91   * To use it, create an AmoebaMultipoleForce object then call addMultipole() once for each atom.  After
92   * an entry has been added, you can modify its force field parameters by calling setMultipoleParameters().
93   * This will have no effect on Contexts that already exist unless you call updateParametersInContext().
94   */
95  public class MultipoleForce extends Force {
96  
97    public MultipoleForce() {
98      super(OpenMM_AmoebaMultipoleForce_create());
99    }
100 
101   /**
102    * Add multipole-related info for a particle
103    *
104    * @param charge              the particle's charge
105    * @param molecularDipole     the particle's molecular dipole (vector of size 3)
106    * @param molecularQuadrupole the particle's molecular quadrupole (vector of size 9)
107    * @param axisType            the particle's axis type
108    * @param multipoleAtomZ      index of first atom used in constructing lab to molecular frames
109    * @param multipoleAtomX      index of second atom used in constructing lab to molecular frames
110    * @param multipoleAtomY      index of second atom used in constructing lab to molecular frames
111    * @param thole               Thole parameter
112    * @param dampingFactor       dampingFactor parameter
113    * @param polarity            polarity parameter
114    * @return the index of the particle that was added
115    */
116   public int addMultipole(double charge, DoubleArray molecularDipole, DoubleArray molecularQuadrupole, int axisType,
117                           int multipoleAtomZ, int multipoleAtomX, int multipoleAtomY, double thole, double dampingFactor, double polarity) {
118     return OpenMM_AmoebaMultipoleForce_addMultipole(pointer, charge, molecularDipole.getPointer(), molecularQuadrupole.getPointer(),
119         axisType, multipoleAtomZ, multipoleAtomX, multipoleAtomY, thole, dampingFactor, polarity);
120   }
121 
122   /**
123    * Destroy the force.
124    */
125   @Override
126   public void destroy() {
127     if (pointer != null) {
128       OpenMM_AmoebaMultipoleForce_destroy(pointer);
129       pointer = null;
130     }
131   }
132 
133   /**
134    * Get the Ewald alpha parameter.  If this is 0 (the default), a value is chosen automatically
135    * based on the Ewald error tolerance.
136    *
137    * @return the Ewald alpha parameter
138    * @deprecated This method exists only for backward compatibility.  Use getPMEParameters() instead.
139    */
140   @Deprecated
141   public double getAEwald() {
142     return OpenMM_AmoebaMultipoleForce_getAEwald(pointer);
143   }
144 
145   /**
146    * Get the covalent map for a given atom index and covalent type.
147    *
148    * @param i            The atom index.
149    * @param covalentType The covalent type.
150    * @return An IntArray representing the covalent map for the specified atom index and covalent type.
151    */
152   public IntArray getCovalentMap(int i, int covalentType) {
153     IntArray covalentMap = new IntArray(0);
154     if (pointer != null) {
155       OpenMM_AmoebaMultipoleForce_getCovalentMap(pointer, i, covalentType, covalentMap.getPointer());
156     }
157     return covalentMap;
158   }
159 
160   /**
161    * Get all covalent maps for a given atom index.
162    *
163    * @param i The atom index.
164    * @return An IntArray containing all covalent maps for the specified atom.
165    */
166   public IntArray getCovalentMaps(int i) {
167     IntArray covalentMaps = new IntArray(0);
168     if (pointer != null) {
169       OpenMM_AmoebaMultipoleForce_getCovalentMaps(pointer, i, covalentMaps.getPointer());
170     }
171     return covalentMaps;
172   }
173 
174   /**
175    * Get the cutoff distance (in nm) being used for nonbonded interactions.  If the NonbondedMethod in use
176    * is NoCutoff, this value will have no effect.
177    *
178    * @return the cutoff distance, measured in nm
179    */
180   public double getCutoffDistance() {
181     return OpenMM_AmoebaMultipoleForce_getCutoffDistance(pointer);
182   }
183 
184   /**
185    * Get the electrostatic potential at specified points.
186    *
187    * @param context The OpenMM context.
188    * @param points  The points at which to calculate the potential.
189    * @return A DoubleArray containing the electrostatic potential at each point.
190    */
191   public DoubleArray getElectrostaticPotential(Context context, DoubleArray points) {
192     DoubleArray potential = new DoubleArray(0);
193     if (context.hasContextPointer() && pointer != null) {
194       OpenMM_AmoebaMultipoleForce_getElectrostaticPotential(pointer, context.getPointer(),
195           points.getPointer(), potential.getPointer());
196     }
197     return potential;
198   }
199 
200   /**
201    * Get the error tolerance for Ewald summation.  This corresponds to the fractional error in the forces
202    * which is acceptable.  This value is used to select the grid dimensions and separation (alpha)
203    * parameter so that the average error level will be less than the tolerance.  There is not a
204    * rigorous guarantee that all forces on all atoms will be less than the tolerance, however.
205    * <p>
206    * This can be overridden by explicitly setting an alpha parameter and grid dimensions to use.
207    */
208   public double getEwaldErrorTolerance() {
209     return OpenMM_AmoebaMultipoleForce_getEwaldErrorTolerance(pointer);
210   }
211 
212   /**
213    * Get the extrapolation coefficients.
214    *
215    * @return A PointerByReference to the extrapolation coefficients.
216    */
217   public PointerByReference getExtrapolationCoefficients() {
218     return OpenMM_AmoebaMultipoleForce_getExtrapolationCoefficients(pointer);
219   }
220 
221   /**
222    * Get the induced dipoles.
223    *
224    * @param context        The OpenMM context.
225    * @param inducedDipoles The induced dipoles (output).
226    */
227   public void getInducedDipoles(Context context, DoubleArray inducedDipoles) {
228     if (context.hasContextPointer() && pointer != null) {
229       OpenMM_AmoebaMultipoleForce_getInducedDipoles(pointer, context.getPointer(), inducedDipoles.getPointer());
230     }
231   }
232 
233   /**
234    * Get the lab frame permanent dipoles.
235    *
236    * @param context The OpenMM context.
237    * @param dipoles The lab frame permanent dipoles (output).
238    */
239   public void getLabFramePermanentDipoles(Context context, DoubleArray dipoles) {
240     if (context.hasContextPointer() && pointer != null) {
241       OpenMM_AmoebaMultipoleForce_getLabFramePermanentDipoles(pointer, context.getPointer(), dipoles.getPointer());
242     }
243   }
244 
245   /**
246    * Get the multipole parameters for a particle.
247    *
248    * @param index               the index of the atom for which to get parameters
249    * @param charge              the particle's charge
250    * @param molecularDipole     the particle's molecular dipole (vector of size 3)
251    * @param molecularQuadrupole the particle's molecular quadrupole (vector of size 9)
252    * @param axisType            the particle's axis type
253    * @param multipoleAtomZ      index of first atom used in constructing lab to molecular frames
254    * @param multipoleAtomX      index of second atom used in constructing lab to molecular frames
255    * @param multipoleAtomY      index of second atom used in constructing lab to molecular frames
256    * @param thole               Thole parameter
257    * @param dampingFactor       dampingFactor parameter
258    * @param polarity            polarity parameter
259    */
260   public void getMultipoleParameters(int index, DoubleByReference charge, PointerByReference molecularDipole,
261                                      PointerByReference molecularQuadrupole, IntByReference axisType,
262                                      IntByReference multipoleAtomZ, IntByReference multipoleAtomX, IntByReference multipoleAtomY,
263                                      DoubleByReference thole, DoubleByReference dampingFactor,
264                                      DoubleByReference polarity) {
265     OpenMM_AmoebaMultipoleForce_getMultipoleParameters(pointer, index, charge, molecularDipole, molecularQuadrupole,
266         axisType, multipoleAtomZ, multipoleAtomX, multipoleAtomY, thole, dampingFactor, polarity);
267   }
268 
269   /**
270    * Get the mutual induced max iterations.
271    *
272    * @return The mutual induced max iterations.
273    */
274   public int getMutualInducedMaxIterations() {
275     return OpenMM_AmoebaMultipoleForce_getMutualInducedMaxIterations(pointer);
276   }
277 
278   /**
279    * Get the mutual induced target epsilon.
280    *
281    * @return The mutual induced target epsilon.
282    */
283   public double getMutualInducedTargetEpsilon() {
284     return OpenMM_AmoebaMultipoleForce_getMutualInducedTargetEpsilon(pointer);
285   }
286 
287   /**
288    * Get the nonbonded method for the multipole force.
289    *
290    * @return The nonbonded method.
291    */
292   public int getNonbondedMethod() {
293     return OpenMM_AmoebaMultipoleForce_getNonbondedMethod(pointer);
294   }
295 
296   /**
297    * Get the number of multipoles.
298    *
299    * @return The number of multipoles.
300    */
301   public int getNumMultipoles() {
302     return OpenMM_AmoebaMultipoleForce_getNumMultipoles(pointer);
303   }
304 
305   /**
306    * Get the parameters to use for PME calculations.  If alpha is 0 (the default), these parameters are
307    * ignored and instead their values are chosen based on the Ewald error tolerance.
308    *
309    * @param alpha the separation parameter
310    * @param nx    the number of grid points along the X axis
311    * @param ny    the number of grid points along the Y axis
312    * @param nz    the number of grid points along the Z axis
313    */
314   public void getPMEParameters(DoubleByReference alpha, IntByReference nx, IntByReference ny, IntByReference nz) {
315     OpenMM_AmoebaMultipoleForce_getPMEParameters(pointer, alpha, nx, ny, nz);
316   }
317 
318   /**
319    * Get the PME parameters in the context.
320    *
321    * @param context The OpenMM context.
322    * @param alpha   The Ewald alpha parameter (output).
323    * @param nx      The PME grid size in x (output).
324    * @param ny      The PME grid size in y (output).
325    * @param nz      The PME grid size in z (output).
326    */
327   public void getPMEParametersInContext(Context context, DoubleByReference alpha,
328                                         IntByReference nx, IntByReference ny, IntByReference nz) {
329     if (context.hasContextPointer() && pointer != null) {
330       OpenMM_AmoebaMultipoleForce_getPMEParametersInContext(pointer, context.getPointer(), alpha, nx, ny, nz);
331     }
332   }
333 
334   /**
335    * Get the PME grid dimensions.
336    *
337    * @return An IntArray containing the PME grid dimensions.
338    */
339   public IntArray getPmeGridDimensions() {
340     IntArray gridDimensions = new IntArray(0);
341     if (pointer != null) {
342       OpenMM_AmoebaMultipoleForce_getPmeGridDimensions(pointer, gridDimensions.getPointer());
343     }
344     return gridDimensions;
345   }
346 
347   /**
348    * Get the PME B-spline order.
349    *
350    * @return The PME B-spline order.
351    */
352   public int getPmeBSplineOrder() {
353     return OpenMM_AmoebaMultipoleForce_getPmeBSplineOrder(pointer);
354   }
355 
356   /**
357    * Get the polarization type.
358    *
359    * @return The polarization type.
360    */
361   public int getPolarizationType() {
362     return OpenMM_AmoebaMultipoleForce_getPolarizationType(pointer);
363   }
364 
365   /**
366    * Get the system multipole moments.
367    *
368    * @param context The OpenMM context.
369    * @param moments The system multipole moments (output).
370    */
371   public void getSystemMultipoleMoments(Context context, DoubleArray moments) {
372     if (context.hasContextPointer() && pointer != null) {
373       OpenMM_AmoebaMultipoleForce_getSystemMultipoleMoments(pointer, context.getPointer(), moments.getPointer());
374     }
375   }
376 
377   /**
378    * Get the total dipoles.
379    *
380    * @param context The OpenMM context.
381    * @param dipoles The total dipoles (output).
382    */
383   public void getTotalDipoles(Context context, DoubleArray dipoles) {
384     if (context.hasContextPointer() && pointer != null) {
385       OpenMM_AmoebaMultipoleForce_getTotalDipoles(pointer, context.getPointer(), dipoles.getPointer());
386     }
387   }
388 
389   /**
390    * Set the Ewald alpha parameter.  If this is 0 (the default), a value is chosen automatically
391    * based on the Ewald error tolerance.
392    *
393    * @param aewald alpha parameter
394    * @deprecated This method exists only for backward compatibility.  Use setPMEParameters() instead.
395    */
396   @Deprecated
397   public void setAEwald(double aewald) {
398     OpenMM_AmoebaMultipoleForce_setAEwald(pointer, aewald);
399   }
400 
401   /**
402    * Set the covalent map.
403    *
404    * @param i            The atom index.
405    * @param covalentType The covalent type.
406    * @param covalentMap  The covalent map.
407    */
408   public void setCovalentMap(int i, int covalentType, IntArray covalentMap) {
409     OpenMM_AmoebaMultipoleForce_setCovalentMap(pointer, i, covalentType, covalentMap.getPointer());
410   }
411 
412   /**
413    * Set the cutoff distance (in nm) being used for nonbonded interactions.  If the NonbondedMethod in use
414    * is NoCutoff, this value will have no effect.
415    *
416    * @param cutoff the cutoff distance, measured in nm
417    */
418   public void setCutoffDistance(double cutoff) {
419     OpenMM_AmoebaMultipoleForce_setCutoffDistance(pointer, cutoff);
420   }
421 
422   /**
423    * Set the error tolerance for Ewald summation.  This corresponds to the fractional error in the forces
424    * which is acceptable.  This value is used to select the grid dimensions and separation (alpha)
425    * parameter so that the average error level will be less than the tolerance.  There is not a
426    * rigorous guarantee that all forces on all atoms will be less than the tolerance, however.
427    * <p>
428    * This can be overridden by explicitly setting an alpha parameter and grid dimensions to use.
429    *
430    * @param ewaldTolerance the error tolerance
431    */
432   public void setEwaldErrorTolerance(double ewaldTolerance) {
433     OpenMM_AmoebaMultipoleForce_setEwaldErrorTolerance(pointer, ewaldTolerance);
434   }
435 
436   /**
437    * Set extrapolation coefficients.
438    *
439    * @param exptCoefficients The extrapolation coefficients.
440    */
441   public void setExtrapolationCoefficients(DoubleArray exptCoefficients) {
442     OpenMM_AmoebaMultipoleForce_setExtrapolationCoefficients(pointer, exptCoefficients.getPointer());
443   }
444 
445   /**
446    * Set the multipole parameters for a particle.
447    *
448    * @param index               the index of the atom for which to set parameters
449    * @param charge              the particle's charge
450    * @param molecularDipole     the particle's molecular dipole (vector of size 3)
451    * @param molecularQuadrupole the particle's molecular quadrupole (vector of size 9)
452    * @param axisType            the particle's axis type
453    * @param multipoleAtomZ      index of first atom used in constructing lab to molecular frames
454    * @param multipoleAtomX      index of second atom used in constructing lab to molecular frames
455    * @param multipoleAtomY      index of second atom used in constructing lab to molecular frames
456    * @param thole               thole parameter
457    * @param dampingFactor       damping factor parameter
458    * @param polarity            polarity parameter
459    */
460   public void setMultipoleParameters(int index, double charge, DoubleArray molecularDipole, DoubleArray molecularQuadrupole,
461                                      int axisType, int multipoleAtomZ, int multipoleAtomX, int multipoleAtomY,
462                                      double thole, double dampingFactor, double polarity) {
463     OpenMM_AmoebaMultipoleForce_setMultipoleParameters(pointer, index, charge,
464         molecularDipole.getPointer(), molecularQuadrupole.getPointer(), axisType, multipoleAtomZ, multipoleAtomX, multipoleAtomY, thole, dampingFactor, polarity);
465   }
466 
467   /**
468    * Set the mutual induced target maximum number of iterations.
469    *
470    * @param iterations The mutual induced max iterations.
471    */
472   public void setMutualInducedMaxIterations(int iterations) {
473     OpenMM_AmoebaMultipoleForce_setMutualInducedMaxIterations(pointer, iterations);
474   }
475 
476   /**
477    * Set the mutual induced target epsilon.
478    *
479    * @param epsilon The mutual induced target epsilon.
480    */
481   public void setMutualInducedTargetEpsilon(double epsilon) {
482     OpenMM_AmoebaMultipoleForce_setMutualInducedTargetEpsilon(pointer, epsilon);
483   }
484 
485   /**
486    * Set the nonbonded method for the multipole force.
487    *
488    * @param method The nonbonded method.
489    */
490   public void setNonbondedMethod(int method) {
491     OpenMM_AmoebaMultipoleForce_setNonbondedMethod(pointer, method);
492   }
493 
494   /**
495    * Set the parameters to use for PME calculations.  If alpha is 0 (the default), these parameters are
496    * ignored and instead their values are chosen based on the Ewald error tolerance.
497    *
498    * @param alpha the separation parameter
499    * @param nx    the number of grid points along the X axis
500    * @param ny    the number of grid points along the Y axis
501    * @param nz    the number of grid points along the Z axis
502    */
503   public void setPMEParameters(double alpha, int nx, int ny, int nz) {
504     OpenMM_AmoebaMultipoleForce_setPMEParameters(pointer, alpha, nx, ny, nz);
505   }
506 
507   /**
508    * Set the PME grid dimensions for the multipole force.
509    *
510    * @param gridDimensions The PME grid dimensions.
511    */
512   public void setPmeGridDimensions(IntArray gridDimensions) {
513     OpenMM_AmoebaMultipoleForce_setPmeGridDimensions(pointer, gridDimensions.getPointer());
514   }
515 
516   /**
517    * Set the polarization method.
518    *
519    * @param method The polarization method.
520    */
521   public void setPolarizationType(int method) {
522     OpenMM_AmoebaMultipoleForce_setPolarizationType(pointer, method);
523   }
524 
525   /**
526    * Update the parameters in the context.
527    *
528    * @param context The OpenMM context.
529    */
530   public void updateParametersInContext(Context context) {
531     if (context.hasContextPointer()) {
532       OpenMM_AmoebaMultipoleForce_updateParametersInContext(pointer, context.getPointer());
533     }
534   }
535 
536   /**
537    * Check if the force uses periodic boundary conditions.
538    *
539    * @return True if the force uses periodic boundary conditions.
540    */
541   @Override
542   public boolean usesPeriodicBoundaryConditions() {
543     int pbc = OpenMM_AmoebaMultipoleForce_usesPeriodicBoundaryConditions(pointer);
544     return pbc == OpenMM_True;
545   }
546 }