1#ifndef RPG_MASK_GEN_FILM_TPP
2#define RPG_MASK_GEN_FILM_TPP
11#include "MaskGenFilm.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>
49 __global__
void _maskDerivative(cudaReal * deriv,
50 cudaReal
const * mask,
51 cudaReal
const interfaceThickness,
52 int const normalVecId,
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;
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;
72 if (normalVecId == 1) d = (cudaReal)ix / (cudaReal)
meshDims.x;
73 if (normalVecId == 0) d = (cudaReal)iy / (cudaReal)
meshDims.y;
75 d = (cudaReal)ix / (cudaReal)
meshDims.x;
80 cudaReal maskVal = mask[id];
83 deriv[id] = maskVal * (maskVal - 1) * 8.0 * (
abs(d - 0.5) - 0.5)
99 __global__
void _generateMask(cudaReal * mask, cudaReal
const nvLength,
100 cudaReal
const interfaceThickness,
101 cudaReal
const excludedThickness,
102 int const normalVecId, dim3
const meshDims)
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;
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;
121 if (normalVecId == 1) d = (cudaReal)ix / (cudaReal)
meshDims.x;
122 if (normalVecId == 0) d = (cudaReal)iy / (cudaReal)
meshDims.y;
124 d = (cudaReal)ix / (cudaReal)
meshDims.x;
135 0.5 * (excludedThickness - nvLength) +
136 abs(d - (nvLength / 2.0))
138 * 4.0 / interfaceThickness
179 system().domain().lattice());
180 if (normalVecParamId != paramId)
return 0.0;
186 double nvLength = system().unitCell().parameter(paramId);
189 double phiTot = system().mask().phiTot();
193 deriv.
allocate(system().mesh().dimensions());
201 _maskDerivative<D><<<gridDims, blockDims>>>
202 (deriv.
cArray(), system().mask().rgrid().cArray(),
207 xi.
allocate(system().mesh().dimensions());
209 if (system().hasExternalFields()) {
210 VecOp::subVV(xi, system().w().rgrid(0), system().h().rgrid(0));
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) {
225 intTerm /= phiTot * deriv.
capacity();
228 if (!sysPtr_->hasFreeEnergy()) {
229 sysPtr_->computeFreeEnergy();
231 double pSys = sysPtr_->pressure();
232 double pTerm = pSys * excludedThickness() /
233 (phiTot * nvLength * nvLength);
235 return pTerm - intTerm;
243 system().domain().lattice());
246 if (nvParamId != paramId)
return stress;
252 UTIL_THROW(
"fBulk must be set before calling modifyStress.");
256 if (!sysPtr_->hasFreeEnergy()) {
257 sysPtr_->computeFreeEnergy();
259 double fSys = sysPtr_->fHelmholtz();
262 double L = system().unitCell().parameter(paramId);
264 double modifiedStress = (stress * (L - excludedThickness()))
267 return modifiedStress;
276 UTIL_CHECK(system().unitCell().isInitialized());
279 system().mask().setFieldIo(system().fieldIo());
282 if (!system().mask().isAllocatedRGrid()) {
283 system().mask().allocateRGrid(system().mesh().dimensions());
285 if (system().iterator().isSymmetric()) {
287 if (!system().mask().isAllocatedBasis()) {
288 system().mask().allocateBasis(system().basis().nBasis());
301 UTIL_CHECK(excludedThickness() > interfaceThickness());
302 UTIL_CHECK(system().mask().isAllocatedRGrid());
303 if (system().iterator().isSymmetric()) {
304 UTIL_CHECK(system().mask().isAllocatedBasis());
309 system().domain().lattice());
310 double nvLength = system().unitCell().parameter(paramId);
314 rGrid.
allocate(system().mesh().dimensions());
322 _generateMask<D><<<gridDims, blockDims>>>
323 (rGrid.
cArray(), nvLength, interfaceThickness(),
327 system().mask().setRGrid(rGrid,system().iterator().isSymmetric());
330 normalVecCurrent_ = systemLatticeVector(normalVecId());
339 if (system().iterator().isFlexible()) {
341 updated = modifyFlexibleParams(system().iterator().flexibleParams(),
342 system().domain().unitCell());
343 sysPtr_->iterator().setFlexibleParams(updated);
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.
~MaskGenFilm()
Destructor.
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.
MaskGenFilm()
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.
A fixed capacity (static) contiguous array with a variable logical size.
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
#define UTIL_THROW(msg)
Macro for throwing an Exception, reporting function, file and line number.
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).
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).
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).
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.
Fields, FFTs, and utilities for periodic boundary conditions (CUDA)
Periodic fields and crystallography.
PSCF package top-level namespace.
Utility classes for scientific computation.