PSCF v1.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/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>
17#include <cmath>
18
19namespace Pscf {
20namespace Rpc
21{
22
23 using namespace Util;
24 using namespace Pscf::Prdc;
25 using namespace Pscf::Prdc::Cpu;
26
27 /*
28 * Default constructor
29 */
30 template <int D>
33 sysPtr_(nullptr)
34 { setClassName("FilmFieldGenExt"); }
35
36 /*
37 * Constructor
38 */
39 template <int D>
42 sysPtr_(&sys)
43 { setClassName("FilmFieldGenExt"); }
44
45 /*
46 * Destructor
47 */
48 template <int D>
51
52 /*
53 * Get contribution to the stress from these external fields
54 */
55 template <int D>
56 double FilmFieldGenExt<D>::stress(int paramId) const
57 {
58 // If walls are athermal then there is no external field, so no
59 // contribution to the stress.
60 if (isAthermal()) return 0.0;
61
62 // If paramId is not normalVecId, there is no stress contribution
63 UTIL_CHECK(normalVecId() >= 0); // normalVecId has been set
65 system().domain().lattice());
66 if (nvParamId != paramId) return 0.0;
67
68 // If this point is reached, calculate the stress contribution
69 // from the external fields.
71 UTIL_CHECK(system().mask().hasData());
72 UTIL_CHECK(system().h().hasData());
73
74 // Setup
75 int nMonomer = system().mixture().nMonomer();
76 int nx = system().domain().mesh().size();
77 RField<D> const & maskRGrid = system().mask().rgrid();
78 RField<D> maskDeriv, hDeriv;
79 maskDeriv.allocate(system().domain().mesh().dimensions());
80 hDeriv.allocate(system().domain().mesh().dimensions());
81 IntVec<3> coords;
82 int x, y, z;
83 int counter = 0;
84 double d, maskVal;
85 double term = 0.0;
86
87 // Create a 3 element vector 'dim' that contains the grid dimensions.
88 // If system is 2D (1D), then the z (y & z) dimensions are set to 1.
89 IntVec<3> dim;
90 for (int ind = 0; ind < 3; ind++) {
91 if (ind < D) {
92 dim[ind] = system().domain().mesh().dimensions()[ind];
93 } else {
94 dim[ind] = 1;
95 }
96 }
97
98 // Compute the derivative of the mask with respect to film thickness
99 for (x = 0; x < dim[0]; x++) {
100 coords[0] = x;
101 for (y = 0; y < dim[1]; y++) {
102 coords[1] = y;
103 for (z = 0; z < dim[2]; z++) {
104 coords[2] = z;
105 maskVal = maskRGrid[counter];
106
107 // Get the distance 'd' traveled along the lattice basis
108 // vector normal to the walls, in reduced coordinates
109 d = (double)coords[normalVecId()] /
110 (double)dim[normalVecId()];
111
112 maskDeriv[counter] = maskVal * (maskVal - 1) * 8.0
113 * (std::abs(d - 0.5) - 0.5)
115 counter++;
116 }
117 }
118 }
119
120 for (int i = 0; i < nMonomer; i++) {
121 // Compute the h field derivative with respect to film thickness
122 counter = 0;
123 for (x = 0; x < dim[0]; x++) {
124 coords[0] = x;
125 for (y = 0; y < dim[1]; y++) {
126 coords[1] = y;
127 for (z = 0; z < dim[2]; z++) {
128 coords[2] = z;
129 if (coords[normalVecId()] < (dim[normalVecId()] / 2)) {
130 hDeriv[counter] = -1.0 * maskDeriv[counter] * chiBottom(i);
131 } else {
132 hDeriv[counter] = -1.0 * maskDeriv[counter] * chiTop(i);
133 }
134 counter++;
135 }
136 }
137 }
138
139 // Get the integral term in the stress
140 RField<D> const & c = system().c().rgrid(i);
141 for (int i = 0; i < nx; i++) {
142 term += c[i] * hDeriv[i];
143 }
144 }
145 term /= (system().mask().phiTot() * nx);
146 return term;
147 }
148
152 template <int D>
154 {
155 // Set chiBottomCurrent_, chiTopCurrent_, and parametersCurrent_
159
160 // If walls are athermal then there is no external field needed.
161 // If external fields already exist in System, clear them, then return.
162 if (isAthermal()) {
163 system().h().clear();
164 return;
165 }
166
167 // If this point is reached, external field must be computed
168 UTIL_CHECK(normalVecId() >= 0);
169 UTIL_CHECK(system().domain().unitCell().isInitialized());
170 UTIL_CHECK(system().mask().hasData());
171 if (system().iterator().isSymmetric()) {
172 UTIL_CHECK(system().domain().basis().isInitialized());
173 }
174
175 int nm = systemNMonomer();
176
177 // Create a 3 element vector 'dim' that contains the grid dimensions.
178 // If system is 2D (1D), then the z (y & z) dimensions are set to 1.
179 IntVec<3> dim;
180 for (int ind = 0; ind < 3; ind++) {
181 if (ind < D) {
182 dim[ind] = system().domain().mesh().dimension(ind);
183 } else {
184 dim[ind] = 1;
185 }
186 }
187
188 // Get pointer to mask RField
189 RField<D> const & maskPtr = system().mask().rgrid();
190
191 // Generate an r-grid representation of the external fields
192 DArray< RField<D> > hRGrid;
193 hRGrid.allocate(nm);
194 for (int i = 0; i < nm; i++) {
195 hRGrid[i].allocate(system().domain().mesh().dimensions());
196 }
197
198 int i, x, y, z;
199 int counter = 0;
200 IntVec<3> coords;
201 double rhoW;
202
203 for (i = 0; i < nm; i++) {
204 for (x = 0; x < dim[0]; x++) {
205 coords[0] = x;
206 for (y = 0; y < dim[1]; y++) {
207 coords[1] = y;
208 for (z = 0; z < dim[2]; z++) {
209 coords[2] = z;
210
211 // Calculate wall volume fraction (rho_w) at gridpoint
212 // (x,y,z)
213 rhoW = maskPtr[counter];
214 if (coords[normalVecId()] < (dim[normalVecId()]/2)) {
215 hRGrid[i][counter++] = (1.0-rhoW) * chiBottom(i);
216 } else {
217 hRGrid[i][counter++] = (1.0-rhoW) * chiTop(i);
218 }
219 }
220 }
221 }
222 counter = 0;
223 }
224
225 // Pass h into the System
226 system().h().setRGrid(hRGrid, system().iterator().isSymmetric());
227 }
228
229}
230}
231
232#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.
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.
RealVec< D > systemLatticeVector(int id) const
Get one of the lattice vectors for this system.
Main class, representing one complete 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.
Definition param_pc.dox:1