Simpatico  v1.10
IntraBondTensorAutoCorr.h
1 #ifndef MCMD_INTRA_BOND_TENSOR_AUTO_CORR_H
2 #define MCMD_INTRA_BOND_TENSOR_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 <util/accumulators/AutoCorrArray.h> // member template
13 #include <util/space/Tensor.h> // member template parameter
14 #include <util/containers/DArray.h> // member template
15 #include <util/archives/serialize.h> // used in method template
16 #include <util/global.h> // used in method template
17 
18 
19 namespace Simp {
20  class Species;
21 }
22 
23 namespace McMd
24 {
25 
26  using namespace Util;
27  using namespace Simp;
28 
51  template <class SystemType>
52  class IntraBondTensorAutoCorr : public SystemAnalyzer<SystemType>
53  {
54 
55  public:
56 
62  IntraBondTensorAutoCorr(SystemType &system);
63 
67  virtual ~IntraBondTensorAutoCorr();
68 
74  virtual void readParameters(std::istream& in);
75 
81  virtual void loadParameters(Serializable::IArchive &ar);
82 
88  virtual void save(Serializable::OArchive &ar);
89 
96  template <class Archive>
97  void serialize(Archive& ar, const unsigned int version);
98 
102  virtual void setup();
103 
109  virtual void sample(long iStep);
110 
114  virtual void output();
115 
116  protected:
117 
127 
128  private:
129 
131  std::ofstream outputFile_;
132 
134  AutoCorrArray<Tensor, double> accumulator_;
135 
137  DArray<Tensor> data_;
138 
140  Species* speciesPtr_;
141 
143  int speciesId_;
144 
146  int nMolecule_;
147 
149  int nBond_;
150 
152  int capacity_;
153 
155  bool isInitialized_;
156 
157  };
158 
159  /*
160  * Serialize to/from an archive.
161  */
162  template <class SystemType>
163  template <class Archive>
165  const unsigned int version)
166  {
167  Analyzer::serialize(ar, version);
168  ar & speciesId_;
169  ar & capacity_;
170 
171  ar & nBond_;
172  ar & nMolecule_;
173  ar & accumulator_;
174  }
175 
176 }
177 #endif
void serialize(Archive &ar, BoundaryEnsemble::Type &data, const unsigned int version)
Serialize a BoundaryEnsemble::Type enum value.
File containing preprocessor macros for error handling.
Classes used by all simpatico molecular simulations.
Saving / output archive for binary ostream.
Auto-correlation function for an ensemble of sequences.
Definition: AutoCorrArray.h:57
Autocorrelation for bond stress of a molecule.
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.