PSCF v1.3
Pscf::Rpc::Block< D > Class Template Reference

Block within a linear or branched block polymer. More...

#include <Block.h>

Inheritance diagram for Pscf::Rpc::Block< D >:
Pscf::BlockTmpl< Propagator< D > > Pscf::Edge

Public Types

using Base = BlockTmpl< Propagator<D> >
 Base class.
using PropagatorT = Propagator<D>
 Propagator type.
using FFTT = FFT<D>
 Fast Fourier Transform type.
using WaveListT = WaveList<D>
 WaveList type.
using FieldIoT = FieldIo<D>
 FieldIo type.
Public Types inherited from Pscf::BlockTmpl< Propagator< D > >
using PropagatorT
 Modified diffusion equation solver (propagator) type.

Public Member Functions

 Block ()
 Constructor.
 ~Block ()
 Destructor.
void associate (Mesh< D > const &mesh, FFT< D > const &fft, UnitCell< D > const &cell, WaveList< D > &wavelist)
 Create permanent associations with related objects.
void allocate (double ds)
 Allocate memory and set contour step size.
void clearUnitCellData ()
 Clear all internal data that depends on the unit cell parameters.
void setLength (double newLength)
 Set or reset block length (only used in thread model).
void setKuhn (double kuhn)
 Set or reset monomer statistical segment length.
void setupSolver (RField< D > const &w)
 Set up the MDE solver for this block.
void stepThread (RField< D > const &qin, RField< D > &qout) const
 Compute one step of solution of MDE for the thread model.
void stepBead (RField< D > const &qin, RField< D > &qout) const
 Compute one step of solution of MDE for the bead model.
void stepFieldBead (RField< D > &q) const
 Apply the exponential field operator for the bead model.
void stepBondBead (RField< D > const &qin, RField< D > &qout) const
 Apply a bond operator for the bead model.
void stepHalfBondBead (RField< D > const &qin, RField< D > &qout) const
 Apply a half-bond operator for the bead model.
void computeConcentrationThread (double prefactor)
 Compute the concentration for this block, for the thread model.
void computeConcentrationBead (double prefactor)
 Compute the concentration for this block, using the bead model.
void computeStressThread (double prefactor)
 Compute stress contribution for this block, using thread model.
void computeStressBead (double prefactor)
 Compute stress contribution for this block, using bead model.
double stress (int n) const
 Get derivative of free energy w/ respect to a unit cell parameter.
double ds () const
 Get contour length step size.
int ns () const
 Get the number of contour grid points, including end points.
Propagator< D > & propagator (int directionId)
 Get a Propagator for a specified direction.
Propagator< D > const & propagator (int directionId) const
 Get a const Propagator for a specified direction.
Propagator< D >::FieldT & cField ()
 Get the associated monomer concentration field.
Propagator< D >::FieldT const & cField () const
 Get the associated const monomer concentration field.
double kuhn () const
 Get monomer statistical segment length.
void setId (int id)
 Set the id for this block.
void setVertexIds (int vertexId0, int vertexId1)
 Set indices of associated vertices.
void setMonomerId (int monomerId)
 Set the monomer type id.
int id () const
 Get the id of this block (unique within the polymer).
int monomerId () const
 Get the monomer type id for this block.
const Pair< int > & vertexIds () const
 Get the pair of associated vertex ids.
int vertexId (int i) const
 Get the id of one associated vertex.
double length () const
 Get the length of this block, in the thread model.
int nBead () const
 Get the number of beads in this block, in the bead model.
Public Member Functions inherited from Pscf::BlockTmpl< Propagator< D > >
 BlockTmpl ()
 Constructor.
virtual ~BlockTmpl ()
 Destructor.
Propagator< D > & propagator (int directionId)
 Get a Propagator for a specified direction.
Propagator< D > const & propagator (int directionId) const
 Get a const Propagator for a specified direction.
Propagator< D >::FieldT & cField ()
 Get the associated monomer concentration field.
Propagator< D >::FieldT const & cField () const
 Get the associated const monomer concentration field.
double kuhn () const
 Get monomer statistical segment length.
Public Member Functions inherited from Pscf::Edge
 Edge ()
 Constructor.
virtual ~Edge ()
 Destructor.
template<class Archive>
void serialize (Archive &ar, unsigned int versionId)
 Serialize to/from archive.
void setId (int id)
 Set the id for this block.
void setMonomerId (int monomerId)
 Set the monomer type id.
