void ReferenceIntegrateDrudeSCFStepKernel::initialize(const System& system, const DrudeSCFIntegrator& integrator, const DrudeForce& force) { // Identify Drude particles. for (int i = 0; i < force.getNumParticles(); i++) { int p, p1, p2, p3, p4; double charge, polarizability, aniso12, aniso34; force.getParticleParameters(i, p, p1, p2, p3, p4, charge, polarizability, aniso12, aniso34); drudeParticles.push_back(p); } // Record particle masses. vector<RealOpenMM> particleMass; for (int i = 0; i < system.getNumParticles(); i++) { double mass = system.getParticleMass(i); particleMass.push_back(mass); particleInvMass.push_back(mass == 0.0 ? 0.0 : 1.0/mass); } // Prepare constraints. int numConstraints = system.getNumConstraints(); if (numConstraints > 0) { vector<pair<int, int> > constraintIndices(numConstraints); vector<RealOpenMM> constraintDistances(numConstraints); for (int i = 0; i < numConstraints; ++i) { int particle1, particle2; double distance; system.getConstraintParameters(i, particle1, particle2, distance); constraintIndices[i].first = particle1; constraintIndices[i].second = particle2; constraintDistances[i] = static_cast<RealOpenMM>(distance); } vector<ReferenceCCMAAlgorithm::AngleInfo> angles; findAnglesForCCMA(system, angles); constraints = new ReferenceCCMAAlgorithm(system.getNumParticles(), numConstraints, constraintIndices, constraintDistances, particleMass, angles, (RealOpenMM)integrator.getConstraintTolerance()); } // Initialize the energy minimizer. minimizerPos = lbfgs_malloc(drudeParticles.size()*3); if (minimizerPos == NULL) throw OpenMMException("DrudeSCFIntegrator: Failed to allocate memory"); lbfgs_parameter_init(&minimizerParams); minimizerParams.linesearch = LBFGS_LINESEARCH_BACKTRACKING_STRONG_WOLFE; if (sizeof(RealOpenMM) < 8) minimizerParams.xtol = 1e-7; }
void ReferenceIntegrateDrudeLangevinStepKernel::initialize(const System& system, const DrudeLangevinIntegrator& integrator, const DrudeForce& force) { SimTKOpenMMUtilities::setRandomNumberSeed((unsigned int) integrator.getRandomNumberSeed()); // Identify particle pairs and ordinary particles. set<int> particles; vector<RealOpenMM> particleMass; for (int i = 0; i < system.getNumParticles(); i++) { particles.insert(i); double mass = system.getParticleMass(i); particleMass.push_back(mass); particleInvMass.push_back(mass == 0.0 ? 0.0 : 1.0/mass); } for (int i = 0; i < force.getNumParticles(); i++) { int p, p1, p2, p3, p4; double charge, polarizability, aniso12, aniso34; force.getParticleParameters(i, p, p1, p2, p3, p4, charge, polarizability, aniso12, aniso34); particles.erase(p); particles.erase(p1); pairParticles.push_back(make_pair(p, p1)); double m1 = system.getParticleMass(p); double m2 = system.getParticleMass(p1); pairInvTotalMass.push_back(1.0/(m1+m2)); pairInvReducedMass.push_back((m1+m2)/(m1*m2)); } normalParticles.insert(normalParticles.begin(), particles.begin(), particles.end()); // Prepare constraints. int numConstraints = system.getNumConstraints(); if (numConstraints > 0) { vector<pair<int, int> > constraintIndices(numConstraints); vector<RealOpenMM> constraintDistances(numConstraints); for (int i = 0; i < numConstraints; ++i) { int particle1, particle2; double distance; system.getConstraintParameters(i, particle1, particle2, distance); constraintIndices[i].first = particle1; constraintIndices[i].second = particle2; constraintDistances[i] = static_cast<RealOpenMM>(distance); } vector<ReferenceCCMAAlgorithm::AngleInfo> angles; findAnglesForCCMA(system, angles); constraints = new ReferenceCCMAAlgorithm(system.getNumParticles(), numConstraints, constraintIndices, constraintDistances, particleMass, angles, (RealOpenMM)integrator.getConstraintTolerance()); } }
void CudaIntegrateDrudeLangevinStepKernel::initialize(const System& system, const DrudeLangevinIntegrator& integrator, const DrudeForce& force) { cu.getPlatformData().initializeContexts(system); cu.getIntegrationUtilities().initRandomNumberGenerator((unsigned int) integrator.getRandomNumberSeed()); // Identify particle pairs and ordinary particles. set<int> particles; vector<int> normalParticleVec; vector<int2> pairParticleVec; for (int i = 0; i < system.getNumParticles(); i++) particles.insert(i); for (int i = 0; i < force.getNumParticles(); i++) { int p, p1, p2, p3, p4; double charge, polarizability, aniso12, aniso34; force.getParticleParameters(i, p, p1, p2, p3, p4, charge, polarizability, aniso12, aniso34); particles.erase(p); particles.erase(p1); pairParticleVec.push_back(make_int2(p, p1)); } normalParticleVec.insert(normalParticleVec.begin(), particles.begin(), particles.end()); normalParticles = CudaArray::create<int>(cu, max((int) normalParticleVec.size(), 1), "drudeNormalParticles"); pairParticles = CudaArray::create<int2>(cu, max((int) pairParticleVec.size(), 1), "drudePairParticles"); if (normalParticleVec.size() > 0) normalParticles->upload(normalParticleVec); if (pairParticleVec.size() > 0) pairParticles->upload(pairParticleVec); // Create kernels. map<string, string> defines; defines["NUM_ATOMS"] = cu.intToString(cu.getNumAtoms()); defines["PADDED_NUM_ATOMS"] = cu.intToString(cu.getPaddedNumAtoms()); defines["NUM_NORMAL_PARTICLES"] = cu.intToString(normalParticleVec.size()); defines["NUM_PAIRS"] = cu.intToString(pairParticleVec.size()); map<string, string> replacements; CUmodule module = cu.createModule(CudaKernelSources::vectorOps+CudaDrudeKernelSources::drudeLangevin, defines, ""); kernel1 = cu.getKernel(module, "integrateDrudeLangevinPart1"); kernel2 = cu.getKernel(module, "integrateDrudeLangevinPart2"); hardwallKernel = cu.getKernel(module, "applyHardWallConstraints"); prevStepSize = -1.0; }
void CudaCalcDrudeForceKernel::initialize(const System& system, const DrudeForce& force) { cu.setAsCurrent(); int numContexts = cu.getPlatformData().contexts.size(); int startParticleIndex = cu.getContextIndex()*force.getNumParticles()/numContexts; int endParticleIndex = (cu.getContextIndex()+1)*force.getNumParticles()/numContexts; int numParticles = endParticleIndex-startParticleIndex; if (numParticles > 0) { // Create the harmonic interaction . vector<vector<int> > atoms(numParticles, vector<int>(5)); particleParams = CudaArray::create<float4>(cu, numParticles, "drudeParticleParams"); vector<float4> paramVector(numParticles); for (int i = 0; i < numParticles; i++) { double charge, polarizability, aniso12, aniso34; force.getParticleParameters(startParticleIndex+i, atoms[i][0], atoms[i][1], atoms[i][2], atoms[i][3], atoms[i][4], charge, polarizability, aniso12, aniso34); double a1 = (atoms[i][2] == -1 ? 1 : aniso12); double a2 = (atoms[i][3] == -1 || atoms[i][4] == -1 ? 1 : aniso34); double a3 = 3-a1-a2; double k3 = ONE_4PI_EPS0*charge*charge/(polarizability*a3); double k1 = ONE_4PI_EPS0*charge*charge/(polarizability*a1) - k3; double k2 = ONE_4PI_EPS0*charge*charge/(polarizability*a2) - k3; if (atoms[i][2] == -1) { atoms[i][2] = 0; k1 = 0; } if (atoms[i][3] == -1 || atoms[i][4] == -1) { atoms[i][3] = 0; atoms[i][4] = 0; k2 = 0; } paramVector[i] = make_float4((float) k1, (float) k2, (float) k3, 0.0f); } particleParams->upload(paramVector); map<string, string> replacements; replacements["PARAMS"] = cu.getBondedUtilities().addArgument(particleParams->getDevicePointer(), "float4"); cu.getBondedUtilities().addInteraction(atoms, cu.replaceStrings(CudaDrudeKernelSources::drudeParticleForce, replacements), force.getForceGroup()); } int startPairIndex = cu.getContextIndex()*force.getNumScreenedPairs()/numContexts; int endPairIndex = (cu.getContextIndex()+1)*force.getNumScreenedPairs()/numContexts; int numPairs = endPairIndex-startPairIndex; if (numPairs > 0) { // Create the screened interaction between dipole pairs. vector<vector<int> > atoms(numPairs, vector<int>(4)); pairParams = CudaArray::create<float2>(cu, numPairs, "drudePairParams"); vector<float2> paramVector(numPairs); for (int i = 0; i < numPairs; i++) { int drude1, drude2; double thole; force.getScreenedPairParameters(startPairIndex+i, drude1, drude2, thole); int p2, p3, p4; double charge1, charge2, polarizability1, polarizability2, aniso12, aniso34; force.getParticleParameters(drude1, atoms[i][0], atoms[i][1], p2, p3, p4, charge1, polarizability1, aniso12, aniso34); force.getParticleParameters(drude2, atoms[i][2], atoms[i][3], p2, p3, p4, charge2, polarizability2, aniso12, aniso34); double screeningScale = thole/pow(polarizability1*polarizability2, 1.0/6.0); double energyScale = ONE_4PI_EPS0*charge1*charge2; paramVector[i] = make_float2((float) screeningScale, (float) energyScale); } pairParams->upload(paramVector); map<string, string> replacements; replacements["PARAMS"] = cu.getBondedUtilities().addArgument(pairParams->getDevicePointer(), "float2"); cu.getBondedUtilities().addInteraction(atoms, cu.replaceStrings(CudaDrudeKernelSources::drudePairForce, replacements), force.getForceGroup()); } cu.addForce(new CudaDrudeForceInfo(force)); }