DuMMForceFieldSubsystem.h

Go to the documentation of this file.
00001 #ifndef SimTK_MOLMODEL_DUMM_FORCE_FIELD_SUBSYSTEM_H_
00002 #define SimTK_MOLMODEL_DUMM_FORCE_FIELD_SUBSYSTEM_H_
00003 
00004 /* -------------------------------------------------------------------------- *
00005  *                      SimTK Core: SimTK Simbody(tm)                         *
00006  * -------------------------------------------------------------------------- *
00007  * This is part of the SimTK Core biosimulation toolkit originating from      *
00008  * Simbios, the NIH National Center for Physics-Based Simulation of           *
00009  * Biological Structures at Stanford, funded under the NIH Roadmap for        *
00010  * Medical Research, grant U54 GM072970. See https://simtk.org.               *
00011  *                                                                            *
00012  * Portions copyright (c) 2007-10 Stanford University and the Authors.        *
00013  * Authors: Michael Sherman                                                   *
00014  * Contributors: Chris Bruns, Peter Eastman, Randy Radmer                     *
00015  *                                                                            *
00016  * Permission is hereby granted, free of charge, to any person obtaining a    *
00017  * copy of this software and associated documentation files (the "Software"), *
00018  * to deal in the Software without restriction, including without limitation  *
00019  * the rights to use, copy, modify, merge, publish, distribute, sublicense,   *
00020  * and/or sell copies of the Software, and to permit persons to whom the      *
00021  * Software is furnished to do so, subject to the following conditions:       *
00022  *                                                                            *
00023  * The above copyright notice and this permission notice shall be included in *
00024  * all copies or substantial portions of the Software.                        *
00025  *                                                                            *
00026  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR *
00027  * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,   *
00028  * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL    *
00029  * THE AUTHORS, CONTRIBUTORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM,    *
00030  * DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR      *
00031  * OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE  *
00032  * USE OR OTHER DEALINGS IN THE SOFTWARE.                                     *
00033  * -------------------------------------------------------------------------- */
00034 
00040 #include "SimTKcommon.h"
00041 #include "simbody/internal/ForceSubsystem.h"
00042 
00043 #include "molmodel/internal/common.h"
00044 #include "molmodel/internal/Biotype.h"
00045 
00046 #include <cassert>
00047 
00048 
00057 namespace SimTK {
00058 
00059 class MolecularMechanicsSystem;
00060 
00066 namespace DuMM {
00067 SimTK_DEFINE_UNIQUE_INDEX_TYPE(AtomIndex);
00068 SimTK_DEFINE_UNIQUE_INDEX_TYPE(AtomClassIndex);
00069 SimTK_DEFINE_UNIQUE_INDEX_TYPE(ChargedAtomTypeIndex);
00070 SimTK_DEFINE_UNIQUE_INDEX_TYPE(BondIndex);
00071 SimTK_DEFINE_UNIQUE_INDEX_TYPE(ClusterIndex);
00072 
00073 
00145 
00153 class CustomBondStretch {
00154 public:
00155     virtual ~CustomBondStretch() {}
00156     virtual Real calcEnergy(Real distance) const = 0;
00157     virtual Real calcForce(Real distance) const = 0;
00158 };
00159 
00167 class CustomBondBend {
00168 public:
00169     virtual ~CustomBondBend() {}
00170     virtual Real calcEnergy(Real bendAngle) const = 0;
00171     virtual Real calcTorque(Real bendAngle) const = 0;
00172 };
00173 
00225 class CustomBondTorsion {
00226 public:
00227     virtual ~CustomBondTorsion() {}
00228     virtual Real calcEnergy(Real dihedralAngle) const = 0;
00229     virtual Real calcTorque(Real dihedralAngle) const = 0;
00230 };
00232 
00256 static const Real Ang2Nm  = (Real)0.1L; 
00257 static const Real Nm2Ang  = (Real)10.L; 
00258 static const Real Deg2Rad = (Real)SimTK_DEGREE_TO_RADIAN; 
00259 static const Real Rad2Deg = (Real)SimTK_RADIAN_TO_DEGREE; 
00260 static const Real KJ2Kcal = (Real)SimTK_KJOULE_TO_KCAL;   
00261 static const Real Kcal2KJ = (Real)SimTK_KCAL_TO_KJOULE;   
00262 
00263 static const Real Sigma2Radius = (Real)std::pow(2.L,  1.L/6.L);
00265 static const Real Radius2Sigma = (Real)std::pow(2.L, -1.L/6.L);
00267 
00268 } // namespace DuMM
00269 
00272 
00287 class SimTK_MOLMODEL_EXPORT DuMMForceFieldSubsystem : public ForceSubsystem {
00288 public:
00290     enum VdwMixingRule {
00291         WaldmanHagler       = 1,    
00292         HalgrenHHG          = 2,    
00293         Jorgensen           = 3,    
00294         LorentzBerthelot    = 4,    
00295         Kong                = 5     
00296     };
00297 
00298     DuMMForceFieldSubsystem();
00299     explicit DuMMForceFieldSubsystem(MolecularMechanicsSystem&);
00300 
00307 
00310     const Vec3& getAtomLocationInGround(const State& state, DuMM::AtomIndex atomIx);
00313     const Vec3& getAtomVelocityInGround(const State& state, DuMM::AtomIndex atomIx);
00314 
00320     const Vec3& getAtomBodyStationInGround(const State& state, DuMM::AtomIndex atomIx);
00321 
00323 
00330 
00334     DuMM::AtomIndex addAtom(DuMM::ChargedAtomTypeIndex chargedAtomTypeIx);
00335 
00338     DuMM::BondIndex addBond(DuMM::AtomIndex atom1Ix, DuMM::AtomIndex atom2Ix);
00339 
00344     DuMM::AtomIndex  getBondAtom(DuMM::BondIndex bond, int which) const;
00345 
00347     int getNumAtoms() const;
00349     int getNumBonds() const;
00350 
00353     Real getAtomMass(DuMM::AtomIndex atomIx) const;
00356     int  getAtomElement(DuMM::AtomIndex atomIx) const;
00359     Real getAtomRadius(DuMM::AtomIndex atomIx) const;
00363     MobilizedBodyIndex getAtomBody(DuMM::AtomIndex atomIx) const;
00366     Vec3 getAtomStationOnBody(DuMM::AtomIndex atomIx) const;
00367 
00371     Vec3 getAtomStationInCluster(DuMM::AtomIndex atomIx, DuMM::ClusterIndex clusterIx) const;
00372 
00376     Vec3 getElementDefaultColor(int atomicNumber) const;
00379     Vec3 getAtomDefaultColor(DuMM::AtomIndex atomIx) const;
00380 
00381 
00383 
00396 
00400     DuMM::ClusterIndex createCluster(const char* clusterName);
00401 
00405     void placeAtomInCluster(DuMM::AtomIndex atomIx, DuMM::ClusterIndex clusterIx, const Vec3& station);
00406 
00410     void placeClusterInCluster(DuMM::ClusterIndex childClusterIndex, DuMM::ClusterIndex parentClusterIndex, 
00411                                const Transform& X_PC);
00412 
00416     MassProperties calcClusterMassProperties(DuMM::ClusterIndex clusterIx, const Transform& X_BC = Transform()) const;
00417 
00421     void attachClusterToBody(DuMM::ClusterIndex clusterIx, MobilizedBodyIndex body, const Transform& X_BC = Transform());
00422 
00424     void attachAtomToBody(DuMM::AtomIndex atomIx, MobilizedBodyIndex body, const Vec3& station = Vec3(0));
00425 
00427     MobilizedBodyIndex getClusterBody(DuMM::ClusterIndex clusterIx) const;
00428 
00431     Transform getClusterPlacementOnBody(DuMM::ClusterIndex clusterIx) const;
00432 
00435     Transform getClusterPlacementInCluster(DuMM::ClusterIndex childClusterIndex, DuMM::ClusterIndex parentClusterIndex) const;
00436 
00438 
00439 
00440         // DEFINE FORCE FIELD PARAMETERS
00441     
00463     void defineAtomClass(DuMM::AtomClassIndex atomClassIx, 
00464                          const char* atomClassName,
00465                          int elementNumber, int expectedValence,
00466                          Real vdwRadiusInNm, Real vdwWellDepthInKJ) {
00467          defineIncompleteAtomClass(atomClassIx, atomClassName, elementNumber, expectedValence);
00468          setAtomClassVdwParameters(atomClassIx, vdwRadiusInNm, vdwWellDepthInKJ);
00469     }
00470 
00473     void defineAtomClass_KA(DuMM::AtomClassIndex atomClassIx, 
00474                             const char* atomClassName,
00475                             int element, int valence,
00476                             Real vdwRadiusInAng, Real vdwWellDepthInKcal)
00477     {
00478         defineAtomClass(atomClassIx, atomClassName, element, valence,
00479             vdwRadiusInAng*DuMM::Ang2Nm, vdwWellDepthInKcal*DuMM::Kcal2KJ);
00480     }
00481     // For backwards compatibility and compactness of expression, permit integer indices
00482     void defineAtomClass_KA(int atomClassIx, const char* atomClassName,
00483                             int element, int valence,
00484                             Real vdwRadiusInAng, Real vdwWellDepthInKcal)
00485     {
00486         defineAtomClass_KA((DuMM::AtomClassIndex)atomClassIx, atomClassName, element, valence, vdwRadiusInAng, vdwWellDepthInKcal);
00487     }
00488 
00489     bool hasAtomClass(DuMM::AtomClassIndex) const; // TODO
00490     bool hasAtomClass(const String& atomClassName) const; // TODO
00491     DuMM::AtomClassIndex getAtomClassIndex(const String& atomClassName) const;
00492     DuMM::AtomClassIndex getNextUnusedAtomClassIndex() const; // TODO
00493     
00494     DuMM::AtomClassIndex getAtomClassIndex(DuMM::AtomIndex atomIx) const;
00495     Real getVdwRadius(DuMM::AtomClassIndex atomClassIx) const;
00496     Real getVdwWellDepth(DuMM::AtomClassIndex atomClassIx) const;
00497 
00499     void defineChargedAtomType(DuMM::ChargedAtomTypeIndex atomTypeIx, const char* atomTypeName,
00500         DuMM::AtomClassIndex atomClassIx, Real partialChargeInE) {
00501             defineIncompleteChargedAtomType(atomTypeIx, atomTypeName, atomClassIx);
00502             setChargedAtomTypeCharge(atomTypeIx, partialChargeInE);
00503     }
00504 
00505     void defineChargedAtomType_KA(DuMM::ChargedAtomTypeIndex atomTypeIx, const char* atomTypeName,
00506                                   DuMM::AtomClassIndex atomClassIx, Real partialChargeInE)
00507     {
00508         defineChargedAtomType(atomTypeIx, atomTypeName, atomClassIx, partialChargeInE); // easy!
00509     }
00510     // For backwards compatibility and compactness of expression, permit integer indices
00511     void defineChargedAtomType_KA(int atomTypeIx, const char* atomTypeName,
00512                                   int atomClassIx, Real partialChargeInE)
00513     {
00514         defineChargedAtomType_KA((DuMM::ChargedAtomTypeIndex) atomTypeIx, atomTypeName, (DuMM::AtomClassIndex)atomClassIx, partialChargeInE);
00515     }
00516 
00517     bool hasChargedAtomType(DuMM::ChargedAtomTypeIndex) const; // TODO
00518     bool hasChargedAtomType(const String& chargedTypeName) const; // TODO
00519     DuMM::ChargedAtomTypeIndex getChargedAtomTypeIndex(const String& chargedTypeName) const; // TODO
00520     DuMM::ChargedAtomTypeIndex getNextUnusedChargedAtomTypeIndex() const; // TODO
00521 
00523 
00531 
00533     const char* getVdwMixingRuleName(VdwMixingRule) const;
00534 
00536     void setVdwMixingRule(VdwMixingRule);
00539     VdwMixingRule getVdwMixingRule() const;
00540 
00541     void setVdw12ScaleFactor(Real); 
00542     void setVdw13ScaleFactor(Real); 
00543     void setVdw14ScaleFactor(Real); 
00544     void setVdw15ScaleFactor(Real); 
00545 
00546     void setCoulomb12ScaleFactor(Real); 
00547     void setCoulomb13ScaleFactor(Real); 
00548     void setCoulomb14ScaleFactor(Real); 
00549     void setCoulomb15ScaleFactor(Real); 
00550 
00551 
00559 
00564     void loadAmber99Parameters();
00565 
00566     // void dumpCForcefieldParameters(std::ostream& os, const String methodName) const;
00567 
00568 
00573     void populateFromTinkerParameterFile(
00574         std::istream& 
00575         );
00576 
00581     void setBiotypeChargedAtomType(
00582         DuMM::ChargedAtomTypeIndex chargedAtomTypeIndex, 
00583         BiotypeIndex biotypeIx 
00584         );
00585 
00589     DuMM::ChargedAtomTypeIndex getBiotypeChargedAtomType(BiotypeIndex biotypeIx) const;
00590 
00594     std::ostream& generateBiotypeChargedAtomTypeSelfCode(std::ostream& os) const;
00595 
00597 
00604 
00624     void defineBondStretch(DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2,
00625                            Real stiffnessInKJperNmSq, Real nominalLengthInNm);
00626 
00638     void defineBondStretch_KA(DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2,
00639                               Real stiffnessInKcalPerAngSq, Real nominalLengthInAng)
00640     {
00641         defineBondStretch(class1, class2, 
00642                           stiffnessInKcalPerAngSq * DuMM::Kcal2KJ/square(DuMM::Ang2Nm),
00643                           nominalLengthInAng      * DuMM::Ang2Nm);
00644     }
00646     void defineBondStretch_KA(int class1, int class2,
00647                               Real stiffnessInKcalPerAngSq, Real nominalLengthInAng)
00648     {   defineBondStretch_KA((DuMM::AtomClassIndex)class1, (DuMM::AtomClassIndex)class2,
00649                              stiffnessInKcalPerAngSq, nominalLengthInAng); }
00650 
00664     void defineCustomBondStretch(DuMM::AtomClassIndex       class1, 
00665                                  DuMM::AtomClassIndex       class2,
00666                                  DuMM::CustomBondStretch*   bondStretchTerm);
00667 
00669 
00676 
00698     void defineBondBend(DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3,
00699                         Real stiffnessInKJPerRadSq, Real nominalAngleInDeg);
00700 
00713     void defineBondBend_KA(DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3,
00714                            Real stiffnessInKcalPerRadSq, Real nominalAngleInDeg) 
00715     {
00716         defineBondBend(class1,class2,class3,
00717                        stiffnessInKcalPerRadSq * DuMM::Kcal2KJ,
00718                        nominalAngleInDeg);
00719     }
00721     void defineBondBend_KA(int class1, int class2, int class3,
00722                            Real stiffnessInKcalPerRadSq, Real nominalAngleInDeg)
00723     {   defineBondBend_KA((DuMM::AtomClassIndex)class1,(DuMM::AtomClassIndex)class2,(DuMM::AtomClassIndex)class3,
00724                           stiffnessInKcalPerRadSq,nominalAngleInDeg); }
00725 
00740     void defineCustomBondBend(DuMM::AtomClassIndex      class1, 
00741                               DuMM::AtomClassIndex      class2, 
00742                               DuMM::AtomClassIndex      class3,
00743                               DuMM::CustomBondBend*     bondBendTerm);
00744 
00746 
00804 
00832     void defineBondTorsion
00833        (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, 
00834         int periodicity, Real ampInKJ, Real phaseInDegrees);
00837     void defineBondTorsion
00838        (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4,  
00839         int periodicity1, Real amp1InKJ, Real phase1InDegrees,
00840         int periodicity2, Real amp2InKJ, Real phase2InDegrees);
00843     void defineBondTorsion
00844        (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, 
00845         int periodicity1, Real amp1InKJ, Real phase1InDegrees,
00846         int periodicity2, Real amp2InKJ, Real phase2InDegrees,
00847         int periodicity3, Real amp3InKJ, Real phase3InDegrees);
00848 
00851     void defineBondTorsion_KA
00852        (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, 
00853         int periodicity1, Real amp1InKcal, Real phase1InDegrees)
00854     { 
00855         defineBondTorsion(class1,class2,class3,class4,
00856                           periodicity1, amp1InKcal * DuMM::Kcal2KJ, phase1InDegrees);
00857     }
00860     void defineBondTorsion_KA
00861        (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4,  
00862         int periodicity1, Real amp1InKcal, Real phase1InDegrees,
00863         int periodicity2, Real amp2InKcal, Real phase2InDegrees)
00864     {
00865         defineBondTorsion(class1,class2,class3,class4,
00866                           periodicity1, amp1InKcal * DuMM::Kcal2KJ, phase1InDegrees,
00867                           periodicity2, amp2InKcal * DuMM::Kcal2KJ, phase2InDegrees);
00868     }
00871     void defineBondTorsion_KA
00872        (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4, 
00873         int periodicity1, Real amp1InKcal, Real phase1InDegrees,
00874         int periodicity2, Real amp2InKcal, Real phase2InDegrees,
00875         int periodicity3, Real amp3InKcal, Real phase3InDegrees)
00876     {
00877         defineBondTorsion(class1,class2,class3,class4,
00878                           periodicity1, amp1InKcal * DuMM::Kcal2KJ, phase1InDegrees,
00879                           periodicity2, amp2InKcal * DuMM::Kcal2KJ, phase2InDegrees,
00880                           periodicity3, amp3InKcal * DuMM::Kcal2KJ, phase3InDegrees);
00881     }
00883     void defineBondTorsion_KA
00884        (int class1, int class2, int class3, int class4, 
00885         int periodicity1, Real amp1InKcal, Real phase1InDegrees)
00886     {
00887         defineBondTorsion_KA
00888                ((DuMM::AtomClassIndex)class1, (DuMM::AtomClassIndex)class2, (DuMM::AtomClassIndex)class3, (DuMM::AtomClassIndex)class4, 
00889                 periodicity1, amp1InKcal, phase1InDegrees);
00890     }
00892     void defineBondTorsion_KA
00893        (int class1, int class2, int class3, int class4, 
00894         int periodicity1, Real amp1InKcal, Real phase1InDegrees,
00895         int periodicity2, Real amp2InKcal, Real phase2InDegrees)
00896     {
00897         defineBondTorsion_KA
00898                ((DuMM::AtomClassIndex)class1, (DuMM::AtomClassIndex)class2, (DuMM::AtomClassIndex)class3, (DuMM::AtomClassIndex)class4, 
00899                 periodicity1, amp1InKcal, phase1InDegrees,
00900                 periodicity2, amp2InKcal, phase2InDegrees);
00901     }
00903     void defineBondTorsion_KA
00904        (int class1, int class2, int class3, int class4, 
00905         int periodicity1, Real amp1InKcal, Real phase1InDegrees,
00906         int periodicity2, Real amp2InKcal, Real phase2InDegrees,
00907         int periodicity3, Real amp3InKcal, Real phase3InDegrees)
00908     {
00909         defineBondTorsion_KA
00910                ((DuMM::AtomClassIndex)class1, (DuMM::AtomClassIndex)class2, (DuMM::AtomClassIndex)class3, (DuMM::AtomClassIndex)class4, 
00911                 periodicity1, amp1InKcal, phase1InDegrees,
00912                 periodicity2, amp2InKcal, phase2InDegrees,
00913                 periodicity3, amp3InKcal, phase3InDegrees);
00914     }
00915 
00916 
00932     void defineCustomBondTorsion(DuMM::AtomClassIndex       class1, 
00933                                  DuMM::AtomClassIndex       class2, 
00934                                  DuMM::AtomClassIndex       class3,
00935                                  DuMM::AtomClassIndex       class4,
00936                                  DuMM::CustomBondTorsion*   bondTorsionTerm);
00937 
00939 
00947 
00948     void defineAmberImproperTorsion
00949        (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4,
00950         int periodicity, Real ampInKJ, Real phaseInDegrees);
00952     void defineAmberImproperTorsion
00953        (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4,
00954         int periodicity1, Real amp1InKJ, Real phase1InDegrees,
00955         int periodicity2, Real amp2InKJ, Real phase2InDegrees);
00957     void defineAmberImproperTorsion
00958        (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4,
00959         int periodicity1, Real amp1InKJ, Real phase1InDegrees,
00960         int periodicity2, Real amp2InKJ, Real phase2InDegrees,
00961         int periodicity3, Real amp3InKJ, Real phase3InDegrees);
00962 
00964     void defineAmberImproperTorsion_KA
00965        (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4,
00966         int periodicity1, Real amp1InKcal, Real phase1InDegrees)
00967     {
00968         defineAmberImproperTorsion(class1,class2,class3,class4,
00969                           periodicity1, amp1InKcal * DuMM::Kcal2KJ, phase1InDegrees);
00970     }
00972     void defineAmberImproperTorsion_KA
00973        (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4,
00974         int periodicity1, Real amp1InKcal, Real phase1InDegrees,
00975         int periodicity2, Real amp2InKcal, Real phase2InDegrees)
00976     {
00977         defineAmberImproperTorsion(class1,class2,class3,class4,
00978                           periodicity1, amp1InKcal * DuMM::Kcal2KJ, phase1InDegrees,
00979                           periodicity2, amp2InKcal * DuMM::Kcal2KJ, phase2InDegrees);
00980     }
00982     void defineAmberImproperTorsion_KA
00983        (DuMM::AtomClassIndex class1, DuMM::AtomClassIndex class2, DuMM::AtomClassIndex class3, DuMM::AtomClassIndex class4,
00984         int periodicity1, Real amp1InKcal, Real phase1InDegrees,
00985         int periodicity2, Real amp2InKcal, Real phase2InDegrees,
00986         int periodicity3, Real amp3InKcal, Real phase3InDegrees)
00987     {
00988         defineAmberImproperTorsion(class1,class2,class3,class4,
00989                           periodicity1, amp1InKcal * DuMM::Kcal2KJ, phase1InDegrees,
00990                           periodicity2, amp2InKcal * DuMM::Kcal2KJ, phase2InDegrees,
00991                           periodicity3, amp3InKcal * DuMM::Kcal2KJ, phase3InDegrees);
00992     }
00994 
00995 
01001     void setSolventDielectric(Real);    // typically 80 for water
01002     void setSoluteDielectric(Real);     // typically 1 or 2 for protein
01003     Real getSolventDielectric() const;
01004     Real getSoluteDielectric() const;
01005 
01006     void setGbsaIncludeAceApproximation(bool);
01007     void setGbsaIncludeAceApproximationOn()  {setGbsaIncludeAceApproximation(true );}
01008     void setGbsaIncludeAceApproximationOff() {setGbsaIncludeAceApproximation(false);}
01010 
01020     void setVdwGlobalScaleFactor(Real);                 
01021     void setCoulombGlobalScaleFactor(Real);             
01022     void setGbsaGlobalScaleFactor(Real);                
01023     void setBondStretchGlobalScaleFactor(Real);         
01024     void setBondBendGlobalScaleFactor(Real);            
01025     void setBondTorsionGlobalScaleFactor(Real);         
01026     void setAmberImproperTorsionGlobalScaleFactor(Real);
01027     void setCustomBondStretchGlobalScaleFactor(Real);   
01028     void setCustomBondBendGlobalScaleFactor(Real);      
01029     void setCustomBondTorsionGlobalScaleFactor(Real);   
01030 
01031     Real getVdwGlobalScaleFactor() const;                   
01032     Real getCoulombGlobalScaleFactor() const;               
01033     Real getGbsaGlobalScaleFactor() const;                  
01034     Real getBondStretchGlobalScaleFactor() const;           
01035     Real getBondBendGlobalScaleFactor() const;              
01036     Real getBondTorsionGlobalScaleFactor() const;           
01037     Real getAmberImproperTorsionGlobalScaleFactor() const;  
01038     Real getCustomBondStretchGlobalScaleFactor() const;     
01039     Real getCustomBondBendGlobalScaleFactor() const;        
01040     Real getCustomBondTorsionGlobalScaleFactor() const;     
01041 
01045     void setAllGlobalScaleFactors(Real s) {
01046         setVdwGlobalScaleFactor(s);
01047         setCoulombGlobalScaleFactor(s);
01048         setGbsaGlobalScaleFactor(s);
01049         setBondStretchGlobalScaleFactor(s);
01050         setBondBendGlobalScaleFactor(s);
01051         setBondTorsionGlobalScaleFactor(s);
01052         setAmberImproperTorsionGlobalScaleFactor(s);
01053         setCustomBondStretchGlobalScaleFactor(s);
01054         setCustomBondBendGlobalScaleFactor(s);
01055         setCustomBondTorsionGlobalScaleFactor(s);
01056     }
01058 
01064 
01068     void setTracing(bool);
01069 
01071     void setUseMultithreadedComputation(bool);
01073     bool getUseMultithreadedComputation() const;
01074 
01081     void setNumThreadsRequested(int);
01086     int getNumThreadsRequested() const;
01087 
01090     bool isUsingMultithreadedComputation() const;
01091 
01094     int getNumThreadsInUse() const;
01095 
01103     void setUseOpenMMAcceleration(bool);
01105     bool getUseOpenMMAcceleration() const;
01106 
01107 
01112     void setAllowOpenMMReference(bool);
01114     bool getAllowOpenMMReference() const;
01115 
01117     bool isUsingOpenMM() const;
01120     std::string getOpenMMPlatformInUse() const;
01121 
01123 
01129 
01131     long long getForceEvaluationCount() const;
01132 
01135     void dump() const;
01136 
01138     void dumpCForceFieldParameters(std::ostream& os, const String& methodName = "loadParameters") const;
01139 
01141     void loadTestMoleculeParameters();
01142 
01144 
01146     void setTraceOpenMM(bool shouldTrace) {setTracing(shouldTrace);}
01147 protected:
01148     SimTK_PIMPL_DOWNCAST(DuMMForceFieldSubsystem, ForceSubsystem);
01149     SimTK_PIMPL_DOWNCAST(DuMMForceFieldSubsystem, Subsystem);
01150 
01151     void defineIncompleteAtomClass(
01152         DuMM::AtomClassIndex classIx, 
01153         const char* name, 
01154         int elementNumber, 
01155         int valence);
01156 
01157     void defineIncompleteAtomClass_KA(
01158         DuMM::AtomClassIndex classIx, 
01159         const char* name, 
01160         int elementNumber, 
01161         int valence) 
01162     {
01163         defineIncompleteAtomClass(
01164             classIx, 
01165             name, 
01166             elementNumber, 
01167             valence
01168             );
01169     }
01170 
01171     void setAtomClassVdwParameters(DuMM::AtomClassIndex atomClassIx, Real vdwRadiusInNm, Real vdwWellDepthInKJPerMol);
01172     void setAtomClassVdwParameters_KA(DuMM::AtomClassIndex atomClassIx, Real radiusInAng, Real wellDepthInKcal) {
01173         setAtomClassVdwParameters(atomClassIx, radiusInAng*DuMM::Ang2Nm, wellDepthInKcal*DuMM::Kcal2KJ);
01174     }
01175 
01176     bool isValidAtomClass(DuMM::AtomClassIndex) const;
01177     
01178     void defineIncompleteChargedAtomType(
01179         DuMM::ChargedAtomTypeIndex typeIx, 
01180         const char* name,
01181         DuMM::AtomClassIndex classIx);
01182 
01183     void defineIncompleteChargedAtomType_KA(
01184         DuMM::ChargedAtomTypeIndex typeIx, 
01185         const char* name,
01186         DuMM::AtomClassIndex classIx) 
01187     {
01188         defineIncompleteChargedAtomType(typeIx, name, classIx);
01189     }
01190 
01191     void setChargedAtomTypeCharge(DuMM::ChargedAtomTypeIndex, Real charge);
01192     void setChargedAtomTypeCharge_KA(DuMM::ChargedAtomTypeIndex chargedAtomTypeIx, Real charge) {
01193         setChargedAtomTypeCharge(chargedAtomTypeIx, charge);
01194     }
01195 
01196 private:
01197     class DuMMForceFieldSubsystemRep& updRep();
01198     const DuMMForceFieldSubsystemRep& getRep() const;
01199 
01200 friend class MolecularMechanicsSystem;
01201 };
01202 
01207 class Amber99ForceSubsystem : public DuMMForceFieldSubsystem {
01208 public:
01209     explicit Amber99ForceSubsystem(MolecularMechanicsSystem& system) 
01210         : DuMMForceFieldSubsystem(system)
01211     {
01212         loadAmber99Parameters();
01213     }
01214 };
01215 
01216 
01217 } // namespace SimTK
01218 
01219 
01221 
01222 #endif // SimTK_MOLMODEL_DUMM_FORCE_FIELD_SUBSYSTEM_H_

Generated on Thu Aug 12 16:37:06 2010 for SimTKcore by  doxygen 1.6.1