Reference

Documentation for Polyorder.jl's public interface.

Contents

Index

General interface

Polyorder.AbstractFieldType
AbstractField{T, N, S<:AbstractArray{T, N}, P} <: AbstractArray{T, N}
const AbstractField1D{T, S, P} = AbstractField{T, 1, S, P}
const AbstractField2D{T, S, P} = AbstractField{T, 2, S, P}
const AbstractField3D{T, S, P} = AbstractField{T, 3, S, P}

An abstract type describing the simulation cell as a discrete grid which can hold either potential fields or density fields data. Note that AbstractField and all its subtypes should conforms to the AbstractArray interface.

Mandate Fields

  • data:S: the actual data.
  • lattice::BravaisLattice: the Bravais lattice of the simulation cell.
  • specie: the polymer specie.
source
Polyorder.AuxiliaryFieldType
AuxiliaryField{T, N, S, P, PF, PB} <: AbstractField{T, N, S, P}
const AuxiliaryField1D{T, S, P, PF, PB} = AuxiliaryField{T, 1, S, P, PF, PB}
const AuxiliaryField2D{T, S, P, PF, PB} = AuxiliaryField{T, 2, S, P, PF, PB}
const AuxiliaryField3D{T, S, P, PF, PB} = AuxiliaryField{T, 3, S, P, PF, PB}

A field representing the potential field of a specie.

Fields

  • data::S: the potential field data, which should be a 1D, 2D or 3D AbstractArray.
  • lattice::BravaisLattice{N, P}: the Bravais lattice of the field.
  • specie::Symbol: the polymer specie of the field.
  • fft::PF: the forward FFT plan.
  • ifft::PB: the backward FFT plan.
source
Polyorder.DensityFieldType
DensityField{T, N, S, P} <: AbstractField{T, N, S, P}
const DensityField1D{T, S, P} = DensityField{T, 1, S, P}
const DensityField2D{T, S, P} = DensityField{T, 2, S, P}
const DensityField3D{T, S, P} = DensityField{T, 3, S, P}

A field representing the density of a specie.

Fields

  • data::S: the density data, which should be a 1D, 2D or 3D AbstractArray.
  • lattice::BravaisLattice{N, P}: the Bravais lattice of the field.
  • specie::Symbol: the polymer specie of the field.
source
Polyorder.AbstractPropagatorType
AbstractPropagator{S<:AbstractArray, P} <: AbstractVector{S}

An abstract type for propagators. A AbstractPropagator object is a vector of data that represents the propagator of a polymer block or a small molecule. The data is usually a vector of S type, where S is a subtype of AbstractArray. The propagator also has a relative length α of the whole molecule, which is important for computing density fields of multi-component systems.

Mandatory fields

  • data::Vector{S}: The data of the propagator.
  • α::P: The relative length of the whole molecule.
source
Polyorder.PropagatorType
Propagator{S, P, B<:PolymerBlock, V} <: AbstractPropagator{S, P}

A Propagator object represents a propagator associated with a polymer block. The common notation for a propagator is q.

Fields

  • data::Vector{S}: The data of the propagator.
  • Ns::Int: The number of contour steps along a block.
  • ds::P: The contour step size.
  • α::P: The relative length of the whole polymer. This is important for computing density fields of multi-component systems.
  • block::B: The corresponding polymer block.
  • direction::Pair{V,V}: The direction of the propagator. It is a pair of nodes (v1 => v2) in the graph of a BlockCopolymerGraph object.
source
Polyorder.PropagatorSmallType
PropagatorSmall{S, P} <: AbstractPropagator{S, P}

A PropagatorSmall object represents a propagator associated with a small molecule.

Fields

  • data::Vector{S}: The data of the propagator, length(data)=1 is assumed.
  • α::P: The relative length of the whole molecule. This is important for computing density fields of multi-component systems.
  • molecule::SmallMolecule: The corresponding small molecule.
source
Polyorder.SimpleFieldModelType
SimpleFieldModel <: AbstractFieldModelType

One speice corresponds to one auxiliary field. For Incompressible system, one additional field to ensure incompressibility.

source
Polyorder.FTolModeType
FTolMode <: AbstractTolMode

Use difference of free energy as tolerance for SCFT convergence.

source
Polyorder.find_field_idxFunction
find_field_idx(sp::Symbol, fields::AbstractVector{T}) where T<:AbstractField

Return the index of the first field with specie sp in fields.

source
Polyorder.find_fieldFunction
find_field(sp::Symbol, fields::AbstractVector{T}) where T<:AbstractField

Return the first field with specie sp in fields.

source
Polyorder.find_propagator_idxFunction
find_propagator_idx(direction, propagators)

Return a single propagator that has a direction as the input direction. direction is a Pair whose key and value, i.e. (key, value), forms an edge in the graph representation of the polymer chain. There are two propagators corresponding to the same edge (key, value), where the one with direction key => value is returned.

source
Polyorder.find_propagatorFunction
find_propagator(direction, propagators)

Find the propagator that has a direction as the input direction from a list of propagators and return it. If no propagator is found, nothing is returned.

source
Polyorder.find_propagator_pair_indicesFunction
find_propagator_pair_indices(block, propagators)

Return a pair of propagators corresonding to the same edge in the graph representation of the block copolymer. It is assumed that the pair of propagators MUST present next to each other.

source
Polyorder.block_dsFunction
block_ds(q::Propagator)

Return the contour step of the propagator q.

source
block_ds(scft::NoncyclicChainSCFT)

Return the contour steps for each block of each component in scft as a Vector{Dict{Pair, Float64}}.

source
Polyorder.block_NsFunction
block_Ns(q::AbstractPropagator)

Return the number of contour steps of the propagator q. If q is a PropagatorSmall, then the number of contour steps is 1.

source
Polyorder.resetFunction
reset(obj::T, args...; kwargs...)

Return a new object of type T based on obj with updates according to args and kwargs.

See also reset!.

source
Polyorder.coordinatesFunction
coordinates(grid::AbstractField1D, space::CrystalSpace=RCSpace())
coordinates(grid::AbstractField2D, space::CrystalSpace=RCSpace())
coordinates(grid::AbstractField3D, space::CrystalSpace=RCSpace())

Return the coordinates of all grid points for each dimension in the particular space given by the space argument.

Arguments

  • grid: an AbstractField object.
  • space::CrystalSpace=RCSpace(): the space of the returned coordinates. See Scattering.CrystalSpace for all supported spaces.
source
Polyorder.num_free_cellsizeFunction
num_free_cellsize(container)::Int

Return the number of free cell sizes of a container that can be treated as free variables for optimization.

source

RPA

Polyorder.RPA.compute_stability_limitFunction
compute_stability_limit(bcp::BlockCopolymer; χN0=1.0, k2_max=100.0)
compute_stability_limit(system::PolymerSystem; χN0=1.0, k2_max=100.0)

Compute the stability limit for a block copolymer or a polymer system. Note that only TwoSpeciesSystem is supported.

source
Polyorder.RPA.form_factorFunction
form_factor(chain::BlockCopolymer, k2::Real)
form_factor(chain::BlockCopolymer, k2::AbstractArray)

Compute the form factor of a block copolymer chain for scattering vector(s) using the RPA approximation.

The form factor is $P(k) = \sum_i \sum_j g_{ij}(k)$ where $i$ and $j$ iterates all blocks of the block copolymer. $g_{ij}(k) = g_D(k, f_i)$ when $i = j$, $g_{ij}(k) = h(k, f_i) h(k, f_j)$ when $i$ and $j$ have common block ends, and $g_{ij}(k) = h(k, f_i) l(k,d_{ij}) h(k, f_j)$ when $i$ and $j$ are at least one block apart. Here $d_{ij}$ is the sum of the lengths of blocks that connect block $i$ and $j$ with shortest path.

source
form_factor(chain::BlockCopolymer, specie1::Symbol, specie2::Symbol k2::Real)
form_factor(chain::BlockCopolymer, specie1::Symbol, specie2::Symbol, k2::AbstractArray)

Return a partial form factor of a block copolymer chain for two specified species. Suppose the specie1 is A, adn specie2 is B, then $P_{AB}(k^2)$ is computed.

source
form_factor(chain::BlockCopolymer, specie::Symbol, k2::Real)
form_factor(chain::BlockCopolymer, specie::Symbol, k2::AbstractArray)

Return a partial form factor of a block copolymer chain for a specified specie. Suppose the specie is A, then $P_{AA}(k^2)$ is computed.

source
Polyorder.RPA.structure_factorFunction
structure_factor(chain::BlockCopolymer, χN, k2::Real)
structure_factor(chain::BlockCopolymer, χN, k2::AbstractVector)

