PSCF v1.1
NrIterator.h
1#ifndef FD1D_NR_ITERATOR_H
2#define FD1D_NR_ITERATOR_H
3
4/*
5* PSCF - Polymer Self-Consistent Field Theory
6*
7* Copyright 2016 - 2022, The Regents of the University of Minnesota
8* Distributed under the terms of the GNU General Public License.
9*/
10
11#include "Iterator.h"
12#include <fd1d/solvers/Mixture.h>
13#include <pscf/math/LuSolver.h>
14#include <util/containers/Array.h>
15#include <util/containers/DArray.h>
16#include <util/containers/DMatrix.h>
17
18namespace Pscf {
19namespace Fd1d
20{
21
22 using namespace Util;
23
29 class NrIterator : public Iterator
30 {
31
32 public:
33
38
43
47 NrIterator();
48
55
59 virtual ~NrIterator();
60
66 void readParameters(std::istream& in);
67
74 int solve(bool isContinuation = false);
75
76 private:
77
79 LuSolver solver_;
80
82 DArray<WField> wFieldsNew_;
83
85 DArray<WField> cFieldsNew_;
86
88 DArray<double> cArray_;
89
91 DArray<double> wArray_;
92
94 DArray<double> residual_;
95
97 DMatrix<double> jacobian_;
98
100 DArray<double> residualNew_;
101
103 DArray<double> dOmega_;
104
106 double epsilon_;
107
109 int maxItr_;
110
112 bool isAllocated_;
113
115 bool newJacobian_;
116
118 bool needsJacobian_;
119
121 bool isCanonical_;
122
123 // Private member functions
124
128 void setup();
129
137 void computeResidual(Array<WField> const & wFields,
138 Array<WField> const & cFields,
139 Array<double>& residual);
140
146 double residualNorm(Array<double> const & residual) const;
147
151 void computeJacobian();
152
160 void incrementWFields(Array<WField> const & wOld,
161 Array<double> const & dW,
162 Array<WField>& wNew);
163
164 };
165
166} // namespace Fd1d
167} // namespace Pscf
168#endif
Base class for iterative solvers for SCF equations.
Newton-Raphson Iterator for SCF equations.
Definition: NrIterator.h:30
NrIterator()
Default constructor.
Definition: NrIterator.cpp:21
int solve(bool isContinuation=false)
Iterate self-consistent field equations to solution.
Definition: NrIterator.cpp:225
Mixture::WField WField
Monomer chemical potential field.
Definition: NrIterator.h:37
void readParameters(std::istream &in)
Read all parameters and initialize.
Definition: NrIterator.cpp:44
Mixture::CField CField
Monomer concentration / volume fraction field.
Definition: NrIterator.h:42
virtual ~NrIterator()
Destructor.
Definition: NrIterator.cpp:41
DArray< System::CField > & cFields()
Get array of all chemical potential fields.
Definition: SystemAccess.h:270
const System & system() const
Get parent System by reference.
Definition: SystemAccess.h:158
DArray< System::WField > & wFields()
Get array of all chemical potential fields.
Definition: SystemAccess.h:250
Main class in SCFT simulation of one system.
Definition: fd1d/System.h:63
Solve Ax=b by LU decomposition of A.
Definition: LuSolver.h:31
Array container class template.
Definition: Array.h:34
Dynamically allocated Matrix.
Definition: DMatrix.h:25
C++ namespace for polymer self-consistent field theory (PSCF).
Utility classes for scientific computation.
Definition: accumulators.mod:1