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 1504 by gezelter, Sat Oct 2 20:41:53 2010 UTC vs.
Revision 1587 by gezelter, Fri Jul 8 20:25:32 2011 UTC

# Line 56 | Line 56 | namespace OpenMD {
56     * being applied to any given pair of atom types.
57     */
58    enum InteractionFamily {
59 <    NO_FAMILY,             /**< No family defined */
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 */
59 >    NO_FAMILY = 0,             /**< No family defined */
60 >    VANDERWAALS_FAMILY = 1,    /**< Long-range dispersion and short-range pauli repulsion */
61 >    ELECTROSTATIC_FAMILY = 2,  /**< Coulombic and point-multipole interactions */
62 >    METALLIC_FAMILY = 3,       /**< Transition metal interactions involving electron density */
63 >    HYDROGENBONDING_FAMILY = 4,/**< Short-range directional interactions */
64 >    N_INTERACTION_FAMILIES = 5
65    };
66  
67 +  typedef Vector<RealType, N_INTERACTION_FAMILIES> potVec;
68 +
69    /**
70     * The InteractionData struct.  
71     *
72 <   * This is used to pass data to specific non-bonded interactions for
73 <   * force calculations.  Not all of the struct members are utilized
74 <   * by any given interaction.
72 >   * This is used to pass pointers to data to specific non-bonded
73 >   * interactions for force calculations.  Not all of the struct
74 >   * members are utilized by any given interaction.
75     */
76    struct InteractionData {
77 <    AtomType* atype1;     /**< pointer to AtomType of first atom */
78 <    AtomType* atype2;     /**< pointer to AtomType of second atom */
79 <    Vector3d d;           /**< interatomic vector (already wrapped into box) */
80 <    RealType rij;         /**< interatomic separation (precomputed) */
81 <    RealType r2;          /**< square of rij (precomputed) */
82 <    RealType rcut;        /**< cutoff radius for this interaction */
83 <    RealType sw;          /**< switching function value at rij (precomputed) */
84 <    RealType vdwMult;     /**< multiplier for van der Waals interactions */
85 <    RealType electroMult; /**< multiplier for electrostatic interactions */
86 <    RealType pot;         /**< total potential */
87 <    RealType vpair;       /**< pair potential */
88 <    Vector3d f1;         /**< force between the two atoms */
89 <    Mat3x3d eFrame1;     /**< pointer to electrostatic frame for first atom */
90 <    Mat3x3d eFrame2;     /**< pointer to electrostatic frame for second atom*/
91 <    RotMat3x3d A1;       /**< pointer to rotation matrix of first atom */
92 <    RotMat3x3d A2;       /**< pointer to rotation matrix of second atom */
93 <    Vector3d t1;         /**< pointer to torque on first atom */
94 <    Vector3d t2;         /**< pointer to torque on second atom */
95 <    RealType rho1;        /**< electron density at first atom */
96 <    RealType rho2;        /**< electron density at second atom */
97 <    RealType dfrho1;      /**< derivative of density functional for atom 1 */
98 <    RealType dfrho2;      /**< derivative of density functional for atom 2 */
99 <    RealType fshift1;     /**< indirect potential contribution from atom 1 */
100 <    RealType fshift2;     /**< indirect potential contribution from atom 2 */
77 >    pair<AtomType*, AtomType*> atypes; /**< pair of atom types interacting */
78 >    Vector3d* d;              /**< interatomic vector (already wrapped into box) */
79 >    RealType* rij;            /**< interatomic separation */
80 >    RealType* r2;             /**< square of rij */
81 >    RealType* rcut;           /**< cutoff radius for this interaction */
82 >    bool shiftedPot;          /**< shift the potential up inside the cutoff? */
83 >    bool shiftedForce;        /**< shifted forces smoothly inside the cutoff? */
84 >    RealType* sw;             /**< switching function value at rij */
85 >    int* topoDist;            /**< topological distance between atoms */
86 >    bool excluded;            /**< is this excluded from *direct* pairwise interactions? (some indirect interactions may still apply) */
87 >    RealType* vdwMult;        /**< multiplier for van der Waals interactions */
88 >    RealType* electroMult;    /**< multiplier for electrostatic interactions */
89 >    potVec* pot;              /**< total potential */
90 >    RealType* vpair;          /**< pair potential */
91 >    RealType* particlePot1;   /**< pointer to particle potential for atom1 */
92 >    RealType* particlePot2;   /**< pointer to particle potential for atom2 */
93 >    Vector3d* f1;             /**< force between the two atoms */
94 >    Mat3x3d* eFrame1;         /**< pointer to electrostatic frame for atom 1 */
95 >    Mat3x3d* eFrame2;         /**< pointer to electrostatic frame for atom 2 */
96 >    RotMat3x3d* A1;           /**< pointer to rotation matrix of first atom */
97 >    RotMat3x3d* A2;           /**< pointer to rotation matrix of second atom */
98 >    Vector3d* t1;             /**< pointer to torque on first atom */
99 >    Vector3d* t2;             /**< pointer to torque on second atom */
100 >    RealType* rho1;           /**< total electron density at first atom */
101 >    RealType* rho2;           /**< total electron density at second atom */
102 >    RealType* frho1;          /**< density functional at first atom */
103 >    RealType* frho2;          /**< density functional at second atom */
104 >    RealType* dfrho1;         /**< derivative of functional for atom 1 */
105 >    RealType* dfrho2;         /**< derivative of functional for atom 2 */
106 >    RealType* skippedCharge1; /**< charge skipped on atom1 in pairwise interaction loop with atom2 */
107 >    RealType* skippedCharge2; /**< charge skipped on atom2 in pairwise interaction loop with atom1 */
108    };
109    
110    /**
111 <   * 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.
111 >   * The SelfData struct.
112     *
113 <   * This is used to pass data for the self-interaction (used by
114 <   * electrostatic methods) that have long-range corrections involving
115 <   * interactions with a medium or a boundary.
113 >   * This is used to pass pointers to data for the self-interaction or
114 >   * derived information on a single atom after a pass through all
115 >   * other interactions.  This is used by electrostatic methods that
116 >   * have long-range corrections involving interactions with a medium
117 >   * or a boundary and also by specific metal interactions for
118 >   * electron density functional calculations.  Not all of the struct
119 >   * members are utilized by any given self interaction.
120     */
121 <  struct SelfCorrectionData {
121 >  struct SelfData {
122      AtomType* atype;        /**< pointer to AtomType of the atom */
123 <    Mat3x3d eFrame;        /**< pointer to electrostatic frame for first atom */
124 <    RealType skippedCharge; /**< charge skipped in normal pairwise interaction loop */
125 <    RealType pot;           /**< total potential contribution */
126 <    Vector3d t;            /**< pointer to resultant torque on atom */
123 >    Mat3x3d* eFrame;        /**< pointer to electrostatic frame for atom */
124 >    RealType* skippedCharge;/**< charge skipped in pairwise interaction loop */
125 >    potVec* pot;            /**< total potential */
126 >    RealType* particlePot;  /**< contribution to potential from this particle */
127 >    Vector3d* t;            /**< pointer to resultant torque on atom */
128 >    RealType* rho;          /**< electron density */
129 >    RealType* frho;         /**< value of density functional for atom */
130 >    RealType* dfrhodrho;    /**< derivative of density functional for atom */
131    };
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   */
132    
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  };
133      
134    /**
135     * The basic interface for non-bonded interactions.  
# Line 178 | Line 138 | namespace OpenMD {
138    public:
139      NonBondedInteraction() {}
140      virtual ~NonBondedInteraction() {}
141 <    virtual void calcForce(InteractionData idat) = 0;
141 >    virtual void calcForce(InteractionData &idat) = 0;
142      virtual InteractionFamily getFamily() = 0;
143 +    virtual RealType getSuggestedCutoffRadius(pair<AtomType*, AtomType*> atypes) = 0;
144      virtual string getName() =  0;
145    };    
146  
# Line 200 | Line 161 | namespace OpenMD {
161    public:
162      ElectrostaticInteraction() : NonBondedInteraction() { }
163      virtual ~ElectrostaticInteraction() {}
164 <    virtual void calcSkipCorrection(SkipCorrectionData skdat) = 0;
165 <    virtual void calcSelfCorrection(SelfCorrectionData scdat) = 0;
205 <    virtual InteractionFamily getFamily() {return ELECTROSTATIC_FAMILY;}
164 >    virtual void calcSelfCorrection(SelfData &sdat) = 0;
165 >    virtual InteractionFamily getFamily() {return ELECTROSTATIC_FAMILY;}    
166    };    
167  
168    /**
# Line 212 | Line 172 | namespace OpenMD {
172    public:
173      MetallicInteraction() : NonBondedInteraction() { }
174      virtual ~MetallicInteraction() {}
175 <    virtual void calcDensity(DensityData ddat) = 0;
176 <    virtual void calcFunctional(FunctionalData fdat) = 0;
175 >    virtual void calcDensity(InteractionData &idat) = 0;
176 >    virtual void calcFunctional(SelfData &sdat) = 0;
177      virtual InteractionFamily getFamily() {return METALLIC_FAMILY;}
178    };
179            

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines