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 1760 by gezelter, Thu Jun 21 19:26:46 2012 UTC

# Line 36 | Line 36
36   * [1]  Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).            
37   * [2]  Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).          
38   * [3]  Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).          
39 < * [4]  Vardeman & Gezelter, in progress (2009).                        
39 > * [4]  Kuang & Gezelter,  J. Chem. Phys. 133, 164101 (2010).
40 > * [5]  Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41   */
42  
43   #ifndef NONBONDED_NONBONDEDINTERACTION_HPP
# Line 56 | Line 57 | namespace OpenMD {
57     * being applied to any given pair of atom types.
58     */
59    enum InteractionFamily {
60 <    NO_FAMILY,             /**< No family defined */
61 <    VANDERWAALS_FAMILY,    /**< Long-range dispersion and short-range pauli repulsion */
62 <    ELECTROSTATIC_FAMILY,  /**< Coulombic and point-multipole interactions */
63 <    METALLIC_FAMILY,       /**< Transition metal interactions involving electron density */
64 <    HYDROGENBONDING_FAMILY /**< Short-range directional interactions */
60 >    NO_FAMILY = 0,             /**< No family defined */
61 >    VANDERWAALS_FAMILY = 1,    /**< Long-range dispersion and short-range pauli repulsion */
62 >    ELECTROSTATIC_FAMILY = 2,  /**< Coulombic and point-multipole interactions */
63 >    METALLIC_FAMILY = 3,       /**< Transition metal interactions involving electron density */
64 >    HYDROGENBONDING_FAMILY = 4,/**< Short-range directional interactions */
65 >    N_INTERACTION_FAMILIES = 5
66    };
67  
68 +  typedef Vector<RealType, N_INTERACTION_FAMILIES> potVec;
69 +
70    /**
71     * The InteractionData struct.  
72     *
73 <   * This is used to pass data to specific non-bonded interactions for
74 <   * force calculations.  Not all of the struct members are utilized
75 <   * by any given interaction.
73 >   * This is used to pass pointers to data to specific non-bonded
74 >   * interactions for force calculations.  Not all of the struct
75 >   * members are utilized by any given interaction.
76     */
77    struct InteractionData {
78 <    AtomType* atype1;     /**< pointer to AtomType of first atom */
79 <    AtomType* atype2;     /**< pointer to AtomType of second atom */
80 <    Vector3d d;           /**< interatomic vector (already wrapped into box) */
81 <    RealType rij;         /**< interatomic separation (precomputed) */
82 <    RealType r2;          /**< square of rij (precomputed) */
83 <    RealType rcut;        /**< cutoff radius for this interaction */
84 <    RealType sw;          /**< switching function value at rij (precomputed) */
85 <    RealType vdwMult;     /**< multiplier for van der Waals interactions */
86 <    RealType electroMult; /**< multiplier for electrostatic interactions */
87 <    RealType pot;         /**< total potential */
88 <    RealType vpair;       /**< pair potential */
89 <    Vector3d f1;         /**< force between the two atoms */
90 <    Mat3x3d eFrame1;     /**< pointer to electrostatic frame for first atom */
91 <    Mat3x3d eFrame2;     /**< pointer to electrostatic frame for second atom*/
92 <    RotMat3x3d A1;       /**< pointer to rotation matrix of first atom */
93 <    RotMat3x3d A2;       /**< pointer to rotation matrix of second atom */
94 <    Vector3d t1;         /**< pointer to torque on first atom */
95 <    Vector3d t2;         /**< pointer to torque on second atom */
96 <    RealType rho1;        /**< electron density at first atom */
97 <    RealType rho2;        /**< electron density at second atom */
98 <    RealType dfrho1;      /**< derivative of density functional for atom 1 */
99 <    RealType dfrho2;      /**< derivative of density functional for atom 2 */
100 <    RealType fshift1;     /**< indirect potential contribution from atom 1 */
101 <    RealType fshift2;     /**< indirect potential contribution from atom 2 */
78 >    pair<AtomType*, AtomType*> atypes; /**< pair of atom types interacting */
79 >    Vector3d* d;              /**< interatomic vector (already wrapped into box) */
80 >    RealType* rij;            /**< interatomic separation */
81 >    RealType* r2;             /**< square of rij */
82 >    RealType* rcut;           /**< cutoff radius for this interaction */
83 >    bool shiftedPot;          /**< shift the potential up inside the cutoff? */
84 >    bool shiftedForce;        /**< shifted forces smoothly inside the cutoff? */
85 >    RealType* sw;             /**< switching function value at rij */
86 >    int* topoDist;            /**< topological distance between atoms */
87 >    bool excluded;            /**< is this excluded from *direct* pairwise interactions? (some indirect interactions may still apply) */
88 >    RealType* vdwMult;        /**< multiplier for van der Waals interactions */
89 >    RealType* electroMult;    /**< multiplier for electrostatic interactions */
90 >    potVec* pot;              /**< total potential */
91 >    potVec* excludedPot;      /**< potential energy excluded from the overall calculation */
92 >    RealType* vpair;          /**< pair potential */
93 >    bool doParticlePot;       /**< should we bother with the particle pot? */
94 >    RealType* particlePot1;   /**< pointer to particle potential for atom1 */
95 >    RealType* particlePot2;   /**< pointer to particle potential for atom2 */
96 >    Vector3d* f1;             /**< force between the two atoms */
97 >    Mat3x3d* eFrame1;         /**< pointer to electrostatic frame for atom 1 */
98 >    Mat3x3d* eFrame2;         /**< pointer to electrostatic frame for atom 2 */
99 >    RotMat3x3d* A1;           /**< pointer to rotation matrix of first atom */
100 >    RotMat3x3d* A2;           /**< pointer to rotation matrix of second atom */
101 >    Vector3d* t1;             /**< pointer to torque on first atom */
102 >    Vector3d* t2;             /**< pointer to torque on second atom */
103 >    RealType* rho1;           /**< total electron density at first atom */
104 >    RealType* rho2;           /**< total electron density at second atom */
105 >    RealType* frho1;          /**< density functional at first atom */
106 >    RealType* frho2;          /**< density functional at second atom */
107 >    RealType* dfrho1;         /**< derivative of functional for atom 1 */
108 >    RealType* dfrho2;         /**< derivative of functional for atom 2 */
109 >    RealType* flucQ1;         /**< fluctuating charge on atom1 */
110 >    RealType* flucQ2;         /**< fluctuating charge on atom2 */
111 >    RealType* dVdFQ1;         /**< fluctuating charge force on atom1 */
112 >    RealType* dVdFQ2;         /**< fluctuating charge force on atom2 */
113 >    Vector3d* eField1;        /**< pointer to electric field on first atom */
114 >    Vector3d* eField2;        /**< pointer to electric field on second atom */
115 >    RealType* skippedCharge1; /**< charge skipped on atom1 in pairwise interaction loop with atom2 */
116 >    RealType* skippedCharge2; /**< charge skipped on atom2 in pairwise interaction loop with atom1 */
117    };
118    
119    /**
120 <   * The SkipCorrectionData struct.
102 <   *
103 <   * This is used to pass data for corrections due to "skipped" pairs
104 <   * of atoms.  These are atoms that are topologically bonded to each
105 <   * other, but which have indirect effects by way of long range
106 <   * interactions.  In the normal pair interaction loop, these atoms
107 <   * are skipped entirely, so a second pass must be made to compute
108 <   * their indirect interactions on each other.
109 <   */
110 <  struct SkipCorrectionData {
111 <    AtomType* atype1;         /**< pointer to AtomType of first atom */
112 <    AtomType* atype2;         /**< pointer to AtomType of second atom */
113 <    Vector3d d;               /**< interatomic vector (already wrapped into box) */
114 <    RealType rij;             /**< interatomic separation (precomputed) */
115 <    RealType skippedCharge1;  /**< charge skipped in normal pairwise interaction loop */
116 <    RealType skippedCharge2;  /**< charge skipped in normal pairwise interaction loop */
117 <    RealType sw;              /**< switching function value at rij (precomputed) */
118 <    RealType electroMult;     /**< multiplier for electrostatic interactions */
119 <    RealType pot;             /**< total potential */
120 <    RealType vpair;           /**< pair potential */
121 <    Vector3d f1;              /**< force correction */
122 <    Mat3x3d eFrame1;         /**< pointer to electrostatic frame for first atom */
123 <    Mat3x3d eFrame2;         /**< pointer to electrostatic frame for second atom*/
124 <    Vector3d t1;             /**< pointer to torque on first atom */
125 <    Vector3d t2;             /**< pointer to torque on second atom */
126 <  };
127 <
128 <  /**
129 <   * The SelfCorrectionData struct.
120 >   * The SelfData struct.
121     *
122 <   * This is used to pass data for the self-interaction (used by
123 <   * electrostatic methods) that have long-range corrections involving
124 <   * interactions with a medium or a boundary.
122 >   * This is used to pass pointers to data for the self-interaction or
123 >   * derived information on a single atom after a pass through all
124 >   * other interactions.  This is used by electrostatic methods that
125 >   * have long-range corrections involving interactions with a medium
126 >   * or a boundary and also by specific metal interactions for
127 >   * electron density functional calculations.  Not all of the struct
128 >   * members are utilized by any given self interaction.
129     */
130 <  struct SelfCorrectionData {
130 >  struct SelfData {
131      AtomType* atype;        /**< pointer to AtomType of the atom */
132 <    Mat3x3d eFrame;        /**< pointer to electrostatic frame for first atom */
133 <    RealType skippedCharge; /**< charge skipped in normal pairwise interaction loop */
134 <    RealType pot;           /**< total potential contribution */
135 <    Vector3d t;            /**< pointer to resultant torque on atom */
132 >    Mat3x3d* eFrame;        /**< pointer to electrostatic frame for atom */
133 >    RealType* skippedCharge;/**< charge skipped in pairwise interaction loop */
134 >    potVec* pot;            /**< total potential */
135 >    bool doParticlePot;     /**< should we bother with the particle pot? */
136 >    RealType* particlePot;  /**< contribution to potential from this particle */
137 >    Vector3d* t;            /**< pointer to resultant torque on atom */
138 >    RealType* rho;          /**< electron density */
139 >    RealType* frho;         /**< value of density functional for atom */
140 >    RealType* dfrhodrho;    /**< derivative of density functional for atom */
141 >    RealType* flucQ;        /**< current value of atom's fluctuating charge */
142 >    RealType* dVdFQ;        /**< fluctuating charge derivative */
143    };
142
143
144
145  /**
146   * The DensityData struct.  
147   *
148   * This is used to pass data to specific metallic interactions for
149   * electron density calculations.  
150   */
144    
152  struct DensityData {
153    AtomType* atype1;     /**< pointer to AtomType of first atom */
154    AtomType* atype2;     /**< pointer to AtomType of second atom */
155    RealType rij;         /**< interatomic separation (precomputed) */
156    RealType rho_i_at_j;  /**< electron density at second atom due to first */
157    RealType rho_j_at_i;  /**< electron density at first atom due to second */
158  };
159  
160  /**
161   * The FunctionalData struct.  
162   *
163   * This is used to pass data to specific metallic interactions for
164   * electron density functional calculations.  
165   */
166  
167  struct FunctionalData {
168    AtomType* atype;     /**< pointer to AtomType of the atom */
169    RealType rho;        /**< electron density (precomputed) */
170    RealType frho;       /**< value of density functional for the atom */
171    RealType dfrhodrho;  /**< derivative of density functional for the atom */
172  };
145      
146    /**
147     * The basic interface for non-bonded interactions.  
# Line 178 | Line 150 | namespace OpenMD {
150    public:
151      NonBondedInteraction() {}
152      virtual ~NonBondedInteraction() {}
153 <    virtual void calcForce(InteractionData idat) = 0;
153 >    virtual void calcForce(InteractionData &idat) = 0;
154      virtual InteractionFamily getFamily() = 0;
155 <    virtual RealType getSuggestedCutoffRadius(AtomType* at1, AtomType* at2) = 0;
155 >    virtual RealType getSuggestedCutoffRadius(pair<AtomType*, AtomType*> atypes) = 0;
156      virtual string getName() =  0;
157    };    
158  
# Line 201 | Line 173 | namespace OpenMD {
173    public:
174      ElectrostaticInteraction() : NonBondedInteraction() { }
175      virtual ~ElectrostaticInteraction() {}
176 <    virtual void calcSkipCorrection(SkipCorrectionData skdat) = 0;
177 <    virtual void calcSelfCorrection(SelfCorrectionData scdat) = 0;
206 <    virtual InteractionFamily getFamily() {return ELECTROSTATIC_FAMILY;}
176 >    virtual void calcSelfCorrection(SelfData &sdat) = 0;
177 >    virtual InteractionFamily getFamily() {return ELECTROSTATIC_FAMILY;}    
178    };    
179  
180    /**
# Line 213 | Line 184 | namespace OpenMD {
184    public:
185      MetallicInteraction() : NonBondedInteraction() { }
186      virtual ~MetallicInteraction() {}
187 <    virtual void calcDensity(DensityData ddat) = 0;
188 <    virtual void calcFunctional(FunctionalData fdat) = 0;
187 >    virtual void calcDensity(InteractionData &idat) = 0;
188 >    virtual void calcFunctional(SelfData &sdat) = 0;
189      virtual InteractionFamily getFamily() {return METALLIC_FAMILY;}
190    };
191            

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines