PSCF v1.1
pspc/solvers/Block.tpp
1#ifndef PSPC_BLOCK_TPP
2#define PSPC_BLOCK_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 "Block.h"
12#include <pscf/mesh/Mesh.h>
13#include <pscf/mesh/MeshIterator.h>
14#include <pscf/crystal/UnitCell.h>
15#include <pscf/crystal/shiftToMinimum.h>
16#include <pscf/math/IntVec.h>
17#include <util/containers/DMatrix.h>
18#include <util/containers/DArray.h>
19#include <util/containers/FArray.h>
20#include <util/containers/FSArray.h>
21
22namespace Pscf {
23namespace Pspc {
24
25 using namespace Util;
26
27 /*
28 * Constructor.
29 */
30 template <int D>
32 : meshPtr_(0),
33 kMeshDimensions_(0),
34 ds_(0.0),
35 dsTarget_(0.0),
36 ns_(0),
37 isAllocated_(false),
38 hasExpKsq_(false)
39 {
40 propagator(0).setBlock(*this);
41 propagator(1).setBlock(*this);
42 }
43
44 /*
45 * Destructor.
46 */
47 template <int D>
49 {}
50
51 template <int D>
52 void Block<D>::setDiscretization(double ds, const Mesh<D>& mesh)
53 {
54 // Preconditions
55 UTIL_CHECK(mesh.size() > 1);
56 UTIL_CHECK(ds > 0.0);
57 UTIL_CHECK(!isAllocated_);
58
59 // Set contour length discretization
60 dsTarget_ = ds;
61 int tempNs;
62 tempNs = floor( length()/(2.0 *ds) + 0.5 );
63 if (tempNs == 0) {
64 tempNs = 1;
65 }
66 ns_ = 2*tempNs + 1;
67 ds_ = length()/double(ns_-1);
68
69 // Set association to mesh
70 meshPtr_ = &mesh;
71
72 fft_.setup(mesh.dimensions());
73
74 // Compute Fourier space kMeshDimensions_
75 for (int i = 0; i < D; ++i) {
76 if (i < D - 1) {
77 kMeshDimensions_[i] = mesh.dimensions()[i];
78 } else {
79 kMeshDimensions_[i] = mesh.dimensions()[i]/2 + 1;
80 }
81 }
82
83 // Set number kSize of points in k-space mesh
84 int kSize = 1;
85 for (int i = 0; i < D; ++i) {
86 kSize *= kMeshDimensions_[i];
87 }
88
89 // Allocate work arrays for MDE solution
90 expKsq_.allocate(kMeshDimensions_);
91 expKsq2_.allocate(kMeshDimensions_);
92 expW_.allocate(mesh.dimensions());
93 expW2_.allocate(mesh.dimensions());
94 qr_.allocate(mesh.dimensions());
95 qk_.allocate(mesh.dimensions());
96 qr2_.allocate(mesh.dimensions());
97 qk2_.allocate(mesh.dimensions());
98
99 // Allocate work array for stress calculation
100 dGsq_.allocate(kSize, 6);
101
102 // Allocate block concentration field
103 cField().allocate(mesh.dimensions());
104
105 // Allocate memory for solutions to MDE (requires ns_)
106 propagator(0).allocate(ns_, mesh);
107 propagator(1).allocate(ns_, mesh);
108
109 isAllocated_ = true;
110 hasExpKsq_ = false;
111 }
112
113 /*
114 * Set or reset the the block length.
115 */
116 template <int D>
117 void Block<D>::setLength(double newLength)
118 {
120
121 if (isAllocated_) { // if setDiscretization has already been called
122 // Reset contour length discretization
123 UTIL_CHECK(dsTarget_ > 0);
124 int oldNs = ns_;
125 int tempNs;
126 tempNs = floor( length()/(2.0 *dsTarget_) + 0.5 );
127 if (tempNs == 0) {
128 tempNs = 1;
129 }
130 ns_ = 2*tempNs + 1;
131 ds_ = length()/double(ns_-1);
132
133 if (oldNs != ns_) {
134 // If propagators are already allocated and ns_ has changed,
135 // reallocate memory for solutions to MDE
136 propagator(0).reallocate(ns_);
137 propagator(1).reallocate(ns_);
138 }
139 }
140
141 hasExpKsq_ = false;
142 }
143
144 /*
145 * Set or reset the the block length.
146 */
147 template <int D>
148 void Block<D>::setKuhn(double kuhn)
149 {
150 BlockTmpl< Propagator<D> >::setKuhn(kuhn);
151 hasExpKsq_ = false;
152 }
153
154 /*
155 * Setup data that depend on the unit cell parameters.
156 */
157 template <int D>
159 {
160 unitCellPtr_ = &unitCell;
161 hasExpKsq_ = false;
162 }
163
164 template <int D>
166 {
167 UTIL_CHECK(isAllocated_);
168 UTIL_CHECK(unitCellPtr_);
169 UTIL_CHECK(unitCellPtr_->isInitialized());
170
171 MeshIterator<D> iter;
172 iter.setDimensions(kMeshDimensions_);
173 IntVec<D> G, Gmin;
174 double Gsq;
175 double factor = -1.0*kuhn()*kuhn()*ds_/6.0;
176 int i;
177 for (iter.begin(); !iter.atEnd(); ++iter) {
178 i = iter.rank();
179 G = iter.position();
180 Gmin = shiftToMinimum(G, mesh().dimensions(), unitCell());
181 Gsq = unitCell().ksq(Gmin);
182 expKsq_[i] = exp(Gsq*factor);
183 expKsq2_[i] = exp(Gsq*factor*0.5);
184 }
185
186 hasExpKsq_ = true;
187 }
188
189 /*
190 * Setup the contour length step algorithm.
191 */
192 template <int D>
193 void
195 {
196 // Preconditions
197 int nx = mesh().size();
198 UTIL_CHECK(nx > 0);
199 UTIL_CHECK(isAllocated_);
200
201 // Compute expW arrays
202 for (int i = 0; i < nx; ++i) {
203
204 // First, check that w[i]*ds_ is not unreasonably large:
205 // (if this condition is not met, solution will have large
206 // error, and user should consider using a smaller ds_)
207 //UTIL_CHECK(std::abs(w[i]*ds_) < 1.0);
208
209 // Calculate values
210 expW_[i] = exp(-0.5*w[i]*ds_);
211 expW2_[i] = exp(-0.5*0.5*w[i]*ds_);
212 }
213
214 // Compute expKsq arrays if necessary
215 if (!hasExpKsq_) {
216 computeExpKsq();
217 }
218
219 }
220
221 /*
222 * Integrate to calculate monomer concentration for this block
223 */
224 template <int D>
225 void Block<D>::computeConcentration(double prefactor)
226 {
227 // Preconditions
228 UTIL_CHECK(isAllocated_);
229 int nx = mesh().size();
230 UTIL_CHECK(nx > 0);
231 UTIL_CHECK(ns_ > 0);
232 UTIL_CHECK(ds_ > 0);
233 UTIL_CHECK(propagator(0).isAllocated());
234 UTIL_CHECK(propagator(1).isAllocated());
235 UTIL_CHECK(cField().capacity() == nx);
236
237 // Initialize cField to zero at all points
238 int i;
239 for (i = 0; i < nx; ++i) {
240 cField()[i] = 0.0;
241 }
242
243 Propagator<D> const & p0 = propagator(0);
244 Propagator<D> const & p1 = propagator(1);
245
246 // Evaluate unnormalized integral
247
248 // Endpoint contributions
249 for (i = 0; i < nx; ++i) {
250 cField()[i] += p0.q(0)[i]*p1.q(ns_ - 1)[i];
251 cField()[i] += p0.q(ns_ -1)[i]*p1.q(0)[i];
252 }
253
254 // Odd indices
255 int j;
256 for (j = 1; j < (ns_ -1); j += 2) {
257 for (i = 0; i < nx; ++i) {
258 cField()[i] += p0.q(j)[i] * p1.q(ns_ - 1 - j)[i] * 4.0;
259 }
260 }
261
262 // Even indices
263 for (j = 2; j < (ns_ -2); j += 2) {
264 for (i = 0; i < nx; ++i) {
265 cField()[i] += p0.q(j)[i] * p1.q(ns_ - 1 - j)[i] * 2.0;
266 }
267 }
268
269 // Normalize the integral
270 prefactor *= ds_ / 3.0;
271 for (i = 0; i < nx; ++i) {
272 cField()[i] *= prefactor;
273 }
274
275 }
276
277 /*
278 * Integrate to Stress exerted by the chain for this block
279 */
280 template <int D>
281 void Block<D>::computeStress(double prefactor)
282 {
283 // Preconditions
284 int nx = mesh().size();
285 UTIL_CHECK(nx > 0);
286 UTIL_CHECK(ns_ > 0);
287 UTIL_CHECK(ds_ > 0);
288 UTIL_CHECK(propagator(0).isAllocated());
289 UTIL_CHECK(propagator(1).isAllocated());
290
291 stress_.clear();
292
293 double dels, normal, increment;
294 int nParam, c, m;
295
296 normal = 3.0*6.0;
297
298 int kSize_ = 1;
299 for (int i = 0; i < D; ++i) {
300 kSize_ *= kMeshDimensions_[i];
301 }
302 nParam = unitCell().nParameter();
303 c = kSize_;
304
306
307 // Initialize work array and stress_ to zero at all points
308 int i;
309 for (i = 0; i < nParam; ++i) {
310 dQ.append(0.0);
311 stress_.append(0.0);
312 }
313
314 computedGsq();
315
316 Propagator<D> const & p0 = propagator(0);
317 Propagator<D> const & p1 = propagator(1);
318
319 // Evaluate unnormalized integral
320 for (int j = 0; j < ns_ ; ++j) {
321
322 qr_ = p0.q(j);
323 fft_.forwardTransform(qr_, qk_);
324
325 qr2_ = p1.q(ns_ - 1 - j);
326 fft_.forwardTransform(qr2_, qk2_);
327
328 dels = ds_;
329
330 if (j != 0 && j != ns_ - 1) {
331 if (j % 2 == 0) {
332 dels = dels*2.0;
333 } else {
334 dels = dels*4.0;
335 }
336 }
337
338 for (int n = 0; n < nParam ; ++n) {
339 increment = 0.0;
340
341 for (m = 0; m < c ; ++m) {
342 double prod = 0;
343 prod = (qk2_[m][0] * qk_[m][0]) + (qk2_[m][1] * qk_[m][1]);
344 prod *= dGsq_(m,n);
345 increment += prod;
346 }
347 increment = (increment * kuhn() * kuhn() * dels)/normal;
348 dQ[n] = dQ[n] - increment;
349 }
350 }
351
352 // Normalize
353 for (i = 0; i < nParam; ++i) {
354 stress_[i] = stress_[i] - (dQ[i] * prefactor);
355 }
356
357 }
358
359 /*
360 * Compute dGsq_ array (derivatives of Gsq for all wavevectors)
361 */
362 template <int D>
364 {
365 IntVec<D> temp;
366 IntVec<D> vec;
367 IntVec<D> Partner;
368 MeshIterator<D> iter;
369 iter.setDimensions(kMeshDimensions_);
370
371 for (int n = 0; n < unitCell().nParameter() ; ++n) {
372 for (iter.begin(); !iter.atEnd(); ++iter) {
373 temp = iter.position();
374 vec = shiftToMinimum(temp, mesh().dimensions(), unitCell());
375 dGsq_(iter.rank(), n) = unitCell().dksq(vec, n);
376 for (int p = 0; p < D; ++p) {
377 if (temp [p] != 0) {
378 Partner[p] = mesh().dimensions()[p] - temp[p];
379 } else {
380 Partner[p] = 0;
381 }
382 }
383 if (Partner[D-1] > kMeshDimensions_[D-1]) {
384 dGsq_(iter.rank(), n) *= 2;
385 }
386 }
387 }
388 }
389
390 /*
391 * Propagate solution by one step.
392 */
393 template <int D>
394 void Block<D>::step(RField<D> const & q, RField<D>& qNew)
395 {
396 // Internal prereconditions
397 UTIL_CHECK(isAllocated_);
398 int nx = mesh().size();
399 int nk = qk_.capacity();
400 UTIL_CHECK(nx > 0);
401 UTIL_CHECK(nk > 0);
402 UTIL_CHECK(qr_.capacity() == nx);
403 UTIL_CHECK(expW_.capacity() == nx);
404 UTIL_CHECK(expKsq_.capacity() == nk);
405 UTIL_CHECK(hasExpKsq_);
406
407 // Preconditions on parameters
409 UTIL_CHECK(q.capacity() == nx);
410 UTIL_CHECK(qNew.isAllocated());
411 UTIL_CHECK(qNew.capacity() == nx);
412
413 // Apply pseudo-spectral algorithm
414
415 // Full step for ds, half-step for ds/2
416 int i;
417 for (i = 0; i < nx; ++i) {
418 qr_[i] = q[i]*expW_[i];
419 qr2_[i] = q[i]*expW2_[i];
420 }
421 fft_.forwardTransform(qr_, qk_);
422 fft_.forwardTransform(qr2_, qk2_);
423 for (i = 0; i < nk; ++i) {
424 qk_[i][0] *= expKsq_[i];
425 qk_[i][1] *= expKsq_[i];
426 qk2_[i][0] *= expKsq2_[i];
427 qk2_[i][1] *= expKsq2_[i];
428 }
429 fft_.inverseTransform(qk_, qr_);
430 fft_.inverseTransform(qk2_, qr2_);
431 for (i = 0; i < nx; ++i) {
432 qr_[i] = qr_[i]*expW_[i];
433 qr2_[i] = qr2_[i]*expW_[i];
434 }
435
436 // Finish second half-step for ds/2
437 fft_.forwardTransform(qr2_, qk2_);
438 for (i = 0; i < nk; ++i) {
439 qk2_[i][0] *= expKsq2_[i];
440 qk2_[i][1] *= expKsq2_[i];
441 }
442 fft_.inverseTransform(qk2_, qr2_);
443 for (i = 0; i < nx; ++i) {
444 qr2_[i] = qr2_[i]*expW2_[i];
445 }
446
447 // Richardson extrapolation
448 for (i = 0; i < nx; ++i) {
449 qNew[i] = (4.0*qr2_[i] - qr_[i])/3.0;
450 }
451 }
452
453}
454}
455#endif
virtual void setLength(double length)
Set the length of this block.
Class template for a block in a block copolymer.
Definition: BlockTmpl.h:89
Propagator< D > & propagator(int directionId)
Get a Propagator for a specified direction.
Definition: BlockTmpl.h:174
An IntVec<D, T> is a D-component vector of elements of integer type T.
Definition: IntVec.h:27
Iterator over points in a Mesh<D>.
Definition: MeshIterator.h:29
int rank() const
Get the rank of current element.
Definition: MeshIterator.h:136
void begin()
Set iterator to the first point in the mesh.
bool atEnd() const
Is this the end (i.e., one past the last point)?
Definition: MeshIterator.h:141
void setDimensions(const IntVec< D > &dimensions)
Set the grid dimensions in all directions.
IntVec< D > position() const
Get current position in the grid, as integer vector.
Definition: MeshIterator.h:122
Description of a regular grid of points in a periodic domain.
Definition: Mesh.h:61
IntVec< D > dimensions() const
Get an IntVec<D> of the grid dimensions.
Definition: Mesh.h:217
int size() const
Get total number of grid points.
Definition: Mesh.h:229
Block within a branched polymer.
void computeStress(double prefactor)
Compute stress contribution for this block.
void step(RField< D > const &q, RField< D > &qNew)
Compute one step of solution of MDE, from step i to i+1.
void setDiscretization(double ds, const Mesh< D > &mesh)
Initialize discretization and allocate required memory.
void computeConcentration(double prefactor)
Compute concentration (volume fraction) for block by integration.
void setLength(double newLength)
Set or reset block length.
void setupUnitCell(const UnitCell< D > &unitCell)
Setup parameters that depend on the unit cell.
void setKuhn(double kuhn)
Set or reset monomer statistical segment length.
void setupSolver(RField< D > const &w)
Set solver for this block.
Field of real double precision values on an FFT mesh.
Definition: RField.h:29
Base template for UnitCell<D> classes, D=1, 2 or 3.
Definition: UnitCell.h:44
int capacity() const
Return allocated size.
Definition: Array.h:159
bool isAllocated() const
Return true if this DArray has been allocated, false otherwise.
Definition: DArray.h:247
A fixed capacity (static) contiguous array with a variable logical size.
Definition: FSArray.h:38
void append(Data const &data)
Append data to the end of the array.
Definition: FSArray.h:258
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
Definition: global.h:68
C++ namespace for polymer self-consistent field theory (PSCF).
Utility classes for scientific computation.
Definition: accumulators.mod:1