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 static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_LangevinMiddleIntegrator_create; 41 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_LangevinMiddleIntegrator_destroy; 42 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_LangevinMiddleIntegrator_getFriction; 43 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_LangevinMiddleIntegrator_getRandomNumberSeed; 44 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_LangevinMiddleIntegrator_getTemperature; 45 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_LangevinMiddleIntegrator_setFriction; 46 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_LangevinMiddleIntegrator_setRandomNumberSeed; 47 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_LangevinMiddleIntegrator_setTemperature; 48 import static edu.uiowa.jopenmm.OpenMMLibrary.OpenMM_LangevinMiddleIntegrator_step; 49 50 /** 51 * This is an Integrator which simulates a System using Langevin dynamics, with 52 * the LFMiddle discretization (J. Phys. Chem. A 2019, 123, 28, 6056-6079). 53 * This method tend to produce more accurate configurational sampling than other 54 * discretizations, such as the one used in LangevinIntegrator. 55 * <p> 56 * The algorithm is closely related to the BAOAB discretization 57 * (Proc. R. Soc. A. 472: 20160138). Both methods produce identical trajectories, 58 * but LFMiddle returns half step (leapfrog) velocities, while BAOAB returns 59 * on-step velocities. The former provide a much more accurate sampling of the 60 * thermal ensemble. 61 */ 62 public class LangevinMiddleIntegrator extends Integrator { 63 64 /** 65 * Create a LangevinMiddleIntegrator. 66 * 67 * @param dt the step size with which to integrate the system (in picoseconds) 68 * @param temp the temperature of the heat bath (in Kelvin) 69 * @param gamma the friction coefficient which couples the system to the heat bath (in inverse picoseconds) 70 */ 71 public LangevinMiddleIntegrator(double dt, double temp, double gamma) { 72 super(OpenMM_LangevinMiddleIntegrator_create(temp, gamma, dt)); 73 } 74 75 /** 76 * Destroy the integrator. 77 */ 78 @Override 79 public void destroy() { 80 if (pointer != null) { 81 OpenMM_LangevinMiddleIntegrator_destroy(pointer); 82 pointer = null; 83 } 84 } 85 86 /** 87 * Get the friction coefficient. 88 * 89 * @return The friction coefficient in inverse picoseconds. 90 */ 91 public double getFriction() { 92 return OpenMM_LangevinMiddleIntegrator_getFriction(pointer); 93 } 94 95 /** 96 * Get the random number seed. 97 * 98 * @return The random number seed. 99 */ 100 public int getRandomNumberSeed() { 101 return OpenMM_LangevinMiddleIntegrator_getRandomNumberSeed(pointer); 102 } 103 104 /** 105 * Get the temperature. 106 * 107 * @return The temperature in Kelvin. 108 */ 109 public double getTemperature() { 110 return OpenMM_LangevinMiddleIntegrator_getTemperature(pointer); 111 } 112 113 /** 114 * Set the friction coefficient. 115 * 116 * @param gamma The friction coefficient in inverse picoseconds. 117 */ 118 public void setFriction(double gamma) { 119 OpenMM_LangevinMiddleIntegrator_setFriction(pointer, gamma); 120 } 121 122 /** 123 * Set the random number seed. 124 * 125 * @param seed The random number seed. 126 */ 127 public void setRandomNumberSeed(int seed) { 128 OpenMM_LangevinMiddleIntegrator_setRandomNumberSeed(pointer, seed); 129 } 130 131 /** 132 * Set the temperature. 133 * 134 * @param temp The temperature in Kelvin. 135 */ 136 public void setTemperature(double temp) { 137 OpenMM_LangevinMiddleIntegrator_setTemperature(pointer, temp); 138 } 139 140 /** 141 * Step the integrator. 142 * 143 * @param steps The number of steps to take. 144 */ 145 @Override 146 public void step(int steps) { 147 OpenMM_LangevinMiddleIntegrator_step(pointer, steps); 148 } 149 }