1#ifndef RPC_EXT_GEN_FILM_TPP
2#define RPC_EXT_GEN_FILM_TPP
11#include "FilmFieldGenExt.h"
12#include <rpc/solvers/Mixture.h>
13#include <rpc/field/Domain.h>
14#include <rpc/field/FieldIo.h>
15#include <rpc/scft/iterator/Iterator.h>
16#include <prdc/cpu/RField.h>
17#include <prdc/crystal/Basis.h>
18#include <prdc/crystal/paramIdConversions.h>
19#include <pscf/math/IntVec.h>
68 system().domain().lattice());
69 if (nvParamId != paramId)
return 0.0;
78 int nMonomer =
system().mixture().nMonomer();
79 int nx =
system().domain().mesh().size();
93 for (
int ind = 0; ind < 3; ind++) {
95 dim[ind] =
system().domain().mesh().dimensions()[ind];
102 for (x = 0; x < dim[0]; x++) {
104 for (y = 0; y < dim[1]; y++) {
106 for (z = 0; z < dim[2]; z++) {
108 maskVal = maskRGrid[counter];
115 maskDeriv[counter] = maskVal * (maskVal - 1) * 8.0
116 * (std::abs(d - 0.5) - 0.5)
123 for (
int i = 0; i < nMonomer; i++) {
126 for (x = 0; x < dim[0]; x++) {
128 for (y = 0; y < dim[1]; y++) {
130 for (z = 0; z < dim[2]; z++) {
133 hDeriv[counter] = -1.0 * maskDeriv[counter] *
chiBottom(i);
135 hDeriv[counter] = -1.0 * maskDeriv[counter] *
chiTop(i);
144 for (
int i = 0; i < nx; i++) {
145 term += c[i] * hDeriv[i];
148 term /= (
system().mask().phiTot() * nx);
174 if (
system().iterator().isSymmetric()) {
183 for (
int ind = 0; ind < 3; ind++) {
185 dim[ind] =
system().domain().mesh().dimension(ind);
197 for (
int i = 0; i < nm; i++) {
206 for (i = 0; i < nm; i++) {
207 for (x = 0; x < dim[0]; x++) {
209 for (y = 0; y < dim[1]; y++) {
211 for (z = 0; z < dim[2]; z++) {
216 rhoW = maskPtr[counter];
218 hRGrid[i][counter++] = (1.0-rhoW) *
chiBottom(i);
220 hRGrid[i][counter++] = (1.0-rhoW) *
chiTop(i);
229 system().h().setRGrid(hRGrid,
system().iterator().isSymmetric());
237 {
return system().domain().groupName(); }
244 {
return system().domain().unitCell().rBasis(
id); }
251 {
return system().mixture().nMonomer(); }
An IntVec<D, T> is a D-component vector of elements of integer type T.
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.
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.
A RealVec<D, T> is D-component vector with elements of floating type T.
double stress(int paramId) const
Get contribution to the stress from the external fields.
System< D > & system()
Get the parent System by non-const reference.
void compute()
Compute the fields and store where the System can access.
void setClassName(const char *className)
Set class name string.
int systemNMonomer() const
Get the number of monomer species for this system.
FilmFieldGenExt()
Default constructor.
std::string systemSpaceGroup() const
Get the space group name for this system.
RealVec< D > systemLatticeVector(int id) const
Get one of the lattice vectors for this system.
~FilmFieldGenExt()
Destructor.
Main class, representing a complete physical system.
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.
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.