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