Reference
Documentation for Polyorder.jl's public interface.
Contents
- Contents
- Index
- General interface
- RPA
- SCFT models
- State of the SCFT model
- Simulation Cell
- Solve the SCFT model
- SCFT Updaters
- MDE solvers
- Cell optimization
- Boundary Conditions
- Configurations
- Visualization
- IO
- Utilities
Index
Polyorder.RPAPolyorder.timerPolyorder.ARDMPolyorder.AbstractAlgorithmPolyorder.AbstractAuxiliaryFieldPolyorder.AbstractBlockModePolyorder.AbstractBoundaryConditionPolyorder.AbstractCellOptAlgorithmPolyorder.AbstractDensityFieldPolyorder.AbstractFieldPolyorder.AbstractFieldModelTypePolyorder.AbstractPropagatorPolyorder.AbstractSCFTPolyorder.AbstractTolModePolyorder.AcceptablePolyorder.AndersonPolyorder.Anderson_S1Polyorder.AuxiliaryFieldPolyorder.BBPolyorder.BoundaryConditionPolyorder.CartesianGridPolyorder.CellOptAlgoConfigPolyorder.CellOptConfigPolyorder.CompressibilityPolyorder.CompressiblePolyorder.ConfigPolyorder.ConfinedAuxiliaryFieldPolyorder.ConfinedDensityFieldPolyorder.ConfinedETDRK4Module.ETDRK4CxPolyorder.ConfinedETDRK4Module.ETDRK4FxCyPolyorder.ConfinedETDRK4Module.ETDRK4PolarPolyorder.ConvergenceStatusPolyorder.CxCyGridPolyorder.CxGridPolyorder.CxGridMDEAlgorithmPolyorder.DataFormatTraitPolyorder.DensityFieldPolyorder.DirichletPolyorder.DisplayTimerPolyorder.DivergePolyorder.EMPECPolyorder.ETDPolyorder.ETDFPolyorder.ETDPECPolyorder.EulerPolyorder.ExchangeFieldModelPolyorder.FTolModePolyorder.FailedPolyorder.FixdsBlockModePolyorder.FixfBlockModePolyorder.FourierFlowsAlgorithmModule.ETDRK4Polyorder.FxCyGridPolyorder.FxCyGridMDEAlgorithmPolyorder.FθCrGridPolyorder.FθCrGridMDEAlgorithmPolyorder.GradientFreeCellOptAlgorithmPolyorder.HDF5FormatPolyorder.IOConfigPolyorder.IncompressiblePolyorder.IntegrationAlgorithmPolyorder.LoggingControlPolyorder.MATFormatPolyorder.MDEAlgorithmPolyorder.MDEConfigPolyorder.MDESmallPolyorder.MDESolverPolyorder.MDESolverConfigPolyorder.NGMRESPolyorder.NGMRES_S1Polyorder.NeedMoreIterationPolyorder.NeedMoreTimePolyorder.NestedSCFTAlgorithmPolyorder.NesterovPolyorder.NeumannPolyorder.NoncyclicChainSCFTPolyorder.NoncyclicChainSCFTPolyorder.OACCELPolyorder.OSFPolyorder.OpenInt4Polyorder.OptimGradientFreeCellOptPolyorder.OptimStressGuidedCellOptPolyorder.OrthogonalTraitPolyorder.OscillatoryPolyorder.OscillatoryProgressPolyorder.POPolyorder.PeriodicPolyorder.PicardPolyorder.PicardIshikawaPolyorder.PicardKhanPolyorder.PicardMannPolyorder.PicardSPolyorder.PolarGridPolyorder.PropagatorPolyorder.PropagatorSmallPolyorder.RQM4Polyorder.ResidualTolModePolyorder.RobinPolyorder.RombergExternalPolyorder.RunningPolyorder.SCFTABPolyorder.SCFTAlgorithmPolyorder.SCFTConfigPolyorder.SCFTModelPolyorder.SCFTSummaryPolyorder.SCFTUpdaterConfigPolyorder.SISPolyorder.SISFPolyorder.SimpleFieldModelPolyorder.SimpsonPolyorder.Simpson1Polyorder.Simpson2Polyorder.SlowProgressPolyorder.SpectralSCFTAlgorithmPolyorder.StressGuidedCellOptAlgorithmPolyorder.StressTensorHelperPolyorder.SuccessfulPolyorder.TerminatedPolyorder.ThresholdObjFunPolyorder.ThresholdStressPolyorder.TwoSpecieFieldModelPolyorder.VariableCellPolyorder.VariableCellOptPolyorder.VersionConfigScattering.AbstractGridPolymer.from_configPolymer.isconfinedPolymer.load_configPolymer.makePolymer.save_configPolymer.to_configPolyorder.FPolyorder.F_DISPolyorder.F_FHPolyorder.F_igPolyorder.F_nigPolyorder.FgPolyorder.Fg_DISPolyorder.Fg_FHPolyorder.HsPolyorder.HwPolyorder.QPolyorder.QsPolyorder.RPA.compute_SmatPolyorder.RPA.compute_Smat_interactingPolyorder.RPA.compute_Smat_non_interactingPolyorder.RPA.compute_stability_limitPolyorder.RPA.form_factorPolyorder.RPA.reciprocal_structure_factorPolyorder.RPA.structure_factorPolyorder._kPolyorder._shift_to_FBZPolyorder.best_N_fftPolyorder.best_contour_discretizationPolyorder.block_NsPolyorder.block_dsPolyorder.block_fPolyorder.cell_solve!Polyorder.cheb_quadrature_clencurtPolyorder.clear_timer!Polyorder.clonePolyorder.compute_density!Polyorder.coordinatesPolyorder.create_MDE_solversPolyorder.create_auxiliaryfieldsPolyorder.create_densityfieldsPolyorder.create_forcesPolyorder.create_mdesolversPolyorder.create_unique_propagatorsPolyorder.create_updaterPolyorder.defaultPolyorder.default_iteration_controlsPolyorder.densityplotPolyorder.dimensionPolyorder.distinct_cell_variablesPolyorder.distinct_kk_tensorPolyorder.distinct_voigtsPolyorder.enthalpyPolyorder.enthalpy_FHPolyorder.entropyPolyorder.entropy_igPolyorder.find_fieldPolyorder.find_field_idxPolyorder.find_propagatorPolyorder.find_propagator_idxPolyorder.find_propagator_indicesPolyorder.find_propagator_pair_indicesPolyorder.force!Polyorder.getβPolyorder.goodPolyorder.gradient_wrt_cellPolyorder.initialize!Polyorder.inner_algosPolyorder.is_nestedPolyorder.is_spectralPolyorder.iscompressiblePolyorder.isincompressiblePolyorder.k2Polyorder.kk_tensorPolyorder.latticePolyorder.list_blockmodesPolyorder.list_cellopt_algorithmsPolyorder.list_dataformatsPolyorder.list_fieldmodelsPolyorder.list_mde_algorithmsPolyorder.list_normsPolyorder.list_scft_algorithmsPolyorder.list_tolmodesPolyorder.mdesolvertypePolyorder.nextfastfftPolyorder.num_distinct_voigtsPolyorder.num_free_cellsizePolyorder.objgradfun!Polyorder.optionsPolyorder.orthogonalPolyorder.precondfun!Polyorder.q!Polyorder.read_densitiesPolyorder.read_fieldsPolyorder.read_trace_cellPolyorder.read_trace_solvePolyorder.record!Polyorder.resetPolyorder.reset!Polyorder.reset!Polyorder.residualPolyorder.save_densitiesPolyorder.save_fieldsPolyorder.save_trace_cellPolyorder.save_trace_solvePolyorder.select_blockmodePolyorder.select_cellopt_algorithmPolyorder.select_dataformatPolyorder.select_fieldmodelPolyorder.select_mde_algorithmPolyorder.select_normPolyorder.select_scft_algorithmPolyorder.select_tolmodePolyorder.set_timer_enabled!Polyorder.solve!Polyorder.solve!Polyorder.starsPolyorder.statusPolyorder.step!Polyorder.stress_tensorPolyorder.stress_tensor!Polyorder.strip_type_paramPolyorder.symmetrize!Polyorder.tensor_to_voigtPolyorder.traceplotPolyorder.update!Polyorder.update_density!Polyorder.update_fft_plans!Polyorder.update_field!Polyorder.update_force!Polyorder.update_propagator!Polyorder.voigt_to_tensorPolyorder.γPolyorder.γ_DISPolyorder.γ_FHPolyorder.γ_igPolyorder.γsPolyorder.γs_DISPolyorder.γs_FHPolyorder.γs_igPolyorder.μPolyorder.μ_DISPolyorder.μ_FHPolyorder.μ_igPolyorder.μsPolyorder.μs_DISPolyorder.μs_FHPolyorder.μs_igPolyorder.μ̃Polyorder.μ̃_DISPolyorder.μ̃_FHPolyorder.μ̃_igPolyorder.μ̃sPolyorder.μ̃s_DISPolyorder.μ̃s_FHPolyorder.μ̃s_igPolyorder.ϕScattering.Crystal.unitcell
General interface
Scattering.AbstractGrid — Type
abstract type AbstractGrid endA grid (lattice) provide a discrete representation of a continuous space. This is the top level abstract type for all kinds of grids (lattices).
Polyorder.CartesianGrid — Type
abstract type CartesianGrid <: AbstractGrid endAbstract type for Cartesian grids.
Polyorder.PolarGrid — Type
abstract type PolarGrid <: AbstractGrid endAbstract type for polar grids.
Polyorder.CxGrid — Type
CxGrid{T, LB, RB} <: CartesianGridA 1D Cartesian grid with Chebyshev nodes.
Fields
edges::SVector{1, T}: The length of the grid.Nx::Int: Number of Chebyshev grid points isNx + 1.x::Vector{T}: The grid points.lbc::LB: Left boundary condition.rbc::RB: Right boundary condition.
Polyorder.FxCyGrid — Type
FxCyGrid{T, LB, RB} <: CartesianGridA 2D Cartesian grid with Fourier grid in the x-direction and Chebyshev grid in the y-direction.
Fields
edges::SVector{2, T}: The lengths of the grid in x and y directions.Nx::Int: Number of Fourier grid points in x direction.Ny::Int: Number of Chebyshev grid points in y direction isNy + 1.x::Vector{T}: Fourier grid points in x direction.y::Vector{T}: Chebyshev grid points in y direction.lbc::LB: Left boundary condition for the y-direction.rbc::RB: Right boundary condition for the y-direction.
Polyorder.FθCrGrid — Type
FθCrGrid{T, LB, RB} <: PolarGridA 2D polar grid with Fourier grid in the θ-direction and Chebyshev grid in the r-direction.
Fields
radius::SVector{1, T}: The radius of the grid.Nθ::Int: Number of Fourier grid points in θ direction.Nr::Int: Number of Chebyshev grid points in r direction.θ::Vector{T}: Fourier grid points in θ direction.r::Vector{T}: Chebyshev grid points in r direction.lbc::LB: Left boundary condition for the r-direction.rbc::RB: Right boundary condition for the r-direction.
Polyorder.CxCyGrid — Type
CxCyGrid{T, LBx, RBx, LBy, RBy} <: CartesianGridA 2D Cartesian grid with Chebyshev nodes in both x and y directions.
Fields
edges::SVector{2, T}: The lengths of the grid in x and y directions.Nx::Int: Number of Chebyshev grid points in x direction isNx + 1.Ny::Int: Number of Chebyshev grid points in y direction isNy + 1.x::Vector{T}: Chebyshev grid points in x direction.y::Vector{T}: Chebyshev grid points in y direction.lbcx::LBx: Left boundary condition for the x-direction.rbcx::RBx: Right boundary condition for the x-direction.lbcy::LBy: Left boundary condition for the y-direction.rbcy::RBy: Right boundary condition for the y-direction.
Polyorder.AbstractField — Type
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.
GPU Support
AbstractField is device-agnostic. The data field can be backed by a CuArray (from CUDA.jl) or any other AbstractArray compatible with KernelAbstractions. When running on GPU, broadcasting and linear algebra operations are automatically dispatched to the appropriate GPU kernels.
Mandate Fields
data:S: the actual data (e.g.ArrayorCuArray).lattice::BravaisLattice: the Bravais lattice of the simulation cell.specie: the polymer specie.
Polyorder.AbstractAuxiliaryField — Type
AbstractAuxiliaryField{T, N, S, P} <: AbstractField{T, N, S, P}Abstract type for auxiliary fields.
Polyorder.AbstractDensityField — Type
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 3DAbstractArray.lattice::BravaisLattice{N, P}: the Bravais lattice of the field.specie::Symbol: the polymer specie of the field.
Polyorder.AuxiliaryField — Type
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.
GPU Support
When data is a CuArray, AuxiliaryField automatically creates and stores CUDA-compatible FFT plans (CuFFT). This enables seamless spectral operations on the GPU.
Fields
data::S: the potential field data, which should be a 1D, 2D or 3DAbstractArray(CPU or GPU).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.
Polyorder.ConfinedAuxiliaryField — Type
ConfinedAuxiliaryField{T, N, S, P, PF, PB, G} <: AbstractAuxiliaryField{T, N, S, P}An auxiliary field for a confined system.
Fields
data::S: The potential field data.grid::G: The grid for the confined system.specie::Symbol: The polymer specie of the field.fft::PF: The forward FFT plan.ifft::PB: The backward FFT plan.
Polyorder.DensityField — Type
DensityField{T, N, S, P} <: AbstractDensityField{T, N, S, P}A field representing the density of a specie.
Fields
data::S: the density data, which should be a 1D, 2D or 3DAbstractArray.lattice::BravaisLattice{N, P}: the Bravais lattice of the field.specie::Symbol: the polymer specie of the field.
Polyorder.ConfinedDensityField — Type
ConfinedDensityField{T, N, S, P, G} <: AbstractDensityField{T, N, S, P}A density field for a confined system.
Fields
data::S: The density data.grid::G: The grid for the confined system.specie::Symbol: The polymer specie of the field.
Polyorder.AbstractPropagator — Type
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::AbstractStorage{S}: The data of the propagator.α::P: The relative length of the whole molecule.
Polyorder.Propagator — Type
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 aBlockCopolymerGraphobject.
Polyorder.PropagatorSmall — Type
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)=1is 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.
Polyorder.AbstractAlgorithm — Type
AbstractAlgorithmAn abstract type for algorithms.
Polyorder.IntegrationAlgorithm — Type
IntegrationAlgorithm <: AbstractAlgorithmAbstract type for integration algorithms.
Polyorder.RombergExternal — Type
RombergExternal <: IntegrationAlgorithmRomberg integration algorithm from external package.
Polyorder.Simpson — Type
Simpson <: IntegrationAlgorithmSimpson integration algorithm.
Polyorder.Simpson1 — Type
Simpson1 <: IntegrationAlgorithmSimpson integration algorithm of first kind.
Polyorder.Simpson2 — Type
Simpson2 <: IntegrationAlgorithmSimpson integration algorithm of second kind.
Polyorder.OpenInt4 — Type
OpenInt4 <: IntegrationAlgorithmAn 4th order integration algorithm for open interval.
Polyorder.AbstractFieldModelType — Type
AbstractFieldModelTypeAbstract type for the type of a field model.
Polyorder.SimpleFieldModel — Type
SimpleFieldModel <: AbstractFieldModelTypeOne speice corresponds to one auxiliary field. For Incompressible system, one additional field to ensure incompressibility.
Polyorder.ExchangeFieldModel — Type
ExchangeFieldModel <: AbstractFieldModelTypeMulti-specie exchange model proposed by Delaney-Fredrickson.
Polyorder.TwoSpecieFieldModel — Type
TwoSpecieFieldModel <: AbstractFieldModelTypeA special exchange model by Fredrickson, see his 2006 Book.
Polyorder.Compressibility — Type
CompressibilityAbstract type for compressibility of a field model.
Polyorder.Compressible — Type
Compressible <: CompressibilityTrait for a compressible field model.
Polyorder.Incompressible — Type
Incompressible <: CompressibilityTrait for an incompressible field model.
Polyorder.AbstractTolMode — Type
AbstractTolModeAbstract type of tolerance mode for testing convergence.
Polyorder.ResidualTolMode — Type
ResidualTolMode <: AbstractTolModeUse residual as tolerance for SCFT convergence.
Polyorder.FTolMode — Type
FTolMode <: AbstractTolModeUse difference of free energy as tolerance for SCFT convergence.
Polyorder.AbstractBlockMode — Type
AbstractBlockModeAbstract type for block mode.
Polyorder.FixfBlockMode — Type
FixfBlockMode <: AbstractBlockModeKeep f of all blocks unchanged. Varying ds or Ns accordingly.
Polyorder.FixdsBlockMode — Type
FixdsBlockMode <: AbstractBlockModeKeep ds of all blocks unchanged. Varying f or Ns accordingly.
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.
Polyorder.update! — Function
update!(object, args...; kwargs...)
update!(object::AbstractSCFT, updater::SCFTAlgorithm, config::Config)In place updating of an object.
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 aSCFTAlgorithm. When it is not provided,scft.updateris used.config=CONFIG: extra configuration for solving the set of SCFT equations.
Polyorder.dimension — Function
dimension(container)::Polymer.SpaceDimensionReturn the dimension of a container.
Polyorder.lattice — Function
lattice(container):AbstractGridReturn the underlying discrete representation of the container. container can be AbstractField, AbstractSCFT.
Scattering.Crystal.unitcell — Function
unitcell(grid::AbstractField)Return the unitcell of the field.
Polyorder.reset — Function
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!.
Polyorder.reset! — Function
Polyorder.clone — Function
clone(obj::T, args...; kwargs...)Return a new object of type T based on obj with updates according to args and kwargs.
Polyorder.default — Function
default(::Type{T}) where T
default(::T) where TReturn a default instance of type T.
Polymer.isconfined — Function
isconfined(scft::AbstractSCFT)::BoolReturn true if the SCFT system is confined.
Polyorder.coordinates — Function
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: anAbstractFieldobject.space::CrystalSpace=RCSpace(): the space of the returned coordinates. SeeScattering.CrystalSpacefor all supported spaces.
Polyorder.num_free_cellsize — Function
num_free_cellsize(container)::IntReturn the number of free cell dimensions of the container that can be treated as free variables for optimization.
Polyorder.iscompressible — Function
iscompressible(c::Compressibility)Check if a field model is compressible.
Polyorder.isincompressible — Function
isincompressible(c::Compressibility)Check if a field model is incompressible.
Polyorder.list_fieldmodels — Function
list_fieldmodels()List available field model types.
Polyorder.select_fieldmodel — Function
select_fieldmodel(n::Symbol)Select a field model type by name.
Polyorder.list_tolmodes — Function
list_tolmodes()List available tolerance modes.
Polyorder.select_tolmode — Function
select_tolmode(n::Symbol)Select a tolerance mode by name.
Polyorder.list_blockmodes — Function
list_blockmodes()List available block modes.
Polyorder.select_blockmode — Function
select_blockmode(n::Symbol)Select a block mode by name.
Polyorder.find_field_idx — Function
find_field_idx(sp::Symbol, fields::AbstractVector{T}) where T<:AbstractFieldReturn the index of the first field with specie sp in fields.
Polyorder.find_field — Function
find_field(sp::Symbol, fields::AbstractVector{T}) where T<:AbstractFieldReturn the first field with specie sp in fields.
Polyorder.find_propagator_idx — Function
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.
Polyorder.find_propagator — Function
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.
Polyorder.find_propagator_indices — Function
find_propagator_indices(specie, propagators)Return a list of indices of propagators whose corresponding specie is the same as the input sp.
Polyorder.find_propagator_pair_indices — Function
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.
Polyorder.block_ds — Function
block_ds(q::Propagator)Return the contour step of the propagator q.
block_ds(scft::NoncyclicChainSCFT)Return the contour steps for each block of each component in scft as a Vector{Dict{Pair, Float64}}.
Polyorder.block_Ns — Function
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.
Polyorder.block_f — Function
block_f(q::Propagator)Return the contour length of the propagator q.
Polyorder._shift_to_FBZ — Function
_shift_to_FBZ(w::AuxiliaryField1D, i, [j], [k])Transform the integral representation of a wave vector in a reciprocal Bravais lattice to the wave vector in a reciprocal Cartesian lattice, and shift it to the first Brillouin zone (FBZ).
Note that the index i, j, k is assumed staring from 1 as consistent to the Julia array indexing. Therefore, in the code we have to subtract 1 first.
For the returned wave vector, 2π is multiplied.
Polyorder._k — Function
_k(w::AuxiliaryField1D, i, [j], [k])Transform the wave vector of the reciprocal grid from QBSpace to QCSpace, and shift it to the first Brillouin zone (FBZ).
For 1D cells and 2D & 3D orthogonal unit cells, it is equivalent to transform the grid origin to the center of grid, which is also the old approach we used to discretize the Laplacian operator.
The free energy calculated via _shift and _shift_to_FBZ is identical for both orthogonal and non-orthogonal unit cell. However, Qin & Morse (PhD thesis) suggested that _shift_to_FBZ is better. It conserves the symmetry of the solution.
RPA
Polyorder.RPA — Module
RPARandom Phase Approximation (RPA) module for analyzing polymer phase behavior and computing structure factors, form factors, and stability limits for polymer systems.
This module provides functions for:
- Computing single chain structure factors and form factors
- Analyzing spinodal stability limits
- Calculating RPA structure factors for multi-component polymer systems
- Computing scattering structure factor matrices
Polyorder.RPA.compute_stability_limit — Function
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.
Polyorder.RPA.form_factor — Function
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.
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.
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.
Polyorder.RPA.reciprocal_structure_factor — Function
reciprocal_structure_factor(chain::BlockCopolymer, χN, k2::Real)
reciprocal_structure_factor(chain::BlockCopolymer, χN, k2::AbstractVector)Return the reciprocal structure factor: $N/S$.
See also RPA.structure_factor.
Polyorder.RPA.structure_factor — Function
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.
Polyorder.RPA.compute_Smat — Function
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.
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.
Polyorder.RPA.compute_Smat_non_interacting — Function
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).
Polyorder.RPA.compute_Smat_interacting — Function
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).
SCFT models
Polyorder.AbstractSCFT — Type
AbstractSCFTAn abstract type for SCFT models.
All concrete subtypes of AbstractSCFT should contain following fields:
system::PolymerSystem: polymer systemfields::AuxiliaryField: potential fieldsϕfields::DensityField: density fieldsupdater::SCFTAlgorithm: algorithm for solving SCFT equations.solvers: algorithm lists for solving MDE equations.
Polyorder.NoncyclicChainSCFT — Type
NoncyclicChainSCFT <: AbstractSCFTAn SCFT model for systems consisting of noncyclic polymer chains w/o small molecules.
Performance Profile
The profile field controls memory/speed trade-off via predefined configurations:
| Profile | storage | symmetry | mde_precompute | Use Case |
|---|---|---|---|---|
:fast | :full | none | true | CPU, small/medium scale |
:balanced | :full | auto | true | Default, balanced choice |
:compact | :full | auto | false | GPU, medium scale |
:minimal | :checkpoint | auto | false | GPU, large scale |
Memory Optimization Fields
profile::Symbol: Performance profile (:fast, :balanced, :compact, :minimal)_mde_buffers: Shared work buffers for periodic MDE solvers (OSF, RQM4, ETDRK4). Contains u, k², û, ak arrays shared across all solver instances. Reduces MDE memory by ~69%. Mutable struct allows selective field reallocation during cell optimization with changeNx=true. Set tonothingfor confined solvers (ETDRK4Cx, etc.) which have their own specialized buffers._expand_buffer_a: Shared buffer for various operations including: - Expandingq[i]from SymmetricStorage without allocation - Cumulative integration intermediate results - Free energy calculation workspace_expand_buffer_b: Buffer for expandingqc[i]from SymmetricStorage without allocation. Used in density calculations with SymmetricStorage.
Polyorder.SCFTAB — Type
SCFTAB{T, AT<:SCFTAlgorithm, MT<:MDEAlgorithm} <: AbstractSCFTThe SCFT model for AB diblock copolymers.
Polyorder.SCFTModel — Type
SCFTModel{<:AbstractSCFT}A thin wrapper to combine AbstractSCFT and Config objects.
Polyorder.NoncyclicChainSCFT — Method
NoncyclicChainSCFT(system, w_or_lattice; kwargs...)
NoncyclicChainSCFT(system, w_or_lattice, ds; kwargs...)
NoncyclicChainSCFT(system, w, mdesolvers; kwargs...)
NoncyclicChainSCFT(config::Config)
NoncyclicChainSCFT(chainscft::NoncyclicChainSCFT)Create a NoncyclicChainSCFT instance in various ways.
Arguments
system::PolymerSystem: aPolymer.PolymerSystemobject which defines the polymer system.w_or_lattice: either anAbstractField(w) or anAbstractGrid(lattice) to specify the simulation cell.w: anAbstractFieldto specify the simulation cell.mdesolvers::Union{T, AbstractDict{Symbol, T}, AbstractDict{Symbol, S}}whereT <: MDESolver, S <: AbstractDict{Symbol, <:MDESolver}: the MDE solvers. It can be a singleMDESolverwhich will be applied to all blocks for all components in thesystem. It can also be a Dict whose element can be a singleMDESolverfor each component, or a Dict whose element isMDESolverfor each block.ds::Union{T, AbstractDict{Symbol, T}, AbstractDict{Symbol, S}}whereT <: Real, S <: AbstractDict{Symbol, Real}: contour steps. The way to specify the contour steps is similar tomdesolvers.config::Config: We can reconstruct aNoncyclicChainSCFTobject from aConfigobject.chainscft::NoncyclicChainSCFT: Reconstruct a brand new NoncyclicChainSCFT instance which is completely independent to the originalchainscft. This means all fields are initialized and stored in a new memory location.
Keywords
NoncyclicChainSCFT(system, w, mdesolvers; kwargs...)
init=randn: method for initializing auxiliary fields. Can be a symbol (:zeros,:ones,:rand,:randn), another SCFT object, a vector of arrays, or a filepath to a.mator.h5file.rng=default_rng(): random number generator for generating random initial fields data.updater::SCFTAlgorithm=SD(0.2): the SCFT updater.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.
NoncyclicChainSCFT(system, w_or_lattice; kwargs...)
ds=0.01: same asdsin Arguments.mde=OSF: same asmdesolversin Arguments.mde_options=(;): aNamedTupleof options passed to the MDE solvers when they are constructed. Ignored ifmdesolversis provided.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.kwargs: additional keywords that can be accepted by the above constructor.
NoncyclicChainSCFT(system, w_or_lattice, ds; kwargs...)
kwargs: keywords that can be accepted by the above constructor.
GPU Support
To run the simulation on a GPU, ensure that the w field passed to the constructor is backed by a CuArray (from CUDA.jl). The NoncyclicChainSCFT object will then automatically handle all operations on the device.
Polyorder.create_auxiliaryfields — Function
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.
Polyorder.create_densityfields — Function
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 ϕ.
Polyorder.create_forces — Function
create_forces(wfields)Return a list of AbstractField objects corresponding to the list of input wfields to hold the forces. This function is used in the NoncyclicChainSCFT constructor.
Polyorder.create_mdesolvers — Function
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.
Polyorder.create_unique_propagators — Function
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: aPolymerSystemobject.grid: aAbstractFieldobject to specify the simulation cell.mdesolverscan be generate viacreate_mdesolversor constructed manually. It should be of typeAbstractDict{Symbol, AbstractDict{Symbol, <:MDESolver}}. Make sure each block of each component insystemhas a correspondingMDESolverobject assigned.
Polyorder.create_MDE_solvers — Function
create_MDE_solvers(mdesolvers, system, block2propagator_list, ws, shared)Create MDE solvers for each component in system with shared work buffers. Returns a Vector{Dict{Pair, MDEAlgorithm}} object.
Solvers with same (algo, b, ds) are shared regardless of low_memory setting.
- low_memory=true: Shared solver uses on-the-fly expd/expw computation
- low_memory=false: Shared solver uses precomputed expd and expw buffer
Arguments
mdesolvers: seecreate_unique_propagators.system: aPolymerSystemobject.block2propagator_list: generated bycreate_unique_propagators.ws: aVector{AbstractField}object to specify the simulation cell.shared:MDEWorkBuffersshared by all periodic MDE solvers, ornothingfor confined solvers.low_memory: If true, use on-the-fly expd computation (slower but less memory).
create_MDE_solvers(mdesolvers, system, block2propagator_list, ws)Backwards-compatible version that auto-creates MDEWorkBuffers for periodic fields. For confined fields (fft=nothing), passes nothing as shared buffers.
Polyorder.mdesolvertype — Function
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.
Polyorder.compute_density! — Function
compute_density!(ϕ::AbstractDensityField, q::Propagator, qc::Propagator, buffer::AbstractArray)Compute the density using cumulative integration, avoiding the need for _qqc storage.
This is the memory-optimized version that uses only a single buffer of size grid_size instead of Ns * grid_size required by the qqc-based version.
Arguments
ϕ: Output density fieldq,qc: Forward and backward propagatorsbuffer: Shared buffer for intermediateq[i] * qc[Ns+1-i]products
Memory Usage
- O(gridsize) for buffer, instead of O(Ns * gridsize) for qqc
compute_density!(ϕ::AbstractDensityField, q::Propagator, qc::Propagator,
buf_a::AbstractArray, buf_b::AbstractArray)Compute the density using zero-allocation expansion from SymmetricStorage.
This is the fully memory-optimized version that:
- Uses cumulative integration (like Phase 1)
- Uses
getindex_into!to expand compressed propagator data into pre-allocated buffers - Avoids ALL temporary memory allocations during density calculation
Arguments
ϕ: Output density fieldq,qc: Forward and backward propagators (may use SymmetricStorage)buf_a: Buffer for expandingq[i]buf_b: Buffer for expandingqc[i]
Memory Usage
- O(2 * grid_size) for buffers, no per-iteration allocations
- Critical for SymmetricStorage which allocates on every
getindexcall
compute_density!(ϕ, component, propagators, shared_buffer; graph=nothing)Memory-optimized version using cumulative integration. Uses shared_buffer (single grid-sized array) instead of qqc (Ns grid-sized arrays).
Polyorder.force! — Function
Fixed point equations are x = g(x) Here, x = wA, wB, η, etc. We define a force as f(x) = g(x) - x which is a function of x.
Polyorder.update_density! — Function
update_density!(scft::AbstractSCFT)Update the density fields given propagators by integration.
Polyorder.update_propagator! — Function
update_propagator!(scft::AbstractSCFT)Update the propagators given potential fields by solving the MDE equations.
Polyorder.update_force! — Function
update_force!(scft::AbstractSCFT)Update the force fields given potential fields and density fields.
Polyorder.update_field! — Function
update_field!(scft::AbstractSCFT)Update the potential fields given force fields and density fields.
Polyorder.q! — Function
q!(chainscft::NoncyclicChainSCFT)q is a function in a fixed point equation by treating SCFT equations as a fixed point equation: $x = q(x)$. The force is $f(x) = q(x) - x$. In this function, $q(x)$ is returned. In addition, forces are also computed and stored in the forces vector.
For a general NoncyclicChainSCFT, we have
\[\mathbf{w} = \mathbf{χ} \mathbf{\tilde{\phi}} + \eta\mathbf{1}\]
where $\mathbf{w}=(w_A, w_B, \dots, w_Z)^T$ is a column vector of auxiliary fields, $\mathbf{χ}$ is the χNMatri, $\mathbf{\tilde{\phi}}=(\phi_A, \phi_B, \dots, \phi_z)^T$ is a column vector of the density (operator) fields, $\eta$ is the incompressible field, which is computed as
\[\eta = \frac{W - 1}{X}\]
with
\[W = \sum_i\sum_j χ_{ij}^{-1} w_i\]
and
\[X = \sum_i\sum_j χ_{ij}^{-1}\]
where $χ_{ij}^{-1}$ is the element of the inverse matrix of $\mathbf{χ}$.
The above equation for $\eta$ is derived by solving the density field as
\[\mathbf{\tilde{\phi}} = \mathbf{χ}^{-1} \mathbf{w} + \eta\mathbf{χ}^{-1}\mathbf{1}\]
and ivoking the incompressible condition:
\[\sum_i \tilde{\phi}_i = 1\]
The notation q is from the reference: 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.
State of the SCFT model
Polyorder.residual — Function
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=norm: operate on each force. Examples:x->norm(x, Inf),x->norm(x, 1)norm2=mean: operate on a vector of forces. This operator is usually an aggregation operation. Examples:mean,maximum.relative=true: whether to compute the relative residual error.
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.
Polyorder.F — Function
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.
Polyorder.F_ig — Function
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.
Polyorder.F_nig — Function
F_nig(scft::AbstractSCFT)Compute the non-ideal gas contribution to the free energy.
This represents the portion of the free energy that arises from non-ideal gas behavior, calculated as the sum of the chain stretch energy Hw(scft) and the incompressibility energy Hs(scft).
Arguments
scft::AbstractSCFT: The SCFT model
Returns
The non-ideal gas free energy contribution.
Polyorder.F_FH — Function
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)
Polyorder.F_DIS — Function
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)
Polyorder.enthalpy — Function
enthalpy(model::AbstractSCFT)::RealCompute the enthalpy of a polymer system from an SCFT model.
Polyorder.enthalpy_FH — Function
enthalpy_FH(chainscft::NoncyclicChainSCFT)The enthalpy according to the Flory-Huggins theory: Σ{ij} (χ{ij}N ϕi ϕj).
Polyorder.entropy — Function
entropy(model::AbstractSCFT, [i])::RealCompute the entropy of a polymer system from an SCFT model. If i is provided, the entropy of the i-th component is computed.
Polyorder.entropy_ig — Function
The ideal gas entropy contribution of the i-th component.
References
- R.B. Thompson; M.W. Matsen, J. Chem. Phys. 2000, 112, 6863-6872.
Polyorder.Hw — Function
Hw(model::AbstractSCFT)::RealCompute the interaction part of the free energy of a polymer system from an SCFT model.
Polyorder.Hs — Function
Hs(model::AbstractSCFT)::RealCompute the Q part of the free energy of a polymer system from an SCFT model.
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)
Polyorder.Q — Function
Q(args...; kwargs...)::RealCompute normalized single chain partition function of a polymer/small molecule component.
Polyorder.Qs — Function
Qs(model::AbstractSCFT)::Vector{Real}Compute normalized single chain partition functions of all components in a polymer system from an SCFT model.
Polyorder.Fg — Function
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}$.
Polyorder.Fg_FH — Function
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.
Polyorder.Fg_DIS — Function
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.
Polyorder.γ — Function
γ(chainscft::NoncyclicChainSCFT, i)Compute the partial derivative of the free energy with respect to the volume fraction of the i-th component in the SCFT model chainscft.
Polyorder.γs — Function
γs(chainscft::NoncyclicChainSCFT)Compute the partial derivative of the free energy with respect to the volume fraction of all components in the SCFT model chainscft.
Polyorder.γ_ig — Function
The analytical solution for the disordered/homogeneous phase of any polymeric systems at the mean-field level is available, which is the same as the Flory-Huggins lattice theory,
\[\tilde{F} = \frac{F}{CV} = \sum_{X,Y}\chi_{XY}N \bar\phi_X \bar\phi_Y + \sum_i \frac{\phi_i}{\alpha_i}\left( \ln\frac{C\phi_i}{\alpha_i} - 1 \right)\]
where $X, Y = \lbrace A, B, C, \dots \rbrace$ and $Y > X$ are the specie types in the polymer system, and $i = 1, 2, 3, \dots$ are the components in the polymer system. The averaging density for the specie X, $\bar\phi_X$, is
\[\bar\phi_X = \sum_i\sum_j \delta_{XS_{ij}} f_{ij} \phi_i\]
where $S_{ij}$ is the specie type in the $j$-th block of the $i$-th molecule, and $\delta_{XY}$ is the Kronecker delta function. Now we can compute its first-order partial derivative
\[\frac{\partial\bar\phi_X}{\partial\phi_i} = \sum_j \delta_{XS_{ij}} f_{ij}\]
Therefore,
\[\gamma_i = \frac{\partial\tilde{F}}{\partial\phi_i} = \sum_{X,Y}\sum_j ( \delta_{XS_{ij}}\bar\phi_Y + \delta_{YS_{ij}}\bar\phi_X)\chi_{XY}N f_{ij} + \frac{1}{\alpha_i}\ln\frac{C\phi_i}{\alpha_i}\]
and the chemical potential can be readily computed
\[\tilde\mu_i = \gamma_i - \gamma_{n_c}\]
for $i=1, 2, \dots, n_c$ with $n_c$ the number of components in the system.
γ_ig(chainscft::NoncyclicChainSCFT, i)Compute the ideal gas contribution to the partial derivative of the free energy with respect to the volume fraction of the i-th component.
Polyorder.γs_ig — Function
γs_ig(system::PolymerSystem)
γs_ig(chainscft::NoncyclicChainSCFT)Compute the ideal gas contribution to the partial derivative of the free energy with respect to the volume fraction of all components.
Polyorder.γ_FH — Function
γ_FH(system::PolymerSystem, i)
γ_FH(chainscft::NoncyclicChainSCFT, i)
γ_DIS(system::PolymerSystem, i)
γ_DIS(chainscft::NoncyclicChainSCFT, i)Compute the partial derivative of the free energy with respect to the volume fraction of the i-th component in a chainscft model according to the Flory Huggins theory, which is only applicable to homogeneous (disordered) phases.
Polyorder.γs_FH — Function
γs_FH(system::PolymerSystem)
γs_FH(chainscft::NoncyclicChainSCFT)
γs_DIS(system::PolymerSystem)
γs_DIS(chainscft::NoncyclicChainSCFT)Compute the partial derivative of the free energy with respect to the volume fraction of all components in a chainscft model according to the Flory-Huggins theory, which is only applicable to homogeneous (disordered) phases.
Polyorder.γ_DIS — Function
γ_FH(system::PolymerSystem, i)
γ_FH(chainscft::NoncyclicChainSCFT, i)
γ_DIS(system::PolymerSystem, i)
γ_DIS(chainscft::NoncyclicChainSCFT, i)Compute the partial derivative of the free energy with respect to the volume fraction of the i-th component in a chainscft model according to the Flory Huggins theory, which is only applicable to homogeneous (disordered) phases.
Polyorder.γs_DIS — Function
γs_FH(system::PolymerSystem)
γs_FH(chainscft::NoncyclicChainSCFT)
γs_DIS(system::PolymerSystem)
γs_DIS(chainscft::NoncyclicChainSCFT)Compute the partial derivative of the free energy with respect to the volume fraction of all components in a chainscft model according to the Flory-Huggins theory, which is only applicable to homogeneous (disordered) phases.
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 inchainscft.
Polyorder.μ̃s — Function
μ̃s(chainscft::NoncyclicChainSCFT)Compute the effective chemical potentials of all components (except the reference component) in the SCFT model chainscft.
Polyorder.μ̃_ig — Function
μ̃_ig(system::PolymerSystem, i)
μ̃_ig(chainscft::NoncyclicChainSCFT, i)Compute the ideal gas contribution to the effective chemical potential of the i-th component in a chainscft model.
Polyorder.μ̃s_ig — Function
μ̃s_ig(system::PolymerSystem)
μ̃s_ig(chainscft::NoncyclicChainSCFT)Compute the ideal gas contribution to the effective chemical potentials of all components (except the reference component) in a chainscft model.
Polyorder.μ̃_FH — Function
μ̃_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 insystemorchainscft.
Polyorder.μ̃s_FH — Function
μ̃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.
Polyorder.μ̃_DIS — Function
μ̃_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 insystemorchainscft.
Polyorder.μ̃s_DIS — Function
μ̃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.
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 inchainscft.
Polyorder.μs — Function
μs(chainscft::NoncyclicChainSCFT)Compute the chemical potentials of all components (except the reference component) in the SCFT model chainscft.
Polyorder.μ_ig — Function
μ_ig(system::PolymerSystem, i)
μ_ig(chainscft::NoncyclicChainSCFT, i)Compute the ideal gas contribution to the chemical potential of the i-th component in the SCFT model chainscft.
Polyorder.μs_ig — Function
μs_ig(system::PolymerSystem)
μs_ig(chainscft::NoncyclicChainSCFT)Compute the ideal gas contribution to the chemical potentials of all components (except the reference component) in the SCFT model chainscft.
Polyorder.μ_FH — Function
μ_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 insystemorchainscft.
Polyorder.μs_FH — Function
μ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.
Polyorder.μ_DIS — Function
μ_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 insystemorchainscft.
Polyorder.μs_DIS — Function
μ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.
Polyorder.gradient_wrt_cell — Function
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.
Polyorder.stress_tensor — Function
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.
stress_tensor(scft::NoncyclicChainSCFT)Optimized stress tensor computation for NoncyclicChainSCFT that reuses _expand_buffer_a and _expand_buffer_b to avoid temporary allocations when propagators use SymmetricStorage.
Phase 3 memory optimization:
- Reuse expand buffers for getindex_into! (avoid SymmetricStorage allocations)
- Reuse pre-allocated complex FFT buffers from StressTensorHelper (avoid FFT temp allocations)
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
stressesmay 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
stresseswith all elements being zero.
Phase 3 optimized version with pre-allocated buffers. Includes both real expand buffers (bufa, bufb) and complex FFT buffers (qk, qck, qt).
stress_tensor!(stresses, cbc::Component{<:SmallMolecule}, propagators::AbstractVector{<:Propagator}, kk)Do nothing because SmallMolecule does not contribute to the stress operator.
Phase 3 fully optimized version: uses getindex_into! to avoid temporary allocations when propagators use SymmetricStorage, AND uses pre-allocated complex FFT buffers.
Buffers:
buf_a,buf_b: Reused from NoncyclicChainSCFT's_expand_buffer_a/bfor getindex_into!qk,qck: Pre-allocated complex FFT buffers from StressTensorHelper
Algorithm: We use qck as a temporary input buffer for FFT operations.
- qck .= buf_a (real → complex)
- mul!(qk, Ti, qck) gives q(-k)
- qck .= buf_b (real → complex, reusing qck as temp)
- mul!(qck, T, qck) is NOT allowed (in-place mul! fails)
Since in-place mul! doesn't work, we need a workaround:
- Use the fact that we can save qk first, then recompute qck
Simulation Cell
Polyorder.StressTensorHelper — Type
StressTensorHelper{K, P, Pi, C}Helper struct for stress tensor calculation.
Fields
distinct_kk_tensor: The distinct kk tensor for stress calculationT: Forward FFT planTi: Inverse FFT plan_fft_buffer_qk: Pre-allocated complex buffer for qk (inverse FFT result)_fft_buffer_qck: Pre-allocated complex buffer for qck (forward FFT result)_fft_buffer_qt: Pre-allocated complex buffer for real-to-complex conversion
Phase 3 optimization: These buffers are allocated once during SCFT construction and reused across all stress_tensor! calls, avoiding repeated allocations.
Polyorder.OrthogonalTrait — Type
OrthogonalTraitAn abstract type to indicate whether the unit cell is orthogonal or not. It has two concrete types: Orthogonal and NonOrthogonal.
Polyorder.orthogonal — Function
orthogonal(::{<:CrystalSystem}) -> OrthogonalTraitDetermine whether the unit cell is orthogonal or not.
See also OrthogonalTrait.
Polyorder.k2 — Function
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.
Polyorder.kk_tensor — Function
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}.
These functions are device-agnostic: they compute on CPU and copy to kk (which can contain CuArrays).
Missing docstring for Polyorder.kk_tensor!. Check Documenter's build log for details.
Polyorder.stars — Function
stars(w::AuxiliaryField)Generate stars from the grid of the field w. The symmetry operations are obtained from the BravaisLattice instance associated with w.
This is an optimized implementation using type-parameterized dimensions and SVector operations for maximum type stability and minimal heap allocations.
It has been tested for space groups:
- 1D: 1, 2
- 2D: 17
- 3D: 223, 225, 229, 230
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.
symmetrize!(w::AuxiliaryField)Symmetrize the field w according to the symmetry operations of the Bravais lattice. The stars is computed on-the-fly.
symmetrize!(w::AuxiliaryField, stars)Symmetrize the field w according to the symmetry operations of the Bravais lattice. The stars is provided by the stars.
Polyorder.distinct_voigts — Function
Line: [e11]
Oblique: [e11, e22, e12]
Rectangular: [e11, e22, 0]
Square, Hexagonal2D, HexRect: [e11, e11, 0]
Triclinic: [e11, e22, e33, e23, e13, e12]
Monoclinic: [e11, e22, e33, 0, e13, 0]
Orthorhombic: [e11, e22, e33, 0, 0, 0]
Tetragonal, Trigonal, Hexagonal, HexOrthorhombic: [e11, e11, e33, 0, 0, 0]
Cubic: [e11, e11, e11, 0, 0, 0]
Polyorder.num_distinct_voigts — Function
num_distinct_voigts(::CrystalSystem) -> IntReturns the number of distinct voigts in the crystal system.
Polyorder.voigt_to_tensor — Function
voigt_to_tensor(::Union{D1,D2,D3}, i) -> TupleConverts voigt notation to tensor notation.
Polyorder.tensor_to_voigt — Function
tensor_to_voigt(::Union{D1,D2,D3}, i, j) -> TupleConverts tensor notation to voigt notation.
Polyorder.distinct_kk_tensor — Function
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.
Polyorder.distinct_cell_variables — Function
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.
Solve the SCFT model
Polyorder.ConvergenceStatus — Type
ConvergenceStatusAbstract type for representing the convergence status of iterative algorithms.
All convergence status types inherit from this base type and indicate whether the algorithm is still running or has terminated with various outcomes.
Polyorder.Terminated — Type
Terminated <: ConvergenceStatusAbstract type for convergence statuses that indicate the algorithm has terminated. This includes both successful and failed termination conditions.
Polyorder.Running — Type
Running <: ConvergenceStatusIndicates that the iterative algorithm is still running and has not yet reached a termination condition.
Polyorder.Failed — Type
Failed <: TerminatedAbstract type for convergence statuses that indicate the algorithm has failed.
Polyorder.Successful — Type
Successful <: TerminatedIndicates that the algorithm has converged successfully and met all convergence criteria.
Polyorder.Acceptable — Type
Acceptable <: TerminatedIndicates that the algorithm has converged to an acceptable solution, though perhaps not meeting the strictest convergence criteria.
Polyorder.Diverge — Type
Diverge <: FailedIndicates that the algorithm has diverged and the solution is moving away from convergence.
Polyorder.Oscillatory — Type
Oscillatory <: FailedIndicates that the algorithm is oscillating and not making progress toward convergence.
Polyorder.NeedMoreIteration — Type
NeedMoreIteration <: FailedIndicates that the algorithm needs more iterations to reach convergence but has reached the maximum iteration limit.
Polyorder.NeedMoreTime — Type
NeedMoreTime <: FailedIndicates that the algorithm needs more time to reach convergence but has reached the time limit.
Polyorder.LoggingControl — Type
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.
Polyorder.OscillatoryProgress — Type
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.
Polyorder.SlowProgress — Type
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.
Polyorder.ThresholdObjFun — Type
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.
Polyorder.ThresholdStress — Type
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
- 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_lossserves as a minimum requirement for residual norm. - When
tolmode == ResidualTolMode: the residual norm is less than a threshold value (tol_loss).
- When
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.
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.
Arguments
scft::AbstractSCFT: the SCFT model to be solved.updater::SCFTAlgorithm: the SCFT algorithm to be used. If not provided, the default updater inscftis used.config::Config: the configuration for the SCFT simulation. If not provided, the default configurationPolyorder.CONFIGis used.
Keywords
controls=default_iteration_controls(updater, config): the iteration controls to control the SCFT simulation. See thedefault_iteration_controlsfunction for a list of default controls. You can usedefault_iteration_controlsfor 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 twoFvalues is smaller thanconfig.scft.tol(1e-8 or lower is recommended), in addtion to the residual norm is smaller thanconfig.scft.max_tol(1e-3 or lower is recommended). This control is recommended for practical usage, since it can significantly reduce the computation time, asFvalue is much easier to converge than the residual norm.Threshold: stop the simulation when the residual norm is less thanconfig.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 thanconfig.cellopt.gtoland one of the same requirement asThresholdObjFunandThreshold(which is determined byconfig.scft.tolmode) are met. This control is supposed to work for theVariableCellupdater.SlowProgress: stop the simulation when the residual norm reduces slower than a pre-given criterion. This control is recommended to work withThresholdbut notThresholdObjFun.OscillatoryProgress: stop the simulation when oscillation in residual norm is detected. this control is recommended to work withThresholdbut notThresholdObjFun.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, andETDare more stable and robust to use.
IO
- A progress meter can be enabled by setting
config.io.progress_solveto betrue.
Polyorder.good — Function
good(status::ConvergenceStatus) -> BoolReturn true if the convergence status indicates a good (successful or acceptable) outcome, false otherwise.
Examples
good(Successful()) # returns true
good(Acceptable()) # returns true
good(Diverge()) # returns falsePolyorder.status — Function
status(control, loss, maxtol, dangertol) -> ConvergenceStatusDetermine the convergence status based on the iteration control type and loss values.
Arguments
control: The iteration control object that determines stopping criterialoss: Current loss/residual valuemaxtol: Maximum tolerance for acceptable convergencedangertol: Danger tolerance threshold indicating divergence
Returns
A ConvergenceStatus indicating the current state of convergence.
Polyorder.default_iteration_controls — Function
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 exceedsconfig.scft.max_iter.TimLimit: stop if time exceedsconfig.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 thanconfig.cellopt.gtoland the condition ofThresholdObjFun(tolmode isFTolMode) orThreshold(tolmode isResidualTolMode) is satisfied.display_stress: display stress.display_unitcell: display unit cell.
Otherwise, following controls are added:
Threshold: if tolmode isResidualTolMode. Stop if loss is lower thanconfig.scft.tol.ThresholdObjFun: if tolmode isFTolMode. Stop if loss is lower thanconfig.scft.tolandconfig.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 byconfig.scft.use_slow_control.OscillatoryProgress: stop if the decrease of loss is ossilatory. Enabled byconfig.scft.use_osc_control.Patience: stop if the loss does not decrease after the number of iteration exceedsconfig.scft.numpatience. Enabled byconfig.scft.use_patience_control.NumberSinceBest: stop if the number of iteration since the best loss exceedsconfig.scft.numbest. Enabled byconfig.scft.use_nsbest_control.LoggingControl: log the loss and free energy into a CSV file. Enabled byconfig.scft.use_log_control.
SCFT Updaters
Polyorder.SCFTAlgorithm — Type
SCFTAlgorithm{T} <: AbstractAlgorithmAbstract 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.
Polyorder.SpectralSCFTAlgorithm — Type
SpectralSCFTAlgorithm{T} <: SCFTAlgorithm{T}Abstract type for spectral SCFT updaters.
Polyorder.NestedSCFTAlgorithm — Type
NestedSCFTAlgorithm{T} <: SCFTAlgorithm{T}Abstract type for nested SCFT updaters.
Polyorder.Picard — Type
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 asPicard{Float64}().PicardMannorSD: 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.
Polyorder.PicardMann — Type
PicardMann{T} <: SCFTAlgorithm{T}Picard-Mann iteration algorithm for SCFT calculations, also known as steepest descent method. A convenient alias is SD.
This algorithm uses a simple mixing scheme where the new field is computed as: w_new = (1-α) * w_old + α * w_picard, where α is the mixing parameter.
Fields
α::T: Mixing parameter (default 0.2)n::Int: Current iteration numberFs::Vector{T}: History of free energy valuesrs::Vector{T}: History of residual valuesevals::Vector{Int}: History of function evaluation counts
Example
updater = PicardMann(0.1) # Create with mixing parameter 0.1Polyorder.PicardKhan — Type
PicardKhan{T} <: SCFTAlgorithm{T}Picard-Khan iteration algorithm for SCFT calculations.
This is a variant of the Picard iteration scheme with different convergence properties. Uses a mixing parameter to blend old and new field values.
Fields
α::T: Mixing parameter (default 0.1)n::Int: Current iteration numberFs::Vector{T}: History of free energy valuesrs::Vector{T}: History of residual valuesevals::Vector{Int}: History of function evaluation counts
Polyorder.PicardIshikawa — Type
PicardIshikawa{T} <: SCFTAlgorithm{T}Picard-Ishikawa iteration algorithm for SCFT calculations.
This is a two-parameter variant of the Picard iteration scheme that can provide improved convergence properties through the use of two mixing parameters.
Fields
α::T: First mixing parameter (default 0.1)β::T: Second mixing parameter (default 0.1)n::Int: Current iteration numberFs::Vector{T}: History of free energy valuesrs::Vector{T}: History of residual valuesevals::Vector{Int}: History of function evaluation counts
Polyorder.PicardS — Type
PicardS{T} <: SCFTAlgorithm{T}Picard-S iteration algorithm for SCFT calculations.
Another variant of the Picard iteration scheme with two mixing parameters, providing an alternative convergence strategy for SCFT calculations.
Fields
α::T: First mixing parameter (default 0.1)β::T: Second mixing parameter (default 0.1)n::Int: Current iteration numberFs::Vector{T}: History of free energy valuesrs::Vector{T}: History of residual valuesevals::Vector{Int}: History of function evaluation counts
GPU Compatibility
Fully GPU-compatible using VectorOfFields to manage field data on the device.
Polyorder.Euler — Type
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.
Polyorder.EMPEC — Type
EMPEC{T} <: SCFTAlgorithm{T}An SCFT updater based on EMPEC method, which is a predictor-corrector Euler method.
Polyorder.SIS — Type
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.
Polyorder.SISF — Type
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.
GPU Compatibility
The implementation is fully GPU-compatible. The potentially expensive calculation of the relaxation matrix κmat is parallelized using custom GPU kernels (perpixel_sis_kernel!), enabling efficient execution on CUDA devices without scalar indexing bottlenecks.
Polyorder.PO — Type
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.
Polyorder.ETD — Type
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.
GPU Compatibility
The implementation is GPU-compatible. It utilizes dedicated GPU kernels (perpixel_etd_kernel!) to compute the relaxation coefficients in parallel, ensuring high performance on modern hardware.
Polyorder.ETDF — Type
ETDF{T} <: SpectralSCFTAlgorithm{T}Exponential time difference SCFT algorithm using a full matrix of relaxation parameters.
Polyorder.ETDPEC — Type
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.
Polyorder.Nesterov — Type
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.
Polyorder.ARDM — Type
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.
Polyorder.NGMRES_S1 — Type
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.
Polyorder.Anderson_S1 — Type
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 β.
Polyorder.NGMRES — Type
NGMRES{T, S<:SCFTAlgorithm{T}} <: NestedSCFTAlgorithm{T}Nonlinear GMRES acceleration NGMRES(m) preconditioned by non-acceleration methods, such SD, SIS, etc. 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.
GPU Compatibility
Like Anderson acceleration, this NGMRES implementation is fully GPU-compatible. It leverages VectorOfFields for device-agnostic field operations and performs heavy linear algebra (QR update) using optimized dot products, keeping large field data on the GPU.
Polyorder.OACCEL — Type
OACCEL{T, S<:SCFTAlgorithm{T}} <: NestedSCFTAlgorithm{T}Objective acceleration (O-ACCEL) is a slight modification to the NGMRES acceleration. It can be preconditioned by a non-acceleration updater.
References
- Riseth, A. N. Objective Acceleration for Unconstrained Optimization. Numer. Linear Algebra Appl. 2019, 26 (1), e2216.
GPU Compatibility
The implementation is fully GPU-compatible. It avoids constructing large matrices X and G by concatenating fields (Density/Potential) which would be inefficient and difficult on GPU. Instead, it computes the matrix A = X'G and vector γ = X'g(y_k) directly using element-wise dot products of VectorOfFields. This ensures that all heavy lifting remains on the device, and only the small (m+1)×(m+1) linear system is solved on the CPU.
Polyorder.Anderson — Type
Anderson{T, S<:SCFTAlgorithm{T}} <: NestedSCFTAlgorithm{T}Anderson acceleration precoditioned by SD (Euler, PicardMann), EMPEC, SIS, ETD, PO, ETDPEC and other non-acceleration methods which implement a step! method.
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.
GPU Compatibility
The implementation is fully GPU-compatible. It uses VectorOfFields to handle field data abstractly, ensuring that operations are dispatched to the correct backend (e.g., CUDA). The QR decomposition (via qradd! and qrdelete!) and least squares solve are handled efficiently by keeping the large field vectors on the device and only operating on scalar projections (dot products) or small m x m matrices on the CPU.
Polyorder.BB — Type
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:
- αk = norm(Δxk)^2 / Δxk^T Δgk)
- α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.lb: lower bound of the step size. Default: 0. When x < lb * x0, the step is reset to x0, which is the initial x.ub: upper bound of the step size. Default: 2. When x > ub * x0, the step is reset to x0, which is the initial x.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_oldandgk_old: for internal use.
References
- Burdakov, O.; Dai, Y.-H.; Huang, N. Stabilized Barzilai-Borwein Method. arXiv [math.OC], 2019.
Polyorder.VariableCell — Type
VariableCell{T, C<:SCFTAlgorithm{T}, W<:SCFTAlgorithm{T}} <: NestedSCFTAlgorithm{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.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_changedand_warmup: for internal use.
Polyorder.list_scft_algorithms — Function
list_scft_algorithms() -> Vector{Symbol}
list_scft_updaters() -> Vector{Symbol}List all available SCFT algorithms.
Polyorder.select_scft_algorithm — Function
select_scft_algorithm(n::Symbol) -> SCFTAlgorithm
select_scft_updater(n::Symbol) -> SCFTAlgorithmSelect the SCFT updater by name.
Polyorder.create_updater — Function
create_updater(config::Config) -> SCFTAlgorithmCreate an SCFT algorithm from a Config instance. Types not in list_scft_algorithms() are not supported and will throw an ArgumentError.
Polyorder.is_spectral — Function
is_spectral(::SCFTAlgorithm) -> BoolCheck if the SCFT algorithm is spectral or not, of type SpectralSCFTAlgorithm.
Polyorder.is_nested — Function
is_nested(::SCFTAlgorithm) -> BoolCheck if the SCFT algorithm is nested or not, of type NestedSCFTAlgorithm.
Polyorder.inner_algos — Function
inner_algos(::SCFTAlgorithm) -> NamedTupleReturn 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.
Polyorder.options — Function
options(::SCFTAlgorithm) -> NamedTupleReturn the options of the SCFT algorithm. No options returns an empty named tuple. Subtypes of SCFTAlgorithm should implement this function if it has options.
Polyorder.record! — Function
record!(updater::SCFTAlgorithm, scft::AbstractSCFT) -> NothingRecord the free energy and residual of the SCFT system after each update step. This function is called by the update! function to keep track of the history of the algorithm.
Polyorder.step! — Function
step!(updater::SCFTAlgorithm, scft::AbstractSCFT, config::Config) -> NothingPerform one step of the SCFT algorithm, updating the state of the SCFT system in place. This function should be implemented by each specific SCFT algorithm.
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.
GPU Compatibility
This function works seamlessly on GPU if scft holds GPU-backed fields. The q! (propagator solver) and force calculations are automatically dispatched to the device.
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.
Polyorder.update_fft_plans! — Function
update_fft_plans!(::SCFTAlgorithm, ::AbstractArray) -> Nothing
update_fft_plans!(::SCFTAlgorithm, ::AbstractField) -> NothingInplace Updating the FFT plans of the SCFT algorithm.
Polyorder.reset! — Method
reset!(algo::SCFTAlgorithm, w::AbstractField) -> NothingInplace resetting the SCFT algorithm to its initial state, including calling update_fft_plans!.
MDE solvers
Polyorder.MDESolver — Type
MDESolverA 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 fromdsstorage::Symbol=:full: storage type for blocks. :full or :checkpoint.k::Int=0: Recomputation interval along a block.shared_cache::Bool=true: Whether to share the cache among blocks.options::Dict{Symbol, Any}=Dict{Symbol, Any}(): kwargs for the algorithm.
Polyorder.MDEAlgorithm — Type
MDEAlgorithm <: AbstractAlgorithmAbstract type for algorithms to solve the MDE equations, i.e. computing propagators from potential fields.
Polyorder.CxGridMDEAlgorithm — Type
CxGridMDEAlgorithm <: MDEAlgorithmAbstract type for MDE algorithms that operate on 1D Chebyshev grids.
Polyorder.FxCyGridMDEAlgorithm — Type
FxCyGridMDEAlgorithm <: MDEAlgorithmAbstract type for MDE algorithms that operate on 2D grids with Fourier grid for x direction and Chebyshev grid for y direction.
Polyorder.FθCrGridMDEAlgorithm — Type
FθCrGridMDEAlgorithm <: MDEAlgorithmAbstract type for MDE algorithms that operate on grids with angular (θ) and radial (r) coordinates.
Polyorder.OSF — Type
OSF{M, K, S} <: MDEAlgorithmOperator-Splitting pseudo-spectral MDE solver with shared work buffers.
Performance Modes
low_memory=true(default in SCFT): On-the-fly computation of expd/expw, minimal memorylow_memory=false(default standalone): Precomputes expd and expw buffer, faster
Fields
shared::M: Reference to MDEWorkBuffers (u, k2, û, ak)expd::K: Precomputed exp(-b²dsk²), or nothing for on-the-fly modeexpw::S: Buffer for exp(-ds/2*w), or nothing for on-the-fly mode
GPU Compatibility
OSF is fully GPU-compatible. It uses AbstractFFTs interface which dispatches to CUDA.CUFFT when fields are CuArrays, and performs all intermediate steps (exp, mul!) on the device.
Memory Efficiency (low_memory=true)
By sharing work buffers across all OSF solvers in an SCFT simulation:
- Before: N solvers × 5 arrays = 5N full-grid arrays
- After: 1 shared × 4 arrays = 4 full-grid arrays
- Savings: ~69% for typical simulations
Performance (low_memory=false)
Precomputes expd once at construction, computes expw once per solve! call:
- Saves 2*(Ns-1) exp() calls per solve! for expw
- Solvers with same (algo, b, ds) can share expd and expw
Missing docstring for Polyorder.OSFCUDA. Check Documenter's build log for details.
Polyorder.RQM4 — Type
RQM4{M, Q, K, S} <: MDEAlgorithmA MDE solver: Ranjan-Qin-Morse 4th order pseudo-spectral method. Uses Richardson extrapolation with two OSF-style solves at different step sizes.
Performance Modes
low_memory=true(default in SCFT): On-the-fly computation of expd/expw, minimal memorylow_memory=false(default standalone): Precomputes expd and expw buffer, faster
Fields
shared::M: Reference to MDEWorkBuffersq2::Q: Auxiliary propagator for half-step computationexpd::K: Precomputed exp(-b²dsk²), or nothing for low_memory modeexpd2::K: Precomputed exp(-b²(ds/2)k²), or nothing for low_memory modeexpw::S: Buffer for exp(-ds/2*w), or nothing for low_memory modeexpw2::S: Buffer for exp(-ds/4*w), or nothing for low_memory mode
GPU Compatibility
RQM4 uses shared MDEWorkBuffers and inherits full GPU compatibility.
Missing docstring for Polyorder.RQM4CUDA. Check Documenter's build log for details.
Polyorder.FourierFlowsAlgorithmModule.ETDRK4 — Type
ETDRK4 <: FourierFlowsAlgorithmA MDE solver based on the ETDRK4 method in FourierFlows.jl.
Polyorder.MDESmall — Type
MDESmall <: MDEAlgorithmThe MDE solver for small molecules.
Polyorder.ConfinedETDRK4Module.ETDRK4Cx — Type
ETDRK4Cx{B1, B2, T, V, M} <: CxGridMDEAlgorithmAn ETDRK4 solver for 1D confined systems with Chebyshev grids.
This solver operates on Chebyshev grids (non-uniform) and is CPU-only. GPU acceleration is not supported for confined systems.
Fields
Lx::T: The domain size.N::Int: The number of grid points.Ns::Int: The number of steps.h::T: The step size.c::T: The scaling parameter.lbc::B1: The left boundary condition.rbc::B2: The right boundary condition.M::Int: The contour grid size.x::V: The Chebyshev grid points.L::M: The linear operator.E::M: The exponential term.E2::M: The second exponential term.f1::M: The first coefficient.f2::M: The second coefficient.f3::M: The third coefficient.f4::M: The fourth coefficient.f5::M: The fifth coefficient.f6::M: The sixth coefficient.a::V: The first temporary variable.b::V: The second temporary variable.c1::V: The third temporary variable.Fa::V: The first force term.Fb::V: The second force term.Fc::V: The third force term.Fq::V: The force term for the propagator.
Polyorder.ConfinedETDRK4Module.ETDRK4FxCy — Type
The PDE is in 2D, du/dt = cLu - wu where u=u(x,y), L=d^2/dx^2 + d^2/dy^2, w=w(x,y), c is a constant. First, do a FFT in x direction to obtain du(kx,y)/dt = c L u(kx,y) - Fx[w(x,y)u(x,y)] where L = D^2 - kx^2, with D^2 the Chebyshev 2nd order differential matrix, and kx^2 the d^2/dx^2 in Fourier space, see detail in the Notebook (page 2013.8.2).
Test: PASSED?
:param:Lx: physical size of the 1D spacial grid. :param:Lx: physical size of the 1D spacial grid. :param:Ns: number of grid points in time. :param:lbc: left boundary condition. :param:rbc: right boundary condition. :param:h: time step.
Polyorder.ConfinedETDRK4Module.ETDRK4Polar — Type
ETDRK4Polar{B1, B2, T, M, V, A, C} <: FθCrGridMDEAlgorithmAn ETDRK4 solver for 2D polar coordinates.
Fields
R::T: The radius of the domain.Nr::Int: The number of grid points in the r direction.Nθ::Int: The number of grid points in the θ direction.Ns::Int: The number of steps.h::T: The step size.c::T: The scaling parameter.lbc::B1: The left boundary condition.rbc::B2: The right boundary condition.M::Int: The contour grid size.r::V: The Chebyshev grid points.L::M: The linear operator.Lk::M: The linear operator in Fourier space.E::A: The exponential term.E2::A: The second exponential term.f1::A: The first coefficient.f2::A: The second coefficient.f3::A: The third coefficient.f4::A: The fourth coefficient.f5::A: The fifth coefficient.f6::A: The sixth coefficient.ak::C: The first temporary variable in Fourier space.bk::C: The second temporary variable in Fourier space.ck::C: The third temporary variable in Fourier space.
Polyorder.cheb_quadrature_clencurt — Function
cheb_quadrature_clencurt(f, w=nothing)Perform Clenshaw-Curtis quadrature on function values f.
Arguments
f: Vector of function values at Chebyshev nodesw: Optional quadrature weights. Ifnothing, computed usingclencurt_weights_fft
Returns
The integral approximation using Clenshaw-Curtis quadrature.
This function implements the Clenshaw-Curtis quadrature rule, which is based on Chebyshev polynomial interpolation and provides high accuracy for smooth functions.
Note: The weights are automatically moved to the same device as f for GPU compatibility.
Polyorder.list_mde_algorithms — Function
list_mde_algorithms()
list_mde_solvers()List all available MDE solvers.
Polyorder.select_mde_algorithm — Function
select_mde_algorithm(n::Symbol)
select_mde_solver(n::Symbol)Select the MDE solver by name.
Polyorder.best_contour_discretization — Function
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).
Cell optimization
Polyorder.AbstractCellOptAlgorithm — Type
AbstractCellOptAlgorithm <: AbstractAlgorithmAbstract type for cell optimization algorithms.
Polyorder.StressGuidedCellOptAlgorithm — Type
StressGuidedCellOptAlgorithm <: AbstractCellOptAlgorithmAbstract type for stress guided cell optimization algorithms.
Polyorder.GradientFreeCellOptAlgorithm — Type
GradientFreeCellOptAlgorithm <: AbstractCellOptAlgorithmAbstract type for gradient free cell optimization algorithms.
Polyorder.VariableCellOpt — Type
VariableCellOptAlgorithm <: AbstractCellOptAlgorithmCell optimization by the Variable Cell method.
Polyorder.OptimStressGuidedCellOpt — Type
OptimCellOptAlgorithm{O} <: AbstractCellOptAlgorithmStress-guided cell optimization by algorithms from the Optim.jl package.
Polyorder.OptimGradientFreeCellOpt — Type
OptimGradientFreeCellOpt{A, O} <: GradientFreeCellOptAlgorithmGradient-free cell optimization by algorithms from the Optim.jl package.
Polyorder.list_cellopt_algorithms — Function
list_cellopt_algorithms()List available cell optimization algorithms.
Polyorder.select_cellopt_algorithm — Function
select_cellopt_algorithm(n::Symbol)Select a cell optimization algorithm by name.
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.
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.
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 aSCFTAlgorithminstance ornothing. When it is not provided or it isnothing,scft.updateris used. Otherwise,updateris used.config: if not provided, a default one which isPolyorder.Config(; scft=SCFTConfig(; tolmode=:F, tol=1e-8), cellopt=CellOptConfig(; algo1=Symbol("Optim.GradientDescent")))is used.
Tips and Tricks
scftis 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_cellto betrue. - 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_wand/orconfig.io.save_ϕto befalse. - 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.
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
scftis 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.
Boundary Conditions
Polyorder.AbstractBoundaryCondition — Type
AbstractBoundaryConditionAbstract type for all boundary conditions used in field theory calculations. This serves as the base type for both general and specific boundary condition types.
Polyorder.BoundaryCondition — Type
BoundaryCondition <: AbstractBoundaryConditionAbstract type for boundary conditions that can be applied to differential equations. This includes Dirichlet, Neumann, and Robin boundary conditions.
Polyorder.Periodic — Type
PeriodicPeriodic boundary condition.
Polyorder.Dirichlet — Type
DirichletDirichlet boundary condition.
Polyorder.Neumann — Type
NeumannNeumann boundary condition.
Polyorder.Robin — Type
Robin{T}Robin boundary condition.
Polyorder.getβ — Function
getβ(::AbstractBoundaryCondition)Get the β value of a boundary condition.
Configurations
Polyorder.Config — Type
ConfigThe configuration object for the SCFT simulation.
Fields
system::PolymerSystemConfig=PolymerSystemConfig(): the configuration of the polymer system.PolymerSystemConfigis defined in Polymer.jl.lattice::LatticeConfig=LatticeConfig(): the configuration of the lattice.LatticeConfigis 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)Polyorder.VersionConfig — Type
VersionConfigThe configuration object for the version information.
Fields
package::String="Polyorder.jl": the name of the package.version::VersionNumber=version(): the version of the package.
Polyorder.MDEConfig — Type
MDEConfigThe 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 bylist_blockmodes().solvers::Dict{<:Pair, MDESolverConfig}=Dict(): the configurations of the MDE solvers. The key is the block (specified as aPair).
Polyorder.SCFTConfig — Type
SCFTConfigThe configuration object for solving SCFT equations.
Fields
compress::Bool=false: compressible or incompressible SCFT model.fieldmodeltype::Symbol=:simple: the type of field model.:simpleis 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 forIterationControl.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 forIterationControl.max_restart::Int: the maximum number of restarting simulation. Default is1. 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:Residualand:F.tol::Float64=1e-5: the target tolerance.maxtol::Float64=1e-3: the maximum allowed target tolerance, should >= tol. Also used byThresholdObjFunastol_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: useLoggingControlcontrol or not.use_slow_control::Bool=false: useSlowProgresscontrol or not.use_osc_control::Bool=false: useOscillatoryProgresscontrol or not.use_nsbest_control::Bool=false: useNumberSinceBestcontrol or not.use_patience_control::Bool=false: usePatiencecontrol or not.k_slow::Int=100: forSlowProgresscontrol.tol_slow::Float64=√eps(1.0): forSlowProgresscontrol.tolF::Float64=√eps(1.0): forSlowProgresscontrol.k_osc::Int=100: forOscillatoryProgresscontrol.tol_osc::Float64=1.0: forOscillatoryProgresscontrol.m_osc::Int=100: forOscillatoryProgresscontrol.numbest::Int=10: forNumberSinceBestcontrol. number of iterations since best error:skip_iter * numbest.numpatience::Int=10: number of iterations allowed when error increases:skip_iter * numpatience.
Polyorder.CellOptConfig — Type
CellOptConfigThe configuration object for the cell optimization.
Fields
type::Symbol=:VariableCellOpt: the type name of a cell optimization method which isAbstractCellOptAlgorithm.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.
Polyorder.IOConfig — Type
IOConfigThe 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 forsolve!.progress_cell::Bool=false: show progress meter forcell_solve!.base_dir::String=".": the base directory for saving files.summary::String="summary": the file name for the final state ofcell_solve!.trace::String="trace": the file name for the trace ofsolve!andcell_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:HDF5and: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 ofcell_solve!or not.save_trace::Bool=true: save the trace ofsolve!andcell_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.
Polyorder.MDESolverConfig — Type
MDESolverConfigConfiguration for MDE (Modified Diffusion Equation) solvers.
Fields
algo::Symbol: Algorithm type for solving MDE equationsds::Float64: Contour step sizeNs::Int: Number of contour stepsstorage::Symbol: Storage type for contour steps: :full or :checkpoint.k::Int: Recomputation interval.shared_cache::Bool: Whether to share the cache among blocks.options::Dict{Symbol, Any}: Additional algorithm-specific options
Polyorder.SCFTUpdaterConfig — Type
SCFTUpdaterConfigConfiguration for SCFT (Self-Consistent Field Theory) updater algorithms.
Fields
type::Symbol: Type of SCFT updater algorithm (default::SD)inner_algos::Dict{Symbol, SCFTUpdaterConfig}: Nested algorithm configurationsoptions::Dict{Symbol, Any}: Algorithm-specific options and parameters
Polyorder.CellOptAlgoConfig — Type
CellOptAlgoConfigConfiguration for unit cell optimization algorithms.
Fields
type::Symbol: Type of cell optimization algorithm (default::VariableCell)inner_algos::Dict{Symbol, SCFTUpdaterConfig}: Inner algorithm configurations (for VariableCell only)options::Dict{Symbol, Any}: Algorithm-specific options and parameters
Polymer.to_config — Function
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.
Polymer.make — Function
make(config::CellOptConfig) -> AbstractCellOptAlgorithm
from_config(config::CellOptConfig) -> AbstractCellOptAlgorithmConstruct an AbstractCellOptAlgorithm instance from a CellOptConfig instance. Types not in list_cellopt_algorithms() are not supported and will throw an ArgumentError.
make(config::Config) -> NoncyclicChainSCFT
from_config(config::Config) -> NoncyclicChainSCFTConstruct an NoncyclicChainSCFT instance from a Config instance.
make(::ObjType{:System}, config)Make an PolymerSystem object.
Note: confinement is not implmented.
Polymer.from_config — Function
from_config is an alias for make
Polymer.save_config — Function
save_config(yamlfile, config)Save a Configurations.OptionType (defined by @option) object to a YAML file.
Polymer.load_config — Function
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").
Visualization
Polyorder.densityplot — Function
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.
Missing docstring for Polyorder.densityplot!. Check Documenter's build log for details.
Polyorder.traceplot — Function
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 applicablevar=:muandvar=:xs.xlog=false: whether to use log scale on the x-axis.ylog=true: whether to use log scale on the y-axis.
Missing docstring for Polyorder.traceplot!. Check Documenter's build log for details.
IO
Polyorder.SCFTSummary — Type
SCFTSummarySummary of an SCFT model including simulation results if any.
Polyorder.DataFormatTrait — Type
DataFormatTraitAn abstract type for data format. Currently, only HDF5 and MAT formats are supported.
Polyorder.HDF5Format — Type
HDF5Format <: DataFormatTraitHDF5 format.
Polyorder.MATFormat — Type
MATFormat <: DataFormatTraitMAT format from Matlab.
Polyorder.list_dataformats — Function
list_dataformats()List the available data formats.
Polyorder.select_dataformat — Function
select_dataformat(name::Symbol)Select a data format by name. Possible names can be listed by list_dataformats.
Polyorder.save_fields — Function
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.
Polyorder.read_fields — Function
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.
Polyorder.save_densities — Function
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.
Polyorder.read_densities — Function
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.
Polyorder.save_trace_solve — Function
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).
Polyorder.read_trace_solve — Function
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.
Polyorder.save_trace_cell — Function
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]".
Polyorder.read_trace_cell — Function
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.
Utilities
Polyorder.DisplayTimer — Type
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()Polyorder.strip_type_param — Function
strip_type_param(T::Type)
strip_type_param(::T) where TObtain the type without parametric types. A Type is returned.
Example:
strip_type_param(Vector{Float64}) returns Array. strip_type_param(OSF{...}) returns OSF.
Polyorder.nextfastfft — Function
nextfastfft(n)Find the next fast number of points for FFT given the number of points n.
Polyorder.best_N_fft — Function
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.
Polyorder.list_norms — Function
list_norms()Return the list of available norms to be passed into residual.
Polyorder.select_norm — Function
select_norm(name::Symbol)Select a norm by name.
Polyorder.timer — Constant
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 timingsPolyorder.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.
Polyorder.clear_timer! — Function
clear_timer!(timer=timer)Reset the timer.