CudaPlatform.cpp 7.16 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
/* -------------------------------------------------------------------------- *
 *                                   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) 2008 Stanford University and the Authors.           *
 * Authors: Peter Eastman                                                     *
 * Contributors:                                                              *
 *                                                                            *
13
14
15
16
 * 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.                                        *
17
 *                                                                            *
18
19
20
21
 * 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.                        *
22
 *                                                                            *
23
24
 * 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/>.      *
25
26
27
28
29
 * -------------------------------------------------------------------------- */

#include "CudaPlatform.h"
#include "CudaKernelFactory.h"
#include "CudaKernels.h"
30
#include "openmm/internal/ContextImpl.h"
31
#include "kernels/gputypes.h"
32
#include "openmm/Context.h"
33
#include "openmm/OpenMMException.h"
34
#include "openmm/System.h"
35
#include <sstream>
36

37
using namespace OpenMM;
38
39
40
using std::map;
using std::string;
using std::stringstream;
41

42
extern "C" OPENMMCUDA_EXPORT void registerPlatforms() {
Peter Eastman's avatar
Peter Eastman committed
43
44
    if (gpuIsAvailable())
        Platform::registerPlatform(new CudaPlatform());
45
46
}

47
48
CudaPlatform::CudaPlatform() {
    CudaKernelFactory* factory = new CudaKernelFactory();
49
    registerKernelFactory(CalcForcesAndEnergyKernel::Name(), factory);
50
    registerKernelFactory(UpdateStateDataKernel::Name(), factory);
51
    registerKernelFactory(ApplyConstraintsKernel::Name(), factory);
52
    registerKernelFactory(VirtualSitesKernel::Name(), factory);
53
    registerKernelFactory(CalcHarmonicBondForceKernel::Name(), factory);
54
    registerKernelFactory(CalcCustomBondForceKernel::Name(), factory);
55
    registerKernelFactory(CalcHarmonicAngleForceKernel::Name(), factory);
56
    registerKernelFactory(CalcCustomAngleForceKernel::Name(), factory);
57
58
    registerKernelFactory(CalcPeriodicTorsionForceKernel::Name(), factory);
    registerKernelFactory(CalcRBTorsionForceKernel::Name(), factory);
59
    registerKernelFactory(CalcCMAPTorsionForceKernel::Name(), factory);
60
    registerKernelFactory(CalcCustomTorsionForceKernel::Name(), factory);
61
    registerKernelFactory(CalcNonbondedForceKernel::Name(), factory);
62
    registerKernelFactory(CalcCustomNonbondedForceKernel::Name(), factory);
63
    registerKernelFactory(CalcGBSAOBCForceKernel::Name(), factory);
64
65
    registerKernelFactory(CalcGBVIForceKernel::Name(), factory);
    registerKernelFactory(CalcCustomExternalForceKernel::Name(), factory);
66
    registerKernelFactory(IntegrateVerletStepKernel::Name(), factory);
67
    registerKernelFactory(IntegrateLangevinStepKernel::Name(), factory);
68
    registerKernelFactory(IntegrateBrownianStepKernel::Name(), factory);
69
    registerKernelFactory(IntegrateVariableVerletStepKernel::Name(), factory);
70
    registerKernelFactory(IntegrateVariableLangevinStepKernel::Name(), factory);
71
    registerKernelFactory(ApplyAndersenThermostatKernel::Name(), factory);
72
    registerKernelFactory(ApplyMonteCarloBarostatKernel::Name(), factory);
73
    registerKernelFactory(CalcKineticEnergyKernel::Name(), factory);
74
    registerKernelFactory(RemoveCMMotionKernel::Name(), factory);
75
    platformProperties.push_back(CudaDevice());
76
    platformProperties.push_back(CudaUseBlockingSync());
77
    setPropertyDefaultValue(CudaDevice(), "0");
78
    setPropertyDefaultValue(CudaUseBlockingSync(), "true");
79
80
81
82
83
84
}

bool CudaPlatform::supportsDoublePrecision() const {
    return false;
}

85
86
87
88
89
90
91
92
93
94
95
96
const string& CudaPlatform::getPropertyValue(const Context& context, const string& property) const {
    const ContextImpl& impl = getContextImpl(context);
    const PlatformData* data = reinterpret_cast<const PlatformData*>(impl.getPlatformData());
    map<string, string>::const_iterator value = data->propertyValues.find(property);
    if (value != data->propertyValues.end())
        return value->second;
    return Platform::getPropertyValue(context, property);
}

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

97
void CudaPlatform::contextCreated(ContextImpl& context, const map<string, string>& properties) const {
98
99
100
101
    System& system = context.getSystem();
    for (int i = 0; i < system.getNumParticles(); i++)
        if (system.isVirtualSite(i))
            throw OpenMMException("CudaPlatform does not support virtual sites");
102
103
104
    for (int i = 0; i < system.getNumForces(); i++)
        if (system.getForce(i).getForceGroup() != 0)
            throw OpenMMException("CudaPlatform does not support force groups");
105
    unsigned int device = 0;
106
107
    const string& devicePropValue = (properties.find(CudaDevice()) == properties.end() ?
            getPropertyDefaultValue(CudaDevice()) : properties.find(CudaDevice())->second);
108
109
    if (devicePropValue.length() > 0)
        stringstream(devicePropValue) >> device;
Peter Eastman's avatar
Peter Eastman committed
110
    int numParticles = context.getSystem().getNumParticles();
111
112
113
    const string& blockingSync = (properties.find(CudaUseBlockingSync()) == properties.end() ?
            getPropertyDefaultValue(CudaUseBlockingSync()) : properties.find(CudaUseBlockingSync())->second);
    _gpuContext* gpu = (_gpuContext*) gpuInit(numParticles, device, blockingSync == "true");
114
    context.setPlatformData(new PlatformData(gpu));
115
116
}

117
void CudaPlatform::contextDestroyed(ContextImpl& context) const {
118
119
120
    PlatformData* data = reinterpret_cast<PlatformData*>(context.getPlatformData());
    gpuShutDown(data->gpu);
    delete data;
121
}
122

123
CudaPlatform::PlatformData::PlatformData(_gpuContext* gpu) : gpu(gpu), removeCM(false), nonbondedMethod(0), customNonbondedMethod(0), hasBonds(false), hasAngles(false),
124
        hasPeriodicTorsions(false), hasRB(false), hasNonbonded(false), hasCustomNonbonded(false), stepCount(0), computeForceCount(0), time(0.0),
Peter Eastman's avatar
Bug fix  
Peter Eastman committed
125
        ewaldSelfEnergy(0.0), dispersionCoefficient(0.0) {
126
127
128
    stringstream device;
    device << gpu->device;
    propertyValues[CudaPlatform::CudaDevice()] = device.str();
129
    propertyValues[CudaPlatform::CudaUseBlockingSync()] = (gpu->useBlockingSync ? "true" : "false");
130
}