PSCF v1.1
NrIterator.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 "NrIterator.h"
9#include <fd1d/System.h>
10#include <pscf/inter/Interaction.h>
11#include <util/misc/Log.h>
12
13#include <math.h>
14
15namespace Pscf {
16namespace Fd1d
17{
18
19 using namespace Util;
20
22 : Iterator(),
23 epsilon_(0.0),
24 maxItr_(0),
25 isAllocated_(false),
26 newJacobian_(false),
27 needsJacobian_(true),
28 isCanonical_(true)
29 { setClassName("NrIterator"); }
30
32 : Iterator(system),
33 epsilon_(0.0),
34 maxItr_(0),
35 isAllocated_(false),
36 newJacobian_(false),
37 needsJacobian_(true),
38 isCanonical_(true)
39 { setClassName("NrIterator"); }
40
42 {}
43
44 void NrIterator::readParameters(std::istream& in)
45 {
46 maxItr_ = 400;
47 read(in, "epsilon", epsilon_);
48 readOptional(in, "maxItr", maxItr_);
49 setup();
50 }
51
52 void NrIterator::setup()
53 {
54 int nm = system().mixture().nMonomer(); // number of monomer types
55 int nx = domain().nx(); // number of grid points
56 UTIL_CHECK(nm > 0);
57 UTIL_CHECK(nx > 0);
58 int nr = nm*nx; // number of residual components
59 if (isAllocated_) {
60 UTIL_CHECK(cArray_.capacity() == nm);
61 UTIL_CHECK(residual_.capacity() == nr);
62 } else {
63 cArray_.allocate(nm);
64 wArray_.allocate(nm);
65 residual_.allocate(nr);
66 jacobian_.allocate(nr, nr);
67 residualNew_.allocate(nr);
68 dOmega_.allocate(nr);
69 wFieldsNew_.allocate(nm);
70 cFieldsNew_.allocate(nm);
71 for (int i = 0; i < nm; ++i) {
72 wFieldsNew_[i].allocate(nx);
73 cFieldsNew_[i].allocate(nx);
74 }
75 solver_.allocate(nr);
76 isAllocated_ = true;
77 }
78 }
79
80 void NrIterator::computeResidual(Array<WField> const & wFields,
81 Array<CField> const & cFields,
82 Array<double>& residual)
83 {
84 int nm = system().mixture().nMonomer(); // number of monomer types
85 int nx = domain().nx(); // number of grid points
86 int i; // grid point index
87 int j; // monomer indices
88 int ir; // residual index
89
90 // Loop over grid points
91 for (i = 0; i < nx; ++i) {
92
93 // Copy volume fractions at grid point i to cArray_
94 for (j = 0; j < nm; ++j) {
95 cArray_[j] = cFields[j][i];
96 }
97
98 // Compute w fields, without Langrange multiplier, from c fields
99 system().interaction().computeW(cArray_, wArray_);
100
101 // Initial residual = wPredicted(from above) - actual w
102 for (j = 0; j < nm; ++j) {
103 ir = j*nx + i;
104 residual[ir] = wArray_[j] - wFields[j][i];
105 }
106
107 // Residuals j = 1, ..., nm-1 are differences from component j=0
108 for (j = 1; j < nm; ++j) {
109 ir = j*nx + i;
110 residual[ir] = residual[ir] - residual[i];
111 }
112
113 // Residual for component j=0 then imposes incompressiblity
114 residual[i] = -1.0;
115 for (j = 0; j < nm; ++j) {
116 residual[i] += cArray_[j];
117 }
118 }
119
120 /*
121 * Note: In canonical ensemble, the spatial integral of the
122 * incompressiblity residual is guaranteed to be zero, as a result of how
123 * volume fractions are computed in SCFT. One of the nx incompressibility
124 * constraints is thus redundant. To avoid this redundancy, replace the
125 * incompressibility residual at the last grid point by a residual that
126 * requires the w field for the last monomer type at the last grid point
127 * to equal zero.
128 */
129
130 if (isCanonical_) {
131 residual[nx-1] = wFields[nm-1][nx-1];
132 }
133
134 }
135
136 /*
137 * Compute Jacobian matrix numerically, by evaluating finite differences.
138 */
139 void NrIterator::computeJacobian()
140 {
141 int nm = system().mixture().nMonomer(); // number of monomer types
142 int nx = domain().nx(); // number of grid points
143 int i; // monomer index
144 int j; // grid point index
145
146 // Copy system().wFields to wFieldsNew.
147 for (i = 0; i < nm; ++i) {
148 for (j = 0; j < nx; ++j) {
149 UTIL_CHECK(nx == wFieldsNew_[i].capacity());
150 UTIL_CHECK(nx == system().wField(i).capacity());
151 wFieldsNew_[i][j] = system().wField(i)[j];
152 }
153 }
154
155 // Compute jacobian, column by column
156 double delta = 0.001;
157 int nr = nm*nx; // number of residual elements
158 int jr; // jacobian row index
159 int jc = 0; // jacobian column index
160 for (i = 0; i < nm; ++i) {
161 for (j = 0; j < nx; ++j) {
162 wFieldsNew_[i][j] += delta;
163 system().mixture().compute(wFieldsNew_, cFieldsNew_);
164 computeResidual(wFieldsNew_, cFieldsNew_, residualNew_);
165 for (jr = 0; jr < nr; ++jr) {
166 jacobian_(jr, jc) =
167 (residualNew_[jr] - residual_[jr])/delta;
168 }
169 wFieldsNew_[i][j] = system().wField(i)[j];
170 ++jc;
171 }
172 }
173
174 // Decompose Jacobian matrix
175 solver_.computeLU(jacobian_);
176
177 }
178
179 void NrIterator::incrementWFields(Array<WField> const & wOld,
180 Array<double> const & dW,
181 Array<WField> & wNew)
182 {
183 int nm = system().mixture().nMonomer(); // number of monomers types
184 int nx = domain().nx(); // number of grid points
185 int i; // monomer index
186 int j; // grid point index
187 int k = 0; // residual element index
188
189 // Add dW
190 for (i = 0; i < nm; ++i) {
191 for (j = 0; j < nx; ++j) {
192 wNew[i][j] = wOld[i][j] - dW[k];
193 ++k;
194 }
195 }
196
197 // If canonical, shift such that last element is exactly zero
198 if (isCanonical_) {
199 double shift = wNew[nm-1][nx-1];
200 for (i = 0; i < nm; ++i) {
201 for (j = 0; j < nx; ++j) {
202 wNew[i][j] -= shift;
203 }
204 }
205 }
206
207 }
208
209 double NrIterator::residualNorm(Array<double> const & residual) const
210 {
211 int nm = system().mixture().nMonomer(); // number of monomer types
212 int nx = domain().nx(); // number of grid points
213 int nr = nm*nx; // number of residual components
214 double value, norm;
215 norm = 0.0;
216 for (int ir = 0; ir < nr; ++ir) {
217 value = fabs(residual[ir]);
218 if (value > norm) {
219 norm = value;
220 }
221 }
222 return norm;
223 }
224
225 int NrIterator::solve(bool isContinuation)
226 {
227 int nm = system().mixture().nMonomer(); // number of monomer types
228 int np = system().mixture().nPolymer(); // number of polymer species
229 int nx = domain().nx(); // number of grid points
230 int nr = nm*nx; // number of residual elements
231
232 // Determine if isCanonical (iff all species ensembles are closed)
233 isCanonical_ = true;
234 Species::Ensemble ensemble;
235 for (int i = 0; i < np; ++i) {
236 ensemble = system().mixture().polymer(i).ensemble();
237 if (ensemble == Species::Unknown) {
238 UTIL_THROW("Unknown species ensemble");
239 }
240 if (ensemble == Species::Open) {
241 isCanonical_ = false;
242 }
243 }
244
245 // If isCanonical, shift so that last element is zero.
246 // Note: This is one of the residuals in this case.
247 if (isCanonical_) {
248 double shift = wFields()[nm-1][nx-1];
249 int i, j;
250 for (i = 0; i < nm; ++i) {
251 for (j = 0; j < nx; ++j) {
252 wFields()[i][j] -= shift;
253 }
254 }
255 }
256
257 // Compute initial residual vector and norm
259 computeResidual(system().wFields(), system().cFields(), residual_);
260 double norm = residualNorm(residual_);
261
262 // Set Jacobian status
263 newJacobian_ = false;
264 if (!isContinuation) {
265 needsJacobian_ = true;
266 }
267
268 // Iterative loop
269 double normNew;
270 int i, j, k;
271 for (i = 0; i < maxItr_; ++i) {
272 Log::file() << "iteration " << i
273 << " , error = " << norm
274 << std::endl;
275
276 if (norm < epsilon_) {
277 Log::file() << "Converged" << std::endl;
279 // Success
280 return 0;
281 }
282
283 if (needsJacobian_) {
284 Log::file() << "Computing jacobian" << std::endl;;
285 computeJacobian();
286 newJacobian_ = true;
287 needsJacobian_ = false;
288 }
289
290 // Compute Newton-Raphson increment dOmega_
291 solver_.solve(residual_, dOmega_);
292
293 // Try full Newton-Raphson update
294 incrementWFields(system().wFields(), dOmega_, wFieldsNew_);
295 system().mixture().compute(wFieldsNew_, cFieldsNew_);
296 computeResidual(wFieldsNew_, cFieldsNew_, residualNew_);
297 normNew = residualNorm(residualNew_);
298
299 // Decrease increment if necessary
300 j = 0;
301 while (normNew > norm && j < 3) {
302 Log::file() << " decreasing increment, error = "
303 << normNew << std::endl;
304 needsJacobian_ = true;
305 for (k = 0; k < nr; ++k) {
306 dOmega_[k] *= 0.66666666;
307 }
308 incrementWFields(system().wFields(), dOmega_, wFieldsNew_);
309 system().mixture().compute(wFieldsNew_, cFieldsNew_);
310 computeResidual(wFieldsNew_, cFieldsNew_, residualNew_);
311 normNew = residualNorm(residualNew_);
312 ++j;
313 }
314
315 // If necessary, try reversing direction
316 if (normNew > norm) {
317 Log::file() << " reversing increment, norm = "
318 << normNew << std::endl;
319 needsJacobian_ = true;
320 for (k = 0; k < nr; ++k) {
321 dOmega_[k] *= -1.000;
322 }
323 incrementWFields(system().wFields(), dOmega_, wFieldsNew_);
324 system().mixture().compute(wFieldsNew_, cFieldsNew_);
325 computeResidual(wFieldsNew_, cFieldsNew_, residualNew_);
326 normNew = residualNorm(residualNew_);
327 }
328
329 // Accept or reject update
330 if (normNew < norm) {
331
332 // Update system fields and residual vector
333 for (j = 0; j < nm; ++j) {
334 for (k = 0; k < nx; ++k) {
335 system().wField(j)[k] = wFieldsNew_[j][k];
336 system().cField(j)[k] = cFieldsNew_[j][k];
337 }
338 }
339 for (j = 0; j < nr; ++j) {
340 residual_[j] = residualNew_[j];
341 }
342 newJacobian_ = false;
343 if (!needsJacobian_) {
344 if (normNew/norm > 0.5) {
345 needsJacobian_ = true;
346 }
347 }
348 norm = normNew;
349 } else {
350 Log::file() << "Iteration failed, norm = "
351 << normNew << std::endl;
352 if (newJacobian_) {
353 return 1;
354 Log::file() << "Unrecoverable failure " << std::endl;
355 } else {
356 Log::file() << "Try rebuilding Jacobian" << std::endl;
357 needsJacobian_ = true;
358 }
359 }
360
361 }
362
363 // Failure
364 return 1;
365 }
366
367} // namespace Fd1d
368} // namespace Pscf
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.
NrIterator()
Default constructor.
Definition: NrIterator.cpp:21
int solve(bool isContinuation=false)
Iterate self-consistent field equations to solution.
Definition: NrIterator.cpp:225
void readParameters(std::istream &in)
Read all parameters and initialize.
Definition: NrIterator.cpp:44
virtual ~NrIterator()
Destructor.
Definition: NrIterator.cpp:41
DArray< System::CField > & cFields()
Get array of all chemical potential fields.
Definition: SystemAccess.h:270
System::WField & wField(int monomerId)
Get chemical potential field for a specific monomer type.
Definition: SystemAccess.h:260
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
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.
Mixture & mixture()
Get Mixture by reference.
Definition: fd1d/System.h:537
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
virtual void computeW(Array< double > const &c, Array< double > &w) const
Compute chemical potential from concentration.
void allocate(int n)
Allocate memory.
Definition: LuSolver.cpp:46
void solve(Array< double > &b, Array< double > &x)
Solve Ax = b for known b to compute x.
Definition: LuSolver.cpp:83
void computeLU(const Matrix< double > &A)
Compute the LU decomposition for later use.
Definition: LuSolver.cpp:63
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
int capacity() const
Return allocated size.
Definition: Array.h:159
void allocate(int capacity)
Allocate the underlying C array.
Definition: DArray.h:199
void allocate(int capacity1, int capacity2)
Allocate memory for a matrix.
Definition: DMatrix.h:170
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.
ScalarParam< Type > & readOptional(std::istream &in, const char *label, Type &value)
Add and read a new optional ScalarParam < Type > object.
#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