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