PSCF v1.3
rpc/environment/FilmFieldGenMask.tpp
1#ifndef RPC_MASK_GEN_FILM_TPP
2#define RPC_MASK_GEN_FILM_TPP
3
4/*
5* PSCF - Polymer Self-Consistent Field
6*
7* Copyright 2015 - 2025, The Regents of the University of Minnesota
8* Distributed under the terms of the GNU General Public License.
9*/
10
11#include "FilmFieldGenMask.h"
12#include <rpc/scft/iterator/Iterator.h>
13#include <rpc/scft/ScftThermo.h>
14#include <rpc/field/FieldIo.h>
15#include <prdc/cpu/RField.h>
16#include <prdc/crystal/UnitCell.h>
17#include <prdc/crystal/paramIdConversions.h>
18#include <pscf/math/IntVec.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_(nullptr)
36 { setClassName("FilmFieldGenMask"); }
37
38 /*
39 * Constructor
40 */
41 template <int D>
44 sysPtr_(&sys)
45 { setClassName("FilmFieldGenMask"); }
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 FilmFieldGenMask<D>::stress(int paramId) const
59 {
60 UTIL_CHECK(sysPtr_);
61 Domain<D> const & domain = system().domain();
62 Mesh<D> const & mesh = domain.mesh();
63
64 int normalVecParamId = convertFullParamIdToReduced<D>(normalVecId(),
65 domain.lattice());
66
67 // If paramId is not normalVecId, there is no stress contribution
68 if (normalVecParamId != paramId) return 0.0;
69
70 // If this point is reached, stress contribution must be calculated
71 UTIL_CHECK(system().mask().hasData());
72
73 // Get the length of the lattice basis vector normal to the walls
74 double nvLength = domain.unitCell().parameter(paramId);
75
76 // Get the volume fraction of the unit cell occupied by polymers
77 double phiTot = system().mask().phiTot();
78
79 // Create a 3 element vector 'dim' that contains the grid dimensions.
80 // If system is 2D (1D), then the z (y & z) dimensions are set to 1.
81 IntVec<3> dim;
82 for (int ind = 0; ind < 3; ind++) {
83 if (ind < D) {
84 dim[ind] = mesh.dimensions()[ind];
85 } else {
86 dim[ind] = 1;
87 }
88 }
89
90 // Create the derivative field in rgrid format
91 RField<D> deriv;
92 IntVec<3> coords;
93 int x, y, z;
94 int counter = 0;
95 double d, maskVal;
96
97 deriv.allocate(mesh.dimensions());
98 RField<D> const & maskRGrid = system().mask().rgrid();
99
100 for (x = 0; x < dim[0]; x++) {
101 coords[0] = x;
102 for (y = 0; y < dim[1]; y++) {
103 coords[1] = y;
104 for (z = 0; z < dim[2]; z++) {
105 coords[2] = z;
106 maskVal = maskRGrid[counter];
107
108 // Get the distance 'd' traveled along the lattice basis
109 // vector normal to the walls, in reduced coordinates
110 d = (double)coords[normalVecId()] /
111 (double)dim[normalVecId()];
112
113 deriv[counter] = maskVal * (maskVal - 1) * 8.0
114 * (std::abs(d - 0.5) - 0.5)
116 counter++;
117 }
118 }
119 }
120
121 // Get xi, the Lagrange multiplier field, in rgrid format
122 int nMonomer = system().mixture().nMonomer();
123 int nx = mesh.size();
124 double chi;
125 RField<D> xi;
126 xi.allocate(mesh.dimensions());
127
128 for (int i = 0; i < nx; i++) {
129 xi[i] = system().w().rgrid(0)[i];
130 }
131
132 for (int in = 0; in < nMonomer; in++) {
133 chi = system().interaction().chi(0,in);
134 if (fabs(chi) > 1e-6) { // if chi is nonzero
135 RField<D> const & cfield = system().c().rgrid(in);
136 for (int i = 0; i < nx; i++) {
137 //xi[i] -= system().c().rgrid(in)[i] * chi;
138 xi[i] -= cfield[i] * chi;
139 }
140 }
141 }
142
143 if (system().h().hasData()) {
144 RField<D> const & hfield = system().h().rgrid(0);
145 for (int i = 0; i < nx; i++) {
146 //xi[i] -= system().h().rgrid(0)[i];
147 xi[i] -= hfield[i];
148 }
149 }
150
151 // Get the integral term in the stress
152 double intTerm = 0.0;
153 for (int i = 0; i < nx; i++) {
154 intTerm += xi[i] * deriv[i];
155 }
156 intTerm /= (phiTot * nx);
157
158 // Get the pressure term in the stress
159 if (!system().scft().hasData()) {
160 sysPtr_->scft().compute();
161 }
162 double pSys = system().scft().pressure();
163 double pTerm = pSys * excludedThickness() /
164 (phiTot * nvLength * nvLength);
165
166 return pTerm - intTerm;
167 }
168
169 template <int D>
170 double FilmFieldGenMask<D>::modifyStress(int paramId, double stress)
171 const
172 {
173 UTIL_CHECK(sysPtr_);
175 system().domain().lattice());
176 if (nvParamId == paramId) {
177 UTIL_CHECK(system().mask().hasData());
178
179 if (!hasFBulk()) {
180 UTIL_THROW("fBulk must be set before calculating this stress.");
181 }
182
183 // Get system free energy
184 if (!system().scft().hasData()) {
185 sysPtr_->scft().compute();
186 }
187 double fSys = system().scft().fHelmholtz();
188
189 // Get the length L of the basis vector normal to the walls
190 double L = system().domain().unitCell().parameter(paramId);
191
192 double modifiedStress = (stress * (L - excludedThickness()))
193 + (fSys - fBulk_);
194
195 return modifiedStress;
196 } else {
197 return stress;
198 }
199 }
200
201 /*
202 * Compute the field and store where the System can access
203 */
204 template <int D>
206 {
207 UTIL_CHECK(sysPtr_);
210 UTIL_CHECK(normalVecId() >= 0);
211 if (system().iterator().isSymmetric()) {
212 UTIL_CHECK(system().domain().basis().isInitialized());
213 }
214
215 // Get the length L of the lattice basis vector normal to the walls
217 system().domain().lattice());
218 double L = system().domain().unitCell().parameter(paramId);
219
220 // Create a 3 element vector 'dim' that contains the grid dimensions.
221 // If system is 2D (1D), then the z (y & z) dimensions are set to 1.
222 IntVec<3> dim;
223 for (int ind = 0; ind < 3; ind++) {
224 if (ind < D) {
225 dim[ind] = system().domain().mesh().dimensions()[ind];
226 } else {
227 dim[ind] = 1;
228 }
229 }
230
231 // Generate an r-grid representation of the walls
232 RField<D> rGrid;
233 rGrid.allocate(system().domain().mesh().dimensions());
234 int x, y, z;
235 int counter = 0;
236 IntVec<3> coords;
237 double d, rhoW;
238
239 for (x = 0; x < dim[0]; x++) {
240 coords[0] = x;
241 for (y = 0; y < dim[1]; y++) {
242 coords[1] = y;
243 for (z = 0; z < dim[2]; z++) {
244 coords[2] = z;
245
246 // Get the distance 'd' traveled along the lattice basis
247 // vector that is normal to the walls
248 d = coords[normalVecId()] * L / dim[normalVecId()];
249
250 // Calculate wall volume fraction (rhoW) at gridpoint (x,y,z)
251 rhoW = 0.5 * (1 + tanh(4 * (((.5 * (excludedThickness()-L)) +
252 fabs(d - (L/2))) / interfaceThickness())));
253 rGrid[counter++] = 1-rhoW;
254 }
255 }
256 }
257
258 // Store this mask in System
259 system().mask().setRGrid(rGrid, system().iterator().isSymmetric());
260
261 // Store lattice vector normal to film used to construct this mask
263
264 }
265
266 /*
267 * Sets flexible lattice parameters to be compatible with the mask.
268 */
269 template <int D>
271 {
272 UTIL_CHECK(sysPtr_);
273 if (system().iterator().isFlexible()) {
274 FSArray<bool,6> updated;
275 updated = modifyFlexibleParams(system().iterator().flexibleParams(),
276 system().domain().unitCell());
277 sysPtr_->iterator().setFlexibleParams(updated);
278 }
279 }
280
281} // namespace Rpc
282} // namespace Pscf
283#endif
An IntVec<D, T> is a D-component vector of elements of integer type T.
Definition IntVec.h:27
Description of a regular grid of points in a periodic domain.
Definition Mesh.h:61
IntVec< D > dimensions() const
Get an IntVec<D> of the grid dimensions.
Definition Mesh.h:217
int size() const
Get total number of grid points.
Definition Mesh.h:229
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.
RealVec< D > normalVecCurrent_
The lattice vector normal to the film used to generate these fields.
int normalVecId() const
Get value of normalVecId.
double fBulk_
Reference free energy used to calculate stress normal to the film.
bool hasFBulk() const
Check whether a value of fBulk was provided.
FSArray< bool, 6 > modifyFlexibleParams(FSArray< bool, 6 > current, UnitCell< D > const &cell) const
Modifies a flexibleParams array to be compatible with this mask.
double interfaceThickness() const
Get value of interfaceThickness.
double excludedThickness() const
Get value of excludedThickness.
Spatial domain for a periodic structure with real fields, on a CPU.
Mesh< D > & mesh()
Get the Mesh by non-const reference.
UnitCell< D > & unitCell()
Get the UnitCell by non-const reference.
UnitCell< D >::LatticeSystem lattice() const
Get the lattice system (enumeration value).
void setFlexibleParams() const
Modifies iterator().flexibleParams_ to be compatible with the mask.
void setClassName(const char *className)
Set class name string.
double modifyStress(int paramId, double stress) const
Modify stress value in direction normal to the film.
RealVec< D > systemLatticeVector(int id) const
Get one of the lattice vectors for this system.
System< D > & system()
Get the parent System by non-const reference.
double stress(int paramId) const
Get contribution to the stress from this mask.
void compute()
Compute the field and store where the System can access.
Main class, representing one complete system.
A fixed capacity (static) contiguous array with a variable logical size.
Definition FSArray.h:38
#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:49
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
Real periodic fields, SCFT and PS-FTS (CPU).
Definition param_pc.dox:2
PSCF package top-level namespace.
Definition param_pc.dox:1