Simpatico  v1.10
mcMd/analyzers/system/StructureFactor.cpp
1 /*
2 * Simpatico - Simulation Package for Polymeric and Molecular Liquids
3 *
4 * Copyright 2010, The Regents of the University of Minnesota
5 * Distributed under the terms of the GNU General Public License.
6 */
7 
8 #include "StructureFactor.h"
9 #include <mcMd/simulation/Simulation.h>
10 #include <mcMd/simulation/McMd_mpi.h>
11 #include <mcMd/chemistry/Molecule.h>
12 #include <mcMd/chemistry/Atom.h>
13 #include <simp/boundary/Boundary.h>
14 #include <util/math/Constants.h>
15 #include <util/space/Dimension.h>
16 #include <util/misc/FileMaster.h>
17 #include <util/archives/Serializable_includes.h>
18 #include <util/misc/ioUtil.h>
19 #include <util/format/Int.h>
20 #include <util/format/Dbl.h>
21 
22 namespace McMd
23 {
24 
25  using namespace Util;
26  using namespace Simp;
27 
28  /*
29  * Constructor.
30  */
32  : SystemAnalyzer<System>(system),
33  isFirstStep_(true),
34  isInitialized_(false)
35  { setClassName("StructureFactor"); }
36 
37  /*
38  * Destructor.
39  */
41  {}
42 
43  /*
44  * Read parameters from file, and allocate memory.
45  */
46  void StructureFactor::readParameters(std::istream& in)
47  {
48  readInterval(in);
50  read<int>(in, "nMode", nMode_);
52  modes_.allocate(nMode_, nAtomType_);
53  readDMatrix<double>(in, "modes", modes_, nMode_, nAtomType_);
54  read<int>(in, "nWave", nWave_);
55  waveIntVectors_.allocate(nWave_);
56  readDArray<IntVector>(in, "waveIntVectors", waveIntVectors_, nWave_);
57 
58  waveVectors_.allocate(nWave_);
59  fourierModes_.allocate(nWave_, nMode_);
61 
62  isInitialized_ = true;
63  }
64 
65  /*
66  * Load state from an archive.
67  */
69  {
71  ar & nAtomType_;
72  loadParameter<int>(ar, "nMode", nMode_);
73  loadDMatrix<double>(ar, "modes", modes_, nMode_, nAtomType_);
74  loadParameter<int>(ar, "nWave", nWave_);
75  loadDArray<IntVector>(ar, "waveIntVectors", waveIntVectors_, nWave_);
76  ar & structureFactors_;
77  ar & nSample_;
78 
79  // Validate
80  if (nAtomType_ != system().simulation().nAtomType()) {
81  UTIL_THROW("Inconsistent values of nAtomType_");
82  }
83  if (modes_.capacity1() != nMode_) {
84  UTIL_THROW("Inconsistent capacity1 for modes array");
85  }
86  if (modes_.capacity2() != nAtomType_) {
87  UTIL_THROW("Inconsistent capacity2 for modes array");
88  }
89  if (waveIntVectors_.capacity() != nWave_) {
90  UTIL_THROW("Inconsistent capacity for waveIntVector");
91  }
92 
93  // Allocate temporary data structures
94  waveVectors_.allocate(nWave_);
95  fourierModes_.allocate(nWave_, nMode_);
96 
97  isInitialized_ = true;
98  }
99 
100  /*
101  * Save state to archive.
102  */
104  { ar & *this; }
105 
106  /*
107  * Clear accumulators.
108  */
110  {
111  if (!isInitialized_) {
112  UTIL_THROW("Error: object is not initialized");
113  }
114  assert (nWave_ > 0);
115  assert (nMode_ > 0);
116 
117  // Clear accumulators
118  int i, j;
119  for (i = 0; i < nWave_; ++i) {
120  for (j = 0; j < nMode_; ++j) {
121  structureFactors_(i, j) = 0.0;
122  }
123  }
124  nSample_ = 0;
125  }
126 
127  /*
128  * Increment structure factors for all wavevectors and modes.
129  */
130  void StructureFactor::sample(long iStep)
131  {
132  if (isAtInterval(iStep)) {
133 
134  std::ios_base::openmode mode = std::ios_base::out;
135  if (!isFirstStep_) {
136  mode = std::ios_base::out | std::ios_base::app;
137  }
139  outputFile_, mode);
140 
141  //fileMaster().openOutputFile(outputFileName("_max.dat"),
142  // outputFile_, !isFirstStep_);
143  isFirstStep_ = false;
144 
145  Vector position;
146  std::complex<double> expFactor;
147  double product;
150  int nSpecies, iSpecies, typeId, i, j;
151 
152  makeWaveVectors();
153 
154  // Set all Fourier modes to zero
155  for (i = 0; i < nWave_; ++i) {
156  for (j = 0; j < nMode_; ++j) {
157  fourierModes_(i, j) = std::complex<double>(0.0, 0.0);
158  }
159  }
160 
161  // Loop over all atoms
162  nSpecies = system().simulation().nSpecies();
163  for (iSpecies = 0; iSpecies < nSpecies; ++iSpecies) {
164  system().begin(iSpecies, molIter);
165  for ( ; molIter.notEnd(); ++molIter) {
166  molIter->begin(atomIter);
167  for ( ; atomIter.notEnd(); ++atomIter) {
168  position = atomIter->position();
169  typeId = atomIter->typeId();
170 
171  // Loop over wavevectors
172  for (i = 0; i < nWave_; ++i) {
173 
174  product = position.dot(waveVectors_[i]);
175  expFactor = exp( product*Constants::Im );
176  for (j = 0; j < nMode_; ++j) {
177  fourierModes_(i, j) += modes_(j, typeId)*expFactor;
178  }
179 
180  }
181 
182  }
183  }
184 
185  }
186 
187  // Increment structure factors
188  double volume = system().boundary().volume();
189  double norm;
190  for (j = 0; j < nMode_; ++j) {
191  double maxValue = 0.0;
192  double maxQ = 0.0;
193  IntVector maxIntVector;
194  for (i = 0; i < nWave_; ++i) {
195  norm = std::norm(fourierModes_(i, j));
196  if (double(norm/volume) >= maxValue) {
197  maxValue = double(norm/volume);
198  maxIntVector = waveIntVectors_[i];
199  maxQ = waveVectors_[i].abs();
200  }
201  structureFactors_(i, j) += norm/volume;
202  }
203 
204  // Output current maximum S(q)
205  outputFile_ << maxIntVector;
206  outputFile_ << Dbl(maxQ, 20, 8);
207  outputFile_ << Dbl(maxValue, 20, 8);
208  outputFile_ << std::endl;
209  }
210 
211  ++nSample_;
212 
213  outputFile_ << std::endl;
214  outputFile_.close();
215  }
216 
217  }
218 
219  /*
220  * Calculate floating point wavevectors, using current boundary.
221  */
223  {
224  Boundary* boundaryPtr = &system().boundary();
225  Vector dWave;
226  int i, j;
227  for (i = 0; i < nWave_; ++i) {
229  for (j = 0; j < Dimension; ++j) {
230  dWave = boundaryPtr->reciprocalBasisVector(j);
231  dWave *= waveIntVectors_[i][j];
232  waveVectors_[i] += dWave;
233  }
234  }
235  }
236 
237  /*
238  * Output final results to output file.
239  */
241  {
242  // Echo parameters to a log file
245  outputFile_.close();
246 
247  // Output structure factors to one file
249  double value;
250  int i, j, k;
251  for (i = 0; i < nWave_; ++i) {
252  for (k = 0; k < Dimension; ++k) {
253  outputFile_ << Int(waveIntVectors_[i][k], 5);
254  }
255  outputFile_ << Dbl(waveVectors_[i].abs(), 20, 8);
256  for (j = 0; j < nMode_; ++j) {
257  value = structureFactors_(i, j)/double(nSample_);
258  outputFile_ << Dbl(value, 18, 8);
259  }
260  outputFile_ << std::endl;
261  }
262  outputFile_.close();
263 
264  #if 0
265  // Output each structure factor to a separate file
266  std::string suffix;
267  int typeId;
268  for (j = 0; j < nMode_; ++j) {
269 
270  // Construct file suffix for this structure factor
271  suffix = std::string(".");
272  suffix += toString(j);
273  suffix += std::string(".dat");
274 
276 
277  // Loop over waves to output structure factor
278  for (i = 0; i < nWave_; ++i) {
279  for (k = 0; k < Dimension; ++k) {
280  outputFile_ << Int(waveIntVectors_[i][k], 5);
281  }
282  outputFile_ << Dbl(waveVectors_[i].abs(), 20, 8);
283  value = structureFactors_(i, j)/double(nSample_);
284  outputFile_ << Dbl(value, 20, 8);
285  outputFile_ << std::endl;
286  }
287 
288  outputFile_.close();
289  }
290  #endif
291 
292  }
293 
294 }
const int Dimension
Dimensionality of space.
Definition: Dimension.h:19
int nSample_
Number of samples thus far.
std::string toString(int n)
Return string representation of an integer.
Definition: ioUtil.cpp:52
A Vector is a Cartesian vector.
Definition: Vector.h:75
float product(float a, float b)
Product for float Data.
Definition: product.h:22
DArray< Vector > waveVectors_
Array of floating point wave vectors (temporary).
void allocate(int capacity1, int capacity2)
Allocate memory for a matrix.
Definition: DMatrix.h:170
void begin(int speciesId, MoleculeIterator &iterator)
Initialize an iterator for molecules of one species in this System.
Definition: System.h:1147
double volume() const
Return unit cell volume.
double dot(const Vector &v) const
Return dot product of this vector and vector v.
Definition: Vector.h:632
An orthorhombic periodic unit cell.
void openOutputFile(const std::string &filename, std::ofstream &out, std::ios_base::openmode mode=std::ios_base::out) const
Open an output file.
Definition: FileMaster.cpp:290
static const Vector Zero
Zero Vector = {0.0, 0.0, 0.0}.
Definition: Vector.h:369
virtual void loadParameters(Serializable::IArchive &ar)
Load parameters from archive.
DMatrix< std::complex< double > > fourierModes_
Instantaneous Fourier amplitudes (temporary)
virtual void output()
Output results to predefined output file.
A set of interacting Molecules enclosed by a Boundary.
Definition: System.h:115
System & system()
Return reference to parent system.
Wrapper for a double precision number, for formatted ostream output.
Definition: Dbl.h:39
Forward iterator for a PArray.
Classes used by all simpatico molecular simulations.
void readOutputFileName(std::istream &in)
Read outputFileName from file.
DMatrix< double > modes_
Array of mode vectors.
Forward const iterator for an Array or a C array.
Saving / output archive for binary ostream.
#define UTIL_THROW(msg)
Macro for throwing an Exception, reporting function, file and line number.
Definition: global.h:51
Simulation & simulation() const
Get the parent Simulation by reference.
Definition: System.h:1055
virtual void writeParam(std::ostream &out)
Write all parameters to an output stream.
const Vector & reciprocalBasisVector(int i) const
Return reciprocal lattice basis vector i.
void readInterval(std::istream &in)
Read interval from file, with error checking.
int nAtomType_
Number of atom types, copied from Simulation::nAtomType().
DMatrix< double > structureFactors_
Structure factor accumulators.
Utility classes for scientific computation.
Definition: accumulators.mod:1
DArray< IntVector > waveIntVectors_
Array of Miller index IntVectors for wavevectors.
Wrapper for an int, for formatted ostream output.
Definition: Int.h:36
virtual void save(Serializable::OArchive &ar)
Save state to archive.
std::ofstream outputFile_
Output file stream.
Template for Analyzer associated with one System.
bool notEnd() const
Is the current pointer not at the end of the array?
Boundary & boundary() const
Get the Boundary by reference.
Definition: System.h:1064
bool isFirstStep_
Is this the first step?
Saving archive for binary istream.
virtual void loadParameters(Serializable::IArchive &ar)
Load state from an archive.
StructureFactor(System &system)
Constructor.
Single-processor Monte Carlo (MC) and molecular dynamics (MD).
int nSpecies() const
Get the number of Species in this Simulation.
virtual void readParameters(std::istream &in)
Read parameters from file.
void setClassName(const char *className)
Set class name string.
An IntVector is an integer Cartesian vector.
Definition: IntVector.h:73
FileMaster & fileMaster()
Get the FileMaster by reference.
bool isAtInterval(long counter) const
Return true iff counter is a multiple of the interval.
const std::string & outputFileName() const
Return outputFileName string.
virtual void sample(long iStep)
Add particles to StructureFactor accumulators.
int nAtomType() const
Get the number of atom types.
bool notEnd() const
Is this not the end of the array?
static const std::complex< double > Im
Square root of -1.
Definition: Constants.h:40
virtual void setup()
Clear accumulators.