Compute the structure factor (scattering factor) of a two-specie block copolymer melt: $S = W / (S^t - 2χW)$ where $S^t$ is the sum of all elemments of matrix $S_{ij}$ where $i,j=A,B$, and $W$ is the determinant of the matrix. Note that here $S/N$ is actually returned.

References

  • Leibler, L. Theory of Microphase Separation in Block Copolymers. Macromolecules 1980, 13 (6), 1602–1617.
  • Ranjan, A. Theoretical Studies of Equilibrium Structure and Linear Response in Diblock Copolymer Melts, University of Minnesota, 2010. https://doi.org/10.1002/9781119990413.ch1.
  • Ranjan, A.; Qin, J.; Morse, D. C. Linear Response and Stability of Ordered Phases of Block Copolymer Melts. Macromolecules 2008, 41 (3), 942–954.
source
Polyorder.RPA.compute_SmatFunction
compute_Smat(bcp::BlockCopolymer, k2::Real)
compute_Smat(bcp::BlockCopolymer, k2::AbstractArray)
compute_Smat!(Smat::AbstractMatrix, bcp::BlockCopolymer, k2::AbstractArray)

Compute the scattering matrix (linear response of an ideal gas of a block copolymer chain) $S$ of a BlockCopolymer instance. $S_{ij}$ is the element of matrix $S$, where $i, j$ are the index of specie type. Take AB diblock copolymer for an example, we have $S = [S_{AA} S_{AB}; S_{BA} S_{BB}]$, where we have $S_{AB} = S_{BA}$. $S_{ij}$ can be readily computed by RPA.form_factor function in the RPA module.

Smat should be a Matrix{<:AbstractArray}, the matrix indices $i, j$ correspond to the species in the bcp. For example, if species(bcp) returns [:A, :B, :C], then Smat[2, 3] => $S_{BC}$.

Note that current implementation is neither memory nor speed efficient because we don't utilize the fact that $S$ matrix is symmetric.

Each element Smat should has the same size of k2, the corresponding element in each element of Smat is the $S$ value for that k2.

source
compute_Smat(system::PolymerSystem, k2::Real)
compute_Smat!(Smat::AbstractMatrix, system::PolymerSystem, k2::Real)
compute_Smat(system::PolymerSystem, k2::AbstractArray)
compute_Smat!(Smat::AbstractMatrix, system::PolymerSystem, k2::AbstractArray)

Compute the $S$ matrix for a polymer system by summing up contributions from all BlockCopolymer components. Note the coefficients of c.ϕ/c.α should be used to scale the contribution from different BlockCopolymer components.

source
Polyorder.RPA.compute_Smat_non_interactingFunction
compute_Smat_non_interacting(system::PolymerSystem, k2::Real)

Compute the Scattering matrix for a non-interacting polymer system, i.e. where $\chi_{ij}=0$ for all $(i,j)$ pairs. It is the $P$ matrix following the notations in my research notes (2022.4.23).

source
Polyorder.RPA.compute_Smat_interactingFunction
compute_Smat_interacting(system::PolymerSystem, k2::Real)

Compute the Scattering matrix for an interacting polymer system. It is the $\tilde{S}$ matrix following the notations in my research notes (2022.4.23).

source

SCFT models

Polyorder.AbstractSCFTType
AbstractSCFT

An abstract type for SCFT models.

All concrete subtypes of AbstractSCFT should contain following fields:

  • system::PolymerSystem: polymer system
  • fields::AuxiliaryField: potential fields
  • ϕfields::DensityField: density fields
  • updater::SCFTAlgorithm: algorithm for solving SCFT equations.
  • solvers: algorithm lists for solving MDE equations.
source
Polyorder.SCFTABType
SCFTAB{T, AT<:SCFTAlgorithm, MT<:MDEAlgorithm} <: AbstractSCFT

The SCFT model for AB diblock copolymers.

source
Polyorder.NoncyclicChainSCFTMethod
NoncyclicChainSCFT(system, w_or_lattice; ds=0.01, mde=OSF, mdesolvers=nothing, kwargs...)
NoncyclicChainSCFT(system, w, mdesolvers; updater=SD(0.2), kwargs...)
NoncyclicChainSCFT(system, lattice, mdesolvers; updater=SD(0.2), spacing=0.15, pow2=false, kwargs...)
NoncyclicChainSCFT(system, w, ds; mde=OSF, updater=SD(0.2), kwargs...)
NoncyclicChainSCFT(system, lattice, ds; mde=OSF, updater=SD(0.2), spacing=0.15, pow2=false, kwargs...)

Create a NoncyclicChainSCFT instance for a polymer system system in a simulation cell w or lattice.

Arguments

  • system::PolymerSystem: a Polymer.PolymerSystem object which defines the polymer system.
  • w::AuxiliaryField: specify the simulation cell.
  • lattice::BravaisLattice: specify the BravaisLattice of the simulation cell.
  • w_or_lattice: either w or lattice.
  • mdesolvers::Union{T, AbstractDict{Symbol, T}, AbstractDict{Symbol, S}} where T <: MDESolver, S <: AbstractDict{Symbol, <:MDESolver}: the MDE solvers. It can be a single MDESolver which will be applied to all blocks for all components in the system. It can also be a Dict whose element can be a single MDESolver for each component, or a Dict whose element is MDESolver for each block.
  • ds::Union{T, AbstractDict{Symbol, T}, AbstractDict{Symbol, S}} where T <: Real, S <: AbstractDict{Symbol, Real}: contour steps. The way to specify the contour steps is similar to mdesolvers.
  • mde::Union{T, AbstractDict{Symbol, T}, AbstractDict{Symbol, S}} where T <: Type{<:MDEAlgorithm}, S <: AbstractDict{Symbol, Type{<:MDEAlgorithm}}: the MDE algorithm types. The way to specify the MDE algorithm types is similar to mdesolvers.

Keyword Arguments

  • updater::SCFTAlgorithm=SD(0.2): the SCFT updater.
  • spacing=0.15: the maximum spacing (maxΔx) for the discretization of the simulation cell in real space.
  • pow2=false: whether the discretization of the simulation cell in real space is a power of 2.
  • fieldmodeltype=SimpleFieldModel(): the field model type.
  • compress=Incompressible(): whether the SCFT model is compressible.
  • symmetrize=false: whether to symmetrize the fields according to the space group of the simulation cell.
source
Polyorder.create_auxiliaryfieldsFunction
create_auxiliaryfields(species, w; fieldmodeltype, compress)

Return a list of AuxiliaryField objects corresponding to the list of input species. The grid and unit cell are given by the input w.

source
Polyorder.create_densityfieldsFunction
create_densityfields(species::AbstractVector{Symbol}, ϕ::AbstractField)

Return a list of DensityField objects corresponding to the list of input species. The grid and unit cell are given by the input ϕ.

source
Polyorder.create_forcesFunction
create_forces(wfields)

Return a list of AuxiliaryField objects corresponding to the list of input wfields to hold the forces. This function is used in the NoncyclicChainSCFT constructor.

source
Polyorder.create_mdesolversFunction
create_mdesolvers(system, mde, ds; kwargs...)
create_mdesolvers(system; mde=OSF, ds=0.01, kwargs...)

