1#ifndef RPC_MASK_GEN_FILM_TPP
2#define RPC_MASK_GEN_FILM_TPP
11#include "FilmFieldGenMask.h"
12#include <rpc/solvers/Mixture.h>
13#include <rpc/scft/iterator/Iterator.h>
14#include <rpc/scft/ScftThermo.h>
15#include <rpc/field/Domain.h>
16#include <rpc/field/FieldIo.h>
17#include <prdc/cpu/RField.h>
18#include <prdc/crystal/UnitCell.h>
19#include <prdc/crystal/Basis.h>
20#include <prdc/crystal/paramIdConversions.h>
21#include <pscf/inter/Interaction.h>
22#include <pscf/math/IntVec.h>
72 if (normalVecParamId != paramId)
return 0.0;
78 double nvLength = domain.
unitCell().parameter(paramId);
81 double phiTot =
system().mask().phiTot();
86 for (
int ind = 0; ind < 3; ind++) {
104 for (x = 0; x < dim[0]; x++) {
106 for (y = 0; y < dim[1]; y++) {
108 for (z = 0; z < dim[2]; z++) {
110 maskVal = maskRGrid[counter];
117 deriv[counter] = maskVal * (maskVal - 1) * 8.0
118 * (std::abs(d - 0.5) - 0.5)
126 int nMonomer =
system().mixture().nMonomer();
127 int nx = mesh.
size();
132 for (
int i = 0; i < nx; i++) {
133 xi[i] =
system().w().rgrid(0)[i];
136 for (
int in = 0; in < nMonomer; in++) {
137 chi =
system().interaction().chi(0,in);
138 if (fabs(chi) > 1e-6) {
140 for (
int i = 0; i < nx; i++) {
142 xi[i] -= cfield[i] * chi;
147 if (
system().h().hasData()) {
149 for (
int i = 0; i < nx; i++) {
156 double intTerm = 0.0;
157 for (
int i = 0; i < nx; i++) {
158 intTerm += xi[i] * deriv[i];
160 intTerm /= (phiTot * nx);
163 if (!
system().scft().hasData()) {
164 sysPtr_->scft().compute();
166 double pSys =
system().scft().pressure();
168 (phiTot * nvLength * nvLength);
170 return pTerm - intTerm;
179 system().domain().lattice());
180 if (nvParamId == paramId) {
184 UTIL_THROW(
"fBulk must be set before calculating this stress.");
188 if (!
system().scft().hasData()) {
189 sysPtr_->scft().compute();
191 double fSys =
system().scft().fHelmholtz();
194 double L =
system().domain().unitCell().parameter(paramId);
199 return modifiedStress;
215 if (
system().iterator().isSymmetric()) {
221 system().domain().lattice());
222 double L =
system().domain().unitCell().parameter(paramId);
227 for (
int ind = 0; ind < 3; ind++) {
229 dim[ind] =
system().domain().mesh().dimensions()[ind];
243 for (x = 0; x < dim[0]; x++) {
245 for (y = 0; y < dim[1]; y++) {
247 for (z = 0; z < dim[2]; z++) {
257 rGrid[counter++] = 1-rhoW;
263 system().mask().setRGrid(rGrid,
system().iterator().isSymmetric());
277 if (
system().iterator().isFlexible()) {
280 system().domain().unitCell());
281 sysPtr_->iterator().setFlexibleParams(updated);
An IntVec<D, T> is a D-component vector of elements of integer type T.
Description of a regular grid of points in a periodic domain.
IntVec< D > dimensions() const
Get an IntVec<D> of the grid dimensions.
int size() const
Get total number of grid points.
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.
Spatial domain for a periodic structure with real fields, on a CPU.
Mesh< D > & mesh()
Get the Mesh by non-const reference.
UnitCell< D > & unitCell()
Get the UnitCell by non-const reference.
UnitCell< D >::LatticeSystem lattice() const
Get the lattice system (enumeration value).
void setFlexibleParams() const
Modifies iterator().flexibleParams_ to be compatible with the mask.
void setClassName(const char *className)
Set class name string.
double modifyStress(int paramId, double stress) const
Modify stress value in direction normal to the film.
RealVec< D > systemLatticeVector(int id) const
Get one of the lattice vectors for this system.
System< D > & system()
Get the parent System by non-const reference.
~FilmFieldGenMask()
Destructor.
FilmFieldGenMask()
Default constructor.
double stress(int paramId) const
Get contribution to the stress from this mask.
void compute()
Compute the field and store where the System can access.
Main class, representing a complete physical system.
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.
Fields and FFTs for periodic boundary conditions (CPU)
Periodic fields and crystallography.
Real periodic fields, SCFT and PS-FTS (CPU).
PSCF package top-level namespace.