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