Construct a MDESolver object for each block of each component in system and return it as a Dict{Symbol, Dict{Symbol, Dict{Symbol, MDESolver}} object. The SmallMolecule object is correctly handled. The return value can be passed to the NoncyclicChainSCFT constructor.

source
Polyorder.create_unique_propagatorsFunction
create_unique_propagators(system, grid, mdesolvers)

Create unique propagators for each component in system and return it as a Vector{Dict{Pair, AbstractPropagator}} object. Note the order of returned list is important, which should correspond to the order of components in system.

Arguments

  • system: a PolymerSystem object.
  • grid: a AbstractField object to specify the simulation cell.
  • mdesolvers can be generate via create_mdesolvers or constructed manually. It should be of type AbstractDict{Symbol, AbstractDict{Symbol, <:MDESolver}}. Make sure each block of each component in system has a corresponding MDESolver object assigned.
source
Polyorder.create_MDE_solversFunction
create_MDE_solvers(mdesolvers, system, block2propagator_list, w)

Create MDE solvers for each component in system and return it as a Vector{Dict{Pair, MDEAlgorithm}} object. Note the order of returned list is important, which should correspond to the order of components in system.

Arguments

source
Polyorder.mdesolvertypeFunction
mdesolvertype(scft::NoncyclicChainSCFT)

Return the types of MDEAlogithm objects in scft as a Vector{Dict{Symbol, Type{<:MDEAlgorithm}}. Note the parameters of the type of each MDEAlgorithm is stripped off.

source
Polyorder.initialize!Function
initialize!(scft::AbstractSCFT, seed)

Initialize the SCFT model with specific seed. The seed may be a list of AxuliaryField objects which are potential fields of an ordered phase. Better initialization leads to correct and fast convergence.

source
Polyorder.compute_density!Function
compute_density!(ϕ::DensityField, q::Propagator, qc::Propagator, qqc::Propagator)

Compute the density of a density field ϕ from the propagator q, qc and qqc.

source

State of the SCFT model

Polyorder.residualFunction
residual(AbstractSCFT; norm1=vecnormInf, relative=false, norm2=vecnormInf)

Return the residual error of an SCFT simulation. norm1 is for evaluation of the residual error for each force. If relative is true, we compute the relative residual error. norm2 is for evalution of the statistics of a list of residual errors from all forces.

Due to performance considerations, we here access the force list of an SCFT simulation using field of a type but not the method. To make sure this work, the AbstractSCFT object must have a forces field which is a list of forces of the SCFT model.

Keyword Arguments

  • norm1=vecnormInf: operate on each force. Examples: x->norm(x, Inf), x->norm(x, 1)
  • norm2=vecnormInf: operate on a vector of forces. This operator is usually an aggregation operation. Examples: mean, maximum.
  • relative=false: whether to compute the relative residual error.
source
residual(scft::AbstractSCFT, config::Config)

Return the residual error of an SCFT simulation. config is a Config object which contains the information of the norm1, norm2 and the relative flag.

source
Polyorder.FFunction
F(system::PolymerSystem)
F(model::AbstractSCFT)

Return the total free energy of a polymer system or an SCFT model. Ideal gas free energy is included.

See also F_ig.

source
Polyorder.F_igFunction
F_ig(system::PolymerSystem)
F_ig(model::AbstractSCFT)

Return the ideal gas free energy of a polymer system or an SCFT model.

See also F.

source
Polyorder.F_FHFunction
F_FH(system::PolymerSystem)
F_FH(scft::AbstractSCFT)
F_DIS(system::PolymerSystem)
F_DIS(scft::AbstractSCFT)

Compute the Flory-Huggins free energy of a polymer system in the SCFT model, which is only applicable to the homogeneous phase (the disordered/DIS phase)

source
Polyorder.enthalpyFunction
enthalpy(model::AbstractSCFT)::Real

Compute the enthalpy of a polymer system from an SCFT model.

source
Polyorder.enthalpy_FHFunction
enthalpy_FH(chainscft::NoncyclicChainSCFT)

The enthalpy according to the Flory-Huggins theory: Σ{ij} (χ{ij}N ϕi ϕj).

source
Polyorder.entropyFunction
entropy(model::AbstractSCFT, [i])::Real

Compute the entropy of a polymer system from an SCFT model. If i is provided, the entropy of the i-th component is computed.

source
Polyorder.entropy_igFunction

The ideal gas entropy contribution of the i-th component.

References

  • R.B. Thompson; M.W. Matsen, J. Chem. Phys. 2000, 112, 6863-6872.
source
Polyorder.HwFunction
Hw(model::AbstractSCFT)::Real

Compute the interaction part of the free energy of a polymer system from an SCFT model.

source
Polyorder.HsFunction
Hs(model::AbstractSCFT)::Real

Compute the Q part of the free energy of a polymer system from an SCFT model.

source
Polyorder.ϕFunction
ϕ(chainscft::NoncyclicChainSCFT, i, sp)

Compute the density operator (i.e. a density field) ϕsp of the i-th component. For example, in a polymer blend AB/A/B, there are two density operators for the AB block copolymer, ϕcA and ϕcB, and one density operator for the A homopolymer, ϕhA, and one density operator for the B homopolymer, ϕ_hB.

  • ϕ_cA = ϕ(chainscft, 1, :A)
  • ϕ_cB = ϕ(chainscft, 1, :B)
  • ϕ_hA = ϕ(chainscft, 2, :A)
  • ϕ_hB = ϕ(chainscft, 3, :B)
source
Polyorder.QFunction
Q(args...; kwargs...)::Real

Compute normalized single chain partition function of a polymer/small molecule component.

source
Polyorder.QsFunction
Qs(model::AbstractSCFT)::Vector{Real}

Compute normalized single chain partition functions of all components in a polymer system from an SCFT model.

source
Polyorder.FgFunction
Fg(chainscft::NoncyclicChainSCFT)

Compute the grand potential of the polymer system in the chainscft model:

\[F_g = F - \sum_{i=1}^{n_c-1} \phi_i \tilde{\mu}_i\]

Note, $F_g$ is also equal to $\mu_{nc}/\alpha_{nc}$.

source
Polyorder.Fg_FHFunction
Fg_FH(system::PolymerSystem)
Fg_FH(chainscft::NoncyclicChainSCFT)
Fg_DIS(system::PolymerSystem)
Fg_DIS(chainscft::NoncyclicChainSCFT)

Compute the grand potential of the polymer system according to the Flory-Huiggins theory, which is only applicable to the disordered/homogeneous phase of the polymer system.

source
Polyorder.μ̃_FHFunction
μ̃_FH(system::PolymerSystem, i)
μ̃_FH(chainscft::NoncyclicChainSCFT, i)
μ̃_DIS(system::PolymerSystem, i)
μ̃_DIS(chainscft::NoncyclicChainSCFT, i)

Compute the effective chemical potential of the i-th component in a chainscft model according to the Flory Huggins theory, which is only applicable to homogeneous (disordered) phases.

Arguments

  • i: 1 to nc-1, where nc is the number of components in system or chainscft.
source
Polyorder.μ̃s_FHFunction
μ̃s_FH(system::PolymerSystem)
μ̃s_FH(chainscft::NoncyclicChainSCFT)
μ̃s_DIS(system::PolymerSystem)
μ̃s_DIS(chainscft::NoncyclicChainSCFT)

Compute the effective chemical potentials of all components (except the reference component) in a chainscft model according to the Flory-Huggins theory, which is only applicable to homogeneous (disordered) phases.

source
Polyorder.μ̃Function
μ̃(chainscft::NoncyclicChainSCFT, i)

Compute the effective chemical potential of the i-th component in the SCFT model chainscft.

Arguments

  • i: 1 to nc-1, where nc is the number of components in chainscft.
source
Polyorder.μ̃sFunction
μ̃s(chainscft::NoncyclicChainSCFT)

Compute the effective chemical potentials of all components (except the reference component) in the SCFT model chainscft.

source
Polyorder.μ_nc_FHFunction
μ_nc_FH(system::PolymerSystem)
μ_nc_FH(chainscft::NoncyclicChainSCFT)
μ_nc_DIS(system::PolymerSystem)
μ_nc_DIS(chainscft::NoncyclicChainSCFT)

Compute the chemical potential of the reference component (the last component) in the chainscft model according to the Flory-Huggins theory, which is only applicable to homogeneous (disordered) phases.

source
Polyorder.μ_FHFunction
μ_FH(system::PolymerSystem, i)
μ_FH(chainscft::NoncyclicChainSCFT, i)
μ_DIS(system::PolymerSystem, i)
μ_DIS(chainscft::NoncyclicChainSCFT, i)

Compute the chemical potential of the i-th component in the SCFT model chainscft according to the Flory-Huggins theory, which is only applicable to homogeneous (disordered) phases.

Arguments

  • i: 1 to nc-1, where nc is the number of components in system or chainscft.
source
Polyorder.μs_FHFunction
μs_FH(system::PolymerSystem)
μs_FH(chainscft::NoncyclicChainSCFT)
μs_DIS(system::PolymerSystem)
μs_DIS(chainscft::NoncyclicChainSCFT)

Compute the chemical potential of all components (except the reference component) in the SCFT model chainscft according to the Flory-Huggins theory, which is only applicable to homogeneous (disordered) phases.

source
Polyorder.μ_ncFunction
μ_nc(chainscft::NoncyclicChainSCFT)

Compute the chemical potential of the reference component (the last component) in the SCFT model chainscft.

source
Polyorder.μFunction
μ(chainscft::NoncyclicChainSCFT, i)

Compute the chemical potential of the i-th component in the SCFT model chainscft.

Arguments

  • i: 1 to nc-1, where nc is the number of components in chainscft.
source
Polyorder.μsFunction
μs(chainscft::NoncyclicChainSCFT)

Compute the chemical potentials of all components (except the reference component) in the SCFT model chainscft.

source
Polyorder.gradient_wrt_cellFunction
gradient_wrt_cell(scft::AbstractSCFT)

Compute the gradient with respect to the cell size and shape, $dH/dh = h\Sigma = \sigma h^{-T}$.

Oblique, Triclinic, Monoclinic, Trigonal: full tensor as a vector of length 9 is returned, the elements are aranged as column vectors. Other crystal system: vector of length numdistinctvoigts is returned.

source
Polyorder.stress_tensorFunction
stress_tensor(scft::AbstractSCFT)

Compute the stress tensor of the input scft instance. The value corresponds to eq. (5.131) in GHF book (2006).

Only distinct elements of the stress tensor are computed according to the symmetry of the crystal system of current simulation lattice.

source
Polyorder.stress_tensor!Function
stress_tensor!(stresses, cbc::Component{<:BlockCopolymer}, block2propagator::Dict{<:Pair, <:PropagatorSmall}, kk, T, Ti)

Compute the stress tensor contributed by the cbc component.

  • The input stresses may carry the stresses of other components. Therefore, we simply added the contribution from current component to it.
  • To compute the contribution of current component, simply pass in stresses with all elements being zero.
source
stress_tensor!(stresses, cbc::Component{<:SmallMolecule}, propagators::AbstractVector{<:Propagator}, kk)

Do nothing because SmallMolecule does not contribute to the stress operator.

source

Simulation Cell

Polyorder.OrthogonalTraitType
OrthogonalTrait

An abstract type to indicate whether the unit cell is orthogonal or not. It has two concrete types: Orthogonal and NonOrthogonal.

source
Polyorder.orthogonalFunction
orthogonal(::{<:CrystalSystem}) -> OrthogonalTrait
orthogonal(w::AuxiliaryField) -> OrthogonalTrait

Determine whether the unit cell is orthogonal or not.

See also OrthogonalTrait.

source
Polyorder.k2Function
k2(w::AuxiliaryField)
k2!(k2, w::AuxiliaryField)

Compute the square of each wave vector associated with the grid of the field w. The wave vector is transformed to the reciprocal Cartesian space and shifted to the first Brillouin zone (FBZ). Both orthogonal and non-orthogonal unit cells are supported.

source
Polyorder.starsFunction
stars(w::AuxiliaryField)

Generate stars from the grid of the field w. The symmetry operations are obtained from the BravaisLattice instance associated with w.

It has been tested for space groups:

  • 1D: 1, 2
  • 2D: 17
  • 3D: 223, 225, 229, 230
source
Polyorder.symmetrize!Function
symmetrize!(scft::AbstractSCFT)

Symmetrize the fields of the SCFT instance scft according to the symmetry operations of the Bravais lattice. The stars is precomuted and stored in the SCFT model.

source
symmetrize!(w::AuxiliaryField)

Symmetrize the field w according to the symmetry operations of the Bravais lattice. The stars is computed on-the-fly.

source
symmetrize!(w::AuxiliaryField, stars)

Symmetrize the field w according to the symmetry operations of the Bravais lattice. The stars is provided by the stars.

source
Polyorder.distinct_voigtsFunction

Line: [e11]

source

Oblique: [e11, e22, e12]

source

Rectangular: [e11, e22, 0]

source

Square, Hexagonal2D, HexRect: [e11, e11, 0]

source

Triclinic: [e11, e22, e33, e23, e13, e12]

source

Monoclinic: [e11, e22, e33, 0, e13, 0]

source

Orthorhombic: [e11, e22, e33, 0, 0, 0]

source

Tetragonal, Trigonal, Hexagonal, HexOrthorhombic: [e11, e11, e33, 0, 0, 0]

source

Cubic: [e11, e11, e11, 0, 0, 0]

source
Polyorder.kk_tensorFunction
kk_tensor(w::AuxiliaryField)
kk_tensor!(kk, w::AuxiliaryField)

Compute the dyad of wave vector k, a tensor whose elements are k_{i}k{j} with i, j = {1, ..., ndim}.

source
Polyorder.distinct_kk_tensorFunction
distinct_kk_tensor(w::AuxiliaryField)
distinct_kk_tensor!(kk, w::AuxiliaryField)

Return the distinct kk tensor for a specific simulation cell.

IMPORTANT: It is assumed that all elements of the input kk are zero.

source
Polyorder.distinct_cell_variablesFunction
distinct_cell_variables(scft::AbstractSCFT)

Return a vector of unit cell variables according to the given crystal system presented in scft. The length of the vector is identical to the gradient_wrt_cell.

source

Solve the SCFT model

Polyorder.solve!Function
solve!(scft::AbstractSCFT, [updater], [config]; kwargs...)

Solve a set of SCFT equations in a fixed unit cell. The input scft is modified in place which will carry the final solution after converging.

Arguments

  • updater: is a SCFTAlgorithm. When it is not provided, scft.updater is used.
  • config=CONFIG: extra configuration for solving the set of SCFT equations.
source
Polyorder.solve!Method
solve!(scft::AbstractSCFT, [updater=scft.updater], [config=CONFIG]; kwargs...)

An example implementation of the interface solve! for SCFT simulations. The function solve! is used to solve the SCFT equations in a fixed unit cell. The input scft is modified in place which will carry the final solution after converging.

Keyword Arguments

  • controls=default_iteration_controls(updater, config): the iteration controls to control the SCFT simulation. See the default_iteration_controls function for a list of default controls. You can use default_iteration_controls for a starting point. Add your own controls to the list or provide a new list of controls as long as these controls follows IterationControl.jl interface.

Convergence Test

IterationControl is used to detect convergence. Several early stopping controls are available:

  • ThresholdObjFun: stop the simulation when the difference of the latest neighboring two F values is smaller than config.scft.tol (1e-8 or lower is recommended), in addtion to the residual norm is smaller than config.scft.max_tol (1e-3 or lower is recommended). This control is recommended for practical usage, since it can significantly reduce the computation time, as F value is much easier to converge than the residual norm.
  • Threshold: stop the simulation when the residual norm is less than config.scft.tol (1e-5 or higher is recommended). Without acceleration, it is very hard to achieve redisual norm down to 1e-6. This control is recommended for benchmark purpose.
  • ThresholdStress: stop the simulation when both the stress norm is less than config.cellopt.gtol and one of the same requirement as ThresholdObjFun and Threshold (which is determined by config.scft.tolmode) are met. This control is supposed to work for the VariableCell updater.
  • SlowProgress: stop the simulation when the residual norm reduces slower than a pre-given criterion. This control is recommended to work with Threshold but not ThresholdObjFun.
  • OscillatoryProgress: stop the simulation when oscillation in residual norm is detected. this control is recommended to work with Threshold but not ThresholdObjFun.
  • NumberSinceBest: stop the simulation when a number of iteration is passed since the best residual norm. Not recommended.
  • Patience: stop the simulation when the number of iterations allowed for increasing residual norm is reached. Not recommended.

Tips and Tricks

  • A convergent solution of the target phase from other system is a good intial guess.
  • If convegent fields from other system are not available, a crude intial auxiliary fields (may be drawn by hand), which assembles the target phase is much better than the fields intialized with random numbers.
  • When even crude intial fields are unavailable and random intialization is used, the updaters SD, SIS, and ETD are more stable and robust to use.

IO

  • A progress meter can be enabled by setting config.io.progress_solve to be true.
source
Polyorder.LoggingControlType
LoggingControl(; dir=".", suffix=".csv", io=open(tempname(dir)*suffix, "w"), nstep=1)

An IterationContorl control for logging the loss and F of the SCFT model during simulation to a file.

source
Polyorder.OscillatoryProgressType
OscillatoryProgress(; k=100, tol=1.0, m=k, stop_message=nothing)

An InterationControl control for detecting oscillating loss during a simulation. If oscillating signal is detected, a terminating signal will be sent.

source
Polyorder.SlowProgressType
SlowProgress(; k=100, tol=√eps(1.0), stop_message=nothing)

An IterationContorl control for detecting a slow progress in reducing lossduring a simulation. If the slow progress is detected, a terminating signal will be sent.

source
Polyorder.ThresholdObjFunType
ThresholdObjFun(; tol_obj=√eps(1.0), tol_loss=1e-4, stop_message=nothing)

An IterationContorl control for terminating a simulation when the tolerance of the objective function value and an acceptable tolerance of the residual have both been reached. If both conditions are matched, a terminating signal will be sent.

source
Polyorder.ThresholdStressType
ThresholdStress(; tol_stress=1e-5, tolmode=FTolMode(), tol_obj=1e-8, tol_loss=1e-4, stop_message=nothing)

An IterationContorl control for terminating a simulation when both following conditions are met:

  • stressnorm is less than a threshold value (`tolstress`)
  • Either
    1. When tolmode == FTolMode(): The difference of the free energies between two neighboring iterations is less than a threshold value (tol_obj). In this case, tol_loss serves as a minimum requirement for residual norm.
    2. When tolmode == ResidualTolMode: the residual norm is less than a threshold value (tol_loss).

This control is designed for variable cell simulations. The SCFT updater must has a field stress_norm which contains a list of stress norms of each iteration.

source
Polyorder.default_iteration_controlsFunction
default_iteration_controls(updater, config)

Create default iteration controls. The controls include:

  • Step: how many steps (config.scft.atom_iter) to perform each iteration.
  • NumberLimit: stop if number of iteration exceeds config.scft.max_iter.
  • TimLimit: stop if time exceeds config.scft.maxtime.
  • InvalidValue: stop if loss is NaN or Inf.
  • number_control: display number of iteration.
  • display_F: display free energy.
  • display_loss: display loss.

When updater is VariableCell, following controls are added:

  • ThresholdStress: stop if stress norm is lower than config.cellopt.gtol and the condition of ThresholdObjFun (tolmode is FTolMode) or Threshold (tolmode is ResidualTolMode) is satisfied.
  • display_stress: display stress.
  • display_unitcell: display unit cell.

Otherwise, following controls are added:

  • Threshold: if tolmode is ResidualTolMode. Stop if loss is lower than config.scft.tol.
  • ThresholdObjFun: if tolmode is FTolMode. Stop if loss is lower than config.scft.tol and config.scft.maxtol.

Following controls are added when the correponding key in config.scft is true:

  • SlowProgress: stop if the decrease of loss is too slow. Enabled by config.scft.use_slow_control.
  • OscillatoryProgress: stop if the decrease of loss is ossilatory. Enabled by config.scft.use_osc_control.
  • Patience: stop if the loss does not decrease after the number of iteration exceeds config.scft.numpatience. Enabled by config.scft.use_patience_control.
  • NumberSinceBest: stop if the number of iteration since the best loss exceeds config.scft.numbest. Enabled by config.scft.use_nsbest_control.
  • LoggingControl: log the loss and free energy into a CSV file. Enabled by config.scft.use_log_control.
source

SCFT Updaters

Polyorder.SCFTAlgorithmType
SCFTAlgorithm{T} <: AbstractAlgorithm

Abstract type for algorithms to solve the SCFT equations, i.e. relaxing auxiliary fields to the equilibrium state in a fixed unit cell.

Following fields are mandated to implement a concrete subtype:

  • n: a counter for iterations.
  • evals: number of SCFT equations evaluations.
  • Fs: free energy at each iteration.
  • rs: residual norm at each iteration.
source
Polyorder.PicardType
Picard{T} <: SCFTAlgorithm{T}
PicardMann{T} <: SCFTAlgorithm{T}
const SD = PicardMann
PicardIshikawa{T} <: SCFTAlgorithm{T}
PicardKhan{T} <: SCFTAlgorithm{T}
PicardS{T} <: SCFTAlgorithm{T}

Picard iteration and its variants for solving fixed point equations.

  • Picard: the simple Picard iteration, which is a PicardMann iteration with α = 1. Should be constructed explicitly with datatype, such as Picard{Float64}().
  • PicardMann or SD: the Picard-Mann iteration, which is equivalent to steepest descent (SD) with constant relaxation parameter α.
  • PicardIshikawa: the Picard Ishikawa iteration.
  • PicardKhan: the Picard Khan iteration.
  • PicardS: the Pciard S-hybrid iteration.

In solving SCFT equations, Picard, PicardKhan and PicardS will diverge due to the one step involed having α=1. PicardIshikawa can converge but less efficient than PicardMann.

Therefore, we recommend PicardMann for solving SCFT equations.

source
Polyorder.EulerType
Euler{T} <: SCFTAlgorithm
EMPEC{T} <: SCFTAlgorithm{T}

Forward Euler method for solving ODE, which is almost the same to the Picard-Mann algorithm for solving fixed point equations. Euler is almost equivalent to PicardMann (SD) except that Euler updates the η field while PicardMann use q! which computes the η field based on current auxiliary potential fields.

EMPEC is a predictor-corrector Euler method.

References

  • Nguyen, N. C.; Fernandez, P.; Freund, R. M.; Peraire, J. Accelerated Residual Methods for the Iterative Solution of Systems of Equations. /SIAM J. Sci. Comput./ 2018, /40/ (5), A3157–A3179.
source
Polyorder.EMPECType
EMPEC{T} <: SCFTAlgorithm{T}

An SCFT updater based on EMPEC method, which is a predictor-corrector Euler method.

source
Polyorder.SISType
SIS{T} <: SpectralSCFTAlgorithm{T}

Semi-implicit SCFT algorithm inspired by Fredrickson, et al.

SIS has two types. Type 1: treat the linear term in the force to be time-implicit, which is generally more stable. Type 2: treat the whole force term time explicit.

source
Polyorder.SISFType
SISF{T} <: SpectralSCFTAlgorithm{T}

Semi-implicit SCFT algorithm with full matrix of relaxation parameters.

SIS has two types. Type 1: treat the linear term in the force to be time-implicit, which is generally more stable. Type 2: treat the whole force term time explicit.

source
Polyorder.POType
PO{T} <: SpectralSCFTAlgorithm{T}

Petersen and Ottinger predictor-corrector SCFT algorithm, see Ref: Düchs, D.; Delaney, K. T.; Fredrickson, G. H. A Multi-Species Exchange Model for Fully Fluctuating Polymer Field Theory Simulations. J. Chem. Phys. 2014, 141 (17), 174103.

source
Polyorder.ETDType
ETD{T} <: SpectralSCFTAlgorithm{T}

Exponential time difference SCFT algorithm, see Ref: Düchs, D.; Delaney, K. T.; Fredrickson, G. H. A Multi-Species Exchange Model for Fully Fluctuating Polymer Field Theory Simulations. J. Chem. Phys. 2014, 141 (17), 174103.

source
Polyorder.ETDFType
ETDF{T} <: SpectralSCFTAlgorithm{T}

Exponential time difference SCFT algorithm using a full matrix of relaxation parameters.

source
Polyorder.ETDPECType
ETDPEC{T} <: SpectralSCFTAlgorithm{T}

Exponential time difference with predictor-corrector SCFT algorithm, see Ref: Düchs, D.; Delaney, K. T.; Fredrickson, G. H. A Multi-Species Exchange Model for Fully Fluctuating Polymer Field Theory Simulations. J. Chem. Phys. 2014, 141 (17), 174103.

source
Polyorder.NesterovType
Nesterov{T} <: SCFTAlgorithm{T}

Nesterov acceleration.

Warmup by a number of PicardMann iterations is required to generate a good initial value for the Nesterov acceleration (otherwise, it may diverge) when starting from a random guess of w fields. When starting from a good guess (for example the previous computed w fields), one can set warmup to 0 or a very small number such as 3.

The normal Nesterov iterations are oscillated in reducing the residual in solving SCFT equations. It is found that this oscillation can be effectively suppressed with a suitable restart mechanism. Normally, gradient restart works the best. However, sometimes the x restart works the best.

When β is fixed at a constant value, Nesterov acceleration is exactly the stationally one-step Anderson acceleration Anderson_S1.

References

  • Mitchell, D.; Ye, N.; De Sterck, H. Nesterov Acceleration of Alternating Least Squares for Canonical Tensor Decomposition: Momentum Step Size Selection and Restart Mechanisms. Numer. Linear Algebra Appl. 2020, 27 (4). https://doi.org/10.1002/nla.2297.
  • Nguyen, N. C.; Fernandez, P.; Freund, R. M.; Peraire, J. Accelerated Residual Methods for the Iterative Solution of Systems of Equations. SIAM J. Sci. Comput. 2018, 40 (5), A3157–A3179.
  • d’Aspremont, A.; Scieur, D.; Taylor, A. Acceleration Methods. arXiv [math.OC], 2021.
source
Polyorder.ARDMType
ARDM{T} <: SCFTAlgorithm{T}

Accelerated residual descent method. ARDM is a small tweaks of Nesterov acceleration precoditioned by SD.

References

    • Nguyen, N. C.; Fernandez, P.; Freund, R. M.; Peraire, J. Accelerated Residual Methods for the Iterative Solution of Systems of Equations. SIAM J. Sci. Comput. 2018, 40 (5), A3157–A3179.
source
Polyorder.Anderson_S1Type
Anderson_S1{T} <: SCFTAlgorithm{T}

Stationary Anderson acceleration with m = 1 history. See the same references in the doc of NGMRESS1. AndersonS1 is essentially the same as the Nesterov acceleration with constant β.

source
Polyorder.AndersonSDType
AndersonSD{T} <: SCFTAlgorithm{T}

Anderson acceleration AA(m) preconditioned by PicardMann iteration (equivalent to steepest descent, SD).

The implementation closely follows Walkers & Ni and inspired by NLsolve.anderson. Functions qrdelete! and qradd! are adapted from NLsolve.jl.

References

  • Walker, H. F. Anderson Acceleration : Algorithms and Implementations; 2011; Vol. 2280, pp 1–8.
  • Walker, H. F.; Ni, P. Anderson Acceleration for Fixed-Point Iterations. SIAM J. Numer. Anal. 2011, 49 (4), 1715–1735.
source
Polyorder.AndersonType
Anderson{T} <: NestedSCFTAlgorithm{T}

Anderson acceleration precoditioned by SD (Euler, PicardMann), EMPEC, SIS, ETD, PO, or ETDPEC method. See AndersonSD for more information.

source
Polyorder.NGMRES_S1Type
NGMRES_S1{T} <: SCFTAlgorithm{T}

Stationary NGMRES with m = 1 history.

References

  • Sterck, H. D.; He, Y. On the Asymptotic Linear Convergence Speed of Anderson Acceleration, Nesterov Acceleration, and Nonlinear GMRES. SIAM J. Sci. Comput. 2021, 43 (5), S21–S46.
source
Polyorder.NGMRES_SDType
NGMRES_SD

Nonlinear GMRES acceleration NGMRES(m) preconditioned by PicardMann iteration (equivalent to steepest descent, SD). NGMRES is very similar to Anderson acceleration.

The implementation closely follows Walkers & Ni and inspired by NLsolve.anderson. Functions qrdelete! and qradd! are adapted from NLsolve.jl.

References

  • Sterck, H. de. Steepest Descent Preconditioning for Nonlinear GMRES Optimization. Numer. Linear Algebra Appl. 2012, 121 (4), 609–635.
  • Sterck, H. de. A Nonlinear GMRES Optimization Algorithm Form Canonical Tensor Decomposition. SIAM J. Sci. Comput. 2012, 34 (3), 1351–1379.
  • Walker, H. F. Anderson Acceleration : Algorithms and Implementations; 2011; Vol. 2280, pp 1–8.
  • Walker, H. F.; Ni, P. Anderson Acceleration for Fixed-Point Iterations. SIAM J. Numer. Anal. 2011, 49 (4), 1715–1735.
source
Polyorder.NGMRESType
NGMRES <: NestedSCFTAlgorithm{T}

NGMRES acceleration precoditioned by SD (Euler, PicardMann), EMPEC, SIS, ETD, PO, or ETDPEC method. See NGMRES_SD for more information.

source
Polyorder.OACCEL_SDType
OACCEL{T} <: NestedSCFTAlgorithm{T}
OACCEL_SD{T} <: SCFTAlgorithm{T}

Objective acceleration (O-ACCEL) is a slight modification to the NGMRES acceleration.

  • OACCEL: OACCEL preconditioned by an adequate SCFTAlgorithm.
  • OACCEL_SD: OACCEL preconditioned by SD in particular.

References

  • Riseth, A. N. Objective Acceleration for Unconstrained Optimization. Numer. Linear Algebra Appl. 2019, 26 (1), e2216.
source
Polyorder.BBType

BB{T} <: SCFTAlgorithm{T}

Barzilai-Borwein (BB) gradient method.

BB is a two-point step size gradient method whose step size is derived from a two-point approximation to the secant equation underlying quasi-Newton methods. It is proved that BB is R-superlinearly convergent for quadratic case.

There are two types of BB step size:

  1. αk = norm(Δxk)^2 / Δxk^T Δgk)
  2. αk = Δxk^T Δgk / norm(Δg_k)^2

where Δxk = xk - x{k-1}, Δgk = gk - g{k-1}.

Fields

  • α: step size only used for the first update. Default: 02.
  • step_max: maximum step size allowed.
  • type: BB1 or BB2. Default is BB2.
  • n: current iteration.
  • Fs: first cell size at each iteration.
  • rs: stress norm for at iteration.
  • αs: BB α at each iteration.
  • evals: number of evaluations till current iteration.
  • xk_old and gk_old: for internal use.

References

  • Burdakov, O.; Dai, Y.-H.; Huang, N. Stabilized Barzilai-Borwein Method. arXiv [math.OC], 2019.
source
Polyorder.VariableCellType

VariableCell{T, C<:SCFTAlgorithm{T}, W<:SCFTAlgorithm{T}} <: SCFTAlgorithm{T}

Variable cell method which updates fields and cell sizes simultaneously. The VariableCell object can be served to both solve! and cell_solve!. However, cell_solve! is recommended because it handles additional details about cell optimization and saving simulation.

Recommended cell updater: BB method.

Fields

  • algoc: algorithm for updating cell sizes, current available: SD, BB.
  • algow: algorithm for updating fields, all SCFTAlgorithm subtypes except BB.
  • block: number of fields updates per cell update.
  • maxΔx: space resolution of the fields.
  • pow2Nx: Grid size must be 2^n if true.
  • changeNx: cell size leads to grid size change if true.
  • n: current iteration.
  • Fs: free energy for each iteration.
  • rs: residual norm of fields for each iteration.
  • stress_norm: stress norm for each iteration.
  • xs: cell sizes for each iteration.
  • μs: chemical potentials at each iteration.
  • evals: cumulative number of SCFT updates at each iteration.
  • _Nx_changed and _warmup: for internal use.
source
Polyorder.is_spectralFunction
is_spectral(::SCFTAlgorithm) -> Bool

Check if the SCFT algorithm is spectral or not, of type SpectralSCFTAlgorithm.

source
Polyorder.is_nestedFunction
is_nested(::SCFTAlgorithm) -> Bool

Check if the SCFT algorithm is nested or not, of type NestedSCFTAlgorithm.

source
Polyorder.inner_algosFunction
inner_algos(::SCFTAlgorithm) -> NamedTuple

Return the inner algorithms of the SCFT algorithm. By default, it returns an empty named tuple for simple algorithms. For nested algorithms, it returns a named tuple of inner algorithms.

source
Polyorder.optionsFunction
options(::SCFTAlgorithm) -> NamedTuple

Return the options of the SCFT algorithm. No options returns an empty named tuple. Subtypes of SCFTAlgorithm should implement this function if it has options.

source
Polyorder.objgradfun!Function
objgradfun!(scft::AbstractSCFT, x)

Evaluate the objective function and its gradient at x for the SCFT Hamiltonian. x is a list of auxiliary fields with length of ns+1 (for incompressible model) and ns (for compressible model), where ns the number of species in the polymer system. Note that number of species is generally different than the number of components in the polymer system.

source
Polyorder.precondfun!Function
precondfun!(scft::AbstractSCFT, x, g, α)

The steepest descent conditioner for acceleration methods. This is essentially the same as PicardMann iteration or the Euler forward methods.

Note the + sign here, which is different from the literature on Nesterov and NGMRES algorithms, but is the same as Picard type iteration methods. The reason here is that in the polymer SCFT theory, the force (gradient) is equal to q(x) - x, while it is x - q(x) in the numerical literature.

source
Polyorder.update_fft_plans!Function
update_fft_plans!(::SCFTAlgorithm, ::AbstractArray) -> Nothing
update_fft_plans!(::SCFTAlgorithm, ::AbstractField) -> Nothing

Inplace Updating the FFT plans of the SCFT algorithm.

source

MDE solvers

Polyorder.MDESolverType
MDESolver

A thin wrapper to combine MDEAlgorithm and contour step size. It simplifies the serialization and deserialization of the configuration of MDE algorithms.

Fields

  • algo::Type{<:MDEAlgorithm}=OSF: algorithm to solve the MDE equations.
  • ds::Real=0.01: contour step size.
  • Ns::Integer=0: number of steps in the contour. 0 stands for actual Ns computed from ds
  • options::Dict{Symbol, Any}=Dict{Symbol, Any}(): kwargs for the algorithm.
source
Polyorder.MDEAlgorithmType
MDEAlgorithm <: AbstractAlgorithm

Abstract type for algorithms to solve the MDE equations, i.e. computing propagators from potential fields.

source
Polyorder.OSFType
OSF <: MDEAlgorithm

A MDE solver: solving the MDE using the Operator-Splitting pseudo-spectral (OSF) method.

source
Polyorder.OSFCUDAType
OSFCUDA <: MDEAlgorithm

A MDE solver: solving the MDE using the Operator-Splitting pseudo-spectral (OSF) method using CUDA.

Note that CUDA.jl should be installed and loaded before using this solver.

source
Polyorder.RQM4Type
RQM4 <: MDEAlgorithm

A MDE solver: Ranjan-Qin_Morse 4th order pseudo-spectral method.

source
Polyorder.RQM4CUDAType
RQM4CUDA <: MDEAlgorithm

RQM4 using CUDA.

Note that CUDA.jl should be installed and loaded before using this solver.

source
Polyorder.best_contour_discretizationFunction
best_contour_discretization(f, ds; Ns_min=MIN_CONTOUR_STEPS)

Return the number of contour steps and the contour step size. The number of contour steps is computed as Ns = floor((f+eps(1.0))/ds) + 1. If Ns < Ns_min, then Ns = Ns_min. The contour step size is computed as f/(Ns-1).

source

Cell optimization

Polyorder.cell_solve!Function
cell_solve!(scft::AbstractSCFT, [updater=scft.updater], [config=CONFIG])

Perform cell optimization for the SCFT model scft. It will pick up the most suitable algorithm for the cell optimization based on the type of updater. In particular, if updater is a VariableCell object, the VariableCellOpt algorithm will be used. Otherwise, the OptimStressGuidedCellOpt algorithm will be used.

See cell_solve!(opt::VariableCellOpt, scft::AbstractSCFT, config::Config=CONFIG), cell_solve!(opt::OptimStressGuidedCellOpt, scft::AbstractSCFT, config::Config=CONFIG), and cell_solve!(opt::OptimGradientFreeCellOpt, scft::AbstractSCFT, config::Config=CONFIG) for more details.

source
cell_solve!(opt::VariableCellOpt, scft::AbstractSCFT, [updater::VariableCell], [config::Config])

Variable cell method which updates fields and unit cell simutaneously.

NOTE: the priority of the VariableCell alogrithm is updater > opt.algo > scft.updater.

It is highly recommended that a converged SCFT instance (e.g. solved by solve!) is used instead of a random intialized SCFT instance. Otherwise, cell_solve! may become highly unstable and even diverging.

source
cell_solve!(opt::OptimStressGuidedCellOpt, scft::AbstractSCFT, [updater=scft.updater], [config=CONFIG])

Using optimization methods to find the zero-stress unit cell with the aid of the gradient with respect to cell size and shape. This is a compromising approach to obtain the stress-free solution of the SCFT equations.

The method is invoked when updater is not a VariableCell object. Do not call cell_solve! with scft.updater being a VariableCell object and updater == nothing!

Arguments

  • updater: either a SCFTAlgorithm instance or nothing. When it is not provided or it is nothing, scft.updater is used. Otherwise, updater is used.
  • config: if not provided, a default one which is Polyorder.Config(; scft=SCFTConfig(; tolmode=:F, tol=1e-8), cellopt=CellOptConfig(; algo1=Symbol("Optim.GradientDescent"))) is used.

Tips and Tricks

  • scft is better to be a converged solution with correct phase in a fixed unit cell. Then an acceleration method, such as Anderson, can be used to drastically reduce the number of iterations.

IO

  • A progress meter can be enabled by setting config.io.progress_cell to be true.
  • Four files: config.yml, summary.yml, tracecell, and tracesolve are saved to the path config.io.base_dir.
  • By default, two additional files: fields.h5, densities.h5 are saved. It can be omitted by setting config.io.save_w and/or config.io.save_ϕ to be false.
  • The data format for fields and densities can be set by config.io.data_format.

Pros

  • 1st order gradients are utilized to accelerate the convergence.
  • Converge very fast when a good initial guess of the cell shape is available.

Cons

  • Have to wait for fields converge which may or may not take extra cost.
source
cell_solve(opt::OptimGradientFreeCellOpt, scft::AbstractSCFT, [updater=scft.updater], [config=CONFIG])

Optimizes the unit cell sizes using gradient-free optimization methods. The input scft instance is intact. An optimized version of a subtype of AbstractSCFT instance is returned.

Tips and Tricks

  • scft is better to be a converged solution with correct phase in a fixed unit cell. Then an acceleration method, such as Anderson, can be used to drastically reduce the number of iterations.

Cons

  • only orthogonal unit cells are supported.
  • 1D optimization using bracketing methods which requires a good guess of initial cell size.
  • Bracketing methods requires a constant interval which makes the optimization unnecessary longer when a good inital value is available, which is the case during the later stage of optimization.
source

Configurations

Polyorder.ConfigType
Config

The configuration object for the SCFT simulation.

Fields

  • system::PolymerSystemConfig=PolymerSystemConfig(): the configuration of the polymer system. PolymerSystemConfig is defined in Polymer.jl.
  • lattice::LatticeConfig=LatticeConfig(): the configuration of the lattice. LatticeConfig is defined in Scattering.jl.
  • io::IOConfig=IOConfig(): the configuration of the input/output.
  • scft::SCFTConfig=SCFTConfig(): the configuration of the SCFT simulation.
  • mde::Dict{Symbol, MDEConfig}=Dict(): the configuration of the MDE solvers for each component of a polymer system.
  • cellopt::CellOptConfig=CellOptConfig(): the configuration of the cell optimization.

A default Config object can be created by Config() or obtained from the global constant CONFIG.

Serialization

Use to_config and from_config to serialize and deserialize a Config or any other Configurations.OptionType object, such as SCFTConfig.

IO

Use Polymer.load_config and Polymer.save_config to load and save a Config object. These functions are defined in Polymer.jl. Config and its fields are all supported.

Example

# save to a YAML file
Polymer.save_config("config.yaml", Polyorder.Config())
Polymer.save_config("scft.yaml", Polyorder.SCFTConfig())
# load from a YAML file
config = Polymer.load_config("config.yaml", Polyorder.Config)
scftconfig = Polymer.load_config("scft.yaml", Polyorder.SCFTConfig)
source
Polyorder.MDEConfigType
MDEConfig

The configuration object for the MDE solvers in a component of a PolymerSystem in a SCFT simulation.

Fields

  • blockmode::Symbol=Fixf: the mode to distribute lengths of all blocks in a BCP. A full list of available modes can be found by list_blockmodes().
  • solvers::Dict{<:Pair, MDESolverConfig}=Dict(): the configurations of the MDE solvers. The key is the block (specified as a Pair).
source
Polyorder.SCFTConfigType
SCFTConfig

The configuration object for solving SCFT equations.

Fields

  • compress::Bool=false: compressible or incompressible SCFT model.
  • fieldmodeltype::Symbol=:simple: the type of field model. :simple is the only supported type at present.
  • symmetrize::Bool=false: symmetrize the fields or not according to the specified space group of the simulation cell.
  • maxΔx::Float64=0.15: the minimum spatial grid resolution, in unit of Rg.
  • pow2Nx::Bool=false: set the number of spatial grid along one dimension a power of 2?
  • atom_iter::Int=1: the number of iterations for IterationControl.Step.
  • min_iter::Int=100: the minimum number of iterations.
  • max_iter::Int=2000: the maximum number of iterations.
  • skip_iter::Int=100: the number of iterations to skip for IterationControl.
  • max_restart::Int: the maximum number of restarting simulation. Default is 1. Not used in the current implementation.
  • maxtime::Float64=0.5: the maximum time allowed for one simulation (in hour unit).
  • relative::Bool=false: is relative residual error?
  • norm::Symbol=:vecnormInf: norm function for computing each force.
  • norm2::Symbol=:vecnormInf: norm function for computing the errors of forces.
  • tolmode::Symbol=:Residual: the tolerance mode. Available options are :Residual and :F.
  • tol::Float64=1e-5: the target tolerance.
  • maxtol::Float64=1e-3: the maximum allowed target tolerance, should >= tol. Also used by ThresholdObjFun as tol_loss.
  • dangertol::Float64=1.0: when the error > dangertol, the simulation does not converge, meaning that the relaxation parameters are too large.
  • use_log_control::Bool=false: use LoggingControl control or not.
  • use_slow_control::Bool=false: use SlowProgress control or not.
  • use_osc_control::Bool=false: use OscillatoryProgress control or not.
  • use_nsbest_control::Bool=false: use NumberSinceBest control or not.
  • use_patience_control::Bool=false: use Patience control or not.
  • k_slow::Int=100: for SlowProgress control.
  • tol_slow::Float64=√eps(1.0): for SlowProgress control.
  • tolF::Float64=√eps(1.0): for SlowProgress control.
  • k_osc::Int=100: for OscillatoryProgress control.
  • tol_osc::Float64=1.0: for OscillatoryProgress control.
  • m_osc::Int=100: for OscillatoryProgress control.
  • numbest::Int=10: for NumberSinceBest control. number of iterations since best error: skip_iter * numbest.
  • numpatience::Int=10: number of iterations allowed when error increases: skip_iter * numpatience.
source
Polyorder.CellOptConfigType
CellOptConfig

The configuration object for the cell optimization.

Fields

  • type::Symbol=:VariableCellOpt: the type name of a cell optimization method which is AbstractCellOptAlgorithm.
  • algo::CellOptAlgoConfig: the configuration of the cell optimization method which can be used to reconstruct the instance of this method.
  • changeNx::Bool=true: change Nx according to lx and maxΔx?
  • tol_stress::Float64=1e-5: the tolerance for stress convergence.
source
Polyorder.IOConfigType
IOConfig

The configuration object for input/output of an SCFT simulation.

Fields

  • verbosity::Int=1: the verbosity level. Set 0 or negative value to suppress all display.
  • progress_solve::Bool=false: show progress meter for solve!.
  • progress_cell::Bool=false: show progress meter for cell_solve!.
  • base_dir::String=".": the base directory for saving files.
  • summary::String="summary": the file name for the final state of cell_solve!.
  • trace::String="trace": the file name for the trace of solve! and cell_solve!.
  • fields::String="fields": the file name for the fields of the SCFT model.
  • densities::String="densities": the file name for the densities of the SCFT model.
  • config::String="config": the file name for the configuration of the SCFT model.
  • data_format::Symbol=:HDF5: the data format for saving files. Available options are :HDF5 and :MAT.
  • save_config::Bool=true: save the configuration of the SCFT model or not.
  • save_w::Bool=true: save the fields of the SCFT model or not.
  • save_ϕ::Bool=true: save the densities of the SCFT model or not.
  • save_summary::Bool=true: save the final state of cell_solve! or not.
  • save_trace::Bool=true: save the trace of solve! and cell_solve! or not.
  • display_interval::Int=100: the number of SCFT iterations for displaying information.
  • record_interval::Int=1000: the number of SCFT iterations for recording information.
  • save_interval::Int=1000: the number of SCFT iterations for saving information.
  • ndigits_F::Int=6: the number of digits for displaying free energy.
  • ndigits_err::Int=0: the number of digits for displaying errors.
source
Polymer.to_configFunction
to_config(scft::AbstractSCFT, [updater], [cellopt], [config])

Serialize an SCFT model scft into a Config instance.

Arguments

  • updater=scft.updater: if not present, default value is used.
  • cellopt=updater isa VariableCell ? VariableCellOpt(updater) : OptimStressGuidedCellOpt(): if not present, default value is used.
  • config=CONFIG: if not present, default value is used.

It is useful when the SCFT model is solved with Config object provided using solve! and cell_solve!. Then the exact configuration can be saved to a YAML file for later use.

source
Polymer.makeFunction
make(config::CellOptConfig) -> AbstractCellOptAlgorithm
from_config(config::CellOptConfig) -> AbstractCellOptAlgorithm

Construct an AbstractCellOptAlgorithm instance from a CellOptConfig instance. Types not in list_cellopt_algorithms() are not supported and will throw an ArgumentError.

source
make(config::Config) -> NoncyclicChainSCFT
from_config(config::Config) -> NoncyclicChainSCFT

Construct an NoncyclicChainSCFT instance from a Config instance.

source
make(::ObjType{:System}, config)

Make an PolymerSystem object.

Note: confinement is not implmented.

source
Polymer.save_configFunction
save_config(yamlfile, config)

Save a Configurations.OptionType (defined by @option) object to a YAML file.

source
Polymer.load_configFunction
load_config(yamlfile, T=Config; top=nothing)

Load a Configurations.OptionType (defined by @option) object from a YAML file. T is the type of the object to be loaded. top points to the top level of a sub configuration in the YAML file to be treated as an object of type T.

For example, if the YAML file is saved from a Config object. Then we can load the CellOptConfig object from the YAML file using Polymer.load_config(yamlfile, CellOptConfig; top="cellopt").

source

Visualization

Polyorder.densityplotFunction
densityplot(scft::NoncyclicChainSCFT; kwargs...)
densityplot!(ax, scft::NoncyclicChainSCFT; kwargs...)

Plot density field(s) for an NoncyclicChainSCFT model scft.

Keyword arguments

  • specie=:A: The specie to plot. For 1D models, this keyword is ignored, and all species are plotted by default.
  • repeat=1: The number of times to repeat the density field. Only applicable to 3D models.
  • isovalue=0.5: The contour level to plot. Only applicable to 3D models.
  • isorange=0.5: The contour range to plot. Only applicable to 3D models.
  • alpha=0.75: The transparency of the contour. Only applicable to 3D models.
source
Polyorder.traceplotFunction
traceplot(scft::NoncyclicChainSCFT; kwargs...)
traceplot!(ax, scft::NoncyclicChainSCFT; kwargs...)
traceplot(updater::SCFTAlgorithm; kwargs...)
traceplot!(ax, updater::SCFTAlgorithm; kwargs...)

Plot convergence trace for SCFTAlgorithm as updater or scft.updater.

Keyword arguments

  • var=:error: The variable to plot. Can be :F, :error, :stress, :mu, or :xs.
  • indices=(1,): when var is a vector, given the indices to plot. Only applicable var=:mu and var=:xs.
  • xlog=false: whether to use log scale on the x-axis.
  • ylog=true: whether to use log scale on the y-axis.
source

IO

Polyorder.save_fieldsFunction
save_fields(scft::AbstractSCFT, config::Config)

Save the potential fields of the SCFT model scft into a file specified by config.io.fields. The file format is determined by config.io.data_format.

source
Polyorder.read_fieldsFunction
read_fields(config::Config)

Read the potential fields output by an SCFT simulation from the file specified by config.io.fields. The file format is determined by config.io.data_format.

source
Polyorder.save_densitiesFunction
save_densities(scft::AbstractSCFT, config::Config)

Save the density fields of the SCFT model scft into a file specified by config.io.densities. The file format is determined by config.io.data_format.

source
Polyorder.read_densitiesFunction
read_densities(config::Config)

Read the density fields output by an SCFT simulation from the file specified by config.io.densities. The file format is determined by config.io.data_format.

source
Polyorder.save_trace_solveFunction
save_trace_solve(trace, config::Config)

Save the trace of each SCFT iteration into a text .csv file whose path and name is given by config. trace is produced by cell_solve!, which is a vector of tuples, each tuple (thus a row in the file) is (nevals, F, residual, [unitcell], [chemical potential]).

The CSV should have at least 3 at most 5 columns. The first column is the number of SCFT equations evaluations, the second column is the free energy, the third column is the residual, the fourth column is the unit cell size (only for variable cell method), the fifth column is the chemical potential (only for variable cell method with multicomponent polymer system).

source
Polyorder.read_trace_solveFunction
read_trace_solve(trace_file)

Read a "tracesolve.csv" file written by `savetracesolve. Typically, acellsolve!will generate a set ofsolve!output. This function splits all thesesolve!output into separate vectors. Therefore, the returned value is a vector of matrix. The length of vector is the number ofsolve!called bycell_solve!`. Each element of the vector is a n x 3 matrix, the 1 - 3 columns correspond to nevals, F, and residual, respectively.

source
Polyorder.save_trace_cellFunction
save_trace_cell(trace, config::Config)

Save the trace of cell_solve! into a text .csv file whose path and name is given by config. trace is produced by cell_solve!, which is a vector of tuples, each tuple (thus a row in the file) is (nevals_cell, nevals_solve, F, residual, stress_norm, unitcell, stress[, chemical_potential]). For Mono-component polymer system, there is no chemical_potential column.

The CSV file has a header, which should be "nevalscell, nevalssolve, F, residual, stressnorm, unitcell, stress[, chemicalpotential]".

source
Polyorder.read_trace_cellFunction
read_trace_cell(trace_file)

Read a "tracecell.csv" file written by `savetracecell. Each line (row) of the file stores the state of the SCFT object after asolve!called bycellsolve!. Note thatunitcell,stressandchemicalpotentialcolumns are vectors. In particular, thechemicalpotential` column does not appear for Mono-component polymer system.

source

Utilities

Polyorder.DisplayTimerType

This type provides a convenient (and casual) way to measure the elapsed time for any line of code. It is different from Base.@time that it only report time and the return type is a Dates.CompoundPeriod, which can be converted to seconds by Polaris.tosecond. When the elapsed time is represented in Dates.CompoundPeriod, it is good for human reading when printing. Thus it is useful for long time running codes which exceeds minutes even hours. The aim of this type is not to profile or seriously measure the perfomance of codes but provide a easy way to sketch how quick codes are running.

Examples

timer = DisplayTimer()
# some codes
elapsed, _ = timer()
# more codes
elapsed_since_init, elapsed_since_last = timer()
source
Polyorder.best_N_fftFunction
best_N_fft(L; maxΔx=0.15, pow2=false)

Find the best number of points for FFT given the length L. The number of points is chosen such that the spacing between points is less than maxΔx. If pow2 is true, the number of points is a power of 2. Otherwise, use nextfastfft to find the best number of points.

source
Polyorder.timerConstant
const timer = TimerOutput()

TimerOutput object used to store Polyorder timings.

For production run, one can disable the timer by Polyorder.set_timer_enabled!(false) and restart Julia process.

See also clear_timer!, set_timer_enabled!.

Examples

Polyorder.clear_timer!()
# run Polyorder codes
Polyorder.timer  # show the timings
source
Polyorder.set_timer_enabled!Function
set_timer_enabled!(state=true)

Enable/Disable Polyorder.timer by writting to LocalPreferences.toml. Note that you should restart your Julia process to activate the new state.

source