PSCF v1.3.2
AmIterator.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 "AmIterator.h"
9#include <r1d/System.h>
10#include <r1d/solvers/Mixture.h>
11#include <pscf/inter/Interaction.h>
12#include <pscf/iterator/NanException.h>
13#include <util/global.h>
14
15namespace Pscf {
16namespace R1d{
17
18 using namespace Util;
19
20 // Public member functions
21
22 /*
23 * Constructor.
24 */
26 : Iterator(system),
27 interaction_()
28 { ParamComposite::setClassName("AmIterator"); }
29
30 /*
31 * Destructor.
32 */
35
36 /*
37 * Read parameters from file.
38 */
39 void AmIterator::readParameters(std::istream& in)
40 {
45
46 const int nMonomer = system().mixture().nMonomer();
47 interaction_.setNMonomer(nMonomer);
48 }
49
50 // Protected virtual function
51
52 /*
53 * Setup before entering iteration loop.
54 */
55 void AmIterator::setup(bool isContinuation)
56 {
58 interaction_.update(system().interaction());
59 }
60
61 // Private virtual functions that interact with parent System
62
63 // Compute and return the number of elements in a field vector
64 int AmIterator::nElements()
65 {
66 const int nm = system().mixture().nMonomer(); // # of monomers
67 const int nx = domain().nx(); // # of grid points
68 return nm*nx;
69 }
70
71 bool AmIterator::hasInitialGuess()
72 {
73 // R1d::System doesn't have an appropriate query function
74 return true;
75 }
76
77 // Get the current fields from the system as a 1D array
78 void AmIterator::getCurrent(DArray<double>& curr)
79 {
80 const int nm = system().mixture().nMonomer(); // # of monomers
81 const int nx = domain().nx(); // # of grid points
82 const DArray< DArray<double> > * currSys = &system().wFields();
83
84 // Straighten out fields into linear arrays
85 for (int i = 0; i < nm; i++) {
86 for (int k = 0; k < nx; k++) {
87 curr[i*nx+k] = (*currSys)[i][k];
88 }
89 }
90 }
91
92 /*
93 * Solve the MDEs for the current system w fields.
94 */
95 void AmIterator::evaluate()
96 { mixture().compute(system().wFields(), system().cFields()); }
97
98 /*
99 * Check whether this is canonical ensemble.
100 */
101 bool AmIterator::isCanonical()
102 {
103 Species::Ensemble ensemble;
104
105 // Check ensemble of all polymers
106 for (int i = 0; i < mixture().nPolymer(); ++i) {
107 ensemble = mixture().polymer(i).ensemble();
108 if (ensemble == Species::Open) {
109 return false;
110 }
111 if (ensemble == Species::Unknown) {
112 UTIL_THROW("Unknown species ensemble");
113 }
114 }
115
116 // Check ensemble of all solvents
117 for (int i = 0; i < mixture().nSolvent(); ++i) {
118 ensemble = mixture().solvent(i).ensemble();
119 if (ensemble == Species::Open) {
120 return false;
121 }
122 if (ensemble == Species::Unknown) {
123 UTIL_THROW("Unknown species ensemble");
124 }
125 }
126
127 // Return true if no species had an open ensemble
128 return true;
129 }
130
131 // Compute the residual for the current system state
132 void AmIterator::getResidual(DArray<double>& resid)
133 {
134 const int nm = system().mixture().nMonomer();
135 const int nx = domain().nx();
136 const int nr = nm*nx;
137
138 // Initialize residuals
139 const double shift = -1.0/interaction_.sumChiInverse();
140 for (int i = 0 ; i < nr; ++i) {
141 resid[i] = shift;
142 }
143
144 // Compute SCF residual vector elements
145 double chi, p;
146 int i, j, k;
147 for (i = 0; i < nm; ++i) {
148 for (j = 0; j < nm; ++j) {
149 DArray<double>& cField = system().cField(j);
150 DArray<double>& wField = system().wField(j);
151 chi = interaction_.chi(i,j);
152 p = interaction_.p(i,j);
153 for (k = 0; k < nx; ++k) {
154 int idx = i*nx + k;
155 resid[idx] += chi*cField[k] - p*wField[k];
156 }
157 }
158 }
159
160 }
161
162 // Update the current system field coordinates
163 void AmIterator::update(DArray<double>& newGuess)
164 {
165 const int nm = mixture().nMonomer(); // # of monomers
166 const int nx = domain().nx(); // # of grid points
167
168 // Set current w fields in system
169 // If canonical, shift to as to set last element to zero
170 const double shift = newGuess[nm*nx - 1];
171 for (int i = 0; i < nm; i++) {
172 DArray<double>& wField = system().wField(i);
173 if (isCanonical()) {
174 for (int k = 0; k < nx; k++) {
175 wField[k] = newGuess[i*nx + k] - shift;
176 }
177 } else {
178 for (int k = 0; k < nx; k++) {
179 wField[k] = newGuess[i*nx + k];
180 }
181 }
182 }
183
184 }
185
186 void AmIterator::outputToLog()
187 {}
188
189}
190}
AmIteratorTmpl< Iterator, DArray< double > > AmTmpl
Alias for base class template.
void readMixingParameters(std::istream &in, bool useLambdaRamp=true)
virtual void readParameters(std::istream &in)
void readParameters(std::istream &in) override
Read all parameters and initialize.
AmIterator(System &system)
Constructor.
~AmIterator()
Destructor.
void setup(bool isContinuation) override
Setup iterator just before entering iteration loop.
Base class for iterative solvers for SCF equations.
Main class in SCFT simulation of one system.
Definition r1d/System.h:65
Ensemble
Statistical ensemble for number of molecules.
Definition Species.h:40
Dynamically allocatable contiguous array template.
Definition DArray.h:32
void setClassName(const char *className)
Set class name string.
File containing preprocessor macros for error handling.
#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.