PSCF v1.3.3
rpc/environment/FilmFieldGenExt.tpp
1#ifndef RPC_EXT_GEN_FILM_TPP
2#define RPC_EXT_GEN_FILM_TPP
3
4/*
5* PSCF - Polymer Self-Consistent Field
6*
7* Copyright 2015 - 2025, The Regents of the University of Minnesota
8* Distributed under the terms of the GNU General Public License.
9*/
10
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>
20#include <cmath>
21
22namespace Pscf {
23namespace Rpc
24{
25
26 using namespace Util;
27 using namespace Pscf::Prdc;
28 using namespace Pscf::Prdc::Cpu;
29
30 /*
31 * Default constructor
32 */
33 template <int D>
36 sysPtr_(nullptr)
37 { setClassName("FilmFieldGenExt"); }
38
39 /*
40 * Constructor
41 */
42 template <int D>
45 sysPtr_(&sys)
46 { setClassName("FilmFieldGenExt"); }
47
48 /*
49 * Destructor
50 */
51 template <int D>
54
55 /*
56 * Get contribution to the stress from these external fields
57 */
58 template <int D>
59 double FilmFieldGenExt<D>::stress(int paramId) const
60 {
61 // If walls are athermal then there is no external field, so no
62 // contribution to the stress.
63 if (isAthermal()) return 0.0;
64
65 // If paramId is not normalVecId, there is no stress contribution
66 UTIL_CHECK(normalVecId() >= 0); // normalVecId has been set
68 system().domain().lattice());
69 if (nvParamId != paramId) return 0.0;
70
71 // If this point is reached, calculate the stress contribution
72 // from the external fields.
74 UTIL_CHECK(system().mask().hasData());
75 UTIL_CHECK(system().h().hasData());
76
77 // Setup
78 int nMonomer = system().mixture().nMonomer();
79 int nx = system().domain().mesh().size();
80 RField<D> const & maskRGrid = system().mask().rgrid();
81 RField<D> maskDeriv, hDeriv;
82 maskDeriv.allocate(system().domain().mesh().dimensions());
83 hDeriv.allocate(system().domain().mesh().dimensions());
84 IntVec<3> coords;
85 int x, y, z;
86 int counter = 0;
87 double d, maskVal;
88 double term = 0.0;
89
90 // Create a 3 element vector 'dim' that contains the grid dimensions.
91 // If system is 2D (1D), then the z (y & z) dimensions are set to 1.
92 IntVec<3> dim;
93 for (int ind = 0; ind < 3; ind++) {
94 if (ind < D) {
95 dim[ind] = system().domain().mesh().dimensions()[ind];
96 } else {
97 dim[ind] = 1;
98 }
99 }
100
101 // Compute the derivative of the mask with respect to film thickness
102 for (x = 0; x < dim[0]; x++) {
103 coords[0] = x;
104 for (y = 0; y < dim[1]; y++) {
105 coords[1] = y;
106 for (z = 0; z < dim[2]; z++) {
107 coords[2] = z;
108 maskVal = maskRGrid[counter];
109
110 // Get the distance 'd' traveled along the lattice basis
111 // vector normal to the walls, in reduced coordinates
112 d = (double)coords[normalVecId()] /
113 (double)dim[normalVecId()];
114
115 maskDeriv[counter] = maskVal * (maskVal - 1) * 8.0
116 * (std::abs(d - 0.5) - 0.5)
118 counter++;
119 }
120 }
121 }
122
123 for (int i = 0; i < nMonomer; i++) {
124 // Compute the h field derivative with respect to film thickness
125 counter = 0;
126 for (x = 0; x < dim[0]; x++) {
127 coords[0] = x;
128 for (y = 0; y < dim[1]; y++) {
129 coords[1] = y;
130 for (z = 0; z < dim[2]; z++) {
131 coords[2] = z;
132 if (coords[normalVecId()] < (dim[normalVecId()] / 2)) {
133 hDeriv[counter] = -1.0 * maskDeriv[counter] * chiBottom(i);
134 } else {
135 hDeriv[counter] = -1.0 * maskDeriv[counter] * chiTop(i);
136 }
137 counter++;
138 }
139 }
140 }
141
142 // Get the integral term in the stress
143 RField<D> const & c = system().c().rgrid(i);
144 for (int i = 0; i < nx; i++) {
145 term += c[i] * hDeriv[i];
146 }
147 }
148 term /= (system().mask().phiTot() * nx);
149 return term;
150 }
151
155 template <int D>
157 {
158 // Set chiBottomCurrent_, chiTopCurrent_, and parametersCurrent_
162
163 // If walls are athermal then there is no external field needed.
164 // If external fields already exist in System, clear them, then return.
165 if (isAthermal()) {
166 system().h().clear();
167 return;
168 }
169
170 // If this point is reached, external field must be computed
171 UTIL_CHECK(normalVecId() >= 0);
172 UTIL_CHECK(system().domain().unitCell().isInitialized());
173 UTIL_CHECK(system().mask().hasData());
174 if (system().iterator().isSymmetric()) {
175 UTIL_CHECK(system().domain().basis().isInitialized());
176 }
177
178 int nm = systemNMonomer();
179
180 // Create a 3 element vector 'dim' that contains the grid dimensions.
181 // If system is 2D (1D), then the z (y & z) dimensions are set to 1.
182 IntVec<3> dim;
183 for (int ind = 0; ind < 3; ind++) {
184 if (ind < D) {
185 dim[ind] = system().domain().mesh().dimension(ind);
186 } else {
187 dim[ind] = 1;
188 }
189 }
190
191 // Get pointer to mask RField
192 RField<D> const & maskPtr = system().mask().rgrid();
193
194 // Generate an r-grid representation of the external fields
195 DArray< RField<D> > hRGrid;
196 hRGrid.allocate(nm);
197 for (int i = 0; i < nm; i++) {
198 hRGrid[i].allocate(system().domain().mesh().dimensions());
199 }
200
201 int i, x, y, z;
202 int counter = 0;
203 IntVec<3> coords;
204 double rhoW;
205
206 for (i = 0; i < nm; i++) {
207 for (x = 0; x < dim[0]; x++) {
208 coords[0] = x;
209 for (y = 0; y < dim[1]; y++) {
210 coords[1] = y;
211 for (z = 0; z < dim[2]; z++) {
212 coords[2] = z;
213
214 // Calculate wall volume fraction (rho_w) at gridpoint
215 // (x,y,z)
216 rhoW = maskPtr[counter];
217 if (coords[normalVecId()] < (dim[normalVecId()]/2)) {
218 hRGrid[i][counter++] = (1.0-rhoW) * chiBottom(i);
219 } else {
220 hRGrid[i][counter++] = (1.0-rhoW) * chiTop(i);
221 }
222 }
223 }
224 }
225 counter = 0;
226 }
227
228 // Pass h into the System
229 system().h().setRGrid(hRGrid, system().iterator().isSymmetric());
230 }
231
232 /*
233 * Get space group name for this system.
234 */
235 template <int D>
237 { return system().domain().groupName(); }
238
239 /*
240 * Get one of the lattice vectors for this system.
241 */
242 template <int D>
244 { return system().domain().unitCell().rBasis(id); }
245
246 /*
247 * Get the number of monomer species for this system.
248 */
249 template <int D>
251 { return system().mixture().nMonomer(); }
252
253}
254}
255
256#endif
An IntVec<D, T> is a D-component vector of elements of integer type T.
Definition IntVec.h:27
Field of real double precision values on an FFT mesh.
Definition cpu/RField.h:29
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.
A RealVec<D, T> is D-component vector with elements of floating type T.
Definition RealVec.h:28
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.
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.
Main class, representing a complete physical system.
Dynamically allocatable contiguous array template.
Definition DArray.h:32
void allocate(int capacity)
Allocate the underlying C array.
Definition DArray.h:199
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
Definition global.h:68
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)
Definition CField.cpp:12
Periodic fields and crystallography.
Definition CField.cpp:11
Real periodic fields, SCFT and PS-FTS (CPU).
Definition param_pc.dox:2
PSCF package top-level namespace.