PSCF v1.1
fd1d/iterator/AmIterator.cpp
1/*
2* PSCF - Polymer Self-Consistent Field Theory
3*
4* Copyright 2016 - 2019, The Regents of the University of Minnesota
5* Distributed under the terms of the GNU General Public License.
6*/
7
8#include "AmIterator.h"
9#include <fd1d/System.h>
10#include <pscf/inter/Interaction.h>
11#include <pscf/iterator/NanException.h>
12#include <util/global.h>
13
14namespace Pscf {
15namespace Fd1d{
16
17 using namespace Util;
18
19 // Constructor
21 : Iterator(system),
22 interaction_()
23 { setClassName("AmIterator"); }
24
25 // Destructor
27 {}
28
29 // Read parameters from file
30 void AmIterator::readParameters(std::istream& in)
31 {
32 // Call parent class readParameters and readErrorType functions
35
36 const int nMonomer = system().mixture().nMonomer();
37 interaction_.setNMonomer(nMonomer);
38 }
39
40 // Setup before entering iteration loop
41 void AmIterator::setup(bool isContinuation)
42 {
44 interaction_.update(system().interaction());
45 }
46
47 // Assign one vector to another: a = b
48 void AmIterator::setEqual(DArray<double>& a, DArray<double> const & b)
49 { a = b; }
50
51 // Compute and return inner product of two vectors
52 double AmIterator::dotProduct(DArray<double> const & a,
53 DArray<double> const & b)
54 {
55 const int n = a.capacity();
56 UTIL_CHECK(n == b.capacity());
57 double product = 0.0;
58 for (int i = 0; i < n; i++) {
59 // if either value is NaN, throw NanException
60 if (std::isnan(a[i]) || std::isnan(b[i])) {
61 throw NanException("AmIterator::dotProduct",__FILE__,__LINE__,0);
62 }
63 product += a[i] * b[i];
64 }
65 return product;
66 }
67
68 // Compute and return maximum element of residual vector.
69 double AmIterator::maxAbs(DArray<double> const & hist)
70 {
71 const int n = hist.capacity();
72 double maxRes = 0.0;
73 double value;
74 for (int i = 0; i < n; i++) {
75 value = hist[i];
76 if (std::isnan(value)) { // if value is NaN, throw NanException
77 throw NanException("AmIterator::dotProduct",__FILE__,__LINE__,0);
78 }
79 if (fabs(value) > maxRes)
80 maxRes = fabs(value);
81 }
82 return maxRes;
83 }
84
85 // Update basis
86 void AmIterator::updateBasis(RingBuffer<DArray<double> > & basis,
87 RingBuffer<DArray<double> > const & hists)
88 {
89 // Make sure at least two histories are stored
90 UTIL_CHECK(hists.size() >= 2);
91
92 const int n = hists[0].capacity();
93 DArray<double> newbasis;
94 newbasis.allocate(n);
95
96 // Basis vector is different of last to vectors
97 for (int i = 0; i < n; i++) {
98 newbasis[i] = hists[0][i] - hists[1][i];
99 }
100
101 basis.append(newbasis);
102 }
103
104 // Add linear combination of basis vector to current trial state
105 void AmIterator::addHistories(DArray<double>& trial,
106 RingBuffer< DArray<double> > const& basis,
107 DArray<double> coeffs,
108 int nHist)
109 {
110 int n = trial.capacity();
111 for (int i = 0; i < nHist; i++) {
112 for (int j = 0; j < n; j++) {
113 trial[j] += coeffs[i] * -1 * basis[i][j];
114 }
115 }
116 }
117
118 void AmIterator::addPredictedError(DArray<double>& fieldTrial,
119 DArray<double> const & resTrial,
120 double lambda)
121 {
122 int n = fieldTrial.capacity();
123 for (int i = 0; i < n; i++) {
124 fieldTrial[i] += lambda * resTrial[i];
125 }
126 }
127
128 bool AmIterator::hasInitialGuess()
129 {
130 // Fd1d::System doesn't hav a hasFields() function
131 return true;
132 }
133
134 // Compute and return the number of elements in a field vector
135 int AmIterator::nElements()
136 {
137 const int nm = system().mixture().nMonomer(); // # of monomers
138 const int nx = domain().nx(); // # of grid points
139 return nm*nx;
140 }
141
142 // Get the current fields from the system as a 1D array
143 void AmIterator::getCurrent(DArray<double>& curr)
144 {
145 const int nm = system().mixture().nMonomer(); // # of monomers
146 const int nx = domain().nx(); // # of grid points
147 const DArray< DArray<double> > * currSys = &system().wFields();
148
149 // Straighten out fields into linear arrays
150 for (int i = 0; i < nm; i++) {
151 for (int k = 0; k < nx; k++) {
152 curr[i*nx+k] = (*currSys)[i][k];
153 }
154 }
155 }
156
157 // Solve the MDEs for the current system w fields
158 void AmIterator::evaluate()
159 {
160 mixture().compute(system().wFields(), system().cFields());
161 }
162
163
164 // Check whether Canonical ensemble
165 bool AmIterator::isCanonical()
166 {
167 Species::Ensemble ensemble;
168
169 // Check ensemble of all polymers
170 for (int i = 0; i < mixture().nPolymer(); ++i) {
171 ensemble = mixture().polymer(i).ensemble();
172 if (ensemble == Species::Open) {
173 return false;
174 }
175 if (ensemble == Species::Unknown) {
176 UTIL_THROW("Unknown species ensemble");
177 }
178 }
179
180 // Check ensemble of all solvents
181 for (int i = 0; i < mixture().nSolvent(); ++i) {
182 ensemble = mixture().solvent(i).ensemble();
183 if (ensemble == Species::Open) {
184 return false;
185 }
186 if (ensemble == Species::Unknown) {
187 UTIL_THROW("Unknown species ensemble");
188 }
189 }
190
191 // Return true if no species had an open ensemble
192 return true;
193 }
194
195 // Compute the residual for the current system state
196 void AmIterator::getResidual(DArray<double>& resid)
197 {
198 const int nm = system().mixture().nMonomer();
199 const int nx = domain().nx();
200 const int nr = nm*nx;
201
202 // Initialize residuals
203 const double shift = -1.0/interaction_.sumChiInverse();
204 for (int i = 0 ; i < nr; ++i) {
205 resid[i] = shift;
206 }
207
208 // Compute SCF residual vector elements
209 double chi, p;
210 int i, j, k;
211 for (i = 0; i < nm; ++i) {
212 for (j = 0; j < nm; ++j) {
213 DArray<double>& cField = system().cField(j);
214 DArray<double>& wField = system().wField(j);
215 chi = interaction_.chi(i,j);
216 p = interaction_.p(i,j);
217 for (k = 0; k < nx; ++k) {
218 int idx = i*nx + k;
219 resid[idx] += chi*cField[k] - p*wField[k];
220 }
221 }
222 }
223
224 }
225
226 // Update the current system field coordinates
227 void AmIterator::update(DArray<double>& newGuess)
228 {
229 const int nm = mixture().nMonomer(); // # of monomers
230 const int nx = domain().nx(); // # of grid points
231
232 // Set current w fields in system
233 // If canonical, shift to as to set last element to zero
234 const double shift = newGuess[nm*nx - 1];
235 for (int i = 0; i < nm; i++) {
236 DArray<double>& wField = system().wField(i);
237 if (isCanonical()) {
238 for (int k = 0; k < nx; k++) {
239 wField[k] = newGuess[i*nx + k] - shift;
240 }
241 } else {
242 for (int k = 0; k < nx; k++) {
243 wField[k] = newGuess[i*nx + k];
244 }
245 }
246 }
247
248 }
249
250 void AmIterator::outputToLog()
251 {}
252
253}
254}
Template for Anderson mixing iterator algorithm.
void readErrorType(std::istream &in)
Read and validate the optional errorType string parameter.
double chi(int i, int j) const
Return one element of the chi matrix.
double sumChiInverse() const
Return sum of elements of the inverse chi matrix.
void setNMonomer(int nMonomer)
Set number of monomers and allocate required memory.
void update(Interaction const &interaction)
Update all computed quantities.
double p(int i, int j) const
Return one element of the potent matrix P.
void readParameters(std::istream &in)
Read all parameters and initialize.
void setup(bool isContinuation)
Setup iterator just before entering iteration loop.
AmIterator(System &system)
Constructor.
Base class for iterative solvers for SCF equations.
Main class in SCFT simulation of one system.
Definition: fd1d/System.h:63
Exception thrown when not-a-number (NaN) is encountered.
Definition: NanException.h:40
Ensemble
Statistical ensemble for number of molecules.
Definition: Species.h:36
int capacity() const
Return allocated size.
Definition: Array.h:159
Dynamically allocatable contiguous array template.
Definition: DArray.h:32
void allocate(int capacity)
Allocate the underlying C array.
Definition: DArray.h:199
Class for storing history of previous values in an array.
Definition: RingBuffer.h:27
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
C++ namespace for polymer self-consistent field theory (PSCF).
Utility classes for scientific computation.
Definition: accumulators.mod:1
float product(float a, float b)
Product for float Data.
Definition: product.h:22