CudaBondedUtilities.cpp 8.11 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
/* -------------------------------------------------------------------------- *
 *                                   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) 2011-2012 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 <http://www.gnu.org/licenses/>.      *
 * -------------------------------------------------------------------------- */

#include "CudaBondedUtilities.h"
#include "CudaExpressionUtilities.h"
29
#include "CudaKernelSources.h"
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
#include "openmm/OpenMMException.h"
#include "CudaNonbondedUtilities.h"
#include <iostream>

using namespace OpenMM;
using namespace std;

CudaBondedUtilities::CudaBondedUtilities(CudaContext& context) : context(context), numForceBuffers(0), maxBonds(0), hasInitializedKernels(false) {
}

CudaBondedUtilities::~CudaBondedUtilities() {
    for (int i = 0; i < (int) atomIndices.size(); i++)
        for (int j = 0; j < (int) atomIndices[i].size(); j++)
            delete atomIndices[i][j];
}

void CudaBondedUtilities::addInteraction(const vector<vector<int> >& atoms, const string& source, int group) {
    if (atoms.size() > 0) {
        forceAtoms.push_back(atoms);
        forceSource.push_back(source);
        forceGroup.push_back(group);
    }
}

std::string CudaBondedUtilities::addArgument(CUdeviceptr data, const string& type) {
    arguments.push_back(data);
    argTypes.push_back(type);
    return "customArg"+context.intToString(arguments.size());
}

void CudaBondedUtilities::addPrefixCode(const string& source) {
61
62
63
    for (int i = 0; i < (int) prefixCode.size(); i++)
        if (prefixCode[i] == source)
            return;
64
65
66
67
68
    prefixCode.push_back(source);
}

void CudaBondedUtilities::initialize(const System& system) {
    int numForces = forceAtoms.size();
69
70
    hasInteractions = (numForces > 0);
    if (!hasInteractions)
71
72
73
74
75
76
77
78
79
80
        return;
    
    // Build the lists of atom indices.
    
    atomIndices.resize(numForces);
    for (int i = 0; i < numForces; i++) {
        int numBonds = forceAtoms[i].size();
        int numAtoms = forceAtoms[i][0].size();
        int startAtom = 0;
        while (startAtom < numAtoms) {
81
            int width = min(numAtoms-startAtom, 4);
82
83
84
85
86
87
88
            if (width == 3)
                width = 4;
            vector<unsigned int> indexVec(width*numBonds);
            for (int bond = 0; bond < numBonds; bond++) {
                for (int atom = 0; atom < width; atom++)
                    indexVec[bond*width+atom] = forceAtoms[i][bond][startAtom+atom];
            }
89
            CudaArray* indices = new CudaArray(context, numBonds, 4*width, "bondedIndices");
90
            indices->upload(&indexVec[0]);
91
92
93
94
95
96
97
98
            atomIndices[i].push_back(indices);
            startAtom += width;
        }
    }

    // Create the kernel.

    stringstream s;
99
    s<<CudaKernelSources::vectorOps;
100
101
    for (int i = 0; i < (int) prefixCode.size(); i++)
        s<<prefixCode[i];
102
    s<<"extern \"C\" __global__ void computeBondedForces(unsigned long long* __restrict__ forceBuffer, real* __restrict__ energyBuffer, const real4* __restrict__ posq, int groups";
103
104
105
    for (int force = 0; force < numForces; force++) {
        for (int i = 0; i < (int) atomIndices[force].size(); i++) {
            int indexWidth = atomIndices[force][i]->getElementSize()/4;
106
            string indexType = "uint"+context.intToString(indexWidth);
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
            s<<", const "<<indexType<<"* __restrict__ atomIndices"<<force<<"_"<<i;
        }
    }
    for (int i = 0; i < (int) arguments.size(); i++)
        s<<", "<<argTypes[i]<<"* customArg"<<(i+1);
    s<<") {\n";
    s<<"real energy = 0;\n";
    for (int force = 0; force < numForces; force++)
        s<<createForceSource(force, forceAtoms[force].size(), forceAtoms[force][0].size(), forceGroup[force], forceSource[force]);
    s<<"energyBuffer[blockIdx.x*blockDim.x+threadIdx.x] += energy;\n";
    s<<"}\n";
    map<string, string> defines;
    defines["PADDED_NUM_ATOMS"] = context.intToString(context.getPaddedNumAtoms());
    CUmodule module = context.createModule(s.str(), defines);
    kernel = context.getKernel(module, "computeBondedForces");
    forceAtoms.clear();
    forceSource.clear();
}

