PSCF v1.1
BinaryRelaxIterator.cpp
1/*
2 * PSCF - Polymer Self-Consistent Field Theory
3 *
4 * Copyright 2016 - 2022, The Regents of the University of Minnesota
5 * Distributed under the terms of the GNU General Public License.
6 */
7
8#include "BinaryRelaxIterator.h"
9#include <fd1d/System.h>
10#include <pscf/inter/Interaction.h>
11#include <util/misc/Log.h>
12#include <util/misc/Timer.h>
13#include <math.h>
14
15namespace Pscf {
16namespace Fd1d
17{
18 using namespace Util;
19
20 // Constructor
22 : Iterator(system),
23 epsilon_(0.0),
24 lambdaPlus_(0.0),
25 lambdaMinus_(0.0),
26 maxItr_(200),
27 isAllocated_(false),
28 isCanonical_(true)
29 { setClassName("BinaryRelaxIterator"); }
30
31 // Destructor
33 {}
34
35 // Read parameter file block and allocate memory
37 {
38 UTIL_CHECK(domain().nx() > 0);
39 read(in, "epsilon", epsilon_);
40 read(in, "maxIter", maxItr_);
41 read(in, "lambdaPlus", lambdaPlus_);
42 read(in, "lambdaMinus", lambdaMinus_);
43
44 allocate();
45 }
46
47 // Allocate required memory (called by readParameters)
48 void BinaryRelaxIterator::allocate()
49 {
50 if (isAllocated_) return;
51 int nm = mixture().nMonomer(); // number of monomer types
52 UTIL_CHECK(nm == 2);
53 int nx = domain().nx(); // number of grid points
54 UTIL_CHECK(nx > 0);
55 cArray_.allocate(nm);
56 wArray_.allocate(nm);
57 dW_.allocate(nm);
58 dWNew_.allocate(nm);
59 wFieldsNew_.allocate(nm);
60 cFieldsNew_.allocate(nm);
61 for (int i = 0; i < nm; ++i) {
62 wFieldsNew_[i].allocate(nx);
63 cFieldsNew_[i].allocate(nx);
64 dW_[i].allocate(nx);
65 dWNew_[i].allocate(nx);
66 }
67 isAllocated_ = true;
68 }
69
70 void BinaryRelaxIterator::computeDW(Array<WField> const & wOld,
71 Array<CField> const & cFields,
72 Array<WField> & dW,
73 double & dWNorm)
74
75 {
76 double dWm, dWp, c0, c1, w0, w1, wm;
77 double dWmNorm = 0.0;
78 double dWpNorm = 0.0;
79 double chi = system().interaction().chi(0,1);
80 int nx = domain().nx(); // number of grid points
81 //For AB diblok
82 for (int i = 0; i < nx; ++i) {
83 c0 = cFields[0][i];
84 c1 = cFields[1][i];
85 w0 = wOld[0][i];
86 w1 = wOld[1][i];
87 wm = w1 - w0;
88 dWp = c1 + c0 - 1.0;
89 dWm = 0.5*( c0 - c1 - wm/chi);
90 dWpNorm += dWp*dWp;
91 dWmNorm += dWm*dWm;
92 dWp *= lambdaPlus_;
93 dWm *= lambdaMinus_;
94 dW[0][i] = dWp - dWm;
95 dW[1][i] = dWp + dWm;
96 }
97 dWNorm = (dWpNorm + dWmNorm)/double(nx);
98 dWNorm = sqrt(dWNorm);
99 }
100
101 void BinaryRelaxIterator::updateWFields(Array<WField> const & wOld,
102 Array<WField> const & dW,
103 Array<WField> & wNew)
104 {
105 int nm = mixture().nMonomer(); // number of monomer types
106 UTIL_CHECK(nm == 2);
107 int nx = domain().nx(); // number of grid points
108 int i; // monomer index
109 int j; // grid point index
110 double w0, w1;
111
112 // Add dW
113 for (j = 0; j < nx; j++){
114 w0 = wOld[0][j];
115 w1 = wOld[1][j];
116 wNew[0][j] = w0 + dW[0][j];
117 wNew[1][j] = w1 + dW[1][j];
118 }
119
120 // If canonical, shift such that last element is exactly zero
121 if (isCanonical_) {
122 double shift = wNew[nm-1][nx-1];
123 for (i = 0; i < nm; ++i) {
124 for (j = 0; j < nx; ++j) {
125 wNew[i][j] -= shift;
126 }
127 }
128 }
129
130 }
131
132 int BinaryRelaxIterator::solve(bool isContinuation)
133 {
134 //Declare Timer
135 Timer timerTotal;
136
137 int nm = mixture().nMonomer(); // number of monomer types
138 UTIL_CHECK(nm == 2);
139 int np = mixture().nPolymer(); // number of polymer species
140 int nx = domain().nx(); // number of grid points
141
142 // Start overall timer
143 timerTotal.start();
144
145 // Determine if isCanonical (iff all species ensembles are closed)
146 isCanonical_ = true;
147 Species::Ensemble ensemble;
148 for (int i = 0; i < np; ++i) {
149 ensemble = mixture().polymer(i).ensemble();
150 if (ensemble == Species::Unknown) {
151 UTIL_THROW("Unknown species ensemble");
152 }
153 if (ensemble == Species::Open) {
154 isCanonical_ = false;
155 }
156 }
157
158 // If isCanonical, shift so that last element is zero.
159 // Note: This is one of the residuals in this case.
160 if (isCanonical_) {
161 double shift = wFields()[nm-1][nx-1];
162 int i, j;
163 for (i = 0; i < nm; ++i) {
164 for (j = 0; j < nx; ++j) {
165 wFields()[i][j] -= shift;
166 }
167 }
168 }
169
170 // Compute initial dWNorm.
172 computeDW(system().wFields(), system().cFields(), dW_, dWNorm_);
173
174 // Iterative loop
175 int i, j, k;
176 for (i = 0; i < maxItr_; ++i) {
177 Log::file() << "iteration " << i
178 << " , error = " << dWNorm_
179 << std::endl;
180
181 if (dWNorm_ < epsilon_) {
182 // Stop timers
183 timerTotal.stop();
184
185 Log::file() << "The epsilon is " << epsilon_<< std::endl;
186 Log::file() << "Converged" << std::endl;
188 // Success
189 Log::file() << "\n\n";
190 // Output timing resultsl;
191 Log::file() << "Total time: "
192 << timerTotal.time() << " s " << std::endl;
193 Log::file() << "Average time cost of each iteration: "
194 << timerTotal.time()/i << " s " << std::endl;
195 Log::file() << "\n\n";
196 return 0;
197 }
198
199 // Try full Fd relaxation update
200 updateWFields(system().wFields(), dW_, wFieldsNew_);
201 mixture().compute(wFieldsNew_, cFieldsNew_);
202 computeDW(wFieldsNew_, cFieldsNew_, dWNew_, dWNormNew_);
203
204 // Decrease increment if necessary
205 j = 0;
206 while (dWNormNew_ > dWNorm_ && j < 3) {
207 //double dWNormDecrease_;
208 Log::file() << " error = "
209 << dWNormNew_ << ", decreasing increment" << std::endl;
210 lambdaPlus_ *= 0.5;
211 lambdaMinus_ *= 0.5;
212 //Print lambdaPlus_ and lambdaMinus_
213 Log::file() << " lambdaPlus = "
214 << lambdaPlus_ << std::endl;
215 Log::file() << " lambdaMinus = "
216 << lambdaMinus_<< std::endl;
217 computeDW(system().wFields(),system().cFields(), dWNew_,
218 dWNormNew_);
219 updateWFields(system().wFields(), dWNew_, wFieldsNew_);
220 mixture().compute(wFieldsNew_, cFieldsNew_);
221 ++j;
222 }
223
224 // Accept or reject update
225 if (dWNormNew_ < dWNorm_) {
226 // Update system fields
227 for (j = 0; j < nm; ++j) {
228 for (k = 0; k < nx; ++k) {
229 system().wField(j)[k] = wFieldsNew_[j][k];
230 system().cField(j)[k] = cFieldsNew_[j][k];
231 dW_[j][k] = dWNew_[j][k];
232 }
233 dWNorm_ = dWNormNew_;
234 }
235 } else {
236 Log::file() << "Iteration failed, norm = "
237 << dWNormNew_ << std::endl;
238 break;
239 }
240
241 }
242 // Failure
243 return 1;
244 }
245
246
247} // namespace Fd1d
248} // namespace Pscf
virtual ~BinaryRelaxIterator()
Destructor.
int solve(bool isContinuation=false)
Iterate self-consistent field equations to solution.
void readParameters(std::istream &in)
Read all parameters and initialize.
BinaryRelaxIterator(System &system)
Constructor.
int nx() const
Get number of spatial grid points, including both endpoints.
Base class for iterative solvers for SCF equations.
void compute(DArray< WField > const &wFields, DArray< CField > &cFields)
Compute concentrations.
DArray< System::CField > & cFields()
Get array of all chemical potential fields.
Definition: SystemAccess.h:270
const Domain & domain() const
Get spatial domain (including grid info) by reference.
Definition: SystemAccess.h:194
const System & system() const
Get parent System by reference.
Definition: SystemAccess.h:158
const Mixture & mixture() const
Get Mixture by reference.
Definition: SystemAccess.h:176
DArray< System::WField > & wFields()
Get array of all chemical potential fields.
Definition: SystemAccess.h:250
Main class in SCFT simulation of one system.
Definition: fd1d/System.h:63
WField & wField(int monomerId)
Get chemical potential field for a specific monomer type.
Definition: fd1d/System.h:603
void computeFreeEnergy()
Compute free energy density and pressure for current fields.
Interaction & interaction()
Get interaction (i.e., excess free energy) by reference.
Definition: fd1d/System.h:549
CField & cField(int monomerId)
Get chemical potential field for a specific monomer type.
Definition: fd1d/System.h:616
double chi(int i, int j) const
Return one element of the chi matrix.
Definition: Interaction.h:163
Polymer & polymer(int id)
Get a polymer object.
Definition: MixtureTmpl.h:220
int nMonomer() const
Get number of monomer types.
Definition: MixtureTmpl.h:190
int nPolymer() const
Get number of polymer species.
Definition: MixtureTmpl.h:194
Ensemble
Statistical ensemble for number of molecules.
Definition: Species.h:36
Array container class template.
Definition: Array.h:34
void allocate(int capacity)
Allocate the underlying C array.
Definition: DArray.h:199
static std::ostream & file()
Get log ostream by reference.
Definition: Log.cpp:57
ScalarParam< Type > & read(std::istream &in, const char *label, Type &value)
Add and read a new required ScalarParam < Type > object.
void setClassName(const char *className)
Set class name string.
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
#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
C++ namespace for polymer self-consistent field theory (PSCF).
Utility classes for scientific computation.
Definition: accumulators.mod:1