12#include "Propagator.h" 
   14#include <prdc/cpu/WaveList.h> 
   15#include <prdc/cpu/FFT.h> 
   16#include <prdc/crystal/UnitCell.h> 
   17#include <prdc/crystal/shiftToMinimum.h> 
   19#include <pscf/chem/Edge.h> 
   20#include <pscf/mesh/Mesh.h> 
   21#include <pscf/mesh/MeshIterator.h> 
   22#include <pscf/math/IntVec.h> 
   23#include <pscf/solvers/BlockTmpl.tpp> 
   25#include <util/containers/DMatrix.h> 
   26#include <util/containers/DArray.h> 
   27#include <util/containers/FSArray.h> 
   44      unitCellPtr_(nullptr),
 
   45      waveListPtr_(nullptr),
 
 
   81      waveListPtr_ = &wavelist;
 
 
   98      UTIL_CHECK(mesh().dimensions() == fft().meshDimensions());
 
  105      expKsq_.allocate(kMeshDimensions_);
 
  106      expKsq2_.allocate(kMeshDimensions_);
 
  107      expW_.allocate(mesh().dimensions());
 
  108      qr_.allocate(mesh().dimensions());
 
  109      qk_.allocate(mesh().dimensions());
 
  110      qr2_.allocate(mesh().dimensions());
 
  111      qk2_.allocate(mesh().dimensions());
 
  113         expW2_.allocate(mesh().dimensions());
 
  115         expWInv_.allocate(mesh().dimensions());
 
  119      cField().allocate(mesh().dimensions());
 
  128         tempNs = floor(
length() / (2.0 *
ds) + 0.5);
 
  133         ds_ = 
length()/double(ns_-1);
 
 
  164         tempNs = floor( 
length()/(2.0 *dsTarget_) + 0.5 );
 
  169         ds_ = 