void setVertexIds (int vertexId0, int vertexId1)
 Set indices of associated vertices.
virtual void setNBead (int nBead)
 Set the number of beads in this block (only valid for bead model).
void setPolymerType (PolymerType::Enum type)
 Set the type of the parent polymer (branched or linear).
int id () const
 Get the id of this block (unique within the polymer).
int monomerId () const
 Get the monomer type id for this block.
const Pair< int > & vertexIds () const
 Get the pair of associated vertex ids.
int vertexId (int i) const
 Get the id of one associated vertex.
double length () const
 Get the length of this block, in the thread model.
int nBead () const
 Get the number of beads in this block, in the bead model.
PolymerType::Enum polymerType () const
 Get the type of the parent polymer (branched or linear).

Detailed Description

template<int D>
class Pscf::Rpc::Block< D >

Block within a linear or branched block polymer.

Derived from BlockTmpl< Propagator<D> >. A BlockTmpl< Propagator<D> > has two Propagator<D> members, and is derived from class Pscf::Edge.

Manual Page

Definition at line 51 of file rpc/solvers/Block.h.

Member Typedef Documentation

◆ Base

template<int D>
using Pscf::Rpc::Block< D >::Base = BlockTmpl< Propagator<D> >

Base class.

Definition at line 59 of file rpc/solvers/Block.h.

◆ PropagatorT

template<int D>
using Pscf::Rpc::Block< D >::PropagatorT = Propagator<D>

Propagator type.

Definition at line 62 of file rpc/solvers/Block.h.

◆ FFTT

template<int D>
using Pscf::Rpc::Block< D >::FFTT = FFT<D>

Fast Fourier Transform type.

Definition at line 65 of file rpc/solvers/Block.h.

◆ WaveListT

template<int D>
using Pscf::Rpc::Block< D >::WaveListT = WaveList<D>

WaveList type.

Definition at line 68 of file rpc/solvers/Block.h.

◆ FieldIoT

template<int D>
using Pscf::Rpc::Block< D >::FieldIoT = FieldIo<D>

FieldIo type.

Definition at line 71 of file rpc/solvers/Block.h.

Constructor & Destructor Documentation

◆ Block()

template<int D>
Pscf::Rpc::Block< D >::Block ( )

Constructor.

Definition at line 39 of file rpc/solvers/Block.tpp.

References propagator().

◆ ~Block()

template<int D>
Pscf::Rpc::Block< D >::~Block ( )

Destructor.

Definition at line 60 of file rpc/solvers/Block.tpp.

Member Function Documentation

◆ associate()

template<int D>
void Pscf::Rpc::Block< D >::associate ( Mesh< D > const & mesh,
FFT< D > const & fft,
UnitCell< D > const & cell,
WaveList< D > & wavelist )

Create permanent associations with related objects.

This function creates associations of this block with the Mesh, FFT, UnitCell and WaveList objects by storing their addresses. It must be called before allocate().

Parameters
meshMesh<D> object, spatial discretization meth
fftFFT<D> object, Fast Fourier Transform
cellUnitCell<D> object, crystallographic unit cell
wavelistWaveList<D>, container for wavevector properties

Definition at line 67 of file rpc/solvers/Block.tpp.

References UTIL_CHECK.

◆ allocate()

template<int D>
void Pscf::Rpc::Block< D >::allocate ( double ds)

Allocate memory and set contour step size.

This function choses a value for the number ns of contour variable grid points for this block, sets the step size, and allocates memory for several private arrays. Spatial grid dimensions are obtained from a pointers to the associated mesh. The associate function must be called before this function.

For the thread model, if PolymerModel::isThread() is true, the value for the number ns of contour variable grid points for this block is chosen to yield a value for the the actual step size length/(ns-1) as close as possible to the input parameter ds (the target step size), consistent with the requirements that ns be an odd integer and ns > 1. These requirements allow use of Simpson's rule for integration with respect to the contour variable s to compute monomer concentration fields and stress contributions.

For the bead model, if PolymerModel::isThread() is true, the value of ns is given by nBead + 2.

The value of input parameter ds is ignored for the bead model.

Parameters
dsdesired (optimal) value for contour length step

Definition at line 88 of file rpc/solvers/Block.tpp.

References cField(), Pscf::Prdc::Cpu::FFT< D >::computeKMesh(), ds(), Pscf::PolymerModel::isBead(), Pscf::PolymerModel::isThread(), length(), nBead(), propagator(), and UTIL_CHECK.

◆ clearUnitCellData()

