Simpatico  v1.10
LinearRouseAutoCorr.h
1 #ifndef MCMD_LINEAR_ROUSE_AUTO_CORR_H
2 #define MCMD_LINEAR_ROUSE_AUTO_CORR_H
3 
4 /*
5 * Simpatico - Simulation Package for Polymeric and Molecular Liquids
6 *
7 * Copyright 2010 - 2017, The Regents of the University of Minnesota
8 * Distributed under the terms of the GNU General Public License.
9 */
10 
11 #include <mcMd/analyzers/SystemAnalyzer.h> // base class template
12 #include <mcMd/simulation/System.h> // base class template parameter
13 #include <util/accumulators/AutoCorrArray.h> // member template
14 #include <util/space/Vector.h> // member template parameter
15 #include <util/containers/DArray.h> // member template
16 
17 #include <util/archives/serialize.h> // used in method template
18 #include <util/global.h> // used in method template
19 
20 namespace Simp {
21  class Species;
22 }
23 
24 namespace McMd
25 {
26 
27  using namespace Util;
28  using namespace Simp;
29 
35  class LinearRouseAutoCorr : public SystemAnalyzer<System>
36  {
37 
38  public:
39 
45  LinearRouseAutoCorr(System &system);
46 
50  virtual ~LinearRouseAutoCorr();
51 
57  virtual void readParameters(std::istream& in);
58 
64  virtual void loadParameters(Serializable::IArchive &ar);
65 
71  virtual void save(Serializable::OArchive &ar);
72 
79  template <class Archive>
80  void serialize(Archive& ar, const unsigned int version);
81 
85  virtual void setup();
86 
92  virtual void sample(long iStep);
93 
97  virtual void output();
98 
99  private:
100 
102  std::ofstream outputFile_;
103 
105  AutoCorrArray<Vector, double> accumulator_;
106 
108  DArray<Vector> data_;
109 
111  DArray<double> projector_;
112 
114  Species* speciesPtr_;
115 
117  int speciesId_;
118 
120  int nMolecule_;
121 
123  int nAtom_;
124 
126  int p_;
127 
129  int capacity_;
130 
132  bool isInitialized_;
133 
134  };
135 
136  /*
137  * Serialize to/from an archive.
138  */
139  template <class Archive>
140  void LinearRouseAutoCorr::serialize(Archive& ar, const unsigned int version)
141  {
142  Analyzer::serialize(ar, version);
143  ar & speciesId_;
144  ar & p_;
145  ar & capacity_;
146  ar & nAtom_;
147  ar & nMolecule_;
148  ar & accumulator_;
149  ar & projector_;
150 
151  }
152 
153 }
154 #endif
void serialize(Archive &ar, BoundaryEnsemble::Type &data, const unsigned int version)
Serialize a BoundaryEnsemble::Type enum value.
A set of interacting Molecules enclosed by a Boundary.
Definition: System.h:115
File containing preprocessor macros for error handling.
Classes used by all simpatico molecular simulations.
Autocorrelation for Rouse mode coefficients of a linear molecule.
Saving / output archive for binary ostream.
Auto-correlation function for an ensemble of sequences.
Definition: AutoCorrArray.h:57
Utility classes for scientific computation.
Definition: accumulators.mod:1
Template for Analyzer associated with one System.
Dynamically allocatable contiguous array template.
Definition: DArray.h:31
Saving archive for binary istream.
Single-processor Monte Carlo (MC) and molecular dynamics (MD).
void serialize(Archive &ar, T &data, const unsigned int version)
Serialize one object of type T.
Definition: serialize.h:29
A Species represents a set of chemically similar molecules.