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 jdk.incubator.vector.DoubleVector; 41 42 import static ffx.numerics.multipole.EwaldTensorGlobal.initEwaldSource; 43 import static ffx.numerics.multipole.EwaldTensorGlobalSIMD.fillEwaldSource; 44 45 /** 46 * The EwaldTensorQI class computes derivatives of erfc(<b>r</b>)/|<b>r</b>| via recursion to 47 * arbitrary order for Cartesian multipoles in a quasi-internal frame. 48 * 49 * @author Michael J. Schnieders 50 * @see <a href="http://doi.org/10.1142/9789812830364_0002" target="_blank"> Matt Challacombe, Eric 51 * Schwegler and Jan Almlof, Modern developments in Hartree-Fock theory: Fast methods for 52 * computing the Coulomb matrix. Computational Chemistry: Review of Current Trends. pp. 53-107, 53 * Ed. J. Leczszynski, World Scientifc, 1996. </a> 54 * @since 1.0 55 */ 56 public class EwaldTensorQISIMD extends CoulombTensorQISIMD { 57 58 /** 59 * These are the "source" terms for the recursion for the screened Coulomb operator erfc(R)/R. 60 */ 61 private final double[] ewaldSource; 62 63 /** 64 * The Ewald convergence parameter. 65 */ 66 private final double beta; 67 68 /** 69 * A work array for generation of source terms that cannot be vectorized (exp and erfc). 70 */ 71 private final double[] work; 72 73 /** 74 * Constructor for EwaldTensorQI. 75 * 76 * @param order Tensor order. 77 * @param beta The Ewald convergence parameter. 78 */ 79 public EwaldTensorQISIMD(int order, double beta) { 80 super(order); 81 this.beta = beta; 82 operator = Operator.SCREENED_COULOMB; 83 84 // Auxiliary terms for screened Coulomb (Sagui et al. Eq. 2.28) 85 ewaldSource = new double[o1]; 86 work = new double[o1]; 87 initEwaldSource(order, beta, ewaldSource); 88 } 89 90 /** 91 * Generate source terms for the Ewald Challacombe et al. recursion. 92 * 93 * @param T000 Location to store the source terms. 94 */ 95 @Override 96 protected void source(DoubleVector[] T000) { 97 // Generate source terms for real space Ewald summation. 98 if (beta > 0.0) { 99 fillEwaldSource(order, beta, ewaldSource, R, T000, work); 100 } else { 101 // For beta = 0, generate tensors for the Coulomb operator. 102 super.source(T000); 103 } 104 } 105 106 }