PSCF v1.2
rpc/fts/perturbation/EinsteinCrystalPerturbation.tpp
1#ifndef RPC_EINSTEIN_CRYSTAL_PERTURBATION_TPP
2#define RPC_EINSTEIN_CRYSTAL_PERTURBATION_TPP
3
4#include "EinsteinCrystalPerturbation.h"
5#include <rpc/fts/simulator/Simulator.h>
6#include <rpc/System.h>
7#include <prdc/cpu/RField.h>
8#include <util/containers/DArray.h>
9#include <util/global.h>
10
11namespace Pscf {
12namespace Rpc {
13
14 using namespace Util;
15
16 /*
17 * Constructor.
18 */
19 template <int D>
21 Simulator<D>& simulator)
22 : Perturbation<D>(simulator),
23 ecHamiltonian_(0.0),
24 unperturbedHamiltonian_(0.0),
25 stateEcHamiltonian_(0.0),
26 stateUnperturbedHamiltonian_(0.0),
27 hasEpsilon_(false)
28 { setClassName("EinsteinCrystal"); }
29
30 /*
31 * Destructor.
32 */
33 template <int D>
36
37 /*
38 * Read parameters from stream, empty default implementation.
39 */
40 template <int D>
42 {
43 read(in, "lambda", lambda_);
44
45 // Allocate and initialize epsilon_ array
46 const int nMonomer = system().mixture().nMonomer();
47 epsilon_.allocate(nMonomer - 1);
48 for (int i = 0; i < nMonomer - 1 ; ++i) {
49 epsilon_[i] = 0.0;
50 }
51
52 // Optionally read the parameters used in Einstein crystal integration
53 hasEpsilon_ =
54 readOptionalDArray(in, "epsilon", epsilon_, nMonomer-1).isActive();
55
56 read(in, "referenceFieldFileName", referenceFieldFileName_);
57 }
58
59 /*
60 * Setup before simulation.
61 */
62 template <int D>
64 {
65 const int nMonomer = system().mixture().nMonomer();
66 const IntVec<D>
67 meshDimensions = system().domain().mesh().dimensions();
68
69 // Allocate memory for reference field
70 w0_.allocate(nMonomer);
71 wc0_.allocate(nMonomer);
72 for (int i = 0; i < nMonomer; ++i) {
73 w0_[i].allocate(meshDimensions);
74 wc0_[i].allocate(meshDimensions);
75 }
76
77 /*
78 * If the user did not input values for epsilon, set
79 * the values of epsilon_ to the nontrivial -1.0 * eigenvalues
80 * of the projected chi matrix by default.
81 */
82 if (!hasEpsilon_){
83 for (int i = 0; i < nMonomer - 1 ; ++i) {
84 epsilon_[i] = -1.0 * simulator().chiEval(i);
85 }
86 }
87
88 // Check that all epsilon values are positive
89 for (int i = 0; i < nMonomer - 1 ; ++i) {
90 UTIL_CHECK(epsilon_[i] > 0.0);
91 }
92
93 // Read in reference field from a file
94 UnitCell<D> tempUnitCell;
95 FieldIo<D> const & fieldIo = system().domain().fieldIo();
96 fieldIo.readFieldsRGrid(referenceFieldFileName_,
97 w0_, tempUnitCell);
98
99 // Compute eigenvector components of the reference field
100 computeWcReference();
101 }
102
103 /*
104 * Compute and return perturbation to Hamiltonian.
105 */
106 template <int D>
107 double
108 EinsteinCrystalPerturbation<D>::hamiltonian(double unperturbedHamiltonian)
109 {
110 // Compute Einstein crystal Hamiltonian
111 const int nMonomer = system().mixture().nMonomer();
112 const int meshSize = system().domain().mesh().size();
113 const double vSystem = system().domain().unitCell().volume();
114 const double vMonomer = system().mixture().vMonomer();
115 const double nMonomerSystem = vSystem / vMonomer;
116 double prefactor, w;
117 ecHamiltonian_ = 0.0;
118
119 for (int j = 0; j < nMonomer - 1; ++j) {
120 RField<D> const & Wc = simulator().wc(j);
121 prefactor = double(nMonomer)/(2.0 * epsilon_[j]);
122 for (int i = 0; i < meshSize; ++i) {
123 w = Wc[i] - wc0_[j][i];
124 ecHamiltonian_ += prefactor*w*w;
125 }
126 }
127
128 // Normalize EC hamiltonian to equal a value per monomer
129 ecHamiltonian_ /= double(meshSize);
130
131 // Compute EC hamiltonian of system
132 ecHamiltonian_ *= nMonomerSystem;
133
134 // Set unperturbedHamiltonian_ member variable
135 unperturbedHamiltonian_ = unperturbedHamiltonian;
136
137 return (1.0 - lambda_)*(ecHamiltonian_ - unperturbedHamiltonian_);
138 }
139
140 /*
141 * Modify functional derivatives, empty default implementation.
142 */
143 template <int D>
144 void
146 {
147 const int meshSize = system().domain().mesh().size();
148 const int nMonomer = system().mixture().nMonomer();
149 const double vMonomer = system().mixture().vMonomer();
150 double DcEC, prefactor;
151
152 // Loop over composition eigenvectors (exclude the last)
153 for (int i = 0; i < nMonomer - 1; ++i) {
154 RField<D>& Dc = dc[i];
155 RField<D> const & Wc = simulator().wc(i);
156 prefactor = double(nMonomer) / (epsilon_[i] * vMonomer);
157
158 // Loop over grid points, compute composite derivative
159 for (int k = 0; k < meshSize; ++k) {
160 DcEC = prefactor * (Wc[k] - wc0_[i][k]);
161 Dc[k] += (1.0 - lambda_) * (DcEC - Dc[k]);
162 }
163
164 }
165 }
166
167 /*
168 * Compute and return derivative of free energy with respect to lambda.
169 */
170 template <int D>
172 return unperturbedHamiltonian_ - ecHamiltonian_;
173 }
174
175 /*
176 * Save any required internal state variables.
177 */
178 template <int D>
180 stateEcHamiltonian_ = ecHamiltonian_;
181 stateUnperturbedHamiltonian_ = unperturbedHamiltonian_;
182 }
183
184 /*
185 * Save any required internal state variables.
186 */
187 template <int D>
189 ecHamiltonian_ = stateEcHamiltonian_;
190 unperturbedHamiltonian_ = stateUnperturbedHamiltonian_;
191 }
192
193 /*
194 * Compute the eigenvector components of the w fields, using the
195 * eigenvectors chiEvecs of the projected chi matrix as a basis.
196 */
197 template <int D>
199 {
200 const int nMonomer = system().mixture().nMonomer();
201 const int meshSize = system().domain().mesh().size();
202 int i, j, k;
203
204 // Loop over eigenvectors (j is an eigenvector index)
205 for (j = 0; j < nMonomer; ++j) {
206
207 // Loop over grid points to zero out field wc_[j]
208 RField<D>& Wc = wc0_[j];
209 for (i = 0; i < meshSize; ++i) {
210 Wc[i] = 0.0;
211 }
212
213 // Loop over monomer types (k is a monomer index)
214 for (k = 0; k < nMonomer; ++k) {
215 double vec = simulator().chiEvecs(j, k)/double(nMonomer);
216
217 // Loop over grid points
218 RField<D> const & Wr = w0_[k];
219 for (i = 0; i < meshSize; ++i) {
220 Wc[i] += vec*Wr[i];
221 }
222 }
223 }
224
225 #if 0
226 Log::file() << "wc " << wc0_.capacity() << "\n";
227 for (i = 0; i < 10; ++i) {
228 Log::file() << "wc_1 " << wc0_[0][i] << "\n";
229 Log::file() << "wc_2 " << wc0_[1][i] << "\n";
230 }
231 #endif
232 }
233
234}
235}
236#endif
An IntVec<D, T> is a D-component vector of elements of integer type T.
Definition IntVec.h:27
Field of real double precision values on an FFT mesh.
Base template for UnitCell<D> classes, D=1, 2 or 3.
Definition rpg/System.h:34
Perturbation for Einstein crystal thermodynamic integration method.
virtual void restoreState()
Restore any required internal state variables.
virtual void readParameters(std::istream &in)
Read parameters from archive.
virtual void saveState()
Save any required internal state variables.
void setClassName(const char *className)
Set class name string.
virtual double df()
Compute and return derivative of free energy.
virtual void incrementDc(DArray< RField< D > > &dc)
Modify the generalized forces to include perturbation.
virtual void setup()
Complete any required initialization.
virtual double hamiltonian(double unperturbedHamiltonian)
Compute and return the perturbation to the Hamiltonian.
File input/output operations and format conversions for fields.
void readFieldsRGrid(std::istream &in, DArray< RField< D > > &fields, UnitCell< D > &unitCell) const override
Read array of RField objects (r-grid fields) from a stream.
Base class for additive perturbations of standard FTS Hamiltonian.
Field theoretic simulator (base class).
Definition rpc/System.h:38
Dynamically allocatable contiguous array template.
static std::ostream & file()
Get log ostream by reference.
Definition Log.cpp:57
File containing preprocessor macros for error handling.
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
Definition global.h:68
PSCF package top-level namespace.
Definition param_pc.dox:1
Utility classes for scientific computation.