template<int D>
void Pscf::Rpc::Block< D >::clearUnitCellData ( )

Clear all internal data that depends on the unit cell parameters.

This function must be called once after every time the unit cell parameters change. The function marks all class member variables that depend on the unit cell parameters as being outdated. All such variables are then recomputed just before they are needed.

Definition at line 193 of file rpc/solvers/Block.tpp.

◆ setLength()

template<int D>
void Pscf::Rpc::Block< D >::setLength ( double newLength)
virtual

Set or reset block length (only used in thread model).

Precondition: PolymerModel::isThread(). An Exception is thrown if this function is called when PolymerModel::isThread() is false.

Parameters
newLengthnew block length

Reimplemented from Pscf::Edge.

Definition at line 150 of file rpc/solvers/Block.tpp.

References Pscf::PolymerModel::isThread(), length(), propagator(), Pscf::Edge::setLength(), and UTIL_CHECK.

◆ setKuhn()

template<int D>
void Pscf::Rpc::Block< D >::setKuhn ( double kuhn)
virtual

Set or reset monomer statistical segment length.

Parameters
kuhnnew monomer statistical segment length.

Reimplemented from Pscf::BlockTmpl< Propagator< D > >.

Definition at line 183 of file rpc/solvers/Block.tpp.

References Pscf::BlockTmpl< Propagator< D > >::BlockTmpl(), kuhn(), and setKuhn().

Referenced by setKuhn().

◆ setupSolver()

template<int D>
void Pscf::Rpc::Block< D >::setupSolver ( RField< D > const & w)

Set up the MDE solver for this block.

This should be called once after every change in w fields, before entering the loop used to solve the MDE for either propagator. This function is called by Polymer<D>::compute.

Parameters
wchemical potential field for this monomer type

Definition at line 239 of file rpc/solvers/Block.tpp.

References Pscf::PolymerModel::isBead(), Pscf::PolymerModel::isThread(), and UTIL_CHECK.

◆ stepThread()

template<int D>
void Pscf::Rpc::Block< D >::stepThread ( RField< D > const & qin,
RField< D > & qout ) const

Compute one step of solution of MDE for the thread model.

This function is called internally by the Propagator::solve function within a loop over steps. It is implemented in the Block class because the same private data structures are needed for the two propagators associated with a Block.

Parameters
qininput slice of q, from step i
qoutoutput slice of q, from step i+1

Definition at line 275 of file rpc/solvers/Block.tpp.

References Util::Array< Data >::capacity(), Pscf::Prdc::Cpu::FftwDArray< Data >::isAllocated(), Pscf::PolymerModel::isThread(), and UTIL_CHECK.

◆ stepBead()

template<int D>
void Pscf::Rpc::Block< D >::stepBead ( RField< D > const & qin,
RField< D > & qout ) const

Compute one step of solution of MDE for the bead model.

This function is called internally by the Propagator::solve function within a loop over steps. It is implemented in the Block class because the same private data structures are needed for the two propagators associated with a Block.

Parameters
qininput slice of q, from step i
qoutoutput slice of q, for step i+1

Definition at line 351 of file rpc/solvers/Block.tpp.

References Pscf::PolymerModel::isBead(), stepBondBead(), stepFieldBead(), and UTIL_CHECK.

◆ stepFieldBead()

template<int D>
void Pscf::Rpc::Block< D >::stepFieldBead ( RField< D > & q) const

Apply the exponential field operator for the bead model.

This function applies exp( -w(r) ), where w(r) is the w-field for the monomer type of this block.

Parameters
qslice of propagator q, modified in place

Definition at line 414 of file rpc/solvers/Block.tpp.

References Util::Array< Data >::capacity(), and UTIL_CHECK.

Referenced by stepBead().

◆ stepBondBead()

template<int D>
void Pscf::Rpc::Block< D >::stepBondBead ( RField< D > const & qin,
RField< D > & qout ) const

Apply a bond operator for the bead model.

This function applies exp( nabla^2 b^2 / 6 ), where nabla^2 denotes a Laplacian operator with eigenvalues given by -G^2 for reciprocal lattice vectors.

Parameters
qininput slice of q, from step i
qoutouptut slice of q, for step i+1

Definition at line 362 of file rpc/solvers/Block.tpp.

References Util::Array< Data >::capacity(), and UTIL_CHECK.

Referenced by stepBead().

◆ stepHalfBondBead()

template<int D>
void Pscf::Rpc::Block< D >::stepHalfBondBead ( RField< D > const & qin,
RField< D > & qout ) const