string CudaBondedUtilities::createForceSource(int forceIndex, int numBonds, int numAtoms, int group, const string& computeForce) {
    maxBonds = max(maxBonds, numBonds);
128
    string suffix[] = {".x", ".y", ".z", ".w"};
129
130
131
132
133
134
    stringstream s;
    s<<"if ((groups&"<<(1<<group)<<") != 0)\n";
    s<<"for (unsigned int index = blockIdx.x*blockDim.x+threadIdx.x; index < "<<numBonds<<"; index += blockDim.x*gridDim.x) {\n";
    int startAtom = 0;
    for (int i = 0; i < (int) atomIndices[forceIndex].size(); i++) {
        int indexWidth = atomIndices[forceIndex][i]->getElementSize()/4;
135
        string indexType = "uint"+context.intToString(indexWidth);
136
        s<<"    "<<indexType<<" atoms"<<i<<" = atomIndices"<<forceIndex<<"_"<<i<<"[index];\n";
137
138
        int atomsToLoad = min(indexWidth, numAtoms-startAtom);
        for (int j = 0; j < atomsToLoad; j++) {
139
            s<<"    unsigned int atom"<<(startAtom+j+1)<<" = atoms"<<i<<suffix[j]<<";\n";
140
            s<<"    real4 pos"<<(startAtom+j+1)<<" = posq[atom"<<(startAtom+j+1)<<"];\n";
141
142
143
144
145
        }
        startAtom += indexWidth;
    }
    s<<computeForce<<"\n";
    for (int i = 0; i < numAtoms; i++) {
146
147
148
149
        s<<"    atomicAdd(&forceBuffer[atom"<<(i+1)<<"], static_cast<unsigned long long>((long long) (force"<<(i+1)<<".x*0xFFFFFFFF)));\n";
        s<<"    atomicAdd(&forceBuffer[atom"<<(i+1)<<"+PADDED_NUM_ATOMS], static_cast<unsigned long long>((long long) (force"<<(i+1)<<".y*0xFFFFFFFF)));\n";
        s<<"    atomicAdd(&forceBuffer[atom"<<(i+1)<<"+PADDED_NUM_ATOMS*2], static_cast<unsigned long long>((long long) (force"<<(i+1)<<".z*0xFFFFFFFF)));\n";
        s<<"    __threadfence_block();\n";
150
151
152
153
154
155
    }
    s<<"}\n";
    return s.str();
}

void CudaBondedUtilities::computeInteractions(int groups) {
156
157
158
159
160
161
162
163
164
165
166
167
    if (!hasInitializedKernels) {
        hasInitializedKernels = true;
        kernelArgs.push_back(&context.getForce().getDevicePointer());
        kernelArgs.push_back(&context.getEnergyBuffer().getDevicePointer());
        kernelArgs.push_back(&context.getPosq().getDevicePointer());
        kernelArgs.push_back(NULL);
        for (int i = 0; i < (int) atomIndices.size(); i++)
            for (int j = 0; j < (int) atomIndices[i].size(); j++)
                kernelArgs.push_back(&atomIndices[i][j]->getDevicePointer());
        for (int i = 0; i < (int) arguments.size(); i++)
            kernelArgs.push_back(&arguments[i]);
    }
168
169
    if (!hasInteractions)
        return;
170
171
    kernelArgs[3] = &groups;
    context.executeKernel(kernel, &kernelArgs[0], maxBonds);
172
}