ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/parallel/ForceDecomposition.cpp
Revision: 1736
Committed: Tue Jun 5 17:51:31 2012 UTC (12 years, 10 months ago) by jmichalk
Original Path: branches/development/src/parallel/ForceDecomposition.cpp
File size: 6590 byte(s)
Log Message:
Fixed fluctuating charge forces

File Contents

# User Rev Content
1 gezelter 1572 /*
2     * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3     *
4     * The University of Notre Dame grants you ("Licensee") a
5     * non-exclusive, royalty free, license to use, modify and
6     * redistribute this software in source and binary code form, provided
7     * that the following conditions are met:
8     *
9     * 1. Redistributions of source code must retain the above copyright
10     * notice, this list of conditions and the following disclaimer.
11     *
12     * 2. Redistributions in binary form must reproduce the above copyright
13     * notice, this list of conditions and the following disclaimer in the
14     * documentation and/or other materials provided with the
15     * distribution.
16     *
17     * This software is provided "AS IS," without a warranty of any
18     * kind. All express or implied conditions, representations and
19     * warranties, including any implied warranty of merchantability,
20     * fitness for a particular purpose or non-infringement, are hereby
21     * excluded. The University of Notre Dame and its licensors shall not
22     * be liable for any damages suffered by licensee as a result of
23     * using, modifying or distributing the software or its
24     * derivatives. In no event will the University of Notre Dame or its
25     * licensors be liable for any lost revenue, profit or data, or for
26     * direct, indirect, special, consequential, incidental or punitive
27     * damages, however caused and regardless of the theory of liability,
28     * arising out of the use of or inability to use software, even if the
29     * University of Notre Dame has been advised of the possibility of
30     * such damages.
31     *
32     * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33     * research, please cite the appropriate papers when you publish your
34     * work. Good starting points are:
35     *
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 gezelter 1665 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40     * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 gezelter 1572 */
42     #include "parallel/ForceDecomposition.hpp"
43     #include "math/Vector3.hpp"
44     #ifdef IS_MPI
45     #include <mpi.h>
46     #endif
47    
48     using namespace std;
49     namespace OpenMD {
50    
51 gezelter 1723 ForceDecomposition::ForceDecomposition(SimInfo* info, InteractionManager* iMan) : info_(info), interactionMan_(iMan), needVelocities_(false) {
52 gezelter 1572 sman_ = info_->getSnapshotManager();
53     storageLayout_ = sman_->getStorageLayout();
54     ff_ = info_->getForceField();
55 gezelter 1576 userChoseCutoff_ = false;
56 gezelter 1572
57     Globals* simParams_ = info_->getSimParams();
58 gezelter 1723 if (simParams_->havePrintHeatFlux()) {
59     if (simParams_->getPrintHeatFlux()) {
60     needVelocities_ = true;
61     }
62     }
63    
64 gezelter 1572 if (simParams_->haveSkinThickness()) {
65     skinThickness_ = simParams_->getSkinThickness();
66     } else {
67     skinThickness_ = 1.0;
68     sprintf(painCave.errMsg,
69     "ForceDecomposition: No value was set for the skinThickness.\n"
70     "\tOpenMD will use a default value of %f Angstroms\n"
71     "\tfor this simulation\n", skinThickness_);
72     painCave.severity = OPENMD_INFO;
73     painCave.isFatal = 0;
74     simError();
75     }
76    
77     // cellOffsets are the partial space for the cell lists used in
78     // constructing the neighbor lists
79     cellOffsets_.clear();
80     cellOffsets_.push_back( Vector3i(0, 0, 0) );
81     cellOffsets_.push_back( Vector3i(1, 0, 0) );
82     cellOffsets_.push_back( Vector3i(1, 1, 0) );
83     cellOffsets_.push_back( Vector3i(0, 1, 0) );
84     cellOffsets_.push_back( Vector3i(-1,1, 0) );
85     cellOffsets_.push_back( Vector3i(0, 0, 1) );
86     cellOffsets_.push_back( Vector3i(1, 0, 1) );
87     cellOffsets_.push_back( Vector3i(1, 1, 1) );
88     cellOffsets_.push_back( Vector3i(0, 1, 1) );
89     cellOffsets_.push_back( Vector3i(-1,1, 1) );
90     cellOffsets_.push_back( Vector3i(-1,0, 1) );
91     cellOffsets_.push_back( Vector3i(-1,-1,1) );
92     cellOffsets_.push_back( Vector3i(0, -1,1) );
93     cellOffsets_.push_back( Vector3i(1, -1,1) );
94     }
95    
96 gezelter 1582 void ForceDecomposition::fillSelfData(SelfData &sdat, int atom1) {
97 gezelter 1583
98 gezelter 1591 sdat.atype = atypesLocal[atom1];
99 gezelter 1613
100     sdat.pot = &embeddingPot;
101 gezelter 1629
102 gezelter 1572 if (storageLayout_ & DataStorage::dslElectroFrame) {
103     sdat.eFrame = &(snap_->atomData.electroFrame[atom1]);
104     }
105    
106     if (storageLayout_ & DataStorage::dslTorque) {
107     sdat.t = &(snap_->atomData.torque[atom1]);
108     }
109    
110     if (storageLayout_ & DataStorage::dslDensity) {
111     sdat.rho = &(snap_->atomData.density[atom1]);
112     }
113    
114     if (storageLayout_ & DataStorage::dslFunctional) {
115     sdat.frho = &(snap_->atomData.functional[atom1]);
116     }
117    
118     if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
119     sdat.dfrhodrho = &(snap_->atomData.functionalDerivative[atom1]);
120     }
121    
122 gezelter 1586 if (storageLayout_ & DataStorage::dslSkippedCharge) {
123     sdat.skippedCharge = &(snap_->atomData.skippedCharge[atom1]);
124     }
125    
126 gezelter 1575 if (storageLayout_ & DataStorage::dslParticlePot) {
127     sdat.particlePot = &(snap_->atomData.particlePot[atom1]);
128     }
129 jmichalk 1736
130     if (storageLayout_ & DataStorage::dslFlucQPosition) {
131     sdat.flucQ = &(snap_->atomData.flucQPos[atom1]);
132     }
133    
134     if (storageLayout_ & DataStorage::dslFlucQForce) {
135     sdat.dVdFQ = &(snap_->atomData.flucQFrc[atom1]);
136     }
137 gezelter 1572 }
138    
139     bool ForceDecomposition::checkNeighborList() {
140    
141     int nGroups = snap_->cgData.position.size();
142 gezelter 1723 if (needVelocities_)
143     snap_->cgData.setStorageLayout(DataStorage::dslPosition | DataStorage::dslVelocity);
144    
145 gezelter 1572 // if we have changed the group identities or haven't set up the
146     // saved positions we automatically will need a neighbor list update:
147    
148 gezelter 1587 if ( saved_CG_positions_.size() != nGroups ) return true;
149 gezelter 1572
150     RealType dispmax = 0.0;
151     Vector3d disp;
152    
153     for (int i = 0; i < nGroups; i++) {
154     disp = snap_->cgData.position[i] - saved_CG_positions_[i];
155     for (int j = 0; j < 3; j++)
156     dispmax = max( abs(disp[j]), dispmax);
157     }
158    
159     #ifdef IS_MPI
160     MPI::COMM_WORLD.Allreduce(&dispmax, &dispmax, 1, MPI::REALTYPE, MPI::MAX);
161     #endif
162    
163     // a conservative test of list skin crossings
164     dispmax = 2.0 * sqrt (3.0 * dispmax * dispmax);
165    
166 gezelter 1582
167 gezelter 1587 if (dispmax > skinThickness_)
168 gezelter 1582 return (dispmax > skinThickness_);
169 gezelter 1587
170 gezelter 1582 return false;
171 gezelter 1572 }
172 gezelter 1723
173     void ForceDecomposition::addToHeatFlux(Vector3d hf) {
174     snap_->frameData.conductiveHeatFlux += hf;
175     }
176     void ForceDecomposition::setHeatFlux(Vector3d hf) {
177     snap_->frameData.conductiveHeatFlux = hf;
178     }
179    
180 gezelter 1572 }

Properties

Name Value
svn:eol-style native