PSCF v1.2
rpc/scft/iterator/ExtGenFilm.tpp
1#ifndef RPC_EXT_GEN_FILM_TPP
2#define RPC_EXT_GEN_FILM_TPP
3
4/*
5* PSCF - Polymer Self-Consistent Field Theory
6*
7* Copyright 2016 - 2022, The Regents of the University of Minnesota
8* Distributed under the terms of the GNU General Public License.
9*/
10
11#include "ExtGenFilm.h"
12#include "Iterator.h"
13#include <rpc/field/FieldIo.h>
14#include <prdc/cpu/RField.h>
15#include <prdc/crystal/paramIdConversions.h>
16#include <pscf/math/IntVec.h>
17#include <util/containers/FArray.h>
18#include <cmath>
19
20namespace Pscf {
21namespace Rpc
22{
23
24 using namespace Util;
25 using namespace Pscf::Prdc;
26 using namespace Pscf::Prdc::Cpu;
27
28 /*
29 * Default constructor
30 */
31 template <int D>
34 sysPtr_(0)
35 { setClassName("ExtGenFilm"); }
36
37 /*
38 * Constructor
39 */
40 template <int D>
43 sysPtr_(&sys)
44 { setClassName("ExtGenFilm"); }
45
46 /*
47 * Destructor
48 */
49 template <int D>
52
53 /*
54 * Get contribution to the stress from these external fields
55 */
56 template <int D>
57 double ExtGenFilm<D>::stressTerm(int paramId) const
58 {
59 // If walls are athermal then there is no external field, so no
60 // contribution to the stress.
61 if (isAthermal()) return 0.0;
62
63 // If paramId is not normalVecId, there is no stress contribution
64 UTIL_CHECK(normalVecId() >= 0); // normalVecId has been set
65 int nvParamId = convertFullParamIdToReduced<D>(normalVecId(),
66 system().domain().lattice());
67 if (nvParamId != paramId) return 0.0;
68
69 // If this point is reached, calculate the stress contribution
70 // from the external fields.
71 UTIL_CHECK(isGenerated());
72 UTIL_CHECK(interfaceThickness() > 0);
73 UTIL_CHECK(system().hasMask());
74 UTIL_CHECK(system().hasExternalFields());
75
76 // Setup
77 int nMonomer = system().mixture().nMonomer();
78 int nx = system().domain().mesh().size();
79 RField<D> const & maskRGrid = system().mask().rgrid();
80 RField<D> maskDeriv, hDeriv;
81 maskDeriv.allocate(system().domain().mesh().dimensions());
82 hDeriv.allocate(system().domain().mesh().dimensions());
83 FArray<int,3> coords;
84 int x, y, z;
85 int counter = 0;
86 double d, maskVal;
87 double term = 0.0;
88
89 // Create a 3 element vector 'dim' that contains the grid dimensions.
90 // If system is 2D (1D), then the z (y & z) dimensions are set to 1.
91 IntVec<3> dim;
92 for (int ind = 0; ind < 3; ind++) {
93 if (ind < D) {
94 dim[ind] = system().domain().mesh().dimensions()[ind];
95 } else {
96 dim[ind] = 1;
97 }
98 }
99
100 // Compute the derivative of the mask with respect to film thickness
101 for (x = 0; x < dim[0]; x++) {
102 coords[0] = x;
103 for (y = 0; y < dim[1]; y++) {
104 coords[1] = y;
105 for (z = 0; z < dim[2]; z++) {
106 coords[2] = z;
107 maskVal = maskRGrid[counter];
108
109 // Get the distance 'd' traveled along the lattice basis
110 // vector normal to the walls, in reduced coordinates
111 d = (double)coords[normalVecId()] /
112 (double)dim[normalVecId()];
113
114 maskDeriv[counter] = maskVal * (maskVal - 1) * 8.0
115 * (std::abs(d - 0.5) - 0.5)
116 / interfaceThickness();
117 counter++;
118 }
119 }
120 }
121
122 for (int i = 0; i < nMonomer; i++) {
123 // Compute the h field derivative with respect to film thickness
124 counter = 0;
125 for (x = 0; x < dim[0]; x++) {
126 coords[0] = x;
127 for (y = 0; y < dim[1]; y++) {
128 coords[1] = y;
129 for (z = 0; z < dim[2]; z++) {
130 coords[2] = z;
131 if (coords[normalVecId()] < (dim[normalVecId()] / 2)) {
132 hDeriv[counter] = -1.0 * maskDeriv[counter] * chiBottom(i);
133 } else {
134 hDeriv[counter] = -1.0 * maskDeriv[counter] * chiTop(i);
135 }
136 counter++;
137 }
138 }
139 }
140
141 // Get the integral term in the stress
142 RField<D> const & c = system().c().rgrid(i);
143 for (int i = 0; i < nx; i++) {
144 term += c[i] * hDeriv[i];
145 }
146 }
147 term /= (system().mask().phiTot() * nx);
148 return term;
149 }
150
151 /*
152 * Allocate container necessary to generate and store field
153 */
154 template <int D>
156 {
157 UTIL_CHECK(system().domain().unitCell().isInitialized());
158
159 // Make sure h field container has access to a fieldIo
160 system().h().setFieldIo(system().domain().fieldIo());
161
162 // Allocate the external field containers if needed
163 if (!isAthermal()) {
164 if (!system().h().isAllocatedRGrid()) {
165 system().h().allocateRGrid(system().domain().mesh().dimensions());
166 }
167 if (system().iterator().isSymmetric()) {
168 UTIL_CHECK(system().domain().basis().isInitialized());
169 if (!system().h().isAllocatedBasis()) {
170 system().h().allocateBasis(system().domain().basis().nBasis());
171 }
172 }
173 }
174 }
175
179 template <int D>
181 {
182 // If walls are athermal then there is no external field needed.
183 // If an external field already exists in the System, we need to
184 // overwrite it with a field of all zeros, otherwise do nothing
185 if ((isAthermal()) && (!isGenerated())) return;
186
187 // If this point is reached, external field must be generated
188 UTIL_CHECK(system().h().isAllocatedRGrid());
189 UTIL_CHECK(system().mask().isAllocatedRGrid());
190 if (system().iterator().isSymmetric()) {
191 UTIL_CHECK(system().h().isAllocatedBasis());
192 UTIL_CHECK(system().mask().isAllocatedBasis());
193 }
194
195 UTIL_CHECK(normalVecId() >= 0);
196
197 // Set chiBottomCurrent_, chiTopCurrent_, and parametersCurrent_
198 chiBottomCurrent_ = chiBottom();
199 chiTopCurrent_ = chiTop();
200 normalVecCurrent_ = systemLatticeVector(normalVecId());
201
202 int nm = systemNMonomer();
203
204 // Create a 3 element vector 'dim' that contains the grid dimensions.
205 // If system is 2D (1D), then the z (y & z) dimensions are set to 1.
206 IntVec<3> dim;
207 for (int ind = 0; ind < 3; ind++) {
208 if (ind < D) {
209 dim[ind] = system().domain().mesh().dimension(ind);
210 } else {
211 dim[ind] = 1;
212 }
213 }
214
215 // Get pointer to mask RField
216 RField<D> const & maskPtr = system().mask().rgrid();
217
218 // Generate an r-grid representation of the external fields
219 DArray< RField<D> > hRGrid;
220 hRGrid.allocate(nm);
221 for (int i = 0; i < nm; i++) {
222 hRGrid[i].allocate(system().domain().mesh().dimensions());
223 }
224
225 int i, x, y, z;
226 int counter = 0;
227 FArray<int,3> coords;
228 double rhoW;
229
230 for (i = 0; i < nm; i++) {
231 for (x = 0; x < dim[0]; x++) {
232 coords[0] = x;
233 for (y = 0; y < dim[1]; y++) {
234 coords[1] = y;
235 for (z = 0; z < dim[2]; z++) {
236 coords[2] = z;
237
238 // Calculate wall volume fraction (rho_w) at gridpoint
239 // (x,y,z)
240 rhoW = maskPtr[counter];
241 if (coords[normalVecId()] < (dim[normalVecId()]/2)) {
242 hRGrid[i][counter++] = (1.0-rhoW) * chiBottom(i);
243 } else {
244 hRGrid[i][counter++] = (1.0-rhoW) * chiTop(i);
245 }
246 }
247 }
248 }
249 counter = 0;
250 }
251
252 // Pass h into the System
253 system().h().setRGrid(hRGrid, system().iterator().isSymmetric());
254 }
255
256}
257}
258
259#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.
void allocate(IntVec< D > const &meshDimensions)
Allocate the underlying C array for an FFT grid.
Base class defining external fields for thin film systems.
void setClassName(const char *className)
Set class name string.
void allocate()
Allocate container necessary to generate and store fields.
void generate()
Generate the fields and store where the Iterator can access.
double stressTerm(int paramId) const
Get contribution to the stress from the external fields.
Main class for SCFT or PS-FTS simulation of one system.
Definition rpc/System.h:100
Dynamically allocatable contiguous array template.
void allocate(int capacity)
Allocate the underlying C array.
Definition DArray.h:199
A fixed size (static) contiguous array template.
Definition FArray.h:47
#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
PSCF package top-level namespace.
Definition param_pc.dox:1
Utility classes for scientific computation.