1#ifndef RPG_MASK_GEN_FILM_TPP
2#define RPG_MASK_GEN_FILM_TPP
11#include "FilmFieldGenMask.h"
12#include <rpg/field/FieldIo.h>
13#include <rpg/scft/iterator/Iterator.h>
14#include <rpg/scft/ScftThermo.h>
15#include <prdc/cpu/RField.h>
16#include <prdc/crystal/UnitCell.h>
17#include <prdc/crystal/paramIdConversions.h>
18#include <prdc/cuda/resources.h>
19#include <pscf/inter/Interaction.h>
51 __global__
void _maskDerivative(cudaReal * deriv,
52 cudaReal
const * mask,
53 cudaReal
const interfaceThickness,
54 int const normalVecId,
58 int ix = blockIdx.x * blockDim.x + threadIdx.x;
59 int iy = blockIdx.y * blockDim.y + threadIdx.y;
60 int iz = blockIdx.z * blockDim.z + threadIdx.z;
70 if (normalVecId == 2) d = (cudaReal)ix / (cudaReal)
meshDims.x;
71 if (normalVecId == 1) d = (cudaReal)iy / (cudaReal)
meshDims.y;
72 if (normalVecId == 0) d = (cudaReal)iz / (cudaReal)
meshDims.z;
74 if (normalVecId == 1) d = (cudaReal)ix / (cudaReal)
meshDims.x;
75 if (normalVecId == 0) d = (cudaReal)iy / (cudaReal)
meshDims.y;
77 d = (cudaReal)ix / (cudaReal)
meshDims.x;
82 cudaReal maskVal = mask[id];
85 deriv[id] = maskVal * (maskVal - 1) * 8.0 * (
abs(d - 0.5) - 0.5)
101 __global__
void _generateMask(cudaReal * mask, cudaReal
const nvLength,
102 cudaReal
const interfaceThickness,
103 cudaReal
const excludedThickness,
104 int const normalVecId, dim3
const meshDims)
107 int ix = blockIdx.x * blockDim.x + threadIdx.x;
108 int iy = blockIdx.y * blockDim.y + threadIdx.y;
109 int iz = blockIdx.z * blockDim.z + threadIdx.z;
119 if (normalVecId == 2) d = (cudaReal)ix / (cudaReal)
meshDims.x;
120 if (normalVecId == 1) d = (cudaReal)iy / (cudaReal)
meshDims.y;
121 if (normalVecId == 0) d = (cudaReal)iz / (cudaReal)
meshDims.z;
123 if (normalVecId == 1) d = (cudaReal)ix / (cudaReal)
meshDims.x;
124 if (normalVecId == 0) d = (cudaReal)iy / (cudaReal)
meshDims.y;
126 d = (cudaReal)ix / (cudaReal)
meshDims.x;
137 0.5 * (excludedThickness - nvLength) +
138 abs(d - (nvLength / 2.0))
140 * 4.0 / interfaceThickness
181 system().domain().lattice());
182 if (normalVecParamId != paramId)
return 0.0;
188 double nvLength =
system().domain().unitCell().parameter(paramId);
191 double phiTot =
system().mask().phiTot();
203 _maskDerivative<D><<<gridDims, blockDims>>>
211 if (
system().h().hasData()) {
217 int nMonomer =
system().mixture().nMonomer();
218 for (
int i = 0; i < nMonomer; i++) {
219 double chi =
system().interaction().chi(0, i);
220 if (fabs(chi) > 1e-6) {
226 double intTerm = Reduce::innerProduct(deriv, xi);
227 intTerm /= phiTot * deriv.
capacity();
230 if (!sysPtr_->scft().hasData()) {
231 sysPtr_->scft().compute();
233 double pSys = sysPtr_->scft().pressure();
235 (phiTot * nvLength * nvLength);
237 return pTerm - intTerm;
245 system().domain().lattice());
248 if (nvParamId != paramId)
return stress;
254 UTIL_THROW(
"fBulk must be set before calling modifyStress.");
258 if (!sysPtr_->scft().hasData()) {
259 sysPtr_->scft().compute();
261 double fSys = sysPtr_->scft().fHelmholtz();
264 double L =
system().domain().unitCell().parameter(paramId);
269 return modifiedStress;
281 if (
system().iterator().isSymmetric()) {
288 system().domain().lattice());
289 double nvLength =
system().domain().unitCell().parameter(paramId);
301 _generateMask<D><<<gridDims, blockDims>>>
306 system().mask().setRGrid(rGrid,
system().iterator().isSymmetric());
318 if (
system().iterator().isFlexible()) {
321 system().domain().unitCell());
322 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.
RealVec< D > normalVecCurrent_
The lattice vector normal to the film used to generate these fields.
int normalVecId() const
Get value of normalVecId.
double fBulk_
Reference free energy used to calculate stress normal to the film.
bool hasFBulk() const
Check whether a value of fBulk was provided.
FSArray< bool, 6 > modifyFlexibleParams(FSArray< bool, 6 > current, UnitCell< D > const &cell) const
Modifies a flexibleParams array to be compatible with this mask.
double interfaceThickness() const
Get value of interfaceThickness.
double excludedThickness() const
Get value of excludedThickness.
FilmFieldGenMaskBase()
Constructor.
void compute()
Compute the field and store where the System can access.
double modifyStress(int paramId, double stress) const
Modify stress value in direction normal to the film.
void setClassName(const char *className)
Set class name string.
System< D > & system()
Get the System associated with this object by reference.
void setFlexibleParams() const
Modifies iterator().flexibleParams_ to be compatible with the mask.
double stress(int paramId) const
Get contribution to the stress from this mask.
~FilmFieldGenMask()
Destructor.
FilmFieldGenMask()
Default constructor.
RealVec< D > systemLatticeVector(int id) const
Get one of the lattice vectors for this system.
Main class, representing one complete 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.
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.
SCFT and PS-FTS with real periodic fields (GPU)
PSCF package top-level namespace.