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