PSCF v1.3
rpg/environment/FilmFieldGenExt.tpp
1#ifndef RPG_EXT_GEN_FILM_TPP
2#define RPG_EXT_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 "FilmFieldGenExt.h"
12#include <rpg/field/FieldIo.h>
13#include <rpg/scft/iterator/Iterator.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("FilmFieldGenExt"); }
170
171 /*
172 * Constructor
173 */
174 template <int D>
177 sysPtr_(&sys)
178 { setClassName("FilmFieldGenExt"); }
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 FilmFieldGenExt<D>::stress(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
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.
206 UTIL_CHECK(system().mask().hasData());
207 UTIL_CHECK(system().h().hasData());
208
209 // Setup
210 int nMonomer = system().mixture().nMonomer();
211
212 DArray< RField<D> > hDerivatives;
213 HostDArray<cudaReal*> hDerivPtrs_h(nMonomer);
214 DeviceArray<cudaReal*> hDerivPtrs(nMonomer);
215
216 hDerivatives.allocate(nMonomer);
217 for (int i = 0; i < nMonomer; i++) {
218 hDerivatives[i].allocate(system().domain().mesh().dimensions());
219 hDerivPtrs_h[i] = hDerivatives[i].cArray();
220 }
221 hDerivPtrs = hDerivPtrs_h;
222
223 HostDArray<cudaReal> chiBottom_h(nMonomer), chiTop_h(nMonomer);
224 DeviceArray<cudaReal> chiBottom_d(nMonomer), chiTop_d(nMonomer);
225 for (int i = 0; i < nMonomer; i++) {
226 chiBottom_h[i] = chiBottom(i);
227 chiTop_h[i] = chiTop(i);
228 }
229 chiBottom_d = chiBottom_h;
230 chiTop_d = chiTop_h;
231
232 // Set D-dimensional GPU configuration
233 ThreadMesh::setConfig(system().domain().mesh().dimensions(), true);
234 dim3 gridDims = ThreadMesh::gridDims();
235 dim3 blockDims = ThreadMesh::blockDims();
236
237 // Compute the derivative of h fields with respect to film thickness
238 _hFieldDerivatives<D><<<gridDims, blockDims>>>
239 (hDerivPtrs.cArray(), system().mask().rgrid().cArray(),
240 chiBottom_d.cArray(), chiTop_d.cArray(), interfaceThickness(),
241 normalVecId(), nMonomer, ThreadMesh::meshDims());
242
243 // Integrate resulting fields to get the stress
244 double stress;
245 for (int i = 0; i < nMonomer; i++) {
246 stress += Reduce::innerProduct(hDerivatives[i], system().c().rgrid(i));
247 }
248 stress /= system().mask().phiTot() * hDerivatives[0].capacity();
249 return stress;
250 }
251
255 template <int D>
257 {
258 // Set chiBottomCurrent_, chiTopCurrent_, and parametersCurrent_
262
263 // If walls are athermal then there is no external field needed.
264 // If external fields already exist in System, clear them, then return.
265 if (isAthermal()) {
266 system().h().clear();
267 return;
268 }
269
270 // If this point is reached, external field must be generated
271 UTIL_CHECK(normalVecId() >= 0);
272 UTIL_CHECK(system().domain().unitCell().isInitialized());
273 UTIL_CHECK(system().mask().hasData());
274 if (system().iterator().isSymmetric()) {
275 UTIL_CHECK(system().domain().basis().isInitialized());
276 }
277
278 // Setup
279 int nMonomer = system().mixture().nMonomer();
280
281 DArray< RField<D> > hFields;
282 HostDArray<cudaReal*> hPtrs_h(nMonomer);
283 DeviceArray<cudaReal*> hPtrs(nMonomer);
284
285 hFields.allocate(nMonomer);
286 for (int i = 0; i < nMonomer; i++) {
287 hFields[i].allocate(system().domain().mesh().dimensions());
288 hPtrs_h[i] = hFields[i].cArray();
289 }
290 hPtrs = hPtrs_h;
291
292 HostDArray<cudaReal> chiBottom_h(nMonomer), chiTop_h(nMonomer);
293 DeviceArray<cudaReal> chiBottom_d(nMonomer), chiTop_d(nMonomer);
294 for (int i = 0; i < nMonomer; i++) {
295 chiBottom_h[i] = chiBottom(i);
296 chiTop_h[i] = chiTop(i);
297 }
298 chiBottom_d = chiBottom_h;
299 chiTop_d = chiTop_h;
300
301 // Set D-dimensional GPU configuration
302 ThreadMesh::setConfig(system().domain().mesh().dimensions(), true);
303 dim3 gridDims = ThreadMesh::gridDims();
304 dim3 blockDims = ThreadMesh::blockDims();
305
306 // Generate the h fields
307 _generateHFields<D><<<gridDims, blockDims>>>
308 (hPtrs.cArray(), system().mask().rgrid().cArray(),
309 chiBottom_d.cArray(), chiTop_d.cArray(), normalVecId(),
310 nMonomer, ThreadMesh::meshDims());
311
312 // Pass h into the System
313 system().h().setRGrid(hFields, system().iterator().isSymmetric());
314 }
315
316}
317}
318
319#endif
Dynamic array on the GPU device with aligned data.
Definition DeviceArray.h:43
Data * cArray()
Return pointer to underlying C array.
Template for dynamic array stored in host CPU memory.
Definition HostDArray.h:43
DArray< double > chiTopCurrent_
The chiTop array used to generate the current external fields.
DArray< double > const & chiTop() const
Get const chiTop array by reference.
DArray< double > chiBottomCurrent_
The chiBottom array used to generate the current external fields.
DArray< double > const & chiBottom() const
Get const chiBottom matrix by reference.
bool isAthermal() const
Are the walls athermal?
double interfaceThickness() const
Get value of interfaceThickness.
RealVec< D > normalVecCurrent_
The lattice vector normal to the film used to generate these fields.
int normalVecId() const
Get value of normalVecId.
void setClassName(const char *className)
Set class name string.
void compute()
Compute the fields and store where the System can access.
RealVec< D > systemLatticeVector(int id) const
Get one of the lattice vectors for this system.
double stress(int paramId) const
Get contribution to the stress from the external fields.
System< D > & system()
Get the System associated with this object by reference.
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
Dynamically allocatable contiguous array template.
Definition DArray.h:32
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.
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