50#include "brains/Register.hpp"
53#include "brains/Thermo.hpp"
55#include "io/DumpWriter.hpp"
58#include "types/FluctuatingChargeAdapter.hpp"
59#include "utils/Constants.hpp"
60#include "utils/simError.h"
66void createMdFile(
const std::string& oldMdFileName,
67 const std::string& newMdFileName, std::vector<int> nMol);
69int main(
int argc,
char* argv[]) {
75 if (cmdline_parser(argc, argv, &args_info) != 0) { exit(1); }
81 strcpy(painCave.errMsg,
"No input file name was specified.\n");
89 strcpy(painCave.errMsg,
"No output file name was specified.\n");
95 double phi = args_info.
rotatePhi_arg * (Constants::PI / 180.0);
97 double psi = args_info.
rotatePsi_arg * (Constants::PI / 180.0);
107 repeatD(0, 0) = repeat.
x();
108 repeatD(1, 1) = repeat.
y();
109 repeatD(2, 2) = repeat.
z();
119 Globals* simParams = oldInfo->getSimParams();
120 std::vector<Component*> components = simParams->getComponents();
121 std::vector<int> nMol;
122 for (vector<Component*>::iterator i = components.begin();
123 i != components.end(); ++i) {
124 int nMolOld = (*i)->getNMol();
125 int nMolNew = nMolOld * repeat.
x() * repeat.
y() * repeat.
z();
126 nMol.push_back(nMolNew);
129 createMdFile(dumpFileName, outFileName, nMol);
138 if (writer == NULL) {
139 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
140 "error in creating DumpWriter");
141 painCave.isFatal = 1;
145 SimInfo::MoleculeIterator miter;
146 Molecule::IntegrableObjectIterator iiter;
147 Molecule::RigidBodyIterator rbIter;
165 for (
int i = 0; i < nframes; i++) {
166 cerr <<
"frame = " << i <<
"\n";
167 dumpReader->readFrame(i);
172 newSnap->setTime(oldSnap->getTime());
175 rotHmat = rotMatrix * oldHmat;
176 newHmat = repeatD * rotHmat;
179 newSnap->setThermostat(oldSnap->getThermostat());
180 newSnap->setBarostat(oldSnap->getBarostat());
187 COM = thermo.getCom();
197 for (sd = mol->beginIntegrableObject(iiter); sd != NULL;
198 sd = mol->nextIntegrableObject(iiter)) {
202 relPos = sd->
getPos() - molCOM;
203 oldPos = molCOM - COM;
214 oldPos = sd->
getPos() - COM;
231 for (
int ii = 0; ii < repeat.
x(); ii++) {
232 for (
int jj = 0; jj < repeat.
y(); jj++) {
233 for (
int kk = 0; kk < repeat.
z(); kk++) {
237 for (sd = mol->beginIntegrableObject(iiter); sd != NULL;
238 sd = mol->nextIntegrableObject(iiter)) {
242 newPos = rotMatrix * oldPos + trans * oldHmat + translate;
250 bodyRotMat = bodyRotMat * rotMatrix.
inverse();
251 sdNew->
setA(bodyRotMat);
252 sdNew->
setJ(rotMatrix * sd->
getJ());
256 atype =
static_cast<Atom*
>(sd)->getAtomType();
258 if (fqa.isFluctuatingCharge()) {
277 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
278 rb = mol->nextRigidBody(rbIter)) {
291void createMdFile(
const std::string& oldMdFileName,
292 const std::string& newMdFileName, std::vector<int> nMol) {
295 const int MAXLEN = 65535;
300 oldMdFile.open(oldMdFileName.c_str());
301 newMdFile.open(newMdFileName.c_str());
303 oldMdFile.getline(buffer, MAXLEN);
306 while (!oldMdFile.eof()) {
308 if (strstr(buffer,
"nMol") != NULL) {
309 if (i < nMol.size()) {
310 snprintf(buffer, MAXLEN,
"\tnMol = %i;", nMol.at(i));
311 newMdFile << buffer << std::endl;
315 newMdFile << buffer << std::endl;
317 oldMdFile.getline(buffer, MAXLEN);
323 if (i != nMol.size()) {
324 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
325 "Couldn't replace the correct number of nMol\n"
326 "\tstatements in component blocks.");
327 painCave.isFatal = 1;
AtomType is what OpenMD looks to for unchanging data about an atom.
int getNFrames()
Returns the number of frames in the dump file.
Vector3d getCom()
Returns the current center of mass position of this molecule.
void updateAtoms()
update the positions of atoms belong to this rigidbody
The only responsibility of SimCreator is to parse the meta-data file and create a SimInfo instance ba...
SimInfo * createSim(const std::string &mdFileName, bool loadInitCoords=true)
Setup Simulation.
One of the heavy-weight classes of OpenMD, SimInfo maintains objects and variables relating to the cu...
Molecule * beginMolecule(MoleculeIterator &i)
Returns the first molecule in this SimInfo and intialize the iterator.
Molecule * nextMolecule(MoleculeIterator &i)
Returns the next avaliable Molecule based on the iterator.
SnapshotManager * getSnapshotManager()
Returns the snapshot manager.
StuntDouble * getIOIndexToIntegrableObject(int index)
return an integral objects by its global index.
The Snapshot class is a repository storing dynamic data during a Simulation.
Mat3x3d getHmat()
Returns the H-Matrix.
void setHmat(const Mat3x3d &m)
Sets the H-Matrix.
void setID(int id)
Sets the id of this Snapshot.
void wrapVector(Vector3d &v)
Wrapping the vector according to periodic boundary condition.
int getID()
Returns the id of this Snapshot.
Snapshot * getCurrentSnapshot()
Returns the pointer of current snapshot.
SquareMatrix3< Real > inverse() const
Sets the value of this matrix to the inverse of itself.
void setupRotMat(const Vector3< Real > &eulerAngles)
Sets this matrix to a rotation matrix by three euler angles @ param euler.
"Don't move, or you're dead! Stand up! Captain, we've got them!"
void setFlucQVel(RealType cvel)
Sets the current charge velocity of this stuntDouble.
RotMat3x3d getA()
Returns the current rotation matrix of this stuntDouble.
Vector3d getVel()
Returns the current velocity of this stuntDouble.
virtual void setA(const RotMat3x3d &a)
Sets the current rotation matrix of this stuntDouble.
void setFlucQPos(RealType charge)
Sets the current fluctuating charge of this stuntDouble.
Vector3d getPos()
Returns the current position of this stuntDouble.
RealType getFlucQPos()
Returns the current fluctuating charge of this stuntDouble.
void setPos(const Vector3d &pos)
Sets the current position of this stuntDouble.
void setVel(const Vector3d &vel)
Sets the current velocity of this stuntDouble.
Vector3d getJ()
Returns the current angular momentum of this stuntDouble (body -fixed).
bool isAtom()
Tests if this stuntDouble is an atom.
bool isDirectional()
Tests if this stuntDouble is a directional one.
void setJ(const Vector3d &angMom)
Sets the current angular momentum of this stuntDouble (body-fixed).
RealType getFlucQVel()
Returns the current charge velocity of this stuntDouble.
Real & z()
Returns reference of the third element of Vector3.
Real & x()
Returns reference of the first element of Vector3.
Real & y()
Returns reference of the second element of Vector3.
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
The header file for the command line option parser generated by GNU Gengetopt version 2....
Where the command line options are stored.
double translateY_arg
translate all y coordinates by some amount (default='0.0').
int repeatY_arg
The number of images to repeat in the y direction (default='0').
unsigned int output_given
Whether output was given.
double translateZ_arg
translate all z coordinates by some amount (default='0.0').
int noCOM_flag
do not use Center of Mass as origin of the box (default=off).
int repeatX_arg
The number of images to repeat in the x direction (default='0').
double rotateTheta_arg
rotate all coordinates Euler angle Theta (default='0.0').
int repairMolecules_arg
rewrap molecules around the molecular center of mass (default='1').
char * output_arg
output file name.
char * input_arg
input dump file.
int repeatZ_arg
The number of images to repeat in the z direction (default='0').
int noWrap_flag
do not rewrap coordinates into the box (default=off).
double rotatePhi_arg
rotate all coordinates Euler angle Phi (default='0.0').
unsigned int input_given
Whether input was given.
double rotatePsi_arg
rotate all coordinates Euler angle Psi (default='0.0').
double translateX_arg
translate all x coordinates by some amount (default='0.0').