10#include <pscf/inter/Interaction.h>
11#include <util/misc/Log.h>
47 read(in,
"epsilon", epsilon_);
52 void NrIterator::setup()
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);
80 void NrIterator::computeResidual(
Array<FieldT> const & wFields,
91 for (i = 0; i < nx; ++i) {
94 for (j = 0; j < nm; ++j) {
102 for (j = 0; j < nm; ++j) {
104 residual[ir] = wArray_[j] -
wFields[j][i];
108 for (j = 1; j < nm; ++j) {
110 residual[ir] = residual[ir] - residual[i];
115 for (j = 0; j < nm; ++j) {
116 residual[i] += cArray_[j];
131 residual[nx-1] =
wFields[nm-1][nx-1];
139 void NrIterator::computeJacobian()
147 for (i = 0; i < nm; ++i) {
148 for (j = 0; j < nx; ++j) {
156 double delta = 0.001;
160 for (i = 0; i < nm; ++i) {
161 for (j = 0; j < nx; ++j) {
162 wFieldsNew_[i][j] += delta;
164 computeResidual(wFieldsNew_, cFieldsNew_, residualNew_);
165 for (jr = 0; jr < nr; ++jr) {
167 (residualNew_[jr] - residual_[jr])/delta;
175 solver_.computeLU(jacobian_);
179 void NrIterator::incrementWFields(
Array<FieldT> const & wOld,
190 for (i = 0; i < nm; ++i) {
191 for (j = 0; j < nx; ++j) {
192 wNew[i][j] = wOld[i][j] - dW[k];
199 double shift = wNew[nm-1][nx-1];
200 for (i = 0; i < nm; ++i) {
201 for (j = 0; j < nx; ++j) {
209 double NrIterator::residualNorm(
Array<double> const & residual)
const
216 for (
int ir = 0; ir < nr; ++ir) {
217 value = fabs(residual[ir]);
235 for (
int i = 0; i < np; ++i) {
237 if (ensemble == Species::Unknown) {
240 if (ensemble == Species::Open) {
241 isCanonical_ =
false;
248 double shift =
wFields()[nm-1][nx-1];
250 for (i = 0; i < nm; ++i) {
251 for (j = 0; j < nx; ++j) {
260 double norm = residualNorm(residual_);
263 newJacobian_ =
false;
264 if (!isContinuation) {
265 needsJacobian_ =
true;
271 for (i = 0; i < maxItr_; ++i) {
273 <<
" , error = " << norm
276 if (norm < epsilon_) {
283 if (needsJacobian_) {
284 Log::file() <<
"Computing jacobian" << std::endl;;
287 needsJacobian_ =
false;
291 solver_.solve(residual_, dOmega_);
296 computeResidual(wFieldsNew_, cFieldsNew_, residualNew_);
297 normNew = residualNorm(residualNew_);
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;
310 computeResidual(wFieldsNew_, cFieldsNew_, residualNew_);
311 normNew = residualNorm(residualNew_);
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;
325 computeResidual(wFieldsNew_, cFieldsNew_, residualNew_);
326 normNew = residualNorm(residualNew_);
330 if (normNew < norm) {
333 for (j = 0; j < nm; ++j) {
334 for (k = 0; k < nx; ++k) {
339 for (j = 0; j < nr; ++j) {
340 residual_[j] = residualNew_[j];
342 newJacobian_ =
false;
343 if (!needsJacobian_) {
344 if (normNew/norm > 0.5) {
345 needsJacobian_ =
true;
350 Log::file() <<
"Iteration failed, norm = "
351 << normNew << std::endl;
354 Log::file() <<
"Unrecoverable failure " << std::endl;
356 Log::file() <<
"Try rebuilding Jacobian" << std::endl;
357 needsJacobian_ =
true;
virtual void computeW(Array< double > const &c, Array< double > &w) const
Compute chemical potential from concentration.
void allocate(int n)
Allocate memory.
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.
Interaction & interaction()
Get interaction (i.e., excess free energy) by reference.
CField & cField(int monomerId)
Get chemical potential field for a specific monomer type.
void computeFreeEnergy()
Compute free energy density and pressure for current fields.
WField & wField(int monomerId)
Get chemical potential field for a specific monomer type.
Mixture & mixture()
Get Mixture by reference.
Ensemble
Statistical ensemble for number of molecules.
Array container class template.
int capacity() const
Return allocated size.
void allocate(int capacity)
Allocate the underlying C array.
void allocate(int capacity1, int capacity2)
Allocate memory for a matrix.
static std::ostream & file()
Get log ostream by reference.
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.
#define UTIL_THROW(msg)
Macro for throwing an Exception, reporting function, file and line number.
SCFT with real 1D fields.
PSCF package top-level namespace.