PSCF v1.2
rpg/scft/iterator/ExtGenFilm.tpp
1#ifndef RPG_EXT_GEN_FILM_TPP
2#define RPG_EXT_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 "ExtGenFilm.h"
12#include "Iterator.h"
13#include <rpg/field/FieldIo.h>
14#include <prdc/cpu/RField.h>
15#include <prdc/crystal/paramIdConversions.h>
16#include <prdc/cuda/resources.h>
17#include <cmath>
18
19namespace Pscf {
20namespace Rpg
21{
22
23 using namespace Util;
24 using namespace Pscf::Prdc;
25 using namespace Pscf::Prdc::Cuda;
26
27 // CUDA kernels:
28 // (defined in anonymous namespace, used only in this file)
29
30 namespace {
31
32 /*
33 * Compute the derivative of h fields with respect to film thickness.
34 *
35 * (note: indices are inverted from their RField representation. So,
36 * for a 3D system, x, y, and z in an RField become z, y, and x in the
37 * CUDA thread grid, respectively. In 2D, x, y becomes y, x. This
38 * allows for memory coalescing; see ThreadMesh::setConfig
39 * documentation for more information.)
40 *
41 * \param hDerivs array of fields containing the h derivatives (output)
42 * \param mask field containing the mask
43 * \param chiBottom array of chiBottom values for each species
44 * \param chiTop array of chiTop values for each species
45 * \param interfaceThickness thickness of the polymer/wall interface
46 * \param normalVecId index of the lattice vector normal to the film
47 * \param nMonomer number of monomer species in the system
48 * \param meshDims grid dimensions of the field in real space
49 */
50 template <int D>
51 __global__ void _hFieldDerivatives(cudaReal ** hDerivs,
52 cudaReal const * mask,
53 cudaReal const * chiBottom,
54 cudaReal const * chiTop,
55 cudaReal const interfaceThickness,
56 int const normalVecId,
57 int const nMonomer,
58 dim3 const meshDims)
59 {
60 // Get position in 3-dimensional grid of threads
61 int ix = blockIdx.x * blockDim.x + threadIdx.x;
62 int iy = blockIdx.y * blockDim.y + threadIdx.y;
63 int iz = blockIdx.z * blockDim.z + threadIdx.z;
64
65 // Get thread index in linearized array
66 int id = ix;
67 if (D > 1) id += meshDims.x * iy;
68 if (D > 2) id += meshDims.x * meshDims.y * iz;
69
70 // Get position d along normalVecId in reduced coordinates
71 cudaReal d;
72 if (D == 3) {
73 if (normalVecId == 2) d = (cudaReal)ix / (cudaReal)meshDims.x;
74 if (normalVecId == 1) d = (cudaReal)iy / (cudaReal)meshDims.y;
75 if (normalVecId == 0) d = (cudaReal)iz / (cudaReal)meshDims.z;
76 } else if (D == 2) {
77 if (normalVecId == 1) d = (cudaReal)ix / (cudaReal)meshDims.x;
78 if (normalVecId == 0) d = (cudaReal)iy / (cudaReal)meshDims.y;
79 } else { // D == 1
80 d = (cudaReal)ix / (cudaReal)meshDims.x;
81 }
82
83 if (ix < meshDims.x && iy < meshDims.y && iz < meshDims.z) {
84 // Get mask value at this gridpoint from global memory
85 cudaReal dPhidL = mask[id];
86
87 // Calculate d\phi_m / dL at this gridpoint
88 dPhidL = dPhidL * (dPhidL - 1) * 8.0 * (abs(d - 0.5) - 0.5)
89 / interfaceThickness;
90
91 // Calculate dh_i / dL at this gridpoint for each species i
92 for (int in = 0; in < nMonomer; in++) {
93 if (d < 0.5) {
94 hDerivs[in][id] = -1.0 * dPhidL * chiBottom[in];
95 } else {
96 hDerivs[in][id] = -1.0 * dPhidL * chiTop[in];
97 }
98 }
99 }
100 }
101
102 /*
103 * Generate the h fields.
104 *
105 * \param hFields array of h fields (output)
106 * \param mask field containing the mask
107 * \param chiBottom array of chiBottom values for each species
108 * \param chiTop array of chiTop values for each species
109 * \param normalVecId index of the lattice vector normal to the film
110 * \param nMonomer number of monomer species in the system
111 * \param meshDims grid dimensions of the field in real space
112 */
113 template <int D>
114 __global__ void _generateHFields(cudaReal ** hFields,
115 cudaReal const * mask,
116 cudaReal const * chiBottom,
117 cudaReal const * chiTop,
118 int const normalVecId,
119 int const nMonomer,
120 dim3 const meshDims)
121 {
122 // Get position in 3-dimensional grid of threads
123 int ix = blockIdx.x * blockDim.x + threadIdx.x;
124 int iy = blockIdx.y * blockDim.y + threadIdx.y;
125 int iz = blockIdx.z * blockDim.z + threadIdx.z;
126
127 // Get thread index in linearized array
128 int id = ix;
129 if (D > 1) id += meshDims.x * iy;
130 if (D > 2) id += meshDims.x * meshDims.y * iz;
131
132 // Determine if this grid point is in the "top" half
133 bool topHalf;
134 if (D == 3) {
135 if (normalVecId == 2) topHalf = (ix > meshDims.x / 2);
136 if (normalVecId == 1) topHalf = (iy > meshDims.y / 2);
137 if (normalVecId == 0) topHalf = (iz > meshDims.z / 2);
138 } else if (D == 2) {
139 if (normalVecId == 1) topHalf = (ix > meshDims.x / 2);
140 if (normalVecId == 0) topHalf = (iy > meshDims.y / 2);
141 } else { // D == 1
142 topHalf = (ix > meshDims.x / 2);
143 }
144
145 if (ix < meshDims.x && iy < meshDims.y && iz < meshDims.z) {
146 // Get mask value at this gridpoint from global memory
147 cudaReal maskVal = mask[id];
148
149 // Calculate h field value at this gridpoint for each species i
150 for (int in = 0; in < nMonomer; in++) {
151 if (topHalf) {
152 hFields[in][id] = (1.0 - maskVal) * chiTop[in];
153 } else {
154 hFields[in][id] = (1.0 - maskVal) * chiBottom[in];
155 }
156 }
157 }
158 }
159
160 }
161
162 /*
163 * Default constructor
164 */
165 template <int D>
168 sysPtr_(0)
169 { setClassName("ExtGenFilm"); }
170
171 /*
172 * Constructor
173 */
174 template <int D>
177 sysPtr_(&sys)
178 { setClassName("ExtGenFilm"); }
179
180 /*
181 * Destructor
182 */
183 template <int D>
186
187 /*
188 * Get contribution to the stress from these external fields
189 */
190 template <int D>
191 double ExtGenFilm<D>::stressTerm(int paramId) const
192 {
193 // If walls are athermal then there is no external field, so no
194 // contribution to the stress.
195 if (isAthermal()) return 0.0;
196
197 // If paramId is not normalVecId, there is no stress contribution
198 UTIL_CHECK(normalVecId() >= 0); // normalVecId has been set
199 int nvParamId = convertFullParamIdToReduced<D>(normalVecId(),
200 system().domain().lattice());
201 if (nvParamId != paramId) return 0.0;
202
203 // If this point is reached, calculate the stress contribution
204 // from the external fields.
205 UTIL_CHECK(isGenerated());
206 UTIL_CHECK(interfaceThickness() > 0);
207 UTIL_CHECK(system().hasMask());
208 UTIL_CHECK(system().hasExternalFields());
209
210 // Setup
211 int nMonomer = system().mixture().nMonomer();
212
213 DArray< RField<D> > hDerivatives;
214 HostDArray<cudaReal*> hDerivPtrs_h(nMonomer);
215 DeviceArray<cudaReal*> hDerivPtrs(nMonomer);
216
217 hDerivatives.allocate(nMonomer);
218 for (int i = 0; i < nMonomer; i++) {
219 hDerivatives[i].allocate(system().mesh().dimensions());
220 hDerivPtrs_h[i] = hDerivatives[i].cArray();
221 }
222 hDerivPtrs = hDerivPtrs_h;
223
224 HostDArray<cudaReal> chiBottom_h(nMonomer), chiTop_h(nMonomer);
225 DeviceArray<cudaReal> chiBottom_d(nMonomer), chiTop_d(nMonomer);
226 for (int i = 0; i < nMonomer; i++) {
227 chiBottom_h[i] = chiBottom(i);
228 chiTop_h[i] = chiTop(i);
229 }
230 chiBottom_d = chiBottom_h;
231 chiTop_d = chiTop_h;
232
233 // Set D-dimensional GPU configuration
234 ThreadMesh::setConfig(system().mesh().dimensions(), true);
235 dim3 gridDims = ThreadMesh::gridDims();
236 dim3 blockDims = ThreadMesh::blockDims();
237
238 // Compute the derivative of h fields with respect to film thickness
239 _hFieldDerivatives<D><<<gridDims, blockDims>>>
240 (hDerivPtrs.cArray(), system().mask().rgrid().cArray(),
241 chiBottom_d.cArray(), chiTop_d.cArray(), interfaceThickness(),
242 normalVecId(), nMonomer, ThreadMesh::meshDims());
243
244 // Integrate resulting fields to get the stress
245 double stress;
246 for (int i = 0; i < nMonomer; i++) {
247 stress += Reduce::innerProduct(hDerivatives[i], system().c().rgrid(i));
248 }
249 stress /= system().mask().phiTot() * hDerivatives[0].capacity();
250 return stress;
251 }
252
253 /*
254 * Allocate container necessary to generate and store field
255 */
256 template <int D>
258 {
259 UTIL_CHECK(system().unitCell().isInitialized());
260
261 // Make sure h field container has access to a fieldIo
262 system().h().setFieldIo(system().domain().fieldIo());
263
264 // Allocate the external field containers if needed
265 if (!isAthermal()) {
266 if (!system().h().isAllocatedRGrid()) {
267 system().h().allocateRGrid(system().mesh().dimensions());
268 }
269 if (system().iterator().isSymmetric()) {
270 UTIL_CHECK(system().basis().isInitialized());
271 if (!system().h().isAllocatedBasis()) {
272 system().h().allocateBasis(system().basis().nBasis());
273 }
274 }
275 }
276 }
277
281 template <int D>
283 {
284 // If walls are athermal then there is no external field needed.
285 // If an external field already exists in the System, we need to
286 // overwrite it with a field of all zeros, otherwise do nothing
287 if ((isAthermal()) && (!isGenerated())) return;
288
289 // If this point is reached, external field must be generated
290 UTIL_CHECK(system().h().isAllocatedRGrid());
291 UTIL_CHECK(system().mask().isAllocatedRGrid());
292 if (system().iterator().isSymmetric()) {
293 UTIL_CHECK(system().h().isAllocatedBasis());
294 UTIL_CHECK(system().mask().isAllocatedBasis());
295 }
296
297 UTIL_CHECK(normalVecId() >= 0);
298
299 // Set chiBottomCurrent_, chiTopCurrent_, and parametersCurrent_
300 chiBottomCurrent_ = chiBottom();
301 chiTopCurrent_ = chiTop();
302 normalVecCurrent_ = systemLatticeVector(normalVecId());
303
304 // Setup
305 int nMonomer = system().mixture().nMonomer();
306
307 DArray< RField<D> > hFields;
308 HostDArray<cudaReal*> hPtrs_h(nMonomer);
309 DeviceArray<cudaReal*> hPtrs(nMonomer);
310
311 hFields.allocate(nMonomer);
312 for (int i = 0; i < nMonomer; i++) {
313 hFields[i].allocate(system().mesh().dimensions());
314 hPtrs_h[i] = hFields[i].cArray();
315 }
316 hPtrs = hPtrs_h;
317
318 HostDArray<cudaReal> chiBottom_h(nMonomer), chiTop_h(nMonomer);
319 DeviceArray<cudaReal> chiBottom_d(nMonomer), chiTop_d(nMonomer);
320 for (int i = 0; i < nMonomer; i++) {
321 chiBottom_h[i] = chiBottom(i);
322 chiTop_h[i] = chiTop(i);
323 }
324 chiBottom_d = chiBottom_h;
325 chiTop_d = chiTop_h;
326
327 // Set D-dimensional GPU configuration
328 ThreadMesh::setConfig(system().mesh().dimensions(), true);
329 dim3 gridDims = ThreadMesh::gridDims();
330 dim3 blockDims = ThreadMesh::blockDims();
331
332 // Generate the h fields
333 _generateHFields<D><<<gridDims, blockDims>>>
334 (hPtrs.cArray(), system().mask().rgrid().cArray(),
335 chiBottom_d.cArray(), chiTop_d.cArray(), normalVecId(),
336 nMonomer, ThreadMesh::meshDims());
337
338 // Pass h into the System
339 system().h().setRGrid(hFields, system().iterator().isSymmetric());
340 }
341
342}
343}
344
345#endif
Dynamic array on the GPU device with aligned data.
Definition rpg/System.h:32
Data * cArray()
Return pointer to underlying C array.
Template for dynamic array stored in host CPU memory.
Definition HostDArray.h:43
Base class defining external fields for thin film systems.
void generate()
Generate the fields and store where the Iterator can access.
void setClassName(const char *className)
Set class name string.
double stressTerm(int paramId) const
Get contribution to the stress from the external fields.
void allocate()
Allocate container necessary to generate and store fields.
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
Dynamically allocatable contiguous array template.
void allocate(int capacity)
Allocate the underlying C array.
Definition DArray.h:199
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
Definition global.h:68
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
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.