PSCF v1.2
rpg/fts/brownian/BdSimulator.tpp
1#ifndef RPG_BD_SIMULATOR_TPP
2#define RPG_BD_SIMULATOR_TPP
3
4/*
5* PSCF - Polymer Self-Consistent Field Theory
6*
7* Copyright 2016 - 2022, The Regents of the University of Minnesota
8* Distributed under the terms of the GNU General Public License.
9*/
10
11#include "BdSimulator.h"
12
13#include <rpg/fts/brownian/BdStep.h>
14#include <rpg/fts/brownian/BdStepFactory.h>
15#include <rpg/fts/analyzer/AnalyzerFactory.h>
16#include <rpg/fts/trajectory/TrajectoryReader.h>
17#include <rpg/fts/trajectory/TrajectoryReaderFactory.h>
18#include <rpg/fts/perturbation/PerturbationFactory.h>
19#include <rpg/fts/perturbation/Perturbation.h>
20#include <rpg/fts/ramp/RampFactory.h>
21#include <rpg/fts/ramp/Ramp.h>
22#include <rpg/fts/compressor/Compressor.h>
23#include <rpg/System.h>
24
25#include <pscf/cuda/CudaRandom.h>
26#include <util/misc/Timer.h>
27#include <util/global.h>
28
29namespace Pscf {
30namespace Rpg {
31
32 using namespace Util;
33
34 /*
35 * Constructor.
36 */
37 template <int D>
39 : Simulator<D>(system),
40 analyzerManager_(*this, system),
41 bdStepPtr_(0),
42 bdStepFactoryPtr_(0),
43 trajectoryReaderFactoryPtr_(0),
44 seed_(0)
45 {
46 setClassName("BdSimulator");
47 bdStepFactoryPtr_ = new BdStepFactory<D>(*this);
48 trajectoryReaderFactoryPtr_
50 }
51
52 /*
53 * Destructor.
54 */
55 template <int D>
57 {
58 if (bdStepFactoryPtr_) {
59 delete bdStepFactoryPtr_;
60 }
61 if (bdStepPtr_) {
62 delete bdStepPtr_;
63 }
64 if (trajectoryReaderFactoryPtr_) {
65 delete trajectoryReaderFactoryPtr_;
66 }
67 }
68
69 /*
70 * Read instructions for creating objects from file.
71 */
72 template <int D>
73 void BdSimulator<D>::readParameters(std::istream &in)
74 {
75 // Optionally read random seed, initialize random number generators
76 readRandomSeed(in);
77
78 // Optionally read and instantiate a BdStep object
79 bool isEnd = false;
80 std::string className;
81 UTIL_CHECK(!hasBdStep());
82 UTIL_CHECK(bdStepFactoryPtr_);
83 bdStepPtr_ = bdStepFactoryPtr_->readObjectOptional(in, *this,
84 className, isEnd);
85
86 // Attempt to read the optional Compressor, Perturbation and Ramp blocks
87 readCompressor(in, isEnd);
88 readPerturbation(in, isEnd);
89 if (hasBdStep()) {
90 UTIL_CHECK(hasCompressor());
91 readRamp(in, isEnd);
92 }
93 // Compressor is required if a BdStep exists
94 // A Ramp is allowed only if a BdStep exists
95
96 // Read AnalyzerManager block
97 Analyzer<D>::baseInterval = 0; // default value
98 readParamCompositeOptional(in, analyzerManager_);
99
100 // Figure out what variables need to be saved in stored state
101 state_.needsCc = false;
102 state_.needsDc = false;
103 state_.needsHamiltonian = false;
104 if (hasBdStep()) {
105 if (stepper().needsCc()){
106 state_.needsCc = true;
107 }
108 if (stepper().needsDc()){
109 state_.needsDc = true;
110 }
111 }
112
113 // Initialize Simulator<D> base class
114 allocate();
115 }
116
117 /*
118 * Setup before the main loop of a simulate or analyze command.
119 */
120 template <int D>
121 void BdSimulator<D>::setup(int nStep)
122 {
123 UTIL_CHECK(system().w().hasData());
124
125 // Eigenanalysis of the projected chi matrix.
126 analyzeChi();
127
128 if (hasPerturbation()) {
129 perturbation().setup();
130 }
131
132 if (hasRamp()) {
133 ramp().setup(nStep);
134 }
135
136 // Solve MDE and compute c-fields for the intial state
137 system().compute();
138
139 // Compress the initial state (adjust pressure-like field)
140 if (hasCompressor()) {
141 compressor().compress();
142 compressor().clearTimers();
143 }
144
145 // Compute field components and Hamiltonian for initial state.
146 computeWc();
147 computeCc();
148 computeDc();
149 computeHamiltonian();
150
151 if (hasBdStep()) {
152 stepper().setup();
153 }
154
155 if (analyzerManager_.size() > 0){
156 analyzerManager_.setup();
157 }
158
159 }
160
161 /*
162 * Perform a field theoretic MC simulation of nStep steps.
163 */
164 template <int D>
166 {
167 UTIL_CHECK(system().w().hasData());
168 UTIL_CHECK(hasBdStep());
169 UTIL_CHECK(hasCompressor());
170
171 // Initial setup
172 setup(nStep);
173
174 // Main simulation loop
175 Timer timer;
176 Timer analyzerTimer;
177 timer.start();
178 iStep_ = 0;
179 if (hasRamp()) {
180 ramp().setParameters(iStep_);
181 }
182
183 // Analysis initial step (if any)
184 analyzerTimer.start();
185 analyzerManager_.sample(iStep_);
186 analyzerTimer.stop();
187
188 for (iTotalStep_ = 0; iTotalStep_ < nStep; ++iTotalStep_) {
189
190 // Take a step (modifies W fields)
191 bool converged;
192 converged = stepper().step();
193
194 if (converged){
195 iStep_++;
196
197 if (hasRamp()) {
198 ramp().setParameters(iStep_);
199 }
200
201 // Analysis (if any)
202 analyzerTimer.start();
203 if (Analyzer<D>::baseInterval != 0) {
204 if (iStep_ % Analyzer<D>::baseInterval == 0) {
205 if (analyzerManager_.size() > 0) {
206 analyzerManager_.sample(iStep_);
207 }
208 }
209 }
210 analyzerTimer.stop();
211
212 } else{
213 Log::file() << "Step: "<< iTotalStep_<< " fail to converge" << "\n";
214 }
215
216 }
217
218 timer.stop();
219 double time = timer.time();
220 double analyzerTime = analyzerTimer.time();
221
222 // Output results analyzers to files
224 analyzerManager_.output();
225 }
226
227 // Output results of ramp
228 if (hasRamp()){
229 Log::file() << std::endl;
230 ramp().output();
231 }
232
233 // Output times for the simulation run
234 Log::file() << std::endl;
235 Log::file() << "nStep " << nStep << std::endl;
236 if (iStep_ != nStep){
237 Log::file() << "nFail Step " << (nStep - iStep_)
238 << std::endl;
239 }
240 Log::file() << "Total run time " << time
241 << " sec" << std::endl;
242 double rStep = double(nStep);
243 Log::file() << "time / nStep " << time / rStep
244 << " sec" << std::endl;
245 Log::file() << "Analyzer run time " << analyzerTime
246 << " sec" << std::endl;
247 Log::file() << std::endl;
248
249 // Output number of times MDE has been solved for the simulation run
250 Log::file() << "MDE counter "
251 << compressor().mdeCounter() << std::endl;
252 Log::file() << std::endl;
253
254 }
255
256 /*
257 * Open, read and analyze a trajectory file
258 */
259 template <int D>
260 void BdSimulator<D>::analyze(int min, int max,
261 std::string classname,
262 std::string filename)
263 {
264 // Preconditions
265 UTIL_CHECK(min >= 0);
266 UTIL_CHECK(max >= min);
268 UTIL_CHECK(analyzerManager_.size() > 0);
269
270 // Construct TrajectoryReader
271 TrajectoryReader<D>* trajectoryReaderPtr;
272 trajectoryReaderPtr = trajectoryReaderFactory().factory(classname);
273 if (!trajectoryReaderPtr) {
274 std::string message;
275 message = "Invalid TrajectoryReader class name " + classname;
276 UTIL_THROW(message.c_str());
277 }
278
279 // Open trajectory file
280 trajectoryReaderPtr->open(filename);
281 trajectoryReaderPtr->readHeader();
282
283 // Main loop over trajectory frames
284 Timer timer;
285 bool hasFrame;
286 timer.start();
287 hasFrame = trajectoryReaderPtr->readFrame();
288
289 for (iStep_ = 0; iStep_ <= max && hasFrame; ++iStep_) {
290 if (hasFrame) {
291 clearData();
292
293 // Initialize analyzers
294 if (iStep_ == min) {
295 setup(iStep_);
296 }
297
298 // Sample property values only for iStep >= min
299 if (iStep_ >= min) {
300 analyzerManager_.sample(iStep_);
301 }
302 }
303
304 hasFrame = trajectoryReaderPtr->readFrame();
305 }
306 timer.stop();
307 Log::file() << "end main loop" << std::endl;
308 int nFrames = iStep_ - min;
309 trajectoryReaderPtr->close();
310 delete trajectoryReaderPtr;
311
312 // Output results of all analyzers to output files
313 analyzerManager_.output();
314
315 // Output number of frames and times
316 Log::file() << std::endl;
317 Log::file() << "# of frames " << nFrames << std::endl;
318 Log::file() << "run time " << timer.time()
319 << " sec" << std::endl;
320 Log::file() << "time / frame " << timer.time()/double(nFrames)
321 << " sec" << std::endl;
322 Log::file() << std::endl;
323
324 }
325
326 /*
327 * Output timer results.
328 */
329 template<int D>
330 void BdSimulator<D>::outputTimers(std::ostream& out)
331 {
332 compressor().outputTimers(out);
333 }
334
335}
336}
337#endif
Abstract base for periodic output and/or analysis actions.
Brownian dynamics simulator.
virtual void outputTimers(std::ostream &out)
Output timing results.
void simulate(int nStep)
Perform a field theoretic Monte-Carlo simulation.
virtual void analyze(int min, int max, std::string classname, std::string filename)
Read and analyze a trajectory file.
void setClassName(const char *className)
Set class name string.
BdSimulator(System< D > &system)
Constructor.
virtual void readParameters(std::istream &in)
Read parameters for a MC simulation.
Factory for subclasses of BdStep.
Field theoretic simulator (base class).
Definition rpg/System.h:41
System< D > & system()
Get parent system by reference.
Main class for calculations that represent one system.
Definition rpg/System.h:107
Factory for subclasses of TrajectoryReader.
virtual bool readFrame()=0
Read a single frame.
virtual void close()=0
Close the trajectory file.
virtual void open(std::string filename)=0
Open trajectory file and read header, if any.
static std::ostream & file()
Get log ostream by reference.
Definition Log.cpp:57
Wall clock timer.
Definition Timer.h:35
void start(TimePoint begin)
Start timing from an externally supplied time.
Definition Timer.cpp:49
double time()
Return the accumulated time, in seconds.
Definition Timer.cpp:37
void stop(TimePoint end)
Stop the clock at an externally supplied time.
Definition Timer.cpp:73
File containing preprocessor macros for error handling.
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
Definition global.h:68
#define UTIL_THROW(msg)
Macro for throwing an Exception, reporting function, file and line number.
Definition global.h:51
PSCF package top-level namespace.
Definition param_pc.dox:1
Utility classes for scientific computation.