Loading [MathJax]/extensions/TeX/AMSsymbols.js
PSCF v1.2
Interaction.cpp
1/*
2* PSCF - Polymer Self-Consistent Field Theory
3*
4* Copyright 2016 - 2022, The Regents of the University of Minnesota
5* Distributed under the terms of the GNU General Public License.
6*/
7
8#include "Interaction.h"
9#include <pscf/math/LuSolver.h>
10
11namespace Pscf {
12
13 using namespace Util;
14
15 /*
16 * Constructor.
17 */
19 : nMonomer_(0)
20 { setClassName("Interaction"); }
21
22 /*
23 * Destructor.
24 */
27
28 /*
29 * Set the number of monomer types.
30 */
31 void Interaction::setNMonomer(int nMonomer)
32 {
33 UTIL_CHECK(nMonomer_ == 0);
35 nMonomer_ = nMonomer;
37 chiInverse_.allocate(nMonomer, nMonomer);
38 }
39
40 /*
41 * Read chi matrix from file.
42 */
43 void Interaction::readParameters(std::istream& in)
44 {
45 UTIL_CHECK(nMonomer() > 0);
46 readDSymmMatrix(in, "chi", chi_, nMonomer());
47
48 // Compute relevant AM iterator quantities.
49 updateMembers();
50 }
51
52 void Interaction::updateMembers()
53 {
54 UTIL_CHECK(chi_.isAllocated());
55 UTIL_CHECK(chiInverse_.isAllocated());
56
57 if (nMonomer() == 2) {
58 double det = chi_(0,0)*chi_(1, 1) - chi_(0,1)*chi_(1,0);
59 double norm = chi_(0,0)*chi_(0, 0) + chi_(1,1)*chi_(1,1)
60 + 2.0*chi_(0,1)*chi_(1,0);
61 if (fabs(det/norm) < 1.0E-8) {
62 UTIL_THROW("Singular chi matrix");
63 }
64 chiInverse_(0,1) = -chi_(0,1)/det;
65 chiInverse_(1,0) = -chi_(1,0)/det;
66 chiInverse_(1,1) = chi_(0,0)/det;
67 chiInverse_(0,0) = chi_(1,1)/det;
68
69 } else {
70 LuSolver solver;
71 solver.allocate(nMonomer());
72 solver.computeLU(chi_);
73 solver.inverse(chiInverse_);
74 }
75
76 }
77
78 void Interaction::setChi(int i, int j, double chi)
79 {
80 chi_(i,j) = chi;
81 if (i != j) {
82 chi_(j,i) = chi;
83 }
84
85 // Compute relevant AM iterator quantities.
86 updateMembers();
87 }
88
89 /*
90 * Compute and return excess Helmholtz free energy per monomer.
91 */
92 double Interaction::fHelmholtz(Array<double> const & c) const
93 {
94 int i, j;
95 double sum = 0.0;
96 for (i = 0; i < nMonomer(); ++i) {
97 for (j = 0; j < nMonomer(); ++j) {
98 sum += chi_(i, j)* c[i]*c[j];
99 }
100 }
101 return 0.5*sum;
102 }
103
104 /*
105 * Compute chemical potential from monomer concentrations
106 */
107 void
109 Array<double>& w) const
110 {
111 int i, j;
112 for (i = 0; i < nMonomer(); ++i) {
113 w[i] = 0.0;
114 for (j = 0; j < nMonomer(); ++j) {
115 w[i] += chi_(i, j)* c[j];
116 }
117 }
118 }
119
120 /*
121 * Compute concentrations and xi from chemical potentials.
122 */
123 void
125 Array<double>& c, double& xi)
126 const
127 {
128 double sum1 = 0.0;
129 double sum2 = 0.0;
130 int i, j;
131 for (i = 0; i < nMonomer(); ++i) {
132 for (j = 0; j < nMonomer(); ++j) {
133 sum1 += chiInverse_(i, j)*w[j];
134 sum2 += chiInverse_(i, j);
135 }
136 }
137 xi = (sum1 - 1.0)/sum2;
138 for (i = 0; i < nMonomer(); ++i) {
139 c[i] = 0;
140 for (j = 0; j < nMonomer(); ++j) {
141 c[i] += chiInverse_(i, j)*( w[j] - xi );
142 }
143 }
144 }
145
146 /*
147 * Compute Langrange multiplier from chemical potentials.
148 */
149 void
151 const
152 {
153 double sum1 = 0.0;
154 double sum2 = 0.0;
155 int i, j;
156 for (i = 0; i < nMonomer(); ++i) {
157 for (j = 0; j < nMonomer(); ++j) {
158 sum1 += chiInverse_(i, j)*w[j];
159 sum2 += chiInverse_(i, j);
160 }
161 }
162 xi = (sum1 - 1.0)/sum2;
163 }
164
165 /*
166 * Return dWdC = chi matrix.
167 */
168 void
170 Matrix<double>& dWdC) const
171 {
172 int i, j;
173 for (i = 0; i < nMonomer(); ++i) {
174 for (j = 0; j < nMonomer(); ++j) {
175 dWdC(i, j) = chi_(i, j);
176 }
177 }
178 }
179
180} // namespace Pscf
virtual double fHelmholtz(Array< double > const &c) const
Compute excess Helmholtz free energy per monomer.
virtual void computeDwDc(Array< double > const &c, Matrix< double > &dWdC) const
Compute second derivatives of free energy.
virtual void readParameters(std::istream &in)
Read chi parameters.
DMatrix< double > const & chi() const
Return the chi matrix by const reference.
virtual void computeW(Array< double > const &c, Array< double > &w) const
Compute chemical potential from concentration.
virtual void computeC(Array< double > const &w, Array< double > &c, double &xi) const
Compute concentration from chemical potential fields.
void setChi(int i, int j, double chi)
Change one element of the chi matrix.
void setNMonomer(int nMonomer)
Set the number of monomer types.
Interaction()
Constructor.
virtual ~Interaction()
Destructor.
int nMonomer() const
Get number of monomer types.
virtual void computeXi(Array< double > const &w, double &xi) const
Compute Langrange multiplier xi from chemical potential fields.
Array container class template.
Definition Array.h:34
void allocate(int capacity1, int capacity2)
Allocate memory for a matrix.
Definition DMatrix.h:170
bool isAllocated() const
Return true if the DMatrix has been allocated, false otherwise.
Definition DMatrix.h:202
Two-dimensional array container template (abstract).
Definition Matrix.h:32
DSymmMatrixParam< Type > & readDSymmMatrix(std::istream &in, const char *label, DMatrix< Type > &matrix, int n)
Add and read a required symmetrix DMatrix.
void setClassName(const char *className)
Set class name string.
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
Definition global.h:68
#define UTIL_THROW(msg)
Macro for throwing an Exception, reporting function, file and line number.
Definition global.h:51
PSCF package top-level namespace.
Definition param_pc.dox:1
Utility classes for scientific computation.