#ifndef OPENMM_OPENCLINTEGRATIONUTILITIES_H_ #define OPENMM_OPENCLINTEGRATIONUTILITIES_H_ /* -------------------------------------------------------------------------- * * OpenMM * * -------------------------------------------------------------------------- * * This is part of the OpenMM molecular simulation toolkit originating from * * Simbios, the NIH National Center for Physics-Based Simulation of * * Biological Structures at Stanford, funded under the NIH Roadmap for * * Medical Research, grant U54 GM072970. See https://simtk.org. * * * * Portions copyright (c) 2009 Stanford University and the Authors. * * Authors: Peter Eastman * * Contributors: * * * * This program is free software: you can redistribute it and/or modify * * it under the terms of the GNU Lesser General Public License as published * * by the Free Software Foundation, either version 3 of the License, or * * (at your option) any later version. * * * * This program is distributed in the hope that it will be useful, * * but WITHOUT ANY WARRANTY; without even the implied warranty of * * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * * GNU Lesser General Public License for more details. * * * * You should have received a copy of the GNU Lesser General Public License * * along with this program. If not, see . * * -------------------------------------------------------------------------- */ #include "openmm/System.h" #include "OpenCLContext.h" #include "openmm/internal/windowsExport.h" namespace OpenMM { /** * This class implements features that are used by many different integrators, including * common workspace arrays, random number generation, and enforcing constraints. */ class OPENMM_EXPORT OpenCLIntegrationUtilities { public: OpenCLIntegrationUtilities(OpenCLContext& context, const System& system); ~OpenCLIntegrationUtilities(); /** * Get the array which contains position deltas. */ OpenCLArray& getPosDelta() { return *posDelta; } /** * Get the array which contains random values. Each element is a float4, whose components * are independent, normally distributed random numbers with mean 0 and variance 1. */ OpenCLArray& getRandom() { return *random; } /** * Get the array which contains the current step size. */ OpenCLArray& getStepSize() { return *stepSize; } /** * Apply constraints to the atom positions. * * @param tol the constraint tolerance */ void applyConstraints(double tol); /** * Initialize the random number generator. */ void initRandomNumberGenerator(unsigned int randomNumberSeed); /** * Ensure that sufficient random numbers are available in the array, and generate new ones if not. * * @param numValues the number of random float4's that will be required * @return the index in the array at which to start reading */ int prepareRandomNumbers(int numValues); private: OpenCLContext& context; cl::Kernel settleKernel; cl::Kernel shakeKernel; cl::Kernel ccmaDirectionsKernel; cl::Kernel ccmaForceKernel; cl::Kernel ccmaMultiplyKernel; cl::Kernel ccmaUpdateKernel; cl::Kernel randomKernel; OpenCLArray* posDelta; OpenCLArray* settleAtoms; OpenCLArray* settleParams; OpenCLArray* shakeAtoms; OpenCLArray* shakeParams; OpenCLArray* random; OpenCLArray* randomSeed; OpenCLArray* stepSize; OpenCLArray* ccmaAtoms; OpenCLArray* ccmaDistance; OpenCLArray* ccmaReducedMass; OpenCLArray* ccmaAtomConstraints; OpenCLArray* ccmaNumAtomConstraints; OpenCLArray* ccmaConstraintMatrixColumn; OpenCLArray* ccmaConstraintMatrixValue; OpenCLArray* ccmaDelta1; OpenCLArray* ccmaDelta2; OpenCLArray* ccmaConverged; cl::Buffer* ccmaConvergedBuffer; cl_int* ccmaConvergedMemory; int randomPos; int lastSeed; bool hasInitializedConstraintKernels; struct ShakeCluster; struct ConstraintOrderer; }; } // namespace OpenMM #endif /*OPENMM_OPENCLINTEGRATIONUTILITIES_H_*/