1#ifndef RPG_EXT_GEN_FILM_TPP
2#define RPG_EXT_GEN_FILM_TPP
11#include "ExtGenFilm.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>
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,
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;
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;
77 if (normalVecId == 1) d = (cudaReal)ix / (cudaReal)
meshDims.x;
78 if (normalVecId == 0) d = (cudaReal)iy / (cudaReal)
meshDims.y;
80 d = (cudaReal)ix / (cudaReal)
meshDims.x;
85 cudaReal dPhidL = mask[id];
88 dPhidL = dPhidL * (dPhidL - 1) * 8.0 * (
abs(d - 0.5) - 0.5)
92 for (
int in = 0; in < nMonomer; in++) {
94 hDerivs[in][id] = -1.0 * dPhidL * chiBottom[in];
96 hDerivs[in][id] = -1.0 * dPhidL * chiTop[in];
114 __global__
void _generateHFields(cudaReal ** hFields,
115 cudaReal
const * mask,
116 cudaReal
const * chiBottom,
117 cudaReal
const * chiTop,
118 int const normalVecId,
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;
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);
139 if (normalVecId == 1) topHalf = (ix >
meshDims.x / 2);
140 if (normalVecId == 0) topHalf = (iy >
meshDims.y / 2);
147 cudaReal maskVal = mask[id];
150 for (
int in = 0; in < nMonomer; in++) {
152 hFields[in][id] = (1.0 - maskVal) * chiTop[in];
154 hFields[in][id] = (1.0 - maskVal) * chiBottom[in];
195 if (isAthermal())
return 0.0;
200 system().domain().lattice());
201 if (nvParamId != paramId)
return 0.0;
211 int nMonomer = system().mixture().nMonomer();
218 for (
int i = 0; i < nMonomer; i++) {
219 hDerivatives[i].
allocate(system().mesh().dimensions());
220 hDerivPtrs_h[i] = hDerivatives[i].
cArray();
222 hDerivPtrs = hDerivPtrs_h;
226 for (
int i = 0; i < nMonomer; i++) {
227 chiBottom_h[i] = chiBottom(i);
228 chiTop_h[i] = chiTop(i);
230 chiBottom_d = chiBottom_h;
239 _hFieldDerivatives<D><<<gridDims, blockDims>>>
240 (hDerivPtrs.
cArray(), system().mask().rgrid().cArray(),
241 chiBottom_d.cArray(), chiTop_d.
cArray(), interfaceThickness(),
246 for (
int i = 0; i < nMonomer; i++) {
249 stress /= system().mask().phiTot() * hDerivatives[0].
capacity();
259 UTIL_CHECK(system().unitCell().isInitialized());
262 system().h().setFieldIo(system().domain().fieldIo());
266 if (!system().h().isAllocatedRGrid()) {
267 system().h().allocateRGrid(system().mesh().dimensions());
269 if (system().iterator().isSymmetric()) {
271 if (!system().h().isAllocatedBasis()) {
272 system().h().allocateBasis(system().basis().nBasis());
287 if ((isAthermal()) && (!isGenerated()))
return;
291 UTIL_CHECK(system().mask().isAllocatedRGrid());
292 if (system().iterator().isSymmetric()) {
294 UTIL_CHECK(system().mask().isAllocatedBasis());
300 chiBottomCurrent_ = chiBottom();
301 chiTopCurrent_ = chiTop();
302 normalVecCurrent_ = systemLatticeVector(normalVecId());
305 int nMonomer = system().mixture().nMonomer();
312 for (
int i = 0; i < nMonomer; i++) {
313 hFields[i].
allocate(system().mesh().dimensions());
314 hPtrs_h[i] = hFields[i].
cArray();
320 for (
int i = 0; i < nMonomer; i++) {
321 chiBottom_h[i] = chiBottom(i);
322 chiTop_h[i] = chiTop(i);
324 chiBottom_d = chiBottom_h;
333 _generateHFields<D><<<gridDims, blockDims>>>
334 (hPtrs.
cArray(), system().mask().rgrid().cArray(),
335 chiBottom_d.cArray(), chiTop_d.
cArray(), normalVecId(),
339 system().h().setRGrid(hFields, system().iterator().isSymmetric());
Dynamic array on the GPU device with aligned data.
Data * cArray()
Return pointer to underlying C array.
Template for dynamic array stored in host CPU memory.
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.
ExtGenFilm()
Default constructor.
Main class for calculations that represent one system.
Data * cArray()
Return a pointer to the underlying C array.
int capacity() const
Return allocated size.
Dynamically allocatable contiguous array template.
void allocate(int capacity)
Allocate the underlying C array.
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
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.
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).
Fields, FFTs, and utilities for periodic boundary conditions (CUDA)
Periodic fields and crystallography.
PSCF package top-level namespace.
Utility classes for scientific computation.