ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/devel_omp/src/parallel/ForceDecomposition.cpp
Revision: 1595
Committed: Tue Jul 19 18:50:04 2011 UTC (13 years, 9 months ago) by chuckv
File size: 7338 byte(s)
Log Message:
Adding initial OpenMP support using new neighbor lists.


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     * [4] Vardeman & Gezelter, in progress (2009).
40     */
41     #include "parallel/ForceDecomposition.hpp"
42     #include "math/Vector3.hpp"
43     #ifdef IS_MPI
44     #include <mpi.h>
45     #endif
46    
47     using namespace std;
48     namespace OpenMD {
49    
50 gezelter 1579 ForceDecomposition::ForceDecomposition(SimInfo* info, InteractionManager* iMan) : info_(info), interactionMan_(iMan) {
51 gezelter 1572 sman_ = info_->getSnapshotManager();
52     storageLayout_ = sman_->getStorageLayout();
53     ff_ = info_->getForceField();
54 gezelter 1576 userChoseCutoff_ = false;
55 gezelter 1572
56     Globals* simParams_ = info_->getSimParams();
57    
58     if (simParams_->haveSkinThickness()) {
59     skinThickness_ = simParams_->getSkinThickness();
60     } else {
61     skinThickness_ = 1.0;
62     sprintf(painCave.errMsg,
63     "ForceDecomposition: No value was set for the skinThickness.\n"
64     "\tOpenMD will use a default value of %f Angstroms\n"
65     "\tfor this simulation\n", skinThickness_);
66     painCave.severity = OPENMD_INFO;
67     painCave.isFatal = 0;
68     simError();
69     }
70    
71     // cellOffsets are the partial space for the cell lists used in
72     // constructing the neighbor lists
73     cellOffsets_.clear();
74     cellOffsets_.push_back( Vector3i(0, 0, 0) );
75     cellOffsets_.push_back( Vector3i(1, 0, 0) );
76     cellOffsets_.push_back( Vector3i(1, 1, 0) );
77     cellOffsets_.push_back( Vector3i(0, 1, 0) );
78     cellOffsets_.push_back( Vector3i(-1,1, 0) );
79     cellOffsets_.push_back( Vector3i(0, 0, 1) );
80     cellOffsets_.push_back( Vector3i(1, 0, 1) );
81     cellOffsets_.push_back( Vector3i(1, 1, 1) );
82     cellOffsets_.push_back( Vector3i(0, 1, 1) );
83     cellOffsets_.push_back( Vector3i(-1,1, 1) );
84     cellOffsets_.push_back( Vector3i(-1,0, 1) );
85     cellOffsets_.push_back( Vector3i(-1,-1,1) );
86     cellOffsets_.push_back( Vector3i(0, -1,1) );
87     cellOffsets_.push_back( Vector3i(1, -1,1) );
88 chuckv 1595
89     // cellAllOffsets are the full space for the cell lists used in
90     // constructing the neighbor lists
91     cellAllOffsets_.clear();
92     cellAllOffsets_.push_back( Vector3i(0, 0, 0) );
93    
94     cellAllOffsets_.push_back( Vector3i(-1, -1, -1) );
95     cellAllOffsets_.push_back( Vector3i(0, -1, -1) );
96     cellAllOffsets_.push_back( Vector3i(1, -1, -1) );
97     cellAllOffsets_.push_back( Vector3i(-1, 0, -1) );
98     cellAllOffsets_.push_back( Vector3i(0, 0, -1) );
99     cellAllOffsets_.push_back( Vector3i(1, 0, -1) );
100     cellAllOffsets_.push_back( Vector3i(-1, 1, -1) );
101     cellAllOffsets_.push_back( Vector3i(0, 1, -1) );
102     cellAllOffsets_.push_back( Vector3i(1, 1, -1) );
103    
104     cellAllOffsets_.push_back( Vector3i(-1, -1, 0) );
105     cellAllOffsets_.push_back( Vector3i(0, -1, 0) );
106     cellAllOffsets_.push_back( Vector3i(1, -1, 0) );
107     cellAllOffsets_.push_back( Vector3i(-1, 0, 0) );
108     // cellAllOffsets_.push_back( Vector3i(0, 0, 0) );
109     cellAllOffsets_.push_back( Vector3i(1, 0, 0) );
110     cellAllOffsets_.push_back( Vector3i(-1, 1, 0) );
111     cellAllOffsets_.push_back( Vector3i(0, 1, 0) );
112     cellAllOffsets_.push_back( Vector3i(1, 1, 0) );
113    
114     cellAllOffsets_.push_back( Vector3i(-1, -1, 1) );
115     cellAllOffsets_.push_back( Vector3i(0, -1, 1) );
116     cellAllOffsets_.push_back( Vector3i(1, -1, 1) );
117     cellAllOffsets_.push_back( Vector3i(-1, 0, 1) );
118     cellAllOffsets_.push_back( Vector3i(0, 0, 1) );
119     cellAllOffsets_.push_back( Vector3i(1, 0, 1) );
120     cellAllOffsets_.push_back( Vector3i(-1, 1, 1) );
121     cellAllOffsets_.push_back( Vector3i(0, 1, 1) );
122     cellAllOffsets_.push_back( Vector3i(1, 1, 1) );
123 gezelter 1572 }
124    
125 gezelter 1582 void ForceDecomposition::fillSelfData(SelfData &sdat, int atom1) {
126 gezelter 1583
127 gezelter 1591 sdat.atype = atypesLocal[atom1];
128 gezelter 1583
129     // Still Missing skippedCharge
130 gezelter 1572 if (storageLayout_ & DataStorage::dslElectroFrame) {
131     sdat.eFrame = &(snap_->atomData.electroFrame[atom1]);
132     }
133    
134     if (storageLayout_ & DataStorage::dslTorque) {
135     sdat.t = &(snap_->atomData.torque[atom1]);
136     }
137    
138     if (storageLayout_ & DataStorage::dslDensity) {
139     sdat.rho = &(snap_->atomData.density[atom1]);
140     }
141    
142     if (storageLayout_ & DataStorage::dslFunctional) {
143     sdat.frho = &(snap_->atomData.functional[atom1]);
144     }
145    
146     if (storageLayout_ & DataStorage::dslFunctionalDerivative) {
147     sdat.dfrhodrho = &(snap_->atomData.functionalDerivative[atom1]);
148     }
149    
150 gezelter 1586 if (storageLayout_ & DataStorage::dslSkippedCharge) {
151     sdat.skippedCharge = &(snap_->atomData.skippedCharge[atom1]);
152     }
153    
154 gezelter 1575 if (storageLayout_ & DataStorage::dslParticlePot) {
155     sdat.particlePot = &(snap_->atomData.particlePot[atom1]);
156     }
157 gezelter 1572 }
158    
159     bool ForceDecomposition::checkNeighborList() {
160    
161     int nGroups = snap_->cgData.position.size();
162    
163     // if we have changed the group identities or haven't set up the
164     // saved positions we automatically will need a neighbor list update:
165    
166 gezelter 1587 if ( saved_CG_positions_.size() != nGroups ) return true;
167 gezelter 1572
168     RealType dispmax = 0.0;
169     Vector3d disp;
170    
171     for (int i = 0; i < nGroups; i++) {
172     disp = snap_->cgData.position[i] - saved_CG_positions_[i];
173     for (int j = 0; j < 3; j++)
174     dispmax = max( abs(disp[j]), dispmax);
175     }
176    
177     #ifdef IS_MPI
178     MPI::COMM_WORLD.Allreduce(&dispmax, &dispmax, 1, MPI::REALTYPE, MPI::MAX);
179     #endif
180    
181     // a conservative test of list skin crossings
182     dispmax = 2.0 * sqrt (3.0 * dispmax * dispmax);
183    
184 gezelter 1582
185 gezelter 1587 if (dispmax > skinThickness_)
186 gezelter 1582 return (dispmax > skinThickness_);
187 gezelter 1587
188 gezelter 1582 return false;
189 gezelter 1572 }
190     }

Properties

Name Value
svn:eol-style native