--- trunk/src/brains/SimInfo.cpp 2006/05/17 21:51:42 963 +++ trunk/src/brains/SimInfo.cpp 2008/04/25 15:14:47 1241 @@ -53,11 +53,13 @@ #include "brains/SimInfo.hpp" #include "math/Vector3.hpp" #include "primitives/Molecule.hpp" +#include "primitives/StuntDouble.hpp" #include "UseTheForce/fCutoffPolicy.h" #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h" #include "UseTheForce/DarkSide/fElectrostaticScreeningMethod.h" #include "UseTheForce/DarkSide/fSwitchingFunctionType.h" #include "UseTheForce/doForces_interface.h" +#include "UseTheForce/DarkSide/neighborLists_interface.h" #include "UseTheForce/DarkSide/electrostatic_interface.h" #include "UseTheForce/DarkSide/switcheroo_interface.h" #include "utils/MemoryUtils.hpp" @@ -66,6 +68,7 @@ #include "io/ForceFieldOptions.hpp" #include "UseTheForce/ForceField.hpp" + #ifdef IS_MPI #include "UseTheForce/mpiComponentPlan.h" #include "UseTheForce/DarkSide/simParallel_interface.h" @@ -89,7 +92,8 @@ namespace oopse { nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0), nAtoms_(0), nBonds_(0), nBends_(0), nTorsions_(0), nRigidBodies_(0), nIntegrableObjects_(0), nCutoffGroups_(0), nConstraints_(0), - sman_(NULL), fortranInitialized_(false) { + sman_(NULL), fortranInitialized_(false), calcBoxDipole_(false), + useAtomicVirial_(true) { MoleculeStamp* molStamp; int nMolWithSameStamp; @@ -152,11 +156,7 @@ namespace oopse { + nGlobalRigidBodies_; nGlobalMols_ = molStampIds_.size(); - -#ifdef IS_MPI molToProcMap_.resize(nGlobalMols_); -#endif - } SimInfo::~SimInfo() { @@ -600,8 +600,11 @@ namespace oopse { /** @deprecate */ int isError = 0; + setupCutoff(); + setupElectrostaticSummationMethod( isError ); setupSwitchingFunction(); + setupAccumulateBoxDipole(); if(isError){ sprintf( painCave.errMsg, @@ -609,9 +612,6 @@ namespace oopse { painCave.isFatal = 1; simError(); } - - - setupCutoff(); calcNdf(); calcNdfRaw(); @@ -661,25 +661,35 @@ namespace oopse { int usePBC = simParams_->getUsePeriodicBoundaryConditions(); int useRF; int useSF; + int useSP; + int useBoxDipole; + std::string myMethod; // set the useRF logical useRF = 0; useSF = 0; + useSP = 0; if (simParams_->haveElectrostaticSummationMethod()) { std::string myMethod = simParams_->getElectrostaticSummationMethod(); toUpper(myMethod); - if (myMethod == "REACTION_FIELD") { - useRF=1; - } else { - if (myMethod == "SHIFTED_FORCE") { - useSF = 1; - } + if (myMethod == "REACTION_FIELD"){ + useRF = 1; + } else if (myMethod == "SHIFTED_FORCE"){ + useSF = 1; + } else if (myMethod == "SHIFTED_POTENTIAL"){ + useSP = 1; } } + + if (simParams_->haveAccumulateBoxDipole()) + if (simParams_->getAccumulateBoxDipole()) + useBoxDipole = 1; + useAtomicVirial_ = simParams_->getUseAtomicVirial(); + //loop over all of the atom types for (i = atomTypes.begin(); i != atomTypes.end(); ++i) { useLennardJones |= (*i)->isLennardJones(); @@ -749,8 +759,17 @@ namespace oopse { MPI_Allreduce(&temp, &useRF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); temp = useSF; - MPI_Allreduce(&temp, &useSF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); + MPI_Allreduce(&temp, &useSF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); + temp = useSP; + MPI_Allreduce(&temp, &useSP, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); + + temp = useBoxDipole; + MPI_Allreduce(&temp, &useBoxDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); + + temp = useAtomicVirial_; + MPI_Allreduce(&temp, &useAtomicVirial_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD); + #endif fInfo_.SIM_uses_PBC = usePBC; @@ -768,21 +787,9 @@ namespace oopse { fInfo_.SIM_uses_FLARB = useFLARB; fInfo_.SIM_uses_RF = useRF; fInfo_.SIM_uses_SF = useSF; - - if( myMethod == "REACTION_FIELD") { - - if (simParams_->haveDielectric()) { - fInfo_.dielect = simParams_->getDielectric(); - } else { - sprintf(painCave.errMsg, - "SimSetup Error: No Dielectric constant was set.\n" - "\tYou are trying to use Reaction Field without" - "\tsetting a dielectric constant!\n"); - painCave.isFatal = 1; - simError(); - } - } - + fInfo_.SIM_uses_SP = useSP; + fInfo_.SIM_uses_BoxDipole = useBoxDipole; + fInfo_.SIM_uses_AtomicVirial = useAtomicVirial_; } void SimInfo::setupFortranSim() { @@ -849,30 +856,38 @@ namespace oopse { int nGlobalExcludes = 0; int* globalExcludes = NULL; int* excludeList = exclude_.getExcludeList(); - setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0], &nExclude, excludeList , - &nGlobalExcludes, globalExcludes, &molMembershipArray[0], - &mfact[0], &nCutoffGroups_, &fortranGlobalGroupMembership[0], &isError); - + setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0], + &nExclude, excludeList , &nGlobalExcludes, globalExcludes, + &molMembershipArray[0], &mfact[0], &nCutoffGroups_, + &fortranGlobalGroupMembership[0], &isError); + if( isError ){ - + sprintf( painCave.errMsg, "There was an error setting the simulation information in fortran.\n" ); painCave.isFatal = 1; painCave.severity = OOPSE_ERROR; simError(); } - -#ifdef IS_MPI + + sprintf( checkPointMsg, "succesfully sent the simulation information to fortran.\n"); - MPIcheckPoint(); -#endif // is_mpi + + errorCheckPoint(); + + // Setup number of neighbors in neighbor list if present + if (simParams_->haveNeighborListNeighbors()) { + int nlistNeighbors = simParams_->getNeighborListNeighbors(); + setNeighbors(&nlistNeighbors); + } + + } -#ifdef IS_MPI void SimInfo::setupFortranParallel() { - +#ifdef IS_MPI //SimInfo is responsible for creating localToGlobalAtomIndex and localToGlobalGroupIndex std::vector localToGlobalAtomIndex(getNAtoms(), 0); std::vector localToGlobalCutoffGroupIndex; @@ -922,13 +937,11 @@ namespace oopse { } sprintf(checkPointMsg, " mpiRefresh successful.\n"); - MPIcheckPoint(); + errorCheckPoint(); - +#endif } -#endif - void SimInfo::setupCutoff() { ForceFieldOptions& forceFieldOptions_ = forceField_->getForceFieldOptions(); @@ -936,6 +949,10 @@ namespace oopse { // Check the cutoff policy int cp = TRADITIONAL_CUTOFF_POLICY; // Set to traditional by default + // Set LJ shifting bools to false + ljsp_ = false; + ljsf_ = false; + std::string myPolicy; if (forceFieldOptions_.haveCutoffPolicy()){ myPolicy = forceFieldOptions_.getCutoffPolicy(); @@ -999,8 +1016,18 @@ namespace oopse { simError(); } } - - notifyFortranCutoffs(&rcut_, &rsw_); + + if (simParams_->haveElectrostaticSummationMethod()) { + std::string myMethod = simParams_->getElectrostaticSummationMethod(); + toUpper(myMethod); + + if (myMethod == "SHIFTED_POTENTIAL") { + ljsp_ = true; + } else if (myMethod == "SHIFTED_FORCE") { + ljsf_ = true; + } + } + notifyFortranCutoffs(&rcut_, &rsw_, &ljsp_, &ljsf_); } else { @@ -1017,7 +1044,15 @@ namespace oopse { if (simParams_->haveElectrostaticSummationMethod()) { std::string myMethod = simParams_->getElectrostaticSummationMethod(); toUpper(myMethod); - if (myMethod == "SHIFTED_POTENTIAL" || myMethod == "SHIFTED_FORCE") { + + // For the time being, we're tethering the LJ shifted behavior to the + // electrostaticSummationMethod keyword options + if (myMethod == "SHIFTED_POTENTIAL") { + ljsp_ = true; + } else if (myMethod == "SHIFTED_FORCE") { + ljsf_ = true; + } + if (myMethod == "SHIFTED_POTENTIAL" || myMethod == "SHIFTED_FORCE") { if (simParams_->haveSwitchingRadius()){ sprintf(painCave.errMsg, "SimInfo Warning: A value was set for the switchingRadius\n" @@ -1040,7 +1075,9 @@ namespace oopse { simError(); rsw_ = 0.85 * rcut_; } - notifyFortranCutoffs(&rcut_, &rsw_); + + notifyFortranCutoffs(&rcut_, &rsw_, &ljsp_, &ljsf_); + } else { // We didn't set rcut explicitly, and we don't have electrostatic atoms, so // We'll punt and let fortran figure out the cutoffs later. @@ -1058,10 +1095,8 @@ namespace oopse { int sm = UNDAMPED; RealType alphaVal; RealType dielectric; - + errorOut = isError; - alphaVal = simParams_->getDampingAlpha(); - dielectric = simParams_->getDielectric(); if (simParams_->haveElectrostaticSummationMethod()) { std::string myMethod = simParams_->getElectrostaticSummationMethod(); @@ -1078,8 +1113,17 @@ namespace oopse { if (myMethod == "SHIFTED_FORCE") { esm = SHIFTED_FORCE; } else { - if (myMethod == "REACTION_FIELD") { + if (myMethod == "REACTION_FIELD") { esm = REACTION_FIELD; + dielectric = simParams_->getDielectric(); + if (!simParams_->haveDielectric()) { + // throw warning + sprintf( painCave.errMsg, + "SimInfo warning: dielectric was not specified in the input file\n\tfor the reaction field correction method.\n" + "\tA default value of %f will be used for the dielectric.\n", dielectric); + painCave.isFatal = 0; + simError(); + } } else { // throw error sprintf( painCave.errMsg, @@ -1106,13 +1150,22 @@ namespace oopse { if (myScreen == "DAMPED") { sm = DAMPED; if (!simParams_->haveDampingAlpha()) { - //throw error + // first set a cutoff dependent alpha value + // we assume alpha depends linearly with rcut from 0 to 20.5 ang + alphaVal = 0.5125 - rcut_* 0.025; + // for values rcut > 20.5, alpha is zero + if (alphaVal < 0) alphaVal = 0; + + // throw warning sprintf( painCave.errMsg, "SimInfo warning: dampingAlpha was not specified in the input file.\n" - "\tA default value of %f (1/ang) will be used.\n", alphaVal); + "\tA default value of %f (1/ang) will be used for the cutoff of\n\t%f (ang).\n", alphaVal, rcut_); painCave.isFatal = 0; simError(); + } else { + alphaVal = simParams_->getDampingAlpha(); } + } else { // throw error sprintf( painCave.errMsg, @@ -1158,6 +1211,17 @@ namespace oopse { // send switching function notification to switcheroo setFunctionType(&ft); + + } + + void SimInfo::setupAccumulateBoxDipole() { + + // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true + if ( simParams_->haveAccumulateBoxDipole() ) + if ( simParams_->getAccumulateBoxDipole() ) { + setAccumulateBoxDipole(); + calcBoxDipole_ = true; + } } @@ -1410,6 +1474,61 @@ namespace oopse { return angularMomentum; } - + StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) { + return IOIndexToIntegrableObject.at(index); + } + + void SimInfo::setIOIndexToIntegrableObject(const std::vector& v) { + IOIndexToIntegrableObject= v; + } + + /* Returns the Volume of the simulation based on a ellipsoid with semi-axes + based on the radius of gyration V=4/3*Pi*R_1*R_2*R_3 + where R_i are related to the principle inertia moments R_i = sqrt(C*I_i/N), this reduces to + V = 4/3*Pi*(C/N)^3/2*sqrt(det(I)). See S.E. Baltazar et. al. Comp. Mat. Sci. 37 (2006) 526-536. + */ + void SimInfo::getGyrationalVolume(RealType &volume){ + Mat3x3d intTensor; + RealType det; + Vector3d dummyAngMom; + RealType sysconstants; + RealType geomCnst; + + geomCnst = 3.0/2.0; + /* Get the inertial tensor and angular momentum for free*/ + getInertiaTensor(intTensor,dummyAngMom); + + det = intTensor.determinant(); + sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_; + volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(det); + return; + } + + void SimInfo::getGyrationalVolume(RealType &volume, RealType &detI){ + Mat3x3d intTensor; + Vector3d dummyAngMom; + RealType sysconstants; + RealType geomCnst; + + geomCnst = 3.0/2.0; + /* Get the inertial tensor and angular momentum for free*/ + getInertiaTensor(intTensor,dummyAngMom); + + detI = intTensor.determinant(); + sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_; + volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(detI); + return; + } +/* + void SimInfo::setStuntDoubleFromGlobalIndex(std::vector v) { + assert( v.size() == nAtoms_ + nRigidBodies_); + sdByGlobalIndex_ = v; + } + + StuntDouble* SimInfo::getStuntDoubleFromGlobalIndex(int index) { + //assert(index < nAtoms_ + nRigidBodies_); + return sdByGlobalIndex_.at(index); + } +*/ }//end namespace oopse