10#include <r1d/solvers/Mixture.h>
11#include <r1d/solvers/Polymer.h>
12#include <pscf/inter/Interaction.h>
13#include <util/misc/Log.h>
49 read(in,
"epsilon", epsilon_);
54 void NrIterator::setup()
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);
82 void NrIterator::computeResidual(
Array<FieldT> const & wFields,
93 for (i = 0; i < nx; ++i) {
96 for (j = 0; j < nm; ++j) {
104 for (j = 0; j < nm; ++j) {
106 residual[ir] = wArray_[j] -
wFields[j][i];
110 for (j = 1; j < nm; ++j) {
112 residual[ir] = residual[ir] - residual[i];
117 for (j = 0; j < nm; ++j) {
118 residual[i] += cArray_[j];
133 residual[nx-1] =
wFields[nm-1][nx-1];
141 void NrIterator::computeJacobian()
149 for (i = 0; i < nm; ++i) {
150 for (j = 0; j < nx; ++j) {
158 double delta = 0.001;
162 for (i = 0; i < nm; ++i) {
163 for (j = 0; j < nx; ++j) {
164 wFieldsNew_[i][j] += delta;
166 computeResidual(wFieldsNew_, cFieldsNew_, residualNew_);
167 for (jr = 0; jr < nr; ++jr) {
169 (residualNew_[jr] - residual_[jr])/delta;
177 solver_.computeLU(jacobian_);
181 void NrIterator::incrementWFields(
Array<FieldT> const & wOld,
192 for (i = 0; i < nm; ++i) {
193 for (j = 0; j < nx; ++j) {
194 wNew[i][j] = wOld[i][j] - dW[k];
201 double shift = wNew[nm-1][nx-1];
202 for (i = 0; i < nm; ++i) {
203 for (j = 0; j < nx; ++j) {
211 double NrIterator::residualNorm(
Array<double> const & residual)
const
218 for (
int ir = 0; ir < nr; ++ir) {
219 value = fabs(residual[ir]);
237 for (
int i = 0; i < np; ++i) {
239 if (ensemble == Species::Unknown) {
242 if (ensemble == Species::Open) {
243 isCanonical_ =
false;
250 double shift =
wFields()[nm-1][nx-1];
252 for (i = 0; i < nm; ++i) {
253 for (j = 0; j < nx; ++j) {
262 double norm = residualNorm(residual_);
265 newJacobian_ =
false;
266 if (!isContinuation) {
267 needsJacobian_ =
true;
273 for (i = 0; i < maxItr_; ++i) {
275 <<
" , error = " << norm
278 if (norm < epsilon_) {
285 if (needsJacobian_) {
286 Log::file() <<
"Computing jacobian" << std::endl;;
289 needsJacobian_ =
false;
293 solver_.solve(residual_, dOmega_);
298 computeResidual(wFieldsNew_, cFieldsNew_, residualNew_);
299 normNew = residualNorm(residualNew_);
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;
312 computeResidual(wFieldsNew_, cFieldsNew_, residualNew_);
313 normNew = residualNorm(residualNew_);
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;
327 computeResidual(wFieldsNew_, cFieldsNew_, residualNew_);
328 normNew = residualNorm(residualNew_);
332 if (normNew < norm) {
335 for (j = 0; j < nm; ++j) {
336 for (k = 0; k < nx; ++k) {
341 for (j = 0; j < nr; ++j) {
342 residual_[j] = residualNew_[j];
344 newJacobian_ =
false;
345 if (!needsJacobian_) {
346 if (normNew/norm > 0.5) {
347 needsJacobian_ =
true;
352 Log::file() <<
"Iteration failed, norm = "
353 << normNew << std::endl;
356 Log::file() <<
"Unrecoverable failure " << std::endl;
358 Log::file() <<
"Try rebuilding Jacobian" << std::endl;
359 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.