1#ifndef RPC_EXT_GEN_FILM_TPP
2#define RPC_EXT_GEN_FILM_TPP
11#include "FilmFieldGenExt.h"
12#include <rpc/field/FieldIo.h>
13#include <rpc/scft/iterator/Iterator.h>
14#include <prdc/cpu/RField.h>
15#include <prdc/crystal/paramIdConversions.h>
16#include <pscf/math/IntVec.h>
65 system().domain().lattice());
66 if (nvParamId != paramId)
return 0.0;
75 int nMonomer =
system().mixture().nMonomer();
76 int nx =
system().domain().mesh().size();
90 for (
int ind = 0; ind < 3; ind++) {
92 dim[ind] =
system().domain().mesh().dimensions()[ind];
99 for (x = 0; x < dim[0]; x++) {
101 for (y = 0; y < dim[1]; y++) {
103 for (z = 0; z < dim[2]; z++) {
105 maskVal = maskRGrid[counter];
112 maskDeriv[counter] = maskVal * (maskVal - 1) * 8.0
113 * (std::abs(d - 0.5) - 0.5)
120 for (
int i = 0; i < nMonomer; i++) {
123 for (x = 0; x < dim[0]; x++) {
125 for (y = 0; y < dim[1]; y++) {
127 for (z = 0; z < dim[2]; z++) {
130 hDeriv[counter] = -1.0 * maskDeriv[counter] *
chiBottom(i);
132 hDeriv[counter] = -1.0 * maskDeriv[counter] *
chiTop(i);
141 for (
int i = 0; i < nx; i++) {
142 term += c[i] * hDeriv[i];
145 term /= (
system().mask().phiTot() * nx);
171 if (
system().iterator().isSymmetric()) {
180 for (
int ind = 0; ind < 3; ind++) {
182 dim[ind] =
system().domain().mesh().dimension(ind);
194 for (
int i = 0; i < nm; i++) {
203 for (i = 0; i < nm; i++) {
204 for (x = 0; x < dim[0]; x++) {
206 for (y = 0; y < dim[1]; y++) {
208 for (z = 0; z < dim[2]; z++) {
213 rhoW = maskPtr[counter];
215 hRGrid[i][counter++] = (1.0-rhoW) *
chiBottom(i);
217 hRGrid[i][counter++] = (1.0-rhoW) *
chiTop(i);
226 system().h().setRGrid(hRGrid,
system().iterator().isSymmetric());
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.
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.
RealVec< D > systemLatticeVector(int id) const
Get one of the lattice vectors for this system.
~FilmFieldGenExt()
Destructor.
Main class, representing one complete 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.