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>
62 tempNs = floor( length()/(2.0 *ds) + 0.5 );
67 ds_ = length()/double(ns_-1);
75 for (
int i = 0; i < D; ++i) {
79 kMeshDimensions_[i] = mesh.
dimensions()[i]/2 + 1;
85 for (
int i = 0; i < D; ++i) {
86 kSize *= kMeshDimensions_[i];
90 expKsq_.allocate(kMeshDimensions_);
91 expKsq2_.allocate(kMeshDimensions_);
100 dGsq_.allocate(kSize, 6);
106 propagator(0).allocate(ns_, mesh);
107 propagator(1).allocate(ns_, mesh);
126 tempNs = floor( length()/(2.0 *dsTarget_) + 0.5 );
131 ds_ = length()/double(ns_-1);
136 propagator(0).reallocate(ns_);
137 propagator(1).reallocate(ns_);
160 unitCellPtr_ = &unitCell;
175 double factor = -1.0*kuhn()*kuhn()*ds_/6.0;
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);
197 int nx = mesh().size();
202 for (
int i = 0; i < nx; ++i) {
210 expW_[i] = exp(-0.5*w[i]*ds_);
211 expW2_[i] = exp(-0.5*0.5*w[i]*ds_);
229 int nx = mesh().size();
239 for (i = 0; i < nx; ++i) {
243 Propagator<D>
const & p0 = propagator(0);
244 Propagator<D>
const & p1 = propagator(1);
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];
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;
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;
270 prefactor *= ds_ / 3.0;
271 for (i = 0; i < nx; ++i) {
272 cField()[i] *= prefactor;
284 int nx = mesh().size();
293 double dels, normal, increment;
299 for (
int i = 0; i < D; ++i) {
300 kSize_ *= kMeshDimensions_[i];
302 nParam = unitCell().nParameter();
309 for (i = 0; i < nParam; ++i) {
316 Propagator<D>
const & p0 = propagator(0);
317 Propagator<D>
const & p1 = propagator(1);
320 for (
int j = 0; j < ns_ ; ++j) {
323 fft_.forwardTransform(qr_, qk_);
325 qr2_ = p1.q(ns_ - 1 - j);
326 fft_.forwardTransform(qr2_, qk2_);
330 if (j != 0 && j != ns_ - 1) {
338 for (
int n = 0; n < nParam ; ++n) {
341 for (m = 0; m < c ; ++m) {
343 prod = (qk2_[m][0] * qk_[m][0]) + (qk2_[m][1] * qk_[m][1]);
347 increment = (increment * kuhn() * kuhn() * dels)/normal;
348 dQ[n] = dQ[n] - increment;
353 for (i = 0; i < nParam; ++i) {
354 stress_[i] = stress_[i] - (dQ[i] * prefactor);
371 for (
int n = 0; n < unitCell().nParameter() ; ++n) {
374 vec = shiftToMinimum(temp, mesh().dimensions(), unitCell());
375 dGsq_(iter.
rank(), n) = unitCell().dksq(vec, n);
376 for (
int p = 0; p < D; ++p) {
378 Partner[p] = mesh().dimensions()[p] - temp[p];
383 if (Partner[D-1] > kMeshDimensions_[D-1]) {
384 dGsq_(iter.
rank(), n) *= 2;
398 int nx = mesh().size();
399 int nk = qk_.capacity();
417 for (i = 0; i < nx; ++i) {
418 qr_[i] = q[i]*expW_[i];
419 qr2_[i] = q[i]*expW2_[i];
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];
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];
437 fft_.forwardTransform(qr2_, qk2_);
438 for (i = 0; i < nk; ++i) {
439 qk2_[i][0] *= expKsq2_[i];
440 qk2_[i][1] *= expKsq2_[i];
442 fft_.inverseTransform(qk2_, qr2_);
443 for (i = 0; i < nx; ++i) {
444 qr2_[i] = qr2_[i]*expW2_[i];
448 for (i = 0; i < nx; ++i) {
449 qNew[i] = (4.0*qr2_[i] - qr_[i])/3.0;
virtual void setLength(double length)
Set the length of this block.
Class template for a block in a block copolymer.
Propagator< D > & propagator(int directionId)
Get a Propagator for a specified direction.
An IntVec<D, T> is a D-component vector of elements of integer type T.
Iterator over points in a Mesh<D>.
int rank() const
Get the rank of current element.
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)?
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.
Description of a regular grid of points in a periodic domain.
IntVec< D > dimensions() const
Get an IntVec<D> of the grid dimensions.
int size() const
Get total number of grid points.
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.
Base template for UnitCell<D> classes, D=1, 2 or 3.
int capacity() const
Return allocated size.
bool isAllocated() const
Return true if this DArray has been allocated, false otherwise.
A fixed capacity (static) contiguous array with a variable logical size.
void append(Data const &data)
Append data to the end of the array.
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
C++ namespace for polymer self-consistent field theory (PSCF).
Utility classes for scientific computation.