OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
equationofstate.cpp
1/*
2 * Copyright (c) 2004-present, The University of Notre Dame. All rights
3 * reserved.
4 *
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 *
8 * 1. Redistributions of source code must retain the above copyright notice,
9 * this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright notice,
12 * this list of conditions and the following disclaimer in the documentation
13 * and/or other materials provided with the distribution.
14 *
15 * 3. Neither the name of the copyright holder nor the names of its
16 * contributors may be used to endorse or promote products derived from
17 * this software without specific prior written permission.
18 *
19 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
20 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
21 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
22 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
23 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
24 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
25 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
26 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
27 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
28 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
29 * POSSIBILITY OF SUCH DAMAGE.
30 *
31 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
32 * research, please cite the appropriate papers when you publish your
33 * work. Good starting points are:
34 *
35 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
36 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
37 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
38 * [4] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
39 * [5] Kuang & Gezelter, Mol. Phys., 110, 691-701 (2012).
40 * [6] Lamichhane, Gezelter & Newman, J. Chem. Phys. 141, 134109 (2014).
41 * [7] Lamichhane, Newman & Gezelter, J. Chem. Phys. 141, 134110 (2014).
42 * [8] Bhattarai, Newman & Gezelter, Phys. Rev. B 99, 094106 (2019).
43 */
44
45#include <cmath>
46#include <fstream>
47#include <iostream>
48#include <string>
49
51#include "brains/Register.hpp"
52#include "brains/SimCreator.hpp"
53#include "brains/SimInfo.hpp"
54#include "brains/Thermo.hpp"
56#include "flucq/FluctuatingChargeDamped.hpp"
57#include "io/DumpReader.hpp"
58#include "io/DumpWriter.hpp"
59#include "types/FluctuatingChargeAdapter.hpp"
60#include "utils/Constants.hpp"
61#include "utils/ProgressBar.hpp"
62#include "utils/simError.h"
63
64using namespace OpenMD;
65using namespace std;
66
67int main(int argc, char* argv[]) {
68 gengetopt_args_info args_info;
69 string omdFileName;
70 string outFileName;
71
72 // parse the command line option
73 if (cmdline_parser(argc, argv, &args_info) != 0) { exit(1); }
74
75 // get the omd file name and meta-data file name
76 if (args_info.input_given) {
77 omdFileName = args_info.input_arg;
78 } else {
79 strcpy(painCave.errMsg, "No input file name was specified.\n");
80 painCave.isFatal = 1;
81 simError();
82 }
83
84 if (args_info.output_given) {
85 outFileName = args_info.output_arg;
86 } else {
87 strcpy(painCave.errMsg, "No output file name was specified.\n");
88 painCave.isFatal = 1;
89 simError();
90 }
91
92 double start_affine = args_info.start_arg;
93 double end_affine = args_info.end_arg;
94 int number = args_info.number_arg;
95
96 RealType affine_step = (end_affine - start_affine) / (number + 1);
97
99
100 SimInfo::MoleculeIterator miter;
101 Molecule::IntegrableObjectIterator iiter;
102 Molecule::RigidBodyIterator rbIter;
103 Molecule* mol;
104 StuntDouble* sd;
105 StuntDouble* sdNew;
106 RigidBody* rb;
107 Mat3x3d oldHmat;
108 Mat3x3d newHmat;
109 Snapshot* oldSnap;
110 Snapshot* newSnap;
111 Vector3d oldPos;
112 Vector3d newPos;
113 AtomType* atype;
114
115 // parse omd file and set up the system
116 SimCreator oldCreator;
117 SimInfo* oldInfo = oldCreator.createSim(omdFileName);
118 oldSnap = oldInfo->getSnapshotManager()->getCurrentSnapshot();
119 oldHmat = oldSnap->getHmat();
120
121 // ProgressBarPtr progressBar {nullptr};
122 // progressBar = std::make_unique<ProgressBar>();
123
124 ofstream eos;
125 eos.open(outFileName.c_str());
126
127 RealType current_affine = start_affine;
128 int countStep = 0;
129 std::cout << "Calculation for EOS started." << std::endl;
130 while (current_affine <= end_affine) {
131 // progressBar->setStatus(countStep,number);
132 // progressBar->update();
133 ++countStep;
134 RealType scaling = std::cbrt(current_affine);
135 Mat3x3d scaleMatrix = Mat3x3d(0.0);
136 scaleMatrix(0, 0) = scaling;
137 scaleMatrix(1, 1) = scaling;
138 scaleMatrix(2, 2) = scaling;
139
140 SimInfo* newInfo = oldCreator.createSim(omdFileName);
141 newSnap = newInfo->getSnapshotManager()->getCurrentSnapshot();
142
143 newHmat = scaleMatrix * oldHmat;
144 newSnap->setHmat(newHmat);
145
146 int newIndex = 0;
147 for (mol = oldInfo->beginMolecule(miter); mol != NULL;
148 mol = oldInfo->nextMolecule(miter)) {
149 for (sd = mol->beginIntegrableObject(iiter); sd != NULL;
150 sd = mol->nextIntegrableObject(iiter)) {
151 oldPos = sd->getPos();
152 oldSnap->wrapVector(oldPos);
153 newPos = scaleMatrix * oldPos;
154 sdNew = newInfo->getIOIndexToIntegrableObject(newIndex);
155 sdNew->setPos(newPos);
156 sdNew->setVel(sd->getVel());
157 if (sd->isAtom()) {
158 atype = static_cast<Atom*>(sd)->getAtomType();
160 if (fqa.isFluctuatingCharge()) {
161 RealType charge = sd->getFlucQPos();
162 sdNew->setFlucQPos(charge);
163 RealType cv = sd->getFlucQVel();
164 sdNew->setFlucQVel(cv);
165 }
166 }
167 }
168
169 newIndex++;
170 }
171
172 for (mol = newInfo->beginMolecule(miter); mol != NULL;
173 mol = newInfo->nextMolecule(miter)) {
174 // change the positions of atoms which belong to the rigidbodies
175 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
176 rb = mol->nextRigidBody(rbIter)) {
177 rb->updateAtoms();
178 rb->updateAtomVel();
179 }
180 }
181
182 ForceManager* fman = new ForceManager(newInfo);
183 fman->initialize();
184
186 flucQ->setForceManager(fman);
187 flucQ->initialize();
188
189 fman->calcForces();
190 Thermo thermo(newInfo);
191 RealType totalEnergy(0);
192 totalEnergy = thermo.getTotalEnergy();
193 eos << current_affine << "\t" << totalEnergy << "\n";
194
195 std::cout << countStep << " data generated." << std::endl;
196
197 current_affine += affine_step;
198 }
199 eos.close();
200}
AtomType is what OpenMD looks to for unchanging data about an atom.
Definition AtomType.hpp:66
abstract class for propagating fluctuating charge variables
ForceManager is responsible for calculating both the short range (bonded) interactions and long range...
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...
Definition SimInfo.hpp:93
Molecule * beginMolecule(MoleculeIterator &i)
Returns the first molecule in this SimInfo and intialize the iterator.
Definition SimInfo.cpp:240
Molecule * nextMolecule(MoleculeIterator &i)
Returns the next avaliable Molecule based on the iterator.
Definition SimInfo.cpp:245
SnapshotManager * getSnapshotManager()
Returns the snapshot manager.
Definition SimInfo.hpp:248
StuntDouble * getIOIndexToIntegrableObject(int index)
return an integral objects by its global index.
Definition SimInfo.cpp:1027
The Snapshot class is a repository storing dynamic data during a Simulation.
Definition Snapshot.hpp:147
Mat3x3d getHmat()
Returns the H-Matrix.
Definition Snapshot.cpp:214
void setHmat(const Mat3x3d &m)
Sets the H-Matrix.
Definition Snapshot.cpp:217
void wrapVector(Vector3d &v)
Wrapping the vector according to periodic boundary condition.
Definition Snapshot.cpp:337
Snapshot * getCurrentSnapshot()
Returns the pointer of current snapshot.
"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.
Vector3d getVel()
Returns the current velocity 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.
bool isAtom()
Tests if this stuntDouble is an atom.
RealType getFlucQVel()
Returns the current charge velocity of this stuntDouble.
The header file for the command line option parser generated by GNU Gengetopt version 2....
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
void registerAll()
register force fields, integrators and optimizers
Definition Register.cpp:140
Where the command line options are stored.
unsigned int output_given
Whether output was given.
char * output_arg
output file name.
int number_arg
number of data points (default='50').
char * input_arg
input dump file.
double start_arg
starting affine scale (default='0.8').
unsigned int input_given
Whether input was given.
double end_arg
ending affine scale (default='1.2').