PSCF v1.3
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/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>
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 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 IntVec<D> meshDimensions = system().domain().mesh().dimensions();
114 const double vSystem = system().domain().unitCell().volume();
115 const double vMonomer = system().mixture().vMonomer();
116 const double nMonomerSystem = vSystem / vMonomer;
117 double prefactor;
118 RField<D> wcs;
119 wcs.allocate(meshDimensions);
120 ecHamiltonian_ = 0.0;
121
122 for (int j = 0; j < nMonomer - 1; ++j) {
123 RField<D> const & Wc = simulator().wc(j);
124 prefactor = double(nMonomer)/(2.0 * epsilon_[j]);
125 VecOp::subVV(wcs, Wc, wc0_[j]); // wcs = Wc - wc0_[j]
126 double wSqure = 0;
127 wSqure = Reduce::innerProduct(wcs, wcs);
128 ecHamiltonian_ += prefactor * wSqure;
129 }
130
131 // Normalize EC hamiltonian to equal a value per monomer
132 ecHamiltonian_ /= double(meshSize);
133
134 // Compute EC hamiltonian of system
135 ecHamiltonian_ *= nMonomerSystem;
136
137 // Set unperturbedHamiltonian_ member variable
138 unperturbedHamiltonian_ = unperturbedHamiltonian;
139
140 return (1.0 - lambda_)*(ecHamiltonian_ - unperturbedHamiltonian_);
141 }
142
143 /*
144 * Modify functional derivatives, empty default implementation.
145 */
146 template <int D>
147 void
149 {
150 const IntVec<D> meshDimensions = system().domain().mesh().dimensions();
151 const int nMonomer = system().mixture().nMonomer();
152 const double vMonomer = system().mixture().vMonomer();
153 double prefactor;
154 RField<D> DcEC;
155 DcEC.allocate(meshDimensions);
156
157 // Loop over composition eigenvectors (exclude the last)
158 for (int i = 0; i < nMonomer - 1; ++i) {
159 RField<D>& Dc = dc[i];
160 RField<D> const & Wc = simulator().wc(i);
161 prefactor = double(nMonomer) / (epsilon_[i] * vMonomer);
162
163 // Compute EC derivative
164 // DcEC = prefactor * (Wc - wc0_);
165 VecOp::subVV(DcEC, Wc, wc0_[i]);
166 VecOp::mulEqS(DcEC, prefactor);
167
168 // Compute composite derivative
169 // Dc = (Dc * lambda_) + (DcEC * (1-lambda_))
170 VecOp::addVcVc(Dc, Dc, lambda_, DcEC, 1 - lambda_);
171 }
172 }
173
174 /*
175 * Compute and return derivative of free energy with respect to lambda.
176 */
177 template <int D>
179 return unperturbedHamiltonian_ - ecHamiltonian_;
180 }
181
182 /*
183 * Save any required internal state variables.
184 */
185 template <int D>
187 stateEcHamiltonian_ = ecHamiltonian_;
188 stateUnperturbedHamiltonian_ = unperturbedHamiltonian_;
189 }
190
191 /*
192 * Save any required internal state variables.
193 */
194 template <int D>
196 ecHamiltonian_ = stateEcHamiltonian_;
197 unperturbedHamiltonian_ = stateUnperturbedHamiltonian_;
198 }
199
200 /*
201 * Compute the eigenvector components of the w fields, using the
202 * eigenvectors chiEvecs of the projected chi matrix as a basis.
203 */
204 template <int D>
205 void EinsteinCrystalPerturbation<D>::computeWcReference()
206 {
207 const int nMonomer = system().mixture().nMonomer();
208 int i, j;
209
210 // Loop over eigenvectors (i is an eigenvector index)
211 for (i = 0; i < nMonomer; ++i) {
212
213 // Loop over grid points to zero out field wc0_[i]
214 RField<D>& Wc = wc0_[i];
215
216 VecOp::eqS(Wc, 0.0); // initialize to 0
217
218 // Loop over monomer types (j is a monomer index)
219 for (j = 0; j < nMonomer; ++j) {
220 cudaReal vec;
221 vec = (cudaReal) simulator().chiEvecs(i, j)/nMonomer;
222
223 // Loop over grid points
224 VecOp::addEqVc(Wc, w0_[j], vec);
225 }
226
227 }
228 }
229
230}
231}
232#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.
Definition cpu/RField.h:29
void allocate(IntVec< D > const &meshDimensions)
Allocate the underlying C array for an FFT grid.
Base template for UnitCell<D> classes, D=1, 2 or 3.
Definition UnitCell.h:56
virtual double df()
Compute and return derivative of free energy.
DArrayParam< Type > & readOptionalDArray(std::istream &in, const char *label, DArray< Type > &array, int n)
Add and read an optional DArray < Type > parameter.
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.
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.
bool readFieldsRGrid(std::istream &in, DArray< RField< D > > &fields, UnitCell< D > &unitCell) const override
Read array of RField objects (r-grid fields) from a stream.
System< D > const & system() const
Get parent System<D> by const reference.
double lambda_
Strength of the perturbation.
Simulator< D > const & simulator() const
Get parent Simulator<D> by const reference.
Perturbation(Simulator< D > &simulator)
Constructor.
Field theoretic simulator (base class).
Dynamically allocatable contiguous array template.
Definition DArray.h:32
File containing preprocessor macros for error handling.
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
Definition global.h:68
void eqS(DeviceArray< cudaReal > &a, const cudaReal b, const int beginIdA, const int n)
Vector assignment, a[i] = b, kernel wrapper (cudaReal).
Definition VecOp.cu:1073
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:1247
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, const cudaReal b, const int beginIdA, const int n)
Vector multiplication in-place, a[i] *= b, kernel wrapper (cudaReal).
Definition VecOp.cu:1918
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
SCFT and PS-FTS with real periodic fields (GPU)
PSCF package top-level namespace.
Definition param_pc.dox:1