HipPlatform.cpp 18.8 KB
Newer Older
1
2
3
/* -------------------------------------------------------------------------- *
 *                                   OpenMM                                   *
 * -------------------------------------------------------------------------- *
Evan Pretti's avatar
Evan Pretti committed
4
5
 * This is part of the OpenMM molecular simulation toolkit.                   *
 * See https://openmm.org/development.                                        *
6
 *                                                                            *
7
 * Portions copyright (c) 2008-2026 Stanford University and the Authors.      *
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
 * Portions copyright (c) 2020 Advanced Micro Devices, Inc.                   *
 * Authors: Peter Eastman, Nicholas Curtis                                    *
 * 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 "HipContext.h"
#include "HipExpressionUtilities.h"
#include "HipPlatform.h"
#include "HipKernelFactory.h"
#include "HipKernels.h"
#include "openmm/Context.h"
#include "openmm/System.h"
#include "openmm/internal/ContextImpl.h"
#include "openmm/internal/hardware.h"
#include <algorithm>
#include <cctype>
#include <sstream>
#include <cstdio>
#ifdef _MSC_VER
40
    #include <Windows.h>
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
#endif
using namespace OpenMM;
using namespace std;

#define CHECK_RESULT(result, prefix) \
    if (result != hipSuccess) { \
        std::stringstream m; \
        m<<prefix<<": "<<HipContext::getErrorString(result)<<" ("<<result<<")"<<" at "<<__FILE__<<":"<<__LINE__; \
        throw OpenMMException(m.str());\
    }


#ifdef OPENMM_COMMON_BUILDING_STATIC_LIBRARY
extern "C" void registerHipPlatform() {
    Platform::registerPlatform(new HipPlatform());
}
#else
extern "C" OPENMM_EXPORT_COMMON void registerPlatforms() {
    Platform::registerPlatform(new HipPlatform());
}
#endif

HipPlatform::HipPlatform() {
    deprecatedPropertyReplacements["HipDeviceIndex"] = HipDeviceIndex();
    deprecatedPropertyReplacements["HipDeviceName"] = HipDeviceName();
    deprecatedPropertyReplacements["HipUseBlockingSync"] = HipUseBlockingSync();
    deprecatedPropertyReplacements["HipPrecision"] = HipPrecision();
    deprecatedPropertyReplacements["HipUseCpuPme"] = HipUseCpuPme();
    deprecatedPropertyReplacements["HipTempDirectory"] = HipTempDirectory();
    deprecatedPropertyReplacements["HipDisablePmeStream"] = HipDisablePmeStream();
    deprecatedPropertyReplacements["HipDeterministicForces"] = HipDeterministicForces();
one's avatar
one committed
72
    deprecatedPropertyReplacements["HipFFTBackend"] = HipFFTBackend();
73
74
75
76
77
    HipKernelFactory* factory = new HipKernelFactory();
    registerKernelFactory(CalcForcesAndEnergyKernel::Name(), factory);
    registerKernelFactory(UpdateStateDataKernel::Name(), factory);
    registerKernelFactory(ApplyConstraintsKernel::Name(), factory);
    registerKernelFactory(VirtualSitesKernel::Name(), factory);
78
    registerKernelFactory(MinimizeKernel::Name(), factory);
79
80
81
82
83
84
85
86
87
    registerKernelFactory(CalcHarmonicBondForceKernel::Name(), factory);
    registerKernelFactory(CalcCustomBondForceKernel::Name(), factory);
    registerKernelFactory(CalcHarmonicAngleForceKernel::Name(), factory);
    registerKernelFactory(CalcCustomAngleForceKernel::Name(), factory);
    registerKernelFactory(CalcPeriodicTorsionForceKernel::Name(), factory);
    registerKernelFactory(CalcRBTorsionForceKernel::Name(), factory);
    registerKernelFactory(CalcCMAPTorsionForceKernel::Name(), factory);
    registerKernelFactory(CalcCustomTorsionForceKernel::Name(), factory);
    registerKernelFactory(CalcNonbondedForceKernel::Name(), factory);
88
    registerKernelFactory(CalcConstantPotentialForceKernel::Name(), factory);
89
90
91
92
93
94
95
    registerKernelFactory(CalcCustomNonbondedForceKernel::Name(), factory);
    registerKernelFactory(CalcGBSAOBCForceKernel::Name(), factory);
    registerKernelFactory(CalcCustomGBForceKernel::Name(), factory);
    registerKernelFactory(CalcCustomExternalForceKernel::Name(), factory);
    registerKernelFactory(CalcCustomHbondForceKernel::Name(), factory);
    registerKernelFactory(CalcCustomCentroidBondForceKernel::Name(), factory);
    registerKernelFactory(CalcCustomCompoundBondForceKernel::Name(), factory);
96
    registerKernelFactory(CalcCustomCPPForceKernel::Name(), factory);
97
    registerKernelFactory(CalcCustomCVForceKernel::Name(), factory);
98
    registerKernelFactory(CalcOrientationRestraintForceKernel::Name(), factory);
Peter Eastman's avatar
Peter Eastman committed
99
    registerKernelFactory(CalcPythonForceKernel::Name(), factory);
100
    registerKernelFactory(CalcRGForceKernel::Name(), factory);
101
102
103
    registerKernelFactory(CalcRMSDForceKernel::Name(), factory);
    registerKernelFactory(CalcCustomManyParticleForceKernel::Name(), factory);
    registerKernelFactory(CalcGayBerneForceKernel::Name(), factory);
Evan Pretti's avatar
Evan Pretti committed
104
    registerKernelFactory(CalcLCPOForceKernel::Name(), factory);
105
106
107
108
109
110
111
    registerKernelFactory(IntegrateVerletStepKernel::Name(), factory);
    registerKernelFactory(IntegrateNoseHooverStepKernel::Name(), factory);
    registerKernelFactory(IntegrateLangevinMiddleStepKernel::Name(), factory);
    registerKernelFactory(IntegrateBrownianStepKernel::Name(), factory);
    registerKernelFactory(IntegrateVariableVerletStepKernel::Name(), factory);
    registerKernelFactory(IntegrateVariableLangevinStepKernel::Name(), factory);
    registerKernelFactory(IntegrateCustomStepKernel::Name(), factory);
Peter Eastman's avatar
Peter Eastman committed
112
    registerKernelFactory(IntegrateDPDStepKernel::Name(), factory);
113
    registerKernelFactory(IntegrateQTBStepKernel::Name(), factory);
114
115
116
    registerKernelFactory(ApplyAndersenThermostatKernel::Name(), factory);
    registerKernelFactory(ApplyMonteCarloBarostatKernel::Name(), factory);
    registerKernelFactory(RemoveCMMotionKernel::Name(), factory);
Anton Gorenko's avatar
Anton Gorenko committed
117
    registerKernelFactory(CalcATMForceKernel::Name(), factory);
118
119
120
121
122
123
124
125
    platformProperties.push_back(HipDeviceIndex());
    platformProperties.push_back(HipDeviceName());
    platformProperties.push_back(HipUseBlockingSync());
    platformProperties.push_back(HipPrecision());
    platformProperties.push_back(HipUseCpuPme());
    platformProperties.push_back(HipTempDirectory());
    platformProperties.push_back(HipDisablePmeStream());
    platformProperties.push_back(HipDeterministicForces());
one's avatar
one committed
126
    platformProperties.push_back(HipFFTBackend());
127
128
129
130
131
132
133
    setPropertyDefaultValue(HipDeviceIndex(), "");
    setPropertyDefaultValue(HipDeviceName(), "");
    setPropertyDefaultValue(HipUseBlockingSync(), "true");
    setPropertyDefaultValue(HipPrecision(), "single");
    setPropertyDefaultValue(HipUseCpuPme(), "false");
    setPropertyDefaultValue(HipDisablePmeStream(), "false");
    setPropertyDefaultValue(HipDeterministicForces(), "false");
one's avatar
one committed
134
    setPropertyDefaultValue(HipFFTBackend(), "vkfft");
135
136
137
#ifdef _MSC_VER
    setPropertyDefaultValue(HipTempDirectory(), string(getenv("TEMP")));
#else
138
139
140
    char* tmpdir = getenv("TMPDIR");
    string tmp = (tmpdir == NULL ? string(P_tmpdir) : string(tmpdir));
    setPropertyDefaultValue(HipTempDirectory(), tmp);
141
#endif
142
143
144
}

double HipPlatform::getSpeed() const {
145
146
147
148
    // Reduce the speed of the HIP platform if there are no HIP devices in the system,
    // so the OpenCL plaform can be selected as default
    int numDevices;
    return hipGetDeviceCount(&numDevices) != hipErrorNoDevice ? 100 : 40;
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
}

bool HipPlatform::supportsDoublePrecision() const {
    return true;
}

const string& HipPlatform::getPropertyValue(const Context& context, const string& property) const {
    const ContextImpl& impl = getContextImpl(context);
    const PlatformData* data = reinterpret_cast<const PlatformData*>(impl.getPlatformData());
    string propertyName = property;
    if (deprecatedPropertyReplacements.find(property) != deprecatedPropertyReplacements.end())
        propertyName = deprecatedPropertyReplacements.find(property)->second;
    map<string, string>::const_iterator value = data->propertyValues.find(propertyName);
    if (value != data->propertyValues.end())
        return value->second;
    return Platform::getPropertyValue(context, property);
}

void HipPlatform::setPropertyValue(Context& context, const string& property, const string& value) const {
}

170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
vector<map<string, string> > HipPlatform::getDevices(const map<string, string>& filters) const {
    // Check for properties that might act as filters.

    int deviceIndex = -1;
    if (filters.find(HipDeviceIndex()) != filters.end())
        stringstream(filters.at(HipDeviceIndex())) >> deviceIndex;
    string deviceName = (filters.find(HipDeviceName()) == filters.end() ? "" : filters.at(HipDeviceName()));

    // Loop over devices.

    vector<map<string, string> > results;
    int numDevices;
    if (hipGetDeviceCount(&numDevices) != hipSuccess)
        numDevices = 0;
    for (int i = 0; i < numDevices; i++) {
        if (deviceIndex != -1 && deviceIndex != i)
            continue;
        char name[1000];
        hipDevice_t device;
        CHECK_RESULT(hipDeviceGet(&device, i), "Error querying device");
        CHECK_RESULT(hipDeviceGetName(name, 1000, device), "Error querying device name");
        stringstream deviceNameStr;
        deviceNameStr << name;
        if (deviceName.size() > 0 && deviceName != deviceNameStr.str())
            continue;
        stringstream deviceIndexStr;
        deviceIndexStr << i;
        map<string, string> properties = {{HipDeviceIndex(), deviceIndexStr.str()},
                                          {HipDeviceName(), deviceNameStr.str()}};
        results.push_back(properties);
    }
    return results;
}

204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
void HipPlatform::contextCreated(ContextImpl& context, const map<string, string>& properties) const {
    const string& devicePropValue = (properties.find(HipDeviceIndex()) == properties.end() ?
            getPropertyDefaultValue(HipDeviceIndex()) : properties.find(HipDeviceIndex())->second);
    string blockingPropValue = (properties.find(HipUseBlockingSync()) == properties.end() ?
            getPropertyDefaultValue(HipUseBlockingSync()) : properties.find(HipUseBlockingSync())->second);
    string precisionPropValue = (properties.find(HipPrecision()) == properties.end() ?
            getPropertyDefaultValue(HipPrecision()) : properties.find(HipPrecision())->second);
    string cpuPmePropValue = (properties.find(HipUseCpuPme()) == properties.end() ?
            getPropertyDefaultValue(HipUseCpuPme()) : properties.find(HipUseCpuPme())->second);
    const string& tempPropValue = (properties.find(HipTempDirectory()) == properties.end() ?
            getPropertyDefaultValue(HipTempDirectory()) : properties.find(HipTempDirectory())->second);
    string pmeStreamPropValue = (properties.find(HipDisablePmeStream()) == properties.end() ?
            getPropertyDefaultValue(HipDisablePmeStream()) : properties.find(HipDisablePmeStream())->second);
    string deterministicForcesValue = (properties.find(HipDeterministicForces()) == properties.end() ?
            getPropertyDefaultValue(HipDeterministicForces()) : properties.find(HipDeterministicForces())->second);
one's avatar
one committed
219
220
    string fftBackendValue = (properties.find(HipFFTBackend()) == properties.end() ?
            getPropertyDefaultValue(HipFFTBackend()) : properties.find(HipFFTBackend())->second);
221
222
223
224
225
    transform(blockingPropValue.begin(), blockingPropValue.end(), blockingPropValue.begin(), ::tolower);
    transform(precisionPropValue.begin(), precisionPropValue.end(), precisionPropValue.begin(), ::tolower);
    transform(cpuPmePropValue.begin(), cpuPmePropValue.end(), cpuPmePropValue.begin(), ::tolower);
    transform(pmeStreamPropValue.begin(), pmeStreamPropValue.end(), pmeStreamPropValue.begin(), ::tolower);
    transform(deterministicForcesValue.begin(), deterministicForcesValue.end(), deterministicForcesValue.begin(), ::tolower);
one's avatar
one committed
226
    transform(fftBackendValue.begin(), fftBackendValue.end(), fftBackendValue.begin(), ::tolower);
227
228
229
230
231
232
233
234
    vector<string> pmeKernelName;
    pmeKernelName.push_back(CalcPmeReciprocalForceKernel::Name());
    if (!supportsKernels(pmeKernelName))
        cpuPmePropValue = "false";
    int threads = getNumProcessors();
    char* threadsEnv = getenv("OPENMM_CPU_THREADS");
    if (threadsEnv != NULL)
        stringstream(threadsEnv) >> threads;
235
    context.setPlatformData(new PlatformData(&context, context.getSystem(), devicePropValue, blockingPropValue, precisionPropValue, cpuPmePropValue, tempPropValue,
one's avatar
one committed
236
            pmeStreamPropValue, deterministicForcesValue, fftBackendValue, threads, NULL));
237
238
239
240
241
242
243
244
245
246
247
}

void HipPlatform::linkedContextCreated(ContextImpl& context, ContextImpl& originalContext) const {
    Platform& platform = originalContext.getPlatform();
    string devicePropValue = platform.getPropertyValue(originalContext.getOwner(), HipDeviceIndex());
    string blockingPropValue = platform.getPropertyValue(originalContext.getOwner(), HipUseBlockingSync());
    string precisionPropValue = platform.getPropertyValue(originalContext.getOwner(), HipPrecision());
    string cpuPmePropValue = platform.getPropertyValue(originalContext.getOwner(), HipUseCpuPme());
    string tempPropValue = platform.getPropertyValue(originalContext.getOwner(), HipTempDirectory());
    string pmeStreamPropValue = platform.getPropertyValue(originalContext.getOwner(), HipDisablePmeStream());
    string deterministicForcesValue = platform.getPropertyValue(originalContext.getOwner(), HipDeterministicForces());
one's avatar
one committed
248
    string fftBackendValue = platform.getPropertyValue(originalContext.getOwner(), HipFFTBackend());
249
    int threads = reinterpret_cast<PlatformData*>(originalContext.getPlatformData())->threads.getNumThreads();
250
    context.setPlatformData(new PlatformData(&context, context.getSystem(), devicePropValue, blockingPropValue, precisionPropValue, cpuPmePropValue, tempPropValue,
one's avatar
one committed
251
            pmeStreamPropValue, deterministicForcesValue, fftBackendValue, threads, &originalContext));
252
253
254
255
256
257
258
259
}

void HipPlatform::contextDestroyed(ContextImpl& context) const {
    PlatformData* data = reinterpret_cast<PlatformData*>(context.getPlatformData());
    delete data;
}

HipPlatform::PlatformData::PlatformData(ContextImpl* context, const System& system, const string& deviceIndexProperty, const string& blockingProperty, const string& precisionProperty,
one's avatar
one committed
260
261
              const string& cpuPmeProperty, const string& tempProperty, const string& pmeStreamProperty,
              const string& deterministicForcesProperty, const string& fftBackendProperty, int numThreads, ContextImpl* originalContext) :
262
263
                context(context), removeCM(false), stepCount(0), computeForceCount(0), time(0.0), hasInitializedContexts(false),
                threads(numThreads) {
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
    bool blocking = (blockingProperty == "true");
    vector<string> devices;
    size_t searchPos = 0, nextPos;
    while ((nextPos = deviceIndexProperty.find_first_of(", ", searchPos)) != string::npos) {
        devices.push_back(deviceIndexProperty.substr(searchPos, nextPos-searchPos));
        searchPos = nextPos+1;
    }
    devices.push_back(deviceIndexProperty.substr(searchPos));
    PlatformData* originalData = NULL;
    if (originalContext != NULL)
        originalData = reinterpret_cast<PlatformData*>(originalContext->getPlatformData());
    try {
        for (int i = 0; i < (int) devices.size(); i++) {
            if (devices[i].length() > 0) {
                int deviceIndex;
                stringstream(devices[i]) >> deviceIndex;
280
                contexts.push_back(new HipContext(system, deviceIndex, blocking, precisionProperty, tempProperty, *this, (originalData == NULL ? NULL : originalData->contexts[i])));
281
282
283
            }
        }
        if (contexts.size() == 0)
284
            contexts.push_back(new HipContext(system, -1, blocking, precisionProperty, tempProperty, *this, (originalData == NULL ? NULL : originalData->contexts[0])));
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
    }
    catch (...) {
        // If an exception was thrown, do our best to clean up memory.

        for (int i = 0; i < (int) contexts.size(); i++)
            delete contexts[i];
        throw;
    }
    stringstream deviceIndex, deviceName;
    for (int i = 0; i < (int) contexts.size(); i++) {
        if (i > 0) {
            deviceIndex << ',';
            deviceName << ',';
        }
        deviceIndex << contexts[i]->getDeviceIndex();
        char name[1000];
        CHECK_RESULT(hipDeviceGetName(name, 1000, contexts[i]->getDevice()), "Error querying device name");
        deviceName << name;
    }

    useCpuPme = (cpuPmeProperty == "true" && !contexts[0]->getUseDoublePrecision());
    disablePmeStream = (pmeStreamProperty == "true");
    deterministicForces = (deterministicForcesProperty == "true");
    propertyValues[HipPlatform::HipDeviceIndex()] = deviceIndex.str();
    propertyValues[HipPlatform::HipDeviceName()] = deviceName.str();
    propertyValues[HipPlatform::HipUseBlockingSync()] = blocking ? "true" : "false";
    propertyValues[HipPlatform::HipPrecision()] = precisionProperty;
    propertyValues[HipPlatform::HipUseCpuPme()] = useCpuPme ? "true" : "false";
    propertyValues[HipPlatform::HipTempDirectory()] = tempProperty;
    propertyValues[HipPlatform::HipDisablePmeStream()] = disablePmeStream ? "true" : "false";
    propertyValues[HipPlatform::HipDeterministicForces()] = deterministicForces ? "true" : "false";
one's avatar
one committed
316
    propertyValues[HipPlatform::HipFFTBackend()] = fftBackendProperty;
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
    contextEnergy.resize(contexts.size());

    // Determine whether peer-to-peer copying is supported, and enable it if so.

    peerAccessSupported = true;
    for (int i = 1; i < contexts.size(); i++) {
        int canAccess;
        hipDeviceCanAccessPeer(&canAccess, contexts[i]->getDevice(), contexts[0]->getDevice());
        if (!canAccess) {
            peerAccessSupported = false;
            break;
        }
    }
}

HipPlatform::PlatformData::~PlatformData() {
    for (int i = 0; i < (int) contexts.size(); i++)
        delete contexts[i];
}

void HipPlatform::PlatformData::initializeContexts(const System& system) {
    if (hasInitializedContexts)
        return;
    for (int i = 0; i < (int) contexts.size(); i++)
        contexts[i]->initialize();
    hasInitializedContexts = true;
}

void HipPlatform::PlatformData::syncContexts() {
    for (int i = 0; i < (int) contexts.size(); i++)
        contexts[i]->getWorkThread().flush();
}