/* -------------------------------------------------------------------------- * * 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: Scott Le Grand, Peter Eastman * * Contributors: * * * * Permission is hereby granted, free of charge, to any person obtaining a * * copy of this software and associated documentation files (the "Software"), * * to deal in the Software without restriction, including without limitation * * the rights to use, copy, modify, merge, publish, distribute, sublicense, * * and/or sell copies of the Software, and to permit persons to whom the * * Software is furnished to do so, subject to the following conditions: * * * * The above copyright notice and this permission notice shall be included in * * all copies or substantial portions of the Software. * * * * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR * * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY, * * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL * * THE AUTHORS, CONTRIBUTORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, * * DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR * * OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE * * USE OR OTHER DEALINGS IN THE SOFTWARE. * * -------------------------------------------------------------------------- */ #include #include #include #include #include #include #include using namespace std; #include "gputypes.h" #include "cudatypes.h" #define UNROLLXX 0 #define UNROLLXY 0 struct Atom { float x; float y; float z; float q; float sig; float eps; float fx; float fy; float fz; }; static __constant__ cudaGmxSimulation cSim; void SetCalculateCDLJForcesSim(gpuContext gpu) { cudaError_t status; status = cudaMemcpyToSymbol(cSim, &gpu->sim, sizeof(cudaGmxSimulation)); RTERROR(status, "cudaMemcpyToSymbol: SetSim copy to cSim failed"); } void GetCalculateCDLJForcesSim(gpuContext gpu) { cudaError_t status; status = cudaMemcpyFromSymbol(&gpu->sim, cSim, sizeof(cudaGmxSimulation)); RTERROR(status, "cudaMemcpyFromSymbol: SetSim copy from cSim failed"); } // Include versions of the kernels for N^2 calculations. #define METHOD_NAME(a, b) a##N2##b #include "kCalculateCDLJForces.h" #define USE_OUTPUT_BUFFER_PER_WARP #undef METHOD_NAME #define METHOD_NAME(a, b) a##N2ByWarp##b #include "kCalculateCDLJForces.h" // Include versions of the kernels with cutoffs. #undef METHOD_NAME #undef USE_OUTPUT_BUFFER_PER_WARP #define USE_CUTOFF #define METHOD_NAME(a, b) a##Cutoff##b #include "kCalculateCDLJForces.h" #include "kFindInteractingBlocks.h" #define USE_OUTPUT_BUFFER_PER_WARP #undef METHOD_NAME #define METHOD_NAME(a, b) a##CutoffByWarp##b #include "kCalculateCDLJForces.h" // Include versions of the kernels with periodic boundary conditions. #undef METHOD_NAME #undef USE_OUTPUT_BUFFER_PER_WARP #define USE_PERIODIC #define METHOD_NAME(a, b) a##Periodic##b #include "kCalculateCDLJForces.h" #include "kFindInteractingBlocks.h" #define USE_OUTPUT_BUFFER_PER_WARP #undef METHOD_NAME #define METHOD_NAME(a, b) a##PeriodicByWarp##b #include "kCalculateCDLJForces.h" __global__ extern void kCalculateCDLJCutoffForces_12_kernel(); void kCalculateCDLJForces(gpuContext gpu) { // printf("kCalculateCDLJCutoffForces\n"); CUDPPResult result; size_t numWithInteractions; switch (gpu->sim.nonbondedMethod) { case NO_CUTOFF: if (gpu->bOutputBufferPerWarp) kCalculateCDLJN2ByWarpForces_kernel<<sim.nonbond_blocks, gpu->sim.nonbond_threads_per_block, sizeof(Atom)*gpu->sim.nonbond_threads_per_block>>>(gpu->sim.pWorkUnit, gpu->sim.workUnits); else kCalculateCDLJN2Forces_kernel<<sim.nonbond_blocks, gpu->sim.nonbond_threads_per_block, sizeof(Atom)*gpu->sim.nonbond_threads_per_block>>>(gpu->sim.pWorkUnit, gpu->sim.workUnits); LAUNCHERROR("kCalculateCDLJN2Forces"); break; case CUTOFF: kFindBlockBoundsCutoff_kernel<<<(gpu->psGridBoundingBox->_length+63)/64, 64>>>(); LAUNCHERROR("kFindBlockBoundsCutoff"); kFindBlocksWithInteractionsCutoff_kernel<<sim.interaction_blocks, gpu->sim.interaction_threads_per_block>>>(); LAUNCHERROR("kFindBlocksWithInteractionsCutoff"); result = cudppCompact(gpu->cudpp, gpu->sim.pInteractingWorkUnit, gpu->sim.pInteractionCount, gpu->sim.pWorkUnit, gpu->sim.pInteractionFlag, gpu->sim.workUnits); if (result != CUDPP_SUCCESS) { printf("Error in cudppCompact: %d\n", result); exit(-1); } gpu->psInteractionCount->Download(); numWithInteractions = gpu->psInteractionCount->_pSysData[0]; if (gpu->bOutputBufferPerWarp) kCalculateCDLJCutoffByWarpForces_kernel<<sim.nonbond_blocks, gpu->sim.nonbond_threads_per_block, sizeof(Atom)*gpu->sim.nonbond_threads_per_block>>>(gpu->sim.pInteractingWorkUnit, numWithInteractions); else kCalculateCDLJCutoffForces_kernel<<sim.nonbond_blocks, gpu->sim.nonbond_threads_per_block, sizeof(Atom)*gpu->sim.nonbond_threads_per_block>>>(gpu->sim.pInteractingWorkUnit, numWithInteractions); LAUNCHERROR("kCalculateCDLJCutoffForces"); break; case PERIODIC: kFindBlockBoundsPeriodic_kernel<<<(gpu->psGridBoundingBox->_length+63)/64, 64>>>(); LAUNCHERROR("kFindBlockBoundsPeriodic"); kFindBlocksWithInteractionsPeriodic_kernel<<sim.interaction_blocks, gpu->sim.interaction_threads_per_block>>>(); LAUNCHERROR("kFindBlocksWithInteractionsPeriodic"); result = cudppCompact(gpu->cudpp, gpu->sim.pInteractingWorkUnit, gpu->sim.pInteractionCount, gpu->sim.pWorkUnit, gpu->sim.pInteractionFlag, gpu->sim.workUnits); if (result != CUDPP_SUCCESS) { printf("Error in cudppCompact: %d\n", result); exit(-1); } gpu->psInteractionCount->Download(); numWithInteractions = gpu->psInteractionCount->_pSysData[0]; if (gpu->bOutputBufferPerWarp) kCalculateCDLJPeriodicByWarpForces_kernel<<sim.nonbond_blocks, gpu->sim.nonbond_threads_per_block, sizeof(Atom)*gpu->sim.nonbond_threads_per_block>>>(gpu->sim.pInteractingWorkUnit, numWithInteractions); else kCalculateCDLJPeriodicForces_kernel<<sim.nonbond_blocks, gpu->sim.nonbond_threads_per_block, sizeof(Atom)*gpu->sim.nonbond_threads_per_block>>>(gpu->sim.pInteractingWorkUnit, numWithInteractions); LAUNCHERROR("kCalculateCDLJPeriodicForces"); } }