/* -------------------------------------------------------------------------- * * 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) 2010-2016 Stanford University and the Authors. * * Authors: 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 "openmm/serialization/CustomNonbondedForceProxy.h" #include "openmm/serialization/SerializationNode.h" #include "openmm/Force.h" #include "openmm/CustomNonbondedForce.h" #include using namespace OpenMM; using namespace std; CustomNonbondedForceProxy::CustomNonbondedForceProxy() : SerializationProxy("CustomNonbondedForce") { } void CustomNonbondedForceProxy::serialize(const void* object, SerializationNode& node) const { node.setIntProperty("version", 2); const CustomNonbondedForce& force = *reinterpret_cast(object); node.setIntProperty("forceGroup", force.getForceGroup()); node.setStringProperty("energy", force.getEnergyFunction()); node.setIntProperty("method", (int) force.getNonbondedMethod()); node.setDoubleProperty("cutoff", force.getCutoffDistance()); node.setBoolProperty("useSwitchingFunction", force.getUseSwitchingFunction()); node.setDoubleProperty("switchingDistance", force.getSwitchingDistance()); node.setBoolProperty("useLongRangeCorrection", force.getUseLongRangeCorrection()); SerializationNode& perParticleParams = node.createChildNode("PerParticleParameters"); for (int i = 0; i < force.getNumPerParticleParameters(); i++) { perParticleParams.createChildNode("Parameter").setStringProperty("name", force.getPerParticleParameterName(i)); } SerializationNode& globalParams = node.createChildNode("GlobalParameters"); for (int i = 0; i < force.getNumGlobalParameters(); i++) { globalParams.createChildNode("Parameter").setStringProperty("name", force.getGlobalParameterName(i)).setDoubleProperty("default", force.getGlobalParameterDefaultValue(i)); } SerializationNode& energyDerivs = node.createChildNode("EnergyParameterDerivatives"); for (int i = 0; i < force.getNumEnergyParameterDerivatives(); i++) { energyDerivs.createChildNode("Parameter").setStringProperty("name", force.getEnergyParameterDerivativeName(i)); } SerializationNode& particles = node.createChildNode("Particles"); for (int i = 0; i < force.getNumParticles(); i++) { vector params; force.getParticleParameters(i, params); SerializationNode& particle = particles.createChildNode("Particle"); for (int j = 0; j < (int) params.size(); j++) { stringstream key; key << "param"; key << j+1; particle.setDoubleProperty(key.str(), params[j]); } } SerializationNode& exclusions = node.createChildNode("Exclusions"); for (int i = 0; i < force.getNumExclusions(); i++) { int particle1, particle2; force.getExclusionParticles(i, particle1, particle2); exclusions.createChildNode("Exclusion").setIntProperty("p1", particle1).setIntProperty("p2", particle2); } SerializationNode& functions = node.createChildNode("Functions"); for (int i = 0; i < force.getNumTabulatedFunctions(); i++) functions.createChildNode("Function", &force.getTabulatedFunction(i)).setStringProperty("name", force.getTabulatedFunctionName(i)); SerializationNode& interactionGroups = node.createChildNode("InteractionGroups"); for (int i = 0; i < force.getNumInteractionGroups(); i++) { SerializationNode& interactionGroup = interactionGroups.createChildNode("InteractionGroup"); std::set set1; std::set set2; force.getInteractionGroupParameters(i, set1, set2); SerializationNode& set1node = interactionGroup.createChildNode("Set1"); for (int p : set1) set1node.createChildNode("Particle").setIntProperty("index", p); SerializationNode& set2node = interactionGroup.createChildNode("Set2"); for (int p : set2) set2node.createChildNode("Particle").setIntProperty("index", p); } } void* CustomNonbondedForceProxy::deserialize(const SerializationNode& node) const { int version = node.getIntProperty("version"); if (version < 1 || version > 2) throw OpenMMException("Unsupported version number"); CustomNonbondedForce* force = NULL; try { CustomNonbondedForce* force = new CustomNonbondedForce(node.getStringProperty("energy")); force->setForceGroup(node.getIntProperty("forceGroup", 0)); force->setNonbondedMethod((CustomNonbondedForce::NonbondedMethod) node.getIntProperty("method")); force->setCutoffDistance(node.getDoubleProperty("cutoff")); force->setUseSwitchingFunction(node.getBoolProperty("useSwitchingFunction", false)); force->setSwitchingDistance(node.getDoubleProperty("switchingDistance", -1.0)); force->setUseLongRangeCorrection(node.getBoolProperty("useLongRangeCorrection", false)); const SerializationNode& perParticleParams = node.getChildNode("PerParticleParameters"); for (auto& parameter : perParticleParams.getChildren()) force->addPerParticleParameter(parameter.getStringProperty("name")); const SerializationNode& globalParams = node.getChildNode("GlobalParameters"); for (auto& parameter : globalParams.getChildren()) force->addGlobalParameter(parameter.getStringProperty("name"), parameter.getDoubleProperty("default")); if (version > 1) { const SerializationNode& energyDerivs = node.getChildNode("EnergyParameterDerivatives"); for (auto& parameter : energyDerivs.getChildren()) force->addEnergyParameterDerivative(parameter.getStringProperty("name")); } const SerializationNode& particles = node.getChildNode("Particles"); vector params(force->getNumPerParticleParameters()); for (auto& particle : particles.getChildren()) { for (int j = 0; j < (int) params.size(); j++) { stringstream key; key << "param"; key << j+1; params[j] = particle.getDoubleProperty(key.str()); } force->addParticle(params); } const SerializationNode& exclusions = node.getChildNode("Exclusions"); for (auto& exclusion : exclusions.getChildren()) force->addExclusion(exclusion.getIntProperty("p1"), exclusion.getIntProperty("p2")); const SerializationNode& functions = node.getChildNode("Functions"); for (auto& function : functions.getChildren()) { if (function.hasProperty("type")) { force->addTabulatedFunction(function.getStringProperty("name"), function.decodeObject()); } else { // This is an old file created before TabulatedFunction existed. const SerializationNode& valuesNode = function.getChildNode("Values"); vector values; for (auto& child : valuesNode.getChildren()) values.push_back(child.getDoubleProperty("v")); force->addTabulatedFunction(function.getStringProperty("name"), new Continuous1DFunction(values, function.getDoubleProperty("min"), function.getDoubleProperty("max"))); } } bool hasInteractionGroups = false; // Older files will be missing this block. for (auto& child : node.getChildren()) if (child.getName() == "InteractionGroups") hasInteractionGroups = true; if (hasInteractionGroups) { const SerializationNode& interactionGroups = node.getChildNode("InteractionGroups"); for (auto& interactionGroup : interactionGroups.getChildren()) { // Get set 1. const SerializationNode& set1node = interactionGroup.getChildNode("Set1"); std::set set1; for (auto& child : set1node.getChildren()) set1.insert(child.getIntProperty("index")); // Get set 2. const SerializationNode& set2node = interactionGroup.getChildNode("Set2"); std::set set2; for (auto& child : set2node.getChildren()) set2.insert(child.getIntProperty("index")); force->addInteractionGroup(set1, set2); } } return force; } catch (...) { if (force != NULL) delete force; throw; } }