Apply a half-bond operator for the bead model.

This function applies exp( nabla^2 b^2 / 12 ), where nabla^2 denotes a Laplacian operator with eigenvalues given by -G^2 for reciprocal lattice vectors. It is used in the Propagator::solve function to deal with half-bonds at block ends.

Parameters
qininput slice of q, from step i
qoutouptut slice of q, for step i+1

Definition at line 388 of file rpc/solvers/Block.tpp.

References Util::Array< Data >::capacity(), and UTIL_CHECK.

◆ computeConcentrationThread()

template<int D>
void Pscf::Rpc::Block< D >::computeConcentrationThread ( double prefactor)

Compute the concentration for this block, for the thread model.

This function is called by Polymer::compute if a thread model is is used.

The "prefactor" parameter must equal \( \phi / (L_{tot} Q) \), where \( \phi \) is the species volume fraction, \( L_{tot} \) is the total length of all blocks in this polymer species and Q is the species partition function.

Upon return, grid point r of the array returned by the member function cField() contains the integal

\[ p \int_{0}^{l} ds q_{0}(r,s) q_{1}(r, L - s) \]

where \( q_{0}(r,s) \) and \( q_{1}(r,s) \) are propagators associated with different directions, \( p \) is the prefactor parameter, and the integral is taken over the length \( L \) of this block. Simpson's rule is used for the integral with respect to s.

Parameters
prefactorconstant multiplying integral over s

Definition at line 432 of file rpc/solvers/Block.tpp.

References cField(), propagator(), Pscf::Rpc::Propagator< D >::q(), and UTIL_CHECK.

◆ computeConcentrationBead()

template<int D>
void Pscf::Rpc::Block< D >::computeConcentrationBead ( double prefactor)

Compute the concentration for this block, using the bead model.

This function is called by Polymer::compute if a bead model is is used.

The "prefactor" parameter must equal \( \phi /(N_{tot} Q) \), where \( \phi \) is the species volume fraction, \( N_{tot} \) is the total number of beads in all blocks of the polymer, and \( Q \) is the species partition function.

Upon return, grid point r of the array returned by member function cField() contains the sum

\[ p \sum_{s} q_{0}(r,s) q_{1}(r, N-s) \exp(W(r)*ds) \]

where \( q_{0}(r,s) \) and \( q_{1}(r, N-s) \) denote complementary propagator slices associated with different
directions but the same bead, and \( p \) is the prefactor parameter. The sum is taken over all beads in this block.

Parameters
prefactorconstant multiplying sum over beads

Definition at line 494 of file rpc/solvers/Block.tpp.

References cField(), propagator(), Pscf::Rpc::Propagator< D >::q(), and UTIL_CHECK.

◆ computeStressThread()

template<int D>
void Pscf::Rpc::Block< D >::computeStressThread ( double prefactor)

Compute stress contribution for this block, using thread model.

This function is called by Polymer<D>::computeStress. The prefactor parameter must be equal to that passed to function computeConcentrationThread(double ).

Parameters
prefactorconstant multiplying integral over s

Definition at line 535 of file rpc/solvers/Block.tpp.

References Util::FSArray< Data, Capacity >::append(), Pscf::Rpc::Propagator< D >::isAllocated(), Pscf::PolymerModel::isThread(), kuhn(), propagator(), Pscf::Rpc::Propagator< D >::q(), and UTIL_CHECK.

◆ computeStressBead()

template<int D>
void Pscf::Rpc::Block< D >::computeStressBead ( double prefactor)

Compute stress contribution for this block, using bead model.

This function is called by Polymer<D>::computeStress. The prefactor parameter must be equal to that passed to function computeConcentrationBead(double ).

Parameters
prefactorconstant multiplying sum over beads

Definition at line 617 of file rpc/solvers/Block.tpp.

References Util::FSArray< Data, Capacity >::append(), Pscf::PolymerModel::isBead(), Pscf::Rpc::Propagator< D >::isHeadEnd(), Pscf::Rpc::Propagator< D >::isSolved(), Pscf::Rpc::Propagator< D >::isTailEnd(), kuhn(), propagator(), Pscf::Rpc::Propagator< D >::q(), and UTIL_CHECK.

◆ stress()

template<int D>
double Pscf::Rpc::Block< D >::stress ( int n) const
inline

Get derivative of free energy w/ respect to a unit cell parameter.

This function returns a value computed by a previous call to the computeStress function.

Parameters
nindex of the unit cell parameter

