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;
39
40 import com.sun.jna.ptr.DoubleByReference;
41 import com.sun.jna.ptr.IntByReference;
42
43 import java.nio.DoubleBuffer;
44 import java.nio.IntBuffer;
45
46 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_Boolean.OpenMM_True;
47 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_HarmonicAngleForce_addAngle;
48 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_HarmonicAngleForce_create;
49 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_HarmonicAngleForce_destroy;
50 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_HarmonicAngleForce_getAngleParameters;
51 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_HarmonicAngleForce_getNumAngles;
52 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_HarmonicAngleForce_setAngleParameters;
53 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_HarmonicAngleForce_setUsesPeriodicBoundaryConditions;
54 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_HarmonicAngleForce_updateParametersInContext;
55 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_HarmonicAngleForce_usesPeriodicBoundaryConditions;
56
57 /**
58 * This class implements an interaction between groups of three particles that varies harmonically with the angle
59 * between them. To use it, create a HarmonicAngleForce object then call addAngle() once for each angle. After
60 * an angle has been added, you can modify its force field parameters by calling setAngleParameters(). This will
61 * have no effect on Contexts that already exist unless you call updateParametersInContext().
62 */
63 public class HarmonicAngleForce extends Force {
64
65 /**
66 * Create a new HarmonicAngleForce.
67 */
68 public HarmonicAngleForce() {
69 super(OpenMM_HarmonicAngleForce_create());
70 }
71
72 /**
73 * Add an angle term to the force field.
74 *
75 * @param particle1 The index of the first particle forming the angle.
76 * @param particle2 The index of the second particle forming the angle.
77 * @param particle3 The index of the third particle forming the angle.
78 * @param angle The equilibrium angle, measured in radians.
79 * @param k The harmonic force constant for the angle, measured in kJ/mol/radianˆ2.
80 * @return The index of the angle that was added.
81 */
82 public int addAngle(int particle1, int particle2, int particle3, double angle, double k) {
83 return OpenMM_HarmonicAngleForce_addAngle(pointer, particle1, particle2, particle3, angle, k);
84 }
85
86 /**
87 * Destroy the force.
88 */
89 @Override
90 public void destroy() {
91 if (pointer != null) {
92 OpenMM_HarmonicAngleForce_destroy(pointer);
93 pointer = null;
94 }
95 }
96
97 /**
98 * Get the force field parameters for an angle term.
99 *
100 * @param index The index of the angle for which to get parameters.
101 * @param particle1 The index of the first particle forming the angle (output).
102 * @param particle2 The index of the second particle forming the angle (output).
103 * @param particle3 The index of the third particle forming the angle (output).
104 * @param angle The equilibrium angle, measured in radians (output).
105 * @param k The harmonic force constant for the angle, measured in kJ/mol/radianˆ2 (output).
106 */
107 public void getAngleParameters(int index, IntByReference particle1, IntByReference particle2,
108 IntByReference particle3, DoubleByReference angle, DoubleByReference k) {
109 OpenMM_HarmonicAngleForce_getAngleParameters(pointer, index, particle1, particle2, particle3, angle, k);
110 }
111
112 /**
113 * Get the force field parameters for an angle term.
114 *
115 * @param index The index of the angle for which to get parameters.
116 * @param particle1 The index of the first particle forming the angle (output).
117 * @param particle2 The index of the second particle forming the angle (output).
118 * @param particle3 The index of the third particle forming the angle (output).
119 * @param angle The equilibrium angle, measured in radians (output).
120 * @param k The harmonic force constant for the angle, measured in kJ/mol/radianˆ2 (output).
121 */
122 public void getAngleParameters(int index, IntBuffer particle1, IntBuffer particle2,
123 IntBuffer particle3, DoubleBuffer angle, DoubleBuffer k) {
124 OpenMM_HarmonicAngleForce_getAngleParameters(pointer, index, particle1, particle2, particle3, angle, k);
125 }
126
127 /**
128 * Get the number of angles in the force.
129 *
130 * @return The number of angles.
131 */
132 public int getNumAngles() {
133 return OpenMM_HarmonicAngleForce_getNumAngles(pointer);
134 }
135
136 /**
137 * Set the force field parameters for an angle term.
138 *
139 * @param index The index of the angle for which to set parameters.
140 * @param particle1 The index of the first particle forming the angle.
141 * @param particle2 The index of the second particle forming the angle.
142 * @param particle3 The index of the third particle forming the angle.
143 * @param angle The equilibrium angle, measured in radians.
144 * @param k The harmonic force constant for the angle, measured in kJ/mol/radianˆ2.
145 */
146 public void setAngleParameters(int index, int particle1, int particle2, int particle3, double angle, double k) {
147 OpenMM_HarmonicAngleForce_setAngleParameters(pointer, index, particle1, particle2, particle3, angle, k);
148 }
149
150 /**
151 * Set whether this force should apply periodic boundary conditions when calculating displacements.
152 *
153 * @param periodic If true, periodic boundary conditions will be applied.
154 */
155 public void setUsesPeriodicBoundaryConditions(boolean periodic) {
156 OpenMM_HarmonicAngleForce_setUsesPeriodicBoundaryConditions(pointer, periodic ? 1 : 0);
157 }
158
159 /**
160 * Update the parameters in the OpenMM Context.
161 *
162 * @param context The OpenMM Context.
163 */
164 public void updateParametersInContext(Context context) {
165 if (context.hasContextPointer()) {
166 OpenMM_HarmonicAngleForce_updateParametersInContext(pointer, context.getPointer());
167 }
168 }
169
170 /**
171 * Check if the force uses periodic boundary conditions.
172 *
173 * @return True if the force uses periodic boundary conditions.
174 */
175 @Override
176 public boolean usesPeriodicBoundaryConditions() {
177 int pbc = OpenMM_HarmonicAngleForce_usesPeriodicBoundaryConditions(pointer);
178 return pbc == OpenMM_True;
179 }
180 }