ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/NonBondedInteraction.hpp
(Generate patch)

Comparing branches/development/src/nonbonded/NonBondedInteraction.hpp (file contents):
Revision 1505 by gezelter, Sun Oct 3 22:18:59 2010 UTC vs.
Revision 1544 by gezelter, Fri Mar 18 19:31:52 2011 UTC

# Line 60 | Line 60 | namespace OpenMD {
60      VANDERWAALS_FAMILY,    /**< Long-range dispersion and short-range pauli repulsion */
61      ELECTROSTATIC_FAMILY,  /**< Coulombic and point-multipole interactions */
62      METALLIC_FAMILY,       /**< Transition metal interactions involving electron density */
63 <    HYDROGENBONDING_FAMILY /**< Short-range directional interactions */
63 >    HYDROGENBONDING_FAMILY,/**< Short-range directional interactions */
64 >    N_INTERACTION_FAMILIES
65    };
66  
67    /**
# Line 80 | Line 81 | namespace OpenMD {
81      RealType sw;          /**< switching function value at rij (precomputed) */
82      RealType vdwMult;     /**< multiplier for van der Waals interactions */
83      RealType electroMult; /**< multiplier for electrostatic interactions */
84 <    RealType pot;         /**< total potential */
85 <    RealType vpair;       /**< pair potential */
86 <    Vector3d f1;         /**< force between the two atoms */
87 <    Mat3x3d eFrame1;     /**< pointer to electrostatic frame for first atom */
88 <    Mat3x3d eFrame2;     /**< pointer to electrostatic frame for second atom*/
89 <    RotMat3x3d A1;       /**< pointer to rotation matrix of first atom */
90 <    RotMat3x3d A2;       /**< pointer to rotation matrix of second atom */
91 <    Vector3d t1;         /**< pointer to torque on first atom */
92 <    Vector3d t2;         /**< pointer to torque on second atom */
84 >    RealType pot[4];      /**< total potential */
85 >    RealType vpair[4];    /**< pair potential */
86 >    Vector3d f1;          /**< force between the two atoms */
87 >    Mat3x3d eFrame1;      /**< pointer to electrostatic frame for first atom */
88 >    Mat3x3d eFrame2;      /**< pointer to electrostatic frame for second atom*/
89 >    RotMat3x3d A1;        /**< pointer to rotation matrix of first atom */
90 >    RotMat3x3d A2;        /**< pointer to rotation matrix of second atom */
91 >    Vector3d t1;          /**< pointer to torque on first atom */
92 >    Vector3d t2;          /**< pointer to torque on second atom */
93      RealType rho1;        /**< electron density at first atom */
94      RealType rho2;        /**< electron density at second atom */
95      RealType dfrho1;      /**< derivative of density functional for atom 1 */
# Line 108 | Line 109 | namespace OpenMD {
109     * their indirect interactions on each other.
110     */
111    struct SkipCorrectionData {
112 <    AtomType* atype1;         /**< pointer to AtomType of first atom */
113 <    AtomType* atype2;         /**< pointer to AtomType of second atom */
114 <    Vector3d d;               /**< interatomic vector (already wrapped into box) */
115 <    RealType rij;             /**< interatomic separation (precomputed) */
116 <    RealType skippedCharge1;  /**< charge skipped in normal pairwise interaction loop */
117 <    RealType skippedCharge2;  /**< charge skipped in normal pairwise interaction loop */
118 <    RealType sw;              /**< switching function value at rij (precomputed) */
119 <    RealType electroMult;     /**< multiplier for electrostatic interactions */
120 <    RealType pot;             /**< total potential */
121 <    RealType vpair;           /**< pair potential */
122 <    Vector3d f1;              /**< force correction */
123 <    Mat3x3d eFrame1;         /**< pointer to electrostatic frame for first atom */
124 <    Mat3x3d eFrame2;         /**< pointer to electrostatic frame for second atom*/
125 <    Vector3d t1;             /**< pointer to torque on first atom */
126 <    Vector3d t2;             /**< pointer to torque on second atom */
112 >    AtomType* atype1;      /**< pointer to AtomType of first atom */
113 >    AtomType* atype2;      /**< pointer to AtomType of second atom */
114 >    Vector3d d;            /**< interatomic vector (already wrapped into box) */
115 >    RealType rij;          /**< interatomic separation (precomputed) */
116 >    RealType skippedCharge1; /**< charge skipped in normal pairwise interaction loop */
117 >    RealType skippedCharge2; /**< charge skipped in normal pairwise interaction loop */
118 >    RealType sw;           /**< switching function value at rij (precomputed) */
119 >    RealType electroMult;  /**< multiplier for electrostatic interactions */
120 >    RealType pot[4];       /**< total potential */
121 >    RealType vpair[4];     /**< pair potential */
122 >    Vector3d f1;           /**< force correction */
123 >    Mat3x3d eFrame1;       /**< pointer to electrostatic frame for first atom */
124 >    Mat3x3d eFrame2;       /**< pointer to electrostatic frame for second atom*/
125 >    Vector3d t1;           /**< pointer to torque on first atom */
126 >    Vector3d t2;           /**< pointer to torque on second atom */
127    };
128  
129    /**
# Line 136 | Line 137 | namespace OpenMD {
137      AtomType* atype;        /**< pointer to AtomType of the atom */
138      Mat3x3d eFrame;        /**< pointer to electrostatic frame for first atom */
139      RealType skippedCharge; /**< charge skipped in normal pairwise interaction loop */
140 <    RealType pot;           /**< total potential contribution */
140 >    RealType pot[4];       /**< total potential contribution */
141      Vector3d t;            /**< pointer to resultant torque on atom */
142    };
143  
# Line 178 | Line 179 | namespace OpenMD {
179    public:
180      NonBondedInteraction() {}
181      virtual ~NonBondedInteraction() {}
182 <    virtual void calcForce(InteractionData idat) = 0;
182 >    virtual void calcForce(InteractionData &idat) = 0;
183      virtual InteractionFamily getFamily() = 0;
184      virtual RealType getSuggestedCutoffRadius(AtomType* at1, AtomType* at2) = 0;
185      virtual string getName() =  0;
# Line 201 | Line 202 | namespace OpenMD {
202    public:
203      ElectrostaticInteraction() : NonBondedInteraction() { }
204      virtual ~ElectrostaticInteraction() {}
205 <    virtual void calcSkipCorrection(SkipCorrectionData skdat) = 0;
206 <    virtual void calcSelfCorrection(SelfCorrectionData scdat) = 0;
207 <    virtual InteractionFamily getFamily() {return ELECTROSTATIC_FAMILY;}
205 >    virtual void calcSkipCorrection(SkipCorrectionData &skdat) = 0;
206 >    virtual void calcSelfCorrection(SelfCorrectionData &scdat) = 0;
207 >    virtual InteractionFamily getFamily() {return ELECTROSTATIC_FAMILY;}    
208    };    
209  
210    /**
# Line 213 | Line 214 | namespace OpenMD {
214    public:
215      MetallicInteraction() : NonBondedInteraction() { }
216      virtual ~MetallicInteraction() {}
217 <    virtual void calcDensity(DensityData ddat) = 0;
218 <    virtual void calcFunctional(FunctionalData fdat) = 0;
217 >    virtual void calcDensity(DensityData &ddat) = 0;
218 >    virtual void calcFunctional(FunctionalData &fdat) = 0;
219      virtual InteractionFamily getFamily() {return METALLIC_FAMILY;}
220    };
221            

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines