// ----------------------------------------------------------------------------- // COPY STATE BACK TO CPU FROM OPENMM // ----------------------------------------------------------------------------- static void myGetOpenMMState(MyOpenMMData* omm, bool wantEnergy, double& timeInPs, double& energyInKcal, MyAtomInfo atoms[]) { int infoMask = 0; infoMask = OpenMM::State::Positions; if (wantEnergy) { infoMask += OpenMM::State::Velocities; // for kinetic energy (cheap) infoMask += OpenMM::State::Energy; // for pot. energy (expensive) } // Forces are also available (and cheap). const OpenMM::State state = omm->context->getState(infoMask); timeInPs = state.getTime(); // OpenMM time is in ps already // Copy OpenMM positions into atoms array and change units from nm to Angstroms. const std::vector<Vec3>& positionsInNm = state.getPositions(); for (int i=0; i < (int)positionsInNm.size(); ++i) for (int j=0; j < 3; ++j) atoms[i].posInAng[j] = positionsInNm[i][j] * OpenMM::AngstromsPerNm; // If energy has been requested, obtain it and convert from kJ to kcal. energyInKcal = 0; if (wantEnergy) energyInKcal = (state.getPotentialEnergy() + state.getKineticEnergy()) * OpenMM::KcalPerKJ; }
void writeTimeSeries(int steps, const OpenMM::State& state, std::ofstream &ots, double qscore) { ots << std::setw(8) << steps; ots << ' ' << std::fixed << std::setw(8) << std::setprecision(3) << qscore; ots << ' ' << std::fixed << std::setw(8) << state.getPotentialEnergy(); ots << ' ' << std::fixed << std::setw(8) << state.getKineticEnergy(); ots << '\n'; }
void OpenMMIntegrator::run( int numTimesteps ) { preStepModify(); bool execute = true; /*#ifdef HAVE_OPENMM_LTMD if( mLTMDParameters.ShouldProtoMolDiagonalize && mLTMDParameters.ShouldForceRediagOnMinFail ){ if( app->eigenInfo.OpenMMMinimize ){ OpenMM::LTMD::Integrator *ltmd = dynamic_cast<OpenMM::LTMD::Integrator*>( integrator ); bool minimizePassed = ltmd->minimize( 50, 2 ); if( minimizePassed || app->eigenInfo.RediagonalizationCount >= 5 ){ if( app->eigenInfo.RediagonalizationCount >= 5 ){ std::cout << "Maximum Rediagonalizations Reached" << std::endl; } execute = true; app->eigenInfo.OpenMMMinimize = false; app->eigenInfo.RediagonalizationCount = 0; std::cout << "Exiting Rediagonalizations" << std::endl; }else{ execute = false; app->eigenInfo.reDiagonalize = true; app->eigenInfo.RediagonalizationCount++; } } } #endif */ integrator->step( numTimesteps ); // Retrive data const OpenMM::State state = context->getState( OpenMM::State::Positions | OpenMM::State::Velocities | OpenMM::State::Forces | OpenMM::State::Energy ); const std::vector<OpenMM::Vec3> positions( state.getPositions() ); const std::vector<OpenMM::Vec3> velocities( state.getVelocities() ); const std::vector<OpenMM::Vec3> forces( state.getForces() ); const unsigned int sz = app->positions.size(); for( unsigned int i = 0; i < sz; ++i ) { for( int j = 0; j < 3; j++ ) { app->positions[i].c[j] = positions[i][j] * Constant::NM_ANGSTROM; //nm to A app->velocities[i].c[j] = velocities[i][j] * Constant::NM_ANGSTROM * Constant::TIMEFACTOR * Constant::FS_PS; //nm/ps to A/fs? ( *myForces )[i].c[j] = forces[i][j] * Constant::INV_NM_ANGSTROM * Constant::KJ_KCAL; //KJ/nm to Kcal/A } } app->energies.clear(); //clear old energies app->energies[ScalarStructure::COULOMB] = app->energies[ScalarStructure::LENNARDJONES] = app->energies[ScalarStructure::BOND] = app->energies[ScalarStructure::ANGLE] = app->energies[ScalarStructure::DIHEDRAL] = app->energies[ScalarStructure::IMPROPER] = 0.0; //save total potential energy app->energies[ScalarStructure::OTHER] = state.getPotentialEnergy() * Constant::KJ_KCAL; //fix time as no forces calculated app->topology->time += numTimesteps * getTimestep(); postStepModify(); }
void simulate(json &molinfo) { // Load any shared libraries containing GPU implementations. OpenMM::Platform::loadPluginsFromDirectory( OpenMM::Platform::getDefaultPluginsDirectory()); // force to use CPU platform that is better than CUDA in my simulation scale OpenMM::Platform& platform = OpenMM::Platform::getPlatformByName("CPU"); // force to use single thread platform.setPropertyDefaultValue("Threads", "1"); // Create a system with forces. OpenMM::System system; const int N_particles = molinfo["resSeq"].size(); std::cout << "# Number of Particles : " << N_particles << '\n'; // set MD paramters const double Temperature = molinfo["parameters"]["Temperature"]; const int SimulationSteps = molinfo["parameters"]["SimulationSteps"]; const double TimePerStepInPs = molinfo["parameters"]["TimePerStepInPs"]; const int NStepSave = molinfo["parameters"]["NStepSave"]; const int RandomSeed = molinfo["parameters"]["RandomSeed"]; // set other parameters const std::string filename = molinfo["output"]["filename"]; // Create Atoms std::vector<OpenMM::Vec3> initPosInNm(N_particles); // map from resSeq to particle index(0..N-1) for (int i=0; i<N_particles; ++i) { int resSeq = molinfo["resSeq"][i]; rs2pi.insert(std::pair<int,int>(resSeq, i)); } double x, y, z; for (int i=0; i<N_particles; ++i) { x = molinfo["position"][i]["position"][0]; y = molinfo["position"][i]["position"][1]; z = molinfo["position"][i]["position"][2]; initPosInNm[i] = OpenMM::Vec3(x, y, z) * OpenMM::NmPerAngstrom; system.addParticle(137.0); // average mass of amino acid residues } // add non-local repulsive force // This pair potential is excluded from bonded pairs OpenMM::CustomNonbondedForce& nonlocalRepulsion = *new OpenMM::CustomNonbondedForce("epsilon_repul*(d_exclusive/r)^12"); system.addForce(&nonlocalRepulsion); nonlocalRepulsion.setNonbondedMethod(OpenMM::CustomNonbondedForce::NonbondedMethod::CutoffNonPeriodic); nonlocalRepulsion.setCutoffDistance(D_ExclusionCutoffInNm); nonlocalRepulsion.addGlobalParameter("d_exclusive", D_ExclusiveInNm); nonlocalRepulsion.addGlobalParameter("epsilon_repul", E_ExclusionPair); for (int i=0; i<N_particles; ++i) { nonlocalRepulsion.addParticle(); } // add non-local Go contact for native contact pairs OpenMM::CustomBondForce& goContactForce = *new OpenMM::CustomBondForce("epsilon_ngo*(5*(r_native/r)^12 - 6*(r_native/r)^10)"); system.addForce(&goContactForce); goContactForce.addGlobalParameter("epsilon_ngo", E_GoContactPair); goContactForce.addPerBondParameter("r_native"); for (int i=0; i<molinfo["contact"].size(); ++i) { int p1 = rs2pi[molinfo["contact"][i]["resSeq"][0]]; int p2 = rs2pi[molinfo["contact"][i]["resSeq"][1]]; double r_native = molinfo["contact"][i]["length"]; std::vector<double> perBondParams = {r_native * OpenMM::NmPerAngstrom}; goContactForce.addBond(p1, p2, perBondParams); nonlocalRepulsion.addExclusion(p1, p2); } // add bonding force or constraints between two particles OpenMM::HarmonicBondForce& bondStretch = *new OpenMM::HarmonicBondForce(); system.addForce(&bondStretch); for (int i=0; i<molinfo["bond"].size(); ++i) { int p1 = rs2pi[molinfo["bond"][i]["resSeq"][0]]; int p2 = rs2pi[molinfo["bond"][i]["resSeq"][1]]; double length = molinfo["bond"][i]["length"]; if (UseConstraints) { system.addConstraint(p1, p2, length * OpenMM::NmPerAngstrom); } else { bondStretch.addBond(p1, p2, length * OpenMM::NmPerAngstrom, K_BondStretch); } nonlocalRepulsion.addExclusion(p1, p2); } // add angle force OpenMM::HarmonicAngleForce& bondBend = *new OpenMM::HarmonicAngleForce(); system.addForce(&bondBend); for (int i=0; i<molinfo["angle"].size(); ++i) { int p1 = rs2pi[molinfo["angle"][i]["resSeq"][0]]; int p2 = rs2pi[molinfo["angle"][i]["resSeq"][1]]; int p3 = rs2pi[molinfo["angle"][i]["resSeq"][2]]; double angle = molinfo["angle"][i]["angle"]; bondBend.addAngle(p1, p2, p3, angle * OpenMM::RadiansPerDegree, K_BondAngle); nonlocalRepulsion.addExclusion(p1, p3); } // add dihedral force OpenMM::PeriodicTorsionForce& bondTorsion = *new OpenMM::PeriodicTorsionForce(); system.addForce(&bondTorsion); for (int i=0; i<molinfo["dihedral"].size(); ++i) { int p1 = rs2pi[molinfo["dihedral"][i]["resSeq"][0]]; int p2 = rs2pi[molinfo["dihedral"][i]["resSeq"][1]]; int p3 = rs2pi[molinfo["dihedral"][i]["resSeq"][2]]; int p4 = rs2pi[molinfo["dihedral"][i]["resSeq"][3]]; double dihedral = molinfo["dihedral"][i]["dihedral"]; bondTorsion.addTorsion(p1, p2, p3, p4, 1, (dihedral+180.0) * OpenMM::RadiansPerDegree, K_BondTorsion1); bondTorsion.addTorsion(p1, p2, p3, p4, 3, (dihedral+180.0) * OpenMM::RadiansPerDegree, K_BondTorsion3); nonlocalRepulsion.addExclusion(p1, p4); } // set langevin integrator OpenMM::LangevinIntegrator integrator(Temperature, LangevinFrictionPerPs, TimePerStepInPs); integrator.setRandomNumberSeed(RandomSeed); // Let OpenMM Context choose best platform OpenMM::Context context(system, integrator, platform); std::cout << "# Using OpenMM Platform: "; std::cout << context.getPlatform().getName().c_str() << std::endl; // Set starting positions of the atoms. Leave time and velocity zero. context.setPositions(initPosInNm); // set initial velocity context.setVelocitiesToTemperature(Temperature, RandomSeed); // set output files std::ofstream opdb, ots; opdb.open(filename + ".pdb", std::ios::out); ots.open(filename + ".ts", std::ios::out); std::cout << "# PDB: " << filename + ".pdb\n"; std::cout << "# TimeSereis: " << filename + ".ts\n"; int infoMask = 0; infoMask = OpenMM::State::Positions; infoMask += OpenMM::State::Energy; for (int frameNum = 0;; ++frameNum) { // Output current state information OpenMM::State state = context.getState(infoMask); int steps = frameNum * NStepSave; double qscore = getQscore(state, molinfo); // write PDB frame writePDBFrame(frameNum, state, molinfo, opdb); writeTimeSeries(steps, state, ots, qscore); // show progress in stdout if (SilentMode == false) { std::cout << '\r'; std::cout << "Q=" << std::setw(5) << std::setprecision(3) << qscore; std::cout << "; T=" << std::setw(5) << integrator.getTemperature(); std::cout << "; E=" << std::setw(5) << state.getPotentialEnergy(); std::cout << std::setw(8) << steps << " steps / "; std::cout << std::setw(8) << SimulationSteps << std::flush; } if (frameNum * NStepSave >= SimulationSteps) break; integrator.step(NStepSave); } std::cout << std::endl; opdb.close(); ots.close(); }