Definition at line 462 of file rpc/solvers/Block.h.

◆ ds()

template<int D>
double Pscf::Rpc::Block< D >::ds ( ) const
inline

Get contour length step size.

Definition at line 457 of file rpc/solvers/Block.h.

Referenced by allocate().

◆ ns()

template<int D>
int Pscf::Rpc::Block< D >::ns ( ) const
inline

Get the number of contour grid points, including end points.

Thread mdoel: For the thread model, ns is always an odd number ns >= 3 chosen to give an even number ns - 1 of steps of with a step length ds = length/(ns - 1) as close as possible to the target value of ds passed to the allocate function.

Bead model: For the bead model, ns is equal to nBead + 2, so as to include slices for all beads and two attached phantom vertices. Beads are indexed 1, ..., ns - 2, vertices have indices 0 and ns - 1.

Definition at line 452 of file rpc/solvers/Block.h.

◆ propagator() [1/2]

template<int D>
Propagator< D > & Pscf::BlockTmpl< Propagator< D > >::propagator ( int directionId)

Get a Propagator for a specified direction.

For a block with v0 = vertexId(0) and v1 = vertexId(1), propagator(0) propagates from vertex v0 to v1, while propagator(1) propagates from vertex v1 to v0.

Parameters
directionIdinteger index for direction (0 or 1)

Referenced by allocate(), Block(), computeConcentrationBead(), computeConcentrationThread(), computeStressBead(), computeStressThread(), and setLength().

◆ propagator() [2/2]

template<int D>
Propagator< D > const & Pscf::BlockTmpl< Propagator< D > >::propagator ( int directionId) const

Get a const Propagator for a specified direction.

See above for number conventions.

Parameters
directionIdinteger index for direction (0 or 1)

◆ cField() [1/2]

template<int D>
Propagator< D >::FieldT & Pscf::BlockTmpl< Propagator< D > >::cField ( )

Get the associated monomer concentration field.

Referenced by allocate(), computeConcentrationBead(), and computeConcentrationThread().

◆ cField() [2/2]

template<int D>
Propagator< D >::FieldT const & Pscf::BlockTmpl< Propagator< D > >::cField ( ) const

Get the associated const monomer concentration field.

◆ kuhn()

template<int D>
double Pscf::BlockTmpl< Propagator< D > >::kuhn ( ) const

Get monomer statistical segment length.

Referenced by computeStressBead(), computeStressThread(), and setKuhn().

◆ setId()

template<int D>
void Pscf::Edge::setId ( int id)

Set the id for this block.

Parameters
idinteger index for this block

Definition at line 89 of file Edge.cpp.

References Pscf::PolymerModel::isThread().

◆ setVertexIds()

template<int D>
void Pscf::Edge::setVertexIds ( int vertexId0,
int vertexId1 )

Set indices of associated vertices.

Parameters
vertexId0integer id of vertex 0
vertexId1integer id of vertex 1

Definition at line 104 of file Edge.cpp.

◆ setMonomerId()

template<int D>
void Pscf::Edge::setMonomerId ( int monomerId)

Set the monomer type id.

Parameters
monomerIdinteger id of monomer type

Definition at line 96 of file Edge.cpp.

◆ id()

template<int D>
int Pscf::Edge::id ( ) const
inline

Get the id of this block (unique within the polymer).

Definition at line 147 of file Edge.h.

◆ monomerId()

template<int D>
int Pscf::Edge::monomerId ( ) const
inline

Get the monomer type id for this block.

Definition at line 152 of file Edge.h.

◆ vertexIds()

template<int D>
const Pair< int > & Pscf::Edge::vertexIds ( ) const
inline

Get the pair of associated vertex ids.

Definition at line 157 of file Edge.h.

◆ vertexId()

template<int D>
int Pscf::Edge::vertexId ( int i) const
inline

Get the id of one associated vertex.

Parameters
iindex of vertex (0 or 1)

Definition at line 164 of file Edge.h.

◆ length()

template<int D>
double Pscf::Edge::length ( ) const
inline

Get the length of this block, in the thread model.

Precondition: PolymerModel::isThread()

Definition at line 171 of file Edge.h.

Referenced by allocate(), and setLength().

◆ nBead()

template<int D>
int Pscf::Edge::nBead ( ) const
inline

Get the number of beads in this block, in the bead model.

Precondition: PolymerModel::isBead()

Definition at line 178 of file Edge.h.

Referenced by allocate().


The documentation for this class was generated from the following files: