PSCF v1.2
rpg/fts/perturbation/EinsteinCrystalPerturbation.tpp
1#ifndef RPG_EINSTEIN_CRYSTAL_PERTURBATION_TPP
2#define RPG_EINSTEIN_CRYSTAL_PERTURBATION_TPP
3
4#include "EinsteinCrystalPerturbation.h"
5#include <rpg/fts/simulator/Simulator.h>
6#include <rpg/fts/VecOpFts.h>
7#include <rpg/System.h>
8#include <prdc/cuda/resources.h>
9#include <util/global.h>
10
11namespace Pscf {
12namespace Rpg {
13
14 using namespace Util;
15 using namespace Prdc::Cuda;
16
17 /*
18 * Constructor.
19 */
20 template <int D>
22 Simulator<D>& simulator)
23 : Perturbation<D>(simulator),
24 ecHamiltonian_(0.0),
25 unperturbedHamiltonian_(0.0),
26 stateEcHamiltonian_(0.0),
27 stateUnperturbedHamiltonian_(0.0),
28 hasEpsilon_(false)
29 { setClassName("EinsteinCrystalPerturbation"); }
30
31 /*
32 * Destructor.
33 */
34 template <int D>
37
38 /*
39 * Read parameters from stream, empty default implementation.
40 */
41 template <int D>
43 {
44 read(in, "lambda", lambda_);
45
46 // Allocate and initialize epsilon_ array
47 const int nMonomer = system().mixture().nMonomer();
48 epsilon_.allocate(nMonomer - 1);
49 for (int i = 0; i < nMonomer - 1 ; ++i) {
50 epsilon_[i] = 0.0;
51 }
52
53 // Optionally read the parameters used in Einstein crystal integration
54 hasEpsilon_ =
55 readOptionalDArray(in, "epsilon", epsilon_, nMonomer-1).isActive();
56
57 read(in, "referenceFieldFileName", referenceFieldFileName_);
58 }
59
60 /*
61 * Setup before simulation.
62 */
63 template <int D>
65 {
66 const int nMonomer = system().mixture().nMonomer();
67 const IntVec<D>
68 meshDimensions = system().domain().mesh().dimensions();
69
70 // Allocate memory for reference field
71 w0_.allocate(nMonomer);
72 wc0_.allocate(nMonomer);
73 for (int i = 0; i < nMonomer; ++i) {
74 w0_[i].allocate(meshDimensions);
75 wc0_[i].allocate(meshDimensions);
76 }
77
78 /*
79 * If the user did not input values for epsilon, set values to -1.0
80 * times the nontrivial eigenvalues of the projected chi matrix.
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 FieldIo<D> const & fieldIo = system().domain().fieldIo();
95 fieldIo.readFieldsRGrid(referenceFieldFileName_,
96 w0_, system().domain().unitCell());
97
98 // Compute eigenvector components of the reference field
99 computeWcReference();
100 }
101
102 /*
103 * Compute and return perturbation to Hamiltonian.
104 */
105 template <int D>
106 double
107 EinsteinCrystalPerturbation<D>::hamiltonian(double unperturbedHamiltonian)
108 {
109 // Compute Einstein crystal Hamiltonian
110 const int nMonomer = system().mixture().nMonomer();
111 const int meshSize = system().domain().mesh().size();
112 const IntVec<D> meshDimensions = system().domain().mesh().dimensions();
113 const double vSystem = system().domain().unitCell().volume();
114 const double vMonomer = system().mixture().vMonomer();
115 const double nMonomerSystem = vSystem / vMonomer;
116 double prefactor;
117 RField<D> wcs;
118 wcs.allocate(meshDimensions);
119 ecHamiltonian_ = 0.0;
120
121 for (int j = 0; j < nMonomer - 1; ++j) {
122 RField<D> const & Wc = simulator().wc(j);
123 prefactor = double(nMonomer)/(2.0 * epsilon_[j]);
124 VecOp::subVV(wcs, Wc, wc0_[j]); // wcs = Wc - wc0_[j]
125 double wSqure = 0;
126 wSqure = Reduce::innerProduct(wcs, wcs);
127 ecHamiltonian_ += prefactor * wSqure;
128 }
129
130 // Normalize EC hamiltonian to equal a value per monomer
131 ecHamiltonian_ /= double(meshSize);
132
133 // Compute EC hamiltonian of system
134 ecHamiltonian_ *= nMonomerSystem;
135
136 // Set unperturbedHamiltonian_ member variable
137 unperturbedHamiltonian_ = unperturbedHamiltonian;
138
139 return (1.0 - lambda_)*(ecHamiltonian_ - unperturbedHamiltonian_);
140 }
141
142 /*
143 * Modify functional derivatives, empty default implementation.
144 */
145 template <int D>
146 void
148 {
149 const IntVec<D> meshDimensions = system().domain().mesh().dimensions();
150 const int nMonomer = system().mixture().nMonomer();
151 const double vMonomer = system().mixture().vMonomer();
152 double prefactor;
153 RField<D> DcEC;
154 DcEC.allocate(meshDimensions);
155
156 // Loop over composition eigenvectors (exclude the last)
157 for (int i = 0; i < nMonomer - 1; ++i) {
158 RField<D>& Dc = dc[i];
159 RField<D> const & Wc = simulator().wc(i);
160 prefactor = double(nMonomer) / (epsilon_[i] * vMonomer);
161
162 // Compute EC derivative
163 // DcEC = prefactor * (Wc - wc0_);
164 VecOp::subVV(DcEC, Wc, wc0_[i]);
165 VecOp::mulEqS(DcEC, prefactor);
166
167 // Compute composite derivative
168 // Dc = (Dc * lambda_) + (DcEC * (1-lambda_))
169 VecOp::addVcVc(Dc, Dc, lambda_, DcEC, 1 - lambda_);
170 }
171 }
172
173 /*
174 * Compute and return derivative of free energy with respect to lambda.
175 */
176 template <int D>
178 return unperturbedHamiltonian_ - ecHamiltonian_;
179 }
180
181 /*
182 * Save any required internal state variables.
183 */
184 template <int D>
186 stateEcHamiltonian_ = ecHamiltonian_;
187 stateUnperturbedHamiltonian_ = unperturbedHamiltonian_;
188 }
189
190 /*
191 * Save any required internal state variables.
192 */
193 template <int D>
195 ecHamiltonian_ = stateEcHamiltonian_;
196 unperturbedHamiltonian_ = stateUnperturbedHamiltonian_;
197 }
198
199 /*
200 * Compute the eigenvector components of the w fields, using the
201 * eigenvectors chiEvecs of the projected chi matrix as a basis.
202 */
203 template <int D>
205 {
206 const int nMonomer = system().mixture().nMonomer();
207 int i, j;
208
209 // Loop over eigenvectors (i is an eigenvector index)
210 for (i = 0; i < nMonomer; ++i) {
211
212 // Loop over grid points to zero out field wc0_[i]
213 RField<D>& Wc = wc0_[i];
214
215 VecOp::eqS(Wc, 0.0); // initialize to 0
216
217 // Loop over monomer types (j is a monomer index)
218 for (j = 0; j < nMonomer; ++j) {
219 cudaReal vec;
220 vec = (cudaReal) simulator().chiEvecs(i, j)/nMonomer;
221
222 // Loop over grid points
223 VecOp::addEqVc(Wc, w0_[j], vec);
224 }
225
226 }
227
228 #if 0
229 // Debugging output
230 std::string filename = "wc";
231 system().fieldIo().writeFieldsRGrid(filename, wc0_,
232 system().domain().unitCell(), false);
233 #endif
234 }
235
236}
237}
238#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.
void allocate(IntVec< D > const &meshDimensions)
Allocate the underlying C array for an FFT grid.
Perturbation for Einstein crystal thermodynamic integration method.
virtual double df()
Compute and return derivative of free energy.
void setClassName(const char *className)
Set class name string.
virtual double hamiltonian(double unperturbedHamiltonian)
Compute and return the perturbation to the Hamiltonian.
virtual void setup()
Complete any required initialization.
virtual void incrementDc(DArray< RField< D > > &dc)
Modify the generalized forces to include perturbation.
virtual void restoreState()
Restore any required internal state variables.
virtual void saveState()
Save any required internal state variables.
virtual void readParameters(std::istream &in)
Read parameters from archive.
File input/output operations and format conversions for fields.
void readFieldsRGrid(std::istream &in, DArray< RField< D > > &fields, UnitCell< D > &unitCell) const
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 rpg/System.h:41
Dynamically allocatable contiguous array template.
File containing preprocessor macros for error handling.
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
Definition global.h:68
cudaReal innerProduct(DeviceArray< cudaReal > const &a, DeviceArray< cudaReal > const &b)
Compute inner product of two real arrays (GPU kernel wrapper).
Definition Reduce.cu:891
void subVV(DeviceArray< cudaReal > &a, DeviceArray< cudaReal > const &b, DeviceArray< cudaReal > const &c, const int beginIdA, const int beginIdB, const int beginIdC, const int n)
Vector subtraction, a[i] = b[i] - c[i], kernel wrapper (cudaReal).
Definition VecOp.cu:1228
void addVcVc(DeviceArray< cudaReal > &a, DeviceArray< cudaReal > const &b, cudaReal const c, DeviceArray< cudaReal > const &d, cudaReal const e)
Vector addition w/ coefficient, a[i] = (b[i]*c) + (d[i]*e), kernel wrapper.
Definition VecOpMisc.cu:304
void mulEqS(DeviceArray< cudaReal > &a, cudaReal const b, const int beginIdA, const int n)
Vector multiplication in-place, a[i] *= b, kernel wrapper (cudaReal).
Definition VecOp.cu:1875
void eqS(DeviceArray< cudaReal > &a, cudaReal const b, const int beginIdA, const int n)
Vector assignment, a[i] = b, kernel wrapper (cudaReal).
Definition VecOp.cu:1054
void addEqVc(DeviceArray< cudaReal > &a, DeviceArray< cudaReal > const &b, cudaReal const c)
Vector addition in-place w/ coefficient, a[i] += b[i] * c, kernel wrapper.
Definition VecOpMisc.cu:343
PSCF package top-level namespace.
Definition param_pc.dox:1
Utility classes for scientific computation.