PSCF v1.3.1
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/inter/Interaction.h>
19#include <pscf/math/IntVec.h>
20#include <cmath>
21
22namespace Pscf {
23namespace Rpc
24{
25
26 using namespace Util;
27 using namespace Pscf::Prdc;
28 using namespace Pscf::Prdc::Cpu;
29
30 /*
31 * Default constructor
32 */
33 template <int D>
36 sysPtr_(nullptr)
37 { setClassName("FilmFieldGenMask"); }
38
39 /*
40 * Constructor
41 */
42 template <int D>
45 sysPtr_(&sys)
46 { setClassName("FilmFieldGenMask"); }
47
48 /*
49 * Destructor
50 */
51 template <int D>
54
55 /*
56 * Get contribution to the stress from this mask
57 */
58 template <int D>
59 double FilmFieldGenMask<D>::stress(int paramId) const
60 {
61 UTIL_CHECK(sysPtr_);
62 Domain<D> const & domain = system().domain();
63 Mesh<D> const & mesh = domain.mesh();
64
65 int normalVecParamId = convertFullParamIdToReduced<D>(normalVecId(),
66 domain.lattice());
67
68 // If paramId is not normalVecId, there is no stress contribution
69 if (normalVecParamId != paramId) return 0.0;
70
71 // If this point is reached, stress contribution must be calculated
72 UTIL_CHECK(system().mask().hasData());
73
74 // Get the length of the lattice basis vector normal to the walls
75 double nvLength = domain.unitCell().parameter(paramId);
76
77 // Get the volume fraction of the unit cell occupied by polymers
78 double phiTot = system().mask().phiTot();
79
80 // Create a 3 element vector 'dim' that contains the grid dimensions.
81 // If system is 2D (1D), then the z (y & z) dimensions are set to 1.
82 IntVec<3> dim;
83 for (int ind = 0; ind < 3; ind++) {
84 if (ind < D) {
85 dim[ind] = mesh.dimensions()[ind];
86 } else {
87 dim[ind] = 1;
88 }
89 }
90
91 // Create the derivative field in rgrid format
92 RField<D> deriv;
93 IntVec<3> coords;
94 int x, y, z;
95 int counter = 0;
96 double d, maskVal;
97
98 deriv.allocate(mesh.dimensions());
99 RField<D> const & maskRGrid = system().mask().rgrid();
100
101 for (x = 0; x < dim[0]; x++) {
102 coords[0] = x;
103 for (y = 0; y < dim[1]; y++) {
104 coords[1] = y;
105 for (z = 0; z < dim[2]; z++) {
106 coords[2] = z;
107 maskVal = maskRGrid[counter];
108
109 // Get the distance 'd' traveled along the lattice basis
110 // vector normal to the walls, in reduced coordinates
111 d = (double)coords[normalVecId()] /
112 (double)dim[normalVecId()];
113
114 deriv[counter] = maskVal * (maskVal - 1) * 8.0
115 * (std::abs(d - 0.5) - 0.5)
117 counter++;
118 }
119 }
120 }
121
122 // Get xi, the Lagrange multiplier field, in rgrid format
123 int nMonomer = system().mixture().nMonomer();
124 int nx = mesh.size();
125 double chi;
126 RField<D> xi;
127 xi.allocate(mesh.dimensions());
128
129 for (int i = 0; i < nx; i++) {
130 xi[i] = system().w().rgrid(0)[i];
131 }
132
133 for (int in = 0; in < nMonomer; in++) {
134 chi = system().interaction().chi(0,in);
135 if (fabs(chi) > 1e-6) { // if chi is nonzero
136 RField<D> const & cfield = system().c().rgrid(in);
137 for (int i = 0; i < nx; i++) {
138 //xi[i] -= system().c().rgrid(in)[i] * chi;
139 xi[i] -= cfield[i] * chi;
140 }
141 }
142 }
143
144 if (system().h().hasData()) {
145 RField<D> const & hfield = system().h().rgrid(0);
146 for (int i = 0; i < nx; i++) {
147 //xi[i] -= system().h().rgrid(0)[i];
148 xi[i] -= hfield[i];
149 }
150 }
151
152 // Get the integral term in the stress
153 double intTerm = 0.0;
154 for (int i = 0; i < nx; i++) {
155 intTerm += xi[i] * deriv[i];
156 }
157 intTerm /= (phiTot * nx);
158
159 // Get the pressure term in the stress
160 if (!system().scft().hasData()) {
161 sysPtr_->scft().compute();
162 }
163 double pSys = system().scft().pressure();
164 double pTerm = pSys * excludedThickness() /
165 (phiTot * nvLength * nvLength);
166
167 return pTerm - intTerm;
168 }
169
170 template <int D>
171 double FilmFieldGenMask<D>::modifyStress(int paramId, double stress)
172 const
173 {
174 UTIL_CHECK(sysPtr_);
176 system().domain().lattice());
177 if (nvParamId == paramId) {
178 UTIL_CHECK(system().mask().hasData());
179
180 if (!hasFBulk()) {
181 UTIL_THROW("fBulk must be set before calculating this stress.");
182 }
183
184 // Get system free energy
185 if (!system().scft().hasData()) {
186 sysPtr_->scft().compute();
187 }
188 double fSys = system().scft().fHelmholtz();
189
190 // Get the length L of the basis vector normal to the walls
191 double L = system().domain().unitCell().parameter(paramId);
192
193 double modifiedStress = (stress * (L - excludedThickness()))
194 + (fSys - fBulk_);
195
196 return modifiedStress;
197 } else {
198 return stress;
199 }
200 }
201
202 /*
203 * Compute the field and store where the System can access
204 */
205 template <int D>
207 {
208 UTIL_CHECK(sysPtr_);
211 UTIL_CHECK(normalVecId() >= 0);
212 if (system().iterator().isSymmetric()) {
213 UTIL_CHECK(system().domain().basis().isInitialized());
214 }
215
216 // Get the length L of the lattice basis vector normal to the walls
218 system().domain().lattice());
219 double L = system().domain().unitCell().parameter(paramId);
220
221 // Create a 3 element vector 'dim' that contains the grid dimensions.
222 // If system is 2D (1D), then the z (y & z) dimensions are set to 1.
223 IntVec<3> dim;
224 for (int ind = 0; ind < 3; ind++) {
225 if (ind < D) {
226 dim[ind] = system().domain().mesh().dimensions()[ind];
227 } else {
228 dim[ind] = 1;
229 }
230 }
231
232 // Generate an r-grid representation of the walls
233 RField<D> rGrid;
234 rGrid.allocate(system().domain().mesh().dimensions());
235 int x, y, z;
236 int counter = 0;
237 IntVec<3> coords;
238 double d, rhoW;
239
240 for (x = 0; x < dim[0]; x++) {
241 coords[0] = x;
242 for (y = 0; y < dim[1]; y++) {
243 coords[1] = y;
244 for (z = 0; z < dim[2]; z++) {
245 coords[2] = z;
246
247 // Get the distance 'd' traveled along the lattice basis
248 // vector that is normal to the walls
249 d = coords[normalVecId()] * L / dim[normalVecId()];
250
251 // Calculate wall volume fraction (rhoW) at gridpoint (x,y,z)
252 rhoW = 0.5 * (1 + tanh(4 * (((.5 * (excludedThickness()-L)) +
253 fabs(d - (L/2))) / interfaceThickness())));
254 rGrid[counter++] = 1-rhoW;
255 }
256 }
257 }
258
259 // Store this mask in System
260 system().mask().setRGrid(rGrid, system().iterator().isSymmetric());
261
262 // Store lattice vector normal to film used to construct this mask
264
265 }
266
267 /*
268 * Sets flexible lattice parameters to be compatible with the mask.
269 */
270 template <int D>
272 {
273 UTIL_CHECK(sysPtr_);
274 if (system().iterator().isFlexible()) {
275 FSArray<bool,6> updated;
276 updated = modifyFlexibleParams(system().iterator().flexibleParams(),
277 system().domain().unitCell());
278 sysPtr_->iterator().setFlexibleParams(updated);
279 }
280 }
281
282} // namespace Rpc
283} // namespace Pscf
284#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 a complete physical 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