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