PSCF v1.2
rpc/scft/iterator/MaskGenFilm.tpp
1#ifndef RPC_MASK_GEN_FILM_TPP
2#define RPC_MASK_GEN_FILM_TPP
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 "MaskGenFilm.h"
12#include "Iterator.h"
13#include <rpc/field/FieldIo.h>
14#include <prdc/cpu/RField.h>
15#include <prdc/crystal/UnitCell.h>
16#include <prdc/crystal/paramIdConversions.h>
17#include <pscf/math/IntVec.h>
18#include <util/containers/FArray.h>
19#include <cmath>
20
21namespace Pscf {
22namespace Rpc
23{
24
25 using namespace Util;
26 using namespace Pscf::Prdc;
27 using namespace Pscf::Prdc::Cpu;
28
29 /*
30 * Default constructor
31 */
32 template <int D>
35 sysPtr_(0)
36 { setClassName("MaskGenFilm"); }
37
38 /*
39 * Constructor
40 */
41 template <int D>
44 sysPtr_(&sys)
45 { setClassName("MaskGenFilm"); }
46
47 /*
48 * Destructor
49 */
50 template <int D>
53
54 /*
55 * Get contribution to the stress from this mask
56 */
57 template <int D>
58 double MaskGenFilm<D>::stressTerm(int paramId) const
59 {
60 int normalVecParamId = convertFullParamIdToReduced<D>(normalVecId(),
61 system().domain().lattice());
62
63 // If paramId is not normalVecId, there is no stress contribution
64 if (normalVecParamId != paramId) return 0.0;
65
66 // If this point is reached, stress contribution must be calculated
67 UTIL_CHECK(system().hasMask());
68
69 // Get the length of the lattice basis vector normal to the walls
70 double nvLength = system().domain().unitCell().parameter(paramId);
71
72 // Get the volume fraction of the unit cell occupied by polymers
73 double phiTot = system().mask().phiTot();
74
75 // Create a 3 element vector 'dim' that contains the grid dimensions.
76 // If system is 2D (1D), then the z (y & z) dimensions are set to 1.
77 IntVec<3> dim;
78 for (int ind = 0; ind < 3; ind++) {
79 if (ind < D) {
80 dim[ind] = system().domain().mesh().dimensions()[ind];
81 } else {
82 dim[ind] = 1;
83 }
84 }
85
86 // Create the derivative field in rgrid format
87 RField<D> deriv;
88 FArray<int,3> coords;
89 int x, y, z;
90 int counter = 0;
91 double d, maskVal;
92
93 deriv.allocate(system().domain().mesh().dimensions());
94 RField<D> const & maskRGrid = system().mask().rgrid();
95
96 for (x = 0; x < dim[0]; x++) {
97 coords[0] = x;
98 for (y = 0; y < dim[1]; y++) {
99 coords[1] = y;
100 for (z = 0; z < dim[2]; z++) {
101 coords[2] = z;
102 maskVal = maskRGrid[counter];
103
104 // Get the distance 'd' traveled along the lattice basis
105 // vector normal to the walls, in reduced coordinates
106 d = (double)coords[normalVecId()] /
107 (double)dim[normalVecId()];
108
109 deriv[counter] = maskVal * (maskVal - 1) * 8.0
110 * (std::abs(d - 0.5) - 0.5)
111 / interfaceThickness();
112 counter++;
113 }
114 }
115 }
116
117 // Get xi, the Lagrange multiplier field, in rgrid format
118 int nMonomer = system().mixture().nMonomer();
119 int nx = system().domain().mesh().size();
120 double chi;
121 RField<D> xi;
122 xi.allocate(system().domain().mesh().dimensions());
123
124 for (int i = 0; i < nx; i++) {
125 xi[i] = system().w().rgrid(0)[i];
126 }
127
128 for (int in = 0; in < nMonomer; in++) {
129 chi = system().interaction().chi(0,in);
130 if (fabs(chi) > 1e-6) { // if chi is nonzero
131 for (int i = 0; i < nx; i++) {
132 xi[i] -= system().c().rgrid(in)[i] * chi;
133 }
134 }
135 }
136
137 if (system().hasExternalFields()) {
138 for (int i = 0; i < nx; i++) {
139 xi[i] -= system().h().rgrid(0)[i];
140 }
141 }
142
143 // Get the integral term in the stress
144 double intTerm = 0.0;
145 for (int i = 0; i < nx; i++) {
146 intTerm += xi[i] * deriv[i];
147 }
148 intTerm /= (phiTot * nx);
149
150 // Get the pressure term in the stress
151 if (!sysPtr_->hasFreeEnergy()) {
152 sysPtr_->computeFreeEnergy();
153 }
154 double pSys = sysPtr_->pressure();
155 double pTerm = pSys * excludedThickness() /
156 (phiTot * nvLength * nvLength);
157
158 return pTerm - intTerm;
159 }
160
161 template <int D>
162 double MaskGenFilm<D>::modifyStress(int paramId, double stress)
163 const
164 {
165 int nvParamId = convertFullParamIdToReduced<D>(normalVecId(),
166 system().domain().lattice());
167 if (nvParamId == paramId) {
168 UTIL_CHECK(system().hasMask());
169
170 if (!hasFBulk()) {
171 UTIL_THROW("fBulk must be set before calculating this stress.");
172 }
173
174 // Get system free energy
175 if (!sysPtr_->hasFreeEnergy()) {
176 sysPtr_->computeFreeEnergy();
177 }
178 double fSys = sysPtr_->fHelmholtz();
179
180 // Get the length L of the basis vector normal to the walls
181 double L = system().domain().unitCell().parameter(paramId);
182
183 double modifiedStress = (stress * (L - excludedThickness()))
184 + (fSys - fBulk_);
185
186 return modifiedStress;
187 } else {
188 return stress;
189 }
190 }
191
192 /*
193 * Allocate container necessary to generate and store field
194 */
195 template <int D>
197 {
198 UTIL_CHECK(system().domain().unitCell().isInitialized());
199
200 // Make sure mask has access to a fieldIo
201 system().mask().setFieldIo(system().domain().fieldIo());
202
203 // Allocate the mask containers if needed
204 if (!system().mask().isAllocatedRGrid()) {
205 system().mask().allocateRGrid(system().domain().mesh().dimensions());
206 }
207 if (system().iterator().isSymmetric()) {
208 UTIL_CHECK(system().domain().basis().isInitialized());
209 if (!system().mask().isAllocatedBasis()) {
210 system().mask().allocateBasis(system().domain().basis().nBasis());
211 }
212 }
213 }
214
215 /*
216 * Generate the field and store where the Iterator can access
217 */
218 template <int D>
220 {
221 UTIL_CHECK(interfaceThickness() > 0);
222 UTIL_CHECK(excludedThickness() > interfaceThickness());
223 UTIL_CHECK(system().mask().isAllocatedRGrid());
224 if (system().iterator().isSymmetric()) {
225 UTIL_CHECK(system().mask().isAllocatedBasis());
226 }
227
228 // Get the length L of the lattice basis vector normal to the walls
229 int paramId = convertFullParamIdToReduced<D>(normalVecId(),
230 system().domain().lattice());
231 double L = system().domain().unitCell().parameter(paramId);
232
233 // Create a 3 element vector 'dim' that contains the grid dimensions.
234 // If system is 2D (1D), then the z (y & z) dimensions are set to 1.
235 IntVec<3> dim;
236 for (int ind = 0; ind < 3; ind++) {
237 if (ind < D) {
238 dim[ind] = system().domain().mesh().dimensions()[ind];
239 } else {
240 dim[ind] = 1;
241 }
242 }
243
244 // Generate an r-grid representation of the walls
245 RField<D> rGrid;
246 rGrid.allocate(system().domain().mesh().dimensions());
247 int x, y, z;
248 int counter = 0;
249 FArray<int,3> coords;
250 double d, rhoW;
251
252 for (x = 0; x < dim[0]; x++) {
253 coords[0] = x;
254 for (y = 0; y < dim[1]; y++) {
255 coords[1] = y;
256 for (z = 0; z < dim[2]; z++) {
257 coords[2] = z;
258
259 // Get the distance 'd' traveled along the lattice basis
260 // vector that is normal to the walls
261 d = coords[normalVecId()] * L / dim[normalVecId()];
262
263 // Calculate wall volume fraction (rhoW) at gridpoint (x,y,z)
264 rhoW = 0.5 * (1 + tanh(4 * (((.5 * (excludedThickness()-L)) +
265 fabs(d - (L/2))) / interfaceThickness())));
266 rGrid[counter++] = 1-rhoW;
267 }
268 }
269 }
270
271 // Store this mask in System
272 system().mask().setRGrid(rGrid, system().iterator().isSymmetric());
273
274 // Store lattice vector normal to film used to construct this mask
275 normalVecCurrent_ = systemLatticeVector(normalVecId());
276
277 }
278
279 /*
280 * Sets flexible lattice parameters to be compatible with the mask.
281 */
282 template <int D>
284 {
285 if (system().iterator().isFlexible()) {
286 FSArray<bool,6> updated;
287 updated = modifyFlexibleParams(system().iterator().flexibleParams(),
288 system().domain().unitCell());
289 sysPtr_->iterator().setFlexibleParams(updated);
290 }
291 }
292
293} // namespace Rpc
294} // namespace Pscf
295#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.
Base class defining mask that imposes thin film confinement.
void setFlexibleParams() const
Modifies iterator().flexibleParams_ to be compatible with the mask.
void setClassName(const char *className)
Set class name string.
double stressTerm(int paramId) const
Get contribution to the stress from this mask.
void generate()
Generate the field and store where the Iterator can access.
void allocate()
Allocate container necessary to generate and store field.
double modifyStress(int paramId, double stress) const
Modify stress value in direction normal to the film.
Main class for SCFT or PS-FTS simulation of one system.
Definition rpc/System.h:100
A fixed size (static) contiguous array template.
Definition FArray.h:47
A fixed capacity (static) contiguous array with a variable logical size.
Definition rpg/System.h:28
#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
int convertFullParamIdToReduced(const int fullId, const typename UnitCell< D >::LatticeSystem lattice)
Convert full lattice parameter index to a reduced index.
Fields and FFTs for periodic boundary conditions (CPU)
Definition CField.cpp:12
Periodic fields and crystallography.
Definition CField.cpp:11
PSCF package top-level namespace.
Definition param_pc.dox:1
Utility classes for scientific computation.