HipPlatform.cpp 16.6 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
 *                                                                            *
Peter Eastman's avatar
Peter Eastman committed
7
 * Portions copyright (c) 2008-2025 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
72
73
74
75
76
77
78
79
80
81
82
83
84
85
#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();
    HipKernelFactory* factory = new HipKernelFactory();
    registerKernelFactory(CalcForcesAndEnergyKernel::Name(), factory);
    registerKernelFactory(UpdateStateDataKernel::Name(), factory);
    registerKernelFactory(ApplyConstraintsKernel::Name(), factory);
    registerKernelFactory(VirtualSitesKernel::Name(), factory);
    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);
86
    registerKernelFactory(CalcConstantPotentialForceKernel::Name(), factory);
87
88
89
90
91
92
93
    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);
94
    registerKernelFactory(CalcCustomCPPForceKernel::Name(), factory);
95
    registerKernelFactory(CalcCustomCVForceKernel::Name(), factory);
96
    registerKernelFactory(CalcOrientationRestraintForceKernel::Name(), factory);
Peter Eastman's avatar
Peter Eastman committed
97
    registerKernelFactory(CalcPythonForceKernel::Name(), factory);
98
    registerKernelFactory(CalcRGForceKernel::Name(), factory);
99
100
101
102
103
104
105
106
107
108
    registerKernelFactory(CalcRMSDForceKernel::Name(), factory);
    registerKernelFactory(CalcCustomManyParticleForceKernel::Name(), factory);
    registerKernelFactory(CalcGayBerneForceKernel::Name(), factory);
    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
109
    registerKernelFactory(IntegrateDPDStepKernel::Name(), factory);
110
    registerKernelFactory(IntegrateQTBStepKernel::Name(), factory);
111
112
113
    registerKernelFactory(ApplyAndersenThermostatKernel::Name(), factory);
    registerKernelFactory(ApplyMonteCarloBarostatKernel::Name(), factory);
    registerKernelFactory(RemoveCMMotionKernel::Name(), factory);
Anton Gorenko's avatar
Anton Gorenko committed
114
    registerKernelFactory(CalcATMForceKernel::Name(), factory);
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
    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());
    setPropertyDefaultValue(HipDeviceIndex(), "");
    setPropertyDefaultValue(HipDeviceName(), "");
    setPropertyDefaultValue(HipUseBlockingSync(), "true");
    setPropertyDefaultValue(HipPrecision(), "single");
    setPropertyDefaultValue(HipUseCpuPme(), "false");
    setPropertyDefaultValue(HipDisablePmeStream(), "false");
    setPropertyDefaultValue(HipDeterministicForces(), "false");
130
131
132
#ifdef _MSC_VER
    setPropertyDefaultValue(HipTempDirectory(), string(getenv("TEMP")));
#else
133
134
135
    char* tmpdir = getenv("TMPDIR");
    string tmp = (tmpdir == NULL ? string(P_tmpdir) : string(tmpdir));
    setPropertyDefaultValue(HipTempDirectory(), tmp);
136
#endif
137
138
139
}

double HipPlatform::getSpeed() const {
140
141
142
143
    // 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;
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
}

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 {
}

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);
    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);
    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;
193
194
    context.setPlatformData(new PlatformData(&context, context.getSystem(), devicePropValue, blockingPropValue, precisionPropValue, cpuPmePropValue, tempPropValue,
            pmeStreamPropValue, deterministicForcesValue, threads, NULL));
195
196
197
198
199
200
201
202
203
204
205
206
}

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());
    int threads = reinterpret_cast<PlatformData*>(originalContext.getPlatformData())->threads.getNumThreads();
207
208
    context.setPlatformData(new PlatformData(&context, context.getSystem(), devicePropValue, blockingPropValue, precisionPropValue, cpuPmePropValue, tempPropValue,
            pmeStreamPropValue, deterministicForcesValue, threads, &originalContext));
209
210
211
212
213
214
215
216
}

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,
217
            const string& cpuPmeProperty, const string& tempProperty, const string& pmeStreamProperty,
218
            const string& deterministicForcesProperty, int numThreads, ContextImpl* originalContext) :
219
220
                context(context), removeCM(false), stepCount(0), computeForceCount(0), time(0.0), hasInitializedContexts(false),
                threads(numThreads) {
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
    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;
237
                contexts.push_back(new HipContext(system, deviceIndex, blocking, precisionProperty, tempProperty, *this, (originalData == NULL ? NULL : originalData->contexts[i])));
238
239
240
            }
        }
        if (contexts.size() == 0)
241
            contexts.push_back(new HipContext(system, -1, blocking, precisionProperty, tempProperty, *this, (originalData == NULL ? NULL : originalData->contexts[0])));
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
    }
    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";
    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();
}