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.numerics.multipole; 39 40 import static ffx.numerics.multipole.TholeTensorGlobal.tholeSource; 41 42 /** 43 * The TholeTensorQI class computes derivatives of Thole damping via recursion to order <= 4 for 44 * Cartesian multipoles in a quasi-internal frame. 45 * 46 * @author Michael J. Schnieders 47 * @see <a href="http://doi.org/10.1142/9789812830364_0002" target="_blank"> Matt Challacombe, Eric 48 * Schwegler and Jan Almlof, Modern developments in Hartree-Fock theory: Fast methods for 49 * computing the Coulomb matrix. Computational Chemistry: Review of Current Trends. pp. 53-107, 50 * Ed. J. Leczszynski, World Scientifc, 1996. </a> 51 * @since 1.0 52 */ 53 public class TholeTensorQI extends CoulombTensorQI { 54 55 /** 56 * Thole damping parameter is set to min(pti,ptk)). 57 */ 58 private double thole; 59 60 /** 61 * AiAk parameter = 1/(alphaI^6*alphaK^6) where alpha is polarizability. 62 */ 63 private double AiAk; 64 65 /** 66 * Constructor for TholeTensorQI. 67 * 68 * @param order Tensor order. 69 * @param thole Thole damping parameter is set to min(pti,ptk)). 70 * @param AiAk parameter = 1/(alphaI^6*alphaK^6) where alpha is polarizability. 71 */ 72 public TholeTensorQI(int order, double thole, double AiAk) { 73 super(order); 74 this.thole = thole; 75 this.AiAk = AiAk; 76 this.operator = Operator.THOLE_FIELD; 77 78 // Source terms are currently defined up to order 4. 79 assert (order <= 4); 80 } 81 82 /** 83 * Set Thole damping parameters 84 * 85 * @param thole a double. 86 * @param AiAk a double. 87 */ 88 public void setThole(double thole, double AiAk) { 89 this.thole = thole; 90 this.AiAk = AiAk; 91 } 92 93 /** 94 * Check if the Thole damping is exponential is greater than zero (or the interaction can be 95 * neglected). 96 * 97 * @param r The separation distance. 98 * @return True if -thole*u^3 is greater than -50.0. 99 */ 100 public boolean checkThole(double r) { 101 return TholeTensorGlobal.checkThole(thole, AiAk, r); 102 } 103 104 /** 105 * Generate source terms for the Challacombe et al. recursion. 106 * 107 * @param T000 Location to store the source terms. 108 */ 109 protected void source(double[] T000) { 110 // Compute the normal Coulomb auxiliary term. 111 super.source(T000); 112 113 // Add the Thole damping terms: edamp = exp(-thole*u^3). 114 tholeSource(thole, AiAk, R, false, T000); 115 } 116 117 }