PSCF v1.3
rpg/environment/FilmFieldGenMask.tpp
1#ifndef RPG_MASK_GEN_FILM_TPP
2#define RPG_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 <rpg/field/FieldIo.h>
13#include <rpg/scft/iterator/Iterator.h>
14#include <rpg/scft/ScftThermo.h>
15#include <prdc/cpu/RField.h>
16#include <prdc/crystal/UnitCell.h>
17#include <prdc/crystal/paramIdConversions.h>
18#include <prdc/cuda/resources.h>
19#include <pscf/inter/Interaction.h>
20#include <cmath>
21
22namespace Pscf {
23namespace Rpg
24{
25
26 using namespace Util;
27 using namespace Pscf::Prdc;
28 using namespace Pscf::Prdc::Cuda;
29
30 // CUDA kernels:
31 // (defined in anonymous namespace, used only in this file)
32
33 namespace {
34
35 /*
36 * Compute the derivative of mask field with respect to film thickness.
37 *
38 * (note: indices are inverted from their RField representation. So,
39 * for a 3D system, x, y, and z in an RField become z, y, and x in the
40 * CUDA thread grid, respectively. In 2D, x, y becomes y, x. This
41 * allows for memory coalescing; see ThreadMesh::setConfig
42 * documentation for more information.)
43 *
44 * \param deriv field containing the mask derivative (output)
45 * \param mask field containing the mask
46 * \param interfaceThickness thickness of the polymer/wall interface
47 * \param normalVecId index of the lattice vector normal to the film
48 * \param meshDims grid dimensions of the field in real space
49 */
50 template <int D>
51 __global__ void _maskDerivative(cudaReal * deriv,
52 cudaReal const * mask,
53 cudaReal const interfaceThickness,
54 int const normalVecId,
55 dim3 const meshDims)
56 {
57 // Get position in 3-dimensional grid of threads
58 int ix = blockIdx.x * blockDim.x + threadIdx.x;
59 int iy = blockIdx.y * blockDim.y + threadIdx.y;
60 int iz = blockIdx.z * blockDim.z + threadIdx.z;
61
62 // Get thread index in linearized array
63 int id = ix;
64 if (D > 1) id += meshDims.x * iy;
65 if (D > 2) id += meshDims.x * meshDims.y * iz;
66
67 // Get position d along normalVecId in reduced coordinates
68 cudaReal d;
69 if (D == 3) {
70 if (normalVecId == 2) d = (cudaReal)ix / (cudaReal)meshDims.x;
71 if (normalVecId == 1) d = (cudaReal)iy / (cudaReal)meshDims.y;
72 if (normalVecId == 0) d = (cudaReal)iz / (cudaReal)meshDims.z;
73 } else if (D == 2) {
74 if (normalVecId == 1) d = (cudaReal)ix / (cudaReal)meshDims.x;
75 if (normalVecId == 0) d = (cudaReal)iy / (cudaReal)meshDims.y;
76 } else { // D == 1
77 d = (cudaReal)ix / (cudaReal)meshDims.x;
78 }
79
80 if (ix < meshDims.x && iy < meshDims.y && iz < meshDims.z) {
81 // Get mask value at this gridpoint from global memory
82 cudaReal maskVal = mask[id];
83
84 // Calculate d\phi_m / dL at this gridpoint
85 deriv[id] = maskVal * (maskVal - 1) * 8.0 * (abs(d - 0.5) - 0.5)
86 / interfaceThickness;
87 }
88 }
89
90 /*
91 * Generate the mask.
92 *
93 * \param mask field containing the mask (output)
94 * \param nvLength length of the lattice vector normal to the film
95 * \param excludedThickness thickness of region occupied by the wall
96 * \param interfaceThickness thickness of the polymer/wall interface
97 * \param normalVecId index of the lattice vector normal to the film
98 * \param meshDims grid dimensions of the field in real space
99 */
100 template <int D>
101 __global__ void _generateMask(cudaReal * mask, cudaReal const nvLength,
102 cudaReal const interfaceThickness,
103 cudaReal const excludedThickness,
104 int const normalVecId, dim3 const meshDims)
105 {
106 // Get position in 3-dimensional grid of threads
107 int ix = blockIdx.x * blockDim.x + threadIdx.x;
108 int iy = blockIdx.y * blockDim.y + threadIdx.y;
109 int iz = blockIdx.z * blockDim.z + threadIdx.z;
110
111 // Get thread index in linearized array
112 int id = ix;
113 if (D > 1) id += meshDims.x * iy;
114 if (D > 2) id += meshDims.x * meshDims.y * iz;
115
116 // Get position d along normalVecId in reduced coordinates
117 cudaReal d;
118 if (D == 3) {
119 if (normalVecId == 2) d = (cudaReal)ix / (cudaReal)meshDims.x;
120 if (normalVecId == 1) d = (cudaReal)iy / (cudaReal)meshDims.y;
121 if (normalVecId == 0) d = (cudaReal)iz / (cudaReal)meshDims.z;
122 } else if (D == 2) {
123 if (normalVecId == 1) d = (cudaReal)ix / (cudaReal)meshDims.x;
124 if (normalVecId == 0) d = (cudaReal)iy / (cudaReal)meshDims.y;
125 } else { // D == 1
126 d = (cudaReal)ix / (cudaReal)meshDims.x;
127 }
128
129 // Convert d to unreduced coordinates
130 d *= nvLength;
131
132 if (ix < meshDims.x && iy < meshDims.y && iz < meshDims.z) {
133 // Calculate wall volume fraction at this gridpoint
134 mask[id] =
135 (0.5 * (1.0 - tanh(
136 (
137 0.5 * (excludedThickness - nvLength) +
138 abs(d - (nvLength / 2.0))
139 )
140 * 4.0 / interfaceThickness
141 )));
142 }
143 }
144
145 }
146
147 /*
148 * Default constructor
149 */
150 template <int D>
153 sysPtr_(0)
154 { setClassName("FilmFieldGenMask"); }
155
156 /*
157 * Constructor
158 */
159 template <int D>
162 sysPtr_(&sys)
163 { setClassName("FilmFieldGenMask"); }
164
165 /*
166 * Destructor
167 */
168 template <int D>
171
172 /*
173 * Get contribution to the stress from this mask
174 */
175 template <int D>
176 double FilmFieldGenMask<D>::stress(int paramId) const
177 {
178 // If paramId is not normalVecId, there is no stress contribution
179 UTIL_CHECK(normalVecId() >= 0);
180 int normalVecParamId = convertFullParamIdToReduced<D>(normalVecId(),
181 system().domain().lattice());
182 if (normalVecParamId != paramId) return 0.0;
183
184 // If this point is reached, stress contribution must be calculated
185 UTIL_CHECK(system().mask().hasData());
186
187 // Get the length of the lattice basis vector normal to the walls
188 double nvLength = system().domain().unitCell().parameter(paramId);
189
190 // Get the volume fraction of the unit cell occupied by polymers
191 double phiTot = system().mask().phiTot();
192
193 // Create the derivative field in rgrid format
194 RField<D> deriv;
195 deriv.allocate(system().domain().mesh().dimensions());
196
197 // Set D-dimensional GPU configuration
198 ThreadMesh::setConfig(system().domain().mesh().dimensions(), true);
199 dim3 gridDims = ThreadMesh::gridDims();
200 dim3 blockDims = ThreadMesh::blockDims();
201
202 // Compute derivative of the mask with respect to film thickness
203 _maskDerivative<D><<<gridDims, blockDims>>>
204 (deriv.cArray(), system().mask().rgrid().cArray(),
206
207 // Get xi, the Lagrange multiplier field, in rgrid format
208 RField<D> xi;
209 xi.allocate(system().domain().mesh().dimensions());
210
211 if (system().h().hasData()) {
212 VecOp::subVV(xi, system().w().rgrid(0), system().h().rgrid(0));
213 } else {
214 VecOp::eqV(xi, system().w().rgrid(0));
215 }
216
217 int nMonomer = system().mixture().nMonomer();
218 for (int i = 0; i < nMonomer; i++) {
219 double chi = system().interaction().chi(0, i);
220 if (fabs(chi) > 1e-6) { // if chi is nonzero
221 VecOp::addEqVc(xi, system().c().rgrid(i), -1.0 * chi);
222 }
223 }
224
225 // Integrate deriv * xi to get the integral term in the stress
226 double intTerm = Reduce::innerProduct(deriv, xi);
227 intTerm /= phiTot * deriv.capacity();
228
229 // Get the pressure term in the stress
230 if (!sysPtr_->scft().hasData()) {
231 sysPtr_->scft().compute();
232 }
233 double pSys = sysPtr_->scft().pressure();
234 double pTerm = pSys * excludedThickness() /
235 (phiTot * nvLength * nvLength);
236
237 return pTerm - intTerm;
238 }
239
240 template <int D>
241 double FilmFieldGenMask<D>::modifyStress(int paramId, double stress)
242 const
243 {
245 system().domain().lattice());
246
247 // If paramId is not normalVecId, there is no stress modification
248 if (nvParamId != paramId) return stress;
249
250 // If this point is reached, stress must be modified
251 UTIL_CHECK(system().mask().hasData());
252
253 if (!hasFBulk()) {
254 UTIL_THROW("fBulk must be set before calling modifyStress.");
255 }
256
257 // Get system free energy
258 if (!sysPtr_->scft().hasData()) {
259 sysPtr_->scft().compute();
260 }
261 double fSys = sysPtr_->scft().fHelmholtz();
262
263 // Get the length L of the basis vector normal to the walls
264 double L = system().domain().unitCell().parameter(paramId);
265
266 double modifiedStress = (stress * (L - excludedThickness()))
267 + (fSys - fBulk_);
268
269 return modifiedStress;
270 }
271
272 /*
273 * Compute the field and store where the System can access
274 */
275 template <int D>
277 {
278 UTIL_CHECK(normalVecId() >= 0);
281 if (system().iterator().isSymmetric()) {
282 UTIL_CHECK(system().domain().basis().isInitialized());
283 }
284
285 // Get the length of the lattice basis vector normal to the walls
286 int paramId;
288 system().domain().lattice());
289 double nvLength = system().domain().unitCell().parameter(paramId);
290
291 // Setup
292 RField<D> rGrid;
293 rGrid.allocate(system().domain().mesh().dimensions());
294
295 // Set D-dimensional GPU configuration
296 ThreadMesh::setConfig(system().domain().mesh().dimensions(), true);
297 dim3 gridDims = ThreadMesh::gridDims();
298 dim3 blockDims = ThreadMesh::blockDims();
299
300 // Generate the mask
301 _generateMask<D><<<gridDims, blockDims>>>
302 (rGrid.cArray(), nvLength, interfaceThickness(),
304
305 // Store this mask in System
306 system().mask().setRGrid(rGrid,system().iterator().isSymmetric());
307
308 // Store lattice vector normal to film used to construct this mask
310 }
311
312 /*
313 * Sets flexible lattice parameters to be compatible with the mask.
314 */
315 template <int D>
317 {
318 if (system().iterator().isFlexible()) {
319 FSArray<bool,6> updated;
320 updated = modifyFlexibleParams(system().iterator().flexibleParams(),
321 system().domain().unitCell());
322 sysPtr_->iterator().setFlexibleParams(updated);
323 }
324 }
325
326} // namespace Rpg
327} // namespace Pscf
328#endif
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.
void compute()
Compute the field and store where the System can access.
double modifyStress(int paramId, double stress) const
Modify stress value in direction normal to the film.
void setClassName(const char *className)
Set class name string.
System< D > & system()
Get the System associated with this object by reference.
void setFlexibleParams() const
Modifies iterator().flexibleParams_ to be compatible with the mask.
double stress(int paramId) const
Get contribution to the stress from this mask.
RealVec< D > systemLatticeVector(int id) const
Get one of the lattice vectors for this system.
Main class, representing one complete system.
Data * cArray()
Return a pointer to the underlying C array.
Definition Array.h:214
int capacity() const
Return allocated size.
Definition Array.h:159
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.
void setConfig(IntVec< D > const &meshDims, bool invert, int blockSize)
Given a multidimensional grid of threads, set execution configuration.
Definition ThreadMesh.cu:88
dim3 gridDims()
Get the multidimensional grid of blocks determined by setConfig.
dim3 meshDims()
Return last requested multidimensional grid of threads.
dim3 blockDims()
Get the dimensions of a single block determined by setConfig.
RT abs(CT const &a)
Return absolute magnitude of a complex number.
void eqV(DeviceArray< cudaReal > &a, DeviceArray< cudaReal > const &b, const int beginIdA, const int beginIdB, const int n)
Vector assignment, a[i] = b[i], kernel wrapper (cudaReal).
Definition VecOp.cu:1039
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 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
Fields, FFTs, and utilities for periodic boundary conditions (CUDA)
Definition Reduce.cpp:14
Periodic fields and crystallography.
Definition CField.cpp:11
SCFT and PS-FTS with real periodic fields (GPU)
PSCF package top-level namespace.
Definition param_pc.dox:1