length()/double(ns_-1);
 
 
  205   void Block<D>::computeExpKsq()
 
  211      double bSqFactor = -1.0*kuhn()*kuhn() / 6.0;
 
  217      if (!waveListPtr_->hasKSq()) {
 
  218         waveListPtr_->computeKSq();
 
  220      RField<D> const & kSq = waveListPtr_->kSq();
 
  228         arg = kSq[i]*bSqFactor;
 
  229         expKsq_[i] = exp(arg);
 
  230         expKsq2_[i] = exp(0.5*arg);
 
  244      int nx = mesh().size();
 
  252         for (
int i = 0; i < nx; ++i) {
 
  255            expW2_[i] = exp(0.5*arg);
 
  259         for (
int i = 0; i < nx; ++i) {
 
  262            expWInv_[i] = 1.0/expW_[i];
 
 
  282      int nx = mesh().size();
 
  285      UTIL_CHECK(mesh().dimensions() == fft().meshDimensions());
 
  289      int nk = qk_.capacity();
 
  308      for (i = 0; i < nx; ++i) {
 
  309         qr_[i] = q[i]*expW_[i];
 
  310         qr2_[i] = q[i]*expW2_[i];
 
  312      fft().forwardTransform(qr_, qk_);
 
  313      fft().forwardTransform(qr2_, qk2_);
 
  314      for (i = 0; i < nk; ++i) {
 
  315         qk_[i][0] *= expKsq_[i];
 
  316         qk_[i][1] *= expKsq_[i];
 
  317         qk2_[i][0] *= expKsq2_[i];
 
  318         qk2_[i][1] *= expKsq2_[i];
 
  320      fft().inverseTransformUnsafe(qk_, qr_); 
 
  321      fft().inverseTransformUnsafe(qk2_, qr2_); 
 
  322      for (i = 0; i < nx; ++i) {
 
  323         qr_[i] = qr_[i]*expW_[i];
 
  324         qr2_[i] = qr2_[i]*expW_[i];
 
  332      fft().forwardTransform(qr2_, qk2_);
 
  333      for (i = 0; i < nk; ++i) {
 
  334         qk2_[i][0] *= expKsq2_[i];
 
  335         qk2_[i][1] *= expKsq2_[i];
 
  337      fft().inverseTransformUnsafe(qk2_, qr2_); 
 
  338      for (i = 0; i < nx; ++i) {
 
  339         qr2_[i] = qr2_[i]*expW2_[i];
 
  343      for (i = 0; i < nx; ++i) {
 
  344         qout[i] = (4.0*qr2_[i] - qr_[i])/3.0;
 
 
  369      int nx = mesh().size();
 
  373      int nk = qk_.capacity();
 
  378      fft().forwardTransform(q, qk_);
 
  379      for (
int i = 0; i < nk; ++i) {
 
  380         qk_[i][0] *= expKsq_[i];
 
  381         qk_[i][1] *= expKsq_[i];
 
  383      fft().inverseTransformUnsafe(qk_, qout); 
 
 
  395      int nx = mesh().size();
 
  399      int nk = qk_.capacity();
 
  404      fft().forwardTransform(q, qk_);
 
  405      for (
int i = 0; i < nk; ++i) {
 
  406         qk_[i][0] *= expKsq2_[i];
 
  407         qk_[i][1] *= expKsq2_[i];
 
  409      fft().inverseTransformUnsafe(qk_, qout); 
 
 
  419      int nx = mesh().size();
 
  425      for (
int i = 0; i < nx; ++i) {
 
 
  438      int nx = mesh().size();
 
  447      for (i = 0; i < nx; ++i) {
 
  458      for (i = 0; i < nx; ++i) {
 
  459         cField()[i] += p0.
q(0)[i]*p1.
q(ns_ - 1)[i];
 
  460         cField()[i] += p0.
q(ns_ -1)[i]*p1.
q(0)[i];
 
  465      for (j = 1; j < (ns_ -1); j += 2) {
 
  468         for (i = 0; i < nx; ++i) {
 
  470            cField()[i] += qf[i] * qr[i] * 4.0;
 
  475      for (j = 2; j < (ns_ -2); j += 2) {
 
  478         for (i = 0; i < nx; ++i) {
 
  480            cField()[i] += qf[i] * qr[i] * 2.0;
 
  485      prefactor *= ds_ / 3.0;
 
  486      for (i = 0; i < nx; ++i) {
 
 
  500      int nx = mesh().size();
 
  509      for (i = 0; i < nx; ++i) {
 
  519      for (j = 1; j < (ns_ -1); ++j) {
 
  522         for (i = 0; i < nx; ++i) {
 
  523            cField()[i] += qf[i] * qr[i] * expWInv_[i];
 
  528      for (i = 0; i < nx; ++i) {
 
 
  541      int nx = mesh().size();
 
  544      UTIL_CHECK(mesh().dimensions() == fft().meshDimensions());
 
  556      if (!waveListPtr_->hasdKSq()) {
 
  557         waveListPtr_->computedKSq();
 
  562      const int nParam = unitCell().nParameter();
 
  563      for (
int i = 0; i < nParam; ++i) {
 
  568      const double bSq = 
kuhn()*
kuhn()/6.0;
 
  569      double dels, prod, increment;
 
  573      for (
int j = 0; j < ns_ ; ++j) {
 
  576         fft().forwardTransform(qr_, qk_);
 
  578         qr2_ = p1.
q(ns_ - 1 - j);
 
  579         fft().forwardTransform(qr2_, qk2_);
 
  583         if (j != 0 && j != ns_ - 1) {
 
  592         for (n = 0; n < nParam ; ++n) {
 
  596            for (m = 0; m < kSize_ ; ++m) {
 
  597               prod = (qk2_[m][0] * qk_[m][0]) + (qk2_[m][1] * qk_[m][1]);
 
  601            increment *= bSq * dels;
 
  602            dQ[n] = dQ[n] - increment;
 
  609      for (
int n = 0; n < nParam; ++n) {
 
  610         stress_.append(-1.0*prefactor*dQ[n]);
 
 
  623      int nx = mesh().size();
 
  626      UTIL_CHECK(mesh().dimensions() == fft().meshDimensions());
 
  638      if (!waveListPtr_->hasdKSq()) {
 
  639         waveListPtr_->computedKSq();
 
  644      int nParam = unitCell().nParameter();
 
  645      for (
int i = 0; i < nParam; ++i) {
 
  649      const double bSq = 
kuhn()*
kuhn()/6.0;
 
  650      double increment, prod;
 
  657         fft().forwardTransform(qr_, qk_);
 
  660         qr2_ = p1.
q(ns_ - 2);
 
  661         fft().forwardTransform(qr2_, qk2_);
 
  664         for (
int n = 0; n < nParam ; ++n) {
 
  668            for (
int m = 0; m < kSize_ ; ++m) {
 
  669               prod = (qk2_[m][0] * qk_[m][0]) + (qk2_[m][1] * qk_[m][1]);
 
  670               prod *= dKSq[m]*expKsq2_[m];
 
  673            increment *= 0.5*bSq;
 
  674            dQ[n] = dQ[n] - increment;
 
  679      for (
int j = 1; j < ns_ - 2 ; ++j) {
 
  683         fft().forwardTransform(qr_, qk_);
 
  686         qr2_ = p1.
q(ns_ - 2 - j);
 
  687         fft().forwardTransform(qr2_, qk2_);
 
  690         for (
int n = 0; n < nParam ; ++n) {
 
  694            for (
int m = 0; m < kSize_ ; ++m) {
 
  695               prod = (qk2_[m][0] * qk_[m][0]) + (qk2_[m][1] * qk_[m][1]);
 
  696               prod *= dKSq[m]*expKsq_[m];
 
  700            dQ[n] = dQ[n] - increment;
 
  710         fft().forwardTransform(qr_, qk_);
 
  714         fft().forwardTransform(qr2_, qk2_);
 
  717         for (
int n = 0; n < nParam ; ++n) {
 
  721            for (
int m = 0; m < kSize_ ; ++m) {
 
  722               prod = (qk2_[m][0] * qk_[m][0]) + (qk2_[m][1] * qk_[m][1]);
 
  723               prod *= dKSq[m]*expKsq2_[m];
 
  726            increment *= 0.5*bSq;
 
  727            dQ[n] = dQ[n] - increment;
 
  734      for (
int i = 0; i < nParam; ++i) {
 
  735         stress_.append(-1.0*prefactor*dQ[i]);
 
 
virtual void setKuhn(double kuhn)
virtual void setLength(double length)
Set the length of this block (only valid for thread model).
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.
Description of a regular grid of points in a periodic domain.
Fourier transform wrapper.
static void computeKMesh(IntVec< D > const &rMeshDimensions, IntVec< D > &kMeshDimensions, int &kSize)
Compute dimensions and size of k-space mesh for DFT of real data.
bool isAllocated() const
Return true if the FftwDArray has been allocated, false otherwise.
Field of real double precision values on an FFT mesh.
Class to calculate and store properties of wavevectors.
Base template for UnitCell<D> classes, D=1, 2 or 3.
void stepBondBead(RField< D > const &qin, RField< D > &qout) const
Apply a bond operator for the bead model.
void setLength(double newLength)
Set or reset block length (only used in thread model).
void computeConcentrationBead(double prefactor)
Compute the concentration for this block, using the bead model.
void setKuhn(double kuhn)
Set or reset monomer statistical segment length.
void clearUnitCellData()
Clear all internal data that depends on the unit cell parameters.
void setupSolver(RField< D > const &w)
Set up the MDE solver for this block.
double kuhn() const
Get monomer statistical segment length.
void allocate(double ds)
Allocate memory and set contour step size.
void computeStressBead(double prefactor)
Compute stress contribution for this block, using bead model.
int nBead() const
Get the number of beads in this block, in the bead model.
void associate(Mesh< D > const &mesh, FFT< D > const &fft, UnitCell< D > const &cell, WaveList< D > &wavelist)
Create permanent associations with related objects.
Propagator< D > & propagator(int directionId)
Get a Propagator for a specified direction.
void computeConcentrationThread(double prefactor)
Compute the concentration for this block, for the thread model.
void stepThread(RField< D > const &qin, RField< D > &qout) const
Compute one step of solution of MDE for the thread model.
void computeStressThread(double prefactor)
Compute stress contribution for this block, using thread model.
double ds() const
Get contour length step size.
void stepFieldBead(RField< D > &q) const
Apply the exponential field operator for the bead model.
void stepHalfBondBead(RField< D > const &qin, RField< D > &qout) const
Apply a half-bond operator for the bead model.
RField< D > & cField()
Get the associated monomer concentration field.
void stepBead(RField< D > const &qin, RField< D > &qout) const
Compute one step of solution of MDE for the bead model.
double length() const
Get the length of this block, in the thread model.
MDE solver for one direction of one block.
bool isHeadEnd() const
Is the head vertex a chain end?
bool isSolved() const
Has the modified diffusion equation been solved?
bool isTailEnd() const
Is the tail vertex a chain end?
bool isAllocated() const
Has memory been allocated for this propagator?
const FieldT & q(int i) const
Return q-field at specified step.
int capacity() const
Return allocated size.
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.
bool isThread()
Is the thread model in use ?
bool isBead()
Is the bead model in use ?
Fields and FFTs for periodic boundary conditions (CPU)
Periodic fields and crystallography.
Real periodic fields, SCFT and PS-FTS (CPU).
PSCF package top-level namespace.