1#ifndef RPG_EXT_GEN_FILM_TPP
2#define RPG_EXT_GEN_FILM_TPP
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>
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];
200 system().domain().lattice());
201 if (nvParamId != paramId)
return 0.0;
210 int nMonomer =
system().mixture().nMonomer();
217 for (
int i = 0; i < nMonomer; i++) {
219 hDerivPtrs_h[i] = hDerivatives[i].
cArray();
221 hDerivPtrs = hDerivPtrs_h;
225 for (
int i = 0; i < nMonomer; i++) {
229 chiBottom_d = chiBottom_h;
238 _hFieldDerivatives<D><<<gridDims, blockDims>>>
245 for (
int i = 0; i < nMonomer; i++) {
246 stress += Reduce::innerProduct(hDerivatives[i],
system().c().rgrid(i));
274 if (
system().iterator().isSymmetric()) {
279 int nMonomer =
system().mixture().nMonomer();
286 for (
int i = 0; i < nMonomer; i++) {
288 hPtrs_h[i] = hFields[i].
cArray();
294 for (
int i = 0; i < nMonomer; i++) {
298 chiBottom_d = chiBottom_h;
307 _generateHFields<D><<<gridDims, blockDims>>>
313 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.
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.
FilmFieldGenExtBase()
Constructor.
FilmFieldGenExt()
Default constructor.
void setClassName(const char *className)
Set class name string.
void compute()
Compute the fields and store where the System can access.
~FilmFieldGenExt()
Destructor.
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.
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.
Fields, FFTs, and utilities for periodic boundary conditions (CUDA)
Periodic fields and crystallography.
SCFT and PS-FTS with real periodic fields (GPU)
PSCF package top-level namespace.