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